1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * GROwing Monsters And Cloning Shrimps
41 #include <catamount/dclock.h>
47 #ifdef HAVE_SYS_TIME_H
60 #include "chargegroup.h"
83 #include "pull_rotation.h"
84 #include "gmx_random.h"
87 #include "gmx_wallcycle.h"
89 #include "nbnxn_atomdata.h"
90 #include "nbnxn_search.h"
91 #include "nbnxn_kernels/nbnxn_kernel_ref.h"
92 #include "nbnxn_kernels/nbnxn_kernel_simd_4xn.h"
93 #include "nbnxn_kernels/nbnxn_kernel_simd_2xnn.h"
94 #include "nbnxn_kernels/nbnxn_kernel_gpu_ref.h"
106 #include "nbnxn_cuda_data_mgmt.h"
107 #include "nbnxn_cuda/nbnxn_cuda.h"
110 typedef struct gmx_timeprint {
115 /* Portable version of ctime_r implemented in src/gmxlib/string2.c, but we do not want it declared in public installed headers */
117 gmx_ctime_r(const time_t *clock, char *buf, int n);
123 #ifdef HAVE_GETTIMEOFDAY
127 gettimeofday(&t, NULL);
129 seconds = (double) t.tv_sec + 1e-6*(double)t.tv_usec;
135 seconds = time(NULL);
142 #define difftime(end, start) ((double)(end)-(double)(start))
144 void print_time(FILE *out, gmx_runtime_t *runtime, gmx_large_int_t step,
145 t_inputrec *ir, t_commrec *cr)
148 char timebuf[STRLEN];
152 #ifndef GMX_THREAD_MPI
158 fprintf(out, "step %s", gmx_step_str(step, buf));
159 if ((step >= ir->nstlist))
161 runtime->last = gmx_gettime();
162 dt = difftime(runtime->last, runtime->real);
163 runtime->time_per_step = dt/(step - ir->init_step + 1);
165 dt = (ir->nsteps + ir->init_step - step)*runtime->time_per_step;
171 finish = (time_t) (runtime->last + dt);
172 gmx_ctime_r(&finish, timebuf, STRLEN);
173 sprintf(buf, "%s", timebuf);
174 buf[strlen(buf)-1] = '\0';
175 fprintf(out, ", will finish %s", buf);
179 fprintf(out, ", remaining runtime: %5d s ", (int)dt);
184 fprintf(out, " performance: %.1f ns/day ",
185 ir->delta_t/1000*24*60*60/runtime->time_per_step);
188 #ifndef GMX_THREAD_MPI
202 static double set_proctime(gmx_runtime_t *runtime)
208 prev = runtime->proc;
209 runtime->proc = dclock();
211 diff = runtime->proc - prev;
215 prev = runtime->proc;
216 runtime->proc = clock();
218 diff = (double)(runtime->proc - prev)/(double)CLOCKS_PER_SEC;
222 /* The counter has probably looped, ignore this data */
229 void runtime_start(gmx_runtime_t *runtime)
231 runtime->real = gmx_gettime();
233 set_proctime(runtime);
234 runtime->realtime = 0;
235 runtime->proctime = 0;
237 runtime->time_per_step = 0;
240 void runtime_end(gmx_runtime_t *runtime)
246 runtime->proctime += set_proctime(runtime);
247 runtime->realtime = now - runtime->real;
251 void runtime_upd_proc(gmx_runtime_t *runtime)
253 runtime->proctime += set_proctime(runtime);
256 void print_date_and_time(FILE *fplog, int nodeid, const char *title,
257 const gmx_runtime_t *runtime)
260 char timebuf[STRLEN];
261 char time_string[STRLEN];
268 tmptime = (time_t) runtime->real;
269 gmx_ctime_r(&tmptime, timebuf, STRLEN);
273 tmptime = (time_t) gmx_gettime();
274 gmx_ctime_r(&tmptime, timebuf, STRLEN);
276 for (i = 0; timebuf[i] >= ' '; i++)
278 time_string[i] = timebuf[i];
280 time_string[i] = '\0';
282 fprintf(fplog, "%s on node %d %s\n", title, nodeid, time_string);
286 static void sum_forces(int start, int end, rvec f[], rvec flr[])
292 pr_rvecs(debug, 0, "fsr", f+start, end-start);
293 pr_rvecs(debug, 0, "flr", flr+start, end-start);
295 for (i = start; (i < end); i++)
297 rvec_inc(f[i], flr[i]);
302 * calc_f_el calculates forces due to an electric field.
304 * force is kJ mol^-1 nm^-1 = e * kJ mol^-1 nm^-1 / e
306 * Et[] contains the parameters for the time dependent
307 * part of the field (not yet used).
308 * Ex[] contains the parameters for
309 * the spatial dependent part of the field. You can have cool periodic
310 * fields in principle, but only a constant field is supported
312 * The function should return the energy due to the electric field
313 * (if any) but for now returns 0.
316 * There can be problems with the virial.
317 * Since the field is not self-consistent this is unavoidable.
318 * For neutral molecules the virial is correct within this approximation.
319 * For neutral systems with many charged molecules the error is small.
320 * But for systems with a net charge or a few charged molecules
321 * the error can be significant when the field is high.
322 * Solution: implement a self-consitent electric field into PME.
324 static void calc_f_el(FILE *fp, int start, int homenr,
325 real charge[], rvec x[], rvec f[],
326 t_cosines Ex[], t_cosines Et[], double t)
332 for (m = 0; (m < DIM); m++)
339 Ext[m] = cos(Et[m].a[0]*(t-t0))*exp(-sqr(t-t0)/(2.0*sqr(Et[m].a[2])));
343 Ext[m] = cos(Et[m].a[0]*t);
352 /* Convert the field strength from V/nm to MD-units */
353 Ext[m] *= Ex[m].a[0]*FIELDFAC;
354 for (i = start; (i < start+homenr); i++)
356 f[i][m] += charge[i]*Ext[m];
366 fprintf(fp, "%10g %10g %10g %10g #FIELD\n", t,
367 Ext[XX]/FIELDFAC, Ext[YY]/FIELDFAC, Ext[ZZ]/FIELDFAC);
371 static void calc_virial(FILE *fplog, int start, int homenr, rvec x[], rvec f[],
372 tensor vir_part, t_graph *graph, matrix box,
373 t_nrnb *nrnb, const t_forcerec *fr, int ePBC)
378 /* The short-range virial from surrounding boxes */
380 calc_vir(fplog, SHIFTS, fr->shift_vec, fr->fshift, vir_part, ePBC == epbcSCREW, box);
381 inc_nrnb(nrnb, eNR_VIRIAL, SHIFTS);
383 /* Calculate partial virial, for local atoms only, based on short range.
384 * Total virial is computed in global_stat, called from do_md
386 f_calc_vir(fplog, start, start+homenr, x, f, vir_part, graph, box);
387 inc_nrnb(nrnb, eNR_VIRIAL, homenr);
389 /* Add position restraint contribution */
390 for (i = 0; i < DIM; i++)
392 vir_part[i][i] += fr->vir_diag_posres[i];
395 /* Add wall contribution */
396 for (i = 0; i < DIM; i++)
398 vir_part[i][ZZ] += fr->vir_wall_z[i];
403 pr_rvecs(debug, 0, "vir_part", vir_part, DIM);
407 static void posres_wrapper(FILE *fplog,
413 matrix box, rvec x[],
415 gmx_enerdata_t *enerd,
423 /* Position restraints always require full pbc */
424 set_pbc(&pbc, ir->ePBC, box);
426 v = posres(top->idef.il[F_POSRES].nr, top->idef.il[F_POSRES].iatoms,
427 top->idef.iparams_posres,
428 (const rvec*)x, fr->f_novirsum, fr->vir_diag_posres,
429 ir->ePBC == epbcNONE ? NULL : &pbc,
430 lambda[efptRESTRAINT], &dvdl,
431 fr->rc_scaling, fr->ePBC, fr->posres_com, fr->posres_comB);
434 fprintf(fplog, sepdvdlformat,
435 interaction_function[F_POSRES].longname, v, dvdl);
437 enerd->term[F_POSRES] += v;
438 /* If just the force constant changes, the FEP term is linear,
439 * but if k changes, it is not.
441 enerd->dvdl_nonlin[efptRESTRAINT] += dvdl;
442 inc_nrnb(nrnb, eNR_POSRES, top->idef.il[F_POSRES].nr/2);
444 if ((ir->fepvals->n_lambda > 0) && (flags & GMX_FORCE_DHDL))
446 for (i = 0; i < enerd->n_lambda; i++)
448 real dvdl_dum, lambda_dum;
450 lambda_dum = (i == 0 ? lambda[efptRESTRAINT] : ir->fepvals->all_lambda[efptRESTRAINT][i-1]);
451 v = posres(top->idef.il[F_POSRES].nr, top->idef.il[F_POSRES].iatoms,
452 top->idef.iparams_posres,
453 (const rvec*)x, NULL, NULL,
454 ir->ePBC == epbcNONE ? NULL : &pbc, lambda_dum, &dvdl,
455 fr->rc_scaling, fr->ePBC, fr->posres_com, fr->posres_comB);
456 enerd->enerpart_lambda[i] += v;
461 static void pull_potential_wrapper(FILE *fplog,
465 matrix box, rvec x[],
469 gmx_enerdata_t *enerd,
476 /* Calculate the center of mass forces, this requires communication,
477 * which is why pull_potential is called close to other communication.
478 * The virial contribution is calculated directly,
479 * which is why we call pull_potential after calc_virial.
481 set_pbc(&pbc, ir->ePBC, box);
483 enerd->term[F_COM_PULL] +=
484 pull_potential(ir->ePull, ir->pull, mdatoms, &pbc,
485 cr, t, lambda[efptRESTRAINT], x, f, vir_force, &dvdl);
488 fprintf(fplog, sepdvdlformat, "Com pull", enerd->term[F_COM_PULL], dvdl);
490 enerd->dvdl_lin[efptRESTRAINT] += dvdl;
493 static void pme_receive_force_ener(FILE *fplog,
496 gmx_wallcycle_t wcycle,
497 gmx_enerdata_t *enerd,
501 float cycles_ppdpme, cycles_seppme;
503 cycles_ppdpme = wallcycle_stop(wcycle, ewcPPDURINGPME);
504 dd_cycles_add(cr->dd, cycles_ppdpme, ddCyclPPduringPME);
506 /* In case of node-splitting, the PP nodes receive the long-range
507 * forces, virial and energy from the PME nodes here.
509 wallcycle_start(wcycle, ewcPP_PMEWAITRECVF);
511 gmx_pme_receive_f(cr, fr->f_novirsum, fr->vir_el_recip, &e, &dvdl,
515 fprintf(fplog, sepdvdlformat, "PME mesh", e, dvdl);
517 enerd->term[F_COUL_RECIP] += e;
518 enerd->dvdl_lin[efptCOUL] += dvdl;
521 dd_cycles_add(cr->dd, cycles_seppme, ddCyclPME);
523 wallcycle_stop(wcycle, ewcPP_PMEWAITRECVF);
526 static void print_large_forces(FILE *fp, t_mdatoms *md, t_commrec *cr,
527 gmx_large_int_t step, real pforce, rvec *x, rvec *f)
531 char buf[STEPSTRSIZE];
534 for (i = md->start; i < md->start+md->homenr; i++)
537 /* We also catch NAN, if the compiler does not optimize this away. */
538 if (fn2 >= pf2 || fn2 != fn2)
540 fprintf(fp, "step %s atom %6d x %8.3f %8.3f %8.3f force %12.5e\n",
541 gmx_step_str(step, buf),
542 ddglatnr(cr->dd, i), x[i][XX], x[i][YY], x[i][ZZ], sqrt(fn2));
547 static void post_process_forces(FILE *fplog,
549 gmx_large_int_t step,
550 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
552 matrix box, rvec x[],
557 t_forcerec *fr, gmx_vsite_t *vsite,
564 /* Spread the mesh force on virtual sites to the other particles...
565 * This is parallellized. MPI communication is performed
566 * if the constructing atoms aren't local.
568 wallcycle_start(wcycle, ewcVSITESPREAD);
569 spread_vsite_f(fplog, vsite, x, fr->f_novirsum, NULL,
570 (flags & GMX_FORCE_VIRIAL), fr->vir_el_recip,
572 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
573 wallcycle_stop(wcycle, ewcVSITESPREAD);
575 if (flags & GMX_FORCE_VIRIAL)
577 /* Now add the forces, this is local */
580 sum_forces(0, fr->f_novirsum_n, f, fr->f_novirsum);
584 sum_forces(mdatoms->start, mdatoms->start+mdatoms->homenr,
587 if (EEL_FULL(fr->eeltype))
589 /* Add the mesh contribution to the virial */
590 m_add(vir_force, fr->vir_el_recip, vir_force);
594 pr_rvecs(debug, 0, "vir_force", vir_force, DIM);
599 if (fr->print_force >= 0)
601 print_large_forces(stderr, mdatoms, cr, step, fr->print_force, x, f);
605 static void do_nb_verlet(t_forcerec *fr,
606 interaction_const_t *ic,
607 gmx_enerdata_t *enerd,
608 int flags, int ilocality,
611 gmx_wallcycle_t wcycle)
613 int nnbl, kernel_type, enr_nbnxn_kernel_ljc, enr_nbnxn_kernel_lj;
615 nonbonded_verlet_group_t *nbvg;
617 if (!(flags & GMX_FORCE_NONBONDED))
619 /* skip non-bonded calculation */
623 nbvg = &fr->nbv->grp[ilocality];
625 /* CUDA kernel launch overhead is already timed separately */
626 if (fr->cutoff_scheme != ecutsVERLET)
628 gmx_incons("Invalid cut-off scheme passed!");
631 if (nbvg->kernel_type != nbnxnk8x8x8_CUDA)
633 wallcycle_sub_start(wcycle, ewcsNONBONDED);
635 switch (nbvg->kernel_type)
637 case nbnxnk4x4_PlainC:
638 nbnxn_kernel_ref(&nbvg->nbl_lists,
644 enerd->grpp.ener[egCOULSR],
646 enerd->grpp.ener[egBHAMSR] :
647 enerd->grpp.ener[egLJSR]);
650 case nbnxnk4xN_SIMD_4xN:
651 nbnxn_kernel_simd_4xn(&nbvg->nbl_lists,
658 enerd->grpp.ener[egCOULSR],
660 enerd->grpp.ener[egBHAMSR] :
661 enerd->grpp.ener[egLJSR]);
663 case nbnxnk4xN_SIMD_2xNN:
664 nbnxn_kernel_simd_2xnn(&nbvg->nbl_lists,
671 enerd->grpp.ener[egCOULSR],
673 enerd->grpp.ener[egBHAMSR] :
674 enerd->grpp.ener[egLJSR]);
677 case nbnxnk8x8x8_CUDA:
678 nbnxn_cuda_launch_kernel(fr->nbv->cu_nbv, nbvg->nbat, flags, ilocality);
681 case nbnxnk8x8x8_PlainC:
682 nbnxn_kernel_gpu_ref(nbvg->nbl_lists.nbl[0],
687 nbvg->nbat->out[0].f,
689 enerd->grpp.ener[egCOULSR],
691 enerd->grpp.ener[egBHAMSR] :
692 enerd->grpp.ener[egLJSR]);
696 gmx_incons("Invalid nonbonded kernel type passed!");
699 if (nbvg->kernel_type != nbnxnk8x8x8_CUDA)
701 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
704 if (EEL_RF(ic->eeltype) || ic->eeltype == eelCUT)
706 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_RF;
708 else if (nbvg->ewald_excl == ewaldexclTable)
710 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_TAB;
714 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_EWALD;
716 enr_nbnxn_kernel_lj = eNR_NBNXN_LJ;
717 if (flags & GMX_FORCE_ENERGY)
719 /* In eNR_??? the nbnxn F+E kernels are always the F kernel + 1 */
720 enr_nbnxn_kernel_ljc += 1;
721 enr_nbnxn_kernel_lj += 1;
724 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc,
725 nbvg->nbl_lists.natpair_ljq);
726 inc_nrnb(nrnb, enr_nbnxn_kernel_lj,
727 nbvg->nbl_lists.natpair_lj);
728 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc-eNR_NBNXN_LJ_RF+eNR_NBNXN_RF,
729 nbvg->nbl_lists.natpair_q);
732 void do_force_cutsVERLET(FILE *fplog, t_commrec *cr,
733 t_inputrec *inputrec,
734 gmx_large_int_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
737 gmx_groups_t *groups,
738 matrix box, rvec x[], history_t *hist,
742 gmx_enerdata_t *enerd, t_fcdata *fcd,
743 real *lambda, t_graph *graph,
744 t_forcerec *fr, interaction_const_t *ic,
745 gmx_vsite_t *vsite, rvec mu_tot,
746 double t, FILE *field, gmx_edsam_t ed,
754 gmx_bool bSepDVDL, bStateChanged, bNS, bFillGrid, bCalcCGCM, bBS;
755 gmx_bool bDoLongRange, bDoForces, bSepLRF, bUseGPU, bUseOrEmulGPU;
756 gmx_bool bDiffKernels = FALSE;
758 rvec vzero, box_diag;
760 float cycles_pme, cycles_force;
761 nonbonded_verlet_t *nbv;
765 nb_kernel_type = fr->nbv->grp[0].kernel_type;
767 start = mdatoms->start;
768 homenr = mdatoms->homenr;
770 bSepDVDL = (fr->bSepDVDL && do_per_step(step, inputrec->nstlog));
772 clear_mat(vir_force);
775 if (DOMAINDECOMP(cr))
777 cg1 = cr->dd->ncg_tot;
788 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
789 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
790 bFillGrid = (bNS && bStateChanged);
791 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
792 bDoLongRange = (fr->bTwinRange && bNS && (flags & GMX_FORCE_DO_LR));
793 bDoForces = (flags & GMX_FORCE_FORCES);
794 bSepLRF = (bDoLongRange && bDoForces && (flags & GMX_FORCE_SEPLRF));
795 bUseGPU = fr->nbv->bUseGPU;
796 bUseOrEmulGPU = bUseGPU || (nbv->grp[0].kernel_type == nbnxnk8x8x8_PlainC);
800 update_forcerec(fplog, fr, box);
802 if (NEED_MUTOT(*inputrec))
804 /* Calculate total (local) dipole moment in a temporary common array.
805 * This makes it possible to sum them over nodes faster.
807 calc_mu(start, homenr,
808 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
813 if (fr->ePBC != epbcNONE)
815 /* Compute shift vectors every step,
816 * because of pressure coupling or box deformation!
818 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
820 calc_shifts(box, fr->shift_vec);
825 put_atoms_in_box_omp(fr->ePBC, box, homenr, x);
826 inc_nrnb(nrnb, eNR_SHIFTX, homenr);
828 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
830 unshift_self(graph, box, x);
834 nbnxn_atomdata_copy_shiftvec(flags & GMX_FORCE_DYNAMICBOX,
835 fr->shift_vec, nbv->grp[0].nbat);
838 if (!(cr->duty & DUTY_PME))
840 /* Send particle coordinates to the pme nodes.
841 * Since this is only implemented for domain decomposition
842 * and domain decomposition does not use the graph,
843 * we do not need to worry about shifting.
846 wallcycle_start(wcycle, ewcPP_PMESENDX);
848 bBS = (inputrec->nwall == 2);
852 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
855 gmx_pme_send_x(cr, bBS ? boxs : box, x,
856 mdatoms->nChargePerturbed, lambda[efptCOUL],
857 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)), step);
859 wallcycle_stop(wcycle, ewcPP_PMESENDX);
863 /* do gridding for pair search */
866 if (graph && bStateChanged)
868 /* Calculate intramolecular shift vectors to make molecules whole */
869 mk_mshift(fplog, graph, fr->ePBC, box, x);
873 box_diag[XX] = box[XX][XX];
874 box_diag[YY] = box[YY][YY];
875 box_diag[ZZ] = box[ZZ][ZZ];
877 wallcycle_start(wcycle, ewcNS);
880 wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
881 nbnxn_put_on_grid(nbv->nbs, fr->ePBC, box,
883 0, mdatoms->homenr, -1, fr->cginfo, x,
885 nbv->grp[eintLocal].kernel_type,
886 nbv->grp[eintLocal].nbat);
887 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
891 wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
892 nbnxn_put_on_grid_nonlocal(nbv->nbs, domdec_zones(cr->dd),
894 nbv->grp[eintNonlocal].kernel_type,
895 nbv->grp[eintNonlocal].nbat);
896 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
899 if (nbv->ngrp == 1 ||
900 nbv->grp[eintNonlocal].nbat == nbv->grp[eintLocal].nbat)
902 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatAll,
903 nbv->nbs, mdatoms, fr->cginfo);
907 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatLocal,
908 nbv->nbs, mdatoms, fr->cginfo);
909 nbnxn_atomdata_set(nbv->grp[eintNonlocal].nbat, eatAll,
910 nbv->nbs, mdatoms, fr->cginfo);
912 wallcycle_stop(wcycle, ewcNS);
915 /* initialize the GPU atom data and copy shift vector */
920 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
921 nbnxn_cuda_init_atomdata(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
922 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
925 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
926 nbnxn_cuda_upload_shiftvec(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
927 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
930 /* do local pair search */
933 wallcycle_start_nocount(wcycle, ewcNS);
934 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
935 nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintLocal].nbat,
938 nbv->min_ci_balanced,
939 &nbv->grp[eintLocal].nbl_lists,
941 nbv->grp[eintLocal].kernel_type,
943 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
947 /* initialize local pair-list on the GPU */
948 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
949 nbv->grp[eintLocal].nbl_lists.nbl[0],
952 wallcycle_stop(wcycle, ewcNS);
956 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
957 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
958 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, FALSE, x,
959 nbv->grp[eintLocal].nbat);
960 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
961 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
966 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
967 /* launch local nonbonded F on GPU */
968 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFNo,
970 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
973 /* Communicate coordinates and sum dipole if necessary +
974 do non-local pair search */
975 if (DOMAINDECOMP(cr))
977 bDiffKernels = (nbv->grp[eintNonlocal].kernel_type !=
978 nbv->grp[eintLocal].kernel_type);
982 /* With GPU+CPU non-bonded calculations we need to copy
983 * the local coordinates to the non-local nbat struct
984 * (in CPU format) as the non-local kernel call also
985 * calculates the local - non-local interactions.
987 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
988 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
989 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, TRUE, x,
990 nbv->grp[eintNonlocal].nbat);
991 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
992 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
997 wallcycle_start_nocount(wcycle, ewcNS);
998 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1002 nbnxn_grid_add_simple(nbv->nbs, nbv->grp[eintNonlocal].nbat);
1005 nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintNonlocal].nbat,
1008 nbv->min_ci_balanced,
1009 &nbv->grp[eintNonlocal].nbl_lists,
1011 nbv->grp[eintNonlocal].kernel_type,
1014 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1016 if (nbv->grp[eintNonlocal].kernel_type == nbnxnk8x8x8_CUDA)
1018 /* initialize non-local pair-list on the GPU */
1019 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
1020 nbv->grp[eintNonlocal].nbl_lists.nbl[0],
1023 wallcycle_stop(wcycle, ewcNS);
1027 wallcycle_start(wcycle, ewcMOVEX);
1028 dd_move_x(cr->dd, box, x);
1030 /* When we don't need the total dipole we sum it in global_stat */
1031 if (bStateChanged && NEED_MUTOT(*inputrec))
1033 gmx_sumd(2*DIM, mu, cr);
1035 wallcycle_stop(wcycle, ewcMOVEX);
1037 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1038 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
1039 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatNonlocal, FALSE, x,
1040 nbv->grp[eintNonlocal].nbat);
1041 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
1042 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1045 if (bUseGPU && !bDiffKernels)
1047 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
1048 /* launch non-local nonbonded F on GPU */
1049 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFNo,
1051 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1057 /* launch D2H copy-back F */
1058 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1059 if (DOMAINDECOMP(cr) && !bDiffKernels)
1061 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintNonlocal].nbat,
1062 flags, eatNonlocal);
1064 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintLocal].nbat,
1066 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1069 if (bStateChanged && NEED_MUTOT(*inputrec))
1073 gmx_sumd(2*DIM, mu, cr);
1076 for (i = 0; i < 2; i++)
1078 for (j = 0; j < DIM; j++)
1080 fr->mu_tot[i][j] = mu[i*DIM + j];
1084 if (fr->efep == efepNO)
1086 copy_rvec(fr->mu_tot[0], mu_tot);
1090 for (j = 0; j < DIM; j++)
1093 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] +
1094 lambda[efptCOUL]*fr->mu_tot[1][j];
1098 /* Reset energies */
1099 reset_enerdata(&(inputrec->opts), fr, bNS, enerd, MASTER(cr));
1100 clear_rvecs(SHIFTS, fr->fshift);
1102 if (DOMAINDECOMP(cr))
1104 if (!(cr->duty & DUTY_PME))
1106 wallcycle_start(wcycle, ewcPPDURINGPME);
1107 dd_force_flop_start(cr->dd, nrnb);
1111 /* Start the force cycle counter.
1112 * This counter is stopped in do_forcelow_level.
1113 * No parallel communication should occur while this counter is running,
1114 * since that will interfere with the dynamic load balancing.
1116 wallcycle_start(wcycle, ewcFORCE);
1119 /* Reset forces for which the virial is calculated separately:
1120 * PME/Ewald forces if necessary */
1121 if (fr->bF_NoVirSum)
1123 if (flags & GMX_FORCE_VIRIAL)
1125 fr->f_novirsum = fr->f_novirsum_alloc;
1128 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1132 clear_rvecs(homenr, fr->f_novirsum+start);
1137 /* We are not calculating the pressure so we do not need
1138 * a separate array for forces that do not contribute
1145 /* Clear the short- and long-range forces */
1146 clear_rvecs(fr->natoms_force_constr, f);
1147 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1149 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1152 clear_rvec(fr->vir_diag_posres);
1154 if (inputrec->ePull == epullCONSTRAINT)
1156 clear_pull_forces(inputrec->pull);
1159 /* update QMMMrec, if necessary */
1162 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1165 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1167 posres_wrapper(fplog, flags, bSepDVDL, inputrec, nrnb, top, box, x,
1168 f, enerd, lambda, fr);
1171 /* Compute the bonded and non-bonded energies and optionally forces */
1172 do_force_lowlevel(fplog, step, fr, inputrec, &(top->idef),
1173 cr, nrnb, wcycle, mdatoms, &(inputrec->opts),
1174 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, mtop, top, fr->born,
1175 &(top->atomtypes), bBornRadii, box,
1176 inputrec->fepvals, lambda, graph, &(top->excls), fr->mu_tot,
1177 flags, &cycles_pme);
1181 if (do_per_step(step, inputrec->nstcalclr))
1183 /* Add the long range forces to the short range forces */
1184 for (i = 0; i < fr->natoms_force_constr; i++)
1186 rvec_add(fr->f_twin[i], f[i], f[i]);
1193 /* Maybe we should move this into do_force_lowlevel */
1194 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFYes,
1199 if (!bUseOrEmulGPU || bDiffKernels)
1203 if (DOMAINDECOMP(cr))
1205 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal,
1206 bDiffKernels ? enbvClearFYes : enbvClearFNo,
1216 aloc = eintNonlocal;
1219 /* Add all the non-bonded force to the normal force array.
1220 * This can be split into a local a non-local part when overlapping
1221 * communication with calculation with domain decomposition.
1223 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1224 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1225 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1226 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatAll, nbv->grp[aloc].nbat, f);
1227 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1228 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1229 wallcycle_start_nocount(wcycle, ewcFORCE);
1231 /* if there are multiple fshift output buffers reduce them */
1232 if ((flags & GMX_FORCE_VIRIAL) &&
1233 nbv->grp[aloc].nbl_lists.nnbl > 1)
1235 nbnxn_atomdata_add_nbat_fshift_to_fshift(nbv->grp[aloc].nbat,
1240 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1244 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1247 if (bUseOrEmulGPU && !bDiffKernels)
1249 /* wait for non-local forces (or calculate in emulation mode) */
1250 if (DOMAINDECOMP(cr))
1254 wallcycle_start(wcycle, ewcWAIT_GPU_NB_NL);
1255 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1256 nbv->grp[eintNonlocal].nbat,
1258 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1260 cycles_force += wallcycle_stop(wcycle, ewcWAIT_GPU_NB_NL);
1264 wallcycle_start_nocount(wcycle, ewcFORCE);
1265 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFYes,
1267 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1269 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1270 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1271 /* skip the reduction if there was no non-local work to do */
1272 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1274 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatNonlocal,
1275 nbv->grp[eintNonlocal].nbat, f);
1277 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1278 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1284 /* Communicate the forces */
1287 wallcycle_start(wcycle, ewcMOVEF);
1288 if (DOMAINDECOMP(cr))
1290 dd_move_f(cr->dd, f, fr->fshift);
1291 /* Do we need to communicate the separate force array
1292 * for terms that do not contribute to the single sum virial?
1293 * Position restraints and electric fields do not introduce
1294 * inter-cg forces, only full electrostatics methods do.
1295 * When we do not calculate the virial, fr->f_novirsum = f,
1296 * so we have already communicated these forces.
1298 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1299 (flags & GMX_FORCE_VIRIAL))
1301 dd_move_f(cr->dd, fr->f_novirsum, NULL);
1305 /* We should not update the shift forces here,
1306 * since f_twin is already included in f.
1308 dd_move_f(cr->dd, fr->f_twin, NULL);
1311 wallcycle_stop(wcycle, ewcMOVEF);
1317 /* wait for local forces (or calculate in emulation mode) */
1320 wallcycle_start(wcycle, ewcWAIT_GPU_NB_L);
1321 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1322 nbv->grp[eintLocal].nbat,
1324 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1326 wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
1328 /* now clear the GPU outputs while we finish the step on the CPU */
1330 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1331 nbnxn_cuda_clear_outputs(nbv->cu_nbv, flags);
1332 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1336 wallcycle_start_nocount(wcycle, ewcFORCE);
1337 do_nb_verlet(fr, ic, enerd, flags, eintLocal,
1338 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
1340 wallcycle_stop(wcycle, ewcFORCE);
1342 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1343 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1344 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1346 /* skip the reduction if there was no non-local work to do */
1347 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatLocal,
1348 nbv->grp[eintLocal].nbat, f);
1350 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1351 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1354 if (DOMAINDECOMP(cr))
1356 dd_force_flop_stop(cr->dd, nrnb);
1359 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1365 if (IR_ELEC_FIELD(*inputrec))
1367 /* Compute forces due to electric field */
1368 calc_f_el(MASTER(cr) ? field : NULL,
1369 start, homenr, mdatoms->chargeA, x, fr->f_novirsum,
1370 inputrec->ex, inputrec->et, t);
1373 /* If we have NoVirSum forces, but we do not calculate the virial,
1374 * we sum fr->f_novirum=f later.
1376 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1378 wallcycle_start(wcycle, ewcVSITESPREAD);
1379 spread_vsite_f(fplog, vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1380 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1381 wallcycle_stop(wcycle, ewcVSITESPREAD);
1385 wallcycle_start(wcycle, ewcVSITESPREAD);
1386 spread_vsite_f(fplog, vsite, x, fr->f_twin, NULL, FALSE, NULL,
1388 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1389 wallcycle_stop(wcycle, ewcVSITESPREAD);
1393 if (flags & GMX_FORCE_VIRIAL)
1395 /* Calculation of the virial must be done after vsites! */
1396 calc_virial(fplog, mdatoms->start, mdatoms->homenr, x, f,
1397 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1401 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1403 pull_potential_wrapper(fplog, bSepDVDL, cr, inputrec, box, x,
1404 f, vir_force, mdatoms, enerd, lambda, t);
1407 if (PAR(cr) && !(cr->duty & DUTY_PME))
1409 /* In case of node-splitting, the PP nodes receive the long-range
1410 * forces, virial and energy from the PME nodes here.
1412 pme_receive_force_ener(fplog, bSepDVDL, cr, wcycle, enerd, fr);
1417 post_process_forces(fplog, cr, step, nrnb, wcycle,
1418 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
1422 /* Sum the potential energy terms from group contributions */
1423 sum_epot(&(inputrec->opts), &(enerd->grpp), enerd->term);
1426 void do_force_cutsGROUP(FILE *fplog, t_commrec *cr,
1427 t_inputrec *inputrec,
1428 gmx_large_int_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1429 gmx_localtop_t *top,
1431 gmx_groups_t *groups,
1432 matrix box, rvec x[], history_t *hist,
1436 gmx_enerdata_t *enerd, t_fcdata *fcd,
1437 real *lambda, t_graph *graph,
1438 t_forcerec *fr, gmx_vsite_t *vsite, rvec mu_tot,
1439 double t, FILE *field, gmx_edsam_t ed,
1440 gmx_bool bBornRadii,
1446 gmx_bool bSepDVDL, bStateChanged, bNS, bFillGrid, bCalcCGCM, bBS;
1447 gmx_bool bDoLongRangeNS, bDoForces, bDoPotential, bSepLRF;
1448 gmx_bool bDoAdressWF;
1450 rvec vzero, box_diag;
1451 real e, v, dvdlambda[efptNR];
1453 float cycles_pme, cycles_force;
1455 start = mdatoms->start;
1456 homenr = mdatoms->homenr;
1458 bSepDVDL = (fr->bSepDVDL && do_per_step(step, inputrec->nstlog));
1460 clear_mat(vir_force);
1464 pd_cg_range(cr, &cg0, &cg1);
1469 if (DOMAINDECOMP(cr))
1471 cg1 = cr->dd->ncg_tot;
1483 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
1484 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
1485 /* Should we update the long-range neighborlists at this step? */
1486 bDoLongRangeNS = fr->bTwinRange && bNS;
1487 /* Should we perform the long-range nonbonded evaluation inside the neighborsearching? */
1488 bFillGrid = (bNS && bStateChanged);
1489 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
1490 bDoForces = (flags & GMX_FORCE_FORCES);
1491 bDoPotential = (flags & GMX_FORCE_ENERGY);
1492 bSepLRF = ((inputrec->nstcalclr > 1) && bDoForces &&
1493 (flags & GMX_FORCE_SEPLRF) && (flags & GMX_FORCE_DO_LR));
1495 /* should probably move this to the forcerec since it doesn't change */
1496 bDoAdressWF = ((fr->adress_type != eAdressOff));
1500 update_forcerec(fplog, fr, box);
1502 if (NEED_MUTOT(*inputrec))
1504 /* Calculate total (local) dipole moment in a temporary common array.
1505 * This makes it possible to sum them over nodes faster.
1507 calc_mu(start, homenr,
1508 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
1513 if (fr->ePBC != epbcNONE)
1515 /* Compute shift vectors every step,
1516 * because of pressure coupling or box deformation!
1518 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
1520 calc_shifts(box, fr->shift_vec);
1525 put_charge_groups_in_box(fplog, cg0, cg1, fr->ePBC, box,
1526 &(top->cgs), x, fr->cg_cm);
1527 inc_nrnb(nrnb, eNR_CGCM, homenr);
1528 inc_nrnb(nrnb, eNR_RESETX, cg1-cg0);
1530 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
1532 unshift_self(graph, box, x);
1537 calc_cgcm(fplog, cg0, cg1, &(top->cgs), x, fr->cg_cm);
1538 inc_nrnb(nrnb, eNR_CGCM, homenr);
1545 move_cgcm(fplog, cr, fr->cg_cm);
1549 pr_rvecs(debug, 0, "cgcm", fr->cg_cm, top->cgs.nr);
1554 if (!(cr->duty & DUTY_PME))
1556 /* Send particle coordinates to the pme nodes.
1557 * Since this is only implemented for domain decomposition
1558 * and domain decomposition does not use the graph,
1559 * we do not need to worry about shifting.
1562 wallcycle_start(wcycle, ewcPP_PMESENDX);
1564 bBS = (inputrec->nwall == 2);
1567 copy_mat(box, boxs);
1568 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
1571 gmx_pme_send_x(cr, bBS ? boxs : box, x,
1572 mdatoms->nChargePerturbed, lambda[efptCOUL],
1573 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)), step);
1575 wallcycle_stop(wcycle, ewcPP_PMESENDX);
1577 #endif /* GMX_MPI */
1579 /* Communicate coordinates and sum dipole if necessary */
1582 wallcycle_start(wcycle, ewcMOVEX);
1583 if (DOMAINDECOMP(cr))
1585 dd_move_x(cr->dd, box, x);
1589 move_x(fplog, cr, GMX_LEFT, GMX_RIGHT, x, nrnb);
1591 wallcycle_stop(wcycle, ewcMOVEX);
1594 /* update adress weight beforehand */
1595 if (bStateChanged && bDoAdressWF)
1597 /* need pbc for adress weight calculation with pbc_dx */
1598 set_pbc(&pbc, inputrec->ePBC, box);
1599 if (fr->adress_site == eAdressSITEcog)
1601 update_adress_weights_cog(top->idef.iparams, top->idef.il, x, fr, mdatoms,
1602 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1604 else if (fr->adress_site == eAdressSITEcom)
1606 update_adress_weights_com(fplog, cg0, cg1, &(top->cgs), x, fr, mdatoms,
1607 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1609 else if (fr->adress_site == eAdressSITEatomatom)
1611 update_adress_weights_atom_per_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1612 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1616 update_adress_weights_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1617 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1621 if (NEED_MUTOT(*inputrec))
1628 gmx_sumd(2*DIM, mu, cr);
1630 for (i = 0; i < 2; i++)
1632 for (j = 0; j < DIM; j++)
1634 fr->mu_tot[i][j] = mu[i*DIM + j];
1638 if (fr->efep == efepNO)
1640 copy_rvec(fr->mu_tot[0], mu_tot);
1644 for (j = 0; j < DIM; j++)
1647 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] + lambda[efptCOUL]*fr->mu_tot[1][j];
1652 /* Reset energies */
1653 reset_enerdata(&(inputrec->opts), fr, bNS, enerd, MASTER(cr));
1654 clear_rvecs(SHIFTS, fr->fshift);
1658 wallcycle_start(wcycle, ewcNS);
1660 if (graph && bStateChanged)
1662 /* Calculate intramolecular shift vectors to make molecules whole */
1663 mk_mshift(fplog, graph, fr->ePBC, box, x);
1666 /* Do the actual neighbour searching and if twin range electrostatics
1667 * also do the calculation of long range forces and energies.
1669 for (i = 0; i < efptNR; i++)
1673 ns(fplog, fr, x, box,
1674 groups, &(inputrec->opts), top, mdatoms,
1675 cr, nrnb, lambda, dvdlambda, &enerd->grpp, bFillGrid,
1679 fprintf(fplog, sepdvdlformat, "LR non-bonded", 0.0, dvdlambda);
1681 enerd->dvdl_lin[efptVDW] += dvdlambda[efptVDW];
1682 enerd->dvdl_lin[efptCOUL] += dvdlambda[efptCOUL];
1684 wallcycle_stop(wcycle, ewcNS);
1687 if (inputrec->implicit_solvent && bNS)
1689 make_gb_nblist(cr, inputrec->gb_algorithm, inputrec->rlist,
1690 x, box, fr, &top->idef, graph, fr->born);
1693 if (DOMAINDECOMP(cr))
1695 if (!(cr->duty & DUTY_PME))
1697 wallcycle_start(wcycle, ewcPPDURINGPME);
1698 dd_force_flop_start(cr->dd, nrnb);
1704 /* Enforced rotation has its own cycle counter that starts after the collective
1705 * coordinates have been communicated. It is added to ddCyclF to allow
1706 * for proper load-balancing */
1707 wallcycle_start(wcycle, ewcROT);
1708 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1709 wallcycle_stop(wcycle, ewcROT);
1712 /* Start the force cycle counter.
1713 * This counter is stopped in do_forcelow_level.
1714 * No parallel communication should occur while this counter is running,
1715 * since that will interfere with the dynamic load balancing.
1717 wallcycle_start(wcycle, ewcFORCE);
1721 /* Reset forces for which the virial is calculated separately:
1722 * PME/Ewald forces if necessary */
1723 if (fr->bF_NoVirSum)
1725 if (flags & GMX_FORCE_VIRIAL)
1727 fr->f_novirsum = fr->f_novirsum_alloc;
1730 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1734 clear_rvecs(homenr, fr->f_novirsum+start);
1739 /* We are not calculating the pressure so we do not need
1740 * a separate array for forces that do not contribute
1747 /* Clear the short- and long-range forces */
1748 clear_rvecs(fr->natoms_force_constr, f);
1749 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1751 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1754 clear_rvec(fr->vir_diag_posres);
1756 if (inputrec->ePull == epullCONSTRAINT)
1758 clear_pull_forces(inputrec->pull);
1761 /* update QMMMrec, if necessary */
1764 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1767 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1769 posres_wrapper(fplog, flags, bSepDVDL, inputrec, nrnb, top, box, x,
1770 f, enerd, lambda, fr);
1773 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_FBPOSRES].nr > 0)
1775 /* Flat-bottomed position restraints always require full pbc */
1776 if (!(bStateChanged && bDoAdressWF))
1778 set_pbc(&pbc, inputrec->ePBC, box);
1780 v = fbposres(top->idef.il[F_FBPOSRES].nr, top->idef.il[F_FBPOSRES].iatoms,
1781 top->idef.iparams_fbposres,
1782 (const rvec*)x, fr->f_novirsum, fr->vir_diag_posres,
1783 inputrec->ePBC == epbcNONE ? NULL : &pbc,
1784 fr->rc_scaling, fr->ePBC, fr->posres_com);
1785 enerd->term[F_FBPOSRES] += v;
1786 inc_nrnb(nrnb, eNR_FBPOSRES, top->idef.il[F_FBPOSRES].nr/2);
1789 /* Compute the bonded and non-bonded energies and optionally forces */
1790 do_force_lowlevel(fplog, step, fr, inputrec, &(top->idef),
1791 cr, nrnb, wcycle, mdatoms, &(inputrec->opts),
1792 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, mtop, top, fr->born,
1793 &(top->atomtypes), bBornRadii, box,
1794 inputrec->fepvals, lambda,
1795 graph, &(top->excls), fr->mu_tot,
1801 if (do_per_step(step, inputrec->nstcalclr))
1803 /* Add the long range forces to the short range forces */
1804 for (i = 0; i < fr->natoms_force_constr; i++)
1806 rvec_add(fr->f_twin[i], f[i], f[i]);
1811 cycles_force = wallcycle_stop(wcycle, ewcFORCE);
1815 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1818 if (DOMAINDECOMP(cr))
1820 dd_force_flop_stop(cr->dd, nrnb);
1823 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1829 if (IR_ELEC_FIELD(*inputrec))
1831 /* Compute forces due to electric field */
1832 calc_f_el(MASTER(cr) ? field : NULL,
1833 start, homenr, mdatoms->chargeA, x, fr->f_novirsum,
1834 inputrec->ex, inputrec->et, t);
1837 if (bDoAdressWF && fr->adress_icor == eAdressICThermoForce)
1839 /* Compute thermodynamic force in hybrid AdResS region */
1840 adress_thermo_force(start, homenr, &(top->cgs), x, fr->f_novirsum, fr, mdatoms,
1841 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1844 /* Communicate the forces */
1847 wallcycle_start(wcycle, ewcMOVEF);
1848 if (DOMAINDECOMP(cr))
1850 dd_move_f(cr->dd, f, fr->fshift);
1851 /* Do we need to communicate the separate force array
1852 * for terms that do not contribute to the single sum virial?
1853 * Position restraints and electric fields do not introduce
1854 * inter-cg forces, only full electrostatics methods do.
1855 * When we do not calculate the virial, fr->f_novirsum = f,
1856 * so we have already communicated these forces.
1858 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1859 (flags & GMX_FORCE_VIRIAL))
1861 dd_move_f(cr->dd, fr->f_novirsum, NULL);
1865 /* We should not update the shift forces here,
1866 * since f_twin is already included in f.
1868 dd_move_f(cr->dd, fr->f_twin, NULL);
1873 pd_move_f(cr, f, nrnb);
1876 pd_move_f(cr, fr->f_twin, nrnb);
1879 wallcycle_stop(wcycle, ewcMOVEF);
1882 /* If we have NoVirSum forces, but we do not calculate the virial,
1883 * we sum fr->f_novirum=f later.
1885 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1887 wallcycle_start(wcycle, ewcVSITESPREAD);
1888 spread_vsite_f(fplog, vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1889 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1890 wallcycle_stop(wcycle, ewcVSITESPREAD);
1894 wallcycle_start(wcycle, ewcVSITESPREAD);
1895 spread_vsite_f(fplog, vsite, x, fr->f_twin, NULL, FALSE, NULL,
1897 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1898 wallcycle_stop(wcycle, ewcVSITESPREAD);
1902 if (flags & GMX_FORCE_VIRIAL)
1904 /* Calculation of the virial must be done after vsites! */
1905 calc_virial(fplog, mdatoms->start, mdatoms->homenr, x, f,
1906 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1910 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1912 pull_potential_wrapper(fplog, bSepDVDL, cr, inputrec, box, x,
1913 f, vir_force, mdatoms, enerd, lambda, t);
1916 /* Add the forces from enforced rotation potentials (if any) */
1919 wallcycle_start(wcycle, ewcROTadd);
1920 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
1921 wallcycle_stop(wcycle, ewcROTadd);
1924 if (PAR(cr) && !(cr->duty & DUTY_PME))
1926 /* In case of node-splitting, the PP nodes receive the long-range
1927 * forces, virial and energy from the PME nodes here.
1929 pme_receive_force_ener(fplog, bSepDVDL, cr, wcycle, enerd, fr);
1934 post_process_forces(fplog, cr, step, nrnb, wcycle,
1935 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
1939 /* Sum the potential energy terms from group contributions */
1940 sum_epot(&(inputrec->opts), &(enerd->grpp), enerd->term);
1943 void do_force(FILE *fplog, t_commrec *cr,
1944 t_inputrec *inputrec,
1945 gmx_large_int_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1946 gmx_localtop_t *top,
1948 gmx_groups_t *groups,
1949 matrix box, rvec x[], history_t *hist,
1953 gmx_enerdata_t *enerd, t_fcdata *fcd,
1954 real *lambda, t_graph *graph,
1956 gmx_vsite_t *vsite, rvec mu_tot,
1957 double t, FILE *field, gmx_edsam_t ed,
1958 gmx_bool bBornRadii,
1961 /* modify force flag if not doing nonbonded */
1962 if (!fr->bNonbonded)
1964 flags &= ~GMX_FORCE_NONBONDED;
1967 switch (inputrec->cutoff_scheme)
1970 do_force_cutsVERLET(fplog, cr, inputrec,
1986 do_force_cutsGROUP(fplog, cr, inputrec,
2001 gmx_incons("Invalid cut-off scheme passed!");
2006 void do_constrain_first(FILE *fplog, gmx_constr_t constr,
2007 t_inputrec *ir, t_mdatoms *md,
2008 t_state *state, rvec *f,
2009 t_graph *graph, t_commrec *cr, t_nrnb *nrnb,
2010 t_forcerec *fr, gmx_localtop_t *top, tensor shake_vir)
2012 int i, m, start, end;
2013 gmx_large_int_t step;
2014 real dt = ir->delta_t;
2018 snew(savex, state->natoms);
2021 end = md->homenr + start;
2025 fprintf(debug, "vcm: start=%d, homenr=%d, end=%d\n",
2026 start, md->homenr, end);
2028 /* Do a first constrain to reset particles... */
2029 step = ir->init_step;
2032 char buf[STEPSTRSIZE];
2033 fprintf(fplog, "\nConstraining the starting coordinates (step %s)\n",
2034 gmx_step_str(step, buf));
2038 /* constrain the current position */
2039 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2040 ir, NULL, cr, step, 0, md,
2041 state->x, state->x, NULL,
2042 fr->bMolPBC, state->box,
2043 state->lambda[efptBONDED], &dvdl_dum,
2044 NULL, NULL, nrnb, econqCoord,
2045 ir->epc == epcMTTK, state->veta, state->veta);
2048 /* constrain the inital velocity, and save it */
2049 /* also may be useful if we need the ekin from the halfstep for velocity verlet */
2050 /* might not yet treat veta correctly */
2051 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2052 ir, NULL, cr, step, 0, md,
2053 state->x, state->v, state->v,
2054 fr->bMolPBC, state->box,
2055 state->lambda[efptBONDED], &dvdl_dum,
2056 NULL, NULL, nrnb, econqVeloc,
2057 ir->epc == epcMTTK, state->veta, state->veta);
2059 /* constrain the inital velocities at t-dt/2 */
2060 if (EI_STATE_VELOCITY(ir->eI) && ir->eI != eiVV)
2062 for (i = start; (i < end); i++)
2064 for (m = 0; (m < DIM); m++)
2066 /* Reverse the velocity */
2067 state->v[i][m] = -state->v[i][m];
2068 /* Store the position at t-dt in buf */
2069 savex[i][m] = state->x[i][m] + dt*state->v[i][m];
2072 /* Shake the positions at t=-dt with the positions at t=0
2073 * as reference coordinates.
2077 char buf[STEPSTRSIZE];
2078 fprintf(fplog, "\nConstraining the coordinates at t0-dt (step %s)\n",
2079 gmx_step_str(step, buf));
2082 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2083 ir, NULL, cr, step, -1, md,
2084 state->x, savex, NULL,
2085 fr->bMolPBC, state->box,
2086 state->lambda[efptBONDED], &dvdl_dum,
2087 state->v, NULL, nrnb, econqCoord,
2088 ir->epc == epcMTTK, state->veta, state->veta);
2090 for (i = start; i < end; i++)
2092 for (m = 0; m < DIM; m++)
2094 /* Re-reverse the velocities */
2095 state->v[i][m] = -state->v[i][m];
2102 void calc_enervirdiff(FILE *fplog, int eDispCorr, t_forcerec *fr)
2104 double eners[2], virs[2], enersum, virsum, y0, f, g, h;
2105 double r0, r1, r, rc3, rc9, ea, eb, ec, pa, pb, pc, pd;
2106 double invscale, invscale2, invscale3;
2107 int ri0, ri1, ri, i, offstart, offset;
2108 real scale, *vdwtab, tabfactor, tmp;
2110 fr->enershiftsix = 0;
2111 fr->enershifttwelve = 0;
2112 fr->enerdiffsix = 0;
2113 fr->enerdifftwelve = 0;
2115 fr->virdifftwelve = 0;
2117 if (eDispCorr != edispcNO)
2119 for (i = 0; i < 2; i++)
2124 if ((fr->vdwtype == evdwSWITCH) || (fr->vdwtype == evdwSHIFT))
2126 if (fr->rvdw_switch == 0)
2129 "With dispersion correction rvdw-switch can not be zero "
2130 "for vdw-type = %s", evdw_names[fr->vdwtype]);
2133 scale = fr->nblists[0].table_elec_vdw.scale;
2134 vdwtab = fr->nblists[0].table_vdw.data;
2136 /* Round the cut-offs to exact table values for precision */
2137 ri0 = floor(fr->rvdw_switch*scale);
2138 ri1 = ceil(fr->rvdw*scale);
2144 if (fr->vdwtype == evdwSHIFT)
2146 /* Determine the constant energy shift below rvdw_switch.
2147 * Table has a scale factor since we have scaled it down to compensate
2148 * for scaling-up c6/c12 with the derivative factors to save flops in analytical kernels.
2150 fr->enershiftsix = (real)(-1.0/(rc3*rc3)) - 6.0*vdwtab[8*ri0];
2151 fr->enershifttwelve = (real)( 1.0/(rc9*rc3)) - 12.0*vdwtab[8*ri0 + 4];
2153 /* Add the constant part from 0 to rvdw_switch.
2154 * This integration from 0 to rvdw_switch overcounts the number
2155 * of interactions by 1, as it also counts the self interaction.
2156 * We will correct for this later.
2158 eners[0] += 4.0*M_PI*fr->enershiftsix*rc3/3.0;
2159 eners[1] += 4.0*M_PI*fr->enershifttwelve*rc3/3.0;
2161 invscale = 1.0/(scale);
2162 invscale2 = invscale*invscale;
2163 invscale3 = invscale*invscale2;
2165 /* following summation derived from cubic spline definition,
2166 Numerical Recipies in C, second edition, p. 113-116. Exact
2167 for the cubic spline. We first calculate the negative of
2168 the energy from rvdw to rvdw_switch, assuming that g(r)=1,
2169 and then add the more standard, abrupt cutoff correction to
2170 that result, yielding the long-range correction for a
2171 switched function. We perform both the pressure and energy
2172 loops at the same time for simplicity, as the computational
2175 for (i = 0; i < 2; i++)
2177 enersum = 0.0; virsum = 0.0;
2181 /* Since the dispersion table has been scaled down a factor 6.0 and the repulsion
2182 * a factor 12.0 to compensate for the c6/c12 parameters inside nbfp[] being scaled
2183 * up (to save flops in kernels), we need to correct for this.
2192 for (ri = ri0; ri < ri1; ri++)
2196 eb = 2.0*invscale2*r;
2200 pb = 3.0*invscale2*r;
2201 pc = 3.0*invscale*r*r;
2204 /* this "8" is from the packing in the vdwtab array - perhaps should be #define'ed? */
2205 offset = 8*ri + offstart;
2206 y0 = vdwtab[offset];
2207 f = vdwtab[offset+1];
2208 g = vdwtab[offset+2];
2209 h = vdwtab[offset+3];
2211 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);
2212 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);
2215 enersum *= 4.0*M_PI*tabfactor;
2216 virsum *= 4.0*M_PI*tabfactor;
2217 eners[i] -= enersum;
2221 /* now add the correction for rvdw_switch to infinity */
2222 eners[0] += -4.0*M_PI/(3.0*rc3);
2223 eners[1] += 4.0*M_PI/(9.0*rc9);
2224 virs[0] += 8.0*M_PI/rc3;
2225 virs[1] += -16.0*M_PI/(3.0*rc9);
2227 else if ((fr->vdwtype == evdwCUT) || (fr->vdwtype == evdwUSER))
2229 if (fr->vdwtype == evdwUSER && fplog)
2232 "WARNING: using dispersion correction with user tables\n");
2234 rc3 = fr->rvdw*fr->rvdw*fr->rvdw;
2236 /* Contribution beyond the cut-off */
2237 eners[0] += -4.0*M_PI/(3.0*rc3);
2238 eners[1] += 4.0*M_PI/(9.0*rc9);
2239 if (fr->vdw_modifier == eintmodPOTSHIFT)
2241 /* Contribution within the cut-off */
2242 eners[0] += -4.0*M_PI/(3.0*rc3);
2243 eners[1] += 4.0*M_PI/(3.0*rc9);
2245 /* Contribution beyond the cut-off */
2246 virs[0] += 8.0*M_PI/rc3;
2247 virs[1] += -16.0*M_PI/(3.0*rc9);
2252 "Dispersion correction is not implemented for vdw-type = %s",
2253 evdw_names[fr->vdwtype]);
2255 fr->enerdiffsix = eners[0];
2256 fr->enerdifftwelve = eners[1];
2257 /* The 0.5 is due to the Gromacs definition of the virial */
2258 fr->virdiffsix = 0.5*virs[0];
2259 fr->virdifftwelve = 0.5*virs[1];
2263 void calc_dispcorr(FILE *fplog, t_inputrec *ir, t_forcerec *fr,
2264 gmx_large_int_t step, int natoms,
2265 matrix box, real lambda, tensor pres, tensor virial,
2266 real *prescorr, real *enercorr, real *dvdlcorr)
2268 gmx_bool bCorrAll, bCorrPres;
2269 real dvdlambda, invvol, dens, ninter, avcsix, avctwelve, enerdiff, svir = 0, spres = 0;
2279 if (ir->eDispCorr != edispcNO)
2281 bCorrAll = (ir->eDispCorr == edispcAllEner ||
2282 ir->eDispCorr == edispcAllEnerPres);
2283 bCorrPres = (ir->eDispCorr == edispcEnerPres ||
2284 ir->eDispCorr == edispcAllEnerPres);
2286 invvol = 1/det(box);
2289 /* Only correct for the interactions with the inserted molecule */
2290 dens = (natoms - fr->n_tpi)*invvol;
2295 dens = natoms*invvol;
2296 ninter = 0.5*natoms;
2299 if (ir->efep == efepNO)
2301 avcsix = fr->avcsix[0];
2302 avctwelve = fr->avctwelve[0];
2306 avcsix = (1 - lambda)*fr->avcsix[0] + lambda*fr->avcsix[1];
2307 avctwelve = (1 - lambda)*fr->avctwelve[0] + lambda*fr->avctwelve[1];
2310 enerdiff = ninter*(dens*fr->enerdiffsix - fr->enershiftsix);
2311 *enercorr += avcsix*enerdiff;
2313 if (ir->efep != efepNO)
2315 dvdlambda += (fr->avcsix[1] - fr->avcsix[0])*enerdiff;
2319 enerdiff = ninter*(dens*fr->enerdifftwelve - fr->enershifttwelve);
2320 *enercorr += avctwelve*enerdiff;
2321 if (fr->efep != efepNO)
2323 dvdlambda += (fr->avctwelve[1] - fr->avctwelve[0])*enerdiff;
2329 svir = ninter*dens*avcsix*fr->virdiffsix/3.0;
2330 if (ir->eDispCorr == edispcAllEnerPres)
2332 svir += ninter*dens*avctwelve*fr->virdifftwelve/3.0;
2334 /* The factor 2 is because of the Gromacs virial definition */
2335 spres = -2.0*invvol*svir*PRESFAC;
2337 for (m = 0; m < DIM; m++)
2339 virial[m][m] += svir;
2340 pres[m][m] += spres;
2345 /* Can't currently control when it prints, for now, just print when degugging */
2350 fprintf(debug, "Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
2356 "Long Range LJ corr.: Epot %10g, Pres: %10g, Vir: %10g\n",
2357 *enercorr, spres, svir);
2361 fprintf(debug, "Long Range LJ corr.: Epot %10g\n", *enercorr);
2365 if (fr->bSepDVDL && do_per_step(step, ir->nstlog))
2367 fprintf(fplog, sepdvdlformat, "Dispersion correction",
2368 *enercorr, dvdlambda);
2370 if (fr->efep != efepNO)
2372 *dvdlcorr += dvdlambda;
2377 void do_pbc_first(FILE *fplog, matrix box, t_forcerec *fr,
2378 t_graph *graph, rvec x[])
2382 fprintf(fplog, "Removing pbc first time\n");
2384 calc_shifts(box, fr->shift_vec);
2387 mk_mshift(fplog, graph, fr->ePBC, box, x);
2390 p_graph(debug, "do_pbc_first 1", graph);
2392 shift_self(graph, box, x);
2393 /* By doing an extra mk_mshift the molecules that are broken
2394 * because they were e.g. imported from another software
2395 * will be made whole again. Such are the healing powers
2398 mk_mshift(fplog, graph, fr->ePBC, box, x);
2401 p_graph(debug, "do_pbc_first 2", graph);
2406 fprintf(fplog, "Done rmpbc\n");
2410 static void low_do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2411 gmx_mtop_t *mtop, rvec x[],
2416 gmx_molblock_t *molb;
2418 if (bFirst && fplog)
2420 fprintf(fplog, "Removing pbc first time\n");
2425 for (mb = 0; mb < mtop->nmolblock; mb++)
2427 molb = &mtop->molblock[mb];
2428 if (molb->natoms_mol == 1 ||
2429 (!bFirst && mtop->moltype[molb->type].cgs.nr == 1))
2431 /* Just one atom or charge group in the molecule, no PBC required */
2432 as += molb->nmol*molb->natoms_mol;
2436 /* Pass NULL iso fplog to avoid graph prints for each molecule type */
2437 mk_graph_ilist(NULL, mtop->moltype[molb->type].ilist,
2438 0, molb->natoms_mol, FALSE, FALSE, graph);
2440 for (mol = 0; mol < molb->nmol; mol++)
2442 mk_mshift(fplog, graph, ePBC, box, x+as);
2444 shift_self(graph, box, x+as);
2445 /* The molecule is whole now.
2446 * We don't need the second mk_mshift call as in do_pbc_first,
2447 * since we no longer need this graph.
2450 as += molb->natoms_mol;
2458 void do_pbc_first_mtop(FILE *fplog, int ePBC, matrix box,
2459 gmx_mtop_t *mtop, rvec x[])
2461 low_do_pbc_mtop(fplog, ePBC, box, mtop, x, TRUE);
2464 void do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2465 gmx_mtop_t *mtop, rvec x[])
2467 low_do_pbc_mtop(fplog, ePBC, box, mtop, x, FALSE);
2470 void finish_run(FILE *fplog, t_commrec *cr, const char *confout,
2471 t_inputrec *inputrec,
2472 t_nrnb nrnb[], gmx_wallcycle_t wcycle,
2473 gmx_runtime_t *runtime,
2474 wallclock_gpu_t *gputimes,
2476 gmx_bool bWriteStat)
2479 t_nrnb *nrnb_tot = NULL;
2483 wallcycle_sum(cr, wcycle);
2489 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM,
2490 cr->mpi_comm_mysim);
2498 #if defined(GMX_MPI) && !defined(GMX_THREAD_MPI)
2501 /* reduce nodetime over all MPI processes in the current simulation */
2503 MPI_Allreduce(&runtime->proctime, &sum, 1, MPI_DOUBLE, MPI_SUM,
2504 cr->mpi_comm_mysim);
2505 runtime->proctime = sum;
2511 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
2518 if ((cr->duty & DUTY_PP) && DOMAINDECOMP(cr))
2520 print_dd_statistics(cr, inputrec, fplog);
2532 snew(nrnb_all, cr->nnodes);
2533 nrnb_all[0] = *nrnb;
2534 for (s = 1; s < cr->nnodes; s++)
2536 MPI_Recv(nrnb_all[s].n, eNRNB, MPI_DOUBLE, s, 0,
2537 cr->mpi_comm_mysim, &stat);
2539 pr_load(fplog, cr, nrnb_all);
2544 MPI_Send(nrnb->n, eNRNB, MPI_DOUBLE, MASTERRANK(cr), 0,
2545 cr->mpi_comm_mysim);
2552 wallcycle_print(fplog, cr->nnodes, cr->npmenodes, runtime->realtime,
2555 if (EI_DYNAMICS(inputrec->eI))
2557 delta_t = inputrec->delta_t;
2566 print_perf(fplog, runtime->proctime, runtime->realtime,
2567 cr->nnodes-cr->npmenodes,
2568 runtime->nsteps_done, delta_t, nbfs, mflop,
2573 print_perf(stderr, runtime->proctime, runtime->realtime,
2574 cr->nnodes-cr->npmenodes,
2575 runtime->nsteps_done, delta_t, nbfs, mflop,
2581 extern void initialize_lambdas(FILE *fplog, t_inputrec *ir, int *fep_state, real *lambda, double *lam0)
2583 /* this function works, but could probably use a logic rewrite to keep all the different
2584 types of efep straight. */
2587 t_lambda *fep = ir->fepvals;
2589 if ((ir->efep == efepNO) && (ir->bSimTemp == FALSE))
2591 for (i = 0; i < efptNR; i++)
2603 *fep_state = fep->init_fep_state; /* this might overwrite the checkpoint
2604 if checkpoint is set -- a kludge is in for now
2606 for (i = 0; i < efptNR; i++)
2608 /* overwrite lambda state with init_lambda for now for backwards compatibility */
2609 if (fep->init_lambda >= 0) /* if it's -1, it was never initializd */
2611 lambda[i] = fep->init_lambda;
2614 lam0[i] = lambda[i];
2619 lambda[i] = fep->all_lambda[i][*fep_state];
2622 lam0[i] = lambda[i];
2628 /* need to rescale control temperatures to match current state */
2629 for (i = 0; i < ir->opts.ngtc; i++)
2631 if (ir->opts.ref_t[i] > 0)
2633 ir->opts.ref_t[i] = ir->simtempvals->temperatures[*fep_state];
2639 /* Send to the log the information on the current lambdas */
2642 fprintf(fplog, "Initial vector of lambda components:[ ");
2643 for (i = 0; i < efptNR; i++)
2645 fprintf(fplog, "%10.4f ", lambda[i]);
2647 fprintf(fplog, "]\n");
2653 void init_md(FILE *fplog,
2654 t_commrec *cr, t_inputrec *ir, const output_env_t oenv,
2655 double *t, double *t0,
2656 real *lambda, int *fep_state, double *lam0,
2657 t_nrnb *nrnb, gmx_mtop_t *mtop,
2659 int nfile, const t_filenm fnm[],
2660 gmx_mdoutf_t **outf, t_mdebin **mdebin,
2661 tensor force_vir, tensor shake_vir, rvec mu_tot,
2662 gmx_bool *bSimAnn, t_vcm **vcm, t_state *state, unsigned long Flags)
2667 /* Initial values */
2668 *t = *t0 = ir->init_t;
2671 for (i = 0; i < ir->opts.ngtc; i++)
2673 /* set bSimAnn if any group is being annealed */
2674 if (ir->opts.annealing[i] != eannNO)
2681 update_annealing_target_temp(&(ir->opts), ir->init_t);
2684 /* Initialize lambda variables */
2685 initialize_lambdas(fplog, ir, fep_state, lambda, lam0);
2689 *upd = init_update(fplog, ir);
2695 *vcm = init_vcm(fplog, &mtop->groups, ir);
2698 if (EI_DYNAMICS(ir->eI) && !(Flags & MD_APPENDFILES))
2700 if (ir->etc == etcBERENDSEN)
2702 please_cite(fplog, "Berendsen84a");
2704 if (ir->etc == etcVRESCALE)
2706 please_cite(fplog, "Bussi2007a");
2714 *outf = init_mdoutf(nfile, fnm, Flags, cr, ir, oenv);
2716 *mdebin = init_mdebin((Flags & MD_APPENDFILES) ? NULL : (*outf)->fp_ene,
2717 mtop, ir, (*outf)->fp_dhdl);
2722 please_cite(fplog, "Fritsch12");
2723 please_cite(fplog, "Junghans10");
2725 /* Initiate variables */
2726 clear_mat(force_vir);
2727 clear_mat(shake_vir);