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"
85 #include "mpelogging.h"
88 #include "gmx_wallcycle.h"
90 #include "nbnxn_atomdata.h"
91 #include "nbnxn_search.h"
92 #include "nbnxn_kernels/nbnxn_kernel_ref.h"
93 #include "nbnxn_kernels/nbnxn_kernel_x86_simd128.h"
94 #include "nbnxn_kernels/nbnxn_kernel_x86_simd256.h"
95 #include "nbnxn_kernels/nbnxn_kernel_gpu_ref.h"
100 #ifdef GMX_THREAD_MPI
107 #include "nbnxn_cuda_data_mgmt.h"
108 #include "nbnxn_cuda/nbnxn_cuda.h"
111 typedef struct gmx_timeprint {
116 /* Portable version of ctime_r implemented in src/gmxlib/string2.c, but we do not want it declared in public installed headers */
118 gmx_ctime_r(const time_t *clock,char *buf, int n);
124 #ifdef HAVE_GETTIMEOFDAY
128 gettimeofday(&t,NULL);
130 seconds = (double) t.tv_sec + 1e-6*(double)t.tv_usec;
136 seconds = time(NULL);
143 #define difftime(end,start) ((double)(end)-(double)(start))
145 void print_time(FILE *out,gmx_runtime_t *runtime,gmx_large_int_t step,
146 t_inputrec *ir, t_commrec *cr)
149 char timebuf[STRLEN];
153 #ifndef GMX_THREAD_MPI
159 fprintf(out,"step %s",gmx_step_str(step,buf));
160 if ((step >= ir->nstlist))
162 runtime->last = gmx_gettime();
163 dt = difftime(runtime->last,runtime->real);
164 runtime->time_per_step = dt/(step - ir->init_step + 1);
166 dt = (ir->nsteps + ir->init_step - step)*runtime->time_per_step;
172 finish = (time_t) (runtime->last + dt);
173 gmx_ctime_r(&finish,timebuf,STRLEN);
174 sprintf(buf,"%s",timebuf);
175 buf[strlen(buf)-1]='\0';
176 fprintf(out,", will finish %s",buf);
179 fprintf(out,", remaining runtime: %5d s ",(int)dt);
183 fprintf(out," performance: %.1f ns/day ",
184 ir->delta_t/1000*24*60*60/runtime->time_per_step);
187 #ifndef GMX_THREAD_MPI
201 static double set_proctime(gmx_runtime_t *runtime)
207 prev = runtime->proc;
208 runtime->proc = dclock();
210 diff = runtime->proc - prev;
214 prev = runtime->proc;
215 runtime->proc = clock();
217 diff = (double)(runtime->proc - prev)/(double)CLOCKS_PER_SEC;
221 /* The counter has probably looped, ignore this data */
228 void runtime_start(gmx_runtime_t *runtime)
230 runtime->real = gmx_gettime();
232 set_proctime(runtime);
233 runtime->realtime = 0;
234 runtime->proctime = 0;
236 runtime->time_per_step = 0;
239 void runtime_end(gmx_runtime_t *runtime)
245 runtime->proctime += set_proctime(runtime);
246 runtime->realtime = now - runtime->real;
250 void runtime_upd_proc(gmx_runtime_t *runtime)
252 runtime->proctime += set_proctime(runtime);
255 void print_date_and_time(FILE *fplog,int nodeid,const char *title,
256 const gmx_runtime_t *runtime)
259 char timebuf[STRLEN];
260 char time_string[STRLEN];
267 tmptime = (time_t) runtime->real;
268 gmx_ctime_r(&tmptime,timebuf,STRLEN);
272 tmptime = (time_t) gmx_gettime();
273 gmx_ctime_r(&tmptime,timebuf,STRLEN);
275 for(i=0; timebuf[i]>=' '; i++)
277 time_string[i]=timebuf[i];
281 fprintf(fplog,"%s on node %d %s\n",title,nodeid,time_string);
285 static void sum_forces(int start,int end,rvec f[],rvec flr[])
290 pr_rvecs(debug,0,"fsr",f+start,end-start);
291 pr_rvecs(debug,0,"flr",flr+start,end-start);
293 for(i=start; (i<end); i++)
294 rvec_inc(f[i],flr[i]);
298 * calc_f_el calculates forces due to an electric field.
300 * force is kJ mol^-1 nm^-1 = e * kJ mol^-1 nm^-1 / e
302 * Et[] contains the parameters for the time dependent
303 * part of the field (not yet used).
304 * Ex[] contains the parameters for
305 * the spatial dependent part of the field. You can have cool periodic
306 * fields in principle, but only a constant field is supported
308 * The function should return the energy due to the electric field
309 * (if any) but for now returns 0.
312 * There can be problems with the virial.
313 * Since the field is not self-consistent this is unavoidable.
314 * For neutral molecules the virial is correct within this approximation.
315 * For neutral systems with many charged molecules the error is small.
316 * But for systems with a net charge or a few charged molecules
317 * the error can be significant when the field is high.
318 * Solution: implement a self-consitent electric field into PME.
320 static void calc_f_el(FILE *fp,int start,int homenr,
321 real charge[],rvec x[],rvec f[],
322 t_cosines Ex[],t_cosines Et[],double t)
328 for(m=0; (m<DIM); m++)
335 Ext[m] = cos(Et[m].a[0]*(t-t0))*exp(-sqr(t-t0)/(2.0*sqr(Et[m].a[2])));
339 Ext[m] = cos(Et[m].a[0]*t);
348 /* Convert the field strength from V/nm to MD-units */
349 Ext[m] *= Ex[m].a[0]*FIELDFAC;
350 for(i=start; (i<start+homenr); i++)
351 f[i][m] += charge[i]*Ext[m];
360 fprintf(fp,"%10g %10g %10g %10g #FIELD\n",t,
361 Ext[XX]/FIELDFAC,Ext[YY]/FIELDFAC,Ext[ZZ]/FIELDFAC);
365 static void calc_virial(FILE *fplog,int start,int homenr,rvec x[],rvec f[],
366 tensor vir_part,t_graph *graph,matrix box,
367 t_nrnb *nrnb,const t_forcerec *fr,int ePBC)
372 /* The short-range virial from surrounding boxes */
374 calc_vir(fplog,SHIFTS,fr->shift_vec,fr->fshift,vir_part,ePBC==epbcSCREW,box);
375 inc_nrnb(nrnb,eNR_VIRIAL,SHIFTS);
377 /* Calculate partial virial, for local atoms only, based on short range.
378 * Total virial is computed in global_stat, called from do_md
380 f_calc_vir(fplog,start,start+homenr,x,f,vir_part,graph,box);
381 inc_nrnb(nrnb,eNR_VIRIAL,homenr);
383 /* Add position restraint contribution */
384 for(i=0; i<DIM; i++) {
385 vir_part[i][i] += fr->vir_diag_posres[i];
388 /* Add wall contribution */
389 for(i=0; i<DIM; i++) {
390 vir_part[i][ZZ] += fr->vir_wall_z[i];
394 pr_rvecs(debug,0,"vir_part",vir_part,DIM);
397 static void posres_wrapper(FILE *fplog,
405 gmx_enerdata_t *enerd,
413 /* Position restraints always require full pbc */
414 set_pbc(&pbc,ir->ePBC,box);
416 v = posres(top->idef.il[F_POSRES].nr,top->idef.il[F_POSRES].iatoms,
417 top->idef.iparams_posres,
418 (const rvec*)x,fr->f_novirsum,fr->vir_diag_posres,
419 ir->ePBC==epbcNONE ? NULL : &pbc,
420 lambda[efptRESTRAINT],&dvdl,
421 fr->rc_scaling,fr->ePBC,fr->posres_com,fr->posres_comB);
424 fprintf(fplog,sepdvdlformat,
425 interaction_function[F_POSRES].longname,v,dvdl);
427 enerd->term[F_POSRES] += v;
428 /* If just the force constant changes, the FEP term is linear,
429 * but if k changes, it is not.
431 enerd->dvdl_nonlin[efptRESTRAINT] += dvdl;
432 inc_nrnb(nrnb,eNR_POSRES,top->idef.il[F_POSRES].nr/2);
434 if ((ir->fepvals->n_lambda > 0) && (flags & GMX_FORCE_DHDL))
436 for(i=0; i<enerd->n_lambda; i++)
438 real dvdl_dum,lambda_dum;
440 lambda_dum = (i==0 ? lambda[efptRESTRAINT] : ir->fepvals->all_lambda[efptRESTRAINT][i-1]);
441 v = posres(top->idef.il[F_POSRES].nr,top->idef.il[F_POSRES].iatoms,
442 top->idef.iparams_posres,
443 (const rvec*)x,NULL,NULL,
444 ir->ePBC==epbcNONE ? NULL : &pbc,lambda_dum,&dvdl,
445 fr->rc_scaling,fr->ePBC,fr->posres_com,fr->posres_comB);
446 enerd->enerpart_lambda[i] += v;
451 static void pull_potential_wrapper(FILE *fplog,
459 gmx_enerdata_t *enerd,
466 /* Calculate the center of mass forces, this requires communication,
467 * which is why pull_potential is called close to other communication.
468 * The virial contribution is calculated directly,
469 * which is why we call pull_potential after calc_virial.
471 set_pbc(&pbc,ir->ePBC,box);
473 enerd->term[F_COM_PULL] +=
474 pull_potential(ir->ePull,ir->pull,mdatoms,&pbc,
475 cr,t,lambda[efptRESTRAINT],x,f,vir_force,&dvdl);
478 fprintf(fplog,sepdvdlformat,"Com pull",enerd->term[F_COM_PULL],dvdl);
480 enerd->dvdl_lin[efptRESTRAINT] += dvdl;
483 static void pme_receive_force_ener(FILE *fplog,
486 gmx_wallcycle_t wcycle,
487 gmx_enerdata_t *enerd,
491 float cycles_ppdpme,cycles_seppme;
493 cycles_ppdpme = wallcycle_stop(wcycle,ewcPPDURINGPME);
494 dd_cycles_add(cr->dd,cycles_ppdpme,ddCyclPPduringPME);
496 /* In case of node-splitting, the PP nodes receive the long-range
497 * forces, virial and energy from the PME nodes here.
499 wallcycle_start(wcycle,ewcPP_PMEWAITRECVF);
501 gmx_pme_receive_f(cr,fr->f_novirsum,fr->vir_el_recip,&e,&dvdl,
505 fprintf(fplog,sepdvdlformat,"PME mesh",e,dvdl);
507 enerd->term[F_COUL_RECIP] += e;
508 enerd->dvdl_lin[efptCOUL] += dvdl;
511 dd_cycles_add(cr->dd,cycles_seppme,ddCyclPME);
513 wallcycle_stop(wcycle,ewcPP_PMEWAITRECVF);
516 static void print_large_forces(FILE *fp,t_mdatoms *md,t_commrec *cr,
517 gmx_large_int_t step,real pforce,rvec *x,rvec *f)
521 char buf[STEPSTRSIZE];
524 for(i=md->start; i<md->start+md->homenr; i++) {
526 /* We also catch NAN, if the compiler does not optimize this away. */
527 if (fn2 >= pf2 || fn2 != fn2) {
528 fprintf(fp,"step %s atom %6d x %8.3f %8.3f %8.3f force %12.5e\n",
529 gmx_step_str(step,buf),
530 ddglatnr(cr->dd,i),x[i][XX],x[i][YY],x[i][ZZ],sqrt(fn2));
535 static void post_process_forces(FILE *fplog,
537 gmx_large_int_t step,
538 t_nrnb *nrnb,gmx_wallcycle_t wcycle,
545 t_forcerec *fr,gmx_vsite_t *vsite,
552 /* Spread the mesh force on virtual sites to the other particles...
553 * This is parallellized. MPI communication is performed
554 * if the constructing atoms aren't local.
556 wallcycle_start(wcycle,ewcVSITESPREAD);
557 spread_vsite_f(fplog,vsite,x,fr->f_novirsum,NULL,
558 (flags & GMX_FORCE_VIRIAL),fr->vir_el_recip,
560 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
561 wallcycle_stop(wcycle,ewcVSITESPREAD);
563 if (flags & GMX_FORCE_VIRIAL)
565 /* Now add the forces, this is local */
568 sum_forces(0,fr->f_novirsum_n,f,fr->f_novirsum);
572 sum_forces(mdatoms->start,mdatoms->start+mdatoms->homenr,
575 if (EEL_FULL(fr->eeltype))
577 /* Add the mesh contribution to the virial */
578 m_add(vir_force,fr->vir_el_recip,vir_force);
582 pr_rvecs(debug,0,"vir_force",vir_force,DIM);
587 if (fr->print_force >= 0)
589 print_large_forces(stderr,mdatoms,cr,step,fr->print_force,x,f);
593 static void do_nb_verlet(t_forcerec *fr,
594 interaction_const_t *ic,
595 gmx_enerdata_t *enerd,
596 int flags, int ilocality,
599 gmx_wallcycle_t wcycle)
601 int nnbl, kernel_type, sh_e;
603 nonbonded_verlet_group_t *nbvg;
605 if (!(flags & GMX_FORCE_NONBONDED))
607 /* skip non-bonded calculation */
611 nbvg = &fr->nbv->grp[ilocality];
613 /* CUDA kernel launch overhead is already timed separately */
614 if (fr->cutoff_scheme != ecutsVERLET)
616 gmx_incons("Invalid cut-off scheme passed!");
619 if (nbvg->kernel_type != nbk8x8x8_CUDA)
621 wallcycle_sub_start(wcycle, ewcsNONBONDED);
623 switch (nbvg->kernel_type)
626 nbnxn_kernel_ref(&nbvg->nbl_lists,
632 enerd->grpp.ener[egCOULSR],
634 enerd->grpp.ener[egBHAMSR] :
635 enerd->grpp.ener[egLJSR]);
638 case nbk4xN_X86_SIMD128:
639 nbnxn_kernel_x86_simd128(&nbvg->nbl_lists,
645 enerd->grpp.ener[egCOULSR],
647 enerd->grpp.ener[egBHAMSR] :
648 enerd->grpp.ener[egLJSR]);
650 case nbk4xN_X86_SIMD256:
651 nbnxn_kernel_x86_simd256(&nbvg->nbl_lists,
657 enerd->grpp.ener[egCOULSR],
659 enerd->grpp.ener[egBHAMSR] :
660 enerd->grpp.ener[egLJSR]);
664 nbnxn_cuda_launch_kernel(fr->nbv->cu_nbv, nbvg->nbat, flags, ilocality);
667 case nbk8x8x8_PlainC:
668 nbnxn_kernel_gpu_ref(nbvg->nbl_lists.nbl[0],
673 nbvg->nbat->out[0].f,
675 enerd->grpp.ener[egCOULSR],
677 enerd->grpp.ener[egBHAMSR] :
678 enerd->grpp.ener[egLJSR]);
682 gmx_incons("Invalid nonbonded kernel type passed!");
685 if (nbvg->kernel_type != nbk8x8x8_CUDA)
687 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
690 /* In eNR_??? the nbnxn F+E kernels are always the F kernel + 1 */
691 sh_e = ((flags & GMX_FORCE_ENERGY) ? 1 : 0);
693 ((EEL_RF(ic->eeltype) || ic->eeltype == eelCUT) ?
694 eNR_NBNXN_LJ_RF : eNR_NBNXN_LJ_TAB) + sh_e,
695 nbvg->nbl_lists.natpair_ljq);
696 inc_nrnb(nrnb,eNR_NBNXN_LJ+sh_e,nbvg->nbl_lists.natpair_lj);
698 ((EEL_RF(ic->eeltype) || ic->eeltype == eelCUT) ?
699 eNR_NBNXN_RF : eNR_NBNXN_TAB)+sh_e,
700 nbvg->nbl_lists.natpair_q);
703 void do_force_cutsVERLET(FILE *fplog,t_commrec *cr,
704 t_inputrec *inputrec,
705 gmx_large_int_t step,t_nrnb *nrnb,gmx_wallcycle_t wcycle,
708 gmx_groups_t *groups,
709 matrix box,rvec x[],history_t *hist,
713 gmx_enerdata_t *enerd,t_fcdata *fcd,
714 real *lambda,t_graph *graph,
715 t_forcerec *fr, interaction_const_t *ic,
716 gmx_vsite_t *vsite,rvec mu_tot,
717 double t,FILE *field,gmx_edsam_t ed,
725 gmx_bool bSepDVDL,bStateChanged,bNS,bFillGrid,bCalcCGCM,bBS;
726 gmx_bool bDoLongRange,bDoForces,bSepLRF,bUseGPU,bUseOrEmulGPU;
727 gmx_bool bDiffKernels=FALSE;
731 float cycles_pme,cycles_force;
732 nonbonded_verlet_t *nbv;
736 nb_kernel_type = fr->nbv->grp[0].kernel_type;
738 start = mdatoms->start;
739 homenr = mdatoms->homenr;
741 bSepDVDL = (fr->bSepDVDL && do_per_step(step,inputrec->nstlog));
743 clear_mat(vir_force);
746 if (DOMAINDECOMP(cr))
748 cg1 = cr->dd->ncg_tot;
759 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
760 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll==FALSE);
761 bFillGrid = (bNS && bStateChanged);
762 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
763 bDoLongRange = (fr->bTwinRange && bNS && (flags & GMX_FORCE_DO_LR));
764 bDoForces = (flags & GMX_FORCE_FORCES);
765 bSepLRF = (bDoLongRange && bDoForces && (flags & GMX_FORCE_SEPLRF));
766 bUseGPU = fr->nbv->bUseGPU;
767 bUseOrEmulGPU = bUseGPU || (nbv->grp[0].kernel_type == nbk8x8x8_PlainC);
771 update_forcerec(fplog,fr,box);
773 if (NEED_MUTOT(*inputrec))
775 /* Calculate total (local) dipole moment in a temporary common array.
776 * This makes it possible to sum them over nodes faster.
778 calc_mu(start,homenr,
779 x,mdatoms->chargeA,mdatoms->chargeB,mdatoms->nChargePerturbed,
784 if (fr->ePBC != epbcNONE) {
785 /* Compute shift vectors every step,
786 * because of pressure coupling or box deformation!
788 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
789 calc_shifts(box,fr->shift_vec);
792 put_atoms_in_box_omp(fr->ePBC,box,homenr,x);
793 inc_nrnb(nrnb,eNR_SHIFTX,homenr);
795 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph) {
796 unshift_self(graph,box,x);
800 nbnxn_atomdata_copy_shiftvec(flags & GMX_FORCE_DYNAMICBOX,
801 fr->shift_vec,nbv->grp[0].nbat);
804 if (!(cr->duty & DUTY_PME)) {
805 /* Send particle coordinates to the pme nodes.
806 * Since this is only implemented for domain decomposition
807 * and domain decomposition does not use the graph,
808 * we do not need to worry about shifting.
811 wallcycle_start(wcycle,ewcPP_PMESENDX);
812 GMX_MPE_LOG(ev_send_coordinates_start);
814 bBS = (inputrec->nwall == 2);
817 svmul(inputrec->wall_ewald_zfac,boxs[ZZ],boxs[ZZ]);
820 gmx_pme_send_x(cr,bBS ? boxs : box,x,
821 mdatoms->nChargePerturbed,lambda[efptCOUL],
822 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),step);
824 GMX_MPE_LOG(ev_send_coordinates_finish);
825 wallcycle_stop(wcycle,ewcPP_PMESENDX);
829 /* do gridding for pair search */
832 if (graph && bStateChanged)
834 /* Calculate intramolecular shift vectors to make molecules whole */
835 mk_mshift(fplog,graph,fr->ePBC,box,x);
839 box_diag[XX] = box[XX][XX];
840 box_diag[YY] = box[YY][YY];
841 box_diag[ZZ] = box[ZZ][ZZ];
843 wallcycle_start(wcycle,ewcNS);
846 wallcycle_sub_start(wcycle,ewcsNBS_GRID_LOCAL);
847 nbnxn_put_on_grid(nbv->nbs,fr->ePBC,box,
849 0,mdatoms->homenr,-1,fr->cginfo,x,
851 nbv->grp[eintLocal].kernel_type,
852 nbv->grp[eintLocal].nbat);
853 wallcycle_sub_stop(wcycle,ewcsNBS_GRID_LOCAL);
857 wallcycle_sub_start(wcycle,ewcsNBS_GRID_NONLOCAL);
858 nbnxn_put_on_grid_nonlocal(nbv->nbs,domdec_zones(cr->dd),
860 nbv->grp[eintNonlocal].kernel_type,
861 nbv->grp[eintNonlocal].nbat);
862 wallcycle_sub_stop(wcycle,ewcsNBS_GRID_NONLOCAL);
865 if (nbv->ngrp == 1 ||
866 nbv->grp[eintNonlocal].nbat == nbv->grp[eintLocal].nbat)
868 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat,eatAll,
869 nbv->nbs,mdatoms,fr->cginfo);
873 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat,eatLocal,
874 nbv->nbs,mdatoms,fr->cginfo);
875 nbnxn_atomdata_set(nbv->grp[eintNonlocal].nbat,eatAll,
876 nbv->nbs,mdatoms,fr->cginfo);
878 wallcycle_stop(wcycle, ewcNS);
881 /* initialize the GPU atom data and copy shift vector */
886 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
887 nbnxn_cuda_init_atomdata(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
888 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
891 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
892 nbnxn_cuda_upload_shiftvec(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
893 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
896 /* do local pair search */
899 wallcycle_start_nocount(wcycle,ewcNS);
900 wallcycle_sub_start(wcycle,ewcsNBS_SEARCH_LOCAL);
901 nbnxn_make_pairlist(nbv->nbs,nbv->grp[eintLocal].nbat,
904 nbv->min_ci_balanced,
905 &nbv->grp[eintLocal].nbl_lists,
907 nbv->grp[eintLocal].kernel_type,
909 wallcycle_sub_stop(wcycle,ewcsNBS_SEARCH_LOCAL);
913 /* initialize local pair-list on the GPU */
914 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
915 nbv->grp[eintLocal].nbl_lists.nbl[0],
918 wallcycle_stop(wcycle, ewcNS);
922 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
923 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
924 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs,eatLocal,FALSE,x,
925 nbv->grp[eintLocal].nbat);
926 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
927 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
932 wallcycle_start(wcycle,ewcLAUNCH_GPU_NB);
933 /* launch local nonbonded F on GPU */
934 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFNo,
936 wallcycle_stop(wcycle,ewcLAUNCH_GPU_NB);
939 /* Communicate coordinates and sum dipole if necessary +
940 do non-local pair search */
941 if (DOMAINDECOMP(cr))
943 bDiffKernels = (nbv->grp[eintNonlocal].kernel_type !=
944 nbv->grp[eintLocal].kernel_type);
948 /* With GPU+CPU non-bonded calculations we need to copy
949 * the local coordinates to the non-local nbat struct
950 * (in CPU format) as the non-local kernel call also
951 * calculates the local - non-local interactions.
953 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
954 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
955 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs,eatLocal,TRUE,x,
956 nbv->grp[eintNonlocal].nbat);
957 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
958 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
963 wallcycle_start_nocount(wcycle,ewcNS);
964 wallcycle_sub_start(wcycle,ewcsNBS_SEARCH_NONLOCAL);
968 nbnxn_grid_add_simple(nbv->nbs,nbv->grp[eintNonlocal].nbat);
971 nbnxn_make_pairlist(nbv->nbs,nbv->grp[eintNonlocal].nbat,
974 nbv->min_ci_balanced,
975 &nbv->grp[eintNonlocal].nbl_lists,
977 nbv->grp[eintNonlocal].kernel_type,
980 wallcycle_sub_stop(wcycle,ewcsNBS_SEARCH_NONLOCAL);
982 if (nbv->grp[eintNonlocal].kernel_type == nbk8x8x8_CUDA)
984 /* initialize non-local pair-list on the GPU */
985 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
986 nbv->grp[eintNonlocal].nbl_lists.nbl[0],
989 wallcycle_stop(wcycle,ewcNS);
993 wallcycle_start(wcycle,ewcMOVEX);
994 dd_move_x(cr->dd,box,x);
996 /* When we don't need the total dipole we sum it in global_stat */
997 if (bStateChanged && NEED_MUTOT(*inputrec))
999 gmx_sumd(2*DIM,mu,cr);
1001 wallcycle_stop(wcycle,ewcMOVEX);
1003 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1004 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
1005 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs,eatNonlocal,FALSE,x,
1006 nbv->grp[eintNonlocal].nbat);
1007 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
1008 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1011 if (bUseGPU && !bDiffKernels)
1013 wallcycle_start(wcycle,ewcLAUNCH_GPU_NB);
1014 /* launch non-local nonbonded F on GPU */
1015 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFNo,
1017 cycles_force += wallcycle_stop(wcycle,ewcLAUNCH_GPU_NB);
1023 /* launch D2H copy-back F */
1024 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1025 if (DOMAINDECOMP(cr) && !bDiffKernels)
1027 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintNonlocal].nbat,
1028 flags, eatNonlocal);
1030 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintLocal].nbat,
1032 cycles_force += wallcycle_stop(wcycle,ewcLAUNCH_GPU_NB);
1035 if (bStateChanged && NEED_MUTOT(*inputrec))
1039 gmx_sumd(2*DIM,mu,cr);
1046 fr->mu_tot[i][j] = mu[i*DIM + j];
1050 if (fr->efep == efepNO)
1052 copy_rvec(fr->mu_tot[0],mu_tot);
1056 for(j=0; j<DIM; j++)
1059 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] +
1060 lambda[efptCOUL]*fr->mu_tot[1][j];
1064 /* Reset energies */
1065 reset_enerdata(&(inputrec->opts),fr,bNS,enerd,MASTER(cr));
1066 clear_rvecs(SHIFTS,fr->fshift);
1068 if (DOMAINDECOMP(cr))
1070 if (!(cr->duty & DUTY_PME))
1072 wallcycle_start(wcycle,ewcPPDURINGPME);
1073 dd_force_flop_start(cr->dd,nrnb);
1077 /* Start the force cycle counter.
1078 * This counter is stopped in do_forcelow_level.
1079 * No parallel communication should occur while this counter is running,
1080 * since that will interfere with the dynamic load balancing.
1082 wallcycle_start(wcycle,ewcFORCE);
1085 /* Reset forces for which the virial is calculated separately:
1086 * PME/Ewald forces if necessary */
1087 if (fr->bF_NoVirSum)
1089 if (flags & GMX_FORCE_VIRIAL)
1091 fr->f_novirsum = fr->f_novirsum_alloc;
1092 GMX_BARRIER(cr->mpi_comm_mygroup);
1095 clear_rvecs(fr->f_novirsum_n,fr->f_novirsum);
1099 clear_rvecs(homenr,fr->f_novirsum+start);
1101 GMX_BARRIER(cr->mpi_comm_mygroup);
1105 /* We are not calculating the pressure so we do not need
1106 * a separate array for forces that do not contribute
1113 /* Clear the short- and long-range forces */
1114 clear_rvecs(fr->natoms_force_constr,f);
1115 if(bSepLRF && do_per_step(step,inputrec->nstcalclr))
1117 clear_rvecs(fr->natoms_force_constr,fr->f_twin);
1120 clear_rvec(fr->vir_diag_posres);
1122 GMX_BARRIER(cr->mpi_comm_mygroup);
1124 if (inputrec->ePull == epullCONSTRAINT)
1126 clear_pull_forces(inputrec->pull);
1129 /* update QMMMrec, if necessary */
1132 update_QMMMrec(cr,fr,x,mdatoms,box,top);
1135 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1137 posres_wrapper(fplog,flags,bSepDVDL,inputrec,nrnb,top,box,x,
1141 /* Compute the bonded and non-bonded energies and optionally forces */
1142 do_force_lowlevel(fplog,step,fr,inputrec,&(top->idef),
1143 cr,nrnb,wcycle,mdatoms,&(inputrec->opts),
1144 x,hist,f, bSepLRF ? fr->f_twin : f,enerd,fcd,mtop,top,fr->born,
1145 &(top->atomtypes),bBornRadii,box,
1146 inputrec->fepvals,lambda,graph,&(top->excls),fr->mu_tot,
1147 flags, &cycles_pme);
1151 if (do_per_step(step,inputrec->nstcalclr))
1153 /* Add the long range forces to the short range forces */
1154 for(i=0; i<fr->natoms_force_constr; i++)
1156 rvec_add(fr->f_twin[i],f[i],f[i]);
1163 /* Maybe we should move this into do_force_lowlevel */
1164 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFYes,
1169 if (!bUseOrEmulGPU || bDiffKernels)
1173 if (DOMAINDECOMP(cr))
1175 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal,
1176 bDiffKernels ? enbvClearFYes : enbvClearFNo,
1186 aloc = eintNonlocal;
1189 /* Add all the non-bonded force to the normal force array.
1190 * This can be split into a local a non-local part when overlapping
1191 * communication with calculation with domain decomposition.
1193 cycles_force += wallcycle_stop(wcycle,ewcFORCE);
1194 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1195 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1196 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs,eatAll,nbv->grp[aloc].nbat,f);
1197 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1198 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1199 wallcycle_start_nocount(wcycle,ewcFORCE);
1201 /* if there are multiple fshift output buffers reduce them */
1202 if ((flags & GMX_FORCE_VIRIAL) &&
1203 nbv->grp[aloc].nbl_lists.nnbl > 1)
1205 nbnxn_atomdata_add_nbat_fshift_to_fshift(nbv->grp[aloc].nbat,
1210 cycles_force += wallcycle_stop(wcycle,ewcFORCE);
1211 GMX_BARRIER(cr->mpi_comm_mygroup);
1215 do_flood(fplog,cr,x,f,ed,box,step,bNS);
1218 if (bUseOrEmulGPU && !bDiffKernels)
1220 /* wait for non-local forces (or calculate in emulation mode) */
1221 if (DOMAINDECOMP(cr))
1225 wallcycle_start(wcycle,ewcWAIT_GPU_NB_NL);
1226 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1227 nbv->grp[eintNonlocal].nbat,
1229 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1231 cycles_force += wallcycle_stop(wcycle,ewcWAIT_GPU_NB_NL);
1235 wallcycle_start_nocount(wcycle,ewcFORCE);
1236 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFYes,
1238 cycles_force += wallcycle_stop(wcycle,ewcFORCE);
1240 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1241 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1242 /* skip the reduction if there was no non-local work to do */
1243 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1245 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs,eatNonlocal,
1246 nbv->grp[eintNonlocal].nbat,f);
1248 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1249 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1255 /* Communicate the forces */
1258 wallcycle_start(wcycle,ewcMOVEF);
1259 if (DOMAINDECOMP(cr))
1261 dd_move_f(cr->dd,f,fr->fshift);
1262 /* Do we need to communicate the separate force array
1263 * for terms that do not contribute to the single sum virial?
1264 * Position restraints and electric fields do not introduce
1265 * inter-cg forces, only full electrostatics methods do.
1266 * When we do not calculate the virial, fr->f_novirsum = f,
1267 * so we have already communicated these forces.
1269 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1270 (flags & GMX_FORCE_VIRIAL))
1272 dd_move_f(cr->dd,fr->f_novirsum,NULL);
1276 /* We should not update the shift forces here,
1277 * since f_twin is already included in f.
1279 dd_move_f(cr->dd,fr->f_twin,NULL);
1282 wallcycle_stop(wcycle,ewcMOVEF);
1288 /* wait for local forces (or calculate in emulation mode) */
1291 wallcycle_start(wcycle,ewcWAIT_GPU_NB_L);
1292 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1293 nbv->grp[eintLocal].nbat,
1295 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1297 wallcycle_stop(wcycle,ewcWAIT_GPU_NB_L);
1299 /* now clear the GPU outputs while we finish the step on the CPU */
1300 nbnxn_cuda_clear_outputs(nbv->cu_nbv, flags);
1304 wallcycle_start_nocount(wcycle,ewcFORCE);
1305 do_nb_verlet(fr, ic, enerd, flags, eintLocal,
1306 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
1308 wallcycle_stop(wcycle,ewcFORCE);
1310 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1311 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1312 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1314 /* skip the reduction if there was no non-local work to do */
1315 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs,eatLocal,
1316 nbv->grp[eintLocal].nbat,f);
1318 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1319 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1322 if (DOMAINDECOMP(cr))
1324 dd_force_flop_stop(cr->dd,nrnb);
1327 dd_cycles_add(cr->dd,cycles_force-cycles_pme,ddCyclF);
1333 if (IR_ELEC_FIELD(*inputrec))
1335 /* Compute forces due to electric field */
1336 calc_f_el(MASTER(cr) ? field : NULL,
1337 start,homenr,mdatoms->chargeA,x,fr->f_novirsum,
1338 inputrec->ex,inputrec->et,t);
1341 /* If we have NoVirSum forces, but we do not calculate the virial,
1342 * we sum fr->f_novirum=f later.
1344 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1346 wallcycle_start(wcycle,ewcVSITESPREAD);
1347 spread_vsite_f(fplog,vsite,x,f,fr->fshift,FALSE,NULL,nrnb,
1348 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
1349 wallcycle_stop(wcycle,ewcVSITESPREAD);
1353 wallcycle_start(wcycle,ewcVSITESPREAD);
1354 spread_vsite_f(fplog,vsite,x,fr->f_twin,NULL,FALSE,NULL,
1356 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
1357 wallcycle_stop(wcycle,ewcVSITESPREAD);
1361 if (flags & GMX_FORCE_VIRIAL)
1363 /* Calculation of the virial must be done after vsites! */
1364 calc_virial(fplog,mdatoms->start,mdatoms->homenr,x,f,
1365 vir_force,graph,box,nrnb,fr,inputrec->ePBC);
1369 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1371 pull_potential_wrapper(fplog,bSepDVDL,cr,inputrec,box,x,
1372 f,vir_force,mdatoms,enerd,lambda,t);
1375 if (PAR(cr) && !(cr->duty & DUTY_PME))
1377 /* In case of node-splitting, the PP nodes receive the long-range
1378 * forces, virial and energy from the PME nodes here.
1380 pme_receive_force_ener(fplog,bSepDVDL,cr,wcycle,enerd,fr);
1385 post_process_forces(fplog,cr,step,nrnb,wcycle,
1386 top,box,x,f,vir_force,mdatoms,graph,fr,vsite,
1390 /* Sum the potential energy terms from group contributions */
1391 sum_epot(&(inputrec->opts),enerd);
1394 void do_force_cutsGROUP(FILE *fplog,t_commrec *cr,
1395 t_inputrec *inputrec,
1396 gmx_large_int_t step,t_nrnb *nrnb,gmx_wallcycle_t wcycle,
1397 gmx_localtop_t *top,
1399 gmx_groups_t *groups,
1400 matrix box,rvec x[],history_t *hist,
1404 gmx_enerdata_t *enerd,t_fcdata *fcd,
1405 real *lambda,t_graph *graph,
1406 t_forcerec *fr,gmx_vsite_t *vsite,rvec mu_tot,
1407 double t,FILE *field,gmx_edsam_t ed,
1408 gmx_bool bBornRadii,
1414 gmx_bool bSepDVDL,bStateChanged,bNS,bFillGrid,bCalcCGCM,bBS;
1415 gmx_bool bDoLongRangeNS,bDoForces,bDoPotential,bSepLRF;
1416 gmx_bool bDoAdressWF;
1418 rvec vzero,box_diag;
1419 real e,v,dvdlambda[efptNR];
1421 float cycles_pme,cycles_force;
1423 start = mdatoms->start;
1424 homenr = mdatoms->homenr;
1426 bSepDVDL = (fr->bSepDVDL && do_per_step(step,inputrec->nstlog));
1428 clear_mat(vir_force);
1432 pd_cg_range(cr,&cg0,&cg1);
1437 if (DOMAINDECOMP(cr))
1439 cg1 = cr->dd->ncg_tot;
1451 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
1452 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll==FALSE);
1453 /* Should we update the long-range neighborlists at this step? */
1454 bDoLongRangeNS = fr->bTwinRange && bNS;
1455 /* Should we perform the long-range nonbonded evaluation inside the neighborsearching? */
1456 bFillGrid = (bNS && bStateChanged);
1457 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
1458 bDoForces = (flags & GMX_FORCE_FORCES);
1459 bDoPotential = (flags & GMX_FORCE_ENERGY);
1460 bSepLRF = ((inputrec->nstcalclr>1) && bDoForces &&
1461 (flags & GMX_FORCE_SEPLRF) && (flags & GMX_FORCE_DO_LR));
1463 /* should probably move this to the forcerec since it doesn't change */
1464 bDoAdressWF = ((fr->adress_type!=eAdressOff));
1468 update_forcerec(fplog,fr,box);
1470 if (NEED_MUTOT(*inputrec))
1472 /* Calculate total (local) dipole moment in a temporary common array.
1473 * This makes it possible to sum them over nodes faster.
1475 calc_mu(start,homenr,
1476 x,mdatoms->chargeA,mdatoms->chargeB,mdatoms->nChargePerturbed,
1481 if (fr->ePBC != epbcNONE) {
1482 /* Compute shift vectors every step,
1483 * because of pressure coupling or box deformation!
1485 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
1486 calc_shifts(box,fr->shift_vec);
1489 put_charge_groups_in_box(fplog,cg0,cg1,fr->ePBC,box,
1490 &(top->cgs),x,fr->cg_cm);
1491 inc_nrnb(nrnb,eNR_CGCM,homenr);
1492 inc_nrnb(nrnb,eNR_RESETX,cg1-cg0);
1494 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph) {
1495 unshift_self(graph,box,x);
1498 else if (bCalcCGCM) {
1499 calc_cgcm(fplog,cg0,cg1,&(top->cgs),x,fr->cg_cm);
1500 inc_nrnb(nrnb,eNR_CGCM,homenr);
1505 move_cgcm(fplog,cr,fr->cg_cm);
1508 pr_rvecs(debug,0,"cgcm",fr->cg_cm,top->cgs.nr);
1512 if (!(cr->duty & DUTY_PME)) {
1513 /* Send particle coordinates to the pme nodes.
1514 * Since this is only implemented for domain decomposition
1515 * and domain decomposition does not use the graph,
1516 * we do not need to worry about shifting.
1519 wallcycle_start(wcycle,ewcPP_PMESENDX);
1520 GMX_MPE_LOG(ev_send_coordinates_start);
1522 bBS = (inputrec->nwall == 2);
1525 svmul(inputrec->wall_ewald_zfac,boxs[ZZ],boxs[ZZ]);
1528 gmx_pme_send_x(cr,bBS ? boxs : box,x,
1529 mdatoms->nChargePerturbed,lambda[efptCOUL],
1530 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),step);
1532 GMX_MPE_LOG(ev_send_coordinates_finish);
1533 wallcycle_stop(wcycle,ewcPP_PMESENDX);
1535 #endif /* GMX_MPI */
1537 /* Communicate coordinates and sum dipole if necessary */
1540 wallcycle_start(wcycle,ewcMOVEX);
1541 if (DOMAINDECOMP(cr))
1543 dd_move_x(cr->dd,box,x);
1547 move_x(fplog,cr,GMX_LEFT,GMX_RIGHT,x,nrnb);
1549 wallcycle_stop(wcycle,ewcMOVEX);
1552 /* update adress weight beforehand */
1553 if(bStateChanged && bDoAdressWF)
1555 /* need pbc for adress weight calculation with pbc_dx */
1556 set_pbc(&pbc,inputrec->ePBC,box);
1557 if(fr->adress_site == eAdressSITEcog)
1559 update_adress_weights_cog(top->idef.iparams,top->idef.il,x,fr,mdatoms,
1560 inputrec->ePBC==epbcNONE ? NULL : &pbc);
1562 else if (fr->adress_site == eAdressSITEcom)
1564 update_adress_weights_com(fplog,cg0,cg1,&(top->cgs),x,fr,mdatoms,
1565 inputrec->ePBC==epbcNONE ? NULL : &pbc);
1567 else if (fr->adress_site == eAdressSITEatomatom){
1568 update_adress_weights_atom_per_atom(cg0,cg1,&(top->cgs),x,fr,mdatoms,
1569 inputrec->ePBC==epbcNONE ? NULL : &pbc);
1573 update_adress_weights_atom(cg0,cg1,&(top->cgs),x,fr,mdatoms,
1574 inputrec->ePBC==epbcNONE ? NULL : &pbc);
1578 if (NEED_MUTOT(*inputrec))
1585 gmx_sumd(2*DIM,mu,cr);
1591 fr->mu_tot[i][j] = mu[i*DIM + j];
1595 if (fr->efep == efepNO)
1597 copy_rvec(fr->mu_tot[0],mu_tot);
1601 for(j=0; j<DIM; j++)
1604 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] + lambda[efptCOUL]*fr->mu_tot[1][j];
1609 /* Reset energies */
1610 reset_enerdata(&(inputrec->opts),fr,bNS,enerd,MASTER(cr));
1611 clear_rvecs(SHIFTS,fr->fshift);
1615 wallcycle_start(wcycle,ewcNS);
1617 if (graph && bStateChanged)
1619 /* Calculate intramolecular shift vectors to make molecules whole */
1620 mk_mshift(fplog,graph,fr->ePBC,box,x);
1623 /* Do the actual neighbour searching and if twin range electrostatics
1624 * also do the calculation of long range forces and energies.
1626 for (i=0;i<efptNR;i++) {dvdlambda[i] = 0;}
1628 groups,&(inputrec->opts),top,mdatoms,
1629 cr,nrnb,lambda,dvdlambda,&enerd->grpp,bFillGrid,
1633 fprintf(fplog,sepdvdlformat,"LR non-bonded",0.0,dvdlambda);
1635 enerd->dvdl_lin[efptVDW] += dvdlambda[efptVDW];
1636 enerd->dvdl_lin[efptCOUL] += dvdlambda[efptCOUL];
1638 wallcycle_stop(wcycle,ewcNS);
1641 if (inputrec->implicit_solvent && bNS)
1643 make_gb_nblist(cr,inputrec->gb_algorithm,inputrec->rlist,
1644 x,box,fr,&top->idef,graph,fr->born);
1647 if (DOMAINDECOMP(cr))
1649 if (!(cr->duty & DUTY_PME))
1651 wallcycle_start(wcycle,ewcPPDURINGPME);
1652 dd_force_flop_start(cr->dd,nrnb);
1658 /* Enforced rotation has its own cycle counter that starts after the collective
1659 * coordinates have been communicated. It is added to ddCyclF to allow
1660 * for proper load-balancing */
1661 wallcycle_start(wcycle,ewcROT);
1662 do_rotation(cr,inputrec,box,x,t,step,wcycle,bNS);
1663 wallcycle_stop(wcycle,ewcROT);
1666 /* Start the force cycle counter.
1667 * This counter is stopped in do_forcelow_level.
1668 * No parallel communication should occur while this counter is running,
1669 * since that will interfere with the dynamic load balancing.
1671 wallcycle_start(wcycle,ewcFORCE);
1675 /* Reset forces for which the virial is calculated separately:
1676 * PME/Ewald forces if necessary */
1677 if (fr->bF_NoVirSum)
1679 if (flags & GMX_FORCE_VIRIAL)
1681 fr->f_novirsum = fr->f_novirsum_alloc;
1682 GMX_BARRIER(cr->mpi_comm_mygroup);
1685 clear_rvecs(fr->f_novirsum_n,fr->f_novirsum);
1689 clear_rvecs(homenr,fr->f_novirsum+start);
1691 GMX_BARRIER(cr->mpi_comm_mygroup);
1695 /* We are not calculating the pressure so we do not need
1696 * a separate array for forces that do not contribute
1703 /* Clear the short- and long-range forces */
1704 clear_rvecs(fr->natoms_force_constr,f);
1705 if(bSepLRF && do_per_step(step,inputrec->nstcalclr))
1707 clear_rvecs(fr->natoms_force_constr,fr->f_twin);
1710 clear_rvec(fr->vir_diag_posres);
1712 GMX_BARRIER(cr->mpi_comm_mygroup);
1714 if (inputrec->ePull == epullCONSTRAINT)
1716 clear_pull_forces(inputrec->pull);
1719 /* update QMMMrec, if necessary */
1722 update_QMMMrec(cr,fr,x,mdatoms,box,top);
1725 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1727 posres_wrapper(fplog,flags,bSepDVDL,inputrec,nrnb,top,box,x,
1731 /* Compute the bonded and non-bonded energies and optionally forces */
1732 do_force_lowlevel(fplog,step,fr,inputrec,&(top->idef),
1733 cr,nrnb,wcycle,mdatoms,&(inputrec->opts),
1734 x,hist,f, bSepLRF ? fr->f_twin : f,enerd,fcd,mtop,top,fr->born,
1735 &(top->atomtypes),bBornRadii,box,
1736 inputrec->fepvals,lambda,
1737 graph,&(top->excls),fr->mu_tot,
1743 if (do_per_step(step,inputrec->nstcalclr))
1745 /* Add the long range forces to the short range forces */
1746 for(i=0; i<fr->natoms_force_constr; i++)
1748 rvec_add(fr->f_twin[i],f[i],f[i]);
1753 cycles_force = wallcycle_stop(wcycle,ewcFORCE);
1754 GMX_BARRIER(cr->mpi_comm_mygroup);
1758 do_flood(fplog,cr,x,f,ed,box,step,bNS);
1761 if (DOMAINDECOMP(cr))
1763 dd_force_flop_stop(cr->dd,nrnb);
1766 dd_cycles_add(cr->dd,cycles_force-cycles_pme,ddCyclF);
1772 if (IR_ELEC_FIELD(*inputrec))
1774 /* Compute forces due to electric field */
1775 calc_f_el(MASTER(cr) ? field : NULL,
1776 start,homenr,mdatoms->chargeA,x,fr->f_novirsum,
1777 inputrec->ex,inputrec->et,t);
1780 if (bDoAdressWF && fr->adress_icor == eAdressICThermoForce)
1782 /* Compute thermodynamic force in hybrid AdResS region */
1783 adress_thermo_force(start,homenr,&(top->cgs),x,fr->f_novirsum,fr,mdatoms,
1784 inputrec->ePBC==epbcNONE ? NULL : &pbc);
1787 /* Communicate the forces */
1790 wallcycle_start(wcycle,ewcMOVEF);
1791 if (DOMAINDECOMP(cr))
1793 dd_move_f(cr->dd,f,fr->fshift);
1794 /* Do we need to communicate the separate force array
1795 * for terms that do not contribute to the single sum virial?
1796 * Position restraints and electric fields do not introduce
1797 * inter-cg forces, only full electrostatics methods do.
1798 * When we do not calculate the virial, fr->f_novirsum = f,
1799 * so we have already communicated these forces.
1801 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1802 (flags & GMX_FORCE_VIRIAL))
1804 dd_move_f(cr->dd,fr->f_novirsum,NULL);
1808 /* We should not update the shift forces here,
1809 * since f_twin is already included in f.
1811 dd_move_f(cr->dd,fr->f_twin,NULL);
1816 pd_move_f(cr,f,nrnb);
1819 pd_move_f(cr,fr->f_twin,nrnb);
1822 wallcycle_stop(wcycle,ewcMOVEF);
1825 /* If we have NoVirSum forces, but we do not calculate the virial,
1826 * we sum fr->f_novirum=f later.
1828 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1830 wallcycle_start(wcycle,ewcVSITESPREAD);
1831 spread_vsite_f(fplog,vsite,x,f,fr->fshift,FALSE,NULL,nrnb,
1832 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
1833 wallcycle_stop(wcycle,ewcVSITESPREAD);
1837 wallcycle_start(wcycle,ewcVSITESPREAD);
1838 spread_vsite_f(fplog,vsite,x,fr->f_twin,NULL,FALSE,NULL,
1840 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
1841 wallcycle_stop(wcycle,ewcVSITESPREAD);
1845 if (flags & GMX_FORCE_VIRIAL)
1847 /* Calculation of the virial must be done after vsites! */
1848 calc_virial(fplog,mdatoms->start,mdatoms->homenr,x,f,
1849 vir_force,graph,box,nrnb,fr,inputrec->ePBC);
1853 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1855 pull_potential_wrapper(fplog,bSepDVDL,cr,inputrec,box,x,
1856 f,vir_force,mdatoms,enerd,lambda,t);
1859 /* Add the forces from enforced rotation potentials (if any) */
1862 wallcycle_start(wcycle,ewcROTadd);
1863 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr,step,t);
1864 wallcycle_stop(wcycle,ewcROTadd);
1867 if (PAR(cr) && !(cr->duty & DUTY_PME))
1869 /* In case of node-splitting, the PP nodes receive the long-range
1870 * forces, virial and energy from the PME nodes here.
1872 pme_receive_force_ener(fplog,bSepDVDL,cr,wcycle,enerd,fr);
1877 post_process_forces(fplog,cr,step,nrnb,wcycle,
1878 top,box,x,f,vir_force,mdatoms,graph,fr,vsite,
1882 /* Sum the potential energy terms from group contributions */
1883 sum_epot(&(inputrec->opts),enerd);
1886 void do_force(FILE *fplog,t_commrec *cr,
1887 t_inputrec *inputrec,
1888 gmx_large_int_t step,t_nrnb *nrnb,gmx_wallcycle_t wcycle,
1889 gmx_localtop_t *top,
1891 gmx_groups_t *groups,
1892 matrix box,rvec x[],history_t *hist,
1896 gmx_enerdata_t *enerd,t_fcdata *fcd,
1897 real *lambda,t_graph *graph,
1899 gmx_vsite_t *vsite,rvec mu_tot,
1900 double t,FILE *field,gmx_edsam_t ed,
1901 gmx_bool bBornRadii,
1904 /* modify force flag if not doing nonbonded */
1905 if (!fr->bNonbonded)
1907 flags &= ~GMX_FORCE_NONBONDED;
1910 switch (inputrec->cutoff_scheme)
1913 do_force_cutsVERLET(fplog, cr, inputrec,
1929 do_force_cutsGROUP(fplog, cr, inputrec,
1944 gmx_incons("Invalid cut-off scheme passed!");
1949 void do_constrain_first(FILE *fplog,gmx_constr_t constr,
1950 t_inputrec *ir,t_mdatoms *md,
1951 t_state *state,rvec *f,
1952 t_graph *graph,t_commrec *cr,t_nrnb *nrnb,
1953 t_forcerec *fr, gmx_localtop_t *top, tensor shake_vir)
1956 gmx_large_int_t step;
1957 real dt=ir->delta_t;
1961 snew(savex,state->natoms);
1964 end = md->homenr + start;
1967 fprintf(debug,"vcm: start=%d, homenr=%d, end=%d\n",
1968 start,md->homenr,end);
1969 /* Do a first constrain to reset particles... */
1970 step = ir->init_step;
1973 char buf[STEPSTRSIZE];
1974 fprintf(fplog,"\nConstraining the starting coordinates (step %s)\n",
1975 gmx_step_str(step,buf));
1979 /* constrain the current position */
1980 constrain(NULL,TRUE,FALSE,constr,&(top->idef),
1981 ir,NULL,cr,step,0,md,
1982 state->x,state->x,NULL,
1983 fr->bMolPBC,state->box,
1984 state->lambda[efptBONDED],&dvdl_dum,
1985 NULL,NULL,nrnb,econqCoord,
1986 ir->epc==epcMTTK,state->veta,state->veta);
1989 /* constrain the inital velocity, and save it */
1990 /* also may be useful if we need the ekin from the halfstep for velocity verlet */
1991 /* might not yet treat veta correctly */
1992 constrain(NULL,TRUE,FALSE,constr,&(top->idef),
1993 ir,NULL,cr,step,0,md,
1994 state->x,state->v,state->v,
1995 fr->bMolPBC,state->box,
1996 state->lambda[efptBONDED],&dvdl_dum,
1997 NULL,NULL,nrnb,econqVeloc,
1998 ir->epc==epcMTTK,state->veta,state->veta);
2000 /* constrain the inital velocities at t-dt/2 */
2001 if (EI_STATE_VELOCITY(ir->eI) && ir->eI!=eiVV)
2003 for(i=start; (i<end); i++)
2005 for(m=0; (m<DIM); m++)
2007 /* Reverse the velocity */
2008 state->v[i][m] = -state->v[i][m];
2009 /* Store the position at t-dt in buf */
2010 savex[i][m] = state->x[i][m] + dt*state->v[i][m];
2013 /* Shake the positions at t=-dt with the positions at t=0
2014 * as reference coordinates.
2018 char buf[STEPSTRSIZE];
2019 fprintf(fplog,"\nConstraining the coordinates at t0-dt (step %s)\n",
2020 gmx_step_str(step,buf));
2023 constrain(NULL,TRUE,FALSE,constr,&(top->idef),
2024 ir,NULL,cr,step,-1,md,
2025 state->x,savex,NULL,
2026 fr->bMolPBC,state->box,
2027 state->lambda[efptBONDED],&dvdl_dum,
2028 state->v,NULL,nrnb,econqCoord,
2029 ir->epc==epcMTTK,state->veta,state->veta);
2031 for(i=start; i<end; i++) {
2032 for(m=0; m<DIM; m++) {
2033 /* Re-reverse the velocities */
2034 state->v[i][m] = -state->v[i][m];
2041 void calc_enervirdiff(FILE *fplog,int eDispCorr,t_forcerec *fr)
2043 double eners[2],virs[2],enersum,virsum,y0,f,g,h;
2044 double r0,r1,r,rc3,rc9,ea,eb,ec,pa,pb,pc,pd;
2045 double invscale,invscale2,invscale3;
2046 int ri0,ri1,ri,i,offstart,offset;
2047 real scale,*vdwtab,tabfactor,tmp;
2049 fr->enershiftsix = 0;
2050 fr->enershifttwelve = 0;
2051 fr->enerdiffsix = 0;
2052 fr->enerdifftwelve = 0;
2054 fr->virdifftwelve = 0;
2056 if (eDispCorr != edispcNO) {
2057 for(i=0; i<2; i++) {
2061 if ((fr->vdwtype == evdwSWITCH) || (fr->vdwtype == evdwSHIFT)) {
2062 if (fr->rvdw_switch == 0)
2064 "With dispersion correction rvdw-switch can not be zero "
2065 "for vdw-type = %s",evdw_names[fr->vdwtype]);
2067 scale = fr->nblists[0].table_elec_vdw.scale;
2068 vdwtab = fr->nblists[0].table_vdw.data;
2070 /* Round the cut-offs to exact table values for precision */
2071 ri0 = floor(fr->rvdw_switch*scale);
2072 ri1 = ceil(fr->rvdw*scale);
2078 if (fr->vdwtype == evdwSHIFT)
2080 /* Determine the constant energy shift below rvdw_switch.
2081 * Table has a scale factor since we have scaled it down to compensate
2082 * for scaling-up c6/c12 with the derivative factors to save flops in analytical kernels.
2084 fr->enershiftsix = (real)(-1.0/(rc3*rc3)) - 6.0*vdwtab[8*ri0];
2085 fr->enershifttwelve = (real)( 1.0/(rc9*rc3)) - 12.0*vdwtab[8*ri0 + 4];
2087 /* Add the constant part from 0 to rvdw_switch.
2088 * This integration from 0 to rvdw_switch overcounts the number
2089 * of interactions by 1, as it also counts the self interaction.
2090 * We will correct for this later.
2092 eners[0] += 4.0*M_PI*fr->enershiftsix*rc3/3.0;
2093 eners[1] += 4.0*M_PI*fr->enershifttwelve*rc3/3.0;
2095 invscale = 1.0/(scale);
2096 invscale2 = invscale*invscale;
2097 invscale3 = invscale*invscale2;
2099 /* following summation derived from cubic spline definition,
2100 Numerical Recipies in C, second edition, p. 113-116. Exact
2101 for the cubic spline. We first calculate the negative of
2102 the energy from rvdw to rvdw_switch, assuming that g(r)=1,
2103 and then add the more standard, abrupt cutoff correction to
2104 that result, yielding the long-range correction for a
2105 switched function. We perform both the pressure and energy
2106 loops at the same time for simplicity, as the computational
2110 enersum = 0.0; virsum = 0.0;
2114 /* Since the dispersion table has been scaled down a factor 6.0 and the repulsion
2115 * a factor 12.0 to compensate for the c6/c12 parameters inside nbfp[] being scaled
2116 * up (to save flops in kernels), we need to correct for this.
2125 for (ri=ri0; ri<ri1; ri++) {
2128 eb = 2.0*invscale2*r;
2132 pb = 3.0*invscale2*r;
2133 pc = 3.0*invscale*r*r;
2136 /* this "8" is from the packing in the vdwtab array - perhaps should be #define'ed? */
2137 offset = 8*ri + offstart;
2138 y0 = vdwtab[offset];
2139 f = vdwtab[offset+1];
2140 g = vdwtab[offset+2];
2141 h = vdwtab[offset+3];
2143 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);
2144 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);
2147 enersum *= 4.0*M_PI*tabfactor;
2148 virsum *= 4.0*M_PI*tabfactor;
2149 eners[i] -= enersum;
2153 /* now add the correction for rvdw_switch to infinity */
2154 eners[0] += -4.0*M_PI/(3.0*rc3);
2155 eners[1] += 4.0*M_PI/(9.0*rc9);
2156 virs[0] += 8.0*M_PI/rc3;
2157 virs[1] += -16.0*M_PI/(3.0*rc9);
2159 else if ((fr->vdwtype == evdwCUT) || (fr->vdwtype == evdwUSER)) {
2160 if (fr->vdwtype == evdwUSER && fplog)
2162 "WARNING: using dispersion correction with user tables\n");
2163 rc3 = fr->rvdw*fr->rvdw*fr->rvdw;
2165 /* Contribution beyond the cut-off */
2166 eners[0] += -4.0*M_PI/(3.0*rc3);
2167 eners[1] += 4.0*M_PI/(9.0*rc9);
2168 if (fr->vdw_modifier==eintmodPOTSHIFT) {
2169 /* Contribution within the cut-off */
2170 eners[0] += -4.0*M_PI/(3.0*rc3);
2171 eners[1] += 4.0*M_PI/(3.0*rc9);
2173 /* Contribution beyond the cut-off */
2174 virs[0] += 8.0*M_PI/rc3;
2175 virs[1] += -16.0*M_PI/(3.0*rc9);
2178 "Dispersion correction is not implemented for vdw-type = %s",
2179 evdw_names[fr->vdwtype]);
2181 fr->enerdiffsix = eners[0];
2182 fr->enerdifftwelve = eners[1];
2183 /* The 0.5 is due to the Gromacs definition of the virial */
2184 fr->virdiffsix = 0.5*virs[0];
2185 fr->virdifftwelve = 0.5*virs[1];
2189 void calc_dispcorr(FILE *fplog,t_inputrec *ir,t_forcerec *fr,
2190 gmx_large_int_t step,int natoms,
2191 matrix box,real lambda,tensor pres,tensor virial,
2192 real *prescorr, real *enercorr, real *dvdlcorr)
2194 gmx_bool bCorrAll,bCorrPres;
2195 real dvdlambda,invvol,dens,ninter,avcsix,avctwelve,enerdiff,svir=0,spres=0;
2205 if (ir->eDispCorr != edispcNO) {
2206 bCorrAll = (ir->eDispCorr == edispcAllEner ||
2207 ir->eDispCorr == edispcAllEnerPres);
2208 bCorrPres = (ir->eDispCorr == edispcEnerPres ||
2209 ir->eDispCorr == edispcAllEnerPres);
2211 invvol = 1/det(box);
2214 /* Only correct for the interactions with the inserted molecule */
2215 dens = (natoms - fr->n_tpi)*invvol;
2220 dens = natoms*invvol;
2221 ninter = 0.5*natoms;
2224 if (ir->efep == efepNO)
2226 avcsix = fr->avcsix[0];
2227 avctwelve = fr->avctwelve[0];
2231 avcsix = (1 - lambda)*fr->avcsix[0] + lambda*fr->avcsix[1];
2232 avctwelve = (1 - lambda)*fr->avctwelve[0] + lambda*fr->avctwelve[1];
2235 enerdiff = ninter*(dens*fr->enerdiffsix - fr->enershiftsix);
2236 *enercorr += avcsix*enerdiff;
2238 if (ir->efep != efepNO)
2240 dvdlambda += (fr->avcsix[1] - fr->avcsix[0])*enerdiff;
2244 enerdiff = ninter*(dens*fr->enerdifftwelve - fr->enershifttwelve);
2245 *enercorr += avctwelve*enerdiff;
2246 if (fr->efep != efepNO)
2248 dvdlambda += (fr->avctwelve[1] - fr->avctwelve[0])*enerdiff;
2254 svir = ninter*dens*avcsix*fr->virdiffsix/3.0;
2255 if (ir->eDispCorr == edispcAllEnerPres)
2257 svir += ninter*dens*avctwelve*fr->virdifftwelve/3.0;
2259 /* The factor 2 is because of the Gromacs virial definition */
2260 spres = -2.0*invvol*svir*PRESFAC;
2262 for(m=0; m<DIM; m++) {
2263 virial[m][m] += svir;
2264 pres[m][m] += spres;
2269 /* Can't currently control when it prints, for now, just print when degugging */
2273 fprintf(debug,"Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
2279 "Long Range LJ corr.: Epot %10g, Pres: %10g, Vir: %10g\n",
2280 *enercorr,spres,svir);
2284 fprintf(debug,"Long Range LJ corr.: Epot %10g\n",*enercorr);
2288 if (fr->bSepDVDL && do_per_step(step,ir->nstlog))
2290 fprintf(fplog,sepdvdlformat,"Dispersion correction",
2291 *enercorr,dvdlambda);
2293 if (fr->efep != efepNO)
2295 *dvdlcorr += dvdlambda;
2300 void do_pbc_first(FILE *fplog,matrix box,t_forcerec *fr,
2301 t_graph *graph,rvec x[])
2304 fprintf(fplog,"Removing pbc first time\n");
2305 calc_shifts(box,fr->shift_vec);
2307 mk_mshift(fplog,graph,fr->ePBC,box,x);
2309 p_graph(debug,"do_pbc_first 1",graph);
2310 shift_self(graph,box,x);
2311 /* By doing an extra mk_mshift the molecules that are broken
2312 * because they were e.g. imported from another software
2313 * will be made whole again. Such are the healing powers
2316 mk_mshift(fplog,graph,fr->ePBC,box,x);
2318 p_graph(debug,"do_pbc_first 2",graph);
2321 fprintf(fplog,"Done rmpbc\n");
2324 static void low_do_pbc_mtop(FILE *fplog,int ePBC,matrix box,
2325 gmx_mtop_t *mtop,rvec x[],
2330 gmx_molblock_t *molb;
2332 if (bFirst && fplog)
2333 fprintf(fplog,"Removing pbc first time\n");
2337 for(mb=0; mb<mtop->nmolblock; mb++) {
2338 molb = &mtop->molblock[mb];
2339 if (molb->natoms_mol == 1 ||
2340 (!bFirst && mtop->moltype[molb->type].cgs.nr == 1)) {
2341 /* Just one atom or charge group in the molecule, no PBC required */
2342 as += molb->nmol*molb->natoms_mol;
2344 /* Pass NULL iso fplog to avoid graph prints for each molecule type */
2345 mk_graph_ilist(NULL,mtop->moltype[molb->type].ilist,
2346 0,molb->natoms_mol,FALSE,FALSE,graph);
2348 for(mol=0; mol<molb->nmol; mol++) {
2349 mk_mshift(fplog,graph,ePBC,box,x+as);
2351 shift_self(graph,box,x+as);
2352 /* The molecule is whole now.
2353 * We don't need the second mk_mshift call as in do_pbc_first,
2354 * since we no longer need this graph.
2357 as += molb->natoms_mol;
2365 void do_pbc_first_mtop(FILE *fplog,int ePBC,matrix box,
2366 gmx_mtop_t *mtop,rvec x[])
2368 low_do_pbc_mtop(fplog,ePBC,box,mtop,x,TRUE);
2371 void do_pbc_mtop(FILE *fplog,int ePBC,matrix box,
2372 gmx_mtop_t *mtop,rvec x[])
2374 low_do_pbc_mtop(fplog,ePBC,box,mtop,x,FALSE);
2377 void finish_run(FILE *fplog,t_commrec *cr,const char *confout,
2378 t_inputrec *inputrec,
2379 t_nrnb nrnb[],gmx_wallcycle_t wcycle,
2380 gmx_runtime_t *runtime,
2381 wallclock_gpu_t *gputimes,
2383 gmx_bool bWriteStat)
2386 t_nrnb *nrnb_tot=NULL;
2390 wallcycle_sum(cr,wcycle);
2396 MPI_Allreduce(nrnb->n,nrnb_tot->n,eNRNB,MPI_DOUBLE,MPI_SUM,
2397 cr->mpi_comm_mysim);
2405 #if defined(GMX_MPI) && !defined(GMX_THREAD_MPI)
2408 /* reduce nodetime over all MPI processes in the current simulation */
2410 MPI_Allreduce(&runtime->proctime,&sum,1,MPI_DOUBLE,MPI_SUM,
2411 cr->mpi_comm_mysim);
2412 runtime->proctime = sum;
2418 print_flop(fplog,nrnb_tot,&nbfs,&mflop);
2425 if ((cr->duty & DUTY_PP) && DOMAINDECOMP(cr))
2427 print_dd_statistics(cr,inputrec,fplog);
2439 snew(nrnb_all,cr->nnodes);
2440 nrnb_all[0] = *nrnb;
2441 for(s=1; s<cr->nnodes; s++)
2443 MPI_Recv(nrnb_all[s].n,eNRNB,MPI_DOUBLE,s,0,
2444 cr->mpi_comm_mysim,&stat);
2446 pr_load(fplog,cr,nrnb_all);
2451 MPI_Send(nrnb->n,eNRNB,MPI_DOUBLE,MASTERRANK(cr),0,
2452 cr->mpi_comm_mysim);
2459 wallcycle_print(fplog,cr->nnodes,cr->npmenodes,runtime->realtime,
2462 if (EI_DYNAMICS(inputrec->eI))
2464 delta_t = inputrec->delta_t;
2473 print_perf(fplog,runtime->proctime,runtime->realtime,
2474 cr->nnodes-cr->npmenodes,
2475 runtime->nsteps_done,delta_t,nbfs,mflop,
2480 print_perf(stderr,runtime->proctime,runtime->realtime,
2481 cr->nnodes-cr->npmenodes,
2482 runtime->nsteps_done,delta_t,nbfs,mflop,
2488 extern void initialize_lambdas(FILE *fplog,t_inputrec *ir,int *fep_state,real *lambda,double *lam0)
2490 /* this function works, but could probably use a logic rewrite to keep all the different
2491 types of efep straight. */
2494 t_lambda *fep = ir->fepvals;
2496 if ((ir->efep==efepNO) && (ir->bSimTemp == FALSE)) {
2497 for (i=0;i<efptNR;i++) {
2506 *fep_state = fep->init_fep_state; /* this might overwrite the checkpoint
2507 if checkpoint is set -- a kludge is in for now
2509 for (i=0;i<efptNR;i++)
2511 /* overwrite lambda state with init_lambda for now for backwards compatibility */
2512 if (fep->init_lambda>=0) /* if it's -1, it was never initializd */
2514 lambda[i] = fep->init_lambda;
2516 lam0[i] = lambda[i];
2521 lambda[i] = fep->all_lambda[i][*fep_state];
2523 lam0[i] = lambda[i];
2528 /* need to rescale control temperatures to match current state */
2529 for (i=0;i<ir->opts.ngtc;i++) {
2530 if (ir->opts.ref_t[i] > 0) {
2531 ir->opts.ref_t[i] = ir->simtempvals->temperatures[*fep_state];
2537 /* Send to the log the information on the current lambdas */
2540 fprintf(fplog,"Initial vector of lambda components:[ ");
2541 for (i=0;i<efptNR;i++)
2543 fprintf(fplog,"%10.4f ",lambda[i]);
2545 fprintf(fplog,"]\n");
2551 void init_md(FILE *fplog,
2552 t_commrec *cr,t_inputrec *ir,const output_env_t oenv,
2553 double *t,double *t0,
2554 real *lambda, int *fep_state, double *lam0,
2555 t_nrnb *nrnb,gmx_mtop_t *mtop,
2557 int nfile,const t_filenm fnm[],
2558 gmx_mdoutf_t **outf,t_mdebin **mdebin,
2559 tensor force_vir,tensor shake_vir,rvec mu_tot,
2560 gmx_bool *bSimAnn,t_vcm **vcm, t_state *state, unsigned long Flags)
2565 /* Initial values */
2566 *t = *t0 = ir->init_t;
2569 for(i=0;i<ir->opts.ngtc;i++)
2571 /* set bSimAnn if any group is being annealed */
2572 if(ir->opts.annealing[i]!=eannNO)
2579 update_annealing_target_temp(&(ir->opts),ir->init_t);
2582 /* Initialize lambda variables */
2583 initialize_lambdas(fplog,ir,fep_state,lambda,lam0);
2587 *upd = init_update(fplog,ir);
2593 *vcm = init_vcm(fplog,&mtop->groups,ir);
2596 if (EI_DYNAMICS(ir->eI) && !(Flags & MD_APPENDFILES))
2598 if (ir->etc == etcBERENDSEN)
2600 please_cite(fplog,"Berendsen84a");
2602 if (ir->etc == etcVRESCALE)
2604 please_cite(fplog,"Bussi2007a");
2612 *outf = init_mdoutf(nfile,fnm,Flags,cr,ir,oenv);
2614 *mdebin = init_mdebin((Flags & MD_APPENDFILES) ? NULL : (*outf)->fp_ene,
2615 mtop,ir, (*outf)->fp_dhdl);
2620 please_cite(fplog,"Fritsch12");
2621 please_cite(fplog,"Junghans10");
2623 /* Initiate variables */
2624 clear_mat(force_vir);
2625 clear_mat(shake_vir);