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