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"
72 #include "gromacs/random/random.h"
75 #include "nbnxn_atomdata.h"
76 #include "nbnxn_search.h"
77 #include "nbnxn_kernels/nbnxn_kernel_ref.h"
78 #include "nbnxn_kernels/simd_4xn/nbnxn_kernel_simd_4xn.h"
79 #include "nbnxn_kernels/simd_2xnn/nbnxn_kernel_simd_2xnn.h"
80 #include "nbnxn_kernels/nbnxn_kernel_gpu_ref.h"
82 #include "gromacs/timing/wallcycle.h"
83 #include "gromacs/timing/walltime_accounting.h"
84 #include "gromacs/utility/gmxmpi.h"
85 #include "gromacs/essentialdynamics/edsam.h"
86 #include "gromacs/pulling/pull.h"
87 #include "gromacs/pulling/pull_rotation.h"
92 #include "nbnxn_cuda_data_mgmt.h"
93 #include "nbnxn_cuda/nbnxn_cuda.h"
95 void print_time(FILE *out,
96 gmx_walltime_accounting_t walltime_accounting,
99 t_commrec gmx_unused *cr)
102 char timebuf[STRLEN];
103 double dt, elapsed_seconds, time_per_step;
106 #ifndef GMX_THREAD_MPI
112 fprintf(out, "step %s", gmx_step_str(step, buf));
113 if ((step >= ir->nstlist))
115 double seconds_since_epoch = gmx_gettime();
116 elapsed_seconds = seconds_since_epoch - walltime_accounting_get_start_time_stamp(walltime_accounting);
117 time_per_step = elapsed_seconds/(step - ir->init_step + 1);
118 dt = (ir->nsteps + ir->init_step - step) * time_per_step;
124 finish = (time_t) (seconds_since_epoch + dt);
125 gmx_ctime_r(&finish, timebuf, STRLEN);
126 sprintf(buf, "%s", timebuf);
127 buf[strlen(buf)-1] = '\0';
128 fprintf(out, ", will finish %s", buf);
132 fprintf(out, ", remaining wall clock time: %5d s ", (int)dt);
137 fprintf(out, " performance: %.1f ns/day ",
138 ir->delta_t/1000*24*60*60/time_per_step);
141 #ifndef GMX_THREAD_MPI
151 void print_date_and_time(FILE *fplog, int nodeid, const char *title,
152 const gmx_walltime_accounting_t walltime_accounting)
155 char timebuf[STRLEN];
156 char time_string[STRLEN];
161 if (walltime_accounting != NULL)
163 tmptime = (time_t) walltime_accounting_get_start_time_stamp(walltime_accounting);
164 gmx_ctime_r(&tmptime, timebuf, STRLEN);
168 tmptime = (time_t) gmx_gettime();
169 gmx_ctime_r(&tmptime, timebuf, STRLEN);
171 for (i = 0; timebuf[i] >= ' '; i++)
173 time_string[i] = timebuf[i];
175 time_string[i] = '\0';
177 fprintf(fplog, "%s on node %d %s\n", title, nodeid, time_string);
181 void print_start(FILE *fplog, t_commrec *cr,
182 gmx_walltime_accounting_t walltime_accounting,
187 sprintf(buf, "Started %s", name);
188 print_date_and_time(fplog, cr->nodeid, buf, walltime_accounting);
191 static void sum_forces(int start, int end, rvec f[], rvec flr[])
197 pr_rvecs(debug, 0, "fsr", f+start, end-start);
198 pr_rvecs(debug, 0, "flr", flr+start, end-start);
200 for (i = start; (i < end); i++)
202 rvec_inc(f[i], flr[i]);
207 * calc_f_el calculates forces due to an electric field.
209 * force is kJ mol^-1 nm^-1 = e * kJ mol^-1 nm^-1 / e
211 * Et[] contains the parameters for the time dependent
212 * part of the field (not yet used).
213 * Ex[] contains the parameters for
214 * the spatial dependent part of the field. You can have cool periodic
215 * fields in principle, but only a constant field is supported
217 * The function should return the energy due to the electric field
218 * (if any) but for now returns 0.
221 * There can be problems with the virial.
222 * Since the field is not self-consistent this is unavoidable.
223 * For neutral molecules the virial is correct within this approximation.
224 * For neutral systems with many charged molecules the error is small.
225 * But for systems with a net charge or a few charged molecules
226 * the error can be significant when the field is high.
227 * Solution: implement a self-consitent electric field into PME.
229 static void calc_f_el(FILE *fp, int start, int homenr,
230 real charge[], rvec f[],
231 t_cosines Ex[], t_cosines Et[], double t)
237 for (m = 0; (m < DIM); m++)
244 Ext[m] = cos(Et[m].a[0]*(t-t0))*exp(-sqr(t-t0)/(2.0*sqr(Et[m].a[2])));
248 Ext[m] = cos(Et[m].a[0]*t);
257 /* Convert the field strength from V/nm to MD-units */
258 Ext[m] *= Ex[m].a[0]*FIELDFAC;
259 for (i = start; (i < start+homenr); i++)
261 f[i][m] += charge[i]*Ext[m];
271 fprintf(fp, "%10g %10g %10g %10g #FIELD\n", t,
272 Ext[XX]/FIELDFAC, Ext[YY]/FIELDFAC, Ext[ZZ]/FIELDFAC);
276 static void calc_virial(int start, int homenr, rvec x[], rvec f[],
277 tensor vir_part, t_graph *graph, matrix box,
278 t_nrnb *nrnb, const t_forcerec *fr, int ePBC)
283 /* The short-range virial from surrounding boxes */
285 calc_vir(SHIFTS, fr->shift_vec, fr->fshift, vir_part, ePBC == epbcSCREW, box);
286 inc_nrnb(nrnb, eNR_VIRIAL, SHIFTS);
288 /* Calculate partial virial, for local atoms only, based on short range.
289 * Total virial is computed in global_stat, called from do_md
291 f_calc_vir(start, start+homenr, x, f, vir_part, graph, box);
292 inc_nrnb(nrnb, eNR_VIRIAL, homenr);
294 /* Add position restraint contribution */
295 for (i = 0; i < DIM; i++)
297 vir_part[i][i] += fr->vir_diag_posres[i];
300 /* Add wall contribution */
301 for (i = 0; i < DIM; i++)
303 vir_part[i][ZZ] += fr->vir_wall_z[i];
308 pr_rvecs(debug, 0, "vir_part", vir_part, DIM);
312 static void posres_wrapper(FILE *fplog,
318 matrix box, rvec x[],
319 gmx_enerdata_t *enerd,
327 /* Position restraints always require full pbc */
328 set_pbc(&pbc, ir->ePBC, box);
330 v = posres(top->idef.il[F_POSRES].nr, top->idef.il[F_POSRES].iatoms,
331 top->idef.iparams_posres,
332 (const rvec*)x, fr->f_novirsum, fr->vir_diag_posres,
333 ir->ePBC == epbcNONE ? NULL : &pbc,
334 lambda[efptRESTRAINT], &dvdl,
335 fr->rc_scaling, fr->ePBC, fr->posres_com, fr->posres_comB);
338 gmx_print_sepdvdl(fplog, interaction_function[F_POSRES].longname, v, dvdl);
340 enerd->term[F_POSRES] += v;
341 /* If just the force constant changes, the FEP term is linear,
342 * but if k changes, it is not.
344 enerd->dvdl_nonlin[efptRESTRAINT] += dvdl;
345 inc_nrnb(nrnb, eNR_POSRES, top->idef.il[F_POSRES].nr/2);
347 if ((ir->fepvals->n_lambda > 0) && (flags & GMX_FORCE_DHDL))
349 for (i = 0; i < enerd->n_lambda; i++)
351 real dvdl_dum, lambda_dum;
353 lambda_dum = (i == 0 ? lambda[efptRESTRAINT] : ir->fepvals->all_lambda[efptRESTRAINT][i-1]);
354 v = posres(top->idef.il[F_POSRES].nr, top->idef.il[F_POSRES].iatoms,
355 top->idef.iparams_posres,
356 (const rvec*)x, NULL, NULL,
357 ir->ePBC == epbcNONE ? NULL : &pbc, lambda_dum, &dvdl,
358 fr->rc_scaling, fr->ePBC, fr->posres_com, fr->posres_comB);
359 enerd->enerpart_lambda[i] += v;
364 static void fbposres_wrapper(t_inputrec *ir,
367 matrix box, rvec x[],
368 gmx_enerdata_t *enerd,
374 /* Flat-bottomed position restraints always require full pbc */
375 set_pbc(&pbc, ir->ePBC, box);
376 v = fbposres(top->idef.il[F_FBPOSRES].nr, top->idef.il[F_FBPOSRES].iatoms,
377 top->idef.iparams_fbposres,
378 (const rvec*)x, fr->f_novirsum, fr->vir_diag_posres,
379 ir->ePBC == epbcNONE ? NULL : &pbc,
380 fr->rc_scaling, fr->ePBC, fr->posres_com);
381 enerd->term[F_FBPOSRES] += v;
382 inc_nrnb(nrnb, eNR_FBPOSRES, top->idef.il[F_FBPOSRES].nr/2);
385 static void pull_potential_wrapper(FILE *fplog,
389 matrix box, rvec x[],
393 gmx_enerdata_t *enerd,
400 /* Calculate the center of mass forces, this requires communication,
401 * which is why pull_potential is called close to other communication.
402 * The virial contribution is calculated directly,
403 * which is why we call pull_potential after calc_virial.
405 set_pbc(&pbc, ir->ePBC, box);
407 enerd->term[F_COM_PULL] +=
408 pull_potential(ir->ePull, ir->pull, mdatoms, &pbc,
409 cr, t, lambda[efptRESTRAINT], x, f, vir_force, &dvdl);
412 gmx_print_sepdvdl(fplog, "Com pull", enerd->term[F_COM_PULL], dvdl);
414 enerd->dvdl_lin[efptRESTRAINT] += dvdl;
417 static void pme_receive_force_ener(FILE *fplog,
420 gmx_wallcycle_t wcycle,
421 gmx_enerdata_t *enerd,
424 real e_q, e_lj, v, dvdl_q, dvdl_lj;
425 float cycles_ppdpme, cycles_seppme;
427 cycles_ppdpme = wallcycle_stop(wcycle, ewcPPDURINGPME);
428 dd_cycles_add(cr->dd, cycles_ppdpme, ddCyclPPduringPME);
430 /* In case of node-splitting, the PP nodes receive the long-range
431 * forces, virial and energy from the PME nodes here.
433 wallcycle_start(wcycle, ewcPP_PMEWAITRECVF);
436 gmx_pme_receive_f(cr, fr->f_novirsum, fr->vir_el_recip, &e_q,
437 fr->vir_lj_recip, &e_lj, &dvdl_q, &dvdl_lj,
441 gmx_print_sepdvdl(fplog, "Electrostatic PME mesh", e_q, dvdl_q);
442 gmx_print_sepdvdl(fplog, "Lennard-Jones PME mesh", e_lj, dvdl_lj);
444 enerd->term[F_COUL_RECIP] += e_q;
445 enerd->term[F_LJ_RECIP] += e_lj;
446 enerd->dvdl_lin[efptCOUL] += dvdl_q;
447 enerd->dvdl_lin[efptVDW] += dvdl_lj;
451 dd_cycles_add(cr->dd, cycles_seppme, ddCyclPME);
453 wallcycle_stop(wcycle, ewcPP_PMEWAITRECVF);
456 static void print_large_forces(FILE *fp, t_mdatoms *md, t_commrec *cr,
457 gmx_int64_t step, real pforce, rvec *x, rvec *f)
461 char buf[STEPSTRSIZE];
464 for (i = 0; i < md->homenr; i++)
467 /* We also catch NAN, if the compiler does not optimize this away. */
468 if (fn2 >= pf2 || fn2 != fn2)
470 fprintf(fp, "step %s atom %6d x %8.3f %8.3f %8.3f force %12.5e\n",
471 gmx_step_str(step, buf),
472 ddglatnr(cr->dd, i), x[i][XX], x[i][YY], x[i][ZZ], sqrt(fn2));
477 static void post_process_forces(t_commrec *cr,
479 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
481 matrix box, rvec x[],
486 t_forcerec *fr, gmx_vsite_t *vsite,
493 /* Spread the mesh force on virtual sites to the other particles...
494 * This is parallellized. MPI communication is performed
495 * if the constructing atoms aren't local.
497 wallcycle_start(wcycle, ewcVSITESPREAD);
498 spread_vsite_f(vsite, x, fr->f_novirsum, NULL,
499 (flags & GMX_FORCE_VIRIAL), fr->vir_el_recip,
501 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
502 wallcycle_stop(wcycle, ewcVSITESPREAD);
504 if (flags & GMX_FORCE_VIRIAL)
506 /* Now add the forces, this is local */
509 sum_forces(0, fr->f_novirsum_n, f, fr->f_novirsum);
513 sum_forces(0, mdatoms->homenr,
516 if (EEL_FULL(fr->eeltype))
518 /* Add the mesh contribution to the virial */
519 m_add(vir_force, fr->vir_el_recip, vir_force);
521 if (EVDW_PME(fr->vdwtype))
523 /* Add the mesh contribution to the virial */
524 m_add(vir_force, fr->vir_lj_recip, vir_force);
528 pr_rvecs(debug, 0, "vir_force", vir_force, DIM);
533 if (fr->print_force >= 0)
535 print_large_forces(stderr, mdatoms, cr, step, fr->print_force, x, f);
539 static void do_nb_verlet(t_forcerec *fr,
540 interaction_const_t *ic,
541 gmx_enerdata_t *enerd,
542 int flags, int ilocality,
545 gmx_wallcycle_t wcycle)
547 int nnbl, kernel_type, enr_nbnxn_kernel_ljc, enr_nbnxn_kernel_lj;
549 nonbonded_verlet_group_t *nbvg;
552 if (!(flags & GMX_FORCE_NONBONDED))
554 /* skip non-bonded calculation */
558 nbvg = &fr->nbv->grp[ilocality];
560 /* CUDA kernel launch overhead is already timed separately */
561 if (fr->cutoff_scheme != ecutsVERLET)
563 gmx_incons("Invalid cut-off scheme passed!");
566 bCUDA = (nbvg->kernel_type == nbnxnk8x8x8_CUDA);
570 wallcycle_sub_start(wcycle, ewcsNONBONDED);
572 switch (nbvg->kernel_type)
574 case nbnxnk4x4_PlainC:
575 nbnxn_kernel_ref(&nbvg->nbl_lists,
581 enerd->grpp.ener[egCOULSR],
583 enerd->grpp.ener[egBHAMSR] :
584 enerd->grpp.ener[egLJSR]);
587 case nbnxnk4xN_SIMD_4xN:
588 nbnxn_kernel_simd_4xn(&nbvg->nbl_lists,
595 enerd->grpp.ener[egCOULSR],
597 enerd->grpp.ener[egBHAMSR] :
598 enerd->grpp.ener[egLJSR]);
600 case nbnxnk4xN_SIMD_2xNN:
601 nbnxn_kernel_simd_2xnn(&nbvg->nbl_lists,
608 enerd->grpp.ener[egCOULSR],
610 enerd->grpp.ener[egBHAMSR] :
611 enerd->grpp.ener[egLJSR]);
614 case nbnxnk8x8x8_CUDA:
615 nbnxn_cuda_launch_kernel(fr->nbv->cu_nbv, nbvg->nbat, flags, ilocality);
618 case nbnxnk8x8x8_PlainC:
619 nbnxn_kernel_gpu_ref(nbvg->nbl_lists.nbl[0],
624 nbvg->nbat->out[0].f,
626 enerd->grpp.ener[egCOULSR],
628 enerd->grpp.ener[egBHAMSR] :
629 enerd->grpp.ener[egLJSR]);
633 gmx_incons("Invalid nonbonded kernel type passed!");
638 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
641 if (EEL_RF(ic->eeltype) || ic->eeltype == eelCUT)
643 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_RF;
645 else if ((!bCUDA && nbvg->ewald_excl == ewaldexclAnalytical) ||
646 (bCUDA && nbnxn_cuda_is_kernel_ewald_analytical(fr->nbv->cu_nbv)))
648 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_EWALD;
652 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_TAB;
654 enr_nbnxn_kernel_lj = eNR_NBNXN_LJ;
655 if (flags & GMX_FORCE_ENERGY)
657 /* In eNR_??? the nbnxn F+E kernels are always the F kernel + 1 */
658 enr_nbnxn_kernel_ljc += 1;
659 enr_nbnxn_kernel_lj += 1;
662 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc,
663 nbvg->nbl_lists.natpair_ljq);
664 inc_nrnb(nrnb, enr_nbnxn_kernel_lj,
665 nbvg->nbl_lists.natpair_lj);
666 /* The Coulomb-only kernels are offset -eNR_NBNXN_LJ_RF+eNR_NBNXN_RF */
667 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc-eNR_NBNXN_LJ_RF+eNR_NBNXN_RF,
668 nbvg->nbl_lists.natpair_q);
670 if (ic->vdw_modifier == eintmodFORCESWITCH)
672 /* We add up the switch cost separately */
673 inc_nrnb(nrnb, eNR_NBNXN_LJ_FSW+((flags & GMX_FORCE_ENERGY) ? 1 : 0),
674 nbvg->nbl_lists.natpair_ljq + nbvg->nbl_lists.natpair_lj);
676 if (ic->vdw_modifier == eintmodPOTSWITCH)
678 /* We add up the switch cost separately */
679 inc_nrnb(nrnb, eNR_NBNXN_LJ_PSW+((flags & GMX_FORCE_ENERGY) ? 1 : 0),
680 nbvg->nbl_lists.natpair_ljq + nbvg->nbl_lists.natpair_lj);
684 void do_force_cutsVERLET(FILE *fplog, t_commrec *cr,
685 t_inputrec *inputrec,
686 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
688 gmx_groups_t gmx_unused *groups,
689 matrix box, rvec x[], history_t *hist,
693 gmx_enerdata_t *enerd, t_fcdata *fcd,
694 real *lambda, t_graph *graph,
695 t_forcerec *fr, interaction_const_t *ic,
696 gmx_vsite_t *vsite, rvec mu_tot,
697 double t, FILE *field, gmx_edsam_t ed,
705 gmx_bool bSepDVDL, bStateChanged, bNS, bFillGrid, bCalcCGCM, bBS;
706 gmx_bool bDoLongRange, bDoForces, bSepLRF, bUseGPU, bUseOrEmulGPU;
707 gmx_bool bDiffKernels = FALSE;
709 rvec vzero, box_diag;
711 float cycles_pme, cycles_force, cycles_wait_gpu;
712 nonbonded_verlet_t *nbv;
717 nb_kernel_type = fr->nbv->grp[0].kernel_type;
720 homenr = mdatoms->homenr;
722 bSepDVDL = (fr->bSepDVDL && do_per_step(step, inputrec->nstlog));
724 clear_mat(vir_force);
727 if (DOMAINDECOMP(cr))
729 cg1 = cr->dd->ncg_tot;
740 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
741 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
742 bFillGrid = (bNS && bStateChanged);
743 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
744 bDoLongRange = (fr->bTwinRange && bNS && (flags & GMX_FORCE_DO_LR));
745 bDoForces = (flags & GMX_FORCE_FORCES);
746 bSepLRF = (bDoLongRange && bDoForces && (flags & GMX_FORCE_SEPLRF));
747 bUseGPU = fr->nbv->bUseGPU;
748 bUseOrEmulGPU = bUseGPU || (nbv->grp[0].kernel_type == nbnxnk8x8x8_PlainC);
752 update_forcerec(fr, box);
754 if (NEED_MUTOT(*inputrec))
756 /* Calculate total (local) dipole moment in a temporary common array.
757 * This makes it possible to sum them over nodes faster.
759 calc_mu(start, homenr,
760 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
765 if (fr->ePBC != epbcNONE)
767 /* Compute shift vectors every step,
768 * because of pressure coupling or box deformation!
770 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
772 calc_shifts(box, fr->shift_vec);
777 put_atoms_in_box_omp(fr->ePBC, box, homenr, x);
778 inc_nrnb(nrnb, eNR_SHIFTX, homenr);
780 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
782 unshift_self(graph, box, x);
786 nbnxn_atomdata_copy_shiftvec(flags & GMX_FORCE_DYNAMICBOX,
787 fr->shift_vec, nbv->grp[0].nbat);
790 if (!(cr->duty & DUTY_PME))
792 /* Send particle coordinates to the pme nodes.
793 * Since this is only implemented for domain decomposition
794 * and domain decomposition does not use the graph,
795 * we do not need to worry about shifting.
800 wallcycle_start(wcycle, ewcPP_PMESENDX);
802 bBS = (inputrec->nwall == 2);
806 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
809 if (EEL_PME(fr->eeltype))
811 pme_flags |= GMX_PME_DO_COULOMB;
814 if (EVDW_PME(fr->vdwtype))
816 pme_flags |= GMX_PME_DO_LJ;
817 if (fr->ljpme_combination_rule == eljpmeLB)
819 pme_flags |= GMX_PME_LJ_LB;
823 gmx_pme_send_coordinates(cr, bBS ? boxs : box, x,
824 mdatoms->nChargePerturbed, mdatoms->nTypePerturbed, lambda[efptCOUL], lambda[efptVDW],
825 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),
828 wallcycle_stop(wcycle, ewcPP_PMESENDX);
832 /* do gridding for pair search */
835 if (graph && bStateChanged)
837 /* Calculate intramolecular shift vectors to make molecules whole */
838 mk_mshift(fplog, graph, fr->ePBC, box, x);
842 box_diag[XX] = box[XX][XX];
843 box_diag[YY] = box[YY][YY];
844 box_diag[ZZ] = box[ZZ][ZZ];
846 wallcycle_start(wcycle, ewcNS);
849 wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
850 nbnxn_put_on_grid(nbv->nbs, fr->ePBC, box,
852 0, mdatoms->homenr, -1, fr->cginfo, x,
854 nbv->grp[eintLocal].kernel_type,
855 nbv->grp[eintLocal].nbat);
856 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
860 wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
861 nbnxn_put_on_grid_nonlocal(nbv->nbs, domdec_zones(cr->dd),
863 nbv->grp[eintNonlocal].kernel_type,
864 nbv->grp[eintNonlocal].nbat);
865 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
868 if (nbv->ngrp == 1 ||
869 nbv->grp[eintNonlocal].nbat == nbv->grp[eintLocal].nbat)
871 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatAll,
872 nbv->nbs, mdatoms, fr->cginfo);
876 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatLocal,
877 nbv->nbs, mdatoms, fr->cginfo);
878 nbnxn_atomdata_set(nbv->grp[eintNonlocal].nbat, eatAll,
879 nbv->nbs, mdatoms, fr->cginfo);
881 wallcycle_stop(wcycle, ewcNS);
884 /* initialize the GPU atom data and copy shift vector */
889 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
890 nbnxn_cuda_init_atomdata(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
891 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
894 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
895 nbnxn_cuda_upload_shiftvec(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
896 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
899 /* do local pair search */
902 wallcycle_start_nocount(wcycle, ewcNS);
903 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
904 nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintLocal].nbat,
907 nbv->min_ci_balanced,
908 &nbv->grp[eintLocal].nbl_lists,
910 nbv->grp[eintLocal].kernel_type,
912 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
916 /* initialize local pair-list on the GPU */
917 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
918 nbv->grp[eintLocal].nbl_lists.nbl[0],
921 wallcycle_stop(wcycle, ewcNS);
925 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
926 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
927 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, FALSE, x,
928 nbv->grp[eintLocal].nbat);
929 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
930 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
935 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
936 /* launch local nonbonded F on GPU */
937 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFNo,
939 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
942 /* Communicate coordinates and sum dipole if necessary +
943 do non-local pair search */
944 if (DOMAINDECOMP(cr))
946 bDiffKernels = (nbv->grp[eintNonlocal].kernel_type !=
947 nbv->grp[eintLocal].kernel_type);
951 /* With GPU+CPU non-bonded calculations we need to copy
952 * the local coordinates to the non-local nbat struct
953 * (in CPU format) as the non-local kernel call also
954 * calculates the local - non-local interactions.
956 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
957 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
958 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, TRUE, x,
959 nbv->grp[eintNonlocal].nbat);
960 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
961 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
966 wallcycle_start_nocount(wcycle, ewcNS);
967 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
971 nbnxn_grid_add_simple(nbv->nbs, nbv->grp[eintNonlocal].nbat);
974 nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintNonlocal].nbat,
977 nbv->min_ci_balanced,
978 &nbv->grp[eintNonlocal].nbl_lists,
980 nbv->grp[eintNonlocal].kernel_type,
983 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
985 if (nbv->grp[eintNonlocal].kernel_type == nbnxnk8x8x8_CUDA)
987 /* initialize non-local pair-list on the GPU */
988 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
989 nbv->grp[eintNonlocal].nbl_lists.nbl[0],
992 wallcycle_stop(wcycle, ewcNS);
996 wallcycle_start(wcycle, ewcMOVEX);
997 dd_move_x(cr->dd, box, x);
999 /* When we don't need the total dipole we sum it in global_stat */
1000 if (bStateChanged && NEED_MUTOT(*inputrec))
1002 gmx_sumd(2*DIM, mu, cr);
1004 wallcycle_stop(wcycle, ewcMOVEX);
1006 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1007 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
1008 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatNonlocal, FALSE, x,
1009 nbv->grp[eintNonlocal].nbat);
1010 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
1011 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1014 if (bUseGPU && !bDiffKernels)
1016 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
1017 /* launch non-local nonbonded F on GPU */
1018 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFNo,
1020 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1026 /* launch D2H copy-back F */
1027 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1028 if (DOMAINDECOMP(cr) && !bDiffKernels)
1030 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintNonlocal].nbat,
1031 flags, eatNonlocal);
1033 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintLocal].nbat,
1035 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1038 if (bStateChanged && NEED_MUTOT(*inputrec))
1042 gmx_sumd(2*DIM, mu, cr);
1045 for (i = 0; i < 2; i++)
1047 for (j = 0; j < DIM; j++)
1049 fr->mu_tot[i][j] = mu[i*DIM + j];
1053 if (fr->efep == efepNO)
1055 copy_rvec(fr->mu_tot[0], mu_tot);
1059 for (j = 0; j < DIM; j++)
1062 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] +
1063 lambda[efptCOUL]*fr->mu_tot[1][j];
1067 /* Reset energies */
1068 reset_enerdata(fr, bNS, enerd, MASTER(cr));
1069 clear_rvecs(SHIFTS, fr->fshift);
1071 if (DOMAINDECOMP(cr) && !(cr->duty & DUTY_PME))
1073 wallcycle_start(wcycle, ewcPPDURINGPME);
1074 dd_force_flop_start(cr->dd, nrnb);
1079 /* Enforced rotation has its own cycle counter that starts after the collective
1080 * coordinates have been communicated. It is added to ddCyclF to allow
1081 * for proper load-balancing */
1082 wallcycle_start(wcycle, ewcROT);
1083 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1084 wallcycle_stop(wcycle, ewcROT);
1087 /* Start the force cycle counter.
1088 * This counter is stopped in do_forcelow_level.
1089 * No parallel communication should occur while this counter is running,
1090 * since that will interfere with the dynamic load balancing.
1092 wallcycle_start(wcycle, ewcFORCE);
1095 /* Reset forces for which the virial is calculated separately:
1096 * PME/Ewald forces if necessary */
1097 if (fr->bF_NoVirSum)
1099 if (flags & GMX_FORCE_VIRIAL)
1101 fr->f_novirsum = fr->f_novirsum_alloc;
1104 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1108 clear_rvecs(homenr, fr->f_novirsum+start);
1113 /* We are not calculating the pressure so we do not need
1114 * a separate array for forces that do not contribute
1121 /* Clear the short- and long-range forces */
1122 clear_rvecs(fr->natoms_force_constr, f);
1123 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1125 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1128 clear_rvec(fr->vir_diag_posres);
1131 if (inputrec->ePull == epullCONSTRAINT)
1133 clear_pull_forces(inputrec->pull);
1136 /* We calculate the non-bonded forces, when done on the CPU, here.
1137 * We do this before calling do_force_lowlevel, as in there bondeds
1138 * forces are calculated before PME, which does communication.
1139 * With this order, non-bonded and bonded force calculation imbalance
1140 * can be balanced out by the domain decomposition load balancing.
1145 /* Maybe we should move this into do_force_lowlevel */
1146 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFYes,
1150 if (!bUseOrEmulGPU || bDiffKernels)
1154 if (DOMAINDECOMP(cr))
1156 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal,
1157 bDiffKernels ? enbvClearFYes : enbvClearFNo,
1167 aloc = eintNonlocal;
1170 /* Add all the non-bonded force to the normal force array.
1171 * This can be split into a local a non-local part when overlapping
1172 * communication with calculation with domain decomposition.
1174 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1175 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1176 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1177 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatAll, nbv->grp[aloc].nbat, f);
1178 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1179 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1180 wallcycle_start_nocount(wcycle, ewcFORCE);
1182 /* if there are multiple fshift output buffers reduce them */
1183 if ((flags & GMX_FORCE_VIRIAL) &&
1184 nbv->grp[aloc].nbl_lists.nnbl > 1)
1186 nbnxn_atomdata_add_nbat_fshift_to_fshift(nbv->grp[aloc].nbat,
1191 /* update QMMMrec, if necessary */
1194 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1197 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1199 posres_wrapper(fplog, flags, bSepDVDL, inputrec, nrnb, top, box, x,
1203 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_FBPOSRES].nr > 0)
1205 fbposres_wrapper(inputrec, nrnb, top, box, x, enerd, fr);
1208 /* Compute the bonded and non-bonded energies and optionally forces */
1209 do_force_lowlevel(fplog, step, fr, inputrec, &(top->idef),
1210 cr, nrnb, wcycle, mdatoms,
1211 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, top, fr->born,
1212 &(top->atomtypes), bBornRadii, box,
1213 inputrec->fepvals, lambda, graph, &(top->excls), fr->mu_tot,
1214 flags, &cycles_pme);
1218 if (do_per_step(step, inputrec->nstcalclr))
1220 /* Add the long range forces to the short range forces */
1221 for (i = 0; i < fr->natoms_force_constr; i++)
1223 rvec_add(fr->f_twin[i], f[i], f[i]);
1228 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1232 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1235 if (bUseOrEmulGPU && !bDiffKernels)
1237 /* wait for non-local forces (or calculate in emulation mode) */
1238 if (DOMAINDECOMP(cr))
1244 wallcycle_start(wcycle, ewcWAIT_GPU_NB_NL);
1245 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1246 nbv->grp[eintNonlocal].nbat,
1248 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1250 cycles_tmp = wallcycle_stop(wcycle, ewcWAIT_GPU_NB_NL);
1251 cycles_wait_gpu += cycles_tmp;
1252 cycles_force += cycles_tmp;
1256 wallcycle_start_nocount(wcycle, ewcFORCE);
1257 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFYes,
1259 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1261 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1262 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1263 /* skip the reduction if there was no non-local work to do */
1264 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1266 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatNonlocal,
1267 nbv->grp[eintNonlocal].nbat, f);
1269 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1270 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1274 if (bDoForces && DOMAINDECOMP(cr))
1276 /* Communicate the forces */
1277 wallcycle_start(wcycle, ewcMOVEF);
1278 dd_move_f(cr->dd, f, fr->fshift);
1279 /* Do we need to communicate the separate force array
1280 * for terms that do not contribute to the single sum virial?
1281 * Position restraints and electric fields do not introduce
1282 * inter-cg forces, only full electrostatics methods do.
1283 * When we do not calculate the virial, fr->f_novirsum = f,
1284 * so we have already communicated these forces.
1286 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1287 (flags & GMX_FORCE_VIRIAL))
1289 dd_move_f(cr->dd, fr->f_novirsum, NULL);
1293 /* We should not update the shift forces here,
1294 * since f_twin is already included in f.
1296 dd_move_f(cr->dd, fr->f_twin, NULL);
1298 wallcycle_stop(wcycle, ewcMOVEF);
1303 /* wait for local forces (or calculate in emulation mode) */
1306 wallcycle_start(wcycle, ewcWAIT_GPU_NB_L);
1307 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1308 nbv->grp[eintLocal].nbat,
1310 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1312 cycles_wait_gpu += wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
1314 /* now clear the GPU outputs while we finish the step on the CPU */
1316 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1317 nbnxn_cuda_clear_outputs(nbv->cu_nbv, flags);
1318 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1322 wallcycle_start_nocount(wcycle, ewcFORCE);
1323 do_nb_verlet(fr, ic, enerd, flags, eintLocal,
1324 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
1326 wallcycle_stop(wcycle, ewcFORCE);
1328 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1329 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1330 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1332 /* skip the reduction if there was no non-local work to do */
1333 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatLocal,
1334 nbv->grp[eintLocal].nbat, f);
1336 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1337 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1340 if (DOMAINDECOMP(cr))
1342 dd_force_flop_stop(cr->dd, nrnb);
1345 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1348 dd_cycles_add(cr->dd, cycles_wait_gpu, ddCyclWaitGPU);
1355 if (IR_ELEC_FIELD(*inputrec))
1357 /* Compute forces due to electric field */
1358 calc_f_el(MASTER(cr) ? field : NULL,
1359 start, homenr, mdatoms->chargeA, fr->f_novirsum,
1360 inputrec->ex, inputrec->et, t);
1363 /* If we have NoVirSum forces, but we do not calculate the virial,
1364 * we sum fr->f_novirum=f later.
1366 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1368 wallcycle_start(wcycle, ewcVSITESPREAD);
1369 spread_vsite_f(vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1370 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1371 wallcycle_stop(wcycle, ewcVSITESPREAD);
1375 wallcycle_start(wcycle, ewcVSITESPREAD);
1376 spread_vsite_f(vsite, x, fr->f_twin, NULL, FALSE, NULL,
1378 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1379 wallcycle_stop(wcycle, ewcVSITESPREAD);
1383 if (flags & GMX_FORCE_VIRIAL)
1385 /* Calculation of the virial must be done after vsites! */
1386 calc_virial(0, mdatoms->homenr, x, f,
1387 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1391 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1393 pull_potential_wrapper(fplog, bSepDVDL, cr, inputrec, box, x,
1394 f, vir_force, mdatoms, enerd, lambda, t);
1397 /* Add the forces from enforced rotation potentials (if any) */
1400 wallcycle_start(wcycle, ewcROTadd);
1401 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
1402 wallcycle_stop(wcycle, ewcROTadd);
1405 if (PAR(cr) && !(cr->duty & DUTY_PME))
1407 /* In case of node-splitting, the PP nodes receive the long-range
1408 * forces, virial and energy from the PME nodes here.
1410 pme_receive_force_ener(fplog, bSepDVDL, cr, wcycle, enerd, fr);
1415 post_process_forces(cr, step, nrnb, wcycle,
1416 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
1420 /* Sum the potential energy terms from group contributions */
1421 sum_epot(&(enerd->grpp), enerd->term);
1424 void do_force_cutsGROUP(FILE *fplog, t_commrec *cr,
1425 t_inputrec *inputrec,
1426 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1427 gmx_localtop_t *top,
1428 gmx_groups_t *groups,
1429 matrix box, rvec x[], history_t *hist,
1433 gmx_enerdata_t *enerd, t_fcdata *fcd,
1434 real *lambda, t_graph *graph,
1435 t_forcerec *fr, gmx_vsite_t *vsite, rvec mu_tot,
1436 double t, FILE *field, gmx_edsam_t ed,
1437 gmx_bool bBornRadii,
1443 gmx_bool bSepDVDL, bStateChanged, bNS, bFillGrid, bCalcCGCM, bBS;
1444 gmx_bool bDoLongRangeNS, bDoForces, bDoPotential, bSepLRF;
1445 gmx_bool bDoAdressWF;
1447 rvec vzero, box_diag;
1448 real e, v, dvdlambda[efptNR];
1450 float cycles_pme, cycles_force;
1453 homenr = mdatoms->homenr;
1455 bSepDVDL = (fr->bSepDVDL && do_per_step(step, inputrec->nstlog));
1457 clear_mat(vir_force);
1460 if (DOMAINDECOMP(cr))
1462 cg1 = cr->dd->ncg_tot;
1473 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
1474 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
1475 /* Should we update the long-range neighborlists at this step? */
1476 bDoLongRangeNS = fr->bTwinRange && bNS;
1477 /* Should we perform the long-range nonbonded evaluation inside the neighborsearching? */
1478 bFillGrid = (bNS && bStateChanged);
1479 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
1480 bDoForces = (flags & GMX_FORCE_FORCES);
1481 bDoPotential = (flags & GMX_FORCE_ENERGY);
1482 bSepLRF = ((inputrec->nstcalclr > 1) && bDoForces &&
1483 (flags & GMX_FORCE_SEPLRF) && (flags & GMX_FORCE_DO_LR));
1485 /* should probably move this to the forcerec since it doesn't change */
1486 bDoAdressWF = ((fr->adress_type != eAdressOff));
1490 update_forcerec(fr, box);
1492 if (NEED_MUTOT(*inputrec))
1494 /* Calculate total (local) dipole moment in a temporary common array.
1495 * This makes it possible to sum them over nodes faster.
1497 calc_mu(start, homenr,
1498 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
1503 if (fr->ePBC != epbcNONE)
1505 /* Compute shift vectors every step,
1506 * because of pressure coupling or box deformation!
1508 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
1510 calc_shifts(box, fr->shift_vec);
1515 put_charge_groups_in_box(fplog, cg0, cg1, fr->ePBC, box,
1516 &(top->cgs), x, fr->cg_cm);
1517 inc_nrnb(nrnb, eNR_CGCM, homenr);
1518 inc_nrnb(nrnb, eNR_RESETX, cg1-cg0);
1520 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
1522 unshift_self(graph, box, x);
1527 calc_cgcm(fplog, cg0, cg1, &(top->cgs), x, fr->cg_cm);
1528 inc_nrnb(nrnb, eNR_CGCM, homenr);
1531 if (bCalcCGCM && gmx_debug_at)
1533 pr_rvecs(debug, 0, "cgcm", fr->cg_cm, top->cgs.nr);
1537 if (!(cr->duty & DUTY_PME))
1539 /* Send particle coordinates to the pme nodes.
1540 * Since this is only implemented for domain decomposition
1541 * and domain decomposition does not use the graph,
1542 * we do not need to worry about shifting.
1547 wallcycle_start(wcycle, ewcPP_PMESENDX);
1549 bBS = (inputrec->nwall == 2);
1552 copy_mat(box, boxs);
1553 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
1556 if (EEL_PME(fr->eeltype))
1558 pme_flags |= GMX_PME_DO_COULOMB;
1561 if (EVDW_PME(fr->vdwtype))
1563 pme_flags |= GMX_PME_DO_LJ;
1564 if (fr->ljpme_combination_rule == eljpmeLB)
1566 pme_flags |= GMX_PME_LJ_LB;
1570 gmx_pme_send_coordinates(cr, bBS ? boxs : box, x,
1571 mdatoms->nChargePerturbed, mdatoms->nTypePerturbed, lambda[efptCOUL], lambda[efptVDW],
1572 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),
1575 wallcycle_stop(wcycle, ewcPP_PMESENDX);
1577 #endif /* GMX_MPI */
1579 /* Communicate coordinates and sum dipole if necessary */
1580 if (DOMAINDECOMP(cr))
1582 wallcycle_start(wcycle, ewcMOVEX);
1583 dd_move_x(cr->dd, box, x);
1584 wallcycle_stop(wcycle, ewcMOVEX);
1587 /* update adress weight beforehand */
1588 if (bStateChanged && bDoAdressWF)
1590 /* need pbc for adress weight calculation with pbc_dx */
1591 set_pbc(&pbc, inputrec->ePBC, box);
1592 if (fr->adress_site == eAdressSITEcog)
1594 update_adress_weights_cog(top->idef.iparams, top->idef.il, x, fr, mdatoms,
1595 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1597 else if (fr->adress_site == eAdressSITEcom)
1599 update_adress_weights_com(fplog, cg0, cg1, &(top->cgs), x, fr, mdatoms,
1600 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1602 else if (fr->adress_site == eAdressSITEatomatom)
1604 update_adress_weights_atom_per_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1605 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1609 update_adress_weights_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1610 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1614 if (NEED_MUTOT(*inputrec))
1621 gmx_sumd(2*DIM, mu, cr);
1623 for (i = 0; i < 2; i++)
1625 for (j = 0; j < DIM; j++)
1627 fr->mu_tot[i][j] = mu[i*DIM + j];
1631 if (fr->efep == efepNO)
1633 copy_rvec(fr->mu_tot[0], mu_tot);
1637 for (j = 0; j < DIM; j++)
1640 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] + lambda[efptCOUL]*fr->mu_tot[1][j];
1645 /* Reset energies */
1646 reset_enerdata(fr, bNS, enerd, MASTER(cr));
1647 clear_rvecs(SHIFTS, fr->fshift);
1651 wallcycle_start(wcycle, ewcNS);
1653 if (graph && bStateChanged)
1655 /* Calculate intramolecular shift vectors to make molecules whole */
1656 mk_mshift(fplog, graph, fr->ePBC, box, x);
1659 /* Do the actual neighbour searching */
1661 groups, top, mdatoms,
1662 cr, nrnb, bFillGrid,
1665 wallcycle_stop(wcycle, ewcNS);
1668 if (inputrec->implicit_solvent && bNS)
1670 make_gb_nblist(cr, inputrec->gb_algorithm,
1671 x, box, fr, &top->idef, graph, fr->born);
1674 if (DOMAINDECOMP(cr) && !(cr->duty & DUTY_PME))
1676 wallcycle_start(wcycle, ewcPPDURINGPME);
1677 dd_force_flop_start(cr->dd, nrnb);
1682 /* Enforced rotation has its own cycle counter that starts after the collective
1683 * coordinates have been communicated. It is added to ddCyclF to allow
1684 * for proper load-balancing */
1685 wallcycle_start(wcycle, ewcROT);
1686 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1687 wallcycle_stop(wcycle, ewcROT);
1690 /* Start the force cycle counter.
1691 * This counter is stopped in do_forcelow_level.
1692 * No parallel communication should occur while this counter is running,
1693 * since that will interfere with the dynamic load balancing.
1695 wallcycle_start(wcycle, ewcFORCE);
1699 /* Reset forces for which the virial is calculated separately:
1700 * PME/Ewald forces if necessary */
1701 if (fr->bF_NoVirSum)
1703 if (flags & GMX_FORCE_VIRIAL)
1705 fr->f_novirsum = fr->f_novirsum_alloc;
1708 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1712 clear_rvecs(homenr, fr->f_novirsum+start);
1717 /* We are not calculating the pressure so we do not need
1718 * a separate array for forces that do not contribute
1725 /* Clear the short- and long-range forces */
1726 clear_rvecs(fr->natoms_force_constr, f);
1727 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1729 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1732 clear_rvec(fr->vir_diag_posres);
1734 if (inputrec->ePull == epullCONSTRAINT)
1736 clear_pull_forces(inputrec->pull);
1739 /* update QMMMrec, if necessary */
1742 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1745 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1747 posres_wrapper(fplog, flags, bSepDVDL, inputrec, nrnb, top, box, x,
1751 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_FBPOSRES].nr > 0)
1753 fbposres_wrapper(inputrec, nrnb, top, box, x, enerd, fr);
1756 /* Compute the bonded and non-bonded energies and optionally forces */
1757 do_force_lowlevel(fplog, step, fr, inputrec, &(top->idef),
1758 cr, nrnb, wcycle, mdatoms,
1759 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, top, fr->born,
1760 &(top->atomtypes), bBornRadii, box,
1761 inputrec->fepvals, lambda,
1762 graph, &(top->excls), fr->mu_tot,
1768 if (do_per_step(step, inputrec->nstcalclr))
1770 /* Add the long range forces to the short range forces */
1771 for (i = 0; i < fr->natoms_force_constr; i++)
1773 rvec_add(fr->f_twin[i], f[i], f[i]);
1778 cycles_force = wallcycle_stop(wcycle, ewcFORCE);
1782 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1785 if (DOMAINDECOMP(cr))
1787 dd_force_flop_stop(cr->dd, nrnb);
1790 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1796 if (IR_ELEC_FIELD(*inputrec))
1798 /* Compute forces due to electric field */
1799 calc_f_el(MASTER(cr) ? field : NULL,
1800 start, homenr, mdatoms->chargeA, fr->f_novirsum,
1801 inputrec->ex, inputrec->et, t);
1804 if (bDoAdressWF && fr->adress_icor == eAdressICThermoForce)
1806 /* Compute thermodynamic force in hybrid AdResS region */
1807 adress_thermo_force(start, homenr, &(top->cgs), x, fr->f_novirsum, fr, mdatoms,
1808 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1811 /* Communicate the forces */
1812 if (DOMAINDECOMP(cr))
1814 wallcycle_start(wcycle, ewcMOVEF);
1815 dd_move_f(cr->dd, f, fr->fshift);
1816 /* Do we need to communicate the separate force array
1817 * for terms that do not contribute to the single sum virial?
1818 * Position restraints and electric fields do not introduce
1819 * inter-cg forces, only full electrostatics methods do.
1820 * When we do not calculate the virial, fr->f_novirsum = f,
1821 * so we have already communicated these forces.
1823 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1824 (flags & GMX_FORCE_VIRIAL))
1826 dd_move_f(cr->dd, fr->f_novirsum, NULL);
1830 /* We should not update the shift forces here,
1831 * since f_twin is already included in f.
1833 dd_move_f(cr->dd, fr->f_twin, NULL);
1835 wallcycle_stop(wcycle, ewcMOVEF);
1838 /* If we have NoVirSum forces, but we do not calculate the virial,
1839 * we sum fr->f_novirum=f later.
1841 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1843 wallcycle_start(wcycle, ewcVSITESPREAD);
1844 spread_vsite_f(vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1845 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1846 wallcycle_stop(wcycle, ewcVSITESPREAD);
1850 wallcycle_start(wcycle, ewcVSITESPREAD);
1851 spread_vsite_f(vsite, x, fr->f_twin, NULL, FALSE, NULL,
1853 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1854 wallcycle_stop(wcycle, ewcVSITESPREAD);
1858 if (flags & GMX_FORCE_VIRIAL)
1860 /* Calculation of the virial must be done after vsites! */
1861 calc_virial(0, mdatoms->homenr, x, f,
1862 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1866 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1868 pull_potential_wrapper(fplog, bSepDVDL, cr, inputrec, box, x,
1869 f, vir_force, mdatoms, enerd, lambda, t);
1872 /* Add the forces from enforced rotation potentials (if any) */
1875 wallcycle_start(wcycle, ewcROTadd);
1876 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
1877 wallcycle_stop(wcycle, ewcROTadd);
1880 if (PAR(cr) && !(cr->duty & DUTY_PME))
1882 /* In case of node-splitting, the PP nodes receive the long-range
1883 * forces, virial and energy from the PME nodes here.
1885 pme_receive_force_ener(fplog, bSepDVDL, cr, wcycle, enerd, fr);
1890 post_process_forces(cr, step, nrnb, wcycle,
1891 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
1895 /* Sum the potential energy terms from group contributions */
1896 sum_epot(&(enerd->grpp), enerd->term);
1899 void do_force(FILE *fplog, t_commrec *cr,
1900 t_inputrec *inputrec,
1901 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1902 gmx_localtop_t *top,
1903 gmx_groups_t *groups,
1904 matrix box, rvec x[], history_t *hist,
1908 gmx_enerdata_t *enerd, t_fcdata *fcd,
1909 real *lambda, t_graph *graph,
1911 gmx_vsite_t *vsite, rvec mu_tot,
1912 double t, FILE *field, gmx_edsam_t ed,
1913 gmx_bool bBornRadii,
1916 /* modify force flag if not doing nonbonded */
1917 if (!fr->bNonbonded)
1919 flags &= ~GMX_FORCE_NONBONDED;
1922 switch (inputrec->cutoff_scheme)
1925 do_force_cutsVERLET(fplog, cr, inputrec,
1941 do_force_cutsGROUP(fplog, cr, inputrec,
1956 gmx_incons("Invalid cut-off scheme passed!");
1961 void do_constrain_first(FILE *fplog, gmx_constr_t constr,
1962 t_inputrec *ir, t_mdatoms *md,
1963 t_state *state, t_commrec *cr, t_nrnb *nrnb,
1964 t_forcerec *fr, gmx_localtop_t *top)
1966 int i, m, start, end;
1968 real dt = ir->delta_t;
1972 snew(savex, state->natoms);
1979 fprintf(debug, "vcm: start=%d, homenr=%d, end=%d\n",
1980 start, md->homenr, end);
1982 /* Do a first constrain to reset particles... */
1983 step = ir->init_step;
1986 char buf[STEPSTRSIZE];
1987 fprintf(fplog, "\nConstraining the starting coordinates (step %s)\n",
1988 gmx_step_str(step, buf));
1992 /* constrain the current position */
1993 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
1994 ir, NULL, cr, step, 0, md,
1995 state->x, state->x, NULL,
1996 fr->bMolPBC, state->box,
1997 state->lambda[efptBONDED], &dvdl_dum,
1998 NULL, NULL, nrnb, econqCoord,
1999 ir->epc == epcMTTK, state->veta, state->veta);
2002 /* constrain the inital velocity, and save it */
2003 /* also may be useful if we need the ekin from the halfstep for velocity verlet */
2004 /* might not yet treat veta correctly */
2005 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2006 ir, NULL, cr, step, 0, md,
2007 state->x, state->v, state->v,
2008 fr->bMolPBC, state->box,
2009 state->lambda[efptBONDED], &dvdl_dum,
2010 NULL, NULL, nrnb, econqVeloc,
2011 ir->epc == epcMTTK, state->veta, state->veta);
2013 /* constrain the inital velocities at t-dt/2 */
2014 if (EI_STATE_VELOCITY(ir->eI) && ir->eI != eiVV)
2016 for (i = start; (i < end); i++)
2018 for (m = 0; (m < DIM); m++)
2020 /* Reverse the velocity */
2021 state->v[i][m] = -state->v[i][m];
2022 /* Store the position at t-dt in buf */
2023 savex[i][m] = state->x[i][m] + dt*state->v[i][m];
2026 /* Shake the positions at t=-dt with the positions at t=0
2027 * as reference coordinates.
2031 char buf[STEPSTRSIZE];
2032 fprintf(fplog, "\nConstraining the coordinates at t0-dt (step %s)\n",
2033 gmx_step_str(step, buf));
2036 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2037 ir, NULL, cr, step, -1, md,
2038 state->x, savex, NULL,
2039 fr->bMolPBC, state->box,
2040 state->lambda[efptBONDED], &dvdl_dum,
2041 state->v, NULL, nrnb, econqCoord,
2042 ir->epc == epcMTTK, state->veta, state->veta);
2044 for (i = start; i < end; i++)
2046 for (m = 0; m < DIM; m++)
2048 /* Re-reverse the velocities */
2049 state->v[i][m] = -state->v[i][m];
2058 integrate_table(real vdwtab[], real scale, int offstart, int rstart, int rend,
2059 double *enerout, double *virout)
2061 double enersum, virsum;
2062 double invscale, invscale2, invscale3;
2063 double r, ea, eb, ec, pa, pb, pc, pd;
2065 int ri, offset, tabfactor;
2067 invscale = 1.0/scale;
2068 invscale2 = invscale*invscale;
2069 invscale3 = invscale*invscale2;
2071 /* Following summation derived from cubic spline definition,
2072 * Numerical Recipies in C, second edition, p. 113-116. Exact for
2073 * the cubic spline. We first calculate the negative of the
2074 * energy from rvdw to rvdw_switch, assuming that g(r)=1, and then
2075 * add the more standard, abrupt cutoff correction to that result,
2076 * yielding the long-range correction for a switched function. We
2077 * perform both the pressure and energy loops at the same time for
2078 * simplicity, as the computational cost is low. */
2082 /* Since the dispersion table has been scaled down a factor
2083 * 6.0 and the repulsion a factor 12.0 to compensate for the
2084 * c6/c12 parameters inside nbfp[] being scaled up (to save
2085 * flops in kernels), we need to correct for this.
2096 for (ri = rstart; ri < rend; ++ri)
2100 eb = 2.0*invscale2*r;
2104 pb = 3.0*invscale2*r;
2105 pc = 3.0*invscale*r*r;
2108 /* this "8" is from the packing in the vdwtab array - perhaps
2109 should be defined? */
2111 offset = 8*ri + offstart;
2112 y0 = vdwtab[offset];
2113 f = vdwtab[offset+1];
2114 g = vdwtab[offset+2];
2115 h = vdwtab[offset+3];
2117 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);
2118 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);
2120 *enerout = 4.0*M_PI*enersum*tabfactor;
2121 *virout = 4.0*M_PI*virsum*tabfactor;
2124 void calc_enervirdiff(FILE *fplog, int eDispCorr, t_forcerec *fr)
2126 double eners[2], virs[2], enersum, virsum, y0, f, g, h;
2127 double r0, r1, r, rc3, rc9, ea, eb, ec, pa, pb, pc, pd;
2128 double invscale, invscale2, invscale3;
2129 int ri0, ri1, ri, i, offstart, offset;
2130 real scale, *vdwtab, tabfactor, tmp;
2132 fr->enershiftsix = 0;
2133 fr->enershifttwelve = 0;
2134 fr->enerdiffsix = 0;
2135 fr->enerdifftwelve = 0;
2137 fr->virdifftwelve = 0;
2139 if (eDispCorr != edispcNO)
2141 for (i = 0; i < 2; i++)
2146 if (fr->vdwtype == evdwSWITCH || fr->vdwtype == evdwSHIFT ||
2147 fr->vdw_modifier == eintmodPOTSWITCH ||
2148 fr->vdw_modifier == eintmodFORCESWITCH)
2150 if (fr->rvdw_switch == 0)
2153 "With dispersion correction rvdw-switch can not be zero "
2154 "for vdw-type = %s", evdw_names[fr->vdwtype]);
2157 scale = fr->nblists[0].table_elec_vdw.scale;
2158 vdwtab = fr->nblists[0].table_vdw.data;
2160 /* Round the cut-offs to exact table values for precision */
2161 ri0 = floor(fr->rvdw_switch*scale);
2162 ri1 = ceil(fr->rvdw*scale);
2168 if (fr->vdwtype == evdwSHIFT ||
2169 fr->vdw_modifier == eintmodFORCESWITCH)
2171 /* Determine the constant energy shift below rvdw_switch.
2172 * Table has a scale factor since we have scaled it down to compensate
2173 * for scaling-up c6/c12 with the derivative factors to save flops in analytical kernels.
2175 fr->enershiftsix = (real)(-1.0/(rc3*rc3)) - 6.0*vdwtab[8*ri0];
2176 fr->enershifttwelve = (real)( 1.0/(rc9*rc3)) - 12.0*vdwtab[8*ri0 + 4];
2178 /* Add the constant part from 0 to rvdw_switch.
2179 * This integration from 0 to rvdw_switch overcounts the number
2180 * of interactions by 1, as it also counts the self interaction.
2181 * We will correct for this later.
2183 eners[0] += 4.0*M_PI*fr->enershiftsix*rc3/3.0;
2184 eners[1] += 4.0*M_PI*fr->enershifttwelve*rc3/3.0;
2185 for (i = 0; i < 2; i++)
2189 integrate_table(vdwtab, scale, (i == 0 ? 0 : 4), ri0, ri1, &enersum, &virsum);
2190 eners[i] -= enersum;
2194 /* now add the correction for rvdw_switch to infinity */
2195 eners[0] += -4.0*M_PI/(3.0*rc3);
2196 eners[1] += 4.0*M_PI/(9.0*rc9);
2197 virs[0] += 8.0*M_PI/rc3;
2198 virs[1] += -16.0*M_PI/(3.0*rc9);
2200 else if (EVDW_PME(fr->vdwtype))
2202 scale = fr->nblists[0].table_vdw.scale;
2203 vdwtab = fr->nblists[0].table_vdw.data;
2205 ri0 = floor(fr->rvdw_switch*scale);
2206 ri1 = ceil(fr->rvdw*scale);
2212 /* Calculate self-interaction coefficient (assuming that
2213 * the reciprocal-space contribution is constant in the
2214 * region that contributes to the self-interaction).
2216 fr->enershiftsix = pow(fr->ewaldcoeff_lj, 6) / 6.0;
2218 /* Add analytical corrections, C6 for the whole range, C12
2219 * from rvdw_switch to infinity.
2222 eners[0] += -pow(sqrt(M_PI)*fr->ewaldcoeff_lj, 3)/3.0;
2223 eners[1] += 4.0*M_PI/(9.0*rc9);
2224 virs[0] += pow(sqrt(M_PI)*fr->ewaldcoeff_lj, 3);
2225 virs[1] += -16.0*M_PI/(3.0*rc9);
2227 else if ((fr->vdwtype == evdwCUT) || (fr->vdwtype == evdwUSER))
2229 if (fr->vdwtype == evdwUSER && fplog)
2232 "WARNING: using dispersion correction with user tables\n");
2234 rc3 = fr->rvdw*fr->rvdw*fr->rvdw;
2236 /* Contribution beyond the cut-off */
2237 eners[0] += -4.0*M_PI/(3.0*rc3);
2238 eners[1] += 4.0*M_PI/(9.0*rc9);
2239 if (fr->vdw_modifier == eintmodPOTSHIFT)
2241 /* Contribution within the cut-off */
2242 eners[0] += -4.0*M_PI/(3.0*rc3);
2243 eners[1] += 4.0*M_PI/(3.0*rc9);
2245 /* Contribution beyond the cut-off */
2246 virs[0] += 8.0*M_PI/rc3;
2247 virs[1] += -16.0*M_PI/(3.0*rc9);
2252 "Dispersion correction is not implemented for vdw-type = %s",
2253 evdw_names[fr->vdwtype]);
2255 fr->enerdiffsix = eners[0];
2256 fr->enerdifftwelve = eners[1];
2257 /* The 0.5 is due to the Gromacs definition of the virial */
2258 fr->virdiffsix = 0.5*virs[0];
2259 fr->virdifftwelve = 0.5*virs[1];
2263 void calc_dispcorr(FILE *fplog, t_inputrec *ir, t_forcerec *fr,
2264 gmx_int64_t step, int natoms,
2265 matrix box, real lambda, tensor pres, tensor virial,
2266 real *prescorr, real *enercorr, real *dvdlcorr)
2268 gmx_bool bCorrAll, bCorrPres;
2269 real dvdlambda, invvol, dens, ninter, avcsix, avctwelve, enerdiff, svir = 0, spres = 0;
2279 if (ir->eDispCorr != edispcNO)
2281 bCorrAll = (ir->eDispCorr == edispcAllEner ||
2282 ir->eDispCorr == edispcAllEnerPres);
2283 bCorrPres = (ir->eDispCorr == edispcEnerPres ||
2284 ir->eDispCorr == edispcAllEnerPres);
2286 invvol = 1/det(box);
2289 /* Only correct for the interactions with the inserted molecule */
2290 dens = (natoms - fr->n_tpi)*invvol;
2295 dens = natoms*invvol;
2296 ninter = 0.5*natoms;
2299 if (ir->efep == efepNO)
2301 avcsix = fr->avcsix[0];
2302 avctwelve = fr->avctwelve[0];
2306 avcsix = (1 - lambda)*fr->avcsix[0] + lambda*fr->avcsix[1];
2307 avctwelve = (1 - lambda)*fr->avctwelve[0] + lambda*fr->avctwelve[1];
2310 enerdiff = ninter*(dens*fr->enerdiffsix - fr->enershiftsix);
2311 *enercorr += avcsix*enerdiff;
2313 if (ir->efep != efepNO)
2315 dvdlambda += (fr->avcsix[1] - fr->avcsix[0])*enerdiff;
2319 enerdiff = ninter*(dens*fr->enerdifftwelve - fr->enershifttwelve);
2320 *enercorr += avctwelve*enerdiff;
2321 if (fr->efep != efepNO)
2323 dvdlambda += (fr->avctwelve[1] - fr->avctwelve[0])*enerdiff;
2329 svir = ninter*dens*avcsix*fr->virdiffsix/3.0;
2330 if (ir->eDispCorr == edispcAllEnerPres)
2332 svir += ninter*dens*avctwelve*fr->virdifftwelve/3.0;
2334 /* The factor 2 is because of the Gromacs virial definition */
2335 spres = -2.0*invvol*svir*PRESFAC;
2337 for (m = 0; m < DIM; m++)
2339 virial[m][m] += svir;
2340 pres[m][m] += spres;
2345 /* Can't currently control when it prints, for now, just print when degugging */
2350 fprintf(debug, "Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
2356 "Long Range LJ corr.: Epot %10g, Pres: %10g, Vir: %10g\n",
2357 *enercorr, spres, svir);
2361 fprintf(debug, "Long Range LJ corr.: Epot %10g\n", *enercorr);
2365 if (fr->bSepDVDL && do_per_step(step, ir->nstlog))
2367 gmx_print_sepdvdl(fplog, "Dispersion correction", *enercorr, dvdlambda);
2369 if (fr->efep != efepNO)
2371 *dvdlcorr += dvdlambda;
2376 void do_pbc_first(FILE *fplog, matrix box, t_forcerec *fr,
2377 t_graph *graph, rvec x[])
2381 fprintf(fplog, "Removing pbc first time\n");
2383 calc_shifts(box, fr->shift_vec);
2386 mk_mshift(fplog, graph, fr->ePBC, box, x);
2389 p_graph(debug, "do_pbc_first 1", graph);
2391 shift_self(graph, box, x);
2392 /* By doing an extra mk_mshift the molecules that are broken
2393 * because they were e.g. imported from another software
2394 * will be made whole again. Such are the healing powers
2397 mk_mshift(fplog, graph, fr->ePBC, box, x);
2400 p_graph(debug, "do_pbc_first 2", graph);
2405 fprintf(fplog, "Done rmpbc\n");
2409 static void low_do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2410 gmx_mtop_t *mtop, rvec x[],
2415 gmx_molblock_t *molb;
2417 if (bFirst && fplog)
2419 fprintf(fplog, "Removing pbc first time\n");
2424 for (mb = 0; mb < mtop->nmolblock; mb++)
2426 molb = &mtop->molblock[mb];
2427 if (molb->natoms_mol == 1 ||
2428 (!bFirst && mtop->moltype[molb->type].cgs.nr == 1))
2430 /* Just one atom or charge group in the molecule, no PBC required */
2431 as += molb->nmol*molb->natoms_mol;
2435 /* Pass NULL iso fplog to avoid graph prints for each molecule type */
2436 mk_graph_ilist(NULL, mtop->moltype[molb->type].ilist,
2437 0, molb->natoms_mol, FALSE, FALSE, graph);
2439 for (mol = 0; mol < molb->nmol; mol++)
2441 mk_mshift(fplog, graph, ePBC, box, x+as);
2443 shift_self(graph, box, x+as);
2444 /* The molecule is whole now.
2445 * We don't need the second mk_mshift call as in do_pbc_first,
2446 * since we no longer need this graph.
2449 as += molb->natoms_mol;
2457 void do_pbc_first_mtop(FILE *fplog, int ePBC, matrix box,
2458 gmx_mtop_t *mtop, rvec x[])
2460 low_do_pbc_mtop(fplog, ePBC, box, mtop, x, TRUE);
2463 void do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2464 gmx_mtop_t *mtop, rvec x[])
2466 low_do_pbc_mtop(fplog, ePBC, box, mtop, x, FALSE);
2469 void finish_run(FILE *fplog, t_commrec *cr,
2470 t_inputrec *inputrec,
2471 t_nrnb nrnb[], gmx_wallcycle_t wcycle,
2472 gmx_walltime_accounting_t walltime_accounting,
2473 wallclock_gpu_t *gputimes,
2474 gmx_bool bWriteStat)
2477 t_nrnb *nrnb_tot = NULL;
2480 double elapsed_time,
2481 elapsed_time_over_all_ranks,
2482 elapsed_time_over_all_threads,
2483 elapsed_time_over_all_threads_over_all_ranks;
2484 wallcycle_sum(cr, wcycle);
2490 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM,
2491 cr->mpi_comm_mysim);
2499 elapsed_time = walltime_accounting_get_elapsed_time(walltime_accounting);
2500 elapsed_time_over_all_ranks = elapsed_time;
2501 elapsed_time_over_all_threads = walltime_accounting_get_elapsed_time_over_all_threads(walltime_accounting);
2502 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
2506 /* reduce elapsed_time over all MPI ranks in the current simulation */
2507 MPI_Allreduce(&elapsed_time,
2508 &elapsed_time_over_all_ranks,
2509 1, MPI_DOUBLE, MPI_SUM,
2510 cr->mpi_comm_mysim);
2511 elapsed_time_over_all_ranks /= cr->nnodes;
2512 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
2513 * current simulation. */
2514 MPI_Allreduce(&elapsed_time_over_all_threads,
2515 &elapsed_time_over_all_threads_over_all_ranks,
2516 1, MPI_DOUBLE, MPI_SUM,
2517 cr->mpi_comm_mysim);
2523 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
2530 if ((cr->duty & DUTY_PP) && DOMAINDECOMP(cr))
2532 print_dd_statistics(cr, inputrec, fplog);
2537 wallcycle_print(fplog, cr->nnodes, cr->npmenodes,
2538 elapsed_time_over_all_ranks,
2541 if (EI_DYNAMICS(inputrec->eI))
2543 delta_t = inputrec->delta_t;
2552 print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks,
2553 elapsed_time_over_all_ranks,
2554 walltime_accounting_get_nsteps_done(walltime_accounting),
2555 delta_t, nbfs, mflop);
2559 print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks,
2560 elapsed_time_over_all_ranks,
2561 walltime_accounting_get_nsteps_done(walltime_accounting),
2562 delta_t, nbfs, mflop);
2567 extern void initialize_lambdas(FILE *fplog, t_inputrec *ir, int *fep_state, real *lambda, double *lam0)
2569 /* this function works, but could probably use a logic rewrite to keep all the different
2570 types of efep straight. */
2573 t_lambda *fep = ir->fepvals;
2575 if ((ir->efep == efepNO) && (ir->bSimTemp == FALSE))
2577 for (i = 0; i < efptNR; i++)
2589 *fep_state = fep->init_fep_state; /* this might overwrite the checkpoint
2590 if checkpoint is set -- a kludge is in for now
2592 for (i = 0; i < efptNR; i++)
2594 /* overwrite lambda state with init_lambda for now for backwards compatibility */
2595 if (fep->init_lambda >= 0) /* if it's -1, it was never initializd */
2597 lambda[i] = fep->init_lambda;
2600 lam0[i] = lambda[i];
2605 lambda[i] = fep->all_lambda[i][*fep_state];
2608 lam0[i] = lambda[i];
2614 /* need to rescale control temperatures to match current state */
2615 for (i = 0; i < ir->opts.ngtc; i++)
2617 if (ir->opts.ref_t[i] > 0)
2619 ir->opts.ref_t[i] = ir->simtempvals->temperatures[*fep_state];
2625 /* Send to the log the information on the current lambdas */
2628 fprintf(fplog, "Initial vector of lambda components:[ ");
2629 for (i = 0; i < efptNR; i++)
2631 fprintf(fplog, "%10.4f ", lambda[i]);
2633 fprintf(fplog, "]\n");
2639 void init_md(FILE *fplog,
2640 t_commrec *cr, t_inputrec *ir, const output_env_t oenv,
2641 double *t, double *t0,
2642 real *lambda, int *fep_state, double *lam0,
2643 t_nrnb *nrnb, gmx_mtop_t *mtop,
2645 int nfile, const t_filenm fnm[],
2646 gmx_mdoutf_t *outf, t_mdebin **mdebin,
2647 tensor force_vir, tensor shake_vir, rvec mu_tot,
2648 gmx_bool *bSimAnn, t_vcm **vcm, unsigned long Flags)
2653 /* Initial values */
2654 *t = *t0 = ir->init_t;
2657 for (i = 0; i < ir->opts.ngtc; i++)
2659 /* set bSimAnn if any group is being annealed */
2660 if (ir->opts.annealing[i] != eannNO)
2667 update_annealing_target_temp(&(ir->opts), ir->init_t);
2670 /* Initialize lambda variables */
2671 initialize_lambdas(fplog, ir, fep_state, lambda, lam0);
2675 *upd = init_update(ir);
2681 *vcm = init_vcm(fplog, &mtop->groups, ir);
2684 if (EI_DYNAMICS(ir->eI) && !(Flags & MD_APPENDFILES))
2686 if (ir->etc == etcBERENDSEN)
2688 please_cite(fplog, "Berendsen84a");
2690 if (ir->etc == etcVRESCALE)
2692 please_cite(fplog, "Bussi2007a");
2700 *outf = init_mdoutf(fplog, nfile, fnm, Flags, cr, ir, mtop, oenv);
2702 *mdebin = init_mdebin((Flags & MD_APPENDFILES) ? NULL : mdoutf_get_fp_ene(*outf),
2703 mtop, ir, mdoutf_get_fp_dhdl(*outf));
2708 please_cite(fplog, "Fritsch12");
2709 please_cite(fplog, "Junghans10");
2711 /* Initiate variables */
2712 clear_mat(force_vir);
2713 clear_mat(shake_vir);