e1c0e5ccb25af9e36e8532c7211b486afe6acbca
[alexxy/gromacs.git] / src / mdlib / sim_util.c
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
3  *
4  *                This source code is part of
5  *
6  *                 G   R   O   M   A   C   S
7  *
8  *          GROningen MAchine for Chemical Simulations
9  *
10  *                        VERSION 3.2.0
11  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13  * Copyright (c) 2001-2004, The GROMACS development team,
14  * check out http://www.gromacs.org for more information.
15
16  * This program is free software; you can redistribute it and/or
17  * modify it under the terms of the GNU General Public License
18  * as published by the Free Software Foundation; either version 2
19  * of the License, or (at your option) any later version.
20  *
21  * If you want to redistribute modifications, please consider that
22  * scientific software is very special. Version control is crucial -
23  * bugs must be traceable. We will be happy to consider code for
24  * inclusion in the official distribution, but derived work must not
25  * be called official GROMACS. Details are found in the README & COPYING
26  * files - if they are missing, get the official version at www.gromacs.org.
27  *
28  * To help us fund GROMACS development, we humbly ask that you cite
29  * the papers on the package - you can find them in the top README file.
30  *
31  * For more info, check our website at http://www.gromacs.org
32  *
33  * And Hey:
34  * GROwing Monsters And Cloning Shrimps
35  */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
39
40 #ifdef GMX_CRAY_XT3
41 #include<catamount/dclock.h>
42 #endif
43
44
45 #include <stdio.h>
46 #include <time.h>
47 #ifdef HAVE_SYS_TIME_H
48 #include <sys/time.h>
49 #endif
50 #include <math.h>
51 #include "typedefs.h"
52 #include "string2.h"
53 #include "gmxfio.h"
54 #include "smalloc.h"
55 #include "names.h"
56 #include "confio.h"
57 #include "mvdata.h"
58 #include "txtdump.h"
59 #include "pbc.h"
60 #include "chargegroup.h"
61 #include "vec.h"
62 #include <time.h>
63 #include "nrnb.h"
64 #include "mshift.h"
65 #include "mdrun.h"
66 #include "sim_util.h"
67 #include "update.h"
68 #include "physics.h"
69 #include "main.h"
70 #include "mdatoms.h"
71 #include "force.h"
72 #include "bondf.h"
73 #include "pme.h"
74 #include "disre.h"
75 #include "orires.h"
76 #include "network.h"
77 #include "calcmu.h"
78 #include "constr.h"
79 #include "xvgr.h"
80 #include "trnio.h"
81 #include "xtcio.h"
82 #include "copyrite.h"
83 #include "pull_rotation.h"
84 #include "gmx_random.h"
85 #include "mpelogging.h"
86 #include "domdec.h"
87 #include "partdec.h"
88 #include "gmx_wallcycle.h"
89 #include "genborn.h"
90 #include "nbnxn_atomdata.h"
91 #include "nbnxn_search.h"
92 #include "nbnxn_kernels/nbnxn_kernel_ref.h"
93 #include "nbnxn_kernels/nbnxn_kernel_x86_simd128.h"
94 #include "nbnxn_kernels/nbnxn_kernel_x86_simd256.h"
95 #include "nbnxn_kernels/nbnxn_kernel_gpu_ref.h"
96
97 #ifdef GMX_LIB_MPI
98 #include <mpi.h>
99 #endif
100 #ifdef GMX_THREAD_MPI
101 #include "tmpi.h"
102 #endif
103
104 #include "adress.h"
105 #include "qmmm.h"
106
107 #include "nbnxn_cuda_data_mgmt.h"
108 #include "nbnxn_cuda/nbnxn_cuda.h"
109
110 #if 0
111 typedef struct gmx_timeprint {
112
113 } t_gmx_timeprint;
114 #endif
115
116 /* Portable version of ctime_r implemented in src/gmxlib/string2.c, but we do not want it declared in public installed headers */
117 char *
118 gmx_ctime_r(const time_t *clock,char *buf, int n);
119
120
121 double
122 gmx_gettime()
123 {
124 #ifdef HAVE_GETTIMEOFDAY
125         struct timeval t;
126         double seconds;
127
128         gettimeofday(&t,NULL);
129
130         seconds = (double) t.tv_sec + 1e-6*(double)t.tv_usec;
131
132         return seconds;
133 #else
134         double  seconds;
135
136         seconds = time(NULL);
137
138         return seconds;
139 #endif
140 }
141
142
143 #define difftime(end,start) ((double)(end)-(double)(start))
144
145 void print_time(FILE *out,gmx_runtime_t *runtime,gmx_large_int_t step,
146                 t_inputrec *ir, t_commrec *cr)
147 {
148     time_t finish;
149     char   timebuf[STRLEN];
150     double dt;
151     char buf[48];
152
153 #ifndef GMX_THREAD_MPI
154     if (!PAR(cr))
155 #endif
156     {
157         fprintf(out,"\r");
158     }
159     fprintf(out,"step %s",gmx_step_str(step,buf));
160     if ((step >= ir->nstlist))
161     {
162         runtime->last = gmx_gettime();
163         dt = difftime(runtime->last,runtime->real);
164         runtime->time_per_step = dt/(step - ir->init_step + 1);
165
166         dt = (ir->nsteps + ir->init_step - step)*runtime->time_per_step;
167
168         if (ir->nsteps >= 0)
169         {
170             if (dt >= 300)
171             {
172                 finish = (time_t) (runtime->last + dt);
173                 gmx_ctime_r(&finish,timebuf,STRLEN);
174                 sprintf(buf,"%s",timebuf);
175                 buf[strlen(buf)-1]='\0';
176                 fprintf(out,", will finish %s",buf);
177             }
178             else
179                 fprintf(out,", remaining runtime: %5d s          ",(int)dt);
180         }
181         else
182         {
183             fprintf(out," performance: %.1f ns/day    ",
184                     ir->delta_t/1000*24*60*60/runtime->time_per_step);
185         }
186     }
187 #ifndef GMX_THREAD_MPI
188     if (PAR(cr))
189     {
190         fprintf(out,"\n");
191     }
192 #endif
193
194     fflush(out);
195 }
196
197 #ifdef NO_CLOCK
198 #define clock() -1
199 #endif
200
201 static double set_proctime(gmx_runtime_t *runtime)
202 {
203     double diff;
204 #ifdef GMX_CRAY_XT3
205     double prev;
206
207     prev = runtime->proc;
208     runtime->proc = dclock();
209
210     diff = runtime->proc - prev;
211 #else
212     clock_t prev;
213
214     prev = runtime->proc;
215     runtime->proc = clock();
216
217     diff = (double)(runtime->proc - prev)/(double)CLOCKS_PER_SEC;
218 #endif
219     if (diff < 0)
220     {
221         /* The counter has probably looped, ignore this data */
222         diff = 0;
223     }
224
225     return diff;
226 }
227
228 void runtime_start(gmx_runtime_t *runtime)
229 {
230     runtime->real = gmx_gettime();
231     runtime->proc          = 0;
232     set_proctime(runtime);
233     runtime->realtime      = 0;
234     runtime->proctime      = 0;
235     runtime->last          = 0;
236     runtime->time_per_step = 0;
237 }
238
239 void runtime_end(gmx_runtime_t *runtime)
240 {
241     double now;
242
243     now = gmx_gettime();
244
245     runtime->proctime += set_proctime(runtime);
246     runtime->realtime  = now - runtime->real;
247     runtime->real      = now;
248 }
249
250 void runtime_upd_proc(gmx_runtime_t *runtime)
251 {
252     runtime->proctime += set_proctime(runtime);
253 }
254
255 void print_date_and_time(FILE *fplog,int nodeid,const char *title,
256                          const gmx_runtime_t *runtime)
257 {
258     int i;
259     char timebuf[STRLEN];
260     char time_string[STRLEN];
261     time_t tmptime;
262
263     if (fplog)
264     {
265         if (runtime != NULL)
266         {
267             tmptime = (time_t) runtime->real;
268             gmx_ctime_r(&tmptime,timebuf,STRLEN);
269         }
270         else
271         {
272             tmptime = (time_t) gmx_gettime();
273             gmx_ctime_r(&tmptime,timebuf,STRLEN);
274         }
275         for(i=0; timebuf[i]>=' '; i++)
276         {
277             time_string[i]=timebuf[i];
278         }
279         time_string[i]='\0';
280
281         fprintf(fplog,"%s on node %d %s\n",title,nodeid,time_string);
282     }
283 }
284
285 static void sum_forces(int start,int end,rvec f[],rvec flr[])
286 {
287   int i;
288
289   if (gmx_debug_at) {
290     pr_rvecs(debug,0,"fsr",f+start,end-start);
291     pr_rvecs(debug,0,"flr",flr+start,end-start);
292   }
293   for(i=start; (i<end); i++)
294     rvec_inc(f[i],flr[i]);
295 }
296
297 /*
298  * calc_f_el calculates forces due to an electric field.
299  *
300  * force is kJ mol^-1 nm^-1 = e * kJ mol^-1 nm^-1 / e
301  *
302  * Et[] contains the parameters for the time dependent
303  * part of the field (not yet used).
304  * Ex[] contains the parameters for
305  * the spatial dependent part of the field. You can have cool periodic
306  * fields in principle, but only a constant field is supported
307  * now.
308  * The function should return the energy due to the electric field
309  * (if any) but for now returns 0.
310  *
311  * WARNING:
312  * There can be problems with the virial.
313  * Since the field is not self-consistent this is unavoidable.
314  * For neutral molecules the virial is correct within this approximation.
315  * For neutral systems with many charged molecules the error is small.
316  * But for systems with a net charge or a few charged molecules
317  * the error can be significant when the field is high.
318  * Solution: implement a self-consitent electric field into PME.
319  */
320 static void calc_f_el(FILE *fp,int  start,int homenr,
321                       real charge[],rvec x[],rvec f[],
322                       t_cosines Ex[],t_cosines Et[],double t)
323 {
324     rvec Ext;
325     real t0;
326     int  i,m;
327
328     for(m=0; (m<DIM); m++)
329     {
330         if (Et[m].n > 0)
331         {
332             if (Et[m].n == 3)
333             {
334                 t0 = Et[m].a[1];
335                 Ext[m] = cos(Et[m].a[0]*(t-t0))*exp(-sqr(t-t0)/(2.0*sqr(Et[m].a[2])));
336             }
337             else
338             {
339                 Ext[m] = cos(Et[m].a[0]*t);
340             }
341         }
342         else
343         {
344             Ext[m] = 1.0;
345         }
346         if (Ex[m].n > 0)
347         {
348             /* Convert the field strength from V/nm to MD-units */
349             Ext[m] *= Ex[m].a[0]*FIELDFAC;
350             for(i=start; (i<start+homenr); i++)
351                 f[i][m] += charge[i]*Ext[m];
352         }
353         else
354         {
355             Ext[m] = 0;
356         }
357     }
358     if (fp != NULL)
359     {
360         fprintf(fp,"%10g  %10g  %10g  %10g #FIELD\n",t,
361                 Ext[XX]/FIELDFAC,Ext[YY]/FIELDFAC,Ext[ZZ]/FIELDFAC);
362     }
363 }
364
365 static void calc_virial(FILE *fplog,int start,int homenr,rvec x[],rvec f[],
366                         tensor vir_part,t_graph *graph,matrix box,
367                         t_nrnb *nrnb,const t_forcerec *fr,int ePBC)
368 {
369   int i,j;
370   tensor virtest;
371
372   /* The short-range virial from surrounding boxes */
373   clear_mat(vir_part);
374   calc_vir(fplog,SHIFTS,fr->shift_vec,fr->fshift,vir_part,ePBC==epbcSCREW,box);
375   inc_nrnb(nrnb,eNR_VIRIAL,SHIFTS);
376
377   /* Calculate partial virial, for local atoms only, based on short range.
378    * Total virial is computed in global_stat, called from do_md
379    */
380   f_calc_vir(fplog,start,start+homenr,x,f,vir_part,graph,box);
381   inc_nrnb(nrnb,eNR_VIRIAL,homenr);
382
383   /* Add position restraint contribution */
384   for(i=0; i<DIM; i++) {
385     vir_part[i][i] += fr->vir_diag_posres[i];
386   }
387
388   /* Add wall contribution */
389   for(i=0; i<DIM; i++) {
390     vir_part[i][ZZ] += fr->vir_wall_z[i];
391   }
392
393   if (debug)
394     pr_rvecs(debug,0,"vir_part",vir_part,DIM);
395 }
396
397 static void posres_wrapper(FILE *fplog,
398                            int flags,
399                            gmx_bool bSepDVDL,
400                            t_inputrec *ir,
401                            t_nrnb *nrnb,
402                            gmx_localtop_t *top,
403                            matrix box,rvec x[],
404                            rvec f[],
405                            gmx_enerdata_t *enerd,
406                            real *lambda,
407                            t_forcerec *fr)
408 {
409     t_pbc pbc;
410     real  v,dvdl;
411     int   i;
412
413     /* Position restraints always require full pbc */
414     set_pbc(&pbc,ir->ePBC,box);
415     dvdl = 0;
416     v = posres(top->idef.il[F_POSRES].nr,top->idef.il[F_POSRES].iatoms,
417                top->idef.iparams_posres,
418                (const rvec*)x,fr->f_novirsum,fr->vir_diag_posres,
419                ir->ePBC==epbcNONE ? NULL : &pbc,
420                lambda[efptRESTRAINT],&dvdl,
421                fr->rc_scaling,fr->ePBC,fr->posres_com,fr->posres_comB);
422     if (bSepDVDL)
423     {
424         fprintf(fplog,sepdvdlformat,
425                 interaction_function[F_POSRES].longname,v,dvdl);
426     }
427     enerd->term[F_POSRES] += v;
428     /* If just the force constant changes, the FEP term is linear,
429      * but if k changes, it is not.
430      */
431     enerd->dvdl_nonlin[efptRESTRAINT] += dvdl;
432     inc_nrnb(nrnb,eNR_POSRES,top->idef.il[F_POSRES].nr/2);
433
434     if ((ir->fepvals->n_lambda > 0) && (flags & GMX_FORCE_DHDL))
435     {
436         for(i=0; i<enerd->n_lambda; i++)
437         {
438             real dvdl_dum,lambda_dum;
439
440             lambda_dum = (i==0 ? lambda[efptRESTRAINT] : ir->fepvals->all_lambda[efptRESTRAINT][i-1]);
441             v = posres(top->idef.il[F_POSRES].nr,top->idef.il[F_POSRES].iatoms,
442                        top->idef.iparams_posres,
443                        (const rvec*)x,NULL,NULL,
444                        ir->ePBC==epbcNONE ? NULL : &pbc,lambda_dum,&dvdl,
445                        fr->rc_scaling,fr->ePBC,fr->posres_com,fr->posres_comB);
446             enerd->enerpart_lambda[i] += v;
447         }
448     }
449 }
450
451 static void pull_potential_wrapper(FILE *fplog,
452                                    gmx_bool bSepDVDL,
453                                    t_commrec *cr,
454                                    t_inputrec *ir,
455                                    matrix box,rvec x[],
456                                    rvec f[],
457                                    tensor vir_force,
458                                    t_mdatoms *mdatoms,
459                                    gmx_enerdata_t *enerd,
460                                    real *lambda,
461                                    double t)
462 {
463     t_pbc  pbc;
464     real   dvdl;
465
466     /* Calculate the center of mass forces, this requires communication,
467      * which is why pull_potential is called close to other communication.
468      * The virial contribution is calculated directly,
469      * which is why we call pull_potential after calc_virial.
470      */
471     set_pbc(&pbc,ir->ePBC,box);
472     dvdl = 0; 
473     enerd->term[F_COM_PULL] +=
474         pull_potential(ir->ePull,ir->pull,mdatoms,&pbc,
475                        cr,t,lambda[efptRESTRAINT],x,f,vir_force,&dvdl);
476     if (bSepDVDL)
477     {
478         fprintf(fplog,sepdvdlformat,"Com pull",enerd->term[F_COM_PULL],dvdl);
479     }
480     enerd->dvdl_lin[efptRESTRAINT] += dvdl;
481 }
482
483 static void pme_receive_force_ener(FILE *fplog,
484                                    gmx_bool bSepDVDL,
485                                    t_commrec *cr,
486                                    gmx_wallcycle_t wcycle,
487                                    gmx_enerdata_t *enerd,
488                                    t_forcerec *fr)
489 {
490     real   e,v,dvdl;    
491     float  cycles_ppdpme,cycles_seppme;
492
493     cycles_ppdpme = wallcycle_stop(wcycle,ewcPPDURINGPME);
494     dd_cycles_add(cr->dd,cycles_ppdpme,ddCyclPPduringPME);
495
496     /* In case of node-splitting, the PP nodes receive the long-range 
497      * forces, virial and energy from the PME nodes here.
498      */    
499     wallcycle_start(wcycle,ewcPP_PMEWAITRECVF);
500     dvdl = 0;
501     gmx_pme_receive_f(cr,fr->f_novirsum,fr->vir_el_recip,&e,&dvdl,
502                       &cycles_seppme);
503     if (bSepDVDL)
504     {
505         fprintf(fplog,sepdvdlformat,"PME mesh",e,dvdl);
506     }
507     enerd->term[F_COUL_RECIP] += e;
508     enerd->dvdl_lin[efptCOUL] += dvdl;
509     if (wcycle)
510     {
511         dd_cycles_add(cr->dd,cycles_seppme,ddCyclPME);
512     }
513     wallcycle_stop(wcycle,ewcPP_PMEWAITRECVF);
514 }
515
516 static void print_large_forces(FILE *fp,t_mdatoms *md,t_commrec *cr,
517                                gmx_large_int_t step,real pforce,rvec *x,rvec *f)
518 {
519   int  i;
520   real pf2,fn2;
521   char buf[STEPSTRSIZE];
522
523   pf2 = sqr(pforce);
524   for(i=md->start; i<md->start+md->homenr; i++) {
525     fn2 = norm2(f[i]);
526     /* We also catch NAN, if the compiler does not optimize this away. */
527     if (fn2 >= pf2 || fn2 != fn2) {
528       fprintf(fp,"step %s  atom %6d  x %8.3f %8.3f %8.3f  force %12.5e\n",
529               gmx_step_str(step,buf),
530               ddglatnr(cr->dd,i),x[i][XX],x[i][YY],x[i][ZZ],sqrt(fn2));
531     }
532   }
533 }
534
535 static void post_process_forces(FILE *fplog,
536                                 t_commrec *cr,
537                                 gmx_large_int_t step,
538                                 t_nrnb *nrnb,gmx_wallcycle_t wcycle,
539                                 gmx_localtop_t *top,
540                                 matrix box,rvec x[],
541                                 rvec f[],
542                                 tensor vir_force,
543                                 t_mdatoms *mdatoms,
544                                 t_graph *graph,
545                                 t_forcerec *fr,gmx_vsite_t *vsite,
546                                 int flags)
547 {
548     if (fr->bF_NoVirSum)
549     {
550         if (vsite)
551         {
552             /* Spread the mesh force on virtual sites to the other particles... 
553              * This is parallellized. MPI communication is performed
554              * if the constructing atoms aren't local.
555              */
556             wallcycle_start(wcycle,ewcVSITESPREAD);
557             spread_vsite_f(fplog,vsite,x,fr->f_novirsum,NULL,
558                            (flags & GMX_FORCE_VIRIAL),fr->vir_el_recip,
559                            nrnb,
560                            &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
561             wallcycle_stop(wcycle,ewcVSITESPREAD);
562         }
563         if (flags & GMX_FORCE_VIRIAL)
564         {
565             /* Now add the forces, this is local */
566             if (fr->bDomDec)
567             {
568                 sum_forces(0,fr->f_novirsum_n,f,fr->f_novirsum);
569             }
570             else
571             {
572                 sum_forces(mdatoms->start,mdatoms->start+mdatoms->homenr,
573                            f,fr->f_novirsum);
574             }
575             if (EEL_FULL(fr->eeltype))
576             {
577                 /* Add the mesh contribution to the virial */
578                 m_add(vir_force,fr->vir_el_recip,vir_force);
579             }
580             if (debug)
581             {
582                 pr_rvecs(debug,0,"vir_force",vir_force,DIM);
583             }
584         }
585     }
586     
587     if (fr->print_force >= 0)
588     {
589         print_large_forces(stderr,mdatoms,cr,step,fr->print_force,x,f);
590     }
591 }
592
593 static void do_nb_verlet(t_forcerec *fr,
594                          interaction_const_t *ic,
595                          gmx_enerdata_t *enerd,
596                          int flags, int ilocality,
597                          int clearF,
598                          t_nrnb *nrnb,
599                          gmx_wallcycle_t wcycle)
600 {
601     int     nnbl, kernel_type, sh_e;
602     char    *env;
603     nonbonded_verlet_group_t  *nbvg;
604
605     if (!(flags & GMX_FORCE_NONBONDED))
606     {
607         /* skip non-bonded calculation */
608         return;
609     }
610
611     nbvg = &fr->nbv->grp[ilocality];
612
613     /* CUDA kernel launch overhead is already timed separately */
614     if (fr->cutoff_scheme != ecutsVERLET)
615     {
616         gmx_incons("Invalid cut-off scheme passed!");
617     }
618
619     if (nbvg->kernel_type != nbk8x8x8_CUDA)
620     {
621         wallcycle_sub_start(wcycle, ewcsNONBONDED);
622     }
623     switch (nbvg->kernel_type)
624     {
625         case nbk4x4_PlainC:
626             nbnxn_kernel_ref(&nbvg->nbl_lists,
627                              nbvg->nbat, ic,
628                              fr->shift_vec,
629                              flags,
630                              clearF,
631                              fr->fshift[0],
632                              enerd->grpp.ener[egCOULSR],
633                              fr->bBHAM ?
634                              enerd->grpp.ener[egBHAMSR] :
635                              enerd->grpp.ener[egLJSR]);
636             break;
637         
638         case nbk4xN_X86_SIMD128:
639             nbnxn_kernel_x86_simd128(&nbvg->nbl_lists,
640                                      nbvg->nbat, ic,
641                                      fr->shift_vec,
642                                      flags,
643                                      clearF,
644                                      fr->fshift[0],
645                                      enerd->grpp.ener[egCOULSR],
646                                      fr->bBHAM ?
647                                      enerd->grpp.ener[egBHAMSR] :
648                                      enerd->grpp.ener[egLJSR]);
649             break;
650         case nbk4xN_X86_SIMD256:
651             nbnxn_kernel_x86_simd256(&nbvg->nbl_lists,
652                                      nbvg->nbat, ic,
653                                      fr->shift_vec,
654                                      flags,
655                                      clearF,
656                                      fr->fshift[0],
657                                      enerd->grpp.ener[egCOULSR],
658                                      fr->bBHAM ?
659                                      enerd->grpp.ener[egBHAMSR] :
660                                      enerd->grpp.ener[egLJSR]);
661             break;
662
663         case nbk8x8x8_CUDA:
664             nbnxn_cuda_launch_kernel(fr->nbv->cu_nbv, nbvg->nbat, flags, ilocality);
665             break;
666
667         case nbk8x8x8_PlainC:
668             nbnxn_kernel_gpu_ref(nbvg->nbl_lists.nbl[0],
669                                  nbvg->nbat, ic,
670                                  fr->shift_vec,
671                                  flags,
672                                  clearF,
673                                  nbvg->nbat->out[0].f,
674                                  fr->fshift[0],
675                                  enerd->grpp.ener[egCOULSR],
676                                  fr->bBHAM ?
677                                  enerd->grpp.ener[egBHAMSR] :
678                                  enerd->grpp.ener[egLJSR]);
679             break;
680
681         default:
682             gmx_incons("Invalid nonbonded kernel type passed!");
683
684     }
685     if (nbvg->kernel_type != nbk8x8x8_CUDA)
686     {
687         wallcycle_sub_stop(wcycle, ewcsNONBONDED);
688     }
689
690     /* In eNR_??? the nbnxn F+E kernels are always the F kernel + 1 */
691     sh_e = ((flags & GMX_FORCE_ENERGY) ? 1 : 0);
692     inc_nrnb(nrnb,
693              ((EEL_RF(ic->eeltype) || ic->eeltype == eelCUT) ?
694               eNR_NBNXN_LJ_RF : eNR_NBNXN_LJ_TAB) + sh_e,
695              nbvg->nbl_lists.natpair_ljq);
696     inc_nrnb(nrnb,eNR_NBNXN_LJ+sh_e,nbvg->nbl_lists.natpair_lj);
697     inc_nrnb(nrnb,
698              ((EEL_RF(ic->eeltype) || ic->eeltype == eelCUT) ?
699               eNR_NBNXN_RF : eNR_NBNXN_TAB)+sh_e,
700              nbvg->nbl_lists.natpair_q);
701 }
702
703 void do_force_cutsVERLET(FILE *fplog,t_commrec *cr,
704               t_inputrec *inputrec,
705               gmx_large_int_t step,t_nrnb *nrnb,gmx_wallcycle_t wcycle,
706               gmx_localtop_t *top,
707               gmx_mtop_t *mtop,
708               gmx_groups_t *groups,
709               matrix box,rvec x[],history_t *hist,
710               rvec f[],
711               tensor vir_force,
712               t_mdatoms *mdatoms,
713               gmx_enerdata_t *enerd,t_fcdata *fcd,
714               real *lambda,t_graph *graph,
715               t_forcerec *fr, interaction_const_t *ic,
716               gmx_vsite_t *vsite,rvec mu_tot,
717               double t,FILE *field,gmx_edsam_t ed,
718               gmx_bool bBornRadii,
719               int flags)
720 {
721     int     cg0,cg1,i,j;
722     int     start,homenr;
723     int     nb_kernel_type;
724     double  mu[2*DIM];
725     gmx_bool   bSepDVDL,bStateChanged,bNS,bFillGrid,bCalcCGCM,bBS;
726     gmx_bool   bDoLongRange,bDoForces,bSepLRF,bUseGPU,bUseOrEmulGPU;
727     gmx_bool   bDiffKernels=FALSE;
728     matrix  boxs;
729     rvec    vzero,box_diag;
730     real    e,v,dvdl;
731     float  cycles_pme,cycles_force;
732     nonbonded_verlet_t *nbv;
733
734     cycles_force = 0;
735     nbv = fr->nbv;
736     nb_kernel_type = fr->nbv->grp[0].kernel_type;
737
738     start  = mdatoms->start;
739     homenr = mdatoms->homenr;
740
741     bSepDVDL = (fr->bSepDVDL && do_per_step(step,inputrec->nstlog));
742
743     clear_mat(vir_force);
744
745     cg0 = 0;
746     if (DOMAINDECOMP(cr))
747     {
748         cg1 = cr->dd->ncg_tot;
749     }
750     else
751     {
752         cg1 = top->cgs.nr;
753     }
754     if (fr->n_tpi > 0)
755     {
756         cg1--;
757     }
758
759     bStateChanged = (flags & GMX_FORCE_STATECHANGED);
760     bNS           = (flags & GMX_FORCE_NS) && (fr->bAllvsAll==FALSE); 
761     bFillGrid     = (bNS && bStateChanged);
762     bCalcCGCM     = (bFillGrid && !DOMAINDECOMP(cr));
763     bDoLongRange  = (fr->bTwinRange && bNS && (flags & GMX_FORCE_DO_LR));
764     bDoForces     = (flags & GMX_FORCE_FORCES);
765     bSepLRF       = (bDoLongRange && bDoForces && (flags & GMX_FORCE_SEPLRF));
766     bUseGPU       = fr->nbv->bUseGPU;
767     bUseOrEmulGPU = bUseGPU || (nbv->grp[0].kernel_type == nbk8x8x8_PlainC);
768
769     if (bStateChanged)
770     {
771         update_forcerec(fplog,fr,box);
772
773         if (NEED_MUTOT(*inputrec))
774         {
775             /* Calculate total (local) dipole moment in a temporary common array.
776              * This makes it possible to sum them over nodes faster.
777              */
778             calc_mu(start,homenr,
779                     x,mdatoms->chargeA,mdatoms->chargeB,mdatoms->nChargePerturbed,
780                     mu,mu+DIM);
781         }
782     }
783
784     if (fr->ePBC != epbcNONE) { 
785         /* Compute shift vectors every step,
786          * because of pressure coupling or box deformation!
787          */
788         if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
789             calc_shifts(box,fr->shift_vec);
790
791         if (bCalcCGCM) { 
792             put_atoms_in_box_omp(fr->ePBC,box,homenr,x);
793             inc_nrnb(nrnb,eNR_SHIFTX,homenr);
794         } 
795         else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph) {
796             unshift_self(graph,box,x);
797         }
798     } 
799
800     nbnxn_atomdata_copy_shiftvec(flags & GMX_FORCE_DYNAMICBOX,
801                                   fr->shift_vec,nbv->grp[0].nbat);
802
803 #ifdef GMX_MPI
804     if (!(cr->duty & DUTY_PME)) {
805         /* Send particle coordinates to the pme nodes.
806          * Since this is only implemented for domain decomposition
807          * and domain decomposition does not use the graph,
808          * we do not need to worry about shifting.
809          */    
810
811         wallcycle_start(wcycle,ewcPP_PMESENDX);
812         GMX_MPE_LOG(ev_send_coordinates_start);
813
814         bBS = (inputrec->nwall == 2);
815         if (bBS) {
816             copy_mat(box,boxs);
817             svmul(inputrec->wall_ewald_zfac,boxs[ZZ],boxs[ZZ]);
818         }
819
820         gmx_pme_send_x(cr,bBS ? boxs : box,x,
821                        mdatoms->nChargePerturbed,lambda[efptCOUL],
822                        (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),step);
823
824         GMX_MPE_LOG(ev_send_coordinates_finish);
825         wallcycle_stop(wcycle,ewcPP_PMESENDX);
826     }
827 #endif /* GMX_MPI */
828
829     /* do gridding for pair search */
830     if (bNS)
831     {
832         if (graph && bStateChanged)
833         {
834             /* Calculate intramolecular shift vectors to make molecules whole */
835             mk_mshift(fplog,graph,fr->ePBC,box,x);
836         }
837
838         clear_rvec(vzero);
839         box_diag[XX] = box[XX][XX];
840         box_diag[YY] = box[YY][YY];
841         box_diag[ZZ] = box[ZZ][ZZ];
842
843         wallcycle_start(wcycle,ewcNS);
844         if (!fr->bDomDec)
845         {
846             wallcycle_sub_start(wcycle,ewcsNBS_GRID_LOCAL);
847             nbnxn_put_on_grid(nbv->nbs,fr->ePBC,box,
848                               0,vzero,box_diag,
849                               0,mdatoms->homenr,-1,fr->cginfo,x,
850                               0,NULL,
851                               nbv->grp[eintLocal].kernel_type,
852                               nbv->grp[eintLocal].nbat);
853             wallcycle_sub_stop(wcycle,ewcsNBS_GRID_LOCAL);
854         }
855         else
856         {
857             wallcycle_sub_start(wcycle,ewcsNBS_GRID_NONLOCAL);
858             nbnxn_put_on_grid_nonlocal(nbv->nbs,domdec_zones(cr->dd),
859                                        fr->cginfo,x,
860                                        nbv->grp[eintNonlocal].kernel_type,
861                                        nbv->grp[eintNonlocal].nbat);
862             wallcycle_sub_stop(wcycle,ewcsNBS_GRID_NONLOCAL);
863         }
864
865         if (nbv->ngrp == 1 ||
866             nbv->grp[eintNonlocal].nbat == nbv->grp[eintLocal].nbat)
867         {
868             nbnxn_atomdata_set(nbv->grp[eintLocal].nbat,eatAll,
869                                 nbv->nbs,mdatoms,fr->cginfo);
870         }
871         else
872         {
873             nbnxn_atomdata_set(nbv->grp[eintLocal].nbat,eatLocal,
874                                 nbv->nbs,mdatoms,fr->cginfo);
875             nbnxn_atomdata_set(nbv->grp[eintNonlocal].nbat,eatAll,
876                                 nbv->nbs,mdatoms,fr->cginfo);
877         }
878         wallcycle_stop(wcycle, ewcNS);
879     }
880
881     /* initialize the GPU atom data and copy shift vector */
882     if (bUseGPU)
883     {
884         if (bNS)
885         {
886             wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
887             nbnxn_cuda_init_atomdata(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
888             wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
889         }
890
891         wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
892         nbnxn_cuda_upload_shiftvec(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
893         wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
894     }
895
896     /* do local pair search */
897     if (bNS)
898     {
899         wallcycle_start_nocount(wcycle,ewcNS);
900         wallcycle_sub_start(wcycle,ewcsNBS_SEARCH_LOCAL);
901         nbnxn_make_pairlist(nbv->nbs,nbv->grp[eintLocal].nbat,
902                             &top->excls,
903                             ic->rlist,
904                             nbv->min_ci_balanced,
905                             &nbv->grp[eintLocal].nbl_lists,
906                             eintLocal,
907                             nbv->grp[eintLocal].kernel_type,
908                             nrnb);
909         wallcycle_sub_stop(wcycle,ewcsNBS_SEARCH_LOCAL);
910
911         if (bUseGPU)
912         {
913             /* initialize local pair-list on the GPU */
914             nbnxn_cuda_init_pairlist(nbv->cu_nbv,
915                                      nbv->grp[eintLocal].nbl_lists.nbl[0],
916                                      eintLocal);
917         }
918         wallcycle_stop(wcycle, ewcNS);
919     }
920     else
921     {
922         wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
923         wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
924         nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs,eatLocal,FALSE,x,
925                                         nbv->grp[eintLocal].nbat);
926         wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
927         wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
928     }
929
930     if (bUseGPU)
931     {
932         wallcycle_start(wcycle,ewcLAUNCH_GPU_NB);
933         /* launch local nonbonded F on GPU */
934         do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFNo,
935                      nrnb, wcycle);
936         wallcycle_stop(wcycle,ewcLAUNCH_GPU_NB);
937     }
938
939     /* Communicate coordinates and sum dipole if necessary + 
940        do non-local pair search */
941     if (DOMAINDECOMP(cr))
942     {
943         bDiffKernels = (nbv->grp[eintNonlocal].kernel_type !=
944                         nbv->grp[eintLocal].kernel_type);
945
946         if (bDiffKernels)
947         {
948             /* With GPU+CPU non-bonded calculations we need to copy
949              * the local coordinates to the non-local nbat struct
950              * (in CPU format) as the non-local kernel call also
951              * calculates the local - non-local interactions.
952              */
953             wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
954             wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
955             nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs,eatLocal,TRUE,x,
956                                              nbv->grp[eintNonlocal].nbat);
957             wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
958             wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
959         }
960
961         if (bNS)
962         {
963             wallcycle_start_nocount(wcycle,ewcNS);
964             wallcycle_sub_start(wcycle,ewcsNBS_SEARCH_NONLOCAL);
965
966             if (bDiffKernels)
967             {
968                 nbnxn_grid_add_simple(nbv->nbs,nbv->grp[eintNonlocal].nbat);
969             }
970
971             nbnxn_make_pairlist(nbv->nbs,nbv->grp[eintNonlocal].nbat,
972                                 &top->excls,
973                                 ic->rlist,
974                                 nbv->min_ci_balanced,
975                                 &nbv->grp[eintNonlocal].nbl_lists,
976                                 eintNonlocal,
977                                 nbv->grp[eintNonlocal].kernel_type,
978                                 nrnb);
979
980             wallcycle_sub_stop(wcycle,ewcsNBS_SEARCH_NONLOCAL);
981
982             if (nbv->grp[eintNonlocal].kernel_type == nbk8x8x8_CUDA)
983             {
984                 /* initialize non-local pair-list on the GPU */
985                 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
986                                          nbv->grp[eintNonlocal].nbl_lists.nbl[0],
987                                          eintNonlocal);
988             }
989             wallcycle_stop(wcycle,ewcNS);
990         } 
991         else
992         {
993             wallcycle_start(wcycle,ewcMOVEX);
994             dd_move_x(cr->dd,box,x);
995
996             /* When we don't need the total dipole we sum it in global_stat */
997             if (bStateChanged && NEED_MUTOT(*inputrec))
998             {
999                 gmx_sumd(2*DIM,mu,cr);
1000             }
1001             wallcycle_stop(wcycle,ewcMOVEX);
1002
1003             wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1004             wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
1005             nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs,eatNonlocal,FALSE,x,
1006                                             nbv->grp[eintNonlocal].nbat);
1007             wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
1008             cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1009         }
1010
1011         if (bUseGPU && !bDiffKernels)
1012         { 
1013             wallcycle_start(wcycle,ewcLAUNCH_GPU_NB);
1014             /* launch non-local nonbonded F on GPU */
1015             do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFNo,
1016                          nrnb, wcycle);
1017             cycles_force += wallcycle_stop(wcycle,ewcLAUNCH_GPU_NB);
1018         }
1019     }
1020
1021     if (bUseGPU)
1022     {
1023         /* launch D2H copy-back F */
1024         wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1025         if (DOMAINDECOMP(cr) && !bDiffKernels)
1026         {
1027             nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintNonlocal].nbat,
1028                                       flags, eatNonlocal);
1029         }
1030         nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintLocal].nbat,
1031                                   flags, eatLocal);
1032         cycles_force += wallcycle_stop(wcycle,ewcLAUNCH_GPU_NB);
1033     }
1034
1035     if (bStateChanged && NEED_MUTOT(*inputrec))
1036     {
1037         if (PAR(cr))
1038         {
1039             gmx_sumd(2*DIM,mu,cr);
1040         } 
1041
1042         for(i=0; i<2; i++)
1043         {
1044             for(j=0;j<DIM;j++)
1045             {
1046                 fr->mu_tot[i][j] = mu[i*DIM + j];
1047             }
1048         }
1049     }
1050     if (fr->efep == efepNO)
1051     {
1052         copy_rvec(fr->mu_tot[0],mu_tot);
1053     }
1054     else
1055     {
1056         for(j=0; j<DIM; j++)
1057         {
1058             mu_tot[j] =
1059                 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] +
1060                 lambda[efptCOUL]*fr->mu_tot[1][j];
1061         }
1062     }
1063
1064     /* Reset energies */
1065     reset_enerdata(&(inputrec->opts),fr,bNS,enerd,MASTER(cr));
1066     clear_rvecs(SHIFTS,fr->fshift);
1067
1068     if (DOMAINDECOMP(cr))
1069     {
1070         if (!(cr->duty & DUTY_PME))
1071         {
1072             wallcycle_start(wcycle,ewcPPDURINGPME);
1073             dd_force_flop_start(cr->dd,nrnb);
1074         }
1075     }
1076     
1077     /* Start the force cycle counter.
1078      * This counter is stopped in do_forcelow_level.
1079      * No parallel communication should occur while this counter is running,
1080      * since that will interfere with the dynamic load balancing.
1081      */
1082     wallcycle_start(wcycle,ewcFORCE);
1083     if (bDoForces)
1084     {
1085         /* Reset forces for which the virial is calculated separately:
1086          * PME/Ewald forces if necessary */
1087         if (fr->bF_NoVirSum) 
1088         {
1089             if (flags & GMX_FORCE_VIRIAL)
1090             {
1091                 fr->f_novirsum = fr->f_novirsum_alloc;
1092                 GMX_BARRIER(cr->mpi_comm_mygroup);
1093                 if (fr->bDomDec)
1094                 {
1095                     clear_rvecs(fr->f_novirsum_n,fr->f_novirsum);
1096                 }
1097                 else
1098                 {
1099                     clear_rvecs(homenr,fr->f_novirsum+start);
1100                 }
1101                 GMX_BARRIER(cr->mpi_comm_mygroup);
1102             }
1103             else
1104             {
1105                 /* We are not calculating the pressure so we do not need
1106                  * a separate array for forces that do not contribute
1107                  * to the pressure.
1108                  */
1109                 fr->f_novirsum = f;
1110             }
1111         }
1112
1113         /* Clear the short- and long-range forces */
1114         clear_rvecs(fr->natoms_force_constr,f);
1115         if(bSepLRF && do_per_step(step,inputrec->nstcalclr))
1116         {
1117             clear_rvecs(fr->natoms_force_constr,fr->f_twin);
1118         }
1119         
1120         clear_rvec(fr->vir_diag_posres);
1121
1122         GMX_BARRIER(cr->mpi_comm_mygroup);
1123     }
1124     if (inputrec->ePull == epullCONSTRAINT)
1125     {
1126         clear_pull_forces(inputrec->pull);
1127     }
1128
1129     /* update QMMMrec, if necessary */
1130     if(fr->bQMMM)
1131     {
1132         update_QMMMrec(cr,fr,x,mdatoms,box,top);
1133     }
1134
1135     if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1136     {
1137         posres_wrapper(fplog,flags,bSepDVDL,inputrec,nrnb,top,box,x,
1138                        f,enerd,lambda,fr);
1139     }
1140
1141     /* Compute the bonded and non-bonded energies and optionally forces */    
1142     /* if we use the GPU turn off the nonbonded */
1143     do_force_lowlevel(fplog,step,fr,inputrec,&(top->idef),
1144                       cr,nrnb,wcycle,mdatoms,&(inputrec->opts),
1145                       x,hist,f, bSepLRF ? fr->f_twin : f,enerd,fcd,mtop,top,fr->born,
1146                       &(top->atomtypes),bBornRadii,box,
1147                       inputrec->fepvals,lambda,graph,&(top->excls),fr->mu_tot,
1148                       ((nb_kernel_type == nbk8x8x8_CUDA || nb_kernel_type == nbk8x8x8_PlainC) 
1149                         ? flags&~GMX_FORCE_NONBONDED : flags),
1150                       &cycles_pme);
1151
1152     if(bSepLRF)
1153     {
1154         if (do_per_step(step,inputrec->nstcalclr))
1155         {
1156             /* Add the long range forces to the short range forces */
1157             for(i=0; i<fr->natoms_force_constr; i++)
1158             {
1159                 rvec_add(fr->f_twin[i],f[i],f[i]);
1160             }
1161         }
1162     }
1163     
1164     if (!bUseOrEmulGPU)
1165     {
1166         /* Maybe we should move this into do_force_lowlevel */
1167         do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFYes,
1168                      nrnb, wcycle);
1169     }
1170         
1171
1172     if (!bUseOrEmulGPU || bDiffKernels)
1173     {
1174         int aloc;
1175
1176         if (DOMAINDECOMP(cr))
1177         {
1178             do_nb_verlet(fr, ic, enerd, flags, eintNonlocal,
1179                          bDiffKernels ? enbvClearFYes : enbvClearFNo,
1180                          nrnb, wcycle);
1181         }
1182
1183         if (!bUseOrEmulGPU)
1184         {
1185             aloc = eintLocal;
1186         }
1187         else
1188         {
1189             aloc = eintNonlocal;
1190         }
1191
1192         /* Add all the non-bonded force to the normal force array.
1193          * This can be split into a local a non-local part when overlapping
1194          * communication with calculation with domain decomposition.
1195          */
1196         cycles_force += wallcycle_stop(wcycle,ewcFORCE);
1197         wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1198         wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1199         nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs,eatAll,nbv->grp[aloc].nbat,f);
1200         wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1201         cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1202         wallcycle_start_nocount(wcycle,ewcFORCE);
1203
1204         /* if there are multiple fshift output buffers reduce them */
1205         if ((flags & GMX_FORCE_VIRIAL) &&
1206             nbv->grp[aloc].nbl_lists.nnbl > 1)
1207         {
1208             nbnxn_atomdata_add_nbat_fshift_to_fshift(nbv->grp[aloc].nbat,
1209                                                       fr->fshift);
1210         }
1211     }
1212     
1213     cycles_force += wallcycle_stop(wcycle,ewcFORCE);
1214     GMX_BARRIER(cr->mpi_comm_mygroup);
1215     
1216     if (ed)
1217     {
1218         do_flood(fplog,cr,x,f,ed,box,step,bNS);
1219     }
1220
1221     if (bUseOrEmulGPU && !bDiffKernels)
1222     {
1223         /* wait for non-local forces (or calculate in emulation mode) */
1224         if (DOMAINDECOMP(cr))
1225         {
1226             if (bUseGPU)
1227             {
1228                 wallcycle_start(wcycle,ewcWAIT_GPU_NB_NL);
1229                 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1230                                     nbv->grp[eintNonlocal].nbat,
1231                                     flags, eatNonlocal,
1232                                     enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1233                                     fr->fshift);
1234                 cycles_force += wallcycle_stop(wcycle,ewcWAIT_GPU_NB_NL);
1235             }
1236             else
1237             {
1238                 wallcycle_start_nocount(wcycle,ewcFORCE);
1239                 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFYes,
1240                              nrnb, wcycle);
1241                 cycles_force += wallcycle_stop(wcycle,ewcFORCE);
1242             }            
1243             wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1244             wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1245             /* skip the reduction if there was no non-local work to do */
1246             if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1247             {
1248                 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs,eatNonlocal,
1249                                                nbv->grp[eintNonlocal].nbat,f);
1250             }
1251             wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1252             cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1253         }
1254     }
1255
1256     if (bDoForces)
1257     {
1258         /* Communicate the forces */
1259         if (PAR(cr))
1260         {
1261             wallcycle_start(wcycle,ewcMOVEF);
1262             if (DOMAINDECOMP(cr))
1263             {
1264                 dd_move_f(cr->dd,f,fr->fshift);
1265                 /* Do we need to communicate the separate force array
1266                  * for terms that do not contribute to the single sum virial?
1267                  * Position restraints and electric fields do not introduce
1268                  * inter-cg forces, only full electrostatics methods do.
1269                  * When we do not calculate the virial, fr->f_novirsum = f,
1270                  * so we have already communicated these forces.
1271                  */
1272                 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1273                     (flags & GMX_FORCE_VIRIAL))
1274                 {
1275                     dd_move_f(cr->dd,fr->f_novirsum,NULL);
1276                 }
1277                 if (bSepLRF)
1278                 {
1279                     /* We should not update the shift forces here,
1280                      * since f_twin is already included in f.
1281                      */
1282                     dd_move_f(cr->dd,fr->f_twin,NULL);
1283                 }
1284             }
1285             wallcycle_stop(wcycle,ewcMOVEF);
1286         }
1287     }
1288  
1289     if (bUseOrEmulGPU)
1290     {
1291         /* wait for local forces (or calculate in emulation mode) */
1292         if (bUseGPU)
1293         {
1294             wallcycle_start(wcycle,ewcWAIT_GPU_NB_L);
1295             nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1296                                 nbv->grp[eintLocal].nbat,
1297                                 flags, eatLocal,
1298                                 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1299                                 fr->fshift);
1300             wallcycle_stop(wcycle,ewcWAIT_GPU_NB_L);
1301
1302             /* now clear the GPU outputs while we finish the step on the CPU */
1303             nbnxn_cuda_clear_outputs(nbv->cu_nbv, flags);
1304         }
1305         else
1306         {            
1307             wallcycle_start_nocount(wcycle,ewcFORCE);
1308             do_nb_verlet(fr, ic, enerd, flags, eintLocal,
1309                          DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
1310                          nrnb, wcycle);
1311             wallcycle_stop(wcycle,ewcFORCE);
1312         }
1313         wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1314         wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1315         if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1316         {
1317             /* skip the reduction if there was no non-local work to do */
1318             nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs,eatLocal,
1319                                            nbv->grp[eintLocal].nbat,f);
1320         }
1321         wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1322         wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1323     }
1324     
1325     if (DOMAINDECOMP(cr))
1326     {
1327         dd_force_flop_stop(cr->dd,nrnb);
1328         if (wcycle)
1329         {
1330             dd_cycles_add(cr->dd,cycles_force-cycles_pme,ddCyclF);
1331         }
1332     }
1333
1334     if (bDoForces)
1335     {
1336         if (IR_ELEC_FIELD(*inputrec))
1337         {
1338             /* Compute forces due to electric field */
1339             calc_f_el(MASTER(cr) ? field : NULL,
1340                       start,homenr,mdatoms->chargeA,x,fr->f_novirsum,
1341                       inputrec->ex,inputrec->et,t);
1342         }
1343
1344         /* If we have NoVirSum forces, but we do not calculate the virial,
1345          * we sum fr->f_novirum=f later.
1346          */
1347         if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1348         {
1349             wallcycle_start(wcycle,ewcVSITESPREAD);
1350             spread_vsite_f(fplog,vsite,x,f,fr->fshift,FALSE,NULL,nrnb,
1351                            &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
1352             wallcycle_stop(wcycle,ewcVSITESPREAD);
1353
1354             if (bSepLRF)
1355             {
1356                 wallcycle_start(wcycle,ewcVSITESPREAD);
1357                 spread_vsite_f(fplog,vsite,x,fr->f_twin,NULL,FALSE,NULL,
1358                                nrnb,
1359                                &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
1360                 wallcycle_stop(wcycle,ewcVSITESPREAD);
1361             }
1362         }
1363
1364         if (flags & GMX_FORCE_VIRIAL)
1365         {
1366             /* Calculation of the virial must be done after vsites! */
1367             calc_virial(fplog,mdatoms->start,mdatoms->homenr,x,f,
1368                         vir_force,graph,box,nrnb,fr,inputrec->ePBC);
1369         }
1370     }
1371
1372     if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1373     {
1374         pull_potential_wrapper(fplog,bSepDVDL,cr,inputrec,box,x,
1375                                f,vir_force,mdatoms,enerd,lambda,t);
1376     }
1377
1378     if (PAR(cr) && !(cr->duty & DUTY_PME))
1379     {
1380         /* In case of node-splitting, the PP nodes receive the long-range 
1381          * forces, virial and energy from the PME nodes here.
1382          */    
1383         pme_receive_force_ener(fplog,bSepDVDL,cr,wcycle,enerd,fr);
1384     }
1385
1386     if (bDoForces)
1387     {
1388         post_process_forces(fplog,cr,step,nrnb,wcycle,
1389                             top,box,x,f,vir_force,mdatoms,graph,fr,vsite,
1390                             flags);
1391     }
1392     
1393     /* Sum the potential energy terms from group contributions */
1394     sum_epot(&(inputrec->opts),enerd);
1395 }
1396
1397 void do_force_cutsGROUP(FILE *fplog,t_commrec *cr,
1398               t_inputrec *inputrec,
1399               gmx_large_int_t step,t_nrnb *nrnb,gmx_wallcycle_t wcycle,
1400               gmx_localtop_t *top,
1401               gmx_mtop_t *mtop,
1402               gmx_groups_t *groups,
1403               matrix box,rvec x[],history_t *hist,
1404               rvec f[],
1405               tensor vir_force,
1406               t_mdatoms *mdatoms,
1407               gmx_enerdata_t *enerd,t_fcdata *fcd,
1408               real *lambda,t_graph *graph,
1409               t_forcerec *fr,gmx_vsite_t *vsite,rvec mu_tot,
1410               double t,FILE *field,gmx_edsam_t ed,
1411               gmx_bool bBornRadii,
1412               int flags)
1413 {
1414     int    cg0,cg1,i,j;
1415     int    start,homenr;
1416     double mu[2*DIM];
1417     gmx_bool   bSepDVDL,bStateChanged,bNS,bFillGrid,bCalcCGCM,bBS;
1418     gmx_bool   bDoLongRangeNS,bDoForces,bDoPotential,bSepLRF;
1419     gmx_bool   bDoAdressWF;
1420     matrix boxs;
1421     rvec   vzero,box_diag;
1422     real   e,v,dvdlambda[efptNR];
1423     t_pbc  pbc;
1424     float  cycles_pme,cycles_force;
1425
1426     start  = mdatoms->start;
1427     homenr = mdatoms->homenr;
1428
1429     bSepDVDL = (fr->bSepDVDL && do_per_step(step,inputrec->nstlog));
1430
1431     clear_mat(vir_force);
1432
1433     if (PARTDECOMP(cr))
1434     {
1435         pd_cg_range(cr,&cg0,&cg1);
1436     }
1437     else
1438     {
1439         cg0 = 0;
1440         if (DOMAINDECOMP(cr))
1441         {
1442             cg1 = cr->dd->ncg_tot;
1443         }
1444         else
1445         {
1446             cg1 = top->cgs.nr;
1447         }
1448         if (fr->n_tpi > 0)
1449         {
1450             cg1--;
1451         }
1452     }
1453
1454     bStateChanged  = (flags & GMX_FORCE_STATECHANGED);
1455     bNS            = (flags & GMX_FORCE_NS) && (fr->bAllvsAll==FALSE);
1456     /* Should we update the long-range neighborlists at this step? */
1457     bDoLongRangeNS = fr->bTwinRange && bNS;
1458     /* Should we perform the long-range nonbonded evaluation inside the neighborsearching? */
1459     bFillGrid      = (bNS && bStateChanged);
1460     bCalcCGCM      = (bFillGrid && !DOMAINDECOMP(cr));
1461     bDoForces      = (flags & GMX_FORCE_FORCES);
1462     bDoPotential   = (flags & GMX_FORCE_ENERGY);
1463     bSepLRF        = ((inputrec->nstcalclr>1) && bDoForces &&
1464                       (flags & GMX_FORCE_SEPLRF) && (flags & GMX_FORCE_DO_LR));
1465
1466     /* should probably move this to the forcerec since it doesn't change */
1467     bDoAdressWF   = ((fr->adress_type!=eAdressOff));
1468
1469     if (bStateChanged)
1470     {
1471         update_forcerec(fplog,fr,box);
1472
1473         if (NEED_MUTOT(*inputrec))
1474         {
1475             /* Calculate total (local) dipole moment in a temporary common array.
1476              * This makes it possible to sum them over nodes faster.
1477              */
1478             calc_mu(start,homenr,
1479                     x,mdatoms->chargeA,mdatoms->chargeB,mdatoms->nChargePerturbed,
1480                     mu,mu+DIM);
1481         }
1482     }
1483
1484     if (fr->ePBC != epbcNONE) { 
1485         /* Compute shift vectors every step,
1486          * because of pressure coupling or box deformation!
1487          */
1488         if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
1489             calc_shifts(box,fr->shift_vec);
1490
1491         if (bCalcCGCM) { 
1492             put_charge_groups_in_box(fplog,cg0,cg1,fr->ePBC,box,
1493                     &(top->cgs),x,fr->cg_cm);
1494             inc_nrnb(nrnb,eNR_CGCM,homenr);
1495             inc_nrnb(nrnb,eNR_RESETX,cg1-cg0);
1496         } 
1497         else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph) {
1498             unshift_self(graph,box,x);
1499         }
1500     } 
1501     else if (bCalcCGCM) {
1502         calc_cgcm(fplog,cg0,cg1,&(top->cgs),x,fr->cg_cm);
1503         inc_nrnb(nrnb,eNR_CGCM,homenr);
1504     }
1505
1506     if (bCalcCGCM) {
1507         if (PAR(cr)) {
1508             move_cgcm(fplog,cr,fr->cg_cm);
1509         }
1510         if (gmx_debug_at)
1511             pr_rvecs(debug,0,"cgcm",fr->cg_cm,top->cgs.nr);
1512     }
1513
1514 #ifdef GMX_MPI
1515     if (!(cr->duty & DUTY_PME)) {
1516         /* Send particle coordinates to the pme nodes.
1517          * Since this is only implemented for domain decomposition
1518          * and domain decomposition does not use the graph,
1519          * we do not need to worry about shifting.
1520          */    
1521
1522         wallcycle_start(wcycle,ewcPP_PMESENDX);
1523         GMX_MPE_LOG(ev_send_coordinates_start);
1524
1525         bBS = (inputrec->nwall == 2);
1526         if (bBS) {
1527             copy_mat(box,boxs);
1528             svmul(inputrec->wall_ewald_zfac,boxs[ZZ],boxs[ZZ]);
1529         }
1530
1531         gmx_pme_send_x(cr,bBS ? boxs : box,x,
1532                        mdatoms->nChargePerturbed,lambda[efptCOUL],
1533                        (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),step);
1534
1535         GMX_MPE_LOG(ev_send_coordinates_finish);
1536         wallcycle_stop(wcycle,ewcPP_PMESENDX);
1537     }
1538 #endif /* GMX_MPI */
1539
1540     /* Communicate coordinates and sum dipole if necessary */
1541     if (PAR(cr))
1542     {
1543         wallcycle_start(wcycle,ewcMOVEX);
1544         if (DOMAINDECOMP(cr))
1545         {
1546             dd_move_x(cr->dd,box,x);
1547         }
1548         else
1549         {
1550             move_x(fplog,cr,GMX_LEFT,GMX_RIGHT,x,nrnb);
1551         }
1552         wallcycle_stop(wcycle,ewcMOVEX);
1553     }
1554
1555     /* update adress weight beforehand */
1556     if(bStateChanged && bDoAdressWF)
1557     {
1558         /* need pbc for adress weight calculation with pbc_dx */
1559         set_pbc(&pbc,inputrec->ePBC,box);
1560         if(fr->adress_site == eAdressSITEcog)
1561         {
1562             update_adress_weights_cog(top->idef.iparams,top->idef.il,x,fr,mdatoms,
1563                                       inputrec->ePBC==epbcNONE ? NULL : &pbc);
1564         }
1565         else if (fr->adress_site == eAdressSITEcom)
1566         {
1567             update_adress_weights_com(fplog,cg0,cg1,&(top->cgs),x,fr,mdatoms,
1568                                       inputrec->ePBC==epbcNONE ? NULL : &pbc);
1569         }
1570         else if (fr->adress_site == eAdressSITEatomatom){
1571             update_adress_weights_atom_per_atom(cg0,cg1,&(top->cgs),x,fr,mdatoms,
1572                                                 inputrec->ePBC==epbcNONE ? NULL : &pbc);
1573         }
1574         else
1575         {
1576             update_adress_weights_atom(cg0,cg1,&(top->cgs),x,fr,mdatoms,
1577                                        inputrec->ePBC==epbcNONE ? NULL : &pbc);
1578         }
1579     }
1580
1581     if (NEED_MUTOT(*inputrec))
1582     {
1583
1584         if (bStateChanged)
1585         {
1586             if (PAR(cr))
1587             {
1588                 gmx_sumd(2*DIM,mu,cr);
1589             }
1590             for(i=0; i<2; i++)
1591             {
1592                 for(j=0;j<DIM;j++)
1593                 {
1594                     fr->mu_tot[i][j] = mu[i*DIM + j];
1595                 }
1596             }
1597         }
1598         if (fr->efep == efepNO)
1599         {
1600             copy_rvec(fr->mu_tot[0],mu_tot);
1601         }
1602         else
1603         {
1604             for(j=0; j<DIM; j++)
1605             {
1606                 mu_tot[j] =
1607                     (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] + lambda[efptCOUL]*fr->mu_tot[1][j];
1608             }
1609         }
1610     }
1611
1612     /* Reset energies */
1613     reset_enerdata(&(inputrec->opts),fr,bNS,enerd,MASTER(cr));
1614     clear_rvecs(SHIFTS,fr->fshift);
1615
1616     if (bNS)
1617     {
1618         wallcycle_start(wcycle,ewcNS);
1619
1620         if (graph && bStateChanged)
1621         {
1622             /* Calculate intramolecular shift vectors to make molecules whole */
1623             mk_mshift(fplog,graph,fr->ePBC,box,x);
1624         }
1625
1626         /* Do the actual neighbour searching and if twin range electrostatics
1627          * also do the calculation of long range forces and energies.
1628          */
1629         for (i=0;i<efptNR;i++) {dvdlambda[i] = 0;}
1630         ns(fplog,fr,x,box,
1631            groups,&(inputrec->opts),top,mdatoms,
1632            cr,nrnb,lambda,dvdlambda,&enerd->grpp,bFillGrid,
1633            bDoLongRangeNS);
1634         if (bSepDVDL)
1635         {
1636             fprintf(fplog,sepdvdlformat,"LR non-bonded",0.0,dvdlambda);
1637         }
1638         enerd->dvdl_lin[efptVDW] += dvdlambda[efptVDW];
1639         enerd->dvdl_lin[efptCOUL] += dvdlambda[efptCOUL];
1640
1641         wallcycle_stop(wcycle,ewcNS);
1642     }
1643
1644     if (inputrec->implicit_solvent && bNS)
1645     {
1646         make_gb_nblist(cr,inputrec->gb_algorithm,inputrec->rlist,
1647                        x,box,fr,&top->idef,graph,fr->born);
1648     }
1649
1650     if (DOMAINDECOMP(cr))
1651     {
1652         if (!(cr->duty & DUTY_PME))
1653         {
1654             wallcycle_start(wcycle,ewcPPDURINGPME);
1655             dd_force_flop_start(cr->dd,nrnb);
1656         }
1657     }
1658
1659     if (inputrec->bRot)
1660     {
1661         /* Enforced rotation has its own cycle counter that starts after the collective
1662          * coordinates have been communicated. It is added to ddCyclF to allow
1663          * for proper load-balancing */
1664         wallcycle_start(wcycle,ewcROT);
1665         do_rotation(cr,inputrec,box,x,t,step,wcycle,bNS);
1666         wallcycle_stop(wcycle,ewcROT);
1667     }
1668
1669     /* Start the force cycle counter.
1670      * This counter is stopped in do_forcelow_level.
1671      * No parallel communication should occur while this counter is running,
1672      * since that will interfere with the dynamic load balancing.
1673      */
1674     wallcycle_start(wcycle,ewcFORCE);
1675     
1676     if (bDoForces)
1677     {
1678         /* Reset forces for which the virial is calculated separately:
1679          * PME/Ewald forces if necessary */
1680         if (fr->bF_NoVirSum)
1681         {
1682             if (flags & GMX_FORCE_VIRIAL)
1683             {
1684                 fr->f_novirsum = fr->f_novirsum_alloc;
1685                 GMX_BARRIER(cr->mpi_comm_mygroup);
1686                 if (fr->bDomDec)
1687                 {
1688                     clear_rvecs(fr->f_novirsum_n,fr->f_novirsum);
1689                 }
1690                 else
1691                 {
1692                     clear_rvecs(homenr,fr->f_novirsum+start);
1693                 }
1694                 GMX_BARRIER(cr->mpi_comm_mygroup);
1695             }
1696             else
1697             {
1698                 /* We are not calculating the pressure so we do not need
1699                  * a separate array for forces that do not contribute
1700                  * to the pressure.
1701                  */
1702                 fr->f_novirsum = f;
1703             }
1704         }
1705
1706         /* Clear the short- and long-range forces */
1707         clear_rvecs(fr->natoms_force_constr,f);
1708         if(bSepLRF && do_per_step(step,inputrec->nstcalclr))
1709         {
1710             clear_rvecs(fr->natoms_force_constr,fr->f_twin);
1711         }
1712         
1713         clear_rvec(fr->vir_diag_posres);
1714
1715         GMX_BARRIER(cr->mpi_comm_mygroup);
1716     }
1717     if (inputrec->ePull == epullCONSTRAINT)
1718     {
1719         clear_pull_forces(inputrec->pull);
1720     }
1721
1722     /* update QMMMrec, if necessary */
1723     if(fr->bQMMM)
1724     {
1725         update_QMMMrec(cr,fr,x,mdatoms,box,top);
1726     }
1727
1728     if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1729     {
1730         posres_wrapper(fplog,flags,bSepDVDL,inputrec,nrnb,top,box,x,
1731                        f,enerd,lambda,fr);
1732     }
1733
1734     /* Compute the bonded and non-bonded energies and optionally forces */
1735     do_force_lowlevel(fplog,step,fr,inputrec,&(top->idef),
1736                       cr,nrnb,wcycle,mdatoms,&(inputrec->opts),
1737                       x,hist,f, bSepLRF ? fr->f_twin : f,enerd,fcd,mtop,top,fr->born,
1738                       &(top->atomtypes),bBornRadii,box,
1739                       inputrec->fepvals,lambda,
1740                       graph,&(top->excls),fr->mu_tot,
1741                       flags,
1742                       &cycles_pme);
1743
1744     if(bSepLRF)
1745     {
1746         if (do_per_step(step,inputrec->nstcalclr))
1747         {
1748             /* Add the long range forces to the short range forces */
1749             for(i=0; i<fr->natoms_force_constr; i++)
1750             {
1751                 rvec_add(fr->f_twin[i],f[i],f[i]);
1752             }
1753         }
1754     }
1755     
1756     cycles_force = wallcycle_stop(wcycle,ewcFORCE);
1757     GMX_BARRIER(cr->mpi_comm_mygroup);
1758
1759     if (ed)
1760     {
1761         do_flood(fplog,cr,x,f,ed,box,step,bNS);
1762     }
1763
1764     if (DOMAINDECOMP(cr))
1765     {
1766         dd_force_flop_stop(cr->dd,nrnb);
1767         if (wcycle)
1768         {
1769             dd_cycles_add(cr->dd,cycles_force-cycles_pme,ddCyclF);
1770         }
1771     }
1772
1773     if (bDoForces)
1774     {
1775         if (IR_ELEC_FIELD(*inputrec))
1776         {
1777             /* Compute forces due to electric field */
1778             calc_f_el(MASTER(cr) ? field : NULL,
1779                       start,homenr,mdatoms->chargeA,x,fr->f_novirsum,
1780                       inputrec->ex,inputrec->et,t);
1781         }
1782
1783         if (bDoAdressWF && fr->adress_icor == eAdressICThermoForce)
1784         {
1785             /* Compute thermodynamic force in hybrid AdResS region */
1786             adress_thermo_force(start,homenr,&(top->cgs),x,fr->f_novirsum,fr,mdatoms,
1787                                 inputrec->ePBC==epbcNONE ? NULL : &pbc);
1788         }
1789
1790         /* Communicate the forces */
1791         if (PAR(cr))
1792         {
1793             wallcycle_start(wcycle,ewcMOVEF);
1794             if (DOMAINDECOMP(cr))
1795             {
1796                 dd_move_f(cr->dd,f,fr->fshift);
1797                 /* Do we need to communicate the separate force array
1798                  * for terms that do not contribute to the single sum virial?
1799                  * Position restraints and electric fields do not introduce
1800                  * inter-cg forces, only full electrostatics methods do.
1801                  * When we do not calculate the virial, fr->f_novirsum = f,
1802                  * so we have already communicated these forces.
1803                  */
1804                 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1805                     (flags & GMX_FORCE_VIRIAL))
1806                 {
1807                     dd_move_f(cr->dd,fr->f_novirsum,NULL);
1808                 }
1809                 if (bSepLRF)
1810                 {
1811                     /* We should not update the shift forces here,
1812                      * since f_twin is already included in f.
1813                      */
1814                     dd_move_f(cr->dd,fr->f_twin,NULL);
1815                 }
1816             }
1817             else
1818             {
1819                 pd_move_f(cr,f,nrnb);
1820                 if (bSepLRF)
1821                 {
1822                     pd_move_f(cr,fr->f_twin,nrnb);
1823                 }
1824             }
1825             wallcycle_stop(wcycle,ewcMOVEF);
1826         }
1827
1828         /* If we have NoVirSum forces, but we do not calculate the virial,
1829          * we sum fr->f_novirum=f later.
1830          */
1831         if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1832         {
1833             wallcycle_start(wcycle,ewcVSITESPREAD);
1834             spread_vsite_f(fplog,vsite,x,f,fr->fshift,FALSE,NULL,nrnb,
1835                            &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
1836             wallcycle_stop(wcycle,ewcVSITESPREAD);
1837
1838             if (bSepLRF)
1839             {
1840                 wallcycle_start(wcycle,ewcVSITESPREAD);
1841                 spread_vsite_f(fplog,vsite,x,fr->f_twin,NULL,FALSE,NULL,
1842                                nrnb,
1843                                &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
1844                 wallcycle_stop(wcycle,ewcVSITESPREAD);
1845             }
1846         }
1847
1848         if (flags & GMX_FORCE_VIRIAL)
1849         {
1850             /* Calculation of the virial must be done after vsites! */
1851             calc_virial(fplog,mdatoms->start,mdatoms->homenr,x,f,
1852                         vir_force,graph,box,nrnb,fr,inputrec->ePBC);
1853         }
1854     }
1855
1856     if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1857     {
1858         pull_potential_wrapper(fplog,bSepDVDL,cr,inputrec,box,x,
1859                                f,vir_force,mdatoms,enerd,lambda,t);
1860     }
1861
1862     /* Add the forces from enforced rotation potentials (if any) */
1863     if (inputrec->bRot)
1864     {
1865         wallcycle_start(wcycle,ewcROTadd);
1866         enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr,step,t);
1867         wallcycle_stop(wcycle,ewcROTadd);
1868     }
1869
1870     if (PAR(cr) && !(cr->duty & DUTY_PME))
1871     {
1872         /* In case of node-splitting, the PP nodes receive the long-range 
1873          * forces, virial and energy from the PME nodes here.
1874          */
1875         pme_receive_force_ener(fplog,bSepDVDL,cr,wcycle,enerd,fr);
1876     }
1877
1878     if (bDoForces)
1879     {
1880         post_process_forces(fplog,cr,step,nrnb,wcycle,
1881                             top,box,x,f,vir_force,mdatoms,graph,fr,vsite,
1882                             flags);
1883     }
1884
1885     /* Sum the potential energy terms from group contributions */
1886     sum_epot(&(inputrec->opts),enerd);
1887 }
1888
1889 void do_force(FILE *fplog,t_commrec *cr,
1890               t_inputrec *inputrec,
1891               gmx_large_int_t step,t_nrnb *nrnb,gmx_wallcycle_t wcycle,
1892               gmx_localtop_t *top,
1893               gmx_mtop_t *mtop,
1894               gmx_groups_t *groups,
1895               matrix box,rvec x[],history_t *hist,
1896               rvec f[],
1897               tensor vir_force,
1898               t_mdatoms *mdatoms,
1899               gmx_enerdata_t *enerd,t_fcdata *fcd,
1900               real *lambda,t_graph *graph,
1901               t_forcerec *fr,
1902               gmx_vsite_t *vsite,rvec mu_tot,
1903               double t,FILE *field,gmx_edsam_t ed,
1904               gmx_bool bBornRadii,
1905               int flags)
1906 {
1907     /* modify force flag if not doing nonbonded */
1908     if (!fr->bNonbonded)
1909     {
1910         flags &= ~GMX_FORCE_NONBONDED;
1911     }
1912
1913     switch (inputrec->cutoff_scheme)
1914     {
1915         case ecutsVERLET:
1916             do_force_cutsVERLET(fplog, cr, inputrec,
1917                                 step, nrnb, wcycle,
1918                                 top, mtop,
1919                                 groups,
1920                                 box, x, hist,
1921                                 f, vir_force,
1922                                 mdatoms,
1923                                 enerd, fcd,
1924                                 lambda, graph,
1925                                 fr, fr->ic, 
1926                                 vsite, mu_tot,
1927                                 t, field, ed,
1928                                 bBornRadii,
1929                                 flags);
1930             break;
1931         case ecutsGROUP:
1932              do_force_cutsGROUP(fplog, cr, inputrec,
1933                                 step, nrnb, wcycle,
1934                                 top, mtop,
1935                                 groups,
1936                                 box, x, hist,
1937                                 f, vir_force,
1938                                 mdatoms,
1939                                 enerd, fcd,
1940                                 lambda, graph,
1941                                 fr, vsite, mu_tot,
1942                                 t, field, ed,
1943                                 bBornRadii,
1944                                 flags);
1945             break;
1946         default:
1947             gmx_incons("Invalid cut-off scheme passed!");
1948     }
1949 }
1950
1951
1952 void do_constrain_first(FILE *fplog,gmx_constr_t constr,
1953                         t_inputrec *ir,t_mdatoms *md,
1954                         t_state *state,rvec *f,
1955                         t_graph *graph,t_commrec *cr,t_nrnb *nrnb,
1956                         t_forcerec *fr, gmx_localtop_t *top, tensor shake_vir)
1957 {
1958     int    i,m,start,end;
1959     gmx_large_int_t step;
1960     real   dt=ir->delta_t;
1961     real   dvdl_dum;
1962     rvec   *savex;
1963
1964     snew(savex,state->natoms);
1965
1966     start = md->start;
1967     end   = md->homenr + start;
1968
1969     if (debug)
1970         fprintf(debug,"vcm: start=%d, homenr=%d, end=%d\n",
1971                 start,md->homenr,end);
1972     /* Do a first constrain to reset particles... */
1973     step = ir->init_step;
1974     if (fplog)
1975     {
1976         char buf[STEPSTRSIZE];
1977         fprintf(fplog,"\nConstraining the starting coordinates (step %s)\n",
1978                 gmx_step_str(step,buf));
1979     }
1980     dvdl_dum = 0;
1981
1982     /* constrain the current position */
1983     constrain(NULL,TRUE,FALSE,constr,&(top->idef),
1984               ir,NULL,cr,step,0,md,
1985               state->x,state->x,NULL,
1986               fr->bMolPBC,state->box,
1987               state->lambda[efptBONDED],&dvdl_dum,
1988               NULL,NULL,nrnb,econqCoord,
1989               ir->epc==epcMTTK,state->veta,state->veta);
1990     if (EI_VV(ir->eI))
1991     {
1992         /* constrain the inital velocity, and save it */
1993         /* also may be useful if we need the ekin from the halfstep for velocity verlet */
1994         /* might not yet treat veta correctly */
1995         constrain(NULL,TRUE,FALSE,constr,&(top->idef),
1996                   ir,NULL,cr,step,0,md,
1997                   state->x,state->v,state->v,
1998                   fr->bMolPBC,state->box,
1999                   state->lambda[efptBONDED],&dvdl_dum,
2000                   NULL,NULL,nrnb,econqVeloc,
2001                   ir->epc==epcMTTK,state->veta,state->veta);
2002     }
2003     /* constrain the inital velocities at t-dt/2 */
2004     if (EI_STATE_VELOCITY(ir->eI) && ir->eI!=eiVV)
2005     {
2006         for(i=start; (i<end); i++)
2007         {
2008             for(m=0; (m<DIM); m++)
2009             {
2010                 /* Reverse the velocity */
2011                 state->v[i][m] = -state->v[i][m];
2012                 /* Store the position at t-dt in buf */
2013                 savex[i][m] = state->x[i][m] + dt*state->v[i][m];
2014             }
2015         }
2016     /* Shake the positions at t=-dt with the positions at t=0
2017      * as reference coordinates.
2018          */
2019         if (fplog)
2020         {
2021             char buf[STEPSTRSIZE];
2022             fprintf(fplog,"\nConstraining the coordinates at t0-dt (step %s)\n",
2023                     gmx_step_str(step,buf));
2024         }
2025         dvdl_dum = 0;
2026         constrain(NULL,TRUE,FALSE,constr,&(top->idef),
2027                   ir,NULL,cr,step,-1,md,
2028                   state->x,savex,NULL,
2029                   fr->bMolPBC,state->box,
2030                   state->lambda[efptBONDED],&dvdl_dum,
2031                   state->v,NULL,nrnb,econqCoord,
2032                   ir->epc==epcMTTK,state->veta,state->veta);
2033         
2034         for(i=start; i<end; i++) {
2035             for(m=0; m<DIM; m++) {
2036                 /* Re-reverse the velocities */
2037                 state->v[i][m] = -state->v[i][m];
2038             }
2039         }
2040     }
2041     sfree(savex);
2042 }
2043
2044 void calc_enervirdiff(FILE *fplog,int eDispCorr,t_forcerec *fr)
2045 {
2046   double eners[2],virs[2],enersum,virsum,y0,f,g,h;
2047   double r0,r1,r,rc3,rc9,ea,eb,ec,pa,pb,pc,pd;
2048   double invscale,invscale2,invscale3;
2049   int    ri0,ri1,ri,i,offstart,offset;
2050   real   scale,*vdwtab,tabfactor,tmp;
2051
2052   fr->enershiftsix = 0;
2053   fr->enershifttwelve = 0;
2054   fr->enerdiffsix = 0;
2055   fr->enerdifftwelve = 0;
2056   fr->virdiffsix = 0;
2057   fr->virdifftwelve = 0;
2058
2059   if (eDispCorr != edispcNO) {
2060     for(i=0; i<2; i++) {
2061       eners[i] = 0;
2062       virs[i]  = 0;
2063     }
2064     if ((fr->vdwtype == evdwSWITCH) || (fr->vdwtype == evdwSHIFT)) {
2065       if (fr->rvdw_switch == 0)
2066         gmx_fatal(FARGS,
2067                   "With dispersion correction rvdw-switch can not be zero "
2068                   "for vdw-type = %s",evdw_names[fr->vdwtype]);
2069
2070       scale  = fr->nblists[0].table_elec_vdw.scale;
2071       vdwtab = fr->nblists[0].table_vdw.data;
2072
2073       /* Round the cut-offs to exact table values for precision */
2074       ri0 = floor(fr->rvdw_switch*scale);
2075       ri1 = ceil(fr->rvdw*scale);
2076       r0  = ri0/scale;
2077       r1  = ri1/scale;
2078       rc3 = r0*r0*r0;
2079       rc9  = rc3*rc3*rc3;
2080
2081       if (fr->vdwtype == evdwSHIFT)
2082       {
2083           /* Determine the constant energy shift below rvdw_switch.
2084            * Table has a scale factor since we have scaled it down to compensate
2085            * for scaling-up c6/c12 with the derivative factors to save flops in analytical kernels.
2086            */
2087           fr->enershiftsix    = (real)(-1.0/(rc3*rc3)) - 6.0*vdwtab[8*ri0];
2088           fr->enershifttwelve = (real)( 1.0/(rc9*rc3)) - 12.0*vdwtab[8*ri0 + 4];
2089       }
2090       /* Add the constant part from 0 to rvdw_switch.
2091        * This integration from 0 to rvdw_switch overcounts the number
2092        * of interactions by 1, as it also counts the self interaction.
2093        * We will correct for this later.
2094        */
2095       eners[0] += 4.0*M_PI*fr->enershiftsix*rc3/3.0;
2096       eners[1] += 4.0*M_PI*fr->enershifttwelve*rc3/3.0;
2097
2098       invscale = 1.0/(scale);
2099       invscale2 = invscale*invscale;
2100       invscale3 = invscale*invscale2;
2101
2102       /* following summation derived from cubic spline definition,
2103         Numerical Recipies in C, second edition, p. 113-116.  Exact
2104         for the cubic spline.  We first calculate the negative of
2105         the energy from rvdw to rvdw_switch, assuming that g(r)=1,
2106         and then add the more standard, abrupt cutoff correction to
2107         that result, yielding the long-range correction for a
2108         switched function.  We perform both the pressure and energy
2109         loops at the same time for simplicity, as the computational
2110         cost is low. */
2111
2112       for (i=0;i<2;i++) {
2113         enersum = 0.0; virsum = 0.0;
2114         if (i==0)
2115         {
2116             offstart = 0;
2117             /* Since the dispersion table has been scaled down a factor 6.0 and the repulsion
2118              * a factor 12.0 to compensate for the c6/c12 parameters inside nbfp[] being scaled
2119              * up (to save flops in kernels), we need to correct for this.
2120              */
2121             tabfactor = 6.0;
2122         }
2123         else
2124         {
2125             offstart = 4;
2126             tabfactor = 12.0;
2127         }
2128         for (ri=ri0; ri<ri1; ri++) {
2129           r = ri*invscale;
2130           ea = invscale3;
2131           eb = 2.0*invscale2*r;
2132           ec = invscale*r*r;
2133
2134           pa = invscale3;
2135           pb = 3.0*invscale2*r;
2136           pc = 3.0*invscale*r*r;
2137           pd = r*r*r;
2138
2139           /* this "8" is from the packing in the vdwtab array - perhaps should be #define'ed? */
2140           offset = 8*ri + offstart;
2141           y0 = vdwtab[offset];
2142           f  = vdwtab[offset+1];
2143           g  = vdwtab[offset+2];
2144           h  = vdwtab[offset+3];
2145
2146           enersum += y0*(ea/3 + eb/2 + ec) + f*(ea/4 + eb/3 + ec/2) + g*(ea/5 + eb/4 + ec/3) + h*(ea/6 + eb/5 + ec/4);
2147           virsum  += f*(pa/4 + pb/3 + pc/2 + pd) + 2*g*(pa/5 + pb/4 + pc/3 + pd/2) + 3*h*(pa/6 + pb/5 + pc/4 + pd/3);
2148         }
2149           
2150         enersum *= 4.0*M_PI*tabfactor;
2151         virsum  *= 4.0*M_PI*tabfactor;
2152         eners[i] -= enersum;
2153         virs[i]  -= virsum;
2154       }
2155
2156       /* now add the correction for rvdw_switch to infinity */
2157       eners[0] += -4.0*M_PI/(3.0*rc3);
2158       eners[1] +=  4.0*M_PI/(9.0*rc9);
2159       virs[0]  +=  8.0*M_PI/rc3;
2160       virs[1]  += -16.0*M_PI/(3.0*rc9);
2161     }
2162     else if ((fr->vdwtype == evdwCUT) || (fr->vdwtype == evdwUSER)) {
2163       if (fr->vdwtype == evdwUSER && fplog)
2164         fprintf(fplog,
2165                 "WARNING: using dispersion correction with user tables\n");
2166       rc3  = fr->rvdw*fr->rvdw*fr->rvdw;
2167       rc9  = rc3*rc3*rc3;
2168       /* Contribution beyond the cut-off */
2169       eners[0] += -4.0*M_PI/(3.0*rc3);
2170       eners[1] +=  4.0*M_PI/(9.0*rc9);
2171       if (fr->vdw_modifier==eintmodPOTSHIFT) {
2172           /* Contribution within the cut-off */
2173           eners[0] += -4.0*M_PI/(3.0*rc3);
2174           eners[1] +=  4.0*M_PI/(3.0*rc9);
2175       }
2176       /* Contribution beyond the cut-off */
2177       virs[0]  +=  8.0*M_PI/rc3;
2178       virs[1]  += -16.0*M_PI/(3.0*rc9);
2179     } else {
2180       gmx_fatal(FARGS,
2181                 "Dispersion correction is not implemented for vdw-type = %s",
2182                 evdw_names[fr->vdwtype]);
2183     }
2184     fr->enerdiffsix    = eners[0];
2185     fr->enerdifftwelve = eners[1];
2186     /* The 0.5 is due to the Gromacs definition of the virial */
2187     fr->virdiffsix     = 0.5*virs[0];
2188     fr->virdifftwelve  = 0.5*virs[1];
2189   }
2190 }
2191
2192 void calc_dispcorr(FILE *fplog,t_inputrec *ir,t_forcerec *fr,
2193                    gmx_large_int_t step,int natoms,
2194                    matrix box,real lambda,tensor pres,tensor virial,
2195                    real *prescorr, real *enercorr, real *dvdlcorr)
2196 {
2197     gmx_bool bCorrAll,bCorrPres;
2198     real dvdlambda,invvol,dens,ninter,avcsix,avctwelve,enerdiff,svir=0,spres=0;
2199     int  m;
2200
2201     *prescorr = 0;
2202     *enercorr = 0;
2203     *dvdlcorr = 0;
2204
2205     clear_mat(virial);
2206     clear_mat(pres);
2207
2208     if (ir->eDispCorr != edispcNO) {
2209         bCorrAll  = (ir->eDispCorr == edispcAllEner ||
2210                      ir->eDispCorr == edispcAllEnerPres);
2211         bCorrPres = (ir->eDispCorr == edispcEnerPres ||
2212                      ir->eDispCorr == edispcAllEnerPres);
2213
2214         invvol = 1/det(box);
2215         if (fr->n_tpi)
2216         {
2217             /* Only correct for the interactions with the inserted molecule */
2218             dens = (natoms - fr->n_tpi)*invvol;
2219             ninter = fr->n_tpi;
2220         }
2221         else
2222         {
2223             dens = natoms*invvol;
2224             ninter = 0.5*natoms;
2225         }
2226
2227         if (ir->efep == efepNO)
2228         {
2229             avcsix    = fr->avcsix[0];
2230             avctwelve = fr->avctwelve[0];
2231         }
2232         else
2233         {
2234             avcsix    = (1 - lambda)*fr->avcsix[0]    + lambda*fr->avcsix[1];
2235             avctwelve = (1 - lambda)*fr->avctwelve[0] + lambda*fr->avctwelve[1];
2236         }
2237
2238         enerdiff = ninter*(dens*fr->enerdiffsix - fr->enershiftsix);
2239         *enercorr += avcsix*enerdiff;
2240         dvdlambda = 0.0;
2241         if (ir->efep != efepNO)
2242         {
2243             dvdlambda += (fr->avcsix[1] - fr->avcsix[0])*enerdiff;
2244         }
2245         if (bCorrAll)
2246         {
2247             enerdiff = ninter*(dens*fr->enerdifftwelve - fr->enershifttwelve);
2248             *enercorr += avctwelve*enerdiff;
2249             if (fr->efep != efepNO)
2250             {
2251                 dvdlambda += (fr->avctwelve[1] - fr->avctwelve[0])*enerdiff;
2252             }
2253         }
2254
2255         if (bCorrPres)
2256         {
2257             svir = ninter*dens*avcsix*fr->virdiffsix/3.0;
2258             if (ir->eDispCorr == edispcAllEnerPres)
2259             {
2260                 svir += ninter*dens*avctwelve*fr->virdifftwelve/3.0;
2261             }
2262             /* The factor 2 is because of the Gromacs virial definition */
2263             spres = -2.0*invvol*svir*PRESFAC;
2264
2265             for(m=0; m<DIM; m++) {
2266                 virial[m][m] += svir;
2267                 pres[m][m] += spres;
2268             }
2269             *prescorr += spres;
2270         }
2271
2272         /* Can't currently control when it prints, for now, just print when degugging */
2273         if (debug)
2274         {
2275             if (bCorrAll) {
2276                 fprintf(debug,"Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
2277                         avcsix,avctwelve);
2278             }
2279             if (bCorrPres)
2280             {
2281                 fprintf(debug,
2282                         "Long Range LJ corr.: Epot %10g, Pres: %10g, Vir: %10g\n",
2283                         *enercorr,spres,svir);
2284             }
2285             else
2286             {
2287                 fprintf(debug,"Long Range LJ corr.: Epot %10g\n",*enercorr);
2288             }
2289         }
2290
2291         if (fr->bSepDVDL && do_per_step(step,ir->nstlog))
2292         {
2293             fprintf(fplog,sepdvdlformat,"Dispersion correction",
2294                     *enercorr,dvdlambda);
2295         }
2296         if (fr->efep != efepNO)
2297         {
2298             *dvdlcorr += dvdlambda;
2299         }
2300     }
2301 }
2302
2303 void do_pbc_first(FILE *fplog,matrix box,t_forcerec *fr,
2304                   t_graph *graph,rvec x[])
2305 {
2306   if (fplog)
2307     fprintf(fplog,"Removing pbc first time\n");
2308   calc_shifts(box,fr->shift_vec);
2309   if (graph) {
2310     mk_mshift(fplog,graph,fr->ePBC,box,x);
2311     if (gmx_debug_at)
2312       p_graph(debug,"do_pbc_first 1",graph);
2313     shift_self(graph,box,x);
2314     /* By doing an extra mk_mshift the molecules that are broken
2315      * because they were e.g. imported from another software
2316      * will be made whole again. Such are the healing powers
2317      * of GROMACS.
2318      */
2319     mk_mshift(fplog,graph,fr->ePBC,box,x);
2320     if (gmx_debug_at)
2321       p_graph(debug,"do_pbc_first 2",graph);
2322   }
2323   if (fplog)
2324     fprintf(fplog,"Done rmpbc\n");
2325 }
2326
2327 static void low_do_pbc_mtop(FILE *fplog,int ePBC,matrix box,
2328                             gmx_mtop_t *mtop,rvec x[],
2329                             gmx_bool bFirst)
2330 {
2331   t_graph *graph;
2332   int mb,as,mol;
2333   gmx_molblock_t *molb;
2334
2335   if (bFirst && fplog)
2336     fprintf(fplog,"Removing pbc first time\n");
2337
2338   snew(graph,1);
2339   as = 0;
2340   for(mb=0; mb<mtop->nmolblock; mb++) {
2341     molb = &mtop->molblock[mb];
2342     if (molb->natoms_mol == 1 ||
2343         (!bFirst && mtop->moltype[molb->type].cgs.nr == 1)) {
2344       /* Just one atom or charge group in the molecule, no PBC required */
2345       as += molb->nmol*molb->natoms_mol;
2346     } else {
2347       /* Pass NULL iso fplog to avoid graph prints for each molecule type */
2348       mk_graph_ilist(NULL,mtop->moltype[molb->type].ilist,
2349                      0,molb->natoms_mol,FALSE,FALSE,graph);
2350
2351       for(mol=0; mol<molb->nmol; mol++) {
2352         mk_mshift(fplog,graph,ePBC,box,x+as);
2353
2354         shift_self(graph,box,x+as);
2355         /* The molecule is whole now.
2356          * We don't need the second mk_mshift call as in do_pbc_first,
2357          * since we no longer need this graph.
2358          */
2359
2360         as += molb->natoms_mol;
2361       }
2362       done_graph(graph);
2363     }
2364   }
2365   sfree(graph);
2366 }
2367
2368 void do_pbc_first_mtop(FILE *fplog,int ePBC,matrix box,
2369                        gmx_mtop_t *mtop,rvec x[])
2370 {
2371   low_do_pbc_mtop(fplog,ePBC,box,mtop,x,TRUE);
2372 }
2373
2374 void do_pbc_mtop(FILE *fplog,int ePBC,matrix box,
2375                  gmx_mtop_t *mtop,rvec x[])
2376 {
2377   low_do_pbc_mtop(fplog,ePBC,box,mtop,x,FALSE);
2378 }
2379
2380 void finish_run(FILE *fplog,t_commrec *cr,const char *confout,
2381                 t_inputrec *inputrec,
2382                 t_nrnb nrnb[],gmx_wallcycle_t wcycle,
2383                 gmx_runtime_t *runtime,
2384                 wallclock_gpu_t *gputimes,
2385                 int omp_nth_pp,
2386                 gmx_bool bWriteStat)
2387 {
2388     int    i,j;
2389     t_nrnb *nrnb_tot=NULL;
2390     real   delta_t;
2391     double nbfs,mflop;
2392
2393     wallcycle_sum(cr,wcycle);
2394
2395     if (cr->nnodes > 1)
2396     {
2397         snew(nrnb_tot,1);
2398 #ifdef GMX_MPI
2399         MPI_Allreduce(nrnb->n,nrnb_tot->n,eNRNB,MPI_DOUBLE,MPI_SUM,
2400                       cr->mpi_comm_mysim);
2401 #endif
2402     }
2403     else
2404     {
2405         nrnb_tot = nrnb;
2406     }
2407
2408 #if defined(GMX_MPI) && !defined(GMX_THREAD_MPI)
2409     if (cr->nnodes > 1)
2410     {
2411         /* reduce nodetime over all MPI processes in the current simulation */
2412         double sum;
2413         MPI_Allreduce(&runtime->proctime,&sum,1,MPI_DOUBLE,MPI_SUM,
2414                       cr->mpi_comm_mysim);
2415         runtime->proctime = sum;
2416     }
2417 #endif
2418
2419     if (SIMMASTER(cr))
2420     {
2421         print_flop(fplog,nrnb_tot,&nbfs,&mflop);
2422     }
2423     if (cr->nnodes > 1)
2424     {
2425         sfree(nrnb_tot);
2426     }
2427
2428     if ((cr->duty & DUTY_PP) && DOMAINDECOMP(cr))
2429     {
2430         print_dd_statistics(cr,inputrec,fplog);
2431     }
2432
2433 #ifdef GMX_MPI
2434     if (PARTDECOMP(cr))
2435     {
2436         if (MASTER(cr))
2437         {
2438             t_nrnb     *nrnb_all;
2439             int        s;
2440             MPI_Status stat;
2441
2442             snew(nrnb_all,cr->nnodes);
2443             nrnb_all[0] = *nrnb;
2444             for(s=1; s<cr->nnodes; s++)
2445             {
2446                 MPI_Recv(nrnb_all[s].n,eNRNB,MPI_DOUBLE,s,0,
2447                          cr->mpi_comm_mysim,&stat);
2448             }
2449             pr_load(fplog,cr,nrnb_all);
2450             sfree(nrnb_all);
2451         }
2452         else
2453         {
2454             MPI_Send(nrnb->n,eNRNB,MPI_DOUBLE,MASTERRANK(cr),0,
2455                      cr->mpi_comm_mysim);
2456         }
2457     }
2458 #endif
2459
2460     if (SIMMASTER(cr))
2461     {
2462         wallcycle_print(fplog,cr->nnodes,cr->npmenodes,runtime->realtime,
2463                         wcycle,gputimes);
2464
2465         if (EI_DYNAMICS(inputrec->eI))
2466         {
2467             delta_t = inputrec->delta_t;
2468         }
2469         else
2470         {
2471             delta_t = 0;
2472         }
2473
2474         if (fplog)
2475         {
2476             print_perf(fplog,runtime->proctime,runtime->realtime,
2477                        cr->nnodes-cr->npmenodes,
2478                        runtime->nsteps_done,delta_t,nbfs,mflop,
2479                        omp_nth_pp);
2480         }
2481         if (bWriteStat)
2482         {
2483             print_perf(stderr,runtime->proctime,runtime->realtime,
2484                        cr->nnodes-cr->npmenodes,
2485                        runtime->nsteps_done,delta_t,nbfs,mflop,
2486                        omp_nth_pp);
2487         }
2488     }
2489 }
2490
2491 extern void initialize_lambdas(FILE *fplog,t_inputrec *ir,int *fep_state,real *lambda,double *lam0)
2492 {
2493     /* this function works, but could probably use a logic rewrite to keep all the different
2494        types of efep straight. */
2495
2496     int i;
2497     t_lambda *fep = ir->fepvals;
2498
2499     if ((ir->efep==efepNO) && (ir->bSimTemp == FALSE)) {
2500         for (i=0;i<efptNR;i++)  {
2501             lambda[i] = 0.0;
2502             if (lam0)
2503             {
2504                 lam0[i] = 0.0;
2505             }
2506         }
2507         return;
2508     } else {
2509         *fep_state = fep->init_fep_state; /* this might overwrite the checkpoint
2510                                              if checkpoint is set -- a kludge is in for now
2511                                              to prevent this.*/
2512         for (i=0;i<efptNR;i++)
2513         {
2514             /* overwrite lambda state with init_lambda for now for backwards compatibility */
2515             if (fep->init_lambda>=0) /* if it's -1, it was never initializd */
2516             {
2517                 lambda[i] = fep->init_lambda;
2518                 if (lam0) {
2519                     lam0[i] = lambda[i];
2520                 }
2521             }
2522             else
2523             {
2524                 lambda[i] = fep->all_lambda[i][*fep_state];
2525                 if (lam0) {
2526                     lam0[i] = lambda[i];
2527                 }
2528             }
2529         }
2530         if (ir->bSimTemp) {
2531             /* need to rescale control temperatures to match current state */
2532             for (i=0;i<ir->opts.ngtc;i++) {
2533                 if (ir->opts.ref_t[i] > 0) {
2534                     ir->opts.ref_t[i] = ir->simtempvals->temperatures[*fep_state];
2535                 }
2536             }
2537         }
2538     }
2539
2540     /* Send to the log the information on the current lambdas */
2541     if (fplog != NULL)
2542     {
2543         fprintf(fplog,"Initial vector of lambda components:[ ");
2544         for (i=0;i<efptNR;i++)
2545         {
2546             fprintf(fplog,"%10.4f ",lambda[i]);
2547         }
2548         fprintf(fplog,"]\n");
2549     }
2550     return;
2551 }
2552
2553
2554 void init_md(FILE *fplog,
2555              t_commrec *cr,t_inputrec *ir,const output_env_t oenv,
2556              double *t,double *t0,
2557              real *lambda, int *fep_state, double *lam0,
2558              t_nrnb *nrnb,gmx_mtop_t *mtop,
2559              gmx_update_t *upd,
2560              int nfile,const t_filenm fnm[],
2561              gmx_mdoutf_t **outf,t_mdebin **mdebin,
2562              tensor force_vir,tensor shake_vir,rvec mu_tot,
2563              gmx_bool *bSimAnn,t_vcm **vcm, t_state *state, unsigned long Flags)
2564 {
2565     int  i,j,n;
2566     real tmpt,mod;
2567
2568     /* Initial values */
2569     *t = *t0       = ir->init_t;
2570
2571     *bSimAnn=FALSE;
2572     for(i=0;i<ir->opts.ngtc;i++)
2573     {
2574         /* set bSimAnn if any group is being annealed */
2575         if(ir->opts.annealing[i]!=eannNO)
2576         {
2577             *bSimAnn = TRUE;
2578         }
2579     }
2580     if (*bSimAnn)
2581     {
2582         update_annealing_target_temp(&(ir->opts),ir->init_t);
2583     }
2584
2585     /* Initialize lambda variables */
2586     initialize_lambdas(fplog,ir,fep_state,lambda,lam0);
2587
2588     if (upd)
2589     {
2590         *upd = init_update(fplog,ir);
2591     }
2592
2593
2594     if (vcm != NULL)
2595     {
2596         *vcm = init_vcm(fplog,&mtop->groups,ir);
2597     }
2598
2599     if (EI_DYNAMICS(ir->eI) && !(Flags & MD_APPENDFILES))
2600     {
2601         if (ir->etc == etcBERENDSEN)
2602         {
2603             please_cite(fplog,"Berendsen84a");
2604         }
2605         if (ir->etc == etcVRESCALE)
2606         {
2607             please_cite(fplog,"Bussi2007a");
2608         }
2609     }
2610
2611     init_nrnb(nrnb);
2612
2613     if (nfile != -1)
2614     {
2615         *outf = init_mdoutf(nfile,fnm,Flags,cr,ir,oenv);
2616
2617         *mdebin = init_mdebin((Flags & MD_APPENDFILES) ? NULL : (*outf)->fp_ene,
2618                               mtop,ir, (*outf)->fp_dhdl);
2619     }
2620
2621     if (ir->bAdress)
2622     {
2623       please_cite(fplog,"Fritsch12");
2624       please_cite(fplog,"Junghans10");
2625     }
2626     /* Initiate variables */
2627     clear_mat(force_vir);
2628     clear_mat(shake_vir);
2629     clear_rvec(mu_tot);
2630
2631     debug_gmx();
2632 }
2633