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