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_ADD_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_ADD_LJ_PSW+((flags & GMX_FORCE_ENERGY) ? 1 : 0),
680 nbvg->nbl_lists.natpair_ljq + nbvg->nbl_lists.natpair_lj);
682 if (ic->vdwtype == evdwPME)
684 /* We add up the LJ Ewald cost separately */
685 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_EWALD+((flags & GMX_FORCE_ENERGY) ? 1 : 0),
686 nbvg->nbl_lists.natpair_ljq + nbvg->nbl_lists.natpair_lj);
690 void do_force_cutsVERLET(FILE *fplog, t_commrec *cr,
691 t_inputrec *inputrec,
692 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
694 gmx_groups_t gmx_unused *groups,
695 matrix box, rvec x[], history_t *hist,
699 gmx_enerdata_t *enerd, t_fcdata *fcd,
700 real *lambda, t_graph *graph,
701 t_forcerec *fr, interaction_const_t *ic,
702 gmx_vsite_t *vsite, rvec mu_tot,
703 double t, FILE *field, gmx_edsam_t ed,
711 gmx_bool bSepDVDL, bStateChanged, bNS, bFillGrid, bCalcCGCM, bBS;
712 gmx_bool bDoLongRange, bDoForces, bSepLRF, bUseGPU, bUseOrEmulGPU;
713 gmx_bool bDiffKernels = FALSE;
715 rvec vzero, box_diag;
717 float cycles_pme, cycles_force, cycles_wait_gpu;
718 nonbonded_verlet_t *nbv;
723 nb_kernel_type = fr->nbv->grp[0].kernel_type;
726 homenr = mdatoms->homenr;
728 bSepDVDL = (fr->bSepDVDL && do_per_step(step, inputrec->nstlog));
730 clear_mat(vir_force);
733 if (DOMAINDECOMP(cr))
735 cg1 = cr->dd->ncg_tot;
746 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
747 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
748 bFillGrid = (bNS && bStateChanged);
749 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
750 bDoLongRange = (fr->bTwinRange && bNS && (flags & GMX_FORCE_DO_LR));
751 bDoForces = (flags & GMX_FORCE_FORCES);
752 bSepLRF = (bDoLongRange && bDoForces && (flags & GMX_FORCE_SEPLRF));
753 bUseGPU = fr->nbv->bUseGPU;
754 bUseOrEmulGPU = bUseGPU || (nbv->grp[0].kernel_type == nbnxnk8x8x8_PlainC);
758 update_forcerec(fr, box);
760 if (NEED_MUTOT(*inputrec))
762 /* Calculate total (local) dipole moment in a temporary common array.
763 * This makes it possible to sum them over nodes faster.
765 calc_mu(start, homenr,
766 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
771 if (fr->ePBC != epbcNONE)
773 /* Compute shift vectors every step,
774 * because of pressure coupling or box deformation!
776 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
778 calc_shifts(box, fr->shift_vec);
783 put_atoms_in_box_omp(fr->ePBC, box, homenr, x);
784 inc_nrnb(nrnb, eNR_SHIFTX, homenr);
786 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
788 unshift_self(graph, box, x);
792 nbnxn_atomdata_copy_shiftvec(flags & GMX_FORCE_DYNAMICBOX,
793 fr->shift_vec, nbv->grp[0].nbat);
796 if (!(cr->duty & DUTY_PME))
798 /* Send particle coordinates to the pme nodes.
799 * Since this is only implemented for domain decomposition
800 * and domain decomposition does not use the graph,
801 * we do not need to worry about shifting.
806 wallcycle_start(wcycle, ewcPP_PMESENDX);
808 bBS = (inputrec->nwall == 2);
812 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
815 if (EEL_PME(fr->eeltype))
817 pme_flags |= GMX_PME_DO_COULOMB;
820 if (EVDW_PME(fr->vdwtype))
822 pme_flags |= GMX_PME_DO_LJ;
823 if (fr->ljpme_combination_rule == eljpmeLB)
825 pme_flags |= GMX_PME_LJ_LB;
829 gmx_pme_send_coordinates(cr, bBS ? boxs : box, x,
830 mdatoms->nChargePerturbed, mdatoms->nTypePerturbed, lambda[efptCOUL], lambda[efptVDW],
831 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),
834 wallcycle_stop(wcycle, ewcPP_PMESENDX);
838 /* do gridding for pair search */
841 if (graph && bStateChanged)
843 /* Calculate intramolecular shift vectors to make molecules whole */
844 mk_mshift(fplog, graph, fr->ePBC, box, x);
848 box_diag[XX] = box[XX][XX];
849 box_diag[YY] = box[YY][YY];
850 box_diag[ZZ] = box[ZZ][ZZ];
852 wallcycle_start(wcycle, ewcNS);
855 wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
856 nbnxn_put_on_grid(nbv->nbs, fr->ePBC, box,
858 0, mdatoms->homenr, -1, fr->cginfo, x,
860 nbv->grp[eintLocal].kernel_type,
861 nbv->grp[eintLocal].nbat);
862 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
866 wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
867 nbnxn_put_on_grid_nonlocal(nbv->nbs, domdec_zones(cr->dd),
869 nbv->grp[eintNonlocal].kernel_type,
870 nbv->grp[eintNonlocal].nbat);
871 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
874 if (nbv->ngrp == 1 ||
875 nbv->grp[eintNonlocal].nbat == nbv->grp[eintLocal].nbat)
877 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatAll,
878 nbv->nbs, mdatoms, fr->cginfo);
882 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatLocal,
883 nbv->nbs, mdatoms, fr->cginfo);
884 nbnxn_atomdata_set(nbv->grp[eintNonlocal].nbat, eatAll,
885 nbv->nbs, mdatoms, fr->cginfo);
887 wallcycle_stop(wcycle, ewcNS);
890 /* initialize the GPU atom data and copy shift vector */
895 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
896 nbnxn_cuda_init_atomdata(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
897 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
900 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
901 nbnxn_cuda_upload_shiftvec(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
902 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
905 /* do local pair search */
908 wallcycle_start_nocount(wcycle, ewcNS);
909 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
910 nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintLocal].nbat,
913 nbv->min_ci_balanced,
914 &nbv->grp[eintLocal].nbl_lists,
916 nbv->grp[eintLocal].kernel_type,
918 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
922 /* initialize local pair-list on the GPU */
923 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
924 nbv->grp[eintLocal].nbl_lists.nbl[0],
927 wallcycle_stop(wcycle, ewcNS);
931 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
932 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
933 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, FALSE, x,
934 nbv->grp[eintLocal].nbat);
935 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
936 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
941 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
942 /* launch local nonbonded F on GPU */
943 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFNo,
945 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
948 /* Communicate coordinates and sum dipole if necessary +
949 do non-local pair search */
950 if (DOMAINDECOMP(cr))
952 bDiffKernels = (nbv->grp[eintNonlocal].kernel_type !=
953 nbv->grp[eintLocal].kernel_type);
957 /* With GPU+CPU non-bonded calculations we need to copy
958 * the local coordinates to the non-local nbat struct
959 * (in CPU format) as the non-local kernel call also
960 * calculates the local - non-local interactions.
962 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
963 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
964 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, TRUE, x,
965 nbv->grp[eintNonlocal].nbat);
966 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
967 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
972 wallcycle_start_nocount(wcycle, ewcNS);
973 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
977 nbnxn_grid_add_simple(nbv->nbs, nbv->grp[eintNonlocal].nbat);
980 nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintNonlocal].nbat,
983 nbv->min_ci_balanced,
984 &nbv->grp[eintNonlocal].nbl_lists,
986 nbv->grp[eintNonlocal].kernel_type,
989 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
991 if (nbv->grp[eintNonlocal].kernel_type == nbnxnk8x8x8_CUDA)
993 /* initialize non-local pair-list on the GPU */
994 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
995 nbv->grp[eintNonlocal].nbl_lists.nbl[0],
998 wallcycle_stop(wcycle, ewcNS);
1002 wallcycle_start(wcycle, ewcMOVEX);
1003 dd_move_x(cr->dd, box, x);
1005 /* When we don't need the total dipole we sum it in global_stat */
1006 if (bStateChanged && NEED_MUTOT(*inputrec))
1008 gmx_sumd(2*DIM, mu, cr);
1010 wallcycle_stop(wcycle, ewcMOVEX);
1012 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1013 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
1014 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatNonlocal, FALSE, x,
1015 nbv->grp[eintNonlocal].nbat);
1016 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
1017 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1020 if (bUseGPU && !bDiffKernels)
1022 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
1023 /* launch non-local nonbonded F on GPU */
1024 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFNo,
1026 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1032 /* launch D2H copy-back F */
1033 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1034 if (DOMAINDECOMP(cr) && !bDiffKernels)
1036 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintNonlocal].nbat,
1037 flags, eatNonlocal);
1039 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintLocal].nbat,
1041 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1044 if (bStateChanged && NEED_MUTOT(*inputrec))
1048 gmx_sumd(2*DIM, mu, cr);
1051 for (i = 0; i < 2; i++)
1053 for (j = 0; j < DIM; j++)
1055 fr->mu_tot[i][j] = mu[i*DIM + j];
1059 if (fr->efep == efepNO)
1061 copy_rvec(fr->mu_tot[0], mu_tot);
1065 for (j = 0; j < DIM; j++)
1068 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] +
1069 lambda[efptCOUL]*fr->mu_tot[1][j];
1073 /* Reset energies */
1074 reset_enerdata(fr, bNS, enerd, MASTER(cr));
1075 clear_rvecs(SHIFTS, fr->fshift);
1077 if (DOMAINDECOMP(cr) && !(cr->duty & DUTY_PME))
1079 wallcycle_start(wcycle, ewcPPDURINGPME);
1080 dd_force_flop_start(cr->dd, nrnb);
1085 /* Enforced rotation has its own cycle counter that starts after the collective
1086 * coordinates have been communicated. It is added to ddCyclF to allow
1087 * for proper load-balancing */
1088 wallcycle_start(wcycle, ewcROT);
1089 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1090 wallcycle_stop(wcycle, ewcROT);
1093 /* Start the force cycle counter.
1094 * This counter is stopped in do_forcelow_level.
1095 * No parallel communication should occur while this counter is running,
1096 * since that will interfere with the dynamic load balancing.
1098 wallcycle_start(wcycle, ewcFORCE);
1101 /* Reset forces for which the virial is calculated separately:
1102 * PME/Ewald forces if necessary */
1103 if (fr->bF_NoVirSum)
1105 if (flags & GMX_FORCE_VIRIAL)
1107 fr->f_novirsum = fr->f_novirsum_alloc;
1110 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1114 clear_rvecs(homenr, fr->f_novirsum+start);
1119 /* We are not calculating the pressure so we do not need
1120 * a separate array for forces that do not contribute
1127 /* Clear the short- and long-range forces */
1128 clear_rvecs(fr->natoms_force_constr, f);
1129 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1131 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1134 clear_rvec(fr->vir_diag_posres);
1137 if (inputrec->ePull == epullCONSTRAINT)
1139 clear_pull_forces(inputrec->pull);
1142 /* We calculate the non-bonded forces, when done on the CPU, here.
1143 * We do this before calling do_force_lowlevel, as in there bondeds
1144 * forces are calculated before PME, which does communication.
1145 * With this order, non-bonded and bonded force calculation imbalance
1146 * can be balanced out by the domain decomposition load balancing.
1151 /* Maybe we should move this into do_force_lowlevel */
1152 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFYes,
1156 if (!bUseOrEmulGPU || bDiffKernels)
1160 if (DOMAINDECOMP(cr))
1162 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal,
1163 bDiffKernels ? enbvClearFYes : enbvClearFNo,
1173 aloc = eintNonlocal;
1176 /* Add all the non-bonded force to the normal force array.
1177 * This can be split into a local a non-local part when overlapping
1178 * communication with calculation with domain decomposition.
1180 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1181 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1182 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1183 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatAll, nbv->grp[aloc].nbat, f);
1184 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1185 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1186 wallcycle_start_nocount(wcycle, ewcFORCE);
1188 /* if there are multiple fshift output buffers reduce them */
1189 if ((flags & GMX_FORCE_VIRIAL) &&
1190 nbv->grp[aloc].nbl_lists.nnbl > 1)
1192 nbnxn_atomdata_add_nbat_fshift_to_fshift(nbv->grp[aloc].nbat,
1197 /* update QMMMrec, if necessary */
1200 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1203 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1205 posres_wrapper(fplog, flags, bSepDVDL, inputrec, nrnb, top, box, x,
1209 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_FBPOSRES].nr > 0)
1211 fbposres_wrapper(inputrec, nrnb, top, box, x, enerd, fr);
1214 /* Compute the bonded and non-bonded energies and optionally forces */
1215 do_force_lowlevel(fplog, step, fr, inputrec, &(top->idef),
1216 cr, nrnb, wcycle, mdatoms,
1217 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, top, fr->born,
1218 &(top->atomtypes), bBornRadii, box,
1219 inputrec->fepvals, lambda, graph, &(top->excls), fr->mu_tot,
1220 flags, &cycles_pme);
1224 if (do_per_step(step, inputrec->nstcalclr))
1226 /* Add the long range forces to the short range forces */
1227 for (i = 0; i < fr->natoms_force_constr; i++)
1229 rvec_add(fr->f_twin[i], f[i], f[i]);
1234 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1238 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1241 if (bUseOrEmulGPU && !bDiffKernels)
1243 /* wait for non-local forces (or calculate in emulation mode) */
1244 if (DOMAINDECOMP(cr))
1250 wallcycle_start(wcycle, ewcWAIT_GPU_NB_NL);
1251 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1252 nbv->grp[eintNonlocal].nbat,
1254 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1256 cycles_tmp = wallcycle_stop(wcycle, ewcWAIT_GPU_NB_NL);
1257 cycles_wait_gpu += cycles_tmp;
1258 cycles_force += cycles_tmp;
1262 wallcycle_start_nocount(wcycle, ewcFORCE);
1263 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFYes,
1265 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1267 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1268 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1269 /* skip the reduction if there was no non-local work to do */
1270 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1272 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatNonlocal,
1273 nbv->grp[eintNonlocal].nbat, f);
1275 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1276 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1280 if (bDoForces && DOMAINDECOMP(cr))
1282 /* Communicate the forces */
1283 wallcycle_start(wcycle, ewcMOVEF);
1284 dd_move_f(cr->dd, f, fr->fshift);
1285 /* Do we need to communicate the separate force array
1286 * for terms that do not contribute to the single sum virial?
1287 * Position restraints and electric fields do not introduce
1288 * inter-cg forces, only full electrostatics methods do.
1289 * When we do not calculate the virial, fr->f_novirsum = f,
1290 * so we have already communicated these forces.
1292 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1293 (flags & GMX_FORCE_VIRIAL))
1295 dd_move_f(cr->dd, fr->f_novirsum, NULL);
1299 /* We should not update the shift forces here,
1300 * since f_twin is already included in f.
1302 dd_move_f(cr->dd, fr->f_twin, NULL);
1304 wallcycle_stop(wcycle, ewcMOVEF);
1309 /* wait for local forces (or calculate in emulation mode) */
1312 wallcycle_start(wcycle, ewcWAIT_GPU_NB_L);
1313 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1314 nbv->grp[eintLocal].nbat,
1316 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1318 cycles_wait_gpu += wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
1320 /* now clear the GPU outputs while we finish the step on the CPU */
1322 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1323 nbnxn_cuda_clear_outputs(nbv->cu_nbv, flags);
1324 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1328 wallcycle_start_nocount(wcycle, ewcFORCE);
1329 do_nb_verlet(fr, ic, enerd, flags, eintLocal,
1330 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
1332 wallcycle_stop(wcycle, ewcFORCE);
1334 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1335 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1336 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1338 /* skip the reduction if there was no non-local work to do */
1339 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatLocal,
1340 nbv->grp[eintLocal].nbat, f);
1342 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1343 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1346 if (DOMAINDECOMP(cr))
1348 dd_force_flop_stop(cr->dd, nrnb);
1351 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1354 dd_cycles_add(cr->dd, cycles_wait_gpu, ddCyclWaitGPU);
1361 if (IR_ELEC_FIELD(*inputrec))
1363 /* Compute forces due to electric field */
1364 calc_f_el(MASTER(cr) ? field : NULL,
1365 start, homenr, mdatoms->chargeA, fr->f_novirsum,
1366 inputrec->ex, inputrec->et, t);
1369 /* If we have NoVirSum forces, but we do not calculate the virial,
1370 * we sum fr->f_novirum=f later.
1372 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1374 wallcycle_start(wcycle, ewcVSITESPREAD);
1375 spread_vsite_f(vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1376 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1377 wallcycle_stop(wcycle, ewcVSITESPREAD);
1381 wallcycle_start(wcycle, ewcVSITESPREAD);
1382 spread_vsite_f(vsite, x, fr->f_twin, NULL, FALSE, NULL,
1384 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1385 wallcycle_stop(wcycle, ewcVSITESPREAD);
1389 if (flags & GMX_FORCE_VIRIAL)
1391 /* Calculation of the virial must be done after vsites! */
1392 calc_virial(0, mdatoms->homenr, x, f,
1393 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1397 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1399 pull_potential_wrapper(fplog, bSepDVDL, cr, inputrec, box, x,
1400 f, vir_force, mdatoms, enerd, lambda, t);
1403 /* Add the forces from enforced rotation potentials (if any) */
1406 wallcycle_start(wcycle, ewcROTadd);
1407 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
1408 wallcycle_stop(wcycle, ewcROTadd);
1411 if (PAR(cr) && !(cr->duty & DUTY_PME))
1413 /* In case of node-splitting, the PP nodes receive the long-range
1414 * forces, virial and energy from the PME nodes here.
1416 pme_receive_force_ener(fplog, bSepDVDL, cr, wcycle, enerd, fr);
1421 post_process_forces(cr, step, nrnb, wcycle,
1422 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
1426 /* Sum the potential energy terms from group contributions */
1427 sum_epot(&(enerd->grpp), enerd->term);
1430 void do_force_cutsGROUP(FILE *fplog, t_commrec *cr,
1431 t_inputrec *inputrec,
1432 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1433 gmx_localtop_t *top,
1434 gmx_groups_t *groups,
1435 matrix box, rvec x[], history_t *hist,
1439 gmx_enerdata_t *enerd, t_fcdata *fcd,
1440 real *lambda, t_graph *graph,
1441 t_forcerec *fr, gmx_vsite_t *vsite, rvec mu_tot,
1442 double t, FILE *field, gmx_edsam_t ed,
1443 gmx_bool bBornRadii,
1449 gmx_bool bSepDVDL, bStateChanged, bNS, bFillGrid, bCalcCGCM, bBS;
1450 gmx_bool bDoLongRangeNS, bDoForces, bDoPotential, bSepLRF;
1451 gmx_bool bDoAdressWF;
1453 rvec vzero, box_diag;
1454 real e, v, dvdlambda[efptNR];
1456 float cycles_pme, cycles_force;
1459 homenr = mdatoms->homenr;
1461 bSepDVDL = (fr->bSepDVDL && do_per_step(step, inputrec->nstlog));
1463 clear_mat(vir_force);
1466 if (DOMAINDECOMP(cr))
1468 cg1 = cr->dd->ncg_tot;
1479 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
1480 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
1481 /* Should we update the long-range neighborlists at this step? */
1482 bDoLongRangeNS = fr->bTwinRange && bNS;
1483 /* Should we perform the long-range nonbonded evaluation inside the neighborsearching? */
1484 bFillGrid = (bNS && bStateChanged);
1485 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
1486 bDoForces = (flags & GMX_FORCE_FORCES);
1487 bDoPotential = (flags & GMX_FORCE_ENERGY);
1488 bSepLRF = ((inputrec->nstcalclr > 1) && bDoForces &&
1489 (flags & GMX_FORCE_SEPLRF) && (flags & GMX_FORCE_DO_LR));
1491 /* should probably move this to the forcerec since it doesn't change */
1492 bDoAdressWF = ((fr->adress_type != eAdressOff));
1496 update_forcerec(fr, box);
1498 if (NEED_MUTOT(*inputrec))
1500 /* Calculate total (local) dipole moment in a temporary common array.
1501 * This makes it possible to sum them over nodes faster.
1503 calc_mu(start, homenr,
1504 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
1509 if (fr->ePBC != epbcNONE)
1511 /* Compute shift vectors every step,
1512 * because of pressure coupling or box deformation!
1514 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
1516 calc_shifts(box, fr->shift_vec);
1521 put_charge_groups_in_box(fplog, cg0, cg1, fr->ePBC, box,
1522 &(top->cgs), x, fr->cg_cm);
1523 inc_nrnb(nrnb, eNR_CGCM, homenr);
1524 inc_nrnb(nrnb, eNR_RESETX, cg1-cg0);
1526 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
1528 unshift_self(graph, box, x);
1533 calc_cgcm(fplog, cg0, cg1, &(top->cgs), x, fr->cg_cm);
1534 inc_nrnb(nrnb, eNR_CGCM, homenr);
1537 if (bCalcCGCM && gmx_debug_at)
1539 pr_rvecs(debug, 0, "cgcm", fr->cg_cm, top->cgs.nr);
1543 if (!(cr->duty & DUTY_PME))
1545 /* Send particle coordinates to the pme nodes.
1546 * Since this is only implemented for domain decomposition
1547 * and domain decomposition does not use the graph,
1548 * we do not need to worry about shifting.
1553 wallcycle_start(wcycle, ewcPP_PMESENDX);
1555 bBS = (inputrec->nwall == 2);
1558 copy_mat(box, boxs);
1559 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
1562 if (EEL_PME(fr->eeltype))
1564 pme_flags |= GMX_PME_DO_COULOMB;
1567 if (EVDW_PME(fr->vdwtype))
1569 pme_flags |= GMX_PME_DO_LJ;
1570 if (fr->ljpme_combination_rule == eljpmeLB)
1572 pme_flags |= GMX_PME_LJ_LB;
1576 gmx_pme_send_coordinates(cr, bBS ? boxs : box, x,
1577 mdatoms->nChargePerturbed, mdatoms->nTypePerturbed, lambda[efptCOUL], lambda[efptVDW],
1578 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),
1581 wallcycle_stop(wcycle, ewcPP_PMESENDX);
1583 #endif /* GMX_MPI */
1585 /* Communicate coordinates and sum dipole if necessary */
1586 if (DOMAINDECOMP(cr))
1588 wallcycle_start(wcycle, ewcMOVEX);
1589 dd_move_x(cr->dd, box, x);
1590 wallcycle_stop(wcycle, ewcMOVEX);
1593 /* update adress weight beforehand */
1594 if (bStateChanged && bDoAdressWF)
1596 /* need pbc for adress weight calculation with pbc_dx */
1597 set_pbc(&pbc, inputrec->ePBC, box);
1598 if (fr->adress_site == eAdressSITEcog)
1600 update_adress_weights_cog(top->idef.iparams, top->idef.il, x, fr, mdatoms,
1601 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1603 else if (fr->adress_site == eAdressSITEcom)
1605 update_adress_weights_com(fplog, cg0, cg1, &(top->cgs), x, fr, mdatoms,
1606 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1608 else if (fr->adress_site == eAdressSITEatomatom)
1610 update_adress_weights_atom_per_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1611 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1615 update_adress_weights_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1616 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1620 if (NEED_MUTOT(*inputrec))
1627 gmx_sumd(2*DIM, mu, cr);
1629 for (i = 0; i < 2; i++)
1631 for (j = 0; j < DIM; j++)
1633 fr->mu_tot[i][j] = mu[i*DIM + j];
1637 if (fr->efep == efepNO)
1639 copy_rvec(fr->mu_tot[0], mu_tot);
1643 for (j = 0; j < DIM; j++)
1646 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] + lambda[efptCOUL]*fr->mu_tot[1][j];
1651 /* Reset energies */
1652 reset_enerdata(fr, bNS, enerd, MASTER(cr));
1653 clear_rvecs(SHIFTS, fr->fshift);
1657 wallcycle_start(wcycle, ewcNS);
1659 if (graph && bStateChanged)
1661 /* Calculate intramolecular shift vectors to make molecules whole */
1662 mk_mshift(fplog, graph, fr->ePBC, box, x);
1665 /* Do the actual neighbour searching */
1667 groups, top, mdatoms,
1668 cr, nrnb, bFillGrid,
1671 wallcycle_stop(wcycle, ewcNS);
1674 if (inputrec->implicit_solvent && bNS)
1676 make_gb_nblist(cr, inputrec->gb_algorithm,
1677 x, box, fr, &top->idef, graph, fr->born);
1680 if (DOMAINDECOMP(cr) && !(cr->duty & DUTY_PME))
1682 wallcycle_start(wcycle, ewcPPDURINGPME);
1683 dd_force_flop_start(cr->dd, nrnb);
1688 /* Enforced rotation has its own cycle counter that starts after the collective
1689 * coordinates have been communicated. It is added to ddCyclF to allow
1690 * for proper load-balancing */
1691 wallcycle_start(wcycle, ewcROT);
1692 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1693 wallcycle_stop(wcycle, ewcROT);
1696 /* Start the force cycle counter.
1697 * This counter is stopped in do_forcelow_level.
1698 * No parallel communication should occur while this counter is running,
1699 * since that will interfere with the dynamic load balancing.
1701 wallcycle_start(wcycle, ewcFORCE);
1705 /* Reset forces for which the virial is calculated separately:
1706 * PME/Ewald forces if necessary */
1707 if (fr->bF_NoVirSum)
1709 if (flags & GMX_FORCE_VIRIAL)
1711 fr->f_novirsum = fr->f_novirsum_alloc;
1714 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1718 clear_rvecs(homenr, fr->f_novirsum+start);
1723 /* We are not calculating the pressure so we do not need
1724 * a separate array for forces that do not contribute
1731 /* Clear the short- and long-range forces */
1732 clear_rvecs(fr->natoms_force_constr, f);
1733 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1735 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1738 clear_rvec(fr->vir_diag_posres);
1740 if (inputrec->ePull == epullCONSTRAINT)
1742 clear_pull_forces(inputrec->pull);
1745 /* update QMMMrec, if necessary */
1748 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1751 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1753 posres_wrapper(fplog, flags, bSepDVDL, inputrec, nrnb, top, box, x,
1757 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_FBPOSRES].nr > 0)
1759 fbposres_wrapper(inputrec, nrnb, top, box, x, enerd, fr);
1762 /* Compute the bonded and non-bonded energies and optionally forces */
1763 do_force_lowlevel(fplog, step, fr, inputrec, &(top->idef),
1764 cr, nrnb, wcycle, mdatoms,
1765 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, top, fr->born,
1766 &(top->atomtypes), bBornRadii, box,
1767 inputrec->fepvals, lambda,
1768 graph, &(top->excls), fr->mu_tot,
1774 if (do_per_step(step, inputrec->nstcalclr))
1776 /* Add the long range forces to the short range forces */
1777 for (i = 0; i < fr->natoms_force_constr; i++)
1779 rvec_add(fr->f_twin[i], f[i], f[i]);
1784 cycles_force = wallcycle_stop(wcycle, ewcFORCE);
1788 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1791 if (DOMAINDECOMP(cr))
1793 dd_force_flop_stop(cr->dd, nrnb);
1796 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1802 if (IR_ELEC_FIELD(*inputrec))
1804 /* Compute forces due to electric field */
1805 calc_f_el(MASTER(cr) ? field : NULL,
1806 start, homenr, mdatoms->chargeA, fr->f_novirsum,
1807 inputrec->ex, inputrec->et, t);
1810 if (bDoAdressWF && fr->adress_icor == eAdressICThermoForce)
1812 /* Compute thermodynamic force in hybrid AdResS region */
1813 adress_thermo_force(start, homenr, &(top->cgs), x, fr->f_novirsum, fr, mdatoms,
1814 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1817 /* Communicate the forces */
1818 if (DOMAINDECOMP(cr))
1820 wallcycle_start(wcycle, ewcMOVEF);
1821 dd_move_f(cr->dd, f, fr->fshift);
1822 /* Do we need to communicate the separate force array
1823 * for terms that do not contribute to the single sum virial?
1824 * Position restraints and electric fields do not introduce
1825 * inter-cg forces, only full electrostatics methods do.
1826 * When we do not calculate the virial, fr->f_novirsum = f,
1827 * so we have already communicated these forces.
1829 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1830 (flags & GMX_FORCE_VIRIAL))
1832 dd_move_f(cr->dd, fr->f_novirsum, NULL);
1836 /* We should not update the shift forces here,
1837 * since f_twin is already included in f.
1839 dd_move_f(cr->dd, fr->f_twin, NULL);
1841 wallcycle_stop(wcycle, ewcMOVEF);
1844 /* If we have NoVirSum forces, but we do not calculate the virial,
1845 * we sum fr->f_novirum=f later.
1847 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1849 wallcycle_start(wcycle, ewcVSITESPREAD);
1850 spread_vsite_f(vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1851 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1852 wallcycle_stop(wcycle, ewcVSITESPREAD);
1856 wallcycle_start(wcycle, ewcVSITESPREAD);
1857 spread_vsite_f(vsite, x, fr->f_twin, NULL, FALSE, NULL,
1859 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1860 wallcycle_stop(wcycle, ewcVSITESPREAD);
1864 if (flags & GMX_FORCE_VIRIAL)
1866 /* Calculation of the virial must be done after vsites! */
1867 calc_virial(0, mdatoms->homenr, x, f,
1868 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1872 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1874 pull_potential_wrapper(fplog, bSepDVDL, cr, inputrec, box, x,
1875 f, vir_force, mdatoms, enerd, lambda, t);
1878 /* Add the forces from enforced rotation potentials (if any) */
1881 wallcycle_start(wcycle, ewcROTadd);
1882 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
1883 wallcycle_stop(wcycle, ewcROTadd);
1886 if (PAR(cr) && !(cr->duty & DUTY_PME))
1888 /* In case of node-splitting, the PP nodes receive the long-range
1889 * forces, virial and energy from the PME nodes here.
1891 pme_receive_force_ener(fplog, bSepDVDL, cr, wcycle, enerd, fr);
1896 post_process_forces(cr, step, nrnb, wcycle,
1897 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
1901 /* Sum the potential energy terms from group contributions */
1902 sum_epot(&(enerd->grpp), enerd->term);
1905 void do_force(FILE *fplog, t_commrec *cr,
1906 t_inputrec *inputrec,
1907 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1908 gmx_localtop_t *top,
1909 gmx_groups_t *groups,
1910 matrix box, rvec x[], history_t *hist,
1914 gmx_enerdata_t *enerd, t_fcdata *fcd,
1915 real *lambda, t_graph *graph,
1917 gmx_vsite_t *vsite, rvec mu_tot,
1918 double t, FILE *field, gmx_edsam_t ed,
1919 gmx_bool bBornRadii,
1922 /* modify force flag if not doing nonbonded */
1923 if (!fr->bNonbonded)
1925 flags &= ~GMX_FORCE_NONBONDED;
1928 switch (inputrec->cutoff_scheme)
1931 do_force_cutsVERLET(fplog, cr, inputrec,
1947 do_force_cutsGROUP(fplog, cr, inputrec,
1962 gmx_incons("Invalid cut-off scheme passed!");
1967 void do_constrain_first(FILE *fplog, gmx_constr_t constr,
1968 t_inputrec *ir, t_mdatoms *md,
1969 t_state *state, t_commrec *cr, t_nrnb *nrnb,
1970 t_forcerec *fr, gmx_localtop_t *top)
1972 int i, m, start, end;
1974 real dt = ir->delta_t;
1978 snew(savex, state->natoms);
1985 fprintf(debug, "vcm: start=%d, homenr=%d, end=%d\n",
1986 start, md->homenr, end);
1988 /* Do a first constrain to reset particles... */
1989 step = ir->init_step;
1992 char buf[STEPSTRSIZE];
1993 fprintf(fplog, "\nConstraining the starting coordinates (step %s)\n",
1994 gmx_step_str(step, buf));
1998 /* constrain the current position */
1999 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2000 ir, NULL, cr, step, 0, md,
2001 state->x, state->x, NULL,
2002 fr->bMolPBC, state->box,
2003 state->lambda[efptBONDED], &dvdl_dum,
2004 NULL, NULL, nrnb, econqCoord,
2005 ir->epc == epcMTTK, state->veta, state->veta);
2008 /* constrain the inital velocity, and save it */
2009 /* also may be useful if we need the ekin from the halfstep for velocity verlet */
2010 /* might not yet treat veta correctly */
2011 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2012 ir, NULL, cr, step, 0, md,
2013 state->x, state->v, state->v,
2014 fr->bMolPBC, state->box,
2015 state->lambda[efptBONDED], &dvdl_dum,
2016 NULL, NULL, nrnb, econqVeloc,
2017 ir->epc == epcMTTK, state->veta, state->veta);
2019 /* constrain the inital velocities at t-dt/2 */
2020 if (EI_STATE_VELOCITY(ir->eI) && ir->eI != eiVV)
2022 for (i = start; (i < end); i++)
2024 for (m = 0; (m < DIM); m++)
2026 /* Reverse the velocity */
2027 state->v[i][m] = -state->v[i][m];
2028 /* Store the position at t-dt in buf */
2029 savex[i][m] = state->x[i][m] + dt*state->v[i][m];
2032 /* Shake the positions at t=-dt with the positions at t=0
2033 * as reference coordinates.
2037 char buf[STEPSTRSIZE];
2038 fprintf(fplog, "\nConstraining the coordinates at t0-dt (step %s)\n",
2039 gmx_step_str(step, buf));
2042 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2043 ir, NULL, cr, step, -1, md,
2044 state->x, savex, NULL,
2045 fr->bMolPBC, state->box,
2046 state->lambda[efptBONDED], &dvdl_dum,
2047 state->v, NULL, nrnb, econqCoord,
2048 ir->epc == epcMTTK, state->veta, state->veta);
2050 for (i = start; i < end; i++)
2052 for (m = 0; m < DIM; m++)
2054 /* Re-reverse the velocities */
2055 state->v[i][m] = -state->v[i][m];
2064 integrate_table(real vdwtab[], real scale, int offstart, int rstart, int rend,
2065 double *enerout, double *virout)
2067 double enersum, virsum;
2068 double invscale, invscale2, invscale3;
2069 double r, ea, eb, ec, pa, pb, pc, pd;
2071 int ri, offset, tabfactor;
2073 invscale = 1.0/scale;
2074 invscale2 = invscale*invscale;
2075 invscale3 = invscale*invscale2;
2077 /* Following summation derived from cubic spline definition,
2078 * Numerical Recipies in C, second edition, p. 113-116. Exact for
2079 * the cubic spline. We first calculate the negative of the
2080 * energy from rvdw to rvdw_switch, assuming that g(r)=1, and then
2081 * add the more standard, abrupt cutoff correction to that result,
2082 * yielding the long-range correction for a switched function. We
2083 * perform both the pressure and energy loops at the same time for
2084 * simplicity, as the computational cost is low. */
2088 /* Since the dispersion table has been scaled down a factor
2089 * 6.0 and the repulsion a factor 12.0 to compensate for the
2090 * c6/c12 parameters inside nbfp[] being scaled up (to save
2091 * flops in kernels), we need to correct for this.
2102 for (ri = rstart; ri < rend; ++ri)
2106 eb = 2.0*invscale2*r;
2110 pb = 3.0*invscale2*r;
2111 pc = 3.0*invscale*r*r;
2114 /* this "8" is from the packing in the vdwtab array - perhaps
2115 should be defined? */
2117 offset = 8*ri + offstart;
2118 y0 = vdwtab[offset];
2119 f = vdwtab[offset+1];
2120 g = vdwtab[offset+2];
2121 h = vdwtab[offset+3];
2123 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);
2124 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);
2126 *enerout = 4.0*M_PI*enersum*tabfactor;
2127 *virout = 4.0*M_PI*virsum*tabfactor;
2130 void calc_enervirdiff(FILE *fplog, int eDispCorr, t_forcerec *fr)
2132 double eners[2], virs[2], enersum, virsum, y0, f, g, h;
2133 double r0, r1, r, rc3, rc9, ea, eb, ec, pa, pb, pc, pd;
2134 double invscale, invscale2, invscale3;
2135 int ri0, ri1, ri, i, offstart, offset;
2136 real scale, *vdwtab, tabfactor, tmp;
2138 fr->enershiftsix = 0;
2139 fr->enershifttwelve = 0;
2140 fr->enerdiffsix = 0;
2141 fr->enerdifftwelve = 0;
2143 fr->virdifftwelve = 0;
2145 if (eDispCorr != edispcNO)
2147 for (i = 0; i < 2; i++)
2152 if (fr->vdwtype == evdwSWITCH || fr->vdwtype == evdwSHIFT ||
2153 fr->vdw_modifier == eintmodPOTSWITCH ||
2154 fr->vdw_modifier == eintmodFORCESWITCH)
2156 if (fr->rvdw_switch == 0)
2159 "With dispersion correction rvdw-switch can not be zero "
2160 "for vdw-type = %s", evdw_names[fr->vdwtype]);
2163 scale = fr->nblists[0].table_elec_vdw.scale;
2164 vdwtab = fr->nblists[0].table_vdw.data;
2166 /* Round the cut-offs to exact table values for precision */
2167 ri0 = floor(fr->rvdw_switch*scale);
2168 ri1 = ceil(fr->rvdw*scale);
2174 if (fr->vdwtype == evdwSHIFT ||
2175 fr->vdw_modifier == eintmodFORCESWITCH)
2177 /* Determine the constant energy shift below rvdw_switch.
2178 * Table has a scale factor since we have scaled it down to compensate
2179 * for scaling-up c6/c12 with the derivative factors to save flops in analytical kernels.
2181 fr->enershiftsix = (real)(-1.0/(rc3*rc3)) - 6.0*vdwtab[8*ri0];
2182 fr->enershifttwelve = (real)( 1.0/(rc9*rc3)) - 12.0*vdwtab[8*ri0 + 4];
2184 /* Add the constant part from 0 to rvdw_switch.
2185 * This integration from 0 to rvdw_switch overcounts the number
2186 * of interactions by 1, as it also counts the self interaction.
2187 * We will correct for this later.
2189 eners[0] += 4.0*M_PI*fr->enershiftsix*rc3/3.0;
2190 eners[1] += 4.0*M_PI*fr->enershifttwelve*rc3/3.0;
2191 for (i = 0; i < 2; i++)
2195 integrate_table(vdwtab, scale, (i == 0 ? 0 : 4), ri0, ri1, &enersum, &virsum);
2196 eners[i] -= enersum;
2200 /* now add the correction for rvdw_switch to infinity */
2201 eners[0] += -4.0*M_PI/(3.0*rc3);
2202 eners[1] += 4.0*M_PI/(9.0*rc9);
2203 virs[0] += 8.0*M_PI/rc3;
2204 virs[1] += -16.0*M_PI/(3.0*rc9);
2206 else if (fr->vdwtype == evdwCUT ||
2207 EVDW_PME(fr->vdwtype) ||
2208 fr->vdwtype == evdwUSER)
2210 if (fr->vdwtype == evdwUSER && fplog)
2213 "WARNING: using dispersion correction with user tables\n");
2216 /* Note that with LJ-PME, the dispersion correction is multiplied
2217 * by the difference between the actual C6 and the value of C6
2218 * that would produce the combination rule.
2219 * This means the normal energy and virial difference formulas
2223 rc3 = fr->rvdw*fr->rvdw*fr->rvdw;
2225 /* Contribution beyond the cut-off */
2226 eners[0] += -4.0*M_PI/(3.0*rc3);
2227 eners[1] += 4.0*M_PI/(9.0*rc9);
2228 if (fr->vdw_modifier == eintmodPOTSHIFT)
2230 /* Contribution within the cut-off */
2231 eners[0] += -4.0*M_PI/(3.0*rc3);
2232 eners[1] += 4.0*M_PI/(3.0*rc9);
2234 /* Contribution beyond the cut-off */
2235 virs[0] += 8.0*M_PI/rc3;
2236 virs[1] += -16.0*M_PI/(3.0*rc9);
2241 "Dispersion correction is not implemented for vdw-type = %s",
2242 evdw_names[fr->vdwtype]);
2245 /* TODO: remove this code once we have group LJ-PME kernels
2246 * that calculate the exact, full LJ param C6/r^6 within the cut-off,
2247 * as the current nbnxn kernels do.
2249 if (fr->vdwtype == evdwPME && fr->cutoff_scheme == ecutsGROUP)
2251 /* Calculate self-interaction coefficient (assuming that
2252 * the reciprocal-space contribution is constant in the
2253 * region that contributes to the self-interaction).
2255 fr->enershiftsix = pow(fr->ewaldcoeff_lj, 6) / 6.0;
2257 eners[0] += -pow(sqrt(M_PI)*fr->ewaldcoeff_lj, 3)/3.0;
2258 virs[0] += pow(sqrt(M_PI)*fr->ewaldcoeff_lj, 3);
2261 fr->enerdiffsix = eners[0];
2262 fr->enerdifftwelve = eners[1];
2263 /* The 0.5 is due to the Gromacs definition of the virial */
2264 fr->virdiffsix = 0.5*virs[0];
2265 fr->virdifftwelve = 0.5*virs[1];
2269 void calc_dispcorr(FILE *fplog, t_inputrec *ir, t_forcerec *fr,
2270 gmx_int64_t step, int natoms,
2271 matrix box, real lambda, tensor pres, tensor virial,
2272 real *prescorr, real *enercorr, real *dvdlcorr)
2274 gmx_bool bCorrAll, bCorrPres;
2275 real dvdlambda, invvol, dens, ninter, avcsix, avctwelve, enerdiff, svir = 0, spres = 0;
2285 if (ir->eDispCorr != edispcNO)
2287 bCorrAll = (ir->eDispCorr == edispcAllEner ||
2288 ir->eDispCorr == edispcAllEnerPres);
2289 bCorrPres = (ir->eDispCorr == edispcEnerPres ||
2290 ir->eDispCorr == edispcAllEnerPres);
2292 invvol = 1/det(box);
2295 /* Only correct for the interactions with the inserted molecule */
2296 dens = (natoms - fr->n_tpi)*invvol;
2301 dens = natoms*invvol;
2302 ninter = 0.5*natoms;
2305 if (ir->efep == efepNO)
2307 avcsix = fr->avcsix[0];
2308 avctwelve = fr->avctwelve[0];
2312 avcsix = (1 - lambda)*fr->avcsix[0] + lambda*fr->avcsix[1];
2313 avctwelve = (1 - lambda)*fr->avctwelve[0] + lambda*fr->avctwelve[1];
2316 enerdiff = ninter*(dens*fr->enerdiffsix - fr->enershiftsix);
2317 *enercorr += avcsix*enerdiff;
2319 if (ir->efep != efepNO)
2321 dvdlambda += (fr->avcsix[1] - fr->avcsix[0])*enerdiff;
2325 enerdiff = ninter*(dens*fr->enerdifftwelve - fr->enershifttwelve);
2326 *enercorr += avctwelve*enerdiff;
2327 if (fr->efep != efepNO)
2329 dvdlambda += (fr->avctwelve[1] - fr->avctwelve[0])*enerdiff;
2335 svir = ninter*dens*avcsix*fr->virdiffsix/3.0;
2336 if (ir->eDispCorr == edispcAllEnerPres)
2338 svir += ninter*dens*avctwelve*fr->virdifftwelve/3.0;
2340 /* The factor 2 is because of the Gromacs virial definition */
2341 spres = -2.0*invvol*svir*PRESFAC;
2343 for (m = 0; m < DIM; m++)
2345 virial[m][m] += svir;
2346 pres[m][m] += spres;
2351 /* Can't currently control when it prints, for now, just print when degugging */
2356 fprintf(debug, "Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
2362 "Long Range LJ corr.: Epot %10g, Pres: %10g, Vir: %10g\n",
2363 *enercorr, spres, svir);
2367 fprintf(debug, "Long Range LJ corr.: Epot %10g\n", *enercorr);
2371 if (fr->bSepDVDL && do_per_step(step, ir->nstlog))
2373 gmx_print_sepdvdl(fplog, "Dispersion correction", *enercorr, dvdlambda);
2375 if (fr->efep != efepNO)
2377 *dvdlcorr += dvdlambda;
2382 void do_pbc_first(FILE *fplog, matrix box, t_forcerec *fr,
2383 t_graph *graph, rvec x[])
2387 fprintf(fplog, "Removing pbc first time\n");
2389 calc_shifts(box, fr->shift_vec);
2392 mk_mshift(fplog, graph, fr->ePBC, box, x);
2395 p_graph(debug, "do_pbc_first 1", graph);
2397 shift_self(graph, box, x);
2398 /* By doing an extra mk_mshift the molecules that are broken
2399 * because they were e.g. imported from another software
2400 * will be made whole again. Such are the healing powers
2403 mk_mshift(fplog, graph, fr->ePBC, box, x);
2406 p_graph(debug, "do_pbc_first 2", graph);
2411 fprintf(fplog, "Done rmpbc\n");
2415 static void low_do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2416 gmx_mtop_t *mtop, rvec x[],
2421 gmx_molblock_t *molb;
2423 if (bFirst && fplog)
2425 fprintf(fplog, "Removing pbc first time\n");
2430 for (mb = 0; mb < mtop->nmolblock; mb++)
2432 molb = &mtop->molblock[mb];
2433 if (molb->natoms_mol == 1 ||
2434 (!bFirst && mtop->moltype[molb->type].cgs.nr == 1))
2436 /* Just one atom or charge group in the molecule, no PBC required */
2437 as += molb->nmol*molb->natoms_mol;
2441 /* Pass NULL iso fplog to avoid graph prints for each molecule type */
2442 mk_graph_ilist(NULL, mtop->moltype[molb->type].ilist,
2443 0, molb->natoms_mol, FALSE, FALSE, graph);
2445 for (mol = 0; mol < molb->nmol; mol++)
2447 mk_mshift(fplog, graph, ePBC, box, x+as);
2449 shift_self(graph, box, x+as);
2450 /* The molecule is whole now.
2451 * We don't need the second mk_mshift call as in do_pbc_first,
2452 * since we no longer need this graph.
2455 as += molb->natoms_mol;
2463 void do_pbc_first_mtop(FILE *fplog, int ePBC, matrix box,
2464 gmx_mtop_t *mtop, rvec x[])
2466 low_do_pbc_mtop(fplog, ePBC, box, mtop, x, TRUE);
2469 void do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2470 gmx_mtop_t *mtop, rvec x[])
2472 low_do_pbc_mtop(fplog, ePBC, box, mtop, x, FALSE);
2475 void finish_run(FILE *fplog, t_commrec *cr,
2476 t_inputrec *inputrec,
2477 t_nrnb nrnb[], gmx_wallcycle_t wcycle,
2478 gmx_walltime_accounting_t walltime_accounting,
2479 wallclock_gpu_t *gputimes,
2480 gmx_bool bWriteStat)
2483 t_nrnb *nrnb_tot = NULL;
2486 double elapsed_time,
2487 elapsed_time_over_all_ranks,
2488 elapsed_time_over_all_threads,
2489 elapsed_time_over_all_threads_over_all_ranks;
2490 wallcycle_sum(cr, wcycle);
2496 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM,
2497 cr->mpi_comm_mysim);
2505 elapsed_time = walltime_accounting_get_elapsed_time(walltime_accounting);
2506 elapsed_time_over_all_ranks = elapsed_time;
2507 elapsed_time_over_all_threads = walltime_accounting_get_elapsed_time_over_all_threads(walltime_accounting);
2508 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
2512 /* reduce elapsed_time over all MPI ranks in the current simulation */
2513 MPI_Allreduce(&elapsed_time,
2514 &elapsed_time_over_all_ranks,
2515 1, MPI_DOUBLE, MPI_SUM,
2516 cr->mpi_comm_mysim);
2517 elapsed_time_over_all_ranks /= cr->nnodes;
2518 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
2519 * current simulation. */
2520 MPI_Allreduce(&elapsed_time_over_all_threads,
2521 &elapsed_time_over_all_threads_over_all_ranks,
2522 1, MPI_DOUBLE, MPI_SUM,
2523 cr->mpi_comm_mysim);
2529 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
2536 if ((cr->duty & DUTY_PP) && DOMAINDECOMP(cr))
2538 print_dd_statistics(cr, inputrec, fplog);
2543 wallcycle_print(fplog, cr->nnodes, cr->npmenodes,
2544 elapsed_time_over_all_ranks,
2547 if (EI_DYNAMICS(inputrec->eI))
2549 delta_t = inputrec->delta_t;
2558 print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks,
2559 elapsed_time_over_all_ranks,
2560 walltime_accounting_get_nsteps_done(walltime_accounting),
2561 delta_t, nbfs, mflop);
2565 print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks,
2566 elapsed_time_over_all_ranks,
2567 walltime_accounting_get_nsteps_done(walltime_accounting),
2568 delta_t, nbfs, mflop);
2573 extern void initialize_lambdas(FILE *fplog, t_inputrec *ir, int *fep_state, real *lambda, double *lam0)
2575 /* this function works, but could probably use a logic rewrite to keep all the different
2576 types of efep straight. */
2579 t_lambda *fep = ir->fepvals;
2581 if ((ir->efep == efepNO) && (ir->bSimTemp == FALSE))
2583 for (i = 0; i < efptNR; i++)
2595 *fep_state = fep->init_fep_state; /* this might overwrite the checkpoint
2596 if checkpoint is set -- a kludge is in for now
2598 for (i = 0; i < efptNR; i++)
2600 /* overwrite lambda state with init_lambda for now for backwards compatibility */
2601 if (fep->init_lambda >= 0) /* if it's -1, it was never initializd */
2603 lambda[i] = fep->init_lambda;
2606 lam0[i] = lambda[i];
2611 lambda[i] = fep->all_lambda[i][*fep_state];
2614 lam0[i] = lambda[i];
2620 /* need to rescale control temperatures to match current state */
2621 for (i = 0; i < ir->opts.ngtc; i++)
2623 if (ir->opts.ref_t[i] > 0)
2625 ir->opts.ref_t[i] = ir->simtempvals->temperatures[*fep_state];
2631 /* Send to the log the information on the current lambdas */
2634 fprintf(fplog, "Initial vector of lambda components:[ ");
2635 for (i = 0; i < efptNR; i++)
2637 fprintf(fplog, "%10.4f ", lambda[i]);
2639 fprintf(fplog, "]\n");
2645 void init_md(FILE *fplog,
2646 t_commrec *cr, t_inputrec *ir, const output_env_t oenv,
2647 double *t, double *t0,
2648 real *lambda, int *fep_state, double *lam0,
2649 t_nrnb *nrnb, gmx_mtop_t *mtop,
2651 int nfile, const t_filenm fnm[],
2652 gmx_mdoutf_t *outf, t_mdebin **mdebin,
2653 tensor force_vir, tensor shake_vir, rvec mu_tot,
2654 gmx_bool *bSimAnn, t_vcm **vcm, unsigned long Flags)
2659 /* Initial values */
2660 *t = *t0 = ir->init_t;
2663 for (i = 0; i < ir->opts.ngtc; i++)
2665 /* set bSimAnn if any group is being annealed */
2666 if (ir->opts.annealing[i] != eannNO)
2673 update_annealing_target_temp(&(ir->opts), ir->init_t);
2676 /* Initialize lambda variables */
2677 initialize_lambdas(fplog, ir, fep_state, lambda, lam0);
2681 *upd = init_update(ir);
2687 *vcm = init_vcm(fplog, &mtop->groups, ir);
2690 if (EI_DYNAMICS(ir->eI) && !(Flags & MD_APPENDFILES))
2692 if (ir->etc == etcBERENDSEN)
2694 please_cite(fplog, "Berendsen84a");
2696 if (ir->etc == etcVRESCALE)
2698 please_cite(fplog, "Bussi2007a");
2706 *outf = init_mdoutf(fplog, nfile, fnm, Flags, cr, ir, mtop, oenv);
2708 *mdebin = init_mdebin((Flags & MD_APPENDFILES) ? NULL : mdoutf_get_fp_ene(*outf),
2709 mtop, ir, mdoutf_get_fp_dhdl(*outf));
2714 please_cite(fplog, "Fritsch12");
2715 please_cite(fplog, "Junghans10");
2717 /* Initiate variables */
2718 clear_mat(force_vir);
2719 clear_mat(shake_vir);