1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
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.
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.
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.
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.
31 * For more info, check our website at http://www.gromacs.org
34 * GROwing Monsters And Cloning Shrimps
41 #include<catamount/dclock.h>
47 #ifdef HAVE_SYS_TIME_H
60 #include "chargegroup.h"
83 #include "pull_rotation.h"
84 #include "gmx_random.h"
87 #include "gmx_wallcycle.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"
106 #include "nbnxn_cuda_data_mgmt.h"
107 #include "nbnxn_cuda/nbnxn_cuda.h"
110 typedef struct gmx_timeprint {
115 /* Portable version of ctime_r implemented in src/gmxlib/string2.c, but we do not want it declared in public installed headers */
117 gmx_ctime_r(const time_t *clock,char *buf, int n);
123 #ifdef HAVE_GETTIMEOFDAY
127 gettimeofday(&t,NULL);
129 seconds = (double) t.tv_sec + 1e-6*(double)t.tv_usec;
135 seconds = time(NULL);
142 #define difftime(end,start) ((double)(end)-(double)(start))
144 void print_time(FILE *out,gmx_runtime_t *runtime,gmx_large_int_t step,
145 t_inputrec *ir, t_commrec *cr)
148 char timebuf[STRLEN];
152 #ifndef GMX_THREAD_MPI
158 fprintf(out,"step %s",gmx_step_str(step,buf));
159 if ((step >= ir->nstlist))
161 runtime->last = gmx_gettime();
162 dt = difftime(runtime->last,runtime->real);
163 runtime->time_per_step = dt/(step - ir->init_step + 1);
165 dt = (ir->nsteps + ir->init_step - step)*runtime->time_per_step;
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);
178 fprintf(out,", remaining runtime: %5d s ",(int)dt);
182 fprintf(out," performance: %.1f ns/day ",
183 ir->delta_t/1000*24*60*60/runtime->time_per_step);
186 #ifndef GMX_THREAD_MPI
200 static double set_proctime(gmx_runtime_t *runtime)
206 prev = runtime->proc;
207 runtime->proc = dclock();
209 diff = runtime->proc - prev;
213 prev = runtime->proc;
214 runtime->proc = clock();
216 diff = (double)(runtime->proc - prev)/(double)CLOCKS_PER_SEC;
220 /* The counter has probably looped, ignore this data */
227 void runtime_start(gmx_runtime_t *runtime)
229 runtime->real = gmx_gettime();
231 set_proctime(runtime);
232 runtime->realtime = 0;
233 runtime->proctime = 0;
235 runtime->time_per_step = 0;
238 void runtime_end(gmx_runtime_t *runtime)
244 runtime->proctime += set_proctime(runtime);
245 runtime->realtime = now - runtime->real;
249 void runtime_upd_proc(gmx_runtime_t *runtime)
251 runtime->proctime += set_proctime(runtime);
254 void print_date_and_time(FILE *fplog,int nodeid,const char *title,
255 const gmx_runtime_t *runtime)
258 char timebuf[STRLEN];
259 char time_string[STRLEN];
266 tmptime = (time_t) runtime->real;
267 gmx_ctime_r(&tmptime,timebuf,STRLEN);
271 tmptime = (time_t) gmx_gettime();
272 gmx_ctime_r(&tmptime,timebuf,STRLEN);
274 for(i=0; timebuf[i]>=' '; i++)
276 time_string[i]=timebuf[i];
280 fprintf(fplog,"%s on node %d %s\n",title,nodeid,time_string);
284 static void sum_forces(int start,int end,rvec f[],rvec flr[])
289 pr_rvecs(debug,0,"fsr",f+start,end-start);
290 pr_rvecs(debug,0,"flr",flr+start,end-start);
292 for(i=start; (i<end); i++)
293 rvec_inc(f[i],flr[i]);
297 * calc_f_el calculates forces due to an electric field.
299 * force is kJ mol^-1 nm^-1 = e * kJ mol^-1 nm^-1 / e
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
307 * The function should return the energy due to the electric field
308 * (if any) but for now returns 0.
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.
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)
327 for(m=0; (m<DIM); m++)
334 Ext[m] = cos(Et[m].a[0]*(t-t0))*exp(-sqr(t-t0)/(2.0*sqr(Et[m].a[2])));
338 Ext[m] = cos(Et[m].a[0]*t);
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];
359 fprintf(fp,"%10g %10g %10g %10g #FIELD\n",t,
360 Ext[XX]/FIELDFAC,Ext[YY]/FIELDFAC,Ext[ZZ]/FIELDFAC);
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)
371 /* The short-range virial from surrounding boxes */
373 calc_vir(fplog,SHIFTS,fr->shift_vec,fr->fshift,vir_part,ePBC==epbcSCREW,box);
374 inc_nrnb(nrnb,eNR_VIRIAL,SHIFTS);
376 /* Calculate partial virial, for local atoms only, based on short range.
377 * Total virial is computed in global_stat, called from do_md
379 f_calc_vir(fplog,start,start+homenr,x,f,vir_part,graph,box);
380 inc_nrnb(nrnb,eNR_VIRIAL,homenr);
382 /* Add position restraint contribution */
383 for(i=0; i<DIM; i++) {
384 vir_part[i][i] += fr->vir_diag_posres[i];
387 /* Add wall contribution */
388 for(i=0; i<DIM; i++) {
389 vir_part[i][ZZ] += fr->vir_wall_z[i];
393 pr_rvecs(debug,0,"vir_part",vir_part,DIM);
396 static void posres_wrapper(FILE *fplog,
404 gmx_enerdata_t *enerd,
412 /* Position restraints always require full pbc */
413 set_pbc(&pbc,ir->ePBC,box);
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);
423 fprintf(fplog,sepdvdlformat,
424 interaction_function[F_POSRES].longname,v,dvdl);
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.
430 enerd->dvdl_nonlin[efptRESTRAINT] += dvdl;
431 inc_nrnb(nrnb,eNR_POSRES,top->idef.il[F_POSRES].nr/2);
433 if ((ir->fepvals->n_lambda > 0) && (flags & GMX_FORCE_DHDL))
435 for(i=0; i<enerd->n_lambda; i++)
437 real dvdl_dum,lambda_dum;
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;
450 static void pull_potential_wrapper(FILE *fplog,
458 gmx_enerdata_t *enerd,
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.
470 set_pbc(&pbc,ir->ePBC,box);
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);
477 fprintf(fplog,sepdvdlformat,"Com pull",enerd->term[F_COM_PULL],dvdl);
479 enerd->dvdl_lin[efptRESTRAINT] += dvdl;
482 static void pme_receive_force_ener(FILE *fplog,
485 gmx_wallcycle_t wcycle,
486 gmx_enerdata_t *enerd,
490 float cycles_ppdpme,cycles_seppme;
492 cycles_ppdpme = wallcycle_stop(wcycle,ewcPPDURINGPME);
493 dd_cycles_add(cr->dd,cycles_ppdpme,ddCyclPPduringPME);
495 /* In case of node-splitting, the PP nodes receive the long-range
496 * forces, virial and energy from the PME nodes here.
498 wallcycle_start(wcycle,ewcPP_PMEWAITRECVF);
500 gmx_pme_receive_f(cr,fr->f_novirsum,fr->vir_el_recip,&e,&dvdl,
504 fprintf(fplog,sepdvdlformat,"PME mesh",e,dvdl);
506 enerd->term[F_COUL_RECIP] += e;
507 enerd->dvdl_lin[efptCOUL] += dvdl;
510 dd_cycles_add(cr->dd,cycles_seppme,ddCyclPME);
512 wallcycle_stop(wcycle,ewcPP_PMEWAITRECVF);
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)
520 char buf[STEPSTRSIZE];
523 for(i=md->start; i<md->start+md->homenr; 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));
534 static void post_process_forces(FILE *fplog,
536 gmx_large_int_t step,
537 t_nrnb *nrnb,gmx_wallcycle_t wcycle,
544 t_forcerec *fr,gmx_vsite_t *vsite,
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.
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,
559 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
560 wallcycle_stop(wcycle,ewcVSITESPREAD);
562 if (flags & GMX_FORCE_VIRIAL)
564 /* Now add the forces, this is local */
567 sum_forces(0,fr->f_novirsum_n,f,fr->f_novirsum);
571 sum_forces(mdatoms->start,mdatoms->start+mdatoms->homenr,
574 if (EEL_FULL(fr->eeltype))
576 /* Add the mesh contribution to the virial */
577 m_add(vir_force,fr->vir_el_recip,vir_force);
581 pr_rvecs(debug,0,"vir_force",vir_force,DIM);
586 if (fr->print_force >= 0)
588 print_large_forces(stderr,mdatoms,cr,step,fr->print_force,x,f);
592 static void do_nb_verlet(t_forcerec *fr,
593 interaction_const_t *ic,
594 gmx_enerdata_t *enerd,
595 int flags, int ilocality,
598 gmx_wallcycle_t wcycle)
600 int nnbl, kernel_type, enr_nbnxn_kernel_ljc, enr_nbnxn_kernel_lj;
602 nonbonded_verlet_group_t *nbvg;
604 if (!(flags & GMX_FORCE_NONBONDED))
606 /* skip non-bonded calculation */
610 nbvg = &fr->nbv->grp[ilocality];
612 /* CUDA kernel launch overhead is already timed separately */
613 if (fr->cutoff_scheme != ecutsVERLET)
615 gmx_incons("Invalid cut-off scheme passed!");
618 if (nbvg->kernel_type != nbnxnk8x8x8_CUDA)
620 wallcycle_sub_start(wcycle, ewcsNONBONDED);
622 switch (nbvg->kernel_type)
624 case nbnxnk4x4_PlainC:
625 nbnxn_kernel_ref(&nbvg->nbl_lists,
631 enerd->grpp.ener[egCOULSR],
633 enerd->grpp.ener[egBHAMSR] :
634 enerd->grpp.ener[egLJSR]);
637 case nbnxnk4xN_SIMD_4xN:
638 nbnxn_kernel_simd_4xn(&nbvg->nbl_lists,
645 enerd->grpp.ener[egCOULSR],
647 enerd->grpp.ener[egBHAMSR] :
648 enerd->grpp.ener[egLJSR]);
650 case nbnxnk4xN_SIMD_2xNN:
651 nbnxn_kernel_simd_2xnn(&nbvg->nbl_lists,
658 enerd->grpp.ener[egCOULSR],
660 enerd->grpp.ener[egBHAMSR] :
661 enerd->grpp.ener[egLJSR]);
664 case nbnxnk8x8x8_CUDA:
665 nbnxn_cuda_launch_kernel(fr->nbv->cu_nbv, nbvg->nbat, flags, ilocality);
668 case nbnxnk8x8x8_PlainC:
669 nbnxn_kernel_gpu_ref(nbvg->nbl_lists.nbl[0],
674 nbvg->nbat->out[0].f,
676 enerd->grpp.ener[egCOULSR],
678 enerd->grpp.ener[egBHAMSR] :
679 enerd->grpp.ener[egLJSR]);
683 gmx_incons("Invalid nonbonded kernel type passed!");
686 if (nbvg->kernel_type != nbnxnk8x8x8_CUDA)
688 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
691 if (EEL_RF(ic->eeltype) || ic->eeltype == eelCUT)
693 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_RF;
695 else if (nbvg->ewald_excl == ewaldexclTable)
697 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_TAB;
701 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_EWALD;
703 enr_nbnxn_kernel_lj = eNR_NBNXN_LJ;
704 if (flags & GMX_FORCE_ENERGY)
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;
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);
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,
724 gmx_groups_t *groups,
725 matrix box,rvec x[],history_t *hist,
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,
741 gmx_bool bSepDVDL,bStateChanged,bNS,bFillGrid,bCalcCGCM,bBS;
742 gmx_bool bDoLongRange,bDoForces,bSepLRF,bUseGPU,bUseOrEmulGPU;
743 gmx_bool bDiffKernels=FALSE;
747 float cycles_pme,cycles_force;
748 nonbonded_verlet_t *nbv;
752 nb_kernel_type = fr->nbv->grp[0].kernel_type;
754 start = mdatoms->start;
755 homenr = mdatoms->homenr;
757 bSepDVDL = (fr->bSepDVDL && do_per_step(step,inputrec->nstlog));
759 clear_mat(vir_force);
762 if (DOMAINDECOMP(cr))
764 cg1 = cr->dd->ncg_tot;
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);
787 update_forcerec(fplog,fr,box);
789 if (NEED_MUTOT(*inputrec))
791 /* Calculate total (local) dipole moment in a temporary common array.
792 * This makes it possible to sum them over nodes faster.
794 calc_mu(start,homenr,
795 x,mdatoms->chargeA,mdatoms->chargeB,mdatoms->nChargePerturbed,
800 if (fr->ePBC != epbcNONE) {
801 /* Compute shift vectors every step,
802 * because of pressure coupling or box deformation!
804 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
805 calc_shifts(box,fr->shift_vec);
808 put_atoms_in_box_omp(fr->ePBC,box,homenr,x);
809 inc_nrnb(nrnb,eNR_SHIFTX,homenr);
811 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph) {
812 unshift_self(graph,box,x);
816 nbnxn_atomdata_copy_shiftvec(flags & GMX_FORCE_DYNAMICBOX,
817 fr->shift_vec,nbv->grp[0].nbat);
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.
827 wallcycle_start(wcycle,ewcPP_PMESENDX);
829 bBS = (inputrec->nwall == 2);
832 svmul(inputrec->wall_ewald_zfac,boxs[ZZ],boxs[ZZ]);
835 gmx_pme_send_x(cr,bBS ? boxs : box,x,
836 mdatoms->nChargePerturbed,lambda[efptCOUL],
837 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),step);
839 wallcycle_stop(wcycle,ewcPP_PMESENDX);
843 /* do gridding for pair search */
846 if (graph && bStateChanged)
848 /* Calculate intramolecular shift vectors to make molecules whole */
849 mk_mshift(fplog,graph,fr->ePBC,box,x);
853 box_diag[XX] = box[XX][XX];
854 box_diag[YY] = box[YY][YY];
855 box_diag[ZZ] = box[ZZ][ZZ];
857 wallcycle_start(wcycle,ewcNS);
860 wallcycle_sub_start(wcycle,ewcsNBS_GRID_LOCAL);
861 nbnxn_put_on_grid(nbv->nbs,fr->ePBC,box,
863 0,mdatoms->homenr,-1,fr->cginfo,x,
865 nbv->grp[eintLocal].kernel_type,
866 nbv->grp[eintLocal].nbat);
867 wallcycle_sub_stop(wcycle,ewcsNBS_GRID_LOCAL);
871 wallcycle_sub_start(wcycle,ewcsNBS_GRID_NONLOCAL);
872 nbnxn_put_on_grid_nonlocal(nbv->nbs,domdec_zones(cr->dd),
874 nbv->grp[eintNonlocal].kernel_type,
875 nbv->grp[eintNonlocal].nbat);
876 wallcycle_sub_stop(wcycle,ewcsNBS_GRID_NONLOCAL);
879 if (nbv->ngrp == 1 ||
880 nbv->grp[eintNonlocal].nbat == nbv->grp[eintLocal].nbat)
882 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat,eatAll,
883 nbv->nbs,mdatoms,fr->cginfo);
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);
892 wallcycle_stop(wcycle, ewcNS);
895 /* initialize the GPU atom data and copy shift vector */
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);
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);
910 /* do local pair search */
913 wallcycle_start_nocount(wcycle,ewcNS);
914 wallcycle_sub_start(wcycle,ewcsNBS_SEARCH_LOCAL);
915 nbnxn_make_pairlist(nbv->nbs,nbv->grp[eintLocal].nbat,
918 nbv->min_ci_balanced,
919 &nbv->grp[eintLocal].nbl_lists,
921 nbv->grp[eintLocal].kernel_type,
923 wallcycle_sub_stop(wcycle,ewcsNBS_SEARCH_LOCAL);
927 /* initialize local pair-list on the GPU */
928 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
929 nbv->grp[eintLocal].nbl_lists.nbl[0],
932 wallcycle_stop(wcycle, ewcNS);
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);
946 wallcycle_start(wcycle,ewcLAUNCH_GPU_NB);
947 /* launch local nonbonded F on GPU */
948 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFNo,
950 wallcycle_stop(wcycle,ewcLAUNCH_GPU_NB);
953 /* Communicate coordinates and sum dipole if necessary +
954 do non-local pair search */
955 if (DOMAINDECOMP(cr))
957 bDiffKernels = (nbv->grp[eintNonlocal].kernel_type !=
958 nbv->grp[eintLocal].kernel_type);
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.
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);
977 wallcycle_start_nocount(wcycle,ewcNS);
978 wallcycle_sub_start(wcycle,ewcsNBS_SEARCH_NONLOCAL);
982 nbnxn_grid_add_simple(nbv->nbs,nbv->grp[eintNonlocal].nbat);
985 nbnxn_make_pairlist(nbv->nbs,nbv->grp[eintNonlocal].nbat,
988 nbv->min_ci_balanced,
989 &nbv->grp[eintNonlocal].nbl_lists,
991 nbv->grp[eintNonlocal].kernel_type,
994 wallcycle_sub_stop(wcycle,ewcsNBS_SEARCH_NONLOCAL);
996 if (nbv->grp[eintNonlocal].kernel_type == nbnxnk8x8x8_CUDA)
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],
1003 wallcycle_stop(wcycle,ewcNS);
1007 wallcycle_start(wcycle,ewcMOVEX);
1008 dd_move_x(cr->dd,box,x);
1010 /* When we don't need the total dipole we sum it in global_stat */
1011 if (bStateChanged && NEED_MUTOT(*inputrec))
1013 gmx_sumd(2*DIM,mu,cr);
1015 wallcycle_stop(wcycle,ewcMOVEX);
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);
1025 if (bUseGPU && !bDiffKernels)
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,
1031 cycles_force += wallcycle_stop(wcycle,ewcLAUNCH_GPU_NB);
1037 /* launch D2H copy-back F */
1038 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1039 if (DOMAINDECOMP(cr) && !bDiffKernels)
1041 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintNonlocal].nbat,
1042 flags, eatNonlocal);
1044 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintLocal].nbat,
1046 cycles_force += wallcycle_stop(wcycle,ewcLAUNCH_GPU_NB);
1049 if (bStateChanged && NEED_MUTOT(*inputrec))
1053 gmx_sumd(2*DIM,mu,cr);
1060 fr->mu_tot[i][j] = mu[i*DIM + j];
1064 if (fr->efep == efepNO)
1066 copy_rvec(fr->mu_tot[0],mu_tot);
1070 for(j=0; j<DIM; j++)
1073 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] +
1074 lambda[efptCOUL]*fr->mu_tot[1][j];
1078 /* Reset energies */
1079 reset_enerdata(&(inputrec->opts),fr,bNS,enerd,MASTER(cr));
1080 clear_rvecs(SHIFTS,fr->fshift);
1082 if (DOMAINDECOMP(cr))
1084 if (!(cr->duty & DUTY_PME))
1086 wallcycle_start(wcycle,ewcPPDURINGPME);
1087 dd_force_flop_start(cr->dd,nrnb);
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.
1096 wallcycle_start(wcycle,ewcFORCE);
1099 /* Reset forces for which the virial is calculated separately:
1100 * PME/Ewald forces if necessary */
1101 if (fr->bF_NoVirSum)
1103 if (flags & GMX_FORCE_VIRIAL)
1105 fr->f_novirsum = fr->f_novirsum_alloc;
1108 clear_rvecs(fr->f_novirsum_n,fr->f_novirsum);
1112 clear_rvecs(homenr,fr->f_novirsum+start);
1117 /* We are not calculating the pressure so we do not need
1118 * a separate array for forces that do not contribute
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))
1129 clear_rvecs(fr->natoms_force_constr,fr->f_twin);
1132 clear_rvec(fr->vir_diag_posres);
1134 if (inputrec->ePull == epullCONSTRAINT)
1136 clear_pull_forces(inputrec->pull);
1139 /* update QMMMrec, if necessary */
1142 update_QMMMrec(cr,fr,x,mdatoms,box,top);
1145 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1147 posres_wrapper(fplog,flags,bSepDVDL,inputrec,nrnb,top,box,x,
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);
1161 if (do_per_step(step,inputrec->nstcalclr))
1163 /* Add the long range forces to the short range forces */
1164 for(i=0; i<fr->natoms_force_constr; i++)
1166 rvec_add(fr->f_twin[i],f[i],f[i]);
1173 /* Maybe we should move this into do_force_lowlevel */
1174 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFYes,
1179 if (!bUseOrEmulGPU || bDiffKernels)
1183 if (DOMAINDECOMP(cr))
1185 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal,
1186 bDiffKernels ? enbvClearFYes : enbvClearFNo,
1196 aloc = eintNonlocal;
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.
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);
1211 /* if there are multiple fshift output buffers reduce them */
1212 if ((flags & GMX_FORCE_VIRIAL) &&
1213 nbv->grp[aloc].nbl_lists.nnbl > 1)
1215 nbnxn_atomdata_add_nbat_fshift_to_fshift(nbv->grp[aloc].nbat,
1220 cycles_force += wallcycle_stop(wcycle,ewcFORCE);
1224 do_flood(cr,inputrec,x,f,ed,box,step,bNS);
1227 if (bUseOrEmulGPU && !bDiffKernels)
1229 /* wait for non-local forces (or calculate in emulation mode) */
1230 if (DOMAINDECOMP(cr))
1234 wallcycle_start(wcycle,ewcWAIT_GPU_NB_NL);
1235 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1236 nbv->grp[eintNonlocal].nbat,
1238 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1240 cycles_force += wallcycle_stop(wcycle,ewcWAIT_GPU_NB_NL);
1244 wallcycle_start_nocount(wcycle,ewcFORCE);
1245 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFYes,
1247 cycles_force += wallcycle_stop(wcycle,ewcFORCE);
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)
1254 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs,eatNonlocal,
1255 nbv->grp[eintNonlocal].nbat,f);
1257 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1258 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1264 /* Communicate the forces */
1267 wallcycle_start(wcycle,ewcMOVEF);
1268 if (DOMAINDECOMP(cr))
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.
1278 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1279 (flags & GMX_FORCE_VIRIAL))
1281 dd_move_f(cr->dd,fr->f_novirsum,NULL);
1285 /* We should not update the shift forces here,
1286 * since f_twin is already included in f.
1288 dd_move_f(cr->dd,fr->f_twin,NULL);
1291 wallcycle_stop(wcycle,ewcMOVEF);
1297 /* wait for local forces (or calculate in emulation mode) */
1300 wallcycle_start(wcycle,ewcWAIT_GPU_NB_L);
1301 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1302 nbv->grp[eintLocal].nbat,
1304 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1306 wallcycle_stop(wcycle,ewcWAIT_GPU_NB_L);
1308 /* now clear the GPU outputs while we finish the step on the CPU */
1310 wallcycle_start_nocount(wcycle,ewcLAUNCH_GPU_NB);
1311 nbnxn_cuda_clear_outputs(nbv->cu_nbv, flags);
1312 wallcycle_stop(wcycle,ewcLAUNCH_GPU_NB);
1316 wallcycle_start_nocount(wcycle,ewcFORCE);
1317 do_nb_verlet(fr, ic, enerd, flags, eintLocal,
1318 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
1320 wallcycle_stop(wcycle,ewcFORCE);
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)
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);
1330 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1331 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1334 if (DOMAINDECOMP(cr))
1336 dd_force_flop_stop(cr->dd,nrnb);
1339 dd_cycles_add(cr->dd,cycles_force-cycles_pme,ddCyclF);
1345 if (IR_ELEC_FIELD(*inputrec))
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);
1353 /* If we have NoVirSum forces, but we do not calculate the virial,
1354 * we sum fr->f_novirum=f later.
1356 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
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);
1365 wallcycle_start(wcycle,ewcVSITESPREAD);
1366 spread_vsite_f(fplog,vsite,x,fr->f_twin,NULL,FALSE,NULL,
1368 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
1369 wallcycle_stop(wcycle,ewcVSITESPREAD);
1373 if (flags & GMX_FORCE_VIRIAL)
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);
1381 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1383 pull_potential_wrapper(fplog,bSepDVDL,cr,inputrec,box,x,
1384 f,vir_force,mdatoms,enerd,lambda,t);
1387 if (PAR(cr) && !(cr->duty & DUTY_PME))
1389 /* In case of node-splitting, the PP nodes receive the long-range
1390 * forces, virial and energy from the PME nodes here.
1392 pme_receive_force_ener(fplog,bSepDVDL,cr,wcycle,enerd,fr);
1397 post_process_forces(fplog,cr,step,nrnb,wcycle,
1398 top,box,x,f,vir_force,mdatoms,graph,fr,vsite,
1402 /* Sum the potential energy terms from group contributions */
1403 sum_epot(&(inputrec->opts),&(enerd->grpp),enerd->term);
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,
1411 gmx_groups_t *groups,
1412 matrix box,rvec x[],history_t *hist,
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,
1426 gmx_bool bSepDVDL,bStateChanged,bNS,bFillGrid,bCalcCGCM,bBS;
1427 gmx_bool bDoLongRangeNS,bDoForces,bDoPotential,bSepLRF;
1428 gmx_bool bDoAdressWF;
1430 rvec vzero,box_diag;
1431 real e,v,dvdlambda[efptNR];
1433 float cycles_pme,cycles_force;
1435 start = mdatoms->start;
1436 homenr = mdatoms->homenr;
1438 bSepDVDL = (fr->bSepDVDL && do_per_step(step,inputrec->nstlog));
1440 clear_mat(vir_force);
1444 pd_cg_range(cr,&cg0,&cg1);
1449 if (DOMAINDECOMP(cr))
1451 cg1 = cr->dd->ncg_tot;
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));
1475 /* should probably move this to the forcerec since it doesn't change */
1476 bDoAdressWF = ((fr->adress_type!=eAdressOff));
1480 update_forcerec(fplog,fr,box);
1482 if (NEED_MUTOT(*inputrec))
1484 /* Calculate total (local) dipole moment in a temporary common array.
1485 * This makes it possible to sum them over nodes faster.
1487 calc_mu(start,homenr,
1488 x,mdatoms->chargeA,mdatoms->chargeB,mdatoms->nChargePerturbed,
1493 if (fr->ePBC != epbcNONE) {
1494 /* Compute shift vectors every step,
1495 * because of pressure coupling or box deformation!
1497 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
1498 calc_shifts(box,fr->shift_vec);
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);
1506 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph) {
1507 unshift_self(graph,box,x);
1510 else if (bCalcCGCM) {
1511 calc_cgcm(fplog,cg0,cg1,&(top->cgs),x,fr->cg_cm);
1512 inc_nrnb(nrnb,eNR_CGCM,homenr);
1517 move_cgcm(fplog,cr,fr->cg_cm);
1520 pr_rvecs(debug,0,"cgcm",fr->cg_cm,top->cgs.nr);
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.
1531 wallcycle_start(wcycle,ewcPP_PMESENDX);
1533 bBS = (inputrec->nwall == 2);
1536 svmul(inputrec->wall_ewald_zfac,boxs[ZZ],boxs[ZZ]);
1539 gmx_pme_send_x(cr,bBS ? boxs : box,x,
1540 mdatoms->nChargePerturbed,lambda[efptCOUL],
1541 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),step);
1543 wallcycle_stop(wcycle,ewcPP_PMESENDX);
1545 #endif /* GMX_MPI */
1547 /* Communicate coordinates and sum dipole if necessary */
1550 wallcycle_start(wcycle,ewcMOVEX);
1551 if (DOMAINDECOMP(cr))
1553 dd_move_x(cr->dd,box,x);
1557 move_x(fplog,cr,GMX_LEFT,GMX_RIGHT,x,nrnb);
1559 wallcycle_stop(wcycle,ewcMOVEX);
1562 /* update adress weight beforehand */
1563 if(bStateChanged && bDoAdressWF)
1565 /* need pbc for adress weight calculation with pbc_dx */
1566 set_pbc(&pbc,inputrec->ePBC,box);
1567 if(fr->adress_site == eAdressSITEcog)
1569 update_adress_weights_cog(top->idef.iparams,top->idef.il,x,fr,mdatoms,
1570 inputrec->ePBC==epbcNONE ? NULL : &pbc);
1572 else if (fr->adress_site == eAdressSITEcom)
1574 update_adress_weights_com(fplog,cg0,cg1,&(top->cgs),x,fr,mdatoms,
1575 inputrec->ePBC==epbcNONE ? NULL : &pbc);
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);
1583 update_adress_weights_atom(cg0,cg1,&(top->cgs),x,fr,mdatoms,
1584 inputrec->ePBC==epbcNONE ? NULL : &pbc);
1588 if (NEED_MUTOT(*inputrec))
1595 gmx_sumd(2*DIM,mu,cr);
1601 fr->mu_tot[i][j] = mu[i*DIM + j];
1605 if (fr->efep == efepNO)
1607 copy_rvec(fr->mu_tot[0],mu_tot);
1611 for(j=0; j<DIM; j++)
1614 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] + lambda[efptCOUL]*fr->mu_tot[1][j];
1619 /* Reset energies */
1620 reset_enerdata(&(inputrec->opts),fr,bNS,enerd,MASTER(cr));
1621 clear_rvecs(SHIFTS,fr->fshift);
1625 wallcycle_start(wcycle,ewcNS);
1627 if (graph && bStateChanged)
1629 /* Calculate intramolecular shift vectors to make molecules whole */
1630 mk_mshift(fplog,graph,fr->ePBC,box,x);
1633 /* Do the actual neighbour searching and if twin range electrostatics
1634 * also do the calculation of long range forces and energies.
1636 for (i=0;i<efptNR;i++) {dvdlambda[i] = 0;}
1638 groups,&(inputrec->opts),top,mdatoms,
1639 cr,nrnb,lambda,dvdlambda,&enerd->grpp,bFillGrid,
1643 fprintf(fplog,sepdvdlformat,"LR non-bonded",0.0,dvdlambda);
1645 enerd->dvdl_lin[efptVDW] += dvdlambda[efptVDW];
1646 enerd->dvdl_lin[efptCOUL] += dvdlambda[efptCOUL];
1648 wallcycle_stop(wcycle,ewcNS);
1651 if (inputrec->implicit_solvent && bNS)
1653 make_gb_nblist(cr,inputrec->gb_algorithm,inputrec->rlist,
1654 x,box,fr,&top->idef,graph,fr->born);
1657 if (DOMAINDECOMP(cr))
1659 if (!(cr->duty & DUTY_PME))
1661 wallcycle_start(wcycle,ewcPPDURINGPME);
1662 dd_force_flop_start(cr->dd,nrnb);
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);
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.
1681 wallcycle_start(wcycle,ewcFORCE);
1685 /* Reset forces for which the virial is calculated separately:
1686 * PME/Ewald forces if necessary */
1687 if (fr->bF_NoVirSum)
1689 if (flags & GMX_FORCE_VIRIAL)
1691 fr->f_novirsum = fr->f_novirsum_alloc;
1694 clear_rvecs(fr->f_novirsum_n,fr->f_novirsum);
1698 clear_rvecs(homenr,fr->f_novirsum+start);
1703 /* We are not calculating the pressure so we do not need
1704 * a separate array for forces that do not contribute
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))
1715 clear_rvecs(fr->natoms_force_constr,fr->f_twin);
1718 clear_rvec(fr->vir_diag_posres);
1720 if (inputrec->ePull == epullCONSTRAINT)
1722 clear_pull_forces(inputrec->pull);
1725 /* update QMMMrec, if necessary */
1728 update_QMMMrec(cr,fr,x,mdatoms,box,top);
1731 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1733 posres_wrapper(fplog,flags,bSepDVDL,inputrec,nrnb,top,box,x,
1737 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_FBPOSRES].nr > 0)
1739 /* Flat-bottomed position restraints always require full pbc */
1740 if(!(bStateChanged && bDoAdressWF))
1742 set_pbc(&pbc,inputrec->ePBC,box);
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);
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,
1765 if (do_per_step(step,inputrec->nstcalclr))
1767 /* Add the long range forces to the short range forces */
1768 for(i=0; i<fr->natoms_force_constr; i++)
1770 rvec_add(fr->f_twin[i],f[i],f[i]);
1775 cycles_force = wallcycle_stop(wcycle,ewcFORCE);
1779 do_flood(cr,inputrec,x,f,ed,box,step,bNS);
1782 if (DOMAINDECOMP(cr))
1784 dd_force_flop_stop(cr->dd,nrnb);
1787 dd_cycles_add(cr->dd,cycles_force-cycles_pme,ddCyclF);
1793 if (IR_ELEC_FIELD(*inputrec))
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);
1801 if (bDoAdressWF && fr->adress_icor == eAdressICThermoForce)
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);
1808 /* Communicate the forces */
1811 wallcycle_start(wcycle,ewcMOVEF);
1812 if (DOMAINDECOMP(cr))
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.
1822 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1823 (flags & GMX_FORCE_VIRIAL))
1825 dd_move_f(cr->dd,fr->f_novirsum,NULL);
1829 /* We should not update the shift forces here,
1830 * since f_twin is already included in f.
1832 dd_move_f(cr->dd,fr->f_twin,NULL);
1837 pd_move_f(cr,f,nrnb);
1840 pd_move_f(cr,fr->f_twin,nrnb);
1843 wallcycle_stop(wcycle,ewcMOVEF);
1846 /* If we have NoVirSum forces, but we do not calculate the virial,
1847 * we sum fr->f_novirum=f later.
1849 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
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);
1858 wallcycle_start(wcycle,ewcVSITESPREAD);
1859 spread_vsite_f(fplog,vsite,x,fr->f_twin,NULL,FALSE,NULL,
1861 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
1862 wallcycle_stop(wcycle,ewcVSITESPREAD);
1866 if (flags & GMX_FORCE_VIRIAL)
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);
1874 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1876 pull_potential_wrapper(fplog,bSepDVDL,cr,inputrec,box,x,
1877 f,vir_force,mdatoms,enerd,lambda,t);
1880 /* Add the forces from enforced rotation potentials (if any) */
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);
1888 if (PAR(cr) && !(cr->duty & DUTY_PME))
1890 /* In case of node-splitting, the PP nodes receive the long-range
1891 * forces, virial and energy from the PME nodes here.
1893 pme_receive_force_ener(fplog,bSepDVDL,cr,wcycle,enerd,fr);
1898 post_process_forces(fplog,cr,step,nrnb,wcycle,
1899 top,box,x,f,vir_force,mdatoms,graph,fr,vsite,
1903 /* Sum the potential energy terms from group contributions */
1904 sum_epot(&(inputrec->opts),&(enerd->grpp),enerd->term);
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,
1912 gmx_groups_t *groups,
1913 matrix box,rvec x[],history_t *hist,
1917 gmx_enerdata_t *enerd,t_fcdata *fcd,
1918 real *lambda,t_graph *graph,
1920 gmx_vsite_t *vsite,rvec mu_tot,
1921 double t,FILE *field,gmx_edsam_t ed,
1922 gmx_bool bBornRadii,
1925 /* modify force flag if not doing nonbonded */
1926 if (!fr->bNonbonded)
1928 flags &= ~GMX_FORCE_NONBONDED;
1931 switch (inputrec->cutoff_scheme)
1934 do_force_cutsVERLET(fplog, cr, inputrec,
1950 do_force_cutsGROUP(fplog, cr, inputrec,
1965 gmx_incons("Invalid cut-off scheme passed!");
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)
1977 gmx_large_int_t step;
1978 real dt=ir->delta_t;
1982 snew(savex,state->natoms);
1985 end = md->homenr + start;
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;
1994 char buf[STEPSTRSIZE];
1995 fprintf(fplog,"\nConstraining the starting coordinates (step %s)\n",
1996 gmx_step_str(step,buf));
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);
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);
2021 /* constrain the inital velocities at t-dt/2 */
2022 if (EI_STATE_VELOCITY(ir->eI) && ir->eI!=eiVV)
2024 for(i=start; (i<end); i++)
2026 for(m=0; (m<DIM); m++)
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];
2034 /* Shake the positions at t=-dt with the positions at t=0
2035 * as reference coordinates.
2039 char buf[STEPSTRSIZE];
2040 fprintf(fplog,"\nConstraining the coordinates at t0-dt (step %s)\n",
2041 gmx_step_str(step,buf));
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);
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];
2062 void calc_enervirdiff(FILE *fplog,int eDispCorr,t_forcerec *fr)
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;
2070 fr->enershiftsix = 0;
2071 fr->enershifttwelve = 0;
2072 fr->enerdiffsix = 0;
2073 fr->enerdifftwelve = 0;
2075 fr->virdifftwelve = 0;
2077 if (eDispCorr != edispcNO) {
2078 for(i=0; i<2; i++) {
2082 if ((fr->vdwtype == evdwSWITCH) || (fr->vdwtype == evdwSHIFT)) {
2083 if (fr->rvdw_switch == 0)
2085 "With dispersion correction rvdw-switch can not be zero "
2086 "for vdw-type = %s",evdw_names[fr->vdwtype]);
2088 scale = fr->nblists[0].table_elec_vdw.scale;
2089 vdwtab = fr->nblists[0].table_vdw.data;
2091 /* Round the cut-offs to exact table values for precision */
2092 ri0 = floor(fr->rvdw_switch*scale);
2093 ri1 = ceil(fr->rvdw*scale);
2099 if (fr->vdwtype == evdwSHIFT)
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.
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];
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.
2113 eners[0] += 4.0*M_PI*fr->enershiftsix*rc3/3.0;
2114 eners[1] += 4.0*M_PI*fr->enershifttwelve*rc3/3.0;
2116 invscale = 1.0/(scale);
2117 invscale2 = invscale*invscale;
2118 invscale3 = invscale*invscale2;
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
2131 enersum = 0.0; virsum = 0.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.
2146 for (ri=ri0; ri<ri1; ri++) {
2149 eb = 2.0*invscale2*r;
2153 pb = 3.0*invscale2*r;
2154 pc = 3.0*invscale*r*r;
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];
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);
2168 enersum *= 4.0*M_PI*tabfactor;
2169 virsum *= 4.0*M_PI*tabfactor;
2170 eners[i] -= enersum;
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);
2180 else if ((fr->vdwtype == evdwCUT) || (fr->vdwtype == evdwUSER)) {
2181 if (fr->vdwtype == evdwUSER && fplog)
2183 "WARNING: using dispersion correction with user tables\n");
2184 rc3 = fr->rvdw*fr->rvdw*fr->rvdw;
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);
2194 /* Contribution beyond the cut-off */
2195 virs[0] += 8.0*M_PI/rc3;
2196 virs[1] += -16.0*M_PI/(3.0*rc9);
2199 "Dispersion correction is not implemented for vdw-type = %s",
2200 evdw_names[fr->vdwtype]);
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];
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)
2215 gmx_bool bCorrAll,bCorrPres;
2216 real dvdlambda,invvol,dens,ninter,avcsix,avctwelve,enerdiff,svir=0,spres=0;
2226 if (ir->eDispCorr != edispcNO) {
2227 bCorrAll = (ir->eDispCorr == edispcAllEner ||
2228 ir->eDispCorr == edispcAllEnerPres);
2229 bCorrPres = (ir->eDispCorr == edispcEnerPres ||
2230 ir->eDispCorr == edispcAllEnerPres);
2232 invvol = 1/det(box);
2235 /* Only correct for the interactions with the inserted molecule */
2236 dens = (natoms - fr->n_tpi)*invvol;
2241 dens = natoms*invvol;
2242 ninter = 0.5*natoms;
2245 if (ir->efep == efepNO)
2247 avcsix = fr->avcsix[0];
2248 avctwelve = fr->avctwelve[0];
2252 avcsix = (1 - lambda)*fr->avcsix[0] + lambda*fr->avcsix[1];
2253 avctwelve = (1 - lambda)*fr->avctwelve[0] + lambda*fr->avctwelve[1];
2256 enerdiff = ninter*(dens*fr->enerdiffsix - fr->enershiftsix);
2257 *enercorr += avcsix*enerdiff;
2259 if (ir->efep != efepNO)
2261 dvdlambda += (fr->avcsix[1] - fr->avcsix[0])*enerdiff;
2265 enerdiff = ninter*(dens*fr->enerdifftwelve - fr->enershifttwelve);
2266 *enercorr += avctwelve*enerdiff;
2267 if (fr->efep != efepNO)
2269 dvdlambda += (fr->avctwelve[1] - fr->avctwelve[0])*enerdiff;
2275 svir = ninter*dens*avcsix*fr->virdiffsix/3.0;
2276 if (ir->eDispCorr == edispcAllEnerPres)
2278 svir += ninter*dens*avctwelve*fr->virdifftwelve/3.0;
2280 /* The factor 2 is because of the Gromacs virial definition */
2281 spres = -2.0*invvol*svir*PRESFAC;
2283 for(m=0; m<DIM; m++) {
2284 virial[m][m] += svir;
2285 pres[m][m] += spres;
2290 /* Can't currently control when it prints, for now, just print when degugging */
2294 fprintf(debug,"Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
2300 "Long Range LJ corr.: Epot %10g, Pres: %10g, Vir: %10g\n",
2301 *enercorr,spres,svir);
2305 fprintf(debug,"Long Range LJ corr.: Epot %10g\n",*enercorr);
2309 if (fr->bSepDVDL && do_per_step(step,ir->nstlog))
2311 fprintf(fplog,sepdvdlformat,"Dispersion correction",
2312 *enercorr,dvdlambda);
2314 if (fr->efep != efepNO)
2316 *dvdlcorr += dvdlambda;
2321 void do_pbc_first(FILE *fplog,matrix box,t_forcerec *fr,
2322 t_graph *graph,rvec x[])
2325 fprintf(fplog,"Removing pbc first time\n");
2326 calc_shifts(box,fr->shift_vec);
2328 mk_mshift(fplog,graph,fr->ePBC,box,x);
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
2337 mk_mshift(fplog,graph,fr->ePBC,box,x);
2339 p_graph(debug,"do_pbc_first 2",graph);
2342 fprintf(fplog,"Done rmpbc\n");
2345 static void low_do_pbc_mtop(FILE *fplog,int ePBC,matrix box,
2346 gmx_mtop_t *mtop,rvec x[],
2351 gmx_molblock_t *molb;
2353 if (bFirst && fplog)
2354 fprintf(fplog,"Removing pbc first time\n");
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;
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);
2369 for(mol=0; mol<molb->nmol; mol++) {
2370 mk_mshift(fplog,graph,ePBC,box,x+as);
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.
2378 as += molb->natoms_mol;
2386 void do_pbc_first_mtop(FILE *fplog,int ePBC,matrix box,
2387 gmx_mtop_t *mtop,rvec x[])
2389 low_do_pbc_mtop(fplog,ePBC,box,mtop,x,TRUE);
2392 void do_pbc_mtop(FILE *fplog,int ePBC,matrix box,
2393 gmx_mtop_t *mtop,rvec x[])
2395 low_do_pbc_mtop(fplog,ePBC,box,mtop,x,FALSE);
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,
2404 gmx_bool bWriteStat)
2407 t_nrnb *nrnb_tot=NULL;
2411 wallcycle_sum(cr,wcycle);
2417 MPI_Allreduce(nrnb->n,nrnb_tot->n,eNRNB,MPI_DOUBLE,MPI_SUM,
2418 cr->mpi_comm_mysim);
2426 #if defined(GMX_MPI) && !defined(GMX_THREAD_MPI)
2429 /* reduce nodetime over all MPI processes in the current simulation */
2431 MPI_Allreduce(&runtime->proctime,&sum,1,MPI_DOUBLE,MPI_SUM,
2432 cr->mpi_comm_mysim);
2433 runtime->proctime = sum;
2439 print_flop(fplog,nrnb_tot,&nbfs,&mflop);
2446 if ((cr->duty & DUTY_PP) && DOMAINDECOMP(cr))
2448 print_dd_statistics(cr,inputrec,fplog);
2460 snew(nrnb_all,cr->nnodes);
2461 nrnb_all[0] = *nrnb;
2462 for(s=1; s<cr->nnodes; s++)
2464 MPI_Recv(nrnb_all[s].n,eNRNB,MPI_DOUBLE,s,0,
2465 cr->mpi_comm_mysim,&stat);
2467 pr_load(fplog,cr,nrnb_all);
2472 MPI_Send(nrnb->n,eNRNB,MPI_DOUBLE,MASTERRANK(cr),0,
2473 cr->mpi_comm_mysim);
2480 wallcycle_print(fplog,cr->nnodes,cr->npmenodes,runtime->realtime,
2483 if (EI_DYNAMICS(inputrec->eI))
2485 delta_t = inputrec->delta_t;
2494 print_perf(fplog,runtime->proctime,runtime->realtime,
2495 cr->nnodes-cr->npmenodes,
2496 runtime->nsteps_done,delta_t,nbfs,mflop,
2501 print_perf(stderr,runtime->proctime,runtime->realtime,
2502 cr->nnodes-cr->npmenodes,
2503 runtime->nsteps_done,delta_t,nbfs,mflop,
2509 extern void initialize_lambdas(FILE *fplog,t_inputrec *ir,int *fep_state,real *lambda,double *lam0)
2511 /* this function works, but could probably use a logic rewrite to keep all the different
2512 types of efep straight. */
2515 t_lambda *fep = ir->fepvals;
2517 if ((ir->efep==efepNO) && (ir->bSimTemp == FALSE)) {
2518 for (i=0;i<efptNR;i++) {
2527 *fep_state = fep->init_fep_state; /* this might overwrite the checkpoint
2528 if checkpoint is set -- a kludge is in for now
2530 for (i=0;i<efptNR;i++)
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 */
2535 lambda[i] = fep->init_lambda;
2537 lam0[i] = lambda[i];
2542 lambda[i] = fep->all_lambda[i][*fep_state];
2544 lam0[i] = lambda[i];
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];
2558 /* Send to the log the information on the current lambdas */
2561 fprintf(fplog,"Initial vector of lambda components:[ ");
2562 for (i=0;i<efptNR;i++)
2564 fprintf(fplog,"%10.4f ",lambda[i]);
2566 fprintf(fplog,"]\n");
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,
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)
2586 /* Initial values */
2587 *t = *t0 = ir->init_t;
2590 for(i=0;i<ir->opts.ngtc;i++)
2592 /* set bSimAnn if any group is being annealed */
2593 if(ir->opts.annealing[i]!=eannNO)
2600 update_annealing_target_temp(&(ir->opts),ir->init_t);
2603 /* Initialize lambda variables */
2604 initialize_lambdas(fplog,ir,fep_state,lambda,lam0);
2608 *upd = init_update(fplog,ir);
2614 *vcm = init_vcm(fplog,&mtop->groups,ir);
2617 if (EI_DYNAMICS(ir->eI) && !(Flags & MD_APPENDFILES))
2619 if (ir->etc == etcBERENDSEN)
2621 please_cite(fplog,"Berendsen84a");
2623 if (ir->etc == etcVRESCALE)
2625 please_cite(fplog,"Bussi2007a");
2633 *outf = init_mdoutf(nfile,fnm,Flags,cr,ir,oenv);
2635 *mdebin = init_mdebin((Flags & MD_APPENDFILES) ? NULL : (*outf)->fp_ene,
2636 mtop,ir, (*outf)->fp_dhdl);
2641 please_cite(fplog,"Fritsch12");
2642 please_cite(fplog,"Junghans10");
2644 /* Initiate variables */
2645 clear_mat(force_vir);
2646 clear_mat(shake_vir);