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,2015, 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.
39 #include "gromacs/legacyheaders/sim_util.h"
48 #include "gromacs/domdec/domdec.h"
49 #include "gromacs/essentialdynamics/edsam.h"
50 #include "gromacs/ewald/pme.h"
51 #include "gromacs/gmxlib/nonbonded/nb_free_energy.h"
52 #include "gromacs/gmxlib/nonbonded/nb_kernel.h"
53 #include "gromacs/imd/imd.h"
54 #include "gromacs/legacyheaders/calcmu.h"
55 #include "gromacs/legacyheaders/chargegroup.h"
56 #include "gromacs/legacyheaders/constr.h"
57 #include "gromacs/legacyheaders/copyrite.h"
58 #include "gromacs/legacyheaders/disre.h"
59 #include "gromacs/legacyheaders/force.h"
60 #include "gromacs/legacyheaders/genborn.h"
61 #include "gromacs/legacyheaders/gmx_omp_nthreads.h"
62 #include "gromacs/legacyheaders/mdatoms.h"
63 #include "gromacs/legacyheaders/mdrun.h"
64 #include "gromacs/legacyheaders/names.h"
65 #include "gromacs/legacyheaders/network.h"
66 #include "gromacs/legacyheaders/nonbonded.h"
67 #include "gromacs/legacyheaders/nrnb.h"
68 #include "gromacs/legacyheaders/orires.h"
69 #include "gromacs/legacyheaders/qmmm.h"
70 #include "gromacs/legacyheaders/txtdump.h"
71 #include "gromacs/legacyheaders/typedefs.h"
72 #include "gromacs/legacyheaders/update.h"
73 #include "gromacs/legacyheaders/types/commrec.h"
74 #include "gromacs/listed-forces/bonded.h"
75 #include "gromacs/math/units.h"
76 #include "gromacs/math/vec.h"
77 #include "gromacs/mdlib/nb_verlet.h"
78 #include "gromacs/mdlib/nbnxn_atomdata.h"
79 #include "gromacs/mdlib/nbnxn_search.h"
80 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda.h"
81 #include "gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_data_mgmt.h"
82 #include "gromacs/mdlib/nbnxn_kernels/nbnxn_kernel_gpu_ref.h"
83 #include "gromacs/mdlib/nbnxn_kernels/nbnxn_kernel_ref.h"
84 #include "gromacs/mdlib/nbnxn_kernels/simd_2xnn/nbnxn_kernel_simd_2xnn.h"
85 #include "gromacs/mdlib/nbnxn_kernels/simd_4xn/nbnxn_kernel_simd_4xn.h"
86 #include "gromacs/pbcutil/ishift.h"
87 #include "gromacs/pbcutil/mshift.h"
88 #include "gromacs/pbcutil/pbc.h"
89 #include "gromacs/pulling/pull.h"
90 #include "gromacs/pulling/pull_rotation.h"
91 #include "gromacs/timing/wallcycle.h"
92 #include "gromacs/timing/walltime_accounting.h"
93 #include "gromacs/utility/cstringutil.h"
94 #include "gromacs/utility/gmxmpi.h"
95 #include "gromacs/utility/smalloc.h"
96 #include "gromacs/utility/sysinfo.h"
100 void print_time(FILE *out,
101 gmx_walltime_accounting_t walltime_accounting,
104 t_commrec gmx_unused *cr)
107 char timebuf[STRLEN];
108 double dt, elapsed_seconds, time_per_step;
111 #ifndef GMX_THREAD_MPI
117 fprintf(out, "step %s", gmx_step_str(step, buf));
118 if ((step >= ir->nstlist))
120 double seconds_since_epoch = gmx_gettime();
121 elapsed_seconds = seconds_since_epoch - walltime_accounting_get_start_time_stamp(walltime_accounting);
122 time_per_step = elapsed_seconds/(step - ir->init_step + 1);
123 dt = (ir->nsteps + ir->init_step - step) * time_per_step;
129 finish = (time_t) (seconds_since_epoch + dt);
130 gmx_ctime_r(&finish, timebuf, STRLEN);
131 sprintf(buf, "%s", timebuf);
132 buf[strlen(buf)-1] = '\0';
133 fprintf(out, ", will finish %s", buf);
137 fprintf(out, ", remaining wall clock time: %5d s ", (int)dt);
142 fprintf(out, " performance: %.1f ns/day ",
143 ir->delta_t/1000*24*60*60/time_per_step);
146 #ifndef GMX_THREAD_MPI
156 void print_date_and_time(FILE *fplog, int nodeid, const char *title,
159 char time_string[STRLEN];
168 char timebuf[STRLEN];
169 time_t temp_time = (time_t) the_time;
171 gmx_ctime_r(&temp_time, timebuf, STRLEN);
172 for (i = 0; timebuf[i] >= ' '; i++)
174 time_string[i] = timebuf[i];
176 time_string[i] = '\0';
179 fprintf(fplog, "%s on rank %d %s\n", title, nodeid, time_string);
182 void print_start(FILE *fplog, t_commrec *cr,
183 gmx_walltime_accounting_t walltime_accounting,
188 sprintf(buf, "Started %s", name);
189 print_date_and_time(fplog, cr->nodeid, buf,
190 walltime_accounting_get_start_time_stamp(walltime_accounting));
193 static void sum_forces(int start, int end, rvec f[], rvec flr[])
199 pr_rvecs(debug, 0, "fsr", f+start, end-start);
200 pr_rvecs(debug, 0, "flr", flr+start, end-start);
202 for (i = start; (i < end); i++)
204 rvec_inc(f[i], flr[i]);
209 * calc_f_el calculates forces due to an electric field.
211 * force is kJ mol^-1 nm^-1 = e * kJ mol^-1 nm^-1 / e
213 * Et[] contains the parameters for the time dependent
215 * Ex[] contains the parameters for
216 * the spatial dependent part of the field.
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-consistent 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)
282 /* The short-range virial from surrounding boxes */
284 calc_vir(SHIFTS, fr->shift_vec, fr->fshift, vir_part, ePBC == epbcSCREW, box);
285 inc_nrnb(nrnb, eNR_VIRIAL, SHIFTS);
287 /* Calculate partial virial, for local atoms only, based on short range.
288 * Total virial is computed in global_stat, called from do_md
290 f_calc_vir(start, start+homenr, x, f, vir_part, graph, box);
291 inc_nrnb(nrnb, eNR_VIRIAL, homenr);
293 /* Add position restraint contribution */
294 for (i = 0; i < DIM; i++)
296 vir_part[i][i] += fr->vir_diag_posres[i];
299 /* Add wall contribution */
300 for (i = 0; i < DIM; i++)
302 vir_part[i][ZZ] += fr->vir_wall_z[i];
307 pr_rvecs(debug, 0, "vir_part", vir_part, DIM);
311 static void pull_potential_wrapper(t_commrec *cr,
313 matrix box, rvec x[],
317 gmx_enerdata_t *enerd,
320 gmx_wallcycle_t wcycle)
325 /* Calculate the center of mass forces, this requires communication,
326 * which is why pull_potential is called close to other communication.
327 * The virial contribution is calculated directly,
328 * which is why we call pull_potential after calc_virial.
330 wallcycle_start(wcycle, ewcPULLPOT);
331 set_pbc(&pbc, ir->ePBC, box);
333 enerd->term[F_COM_PULL] +=
334 pull_potential(ir->ePull, ir->pull, mdatoms, &pbc,
335 cr, t, lambda[efptRESTRAINT], x, f, vir_force, &dvdl);
336 enerd->dvdl_lin[efptRESTRAINT] += dvdl;
337 wallcycle_stop(wcycle, ewcPULLPOT);
340 static void pme_receive_force_ener(t_commrec *cr,
341 gmx_wallcycle_t wcycle,
342 gmx_enerdata_t *enerd,
345 real e_q, e_lj, dvdl_q, dvdl_lj;
346 float cycles_ppdpme, cycles_seppme;
348 cycles_ppdpme = wallcycle_stop(wcycle, ewcPPDURINGPME);
349 dd_cycles_add(cr->dd, cycles_ppdpme, ddCyclPPduringPME);
351 /* In case of node-splitting, the PP nodes receive the long-range
352 * forces, virial and energy from the PME nodes here.
354 wallcycle_start(wcycle, ewcPP_PMEWAITRECVF);
357 gmx_pme_receive_f(cr, fr->f_novirsum, fr->vir_el_recip, &e_q,
358 fr->vir_lj_recip, &e_lj, &dvdl_q, &dvdl_lj,
360 enerd->term[F_COUL_RECIP] += e_q;
361 enerd->term[F_LJ_RECIP] += e_lj;
362 enerd->dvdl_lin[efptCOUL] += dvdl_q;
363 enerd->dvdl_lin[efptVDW] += dvdl_lj;
367 dd_cycles_add(cr->dd, cycles_seppme, ddCyclPME);
369 wallcycle_stop(wcycle, ewcPP_PMEWAITRECVF);
372 static void print_large_forces(FILE *fp, t_mdatoms *md, t_commrec *cr,
373 gmx_int64_t step, real pforce, rvec *x, rvec *f)
377 char buf[STEPSTRSIZE];
380 for (i = 0; i < md->homenr; i++)
383 /* We also catch NAN, if the compiler does not optimize this away. */
384 if (fn2 >= pf2 || fn2 != fn2)
386 fprintf(fp, "step %s atom %6d x %8.3f %8.3f %8.3f force %12.5e\n",
387 gmx_step_str(step, buf),
388 ddglatnr(cr->dd, i), x[i][XX], x[i][YY], x[i][ZZ], sqrt(fn2));
393 static void post_process_forces(t_commrec *cr,
395 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
397 matrix box, rvec x[],
402 t_forcerec *fr, gmx_vsite_t *vsite,
409 /* Spread the mesh force on virtual sites to the other particles...
410 * This is parallellized. MPI communication is performed
411 * if the constructing atoms aren't local.
413 wallcycle_start(wcycle, ewcVSITESPREAD);
414 spread_vsite_f(vsite, x, fr->f_novirsum, NULL,
415 (flags & GMX_FORCE_VIRIAL), fr->vir_el_recip,
417 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
418 wallcycle_stop(wcycle, ewcVSITESPREAD);
420 if (flags & GMX_FORCE_VIRIAL)
422 /* Now add the forces, this is local */
425 sum_forces(0, fr->f_novirsum_n, f, fr->f_novirsum);
429 sum_forces(0, mdatoms->homenr,
432 if (EEL_FULL(fr->eeltype))
434 /* Add the mesh contribution to the virial */
435 m_add(vir_force, fr->vir_el_recip, vir_force);
437 if (EVDW_PME(fr->vdwtype))
439 /* Add the mesh contribution to the virial */
440 m_add(vir_force, fr->vir_lj_recip, vir_force);
444 pr_rvecs(debug, 0, "vir_force", vir_force, DIM);
449 if (fr->print_force >= 0)
451 print_large_forces(stderr, mdatoms, cr, step, fr->print_force, x, f);
455 static void do_nb_verlet(t_forcerec *fr,
456 interaction_const_t *ic,
457 gmx_enerdata_t *enerd,
458 int flags, int ilocality,
461 gmx_wallcycle_t wcycle)
463 int enr_nbnxn_kernel_ljc, enr_nbnxn_kernel_lj;
464 nonbonded_verlet_group_t *nbvg;
467 if (!(flags & GMX_FORCE_NONBONDED))
469 /* skip non-bonded calculation */
473 nbvg = &fr->nbv->grp[ilocality];
475 /* CUDA kernel launch overhead is already timed separately */
476 if (fr->cutoff_scheme != ecutsVERLET)
478 gmx_incons("Invalid cut-off scheme passed!");
481 bCUDA = (nbvg->kernel_type == nbnxnk8x8x8_CUDA);
485 wallcycle_sub_start(wcycle, ewcsNONBONDED);
487 switch (nbvg->kernel_type)
489 case nbnxnk4x4_PlainC:
490 nbnxn_kernel_ref(&nbvg->nbl_lists,
496 enerd->grpp.ener[egCOULSR],
498 enerd->grpp.ener[egBHAMSR] :
499 enerd->grpp.ener[egLJSR]);
502 case nbnxnk4xN_SIMD_4xN:
503 nbnxn_kernel_simd_4xn(&nbvg->nbl_lists,
510 enerd->grpp.ener[egCOULSR],
512 enerd->grpp.ener[egBHAMSR] :
513 enerd->grpp.ener[egLJSR]);
515 case nbnxnk4xN_SIMD_2xNN:
516 nbnxn_kernel_simd_2xnn(&nbvg->nbl_lists,
523 enerd->grpp.ener[egCOULSR],
525 enerd->grpp.ener[egBHAMSR] :
526 enerd->grpp.ener[egLJSR]);
529 case nbnxnk8x8x8_CUDA:
530 nbnxn_cuda_launch_kernel(fr->nbv->cu_nbv, nbvg->nbat, flags, ilocality);
533 case nbnxnk8x8x8_PlainC:
534 nbnxn_kernel_gpu_ref(nbvg->nbl_lists.nbl[0],
539 nbvg->nbat->out[0].f,
541 enerd->grpp.ener[egCOULSR],
543 enerd->grpp.ener[egBHAMSR] :
544 enerd->grpp.ener[egLJSR]);
548 gmx_incons("Invalid nonbonded kernel type passed!");
553 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
556 if (EEL_RF(ic->eeltype) || ic->eeltype == eelCUT)
558 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_RF;
560 else if ((!bCUDA && nbvg->ewald_excl == ewaldexclAnalytical) ||
561 (bCUDA && nbnxn_cuda_is_kernel_ewald_analytical(fr->nbv->cu_nbv)))
563 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_EWALD;
567 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_TAB;
569 enr_nbnxn_kernel_lj = eNR_NBNXN_LJ;
570 if (flags & GMX_FORCE_ENERGY)
572 /* In eNR_??? the nbnxn F+E kernels are always the F kernel + 1 */
573 enr_nbnxn_kernel_ljc += 1;
574 enr_nbnxn_kernel_lj += 1;
577 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc,
578 nbvg->nbl_lists.natpair_ljq);
579 inc_nrnb(nrnb, enr_nbnxn_kernel_lj,
580 nbvg->nbl_lists.natpair_lj);
581 /* The Coulomb-only kernels are offset -eNR_NBNXN_LJ_RF+eNR_NBNXN_RF */
582 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc-eNR_NBNXN_LJ_RF+eNR_NBNXN_RF,
583 nbvg->nbl_lists.natpair_q);
585 if (ic->vdw_modifier == eintmodFORCESWITCH)
587 /* We add up the switch cost separately */
588 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_FSW+((flags & GMX_FORCE_ENERGY) ? 1 : 0),
589 nbvg->nbl_lists.natpair_ljq + nbvg->nbl_lists.natpair_lj);
591 if (ic->vdw_modifier == eintmodPOTSWITCH)
593 /* We add up the switch cost separately */
594 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_PSW+((flags & GMX_FORCE_ENERGY) ? 1 : 0),
595 nbvg->nbl_lists.natpair_ljq + nbvg->nbl_lists.natpair_lj);
597 if (ic->vdwtype == evdwPME)
599 /* We add up the LJ Ewald cost separately */
600 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_EWALD+((flags & GMX_FORCE_ENERGY) ? 1 : 0),
601 nbvg->nbl_lists.natpair_ljq + nbvg->nbl_lists.natpair_lj);
605 static void do_nb_verlet_fep(nbnxn_pairlist_set_t *nbl_lists,
612 gmx_enerdata_t *enerd,
615 gmx_wallcycle_t wcycle)
618 nb_kernel_data_t kernel_data;
620 real dvdl_nb[efptNR];
625 /* Add short-range interactions */
626 donb_flags |= GMX_NONBONDED_DO_SR;
628 /* Currently all group scheme kernels always calculate (shift-)forces */
629 if (flags & GMX_FORCE_FORCES)
631 donb_flags |= GMX_NONBONDED_DO_FORCE;
633 if (flags & GMX_FORCE_VIRIAL)
635 donb_flags |= GMX_NONBONDED_DO_SHIFTFORCE;
637 if (flags & GMX_FORCE_ENERGY)
639 donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
641 if (flags & GMX_FORCE_DO_LR)
643 donb_flags |= GMX_NONBONDED_DO_LR;
646 kernel_data.flags = donb_flags;
647 kernel_data.lambda = lambda;
648 kernel_data.dvdl = dvdl_nb;
650 kernel_data.energygrp_elec = enerd->grpp.ener[egCOULSR];
651 kernel_data.energygrp_vdw = enerd->grpp.ener[egLJSR];
653 /* reset free energy components */
654 for (i = 0; i < efptNR; i++)
659 assert(gmx_omp_nthreads_get(emntNonbonded) == nbl_lists->nnbl);
661 wallcycle_sub_start(wcycle, ewcsNONBONDED);
662 #pragma omp parallel for schedule(static) num_threads(nbl_lists->nnbl)
663 for (th = 0; th < nbl_lists->nnbl; th++)
665 gmx_nb_free_energy_kernel(nbl_lists->nbl_fep[th],
666 x, f, fr, mdatoms, &kernel_data, nrnb);
669 if (fepvals->sc_alpha != 0)
671 enerd->dvdl_nonlin[efptVDW] += dvdl_nb[efptVDW];
672 enerd->dvdl_nonlin[efptCOUL] += dvdl_nb[efptCOUL];
676 enerd->dvdl_lin[efptVDW] += dvdl_nb[efptVDW];
677 enerd->dvdl_lin[efptCOUL] += dvdl_nb[efptCOUL];
680 /* If we do foreign lambda and we have soft-core interactions
681 * we have to recalculate the (non-linear) energies contributions.
683 if (fepvals->n_lambda > 0 && (flags & GMX_FORCE_DHDL) && fepvals->sc_alpha != 0)
685 kernel_data.flags = (donb_flags & ~(GMX_NONBONDED_DO_FORCE | GMX_NONBONDED_DO_SHIFTFORCE)) | GMX_NONBONDED_DO_FOREIGNLAMBDA;
686 kernel_data.lambda = lam_i;
687 kernel_data.energygrp_elec = enerd->foreign_grpp.ener[egCOULSR];
688 kernel_data.energygrp_vdw = enerd->foreign_grpp.ener[egLJSR];
689 /* Note that we add to kernel_data.dvdl, but ignore the result */
691 for (i = 0; i < enerd->n_lambda; i++)
693 for (j = 0; j < efptNR; j++)
695 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
697 reset_foreign_enerdata(enerd);
698 #pragma omp parallel for schedule(static) num_threads(nbl_lists->nnbl)
699 for (th = 0; th < nbl_lists->nnbl; th++)
701 gmx_nb_free_energy_kernel(nbl_lists->nbl_fep[th],
702 x, f, fr, mdatoms, &kernel_data, nrnb);
705 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
706 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
710 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
713 gmx_bool use_GPU(const nonbonded_verlet_t *nbv)
715 return nbv != NULL && nbv->bUseGPU;
718 void do_force_cutsVERLET(FILE *fplog, t_commrec *cr,
719 t_inputrec *inputrec,
720 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
722 gmx_groups_t gmx_unused *groups,
723 matrix box, rvec x[], history_t *hist,
727 gmx_enerdata_t *enerd, t_fcdata *fcd,
728 real *lambda, t_graph *graph,
729 t_forcerec *fr, interaction_const_t *ic,
730 gmx_vsite_t *vsite, rvec mu_tot,
731 double t, FILE *field, gmx_edsam_t ed,
738 gmx_bool bStateChanged, bNS, bFillGrid, bCalcCGCM;
739 gmx_bool bDoLongRange, bDoForces, bSepLRF, bUseGPU, bUseOrEmulGPU;
740 gmx_bool bDiffKernels = FALSE;
741 rvec vzero, box_diag;
742 float cycles_pme, cycles_force, cycles_wait_gpu;
743 nonbonded_verlet_t *nbv;
750 homenr = mdatoms->homenr;
752 clear_mat(vir_force);
754 if (DOMAINDECOMP(cr))
756 cg1 = cr->dd->ncg_tot;
767 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
768 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
769 bFillGrid = (bNS && bStateChanged);
770 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
771 bDoLongRange = (fr->bTwinRange && bNS && (flags & GMX_FORCE_DO_LR));
772 bDoForces = (flags & GMX_FORCE_FORCES);
773 bSepLRF = (bDoLongRange && bDoForces && (flags & GMX_FORCE_SEPLRF));
774 bUseGPU = fr->nbv->bUseGPU;
775 bUseOrEmulGPU = bUseGPU || (nbv->grp[0].kernel_type == nbnxnk8x8x8_PlainC);
779 update_forcerec(fr, box);
781 if (NEED_MUTOT(*inputrec))
783 /* Calculate total (local) dipole moment in a temporary common array.
784 * This makes it possible to sum them over nodes faster.
786 calc_mu(start, homenr,
787 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
792 if (fr->ePBC != epbcNONE)
794 /* Compute shift vectors every step,
795 * because of pressure coupling or box deformation!
797 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
799 calc_shifts(box, fr->shift_vec);
804 put_atoms_in_box_omp(fr->ePBC, box, homenr, x);
805 inc_nrnb(nrnb, eNR_SHIFTX, homenr);
807 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
809 unshift_self(graph, box, x);
813 nbnxn_atomdata_copy_shiftvec(flags & GMX_FORCE_DYNAMICBOX,
814 fr->shift_vec, nbv->grp[0].nbat);
817 if (!(cr->duty & DUTY_PME))
822 /* Send particle coordinates to the pme nodes.
823 * Since this is only implemented for domain decomposition
824 * and domain decomposition does not use the graph,
825 * we do not need to worry about shifting.
830 wallcycle_start(wcycle, ewcPP_PMESENDX);
832 bBS = (inputrec->nwall == 2);
836 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
839 if (EEL_PME(fr->eeltype))
841 pme_flags |= GMX_PME_DO_COULOMB;
844 if (EVDW_PME(fr->vdwtype))
846 pme_flags |= GMX_PME_DO_LJ;
849 gmx_pme_send_coordinates(cr, bBS ? boxs : box, x,
850 mdatoms->nChargePerturbed, mdatoms->nTypePerturbed, lambda[efptCOUL], lambda[efptVDW],
851 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),
854 wallcycle_stop(wcycle, ewcPP_PMESENDX);
858 /* do gridding for pair search */
861 if (graph && bStateChanged)
863 /* Calculate intramolecular shift vectors to make molecules whole */
864 mk_mshift(fplog, graph, fr->ePBC, box, x);
868 box_diag[XX] = box[XX][XX];
869 box_diag[YY] = box[YY][YY];
870 box_diag[ZZ] = box[ZZ][ZZ];
872 wallcycle_start(wcycle, ewcNS);
875 wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
876 nbnxn_put_on_grid(nbv->nbs, fr->ePBC, box,
878 0, mdatoms->homenr, -1, fr->cginfo, x,
880 nbv->grp[eintLocal].kernel_type,
881 nbv->grp[eintLocal].nbat);
882 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
886 wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
887 nbnxn_put_on_grid_nonlocal(nbv->nbs, domdec_zones(cr->dd),
889 nbv->grp[eintNonlocal].kernel_type,
890 nbv->grp[eintNonlocal].nbat);
891 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
894 if (nbv->ngrp == 1 ||
895 nbv->grp[eintNonlocal].nbat == nbv->grp[eintLocal].nbat)
897 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatAll,
898 nbv->nbs, mdatoms, fr->cginfo);
902 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatLocal,
903 nbv->nbs, mdatoms, fr->cginfo);
904 nbnxn_atomdata_set(nbv->grp[eintNonlocal].nbat, eatAll,
905 nbv->nbs, mdatoms, fr->cginfo);
907 wallcycle_stop(wcycle, ewcNS);
910 /* initialize the GPU atom data and copy shift vector */
915 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
916 nbnxn_cuda_init_atomdata(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
917 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
920 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
921 nbnxn_cuda_upload_shiftvec(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
922 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
925 /* do local pair search */
928 wallcycle_start_nocount(wcycle, ewcNS);
929 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
930 nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintLocal].nbat,
933 nbv->min_ci_balanced,
934 &nbv->grp[eintLocal].nbl_lists,
936 nbv->grp[eintLocal].kernel_type,
938 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
942 /* initialize local pair-list on the GPU */
943 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
944 nbv->grp[eintLocal].nbl_lists.nbl[0],
947 wallcycle_stop(wcycle, ewcNS);
951 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
952 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
953 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, FALSE, x,
954 nbv->grp[eintLocal].nbat);
955 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
956 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
961 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
962 /* launch local nonbonded F on GPU */
963 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFNo,
965 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
968 /* Communicate coordinates and sum dipole if necessary +
969 do non-local pair search */
970 if (DOMAINDECOMP(cr))
972 bDiffKernels = (nbv->grp[eintNonlocal].kernel_type !=
973 nbv->grp[eintLocal].kernel_type);
977 /* With GPU+CPU non-bonded calculations we need to copy
978 * the local coordinates to the non-local nbat struct
979 * (in CPU format) as the non-local kernel call also
980 * calculates the local - non-local interactions.
982 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
983 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
984 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, TRUE, x,
985 nbv->grp[eintNonlocal].nbat);
986 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
987 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
992 wallcycle_start_nocount(wcycle, ewcNS);
993 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
997 nbnxn_grid_add_simple(nbv->nbs, nbv->grp[eintNonlocal].nbat);
1000 nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintNonlocal].nbat,
1003 nbv->min_ci_balanced,
1004 &nbv->grp[eintNonlocal].nbl_lists,
1006 nbv->grp[eintNonlocal].kernel_type,
1009 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1011 if (nbv->grp[eintNonlocal].kernel_type == nbnxnk8x8x8_CUDA)
1013 /* initialize non-local pair-list on the GPU */
1014 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
1015 nbv->grp[eintNonlocal].nbl_lists.nbl[0],
1018 wallcycle_stop(wcycle, ewcNS);
1022 wallcycle_start(wcycle, ewcMOVEX);
1023 dd_move_x(cr->dd, box, x);
1025 /* When we don't need the total dipole we sum it in global_stat */
1026 if (bStateChanged && NEED_MUTOT(*inputrec))
1028 gmx_sumd(2*DIM, mu, cr);
1030 wallcycle_stop(wcycle, ewcMOVEX);
1032 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1033 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
1034 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatNonlocal, FALSE, x,
1035 nbv->grp[eintNonlocal].nbat);
1036 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
1037 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1040 if (bUseGPU && !bDiffKernels)
1042 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
1043 /* launch non-local nonbonded F on GPU */
1044 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFNo,
1046 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1052 /* launch D2H copy-back F */
1053 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1054 if (DOMAINDECOMP(cr) && !bDiffKernels)
1056 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintNonlocal].nbat,
1057 flags, eatNonlocal);
1059 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintLocal].nbat,
1061 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1064 if (bStateChanged && NEED_MUTOT(*inputrec))
1068 gmx_sumd(2*DIM, mu, cr);
1071 for (i = 0; i < 2; i++)
1073 for (j = 0; j < DIM; j++)
1075 fr->mu_tot[i][j] = mu[i*DIM + j];
1079 if (fr->efep == efepNO)
1081 copy_rvec(fr->mu_tot[0], mu_tot);
1085 for (j = 0; j < DIM; j++)
1088 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] +
1089 lambda[efptCOUL]*fr->mu_tot[1][j];
1093 /* Reset energies */
1094 reset_enerdata(fr, bNS, enerd, MASTER(cr));
1095 clear_rvecs(SHIFTS, fr->fshift);
1097 if (DOMAINDECOMP(cr) && !(cr->duty & DUTY_PME))
1099 wallcycle_start(wcycle, ewcPPDURINGPME);
1100 dd_force_flop_start(cr->dd, nrnb);
1105 /* Enforced rotation has its own cycle counter that starts after the collective
1106 * coordinates have been communicated. It is added to ddCyclF to allow
1107 * for proper load-balancing */
1108 wallcycle_start(wcycle, ewcROT);
1109 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1110 wallcycle_stop(wcycle, ewcROT);
1113 /* Start the force cycle counter.
1114 * This counter is stopped in do_forcelow_level.
1115 * No parallel communication should occur while this counter is running,
1116 * since that will interfere with the dynamic load balancing.
1118 wallcycle_start(wcycle, ewcFORCE);
1121 /* Reset forces for which the virial is calculated separately:
1122 * PME/Ewald forces if necessary */
1123 if (fr->bF_NoVirSum)
1125 if (flags & GMX_FORCE_VIRIAL)
1127 fr->f_novirsum = fr->f_novirsum_alloc;
1130 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1134 clear_rvecs(homenr, fr->f_novirsum+start);
1139 /* We are not calculating the pressure so we do not need
1140 * a separate array for forces that do not contribute
1147 /* Clear the short- and long-range forces */
1148 clear_rvecs(fr->natoms_force_constr, f);
1149 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1151 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1154 clear_rvec(fr->vir_diag_posres);
1157 if (inputrec->ePull == epullCONSTRAINT)
1159 clear_pull_forces(inputrec->pull);
1162 /* We calculate the non-bonded forces, when done on the CPU, here.
1163 * We do this before calling do_force_lowlevel, because in that
1164 * function, the listed forces are calculated before PME, which
1165 * does communication. With this order, non-bonded and listed
1166 * force calculation imbalance can be balanced out by the domain
1167 * decomposition load balancing.
1172 /* Maybe we should move this into do_force_lowlevel */
1173 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFYes,
1177 if (fr->efep != efepNO)
1179 /* Calculate the local and non-local free energy interactions here.
1180 * Happens here on the CPU both with and without GPU.
1182 if (fr->nbv->grp[eintLocal].nbl_lists.nbl_fep[0]->nrj > 0)
1184 do_nb_verlet_fep(&fr->nbv->grp[eintLocal].nbl_lists,
1186 inputrec->fepvals, lambda,
1187 enerd, flags, nrnb, wcycle);
1190 if (DOMAINDECOMP(cr) &&
1191 fr->nbv->grp[eintNonlocal].nbl_lists.nbl_fep[0]->nrj > 0)
1193 do_nb_verlet_fep(&fr->nbv->grp[eintNonlocal].nbl_lists,
1195 inputrec->fepvals, lambda,
1196 enerd, flags, nrnb, wcycle);
1200 if (!bUseOrEmulGPU || bDiffKernels)
1204 if (DOMAINDECOMP(cr))
1206 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal,
1207 bDiffKernels ? enbvClearFYes : enbvClearFNo,
1217 aloc = eintNonlocal;
1220 /* Add all the non-bonded force to the normal force array.
1221 * This can be split into a local a non-local part when overlapping
1222 * communication with calculation with domain decomposition.
1224 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1225 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1226 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1227 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatAll, nbv->grp[aloc].nbat, f);
1228 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1229 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1230 wallcycle_start_nocount(wcycle, ewcFORCE);
1232 /* if there are multiple fshift output buffers reduce them */
1233 if ((flags & GMX_FORCE_VIRIAL) &&
1234 nbv->grp[aloc].nbl_lists.nnbl > 1)
1236 nbnxn_atomdata_add_nbat_fshift_to_fshift(nbv->grp[aloc].nbat,
1241 /* update QMMMrec, if necessary */
1244 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1247 /* Compute the bonded and non-bonded energies and optionally forces */
1248 do_force_lowlevel(fr, inputrec, &(top->idef),
1249 cr, nrnb, wcycle, mdatoms,
1250 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, top, fr->born,
1252 inputrec->fepvals, lambda, graph, &(top->excls), fr->mu_tot,
1253 flags, &cycles_pme);
1257 if (do_per_step(step, inputrec->nstcalclr))
1259 /* Add the long range forces to the short range forces */
1260 for (i = 0; i < fr->natoms_force_constr; i++)
1262 rvec_add(fr->f_twin[i], f[i], f[i]);
1267 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1271 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1274 if (bUseOrEmulGPU && !bDiffKernels)
1276 /* wait for non-local forces (or calculate in emulation mode) */
1277 if (DOMAINDECOMP(cr))
1283 wallcycle_start(wcycle, ewcWAIT_GPU_NB_NL);
1284 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1285 nbv->grp[eintNonlocal].nbat,
1287 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1289 cycles_tmp = wallcycle_stop(wcycle, ewcWAIT_GPU_NB_NL);
1290 cycles_wait_gpu += cycles_tmp;
1291 cycles_force += cycles_tmp;
1295 wallcycle_start_nocount(wcycle, ewcFORCE);
1296 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFYes,
1298 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1300 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1301 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1302 /* skip the reduction if there was no non-local work to do */
1303 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1305 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatNonlocal,
1306 nbv->grp[eintNonlocal].nbat, f);
1308 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1309 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1313 if (bDoForces && DOMAINDECOMP(cr))
1317 /* We are done with the CPU compute, but the GPU local non-bonded
1318 * kernel can still be running while we communicate the forces.
1319 * We start a counter here, so we can, hopefully, time the rest
1320 * of the GPU kernel execution and data transfer.
1322 wallcycle_start(wcycle, ewcWAIT_GPU_NB_L_EST);
1325 /* Communicate the forces */
1326 wallcycle_start(wcycle, ewcMOVEF);
1327 dd_move_f(cr->dd, f, fr->fshift);
1330 /* We should not update the shift forces here,
1331 * since f_twin is already included in f.
1333 dd_move_f(cr->dd, fr->f_twin, NULL);
1335 wallcycle_stop(wcycle, ewcMOVEF);
1340 /* wait for local forces (or calculate in emulation mode) */
1343 float cycles_tmp, cycles_wait_est;
1344 const float cuda_api_overhead_margin = 50000.0f; /* cycles */
1346 wallcycle_start(wcycle, ewcWAIT_GPU_NB_L);
1347 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1348 nbv->grp[eintLocal].nbat,
1350 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1352 cycles_tmp = wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
1354 if (bDoForces && DOMAINDECOMP(cr))
1356 cycles_wait_est = wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L_EST);
1358 if (cycles_tmp < cuda_api_overhead_margin)
1360 /* We measured few cycles, it could be that the kernel
1361 * and transfer finished earlier and there was no actual
1362 * wait time, only API call overhead.
1363 * Then the actual time could be anywhere between 0 and
1364 * cycles_wait_est. As a compromise, we use half the time.
1366 cycles_wait_est *= 0.5f;
1371 /* No force communication so we actually timed the wait */
1372 cycles_wait_est = cycles_tmp;
1374 /* Even though this is after dd_move_f, the actual task we are
1375 * waiting for runs asynchronously with dd_move_f and we usually
1376 * have nothing to balance it with, so we can and should add
1377 * the time to the force time for load balancing.
1379 cycles_force += cycles_wait_est;
1380 cycles_wait_gpu += cycles_wait_est;
1382 /* now clear the GPU outputs while we finish the step on the CPU */
1384 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1385 nbnxn_cuda_clear_outputs(nbv->cu_nbv, flags);
1386 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1390 wallcycle_start_nocount(wcycle, ewcFORCE);
1391 do_nb_verlet(fr, ic, enerd, flags, eintLocal,
1392 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
1394 wallcycle_stop(wcycle, ewcFORCE);
1396 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1397 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1398 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1400 /* skip the reduction if there was no non-local work to do */
1401 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatLocal,
1402 nbv->grp[eintLocal].nbat, f);
1404 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1405 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1408 if (DOMAINDECOMP(cr))
1410 dd_force_flop_stop(cr->dd, nrnb);
1413 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1416 dd_cycles_add(cr->dd, cycles_wait_gpu, ddCyclWaitGPU);
1423 if (IR_ELEC_FIELD(*inputrec))
1425 /* Compute forces due to electric field */
1426 calc_f_el(MASTER(cr) ? field : NULL,
1427 start, homenr, mdatoms->chargeA, fr->f_novirsum,
1428 inputrec->ex, inputrec->et, t);
1431 /* If we have NoVirSum forces, but we do not calculate the virial,
1432 * we sum fr->f_novirsum=f later.
1434 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1436 wallcycle_start(wcycle, ewcVSITESPREAD);
1437 spread_vsite_f(vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1438 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1439 wallcycle_stop(wcycle, ewcVSITESPREAD);
1443 wallcycle_start(wcycle, ewcVSITESPREAD);
1444 spread_vsite_f(vsite, x, fr->f_twin, NULL, FALSE, NULL,
1446 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1447 wallcycle_stop(wcycle, ewcVSITESPREAD);
1451 if (flags & GMX_FORCE_VIRIAL)
1453 /* Calculation of the virial must be done after vsites! */
1454 calc_virial(0, mdatoms->homenr, x, f,
1455 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1459 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1461 /* Since the COM pulling is always done mass-weighted, no forces are
1462 * applied to vsites and this call can be done after vsite spreading.
1464 pull_potential_wrapper(cr, inputrec, box, x,
1465 f, vir_force, mdatoms, enerd, lambda, t,
1469 /* Add the forces from enforced rotation potentials (if any) */
1472 wallcycle_start(wcycle, ewcROTadd);
1473 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
1474 wallcycle_stop(wcycle, ewcROTadd);
1477 /* Add forces from interactive molecular dynamics (IMD), if bIMD == TRUE. */
1478 IMD_apply_forces(inputrec->bIMD, inputrec->imd, cr, f, wcycle);
1480 if (PAR(cr) && !(cr->duty & DUTY_PME))
1482 /* In case of node-splitting, the PP nodes receive the long-range
1483 * forces, virial and energy from the PME nodes here.
1485 pme_receive_force_ener(cr, wcycle, enerd, fr);
1490 post_process_forces(cr, step, nrnb, wcycle,
1491 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
1495 /* Sum the potential energy terms from group contributions */
1496 sum_epot(&(enerd->grpp), enerd->term);
1499 void do_force_cutsGROUP(FILE *fplog, t_commrec *cr,
1500 t_inputrec *inputrec,
1501 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1502 gmx_localtop_t *top,
1503 gmx_groups_t *groups,
1504 matrix box, rvec x[], history_t *hist,
1508 gmx_enerdata_t *enerd, t_fcdata *fcd,
1509 real *lambda, t_graph *graph,
1510 t_forcerec *fr, gmx_vsite_t *vsite, rvec mu_tot,
1511 double t, FILE *field, gmx_edsam_t ed,
1512 gmx_bool bBornRadii,
1518 gmx_bool bStateChanged, bNS, bFillGrid, bCalcCGCM;
1519 gmx_bool bDoLongRangeNS, bDoForces, bSepLRF;
1520 gmx_bool bDoAdressWF;
1522 float cycles_pme, cycles_force;
1525 homenr = mdatoms->homenr;
1527 clear_mat(vir_force);
1530 if (DOMAINDECOMP(cr))
1532 cg1 = cr->dd->ncg_tot;
1543 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
1544 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
1545 /* Should we update the long-range neighborlists at this step? */
1546 bDoLongRangeNS = fr->bTwinRange && bNS;
1547 /* Should we perform the long-range nonbonded evaluation inside the neighborsearching? */
1548 bFillGrid = (bNS && bStateChanged);
1549 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
1550 bDoForces = (flags & GMX_FORCE_FORCES);
1551 bSepLRF = ((inputrec->nstcalclr > 1) && bDoForces &&
1552 (flags & GMX_FORCE_SEPLRF) && (flags & GMX_FORCE_DO_LR));
1554 /* should probably move this to the forcerec since it doesn't change */
1555 bDoAdressWF = ((fr->adress_type != eAdressOff));
1559 update_forcerec(fr, box);
1561 if (NEED_MUTOT(*inputrec))
1563 /* Calculate total (local) dipole moment in a temporary common array.
1564 * This makes it possible to sum them over nodes faster.
1566 calc_mu(start, homenr,
1567 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
1572 if (fr->ePBC != epbcNONE)
1574 /* Compute shift vectors every step,
1575 * because of pressure coupling or box deformation!
1577 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
1579 calc_shifts(box, fr->shift_vec);
1584 put_charge_groups_in_box(fplog, cg0, cg1, fr->ePBC, box,
1585 &(top->cgs), x, fr->cg_cm);
1586 inc_nrnb(nrnb, eNR_CGCM, homenr);
1587 inc_nrnb(nrnb, eNR_RESETX, cg1-cg0);
1589 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
1591 unshift_self(graph, box, x);
1596 calc_cgcm(fplog, cg0, cg1, &(top->cgs), x, fr->cg_cm);
1597 inc_nrnb(nrnb, eNR_CGCM, homenr);
1600 if (bCalcCGCM && gmx_debug_at)
1602 pr_rvecs(debug, 0, "cgcm", fr->cg_cm, top->cgs.nr);
1606 if (!(cr->duty & DUTY_PME))
1611 /* Send particle coordinates to the pme nodes.
1612 * Since this is only implemented for domain decomposition
1613 * and domain decomposition does not use the graph,
1614 * we do not need to worry about shifting.
1619 wallcycle_start(wcycle, ewcPP_PMESENDX);
1621 bBS = (inputrec->nwall == 2);
1624 copy_mat(box, boxs);
1625 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
1628 if (EEL_PME(fr->eeltype))
1630 pme_flags |= GMX_PME_DO_COULOMB;
1633 if (EVDW_PME(fr->vdwtype))
1635 pme_flags |= GMX_PME_DO_LJ;
1638 gmx_pme_send_coordinates(cr, bBS ? boxs : box, x,
1639 mdatoms->nChargePerturbed, mdatoms->nTypePerturbed, lambda[efptCOUL], lambda[efptVDW],
1640 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),
1643 wallcycle_stop(wcycle, ewcPP_PMESENDX);
1645 #endif /* GMX_MPI */
1647 /* Communicate coordinates and sum dipole if necessary */
1648 if (DOMAINDECOMP(cr))
1650 wallcycle_start(wcycle, ewcMOVEX);
1651 dd_move_x(cr->dd, box, x);
1652 wallcycle_stop(wcycle, ewcMOVEX);
1655 /* update adress weight beforehand */
1656 if (bStateChanged && bDoAdressWF)
1658 /* need pbc for adress weight calculation with pbc_dx */
1659 set_pbc(&pbc, inputrec->ePBC, box);
1660 if (fr->adress_site == eAdressSITEcog)
1662 update_adress_weights_cog(top->idef.iparams, top->idef.il, x, fr, mdatoms,
1663 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1665 else if (fr->adress_site == eAdressSITEcom)
1667 update_adress_weights_com(fplog, cg0, cg1, &(top->cgs), x, fr, mdatoms,
1668 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1670 else if (fr->adress_site == eAdressSITEatomatom)
1672 update_adress_weights_atom_per_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1673 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1677 update_adress_weights_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1678 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1682 if (NEED_MUTOT(*inputrec))
1689 gmx_sumd(2*DIM, mu, cr);
1691 for (i = 0; i < 2; i++)
1693 for (j = 0; j < DIM; j++)
1695 fr->mu_tot[i][j] = mu[i*DIM + j];
1699 if (fr->efep == efepNO)
1701 copy_rvec(fr->mu_tot[0], mu_tot);
1705 for (j = 0; j < DIM; j++)
1708 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] + lambda[efptCOUL]*fr->mu_tot[1][j];
1713 /* Reset energies */
1714 reset_enerdata(fr, bNS, enerd, MASTER(cr));
1715 clear_rvecs(SHIFTS, fr->fshift);
1719 wallcycle_start(wcycle, ewcNS);
1721 if (graph && bStateChanged)
1723 /* Calculate intramolecular shift vectors to make molecules whole */
1724 mk_mshift(fplog, graph, fr->ePBC, box, x);
1727 /* Do the actual neighbour searching */
1729 groups, top, mdatoms,
1730 cr, nrnb, bFillGrid,
1733 wallcycle_stop(wcycle, ewcNS);
1736 if (inputrec->implicit_solvent && bNS)
1738 make_gb_nblist(cr, inputrec->gb_algorithm,
1739 x, box, fr, &top->idef, graph, fr->born);
1742 if (DOMAINDECOMP(cr) && !(cr->duty & DUTY_PME))
1744 wallcycle_start(wcycle, ewcPPDURINGPME);
1745 dd_force_flop_start(cr->dd, nrnb);
1750 /* Enforced rotation has its own cycle counter that starts after the collective
1751 * coordinates have been communicated. It is added to ddCyclF to allow
1752 * for proper load-balancing */
1753 wallcycle_start(wcycle, ewcROT);
1754 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1755 wallcycle_stop(wcycle, ewcROT);
1758 /* Start the force cycle counter.
1759 * This counter is stopped in do_forcelow_level.
1760 * No parallel communication should occur while this counter is running,
1761 * since that will interfere with the dynamic load balancing.
1763 wallcycle_start(wcycle, ewcFORCE);
1767 /* Reset forces for which the virial is calculated separately:
1768 * PME/Ewald forces if necessary */
1769 if (fr->bF_NoVirSum)
1771 if (flags & GMX_FORCE_VIRIAL)
1773 fr->f_novirsum = fr->f_novirsum_alloc;
1776 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1780 clear_rvecs(homenr, fr->f_novirsum+start);
1785 /* We are not calculating the pressure so we do not need
1786 * a separate array for forces that do not contribute
1793 /* Clear the short- and long-range forces */
1794 clear_rvecs(fr->natoms_force_constr, f);
1795 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1797 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1800 clear_rvec(fr->vir_diag_posres);
1802 if (inputrec->ePull == epullCONSTRAINT)
1804 clear_pull_forces(inputrec->pull);
1807 /* update QMMMrec, if necessary */
1810 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1813 /* Compute the bonded and non-bonded energies and optionally forces */
1814 do_force_lowlevel(fr, inputrec, &(top->idef),
1815 cr, nrnb, wcycle, mdatoms,
1816 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, top, fr->born,
1818 inputrec->fepvals, lambda,
1819 graph, &(top->excls), fr->mu_tot,
1825 if (do_per_step(step, inputrec->nstcalclr))
1827 /* Add the long range forces to the short range forces */
1828 for (i = 0; i < fr->natoms_force_constr; i++)
1830 rvec_add(fr->f_twin[i], f[i], f[i]);
1835 cycles_force = wallcycle_stop(wcycle, ewcFORCE);
1839 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1842 if (DOMAINDECOMP(cr))
1844 dd_force_flop_stop(cr->dd, nrnb);
1847 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1853 if (IR_ELEC_FIELD(*inputrec))
1855 /* Compute forces due to electric field */
1856 calc_f_el(MASTER(cr) ? field : NULL,
1857 start, homenr, mdatoms->chargeA, fr->f_novirsum,
1858 inputrec->ex, inputrec->et, t);
1861 if (bDoAdressWF && fr->adress_icor == eAdressICThermoForce)
1863 /* Compute thermodynamic force in hybrid AdResS region */
1864 adress_thermo_force(start, homenr, &(top->cgs), x, fr->f_novirsum, fr, mdatoms,
1865 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1868 /* Communicate the forces */
1869 if (DOMAINDECOMP(cr))
1871 wallcycle_start(wcycle, ewcMOVEF);
1872 dd_move_f(cr->dd, f, fr->fshift);
1873 /* Do we need to communicate the separate force array
1874 * for terms that do not contribute to the single sum virial?
1875 * Position restraints and electric fields do not introduce
1876 * inter-cg forces, only full electrostatics methods do.
1877 * When we do not calculate the virial, fr->f_novirsum = f,
1878 * so we have already communicated these forces.
1880 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1881 (flags & GMX_FORCE_VIRIAL))
1883 dd_move_f(cr->dd, fr->f_novirsum, NULL);
1887 /* We should not update the shift forces here,
1888 * since f_twin is already included in f.
1890 dd_move_f(cr->dd, fr->f_twin, NULL);
1892 wallcycle_stop(wcycle, ewcMOVEF);
1895 /* If we have NoVirSum forces, but we do not calculate the virial,
1896 * we sum fr->f_novirsum=f later.
1898 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1900 wallcycle_start(wcycle, ewcVSITESPREAD);
1901 spread_vsite_f(vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1902 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1903 wallcycle_stop(wcycle, ewcVSITESPREAD);
1907 wallcycle_start(wcycle, ewcVSITESPREAD);
1908 spread_vsite_f(vsite, x, fr->f_twin, NULL, FALSE, NULL,
1910 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1911 wallcycle_stop(wcycle, ewcVSITESPREAD);
1915 if (flags & GMX_FORCE_VIRIAL)
1917 /* Calculation of the virial must be done after vsites! */
1918 calc_virial(0, mdatoms->homenr, x, f,
1919 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1923 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1925 pull_potential_wrapper(cr, inputrec, box, x,
1926 f, vir_force, mdatoms, enerd, lambda, t,
1930 /* Add the forces from enforced rotation potentials (if any) */
1933 wallcycle_start(wcycle, ewcROTadd);
1934 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
1935 wallcycle_stop(wcycle, ewcROTadd);
1938 /* Add forces from interactive molecular dynamics (IMD), if bIMD == TRUE. */
1939 IMD_apply_forces(inputrec->bIMD, inputrec->imd, cr, f, wcycle);
1941 if (PAR(cr) && !(cr->duty & DUTY_PME))
1943 /* In case of node-splitting, the PP nodes receive the long-range
1944 * forces, virial and energy from the PME nodes here.
1946 pme_receive_force_ener(cr, wcycle, enerd, fr);
1951 post_process_forces(cr, step, nrnb, wcycle,
1952 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
1956 /* Sum the potential energy terms from group contributions */
1957 sum_epot(&(enerd->grpp), enerd->term);
1960 void do_force(FILE *fplog, t_commrec *cr,
1961 t_inputrec *inputrec,
1962 gmx_int64_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1963 gmx_localtop_t *top,
1964 gmx_groups_t *groups,
1965 matrix box, rvec x[], history_t *hist,
1969 gmx_enerdata_t *enerd, t_fcdata *fcd,
1970 real *lambda, t_graph *graph,
1972 gmx_vsite_t *vsite, rvec mu_tot,
1973 double t, FILE *field, gmx_edsam_t ed,
1974 gmx_bool bBornRadii,
1977 /* modify force flag if not doing nonbonded */
1978 if (!fr->bNonbonded)
1980 flags &= ~GMX_FORCE_NONBONDED;
1983 switch (inputrec->cutoff_scheme)
1986 do_force_cutsVERLET(fplog, cr, inputrec,
2002 do_force_cutsGROUP(fplog, cr, inputrec,
2017 gmx_incons("Invalid cut-off scheme passed!");
2022 void do_constrain_first(FILE *fplog, gmx_constr_t constr,
2023 t_inputrec *ir, t_mdatoms *md,
2024 t_state *state, t_commrec *cr, t_nrnb *nrnb,
2025 t_forcerec *fr, gmx_localtop_t *top)
2027 int i, m, start, end;
2029 real dt = ir->delta_t;
2033 snew(savex, state->natoms);
2040 fprintf(debug, "vcm: start=%d, homenr=%d, end=%d\n",
2041 start, md->homenr, end);
2043 /* Do a first constrain to reset particles... */
2044 step = ir->init_step;
2047 char buf[STEPSTRSIZE];
2048 fprintf(fplog, "\nConstraining the starting coordinates (step %s)\n",
2049 gmx_step_str(step, buf));
2053 /* constrain the current position */
2054 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2055 ir, cr, step, 0, 1.0, md,
2056 state->x, state->x, NULL,
2057 fr->bMolPBC, state->box,
2058 state->lambda[efptBONDED], &dvdl_dum,
2059 NULL, NULL, nrnb, econqCoord);
2062 /* constrain the inital velocity, and save it */
2063 /* also may be useful if we need the ekin from the halfstep for velocity verlet */
2064 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2065 ir, cr, step, 0, 1.0, md,
2066 state->x, state->v, state->v,
2067 fr->bMolPBC, state->box,
2068 state->lambda[efptBONDED], &dvdl_dum,
2069 NULL, NULL, nrnb, econqVeloc);
2071 /* constrain the inital velocities at t-dt/2 */
2072 if (EI_STATE_VELOCITY(ir->eI) && ir->eI != eiVV)
2074 for (i = start; (i < end); i++)
2076 for (m = 0; (m < DIM); m++)
2078 /* Reverse the velocity */
2079 state->v[i][m] = -state->v[i][m];
2080 /* Store the position at t-dt in buf */
2081 savex[i][m] = state->x[i][m] + dt*state->v[i][m];
2084 /* Shake the positions at t=-dt with the positions at t=0
2085 * as reference coordinates.
2089 char buf[STEPSTRSIZE];
2090 fprintf(fplog, "\nConstraining the coordinates at t0-dt (step %s)\n",
2091 gmx_step_str(step, buf));
2094 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2095 ir, cr, step, -1, 1.0, md,
2096 state->x, savex, NULL,
2097 fr->bMolPBC, state->box,
2098 state->lambda[efptBONDED], &dvdl_dum,
2099 state->v, NULL, nrnb, econqCoord);
2101 for (i = start; i < end; i++)
2103 for (m = 0; m < DIM; m++)
2105 /* Re-reverse the velocities */
2106 state->v[i][m] = -state->v[i][m];
2115 integrate_table(real vdwtab[], real scale, int offstart, int rstart, int rend,
2116 double *enerout, double *virout)
2118 double enersum, virsum;
2119 double invscale, invscale2, invscale3;
2120 double r, ea, eb, ec, pa, pb, pc, pd;
2125 invscale = 1.0/scale;
2126 invscale2 = invscale*invscale;
2127 invscale3 = invscale*invscale2;
2129 /* Following summation derived from cubic spline definition,
2130 * Numerical Recipies in C, second edition, p. 113-116. Exact for
2131 * the cubic spline. We first calculate the negative of the
2132 * energy from rvdw to rvdw_switch, assuming that g(r)=1, and then
2133 * add the more standard, abrupt cutoff correction to that result,
2134 * yielding the long-range correction for a switched function. We
2135 * perform both the pressure and energy loops at the same time for
2136 * simplicity, as the computational cost is low. */
2140 /* Since the dispersion table has been scaled down a factor
2141 * 6.0 and the repulsion a factor 12.0 to compensate for the
2142 * c6/c12 parameters inside nbfp[] being scaled up (to save
2143 * flops in kernels), we need to correct for this.
2154 for (ri = rstart; ri < rend; ++ri)
2158 eb = 2.0*invscale2*r;
2162 pb = 3.0*invscale2*r;
2163 pc = 3.0*invscale*r*r;
2166 /* this "8" is from the packing in the vdwtab array - perhaps
2167 should be defined? */
2169 offset = 8*ri + offstart;
2170 y0 = vdwtab[offset];
2171 f = vdwtab[offset+1];
2172 g = vdwtab[offset+2];
2173 h = vdwtab[offset+3];
2175 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);
2176 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);
2178 *enerout = 4.0*M_PI*enersum*tabfactor;
2179 *virout = 4.0*M_PI*virsum*tabfactor;
2182 void calc_enervirdiff(FILE *fplog, int eDispCorr, t_forcerec *fr)
2184 double eners[2], virs[2], enersum, virsum;
2185 double r0, rc3, rc9;
2187 real scale, *vdwtab;
2189 fr->enershiftsix = 0;
2190 fr->enershifttwelve = 0;
2191 fr->enerdiffsix = 0;
2192 fr->enerdifftwelve = 0;
2194 fr->virdifftwelve = 0;
2196 if (eDispCorr != edispcNO)
2198 for (i = 0; i < 2; i++)
2203 if ((fr->vdw_modifier == eintmodPOTSHIFT) ||
2204 (fr->vdw_modifier == eintmodPOTSWITCH) ||
2205 (fr->vdw_modifier == eintmodFORCESWITCH) ||
2206 (fr->vdwtype == evdwSHIFT) ||
2207 (fr->vdwtype == evdwSWITCH))
2209 if (((fr->vdw_modifier == eintmodPOTSWITCH) ||
2210 (fr->vdw_modifier == eintmodFORCESWITCH) ||
2211 (fr->vdwtype == evdwSWITCH)) && fr->rvdw_switch == 0)
2214 "With dispersion correction rvdw-switch can not be zero "
2215 "for vdw-type = %s", evdw_names[fr->vdwtype]);
2218 scale = fr->nblists[0].table_vdw.scale;
2219 vdwtab = fr->nblists[0].table_vdw.data;
2221 /* Round the cut-offs to exact table values for precision */
2222 ri0 = static_cast<int>(floor(fr->rvdw_switch*scale));
2223 ri1 = static_cast<int>(ceil(fr->rvdw*scale));
2225 /* The code below has some support for handling force-switching, i.e.
2226 * when the force (instead of potential) is switched over a limited
2227 * region. This leads to a constant shift in the potential inside the
2228 * switching region, which we can handle by adding a constant energy
2229 * term in the force-switch case just like when we do potential-shift.
2231 * For now this is not enabled, but to keep the functionality in the
2232 * code we check separately for switch and shift. When we do force-switch
2233 * the shifting point is rvdw_switch, while it is the cutoff when we
2234 * have a classical potential-shift.
2236 * For a pure potential-shift the potential has a constant shift
2237 * all the way out to the cutoff, and that is it. For other forms
2238 * we need to calculate the constant shift up to the point where we
2239 * start modifying the potential.
2241 ri0 = (fr->vdw_modifier == eintmodPOTSHIFT) ? ri1 : ri0;
2247 if ((fr->vdw_modifier == eintmodFORCESWITCH) ||
2248 (fr->vdwtype == evdwSHIFT))
2250 /* Determine the constant energy shift below rvdw_switch.
2251 * Table has a scale factor since we have scaled it down to compensate
2252 * for scaling-up c6/c12 with the derivative factors to save flops in analytical kernels.
2254 fr->enershiftsix = (real)(-1.0/(rc3*rc3)) - 6.0*vdwtab[8*ri0];
2255 fr->enershifttwelve = (real)( 1.0/(rc9*rc3)) - 12.0*vdwtab[8*ri0 + 4];
2257 else if (fr->vdw_modifier == eintmodPOTSHIFT)
2259 fr->enershiftsix = (real)(-1.0/(rc3*rc3));
2260 fr->enershifttwelve = (real)( 1.0/(rc9*rc3));
2263 /* Add the constant part from 0 to rvdw_switch.
2264 * This integration from 0 to rvdw_switch overcounts the number
2265 * of interactions by 1, as it also counts the self interaction.
2266 * We will correct for this later.
2268 eners[0] += 4.0*M_PI*fr->enershiftsix*rc3/3.0;
2269 eners[1] += 4.0*M_PI*fr->enershifttwelve*rc3/3.0;
2271 /* Calculate the contribution in the range [r0,r1] where we
2272 * modify the potential. For a pure potential-shift modifier we will
2273 * have ri0==ri1, and there will not be any contribution here.
2275 for (i = 0; i < 2; i++)
2279 integrate_table(vdwtab, scale, (i == 0 ? 0 : 4), ri0, ri1, &enersum, &virsum);
2280 eners[i] -= enersum;
2284 /* Alright: Above we compensated by REMOVING the parts outside r0
2285 * corresponding to the ideal VdW 1/r6 and /r12 potentials.
2287 * Regardless of whether r0 is the point where we start switching,
2288 * or the cutoff where we calculated the constant shift, we include
2289 * all the parts we are missing out to infinity from r0 by
2290 * calculating the analytical dispersion correction.
2292 eners[0] += -4.0*M_PI/(3.0*rc3);
2293 eners[1] += 4.0*M_PI/(9.0*rc9);
2294 virs[0] += 8.0*M_PI/rc3;
2295 virs[1] += -16.0*M_PI/(3.0*rc9);
2297 else if (fr->vdwtype == evdwCUT ||
2298 EVDW_PME(fr->vdwtype) ||
2299 fr->vdwtype == evdwUSER)
2301 if (fr->vdwtype == evdwUSER && fplog)
2304 "WARNING: using dispersion correction with user tables\n");
2307 /* Note that with LJ-PME, the dispersion correction is multiplied
2308 * by the difference between the actual C6 and the value of C6
2309 * that would produce the combination rule.
2310 * This means the normal energy and virial difference formulas
2314 rc3 = fr->rvdw*fr->rvdw*fr->rvdw;
2316 /* Contribution beyond the cut-off */
2317 eners[0] += -4.0*M_PI/(3.0*rc3);
2318 eners[1] += 4.0*M_PI/(9.0*rc9);
2319 if (fr->vdw_modifier == eintmodPOTSHIFT)
2321 /* Contribution within the cut-off */
2322 eners[0] += -4.0*M_PI/(3.0*rc3);
2323 eners[1] += 4.0*M_PI/(3.0*rc9);
2325 /* Contribution beyond the cut-off */
2326 virs[0] += 8.0*M_PI/rc3;
2327 virs[1] += -16.0*M_PI/(3.0*rc9);
2332 "Dispersion correction is not implemented for vdw-type = %s",
2333 evdw_names[fr->vdwtype]);
2336 /* When we deprecate the group kernels the code below can go too */
2337 if (fr->vdwtype == evdwPME && fr->cutoff_scheme == ecutsGROUP)
2339 /* Calculate self-interaction coefficient (assuming that
2340 * the reciprocal-space contribution is constant in the
2341 * region that contributes to the self-interaction).
2343 fr->enershiftsix = pow(fr->ewaldcoeff_lj, 6) / 6.0;
2345 eners[0] += -pow(sqrt(M_PI)*fr->ewaldcoeff_lj, 3)/3.0;
2346 virs[0] += pow(sqrt(M_PI)*fr->ewaldcoeff_lj, 3);
2349 fr->enerdiffsix = eners[0];
2350 fr->enerdifftwelve = eners[1];
2351 /* The 0.5 is due to the Gromacs definition of the virial */
2352 fr->virdiffsix = 0.5*virs[0];
2353 fr->virdifftwelve = 0.5*virs[1];
2357 void calc_dispcorr(t_inputrec *ir, t_forcerec *fr,
2359 matrix box, real lambda, tensor pres, tensor virial,
2360 real *prescorr, real *enercorr, real *dvdlcorr)
2362 gmx_bool bCorrAll, bCorrPres;
2363 real dvdlambda, invvol, dens, ninter, avcsix, avctwelve, enerdiff, svir = 0, spres = 0;
2373 if (ir->eDispCorr != edispcNO)
2375 bCorrAll = (ir->eDispCorr == edispcAllEner ||
2376 ir->eDispCorr == edispcAllEnerPres);
2377 bCorrPres = (ir->eDispCorr == edispcEnerPres ||
2378 ir->eDispCorr == edispcAllEnerPres);
2380 invvol = 1/det(box);
2383 /* Only correct for the interactions with the inserted molecule */
2384 dens = (natoms - fr->n_tpi)*invvol;
2389 dens = natoms*invvol;
2390 ninter = 0.5*natoms;
2393 if (ir->efep == efepNO)
2395 avcsix = fr->avcsix[0];
2396 avctwelve = fr->avctwelve[0];
2400 avcsix = (1 - lambda)*fr->avcsix[0] + lambda*fr->avcsix[1];
2401 avctwelve = (1 - lambda)*fr->avctwelve[0] + lambda*fr->avctwelve[1];
2404 enerdiff = ninter*(dens*fr->enerdiffsix - fr->enershiftsix);
2405 *enercorr += avcsix*enerdiff;
2407 if (ir->efep != efepNO)
2409 dvdlambda += (fr->avcsix[1] - fr->avcsix[0])*enerdiff;
2413 enerdiff = ninter*(dens*fr->enerdifftwelve - fr->enershifttwelve);
2414 *enercorr += avctwelve*enerdiff;
2415 if (fr->efep != efepNO)
2417 dvdlambda += (fr->avctwelve[1] - fr->avctwelve[0])*enerdiff;
2423 svir = ninter*dens*avcsix*fr->virdiffsix/3.0;
2424 if (ir->eDispCorr == edispcAllEnerPres)
2426 svir += ninter*dens*avctwelve*fr->virdifftwelve/3.0;
2428 /* The factor 2 is because of the Gromacs virial definition */
2429 spres = -2.0*invvol*svir*PRESFAC;
2431 for (m = 0; m < DIM; m++)
2433 virial[m][m] += svir;
2434 pres[m][m] += spres;
2439 /* Can't currently control when it prints, for now, just print when degugging */
2444 fprintf(debug, "Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
2450 "Long Range LJ corr.: Epot %10g, Pres: %10g, Vir: %10g\n",
2451 *enercorr, spres, svir);
2455 fprintf(debug, "Long Range LJ corr.: Epot %10g\n", *enercorr);
2459 if (fr->efep != efepNO)
2461 *dvdlcorr += dvdlambda;
2466 void do_pbc_first(FILE *fplog, matrix box, t_forcerec *fr,
2467 t_graph *graph, rvec x[])
2471 fprintf(fplog, "Removing pbc first time\n");
2473 calc_shifts(box, fr->shift_vec);
2476 mk_mshift(fplog, graph, fr->ePBC, box, x);
2479 p_graph(debug, "do_pbc_first 1", graph);
2481 shift_self(graph, box, x);
2482 /* By doing an extra mk_mshift the molecules that are broken
2483 * because they were e.g. imported from another software
2484 * will be made whole again. Such are the healing powers
2487 mk_mshift(fplog, graph, fr->ePBC, box, x);
2490 p_graph(debug, "do_pbc_first 2", graph);
2495 fprintf(fplog, "Done rmpbc\n");
2499 static void low_do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2500 gmx_mtop_t *mtop, rvec x[],
2505 gmx_molblock_t *molb;
2507 if (bFirst && fplog)
2509 fprintf(fplog, "Removing pbc first time\n");
2514 for (mb = 0; mb < mtop->nmolblock; mb++)
2516 molb = &mtop->molblock[mb];
2517 if (molb->natoms_mol == 1 ||
2518 (!bFirst && mtop->moltype[molb->type].cgs.nr == 1))
2520 /* Just one atom or charge group in the molecule, no PBC required */
2521 as += molb->nmol*molb->natoms_mol;
2525 /* Pass NULL iso fplog to avoid graph prints for each molecule type */
2526 mk_graph_ilist(NULL, mtop->moltype[molb->type].ilist,
2527 0, molb->natoms_mol, FALSE, FALSE, graph);
2529 for (mol = 0; mol < molb->nmol; mol++)
2531 mk_mshift(fplog, graph, ePBC, box, x+as);
2533 shift_self(graph, box, x+as);
2534 /* The molecule is whole now.
2535 * We don't need the second mk_mshift call as in do_pbc_first,
2536 * since we no longer need this graph.
2539 as += molb->natoms_mol;
2547 void do_pbc_first_mtop(FILE *fplog, int ePBC, matrix box,
2548 gmx_mtop_t *mtop, rvec x[])
2550 low_do_pbc_mtop(fplog, ePBC, box, mtop, x, TRUE);
2553 void do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2554 gmx_mtop_t *mtop, rvec x[])
2556 low_do_pbc_mtop(fplog, ePBC, box, mtop, x, FALSE);
2559 void finish_run(FILE *fplog, t_commrec *cr,
2560 t_inputrec *inputrec,
2561 t_nrnb nrnb[], gmx_wallcycle_t wcycle,
2562 gmx_walltime_accounting_t walltime_accounting,
2563 nonbonded_verlet_t *nbv,
2564 gmx_bool bWriteStat)
2566 t_nrnb *nrnb_tot = NULL;
2568 double nbfs = 0, mflop = 0;
2569 double elapsed_time,
2570 elapsed_time_over_all_ranks,
2571 elapsed_time_over_all_threads,
2572 elapsed_time_over_all_threads_over_all_ranks;
2573 wallcycle_sum(cr, wcycle);
2579 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM,
2580 cr->mpi_comm_mysim);
2588 elapsed_time = walltime_accounting_get_elapsed_time(walltime_accounting);
2589 elapsed_time_over_all_ranks = elapsed_time;
2590 elapsed_time_over_all_threads = walltime_accounting_get_elapsed_time_over_all_threads(walltime_accounting);
2591 elapsed_time_over_all_threads_over_all_ranks = elapsed_time_over_all_threads;
2595 /* reduce elapsed_time over all MPI ranks in the current simulation */
2596 MPI_Allreduce(&elapsed_time,
2597 &elapsed_time_over_all_ranks,
2598 1, MPI_DOUBLE, MPI_SUM,
2599 cr->mpi_comm_mysim);
2600 elapsed_time_over_all_ranks /= cr->nnodes;
2601 /* Reduce elapsed_time_over_all_threads over all MPI ranks in the
2602 * current simulation. */
2603 MPI_Allreduce(&elapsed_time_over_all_threads,
2604 &elapsed_time_over_all_threads_over_all_ranks,
2605 1, MPI_DOUBLE, MPI_SUM,
2606 cr->mpi_comm_mysim);
2612 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
2619 if ((cr->duty & DUTY_PP) && DOMAINDECOMP(cr))
2621 print_dd_statistics(cr, inputrec, fplog);
2626 wallclock_gpu_t* gputimes = use_GPU(nbv) ?
2627 nbnxn_cuda_get_timings(nbv->cu_nbv) : NULL;
2628 wallcycle_print(fplog, cr->nnodes, cr->npmenodes,
2629 elapsed_time_over_all_ranks,
2632 if (EI_DYNAMICS(inputrec->eI))
2634 delta_t = inputrec->delta_t;
2639 print_perf(fplog, elapsed_time_over_all_threads_over_all_ranks,
2640 elapsed_time_over_all_ranks,
2641 walltime_accounting_get_nsteps_done(walltime_accounting),
2642 delta_t, nbfs, mflop);
2646 print_perf(stderr, elapsed_time_over_all_threads_over_all_ranks,
2647 elapsed_time_over_all_ranks,
2648 walltime_accounting_get_nsteps_done(walltime_accounting),
2649 delta_t, nbfs, mflop);
2654 extern void initialize_lambdas(FILE *fplog, t_inputrec *ir, int *fep_state, real *lambda, double *lam0)
2656 /* this function works, but could probably use a logic rewrite to keep all the different
2657 types of efep straight. */
2660 t_lambda *fep = ir->fepvals;
2662 if ((ir->efep == efepNO) && (ir->bSimTemp == FALSE))
2664 for (i = 0; i < efptNR; i++)
2676 *fep_state = fep->init_fep_state; /* this might overwrite the checkpoint
2677 if checkpoint is set -- a kludge is in for now
2679 for (i = 0; i < efptNR; i++)
2681 /* overwrite lambda state with init_lambda for now for backwards compatibility */
2682 if (fep->init_lambda >= 0) /* if it's -1, it was never initializd */
2684 lambda[i] = fep->init_lambda;
2687 lam0[i] = lambda[i];
2692 lambda[i] = fep->all_lambda[i][*fep_state];
2695 lam0[i] = lambda[i];
2701 /* need to rescale control temperatures to match current state */
2702 for (i = 0; i < ir->opts.ngtc; i++)
2704 if (ir->opts.ref_t[i] > 0)
2706 ir->opts.ref_t[i] = ir->simtempvals->temperatures[*fep_state];
2712 /* Send to the log the information on the current lambdas */
2715 fprintf(fplog, "Initial vector of lambda components:[ ");
2716 for (i = 0; i < efptNR; i++)
2718 fprintf(fplog, "%10.4f ", lambda[i]);
2720 fprintf(fplog, "]\n");
2726 void init_md(FILE *fplog,
2727 t_commrec *cr, t_inputrec *ir, const output_env_t oenv,
2728 double *t, double *t0,
2729 real *lambda, int *fep_state, double *lam0,
2730 t_nrnb *nrnb, gmx_mtop_t *mtop,
2732 int nfile, const t_filenm fnm[],
2733 gmx_mdoutf_t *outf, t_mdebin **mdebin,
2734 tensor force_vir, tensor shake_vir, rvec mu_tot,
2735 gmx_bool *bSimAnn, t_vcm **vcm, unsigned long Flags,
2736 gmx_wallcycle_t wcycle)
2740 /* Initial values */
2741 *t = *t0 = ir->init_t;
2744 for (i = 0; i < ir->opts.ngtc; i++)
2746 /* set bSimAnn if any group is being annealed */
2747 if (ir->opts.annealing[i] != eannNO)
2754 update_annealing_target_temp(&(ir->opts), ir->init_t);
2757 /* Initialize lambda variables */
2758 initialize_lambdas(fplog, ir, fep_state, lambda, lam0);
2762 *upd = init_update(ir);
2768 *vcm = init_vcm(fplog, &mtop->groups, ir);
2771 if (EI_DYNAMICS(ir->eI) && !(Flags & MD_APPENDFILES))
2773 if (ir->etc == etcBERENDSEN)
2775 please_cite(fplog, "Berendsen84a");
2777 if (ir->etc == etcVRESCALE)
2779 please_cite(fplog, "Bussi2007a");
2781 if (ir->eI == eiSD1)
2783 please_cite(fplog, "Goga2012");
2786 if ((ir->et[XX].n > 0) || (ir->et[YY].n > 0) || (ir->et[ZZ].n > 0))
2788 please_cite(fplog, "Caleman2008a");
2794 *outf = init_mdoutf(fplog, nfile, fnm, Flags, cr, ir, mtop, oenv, wcycle);
2796 *mdebin = init_mdebin((Flags & MD_APPENDFILES) ? NULL : mdoutf_get_fp_ene(*outf),
2797 mtop, ir, mdoutf_get_fp_dhdl(*outf));
2802 please_cite(fplog, "Fritsch12");
2803 please_cite(fplog, "Junghans10");
2805 /* Initiate variables */
2806 clear_mat(force_vir);
2807 clear_mat(shake_vir);