2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
42 #ifdef HAVE_SYS_TIME_H
53 #include "chargegroup.h"
73 #include "pull_rotation.h"
74 #include "gmx_random.h"
78 #include "nbnxn_atomdata.h"
79 #include "nbnxn_search.h"
80 #include "nbnxn_kernels/nbnxn_kernel_ref.h"
81 #include "nbnxn_kernels/simd_4xn/nbnxn_kernel_simd_4xn.h"
82 #include "nbnxn_kernels/simd_2xnn/nbnxn_kernel_simd_2xnn.h"
83 #include "nbnxn_kernels/nbnxn_kernel_gpu_ref.h"
85 #include "gromacs/timing/wallcycle.h"
86 #include "gromacs/timing/walltime_accounting.h"
87 #include "gromacs/utility/gmxmpi.h"
92 #include "nbnxn_cuda_data_mgmt.h"
93 #include "nbnxn_cuda/nbnxn_cuda.h"
95 void print_time(FILE *out,
96 gmx_walltime_accounting_t walltime_accounting,
99 t_commrec gmx_unused *cr)
102 char timebuf[STRLEN];
103 double dt, elapsed_seconds, time_per_step;
106 #ifndef GMX_THREAD_MPI
112 fprintf(out, "step %s", gmx_step_str(step, buf));
113 if ((step >= ir->nstlist))
115 double seconds_since_epoch = gmx_gettime();
116 elapsed_seconds = seconds_since_epoch - walltime_accounting_get_start_time_stamp(walltime_accounting);
117 time_per_step = elapsed_seconds/(step - ir->init_step + 1);
118 dt = (ir->nsteps + ir->init_step - step) * time_per_step;
124 finish = (time_t) (seconds_since_epoch + dt);
125 gmx_ctime_r(&finish, timebuf, STRLEN);
126 sprintf(buf, "%s", timebuf);
127 buf[strlen(buf)-1] = '\0';
128 fprintf(out, ", will finish %s", buf);
132 fprintf(out, ", remaining wall clock time: %5d s ", (int)dt);
137 fprintf(out, " performance: %.1f ns/day ",
138 ir->delta_t/1000*24*60*60/time_per_step);
141 #ifndef GMX_THREAD_MPI
151 void print_date_and_time(FILE *fplog, int nodeid, const char *title,
152 const gmx_walltime_accounting_t walltime_accounting)
155 char timebuf[STRLEN];
156 char time_string[STRLEN];
161 if (walltime_accounting != NULL)
163 tmptime = (time_t) walltime_accounting_get_start_time_stamp(walltime_accounting);
164 gmx_ctime_r(&tmptime, timebuf, STRLEN);
168 tmptime = (time_t) gmx_gettime();
169 gmx_ctime_r(&tmptime, timebuf, STRLEN);
171 for (i = 0; timebuf[i] >= ' '; i++)
173 time_string[i] = timebuf[i];
175 time_string[i] = '\0';
177 fprintf(fplog, "%s on node %d %s\n", title, nodeid, time_string);
181 static void sum_forces(int start, int end, rvec f[], rvec flr[])
187 pr_rvecs(debug, 0, "fsr", f+start, end-start);
188 pr_rvecs(debug, 0, "flr", flr+start, end-start);
190 for (i = start; (i < end); i++)
192 rvec_inc(f[i], flr[i]);
197 * calc_f_el calculates forces due to an electric field.
199 * force is kJ mol^-1 nm^-1 = e * kJ mol^-1 nm^-1 / e
201 * Et[] contains the parameters for the time dependent
202 * part of the field (not yet used).
203 * Ex[] contains the parameters for
204 * the spatial dependent part of the field. You can have cool periodic
205 * fields in principle, but only a constant field is supported
207 * The function should return the energy due to the electric field
208 * (if any) but for now returns 0.
211 * There can be problems with the virial.
212 * Since the field is not self-consistent this is unavoidable.
213 * For neutral molecules the virial is correct within this approximation.
214 * For neutral systems with many charged molecules the error is small.
215 * But for systems with a net charge or a few charged molecules
216 * the error can be significant when the field is high.
217 * Solution: implement a self-consitent electric field into PME.
219 static void calc_f_el(FILE *fp, int start, int homenr,
220 real charge[], rvec f[],
221 t_cosines Ex[], t_cosines Et[], double t)
227 for (m = 0; (m < DIM); m++)
234 Ext[m] = cos(Et[m].a[0]*(t-t0))*exp(-sqr(t-t0)/(2.0*sqr(Et[m].a[2])));
238 Ext[m] = cos(Et[m].a[0]*t);
247 /* Convert the field strength from V/nm to MD-units */
248 Ext[m] *= Ex[m].a[0]*FIELDFAC;
249 for (i = start; (i < start+homenr); i++)
251 f[i][m] += charge[i]*Ext[m];
261 fprintf(fp, "%10g %10g %10g %10g #FIELD\n", t,
262 Ext[XX]/FIELDFAC, Ext[YY]/FIELDFAC, Ext[ZZ]/FIELDFAC);
266 static void calc_virial(int start, int homenr, rvec x[], rvec f[],
267 tensor vir_part, t_graph *graph, matrix box,
268 t_nrnb *nrnb, const t_forcerec *fr, int ePBC)
273 /* The short-range virial from surrounding boxes */
275 calc_vir(SHIFTS, fr->shift_vec, fr->fshift, vir_part, ePBC == epbcSCREW, box);
276 inc_nrnb(nrnb, eNR_VIRIAL, SHIFTS);
278 /* Calculate partial virial, for local atoms only, based on short range.
279 * Total virial is computed in global_stat, called from do_md
281 f_calc_vir(start, start+homenr, x, f, vir_part, graph, box);
282 inc_nrnb(nrnb, eNR_VIRIAL, homenr);
284 /* Add position restraint contribution */
285 for (i = 0; i < DIM; i++)
287 vir_part[i][i] += fr->vir_diag_posres[i];
290 /* Add wall contribution */
291 for (i = 0; i < DIM; i++)
293 vir_part[i][ZZ] += fr->vir_wall_z[i];
298 pr_rvecs(debug, 0, "vir_part", vir_part, DIM);
302 static void posres_wrapper(FILE *fplog,
308 matrix box, rvec x[],
309 gmx_enerdata_t *enerd,
317 /* Position restraints always require full pbc */
318 set_pbc(&pbc, ir->ePBC, box);
320 v = posres(top->idef.il[F_POSRES].nr, top->idef.il[F_POSRES].iatoms,
321 top->idef.iparams_posres,
322 (const rvec*)x, fr->f_novirsum, fr->vir_diag_posres,
323 ir->ePBC == epbcNONE ? NULL : &pbc,
324 lambda[efptRESTRAINT], &dvdl,
325 fr->rc_scaling, fr->ePBC, fr->posres_com, fr->posres_comB);
328 gmx_print_sepdvdl(fplog, interaction_function[F_POSRES].longname, v, dvdl);
330 enerd->term[F_POSRES] += v;
331 /* If just the force constant changes, the FEP term is linear,
332 * but if k changes, it is not.
334 enerd->dvdl_nonlin[efptRESTRAINT] += dvdl;
335 inc_nrnb(nrnb, eNR_POSRES, top->idef.il[F_POSRES].nr/2);
337 if ((ir->fepvals->n_lambda > 0) && (flags & GMX_FORCE_DHDL))
339 for (i = 0; i < enerd->n_lambda; i++)
341 real dvdl_dum, lambda_dum;
343 lambda_dum = (i == 0 ? lambda[efptRESTRAINT] : ir->fepvals->all_lambda[efptRESTRAINT][i-1]);
344 v = posres(top->idef.il[F_POSRES].nr, top->idef.il[F_POSRES].iatoms,
345 top->idef.iparams_posres,
346 (const rvec*)x, NULL, NULL,
347 ir->ePBC == epbcNONE ? NULL : &pbc, lambda_dum, &dvdl,
348 fr->rc_scaling, fr->ePBC, fr->posres_com, fr->posres_comB);
349 enerd->enerpart_lambda[i] += v;
354 static void fbposres_wrapper(t_inputrec *ir,
357 matrix box, rvec x[],
358 gmx_enerdata_t *enerd,
364 /* Flat-bottomed position restraints always require full pbc */
365 set_pbc(&pbc, ir->ePBC, box);
366 v = fbposres(top->idef.il[F_FBPOSRES].nr, top->idef.il[F_FBPOSRES].iatoms,
367 top->idef.iparams_fbposres,
368 (const rvec*)x, fr->f_novirsum, fr->vir_diag_posres,
369 ir->ePBC == epbcNONE ? NULL : &pbc,
370 fr->rc_scaling, fr->ePBC, fr->posres_com);
371 enerd->term[F_FBPOSRES] += v;
372 inc_nrnb(nrnb, eNR_FBPOSRES, top->idef.il[F_FBPOSRES].nr/2);
375 static void pull_potential_wrapper(FILE *fplog,
379 matrix box, rvec x[],
383 gmx_enerdata_t *enerd,
390 /* Calculate the center of mass forces, this requires communication,
391 * which is why pull_potential is called close to other communication.
392 * The virial contribution is calculated directly,
393 * which is why we call pull_potential after calc_virial.
395 set_pbc(&pbc, ir->ePBC, box);
397 enerd->term[F_COM_PULL] +=
398 pull_potential(ir->ePull, ir->pull, mdatoms, &pbc,
399 cr, t, lambda[efptRESTRAINT], x, f, vir_force, &dvdl);
402 gmx_print_sepdvdl(fplog, "Com pull", enerd->term[F_COM_PULL], dvdl);
404 enerd->dvdl_lin[efptRESTRAINT] += dvdl;
407 static void pme_receive_force_ener(FILE *fplog,
410 gmx_wallcycle_t wcycle,
411 gmx_enerdata_t *enerd,
414 real e_q, e_lj, v, dvdl_q, dvdl_lj;
415 float cycles_ppdpme, cycles_seppme;
417 cycles_ppdpme = wallcycle_stop(wcycle, ewcPPDURINGPME);
418 dd_cycles_add(cr->dd, cycles_ppdpme, ddCyclPPduringPME);
420 /* In case of node-splitting, the PP nodes receive the long-range
421 * forces, virial and energy from the PME nodes here.
423 wallcycle_start(wcycle, ewcPP_PMEWAITRECVF);
426 gmx_pme_receive_f(cr, fr->f_novirsum, fr->vir_el_recip, &e_q,
427 fr->vir_lj_recip, &e_lj, &dvdl_q, &dvdl_lj,
431 gmx_print_sepdvdl(fplog, "Electrostatic PME mesh", e_q, dvdl_q);
432 gmx_print_sepdvdl(fplog, "Lennard-Jones PME mesh", e_lj, dvdl_lj);
434 enerd->term[F_COUL_RECIP] += e_q;
435 enerd->term[F_LJ_RECIP] += e_lj;
436 enerd->dvdl_lin[efptCOUL] += dvdl_q;
437 enerd->dvdl_lin[efptVDW] += dvdl_lj;
441 dd_cycles_add(cr->dd, cycles_seppme, ddCyclPME);
443 wallcycle_stop(wcycle, ewcPP_PMEWAITRECVF);
446 static void print_large_forces(FILE *fp, t_mdatoms *md, t_commrec *cr,
447 gmx_int64_t step, real pforce, rvec *x, rvec *f)
451 char buf[STEPSTRSIZE];
454 for (i = md->start; i < md->start+md->homenr; i++)
457 /* We also catch NAN, if the compiler does not optimize this away. */
458 if (fn2 >= pf2 || fn2 != fn2)
460 fprintf(fp, "step %s atom %6d x %8.3f %8.3f %8.3f force %12.5e\n",
461 gmx_step_str(step, buf),
462 ddglatnr(cr->dd, i), x[i][XX], x[i][YY], x[i][ZZ], sqrt(fn2));
467 static void post_process_forces(t_commrec *cr,
469 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
471 matrix box, rvec x[],
476 t_forcerec *fr, gmx_vsite_t *vsite,
483 /* Spread the mesh force on virtual sites to the other particles...
484 * This is parallellized. MPI communication is performed
485 * if the constructing atoms aren't local.
487 wallcycle_start(wcycle, ewcVSITESPREAD);
488 spread_vsite_f(vsite, x, fr->f_novirsum, NULL,
489 (flags & GMX_FORCE_VIRIAL), fr->vir_el_recip,
491 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
492 wallcycle_stop(wcycle, ewcVSITESPREAD);
494 if (flags & GMX_FORCE_VIRIAL)
496 /* Now add the forces, this is local */
499 sum_forces(0, fr->f_novirsum_n, f, fr->f_novirsum);
503 sum_forces(mdatoms->start, mdatoms->start+mdatoms->homenr,
506 if (EEL_FULL(fr->eeltype))
508 /* Add the mesh contribution to the virial */
509 m_add(vir_force, fr->vir_el_recip, vir_force);
511 if (EVDW_PME(fr->vdwtype))
513 /* Add the mesh contribution to the virial */
514 m_add(vir_force, fr->vir_lj_recip, vir_force);
518 pr_rvecs(debug, 0, "vir_force", vir_force, DIM);
523 if (fr->print_force >= 0)
525 print_large_forces(stderr, mdatoms, cr, step, fr->print_force, x, f);
529 static void do_nb_verlet(t_forcerec *fr,
530 interaction_const_t *ic,
531 gmx_enerdata_t *enerd,
532 int flags, int ilocality,
535 gmx_wallcycle_t wcycle)
537 int nnbl, kernel_type, enr_nbnxn_kernel_ljc, enr_nbnxn_kernel_lj;
539 nonbonded_verlet_group_t *nbvg;
542 if (!(flags & GMX_FORCE_NONBONDED))
544 /* skip non-bonded calculation */
548 nbvg = &fr->nbv->grp[ilocality];
550 /* CUDA kernel launch overhead is already timed separately */
551 if (fr->cutoff_scheme != ecutsVERLET)
553 gmx_incons("Invalid cut-off scheme passed!");
556 bCUDA = (nbvg->kernel_type == nbnxnk8x8x8_CUDA);
560 wallcycle_sub_start(wcycle, ewcsNONBONDED);
562 switch (nbvg->kernel_type)
564 case nbnxnk4x4_PlainC:
565 nbnxn_kernel_ref(&nbvg->nbl_lists,
571 enerd->grpp.ener[egCOULSR],
573 enerd->grpp.ener[egBHAMSR] :
574 enerd->grpp.ener[egLJSR]);
577 case nbnxnk4xN_SIMD_4xN:
578 nbnxn_kernel_simd_4xn(&nbvg->nbl_lists,
585 enerd->grpp.ener[egCOULSR],
587 enerd->grpp.ener[egBHAMSR] :
588 enerd->grpp.ener[egLJSR]);
590 case nbnxnk4xN_SIMD_2xNN:
591 nbnxn_kernel_simd_2xnn(&nbvg->nbl_lists,
598 enerd->grpp.ener[egCOULSR],
600 enerd->grpp.ener[egBHAMSR] :
601 enerd->grpp.ener[egLJSR]);
604 case nbnxnk8x8x8_CUDA:
605 nbnxn_cuda_launch_kernel(fr->nbv->cu_nbv, nbvg->nbat, flags, ilocality);
608 case nbnxnk8x8x8_PlainC:
609 nbnxn_kernel_gpu_ref(nbvg->nbl_lists.nbl[0],
614 nbvg->nbat->out[0].f,
616 enerd->grpp.ener[egCOULSR],
618 enerd->grpp.ener[egBHAMSR] :
619 enerd->grpp.ener[egLJSR]);
623 gmx_incons("Invalid nonbonded kernel type passed!");
628 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
631 if (EEL_RF(ic->eeltype) || ic->eeltype == eelCUT)
633 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_RF;
635 else if ((!bCUDA && nbvg->ewald_excl == ewaldexclAnalytical) ||
636 (bCUDA && nbnxn_cuda_is_kernel_ewald_analytical(fr->nbv->cu_nbv)))
638 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_EWALD;
642 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_TAB;
644 enr_nbnxn_kernel_lj = eNR_NBNXN_LJ;
645 if (flags & GMX_FORCE_ENERGY)
647 /* In eNR_??? the nbnxn F+E kernels are always the F kernel + 1 */
648 enr_nbnxn_kernel_ljc += 1;
649 enr_nbnxn_kernel_lj += 1;
652 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc,
653 nbvg->nbl_lists.natpair_ljq);
654 inc_nrnb(nrnb, enr_nbnxn_kernel_lj,
655 nbvg->nbl_lists.natpair_lj);
656 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc-eNR_NBNXN_LJ_RF+eNR_NBNXN_RF,
657 nbvg->nbl_lists.natpair_q);
660 void do_force_cutsVERLET(FILE *fplog, t_commrec *cr,
661 t_inputrec *inputrec,
662 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
664 gmx_groups_t gmx_unused *groups,
665 matrix box, rvec x[], history_t *hist,
669 gmx_enerdata_t *enerd, t_fcdata *fcd,
670 real *lambda, t_graph *graph,
671 t_forcerec *fr, interaction_const_t *ic,
672 gmx_vsite_t *vsite, rvec mu_tot,
673 double t, FILE *field, gmx_edsam_t ed,
681 gmx_bool bSepDVDL, bStateChanged, bNS, bFillGrid, bCalcCGCM, bBS;
682 gmx_bool bDoLongRange, bDoForces, bSepLRF, bUseGPU, bUseOrEmulGPU;
683 gmx_bool bDiffKernels = FALSE;
685 rvec vzero, box_diag;
687 float cycles_pme, cycles_force, cycles_wait_gpu;
688 nonbonded_verlet_t *nbv;
693 nb_kernel_type = fr->nbv->grp[0].kernel_type;
695 start = mdatoms->start;
696 homenr = mdatoms->homenr;
698 bSepDVDL = (fr->bSepDVDL && do_per_step(step, inputrec->nstlog));
700 clear_mat(vir_force);
703 if (DOMAINDECOMP(cr))
705 cg1 = cr->dd->ncg_tot;
716 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
717 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
718 bFillGrid = (bNS && bStateChanged);
719 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
720 bDoLongRange = (fr->bTwinRange && bNS && (flags & GMX_FORCE_DO_LR));
721 bDoForces = (flags & GMX_FORCE_FORCES);
722 bSepLRF = (bDoLongRange && bDoForces && (flags & GMX_FORCE_SEPLRF));
723 bUseGPU = fr->nbv->bUseGPU;
724 bUseOrEmulGPU = bUseGPU || (nbv->grp[0].kernel_type == nbnxnk8x8x8_PlainC);
728 update_forcerec(fr, box);
730 if (NEED_MUTOT(*inputrec))
732 /* Calculate total (local) dipole moment in a temporary common array.
733 * This makes it possible to sum them over nodes faster.
735 calc_mu(start, homenr,
736 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
741 if (fr->ePBC != epbcNONE)
743 /* Compute shift vectors every step,
744 * because of pressure coupling or box deformation!
746 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
748 calc_shifts(box, fr->shift_vec);
753 put_atoms_in_box_omp(fr->ePBC, box, homenr, x);
754 inc_nrnb(nrnb, eNR_SHIFTX, homenr);
756 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
758 unshift_self(graph, box, x);
762 nbnxn_atomdata_copy_shiftvec(flags & GMX_FORCE_DYNAMICBOX,
763 fr->shift_vec, nbv->grp[0].nbat);
766 if (!(cr->duty & DUTY_PME))
768 /* Send particle coordinates to the pme nodes.
769 * Since this is only implemented for domain decomposition
770 * and domain decomposition does not use the graph,
771 * we do not need to worry about shifting.
776 wallcycle_start(wcycle, ewcPP_PMESENDX);
778 bBS = (inputrec->nwall == 2);
782 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
785 if (EEL_PME(fr->eeltype))
787 pme_flags |= GMX_PME_DO_COULOMB;
790 if (EVDW_PME(fr->vdwtype))
792 pme_flags |= GMX_PME_DO_LJ;
793 if (fr->ljpme_combination_rule == eljpmeLB)
795 pme_flags |= GMX_PME_LJ_LB;
799 gmx_pme_send_coordinates(cr, bBS ? boxs : box, x,
800 mdatoms->nChargePerturbed, mdatoms->nTypePerturbed, lambda[efptCOUL], lambda[efptVDW],
801 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),
804 wallcycle_stop(wcycle, ewcPP_PMESENDX);
808 /* do gridding for pair search */
811 if (graph && bStateChanged)
813 /* Calculate intramolecular shift vectors to make molecules whole */
814 mk_mshift(fplog, graph, fr->ePBC, box, x);
818 box_diag[XX] = box[XX][XX];
819 box_diag[YY] = box[YY][YY];
820 box_diag[ZZ] = box[ZZ][ZZ];
822 wallcycle_start(wcycle, ewcNS);
825 wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
826 nbnxn_put_on_grid(nbv->nbs, fr->ePBC, box,
828 0, mdatoms->homenr, -1, fr->cginfo, x,
830 nbv->grp[eintLocal].kernel_type,
831 nbv->grp[eintLocal].nbat);
832 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
836 wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
837 nbnxn_put_on_grid_nonlocal(nbv->nbs, domdec_zones(cr->dd),
839 nbv->grp[eintNonlocal].kernel_type,
840 nbv->grp[eintNonlocal].nbat);
841 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
844 if (nbv->ngrp == 1 ||
845 nbv->grp[eintNonlocal].nbat == nbv->grp[eintLocal].nbat)
847 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatAll,
848 nbv->nbs, mdatoms, fr->cginfo);
852 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatLocal,
853 nbv->nbs, mdatoms, fr->cginfo);
854 nbnxn_atomdata_set(nbv->grp[eintNonlocal].nbat, eatAll,
855 nbv->nbs, mdatoms, fr->cginfo);
857 wallcycle_stop(wcycle, ewcNS);
860 /* initialize the GPU atom data and copy shift vector */
865 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
866 nbnxn_cuda_init_atomdata(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
867 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
870 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
871 nbnxn_cuda_upload_shiftvec(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
872 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
875 /* do local pair search */
878 wallcycle_start_nocount(wcycle, ewcNS);
879 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
880 nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintLocal].nbat,
883 nbv->min_ci_balanced,
884 &nbv->grp[eintLocal].nbl_lists,
886 nbv->grp[eintLocal].kernel_type,
888 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
892 /* initialize local pair-list on the GPU */
893 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
894 nbv->grp[eintLocal].nbl_lists.nbl[0],
897 wallcycle_stop(wcycle, ewcNS);
901 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
902 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
903 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, FALSE, x,
904 nbv->grp[eintLocal].nbat);
905 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
906 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
911 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
912 /* launch local nonbonded F on GPU */
913 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFNo,
915 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
918 /* Communicate coordinates and sum dipole if necessary +
919 do non-local pair search */
920 if (DOMAINDECOMP(cr))
922 bDiffKernels = (nbv->grp[eintNonlocal].kernel_type !=
923 nbv->grp[eintLocal].kernel_type);
927 /* With GPU+CPU non-bonded calculations we need to copy
928 * the local coordinates to the non-local nbat struct
929 * (in CPU format) as the non-local kernel call also
930 * calculates the local - non-local interactions.
932 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
933 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
934 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, TRUE, x,
935 nbv->grp[eintNonlocal].nbat);
936 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
937 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
942 wallcycle_start_nocount(wcycle, ewcNS);
943 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
947 nbnxn_grid_add_simple(nbv->nbs, nbv->grp[eintNonlocal].nbat);
950 nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintNonlocal].nbat,
953 nbv->min_ci_balanced,
954 &nbv->grp[eintNonlocal].nbl_lists,
956 nbv->grp[eintNonlocal].kernel_type,
959 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
961 if (nbv->grp[eintNonlocal].kernel_type == nbnxnk8x8x8_CUDA)
963 /* initialize non-local pair-list on the GPU */
964 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
965 nbv->grp[eintNonlocal].nbl_lists.nbl[0],
968 wallcycle_stop(wcycle, ewcNS);
972 wallcycle_start(wcycle, ewcMOVEX);
973 dd_move_x(cr->dd, box, x);
975 /* When we don't need the total dipole we sum it in global_stat */
976 if (bStateChanged && NEED_MUTOT(*inputrec))
978 gmx_sumd(2*DIM, mu, cr);
980 wallcycle_stop(wcycle, ewcMOVEX);
982 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
983 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
984 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatNonlocal, FALSE, x,
985 nbv->grp[eintNonlocal].nbat);
986 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
987 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
990 if (bUseGPU && !bDiffKernels)
992 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
993 /* launch non-local nonbonded F on GPU */
994 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFNo,
996 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1002 /* launch D2H copy-back F */
1003 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1004 if (DOMAINDECOMP(cr) && !bDiffKernels)
1006 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintNonlocal].nbat,
1007 flags, eatNonlocal);
1009 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintLocal].nbat,
1011 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1014 if (bStateChanged && NEED_MUTOT(*inputrec))
1018 gmx_sumd(2*DIM, mu, cr);
1021 for (i = 0; i < 2; i++)
1023 for (j = 0; j < DIM; j++)
1025 fr->mu_tot[i][j] = mu[i*DIM + j];
1029 if (fr->efep == efepNO)
1031 copy_rvec(fr->mu_tot[0], mu_tot);
1035 for (j = 0; j < DIM; j++)
1038 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] +
1039 lambda[efptCOUL]*fr->mu_tot[1][j];
1043 /* Reset energies */
1044 reset_enerdata(fr, bNS, enerd, MASTER(cr));
1045 clear_rvecs(SHIFTS, fr->fshift);
1047 if (DOMAINDECOMP(cr))
1049 if (!(cr->duty & DUTY_PME))
1051 wallcycle_start(wcycle, ewcPPDURINGPME);
1052 dd_force_flop_start(cr->dd, nrnb);
1058 /* Enforced rotation has its own cycle counter that starts after the collective
1059 * coordinates have been communicated. It is added to ddCyclF to allow
1060 * for proper load-balancing */
1061 wallcycle_start(wcycle, ewcROT);
1062 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1063 wallcycle_stop(wcycle, ewcROT);
1066 /* Start the force cycle counter.
1067 * This counter is stopped in do_forcelow_level.
1068 * No parallel communication should occur while this counter is running,
1069 * since that will interfere with the dynamic load balancing.
1071 wallcycle_start(wcycle, ewcFORCE);
1074 /* Reset forces for which the virial is calculated separately:
1075 * PME/Ewald forces if necessary */
1076 if (fr->bF_NoVirSum)
1078 if (flags & GMX_FORCE_VIRIAL)
1080 fr->f_novirsum = fr->f_novirsum_alloc;
1083 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1087 clear_rvecs(homenr, fr->f_novirsum+start);
1092 /* We are not calculating the pressure so we do not need
1093 * a separate array for forces that do not contribute
1100 /* Clear the short- and long-range forces */
1101 clear_rvecs(fr->natoms_force_constr, f);
1102 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1104 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1107 clear_rvec(fr->vir_diag_posres);
1110 if (inputrec->ePull == epullCONSTRAINT)
1112 clear_pull_forces(inputrec->pull);
1115 /* We calculate the non-bonded forces, when done on the CPU, here.
1116 * We do this before calling do_force_lowlevel, as in there bondeds
1117 * forces are calculated before PME, which does communication.
1118 * With this order, non-bonded and bonded force calculation imbalance
1119 * can be balanced out by the domain decomposition load balancing.
1124 /* Maybe we should move this into do_force_lowlevel */
1125 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFYes,
1129 if (!bUseOrEmulGPU || bDiffKernels)
1133 if (DOMAINDECOMP(cr))
1135 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal,
1136 bDiffKernels ? enbvClearFYes : enbvClearFNo,
1146 aloc = eintNonlocal;
1149 /* Add all the non-bonded force to the normal force array.
1150 * This can be split into a local a non-local part when overlapping
1151 * communication with calculation with domain decomposition.
1153 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1154 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1155 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1156 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatAll, nbv->grp[aloc].nbat, f);
1157 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1158 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1159 wallcycle_start_nocount(wcycle, ewcFORCE);
1161 /* if there are multiple fshift output buffers reduce them */
1162 if ((flags & GMX_FORCE_VIRIAL) &&
1163 nbv->grp[aloc].nbl_lists.nnbl > 1)
1165 nbnxn_atomdata_add_nbat_fshift_to_fshift(nbv->grp[aloc].nbat,
1170 /* update QMMMrec, if necessary */
1173 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1176 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1178 posres_wrapper(fplog, flags, bSepDVDL, inputrec, nrnb, top, box, x,
1182 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_FBPOSRES].nr > 0)
1184 fbposres_wrapper(inputrec, nrnb, top, box, x, enerd, fr);
1187 /* Compute the bonded and non-bonded energies and optionally forces */
1188 do_force_lowlevel(fplog, step, fr, inputrec, &(top->idef),
1189 cr, nrnb, wcycle, mdatoms,
1190 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, top, fr->born,
1191 &(top->atomtypes), bBornRadii, box,
1192 inputrec->fepvals, lambda, graph, &(top->excls), fr->mu_tot,
1193 flags, &cycles_pme);
1197 if (do_per_step(step, inputrec->nstcalclr))
1199 /* Add the long range forces to the short range forces */
1200 for (i = 0; i < fr->natoms_force_constr; i++)
1202 rvec_add(fr->f_twin[i], f[i], f[i]);
1207 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1211 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1214 if (bUseOrEmulGPU && !bDiffKernels)
1216 /* wait for non-local forces (or calculate in emulation mode) */
1217 if (DOMAINDECOMP(cr))
1223 wallcycle_start(wcycle, ewcWAIT_GPU_NB_NL);
1224 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1225 nbv->grp[eintNonlocal].nbat,
1227 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1229 cycles_tmp = wallcycle_stop(wcycle, ewcWAIT_GPU_NB_NL);
1230 cycles_wait_gpu += cycles_tmp;
1231 cycles_force += cycles_tmp;
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 cycles_wait_gpu += wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
1299 /* now clear the GPU outputs while we finish the step on the CPU */
1301 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1302 nbnxn_cuda_clear_outputs(nbv->cu_nbv, flags);
1303 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1307 wallcycle_start_nocount(wcycle, ewcFORCE);
1308 do_nb_verlet(fr, ic, enerd, flags, eintLocal,
1309 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
1311 wallcycle_stop(wcycle, ewcFORCE);
1313 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1314 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1315 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1317 /* skip the reduction if there was no non-local work to do */
1318 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatLocal,
1319 nbv->grp[eintLocal].nbat, f);
1321 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1322 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1325 if (DOMAINDECOMP(cr))
1327 dd_force_flop_stop(cr->dd, nrnb);
1330 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1333 dd_cycles_add(cr->dd, cycles_wait_gpu, ddCyclWaitGPU);
1340 if (IR_ELEC_FIELD(*inputrec))
1342 /* Compute forces due to electric field */
1343 calc_f_el(MASTER(cr) ? field : NULL,
1344 start, homenr, mdatoms->chargeA, fr->f_novirsum,
1345 inputrec->ex, inputrec->et, t);
1348 /* If we have NoVirSum forces, but we do not calculate the virial,
1349 * we sum fr->f_novirum=f later.
1351 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1353 wallcycle_start(wcycle, ewcVSITESPREAD);
1354 spread_vsite_f(vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1355 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1356 wallcycle_stop(wcycle, ewcVSITESPREAD);
1360 wallcycle_start(wcycle, ewcVSITESPREAD);
1361 spread_vsite_f(vsite, x, fr->f_twin, NULL, FALSE, NULL,
1363 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1364 wallcycle_stop(wcycle, ewcVSITESPREAD);
1368 if (flags & GMX_FORCE_VIRIAL)
1370 /* Calculation of the virial must be done after vsites! */
1371 calc_virial(mdatoms->start, mdatoms->homenr, x, f,
1372 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1376 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1378 pull_potential_wrapper(fplog, bSepDVDL, cr, inputrec, box, x,
1379 f, vir_force, mdatoms, enerd, lambda, t);
1382 /* Add the forces from enforced rotation potentials (if any) */
1385 wallcycle_start(wcycle, ewcROTadd);
1386 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
1387 wallcycle_stop(wcycle, ewcROTadd);
1390 if (PAR(cr) && !(cr->duty & DUTY_PME))
1392 /* In case of node-splitting, the PP nodes receive the long-range
1393 * forces, virial and energy from the PME nodes here.
1395 pme_receive_force_ener(fplog, bSepDVDL, cr, wcycle, enerd, fr);
1400 post_process_forces(cr, step, nrnb, wcycle,
1401 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
1405 /* Sum the potential energy terms from group contributions */
1406 sum_epot(&(enerd->grpp), enerd->term);
1409 void do_force_cutsGROUP(FILE *fplog, t_commrec *cr,
1410 t_inputrec *inputrec,
1411 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1412 gmx_localtop_t *top,
1413 gmx_groups_t *groups,
1414 matrix box, rvec x[], history_t *hist,
1418 gmx_enerdata_t *enerd, t_fcdata *fcd,
1419 real *lambda, t_graph *graph,
1420 t_forcerec *fr, gmx_vsite_t *vsite, rvec mu_tot,
1421 double t, FILE *field, gmx_edsam_t ed,
1422 gmx_bool bBornRadii,
1428 gmx_bool bSepDVDL, bStateChanged, bNS, bFillGrid, bCalcCGCM, bBS;
1429 gmx_bool bDoLongRangeNS, bDoForces, bDoPotential, bSepLRF;
1430 gmx_bool bDoAdressWF;
1432 rvec vzero, box_diag;
1433 real e, v, dvdlambda[efptNR];
1435 float cycles_pme, cycles_force;
1437 start = mdatoms->start;
1438 homenr = mdatoms->homenr;
1440 bSepDVDL = (fr->bSepDVDL && do_per_step(step, inputrec->nstlog));
1442 clear_mat(vir_force);
1446 pd_cg_range(cr, &cg0, &cg1);
1451 if (DOMAINDECOMP(cr))
1453 cg1 = cr->dd->ncg_tot;
1465 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
1466 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
1467 /* Should we update the long-range neighborlists at this step? */
1468 bDoLongRangeNS = fr->bTwinRange && bNS;
1469 /* Should we perform the long-range nonbonded evaluation inside the neighborsearching? */
1470 bFillGrid = (bNS && bStateChanged);
1471 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
1472 bDoForces = (flags & GMX_FORCE_FORCES);
1473 bDoPotential = (flags & GMX_FORCE_ENERGY);
1474 bSepLRF = ((inputrec->nstcalclr > 1) && bDoForces &&
1475 (flags & GMX_FORCE_SEPLRF) && (flags & GMX_FORCE_DO_LR));
1477 /* should probably move this to the forcerec since it doesn't change */
1478 bDoAdressWF = ((fr->adress_type != eAdressOff));
1482 update_forcerec(fr, box);
1484 if (NEED_MUTOT(*inputrec))
1486 /* Calculate total (local) dipole moment in a temporary common array.
1487 * This makes it possible to sum them over nodes faster.
1489 calc_mu(start, homenr,
1490 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
1495 if (fr->ePBC != epbcNONE)
1497 /* Compute shift vectors every step,
1498 * because of pressure coupling or box deformation!
1500 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
1502 calc_shifts(box, fr->shift_vec);
1507 put_charge_groups_in_box(fplog, cg0, cg1, fr->ePBC, box,
1508 &(top->cgs), x, fr->cg_cm);
1509 inc_nrnb(nrnb, eNR_CGCM, homenr);
1510 inc_nrnb(nrnb, eNR_RESETX, cg1-cg0);
1512 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
1514 unshift_self(graph, box, x);
1519 calc_cgcm(fplog, cg0, cg1, &(top->cgs), x, fr->cg_cm);
1520 inc_nrnb(nrnb, eNR_CGCM, homenr);
1527 move_cgcm(fplog, cr, fr->cg_cm);
1531 pr_rvecs(debug, 0, "cgcm", fr->cg_cm, top->cgs.nr);
1536 if (!(cr->duty & DUTY_PME))
1538 /* Send particle coordinates to the pme nodes.
1539 * Since this is only implemented for domain decomposition
1540 * and domain decomposition does not use the graph,
1541 * we do not need to worry about shifting.
1546 wallcycle_start(wcycle, ewcPP_PMESENDX);
1548 bBS = (inputrec->nwall == 2);
1551 copy_mat(box, boxs);
1552 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
1555 if (EEL_PME(fr->eeltype))
1557 pme_flags |= GMX_PME_DO_COULOMB;
1560 if (EVDW_PME(fr->vdwtype))
1562 pme_flags |= GMX_PME_DO_LJ;
1563 if (fr->ljpme_combination_rule == eljpmeLB)
1565 pme_flags |= GMX_PME_LJ_LB;
1569 gmx_pme_send_coordinates(cr, bBS ? boxs : box, x,
1570 mdatoms->nChargePerturbed, mdatoms->nTypePerturbed, lambda[efptCOUL], lambda[efptVDW],
1571 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),
1574 wallcycle_stop(wcycle, ewcPP_PMESENDX);
1576 #endif /* GMX_MPI */
1578 /* Communicate coordinates and sum dipole if necessary */
1581 wallcycle_start(wcycle, ewcMOVEX);
1582 if (DOMAINDECOMP(cr))
1584 dd_move_x(cr->dd, box, x);
1588 move_x(cr, x, nrnb);
1590 wallcycle_stop(wcycle, ewcMOVEX);
1593 /* update adress weight beforehand */
1594 if (bStateChanged && bDoAdressWF)
1596 /* need pbc for adress weight calculation with pbc_dx */
1597 set_pbc(&pbc, inputrec->ePBC, box);
1598 if (fr->adress_site == eAdressSITEcog)
1600 update_adress_weights_cog(top->idef.iparams, top->idef.il, x, fr, mdatoms,
1601 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1603 else if (fr->adress_site == eAdressSITEcom)
1605 update_adress_weights_com(fplog, cg0, cg1, &(top->cgs), x, fr, mdatoms,
1606 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1608 else if (fr->adress_site == eAdressSITEatomatom)
1610 update_adress_weights_atom_per_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1611 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1615 update_adress_weights_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1616 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1620 if (NEED_MUTOT(*inputrec))
1627 gmx_sumd(2*DIM, mu, cr);
1629 for (i = 0; i < 2; i++)
1631 for (j = 0; j < DIM; j++)
1633 fr->mu_tot[i][j] = mu[i*DIM + j];
1637 if (fr->efep == efepNO)
1639 copy_rvec(fr->mu_tot[0], mu_tot);
1643 for (j = 0; j < DIM; j++)
1646 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] + lambda[efptCOUL]*fr->mu_tot[1][j];
1651 /* Reset energies */
1652 reset_enerdata(fr, bNS, enerd, MASTER(cr));
1653 clear_rvecs(SHIFTS, fr->fshift);
1657 wallcycle_start(wcycle, ewcNS);
1659 if (graph && bStateChanged)
1661 /* Calculate intramolecular shift vectors to make molecules whole */
1662 mk_mshift(fplog, graph, fr->ePBC, box, x);
1665 /* Do the actual neighbour searching */
1667 groups, top, mdatoms,
1668 cr, nrnb, bFillGrid,
1671 wallcycle_stop(wcycle, ewcNS);
1674 if (inputrec->implicit_solvent && bNS)
1676 make_gb_nblist(cr, inputrec->gb_algorithm,
1677 x, box, fr, &top->idef, graph, fr->born);
1680 if (DOMAINDECOMP(cr))
1682 if (!(cr->duty & DUTY_PME))
1684 wallcycle_start(wcycle, ewcPPDURINGPME);
1685 dd_force_flop_start(cr->dd, nrnb);
1691 /* Enforced rotation has its own cycle counter that starts after the collective
1692 * coordinates have been communicated. It is added to ddCyclF to allow
1693 * for proper load-balancing */
1694 wallcycle_start(wcycle, ewcROT);
1695 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1696 wallcycle_stop(wcycle, ewcROT);
1699 /* Start the force cycle counter.
1700 * This counter is stopped in do_forcelow_level.
1701 * No parallel communication should occur while this counter is running,
1702 * since that will interfere with the dynamic load balancing.
1704 wallcycle_start(wcycle, ewcFORCE);
1708 /* Reset forces for which the virial is calculated separately:
1709 * PME/Ewald forces if necessary */
1710 if (fr->bF_NoVirSum)
1712 if (flags & GMX_FORCE_VIRIAL)
1714 fr->f_novirsum = fr->f_novirsum_alloc;
1717 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1721 clear_rvecs(homenr, fr->f_novirsum+start);
1726 /* We are not calculating the pressure so we do not need
1727 * a separate array for forces that do not contribute
1734 /* Clear the short- and long-range forces */
1735 clear_rvecs(fr->natoms_force_constr, f);
1736 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1738 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1741 clear_rvec(fr->vir_diag_posres);
1743 if (inputrec->ePull == epullCONSTRAINT)
1745 clear_pull_forces(inputrec->pull);
1748 /* update QMMMrec, if necessary */
1751 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1754 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1756 posres_wrapper(fplog, flags, bSepDVDL, inputrec, nrnb, top, box, x,
1760 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_FBPOSRES].nr > 0)
1762 fbposres_wrapper(inputrec, nrnb, top, box, x, enerd, fr);
1765 /* Compute the bonded and non-bonded energies and optionally forces */
1766 do_force_lowlevel(fplog, step, fr, inputrec, &(top->idef),
1767 cr, nrnb, wcycle, mdatoms,
1768 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, top, fr->born,
1769 &(top->atomtypes), bBornRadii, box,
1770 inputrec->fepvals, lambda,
1771 graph, &(top->excls), fr->mu_tot,
1777 if (do_per_step(step, inputrec->nstcalclr))
1779 /* Add the long range forces to the short range forces */
1780 for (i = 0; i < fr->natoms_force_constr; i++)
1782 rvec_add(fr->f_twin[i], f[i], f[i]);
1787 cycles_force = wallcycle_stop(wcycle, ewcFORCE);
1791 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1794 if (DOMAINDECOMP(cr))
1796 dd_force_flop_stop(cr->dd, nrnb);
1799 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1805 if (IR_ELEC_FIELD(*inputrec))
1807 /* Compute forces due to electric field */
1808 calc_f_el(MASTER(cr) ? field : NULL,
1809 start, homenr, mdatoms->chargeA, fr->f_novirsum,
1810 inputrec->ex, inputrec->et, t);
1813 if (bDoAdressWF && fr->adress_icor == eAdressICThermoForce)
1815 /* Compute thermodynamic force in hybrid AdResS region */
1816 adress_thermo_force(start, homenr, &(top->cgs), x, fr->f_novirsum, fr, mdatoms,
1817 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1820 /* Communicate the forces */
1823 wallcycle_start(wcycle, ewcMOVEF);
1824 if (DOMAINDECOMP(cr))
1826 dd_move_f(cr->dd, f, fr->fshift);
1827 /* Do we need to communicate the separate force array
1828 * for terms that do not contribute to the single sum virial?
1829 * Position restraints and electric fields do not introduce
1830 * inter-cg forces, only full electrostatics methods do.
1831 * When we do not calculate the virial, fr->f_novirsum = f,
1832 * so we have already communicated these forces.
1834 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1835 (flags & GMX_FORCE_VIRIAL))
1837 dd_move_f(cr->dd, fr->f_novirsum, NULL);
1841 /* We should not update the shift forces here,
1842 * since f_twin is already included in f.
1844 dd_move_f(cr->dd, fr->f_twin, NULL);
1849 pd_move_f(cr, f, nrnb);
1852 pd_move_f(cr, fr->f_twin, nrnb);
1855 wallcycle_stop(wcycle, ewcMOVEF);
1858 /* If we have NoVirSum forces, but we do not calculate the virial,
1859 * we sum fr->f_novirum=f later.
1861 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1863 wallcycle_start(wcycle, ewcVSITESPREAD);
1864 spread_vsite_f(vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1865 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1866 wallcycle_stop(wcycle, ewcVSITESPREAD);
1870 wallcycle_start(wcycle, ewcVSITESPREAD);
1871 spread_vsite_f(vsite, x, fr->f_twin, NULL, FALSE, NULL,
1873 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1874 wallcycle_stop(wcycle, ewcVSITESPREAD);
1878 if (flags & GMX_FORCE_VIRIAL)
1880 /* Calculation of the virial must be done after vsites! */
1881 calc_virial(mdatoms->start, mdatoms->homenr, x, f,
1882 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1886 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1888 pull_potential_wrapper(fplog, bSepDVDL, cr, inputrec, box, x,
1889 f, vir_force, mdatoms, enerd, lambda, t);
1892 /* Add the forces from enforced rotation potentials (if any) */
1895 wallcycle_start(wcycle, ewcROTadd);
1896 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
1897 wallcycle_stop(wcycle, ewcROTadd);
1900 if (PAR(cr) && !(cr->duty & DUTY_PME))
1902 /* In case of node-splitting, the PP nodes receive the long-range
1903 * forces, virial and energy from the PME nodes here.
1905 pme_receive_force_ener(fplog, bSepDVDL, cr, wcycle, enerd, fr);
1910 post_process_forces(cr, step, nrnb, wcycle,
1911 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
1915 /* Sum the potential energy terms from group contributions */
1916 sum_epot(&(enerd->grpp), enerd->term);
1919 void do_force(FILE *fplog, t_commrec *cr,
1920 t_inputrec *inputrec,
1921 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1922 gmx_localtop_t *top,
1923 gmx_groups_t *groups,
1924 matrix box, rvec x[], history_t *hist,
1928 gmx_enerdata_t *enerd, t_fcdata *fcd,
1929 real *lambda, t_graph *graph,
1931 gmx_vsite_t *vsite, rvec mu_tot,
1932 double t, FILE *field, gmx_edsam_t ed,
1933 gmx_bool bBornRadii,
1936 /* modify force flag if not doing nonbonded */
1937 if (!fr->bNonbonded)
1939 flags &= ~GMX_FORCE_NONBONDED;
1942 switch (inputrec->cutoff_scheme)
1945 do_force_cutsVERLET(fplog, cr, inputrec,
1961 do_force_cutsGROUP(fplog, cr, inputrec,
1976 gmx_incons("Invalid cut-off scheme passed!");
1981 void do_constrain_first(FILE *fplog, gmx_constr_t constr,
1982 t_inputrec *ir, t_mdatoms *md,
1983 t_state *state, t_commrec *cr, t_nrnb *nrnb,
1984 t_forcerec *fr, gmx_localtop_t *top)
1986 int i, m, start, end;
1988 real dt = ir->delta_t;
1992 snew(savex, state->natoms);
1995 end = md->homenr + start;
1999 fprintf(debug, "vcm: start=%d, homenr=%d, end=%d\n",
2000 start, md->homenr, end);
2002 /* Do a first constrain to reset particles... */
2003 step = ir->init_step;
2006 char buf[STEPSTRSIZE];
2007 fprintf(fplog, "\nConstraining the starting coordinates (step %s)\n",
2008 gmx_step_str(step, buf));
2012 /* constrain the current position */
2013 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2014 ir, NULL, cr, step, 0, md,
2015 state->x, state->x, NULL,
2016 fr->bMolPBC, state->box,
2017 state->lambda[efptBONDED], &dvdl_dum,
2018 NULL, NULL, nrnb, econqCoord,
2019 ir->epc == epcMTTK, state->veta, state->veta);
2022 /* constrain the inital velocity, and save it */
2023 /* also may be useful if we need the ekin from the halfstep for velocity verlet */
2024 /* might not yet treat veta correctly */
2025 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2026 ir, NULL, cr, step, 0, md,
2027 state->x, state->v, state->v,
2028 fr->bMolPBC, state->box,
2029 state->lambda[efptBONDED], &dvdl_dum,
2030 NULL, NULL, nrnb, econqVeloc,
2031 ir->epc == epcMTTK, state->veta, state->veta);
2033 /* constrain the inital velocities at t-dt/2 */
2034 if (EI_STATE_VELOCITY(ir->eI) && ir->eI != eiVV)
2036 for (i = start; (i < end); i++)
2038 for (m = 0; (m < DIM); m++)
2040 /* Reverse the velocity */
2041 state->v[i][m] = -state->v[i][m];
2042 /* Store the position at t-dt in buf */
2043 savex[i][m] = state->x[i][m] + dt*state->v[i][m];
2046 /* Shake the positions at t=-dt with the positions at t=0
2047 * as reference coordinates.
2051 char buf[STEPSTRSIZE];
2052 fprintf(fplog, "\nConstraining the coordinates at t0-dt (step %s)\n",
2053 gmx_step_str(step, buf));
2056 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2057 ir, NULL, cr, step, -1, md,
2058 state->x, savex, NULL,
2059 fr->bMolPBC, state->box,
2060 state->lambda[efptBONDED], &dvdl_dum,
2061 state->v, NULL, nrnb, econqCoord,
2062 ir->epc == epcMTTK, state->veta, state->veta);
2064 for (i = start; i < end; i++)
2066 for (m = 0; m < DIM; m++)
2068 /* Re-reverse the velocities */
2069 state->v[i][m] = -state->v[i][m];
2078 integrate_table(real vdwtab[], real scale, int offstart, int rstart, int rend,
2079 double *enerout, double *virout)
2081 double enersum, virsum;
2082 double invscale, invscale2, invscale3;
2083 double r, ea, eb, ec, pa, pb, pc, pd;
2085 int ri, offset, tabfactor;
2087 invscale = 1.0/scale;
2088 invscale2 = invscale*invscale;
2089 invscale3 = invscale*invscale2;
2091 /* Following summation derived from cubic spline definition,
2092 * Numerical Recipies in C, second edition, p. 113-116. Exact for
2093 * the cubic spline. We first calculate the negative of the
2094 * energy from rvdw to rvdw_switch, assuming that g(r)=1, and then
2095 * add the more standard, abrupt cutoff correction to that result,
2096 * yielding the long-range correction for a switched function. We
2097 * perform both the pressure and energy loops at the same time for
2098 * simplicity, as the computational cost is low. */
2102 /* Since the dispersion table has been scaled down a factor
2103 * 6.0 and the repulsion a factor 12.0 to compensate for the
2104 * c6/c12 parameters inside nbfp[] being scaled up (to save
2105 * flops in kernels), we need to correct for this.
2116 for (ri = rstart; ri < rend; ++ri)
2120 eb = 2.0*invscale2*r;
2124 pb = 3.0*invscale2*r;
2125 pc = 3.0*invscale*r*r;
2128 /* this "8" is from the packing in the vdwtab array - perhaps
2129 should be defined? */
2131 offset = 8*ri + offstart;
2132 y0 = vdwtab[offset];
2133 f = vdwtab[offset+1];
2134 g = vdwtab[offset+2];
2135 h = vdwtab[offset+3];
2137 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);
2138 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);
2140 *enerout = 4.0*M_PI*enersum*tabfactor;
2141 *virout = 4.0*M_PI*virsum*tabfactor;
2144 void calc_enervirdiff(FILE *fplog, int eDispCorr, t_forcerec *fr)
2146 double eners[2], virs[2], enersum, virsum, y0, f, g, h;
2147 double r0, r1, r, rc3, rc9, ea, eb, ec, pa, pb, pc, pd;
2148 double invscale, invscale2, invscale3;
2149 int ri0, ri1, ri, i, offstart, offset;
2150 real scale, *vdwtab, tabfactor, tmp;
2152 fr->enershiftsix = 0;
2153 fr->enershifttwelve = 0;
2154 fr->enerdiffsix = 0;
2155 fr->enerdifftwelve = 0;
2157 fr->virdifftwelve = 0;
2159 if (eDispCorr != edispcNO)
2161 for (i = 0; i < 2; i++)
2166 if ((fr->vdwtype == evdwSWITCH) || (fr->vdwtype == evdwSHIFT))
2168 if (fr->rvdw_switch == 0)
2171 "With dispersion correction rvdw-switch can not be zero "
2172 "for vdw-type = %s", evdw_names[fr->vdwtype]);
2175 scale = fr->nblists[0].table_elec_vdw.scale;
2176 vdwtab = fr->nblists[0].table_vdw.data;
2178 /* Round the cut-offs to exact table values for precision */
2179 ri0 = floor(fr->rvdw_switch*scale);
2180 ri1 = ceil(fr->rvdw*scale);
2186 if (fr->vdwtype == evdwSHIFT)
2188 /* Determine the constant energy shift below rvdw_switch.
2189 * Table has a scale factor since we have scaled it down to compensate
2190 * for scaling-up c6/c12 with the derivative factors to save flops in analytical kernels.
2192 fr->enershiftsix = (real)(-1.0/(rc3*rc3)) - 6.0*vdwtab[8*ri0];
2193 fr->enershifttwelve = (real)( 1.0/(rc9*rc3)) - 12.0*vdwtab[8*ri0 + 4];
2195 /* Add the constant part from 0 to rvdw_switch.
2196 * This integration from 0 to rvdw_switch overcounts the number
2197 * of interactions by 1, as it also counts the self interaction.
2198 * We will correct for this later.
2200 eners[0] += 4.0*M_PI*fr->enershiftsix*rc3/3.0;
2201 eners[1] += 4.0*M_PI*fr->enershifttwelve*rc3/3.0;
2202 for (i = 0; i < 2; i++)
2206 integrate_table(vdwtab, scale, (i == 0 ? 0 : 4), ri0, ri1, &enersum, &virsum);
2207 eners[i] -= enersum;
2211 /* now add the correction for rvdw_switch to infinity */
2212 eners[0] += -4.0*M_PI/(3.0*rc3);
2213 eners[1] += 4.0*M_PI/(9.0*rc9);
2214 virs[0] += 8.0*M_PI/rc3;
2215 virs[1] += -16.0*M_PI/(3.0*rc9);
2217 else if (EVDW_PME(fr->vdwtype))
2219 if (EVDW_SWITCHED(fr->vdwtype) && fr->rvdw_switch == 0)
2222 "With dispersion correction rvdw-switch can not be zero "
2223 "for vdw-type = %s", evdw_names[fr->vdwtype]);
2226 scale = fr->nblists[0].table_vdw.scale;
2227 vdwtab = fr->nblists[0].table_vdw.data;
2229 ri0 = floor(fr->rvdw_switch*scale);
2230 ri1 = ceil(fr->rvdw*scale);
2236 /* Calculate self-interaction coefficient (assuming that
2237 * the reciprocal-space contribution is constant in the
2238 * region that contributes to the self-interaction).
2240 fr->enershiftsix = pow(fr->ewaldcoeff_lj, 6) / 6.0;
2242 /* Calculate C12 values as without PME. */
2243 if (EVDW_SWITCHED(fr->vdwtype))
2247 integrate_table(vdwtab, scale, 4, ri0, ri1, &enersum, &virsum);
2248 eners[1] -= enersum;
2251 /* Add analytical corrections, C6 for the whole range, C12
2252 * from rvdw_switch to infinity.
2255 eners[0] += -pow(sqrt(M_PI)*fr->ewaldcoeff_lj, 3)/3.0;
2256 eners[1] += 4.0*M_PI/(9.0*rc9);
2257 virs[0] += pow(sqrt(M_PI)*fr->ewaldcoeff_lj, 3);
2258 virs[1] += -16.0*M_PI/(3.0*rc9);
2260 else if ((fr->vdwtype == evdwCUT) || (fr->vdwtype == evdwUSER))
2262 if (fr->vdwtype == evdwUSER && fplog)
2265 "WARNING: using dispersion correction with user tables\n");
2267 rc3 = fr->rvdw*fr->rvdw*fr->rvdw;
2269 /* Contribution beyond the cut-off */
2270 eners[0] += -4.0*M_PI/(3.0*rc3);
2271 eners[1] += 4.0*M_PI/(9.0*rc9);
2272 if (fr->vdw_modifier == eintmodPOTSHIFT)
2274 /* Contribution within the cut-off */
2275 eners[0] += -4.0*M_PI/(3.0*rc3);
2276 eners[1] += 4.0*M_PI/(3.0*rc9);
2278 /* Contribution beyond the cut-off */
2279 virs[0] += 8.0*M_PI/rc3;
2280 virs[1] += -16.0*M_PI/(3.0*rc9);
2285 "Dispersion correction is not implemented for vdw-type = %s",
2286 evdw_names[fr->vdwtype]);
2288 fr->enerdiffsix = eners[0];
2289 fr->enerdifftwelve = eners[1];
2290 /* The 0.5 is due to the Gromacs definition of the virial */
2291 fr->virdiffsix = 0.5*virs[0];
2292 fr->virdifftwelve = 0.5*virs[1];
2296 void calc_dispcorr(FILE *fplog, t_inputrec *ir, t_forcerec *fr,
2297 gmx_int64_t step, int natoms,
2298 matrix box, real lambda, tensor pres, tensor virial,
2299 real *prescorr, real *enercorr, real *dvdlcorr)
2301 gmx_bool bCorrAll, bCorrPres;
2302 real dvdlambda, invvol, dens, ninter, avcsix, avctwelve, enerdiff, svir = 0, spres = 0;
2312 if (ir->eDispCorr != edispcNO)
2314 bCorrAll = (ir->eDispCorr == edispcAllEner ||
2315 ir->eDispCorr == edispcAllEnerPres);
2316 bCorrPres = (ir->eDispCorr == edispcEnerPres ||
2317 ir->eDispCorr == edispcAllEnerPres);
2319 invvol = 1/det(box);
2322 /* Only correct for the interactions with the inserted molecule */
2323 dens = (natoms - fr->n_tpi)*invvol;
2328 dens = natoms*invvol;
2329 ninter = 0.5*natoms;
2332 if (ir->efep == efepNO)
2334 avcsix = fr->avcsix[0];
2335 avctwelve = fr->avctwelve[0];
2339 avcsix = (1 - lambda)*fr->avcsix[0] + lambda*fr->avcsix[1];
2340 avctwelve = (1 - lambda)*fr->avctwelve[0] + lambda*fr->avctwelve[1];
2343 enerdiff = ninter*(dens*fr->enerdiffsix - fr->enershiftsix);
2344 *enercorr += avcsix*enerdiff;
2346 if (ir->efep != efepNO)
2348 dvdlambda += (fr->avcsix[1] - fr->avcsix[0])*enerdiff;
2352 enerdiff = ninter*(dens*fr->enerdifftwelve - fr->enershifttwelve);
2353 *enercorr += avctwelve*enerdiff;
2354 if (fr->efep != efepNO)
2356 dvdlambda += (fr->avctwelve[1] - fr->avctwelve[0])*enerdiff;
2362 svir = ninter*dens*avcsix*fr->virdiffsix/3.0;
2363 if (ir->eDispCorr == edispcAllEnerPres)
2365 svir += ninter*dens*avctwelve*fr->virdifftwelve/3.0;
2367 /* The factor 2 is because of the Gromacs virial definition */
2368 spres = -2.0*invvol*svir*PRESFAC;
2370 for (m = 0; m < DIM; m++)
2372 virial[m][m] += svir;
2373 pres[m][m] += spres;
2378 /* Can't currently control when it prints, for now, just print when degugging */
2383 fprintf(debug, "Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
2389 "Long Range LJ corr.: Epot %10g, Pres: %10g, Vir: %10g\n",
2390 *enercorr, spres, svir);
2394 fprintf(debug, "Long Range LJ corr.: Epot %10g\n", *enercorr);
2398 if (fr->bSepDVDL && do_per_step(step, ir->nstlog))
2400 gmx_print_sepdvdl(fplog, "Dispersion correction", *enercorr, dvdlambda);
2402 if (fr->efep != efepNO)
2404 *dvdlcorr += dvdlambda;
2409 void do_pbc_first(FILE *fplog, matrix box, t_forcerec *fr,
2410 t_graph *graph, rvec x[])
2414 fprintf(fplog, "Removing pbc first time\n");
2416 calc_shifts(box, fr->shift_vec);
2419 mk_mshift(fplog, graph, fr->ePBC, box, x);
2422 p_graph(debug, "do_pbc_first 1", graph);
2424 shift_self(graph, box, x);
2425 /* By doing an extra mk_mshift the molecules that are broken
2426 * because they were e.g. imported from another software
2427 * will be made whole again. Such are the healing powers
2430 mk_mshift(fplog, graph, fr->ePBC, box, x);
2433 p_graph(debug, "do_pbc_first 2", graph);
2438 fprintf(fplog, "Done rmpbc\n");
2442 static void low_do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2443 gmx_mtop_t *mtop, rvec x[],
2448 gmx_molblock_t *molb;
2450 if (bFirst && fplog)
2452 fprintf(fplog, "Removing pbc first time\n");
2457 for (mb = 0; mb < mtop->nmolblock; mb++)
2459 molb = &mtop->molblock[mb];
2460 if (molb->natoms_mol == 1 ||
2461 (!bFirst && mtop->moltype[molb->type].cgs.nr == 1))
2463 /* Just one atom or charge group in the molecule, no PBC required */
2464 as += molb->nmol*molb->natoms_mol;
2468 /* Pass NULL iso fplog to avoid graph prints for each molecule type */
2469 mk_graph_ilist(NULL, mtop->moltype[molb->type].ilist,
2470 0, molb->natoms_mol, FALSE, FALSE, graph);
2472 for (mol = 0; mol < molb->nmol; mol++)
2474 mk_mshift(fplog, graph, ePBC, box, x+as);
2476 shift_self(graph, box, x+as);
2477 /* The molecule is whole now.
2478 * We don't need the second mk_mshift call as in do_pbc_first,
2479 * since we no longer need this graph.
2482 as += molb->natoms_mol;
2490 void do_pbc_first_mtop(FILE *fplog, int ePBC, matrix box,
2491 gmx_mtop_t *mtop, rvec x[])
2493 low_do_pbc_mtop(fplog, ePBC, box, mtop, x, TRUE);
2496 void do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2497 gmx_mtop_t *mtop, rvec x[])
2499 low_do_pbc_mtop(fplog, ePBC, box, mtop, x, FALSE);
2502 void finish_run(FILE *fplog, t_commrec *cr,
2503 t_inputrec *inputrec,
2504 t_nrnb nrnb[], gmx_wallcycle_t wcycle,
2505 gmx_walltime_accounting_t walltime_accounting,
2506 wallclock_gpu_t *gputimes,
2507 gmx_bool bWriteStat)
2510 t_nrnb *nrnb_tot = NULL;
2513 double elapsed_time,
2514 elapsed_time_over_all_ranks,
2515 elapsed_time_over_all_threads,
2516 elapsed_time_over_all_threads_over_all_ranks;
2517 wallcycle_sum(cr, wcycle);
2523 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM,
2524 cr->mpi_comm_mysim);
2532 elapsed_time = walltime_accounting_get_elapsed_time(walltime_accounting);
2533 elapsed_time_over_all_ranks = elapsed_time;
2534 elapsed_time_over_all_threads = walltime_accounting_get_elapsed_time_over_all_threads(walltime_accounting);
2535 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
2539 /* reduce elapsed_time over all MPI ranks in the current simulation */
2540 MPI_Allreduce(&elapsed_time,
2541 &elapsed_time_over_all_ranks,
2542 1, MPI_DOUBLE, MPI_SUM,
2543 cr->mpi_comm_mysim);
2544 elapsed_time_over_all_ranks /= cr->nnodes;
2545 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
2546 * current simulation. */
2547 MPI_Allreduce(&elapsed_time_over_all_threads,
2548 &elapsed_time_over_all_threads_over_all_ranks,
2549 1, MPI_DOUBLE, MPI_SUM,
2550 cr->mpi_comm_mysim);
2556 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
2563 if ((cr->duty & DUTY_PP) && DOMAINDECOMP(cr))
2565 print_dd_statistics(cr, inputrec, fplog);
2577 snew(nrnb_all, cr->nnodes);
2578 nrnb_all[0] = *nrnb;
2579 for (s = 1; s < cr->nnodes; s++)
2581 MPI_Recv(nrnb_all[s].n, eNRNB, MPI_DOUBLE, s, 0,
2582 cr->mpi_comm_mysim, &stat);
2584 pr_load(fplog, cr, nrnb_all);
2589 MPI_Send(nrnb->n, eNRNB, MPI_DOUBLE, MASTERRANK(cr), 0,
2590 cr->mpi_comm_mysim);
2597 wallcycle_print(fplog, cr->nnodes, cr->npmenodes,
2598 elapsed_time_over_all_ranks,
2601 if (EI_DYNAMICS(inputrec->eI))
2603 delta_t = inputrec->delta_t;
2612 print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks,
2613 elapsed_time_over_all_ranks,
2614 walltime_accounting_get_nsteps_done(walltime_accounting),
2615 delta_t, nbfs, mflop);
2619 print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks,
2620 elapsed_time_over_all_ranks,
2621 walltime_accounting_get_nsteps_done(walltime_accounting),
2622 delta_t, nbfs, mflop);
2627 extern void initialize_lambdas(FILE *fplog, t_inputrec *ir, int *fep_state, real *lambda, double *lam0)
2629 /* this function works, but could probably use a logic rewrite to keep all the different
2630 types of efep straight. */
2633 t_lambda *fep = ir->fepvals;
2635 if ((ir->efep == efepNO) && (ir->bSimTemp == FALSE))
2637 for (i = 0; i < efptNR; i++)
2649 *fep_state = fep->init_fep_state; /* this might overwrite the checkpoint
2650 if checkpoint is set -- a kludge is in for now
2652 for (i = 0; i < efptNR; i++)
2654 /* overwrite lambda state with init_lambda for now for backwards compatibility */
2655 if (fep->init_lambda >= 0) /* if it's -1, it was never initializd */
2657 lambda[i] = fep->init_lambda;
2660 lam0[i] = lambda[i];
2665 lambda[i] = fep->all_lambda[i][*fep_state];
2668 lam0[i] = lambda[i];
2674 /* need to rescale control temperatures to match current state */
2675 for (i = 0; i < ir->opts.ngtc; i++)
2677 if (ir->opts.ref_t[i] > 0)
2679 ir->opts.ref_t[i] = ir->simtempvals->temperatures[*fep_state];
2685 /* Send to the log the information on the current lambdas */
2688 fprintf(fplog, "Initial vector of lambda components:[ ");
2689 for (i = 0; i < efptNR; i++)
2691 fprintf(fplog, "%10.4f ", lambda[i]);
2693 fprintf(fplog, "]\n");
2699 void init_md(FILE *fplog,
2700 t_commrec *cr, t_inputrec *ir, const output_env_t oenv,
2701 double *t, double *t0,
2702 real *lambda, int *fep_state, double *lam0,
2703 t_nrnb *nrnb, gmx_mtop_t *mtop,
2705 int nfile, const t_filenm fnm[],
2706 gmx_mdoutf_t **outf, t_mdebin **mdebin,
2707 tensor force_vir, tensor shake_vir, rvec mu_tot,
2708 gmx_bool *bSimAnn, t_vcm **vcm, unsigned long Flags)
2713 /* Initial values */
2714 *t = *t0 = ir->init_t;
2717 for (i = 0; i < ir->opts.ngtc; i++)
2719 /* set bSimAnn if any group is being annealed */
2720 if (ir->opts.annealing[i] != eannNO)
2727 update_annealing_target_temp(&(ir->opts), ir->init_t);
2730 /* Initialize lambda variables */
2731 initialize_lambdas(fplog, ir, fep_state, lambda, lam0);
2735 *upd = init_update(ir);
2741 *vcm = init_vcm(fplog, &mtop->groups, ir);
2744 if (EI_DYNAMICS(ir->eI) && !(Flags & MD_APPENDFILES))
2746 if (ir->etc == etcBERENDSEN)
2748 please_cite(fplog, "Berendsen84a");
2750 if (ir->etc == etcVRESCALE)
2752 please_cite(fplog, "Bussi2007a");
2760 *outf = init_mdoutf(nfile, fnm, Flags, cr, ir, oenv);
2762 *mdebin = init_mdebin((Flags & MD_APPENDFILES) ? NULL : (*outf)->fp_ene,
2763 mtop, ir, (*outf)->fp_dhdl);
2768 please_cite(fplog, "Fritsch12");
2769 please_cite(fplog, "Junghans10");
2771 /* Initiate variables */
2772 clear_mat(force_vir);
2773 clear_mat(shake_vir);