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_x86_simd128.h"
93 #include "nbnxn_kernels/nbnxn_kernel_x86_simd256.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, sh_e;
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 != nbk8x8x8_CUDA)
620 wallcycle_sub_start(wcycle, ewcsNONBONDED);
622 switch (nbvg->kernel_type)
625 nbnxn_kernel_ref(&nbvg->nbl_lists,
631 enerd->grpp.ener[egCOULSR],
633 enerd->grpp.ener[egBHAMSR] :
634 enerd->grpp.ener[egLJSR]);
637 case nbk4xN_X86_SIMD128:
638 nbnxn_kernel_x86_simd128(&nbvg->nbl_lists,
644 enerd->grpp.ener[egCOULSR],
646 enerd->grpp.ener[egBHAMSR] :
647 enerd->grpp.ener[egLJSR]);
649 case nbk4xN_X86_SIMD256:
650 nbnxn_kernel_x86_simd256(&nbvg->nbl_lists,
656 enerd->grpp.ener[egCOULSR],
658 enerd->grpp.ener[egBHAMSR] :
659 enerd->grpp.ener[egLJSR]);
663 nbnxn_cuda_launch_kernel(fr->nbv->cu_nbv, nbvg->nbat, flags, ilocality);
666 case nbk8x8x8_PlainC:
667 nbnxn_kernel_gpu_ref(nbvg->nbl_lists.nbl[0],
672 nbvg->nbat->out[0].f,
674 enerd->grpp.ener[egCOULSR],
676 enerd->grpp.ener[egBHAMSR] :
677 enerd->grpp.ener[egLJSR]);
681 gmx_incons("Invalid nonbonded kernel type passed!");
684 if (nbvg->kernel_type != nbk8x8x8_CUDA)
686 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
689 /* In eNR_??? the nbnxn F+E kernels are always the F kernel + 1 */
690 sh_e = ((flags & GMX_FORCE_ENERGY) ? 1 : 0);
692 ((EEL_RF(ic->eeltype) || ic->eeltype == eelCUT) ?
693 eNR_NBNXN_LJ_RF : eNR_NBNXN_LJ_TAB) + sh_e,
694 nbvg->nbl_lists.natpair_ljq);
695 inc_nrnb(nrnb,eNR_NBNXN_LJ+sh_e,nbvg->nbl_lists.natpair_lj);
697 ((EEL_RF(ic->eeltype) || ic->eeltype == eelCUT) ?
698 eNR_NBNXN_RF : eNR_NBNXN_TAB)+sh_e,
699 nbvg->nbl_lists.natpair_q);
702 void do_force_cutsVERLET(FILE *fplog,t_commrec *cr,
703 t_inputrec *inputrec,
704 gmx_large_int_t step,t_nrnb *nrnb,gmx_wallcycle_t wcycle,
707 gmx_groups_t *groups,
708 matrix box,rvec x[],history_t *hist,
712 gmx_enerdata_t *enerd,t_fcdata *fcd,
713 real *lambda,t_graph *graph,
714 t_forcerec *fr, interaction_const_t *ic,
715 gmx_vsite_t *vsite,rvec mu_tot,
716 double t,FILE *field,gmx_edsam_t ed,
724 gmx_bool bSepDVDL,bStateChanged,bNS,bFillGrid,bCalcCGCM,bBS;
725 gmx_bool bDoLongRange,bDoForces,bSepLRF,bUseGPU,bUseOrEmulGPU;
726 gmx_bool bDiffKernels=FALSE;
730 float cycles_pme,cycles_force;
731 nonbonded_verlet_t *nbv;
735 nb_kernel_type = fr->nbv->grp[0].kernel_type;
737 start = mdatoms->start;
738 homenr = mdatoms->homenr;
740 bSepDVDL = (fr->bSepDVDL && do_per_step(step,inputrec->nstlog));
742 clear_mat(vir_force);
745 if (DOMAINDECOMP(cr))
747 cg1 = cr->dd->ncg_tot;
758 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
759 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll==FALSE);
760 bFillGrid = (bNS && bStateChanged);
761 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
762 bDoLongRange = (fr->bTwinRange && bNS && (flags & GMX_FORCE_DO_LR));
763 bDoForces = (flags & GMX_FORCE_FORCES);
764 bSepLRF = (bDoLongRange && bDoForces && (flags & GMX_FORCE_SEPLRF));
765 bUseGPU = fr->nbv->bUseGPU;
766 bUseOrEmulGPU = bUseGPU || (nbv->grp[0].kernel_type == nbk8x8x8_PlainC);
770 update_forcerec(fplog,fr,box);
772 if (NEED_MUTOT(*inputrec))
774 /* Calculate total (local) dipole moment in a temporary common array.
775 * This makes it possible to sum them over nodes faster.
777 calc_mu(start,homenr,
778 x,mdatoms->chargeA,mdatoms->chargeB,mdatoms->nChargePerturbed,
783 if (fr->ePBC != epbcNONE) {
784 /* Compute shift vectors every step,
785 * because of pressure coupling or box deformation!
787 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
788 calc_shifts(box,fr->shift_vec);
791 put_atoms_in_box_omp(fr->ePBC,box,homenr,x);
792 inc_nrnb(nrnb,eNR_SHIFTX,homenr);
794 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph) {
795 unshift_self(graph,box,x);
799 nbnxn_atomdata_copy_shiftvec(flags & GMX_FORCE_DYNAMICBOX,
800 fr->shift_vec,nbv->grp[0].nbat);
803 if (!(cr->duty & DUTY_PME)) {
804 /* Send particle coordinates to the pme nodes.
805 * Since this is only implemented for domain decomposition
806 * and domain decomposition does not use the graph,
807 * we do not need to worry about shifting.
810 wallcycle_start(wcycle,ewcPP_PMESENDX);
812 bBS = (inputrec->nwall == 2);
815 svmul(inputrec->wall_ewald_zfac,boxs[ZZ],boxs[ZZ]);
818 gmx_pme_send_x(cr,bBS ? boxs : box,x,
819 mdatoms->nChargePerturbed,lambda[efptCOUL],
820 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),step);
822 wallcycle_stop(wcycle,ewcPP_PMESENDX);
826 /* do gridding for pair search */
829 if (graph && bStateChanged)
831 /* Calculate intramolecular shift vectors to make molecules whole */
832 mk_mshift(fplog,graph,fr->ePBC,box,x);
836 box_diag[XX] = box[XX][XX];
837 box_diag[YY] = box[YY][YY];
838 box_diag[ZZ] = box[ZZ][ZZ];
840 wallcycle_start(wcycle,ewcNS);
843 wallcycle_sub_start(wcycle,ewcsNBS_GRID_LOCAL);
844 nbnxn_put_on_grid(nbv->nbs,fr->ePBC,box,
846 0,mdatoms->homenr,-1,fr->cginfo,x,
848 nbv->grp[eintLocal].kernel_type,
849 nbv->grp[eintLocal].nbat);
850 wallcycle_sub_stop(wcycle,ewcsNBS_GRID_LOCAL);
854 wallcycle_sub_start(wcycle,ewcsNBS_GRID_NONLOCAL);
855 nbnxn_put_on_grid_nonlocal(nbv->nbs,domdec_zones(cr->dd),
857 nbv->grp[eintNonlocal].kernel_type,
858 nbv->grp[eintNonlocal].nbat);
859 wallcycle_sub_stop(wcycle,ewcsNBS_GRID_NONLOCAL);
862 if (nbv->ngrp == 1 ||
863 nbv->grp[eintNonlocal].nbat == nbv->grp[eintLocal].nbat)
865 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat,eatAll,
866 nbv->nbs,mdatoms,fr->cginfo);
870 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat,eatLocal,
871 nbv->nbs,mdatoms,fr->cginfo);
872 nbnxn_atomdata_set(nbv->grp[eintNonlocal].nbat,eatAll,
873 nbv->nbs,mdatoms,fr->cginfo);
875 wallcycle_stop(wcycle, ewcNS);
878 /* initialize the GPU atom data and copy shift vector */
883 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
884 nbnxn_cuda_init_atomdata(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
885 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
888 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
889 nbnxn_cuda_upload_shiftvec(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
890 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
893 /* do local pair search */
896 wallcycle_start_nocount(wcycle,ewcNS);
897 wallcycle_sub_start(wcycle,ewcsNBS_SEARCH_LOCAL);
898 nbnxn_make_pairlist(nbv->nbs,nbv->grp[eintLocal].nbat,
901 nbv->min_ci_balanced,
902 &nbv->grp[eintLocal].nbl_lists,
904 nbv->grp[eintLocal].kernel_type,
906 wallcycle_sub_stop(wcycle,ewcsNBS_SEARCH_LOCAL);
910 /* initialize local pair-list on the GPU */
911 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
912 nbv->grp[eintLocal].nbl_lists.nbl[0],
915 wallcycle_stop(wcycle, ewcNS);
919 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
920 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
921 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs,eatLocal,FALSE,x,
922 nbv->grp[eintLocal].nbat);
923 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
924 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
929 wallcycle_start(wcycle,ewcLAUNCH_GPU_NB);
930 /* launch local nonbonded F on GPU */
931 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFNo,
933 wallcycle_stop(wcycle,ewcLAUNCH_GPU_NB);
936 /* Communicate coordinates and sum dipole if necessary +
937 do non-local pair search */
938 if (DOMAINDECOMP(cr))
940 bDiffKernels = (nbv->grp[eintNonlocal].kernel_type !=
941 nbv->grp[eintLocal].kernel_type);
945 /* With GPU+CPU non-bonded calculations we need to copy
946 * the local coordinates to the non-local nbat struct
947 * (in CPU format) as the non-local kernel call also
948 * calculates the local - non-local interactions.
950 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
951 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
952 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs,eatLocal,TRUE,x,
953 nbv->grp[eintNonlocal].nbat);
954 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
955 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
960 wallcycle_start_nocount(wcycle,ewcNS);
961 wallcycle_sub_start(wcycle,ewcsNBS_SEARCH_NONLOCAL);
965 nbnxn_grid_add_simple(nbv->nbs,nbv->grp[eintNonlocal].nbat);
968 nbnxn_make_pairlist(nbv->nbs,nbv->grp[eintNonlocal].nbat,
971 nbv->min_ci_balanced,
972 &nbv->grp[eintNonlocal].nbl_lists,
974 nbv->grp[eintNonlocal].kernel_type,
977 wallcycle_sub_stop(wcycle,ewcsNBS_SEARCH_NONLOCAL);
979 if (nbv->grp[eintNonlocal].kernel_type == nbk8x8x8_CUDA)
981 /* initialize non-local pair-list on the GPU */
982 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
983 nbv->grp[eintNonlocal].nbl_lists.nbl[0],
986 wallcycle_stop(wcycle,ewcNS);
990 wallcycle_start(wcycle,ewcMOVEX);
991 dd_move_x(cr->dd,box,x);
993 /* When we don't need the total dipole we sum it in global_stat */
994 if (bStateChanged && NEED_MUTOT(*inputrec))
996 gmx_sumd(2*DIM,mu,cr);
998 wallcycle_stop(wcycle,ewcMOVEX);
1000 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1001 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
1002 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs,eatNonlocal,FALSE,x,
1003 nbv->grp[eintNonlocal].nbat);
1004 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
1005 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1008 if (bUseGPU && !bDiffKernels)
1010 wallcycle_start(wcycle,ewcLAUNCH_GPU_NB);
1011 /* launch non-local nonbonded F on GPU */
1012 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFNo,
1014 cycles_force += wallcycle_stop(wcycle,ewcLAUNCH_GPU_NB);
1020 /* launch D2H copy-back F */
1021 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1022 if (DOMAINDECOMP(cr) && !bDiffKernels)
1024 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintNonlocal].nbat,
1025 flags, eatNonlocal);
1027 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintLocal].nbat,
1029 cycles_force += wallcycle_stop(wcycle,ewcLAUNCH_GPU_NB);
1032 if (bStateChanged && NEED_MUTOT(*inputrec))
1036 gmx_sumd(2*DIM,mu,cr);
1043 fr->mu_tot[i][j] = mu[i*DIM + j];
1047 if (fr->efep == efepNO)
1049 copy_rvec(fr->mu_tot[0],mu_tot);
1053 for(j=0; j<DIM; j++)
1056 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] +
1057 lambda[efptCOUL]*fr->mu_tot[1][j];
1061 /* Reset energies */
1062 reset_enerdata(&(inputrec->opts),fr,bNS,enerd,MASTER(cr));
1063 clear_rvecs(SHIFTS,fr->fshift);
1065 if (DOMAINDECOMP(cr))
1067 if (!(cr->duty & DUTY_PME))
1069 wallcycle_start(wcycle,ewcPPDURINGPME);
1070 dd_force_flop_start(cr->dd,nrnb);
1074 /* Start the force cycle counter.
1075 * This counter is stopped in do_forcelow_level.
1076 * No parallel communication should occur while this counter is running,
1077 * since that will interfere with the dynamic load balancing.
1079 wallcycle_start(wcycle,ewcFORCE);
1082 /* Reset forces for which the virial is calculated separately:
1083 * PME/Ewald forces if necessary */
1084 if (fr->bF_NoVirSum)
1086 if (flags & GMX_FORCE_VIRIAL)
1088 fr->f_novirsum = fr->f_novirsum_alloc;
1091 clear_rvecs(fr->f_novirsum_n,fr->f_novirsum);
1095 clear_rvecs(homenr,fr->f_novirsum+start);
1100 /* We are not calculating the pressure so we do not need
1101 * a separate array for forces that do not contribute
1108 /* Clear the short- and long-range forces */
1109 clear_rvecs(fr->natoms_force_constr,f);
1110 if(bSepLRF && do_per_step(step,inputrec->nstcalclr))
1112 clear_rvecs(fr->natoms_force_constr,fr->f_twin);
1115 clear_rvec(fr->vir_diag_posres);
1117 if (inputrec->ePull == epullCONSTRAINT)
1119 clear_pull_forces(inputrec->pull);
1122 /* update QMMMrec, if necessary */
1125 update_QMMMrec(cr,fr,x,mdatoms,box,top);
1128 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1130 posres_wrapper(fplog,flags,bSepDVDL,inputrec,nrnb,top,box,x,
1134 /* Compute the bonded and non-bonded energies and optionally forces */
1135 do_force_lowlevel(fplog,step,fr,inputrec,&(top->idef),
1136 cr,nrnb,wcycle,mdatoms,&(inputrec->opts),
1137 x,hist,f, bSepLRF ? fr->f_twin : f,enerd,fcd,mtop,top,fr->born,
1138 &(top->atomtypes),bBornRadii,box,
1139 inputrec->fepvals,lambda,graph,&(top->excls),fr->mu_tot,
1140 flags, &cycles_pme);
1144 if (do_per_step(step,inputrec->nstcalclr))
1146 /* Add the long range forces to the short range forces */
1147 for(i=0; i<fr->natoms_force_constr; i++)
1149 rvec_add(fr->f_twin[i],f[i],f[i]);
1156 /* Maybe we should move this into do_force_lowlevel */
1157 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFYes,
1162 if (!bUseOrEmulGPU || bDiffKernels)
1166 if (DOMAINDECOMP(cr))
1168 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal,
1169 bDiffKernels ? enbvClearFYes : enbvClearFNo,
1179 aloc = eintNonlocal;
1182 /* Add all the non-bonded force to the normal force array.
1183 * This can be split into a local a non-local part when overlapping
1184 * communication with calculation with domain decomposition.
1186 cycles_force += wallcycle_stop(wcycle,ewcFORCE);
1187 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1188 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1189 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs,eatAll,nbv->grp[aloc].nbat,f);
1190 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1191 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1192 wallcycle_start_nocount(wcycle,ewcFORCE);
1194 /* if there are multiple fshift output buffers reduce them */
1195 if ((flags & GMX_FORCE_VIRIAL) &&
1196 nbv->grp[aloc].nbl_lists.nnbl > 1)
1198 nbnxn_atomdata_add_nbat_fshift_to_fshift(nbv->grp[aloc].nbat,
1203 cycles_force += wallcycle_stop(wcycle,ewcFORCE);
1207 do_flood(fplog,cr,x,f,ed,box,step,bNS);
1210 if (bUseOrEmulGPU && !bDiffKernels)
1212 /* wait for non-local forces (or calculate in emulation mode) */
1213 if (DOMAINDECOMP(cr))
1217 wallcycle_start(wcycle,ewcWAIT_GPU_NB_NL);
1218 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1219 nbv->grp[eintNonlocal].nbat,
1221 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1223 cycles_force += wallcycle_stop(wcycle,ewcWAIT_GPU_NB_NL);
1227 wallcycle_start_nocount(wcycle,ewcFORCE);
1228 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFYes,
1230 cycles_force += wallcycle_stop(wcycle,ewcFORCE);
1232 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1233 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1234 /* skip the reduction if there was no non-local work to do */
1235 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1237 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs,eatNonlocal,
1238 nbv->grp[eintNonlocal].nbat,f);
1240 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1241 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1247 /* Communicate the forces */
1250 wallcycle_start(wcycle,ewcMOVEF);
1251 if (DOMAINDECOMP(cr))
1253 dd_move_f(cr->dd,f,fr->fshift);
1254 /* Do we need to communicate the separate force array
1255 * for terms that do not contribute to the single sum virial?
1256 * Position restraints and electric fields do not introduce
1257 * inter-cg forces, only full electrostatics methods do.
1258 * When we do not calculate the virial, fr->f_novirsum = f,
1259 * so we have already communicated these forces.
1261 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1262 (flags & GMX_FORCE_VIRIAL))
1264 dd_move_f(cr->dd,fr->f_novirsum,NULL);
1268 /* We should not update the shift forces here,
1269 * since f_twin is already included in f.
1271 dd_move_f(cr->dd,fr->f_twin,NULL);
1274 wallcycle_stop(wcycle,ewcMOVEF);
1280 /* wait for local forces (or calculate in emulation mode) */
1283 wallcycle_start(wcycle,ewcWAIT_GPU_NB_L);
1284 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1285 nbv->grp[eintLocal].nbat,
1287 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1289 wallcycle_stop(wcycle,ewcWAIT_GPU_NB_L);
1291 /* now clear the GPU outputs while we finish the step on the CPU */
1292 nbnxn_cuda_clear_outputs(nbv->cu_nbv, flags);
1296 wallcycle_start_nocount(wcycle,ewcFORCE);
1297 do_nb_verlet(fr, ic, enerd, flags, eintLocal,
1298 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
1300 wallcycle_stop(wcycle,ewcFORCE);
1302 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1303 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1304 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1306 /* skip the reduction if there was no non-local work to do */
1307 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs,eatLocal,
1308 nbv->grp[eintLocal].nbat,f);
1310 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1311 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1314 if (DOMAINDECOMP(cr))
1316 dd_force_flop_stop(cr->dd,nrnb);
1319 dd_cycles_add(cr->dd,cycles_force-cycles_pme,ddCyclF);
1325 if (IR_ELEC_FIELD(*inputrec))
1327 /* Compute forces due to electric field */
1328 calc_f_el(MASTER(cr) ? field : NULL,
1329 start,homenr,mdatoms->chargeA,x,fr->f_novirsum,
1330 inputrec->ex,inputrec->et,t);
1333 /* If we have NoVirSum forces, but we do not calculate the virial,
1334 * we sum fr->f_novirum=f later.
1336 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1338 wallcycle_start(wcycle,ewcVSITESPREAD);
1339 spread_vsite_f(fplog,vsite,x,f,fr->fshift,FALSE,NULL,nrnb,
1340 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
1341 wallcycle_stop(wcycle,ewcVSITESPREAD);
1345 wallcycle_start(wcycle,ewcVSITESPREAD);
1346 spread_vsite_f(fplog,vsite,x,fr->f_twin,NULL,FALSE,NULL,
1348 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
1349 wallcycle_stop(wcycle,ewcVSITESPREAD);
1353 if (flags & GMX_FORCE_VIRIAL)
1355 /* Calculation of the virial must be done after vsites! */
1356 calc_virial(fplog,mdatoms->start,mdatoms->homenr,x,f,
1357 vir_force,graph,box,nrnb,fr,inputrec->ePBC);
1361 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1363 pull_potential_wrapper(fplog,bSepDVDL,cr,inputrec,box,x,
1364 f,vir_force,mdatoms,enerd,lambda,t);
1367 if (PAR(cr) && !(cr->duty & DUTY_PME))
1369 /* In case of node-splitting, the PP nodes receive the long-range
1370 * forces, virial and energy from the PME nodes here.
1372 pme_receive_force_ener(fplog,bSepDVDL,cr,wcycle,enerd,fr);
1377 post_process_forces(fplog,cr,step,nrnb,wcycle,
1378 top,box,x,f,vir_force,mdatoms,graph,fr,vsite,
1382 /* Sum the potential energy terms from group contributions */
1383 sum_epot(&(inputrec->opts),enerd);
1386 void do_force_cutsGROUP(FILE *fplog,t_commrec *cr,
1387 t_inputrec *inputrec,
1388 gmx_large_int_t step,t_nrnb *nrnb,gmx_wallcycle_t wcycle,
1389 gmx_localtop_t *top,
1391 gmx_groups_t *groups,
1392 matrix box,rvec x[],history_t *hist,
1396 gmx_enerdata_t *enerd,t_fcdata *fcd,
1397 real *lambda,t_graph *graph,
1398 t_forcerec *fr,gmx_vsite_t *vsite,rvec mu_tot,
1399 double t,FILE *field,gmx_edsam_t ed,
1400 gmx_bool bBornRadii,
1406 gmx_bool bSepDVDL,bStateChanged,bNS,bFillGrid,bCalcCGCM,bBS;
1407 gmx_bool bDoLongRangeNS,bDoForces,bDoPotential,bSepLRF;
1408 gmx_bool bDoAdressWF;
1410 rvec vzero,box_diag;
1411 real e,v,dvdlambda[efptNR];
1413 float cycles_pme,cycles_force;
1415 start = mdatoms->start;
1416 homenr = mdatoms->homenr;
1418 bSepDVDL = (fr->bSepDVDL && do_per_step(step,inputrec->nstlog));
1420 clear_mat(vir_force);
1424 pd_cg_range(cr,&cg0,&cg1);
1429 if (DOMAINDECOMP(cr))
1431 cg1 = cr->dd->ncg_tot;
1443 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
1444 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll==FALSE);
1445 /* Should we update the long-range neighborlists at this step? */
1446 bDoLongRangeNS = fr->bTwinRange && bNS;
1447 /* Should we perform the long-range nonbonded evaluation inside the neighborsearching? */
1448 bFillGrid = (bNS && bStateChanged);
1449 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
1450 bDoForces = (flags & GMX_FORCE_FORCES);
1451 bDoPotential = (flags & GMX_FORCE_ENERGY);
1452 bSepLRF = ((inputrec->nstcalclr>1) && bDoForces &&
1453 (flags & GMX_FORCE_SEPLRF) && (flags & GMX_FORCE_DO_LR));
1455 /* should probably move this to the forcerec since it doesn't change */
1456 bDoAdressWF = ((fr->adress_type!=eAdressOff));
1460 update_forcerec(fplog,fr,box);
1462 if (NEED_MUTOT(*inputrec))
1464 /* Calculate total (local) dipole moment in a temporary common array.
1465 * This makes it possible to sum them over nodes faster.
1467 calc_mu(start,homenr,
1468 x,mdatoms->chargeA,mdatoms->chargeB,mdatoms->nChargePerturbed,
1473 if (fr->ePBC != epbcNONE) {
1474 /* Compute shift vectors every step,
1475 * because of pressure coupling or box deformation!
1477 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
1478 calc_shifts(box,fr->shift_vec);
1481 put_charge_groups_in_box(fplog,cg0,cg1,fr->ePBC,box,
1482 &(top->cgs),x,fr->cg_cm);
1483 inc_nrnb(nrnb,eNR_CGCM,homenr);
1484 inc_nrnb(nrnb,eNR_RESETX,cg1-cg0);
1486 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph) {
1487 unshift_self(graph,box,x);
1490 else if (bCalcCGCM) {
1491 calc_cgcm(fplog,cg0,cg1,&(top->cgs),x,fr->cg_cm);
1492 inc_nrnb(nrnb,eNR_CGCM,homenr);
1497 move_cgcm(fplog,cr,fr->cg_cm);
1500 pr_rvecs(debug,0,"cgcm",fr->cg_cm,top->cgs.nr);
1504 if (!(cr->duty & DUTY_PME)) {
1505 /* Send particle coordinates to the pme nodes.
1506 * Since this is only implemented for domain decomposition
1507 * and domain decomposition does not use the graph,
1508 * we do not need to worry about shifting.
1511 wallcycle_start(wcycle,ewcPP_PMESENDX);
1513 bBS = (inputrec->nwall == 2);
1516 svmul(inputrec->wall_ewald_zfac,boxs[ZZ],boxs[ZZ]);
1519 gmx_pme_send_x(cr,bBS ? boxs : box,x,
1520 mdatoms->nChargePerturbed,lambda[efptCOUL],
1521 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),step);
1523 wallcycle_stop(wcycle,ewcPP_PMESENDX);
1525 #endif /* GMX_MPI */
1527 /* Communicate coordinates and sum dipole if necessary */
1530 wallcycle_start(wcycle,ewcMOVEX);
1531 if (DOMAINDECOMP(cr))
1533 dd_move_x(cr->dd,box,x);
1537 move_x(fplog,cr,GMX_LEFT,GMX_RIGHT,x,nrnb);
1539 wallcycle_stop(wcycle,ewcMOVEX);
1542 /* update adress weight beforehand */
1543 if(bStateChanged && bDoAdressWF)
1545 /* need pbc for adress weight calculation with pbc_dx */
1546 set_pbc(&pbc,inputrec->ePBC,box);
1547 if(fr->adress_site == eAdressSITEcog)
1549 update_adress_weights_cog(top->idef.iparams,top->idef.il,x,fr,mdatoms,
1550 inputrec->ePBC==epbcNONE ? NULL : &pbc);
1552 else if (fr->adress_site == eAdressSITEcom)
1554 update_adress_weights_com(fplog,cg0,cg1,&(top->cgs),x,fr,mdatoms,
1555 inputrec->ePBC==epbcNONE ? NULL : &pbc);
1557 else if (fr->adress_site == eAdressSITEatomatom){
1558 update_adress_weights_atom_per_atom(cg0,cg1,&(top->cgs),x,fr,mdatoms,
1559 inputrec->ePBC==epbcNONE ? NULL : &pbc);
1563 update_adress_weights_atom(cg0,cg1,&(top->cgs),x,fr,mdatoms,
1564 inputrec->ePBC==epbcNONE ? NULL : &pbc);
1568 if (NEED_MUTOT(*inputrec))
1575 gmx_sumd(2*DIM,mu,cr);
1581 fr->mu_tot[i][j] = mu[i*DIM + j];
1585 if (fr->efep == efepNO)
1587 copy_rvec(fr->mu_tot[0],mu_tot);
1591 for(j=0; j<DIM; j++)
1594 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] + lambda[efptCOUL]*fr->mu_tot[1][j];
1599 /* Reset energies */
1600 reset_enerdata(&(inputrec->opts),fr,bNS,enerd,MASTER(cr));
1601 clear_rvecs(SHIFTS,fr->fshift);
1605 wallcycle_start(wcycle,ewcNS);
1607 if (graph && bStateChanged)
1609 /* Calculate intramolecular shift vectors to make molecules whole */
1610 mk_mshift(fplog,graph,fr->ePBC,box,x);
1613 /* Do the actual neighbour searching and if twin range electrostatics
1614 * also do the calculation of long range forces and energies.
1616 for (i=0;i<efptNR;i++) {dvdlambda[i] = 0;}
1618 groups,&(inputrec->opts),top,mdatoms,
1619 cr,nrnb,lambda,dvdlambda,&enerd->grpp,bFillGrid,
1623 fprintf(fplog,sepdvdlformat,"LR non-bonded",0.0,dvdlambda);
1625 enerd->dvdl_lin[efptVDW] += dvdlambda[efptVDW];
1626 enerd->dvdl_lin[efptCOUL] += dvdlambda[efptCOUL];
1628 wallcycle_stop(wcycle,ewcNS);
1631 if (inputrec->implicit_solvent && bNS)
1633 make_gb_nblist(cr,inputrec->gb_algorithm,inputrec->rlist,
1634 x,box,fr,&top->idef,graph,fr->born);
1637 if (DOMAINDECOMP(cr))
1639 if (!(cr->duty & DUTY_PME))
1641 wallcycle_start(wcycle,ewcPPDURINGPME);
1642 dd_force_flop_start(cr->dd,nrnb);
1648 /* Enforced rotation has its own cycle counter that starts after the collective
1649 * coordinates have been communicated. It is added to ddCyclF to allow
1650 * for proper load-balancing */
1651 wallcycle_start(wcycle,ewcROT);
1652 do_rotation(cr,inputrec,box,x,t,step,wcycle,bNS);
1653 wallcycle_stop(wcycle,ewcROT);
1656 /* Start the force cycle counter.
1657 * This counter is stopped in do_forcelow_level.
1658 * No parallel communication should occur while this counter is running,
1659 * since that will interfere with the dynamic load balancing.
1661 wallcycle_start(wcycle,ewcFORCE);
1665 /* Reset forces for which the virial is calculated separately:
1666 * PME/Ewald forces if necessary */
1667 if (fr->bF_NoVirSum)
1669 if (flags & GMX_FORCE_VIRIAL)
1671 fr->f_novirsum = fr->f_novirsum_alloc;
1674 clear_rvecs(fr->f_novirsum_n,fr->f_novirsum);
1678 clear_rvecs(homenr,fr->f_novirsum+start);
1683 /* We are not calculating the pressure so we do not need
1684 * a separate array for forces that do not contribute
1691 /* Clear the short- and long-range forces */
1692 clear_rvecs(fr->natoms_force_constr,f);
1693 if(bSepLRF && do_per_step(step,inputrec->nstcalclr))
1695 clear_rvecs(fr->natoms_force_constr,fr->f_twin);
1698 clear_rvec(fr->vir_diag_posres);
1700 if (inputrec->ePull == epullCONSTRAINT)
1702 clear_pull_forces(inputrec->pull);
1705 /* update QMMMrec, if necessary */
1708 update_QMMMrec(cr,fr,x,mdatoms,box,top);
1711 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1713 posres_wrapper(fplog,flags,bSepDVDL,inputrec,nrnb,top,box,x,
1717 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_FBPOSRES].nr > 0)
1719 /* Flat-bottomed position restraints always require full pbc */
1720 if(!(bStateChanged && bDoAdressWF))
1722 set_pbc(&pbc,inputrec->ePBC,box);
1724 v = fbposres(top->idef.il[F_FBPOSRES].nr,top->idef.il[F_FBPOSRES].iatoms,
1725 top->idef.iparams_fbposres,
1726 (const rvec*)x,fr->f_novirsum,fr->vir_diag_posres,
1727 inputrec->ePBC==epbcNONE ? NULL : &pbc,
1728 fr->rc_scaling,fr->ePBC,fr->posres_com);
1729 enerd->term[F_FBPOSRES] += v;
1730 inc_nrnb(nrnb,eNR_FBPOSRES,top->idef.il[F_FBPOSRES].nr/2);
1733 /* Compute the bonded and non-bonded energies and optionally forces */
1734 do_force_lowlevel(fplog,step,fr,inputrec,&(top->idef),
1735 cr,nrnb,wcycle,mdatoms,&(inputrec->opts),
1736 x,hist,f, bSepLRF ? fr->f_twin : f,enerd,fcd,mtop,top,fr->born,
1737 &(top->atomtypes),bBornRadii,box,
1738 inputrec->fepvals,lambda,
1739 graph,&(top->excls),fr->mu_tot,
1745 if (do_per_step(step,inputrec->nstcalclr))
1747 /* Add the long range forces to the short range forces */
1748 for(i=0; i<fr->natoms_force_constr; i++)
1750 rvec_add(fr->f_twin[i],f[i],f[i]);
1755 cycles_force = wallcycle_stop(wcycle,ewcFORCE);
1759 do_flood(fplog,cr,x,f,ed,box,step,bNS);
1762 if (DOMAINDECOMP(cr))
1764 dd_force_flop_stop(cr->dd,nrnb);
1767 dd_cycles_add(cr->dd,cycles_force-cycles_pme,ddCyclF);
1773 if (IR_ELEC_FIELD(*inputrec))
1775 /* Compute forces due to electric field */
1776 calc_f_el(MASTER(cr) ? field : NULL,
1777 start,homenr,mdatoms->chargeA,x,fr->f_novirsum,
1778 inputrec->ex,inputrec->et,t);
1781 if (bDoAdressWF && fr->adress_icor == eAdressICThermoForce)
1783 /* Compute thermodynamic force in hybrid AdResS region */
1784 adress_thermo_force(start,homenr,&(top->cgs),x,fr->f_novirsum,fr,mdatoms,
1785 inputrec->ePBC==epbcNONE ? NULL : &pbc);
1788 /* Communicate the forces */
1791 wallcycle_start(wcycle,ewcMOVEF);
1792 if (DOMAINDECOMP(cr))
1794 dd_move_f(cr->dd,f,fr->fshift);
1795 /* Do we need to communicate the separate force array
1796 * for terms that do not contribute to the single sum virial?
1797 * Position restraints and electric fields do not introduce
1798 * inter-cg forces, only full electrostatics methods do.
1799 * When we do not calculate the virial, fr->f_novirsum = f,
1800 * so we have already communicated these forces.
1802 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1803 (flags & GMX_FORCE_VIRIAL))
1805 dd_move_f(cr->dd,fr->f_novirsum,NULL);
1809 /* We should not update the shift forces here,
1810 * since f_twin is already included in f.
1812 dd_move_f(cr->dd,fr->f_twin,NULL);
1817 pd_move_f(cr,f,nrnb);
1820 pd_move_f(cr,fr->f_twin,nrnb);
1823 wallcycle_stop(wcycle,ewcMOVEF);
1826 /* If we have NoVirSum forces, but we do not calculate the virial,
1827 * we sum fr->f_novirum=f later.
1829 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1831 wallcycle_start(wcycle,ewcVSITESPREAD);
1832 spread_vsite_f(fplog,vsite,x,f,fr->fshift,FALSE,NULL,nrnb,
1833 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
1834 wallcycle_stop(wcycle,ewcVSITESPREAD);
1838 wallcycle_start(wcycle,ewcVSITESPREAD);
1839 spread_vsite_f(fplog,vsite,x,fr->f_twin,NULL,FALSE,NULL,
1841 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
1842 wallcycle_stop(wcycle,ewcVSITESPREAD);
1846 if (flags & GMX_FORCE_VIRIAL)
1848 /* Calculation of the virial must be done after vsites! */
1849 calc_virial(fplog,mdatoms->start,mdatoms->homenr,x,f,
1850 vir_force,graph,box,nrnb,fr,inputrec->ePBC);
1854 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1856 pull_potential_wrapper(fplog,bSepDVDL,cr,inputrec,box,x,
1857 f,vir_force,mdatoms,enerd,lambda,t);
1860 /* Add the forces from enforced rotation potentials (if any) */
1863 wallcycle_start(wcycle,ewcROTadd);
1864 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr,step,t);
1865 wallcycle_stop(wcycle,ewcROTadd);
1868 if (PAR(cr) && !(cr->duty & DUTY_PME))
1870 /* In case of node-splitting, the PP nodes receive the long-range
1871 * forces, virial and energy from the PME nodes here.
1873 pme_receive_force_ener(fplog,bSepDVDL,cr,wcycle,enerd,fr);
1878 post_process_forces(fplog,cr,step,nrnb,wcycle,
1879 top,box,x,f,vir_force,mdatoms,graph,fr,vsite,
1883 /* Sum the potential energy terms from group contributions */
1884 sum_epot(&(inputrec->opts),enerd);
1887 void do_force(FILE *fplog,t_commrec *cr,
1888 t_inputrec *inputrec,
1889 gmx_large_int_t step,t_nrnb *nrnb,gmx_wallcycle_t wcycle,
1890 gmx_localtop_t *top,
1892 gmx_groups_t *groups,
1893 matrix box,rvec x[],history_t *hist,
1897 gmx_enerdata_t *enerd,t_fcdata *fcd,
1898 real *lambda,t_graph *graph,
1900 gmx_vsite_t *vsite,rvec mu_tot,
1901 double t,FILE *field,gmx_edsam_t ed,
1902 gmx_bool bBornRadii,
1905 /* modify force flag if not doing nonbonded */
1906 if (!fr->bNonbonded)
1908 flags &= ~GMX_FORCE_NONBONDED;
1911 switch (inputrec->cutoff_scheme)
1914 do_force_cutsVERLET(fplog, cr, inputrec,
1930 do_force_cutsGROUP(fplog, cr, inputrec,
1945 gmx_incons("Invalid cut-off scheme passed!");
1950 void do_constrain_first(FILE *fplog,gmx_constr_t constr,
1951 t_inputrec *ir,t_mdatoms *md,
1952 t_state *state,rvec *f,
1953 t_graph *graph,t_commrec *cr,t_nrnb *nrnb,
1954 t_forcerec *fr, gmx_localtop_t *top, tensor shake_vir)
1957 gmx_large_int_t step;
1958 real dt=ir->delta_t;
1962 snew(savex,state->natoms);
1965 end = md->homenr + start;
1968 fprintf(debug,"vcm: start=%d, homenr=%d, end=%d\n",
1969 start,md->homenr,end);
1970 /* Do a first constrain to reset particles... */
1971 step = ir->init_step;
1974 char buf[STEPSTRSIZE];
1975 fprintf(fplog,"\nConstraining the starting coordinates (step %s)\n",
1976 gmx_step_str(step,buf));
1980 /* constrain the current position */
1981 constrain(NULL,TRUE,FALSE,constr,&(top->idef),
1982 ir,NULL,cr,step,0,md,
1983 state->x,state->x,NULL,
1984 fr->bMolPBC,state->box,
1985 state->lambda[efptBONDED],&dvdl_dum,
1986 NULL,NULL,nrnb,econqCoord,
1987 ir->epc==epcMTTK,state->veta,state->veta);
1990 /* constrain the inital velocity, and save it */
1991 /* also may be useful if we need the ekin from the halfstep for velocity verlet */
1992 /* might not yet treat veta correctly */
1993 constrain(NULL,TRUE,FALSE,constr,&(top->idef),
1994 ir,NULL,cr,step,0,md,
1995 state->x,state->v,state->v,
1996 fr->bMolPBC,state->box,
1997 state->lambda[efptBONDED],&dvdl_dum,
1998 NULL,NULL,nrnb,econqVeloc,
1999 ir->epc==epcMTTK,state->veta,state->veta);
2001 /* constrain the inital velocities at t-dt/2 */
2002 if (EI_STATE_VELOCITY(ir->eI) && ir->eI!=eiVV)
2004 for(i=start; (i<end); i++)
2006 for(m=0; (m<DIM); m++)
2008 /* Reverse the velocity */
2009 state->v[i][m] = -state->v[i][m];
2010 /* Store the position at t-dt in buf */
2011 savex[i][m] = state->x[i][m] + dt*state->v[i][m];
2014 /* Shake the positions at t=-dt with the positions at t=0
2015 * as reference coordinates.
2019 char buf[STEPSTRSIZE];
2020 fprintf(fplog,"\nConstraining the coordinates at t0-dt (step %s)\n",
2021 gmx_step_str(step,buf));
2024 constrain(NULL,TRUE,FALSE,constr,&(top->idef),
2025 ir,NULL,cr,step,-1,md,
2026 state->x,savex,NULL,
2027 fr->bMolPBC,state->box,
2028 state->lambda[efptBONDED],&dvdl_dum,
2029 state->v,NULL,nrnb,econqCoord,
2030 ir->epc==epcMTTK,state->veta,state->veta);
2032 for(i=start; i<end; i++) {
2033 for(m=0; m<DIM; m++) {
2034 /* Re-reverse the velocities */
2035 state->v[i][m] = -state->v[i][m];
2042 void calc_enervirdiff(FILE *fplog,int eDispCorr,t_forcerec *fr)
2044 double eners[2],virs[2],enersum,virsum,y0,f,g,h;
2045 double r0,r1,r,rc3,rc9,ea,eb,ec,pa,pb,pc,pd;
2046 double invscale,invscale2,invscale3;
2047 int ri0,ri1,ri,i,offstart,offset;
2048 real scale,*vdwtab,tabfactor,tmp;
2050 fr->enershiftsix = 0;
2051 fr->enershifttwelve = 0;
2052 fr->enerdiffsix = 0;
2053 fr->enerdifftwelve = 0;
2055 fr->virdifftwelve = 0;
2057 if (eDispCorr != edispcNO) {
2058 for(i=0; i<2; i++) {
2062 if ((fr->vdwtype == evdwSWITCH) || (fr->vdwtype == evdwSHIFT)) {
2063 if (fr->rvdw_switch == 0)
2065 "With dispersion correction rvdw-switch can not be zero "
2066 "for vdw-type = %s",evdw_names[fr->vdwtype]);
2068 scale = fr->nblists[0].table_elec_vdw.scale;
2069 vdwtab = fr->nblists[0].table_vdw.data;
2071 /* Round the cut-offs to exact table values for precision */
2072 ri0 = floor(fr->rvdw_switch*scale);
2073 ri1 = ceil(fr->rvdw*scale);
2079 if (fr->vdwtype == evdwSHIFT)
2081 /* Determine the constant energy shift below rvdw_switch.
2082 * Table has a scale factor since we have scaled it down to compensate
2083 * for scaling-up c6/c12 with the derivative factors to save flops in analytical kernels.
2085 fr->enershiftsix = (real)(-1.0/(rc3*rc3)) - 6.0*vdwtab[8*ri0];
2086 fr->enershifttwelve = (real)( 1.0/(rc9*rc3)) - 12.0*vdwtab[8*ri0 + 4];
2088 /* Add the constant part from 0 to rvdw_switch.
2089 * This integration from 0 to rvdw_switch overcounts the number
2090 * of interactions by 1, as it also counts the self interaction.
2091 * We will correct for this later.
2093 eners[0] += 4.0*M_PI*fr->enershiftsix*rc3/3.0;
2094 eners[1] += 4.0*M_PI*fr->enershifttwelve*rc3/3.0;
2096 invscale = 1.0/(scale);
2097 invscale2 = invscale*invscale;
2098 invscale3 = invscale*invscale2;
2100 /* following summation derived from cubic spline definition,
2101 Numerical Recipies in C, second edition, p. 113-116. Exact
2102 for the cubic spline. We first calculate the negative of
2103 the energy from rvdw to rvdw_switch, assuming that g(r)=1,
2104 and then add the more standard, abrupt cutoff correction to
2105 that result, yielding the long-range correction for a
2106 switched function. We perform both the pressure and energy
2107 loops at the same time for simplicity, as the computational
2111 enersum = 0.0; virsum = 0.0;
2115 /* Since the dispersion table has been scaled down a factor 6.0 and the repulsion
2116 * a factor 12.0 to compensate for the c6/c12 parameters inside nbfp[] being scaled
2117 * up (to save flops in kernels), we need to correct for this.
2126 for (ri=ri0; ri<ri1; ri++) {
2129 eb = 2.0*invscale2*r;
2133 pb = 3.0*invscale2*r;
2134 pc = 3.0*invscale*r*r;
2137 /* this "8" is from the packing in the vdwtab array - perhaps should be #define'ed? */
2138 offset = 8*ri + offstart;
2139 y0 = vdwtab[offset];
2140 f = vdwtab[offset+1];
2141 g = vdwtab[offset+2];
2142 h = vdwtab[offset+3];
2144 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);
2145 virsum += f*(pa/4 + pb/3 + pc/2 + pd) + 2*g*(pa/5 + pb/4 + pc/3 + pd/2) + 3*h*(pa/6 + pb/5 + pc/4 + pd/3);
2148 enersum *= 4.0*M_PI*tabfactor;
2149 virsum *= 4.0*M_PI*tabfactor;
2150 eners[i] -= enersum;
2154 /* now add the correction for rvdw_switch to infinity */
2155 eners[0] += -4.0*M_PI/(3.0*rc3);
2156 eners[1] += 4.0*M_PI/(9.0*rc9);
2157 virs[0] += 8.0*M_PI/rc3;
2158 virs[1] += -16.0*M_PI/(3.0*rc9);
2160 else if ((fr->vdwtype == evdwCUT) || (fr->vdwtype == evdwUSER)) {
2161 if (fr->vdwtype == evdwUSER && fplog)
2163 "WARNING: using dispersion correction with user tables\n");
2164 rc3 = fr->rvdw*fr->rvdw*fr->rvdw;
2166 /* Contribution beyond the cut-off */
2167 eners[0] += -4.0*M_PI/(3.0*rc3);
2168 eners[1] += 4.0*M_PI/(9.0*rc9);
2169 if (fr->vdw_modifier==eintmodPOTSHIFT) {
2170 /* Contribution within the cut-off */
2171 eners[0] += -4.0*M_PI/(3.0*rc3);
2172 eners[1] += 4.0*M_PI/(3.0*rc9);
2174 /* Contribution beyond the cut-off */
2175 virs[0] += 8.0*M_PI/rc3;
2176 virs[1] += -16.0*M_PI/(3.0*rc9);
2179 "Dispersion correction is not implemented for vdw-type = %s",
2180 evdw_names[fr->vdwtype]);
2182 fr->enerdiffsix = eners[0];
2183 fr->enerdifftwelve = eners[1];
2184 /* The 0.5 is due to the Gromacs definition of the virial */
2185 fr->virdiffsix = 0.5*virs[0];
2186 fr->virdifftwelve = 0.5*virs[1];
2190 void calc_dispcorr(FILE *fplog,t_inputrec *ir,t_forcerec *fr,
2191 gmx_large_int_t step,int natoms,
2192 matrix box,real lambda,tensor pres,tensor virial,
2193 real *prescorr, real *enercorr, real *dvdlcorr)
2195 gmx_bool bCorrAll,bCorrPres;
2196 real dvdlambda,invvol,dens,ninter,avcsix,avctwelve,enerdiff,svir=0,spres=0;
2206 if (ir->eDispCorr != edispcNO) {
2207 bCorrAll = (ir->eDispCorr == edispcAllEner ||
2208 ir->eDispCorr == edispcAllEnerPres);
2209 bCorrPres = (ir->eDispCorr == edispcEnerPres ||
2210 ir->eDispCorr == edispcAllEnerPres);
2212 invvol = 1/det(box);
2215 /* Only correct for the interactions with the inserted molecule */
2216 dens = (natoms - fr->n_tpi)*invvol;
2221 dens = natoms*invvol;
2222 ninter = 0.5*natoms;
2225 if (ir->efep == efepNO)
2227 avcsix = fr->avcsix[0];
2228 avctwelve = fr->avctwelve[0];
2232 avcsix = (1 - lambda)*fr->avcsix[0] + lambda*fr->avcsix[1];
2233 avctwelve = (1 - lambda)*fr->avctwelve[0] + lambda*fr->avctwelve[1];
2236 enerdiff = ninter*(dens*fr->enerdiffsix - fr->enershiftsix);
2237 *enercorr += avcsix*enerdiff;
2239 if (ir->efep != efepNO)
2241 dvdlambda += (fr->avcsix[1] - fr->avcsix[0])*enerdiff;
2245 enerdiff = ninter*(dens*fr->enerdifftwelve - fr->enershifttwelve);
2246 *enercorr += avctwelve*enerdiff;
2247 if (fr->efep != efepNO)
2249 dvdlambda += (fr->avctwelve[1] - fr->avctwelve[0])*enerdiff;
2255 svir = ninter*dens*avcsix*fr->virdiffsix/3.0;
2256 if (ir->eDispCorr == edispcAllEnerPres)
2258 svir += ninter*dens*avctwelve*fr->virdifftwelve/3.0;
2260 /* The factor 2 is because of the Gromacs virial definition */
2261 spres = -2.0*invvol*svir*PRESFAC;
2263 for(m=0; m<DIM; m++) {
2264 virial[m][m] += svir;
2265 pres[m][m] += spres;
2270 /* Can't currently control when it prints, for now, just print when degugging */
2274 fprintf(debug,"Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
2280 "Long Range LJ corr.: Epot %10g, Pres: %10g, Vir: %10g\n",
2281 *enercorr,spres,svir);
2285 fprintf(debug,"Long Range LJ corr.: Epot %10g\n",*enercorr);
2289 if (fr->bSepDVDL && do_per_step(step,ir->nstlog))
2291 fprintf(fplog,sepdvdlformat,"Dispersion correction",
2292 *enercorr,dvdlambda);
2294 if (fr->efep != efepNO)
2296 *dvdlcorr += dvdlambda;
2301 void do_pbc_first(FILE *fplog,matrix box,t_forcerec *fr,
2302 t_graph *graph,rvec x[])
2305 fprintf(fplog,"Removing pbc first time\n");
2306 calc_shifts(box,fr->shift_vec);
2308 mk_mshift(fplog,graph,fr->ePBC,box,x);
2310 p_graph(debug,"do_pbc_first 1",graph);
2311 shift_self(graph,box,x);
2312 /* By doing an extra mk_mshift the molecules that are broken
2313 * because they were e.g. imported from another software
2314 * will be made whole again. Such are the healing powers
2317 mk_mshift(fplog,graph,fr->ePBC,box,x);
2319 p_graph(debug,"do_pbc_first 2",graph);
2322 fprintf(fplog,"Done rmpbc\n");
2325 static void low_do_pbc_mtop(FILE *fplog,int ePBC,matrix box,
2326 gmx_mtop_t *mtop,rvec x[],
2331 gmx_molblock_t *molb;
2333 if (bFirst && fplog)
2334 fprintf(fplog,"Removing pbc first time\n");
2338 for(mb=0; mb<mtop->nmolblock; mb++) {
2339 molb = &mtop->molblock[mb];
2340 if (molb->natoms_mol == 1 ||
2341 (!bFirst && mtop->moltype[molb->type].cgs.nr == 1)) {
2342 /* Just one atom or charge group in the molecule, no PBC required */
2343 as += molb->nmol*molb->natoms_mol;
2345 /* Pass NULL iso fplog to avoid graph prints for each molecule type */
2346 mk_graph_ilist(NULL,mtop->moltype[molb->type].ilist,
2347 0,molb->natoms_mol,FALSE,FALSE,graph);
2349 for(mol=0; mol<molb->nmol; mol++) {
2350 mk_mshift(fplog,graph,ePBC,box,x+as);
2352 shift_self(graph,box,x+as);
2353 /* The molecule is whole now.
2354 * We don't need the second mk_mshift call as in do_pbc_first,
2355 * since we no longer need this graph.
2358 as += molb->natoms_mol;
2366 void do_pbc_first_mtop(FILE *fplog,int ePBC,matrix box,
2367 gmx_mtop_t *mtop,rvec x[])
2369 low_do_pbc_mtop(fplog,ePBC,box,mtop,x,TRUE);
2372 void do_pbc_mtop(FILE *fplog,int ePBC,matrix box,
2373 gmx_mtop_t *mtop,rvec x[])
2375 low_do_pbc_mtop(fplog,ePBC,box,mtop,x,FALSE);
2378 void finish_run(FILE *fplog,t_commrec *cr,const char *confout,
2379 t_inputrec *inputrec,
2380 t_nrnb nrnb[],gmx_wallcycle_t wcycle,
2381 gmx_runtime_t *runtime,
2382 wallclock_gpu_t *gputimes,
2384 gmx_bool bWriteStat)
2387 t_nrnb *nrnb_tot=NULL;
2391 wallcycle_sum(cr,wcycle);
2397 MPI_Allreduce(nrnb->n,nrnb_tot->n,eNRNB,MPI_DOUBLE,MPI_SUM,
2398 cr->mpi_comm_mysim);
2406 #if defined(GMX_MPI) && !defined(GMX_THREAD_MPI)
2409 /* reduce nodetime over all MPI processes in the current simulation */
2411 MPI_Allreduce(&runtime->proctime,&sum,1,MPI_DOUBLE,MPI_SUM,
2412 cr->mpi_comm_mysim);
2413 runtime->proctime = sum;
2419 print_flop(fplog,nrnb_tot,&nbfs,&mflop);
2426 if ((cr->duty & DUTY_PP) && DOMAINDECOMP(cr))
2428 print_dd_statistics(cr,inputrec,fplog);
2440 snew(nrnb_all,cr->nnodes);
2441 nrnb_all[0] = *nrnb;
2442 for(s=1; s<cr->nnodes; s++)
2444 MPI_Recv(nrnb_all[s].n,eNRNB,MPI_DOUBLE,s,0,
2445 cr->mpi_comm_mysim,&stat);
2447 pr_load(fplog,cr,nrnb_all);
2452 MPI_Send(nrnb->n,eNRNB,MPI_DOUBLE,MASTERRANK(cr),0,
2453 cr->mpi_comm_mysim);
2460 wallcycle_print(fplog,cr->nnodes,cr->npmenodes,runtime->realtime,
2463 if (EI_DYNAMICS(inputrec->eI))
2465 delta_t = inputrec->delta_t;
2474 print_perf(fplog,runtime->proctime,runtime->realtime,
2475 cr->nnodes-cr->npmenodes,
2476 runtime->nsteps_done,delta_t,nbfs,mflop,
2481 print_perf(stderr,runtime->proctime,runtime->realtime,
2482 cr->nnodes-cr->npmenodes,
2483 runtime->nsteps_done,delta_t,nbfs,mflop,
2489 extern void initialize_lambdas(FILE *fplog,t_inputrec *ir,int *fep_state,real *lambda,double *lam0)
2491 /* this function works, but could probably use a logic rewrite to keep all the different
2492 types of efep straight. */
2495 t_lambda *fep = ir->fepvals;
2497 if ((ir->efep==efepNO) && (ir->bSimTemp == FALSE)) {
2498 for (i=0;i<efptNR;i++) {
2507 *fep_state = fep->init_fep_state; /* this might overwrite the checkpoint
2508 if checkpoint is set -- a kludge is in for now
2510 for (i=0;i<efptNR;i++)
2512 /* overwrite lambda state with init_lambda for now for backwards compatibility */
2513 if (fep->init_lambda>=0) /* if it's -1, it was never initializd */
2515 lambda[i] = fep->init_lambda;
2517 lam0[i] = lambda[i];
2522 lambda[i] = fep->all_lambda[i][*fep_state];
2524 lam0[i] = lambda[i];
2529 /* need to rescale control temperatures to match current state */
2530 for (i=0;i<ir->opts.ngtc;i++) {
2531 if (ir->opts.ref_t[i] > 0) {
2532 ir->opts.ref_t[i] = ir->simtempvals->temperatures[*fep_state];
2538 /* Send to the log the information on the current lambdas */
2541 fprintf(fplog,"Initial vector of lambda components:[ ");
2542 for (i=0;i<efptNR;i++)
2544 fprintf(fplog,"%10.4f ",lambda[i]);
2546 fprintf(fplog,"]\n");
2552 void init_md(FILE *fplog,
2553 t_commrec *cr,t_inputrec *ir,const output_env_t oenv,
2554 double *t,double *t0,
2555 real *lambda, int *fep_state, double *lam0,
2556 t_nrnb *nrnb,gmx_mtop_t *mtop,
2558 int nfile,const t_filenm fnm[],
2559 gmx_mdoutf_t **outf,t_mdebin **mdebin,
2560 tensor force_vir,tensor shake_vir,rvec mu_tot,
2561 gmx_bool *bSimAnn,t_vcm **vcm, t_state *state, unsigned long Flags)
2566 /* Initial values */
2567 *t = *t0 = ir->init_t;
2570 for(i=0;i<ir->opts.ngtc;i++)
2572 /* set bSimAnn if any group is being annealed */
2573 if(ir->opts.annealing[i]!=eannNO)
2580 update_annealing_target_temp(&(ir->opts),ir->init_t);
2583 /* Initialize lambda variables */
2584 initialize_lambdas(fplog,ir,fep_state,lambda,lam0);
2588 *upd = init_update(fplog,ir);
2594 *vcm = init_vcm(fplog,&mtop->groups,ir);
2597 if (EI_DYNAMICS(ir->eI) && !(Flags & MD_APPENDFILES))
2599 if (ir->etc == etcBERENDSEN)
2601 please_cite(fplog,"Berendsen84a");
2603 if (ir->etc == etcVRESCALE)
2605 please_cite(fplog,"Bussi2007a");
2613 *outf = init_mdoutf(nfile,fnm,Flags,cr,ir,oenv);
2615 *mdebin = init_mdebin((Flags & MD_APPENDFILES) ? NULL : (*outf)->fp_ene,
2616 mtop,ir, (*outf)->fp_dhdl);
2621 please_cite(fplog,"Fritsch12");
2622 please_cite(fplog,"Junghans10");
2624 /* Initiate variables */
2625 clear_mat(force_vir);
2626 clear_mat(shake_vir);