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
52 #include "chargegroup.h"
74 #include "nbnxn_atomdata.h"
75 #include "nbnxn_search.h"
76 #include "nbnxn_kernels/nbnxn_kernel_ref.h"
77 #include "nbnxn_kernels/simd_4xn/nbnxn_kernel_simd_4xn.h"
78 #include "nbnxn_kernels/simd_2xnn/nbnxn_kernel_simd_2xnn.h"
79 #include "nbnxn_kernels/nbnxn_kernel_gpu_ref.h"
81 #include "gromacs/timing/wallcycle.h"
82 #include "gromacs/timing/walltime_accounting.h"
83 #include "gromacs/utility/gmxmpi.h"
84 #include "gromacs/essentialdynamics/edsam.h"
85 #include "gromacs/pulling/pull.h"
86 #include "gromacs/pulling/pull_rotation.h"
91 #include "nbnxn_cuda_data_mgmt.h"
92 #include "nbnxn_cuda/nbnxn_cuda.h"
94 void print_time(FILE *out,
95 gmx_walltime_accounting_t walltime_accounting,
98 t_commrec gmx_unused *cr)
101 char timebuf[STRLEN];
102 double dt, elapsed_seconds, time_per_step;
105 #ifndef GMX_THREAD_MPI
111 fprintf(out, "step %s", gmx_step_str(step, buf));
112 if ((step >= ir->nstlist))
114 double seconds_since_epoch = gmx_gettime();
115 elapsed_seconds = seconds_since_epoch - walltime_accounting_get_start_time_stamp(walltime_accounting);
116 time_per_step = elapsed_seconds/(step - ir->init_step + 1);
117 dt = (ir->nsteps + ir->init_step - step) * time_per_step;
123 finish = (time_t) (seconds_since_epoch + dt);
124 gmx_ctime_r(&finish, timebuf, STRLEN);
125 sprintf(buf, "%s", timebuf);
126 buf[strlen(buf)-1] = '\0';
127 fprintf(out, ", will finish %s", buf);
131 fprintf(out, ", remaining wall clock time: %5d s ", (int)dt);
136 fprintf(out, " performance: %.1f ns/day ",
137 ir->delta_t/1000*24*60*60/time_per_step);
140 #ifndef GMX_THREAD_MPI
150 void print_date_and_time(FILE *fplog, int nodeid, const char *title,
151 const gmx_walltime_accounting_t walltime_accounting)
154 char timebuf[STRLEN];
155 char time_string[STRLEN];
160 if (walltime_accounting != NULL)
162 tmptime = (time_t) walltime_accounting_get_start_time_stamp(walltime_accounting);
163 gmx_ctime_r(&tmptime, timebuf, STRLEN);
167 tmptime = (time_t) gmx_gettime();
168 gmx_ctime_r(&tmptime, timebuf, STRLEN);
170 for (i = 0; timebuf[i] >= ' '; i++)
172 time_string[i] = timebuf[i];
174 time_string[i] = '\0';
176 fprintf(fplog, "%s on node %d %s\n", title, nodeid, time_string);
180 void print_start(FILE *fplog, t_commrec *cr,
181 gmx_walltime_accounting_t walltime_accounting,
186 sprintf(buf, "Started %s", name);
187 print_date_and_time(fplog, cr->nodeid, buf, walltime_accounting);
190 static void sum_forces(int start, int end, rvec f[], rvec flr[])
196 pr_rvecs(debug, 0, "fsr", f+start, end-start);
197 pr_rvecs(debug, 0, "flr", flr+start, end-start);
199 for (i = start; (i < end); i++)
201 rvec_inc(f[i], flr[i]);
206 * calc_f_el calculates forces due to an electric field.
208 * force is kJ mol^-1 nm^-1 = e * kJ mol^-1 nm^-1 / e
210 * Et[] contains the parameters for the time dependent
211 * part of the field (not yet used).
212 * Ex[] contains the parameters for
213 * the spatial dependent part of the field. You can have cool periodic
214 * fields in principle, but only a constant field is supported
216 * The function should return the energy due to the electric field
217 * (if any) but for now returns 0.
220 * There can be problems with the virial.
221 * Since the field is not self-consistent this is unavoidable.
222 * For neutral molecules the virial is correct within this approximation.
223 * For neutral systems with many charged molecules the error is small.
224 * But for systems with a net charge or a few charged molecules
225 * the error can be significant when the field is high.
226 * Solution: implement a self-consitent electric field into PME.
228 static void calc_f_el(FILE *fp, int start, int homenr,
229 real charge[], rvec f[],
230 t_cosines Ex[], t_cosines Et[], double t)
236 for (m = 0; (m < DIM); m++)
243 Ext[m] = cos(Et[m].a[0]*(t-t0))*exp(-sqr(t-t0)/(2.0*sqr(Et[m].a[2])));
247 Ext[m] = cos(Et[m].a[0]*t);
256 /* Convert the field strength from V/nm to MD-units */
257 Ext[m] *= Ex[m].a[0]*FIELDFAC;
258 for (i = start; (i < start+homenr); i++)
260 f[i][m] += charge[i]*Ext[m];
270 fprintf(fp, "%10g %10g %10g %10g #FIELD\n", t,
271 Ext[XX]/FIELDFAC, Ext[YY]/FIELDFAC, Ext[ZZ]/FIELDFAC);
275 static void calc_virial(int start, int homenr, rvec x[], rvec f[],
276 tensor vir_part, t_graph *graph, matrix box,
277 t_nrnb *nrnb, const t_forcerec *fr, int ePBC)
282 /* The short-range virial from surrounding boxes */
284 calc_vir(SHIFTS, fr->shift_vec, fr->fshift, vir_part, ePBC == epbcSCREW, box);
285 inc_nrnb(nrnb, eNR_VIRIAL, SHIFTS);
287 /* Calculate partial virial, for local atoms only, based on short range.
288 * Total virial is computed in global_stat, called from do_md
290 f_calc_vir(start, start+homenr, x, f, vir_part, graph, box);
291 inc_nrnb(nrnb, eNR_VIRIAL, homenr);
293 /* Add position restraint contribution */
294 for (i = 0; i < DIM; i++)
296 vir_part[i][i] += fr->vir_diag_posres[i];
299 /* Add wall contribution */
300 for (i = 0; i < DIM; i++)
302 vir_part[i][ZZ] += fr->vir_wall_z[i];
307 pr_rvecs(debug, 0, "vir_part", vir_part, DIM);
311 static void posres_wrapper(FILE *fplog,
317 matrix box, rvec x[],
318 gmx_enerdata_t *enerd,
326 /* Position restraints always require full pbc */
327 set_pbc(&pbc, ir->ePBC, box);
329 v = posres(top->idef.il[F_POSRES].nr, top->idef.il[F_POSRES].iatoms,
330 top->idef.iparams_posres,
331 (const rvec*)x, fr->f_novirsum, fr->vir_diag_posres,
332 ir->ePBC == epbcNONE ? NULL : &pbc,
333 lambda[efptRESTRAINT], &dvdl,
334 fr->rc_scaling, fr->ePBC, fr->posres_com, fr->posres_comB);
337 gmx_print_sepdvdl(fplog, interaction_function[F_POSRES].longname, v, dvdl);
339 enerd->term[F_POSRES] += v;
340 /* If just the force constant changes, the FEP term is linear,
341 * but if k changes, it is not.
343 enerd->dvdl_nonlin[efptRESTRAINT] += dvdl;
344 inc_nrnb(nrnb, eNR_POSRES, top->idef.il[F_POSRES].nr/2);
346 if ((ir->fepvals->n_lambda > 0) && (flags & GMX_FORCE_DHDL))
348 for (i = 0; i < enerd->n_lambda; i++)
350 real dvdl_dum, lambda_dum;
352 lambda_dum = (i == 0 ? lambda[efptRESTRAINT] : ir->fepvals->all_lambda[efptRESTRAINT][i-1]);
353 v = posres(top->idef.il[F_POSRES].nr, top->idef.il[F_POSRES].iatoms,
354 top->idef.iparams_posres,
355 (const rvec*)x, NULL, NULL,
356 ir->ePBC == epbcNONE ? NULL : &pbc, lambda_dum, &dvdl,
357 fr->rc_scaling, fr->ePBC, fr->posres_com, fr->posres_comB);
358 enerd->enerpart_lambda[i] += v;
363 static void fbposres_wrapper(t_inputrec *ir,
366 matrix box, rvec x[],
367 gmx_enerdata_t *enerd,
373 /* Flat-bottomed position restraints always require full pbc */
374 set_pbc(&pbc, ir->ePBC, box);
375 v = fbposres(top->idef.il[F_FBPOSRES].nr, top->idef.il[F_FBPOSRES].iatoms,
376 top->idef.iparams_fbposres,
377 (const rvec*)x, fr->f_novirsum, fr->vir_diag_posres,
378 ir->ePBC == epbcNONE ? NULL : &pbc,
379 fr->rc_scaling, fr->ePBC, fr->posres_com);
380 enerd->term[F_FBPOSRES] += v;
381 inc_nrnb(nrnb, eNR_FBPOSRES, top->idef.il[F_FBPOSRES].nr/2);
384 static void pull_potential_wrapper(FILE *fplog,
388 matrix box, rvec x[],
392 gmx_enerdata_t *enerd,
399 /* Calculate the center of mass forces, this requires communication,
400 * which is why pull_potential is called close to other communication.
401 * The virial contribution is calculated directly,
402 * which is why we call pull_potential after calc_virial.
404 set_pbc(&pbc, ir->ePBC, box);
406 enerd->term[F_COM_PULL] +=
407 pull_potential(ir->ePull, ir->pull, mdatoms, &pbc,
408 cr, t, lambda[efptRESTRAINT], x, f, vir_force, &dvdl);
411 gmx_print_sepdvdl(fplog, "Com pull", enerd->term[F_COM_PULL], dvdl);
413 enerd->dvdl_lin[efptRESTRAINT] += dvdl;
416 static void pme_receive_force_ener(FILE *fplog,
419 gmx_wallcycle_t wcycle,
420 gmx_enerdata_t *enerd,
423 real e_q, e_lj, v, dvdl_q, dvdl_lj;
424 float cycles_ppdpme, cycles_seppme;
426 cycles_ppdpme = wallcycle_stop(wcycle, ewcPPDURINGPME);
427 dd_cycles_add(cr->dd, cycles_ppdpme, ddCyclPPduringPME);
429 /* In case of node-splitting, the PP nodes receive the long-range
430 * forces, virial and energy from the PME nodes here.
432 wallcycle_start(wcycle, ewcPP_PMEWAITRECVF);
435 gmx_pme_receive_f(cr, fr->f_novirsum, fr->vir_el_recip, &e_q,
436 fr->vir_lj_recip, &e_lj, &dvdl_q, &dvdl_lj,
440 gmx_print_sepdvdl(fplog, "Electrostatic PME mesh", e_q, dvdl_q);
441 gmx_print_sepdvdl(fplog, "Lennard-Jones PME mesh", e_lj, dvdl_lj);
443 enerd->term[F_COUL_RECIP] += e_q;
444 enerd->term[F_LJ_RECIP] += e_lj;
445 enerd->dvdl_lin[efptCOUL] += dvdl_q;
446 enerd->dvdl_lin[efptVDW] += dvdl_lj;
450 dd_cycles_add(cr->dd, cycles_seppme, ddCyclPME);
452 wallcycle_stop(wcycle, ewcPP_PMEWAITRECVF);
455 static void print_large_forces(FILE *fp, t_mdatoms *md, t_commrec *cr,
456 gmx_int64_t step, real pforce, rvec *x, rvec *f)
460 char buf[STEPSTRSIZE];
463 for (i = 0; i < md->homenr; i++)
466 /* We also catch NAN, if the compiler does not optimize this away. */
467 if (fn2 >= pf2 || fn2 != fn2)
469 fprintf(fp, "step %s atom %6d x %8.3f %8.3f %8.3f force %12.5e\n",
470 gmx_step_str(step, buf),
471 ddglatnr(cr->dd, i), x[i][XX], x[i][YY], x[i][ZZ], sqrt(fn2));
476 static void post_process_forces(t_commrec *cr,
478 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
480 matrix box, rvec x[],
485 t_forcerec *fr, gmx_vsite_t *vsite,
492 /* Spread the mesh force on virtual sites to the other particles...
493 * This is parallellized. MPI communication is performed
494 * if the constructing atoms aren't local.
496 wallcycle_start(wcycle, ewcVSITESPREAD);
497 spread_vsite_f(vsite, x, fr->f_novirsum, NULL,
498 (flags & GMX_FORCE_VIRIAL), fr->vir_el_recip,
500 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
501 wallcycle_stop(wcycle, ewcVSITESPREAD);
503 if (flags & GMX_FORCE_VIRIAL)
505 /* Now add the forces, this is local */
508 sum_forces(0, fr->f_novirsum_n, f, fr->f_novirsum);
512 sum_forces(0, mdatoms->homenr,
515 if (EEL_FULL(fr->eeltype))
517 /* Add the mesh contribution to the virial */
518 m_add(vir_force, fr->vir_el_recip, vir_force);
520 if (EVDW_PME(fr->vdwtype))
522 /* Add the mesh contribution to the virial */
523 m_add(vir_force, fr->vir_lj_recip, vir_force);
527 pr_rvecs(debug, 0, "vir_force", vir_force, DIM);
532 if (fr->print_force >= 0)
534 print_large_forces(stderr, mdatoms, cr, step, fr->print_force, x, f);
538 static void do_nb_verlet(t_forcerec *fr,
539 interaction_const_t *ic,
540 gmx_enerdata_t *enerd,
541 int flags, int ilocality,
544 gmx_wallcycle_t wcycle)
546 int nnbl, kernel_type, enr_nbnxn_kernel_ljc, enr_nbnxn_kernel_lj;
548 nonbonded_verlet_group_t *nbvg;
551 if (!(flags & GMX_FORCE_NONBONDED))
553 /* skip non-bonded calculation */
557 nbvg = &fr->nbv->grp[ilocality];
559 /* CUDA kernel launch overhead is already timed separately */
560 if (fr->cutoff_scheme != ecutsVERLET)
562 gmx_incons("Invalid cut-off scheme passed!");
565 bCUDA = (nbvg->kernel_type == nbnxnk8x8x8_CUDA);
569 wallcycle_sub_start(wcycle, ewcsNONBONDED);
571 switch (nbvg->kernel_type)
573 case nbnxnk4x4_PlainC:
574 nbnxn_kernel_ref(&nbvg->nbl_lists,
580 enerd->grpp.ener[egCOULSR],
582 enerd->grpp.ener[egBHAMSR] :
583 enerd->grpp.ener[egLJSR]);
586 case nbnxnk4xN_SIMD_4xN:
587 nbnxn_kernel_simd_4xn(&nbvg->nbl_lists,
594 enerd->grpp.ener[egCOULSR],
596 enerd->grpp.ener[egBHAMSR] :
597 enerd->grpp.ener[egLJSR]);
599 case nbnxnk4xN_SIMD_2xNN:
600 nbnxn_kernel_simd_2xnn(&nbvg->nbl_lists,
607 enerd->grpp.ener[egCOULSR],
609 enerd->grpp.ener[egBHAMSR] :
610 enerd->grpp.ener[egLJSR]);
613 case nbnxnk8x8x8_CUDA:
614 nbnxn_cuda_launch_kernel(fr->nbv->cu_nbv, nbvg->nbat, flags, ilocality);
617 case nbnxnk8x8x8_PlainC:
618 nbnxn_kernel_gpu_ref(nbvg->nbl_lists.nbl[0],
623 nbvg->nbat->out[0].f,
625 enerd->grpp.ener[egCOULSR],
627 enerd->grpp.ener[egBHAMSR] :
628 enerd->grpp.ener[egLJSR]);
632 gmx_incons("Invalid nonbonded kernel type passed!");
637 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
640 if (EEL_RF(ic->eeltype) || ic->eeltype == eelCUT)
642 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_RF;
644 else if ((!bCUDA && nbvg->ewald_excl == ewaldexclAnalytical) ||
645 (bCUDA && nbnxn_cuda_is_kernel_ewald_analytical(fr->nbv->cu_nbv)))
647 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_EWALD;
651 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_TAB;
653 enr_nbnxn_kernel_lj = eNR_NBNXN_LJ;
654 if (flags & GMX_FORCE_ENERGY)
656 /* In eNR_??? the nbnxn F+E kernels are always the F kernel + 1 */
657 enr_nbnxn_kernel_ljc += 1;
658 enr_nbnxn_kernel_lj += 1;
661 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc,
662 nbvg->nbl_lists.natpair_ljq);
663 inc_nrnb(nrnb, enr_nbnxn_kernel_lj,
664 nbvg->nbl_lists.natpair_lj);
665 /* The Coulomb-only kernels are offset -eNR_NBNXN_LJ_RF+eNR_NBNXN_RF */
666 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc-eNR_NBNXN_LJ_RF+eNR_NBNXN_RF,
667 nbvg->nbl_lists.natpair_q);
669 if (ic->vdw_modifier == eintmodFORCESWITCH)
671 /* We add up the switch cost separately */
672 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_FSW+((flags & GMX_FORCE_ENERGY) ? 1 : 0),
673 nbvg->nbl_lists.natpair_ljq + nbvg->nbl_lists.natpair_lj);
675 if (ic->vdw_modifier == eintmodPOTSWITCH)
677 /* We add up the switch cost separately */
678 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_PSW+((flags & GMX_FORCE_ENERGY) ? 1 : 0),
679 nbvg->nbl_lists.natpair_ljq + nbvg->nbl_lists.natpair_lj);
681 if (ic->vdwtype == evdwPME)
683 /* We add up the LJ Ewald cost separately */
684 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_EWALD+((flags & GMX_FORCE_ENERGY) ? 1 : 0),
685 nbvg->nbl_lists.natpair_ljq + nbvg->nbl_lists.natpair_lj);
689 void do_force_cutsVERLET(FILE *fplog, t_commrec *cr,
690 t_inputrec *inputrec,
691 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
693 gmx_groups_t gmx_unused *groups,
694 matrix box, rvec x[], history_t *hist,
698 gmx_enerdata_t *enerd, t_fcdata *fcd,
699 real *lambda, t_graph *graph,
700 t_forcerec *fr, interaction_const_t *ic,
701 gmx_vsite_t *vsite, rvec mu_tot,
702 double t, FILE *field, gmx_edsam_t ed,
710 gmx_bool bSepDVDL, bStateChanged, bNS, bFillGrid, bCalcCGCM, bBS;
711 gmx_bool bDoLongRange, bDoForces, bSepLRF, bUseGPU, bUseOrEmulGPU;
712 gmx_bool bDiffKernels = FALSE;
714 rvec vzero, box_diag;
716 float cycles_pme, cycles_force, cycles_wait_gpu;
717 nonbonded_verlet_t *nbv;
722 nb_kernel_type = fr->nbv->grp[0].kernel_type;
725 homenr = mdatoms->homenr;
727 bSepDVDL = (fr->bSepDVDL && do_per_step(step, inputrec->nstlog));
729 clear_mat(vir_force);
732 if (DOMAINDECOMP(cr))
734 cg1 = cr->dd->ncg_tot;
745 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
746 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
747 bFillGrid = (bNS && bStateChanged);
748 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
749 bDoLongRange = (fr->bTwinRange && bNS && (flags & GMX_FORCE_DO_LR));
750 bDoForces = (flags & GMX_FORCE_FORCES);
751 bSepLRF = (bDoLongRange && bDoForces && (flags & GMX_FORCE_SEPLRF));
752 bUseGPU = fr->nbv->bUseGPU;
753 bUseOrEmulGPU = bUseGPU || (nbv->grp[0].kernel_type == nbnxnk8x8x8_PlainC);
757 update_forcerec(fr, box);
759 if (NEED_MUTOT(*inputrec))
761 /* Calculate total (local) dipole moment in a temporary common array.
762 * This makes it possible to sum them over nodes faster.
764 calc_mu(start, homenr,
765 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
770 if (fr->ePBC != epbcNONE)
772 /* Compute shift vectors every step,
773 * because of pressure coupling or box deformation!
775 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
777 calc_shifts(box, fr->shift_vec);
782 put_atoms_in_box_omp(fr->ePBC, box, homenr, x);
783 inc_nrnb(nrnb, eNR_SHIFTX, homenr);
785 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
787 unshift_self(graph, box, x);
791 nbnxn_atomdata_copy_shiftvec(flags & GMX_FORCE_DYNAMICBOX,
792 fr->shift_vec, nbv->grp[0].nbat);
795 if (!(cr->duty & DUTY_PME))
797 /* Send particle coordinates to the pme nodes.
798 * Since this is only implemented for domain decomposition
799 * and domain decomposition does not use the graph,
800 * we do not need to worry about shifting.
805 wallcycle_start(wcycle, ewcPP_PMESENDX);
807 bBS = (inputrec->nwall == 2);
811 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
814 if (EEL_PME(fr->eeltype))
816 pme_flags |= GMX_PME_DO_COULOMB;
819 if (EVDW_PME(fr->vdwtype))
821 pme_flags |= GMX_PME_DO_LJ;
822 if (fr->ljpme_combination_rule == eljpmeLB)
824 pme_flags |= GMX_PME_LJ_LB;
828 gmx_pme_send_coordinates(cr, bBS ? boxs : box, x,
829 mdatoms->nChargePerturbed, mdatoms->nTypePerturbed, lambda[efptCOUL], lambda[efptVDW],
830 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),
833 wallcycle_stop(wcycle, ewcPP_PMESENDX);
837 /* do gridding for pair search */
840 if (graph && bStateChanged)
842 /* Calculate intramolecular shift vectors to make molecules whole */
843 mk_mshift(fplog, graph, fr->ePBC, box, x);
847 box_diag[XX] = box[XX][XX];
848 box_diag[YY] = box[YY][YY];
849 box_diag[ZZ] = box[ZZ][ZZ];
851 wallcycle_start(wcycle, ewcNS);
854 wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
855 nbnxn_put_on_grid(nbv->nbs, fr->ePBC, box,
857 0, mdatoms->homenr, -1, fr->cginfo, x,
859 nbv->grp[eintLocal].kernel_type,
860 nbv->grp[eintLocal].nbat);
861 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
865 wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
866 nbnxn_put_on_grid_nonlocal(nbv->nbs, domdec_zones(cr->dd),
868 nbv->grp[eintNonlocal].kernel_type,
869 nbv->grp[eintNonlocal].nbat);
870 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
873 if (nbv->ngrp == 1 ||
874 nbv->grp[eintNonlocal].nbat == nbv->grp[eintLocal].nbat)
876 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatAll,
877 nbv->nbs, mdatoms, fr->cginfo);
881 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatLocal,
882 nbv->nbs, mdatoms, fr->cginfo);
883 nbnxn_atomdata_set(nbv->grp[eintNonlocal].nbat, eatAll,
884 nbv->nbs, mdatoms, fr->cginfo);
886 wallcycle_stop(wcycle, ewcNS);
889 /* initialize the GPU atom data and copy shift vector */
894 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
895 nbnxn_cuda_init_atomdata(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
896 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
899 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
900 nbnxn_cuda_upload_shiftvec(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
901 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
904 /* do local pair search */
907 wallcycle_start_nocount(wcycle, ewcNS);
908 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
909 nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintLocal].nbat,
912 nbv->min_ci_balanced,
913 &nbv->grp[eintLocal].nbl_lists,
915 nbv->grp[eintLocal].kernel_type,
917 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
921 /* initialize local pair-list on the GPU */
922 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
923 nbv->grp[eintLocal].nbl_lists.nbl[0],
926 wallcycle_stop(wcycle, ewcNS);
930 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
931 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
932 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, FALSE, x,
933 nbv->grp[eintLocal].nbat);
934 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
935 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
940 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
941 /* launch local nonbonded F on GPU */
942 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFNo,
944 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
947 /* Communicate coordinates and sum dipole if necessary +
948 do non-local pair search */
949 if (DOMAINDECOMP(cr))
951 bDiffKernels = (nbv->grp[eintNonlocal].kernel_type !=
952 nbv->grp[eintLocal].kernel_type);
956 /* With GPU+CPU non-bonded calculations we need to copy
957 * the local coordinates to the non-local nbat struct
958 * (in CPU format) as the non-local kernel call also
959 * calculates the local - non-local interactions.
961 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
962 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
963 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, TRUE, x,
964 nbv->grp[eintNonlocal].nbat);
965 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
966 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
971 wallcycle_start_nocount(wcycle, ewcNS);
972 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
976 nbnxn_grid_add_simple(nbv->nbs, nbv->grp[eintNonlocal].nbat);
979 nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintNonlocal].nbat,
982 nbv->min_ci_balanced,
983 &nbv->grp[eintNonlocal].nbl_lists,
985 nbv->grp[eintNonlocal].kernel_type,
988 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
990 if (nbv->grp[eintNonlocal].kernel_type == nbnxnk8x8x8_CUDA)
992 /* initialize non-local pair-list on the GPU */
993 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
994 nbv->grp[eintNonlocal].nbl_lists.nbl[0],
997 wallcycle_stop(wcycle, ewcNS);
1001 wallcycle_start(wcycle, ewcMOVEX);
1002 dd_move_x(cr->dd, box, x);
1004 /* When we don't need the total dipole we sum it in global_stat */
1005 if (bStateChanged && NEED_MUTOT(*inputrec))
1007 gmx_sumd(2*DIM, mu, cr);
1009 wallcycle_stop(wcycle, ewcMOVEX);
1011 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1012 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
1013 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatNonlocal, FALSE, x,
1014 nbv->grp[eintNonlocal].nbat);
1015 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
1016 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1019 if (bUseGPU && !bDiffKernels)
1021 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
1022 /* launch non-local nonbonded F on GPU */
1023 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFNo,
1025 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1031 /* launch D2H copy-back F */
1032 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1033 if (DOMAINDECOMP(cr) && !bDiffKernels)
1035 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintNonlocal].nbat,
1036 flags, eatNonlocal);
1038 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintLocal].nbat,
1040 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1043 if (bStateChanged && NEED_MUTOT(*inputrec))
1047 gmx_sumd(2*DIM, mu, cr);
1050 for (i = 0; i < 2; i++)
1052 for (j = 0; j < DIM; j++)
1054 fr->mu_tot[i][j] = mu[i*DIM + j];
1058 if (fr->efep == efepNO)
1060 copy_rvec(fr->mu_tot[0], mu_tot);
1064 for (j = 0; j < DIM; j++)
1067 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] +
1068 lambda[efptCOUL]*fr->mu_tot[1][j];
1072 /* Reset energies */
1073 reset_enerdata(fr, bNS, enerd, MASTER(cr));
1074 clear_rvecs(SHIFTS, fr->fshift);
1076 if (DOMAINDECOMP(cr) && !(cr->duty & DUTY_PME))
1078 wallcycle_start(wcycle, ewcPPDURINGPME);
1079 dd_force_flop_start(cr->dd, nrnb);
1084 /* Enforced rotation has its own cycle counter that starts after the collective
1085 * coordinates have been communicated. It is added to ddCyclF to allow
1086 * for proper load-balancing */
1087 wallcycle_start(wcycle, ewcROT);
1088 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1089 wallcycle_stop(wcycle, ewcROT);
1092 /* Start the force cycle counter.
1093 * This counter is stopped in do_forcelow_level.
1094 * No parallel communication should occur while this counter is running,
1095 * since that will interfere with the dynamic load balancing.
1097 wallcycle_start(wcycle, ewcFORCE);
1100 /* Reset forces for which the virial is calculated separately:
1101 * PME/Ewald forces if necessary */
1102 if (fr->bF_NoVirSum)
1104 if (flags & GMX_FORCE_VIRIAL)
1106 fr->f_novirsum = fr->f_novirsum_alloc;
1109 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1113 clear_rvecs(homenr, fr->f_novirsum+start);
1118 /* We are not calculating the pressure so we do not need
1119 * a separate array for forces that do not contribute
1126 /* Clear the short- and long-range forces */
1127 clear_rvecs(fr->natoms_force_constr, f);
1128 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1130 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1133 clear_rvec(fr->vir_diag_posres);
1136 if (inputrec->ePull == epullCONSTRAINT)
1138 clear_pull_forces(inputrec->pull);
1141 /* We calculate the non-bonded forces, when done on the CPU, here.
1142 * We do this before calling do_force_lowlevel, as in there bondeds
1143 * forces are calculated before PME, which does communication.
1144 * With this order, non-bonded and bonded force calculation imbalance
1145 * can be balanced out by the domain decomposition load balancing.
1150 /* Maybe we should move this into do_force_lowlevel */
1151 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFYes,
1155 if (!bUseOrEmulGPU || bDiffKernels)
1159 if (DOMAINDECOMP(cr))
1161 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal,
1162 bDiffKernels ? enbvClearFYes : enbvClearFNo,
1172 aloc = eintNonlocal;
1175 /* Add all the non-bonded force to the normal force array.
1176 * This can be split into a local a non-local part when overlapping
1177 * communication with calculation with domain decomposition.
1179 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1180 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1181 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1182 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatAll, nbv->grp[aloc].nbat, f);
1183 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1184 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1185 wallcycle_start_nocount(wcycle, ewcFORCE);
1187 /* if there are multiple fshift output buffers reduce them */
1188 if ((flags & GMX_FORCE_VIRIAL) &&
1189 nbv->grp[aloc].nbl_lists.nnbl > 1)
1191 nbnxn_atomdata_add_nbat_fshift_to_fshift(nbv->grp[aloc].nbat,
1196 /* update QMMMrec, if necessary */
1199 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1202 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1204 posres_wrapper(fplog, flags, bSepDVDL, inputrec, nrnb, top, box, x,
1208 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_FBPOSRES].nr > 0)
1210 fbposres_wrapper(inputrec, nrnb, top, box, x, enerd, fr);
1213 /* Compute the bonded and non-bonded energies and optionally forces */
1214 do_force_lowlevel(fplog, step, fr, inputrec, &(top->idef),
1215 cr, nrnb, wcycle, mdatoms,
1216 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, top, fr->born,
1217 &(top->atomtypes), bBornRadii, box,
1218 inputrec->fepvals, lambda, graph, &(top->excls), fr->mu_tot,
1219 flags, &cycles_pme);
1223 if (do_per_step(step, inputrec->nstcalclr))
1225 /* Add the long range forces to the short range forces */
1226 for (i = 0; i < fr->natoms_force_constr; i++)
1228 rvec_add(fr->f_twin[i], f[i], f[i]);
1233 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1237 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1240 if (bUseOrEmulGPU && !bDiffKernels)
1242 /* wait for non-local forces (or calculate in emulation mode) */
1243 if (DOMAINDECOMP(cr))
1249 wallcycle_start(wcycle, ewcWAIT_GPU_NB_NL);
1250 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1251 nbv->grp[eintNonlocal].nbat,
1253 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1255 cycles_tmp = wallcycle_stop(wcycle, ewcWAIT_GPU_NB_NL);
1256 cycles_wait_gpu += cycles_tmp;
1257 cycles_force += cycles_tmp;
1261 wallcycle_start_nocount(wcycle, ewcFORCE);
1262 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFYes,
1264 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1266 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1267 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1268 /* skip the reduction if there was no non-local work to do */
1269 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1271 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatNonlocal,
1272 nbv->grp[eintNonlocal].nbat, f);
1274 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1275 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1279 if (bDoForces && DOMAINDECOMP(cr))
1281 /* Communicate the forces */
1282 wallcycle_start(wcycle, ewcMOVEF);
1283 dd_move_f(cr->dd, f, fr->fshift);
1284 /* Do we need to communicate the separate force array
1285 * for terms that do not contribute to the single sum virial?
1286 * Position restraints and electric fields do not introduce
1287 * inter-cg forces, only full electrostatics methods do.
1288 * When we do not calculate the virial, fr->f_novirsum = f,
1289 * so we have already communicated these forces.
1291 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1292 (flags & GMX_FORCE_VIRIAL))
1294 dd_move_f(cr->dd, fr->f_novirsum, NULL);
1298 /* We should not update the shift forces here,
1299 * since f_twin is already included in f.
1301 dd_move_f(cr->dd, fr->f_twin, NULL);
1303 wallcycle_stop(wcycle, ewcMOVEF);
1308 /* wait for local forces (or calculate in emulation mode) */
1311 wallcycle_start(wcycle, ewcWAIT_GPU_NB_L);
1312 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1313 nbv->grp[eintLocal].nbat,
1315 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1317 cycles_wait_gpu += wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
1319 /* now clear the GPU outputs while we finish the step on the CPU */
1321 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1322 nbnxn_cuda_clear_outputs(nbv->cu_nbv, flags);
1323 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1327 wallcycle_start_nocount(wcycle, ewcFORCE);
1328 do_nb_verlet(fr, ic, enerd, flags, eintLocal,
1329 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
1331 wallcycle_stop(wcycle, ewcFORCE);
1333 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1334 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1335 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1337 /* skip the reduction if there was no non-local work to do */
1338 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatLocal,
1339 nbv->grp[eintLocal].nbat, f);
1341 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1342 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1345 if (DOMAINDECOMP(cr))
1347 dd_force_flop_stop(cr->dd, nrnb);
1350 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1353 dd_cycles_add(cr->dd, cycles_wait_gpu, ddCyclWaitGPU);
1360 if (IR_ELEC_FIELD(*inputrec))
1362 /* Compute forces due to electric field */
1363 calc_f_el(MASTER(cr) ? field : NULL,
1364 start, homenr, mdatoms->chargeA, fr->f_novirsum,
1365 inputrec->ex, inputrec->et, t);
1368 /* If we have NoVirSum forces, but we do not calculate the virial,
1369 * we sum fr->f_novirum=f later.
1371 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1373 wallcycle_start(wcycle, ewcVSITESPREAD);
1374 spread_vsite_f(vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1375 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1376 wallcycle_stop(wcycle, ewcVSITESPREAD);
1380 wallcycle_start(wcycle, ewcVSITESPREAD);
1381 spread_vsite_f(vsite, x, fr->f_twin, NULL, FALSE, NULL,
1383 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1384 wallcycle_stop(wcycle, ewcVSITESPREAD);
1388 if (flags & GMX_FORCE_VIRIAL)
1390 /* Calculation of the virial must be done after vsites! */
1391 calc_virial(0, mdatoms->homenr, x, f,
1392 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1396 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1398 pull_potential_wrapper(fplog, bSepDVDL, cr, inputrec, box, x,
1399 f, vir_force, mdatoms, enerd, lambda, t);
1402 /* Add the forces from enforced rotation potentials (if any) */
1405 wallcycle_start(wcycle, ewcROTadd);
1406 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
1407 wallcycle_stop(wcycle, ewcROTadd);
1410 if (PAR(cr) && !(cr->duty & DUTY_PME))
1412 /* In case of node-splitting, the PP nodes receive the long-range
1413 * forces, virial and energy from the PME nodes here.
1415 pme_receive_force_ener(fplog, bSepDVDL, cr, wcycle, enerd, fr);
1420 post_process_forces(cr, step, nrnb, wcycle,
1421 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
1425 /* Sum the potential energy terms from group contributions */
1426 sum_epot(&(enerd->grpp), enerd->term);
1429 void do_force_cutsGROUP(FILE *fplog, t_commrec *cr,
1430 t_inputrec *inputrec,
1431 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1432 gmx_localtop_t *top,
1433 gmx_groups_t *groups,
1434 matrix box, rvec x[], history_t *hist,
1438 gmx_enerdata_t *enerd, t_fcdata *fcd,
1439 real *lambda, t_graph *graph,
1440 t_forcerec *fr, gmx_vsite_t *vsite, rvec mu_tot,
1441 double t, FILE *field, gmx_edsam_t ed,
1442 gmx_bool bBornRadii,
1448 gmx_bool bSepDVDL, bStateChanged, bNS, bFillGrid, bCalcCGCM, bBS;
1449 gmx_bool bDoLongRangeNS, bDoForces, bDoPotential, bSepLRF;
1450 gmx_bool bDoAdressWF;
1452 rvec vzero, box_diag;
1453 real e, v, dvdlambda[efptNR];
1455 float cycles_pme, cycles_force;
1458 homenr = mdatoms->homenr;
1460 bSepDVDL = (fr->bSepDVDL && do_per_step(step, inputrec->nstlog));
1462 clear_mat(vir_force);
1465 if (DOMAINDECOMP(cr))
1467 cg1 = cr->dd->ncg_tot;
1478 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
1479 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
1480 /* Should we update the long-range neighborlists at this step? */
1481 bDoLongRangeNS = fr->bTwinRange && bNS;
1482 /* Should we perform the long-range nonbonded evaluation inside the neighborsearching? */
1483 bFillGrid = (bNS && bStateChanged);
1484 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
1485 bDoForces = (flags & GMX_FORCE_FORCES);
1486 bDoPotential = (flags & GMX_FORCE_ENERGY);
1487 bSepLRF = ((inputrec->nstcalclr > 1) && bDoForces &&
1488 (flags & GMX_FORCE_SEPLRF) && (flags & GMX_FORCE_DO_LR));
1490 /* should probably move this to the forcerec since it doesn't change */
1491 bDoAdressWF = ((fr->adress_type != eAdressOff));
1495 update_forcerec(fr, box);
1497 if (NEED_MUTOT(*inputrec))
1499 /* Calculate total (local) dipole moment in a temporary common array.
1500 * This makes it possible to sum them over nodes faster.
1502 calc_mu(start, homenr,
1503 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
1508 if (fr->ePBC != epbcNONE)
1510 /* Compute shift vectors every step,
1511 * because of pressure coupling or box deformation!
1513 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
1515 calc_shifts(box, fr->shift_vec);
1520 put_charge_groups_in_box(fplog, cg0, cg1, fr->ePBC, box,
1521 &(top->cgs), x, fr->cg_cm);
1522 inc_nrnb(nrnb, eNR_CGCM, homenr);
1523 inc_nrnb(nrnb, eNR_RESETX, cg1-cg0);
1525 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
1527 unshift_self(graph, box, x);
1532 calc_cgcm(fplog, cg0, cg1, &(top->cgs), x, fr->cg_cm);
1533 inc_nrnb(nrnb, eNR_CGCM, homenr);
1536 if (bCalcCGCM && gmx_debug_at)
1538 pr_rvecs(debug, 0, "cgcm", fr->cg_cm, top->cgs.nr);
1542 if (!(cr->duty & DUTY_PME))
1544 /* Send particle coordinates to the pme nodes.
1545 * Since this is only implemented for domain decomposition
1546 * and domain decomposition does not use the graph,
1547 * we do not need to worry about shifting.
1552 wallcycle_start(wcycle, ewcPP_PMESENDX);
1554 bBS = (inputrec->nwall == 2);
1557 copy_mat(box, boxs);
1558 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
1561 if (EEL_PME(fr->eeltype))
1563 pme_flags |= GMX_PME_DO_COULOMB;
1566 if (EVDW_PME(fr->vdwtype))
1568 pme_flags |= GMX_PME_DO_LJ;
1569 if (fr->ljpme_combination_rule == eljpmeLB)
1571 pme_flags |= GMX_PME_LJ_LB;
1575 gmx_pme_send_coordinates(cr, bBS ? boxs : box, x,
1576 mdatoms->nChargePerturbed, mdatoms->nTypePerturbed, lambda[efptCOUL], lambda[efptVDW],
1577 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),
1580 wallcycle_stop(wcycle, ewcPP_PMESENDX);
1582 #endif /* GMX_MPI */
1584 /* Communicate coordinates and sum dipole if necessary */
1585 if (DOMAINDECOMP(cr))
1587 wallcycle_start(wcycle, ewcMOVEX);
1588 dd_move_x(cr->dd, box, x);
1589 wallcycle_stop(wcycle, ewcMOVEX);
1592 /* update adress weight beforehand */
1593 if (bStateChanged && bDoAdressWF)
1595 /* need pbc for adress weight calculation with pbc_dx */
1596 set_pbc(&pbc, inputrec->ePBC, box);
1597 if (fr->adress_site == eAdressSITEcog)
1599 update_adress_weights_cog(top->idef.iparams, top->idef.il, x, fr, mdatoms,
1600 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1602 else if (fr->adress_site == eAdressSITEcom)
1604 update_adress_weights_com(fplog, cg0, cg1, &(top->cgs), x, fr, mdatoms,
1605 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1607 else if (fr->adress_site == eAdressSITEatomatom)
1609 update_adress_weights_atom_per_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1610 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1614 update_adress_weights_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1615 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1619 if (NEED_MUTOT(*inputrec))
1626 gmx_sumd(2*DIM, mu, cr);
1628 for (i = 0; i < 2; i++)
1630 for (j = 0; j < DIM; j++)
1632 fr->mu_tot[i][j] = mu[i*DIM + j];
1636 if (fr->efep == efepNO)
1638 copy_rvec(fr->mu_tot[0], mu_tot);
1642 for (j = 0; j < DIM; j++)
1645 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] + lambda[efptCOUL]*fr->mu_tot[1][j];
1650 /* Reset energies */
1651 reset_enerdata(fr, bNS, enerd, MASTER(cr));
1652 clear_rvecs(SHIFTS, fr->fshift);
1656 wallcycle_start(wcycle, ewcNS);
1658 if (graph && bStateChanged)
1660 /* Calculate intramolecular shift vectors to make molecules whole */
1661 mk_mshift(fplog, graph, fr->ePBC, box, x);
1664 /* Do the actual neighbour searching */
1666 groups, top, mdatoms,
1667 cr, nrnb, bFillGrid,
1670 wallcycle_stop(wcycle, ewcNS);
1673 if (inputrec->implicit_solvent && bNS)
1675 make_gb_nblist(cr, inputrec->gb_algorithm,
1676 x, box, fr, &top->idef, graph, fr->born);
1679 if (DOMAINDECOMP(cr) && !(cr->duty & DUTY_PME))
1681 wallcycle_start(wcycle, ewcPPDURINGPME);
1682 dd_force_flop_start(cr->dd, nrnb);
1687 /* Enforced rotation has its own cycle counter that starts after the collective
1688 * coordinates have been communicated. It is added to ddCyclF to allow
1689 * for proper load-balancing */
1690 wallcycle_start(wcycle, ewcROT);
1691 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1692 wallcycle_stop(wcycle, ewcROT);
1695 /* Start the force cycle counter.
1696 * This counter is stopped in do_forcelow_level.
1697 * No parallel communication should occur while this counter is running,
1698 * since that will interfere with the dynamic load balancing.
1700 wallcycle_start(wcycle, ewcFORCE);
1704 /* Reset forces for which the virial is calculated separately:
1705 * PME/Ewald forces if necessary */
1706 if (fr->bF_NoVirSum)
1708 if (flags & GMX_FORCE_VIRIAL)
1710 fr->f_novirsum = fr->f_novirsum_alloc;
1713 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1717 clear_rvecs(homenr, fr->f_novirsum+start);
1722 /* We are not calculating the pressure so we do not need
1723 * a separate array for forces that do not contribute
1730 /* Clear the short- and long-range forces */
1731 clear_rvecs(fr->natoms_force_constr, f);
1732 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1734 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1737 clear_rvec(fr->vir_diag_posres);
1739 if (inputrec->ePull == epullCONSTRAINT)
1741 clear_pull_forces(inputrec->pull);
1744 /* update QMMMrec, if necessary */
1747 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1750 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1752 posres_wrapper(fplog, flags, bSepDVDL, inputrec, nrnb, top, box, x,
1756 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_FBPOSRES].nr > 0)
1758 fbposres_wrapper(inputrec, nrnb, top, box, x, enerd, fr);
1761 /* Compute the bonded and non-bonded energies and optionally forces */
1762 do_force_lowlevel(fplog, step, fr, inputrec, &(top->idef),
1763 cr, nrnb, wcycle, mdatoms,
1764 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, top, fr->born,
1765 &(top->atomtypes), bBornRadii, box,
1766 inputrec->fepvals, lambda,
1767 graph, &(top->excls), fr->mu_tot,
1773 if (do_per_step(step, inputrec->nstcalclr))
1775 /* Add the long range forces to the short range forces */
1776 for (i = 0; i < fr->natoms_force_constr; i++)
1778 rvec_add(fr->f_twin[i], f[i], f[i]);
1783 cycles_force = wallcycle_stop(wcycle, ewcFORCE);
1787 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1790 if (DOMAINDECOMP(cr))
1792 dd_force_flop_stop(cr->dd, nrnb);
1795 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1801 if (IR_ELEC_FIELD(*inputrec))
1803 /* Compute forces due to electric field */
1804 calc_f_el(MASTER(cr) ? field : NULL,
1805 start, homenr, mdatoms->chargeA, fr->f_novirsum,
1806 inputrec->ex, inputrec->et, t);
1809 if (bDoAdressWF && fr->adress_icor == eAdressICThermoForce)
1811 /* Compute thermodynamic force in hybrid AdResS region */
1812 adress_thermo_force(start, homenr, &(top->cgs), x, fr->f_novirsum, fr, mdatoms,
1813 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1816 /* Communicate the forces */
1817 if (DOMAINDECOMP(cr))
1819 wallcycle_start(wcycle, ewcMOVEF);
1820 dd_move_f(cr->dd, f, fr->fshift);
1821 /* Do we need to communicate the separate force array
1822 * for terms that do not contribute to the single sum virial?
1823 * Position restraints and electric fields do not introduce
1824 * inter-cg forces, only full electrostatics methods do.
1825 * When we do not calculate the virial, fr->f_novirsum = f,
1826 * so we have already communicated these forces.
1828 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1829 (flags & GMX_FORCE_VIRIAL))
1831 dd_move_f(cr->dd, fr->f_novirsum, NULL);
1835 /* We should not update the shift forces here,
1836 * since f_twin is already included in f.
1838 dd_move_f(cr->dd, fr->f_twin, NULL);
1840 wallcycle_stop(wcycle, ewcMOVEF);
1843 /* If we have NoVirSum forces, but we do not calculate the virial,
1844 * we sum fr->f_novirum=f later.
1846 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1848 wallcycle_start(wcycle, ewcVSITESPREAD);
1849 spread_vsite_f(vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1850 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1851 wallcycle_stop(wcycle, ewcVSITESPREAD);
1855 wallcycle_start(wcycle, ewcVSITESPREAD);
1856 spread_vsite_f(vsite, x, fr->f_twin, NULL, FALSE, NULL,
1858 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1859 wallcycle_stop(wcycle, ewcVSITESPREAD);
1863 if (flags & GMX_FORCE_VIRIAL)
1865 /* Calculation of the virial must be done after vsites! */
1866 calc_virial(0, mdatoms->homenr, x, f,
1867 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1871 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1873 pull_potential_wrapper(fplog, bSepDVDL, cr, inputrec, box, x,
1874 f, vir_force, mdatoms, enerd, lambda, t);
1877 /* Add the forces from enforced rotation potentials (if any) */
1880 wallcycle_start(wcycle, ewcROTadd);
1881 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
1882 wallcycle_stop(wcycle, ewcROTadd);
1885 if (PAR(cr) && !(cr->duty & DUTY_PME))
1887 /* In case of node-splitting, the PP nodes receive the long-range
1888 * forces, virial and energy from the PME nodes here.
1890 pme_receive_force_ener(fplog, bSepDVDL, cr, wcycle, enerd, fr);
1895 post_process_forces(cr, step, nrnb, wcycle,
1896 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
1900 /* Sum the potential energy terms from group contributions */
1901 sum_epot(&(enerd->grpp), enerd->term);
1904 void do_force(FILE *fplog, t_commrec *cr,
1905 t_inputrec *inputrec,
1906 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1907 gmx_localtop_t *top,
1908 gmx_groups_t *groups,
1909 matrix box, rvec x[], history_t *hist,
1913 gmx_enerdata_t *enerd, t_fcdata *fcd,
1914 real *lambda, t_graph *graph,
1916 gmx_vsite_t *vsite, rvec mu_tot,
1917 double t, FILE *field, gmx_edsam_t ed,
1918 gmx_bool bBornRadii,
1921 /* modify force flag if not doing nonbonded */
1922 if (!fr->bNonbonded)
1924 flags &= ~GMX_FORCE_NONBONDED;
1927 switch (inputrec->cutoff_scheme)
1930 do_force_cutsVERLET(fplog, cr, inputrec,
1946 do_force_cutsGROUP(fplog, cr, inputrec,
1961 gmx_incons("Invalid cut-off scheme passed!");
1966 void do_constrain_first(FILE *fplog, gmx_constr_t constr,
1967 t_inputrec *ir, t_mdatoms *md,
1968 t_state *state, t_commrec *cr, t_nrnb *nrnb,
1969 t_forcerec *fr, gmx_localtop_t *top)
1971 int i, m, start, end;
1973 real dt = ir->delta_t;
1977 snew(savex, state->natoms);
1984 fprintf(debug, "vcm: start=%d, homenr=%d, end=%d\n",
1985 start, md->homenr, end);
1987 /* Do a first constrain to reset particles... */
1988 step = ir->init_step;
1991 char buf[STEPSTRSIZE];
1992 fprintf(fplog, "\nConstraining the starting coordinates (step %s)\n",
1993 gmx_step_str(step, buf));
1997 /* constrain the current position */
1998 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
1999 ir, NULL, cr, step, 0, md,
2000 state->x, state->x, NULL,
2001 fr->bMolPBC, state->box,
2002 state->lambda[efptBONDED], &dvdl_dum,
2003 NULL, NULL, nrnb, econqCoord,
2004 ir->epc == epcMTTK, state->veta, state->veta);
2007 /* constrain the inital velocity, and save it */
2008 /* also may be useful if we need the ekin from the halfstep for velocity verlet */
2009 /* might not yet treat veta correctly */
2010 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2011 ir, NULL, cr, step, 0, md,
2012 state->x, state->v, state->v,
2013 fr->bMolPBC, state->box,
2014 state->lambda[efptBONDED], &dvdl_dum,
2015 NULL, NULL, nrnb, econqVeloc,
2016 ir->epc == epcMTTK, state->veta, state->veta);
2018 /* constrain the inital velocities at t-dt/2 */
2019 if (EI_STATE_VELOCITY(ir->eI) && ir->eI != eiVV)
2021 for (i = start; (i < end); i++)
2023 for (m = 0; (m < DIM); m++)
2025 /* Reverse the velocity */
2026 state->v[i][m] = -state->v[i][m];
2027 /* Store the position at t-dt in buf */
2028 savex[i][m] = state->x[i][m] + dt*state->v[i][m];
2031 /* Shake the positions at t=-dt with the positions at t=0
2032 * as reference coordinates.
2036 char buf[STEPSTRSIZE];
2037 fprintf(fplog, "\nConstraining the coordinates at t0-dt (step %s)\n",
2038 gmx_step_str(step, buf));
2041 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2042 ir, NULL, cr, step, -1, md,
2043 state->x, savex, NULL,
2044 fr->bMolPBC, state->box,
2045 state->lambda[efptBONDED], &dvdl_dum,
2046 state->v, NULL, nrnb, econqCoord,
2047 ir->epc == epcMTTK, state->veta, state->veta);
2049 for (i = start; i < end; i++)
2051 for (m = 0; m < DIM; m++)
2053 /* Re-reverse the velocities */
2054 state->v[i][m] = -state->v[i][m];
2063 integrate_table(real vdwtab[], real scale, int offstart, int rstart, int rend,
2064 double *enerout, double *virout)
2066 double enersum, virsum;
2067 double invscale, invscale2, invscale3;
2068 double r, ea, eb, ec, pa, pb, pc, pd;
2070 int ri, offset, tabfactor;
2072 invscale = 1.0/scale;
2073 invscale2 = invscale*invscale;
2074 invscale3 = invscale*invscale2;
2076 /* Following summation derived from cubic spline definition,
2077 * Numerical Recipies in C, second edition, p. 113-116. Exact for
2078 * the cubic spline. We first calculate the negative of the
2079 * energy from rvdw to rvdw_switch, assuming that g(r)=1, and then
2080 * add the more standard, abrupt cutoff correction to that result,
2081 * yielding the long-range correction for a switched function. We
2082 * perform both the pressure and energy loops at the same time for
2083 * simplicity, as the computational cost is low. */
2087 /* Since the dispersion table has been scaled down a factor
2088 * 6.0 and the repulsion a factor 12.0 to compensate for the
2089 * c6/c12 parameters inside nbfp[] being scaled up (to save
2090 * flops in kernels), we need to correct for this.
2101 for (ri = rstart; ri < rend; ++ri)
2105 eb = 2.0*invscale2*r;
2109 pb = 3.0*invscale2*r;
2110 pc = 3.0*invscale*r*r;
2113 /* this "8" is from the packing in the vdwtab array - perhaps
2114 should be defined? */
2116 offset = 8*ri + offstart;
2117 y0 = vdwtab[offset];
2118 f = vdwtab[offset+1];
2119 g = vdwtab[offset+2];
2120 h = vdwtab[offset+3];
2122 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);
2123 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);
2125 *enerout = 4.0*M_PI*enersum*tabfactor;
2126 *virout = 4.0*M_PI*virsum*tabfactor;
2129 void calc_enervirdiff(FILE *fplog, int eDispCorr, t_forcerec *fr)
2131 double eners[2], virs[2], enersum, virsum, y0, f, g, h;
2132 double r0, r1, r, rc3, rc9, ea, eb, ec, pa, pb, pc, pd;
2133 double invscale, invscale2, invscale3;
2134 int ri0, ri1, ri, i, offstart, offset;
2135 real scale, *vdwtab, tabfactor, tmp;
2137 fr->enershiftsix = 0;
2138 fr->enershifttwelve = 0;
2139 fr->enerdiffsix = 0;
2140 fr->enerdifftwelve = 0;
2142 fr->virdifftwelve = 0;
2144 if (eDispCorr != edispcNO)
2146 for (i = 0; i < 2; i++)
2151 if (fr->vdwtype == evdwSWITCH || fr->vdwtype == evdwSHIFT ||
2152 fr->vdw_modifier == eintmodPOTSWITCH ||
2153 fr->vdw_modifier == eintmodFORCESWITCH)
2155 if (fr->rvdw_switch == 0)
2158 "With dispersion correction rvdw-switch can not be zero "
2159 "for vdw-type = %s", evdw_names[fr->vdwtype]);
2162 scale = fr->nblists[0].table_elec_vdw.scale;
2163 vdwtab = fr->nblists[0].table_vdw.data;
2165 /* Round the cut-offs to exact table values for precision */
2166 ri0 = floor(fr->rvdw_switch*scale);
2167 ri1 = ceil(fr->rvdw*scale);
2173 if (fr->vdwtype == evdwSHIFT ||
2174 fr->vdw_modifier == eintmodFORCESWITCH)
2176 /* Determine the constant energy shift below rvdw_switch.
2177 * Table has a scale factor since we have scaled it down to compensate
2178 * for scaling-up c6/c12 with the derivative factors to save flops in analytical kernels.
2180 fr->enershiftsix = (real)(-1.0/(rc3*rc3)) - 6.0*vdwtab[8*ri0];
2181 fr->enershifttwelve = (real)( 1.0/(rc9*rc3)) - 12.0*vdwtab[8*ri0 + 4];
2183 /* Add the constant part from 0 to rvdw_switch.
2184 * This integration from 0 to rvdw_switch overcounts the number
2185 * of interactions by 1, as it also counts the self interaction.
2186 * We will correct for this later.
2188 eners[0] += 4.0*M_PI*fr->enershiftsix*rc3/3.0;
2189 eners[1] += 4.0*M_PI*fr->enershifttwelve*rc3/3.0;
2190 for (i = 0; i < 2; i++)
2194 integrate_table(vdwtab, scale, (i == 0 ? 0 : 4), ri0, ri1, &enersum, &virsum);
2195 eners[i] -= enersum;
2199 /* now add the correction for rvdw_switch to infinity */
2200 eners[0] += -4.0*M_PI/(3.0*rc3);
2201 eners[1] += 4.0*M_PI/(9.0*rc9);
2202 virs[0] += 8.0*M_PI/rc3;
2203 virs[1] += -16.0*M_PI/(3.0*rc9);
2205 else if (fr->vdwtype == evdwCUT ||
2206 EVDW_PME(fr->vdwtype) ||
2207 fr->vdwtype == evdwUSER)
2209 if (fr->vdwtype == evdwUSER && fplog)
2212 "WARNING: using dispersion correction with user tables\n");
2215 /* Note that with LJ-PME, the dispersion correction is multiplied
2216 * by the difference between the actual C6 and the value of C6
2217 * that would produce the combination rule.
2218 * This means the normal energy and virial difference formulas
2222 rc3 = fr->rvdw*fr->rvdw*fr->rvdw;
2224 /* Contribution beyond the cut-off */
2225 eners[0] += -4.0*M_PI/(3.0*rc3);
2226 eners[1] += 4.0*M_PI/(9.0*rc9);
2227 if (fr->vdw_modifier == eintmodPOTSHIFT)
2229 /* Contribution within the cut-off */
2230 eners[0] += -4.0*M_PI/(3.0*rc3);
2231 eners[1] += 4.0*M_PI/(3.0*rc9);
2233 /* Contribution beyond the cut-off */
2234 virs[0] += 8.0*M_PI/rc3;
2235 virs[1] += -16.0*M_PI/(3.0*rc9);
2240 "Dispersion correction is not implemented for vdw-type = %s",
2241 evdw_names[fr->vdwtype]);
2244 /* TODO: remove this code once we have group LJ-PME kernels
2245 * that calculate the exact, full LJ param C6/r^6 within the cut-off,
2246 * as the current nbnxn kernels do.
2248 if (fr->vdwtype == evdwPME && fr->cutoff_scheme == ecutsGROUP)
2250 /* Calculate self-interaction coefficient (assuming that
2251 * the reciprocal-space contribution is constant in the
2252 * region that contributes to the self-interaction).
2254 fr->enershiftsix = pow(fr->ewaldcoeff_lj, 6) / 6.0;
2256 eners[0] += -pow(sqrt(M_PI)*fr->ewaldcoeff_lj, 3)/3.0;
2257 virs[0] += pow(sqrt(M_PI)*fr->ewaldcoeff_lj, 3);
2260 fr->enerdiffsix = eners[0];
2261 fr->enerdifftwelve = eners[1];
2262 /* The 0.5 is due to the Gromacs definition of the virial */
2263 fr->virdiffsix = 0.5*virs[0];
2264 fr->virdifftwelve = 0.5*virs[1];
2268 void calc_dispcorr(FILE *fplog, t_inputrec *ir, t_forcerec *fr,
2269 gmx_int64_t step, int natoms,
2270 matrix box, real lambda, tensor pres, tensor virial,
2271 real *prescorr, real *enercorr, real *dvdlcorr)
2273 gmx_bool bCorrAll, bCorrPres;
2274 real dvdlambda, invvol, dens, ninter, avcsix, avctwelve, enerdiff, svir = 0, spres = 0;
2284 if (ir->eDispCorr != edispcNO)
2286 bCorrAll = (ir->eDispCorr == edispcAllEner ||
2287 ir->eDispCorr == edispcAllEnerPres);
2288 bCorrPres = (ir->eDispCorr == edispcEnerPres ||
2289 ir->eDispCorr == edispcAllEnerPres);
2291 invvol = 1/det(box);
2294 /* Only correct for the interactions with the inserted molecule */
2295 dens = (natoms - fr->n_tpi)*invvol;
2300 dens = natoms*invvol;
2301 ninter = 0.5*natoms;
2304 if (ir->efep == efepNO)
2306 avcsix = fr->avcsix[0];
2307 avctwelve = fr->avctwelve[0];
2311 avcsix = (1 - lambda)*fr->avcsix[0] + lambda*fr->avcsix[1];
2312 avctwelve = (1 - lambda)*fr->avctwelve[0] + lambda*fr->avctwelve[1];
2315 enerdiff = ninter*(dens*fr->enerdiffsix - fr->enershiftsix);
2316 *enercorr += avcsix*enerdiff;
2318 if (ir->efep != efepNO)
2320 dvdlambda += (fr->avcsix[1] - fr->avcsix[0])*enerdiff;
2324 enerdiff = ninter*(dens*fr->enerdifftwelve - fr->enershifttwelve);
2325 *enercorr += avctwelve*enerdiff;
2326 if (fr->efep != efepNO)
2328 dvdlambda += (fr->avctwelve[1] - fr->avctwelve[0])*enerdiff;
2334 svir = ninter*dens*avcsix*fr->virdiffsix/3.0;
2335 if (ir->eDispCorr == edispcAllEnerPres)
2337 svir += ninter*dens*avctwelve*fr->virdifftwelve/3.0;
2339 /* The factor 2 is because of the Gromacs virial definition */
2340 spres = -2.0*invvol*svir*PRESFAC;
2342 for (m = 0; m < DIM; m++)
2344 virial[m][m] += svir;
2345 pres[m][m] += spres;
2350 /* Can't currently control when it prints, for now, just print when degugging */
2355 fprintf(debug, "Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
2361 "Long Range LJ corr.: Epot %10g, Pres: %10g, Vir: %10g\n",
2362 *enercorr, spres, svir);
2366 fprintf(debug, "Long Range LJ corr.: Epot %10g\n", *enercorr);
2370 if (fr->bSepDVDL && do_per_step(step, ir->nstlog))
2372 gmx_print_sepdvdl(fplog, "Dispersion correction", *enercorr, dvdlambda);
2374 if (fr->efep != efepNO)
2376 *dvdlcorr += dvdlambda;
2381 void do_pbc_first(FILE *fplog, matrix box, t_forcerec *fr,
2382 t_graph *graph, rvec x[])
2386 fprintf(fplog, "Removing pbc first time\n");
2388 calc_shifts(box, fr->shift_vec);
2391 mk_mshift(fplog, graph, fr->ePBC, box, x);
2394 p_graph(debug, "do_pbc_first 1", graph);
2396 shift_self(graph, box, x);
2397 /* By doing an extra mk_mshift the molecules that are broken
2398 * because they were e.g. imported from another software
2399 * will be made whole again. Such are the healing powers
2402 mk_mshift(fplog, graph, fr->ePBC, box, x);
2405 p_graph(debug, "do_pbc_first 2", graph);
2410 fprintf(fplog, "Done rmpbc\n");
2414 static void low_do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2415 gmx_mtop_t *mtop, rvec x[],
2420 gmx_molblock_t *molb;
2422 if (bFirst && fplog)
2424 fprintf(fplog, "Removing pbc first time\n");
2429 for (mb = 0; mb < mtop->nmolblock; mb++)
2431 molb = &mtop->molblock[mb];
2432 if (molb->natoms_mol == 1 ||
2433 (!bFirst && mtop->moltype[molb->type].cgs.nr == 1))
2435 /* Just one atom or charge group in the molecule, no PBC required */
2436 as += molb->nmol*molb->natoms_mol;
2440 /* Pass NULL iso fplog to avoid graph prints for each molecule type */
2441 mk_graph_ilist(NULL, mtop->moltype[molb->type].ilist,
2442 0, molb->natoms_mol, FALSE, FALSE, graph);
2444 for (mol = 0; mol < molb->nmol; mol++)
2446 mk_mshift(fplog, graph, ePBC, box, x+as);
2448 shift_self(graph, box, x+as);
2449 /* The molecule is whole now.
2450 * We don't need the second mk_mshift call as in do_pbc_first,
2451 * since we no longer need this graph.
2454 as += molb->natoms_mol;
2462 void do_pbc_first_mtop(FILE *fplog, int ePBC, matrix box,
2463 gmx_mtop_t *mtop, rvec x[])
2465 low_do_pbc_mtop(fplog, ePBC, box, mtop, x, TRUE);
2468 void do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2469 gmx_mtop_t *mtop, rvec x[])
2471 low_do_pbc_mtop(fplog, ePBC, box, mtop, x, FALSE);
2474 void finish_run(FILE *fplog, t_commrec *cr,
2475 t_inputrec *inputrec,
2476 t_nrnb nrnb[], gmx_wallcycle_t wcycle,
2477 gmx_walltime_accounting_t walltime_accounting,
2478 wallclock_gpu_t *gputimes,
2479 gmx_bool bWriteStat)
2482 t_nrnb *nrnb_tot = NULL;
2485 double elapsed_time,
2486 elapsed_time_over_all_ranks,
2487 elapsed_time_over_all_threads,
2488 elapsed_time_over_all_threads_over_all_ranks;
2489 wallcycle_sum(cr, wcycle);
2495 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM,
2496 cr->mpi_comm_mysim);
2504 elapsed_time = walltime_accounting_get_elapsed_time(walltime_accounting);
2505 elapsed_time_over_all_ranks = elapsed_time;
2506 elapsed_time_over_all_threads = walltime_accounting_get_elapsed_time_over_all_threads(walltime_accounting);
2507 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
2511 /* reduce elapsed_time over all MPI ranks in the current simulation */
2512 MPI_Allreduce(&elapsed_time,
2513 &elapsed_time_over_all_ranks,
2514 1, MPI_DOUBLE, MPI_SUM,
2515 cr->mpi_comm_mysim);
2516 elapsed_time_over_all_ranks /= cr->nnodes;
2517 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
2518 * current simulation. */
2519 MPI_Allreduce(&elapsed_time_over_all_threads,
2520 &elapsed_time_over_all_threads_over_all_ranks,
2521 1, MPI_DOUBLE, MPI_SUM,
2522 cr->mpi_comm_mysim);
2528 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
2535 if ((cr->duty & DUTY_PP) && DOMAINDECOMP(cr))
2537 print_dd_statistics(cr, inputrec, fplog);
2542 wallcycle_print(fplog, cr->nnodes, cr->npmenodes,
2543 elapsed_time_over_all_ranks,
2546 if (EI_DYNAMICS(inputrec->eI))
2548 delta_t = inputrec->delta_t;
2557 print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks,
2558 elapsed_time_over_all_ranks,
2559 walltime_accounting_get_nsteps_done(walltime_accounting),
2560 delta_t, nbfs, mflop);
2564 print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks,
2565 elapsed_time_over_all_ranks,
2566 walltime_accounting_get_nsteps_done(walltime_accounting),
2567 delta_t, nbfs, mflop);
2572 extern void initialize_lambdas(FILE *fplog, t_inputrec *ir, int *fep_state, real *lambda, double *lam0)
2574 /* this function works, but could probably use a logic rewrite to keep all the different
2575 types of efep straight. */
2578 t_lambda *fep = ir->fepvals;
2580 if ((ir->efep == efepNO) && (ir->bSimTemp == FALSE))
2582 for (i = 0; i < efptNR; i++)
2594 *fep_state = fep->init_fep_state; /* this might overwrite the checkpoint
2595 if checkpoint is set -- a kludge is in for now
2597 for (i = 0; i < efptNR; i++)
2599 /* overwrite lambda state with init_lambda for now for backwards compatibility */
2600 if (fep->init_lambda >= 0) /* if it's -1, it was never initializd */
2602 lambda[i] = fep->init_lambda;
2605 lam0[i] = lambda[i];
2610 lambda[i] = fep->all_lambda[i][*fep_state];
2613 lam0[i] = lambda[i];
2619 /* need to rescale control temperatures to match current state */
2620 for (i = 0; i < ir->opts.ngtc; i++)
2622 if (ir->opts.ref_t[i] > 0)
2624 ir->opts.ref_t[i] = ir->simtempvals->temperatures[*fep_state];
2630 /* Send to the log the information on the current lambdas */
2633 fprintf(fplog, "Initial vector of lambda components:[ ");
2634 for (i = 0; i < efptNR; i++)
2636 fprintf(fplog, "%10.4f ", lambda[i]);
2638 fprintf(fplog, "]\n");
2644 void init_md(FILE *fplog,
2645 t_commrec *cr, t_inputrec *ir, const output_env_t oenv,
2646 double *t, double *t0,
2647 real *lambda, int *fep_state, double *lam0,
2648 t_nrnb *nrnb, gmx_mtop_t *mtop,
2650 int nfile, const t_filenm fnm[],
2651 gmx_mdoutf_t *outf, t_mdebin **mdebin,
2652 tensor force_vir, tensor shake_vir, rvec mu_tot,
2653 gmx_bool *bSimAnn, t_vcm **vcm, unsigned long Flags)
2658 /* Initial values */
2659 *t = *t0 = ir->init_t;
2662 for (i = 0; i < ir->opts.ngtc; i++)
2664 /* set bSimAnn if any group is being annealed */
2665 if (ir->opts.annealing[i] != eannNO)
2672 update_annealing_target_temp(&(ir->opts), ir->init_t);
2675 /* Initialize lambda variables */
2676 initialize_lambdas(fplog, ir, fep_state, lambda, lam0);
2680 *upd = init_update(ir);
2686 *vcm = init_vcm(fplog, &mtop->groups, ir);
2689 if (EI_DYNAMICS(ir->eI) && !(Flags & MD_APPENDFILES))
2691 if (ir->etc == etcBERENDSEN)
2693 please_cite(fplog, "Berendsen84a");
2695 if (ir->etc == etcVRESCALE)
2697 please_cite(fplog, "Bussi2007a");
2705 *outf = init_mdoutf(fplog, nfile, fnm, Flags, cr, ir, mtop, oenv);
2707 *mdebin = init_mdebin((Flags & MD_APPENDFILES) ? NULL : mdoutf_get_fp_ene(*outf),
2708 mtop, ir, mdoutf_get_fp_dhdl(*outf));
2713 please_cite(fplog, "Fritsch12");
2714 please_cite(fplog, "Junghans10");
2716 /* Initiate variables */
2717 clear_mat(force_vir);
2718 clear_mat(shake_vir);