Merge 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_simd_4xn.h"
93 #include "nbnxn_kernels/nbnxn_kernel_simd_2xnn.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, enr_nbnxn_kernel_ljc, enr_nbnxn_kernel_lj;
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 != nbnxnk8x8x8_CUDA)
619     {
620         wallcycle_sub_start(wcycle, ewcsNONBONDED);
621     }
622     switch (nbvg->kernel_type)
623     {
624         case nbnxnk4x4_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 nbnxnk4xN_SIMD_4xN:
638             nbnxn_kernel_simd_4xn(&nbvg->nbl_lists,
639                                   nbvg->nbat, ic,
640                                   nbvg->ewald_excl,
641                                   fr->shift_vec,
642                                   flags,
643                                   clearF,
644                                   fr->fshift[0],
645                                   enerd->grpp.ener[egCOULSR],
646                                   fr->bBHAM ?
647                                   enerd->grpp.ener[egBHAMSR] :
648                                   enerd->grpp.ener[egLJSR]);
649             break;
650         case nbnxnk4xN_SIMD_2xNN:
651             nbnxn_kernel_simd_2xnn(&nbvg->nbl_lists,
652                                    nbvg->nbat, ic,
653                                    nbvg->ewald_excl,
654                                    fr->shift_vec,
655                                    flags,
656                                    clearF,
657                                    fr->fshift[0],
658                                    enerd->grpp.ener[egCOULSR],
659                                    fr->bBHAM ?
660                                    enerd->grpp.ener[egBHAMSR] :
661                                    enerd->grpp.ener[egLJSR]);
662             break;
663
664         case nbnxnk8x8x8_CUDA:
665             nbnxn_cuda_launch_kernel(fr->nbv->cu_nbv, nbvg->nbat, flags, ilocality);
666             break;
667
668         case nbnxnk8x8x8_PlainC:
669             nbnxn_kernel_gpu_ref(nbvg->nbl_lists.nbl[0],
670                                  nbvg->nbat, ic,
671                                  fr->shift_vec,
672                                  flags,
673                                  clearF,
674                                  nbvg->nbat->out[0].f,
675                                  fr->fshift[0],
676                                  enerd->grpp.ener[egCOULSR],
677                                  fr->bBHAM ?
678                                  enerd->grpp.ener[egBHAMSR] :
679                                  enerd->grpp.ener[egLJSR]);
680             break;
681
682         default:
683             gmx_incons("Invalid nonbonded kernel type passed!");
684
685     }
686     if (nbvg->kernel_type != nbnxnk8x8x8_CUDA)
687     {
688         wallcycle_sub_stop(wcycle, ewcsNONBONDED);
689     }
690
691     if (EEL_RF(ic->eeltype) || ic->eeltype == eelCUT)
692     {
693         enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_RF;
694     }
695     else if (nbvg->ewald_excl == ewaldexclTable)
696     {
697         enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_TAB;
698     }
699     else
700     {
701         enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_EWALD;
702     }
703     enr_nbnxn_kernel_lj = eNR_NBNXN_LJ;
704     if (flags & GMX_FORCE_ENERGY)
705     {
706         /* In eNR_??? the nbnxn F+E kernels are always the F kernel + 1 */
707         enr_nbnxn_kernel_ljc += 1;
708         enr_nbnxn_kernel_lj  += 1;
709     }
710
711     inc_nrnb(nrnb,enr_nbnxn_kernel_ljc,
712              nbvg->nbl_lists.natpair_ljq);
713     inc_nrnb(nrnb,enr_nbnxn_kernel_lj,
714              nbvg->nbl_lists.natpair_lj);
715     inc_nrnb(nrnb,enr_nbnxn_kernel_ljc-eNR_NBNXN_LJ_RF+eNR_NBNXN_RF,
716              nbvg->nbl_lists.natpair_q);
717 }
718
719 void do_force_cutsVERLET(FILE *fplog,t_commrec *cr,
720               t_inputrec *inputrec,
721               gmx_large_int_t step,t_nrnb *nrnb,gmx_wallcycle_t wcycle,
722               gmx_localtop_t *top,
723               gmx_mtop_t *mtop,
724               gmx_groups_t *groups,
725               matrix box,rvec x[],history_t *hist,
726               rvec f[],
727               tensor vir_force,
728               t_mdatoms *mdatoms,
729               gmx_enerdata_t *enerd,t_fcdata *fcd,
730               real *lambda,t_graph *graph,
731               t_forcerec *fr, interaction_const_t *ic,
732               gmx_vsite_t *vsite,rvec mu_tot,
733               double t,FILE *field,gmx_edsam_t ed,
734               gmx_bool bBornRadii,
735               int flags)
736 {
737     int     cg0,cg1,i,j;
738     int     start,homenr;
739     int     nb_kernel_type;
740     double  mu[2*DIM];
741     gmx_bool   bSepDVDL,bStateChanged,bNS,bFillGrid,bCalcCGCM,bBS;
742     gmx_bool   bDoLongRange,bDoForces,bSepLRF,bUseGPU,bUseOrEmulGPU;
743     gmx_bool   bDiffKernels=FALSE;
744     matrix  boxs;
745     rvec    vzero,box_diag;
746     real    e,v,dvdl;
747     float  cycles_pme,cycles_force;
748     nonbonded_verlet_t *nbv;
749
750     cycles_force = 0;
751     nbv = fr->nbv;
752     nb_kernel_type = fr->nbv->grp[0].kernel_type;
753
754     start  = mdatoms->start;
755     homenr = mdatoms->homenr;
756
757     bSepDVDL = (fr->bSepDVDL && do_per_step(step,inputrec->nstlog));
758
759     clear_mat(vir_force);
760
761     cg0 = 0;
762     if (DOMAINDECOMP(cr))
763     {
764         cg1 = cr->dd->ncg_tot;
765     }
766     else
767     {
768         cg1 = top->cgs.nr;
769     }
770     if (fr->n_tpi > 0)
771     {
772         cg1--;
773     }
774
775     bStateChanged = (flags & GMX_FORCE_STATECHANGED);
776     bNS           = (flags & GMX_FORCE_NS) && (fr->bAllvsAll==FALSE); 
777     bFillGrid     = (bNS && bStateChanged);
778     bCalcCGCM     = (bFillGrid && !DOMAINDECOMP(cr));
779     bDoLongRange  = (fr->bTwinRange && bNS && (flags & GMX_FORCE_DO_LR));
780     bDoForces     = (flags & GMX_FORCE_FORCES);
781     bSepLRF       = (bDoLongRange && bDoForces && (flags & GMX_FORCE_SEPLRF));
782     bUseGPU       = fr->nbv->bUseGPU;
783     bUseOrEmulGPU = bUseGPU || (nbv->grp[0].kernel_type == nbnxnk8x8x8_PlainC);
784
785     if (bStateChanged)
786     {
787         update_forcerec(fplog,fr,box);
788
789         if (NEED_MUTOT(*inputrec))
790         {
791             /* Calculate total (local) dipole moment in a temporary common array.
792              * This makes it possible to sum them over nodes faster.
793              */
794             calc_mu(start,homenr,
795                     x,mdatoms->chargeA,mdatoms->chargeB,mdatoms->nChargePerturbed,
796                     mu,mu+DIM);
797         }
798     }
799
800     if (fr->ePBC != epbcNONE) { 
801         /* Compute shift vectors every step,
802          * because of pressure coupling or box deformation!
803          */
804         if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
805             calc_shifts(box,fr->shift_vec);
806
807         if (bCalcCGCM) { 
808             put_atoms_in_box_omp(fr->ePBC,box,homenr,x);
809             inc_nrnb(nrnb,eNR_SHIFTX,homenr);
810         } 
811         else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph) {
812             unshift_self(graph,box,x);
813         }
814     } 
815
816     nbnxn_atomdata_copy_shiftvec(flags & GMX_FORCE_DYNAMICBOX,
817                                   fr->shift_vec,nbv->grp[0].nbat);
818
819 #ifdef GMX_MPI
820     if (!(cr->duty & DUTY_PME)) {
821         /* Send particle coordinates to the pme nodes.
822          * Since this is only implemented for domain decomposition
823          * and domain decomposition does not use the graph,
824          * we do not need to worry about shifting.
825          */    
826
827         wallcycle_start(wcycle,ewcPP_PMESENDX);
828
829         bBS = (inputrec->nwall == 2);
830         if (bBS) {
831             copy_mat(box,boxs);
832             svmul(inputrec->wall_ewald_zfac,boxs[ZZ],boxs[ZZ]);
833         }
834
835         gmx_pme_send_x(cr,bBS ? boxs : box,x,
836                        mdatoms->nChargePerturbed,lambda[efptCOUL],
837                        (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),step);
838
839         wallcycle_stop(wcycle,ewcPP_PMESENDX);
840     }
841 #endif /* GMX_MPI */
842
843     /* do gridding for pair search */
844     if (bNS)
845     {
846         if (graph && bStateChanged)
847         {
848             /* Calculate intramolecular shift vectors to make molecules whole */
849             mk_mshift(fplog,graph,fr->ePBC,box,x);
850         }
851
852         clear_rvec(vzero);
853         box_diag[XX] = box[XX][XX];
854         box_diag[YY] = box[YY][YY];
855         box_diag[ZZ] = box[ZZ][ZZ];
856
857         wallcycle_start(wcycle,ewcNS);
858         if (!fr->bDomDec)
859         {
860             wallcycle_sub_start(wcycle,ewcsNBS_GRID_LOCAL);
861             nbnxn_put_on_grid(nbv->nbs,fr->ePBC,box,
862                               0,vzero,box_diag,
863                               0,mdatoms->homenr,-1,fr->cginfo,x,
864                               0,NULL,
865                               nbv->grp[eintLocal].kernel_type,
866                               nbv->grp[eintLocal].nbat);
867             wallcycle_sub_stop(wcycle,ewcsNBS_GRID_LOCAL);
868         }
869         else
870         {
871             wallcycle_sub_start(wcycle,ewcsNBS_GRID_NONLOCAL);
872             nbnxn_put_on_grid_nonlocal(nbv->nbs,domdec_zones(cr->dd),
873                                        fr->cginfo,x,
874                                        nbv->grp[eintNonlocal].kernel_type,
875                                        nbv->grp[eintNonlocal].nbat);
876             wallcycle_sub_stop(wcycle,ewcsNBS_GRID_NONLOCAL);
877         }
878
879         if (nbv->ngrp == 1 ||
880             nbv->grp[eintNonlocal].nbat == nbv->grp[eintLocal].nbat)
881         {
882             nbnxn_atomdata_set(nbv->grp[eintLocal].nbat,eatAll,
883                                 nbv->nbs,mdatoms,fr->cginfo);
884         }
885         else
886         {
887             nbnxn_atomdata_set(nbv->grp[eintLocal].nbat,eatLocal,
888                                 nbv->nbs,mdatoms,fr->cginfo);
889             nbnxn_atomdata_set(nbv->grp[eintNonlocal].nbat,eatAll,
890                                 nbv->nbs,mdatoms,fr->cginfo);
891         }
892         wallcycle_stop(wcycle, ewcNS);
893     }
894
895     /* initialize the GPU atom data and copy shift vector */
896     if (bUseGPU)
897     {
898         if (bNS)
899         {
900             wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
901             nbnxn_cuda_init_atomdata(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
902             wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
903         }
904
905         wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
906         nbnxn_cuda_upload_shiftvec(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
907         wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
908     }
909
910     /* do local pair search */
911     if (bNS)
912     {
913         wallcycle_start_nocount(wcycle,ewcNS);
914         wallcycle_sub_start(wcycle,ewcsNBS_SEARCH_LOCAL);
915         nbnxn_make_pairlist(nbv->nbs,nbv->grp[eintLocal].nbat,
916                             &top->excls,
917                             ic->rlist,
918                             nbv->min_ci_balanced,
919                             &nbv->grp[eintLocal].nbl_lists,
920                             eintLocal,
921                             nbv->grp[eintLocal].kernel_type,
922                             nrnb);
923         wallcycle_sub_stop(wcycle,ewcsNBS_SEARCH_LOCAL);
924
925         if (bUseGPU)
926         {
927             /* initialize local pair-list on the GPU */
928             nbnxn_cuda_init_pairlist(nbv->cu_nbv,
929                                      nbv->grp[eintLocal].nbl_lists.nbl[0],
930                                      eintLocal);
931         }
932         wallcycle_stop(wcycle, ewcNS);
933     }
934     else
935     {
936         wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
937         wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
938         nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs,eatLocal,FALSE,x,
939                                         nbv->grp[eintLocal].nbat);
940         wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
941         wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
942     }
943
944     if (bUseGPU)
945     {
946         wallcycle_start(wcycle,ewcLAUNCH_GPU_NB);
947         /* launch local nonbonded F on GPU */
948         do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFNo,
949                      nrnb, wcycle);
950         wallcycle_stop(wcycle,ewcLAUNCH_GPU_NB);
951     }
952
953     /* Communicate coordinates and sum dipole if necessary + 
954        do non-local pair search */
955     if (DOMAINDECOMP(cr))
956     {
957         bDiffKernels = (nbv->grp[eintNonlocal].kernel_type !=
958                         nbv->grp[eintLocal].kernel_type);
959
960         if (bDiffKernels)
961         {
962             /* With GPU+CPU non-bonded calculations we need to copy
963              * the local coordinates to the non-local nbat struct
964              * (in CPU format) as the non-local kernel call also
965              * calculates the local - non-local interactions.
966              */
967             wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
968             wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
969             nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs,eatLocal,TRUE,x,
970                                              nbv->grp[eintNonlocal].nbat);
971             wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
972             wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
973         }
974
975         if (bNS)
976         {
977             wallcycle_start_nocount(wcycle,ewcNS);
978             wallcycle_sub_start(wcycle,ewcsNBS_SEARCH_NONLOCAL);
979
980             if (bDiffKernels)
981             {
982                 nbnxn_grid_add_simple(nbv->nbs,nbv->grp[eintNonlocal].nbat);
983             }
984
985             nbnxn_make_pairlist(nbv->nbs,nbv->grp[eintNonlocal].nbat,
986                                 &top->excls,
987                                 ic->rlist,
988                                 nbv->min_ci_balanced,
989                                 &nbv->grp[eintNonlocal].nbl_lists,
990                                 eintNonlocal,
991                                 nbv->grp[eintNonlocal].kernel_type,
992                                 nrnb);
993
994             wallcycle_sub_stop(wcycle,ewcsNBS_SEARCH_NONLOCAL);
995
996             if (nbv->grp[eintNonlocal].kernel_type == nbnxnk8x8x8_CUDA)
997             {
998                 /* initialize non-local pair-list on the GPU */
999                 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
1000                                          nbv->grp[eintNonlocal].nbl_lists.nbl[0],
1001                                          eintNonlocal);
1002             }
1003             wallcycle_stop(wcycle,ewcNS);
1004         } 
1005         else
1006         {
1007             wallcycle_start(wcycle,ewcMOVEX);
1008             dd_move_x(cr->dd,box,x);
1009
1010             /* When we don't need the total dipole we sum it in global_stat */
1011             if (bStateChanged && NEED_MUTOT(*inputrec))
1012             {
1013                 gmx_sumd(2*DIM,mu,cr);
1014             }
1015             wallcycle_stop(wcycle,ewcMOVEX);
1016
1017             wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1018             wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
1019             nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs,eatNonlocal,FALSE,x,
1020                                             nbv->grp[eintNonlocal].nbat);
1021             wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
1022             cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1023         }
1024
1025         if (bUseGPU && !bDiffKernels)
1026         { 
1027             wallcycle_start(wcycle,ewcLAUNCH_GPU_NB);
1028             /* launch non-local nonbonded F on GPU */
1029             do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFNo,
1030                          nrnb, wcycle);
1031             cycles_force += wallcycle_stop(wcycle,ewcLAUNCH_GPU_NB);
1032         }
1033     }
1034
1035     if (bUseGPU)
1036     {
1037         /* launch D2H copy-back F */
1038         wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1039         if (DOMAINDECOMP(cr) && !bDiffKernels)
1040         {
1041             nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintNonlocal].nbat,
1042                                       flags, eatNonlocal);
1043         }
1044         nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintLocal].nbat,
1045                                   flags, eatLocal);
1046         cycles_force += wallcycle_stop(wcycle,ewcLAUNCH_GPU_NB);
1047     }
1048
1049     if (bStateChanged && NEED_MUTOT(*inputrec))
1050     {
1051         if (PAR(cr))
1052         {
1053             gmx_sumd(2*DIM,mu,cr);
1054         } 
1055
1056         for(i=0; i<2; i++)
1057         {
1058             for(j=0;j<DIM;j++)
1059             {
1060                 fr->mu_tot[i][j] = mu[i*DIM + j];
1061             }
1062         }
1063     }
1064     if (fr->efep == efepNO)
1065     {
1066         copy_rvec(fr->mu_tot[0],mu_tot);
1067     }
1068     else
1069     {
1070         for(j=0; j<DIM; j++)
1071         {
1072             mu_tot[j] =
1073                 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] +
1074                 lambda[efptCOUL]*fr->mu_tot[1][j];
1075         }
1076     }
1077
1078     /* Reset energies */
1079     reset_enerdata(&(inputrec->opts),fr,bNS,enerd,MASTER(cr));
1080     clear_rvecs(SHIFTS,fr->fshift);
1081
1082     if (DOMAINDECOMP(cr))
1083     {
1084         if (!(cr->duty & DUTY_PME))
1085         {
1086             wallcycle_start(wcycle,ewcPPDURINGPME);
1087             dd_force_flop_start(cr->dd,nrnb);
1088         }
1089     }
1090     
1091     /* Start the force cycle counter.
1092      * This counter is stopped in do_forcelow_level.
1093      * No parallel communication should occur while this counter is running,
1094      * since that will interfere with the dynamic load balancing.
1095      */
1096     wallcycle_start(wcycle,ewcFORCE);
1097     if (bDoForces)
1098     {
1099         /* Reset forces for which the virial is calculated separately:
1100          * PME/Ewald forces if necessary */
1101         if (fr->bF_NoVirSum) 
1102         {
1103             if (flags & GMX_FORCE_VIRIAL)
1104             {
1105                 fr->f_novirsum = fr->f_novirsum_alloc;
1106                 if (fr->bDomDec)
1107                 {
1108                     clear_rvecs(fr->f_novirsum_n,fr->f_novirsum);
1109                 }
1110                 else
1111                 {
1112                     clear_rvecs(homenr,fr->f_novirsum+start);
1113                 }
1114             }
1115             else
1116             {
1117                 /* We are not calculating the pressure so we do not need
1118                  * a separate array for forces that do not contribute
1119                  * to the pressure.
1120                  */
1121                 fr->f_novirsum = f;
1122             }
1123         }
1124
1125         /* Clear the short- and long-range forces */
1126         clear_rvecs(fr->natoms_force_constr,f);
1127         if(bSepLRF && do_per_step(step,inputrec->nstcalclr))
1128         {
1129             clear_rvecs(fr->natoms_force_constr,fr->f_twin);
1130         }
1131         
1132         clear_rvec(fr->vir_diag_posres);
1133     }
1134     if (inputrec->ePull == epullCONSTRAINT)
1135     {
1136         clear_pull_forces(inputrec->pull);
1137     }
1138
1139     /* update QMMMrec, if necessary */
1140     if(fr->bQMMM)
1141     {
1142         update_QMMMrec(cr,fr,x,mdatoms,box,top);
1143     }
1144
1145     if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1146     {
1147         posres_wrapper(fplog,flags,bSepDVDL,inputrec,nrnb,top,box,x,
1148                        f,enerd,lambda,fr);
1149     }
1150
1151     /* Compute the bonded and non-bonded energies and optionally forces */    
1152     do_force_lowlevel(fplog,step,fr,inputrec,&(top->idef),
1153                       cr,nrnb,wcycle,mdatoms,&(inputrec->opts),
1154                       x,hist,f, bSepLRF ? fr->f_twin : f,enerd,fcd,mtop,top,fr->born,
1155                       &(top->atomtypes),bBornRadii,box,
1156                       inputrec->fepvals,lambda,graph,&(top->excls),fr->mu_tot,
1157                       flags, &cycles_pme);
1158
1159     if(bSepLRF)
1160     {
1161         if (do_per_step(step,inputrec->nstcalclr))
1162         {
1163             /* Add the long range forces to the short range forces */
1164             for(i=0; i<fr->natoms_force_constr; i++)
1165             {
1166                 rvec_add(fr->f_twin[i],f[i],f[i]);
1167             }
1168         }
1169     }
1170     
1171     if (!bUseOrEmulGPU)
1172     {
1173         /* Maybe we should move this into do_force_lowlevel */
1174         do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFYes,
1175                      nrnb, wcycle);
1176     }
1177         
1178
1179     if (!bUseOrEmulGPU || bDiffKernels)
1180     {
1181         int aloc;
1182
1183         if (DOMAINDECOMP(cr))
1184         {
1185             do_nb_verlet(fr, ic, enerd, flags, eintNonlocal,
1186                          bDiffKernels ? enbvClearFYes : enbvClearFNo,
1187                          nrnb, wcycle);
1188         }
1189
1190         if (!bUseOrEmulGPU)
1191         {
1192             aloc = eintLocal;
1193         }
1194         else
1195         {
1196             aloc = eintNonlocal;
1197         }
1198
1199         /* Add all the non-bonded force to the normal force array.
1200          * This can be split into a local a non-local part when overlapping
1201          * communication with calculation with domain decomposition.
1202          */
1203         cycles_force += wallcycle_stop(wcycle,ewcFORCE);
1204         wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1205         wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1206         nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs,eatAll,nbv->grp[aloc].nbat,f);
1207         wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1208         cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1209         wallcycle_start_nocount(wcycle,ewcFORCE);
1210
1211         /* if there are multiple fshift output buffers reduce them */
1212         if ((flags & GMX_FORCE_VIRIAL) &&
1213             nbv->grp[aloc].nbl_lists.nnbl > 1)
1214         {
1215             nbnxn_atomdata_add_nbat_fshift_to_fshift(nbv->grp[aloc].nbat,
1216                                                       fr->fshift);
1217         }
1218     }
1219     
1220     cycles_force += wallcycle_stop(wcycle,ewcFORCE);
1221     
1222     if (ed)
1223     {
1224         do_flood(cr,inputrec,x,f,ed,box,step,bNS);
1225     }
1226
1227     if (bUseOrEmulGPU && !bDiffKernels)
1228     {
1229         /* wait for non-local forces (or calculate in emulation mode) */
1230         if (DOMAINDECOMP(cr))
1231         {
1232             if (bUseGPU)
1233             {
1234                 wallcycle_start(wcycle,ewcWAIT_GPU_NB_NL);
1235                 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1236                                     nbv->grp[eintNonlocal].nbat,
1237                                     flags, eatNonlocal,
1238                                     enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1239                                     fr->fshift);
1240                 cycles_force += wallcycle_stop(wcycle,ewcWAIT_GPU_NB_NL);
1241             }
1242             else
1243             {
1244                 wallcycle_start_nocount(wcycle,ewcFORCE);
1245                 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFYes,
1246                              nrnb, wcycle);
1247                 cycles_force += wallcycle_stop(wcycle,ewcFORCE);
1248             }            
1249             wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1250             wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1251             /* skip the reduction if there was no non-local work to do */
1252             if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1253             {
1254                 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs,eatNonlocal,
1255                                                nbv->grp[eintNonlocal].nbat,f);
1256             }
1257             wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1258             cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1259         }
1260     }
1261
1262     if (bDoForces)
1263     {
1264         /* Communicate the forces */
1265         if (PAR(cr))
1266         {
1267             wallcycle_start(wcycle,ewcMOVEF);
1268             if (DOMAINDECOMP(cr))
1269             {
1270                 dd_move_f(cr->dd,f,fr->fshift);
1271                 /* Do we need to communicate the separate force array
1272                  * for terms that do not contribute to the single sum virial?
1273                  * Position restraints and electric fields do not introduce
1274                  * inter-cg forces, only full electrostatics methods do.
1275                  * When we do not calculate the virial, fr->f_novirsum = f,
1276                  * so we have already communicated these forces.
1277                  */
1278                 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1279                     (flags & GMX_FORCE_VIRIAL))
1280                 {
1281                     dd_move_f(cr->dd,fr->f_novirsum,NULL);
1282                 }
1283                 if (bSepLRF)
1284                 {
1285                     /* We should not update the shift forces here,
1286                      * since f_twin is already included in f.
1287                      */
1288                     dd_move_f(cr->dd,fr->f_twin,NULL);
1289                 }
1290             }
1291             wallcycle_stop(wcycle,ewcMOVEF);
1292         }
1293     }
1294  
1295     if (bUseOrEmulGPU)
1296     {
1297         /* wait for local forces (or calculate in emulation mode) */
1298         if (bUseGPU)
1299         {
1300             wallcycle_start(wcycle,ewcWAIT_GPU_NB_L);
1301             nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1302                                 nbv->grp[eintLocal].nbat,
1303                                 flags, eatLocal,
1304                                 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1305                                 fr->fshift);
1306             wallcycle_stop(wcycle,ewcWAIT_GPU_NB_L);
1307
1308             /* now clear the GPU outputs while we finish the step on the CPU */
1309
1310             wallcycle_start_nocount(wcycle,ewcLAUNCH_GPU_NB);
1311             nbnxn_cuda_clear_outputs(nbv->cu_nbv, flags);
1312             wallcycle_stop(wcycle,ewcLAUNCH_GPU_NB);
1313         }
1314         else
1315         {            
1316             wallcycle_start_nocount(wcycle,ewcFORCE);
1317             do_nb_verlet(fr, ic, enerd, flags, eintLocal,
1318                          DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
1319                          nrnb, wcycle);
1320             wallcycle_stop(wcycle,ewcFORCE);
1321         }
1322         wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1323         wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1324         if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1325         {
1326             /* skip the reduction if there was no non-local work to do */
1327             nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs,eatLocal,
1328                                            nbv->grp[eintLocal].nbat,f);
1329         }
1330         wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1331         wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1332     }
1333     
1334     if (DOMAINDECOMP(cr))
1335     {
1336         dd_force_flop_stop(cr->dd,nrnb);
1337         if (wcycle)
1338         {
1339             dd_cycles_add(cr->dd,cycles_force-cycles_pme,ddCyclF);
1340         }
1341     }
1342
1343     if (bDoForces)
1344     {
1345         if (IR_ELEC_FIELD(*inputrec))
1346         {
1347             /* Compute forces due to electric field */
1348             calc_f_el(MASTER(cr) ? field : NULL,
1349                       start,homenr,mdatoms->chargeA,x,fr->f_novirsum,
1350                       inputrec->ex,inputrec->et,t);
1351         }
1352
1353         /* If we have NoVirSum forces, but we do not calculate the virial,
1354          * we sum fr->f_novirum=f later.
1355          */
1356         if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1357         {
1358             wallcycle_start(wcycle,ewcVSITESPREAD);
1359             spread_vsite_f(fplog,vsite,x,f,fr->fshift,FALSE,NULL,nrnb,
1360                            &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
1361             wallcycle_stop(wcycle,ewcVSITESPREAD);
1362
1363             if (bSepLRF)
1364             {
1365                 wallcycle_start(wcycle,ewcVSITESPREAD);
1366                 spread_vsite_f(fplog,vsite,x,fr->f_twin,NULL,FALSE,NULL,
1367                                nrnb,
1368                                &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
1369                 wallcycle_stop(wcycle,ewcVSITESPREAD);
1370             }
1371         }
1372
1373         if (flags & GMX_FORCE_VIRIAL)
1374         {
1375             /* Calculation of the virial must be done after vsites! */
1376             calc_virial(fplog,mdatoms->start,mdatoms->homenr,x,f,
1377                         vir_force,graph,box,nrnb,fr,inputrec->ePBC);
1378         }
1379     }
1380
1381     if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1382     {
1383         pull_potential_wrapper(fplog,bSepDVDL,cr,inputrec,box,x,
1384                                f,vir_force,mdatoms,enerd,lambda,t);
1385     }
1386
1387     if (PAR(cr) && !(cr->duty & DUTY_PME))
1388     {
1389         /* In case of node-splitting, the PP nodes receive the long-range 
1390          * forces, virial and energy from the PME nodes here.
1391          */    
1392         pme_receive_force_ener(fplog,bSepDVDL,cr,wcycle,enerd,fr);
1393     }
1394
1395     if (bDoForces)
1396     {
1397         post_process_forces(fplog,cr,step,nrnb,wcycle,
1398                             top,box,x,f,vir_force,mdatoms,graph,fr,vsite,
1399                             flags);
1400     }
1401     
1402     /* Sum the potential energy terms from group contributions */
1403     sum_epot(&(inputrec->opts),&(enerd->grpp),enerd->term);
1404 }
1405
1406 void do_force_cutsGROUP(FILE *fplog,t_commrec *cr,
1407               t_inputrec *inputrec,
1408               gmx_large_int_t step,t_nrnb *nrnb,gmx_wallcycle_t wcycle,
1409               gmx_localtop_t *top,
1410               gmx_mtop_t *mtop,
1411               gmx_groups_t *groups,
1412               matrix box,rvec x[],history_t *hist,
1413               rvec f[],
1414               tensor vir_force,
1415               t_mdatoms *mdatoms,
1416               gmx_enerdata_t *enerd,t_fcdata *fcd,
1417               real *lambda,t_graph *graph,
1418               t_forcerec *fr,gmx_vsite_t *vsite,rvec mu_tot,
1419               double t,FILE *field,gmx_edsam_t ed,
1420               gmx_bool bBornRadii,
1421               int flags)
1422 {
1423     int    cg0,cg1,i,j;
1424     int    start,homenr;
1425     double mu[2*DIM];
1426     gmx_bool   bSepDVDL,bStateChanged,bNS,bFillGrid,bCalcCGCM,bBS;
1427     gmx_bool   bDoLongRangeNS,bDoForces,bDoPotential,bSepLRF;
1428     gmx_bool   bDoAdressWF;
1429     matrix boxs;
1430     rvec   vzero,box_diag;
1431     real   e,v,dvdlambda[efptNR];
1432     t_pbc  pbc;
1433     float  cycles_pme,cycles_force;
1434
1435     start  = mdatoms->start;
1436     homenr = mdatoms->homenr;
1437
1438     bSepDVDL = (fr->bSepDVDL && do_per_step(step,inputrec->nstlog));
1439
1440     clear_mat(vir_force);
1441
1442     if (PARTDECOMP(cr))
1443     {
1444         pd_cg_range(cr,&cg0,&cg1);
1445     }
1446     else
1447     {
1448         cg0 = 0;
1449         if (DOMAINDECOMP(cr))
1450         {
1451             cg1 = cr->dd->ncg_tot;
1452         }
1453         else
1454         {
1455             cg1 = top->cgs.nr;
1456         }
1457         if (fr->n_tpi > 0)
1458         {
1459             cg1--;
1460         }
1461     }
1462
1463     bStateChanged  = (flags & GMX_FORCE_STATECHANGED);
1464     bNS            = (flags & GMX_FORCE_NS) && (fr->bAllvsAll==FALSE);
1465     /* Should we update the long-range neighborlists at this step? */
1466     bDoLongRangeNS = fr->bTwinRange && bNS;
1467     /* Should we perform the long-range nonbonded evaluation inside the neighborsearching? */
1468     bFillGrid      = (bNS && bStateChanged);
1469     bCalcCGCM      = (bFillGrid && !DOMAINDECOMP(cr));
1470     bDoForces      = (flags & GMX_FORCE_FORCES);
1471     bDoPotential   = (flags & GMX_FORCE_ENERGY);
1472     bSepLRF        = ((inputrec->nstcalclr>1) && bDoForces &&
1473                       (flags & GMX_FORCE_SEPLRF) && (flags & GMX_FORCE_DO_LR));
1474
1475     /* should probably move this to the forcerec since it doesn't change */
1476     bDoAdressWF   = ((fr->adress_type!=eAdressOff));
1477
1478     if (bStateChanged)
1479     {
1480         update_forcerec(fplog,fr,box);
1481
1482         if (NEED_MUTOT(*inputrec))
1483         {
1484             /* Calculate total (local) dipole moment in a temporary common array.
1485              * This makes it possible to sum them over nodes faster.
1486              */
1487             calc_mu(start,homenr,
1488                     x,mdatoms->chargeA,mdatoms->chargeB,mdatoms->nChargePerturbed,
1489                     mu,mu+DIM);
1490         }
1491     }
1492
1493     if (fr->ePBC != epbcNONE) { 
1494         /* Compute shift vectors every step,
1495          * because of pressure coupling or box deformation!
1496          */
1497         if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
1498             calc_shifts(box,fr->shift_vec);
1499
1500         if (bCalcCGCM) { 
1501             put_charge_groups_in_box(fplog,cg0,cg1,fr->ePBC,box,
1502                     &(top->cgs),x,fr->cg_cm);
1503             inc_nrnb(nrnb,eNR_CGCM,homenr);
1504             inc_nrnb(nrnb,eNR_RESETX,cg1-cg0);
1505         } 
1506         else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph) {
1507             unshift_self(graph,box,x);
1508         }
1509     } 
1510     else if (bCalcCGCM) {
1511         calc_cgcm(fplog,cg0,cg1,&(top->cgs),x,fr->cg_cm);
1512         inc_nrnb(nrnb,eNR_CGCM,homenr);
1513     }
1514
1515     if (bCalcCGCM) {
1516         if (PAR(cr)) {
1517             move_cgcm(fplog,cr,fr->cg_cm);
1518         }
1519         if (gmx_debug_at)
1520             pr_rvecs(debug,0,"cgcm",fr->cg_cm,top->cgs.nr);
1521     }
1522
1523 #ifdef GMX_MPI
1524     if (!(cr->duty & DUTY_PME)) {
1525         /* Send particle coordinates to the pme nodes.
1526          * Since this is only implemented for domain decomposition
1527          * and domain decomposition does not use the graph,
1528          * we do not need to worry about shifting.
1529          */    
1530
1531         wallcycle_start(wcycle,ewcPP_PMESENDX);
1532
1533         bBS = (inputrec->nwall == 2);
1534         if (bBS) {
1535             copy_mat(box,boxs);
1536             svmul(inputrec->wall_ewald_zfac,boxs[ZZ],boxs[ZZ]);
1537         }
1538
1539         gmx_pme_send_x(cr,bBS ? boxs : box,x,
1540                        mdatoms->nChargePerturbed,lambda[efptCOUL],
1541                        (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),step);
1542
1543         wallcycle_stop(wcycle,ewcPP_PMESENDX);
1544     }
1545 #endif /* GMX_MPI */
1546
1547     /* Communicate coordinates and sum dipole if necessary */
1548     if (PAR(cr))
1549     {
1550         wallcycle_start(wcycle,ewcMOVEX);
1551         if (DOMAINDECOMP(cr))
1552         {
1553             dd_move_x(cr->dd,box,x);
1554         }
1555         else
1556         {
1557             move_x(fplog,cr,GMX_LEFT,GMX_RIGHT,x,nrnb);
1558         }
1559         wallcycle_stop(wcycle,ewcMOVEX);
1560     }
1561
1562     /* update adress weight beforehand */
1563     if(bStateChanged && bDoAdressWF)
1564     {
1565         /* need pbc for adress weight calculation with pbc_dx */
1566         set_pbc(&pbc,inputrec->ePBC,box);
1567         if(fr->adress_site == eAdressSITEcog)
1568         {
1569             update_adress_weights_cog(top->idef.iparams,top->idef.il,x,fr,mdatoms,
1570                                       inputrec->ePBC==epbcNONE ? NULL : &pbc);
1571         }
1572         else if (fr->adress_site == eAdressSITEcom)
1573         {
1574             update_adress_weights_com(fplog,cg0,cg1,&(top->cgs),x,fr,mdatoms,
1575                                       inputrec->ePBC==epbcNONE ? NULL : &pbc);
1576         }
1577         else if (fr->adress_site == eAdressSITEatomatom){
1578             update_adress_weights_atom_per_atom(cg0,cg1,&(top->cgs),x,fr,mdatoms,
1579                                                 inputrec->ePBC==epbcNONE ? NULL : &pbc);
1580         }
1581         else
1582         {
1583             update_adress_weights_atom(cg0,cg1,&(top->cgs),x,fr,mdatoms,
1584                                        inputrec->ePBC==epbcNONE ? NULL : &pbc);
1585         }
1586     }
1587
1588     if (NEED_MUTOT(*inputrec))
1589     {
1590
1591         if (bStateChanged)
1592         {
1593             if (PAR(cr))
1594             {
1595                 gmx_sumd(2*DIM,mu,cr);
1596             }
1597             for(i=0; i<2; i++)
1598             {
1599                 for(j=0;j<DIM;j++)
1600                 {
1601                     fr->mu_tot[i][j] = mu[i*DIM + j];
1602                 }
1603             }
1604         }
1605         if (fr->efep == efepNO)
1606         {
1607             copy_rvec(fr->mu_tot[0],mu_tot);
1608         }
1609         else
1610         {
1611             for(j=0; j<DIM; j++)
1612             {
1613                 mu_tot[j] =
1614                     (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] + lambda[efptCOUL]*fr->mu_tot[1][j];
1615             }
1616         }
1617     }
1618
1619     /* Reset energies */
1620     reset_enerdata(&(inputrec->opts),fr,bNS,enerd,MASTER(cr));
1621     clear_rvecs(SHIFTS,fr->fshift);
1622
1623     if (bNS)
1624     {
1625         wallcycle_start(wcycle,ewcNS);
1626
1627         if (graph && bStateChanged)
1628         {
1629             /* Calculate intramolecular shift vectors to make molecules whole */
1630             mk_mshift(fplog,graph,fr->ePBC,box,x);
1631         }
1632
1633         /* Do the actual neighbour searching and if twin range electrostatics
1634          * also do the calculation of long range forces and energies.
1635          */
1636         for (i=0;i<efptNR;i++) {dvdlambda[i] = 0;}
1637         ns(fplog,fr,x,box,
1638            groups,&(inputrec->opts),top,mdatoms,
1639            cr,nrnb,lambda,dvdlambda,&enerd->grpp,bFillGrid,
1640            bDoLongRangeNS);
1641         if (bSepDVDL)
1642         {
1643             fprintf(fplog,sepdvdlformat,"LR non-bonded",0.0,dvdlambda);
1644         }
1645         enerd->dvdl_lin[efptVDW] += dvdlambda[efptVDW];
1646         enerd->dvdl_lin[efptCOUL] += dvdlambda[efptCOUL];
1647
1648         wallcycle_stop(wcycle,ewcNS);
1649     }
1650
1651     if (inputrec->implicit_solvent && bNS)
1652     {
1653         make_gb_nblist(cr,inputrec->gb_algorithm,inputrec->rlist,
1654                        x,box,fr,&top->idef,graph,fr->born);
1655     }
1656
1657     if (DOMAINDECOMP(cr))
1658     {
1659         if (!(cr->duty & DUTY_PME))
1660         {
1661             wallcycle_start(wcycle,ewcPPDURINGPME);
1662             dd_force_flop_start(cr->dd,nrnb);
1663         }
1664     }
1665
1666     if (inputrec->bRot)
1667     {
1668         /* Enforced rotation has its own cycle counter that starts after the collective
1669          * coordinates have been communicated. It is added to ddCyclF to allow
1670          * for proper load-balancing */
1671         wallcycle_start(wcycle,ewcROT);
1672         do_rotation(cr,inputrec,box,x,t,step,wcycle,bNS);
1673         wallcycle_stop(wcycle,ewcROT);
1674     }
1675
1676     /* Start the force cycle counter.
1677      * This counter is stopped in do_forcelow_level.
1678      * No parallel communication should occur while this counter is running,
1679      * since that will interfere with the dynamic load balancing.
1680      */
1681     wallcycle_start(wcycle,ewcFORCE);
1682     
1683     if (bDoForces)
1684     {
1685         /* Reset forces for which the virial is calculated separately:
1686          * PME/Ewald forces if necessary */
1687         if (fr->bF_NoVirSum)
1688         {
1689             if (flags & GMX_FORCE_VIRIAL)
1690             {
1691                 fr->f_novirsum = fr->f_novirsum_alloc;
1692                 if (fr->bDomDec)
1693                 {
1694                     clear_rvecs(fr->f_novirsum_n,fr->f_novirsum);
1695                 }
1696                 else
1697                 {
1698                     clear_rvecs(homenr,fr->f_novirsum+start);
1699                 }
1700             }
1701             else
1702             {
1703                 /* We are not calculating the pressure so we do not need
1704                  * a separate array for forces that do not contribute
1705                  * to the pressure.
1706                  */
1707                 fr->f_novirsum = f;
1708             }
1709         }
1710
1711         /* Clear the short- and long-range forces */
1712         clear_rvecs(fr->natoms_force_constr,f);
1713         if(bSepLRF && do_per_step(step,inputrec->nstcalclr))
1714         {
1715             clear_rvecs(fr->natoms_force_constr,fr->f_twin);
1716         }
1717         
1718         clear_rvec(fr->vir_diag_posres);
1719     }
1720     if (inputrec->ePull == epullCONSTRAINT)
1721     {
1722         clear_pull_forces(inputrec->pull);
1723     }
1724
1725     /* update QMMMrec, if necessary */
1726     if(fr->bQMMM)
1727     {
1728         update_QMMMrec(cr,fr,x,mdatoms,box,top);
1729     }
1730
1731     if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1732     {
1733         posres_wrapper(fplog,flags,bSepDVDL,inputrec,nrnb,top,box,x,
1734                        f,enerd,lambda,fr);
1735     }
1736
1737     if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_FBPOSRES].nr > 0)
1738     {
1739         /* Flat-bottomed position restraints always require full pbc */
1740         if(!(bStateChanged && bDoAdressWF))
1741         {
1742             set_pbc(&pbc,inputrec->ePBC,box);
1743         }
1744         v = fbposres(top->idef.il[F_FBPOSRES].nr,top->idef.il[F_FBPOSRES].iatoms,
1745                      top->idef.iparams_fbposres,
1746                      (const rvec*)x,fr->f_novirsum,fr->vir_diag_posres,
1747                      inputrec->ePBC==epbcNONE ? NULL : &pbc,
1748                      fr->rc_scaling,fr->ePBC,fr->posres_com);
1749         enerd->term[F_FBPOSRES] += v;
1750         inc_nrnb(nrnb,eNR_FBPOSRES,top->idef.il[F_FBPOSRES].nr/2);
1751     }
1752
1753     /* Compute the bonded and non-bonded energies and optionally forces */
1754     do_force_lowlevel(fplog,step,fr,inputrec,&(top->idef),
1755                       cr,nrnb,wcycle,mdatoms,&(inputrec->opts),
1756                       x,hist,f, bSepLRF ? fr->f_twin : f,enerd,fcd,mtop,top,fr->born,
1757                       &(top->atomtypes),bBornRadii,box,
1758                       inputrec->fepvals,lambda,
1759                       graph,&(top->excls),fr->mu_tot,
1760                       flags,
1761                       &cycles_pme);
1762
1763     if(bSepLRF)
1764     {
1765         if (do_per_step(step,inputrec->nstcalclr))
1766         {
1767             /* Add the long range forces to the short range forces */
1768             for(i=0; i<fr->natoms_force_constr; i++)
1769             {
1770                 rvec_add(fr->f_twin[i],f[i],f[i]);
1771             }
1772         }
1773     }
1774     
1775     cycles_force = wallcycle_stop(wcycle,ewcFORCE);
1776
1777     if (ed)
1778     {
1779         do_flood(cr,inputrec,x,f,ed,box,step,bNS);
1780     }
1781
1782     if (DOMAINDECOMP(cr))
1783     {
1784         dd_force_flop_stop(cr->dd,nrnb);
1785         if (wcycle)
1786         {
1787             dd_cycles_add(cr->dd,cycles_force-cycles_pme,ddCyclF);
1788         }
1789     }
1790
1791     if (bDoForces)
1792     {
1793         if (IR_ELEC_FIELD(*inputrec))
1794         {
1795             /* Compute forces due to electric field */
1796             calc_f_el(MASTER(cr) ? field : NULL,
1797                       start,homenr,mdatoms->chargeA,x,fr->f_novirsum,
1798                       inputrec->ex,inputrec->et,t);
1799         }
1800
1801         if (bDoAdressWF && fr->adress_icor == eAdressICThermoForce)
1802         {
1803             /* Compute thermodynamic force in hybrid AdResS region */
1804             adress_thermo_force(start,homenr,&(top->cgs),x,fr->f_novirsum,fr,mdatoms,
1805                                 inputrec->ePBC==epbcNONE ? NULL : &pbc);
1806         }
1807
1808         /* Communicate the forces */
1809         if (PAR(cr))
1810         {
1811             wallcycle_start(wcycle,ewcMOVEF);
1812             if (DOMAINDECOMP(cr))
1813             {
1814                 dd_move_f(cr->dd,f,fr->fshift);
1815                 /* Do we need to communicate the separate force array
1816                  * for terms that do not contribute to the single sum virial?
1817                  * Position restraints and electric fields do not introduce
1818                  * inter-cg forces, only full electrostatics methods do.
1819                  * When we do not calculate the virial, fr->f_novirsum = f,
1820                  * so we have already communicated these forces.
1821                  */
1822                 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1823                     (flags & GMX_FORCE_VIRIAL))
1824                 {
1825                     dd_move_f(cr->dd,fr->f_novirsum,NULL);
1826                 }
1827                 if (bSepLRF)
1828                 {
1829                     /* We should not update the shift forces here,
1830                      * since f_twin is already included in f.
1831                      */
1832                     dd_move_f(cr->dd,fr->f_twin,NULL);
1833                 }
1834             }
1835             else
1836             {
1837                 pd_move_f(cr,f,nrnb);
1838                 if (bSepLRF)
1839                 {
1840                     pd_move_f(cr,fr->f_twin,nrnb);
1841                 }
1842             }
1843             wallcycle_stop(wcycle,ewcMOVEF);
1844         }
1845
1846         /* If we have NoVirSum forces, but we do not calculate the virial,
1847          * we sum fr->f_novirum=f later.
1848          */
1849         if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1850         {
1851             wallcycle_start(wcycle,ewcVSITESPREAD);
1852             spread_vsite_f(fplog,vsite,x,f,fr->fshift,FALSE,NULL,nrnb,
1853                            &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
1854             wallcycle_stop(wcycle,ewcVSITESPREAD);
1855
1856             if (bSepLRF)
1857             {
1858                 wallcycle_start(wcycle,ewcVSITESPREAD);
1859                 spread_vsite_f(fplog,vsite,x,fr->f_twin,NULL,FALSE,NULL,
1860                                nrnb,
1861                                &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
1862                 wallcycle_stop(wcycle,ewcVSITESPREAD);
1863             }
1864         }
1865
1866         if (flags & GMX_FORCE_VIRIAL)
1867         {
1868             /* Calculation of the virial must be done after vsites! */
1869             calc_virial(fplog,mdatoms->start,mdatoms->homenr,x,f,
1870                         vir_force,graph,box,nrnb,fr,inputrec->ePBC);
1871         }
1872     }
1873
1874     if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1875     {
1876         pull_potential_wrapper(fplog,bSepDVDL,cr,inputrec,box,x,
1877                                f,vir_force,mdatoms,enerd,lambda,t);
1878     }
1879
1880     /* Add the forces from enforced rotation potentials (if any) */
1881     if (inputrec->bRot)
1882     {
1883         wallcycle_start(wcycle,ewcROTadd);
1884         enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr,step,t);
1885         wallcycle_stop(wcycle,ewcROTadd);
1886     }
1887
1888     if (PAR(cr) && !(cr->duty & DUTY_PME))
1889     {
1890         /* In case of node-splitting, the PP nodes receive the long-range 
1891          * forces, virial and energy from the PME nodes here.
1892          */
1893         pme_receive_force_ener(fplog,bSepDVDL,cr,wcycle,enerd,fr);
1894     }
1895
1896     if (bDoForces)
1897     {
1898         post_process_forces(fplog,cr,step,nrnb,wcycle,
1899                             top,box,x,f,vir_force,mdatoms,graph,fr,vsite,
1900                             flags);
1901     }
1902
1903     /* Sum the potential energy terms from group contributions */
1904     sum_epot(&(inputrec->opts),&(enerd->grpp),enerd->term);
1905 }
1906
1907 void do_force(FILE *fplog,t_commrec *cr,
1908               t_inputrec *inputrec,
1909               gmx_large_int_t step,t_nrnb *nrnb,gmx_wallcycle_t wcycle,
1910               gmx_localtop_t *top,
1911               gmx_mtop_t *mtop,
1912               gmx_groups_t *groups,
1913               matrix box,rvec x[],history_t *hist,
1914               rvec f[],
1915               tensor vir_force,
1916               t_mdatoms *mdatoms,
1917               gmx_enerdata_t *enerd,t_fcdata *fcd,
1918               real *lambda,t_graph *graph,
1919               t_forcerec *fr,
1920               gmx_vsite_t *vsite,rvec mu_tot,
1921               double t,FILE *field,gmx_edsam_t ed,
1922               gmx_bool bBornRadii,
1923               int flags)
1924 {
1925     /* modify force flag if not doing nonbonded */
1926     if (!fr->bNonbonded)
1927     {
1928         flags &= ~GMX_FORCE_NONBONDED;
1929     }
1930
1931     switch (inputrec->cutoff_scheme)
1932     {
1933         case ecutsVERLET:
1934             do_force_cutsVERLET(fplog, cr, inputrec,
1935                                 step, nrnb, wcycle,
1936                                 top, mtop,
1937                                 groups,
1938                                 box, x, hist,
1939                                 f, vir_force,
1940                                 mdatoms,
1941                                 enerd, fcd,
1942                                 lambda, graph,
1943                                 fr, fr->ic, 
1944                                 vsite, mu_tot,
1945                                 t, field, ed,
1946                                 bBornRadii,
1947                                 flags);
1948             break;
1949         case ecutsGROUP:
1950              do_force_cutsGROUP(fplog, cr, inputrec,
1951                                 step, nrnb, wcycle,
1952                                 top, mtop,
1953                                 groups,
1954                                 box, x, hist,
1955                                 f, vir_force,
1956                                 mdatoms,
1957                                 enerd, fcd,
1958                                 lambda, graph,
1959                                 fr, vsite, mu_tot,
1960                                 t, field, ed,
1961                                 bBornRadii,
1962                                 flags);
1963             break;
1964         default:
1965             gmx_incons("Invalid cut-off scheme passed!");
1966     }
1967 }
1968
1969
1970 void do_constrain_first(FILE *fplog,gmx_constr_t constr,
1971                         t_inputrec *ir,t_mdatoms *md,
1972                         t_state *state,rvec *f,
1973                         t_graph *graph,t_commrec *cr,t_nrnb *nrnb,
1974                         t_forcerec *fr, gmx_localtop_t *top, tensor shake_vir)
1975 {
1976     int    i,m,start,end;
1977     gmx_large_int_t step;
1978     real   dt=ir->delta_t;
1979     real   dvdl_dum;
1980     rvec   *savex;
1981
1982     snew(savex,state->natoms);
1983
1984     start = md->start;
1985     end   = md->homenr + start;
1986
1987     if (debug)
1988         fprintf(debug,"vcm: start=%d, homenr=%d, end=%d\n",
1989                 start,md->homenr,end);
1990     /* Do a first constrain to reset particles... */
1991     step = ir->init_step;
1992     if (fplog)
1993     {
1994         char buf[STEPSTRSIZE];
1995         fprintf(fplog,"\nConstraining the starting coordinates (step %s)\n",
1996                 gmx_step_str(step,buf));
1997     }
1998     dvdl_dum = 0;
1999
2000     /* constrain the current position */
2001     constrain(NULL,TRUE,FALSE,constr,&(top->idef),
2002               ir,NULL,cr,step,0,md,
2003               state->x,state->x,NULL,
2004               fr->bMolPBC,state->box,
2005               state->lambda[efptBONDED],&dvdl_dum,
2006               NULL,NULL,nrnb,econqCoord,
2007               ir->epc==epcMTTK,state->veta,state->veta);
2008     if (EI_VV(ir->eI))
2009     {
2010         /* constrain the inital velocity, and save it */
2011         /* also may be useful if we need the ekin from the halfstep for velocity verlet */
2012         /* might not yet treat veta correctly */
2013         constrain(NULL,TRUE,FALSE,constr,&(top->idef),
2014                   ir,NULL,cr,step,0,md,
2015                   state->x,state->v,state->v,
2016                   fr->bMolPBC,state->box,
2017                   state->lambda[efptBONDED],&dvdl_dum,
2018                   NULL,NULL,nrnb,econqVeloc,
2019                   ir->epc==epcMTTK,state->veta,state->veta);
2020     }
2021     /* constrain the inital velocities at t-dt/2 */
2022     if (EI_STATE_VELOCITY(ir->eI) && ir->eI!=eiVV)
2023     {
2024         for(i=start; (i<end); i++)
2025         {
2026             for(m=0; (m<DIM); m++)
2027             {
2028                 /* Reverse the velocity */
2029                 state->v[i][m] = -state->v[i][m];
2030                 /* Store the position at t-dt in buf */
2031                 savex[i][m] = state->x[i][m] + dt*state->v[i][m];
2032             }
2033         }
2034     /* Shake the positions at t=-dt with the positions at t=0
2035      * as reference coordinates.
2036          */
2037         if (fplog)
2038         {
2039             char buf[STEPSTRSIZE];
2040             fprintf(fplog,"\nConstraining the coordinates at t0-dt (step %s)\n",
2041                     gmx_step_str(step,buf));
2042         }
2043         dvdl_dum = 0;
2044         constrain(NULL,TRUE,FALSE,constr,&(top->idef),
2045                   ir,NULL,cr,step,-1,md,
2046                   state->x,savex,NULL,
2047                   fr->bMolPBC,state->box,
2048                   state->lambda[efptBONDED],&dvdl_dum,
2049                   state->v,NULL,nrnb,econqCoord,
2050                   ir->epc==epcMTTK,state->veta,state->veta);
2051         
2052         for(i=start; i<end; i++) {
2053             for(m=0; m<DIM; m++) {
2054                 /* Re-reverse the velocities */
2055                 state->v[i][m] = -state->v[i][m];
2056             }
2057         }
2058     }
2059     sfree(savex);
2060 }
2061
2062 void calc_enervirdiff(FILE *fplog,int eDispCorr,t_forcerec *fr)
2063 {
2064   double eners[2],virs[2],enersum,virsum,y0,f,g,h;
2065   double r0,r1,r,rc3,rc9,ea,eb,ec,pa,pb,pc,pd;
2066   double invscale,invscale2,invscale3;
2067   int    ri0,ri1,ri,i,offstart,offset;
2068   real   scale,*vdwtab,tabfactor,tmp;
2069
2070   fr->enershiftsix = 0;
2071   fr->enershifttwelve = 0;
2072   fr->enerdiffsix = 0;
2073   fr->enerdifftwelve = 0;
2074   fr->virdiffsix = 0;
2075   fr->virdifftwelve = 0;
2076
2077   if (eDispCorr != edispcNO) {
2078     for(i=0; i<2; i++) {
2079       eners[i] = 0;
2080       virs[i]  = 0;
2081     }
2082     if ((fr->vdwtype == evdwSWITCH) || (fr->vdwtype == evdwSHIFT)) {
2083       if (fr->rvdw_switch == 0)
2084         gmx_fatal(FARGS,
2085                   "With dispersion correction rvdw-switch can not be zero "
2086                   "for vdw-type = %s",evdw_names[fr->vdwtype]);
2087
2088       scale  = fr->nblists[0].table_elec_vdw.scale;
2089       vdwtab = fr->nblists[0].table_vdw.data;
2090
2091       /* Round the cut-offs to exact table values for precision */
2092       ri0 = floor(fr->rvdw_switch*scale);
2093       ri1 = ceil(fr->rvdw*scale);
2094       r0  = ri0/scale;
2095       r1  = ri1/scale;
2096       rc3 = r0*r0*r0;
2097       rc9  = rc3*rc3*rc3;
2098
2099       if (fr->vdwtype == evdwSHIFT)
2100       {
2101           /* Determine the constant energy shift below rvdw_switch.
2102            * Table has a scale factor since we have scaled it down to compensate
2103            * for scaling-up c6/c12 with the derivative factors to save flops in analytical kernels.
2104            */
2105           fr->enershiftsix    = (real)(-1.0/(rc3*rc3)) - 6.0*vdwtab[8*ri0];
2106           fr->enershifttwelve = (real)( 1.0/(rc9*rc3)) - 12.0*vdwtab[8*ri0 + 4];
2107       }
2108       /* Add the constant part from 0 to rvdw_switch.
2109        * This integration from 0 to rvdw_switch overcounts the number
2110        * of interactions by 1, as it also counts the self interaction.
2111        * We will correct for this later.
2112        */
2113       eners[0] += 4.0*M_PI*fr->enershiftsix*rc3/3.0;
2114       eners[1] += 4.0*M_PI*fr->enershifttwelve*rc3/3.0;
2115
2116       invscale = 1.0/(scale);
2117       invscale2 = invscale*invscale;
2118       invscale3 = invscale*invscale2;
2119
2120       /* following summation derived from cubic spline definition,
2121         Numerical Recipies in C, second edition, p. 113-116.  Exact
2122         for the cubic spline.  We first calculate the negative of
2123         the energy from rvdw to rvdw_switch, assuming that g(r)=1,
2124         and then add the more standard, abrupt cutoff correction to
2125         that result, yielding the long-range correction for a
2126         switched function.  We perform both the pressure and energy
2127         loops at the same time for simplicity, as the computational
2128         cost is low. */
2129
2130       for (i=0;i<2;i++) {
2131         enersum = 0.0; virsum = 0.0;
2132         if (i==0)
2133         {
2134             offstart = 0;
2135             /* Since the dispersion table has been scaled down a factor 6.0 and the repulsion
2136              * a factor 12.0 to compensate for the c6/c12 parameters inside nbfp[] being scaled
2137              * up (to save flops in kernels), we need to correct for this.
2138              */
2139             tabfactor = 6.0;
2140         }
2141         else
2142         {
2143             offstart = 4;
2144             tabfactor = 12.0;
2145         }
2146         for (ri=ri0; ri<ri1; ri++) {
2147           r = ri*invscale;
2148           ea = invscale3;
2149           eb = 2.0*invscale2*r;
2150           ec = invscale*r*r;
2151
2152           pa = invscale3;
2153           pb = 3.0*invscale2*r;
2154           pc = 3.0*invscale*r*r;
2155           pd = r*r*r;
2156
2157           /* this "8" is from the packing in the vdwtab array - perhaps should be #define'ed? */
2158           offset = 8*ri + offstart;
2159           y0 = vdwtab[offset];
2160           f  = vdwtab[offset+1];
2161           g  = vdwtab[offset+2];
2162           h  = vdwtab[offset+3];
2163
2164           enersum += y0*(ea/3 + eb/2 + ec) + f*(ea/4 + eb/3 + ec/2) + g*(ea/5 + eb/4 + ec/3) + h*(ea/6 + eb/5 + ec/4);
2165           virsum  += f*(pa/4 + pb/3 + pc/2 + pd) + 2*g*(pa/5 + pb/4 + pc/3 + pd/2) + 3*h*(pa/6 + pb/5 + pc/4 + pd/3);
2166         }
2167           
2168         enersum *= 4.0*M_PI*tabfactor;
2169         virsum  *= 4.0*M_PI*tabfactor;
2170         eners[i] -= enersum;
2171         virs[i]  -= virsum;
2172       }
2173
2174       /* now add the correction for rvdw_switch to infinity */
2175       eners[0] += -4.0*M_PI/(3.0*rc3);
2176       eners[1] +=  4.0*M_PI/(9.0*rc9);
2177       virs[0]  +=  8.0*M_PI/rc3;
2178       virs[1]  += -16.0*M_PI/(3.0*rc9);
2179     }
2180     else if ((fr->vdwtype == evdwCUT) || (fr->vdwtype == evdwUSER)) {
2181       if (fr->vdwtype == evdwUSER && fplog)
2182         fprintf(fplog,
2183                 "WARNING: using dispersion correction with user tables\n");
2184       rc3  = fr->rvdw*fr->rvdw*fr->rvdw;
2185       rc9  = rc3*rc3*rc3;
2186       /* Contribution beyond the cut-off */
2187       eners[0] += -4.0*M_PI/(3.0*rc3);
2188       eners[1] +=  4.0*M_PI/(9.0*rc9);
2189       if (fr->vdw_modifier==eintmodPOTSHIFT) {
2190           /* Contribution within the cut-off */
2191           eners[0] += -4.0*M_PI/(3.0*rc3);
2192           eners[1] +=  4.0*M_PI/(3.0*rc9);
2193       }
2194       /* Contribution beyond the cut-off */
2195       virs[0]  +=  8.0*M_PI/rc3;
2196       virs[1]  += -16.0*M_PI/(3.0*rc9);
2197     } else {
2198       gmx_fatal(FARGS,
2199                 "Dispersion correction is not implemented for vdw-type = %s",
2200                 evdw_names[fr->vdwtype]);
2201     }
2202     fr->enerdiffsix    = eners[0];
2203     fr->enerdifftwelve = eners[1];
2204     /* The 0.5 is due to the Gromacs definition of the virial */
2205     fr->virdiffsix     = 0.5*virs[0];
2206     fr->virdifftwelve  = 0.5*virs[1];
2207   }
2208 }
2209
2210 void calc_dispcorr(FILE *fplog,t_inputrec *ir,t_forcerec *fr,
2211                    gmx_large_int_t step,int natoms,
2212                    matrix box,real lambda,tensor pres,tensor virial,
2213                    real *prescorr, real *enercorr, real *dvdlcorr)
2214 {
2215     gmx_bool bCorrAll,bCorrPres;
2216     real dvdlambda,invvol,dens,ninter,avcsix,avctwelve,enerdiff,svir=0,spres=0;
2217     int  m;
2218
2219     *prescorr = 0;
2220     *enercorr = 0;
2221     *dvdlcorr = 0;
2222
2223     clear_mat(virial);
2224     clear_mat(pres);
2225
2226     if (ir->eDispCorr != edispcNO) {
2227         bCorrAll  = (ir->eDispCorr == edispcAllEner ||
2228                      ir->eDispCorr == edispcAllEnerPres);
2229         bCorrPres = (ir->eDispCorr == edispcEnerPres ||
2230                      ir->eDispCorr == edispcAllEnerPres);
2231
2232         invvol = 1/det(box);
2233         if (fr->n_tpi)
2234         {
2235             /* Only correct for the interactions with the inserted molecule */
2236             dens = (natoms - fr->n_tpi)*invvol;
2237             ninter = fr->n_tpi;
2238         }
2239         else
2240         {
2241             dens = natoms*invvol;
2242             ninter = 0.5*natoms;
2243         }
2244
2245         if (ir->efep == efepNO)
2246         {
2247             avcsix    = fr->avcsix[0];
2248             avctwelve = fr->avctwelve[0];
2249         }
2250         else
2251         {
2252             avcsix    = (1 - lambda)*fr->avcsix[0]    + lambda*fr->avcsix[1];
2253             avctwelve = (1 - lambda)*fr->avctwelve[0] + lambda*fr->avctwelve[1];
2254         }
2255
2256         enerdiff = ninter*(dens*fr->enerdiffsix - fr->enershiftsix);
2257         *enercorr += avcsix*enerdiff;
2258         dvdlambda = 0.0;
2259         if (ir->efep != efepNO)
2260         {
2261             dvdlambda += (fr->avcsix[1] - fr->avcsix[0])*enerdiff;
2262         }
2263         if (bCorrAll)
2264         {
2265             enerdiff = ninter*(dens*fr->enerdifftwelve - fr->enershifttwelve);
2266             *enercorr += avctwelve*enerdiff;
2267             if (fr->efep != efepNO)
2268             {
2269                 dvdlambda += (fr->avctwelve[1] - fr->avctwelve[0])*enerdiff;
2270             }
2271         }
2272
2273         if (bCorrPres)
2274         {
2275             svir = ninter*dens*avcsix*fr->virdiffsix/3.0;
2276             if (ir->eDispCorr == edispcAllEnerPres)
2277             {
2278                 svir += ninter*dens*avctwelve*fr->virdifftwelve/3.0;
2279             }
2280             /* The factor 2 is because of the Gromacs virial definition */
2281             spres = -2.0*invvol*svir*PRESFAC;
2282
2283             for(m=0; m<DIM; m++) {
2284                 virial[m][m] += svir;
2285                 pres[m][m] += spres;
2286             }
2287             *prescorr += spres;
2288         }
2289
2290         /* Can't currently control when it prints, for now, just print when degugging */
2291         if (debug)
2292         {
2293             if (bCorrAll) {
2294                 fprintf(debug,"Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
2295                         avcsix,avctwelve);
2296             }
2297             if (bCorrPres)
2298             {
2299                 fprintf(debug,
2300                         "Long Range LJ corr.: Epot %10g, Pres: %10g, Vir: %10g\n",
2301                         *enercorr,spres,svir);
2302             }
2303             else
2304             {
2305                 fprintf(debug,"Long Range LJ corr.: Epot %10g\n",*enercorr);
2306             }
2307         }
2308
2309         if (fr->bSepDVDL && do_per_step(step,ir->nstlog))
2310         {
2311             fprintf(fplog,sepdvdlformat,"Dispersion correction",
2312                     *enercorr,dvdlambda);
2313         }
2314         if (fr->efep != efepNO)
2315         {
2316             *dvdlcorr += dvdlambda;
2317         }
2318     }
2319 }
2320
2321 void do_pbc_first(FILE *fplog,matrix box,t_forcerec *fr,
2322                   t_graph *graph,rvec x[])
2323 {
2324   if (fplog)
2325     fprintf(fplog,"Removing pbc first time\n");
2326   calc_shifts(box,fr->shift_vec);
2327   if (graph) {
2328     mk_mshift(fplog,graph,fr->ePBC,box,x);
2329     if (gmx_debug_at)
2330       p_graph(debug,"do_pbc_first 1",graph);
2331     shift_self(graph,box,x);
2332     /* By doing an extra mk_mshift the molecules that are broken
2333      * because they were e.g. imported from another software
2334      * will be made whole again. Such are the healing powers
2335      * of GROMACS.
2336      */
2337     mk_mshift(fplog,graph,fr->ePBC,box,x);
2338     if (gmx_debug_at)
2339       p_graph(debug,"do_pbc_first 2",graph);
2340   }
2341   if (fplog)
2342     fprintf(fplog,"Done rmpbc\n");
2343 }
2344
2345 static void low_do_pbc_mtop(FILE *fplog,int ePBC,matrix box,
2346                             gmx_mtop_t *mtop,rvec x[],
2347                             gmx_bool bFirst)
2348 {
2349   t_graph *graph;
2350   int mb,as,mol;
2351   gmx_molblock_t *molb;
2352
2353   if (bFirst && fplog)
2354     fprintf(fplog,"Removing pbc first time\n");
2355
2356   snew(graph,1);
2357   as = 0;
2358   for(mb=0; mb<mtop->nmolblock; mb++) {
2359     molb = &mtop->molblock[mb];
2360     if (molb->natoms_mol == 1 ||
2361         (!bFirst && mtop->moltype[molb->type].cgs.nr == 1)) {
2362       /* Just one atom or charge group in the molecule, no PBC required */
2363       as += molb->nmol*molb->natoms_mol;
2364     } else {
2365       /* Pass NULL iso fplog to avoid graph prints for each molecule type */
2366       mk_graph_ilist(NULL,mtop->moltype[molb->type].ilist,
2367                      0,molb->natoms_mol,FALSE,FALSE,graph);
2368
2369       for(mol=0; mol<molb->nmol; mol++) {
2370         mk_mshift(fplog,graph,ePBC,box,x+as);
2371
2372         shift_self(graph,box,x+as);
2373         /* The molecule is whole now.
2374          * We don't need the second mk_mshift call as in do_pbc_first,
2375          * since we no longer need this graph.
2376          */
2377
2378         as += molb->natoms_mol;
2379       }
2380       done_graph(graph);
2381     }
2382   }
2383   sfree(graph);
2384 }
2385
2386 void do_pbc_first_mtop(FILE *fplog,int ePBC,matrix box,
2387                        gmx_mtop_t *mtop,rvec x[])
2388 {
2389   low_do_pbc_mtop(fplog,ePBC,box,mtop,x,TRUE);
2390 }
2391
2392 void do_pbc_mtop(FILE *fplog,int ePBC,matrix box,
2393                  gmx_mtop_t *mtop,rvec x[])
2394 {
2395   low_do_pbc_mtop(fplog,ePBC,box,mtop,x,FALSE);
2396 }
2397
2398 void finish_run(FILE *fplog,t_commrec *cr,const char *confout,
2399                 t_inputrec *inputrec,
2400                 t_nrnb nrnb[],gmx_wallcycle_t wcycle,
2401                 gmx_runtime_t *runtime,
2402                 wallclock_gpu_t *gputimes,
2403                 int omp_nth_pp,
2404                 gmx_bool bWriteStat)
2405 {
2406     int    i,j;
2407     t_nrnb *nrnb_tot=NULL;
2408     real   delta_t;
2409     double nbfs,mflop;
2410
2411     wallcycle_sum(cr,wcycle);
2412
2413     if (cr->nnodes > 1)
2414     {
2415         snew(nrnb_tot,1);
2416 #ifdef GMX_MPI
2417         MPI_Allreduce(nrnb->n,nrnb_tot->n,eNRNB,MPI_DOUBLE,MPI_SUM,
2418                       cr->mpi_comm_mysim);
2419 #endif
2420     }
2421     else
2422     {
2423         nrnb_tot = nrnb;
2424     }
2425
2426 #if defined(GMX_MPI) && !defined(GMX_THREAD_MPI)
2427     if (cr->nnodes > 1)
2428     {
2429         /* reduce nodetime over all MPI processes in the current simulation */
2430         double sum;
2431         MPI_Allreduce(&runtime->proctime,&sum,1,MPI_DOUBLE,MPI_SUM,
2432                       cr->mpi_comm_mysim);
2433         runtime->proctime = sum;
2434     }
2435 #endif
2436
2437     if (SIMMASTER(cr))
2438     {
2439         print_flop(fplog,nrnb_tot,&nbfs,&mflop);
2440     }
2441     if (cr->nnodes > 1)
2442     {
2443         sfree(nrnb_tot);
2444     }
2445
2446     if ((cr->duty & DUTY_PP) && DOMAINDECOMP(cr))
2447     {
2448         print_dd_statistics(cr,inputrec,fplog);
2449     }
2450
2451 #ifdef GMX_MPI
2452     if (PARTDECOMP(cr))
2453     {
2454         if (MASTER(cr))
2455         {
2456             t_nrnb     *nrnb_all;
2457             int        s;
2458             MPI_Status stat;
2459
2460             snew(nrnb_all,cr->nnodes);
2461             nrnb_all[0] = *nrnb;
2462             for(s=1; s<cr->nnodes; s++)
2463             {
2464                 MPI_Recv(nrnb_all[s].n,eNRNB,MPI_DOUBLE,s,0,
2465                          cr->mpi_comm_mysim,&stat);
2466             }
2467             pr_load(fplog,cr,nrnb_all);
2468             sfree(nrnb_all);
2469         }
2470         else
2471         {
2472             MPI_Send(nrnb->n,eNRNB,MPI_DOUBLE,MASTERRANK(cr),0,
2473                      cr->mpi_comm_mysim);
2474         }
2475     }
2476 #endif
2477
2478     if (SIMMASTER(cr))
2479     {
2480         wallcycle_print(fplog,cr->nnodes,cr->npmenodes,runtime->realtime,
2481                         wcycle,gputimes);
2482
2483         if (EI_DYNAMICS(inputrec->eI))
2484         {
2485             delta_t = inputrec->delta_t;
2486         }
2487         else
2488         {
2489             delta_t = 0;
2490         }
2491
2492         if (fplog)
2493         {
2494             print_perf(fplog,runtime->proctime,runtime->realtime,
2495                        cr->nnodes-cr->npmenodes,
2496                        runtime->nsteps_done,delta_t,nbfs,mflop,
2497                        omp_nth_pp);
2498         }
2499         if (bWriteStat)
2500         {
2501             print_perf(stderr,runtime->proctime,runtime->realtime,
2502                        cr->nnodes-cr->npmenodes,
2503                        runtime->nsteps_done,delta_t,nbfs,mflop,
2504                        omp_nth_pp);
2505         }
2506     }
2507 }
2508
2509 extern void initialize_lambdas(FILE *fplog,t_inputrec *ir,int *fep_state,real *lambda,double *lam0)
2510 {
2511     /* this function works, but could probably use a logic rewrite to keep all the different
2512        types of efep straight. */
2513
2514     int i;
2515     t_lambda *fep = ir->fepvals;
2516
2517     if ((ir->efep==efepNO) && (ir->bSimTemp == FALSE)) {
2518         for (i=0;i<efptNR;i++)  {
2519             lambda[i] = 0.0;
2520             if (lam0)
2521             {
2522                 lam0[i] = 0.0;
2523             }
2524         }
2525         return;
2526     } else {
2527         *fep_state = fep->init_fep_state; /* this might overwrite the checkpoint
2528                                              if checkpoint is set -- a kludge is in for now
2529                                              to prevent this.*/
2530         for (i=0;i<efptNR;i++)
2531         {
2532             /* overwrite lambda state with init_lambda for now for backwards compatibility */
2533             if (fep->init_lambda>=0) /* if it's -1, it was never initializd */
2534             {
2535                 lambda[i] = fep->init_lambda;
2536                 if (lam0) {
2537                     lam0[i] = lambda[i];
2538                 }
2539             }
2540             else
2541             {
2542                 lambda[i] = fep->all_lambda[i][*fep_state];
2543                 if (lam0) {
2544                     lam0[i] = lambda[i];
2545                 }
2546             }
2547         }
2548         if (ir->bSimTemp) {
2549             /* need to rescale control temperatures to match current state */
2550             for (i=0;i<ir->opts.ngtc;i++) {
2551                 if (ir->opts.ref_t[i] > 0) {
2552                     ir->opts.ref_t[i] = ir->simtempvals->temperatures[*fep_state];
2553                 }
2554             }
2555         }
2556     }
2557
2558     /* Send to the log the information on the current lambdas */
2559     if (fplog != NULL)
2560     {
2561         fprintf(fplog,"Initial vector of lambda components:[ ");
2562         for (i=0;i<efptNR;i++)
2563         {
2564             fprintf(fplog,"%10.4f ",lambda[i]);
2565         }
2566         fprintf(fplog,"]\n");
2567     }
2568     return;
2569 }
2570
2571
2572 void init_md(FILE *fplog,
2573              t_commrec *cr,t_inputrec *ir,const output_env_t oenv,
2574              double *t,double *t0,
2575              real *lambda, int *fep_state, double *lam0,
2576              t_nrnb *nrnb,gmx_mtop_t *mtop,
2577              gmx_update_t *upd,
2578              int nfile,const t_filenm fnm[],
2579              gmx_mdoutf_t **outf,t_mdebin **mdebin,
2580              tensor force_vir,tensor shake_vir,rvec mu_tot,
2581              gmx_bool *bSimAnn,t_vcm **vcm, t_state *state, unsigned long Flags)
2582 {
2583     int  i,j,n;
2584     real tmpt,mod;
2585
2586     /* Initial values */
2587     *t = *t0       = ir->init_t;
2588
2589     *bSimAnn=FALSE;
2590     for(i=0;i<ir->opts.ngtc;i++)
2591     {
2592         /* set bSimAnn if any group is being annealed */
2593         if(ir->opts.annealing[i]!=eannNO)
2594         {
2595             *bSimAnn = TRUE;
2596         }
2597     }
2598     if (*bSimAnn)
2599     {
2600         update_annealing_target_temp(&(ir->opts),ir->init_t);
2601     }
2602
2603     /* Initialize lambda variables */
2604     initialize_lambdas(fplog,ir,fep_state,lambda,lam0);
2605
2606     if (upd)
2607     {
2608         *upd = init_update(fplog,ir);
2609     }
2610
2611
2612     if (vcm != NULL)
2613     {
2614         *vcm = init_vcm(fplog,&mtop->groups,ir);
2615     }
2616
2617     if (EI_DYNAMICS(ir->eI) && !(Flags & MD_APPENDFILES))
2618     {
2619         if (ir->etc == etcBERENDSEN)
2620         {
2621             please_cite(fplog,"Berendsen84a");
2622         }
2623         if (ir->etc == etcVRESCALE)
2624         {
2625             please_cite(fplog,"Bussi2007a");
2626         }
2627     }
2628
2629     init_nrnb(nrnb);
2630
2631     if (nfile != -1)
2632     {
2633         *outf = init_mdoutf(nfile,fnm,Flags,cr,ir,oenv);
2634
2635         *mdebin = init_mdebin((Flags & MD_APPENDFILES) ? NULL : (*outf)->fp_ene,
2636                               mtop,ir, (*outf)->fp_dhdl);
2637     }
2638
2639     if (ir->bAdress)
2640     {
2641       please_cite(fplog,"Fritsch12");
2642       please_cite(fplog,"Junghans10");
2643     }
2644     /* Initiate variables */
2645     clear_mat(force_vir);
2646     clear_mat(shake_vir);
2647     clear_rvec(mu_tot);
2648
2649     debug_gmx();
2650 }
2651