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,2014, 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 "gromacs/random/random.h"
77 #include "nbnxn_atomdata.h"
78 #include "nbnxn_search.h"
79 #include "nbnxn_kernels/nbnxn_kernel_ref.h"
80 #include "nbnxn_kernels/simd_4xn/nbnxn_kernel_simd_4xn.h"
81 #include "nbnxn_kernels/simd_2xnn/nbnxn_kernel_simd_2xnn.h"
82 #include "nbnxn_kernels/nbnxn_kernel_gpu_ref.h"
84 #include "gromacs/timing/wallcycle.h"
85 #include "gromacs/timing/walltime_accounting.h"
86 #include "gromacs/utility/gmxmpi.h"
87 #include "gromacs/essentialdynamics/edsam.h"
88 #include "gromacs/pulling/pull.h"
89 #include "gromacs/pulling/pull_rotation.h"
94 #include "nbnxn_cuda_data_mgmt.h"
95 #include "nbnxn_cuda/nbnxn_cuda.h"
97 void print_time(FILE *out,
98 gmx_walltime_accounting_t walltime_accounting,
101 t_commrec gmx_unused *cr)
104 char timebuf[STRLEN];
105 double dt, elapsed_seconds, time_per_step;
108 #ifndef GMX_THREAD_MPI
114 fprintf(out, "step %s", gmx_step_str(step, buf));
115 if ((step >= ir->nstlist))
117 double seconds_since_epoch = gmx_gettime();
118 elapsed_seconds = seconds_since_epoch - walltime_accounting_get_start_time_stamp(walltime_accounting);
119 time_per_step = elapsed_seconds/(step - ir->init_step + 1);
120 dt = (ir->nsteps + ir->init_step - step) * time_per_step;
126 finish = (time_t) (seconds_since_epoch + dt);
127 gmx_ctime_r(&finish, timebuf, STRLEN);
128 sprintf(buf, "%s", timebuf);
129 buf[strlen(buf)-1] = '\0';
130 fprintf(out, ", will finish %s", buf);
134 fprintf(out, ", remaining wall clock time: %5d s ", (int)dt);
139 fprintf(out, " performance: %.1f ns/day ",
140 ir->delta_t/1000*24*60*60/time_per_step);
143 #ifndef GMX_THREAD_MPI
153 void print_date_and_time(FILE *fplog, int nodeid, const char *title,
154 const gmx_walltime_accounting_t walltime_accounting)
157 char timebuf[STRLEN];
158 char time_string[STRLEN];
163 if (walltime_accounting != NULL)
165 tmptime = (time_t) walltime_accounting_get_start_time_stamp(walltime_accounting);
166 gmx_ctime_r(&tmptime, timebuf, STRLEN);
170 tmptime = (time_t) gmx_gettime();
171 gmx_ctime_r(&tmptime, timebuf, STRLEN);
173 for (i = 0; timebuf[i] >= ' '; i++)
175 time_string[i] = timebuf[i];
177 time_string[i] = '\0';
179 fprintf(fplog, "%s on node %d %s\n", title, nodeid, time_string);
183 void print_start(FILE *fplog, t_commrec *cr,
184 gmx_walltime_accounting_t walltime_accounting,
189 sprintf(buf, "Started %s", name);
190 print_date_and_time(fplog, cr->nodeid, buf, walltime_accounting);
193 static void sum_forces(int start, int end, rvec f[], rvec flr[])
199 pr_rvecs(debug, 0, "fsr", f+start, end-start);
200 pr_rvecs(debug, 0, "flr", flr+start, end-start);
202 for (i = start; (i < end); i++)
204 rvec_inc(f[i], flr[i]);
209 * calc_f_el calculates forces due to an electric field.
211 * force is kJ mol^-1 nm^-1 = e * kJ mol^-1 nm^-1 / e
213 * Et[] contains the parameters for the time dependent
214 * part of the field (not yet used).
215 * Ex[] contains the parameters for
216 * the spatial dependent part of the field. You can have cool periodic
217 * fields in principle, but only a constant field is supported
219 * The function should return the energy due to the electric field
220 * (if any) but for now returns 0.
223 * There can be problems with the virial.
224 * Since the field is not self-consistent this is unavoidable.
225 * For neutral molecules the virial is correct within this approximation.
226 * For neutral systems with many charged molecules the error is small.
227 * But for systems with a net charge or a few charged molecules
228 * the error can be significant when the field is high.
229 * Solution: implement a self-consitent electric field into PME.
231 static void calc_f_el(FILE *fp, int start, int homenr,
232 real charge[], rvec f[],
233 t_cosines Ex[], t_cosines Et[], double t)
239 for (m = 0; (m < DIM); m++)
246 Ext[m] = cos(Et[m].a[0]*(t-t0))*exp(-sqr(t-t0)/(2.0*sqr(Et[m].a[2])));
250 Ext[m] = cos(Et[m].a[0]*t);
259 /* Convert the field strength from V/nm to MD-units */
260 Ext[m] *= Ex[m].a[0]*FIELDFAC;
261 for (i = start; (i < start+homenr); i++)
263 f[i][m] += charge[i]*Ext[m];
273 fprintf(fp, "%10g %10g %10g %10g #FIELD\n", t,
274 Ext[XX]/FIELDFAC, Ext[YY]/FIELDFAC, Ext[ZZ]/FIELDFAC);
278 static void calc_virial(int start, int homenr, rvec x[], rvec f[],
279 tensor vir_part, t_graph *graph, matrix box,
280 t_nrnb *nrnb, const t_forcerec *fr, int ePBC)
285 /* The short-range virial from surrounding boxes */
287 calc_vir(SHIFTS, fr->shift_vec, fr->fshift, vir_part, ePBC == epbcSCREW, box);
288 inc_nrnb(nrnb, eNR_VIRIAL, SHIFTS);
290 /* Calculate partial virial, for local atoms only, based on short range.
291 * Total virial is computed in global_stat, called from do_md
293 f_calc_vir(start, start+homenr, x, f, vir_part, graph, box);
294 inc_nrnb(nrnb, eNR_VIRIAL, homenr);
296 /* Add position restraint contribution */
297 for (i = 0; i < DIM; i++)
299 vir_part[i][i] += fr->vir_diag_posres[i];
302 /* Add wall contribution */
303 for (i = 0; i < DIM; i++)
305 vir_part[i][ZZ] += fr->vir_wall_z[i];
310 pr_rvecs(debug, 0, "vir_part", vir_part, DIM);
314 static void posres_wrapper(FILE *fplog,
320 matrix box, rvec x[],
321 gmx_enerdata_t *enerd,
329 /* Position restraints always require full pbc */
330 set_pbc(&pbc, ir->ePBC, box);
332 v = posres(top->idef.il[F_POSRES].nr, top->idef.il[F_POSRES].iatoms,
333 top->idef.iparams_posres,
334 (const rvec*)x, fr->f_novirsum, fr->vir_diag_posres,
335 ir->ePBC == epbcNONE ? NULL : &pbc,
336 lambda[efptRESTRAINT], &dvdl,
337 fr->rc_scaling, fr->ePBC, fr->posres_com, fr->posres_comB);
340 gmx_print_sepdvdl(fplog, interaction_function[F_POSRES].longname, v, dvdl);
342 enerd->term[F_POSRES] += v;
343 /* If just the force constant changes, the FEP term is linear,
344 * but if k changes, it is not.
346 enerd->dvdl_nonlin[efptRESTRAINT] += dvdl;
347 inc_nrnb(nrnb, eNR_POSRES, top->idef.il[F_POSRES].nr/2);
349 if ((ir->fepvals->n_lambda > 0) && (flags & GMX_FORCE_DHDL))
351 for (i = 0; i < enerd->n_lambda; i++)
353 real dvdl_dum, lambda_dum;
355 lambda_dum = (i == 0 ? lambda[efptRESTRAINT] : ir->fepvals->all_lambda[efptRESTRAINT][i-1]);
356 v = posres(top->idef.il[F_POSRES].nr, top->idef.il[F_POSRES].iatoms,
357 top->idef.iparams_posres,
358 (const rvec*)x, NULL, NULL,
359 ir->ePBC == epbcNONE ? NULL : &pbc, lambda_dum, &dvdl,
360 fr->rc_scaling, fr->ePBC, fr->posres_com, fr->posres_comB);
361 enerd->enerpart_lambda[i] += v;
366 static void fbposres_wrapper(t_inputrec *ir,
369 matrix box, rvec x[],
370 gmx_enerdata_t *enerd,
376 /* Flat-bottomed position restraints always require full pbc */
377 set_pbc(&pbc, ir->ePBC, box);
378 v = fbposres(top->idef.il[F_FBPOSRES].nr, top->idef.il[F_FBPOSRES].iatoms,
379 top->idef.iparams_fbposres,
380 (const rvec*)x, fr->f_novirsum, fr->vir_diag_posres,
381 ir->ePBC == epbcNONE ? NULL : &pbc,
382 fr->rc_scaling, fr->ePBC, fr->posres_com);
383 enerd->term[F_FBPOSRES] += v;
384 inc_nrnb(nrnb, eNR_FBPOSRES, top->idef.il[F_FBPOSRES].nr/2);
387 static void pull_potential_wrapper(FILE *fplog,
391 matrix box, rvec x[],
395 gmx_enerdata_t *enerd,
402 /* Calculate the center of mass forces, this requires communication,
403 * which is why pull_potential is called close to other communication.
404 * The virial contribution is calculated directly,
405 * which is why we call pull_potential after calc_virial.
407 set_pbc(&pbc, ir->ePBC, box);
409 enerd->term[F_COM_PULL] +=
410 pull_potential(ir->ePull, ir->pull, mdatoms, &pbc,
411 cr, t, lambda[efptRESTRAINT], x, f, vir_force, &dvdl);
414 gmx_print_sepdvdl(fplog, "Com pull", enerd->term[F_COM_PULL], dvdl);
416 enerd->dvdl_lin[efptRESTRAINT] += dvdl;
419 static void pme_receive_force_ener(FILE *fplog,
422 gmx_wallcycle_t wcycle,
423 gmx_enerdata_t *enerd,
426 real e_q, e_lj, v, dvdl_q, dvdl_lj;
427 float cycles_ppdpme, cycles_seppme;
429 cycles_ppdpme = wallcycle_stop(wcycle, ewcPPDURINGPME);
430 dd_cycles_add(cr->dd, cycles_ppdpme, ddCyclPPduringPME);
432 /* In case of node-splitting, the PP nodes receive the long-range
433 * forces, virial and energy from the PME nodes here.
435 wallcycle_start(wcycle, ewcPP_PMEWAITRECVF);
438 gmx_pme_receive_f(cr, fr->f_novirsum, fr->vir_el_recip, &e_q,
439 fr->vir_lj_recip, &e_lj, &dvdl_q, &dvdl_lj,
443 gmx_print_sepdvdl(fplog, "Electrostatic PME mesh", e_q, dvdl_q);
444 gmx_print_sepdvdl(fplog, "Lennard-Jones PME mesh", e_lj, dvdl_lj);
446 enerd->term[F_COUL_RECIP] += e_q;
447 enerd->term[F_LJ_RECIP] += e_lj;
448 enerd->dvdl_lin[efptCOUL] += dvdl_q;
449 enerd->dvdl_lin[efptVDW] += dvdl_lj;
453 dd_cycles_add(cr->dd, cycles_seppme, ddCyclPME);
455 wallcycle_stop(wcycle, ewcPP_PMEWAITRECVF);
458 static void print_large_forces(FILE *fp, t_mdatoms *md, t_commrec *cr,
459 gmx_int64_t step, real pforce, rvec *x, rvec *f)
463 char buf[STEPSTRSIZE];
466 for (i = md->start; i < md->start+md->homenr; i++)
469 /* We also catch NAN, if the compiler does not optimize this away. */
470 if (fn2 >= pf2 || fn2 != fn2)
472 fprintf(fp, "step %s atom %6d x %8.3f %8.3f %8.3f force %12.5e\n",
473 gmx_step_str(step, buf),
474 ddglatnr(cr->dd, i), x[i][XX], x[i][YY], x[i][ZZ], sqrt(fn2));
479 static void post_process_forces(t_commrec *cr,
481 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
483 matrix box, rvec x[],
488 t_forcerec *fr, gmx_vsite_t *vsite,
495 /* Spread the mesh force on virtual sites to the other particles...
496 * This is parallellized. MPI communication is performed
497 * if the constructing atoms aren't local.
499 wallcycle_start(wcycle, ewcVSITESPREAD);
500 spread_vsite_f(vsite, x, fr->f_novirsum, NULL,
501 (flags & GMX_FORCE_VIRIAL), fr->vir_el_recip,
503 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
504 wallcycle_stop(wcycle, ewcVSITESPREAD);
506 if (flags & GMX_FORCE_VIRIAL)
508 /* Now add the forces, this is local */
511 sum_forces(0, fr->f_novirsum_n, f, fr->f_novirsum);
515 sum_forces(mdatoms->start, mdatoms->start+mdatoms->homenr,
518 if (EEL_FULL(fr->eeltype))
520 /* Add the mesh contribution to the virial */
521 m_add(vir_force, fr->vir_el_recip, vir_force);
523 if (EVDW_PME(fr->vdwtype))
525 /* Add the mesh contribution to the virial */
526 m_add(vir_force, fr->vir_lj_recip, vir_force);
530 pr_rvecs(debug, 0, "vir_force", vir_force, DIM);
535 if (fr->print_force >= 0)
537 print_large_forces(stderr, mdatoms, cr, step, fr->print_force, x, f);
541 static void do_nb_verlet(t_forcerec *fr,
542 interaction_const_t *ic,
543 gmx_enerdata_t *enerd,
544 int flags, int ilocality,
547 gmx_wallcycle_t wcycle)
549 int nnbl, kernel_type, enr_nbnxn_kernel_ljc, enr_nbnxn_kernel_lj;
551 nonbonded_verlet_group_t *nbvg;
554 if (!(flags & GMX_FORCE_NONBONDED))
556 /* skip non-bonded calculation */
560 nbvg = &fr->nbv->grp[ilocality];
562 /* CUDA kernel launch overhead is already timed separately */
563 if (fr->cutoff_scheme != ecutsVERLET)
565 gmx_incons("Invalid cut-off scheme passed!");
568 bCUDA = (nbvg->kernel_type == nbnxnk8x8x8_CUDA);
572 wallcycle_sub_start(wcycle, ewcsNONBONDED);
574 switch (nbvg->kernel_type)
576 case nbnxnk4x4_PlainC:
577 nbnxn_kernel_ref(&nbvg->nbl_lists,
583 enerd->grpp.ener[egCOULSR],
585 enerd->grpp.ener[egBHAMSR] :
586 enerd->grpp.ener[egLJSR]);
589 case nbnxnk4xN_SIMD_4xN:
590 nbnxn_kernel_simd_4xn(&nbvg->nbl_lists,
597 enerd->grpp.ener[egCOULSR],
599 enerd->grpp.ener[egBHAMSR] :
600 enerd->grpp.ener[egLJSR]);
602 case nbnxnk4xN_SIMD_2xNN:
603 nbnxn_kernel_simd_2xnn(&nbvg->nbl_lists,
610 enerd->grpp.ener[egCOULSR],
612 enerd->grpp.ener[egBHAMSR] :
613 enerd->grpp.ener[egLJSR]);
616 case nbnxnk8x8x8_CUDA:
617 nbnxn_cuda_launch_kernel(fr->nbv->cu_nbv, nbvg->nbat, flags, ilocality);
620 case nbnxnk8x8x8_PlainC:
621 nbnxn_kernel_gpu_ref(nbvg->nbl_lists.nbl[0],
626 nbvg->nbat->out[0].f,
628 enerd->grpp.ener[egCOULSR],
630 enerd->grpp.ener[egBHAMSR] :
631 enerd->grpp.ener[egLJSR]);
635 gmx_incons("Invalid nonbonded kernel type passed!");
640 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
643 if (EEL_RF(ic->eeltype) || ic->eeltype == eelCUT)
645 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_RF;
647 else if ((!bCUDA && nbvg->ewald_excl == ewaldexclAnalytical) ||
648 (bCUDA && nbnxn_cuda_is_kernel_ewald_analytical(fr->nbv->cu_nbv)))
650 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_EWALD;
654 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_TAB;
656 enr_nbnxn_kernel_lj = eNR_NBNXN_LJ;
657 if (flags & GMX_FORCE_ENERGY)
659 /* In eNR_??? the nbnxn F+E kernels are always the F kernel + 1 */
660 enr_nbnxn_kernel_ljc += 1;
661 enr_nbnxn_kernel_lj += 1;
664 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc,
665 nbvg->nbl_lists.natpair_ljq);
666 inc_nrnb(nrnb, enr_nbnxn_kernel_lj,
667 nbvg->nbl_lists.natpair_lj);
668 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc-eNR_NBNXN_LJ_RF+eNR_NBNXN_RF,
669 nbvg->nbl_lists.natpair_q);
672 void do_force_cutsVERLET(FILE *fplog, t_commrec *cr,
673 t_inputrec *inputrec,
674 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
676 gmx_groups_t gmx_unused *groups,
677 matrix box, rvec x[], history_t *hist,
681 gmx_enerdata_t *enerd, t_fcdata *fcd,
682 real *lambda, t_graph *graph,
683 t_forcerec *fr, interaction_const_t *ic,
684 gmx_vsite_t *vsite, rvec mu_tot,
685 double t, FILE *field, gmx_edsam_t ed,
693 gmx_bool bSepDVDL, bStateChanged, bNS, bFillGrid, bCalcCGCM, bBS;
694 gmx_bool bDoLongRange, bDoForces, bSepLRF, bUseGPU, bUseOrEmulGPU;
695 gmx_bool bDiffKernels = FALSE;
697 rvec vzero, box_diag;
699 float cycles_pme, cycles_force, cycles_wait_gpu;
700 nonbonded_verlet_t *nbv;
705 nb_kernel_type = fr->nbv->grp[0].kernel_type;
707 start = mdatoms->start;
708 homenr = mdatoms->homenr;
710 bSepDVDL = (fr->bSepDVDL && do_per_step(step, inputrec->nstlog));
712 clear_mat(vir_force);
715 if (DOMAINDECOMP(cr))
717 cg1 = cr->dd->ncg_tot;
728 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
729 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
730 bFillGrid = (bNS && bStateChanged);
731 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
732 bDoLongRange = (fr->bTwinRange && bNS && (flags & GMX_FORCE_DO_LR));
733 bDoForces = (flags & GMX_FORCE_FORCES);
734 bSepLRF = (bDoLongRange && bDoForces && (flags & GMX_FORCE_SEPLRF));
735 bUseGPU = fr->nbv->bUseGPU;
736 bUseOrEmulGPU = bUseGPU || (nbv->grp[0].kernel_type == nbnxnk8x8x8_PlainC);
740 update_forcerec(fr, box);
742 if (NEED_MUTOT(*inputrec))
744 /* Calculate total (local) dipole moment in a temporary common array.
745 * This makes it possible to sum them over nodes faster.
747 calc_mu(start, homenr,
748 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
753 if (fr->ePBC != epbcNONE)
755 /* Compute shift vectors every step,
756 * because of pressure coupling or box deformation!
758 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
760 calc_shifts(box, fr->shift_vec);
765 put_atoms_in_box_omp(fr->ePBC, box, homenr, x);
766 inc_nrnb(nrnb, eNR_SHIFTX, homenr);
768 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
770 unshift_self(graph, box, x);
774 nbnxn_atomdata_copy_shiftvec(flags & GMX_FORCE_DYNAMICBOX,
775 fr->shift_vec, nbv->grp[0].nbat);
778 if (!(cr->duty & DUTY_PME))
780 /* Send particle coordinates to the pme nodes.
781 * Since this is only implemented for domain decomposition
782 * and domain decomposition does not use the graph,
783 * we do not need to worry about shifting.
788 wallcycle_start(wcycle, ewcPP_PMESENDX);
790 bBS = (inputrec->nwall == 2);
794 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
797 if (EEL_PME(fr->eeltype))
799 pme_flags |= GMX_PME_DO_COULOMB;
802 if (EVDW_PME(fr->vdwtype))
804 pme_flags |= GMX_PME_DO_LJ;
805 if (fr->ljpme_combination_rule == eljpmeLB)
807 pme_flags |= GMX_PME_LJ_LB;
811 gmx_pme_send_coordinates(cr, bBS ? boxs : box, x,
812 mdatoms->nChargePerturbed, mdatoms->nTypePerturbed, lambda[efptCOUL], lambda[efptVDW],
813 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),
816 wallcycle_stop(wcycle, ewcPP_PMESENDX);
820 /* do gridding for pair search */
823 if (graph && bStateChanged)
825 /* Calculate intramolecular shift vectors to make molecules whole */
826 mk_mshift(fplog, graph, fr->ePBC, box, x);
830 box_diag[XX] = box[XX][XX];
831 box_diag[YY] = box[YY][YY];
832 box_diag[ZZ] = box[ZZ][ZZ];
834 wallcycle_start(wcycle, ewcNS);
837 wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
838 nbnxn_put_on_grid(nbv->nbs, fr->ePBC, box,
840 0, mdatoms->homenr, -1, fr->cginfo, x,
842 nbv->grp[eintLocal].kernel_type,
843 nbv->grp[eintLocal].nbat);
844 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
848 wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
849 nbnxn_put_on_grid_nonlocal(nbv->nbs, domdec_zones(cr->dd),
851 nbv->grp[eintNonlocal].kernel_type,
852 nbv->grp[eintNonlocal].nbat);
853 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
856 if (nbv->ngrp == 1 ||
857 nbv->grp[eintNonlocal].nbat == nbv->grp[eintLocal].nbat)
859 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatAll,
860 nbv->nbs, mdatoms, fr->cginfo);
864 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatLocal,
865 nbv->nbs, mdatoms, fr->cginfo);
866 nbnxn_atomdata_set(nbv->grp[eintNonlocal].nbat, eatAll,
867 nbv->nbs, mdatoms, fr->cginfo);
869 wallcycle_stop(wcycle, ewcNS);
872 /* initialize the GPU atom data and copy shift vector */
877 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
878 nbnxn_cuda_init_atomdata(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
879 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
882 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
883 nbnxn_cuda_upload_shiftvec(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
884 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
887 /* do local pair search */
890 wallcycle_start_nocount(wcycle, ewcNS);
891 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
892 nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintLocal].nbat,
895 nbv->min_ci_balanced,
896 &nbv->grp[eintLocal].nbl_lists,
898 nbv->grp[eintLocal].kernel_type,
900 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
904 /* initialize local pair-list on the GPU */
905 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
906 nbv->grp[eintLocal].nbl_lists.nbl[0],
909 wallcycle_stop(wcycle, ewcNS);
913 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
914 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
915 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, FALSE, x,
916 nbv->grp[eintLocal].nbat);
917 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
918 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
923 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
924 /* launch local nonbonded F on GPU */
925 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFNo,
927 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
930 /* Communicate coordinates and sum dipole if necessary +
931 do non-local pair search */
932 if (DOMAINDECOMP(cr))
934 bDiffKernels = (nbv->grp[eintNonlocal].kernel_type !=
935 nbv->grp[eintLocal].kernel_type);
939 /* With GPU+CPU non-bonded calculations we need to copy
940 * the local coordinates to the non-local nbat struct
941 * (in CPU format) as the non-local kernel call also
942 * calculates the local - non-local interactions.
944 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
945 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
946 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, TRUE, x,
947 nbv->grp[eintNonlocal].nbat);
948 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
949 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
954 wallcycle_start_nocount(wcycle, ewcNS);
955 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
959 nbnxn_grid_add_simple(nbv->nbs, nbv->grp[eintNonlocal].nbat);
962 nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintNonlocal].nbat,
965 nbv->min_ci_balanced,
966 &nbv->grp[eintNonlocal].nbl_lists,
968 nbv->grp[eintNonlocal].kernel_type,
971 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
973 if (nbv->grp[eintNonlocal].kernel_type == nbnxnk8x8x8_CUDA)
975 /* initialize non-local pair-list on the GPU */
976 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
977 nbv->grp[eintNonlocal].nbl_lists.nbl[0],
980 wallcycle_stop(wcycle, ewcNS);
984 wallcycle_start(wcycle, ewcMOVEX);
985 dd_move_x(cr->dd, box, x);
987 /* When we don't need the total dipole we sum it in global_stat */
988 if (bStateChanged && NEED_MUTOT(*inputrec))
990 gmx_sumd(2*DIM, mu, cr);
992 wallcycle_stop(wcycle, ewcMOVEX);
994 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
995 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
996 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatNonlocal, FALSE, x,
997 nbv->grp[eintNonlocal].nbat);
998 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
999 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1002 if (bUseGPU && !bDiffKernels)
1004 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
1005 /* launch non-local nonbonded F on GPU */
1006 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFNo,
1008 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1014 /* launch D2H copy-back F */
1015 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1016 if (DOMAINDECOMP(cr) && !bDiffKernels)
1018 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintNonlocal].nbat,
1019 flags, eatNonlocal);
1021 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintLocal].nbat,
1023 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1026 if (bStateChanged && NEED_MUTOT(*inputrec))
1030 gmx_sumd(2*DIM, mu, cr);
1033 for (i = 0; i < 2; i++)
1035 for (j = 0; j < DIM; j++)
1037 fr->mu_tot[i][j] = mu[i*DIM + j];
1041 if (fr->efep == efepNO)
1043 copy_rvec(fr->mu_tot[0], mu_tot);
1047 for (j = 0; j < DIM; j++)
1050 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] +
1051 lambda[efptCOUL]*fr->mu_tot[1][j];
1055 /* Reset energies */
1056 reset_enerdata(fr, bNS, enerd, MASTER(cr));
1057 clear_rvecs(SHIFTS, fr->fshift);
1059 if (DOMAINDECOMP(cr))
1061 if (!(cr->duty & DUTY_PME))
1063 wallcycle_start(wcycle, ewcPPDURINGPME);
1064 dd_force_flop_start(cr->dd, nrnb);
1070 /* Enforced rotation has its own cycle counter that starts after the collective
1071 * coordinates have been communicated. It is added to ddCyclF to allow
1072 * for proper load-balancing */
1073 wallcycle_start(wcycle, ewcROT);
1074 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1075 wallcycle_stop(wcycle, ewcROT);
1078 /* Start the force cycle counter.
1079 * This counter is stopped in do_forcelow_level.
1080 * No parallel communication should occur while this counter is running,
1081 * since that will interfere with the dynamic load balancing.
1083 wallcycle_start(wcycle, ewcFORCE);
1086 /* Reset forces for which the virial is calculated separately:
1087 * PME/Ewald forces if necessary */
1088 if (fr->bF_NoVirSum)
1090 if (flags & GMX_FORCE_VIRIAL)
1092 fr->f_novirsum = fr->f_novirsum_alloc;
1095 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1099 clear_rvecs(homenr, fr->f_novirsum+start);
1104 /* We are not calculating the pressure so we do not need
1105 * a separate array for forces that do not contribute
1112 /* Clear the short- and long-range forces */
1113 clear_rvecs(fr->natoms_force_constr, f);
1114 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1116 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1119 clear_rvec(fr->vir_diag_posres);
1122 if (inputrec->ePull == epullCONSTRAINT)
1124 clear_pull_forces(inputrec->pull);
1127 /* We calculate the non-bonded forces, when done on the CPU, here.
1128 * We do this before calling do_force_lowlevel, as in there bondeds
1129 * forces are calculated before PME, which does communication.
1130 * With this order, non-bonded and bonded force calculation imbalance
1131 * can be balanced out by the domain decomposition load balancing.
1136 /* Maybe we should move this into do_force_lowlevel */
1137 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFYes,
1141 if (!bUseOrEmulGPU || bDiffKernels)
1145 if (DOMAINDECOMP(cr))
1147 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal,
1148 bDiffKernels ? enbvClearFYes : enbvClearFNo,
1158 aloc = eintNonlocal;
1161 /* Add all the non-bonded force to the normal force array.
1162 * This can be split into a local a non-local part when overlapping
1163 * communication with calculation with domain decomposition.
1165 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1166 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1167 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1168 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatAll, nbv->grp[aloc].nbat, f);
1169 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1170 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1171 wallcycle_start_nocount(wcycle, ewcFORCE);
1173 /* if there are multiple fshift output buffers reduce them */
1174 if ((flags & GMX_FORCE_VIRIAL) &&
1175 nbv->grp[aloc].nbl_lists.nnbl > 1)
1177 nbnxn_atomdata_add_nbat_fshift_to_fshift(nbv->grp[aloc].nbat,
1182 /* update QMMMrec, if necessary */
1185 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1188 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1190 posres_wrapper(fplog, flags, bSepDVDL, inputrec, nrnb, top, box, x,
1194 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_FBPOSRES].nr > 0)
1196 fbposres_wrapper(inputrec, nrnb, top, box, x, enerd, fr);
1199 /* Compute the bonded and non-bonded energies and optionally forces */
1200 do_force_lowlevel(fplog, step, fr, inputrec, &(top->idef),
1201 cr, nrnb, wcycle, mdatoms,
1202 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, top, fr->born,
1203 &(top->atomtypes), bBornRadii, box,
1204 inputrec->fepvals, lambda, graph, &(top->excls), fr->mu_tot,
1205 flags, &cycles_pme);
1209 if (do_per_step(step, inputrec->nstcalclr))
1211 /* Add the long range forces to the short range forces */
1212 for (i = 0; i < fr->natoms_force_constr; i++)
1214 rvec_add(fr->f_twin[i], f[i], f[i]);
1219 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1223 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1226 if (bUseOrEmulGPU && !bDiffKernels)
1228 /* wait for non-local forces (or calculate in emulation mode) */
1229 if (DOMAINDECOMP(cr))
1235 wallcycle_start(wcycle, ewcWAIT_GPU_NB_NL);
1236 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1237 nbv->grp[eintNonlocal].nbat,
1239 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1241 cycles_tmp = wallcycle_stop(wcycle, ewcWAIT_GPU_NB_NL);
1242 cycles_wait_gpu += cycles_tmp;
1243 cycles_force += cycles_tmp;
1247 wallcycle_start_nocount(wcycle, ewcFORCE);
1248 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFYes,
1250 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1252 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1253 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1254 /* skip the reduction if there was no non-local work to do */
1255 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1257 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatNonlocal,
1258 nbv->grp[eintNonlocal].nbat, f);
1260 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1261 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1267 /* Communicate the forces */
1270 wallcycle_start(wcycle, ewcMOVEF);
1271 if (DOMAINDECOMP(cr))
1273 dd_move_f(cr->dd, f, fr->fshift);
1274 /* Do we need to communicate the separate force array
1275 * for terms that do not contribute to the single sum virial?
1276 * Position restraints and electric fields do not introduce
1277 * inter-cg forces, only full electrostatics methods do.
1278 * When we do not calculate the virial, fr->f_novirsum = f,
1279 * so we have already communicated these forces.
1281 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1282 (flags & GMX_FORCE_VIRIAL))
1284 dd_move_f(cr->dd, fr->f_novirsum, NULL);
1288 /* We should not update the shift forces here,
1289 * since f_twin is already included in f.
1291 dd_move_f(cr->dd, fr->f_twin, NULL);
1294 wallcycle_stop(wcycle, ewcMOVEF);
1300 /* wait for local forces (or calculate in emulation mode) */
1303 wallcycle_start(wcycle, ewcWAIT_GPU_NB_L);
1304 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1305 nbv->grp[eintLocal].nbat,
1307 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1309 cycles_wait_gpu += wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
1311 /* now clear the GPU outputs while we finish the step on the CPU */
1313 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1314 nbnxn_cuda_clear_outputs(nbv->cu_nbv, flags);
1315 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1319 wallcycle_start_nocount(wcycle, ewcFORCE);
1320 do_nb_verlet(fr, ic, enerd, flags, eintLocal,
1321 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
1323 wallcycle_stop(wcycle, ewcFORCE);
1325 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1326 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1327 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1329 /* skip the reduction if there was no non-local work to do */
1330 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatLocal,
1331 nbv->grp[eintLocal].nbat, f);
1333 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1334 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1337 if (DOMAINDECOMP(cr))
1339 dd_force_flop_stop(cr->dd, nrnb);
1342 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1345 dd_cycles_add(cr->dd, cycles_wait_gpu, ddCyclWaitGPU);
1352 if (IR_ELEC_FIELD(*inputrec))
1354 /* Compute forces due to electric field */
1355 calc_f_el(MASTER(cr) ? field : NULL,
1356 start, homenr, mdatoms->chargeA, fr->f_novirsum,
1357 inputrec->ex, inputrec->et, t);
1360 /* If we have NoVirSum forces, but we do not calculate the virial,
1361 * we sum fr->f_novirum=f later.
1363 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1365 wallcycle_start(wcycle, ewcVSITESPREAD);
1366 spread_vsite_f(vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1367 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1368 wallcycle_stop(wcycle, ewcVSITESPREAD);
1372 wallcycle_start(wcycle, ewcVSITESPREAD);
1373 spread_vsite_f(vsite, x, fr->f_twin, NULL, FALSE, NULL,
1375 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1376 wallcycle_stop(wcycle, ewcVSITESPREAD);
1380 if (flags & GMX_FORCE_VIRIAL)
1382 /* Calculation of the virial must be done after vsites! */
1383 calc_virial(mdatoms->start, mdatoms->homenr, x, f,
1384 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1388 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1390 pull_potential_wrapper(fplog, bSepDVDL, cr, inputrec, box, x,
1391 f, vir_force, mdatoms, enerd, lambda, t);
1394 /* Add the forces from enforced rotation potentials (if any) */
1397 wallcycle_start(wcycle, ewcROTadd);
1398 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
1399 wallcycle_stop(wcycle, ewcROTadd);
1402 if (PAR(cr) && !(cr->duty & DUTY_PME))
1404 /* In case of node-splitting, the PP nodes receive the long-range
1405 * forces, virial and energy from the PME nodes here.
1407 pme_receive_force_ener(fplog, bSepDVDL, cr, wcycle, enerd, fr);
1412 post_process_forces(cr, step, nrnb, wcycle,
1413 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
1417 /* Sum the potential energy terms from group contributions */
1418 sum_epot(&(enerd->grpp), enerd->term);
1421 void do_force_cutsGROUP(FILE *fplog, t_commrec *cr,
1422 t_inputrec *inputrec,
1423 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1424 gmx_localtop_t *top,
1425 gmx_groups_t *groups,
1426 matrix box, rvec x[], history_t *hist,
1430 gmx_enerdata_t *enerd, t_fcdata *fcd,
1431 real *lambda, t_graph *graph,
1432 t_forcerec *fr, gmx_vsite_t *vsite, rvec mu_tot,
1433 double t, FILE *field, gmx_edsam_t ed,
1434 gmx_bool bBornRadii,
1440 gmx_bool bSepDVDL, bStateChanged, bNS, bFillGrid, bCalcCGCM, bBS;
1441 gmx_bool bDoLongRangeNS, bDoForces, bDoPotential, bSepLRF;
1442 gmx_bool bDoAdressWF;
1444 rvec vzero, box_diag;
1445 real e, v, dvdlambda[efptNR];
1447 float cycles_pme, cycles_force;
1449 start = mdatoms->start;
1450 homenr = mdatoms->homenr;
1452 bSepDVDL = (fr->bSepDVDL && do_per_step(step, inputrec->nstlog));
1454 clear_mat(vir_force);
1458 pd_cg_range(cr, &cg0, &cg1);
1463 if (DOMAINDECOMP(cr))
1465 cg1 = cr->dd->ncg_tot;
1477 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
1478 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
1479 /* Should we update the long-range neighborlists at this step? */
1480 bDoLongRangeNS = fr->bTwinRange && bNS;
1481 /* Should we perform the long-range nonbonded evaluation inside the neighborsearching? */
1482 bFillGrid = (bNS && bStateChanged);
1483 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
1484 bDoForces = (flags & GMX_FORCE_FORCES);
1485 bDoPotential = (flags & GMX_FORCE_ENERGY);
1486 bSepLRF = ((inputrec->nstcalclr > 1) && bDoForces &&
1487 (flags & GMX_FORCE_SEPLRF) && (flags & GMX_FORCE_DO_LR));
1489 /* should probably move this to the forcerec since it doesn't change */
1490 bDoAdressWF = ((fr->adress_type != eAdressOff));
1494 update_forcerec(fr, box);
1496 if (NEED_MUTOT(*inputrec))
1498 /* Calculate total (local) dipole moment in a temporary common array.
1499 * This makes it possible to sum them over nodes faster.
1501 calc_mu(start, homenr,
1502 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
1507 if (fr->ePBC != epbcNONE)
1509 /* Compute shift vectors every step,
1510 * because of pressure coupling or box deformation!
1512 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
1514 calc_shifts(box, fr->shift_vec);
1519 put_charge_groups_in_box(fplog, cg0, cg1, fr->ePBC, box,
1520 &(top->cgs), x, fr->cg_cm);
1521 inc_nrnb(nrnb, eNR_CGCM, homenr);
1522 inc_nrnb(nrnb, eNR_RESETX, cg1-cg0);
1524 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
1526 unshift_self(graph, box, x);
1531 calc_cgcm(fplog, cg0, cg1, &(top->cgs), x, fr->cg_cm);
1532 inc_nrnb(nrnb, eNR_CGCM, homenr);
1539 move_cgcm(fplog, cr, fr->cg_cm);
1543 pr_rvecs(debug, 0, "cgcm", fr->cg_cm, top->cgs.nr);
1548 if (!(cr->duty & DUTY_PME))
1550 /* Send particle coordinates to the pme nodes.
1551 * Since this is only implemented for domain decomposition
1552 * and domain decomposition does not use the graph,
1553 * we do not need to worry about shifting.
1558 wallcycle_start(wcycle, ewcPP_PMESENDX);
1560 bBS = (inputrec->nwall == 2);
1563 copy_mat(box, boxs);
1564 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
1567 if (EEL_PME(fr->eeltype))
1569 pme_flags |= GMX_PME_DO_COULOMB;
1572 if (EVDW_PME(fr->vdwtype))
1574 pme_flags |= GMX_PME_DO_LJ;
1575 if (fr->ljpme_combination_rule == eljpmeLB)
1577 pme_flags |= GMX_PME_LJ_LB;
1581 gmx_pme_send_coordinates(cr, bBS ? boxs : box, x,
1582 mdatoms->nChargePerturbed, mdatoms->nTypePerturbed, lambda[efptCOUL], lambda[efptVDW],
1583 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),
1586 wallcycle_stop(wcycle, ewcPP_PMESENDX);
1588 #endif /* GMX_MPI */
1590 /* Communicate coordinates and sum dipole if necessary */
1593 wallcycle_start(wcycle, ewcMOVEX);
1594 if (DOMAINDECOMP(cr))
1596 dd_move_x(cr->dd, box, x);
1600 move_x(cr, x, nrnb);
1602 wallcycle_stop(wcycle, ewcMOVEX);
1605 /* update adress weight beforehand */
1606 if (bStateChanged && bDoAdressWF)
1608 /* need pbc for adress weight calculation with pbc_dx */
1609 set_pbc(&pbc, inputrec->ePBC, box);
1610 if (fr->adress_site == eAdressSITEcog)
1612 update_adress_weights_cog(top->idef.iparams, top->idef.il, x, fr, mdatoms,
1613 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1615 else if (fr->adress_site == eAdressSITEcom)
1617 update_adress_weights_com(fplog, cg0, cg1, &(top->cgs), x, fr, mdatoms,
1618 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1620 else if (fr->adress_site == eAdressSITEatomatom)
1622 update_adress_weights_atom_per_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1623 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1627 update_adress_weights_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1628 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1632 if (NEED_MUTOT(*inputrec))
1639 gmx_sumd(2*DIM, mu, cr);
1641 for (i = 0; i < 2; i++)
1643 for (j = 0; j < DIM; j++)
1645 fr->mu_tot[i][j] = mu[i*DIM + j];
1649 if (fr->efep == efepNO)
1651 copy_rvec(fr->mu_tot[0], mu_tot);
1655 for (j = 0; j < DIM; j++)
1658 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] + lambda[efptCOUL]*fr->mu_tot[1][j];
1663 /* Reset energies */
1664 reset_enerdata(fr, bNS, enerd, MASTER(cr));
1665 clear_rvecs(SHIFTS, fr->fshift);
1669 wallcycle_start(wcycle, ewcNS);
1671 if (graph && bStateChanged)
1673 /* Calculate intramolecular shift vectors to make molecules whole */
1674 mk_mshift(fplog, graph, fr->ePBC, box, x);
1677 /* Do the actual neighbour searching */
1679 groups, top, mdatoms,
1680 cr, nrnb, bFillGrid,
1683 wallcycle_stop(wcycle, ewcNS);
1686 if (inputrec->implicit_solvent && bNS)
1688 make_gb_nblist(cr, inputrec->gb_algorithm,
1689 x, box, fr, &top->idef, graph, fr->born);
1692 if (DOMAINDECOMP(cr))
1694 if (!(cr->duty & DUTY_PME))
1696 wallcycle_start(wcycle, ewcPPDURINGPME);
1697 dd_force_flop_start(cr->dd, nrnb);
1703 /* Enforced rotation has its own cycle counter that starts after the collective
1704 * coordinates have been communicated. It is added to ddCyclF to allow
1705 * for proper load-balancing */
1706 wallcycle_start(wcycle, ewcROT);
1707 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1708 wallcycle_stop(wcycle, ewcROT);
1711 /* Start the force cycle counter.
1712 * This counter is stopped in do_forcelow_level.
1713 * No parallel communication should occur while this counter is running,
1714 * since that will interfere with the dynamic load balancing.
1716 wallcycle_start(wcycle, ewcFORCE);
1720 /* Reset forces for which the virial is calculated separately:
1721 * PME/Ewald forces if necessary */
1722 if (fr->bF_NoVirSum)
1724 if (flags & GMX_FORCE_VIRIAL)
1726 fr->f_novirsum = fr->f_novirsum_alloc;
1729 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1733 clear_rvecs(homenr, fr->f_novirsum+start);
1738 /* We are not calculating the pressure so we do not need
1739 * a separate array for forces that do not contribute
1746 /* Clear the short- and long-range forces */
1747 clear_rvecs(fr->natoms_force_constr, f);
1748 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1750 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1753 clear_rvec(fr->vir_diag_posres);
1755 if (inputrec->ePull == epullCONSTRAINT)
1757 clear_pull_forces(inputrec->pull);
1760 /* update QMMMrec, if necessary */
1763 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1766 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1768 posres_wrapper(fplog, flags, bSepDVDL, inputrec, nrnb, top, box, x,
1772 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_FBPOSRES].nr > 0)
1774 fbposres_wrapper(inputrec, nrnb, top, box, x, enerd, fr);
1777 /* Compute the bonded and non-bonded energies and optionally forces */
1778 do_force_lowlevel(fplog, step, fr, inputrec, &(top->idef),
1779 cr, nrnb, wcycle, mdatoms,
1780 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, top, fr->born,
1781 &(top->atomtypes), bBornRadii, box,
1782 inputrec->fepvals, lambda,
1783 graph, &(top->excls), fr->mu_tot,
1789 if (do_per_step(step, inputrec->nstcalclr))
1791 /* Add the long range forces to the short range forces */
1792 for (i = 0; i < fr->natoms_force_constr; i++)
1794 rvec_add(fr->f_twin[i], f[i], f[i]);
1799 cycles_force = wallcycle_stop(wcycle, ewcFORCE);
1803 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1806 if (DOMAINDECOMP(cr))
1808 dd_force_flop_stop(cr->dd, nrnb);
1811 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1817 if (IR_ELEC_FIELD(*inputrec))
1819 /* Compute forces due to electric field */
1820 calc_f_el(MASTER(cr) ? field : NULL,
1821 start, homenr, mdatoms->chargeA, fr->f_novirsum,
1822 inputrec->ex, inputrec->et, t);
1825 if (bDoAdressWF && fr->adress_icor == eAdressICThermoForce)
1827 /* Compute thermodynamic force in hybrid AdResS region */
1828 adress_thermo_force(start, homenr, &(top->cgs), x, fr->f_novirsum, fr, mdatoms,
1829 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1832 /* Communicate the forces */
1835 wallcycle_start(wcycle, ewcMOVEF);
1836 if (DOMAINDECOMP(cr))
1838 dd_move_f(cr->dd, f, fr->fshift);
1839 /* Do we need to communicate the separate force array
1840 * for terms that do not contribute to the single sum virial?
1841 * Position restraints and electric fields do not introduce
1842 * inter-cg forces, only full electrostatics methods do.
1843 * When we do not calculate the virial, fr->f_novirsum = f,
1844 * so we have already communicated these forces.
1846 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1847 (flags & GMX_FORCE_VIRIAL))
1849 dd_move_f(cr->dd, fr->f_novirsum, NULL);
1853 /* We should not update the shift forces here,
1854 * since f_twin is already included in f.
1856 dd_move_f(cr->dd, fr->f_twin, NULL);
1861 pd_move_f(cr, f, nrnb);
1864 pd_move_f(cr, fr->f_twin, nrnb);
1867 wallcycle_stop(wcycle, ewcMOVEF);
1870 /* If we have NoVirSum forces, but we do not calculate the virial,
1871 * we sum fr->f_novirum=f later.
1873 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1875 wallcycle_start(wcycle, ewcVSITESPREAD);
1876 spread_vsite_f(vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1877 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1878 wallcycle_stop(wcycle, ewcVSITESPREAD);
1882 wallcycle_start(wcycle, ewcVSITESPREAD);
1883 spread_vsite_f(vsite, x, fr->f_twin, NULL, FALSE, NULL,
1885 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1886 wallcycle_stop(wcycle, ewcVSITESPREAD);
1890 if (flags & GMX_FORCE_VIRIAL)
1892 /* Calculation of the virial must be done after vsites! */
1893 calc_virial(mdatoms->start, mdatoms->homenr, x, f,
1894 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1898 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1900 pull_potential_wrapper(fplog, bSepDVDL, cr, inputrec, box, x,
1901 f, vir_force, mdatoms, enerd, lambda, t);
1904 /* Add the forces from enforced rotation potentials (if any) */
1907 wallcycle_start(wcycle, ewcROTadd);
1908 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
1909 wallcycle_stop(wcycle, ewcROTadd);
1912 if (PAR(cr) && !(cr->duty & DUTY_PME))
1914 /* In case of node-splitting, the PP nodes receive the long-range
1915 * forces, virial and energy from the PME nodes here.
1917 pme_receive_force_ener(fplog, bSepDVDL, cr, wcycle, enerd, fr);
1922 post_process_forces(cr, step, nrnb, wcycle,
1923 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
1927 /* Sum the potential energy terms from group contributions */
1928 sum_epot(&(enerd->grpp), enerd->term);
1931 void do_force(FILE *fplog, t_commrec *cr,
1932 t_inputrec *inputrec,
1933 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1934 gmx_localtop_t *top,
1935 gmx_groups_t *groups,
1936 matrix box, rvec x[], history_t *hist,
1940 gmx_enerdata_t *enerd, t_fcdata *fcd,
1941 real *lambda, t_graph *graph,
1943 gmx_vsite_t *vsite, rvec mu_tot,
1944 double t, FILE *field, gmx_edsam_t ed,
1945 gmx_bool bBornRadii,
1948 /* modify force flag if not doing nonbonded */
1949 if (!fr->bNonbonded)
1951 flags &= ~GMX_FORCE_NONBONDED;
1954 switch (inputrec->cutoff_scheme)
1957 do_force_cutsVERLET(fplog, cr, inputrec,
1973 do_force_cutsGROUP(fplog, cr, inputrec,
1988 gmx_incons("Invalid cut-off scheme passed!");
1993 void do_constrain_first(FILE *fplog, gmx_constr_t constr,
1994 t_inputrec *ir, t_mdatoms *md,
1995 t_state *state, t_commrec *cr, t_nrnb *nrnb,
1996 t_forcerec *fr, gmx_localtop_t *top)
1998 int i, m, start, end;
2000 real dt = ir->delta_t;
2004 snew(savex, state->natoms);
2007 end = md->homenr + start;
2011 fprintf(debug, "vcm: start=%d, homenr=%d, end=%d\n",
2012 start, md->homenr, end);
2014 /* Do a first constrain to reset particles... */
2015 step = ir->init_step;
2018 char buf[STEPSTRSIZE];
2019 fprintf(fplog, "\nConstraining the starting coordinates (step %s)\n",
2020 gmx_step_str(step, buf));
2024 /* constrain the current position */
2025 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2026 ir, NULL, cr, step, 0, md,
2027 state->x, state->x, NULL,
2028 fr->bMolPBC, state->box,
2029 state->lambda[efptBONDED], &dvdl_dum,
2030 NULL, NULL, nrnb, econqCoord,
2031 ir->epc == epcMTTK, state->veta, state->veta);
2034 /* constrain the inital velocity, and save it */
2035 /* also may be useful if we need the ekin from the halfstep for velocity verlet */
2036 /* might not yet treat veta correctly */
2037 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2038 ir, NULL, cr, step, 0, md,
2039 state->x, state->v, state->v,
2040 fr->bMolPBC, state->box,
2041 state->lambda[efptBONDED], &dvdl_dum,
2042 NULL, NULL, nrnb, econqVeloc,
2043 ir->epc == epcMTTK, state->veta, state->veta);
2045 /* constrain the inital velocities at t-dt/2 */
2046 if (EI_STATE_VELOCITY(ir->eI) && ir->eI != eiVV)
2048 for (i = start; (i < end); i++)
2050 for (m = 0; (m < DIM); m++)
2052 /* Reverse the velocity */
2053 state->v[i][m] = -state->v[i][m];
2054 /* Store the position at t-dt in buf */
2055 savex[i][m] = state->x[i][m] + dt*state->v[i][m];
2058 /* Shake the positions at t=-dt with the positions at t=0
2059 * as reference coordinates.
2063 char buf[STEPSTRSIZE];
2064 fprintf(fplog, "\nConstraining the coordinates at t0-dt (step %s)\n",
2065 gmx_step_str(step, buf));
2068 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2069 ir, NULL, cr, step, -1, md,
2070 state->x, savex, NULL,
2071 fr->bMolPBC, state->box,
2072 state->lambda[efptBONDED], &dvdl_dum,
2073 state->v, NULL, nrnb, econqCoord,
2074 ir->epc == epcMTTK, state->veta, state->veta);
2076 for (i = start; i < end; i++)
2078 for (m = 0; m < DIM; m++)
2080 /* Re-reverse the velocities */
2081 state->v[i][m] = -state->v[i][m];
2090 integrate_table(real vdwtab[], real scale, int offstart, int rstart, int rend,
2091 double *enerout, double *virout)
2093 double enersum, virsum;
2094 double invscale, invscale2, invscale3;
2095 double r, ea, eb, ec, pa, pb, pc, pd;
2097 int ri, offset, tabfactor;
2099 invscale = 1.0/scale;
2100 invscale2 = invscale*invscale;
2101 invscale3 = invscale*invscale2;
2103 /* Following summation derived from cubic spline definition,
2104 * Numerical Recipies in C, second edition, p. 113-116. Exact for
2105 * the cubic spline. We first calculate the negative of the
2106 * energy from rvdw to rvdw_switch, assuming that g(r)=1, and then
2107 * add the more standard, abrupt cutoff correction to that result,
2108 * yielding the long-range correction for a switched function. We
2109 * perform both the pressure and energy loops at the same time for
2110 * simplicity, as the computational cost is low. */
2114 /* Since the dispersion table has been scaled down a factor
2115 * 6.0 and the repulsion a factor 12.0 to compensate for the
2116 * c6/c12 parameters inside nbfp[] being scaled up (to save
2117 * flops in kernels), we need to correct for this.
2128 for (ri = rstart; ri < rend; ++ri)
2132 eb = 2.0*invscale2*r;
2136 pb = 3.0*invscale2*r;
2137 pc = 3.0*invscale*r*r;
2140 /* this "8" is from the packing in the vdwtab array - perhaps
2141 should be defined? */
2143 offset = 8*ri + offstart;
2144 y0 = vdwtab[offset];
2145 f = vdwtab[offset+1];
2146 g = vdwtab[offset+2];
2147 h = vdwtab[offset+3];
2149 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);
2150 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);
2152 *enerout = 4.0*M_PI*enersum*tabfactor;
2153 *virout = 4.0*M_PI*virsum*tabfactor;
2156 void calc_enervirdiff(FILE *fplog, int eDispCorr, t_forcerec *fr)
2158 double eners[2], virs[2], enersum, virsum, y0, f, g, h;
2159 double r0, r1, r, rc3, rc9, ea, eb, ec, pa, pb, pc, pd;
2160 double invscale, invscale2, invscale3;
2161 int ri0, ri1, ri, i, offstart, offset;
2162 real scale, *vdwtab, tabfactor, tmp;
2164 fr->enershiftsix = 0;
2165 fr->enershifttwelve = 0;
2166 fr->enerdiffsix = 0;
2167 fr->enerdifftwelve = 0;
2169 fr->virdifftwelve = 0;
2171 if (eDispCorr != edispcNO)
2173 for (i = 0; i < 2; i++)
2178 if ((fr->vdwtype == evdwSWITCH) || (fr->vdwtype == evdwSHIFT))
2180 if (fr->rvdw_switch == 0)
2183 "With dispersion correction rvdw-switch can not be zero "
2184 "for vdw-type = %s", evdw_names[fr->vdwtype]);
2187 scale = fr->nblists[0].table_elec_vdw.scale;
2188 vdwtab = fr->nblists[0].table_vdw.data;
2190 /* Round the cut-offs to exact table values for precision */
2191 ri0 = floor(fr->rvdw_switch*scale);
2192 ri1 = ceil(fr->rvdw*scale);
2198 if (fr->vdwtype == evdwSHIFT)
2200 /* Determine the constant energy shift below rvdw_switch.
2201 * Table has a scale factor since we have scaled it down to compensate
2202 * for scaling-up c6/c12 with the derivative factors to save flops in analytical kernels.
2204 fr->enershiftsix = (real)(-1.0/(rc3*rc3)) - 6.0*vdwtab[8*ri0];
2205 fr->enershifttwelve = (real)( 1.0/(rc9*rc3)) - 12.0*vdwtab[8*ri0 + 4];
2207 /* Add the constant part from 0 to rvdw_switch.
2208 * This integration from 0 to rvdw_switch overcounts the number
2209 * of interactions by 1, as it also counts the self interaction.
2210 * We will correct for this later.
2212 eners[0] += 4.0*M_PI*fr->enershiftsix*rc3/3.0;
2213 eners[1] += 4.0*M_PI*fr->enershifttwelve*rc3/3.0;
2214 for (i = 0; i < 2; i++)
2218 integrate_table(vdwtab, scale, (i == 0 ? 0 : 4), ri0, ri1, &enersum, &virsum);
2219 eners[i] -= enersum;
2223 /* now add the correction for rvdw_switch to infinity */
2224 eners[0] += -4.0*M_PI/(3.0*rc3);
2225 eners[1] += 4.0*M_PI/(9.0*rc9);
2226 virs[0] += 8.0*M_PI/rc3;
2227 virs[1] += -16.0*M_PI/(3.0*rc9);
2229 else if (EVDW_PME(fr->vdwtype))
2231 if (EVDW_SWITCHED(fr->vdwtype) && fr->rvdw_switch == 0)
2234 "With dispersion correction rvdw-switch can not be zero "
2235 "for vdw-type = %s", evdw_names[fr->vdwtype]);
2238 scale = fr->nblists[0].table_vdw.scale;
2239 vdwtab = fr->nblists[0].table_vdw.data;
2241 ri0 = floor(fr->rvdw_switch*scale);
2242 ri1 = ceil(fr->rvdw*scale);
2248 /* Calculate self-interaction coefficient (assuming that
2249 * the reciprocal-space contribution is constant in the
2250 * region that contributes to the self-interaction).
2252 fr->enershiftsix = pow(fr->ewaldcoeff_lj, 6) / 6.0;
2254 /* Calculate C12 values as without PME. */
2255 if (EVDW_SWITCHED(fr->vdwtype))
2259 integrate_table(vdwtab, scale, 4, ri0, ri1, &enersum, &virsum);
2260 eners[1] -= enersum;
2263 /* Add analytical corrections, C6 for the whole range, C12
2264 * from rvdw_switch to infinity.
2267 eners[0] += -pow(sqrt(M_PI)*fr->ewaldcoeff_lj, 3)/3.0;
2268 eners[1] += 4.0*M_PI/(9.0*rc9);
2269 virs[0] += pow(sqrt(M_PI)*fr->ewaldcoeff_lj, 3);
2270 virs[1] += -16.0*M_PI/(3.0*rc9);
2272 else if ((fr->vdwtype == evdwCUT) || (fr->vdwtype == evdwUSER))
2274 if (fr->vdwtype == evdwUSER && fplog)
2277 "WARNING: using dispersion correction with user tables\n");
2279 rc3 = fr->rvdw*fr->rvdw*fr->rvdw;
2281 /* Contribution beyond the cut-off */
2282 eners[0] += -4.0*M_PI/(3.0*rc3);
2283 eners[1] += 4.0*M_PI/(9.0*rc9);
2284 if (fr->vdw_modifier == eintmodPOTSHIFT)
2286 /* Contribution within the cut-off */
2287 eners[0] += -4.0*M_PI/(3.0*rc3);
2288 eners[1] += 4.0*M_PI/(3.0*rc9);
2290 /* Contribution beyond the cut-off */
2291 virs[0] += 8.0*M_PI/rc3;
2292 virs[1] += -16.0*M_PI/(3.0*rc9);
2297 "Dispersion correction is not implemented for vdw-type = %s",
2298 evdw_names[fr->vdwtype]);
2300 fr->enerdiffsix = eners[0];
2301 fr->enerdifftwelve = eners[1];
2302 /* The 0.5 is due to the Gromacs definition of the virial */
2303 fr->virdiffsix = 0.5*virs[0];
2304 fr->virdifftwelve = 0.5*virs[1];
2308 void calc_dispcorr(FILE *fplog, t_inputrec *ir, t_forcerec *fr,
2309 gmx_int64_t step, int natoms,
2310 matrix box, real lambda, tensor pres, tensor virial,
2311 real *prescorr, real *enercorr, real *dvdlcorr)
2313 gmx_bool bCorrAll, bCorrPres;
2314 real dvdlambda, invvol, dens, ninter, avcsix, avctwelve, enerdiff, svir = 0, spres = 0;
2324 if (ir->eDispCorr != edispcNO)
2326 bCorrAll = (ir->eDispCorr == edispcAllEner ||
2327 ir->eDispCorr == edispcAllEnerPres);
2328 bCorrPres = (ir->eDispCorr == edispcEnerPres ||
2329 ir->eDispCorr == edispcAllEnerPres);
2331 invvol = 1/det(box);
2334 /* Only correct for the interactions with the inserted molecule */
2335 dens = (natoms - fr->n_tpi)*invvol;
2340 dens = natoms*invvol;
2341 ninter = 0.5*natoms;
2344 if (ir->efep == efepNO)
2346 avcsix = fr->avcsix[0];
2347 avctwelve = fr->avctwelve[0];
2351 avcsix = (1 - lambda)*fr->avcsix[0] + lambda*fr->avcsix[1];
2352 avctwelve = (1 - lambda)*fr->avctwelve[0] + lambda*fr->avctwelve[1];
2355 enerdiff = ninter*(dens*fr->enerdiffsix - fr->enershiftsix);
2356 *enercorr += avcsix*enerdiff;
2358 if (ir->efep != efepNO)
2360 dvdlambda += (fr->avcsix[1] - fr->avcsix[0])*enerdiff;
2364 enerdiff = ninter*(dens*fr->enerdifftwelve - fr->enershifttwelve);
2365 *enercorr += avctwelve*enerdiff;
2366 if (fr->efep != efepNO)
2368 dvdlambda += (fr->avctwelve[1] - fr->avctwelve[0])*enerdiff;
2374 svir = ninter*dens*avcsix*fr->virdiffsix/3.0;
2375 if (ir->eDispCorr == edispcAllEnerPres)
2377 svir += ninter*dens*avctwelve*fr->virdifftwelve/3.0;
2379 /* The factor 2 is because of the Gromacs virial definition */
2380 spres = -2.0*invvol*svir*PRESFAC;
2382 for (m = 0; m < DIM; m++)
2384 virial[m][m] += svir;
2385 pres[m][m] += spres;
2390 /* Can't currently control when it prints, for now, just print when degugging */
2395 fprintf(debug, "Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
2401 "Long Range LJ corr.: Epot %10g, Pres: %10g, Vir: %10g\n",
2402 *enercorr, spres, svir);
2406 fprintf(debug, "Long Range LJ corr.: Epot %10g\n", *enercorr);
2410 if (fr->bSepDVDL && do_per_step(step, ir->nstlog))
2412 gmx_print_sepdvdl(fplog, "Dispersion correction", *enercorr, dvdlambda);
2414 if (fr->efep != efepNO)
2416 *dvdlcorr += dvdlambda;
2421 void do_pbc_first(FILE *fplog, matrix box, t_forcerec *fr,
2422 t_graph *graph, rvec x[])
2426 fprintf(fplog, "Removing pbc first time\n");
2428 calc_shifts(box, fr->shift_vec);
2431 mk_mshift(fplog, graph, fr->ePBC, box, x);
2434 p_graph(debug, "do_pbc_first 1", graph);
2436 shift_self(graph, box, x);
2437 /* By doing an extra mk_mshift the molecules that are broken
2438 * because they were e.g. imported from another software
2439 * will be made whole again. Such are the healing powers
2442 mk_mshift(fplog, graph, fr->ePBC, box, x);
2445 p_graph(debug, "do_pbc_first 2", graph);
2450 fprintf(fplog, "Done rmpbc\n");
2454 static void low_do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2455 gmx_mtop_t *mtop, rvec x[],
2460 gmx_molblock_t *molb;
2462 if (bFirst && fplog)
2464 fprintf(fplog, "Removing pbc first time\n");
2469 for (mb = 0; mb < mtop->nmolblock; mb++)
2471 molb = &mtop->molblock[mb];
2472 if (molb->natoms_mol == 1 ||
2473 (!bFirst && mtop->moltype[molb->type].cgs.nr == 1))
2475 /* Just one atom or charge group in the molecule, no PBC required */
2476 as += molb->nmol*molb->natoms_mol;
2480 /* Pass NULL iso fplog to avoid graph prints for each molecule type */
2481 mk_graph_ilist(NULL, mtop->moltype[molb->type].ilist,
2482 0, molb->natoms_mol, FALSE, FALSE, graph);
2484 for (mol = 0; mol < molb->nmol; mol++)
2486 mk_mshift(fplog, graph, ePBC, box, x+as);
2488 shift_self(graph, box, x+as);
2489 /* The molecule is whole now.
2490 * We don't need the second mk_mshift call as in do_pbc_first,
2491 * since we no longer need this graph.
2494 as += molb->natoms_mol;
2502 void do_pbc_first_mtop(FILE *fplog, int ePBC, matrix box,
2503 gmx_mtop_t *mtop, rvec x[])
2505 low_do_pbc_mtop(fplog, ePBC, box, mtop, x, TRUE);
2508 void do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2509 gmx_mtop_t *mtop, rvec x[])
2511 low_do_pbc_mtop(fplog, ePBC, box, mtop, x, FALSE);
2514 void finish_run(FILE *fplog, t_commrec *cr,
2515 t_inputrec *inputrec,
2516 t_nrnb nrnb[], gmx_wallcycle_t wcycle,
2517 gmx_walltime_accounting_t walltime_accounting,
2518 wallclock_gpu_t *gputimes,
2519 gmx_bool bWriteStat)
2522 t_nrnb *nrnb_tot = NULL;
2525 double elapsed_time,
2526 elapsed_time_over_all_ranks,
2527 elapsed_time_over_all_threads,
2528 elapsed_time_over_all_threads_over_all_ranks;
2529 wallcycle_sum(cr, wcycle);
2535 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM,
2536 cr->mpi_comm_mysim);
2544 elapsed_time = walltime_accounting_get_elapsed_time(walltime_accounting);
2545 elapsed_time_over_all_ranks = elapsed_time;
2546 elapsed_time_over_all_threads = walltime_accounting_get_elapsed_time_over_all_threads(walltime_accounting);
2547 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
2551 /* reduce elapsed_time over all MPI ranks in the current simulation */
2552 MPI_Allreduce(&elapsed_time,
2553 &elapsed_time_over_all_ranks,
2554 1, MPI_DOUBLE, MPI_SUM,
2555 cr->mpi_comm_mysim);
2556 elapsed_time_over_all_ranks /= cr->nnodes;
2557 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
2558 * current simulation. */
2559 MPI_Allreduce(&elapsed_time_over_all_threads,
2560 &elapsed_time_over_all_threads_over_all_ranks,
2561 1, MPI_DOUBLE, MPI_SUM,
2562 cr->mpi_comm_mysim);
2568 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
2575 if ((cr->duty & DUTY_PP) && DOMAINDECOMP(cr))
2577 print_dd_statistics(cr, inputrec, fplog);
2589 snew(nrnb_all, cr->nnodes);
2590 nrnb_all[0] = *nrnb;
2591 for (s = 1; s < cr->nnodes; s++)
2593 MPI_Recv(nrnb_all[s].n, eNRNB, MPI_DOUBLE, s, 0,
2594 cr->mpi_comm_mysim, &stat);
2596 pr_load(fplog, cr, nrnb_all);
2601 MPI_Send(nrnb->n, eNRNB, MPI_DOUBLE, MASTERRANK(cr), 0,
2602 cr->mpi_comm_mysim);
2609 wallcycle_print(fplog, cr->nnodes, cr->npmenodes,
2610 elapsed_time_over_all_ranks,
2613 if (EI_DYNAMICS(inputrec->eI))
2615 delta_t = inputrec->delta_t;
2624 print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks,
2625 elapsed_time_over_all_ranks,
2626 walltime_accounting_get_nsteps_done(walltime_accounting),
2627 delta_t, nbfs, mflop);
2631 print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks,
2632 elapsed_time_over_all_ranks,
2633 walltime_accounting_get_nsteps_done(walltime_accounting),
2634 delta_t, nbfs, mflop);
2639 extern void initialize_lambdas(FILE *fplog, t_inputrec *ir, int *fep_state, real *lambda, double *lam0)
2641 /* this function works, but could probably use a logic rewrite to keep all the different
2642 types of efep straight. */
2645 t_lambda *fep = ir->fepvals;
2647 if ((ir->efep == efepNO) && (ir->bSimTemp == FALSE))
2649 for (i = 0; i < efptNR; i++)
2661 *fep_state = fep->init_fep_state; /* this might overwrite the checkpoint
2662 if checkpoint is set -- a kludge is in for now
2664 for (i = 0; i < efptNR; i++)
2666 /* overwrite lambda state with init_lambda for now for backwards compatibility */
2667 if (fep->init_lambda >= 0) /* if it's -1, it was never initializd */
2669 lambda[i] = fep->init_lambda;
2672 lam0[i] = lambda[i];
2677 lambda[i] = fep->all_lambda[i][*fep_state];
2680 lam0[i] = lambda[i];
2686 /* need to rescale control temperatures to match current state */
2687 for (i = 0; i < ir->opts.ngtc; i++)
2689 if (ir->opts.ref_t[i] > 0)
2691 ir->opts.ref_t[i] = ir->simtempvals->temperatures[*fep_state];
2697 /* Send to the log the information on the current lambdas */
2700 fprintf(fplog, "Initial vector of lambda components:[ ");
2701 for (i = 0; i < efptNR; i++)
2703 fprintf(fplog, "%10.4f ", lambda[i]);
2705 fprintf(fplog, "]\n");
2711 void init_md(FILE *fplog,
2712 t_commrec *cr, t_inputrec *ir, const output_env_t oenv,
2713 double *t, double *t0,
2714 real *lambda, int *fep_state, double *lam0,
2715 t_nrnb *nrnb, gmx_mtop_t *mtop,
2717 int nfile, const t_filenm fnm[],
2718 gmx_mdoutf_t *outf, t_mdebin **mdebin,
2719 tensor force_vir, tensor shake_vir, rvec mu_tot,
2720 gmx_bool *bSimAnn, t_vcm **vcm, unsigned long Flags)
2725 /* Initial values */
2726 *t = *t0 = ir->init_t;
2729 for (i = 0; i < ir->opts.ngtc; i++)
2731 /* set bSimAnn if any group is being annealed */
2732 if (ir->opts.annealing[i] != eannNO)
2739 update_annealing_target_temp(&(ir->opts), ir->init_t);
2742 /* Initialize lambda variables */
2743 initialize_lambdas(fplog, ir, fep_state, lambda, lam0);
2747 *upd = init_update(ir);
2753 *vcm = init_vcm(fplog, &mtop->groups, ir);
2756 if (EI_DYNAMICS(ir->eI) && !(Flags & MD_APPENDFILES))
2758 if (ir->etc == etcBERENDSEN)
2760 please_cite(fplog, "Berendsen84a");
2762 if (ir->etc == etcVRESCALE)
2764 please_cite(fplog, "Bussi2007a");
2772 *outf = init_mdoutf(nfile, fnm, Flags, cr, ir, mtop, oenv);
2774 *mdebin = init_mdebin((Flags & MD_APPENDFILES) ? NULL : mdoutf_get_fp_ene(*outf),
2775 mtop, ir, mdoutf_get_fp_dhdl(*outf));
2780 please_cite(fplog, "Fritsch12");
2781 please_cite(fplog, "Junghans10");
2783 /* Initiate variables */
2784 clear_mat(force_vir);
2785 clear_mat(shake_vir);