1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * GROwing Monsters And Cloning Shrimps
41 #include <catamount/dclock.h>
47 #ifdef HAVE_SYS_TIME_H
60 #include "chargegroup.h"
83 #include "pull_rotation.h"
84 #include "gmx_random.h"
87 #include "gmx_wallcycle.h"
89 #include "nbnxn_atomdata.h"
90 #include "nbnxn_search.h"
91 #include "nbnxn_kernels/nbnxn_kernel_ref.h"
92 #include "nbnxn_kernels/nbnxn_kernel_simd_4xn.h"
93 #include "nbnxn_kernels/nbnxn_kernel_simd_2xnn.h"
94 #include "nbnxn_kernels/nbnxn_kernel_gpu_ref.h"
96 #include "gromacs/utility/gmxmpi.h"
101 #include "nbnxn_cuda_data_mgmt.h"
102 #include "nbnxn_cuda/nbnxn_cuda.h"
107 #ifdef HAVE_GETTIMEOFDAY
111 gettimeofday(&t, NULL);
113 seconds = (double) t.tv_sec + 1e-6*(double)t.tv_usec;
119 seconds = time(NULL);
126 #define difftime(end, start) ((double)(end)-(double)(start))
128 void print_time(FILE *out, gmx_runtime_t *runtime, gmx_large_int_t step,
129 t_inputrec *ir, t_commrec *cr)
132 char timebuf[STRLEN];
136 #ifndef GMX_THREAD_MPI
142 fprintf(out, "step %s", gmx_step_str(step, buf));
143 if ((step >= ir->nstlist))
145 runtime->last = gmx_gettime();
146 dt = difftime(runtime->last, runtime->real);
147 runtime->time_per_step = dt/(step - ir->init_step + 1);
149 dt = (ir->nsteps + ir->init_step - step)*runtime->time_per_step;
155 finish = (time_t) (runtime->last + dt);
156 gmx_ctime_r(&finish, timebuf, STRLEN);
157 sprintf(buf, "%s", timebuf);
158 buf[strlen(buf)-1] = '\0';
159 fprintf(out, ", will finish %s", buf);
163 fprintf(out, ", remaining runtime: %5d s ", (int)dt);
168 fprintf(out, " performance: %.1f ns/day ",
169 ir->delta_t/1000*24*60*60/runtime->time_per_step);
172 #ifndef GMX_THREAD_MPI
186 static double set_proctime(gmx_runtime_t *runtime)
192 prev = runtime->proc;
193 runtime->proc = dclock();
195 diff = runtime->proc - prev;
199 prev = runtime->proc;
200 runtime->proc = clock();
202 diff = (double)(runtime->proc - prev)/(double)CLOCKS_PER_SEC;
206 /* The counter has probably looped, ignore this data */
213 void runtime_start(gmx_runtime_t *runtime)
215 runtime->real = gmx_gettime();
217 set_proctime(runtime);
218 runtime->realtime = 0;
219 runtime->proctime = 0;
221 runtime->time_per_step = 0;
224 void runtime_end(gmx_runtime_t *runtime)
230 runtime->proctime += set_proctime(runtime);
231 runtime->realtime = now - runtime->real;
235 void runtime_upd_proc(gmx_runtime_t *runtime)
237 runtime->proctime += set_proctime(runtime);
240 void print_date_and_time(FILE *fplog, int nodeid, const char *title,
241 const gmx_runtime_t *runtime)
244 char timebuf[STRLEN];
245 char time_string[STRLEN];
252 tmptime = (time_t) runtime->real;
253 gmx_ctime_r(&tmptime, timebuf, STRLEN);
257 tmptime = (time_t) gmx_gettime();
258 gmx_ctime_r(&tmptime, timebuf, STRLEN);
260 for (i = 0; timebuf[i] >= ' '; i++)
262 time_string[i] = timebuf[i];
264 time_string[i] = '\0';
266 fprintf(fplog, "%s on node %d %s\n", title, nodeid, time_string);
270 static void sum_forces(int start, int end, rvec f[], rvec flr[])
276 pr_rvecs(debug, 0, "fsr", f+start, end-start);
277 pr_rvecs(debug, 0, "flr", flr+start, end-start);
279 for (i = start; (i < end); i++)
281 rvec_inc(f[i], flr[i]);
286 * calc_f_el calculates forces due to an electric field.
288 * force is kJ mol^-1 nm^-1 = e * kJ mol^-1 nm^-1 / e
290 * Et[] contains the parameters for the time dependent
291 * part of the field (not yet used).
292 * Ex[] contains the parameters for
293 * the spatial dependent part of the field. You can have cool periodic
294 * fields in principle, but only a constant field is supported
296 * The function should return the energy due to the electric field
297 * (if any) but for now returns 0.
300 * There can be problems with the virial.
301 * Since the field is not self-consistent this is unavoidable.
302 * For neutral molecules the virial is correct within this approximation.
303 * For neutral systems with many charged molecules the error is small.
304 * But for systems with a net charge or a few charged molecules
305 * the error can be significant when the field is high.
306 * Solution: implement a self-consitent electric field into PME.
308 static void calc_f_el(FILE *fp, int start, int homenr,
309 real charge[], rvec x[], rvec f[],
310 t_cosines Ex[], t_cosines Et[], double t)
316 for (m = 0; (m < DIM); m++)
323 Ext[m] = cos(Et[m].a[0]*(t-t0))*exp(-sqr(t-t0)/(2.0*sqr(Et[m].a[2])));
327 Ext[m] = cos(Et[m].a[0]*t);
336 /* Convert the field strength from V/nm to MD-units */
337 Ext[m] *= Ex[m].a[0]*FIELDFAC;
338 for (i = start; (i < start+homenr); i++)
340 f[i][m] += charge[i]*Ext[m];
350 fprintf(fp, "%10g %10g %10g %10g #FIELD\n", t,
351 Ext[XX]/FIELDFAC, Ext[YY]/FIELDFAC, Ext[ZZ]/FIELDFAC);
355 static void calc_virial(FILE *fplog, int start, int homenr, rvec x[], rvec f[],
356 tensor vir_part, t_graph *graph, matrix box,
357 t_nrnb *nrnb, const t_forcerec *fr, int ePBC)
362 /* The short-range virial from surrounding boxes */
364 calc_vir(fplog, SHIFTS, fr->shift_vec, fr->fshift, vir_part, ePBC == epbcSCREW, box);
365 inc_nrnb(nrnb, eNR_VIRIAL, SHIFTS);
367 /* Calculate partial virial, for local atoms only, based on short range.
368 * Total virial is computed in global_stat, called from do_md
370 f_calc_vir(fplog, start, start+homenr, x, f, vir_part, graph, box);
371 inc_nrnb(nrnb, eNR_VIRIAL, homenr);
373 /* Add position restraint contribution */
374 for (i = 0; i < DIM; i++)
376 vir_part[i][i] += fr->vir_diag_posres[i];
379 /* Add wall contribution */
380 for (i = 0; i < DIM; i++)
382 vir_part[i][ZZ] += fr->vir_wall_z[i];
387 pr_rvecs(debug, 0, "vir_part", vir_part, DIM);
391 static void posres_wrapper(FILE *fplog,
397 matrix box, rvec x[],
399 gmx_enerdata_t *enerd,
407 /* Position restraints always require full pbc */
408 set_pbc(&pbc, ir->ePBC, box);
410 v = posres(top->idef.il[F_POSRES].nr, top->idef.il[F_POSRES].iatoms,
411 top->idef.iparams_posres,
412 (const rvec*)x, fr->f_novirsum, fr->vir_diag_posres,
413 ir->ePBC == epbcNONE ? NULL : &pbc,
414 lambda[efptRESTRAINT], &dvdl,
415 fr->rc_scaling, fr->ePBC, fr->posres_com, fr->posres_comB);
418 gmx_print_sepdvdl(fplog, interaction_function[F_POSRES].longname, v, dvdl);
420 enerd->term[F_POSRES] += v;
421 /* If just the force constant changes, the FEP term is linear,
422 * but if k changes, it is not.
424 enerd->dvdl_nonlin[efptRESTRAINT] += dvdl;
425 inc_nrnb(nrnb, eNR_POSRES, top->idef.il[F_POSRES].nr/2);
427 if ((ir->fepvals->n_lambda > 0) && (flags & GMX_FORCE_DHDL))
429 for (i = 0; i < enerd->n_lambda; i++)
431 real dvdl_dum, lambda_dum;
433 lambda_dum = (i == 0 ? lambda[efptRESTRAINT] : ir->fepvals->all_lambda[efptRESTRAINT][i-1]);
434 v = posres(top->idef.il[F_POSRES].nr, top->idef.il[F_POSRES].iatoms,
435 top->idef.iparams_posres,
436 (const rvec*)x, NULL, NULL,
437 ir->ePBC == epbcNONE ? NULL : &pbc, lambda_dum, &dvdl,
438 fr->rc_scaling, fr->ePBC, fr->posres_com, fr->posres_comB);
439 enerd->enerpart_lambda[i] += v;
444 static void pull_potential_wrapper(FILE *fplog,
448 matrix box, rvec x[],
452 gmx_enerdata_t *enerd,
459 /* Calculate the center of mass forces, this requires communication,
460 * which is why pull_potential is called close to other communication.
461 * The virial contribution is calculated directly,
462 * which is why we call pull_potential after calc_virial.
464 set_pbc(&pbc, ir->ePBC, box);
466 enerd->term[F_COM_PULL] +=
467 pull_potential(ir->ePull, ir->pull, mdatoms, &pbc,
468 cr, t, lambda[efptRESTRAINT], x, f, vir_force, &dvdl);
471 gmx_print_sepdvdl(fplog, "Com pull", enerd->term[F_COM_PULL], dvdl);
473 enerd->dvdl_lin[efptRESTRAINT] += dvdl;
476 static void pme_receive_force_ener(FILE *fplog,
479 gmx_wallcycle_t wcycle,
480 gmx_enerdata_t *enerd,
484 float cycles_ppdpme, cycles_seppme;
486 cycles_ppdpme = wallcycle_stop(wcycle, ewcPPDURINGPME);
487 dd_cycles_add(cr->dd, cycles_ppdpme, ddCyclPPduringPME);
489 /* In case of node-splitting, the PP nodes receive the long-range
490 * forces, virial and energy from the PME nodes here.
492 wallcycle_start(wcycle, ewcPP_PMEWAITRECVF);
494 gmx_pme_receive_f(cr, fr->f_novirsum, fr->vir_el_recip, &e, &dvdl,
498 gmx_print_sepdvdl(fplog, "PME mesh", e, dvdl);
500 enerd->term[F_COUL_RECIP] += e;
501 enerd->dvdl_lin[efptCOUL] += dvdl;
504 dd_cycles_add(cr->dd, cycles_seppme, ddCyclPME);
506 wallcycle_stop(wcycle, ewcPP_PMEWAITRECVF);
509 static void print_large_forces(FILE *fp, t_mdatoms *md, t_commrec *cr,
510 gmx_large_int_t step, real pforce, rvec *x, rvec *f)
514 char buf[STEPSTRSIZE];
517 for (i = md->start; i < md->start+md->homenr; i++)
520 /* We also catch NAN, if the compiler does not optimize this away. */
521 if (fn2 >= pf2 || fn2 != fn2)
523 fprintf(fp, "step %s atom %6d x %8.3f %8.3f %8.3f force %12.5e\n",
524 gmx_step_str(step, buf),
525 ddglatnr(cr->dd, i), x[i][XX], x[i][YY], x[i][ZZ], sqrt(fn2));
530 static void post_process_forces(FILE *fplog,
532 gmx_large_int_t step,
533 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
535 matrix box, rvec x[],
540 t_forcerec *fr, gmx_vsite_t *vsite,
547 /* Spread the mesh force on virtual sites to the other particles...
548 * This is parallellized. MPI communication is performed
549 * if the constructing atoms aren't local.
551 wallcycle_start(wcycle, ewcVSITESPREAD);
552 spread_vsite_f(fplog, vsite, x, fr->f_novirsum, NULL,
553 (flags & GMX_FORCE_VIRIAL), fr->vir_el_recip,
555 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
556 wallcycle_stop(wcycle, ewcVSITESPREAD);
558 if (flags & GMX_FORCE_VIRIAL)
560 /* Now add the forces, this is local */
563 sum_forces(0, fr->f_novirsum_n, f, fr->f_novirsum);
567 sum_forces(mdatoms->start, mdatoms->start+mdatoms->homenr,
570 if (EEL_FULL(fr->eeltype))
572 /* Add the mesh contribution to the virial */
573 m_add(vir_force, fr->vir_el_recip, vir_force);
577 pr_rvecs(debug, 0, "vir_force", vir_force, DIM);
582 if (fr->print_force >= 0)
584 print_large_forces(stderr, mdatoms, cr, step, fr->print_force, x, f);
588 static void do_nb_verlet(t_forcerec *fr,
589 interaction_const_t *ic,
590 gmx_enerdata_t *enerd,
591 int flags, int ilocality,
594 gmx_wallcycle_t wcycle)
596 int nnbl, kernel_type, enr_nbnxn_kernel_ljc, enr_nbnxn_kernel_lj;
598 nonbonded_verlet_group_t *nbvg;
601 if (!(flags & GMX_FORCE_NONBONDED))
603 /* skip non-bonded calculation */
607 nbvg = &fr->nbv->grp[ilocality];
609 /* CUDA kernel launch overhead is already timed separately */
610 if (fr->cutoff_scheme != ecutsVERLET)
612 gmx_incons("Invalid cut-off scheme passed!");
615 bCUDA = (nbvg->kernel_type == nbnxnk8x8x8_CUDA);
619 wallcycle_sub_start(wcycle, ewcsNONBONDED);
621 switch (nbvg->kernel_type)
623 case nbnxnk4x4_PlainC:
624 nbnxn_kernel_ref(&nbvg->nbl_lists,
630 enerd->grpp.ener[egCOULSR],
632 enerd->grpp.ener[egBHAMSR] :
633 enerd->grpp.ener[egLJSR]);
636 case nbnxnk4xN_SIMD_4xN:
637 nbnxn_kernel_simd_4xn(&nbvg->nbl_lists,
644 enerd->grpp.ener[egCOULSR],
646 enerd->grpp.ener[egBHAMSR] :
647 enerd->grpp.ener[egLJSR]);
649 case nbnxnk4xN_SIMD_2xNN:
650 nbnxn_kernel_simd_2xnn(&nbvg->nbl_lists,
657 enerd->grpp.ener[egCOULSR],
659 enerd->grpp.ener[egBHAMSR] :
660 enerd->grpp.ener[egLJSR]);
663 case nbnxnk8x8x8_CUDA:
664 nbnxn_cuda_launch_kernel(fr->nbv->cu_nbv, nbvg->nbat, flags, ilocality);
667 case nbnxnk8x8x8_PlainC:
668 nbnxn_kernel_gpu_ref(nbvg->nbl_lists.nbl[0],
673 nbvg->nbat->out[0].f,
675 enerd->grpp.ener[egCOULSR],
677 enerd->grpp.ener[egBHAMSR] :
678 enerd->grpp.ener[egLJSR]);
682 gmx_incons("Invalid nonbonded kernel type passed!");
687 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
690 if (EEL_RF(ic->eeltype) || ic->eeltype == eelCUT)
692 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_RF;
694 else if ((!bCUDA && nbvg->ewald_excl == ewaldexclAnalytical) ||
695 (bCUDA && nbnxn_cuda_is_kernel_ewald_analytical(fr->nbv->cu_nbv)))
697 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_EWALD;
701 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_TAB;
703 enr_nbnxn_kernel_lj = eNR_NBNXN_LJ;
704 if (flags & GMX_FORCE_ENERGY)
706 /* In eNR_??? the nbnxn F+E kernels are always the F kernel + 1 */
707 enr_nbnxn_kernel_ljc += 1;
708 enr_nbnxn_kernel_lj += 1;
711 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc,
712 nbvg->nbl_lists.natpair_ljq);
713 inc_nrnb(nrnb, enr_nbnxn_kernel_lj,
714 nbvg->nbl_lists.natpair_lj);
715 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc-eNR_NBNXN_LJ_RF+eNR_NBNXN_RF,
716 nbvg->nbl_lists.natpair_q);
719 void do_force_cutsVERLET(FILE *fplog, t_commrec *cr,
720 t_inputrec *inputrec,
721 gmx_large_int_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
724 gmx_groups_t *groups,
725 matrix box, rvec x[], history_t *hist,
729 gmx_enerdata_t *enerd, t_fcdata *fcd,
730 real *lambda, t_graph *graph,
731 t_forcerec *fr, interaction_const_t *ic,
732 gmx_vsite_t *vsite, rvec mu_tot,
733 double t, FILE *field, gmx_edsam_t ed,
741 gmx_bool bSepDVDL, bStateChanged, bNS, bFillGrid, bCalcCGCM, bBS;
742 gmx_bool bDoLongRange, bDoForces, bSepLRF, bUseGPU, bUseOrEmulGPU;
743 gmx_bool bDiffKernels = FALSE;
745 rvec vzero, box_diag;
747 float cycles_pme, cycles_force;
748 nonbonded_verlet_t *nbv;
752 nb_kernel_type = fr->nbv->grp[0].kernel_type;
754 start = mdatoms->start;
755 homenr = mdatoms->homenr;
757 bSepDVDL = (fr->bSepDVDL && do_per_step(step, inputrec->nstlog));
759 clear_mat(vir_force);
762 if (DOMAINDECOMP(cr))
764 cg1 = cr->dd->ncg_tot;
775 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
776 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
777 bFillGrid = (bNS && bStateChanged);
778 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
779 bDoLongRange = (fr->bTwinRange && bNS && (flags & GMX_FORCE_DO_LR));
780 bDoForces = (flags & GMX_FORCE_FORCES);
781 bSepLRF = (bDoLongRange && bDoForces && (flags & GMX_FORCE_SEPLRF));
782 bUseGPU = fr->nbv->bUseGPU;
783 bUseOrEmulGPU = bUseGPU || (nbv->grp[0].kernel_type == nbnxnk8x8x8_PlainC);
787 update_forcerec(fplog, fr, box);
789 if (NEED_MUTOT(*inputrec))
791 /* Calculate total (local) dipole moment in a temporary common array.
792 * This makes it possible to sum them over nodes faster.
794 calc_mu(start, homenr,
795 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
800 if (fr->ePBC != epbcNONE)
802 /* Compute shift vectors every step,
803 * because of pressure coupling or box deformation!
805 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
807 calc_shifts(box, fr->shift_vec);
812 put_atoms_in_box_omp(fr->ePBC, box, homenr, x);
813 inc_nrnb(nrnb, eNR_SHIFTX, homenr);
815 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
817 unshift_self(graph, box, x);
821 nbnxn_atomdata_copy_shiftvec(flags & GMX_FORCE_DYNAMICBOX,
822 fr->shift_vec, nbv->grp[0].nbat);
825 if (!(cr->duty & DUTY_PME))
827 /* Send particle coordinates to the pme nodes.
828 * Since this is only implemented for domain decomposition
829 * and domain decomposition does not use the graph,
830 * we do not need to worry about shifting.
833 wallcycle_start(wcycle, ewcPP_PMESENDX);
835 bBS = (inputrec->nwall == 2);
839 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
842 gmx_pme_send_x(cr, bBS ? boxs : box, x,
843 mdatoms->nChargePerturbed, lambda[efptCOUL],
844 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)), step);
846 wallcycle_stop(wcycle, ewcPP_PMESENDX);
850 /* do gridding for pair search */
853 if (graph && bStateChanged)
855 /* Calculate intramolecular shift vectors to make molecules whole */
856 mk_mshift(fplog, graph, fr->ePBC, box, x);
860 box_diag[XX] = box[XX][XX];
861 box_diag[YY] = box[YY][YY];
862 box_diag[ZZ] = box[ZZ][ZZ];
864 wallcycle_start(wcycle, ewcNS);
867 wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
868 nbnxn_put_on_grid(nbv->nbs, fr->ePBC, box,
870 0, mdatoms->homenr, -1, fr->cginfo, x,
872 nbv->grp[eintLocal].kernel_type,
873 nbv->grp[eintLocal].nbat);
874 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
878 wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
879 nbnxn_put_on_grid_nonlocal(nbv->nbs, domdec_zones(cr->dd),
881 nbv->grp[eintNonlocal].kernel_type,
882 nbv->grp[eintNonlocal].nbat);
883 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
886 if (nbv->ngrp == 1 ||
887 nbv->grp[eintNonlocal].nbat == nbv->grp[eintLocal].nbat)
889 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatAll,
890 nbv->nbs, mdatoms, fr->cginfo);
894 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatLocal,
895 nbv->nbs, mdatoms, fr->cginfo);
896 nbnxn_atomdata_set(nbv->grp[eintNonlocal].nbat, eatAll,
897 nbv->nbs, mdatoms, fr->cginfo);
899 wallcycle_stop(wcycle, ewcNS);
902 /* initialize the GPU atom data and copy shift vector */
907 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
908 nbnxn_cuda_init_atomdata(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
909 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
912 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
913 nbnxn_cuda_upload_shiftvec(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
914 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
917 /* do local pair search */
920 wallcycle_start_nocount(wcycle, ewcNS);
921 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
922 nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintLocal].nbat,
925 nbv->min_ci_balanced,
926 &nbv->grp[eintLocal].nbl_lists,
928 nbv->grp[eintLocal].kernel_type,
930 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
934 /* initialize local pair-list on the GPU */
935 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
936 nbv->grp[eintLocal].nbl_lists.nbl[0],
939 wallcycle_stop(wcycle, ewcNS);
943 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
944 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
945 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, FALSE, x,
946 nbv->grp[eintLocal].nbat);
947 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
948 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
953 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
954 /* launch local nonbonded F on GPU */
955 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFNo,
957 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
960 /* Communicate coordinates and sum dipole if necessary +
961 do non-local pair search */
962 if (DOMAINDECOMP(cr))
964 bDiffKernels = (nbv->grp[eintNonlocal].kernel_type !=
965 nbv->grp[eintLocal].kernel_type);
969 /* With GPU+CPU non-bonded calculations we need to copy
970 * the local coordinates to the non-local nbat struct
971 * (in CPU format) as the non-local kernel call also
972 * calculates the local - non-local interactions.
974 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
975 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
976 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, TRUE, x,
977 nbv->grp[eintNonlocal].nbat);
978 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
979 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
984 wallcycle_start_nocount(wcycle, ewcNS);
985 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
989 nbnxn_grid_add_simple(nbv->nbs, nbv->grp[eintNonlocal].nbat);
992 nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintNonlocal].nbat,
995 nbv->min_ci_balanced,
996 &nbv->grp[eintNonlocal].nbl_lists,
998 nbv->grp[eintNonlocal].kernel_type,
1001 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1003 if (nbv->grp[eintNonlocal].kernel_type == nbnxnk8x8x8_CUDA)
1005 /* initialize non-local pair-list on the GPU */
1006 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
1007 nbv->grp[eintNonlocal].nbl_lists.nbl[0],
1010 wallcycle_stop(wcycle, ewcNS);
1014 wallcycle_start(wcycle, ewcMOVEX);
1015 dd_move_x(cr->dd, box, x);
1017 /* When we don't need the total dipole we sum it in global_stat */
1018 if (bStateChanged && NEED_MUTOT(*inputrec))
1020 gmx_sumd(2*DIM, mu, cr);
1022 wallcycle_stop(wcycle, ewcMOVEX);
1024 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1025 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
1026 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatNonlocal, FALSE, x,
1027 nbv->grp[eintNonlocal].nbat);
1028 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
1029 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1032 if (bUseGPU && !bDiffKernels)
1034 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
1035 /* launch non-local nonbonded F on GPU */
1036 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFNo,
1038 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1044 /* launch D2H copy-back F */
1045 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1046 if (DOMAINDECOMP(cr) && !bDiffKernels)
1048 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintNonlocal].nbat,
1049 flags, eatNonlocal);
1051 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintLocal].nbat,
1053 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1056 if (bStateChanged && NEED_MUTOT(*inputrec))
1060 gmx_sumd(2*DIM, mu, cr);
1063 for (i = 0; i < 2; i++)
1065 for (j = 0; j < DIM; j++)
1067 fr->mu_tot[i][j] = mu[i*DIM + j];
1071 if (fr->efep == efepNO)
1073 copy_rvec(fr->mu_tot[0], mu_tot);
1077 for (j = 0; j < DIM; j++)
1080 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] +
1081 lambda[efptCOUL]*fr->mu_tot[1][j];
1085 /* Reset energies */
1086 reset_enerdata(&(inputrec->opts), fr, bNS, enerd, MASTER(cr));
1087 clear_rvecs(SHIFTS, fr->fshift);
1089 if (DOMAINDECOMP(cr))
1091 if (!(cr->duty & DUTY_PME))
1093 wallcycle_start(wcycle, ewcPPDURINGPME);
1094 dd_force_flop_start(cr->dd, nrnb);
1100 /* Enforced rotation has its own cycle counter that starts after the collective
1101 * coordinates have been communicated. It is added to ddCyclF to allow
1102 * for proper load-balancing */
1103 wallcycle_start(wcycle, ewcROT);
1104 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1105 wallcycle_stop(wcycle, ewcROT);
1108 /* Start the force cycle counter.
1109 * This counter is stopped in do_forcelow_level.
1110 * No parallel communication should occur while this counter is running,
1111 * since that will interfere with the dynamic load balancing.
1113 wallcycle_start(wcycle, ewcFORCE);
1116 /* Reset forces for which the virial is calculated separately:
1117 * PME/Ewald forces if necessary */
1118 if (fr->bF_NoVirSum)
1120 if (flags & GMX_FORCE_VIRIAL)
1122 fr->f_novirsum = fr->f_novirsum_alloc;
1125 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1129 clear_rvecs(homenr, fr->f_novirsum+start);
1134 /* We are not calculating the pressure so we do not need
1135 * a separate array for forces that do not contribute
1142 /* Clear the short- and long-range forces */
1143 clear_rvecs(fr->natoms_force_constr, f);
1144 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1146 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1149 clear_rvec(fr->vir_diag_posres);
1152 if (inputrec->ePull == epullCONSTRAINT)
1154 clear_pull_forces(inputrec->pull);
1157 /* We calculate the non-bonded forces, when done on the CPU, here.
1158 * We do this before calling do_force_lowlevel, as in there bondeds
1159 * forces are calculated before PME, which does communication.
1160 * With this order, non-bonded and bonded force calculation imbalance
1161 * can be balanced out by the domain decomposition load balancing.
1166 /* Maybe we should move this into do_force_lowlevel */
1167 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFYes,
1171 if (!bUseOrEmulGPU || bDiffKernels)
1175 if (DOMAINDECOMP(cr))
1177 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal,
1178 bDiffKernels ? enbvClearFYes : enbvClearFNo,
1188 aloc = eintNonlocal;
1191 /* Add all the non-bonded force to the normal force array.
1192 * This can be split into a local a non-local part when overlapping
1193 * communication with calculation with domain decomposition.
1195 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1196 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1197 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1198 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatAll, nbv->grp[aloc].nbat, f);
1199 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1200 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1201 wallcycle_start_nocount(wcycle, ewcFORCE);
1203 /* if there are multiple fshift output buffers reduce them */
1204 if ((flags & GMX_FORCE_VIRIAL) &&
1205 nbv->grp[aloc].nbl_lists.nnbl > 1)
1207 nbnxn_atomdata_add_nbat_fshift_to_fshift(nbv->grp[aloc].nbat,
1212 /* update QMMMrec, if necessary */
1215 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1218 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1220 posres_wrapper(fplog, flags, bSepDVDL, inputrec, nrnb, top, box, x,
1221 f, enerd, lambda, fr);
1224 /* Compute the bonded and non-bonded energies and optionally forces */
1225 do_force_lowlevel(fplog, step, fr, inputrec, &(top->idef),
1226 cr, nrnb, wcycle, mdatoms, &(inputrec->opts),
1227 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, mtop, top, fr->born,
1228 &(top->atomtypes), bBornRadii, box,
1229 inputrec->fepvals, lambda, graph, &(top->excls), fr->mu_tot,
1230 flags, &cycles_pme);
1234 if (do_per_step(step, inputrec->nstcalclr))
1236 /* Add the long range forces to the short range forces */
1237 for (i = 0; i < fr->natoms_force_constr; i++)
1239 rvec_add(fr->f_twin[i], f[i], f[i]);
1244 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1248 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1251 if (bUseOrEmulGPU && !bDiffKernels)
1253 /* wait for non-local forces (or calculate in emulation mode) */
1254 if (DOMAINDECOMP(cr))
1258 wallcycle_start(wcycle, ewcWAIT_GPU_NB_NL);
1259 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1260 nbv->grp[eintNonlocal].nbat,
1262 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1264 cycles_force += wallcycle_stop(wcycle, ewcWAIT_GPU_NB_NL);
1268 wallcycle_start_nocount(wcycle, ewcFORCE);
1269 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFYes,
1271 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1273 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1274 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1275 /* skip the reduction if there was no non-local work to do */
1276 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1278 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatNonlocal,
1279 nbv->grp[eintNonlocal].nbat, f);
1281 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1282 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1288 /* Communicate the forces */
1291 wallcycle_start(wcycle, ewcMOVEF);
1292 if (DOMAINDECOMP(cr))
1294 dd_move_f(cr->dd, f, fr->fshift);
1295 /* Do we need to communicate the separate force array
1296 * for terms that do not contribute to the single sum virial?
1297 * Position restraints and electric fields do not introduce
1298 * inter-cg forces, only full electrostatics methods do.
1299 * When we do not calculate the virial, fr->f_novirsum = f,
1300 * so we have already communicated these forces.
1302 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1303 (flags & GMX_FORCE_VIRIAL))
1305 dd_move_f(cr->dd, fr->f_novirsum, NULL);
1309 /* We should not update the shift forces here,
1310 * since f_twin is already included in f.
1312 dd_move_f(cr->dd, fr->f_twin, NULL);
1315 wallcycle_stop(wcycle, ewcMOVEF);
1321 /* wait for local forces (or calculate in emulation mode) */
1324 wallcycle_start(wcycle, ewcWAIT_GPU_NB_L);
1325 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1326 nbv->grp[eintLocal].nbat,
1328 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1330 wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
1332 /* now clear the GPU outputs while we finish the step on the CPU */
1334 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1335 nbnxn_cuda_clear_outputs(nbv->cu_nbv, flags);
1336 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1340 wallcycle_start_nocount(wcycle, ewcFORCE);
1341 do_nb_verlet(fr, ic, enerd, flags, eintLocal,
1342 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
1344 wallcycle_stop(wcycle, ewcFORCE);
1346 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1347 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1348 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1350 /* skip the reduction if there was no non-local work to do */
1351 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatLocal,
1352 nbv->grp[eintLocal].nbat, f);
1354 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1355 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1358 if (DOMAINDECOMP(cr))
1360 dd_force_flop_stop(cr->dd, nrnb);
1363 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1369 if (IR_ELEC_FIELD(*inputrec))
1371 /* Compute forces due to electric field */
1372 calc_f_el(MASTER(cr) ? field : NULL,
1373 start, homenr, mdatoms->chargeA, x, fr->f_novirsum,
1374 inputrec->ex, inputrec->et, t);
1377 /* If we have NoVirSum forces, but we do not calculate the virial,
1378 * we sum fr->f_novirum=f later.
1380 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1382 wallcycle_start(wcycle, ewcVSITESPREAD);
1383 spread_vsite_f(fplog, vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1384 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1385 wallcycle_stop(wcycle, ewcVSITESPREAD);
1389 wallcycle_start(wcycle, ewcVSITESPREAD);
1390 spread_vsite_f(fplog, vsite, x, fr->f_twin, NULL, FALSE, NULL,
1392 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1393 wallcycle_stop(wcycle, ewcVSITESPREAD);
1397 if (flags & GMX_FORCE_VIRIAL)
1399 /* Calculation of the virial must be done after vsites! */
1400 calc_virial(fplog, mdatoms->start, mdatoms->homenr, x, f,
1401 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1405 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1407 pull_potential_wrapper(fplog, bSepDVDL, cr, inputrec, box, x,
1408 f, vir_force, mdatoms, enerd, lambda, t);
1411 /* Add the forces from enforced rotation potentials (if any) */
1414 wallcycle_start(wcycle, ewcROTadd);
1415 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
1416 wallcycle_stop(wcycle, ewcROTadd);
1419 if (PAR(cr) && !(cr->duty & DUTY_PME))
1421 /* In case of node-splitting, the PP nodes receive the long-range
1422 * forces, virial and energy from the PME nodes here.
1424 pme_receive_force_ener(fplog, bSepDVDL, cr, wcycle, enerd, fr);
1429 post_process_forces(fplog, cr, step, nrnb, wcycle,
1430 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
1434 /* Sum the potential energy terms from group contributions */
1435 sum_epot(&(inputrec->opts), &(enerd->grpp), enerd->term);
1438 void do_force_cutsGROUP(FILE *fplog, t_commrec *cr,
1439 t_inputrec *inputrec,
1440 gmx_large_int_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1441 gmx_localtop_t *top,
1443 gmx_groups_t *groups,
1444 matrix box, rvec x[], history_t *hist,
1448 gmx_enerdata_t *enerd, t_fcdata *fcd,
1449 real *lambda, t_graph *graph,
1450 t_forcerec *fr, gmx_vsite_t *vsite, rvec mu_tot,
1451 double t, FILE *field, gmx_edsam_t ed,
1452 gmx_bool bBornRadii,
1458 gmx_bool bSepDVDL, bStateChanged, bNS, bFillGrid, bCalcCGCM, bBS;
1459 gmx_bool bDoLongRangeNS, bDoForces, bDoPotential, bSepLRF;
1460 gmx_bool bDoAdressWF;
1462 rvec vzero, box_diag;
1463 real e, v, dvdlambda[efptNR];
1465 float cycles_pme, cycles_force;
1467 start = mdatoms->start;
1468 homenr = mdatoms->homenr;
1470 bSepDVDL = (fr->bSepDVDL && do_per_step(step, inputrec->nstlog));
1472 clear_mat(vir_force);
1476 pd_cg_range(cr, &cg0, &cg1);
1481 if (DOMAINDECOMP(cr))
1483 cg1 = cr->dd->ncg_tot;
1495 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
1496 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
1497 /* Should we update the long-range neighborlists at this step? */
1498 bDoLongRangeNS = fr->bTwinRange && bNS;
1499 /* Should we perform the long-range nonbonded evaluation inside the neighborsearching? */
1500 bFillGrid = (bNS && bStateChanged);
1501 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
1502 bDoForces = (flags & GMX_FORCE_FORCES);
1503 bDoPotential = (flags & GMX_FORCE_ENERGY);
1504 bSepLRF = ((inputrec->nstcalclr > 1) && bDoForces &&
1505 (flags & GMX_FORCE_SEPLRF) && (flags & GMX_FORCE_DO_LR));
1507 /* should probably move this to the forcerec since it doesn't change */
1508 bDoAdressWF = ((fr->adress_type != eAdressOff));
1512 update_forcerec(fplog, fr, box);
1514 if (NEED_MUTOT(*inputrec))
1516 /* Calculate total (local) dipole moment in a temporary common array.
1517 * This makes it possible to sum them over nodes faster.
1519 calc_mu(start, homenr,
1520 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
1525 if (fr->ePBC != epbcNONE)
1527 /* Compute shift vectors every step,
1528 * because of pressure coupling or box deformation!
1530 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
1532 calc_shifts(box, fr->shift_vec);
1537 put_charge_groups_in_box(fplog, cg0, cg1, fr->ePBC, box,
1538 &(top->cgs), x, fr->cg_cm);
1539 inc_nrnb(nrnb, eNR_CGCM, homenr);
1540 inc_nrnb(nrnb, eNR_RESETX, cg1-cg0);
1542 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
1544 unshift_self(graph, box, x);
1549 calc_cgcm(fplog, cg0, cg1, &(top->cgs), x, fr->cg_cm);
1550 inc_nrnb(nrnb, eNR_CGCM, homenr);
1557 move_cgcm(fplog, cr, fr->cg_cm);
1561 pr_rvecs(debug, 0, "cgcm", fr->cg_cm, top->cgs.nr);
1566 if (!(cr->duty & DUTY_PME))
1568 /* Send particle coordinates to the pme nodes.
1569 * Since this is only implemented for domain decomposition
1570 * and domain decomposition does not use the graph,
1571 * we do not need to worry about shifting.
1574 wallcycle_start(wcycle, ewcPP_PMESENDX);
1576 bBS = (inputrec->nwall == 2);
1579 copy_mat(box, boxs);
1580 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
1583 gmx_pme_send_x(cr, bBS ? boxs : box, x,
1584 mdatoms->nChargePerturbed, lambda[efptCOUL],
1585 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)), step);
1587 wallcycle_stop(wcycle, ewcPP_PMESENDX);
1589 #endif /* GMX_MPI */
1591 /* Communicate coordinates and sum dipole if necessary */
1594 wallcycle_start(wcycle, ewcMOVEX);
1595 if (DOMAINDECOMP(cr))
1597 dd_move_x(cr->dd, box, x);
1601 move_x(fplog, cr, GMX_LEFT, GMX_RIGHT, x, nrnb);
1603 wallcycle_stop(wcycle, ewcMOVEX);
1606 /* update adress weight beforehand */
1607 if (bStateChanged && bDoAdressWF)
1609 /* need pbc for adress weight calculation with pbc_dx */
1610 set_pbc(&pbc, inputrec->ePBC, box);
1611 if (fr->adress_site == eAdressSITEcog)
1613 update_adress_weights_cog(top->idef.iparams, top->idef.il, x, fr, mdatoms,
1614 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1616 else if (fr->adress_site == eAdressSITEcom)
1618 update_adress_weights_com(fplog, cg0, cg1, &(top->cgs), x, fr, mdatoms,
1619 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1621 else if (fr->adress_site == eAdressSITEatomatom)
1623 update_adress_weights_atom_per_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1624 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1628 update_adress_weights_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1629 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1633 if (NEED_MUTOT(*inputrec))
1640 gmx_sumd(2*DIM, mu, cr);
1642 for (i = 0; i < 2; i++)
1644 for (j = 0; j < DIM; j++)
1646 fr->mu_tot[i][j] = mu[i*DIM + j];
1650 if (fr->efep == efepNO)
1652 copy_rvec(fr->mu_tot[0], mu_tot);
1656 for (j = 0; j < DIM; j++)
1659 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] + lambda[efptCOUL]*fr->mu_tot[1][j];
1664 /* Reset energies */
1665 reset_enerdata(&(inputrec->opts), fr, bNS, enerd, MASTER(cr));
1666 clear_rvecs(SHIFTS, fr->fshift);
1670 wallcycle_start(wcycle, ewcNS);
1672 if (graph && bStateChanged)
1674 /* Calculate intramolecular shift vectors to make molecules whole */
1675 mk_mshift(fplog, graph, fr->ePBC, box, x);
1678 /* Do the actual neighbour searching and if twin range electrostatics
1679 * also do the calculation of long range forces and energies.
1681 for (i = 0; i < efptNR; i++)
1685 ns(fplog, fr, x, box,
1686 groups, &(inputrec->opts), top, mdatoms,
1687 cr, nrnb, lambda, dvdlambda, &enerd->grpp, bFillGrid,
1690 wallcycle_stop(wcycle, ewcNS);
1693 if (inputrec->implicit_solvent && bNS)
1695 make_gb_nblist(cr, inputrec->gb_algorithm, inputrec->rlist,
1696 x, box, fr, &top->idef, graph, fr->born);
1699 if (DOMAINDECOMP(cr))
1701 if (!(cr->duty & DUTY_PME))
1703 wallcycle_start(wcycle, ewcPPDURINGPME);
1704 dd_force_flop_start(cr->dd, nrnb);
1710 /* Enforced rotation has its own cycle counter that starts after the collective
1711 * coordinates have been communicated. It is added to ddCyclF to allow
1712 * for proper load-balancing */
1713 wallcycle_start(wcycle, ewcROT);
1714 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1715 wallcycle_stop(wcycle, ewcROT);
1718 /* Start the force cycle counter.
1719 * This counter is stopped in do_forcelow_level.
1720 * No parallel communication should occur while this counter is running,
1721 * since that will interfere with the dynamic load balancing.
1723 wallcycle_start(wcycle, ewcFORCE);
1727 /* Reset forces for which the virial is calculated separately:
1728 * PME/Ewald forces if necessary */
1729 if (fr->bF_NoVirSum)
1731 if (flags & GMX_FORCE_VIRIAL)
1733 fr->f_novirsum = fr->f_novirsum_alloc;
1736 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1740 clear_rvecs(homenr, fr->f_novirsum+start);
1745 /* We are not calculating the pressure so we do not need
1746 * a separate array for forces that do not contribute
1753 /* Clear the short- and long-range forces */
1754 clear_rvecs(fr->natoms_force_constr, f);
1755 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1757 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1760 clear_rvec(fr->vir_diag_posres);
1762 if (inputrec->ePull == epullCONSTRAINT)
1764 clear_pull_forces(inputrec->pull);
1767 /* update QMMMrec, if necessary */
1770 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1773 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1775 posres_wrapper(fplog, flags, bSepDVDL, inputrec, nrnb, top, box, x,
1776 f, enerd, lambda, fr);
1779 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_FBPOSRES].nr > 0)
1781 /* Flat-bottomed position restraints always require full pbc */
1782 if (!(bStateChanged && bDoAdressWF))
1784 set_pbc(&pbc, inputrec->ePBC, box);
1786 v = fbposres(top->idef.il[F_FBPOSRES].nr, top->idef.il[F_FBPOSRES].iatoms,
1787 top->idef.iparams_fbposres,
1788 (const rvec*)x, fr->f_novirsum, fr->vir_diag_posres,
1789 inputrec->ePBC == epbcNONE ? NULL : &pbc,
1790 fr->rc_scaling, fr->ePBC, fr->posres_com);
1791 enerd->term[F_FBPOSRES] += v;
1792 inc_nrnb(nrnb, eNR_FBPOSRES, top->idef.il[F_FBPOSRES].nr/2);
1795 /* Compute the bonded and non-bonded energies and optionally forces */
1796 do_force_lowlevel(fplog, step, fr, inputrec, &(top->idef),
1797 cr, nrnb, wcycle, mdatoms, &(inputrec->opts),
1798 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, mtop, top, fr->born,
1799 &(top->atomtypes), bBornRadii, box,
1800 inputrec->fepvals, lambda,
1801 graph, &(top->excls), fr->mu_tot,
1807 if (do_per_step(step, inputrec->nstcalclr))
1809 /* Add the long range forces to the short range forces */
1810 for (i = 0; i < fr->natoms_force_constr; i++)
1812 rvec_add(fr->f_twin[i], f[i], f[i]);
1817 cycles_force = wallcycle_stop(wcycle, ewcFORCE);
1821 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1824 if (DOMAINDECOMP(cr))
1826 dd_force_flop_stop(cr->dd, nrnb);
1829 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1835 if (IR_ELEC_FIELD(*inputrec))
1837 /* Compute forces due to electric field */
1838 calc_f_el(MASTER(cr) ? field : NULL,
1839 start, homenr, mdatoms->chargeA, x, fr->f_novirsum,
1840 inputrec->ex, inputrec->et, t);
1843 if (bDoAdressWF && fr->adress_icor == eAdressICThermoForce)
1845 /* Compute thermodynamic force in hybrid AdResS region */
1846 adress_thermo_force(start, homenr, &(top->cgs), x, fr->f_novirsum, fr, mdatoms,
1847 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1850 /* Communicate the forces */
1853 wallcycle_start(wcycle, ewcMOVEF);
1854 if (DOMAINDECOMP(cr))
1856 dd_move_f(cr->dd, f, fr->fshift);
1857 /* Do we need to communicate the separate force array
1858 * for terms that do not contribute to the single sum virial?
1859 * Position restraints and electric fields do not introduce
1860 * inter-cg forces, only full electrostatics methods do.
1861 * When we do not calculate the virial, fr->f_novirsum = f,
1862 * so we have already communicated these forces.
1864 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1865 (flags & GMX_FORCE_VIRIAL))
1867 dd_move_f(cr->dd, fr->f_novirsum, NULL);
1871 /* We should not update the shift forces here,
1872 * since f_twin is already included in f.
1874 dd_move_f(cr->dd, fr->f_twin, NULL);
1879 pd_move_f(cr, f, nrnb);
1882 pd_move_f(cr, fr->f_twin, nrnb);
1885 wallcycle_stop(wcycle, ewcMOVEF);
1888 /* If we have NoVirSum forces, but we do not calculate the virial,
1889 * we sum fr->f_novirum=f later.
1891 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1893 wallcycle_start(wcycle, ewcVSITESPREAD);
1894 spread_vsite_f(fplog, vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1895 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1896 wallcycle_stop(wcycle, ewcVSITESPREAD);
1900 wallcycle_start(wcycle, ewcVSITESPREAD);
1901 spread_vsite_f(fplog, vsite, x, fr->f_twin, NULL, FALSE, NULL,
1903 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1904 wallcycle_stop(wcycle, ewcVSITESPREAD);
1908 if (flags & GMX_FORCE_VIRIAL)
1910 /* Calculation of the virial must be done after vsites! */
1911 calc_virial(fplog, mdatoms->start, mdatoms->homenr, x, f,
1912 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1916 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1918 pull_potential_wrapper(fplog, bSepDVDL, cr, inputrec, box, x,
1919 f, vir_force, mdatoms, enerd, lambda, t);
1922 /* Add the forces from enforced rotation potentials (if any) */
1925 wallcycle_start(wcycle, ewcROTadd);
1926 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
1927 wallcycle_stop(wcycle, ewcROTadd);
1930 if (PAR(cr) && !(cr->duty & DUTY_PME))
1932 /* In case of node-splitting, the PP nodes receive the long-range
1933 * forces, virial and energy from the PME nodes here.
1935 pme_receive_force_ener(fplog, bSepDVDL, cr, wcycle, enerd, fr);
1940 post_process_forces(fplog, cr, step, nrnb, wcycle,
1941 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
1945 /* Sum the potential energy terms from group contributions */
1946 sum_epot(&(inputrec->opts), &(enerd->grpp), enerd->term);
1949 void do_force(FILE *fplog, t_commrec *cr,
1950 t_inputrec *inputrec,
1951 gmx_large_int_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1952 gmx_localtop_t *top,
1954 gmx_groups_t *groups,
1955 matrix box, rvec x[], history_t *hist,
1959 gmx_enerdata_t *enerd, t_fcdata *fcd,
1960 real *lambda, t_graph *graph,
1962 gmx_vsite_t *vsite, rvec mu_tot,
1963 double t, FILE *field, gmx_edsam_t ed,
1964 gmx_bool bBornRadii,
1967 /* modify force flag if not doing nonbonded */
1968 if (!fr->bNonbonded)
1970 flags &= ~GMX_FORCE_NONBONDED;
1973 switch (inputrec->cutoff_scheme)
1976 do_force_cutsVERLET(fplog, cr, inputrec,
1992 do_force_cutsGROUP(fplog, cr, inputrec,
2007 gmx_incons("Invalid cut-off scheme passed!");
2012 void do_constrain_first(FILE *fplog, gmx_constr_t constr,
2013 t_inputrec *ir, t_mdatoms *md,
2014 t_state *state, rvec *f,
2015 t_graph *graph, t_commrec *cr, t_nrnb *nrnb,
2016 t_forcerec *fr, gmx_localtop_t *top, tensor shake_vir)
2018 int i, m, start, end;
2019 gmx_large_int_t step;
2020 real dt = ir->delta_t;
2024 snew(savex, state->natoms);
2027 end = md->homenr + start;
2031 fprintf(debug, "vcm: start=%d, homenr=%d, end=%d\n",
2032 start, md->homenr, end);
2034 /* Do a first constrain to reset particles... */
2035 step = ir->init_step;
2038 char buf[STEPSTRSIZE];
2039 fprintf(fplog, "\nConstraining the starting coordinates (step %s)\n",
2040 gmx_step_str(step, buf));
2044 /* constrain the current position */
2045 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2046 ir, NULL, cr, step, 0, md,
2047 state->x, state->x, NULL,
2048 fr->bMolPBC, state->box,
2049 state->lambda[efptBONDED], &dvdl_dum,
2050 NULL, NULL, nrnb, econqCoord,
2051 ir->epc == epcMTTK, state->veta, state->veta);
2054 /* constrain the inital velocity, and save it */
2055 /* also may be useful if we need the ekin from the halfstep for velocity verlet */
2056 /* might not yet treat veta correctly */
2057 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2058 ir, NULL, cr, step, 0, md,
2059 state->x, state->v, state->v,
2060 fr->bMolPBC, state->box,
2061 state->lambda[efptBONDED], &dvdl_dum,
2062 NULL, NULL, nrnb, econqVeloc,
2063 ir->epc == epcMTTK, state->veta, state->veta);
2065 /* constrain the inital velocities at t-dt/2 */
2066 if (EI_STATE_VELOCITY(ir->eI) && ir->eI != eiVV)
2068 for (i = start; (i < end); i++)
2070 for (m = 0; (m < DIM); m++)
2072 /* Reverse the velocity */
2073 state->v[i][m] = -state->v[i][m];
2074 /* Store the position at t-dt in buf */
2075 savex[i][m] = state->x[i][m] + dt*state->v[i][m];
2078 /* Shake the positions at t=-dt with the positions at t=0
2079 * as reference coordinates.
2083 char buf[STEPSTRSIZE];
2084 fprintf(fplog, "\nConstraining the coordinates at t0-dt (step %s)\n",
2085 gmx_step_str(step, buf));
2088 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2089 ir, NULL, cr, step, -1, md,
2090 state->x, savex, NULL,
2091 fr->bMolPBC, state->box,
2092 state->lambda[efptBONDED], &dvdl_dum,
2093 state->v, NULL, nrnb, econqCoord,
2094 ir->epc == epcMTTK, state->veta, state->veta);
2096 for (i = start; i < end; i++)
2098 for (m = 0; m < DIM; m++)
2100 /* Re-reverse the velocities */
2101 state->v[i][m] = -state->v[i][m];
2108 void calc_enervirdiff(FILE *fplog, int eDispCorr, t_forcerec *fr)
2110 double eners[2], virs[2], enersum, virsum, y0, f, g, h;
2111 double r0, r1, r, rc3, rc9, ea, eb, ec, pa, pb, pc, pd;
2112 double invscale, invscale2, invscale3;
2113 int ri0, ri1, ri, i, offstart, offset;
2114 real scale, *vdwtab, tabfactor, tmp;
2116 fr->enershiftsix = 0;
2117 fr->enershifttwelve = 0;
2118 fr->enerdiffsix = 0;
2119 fr->enerdifftwelve = 0;
2121 fr->virdifftwelve = 0;
2123 if (eDispCorr != edispcNO)
2125 for (i = 0; i < 2; i++)
2130 if ((fr->vdwtype == evdwSWITCH) || (fr->vdwtype == evdwSHIFT))
2132 if (fr->rvdw_switch == 0)
2135 "With dispersion correction rvdw-switch can not be zero "
2136 "for vdw-type = %s", evdw_names[fr->vdwtype]);
2139 scale = fr->nblists[0].table_elec_vdw.scale;
2140 vdwtab = fr->nblists[0].table_vdw.data;
2142 /* Round the cut-offs to exact table values for precision */
2143 ri0 = floor(fr->rvdw_switch*scale);
2144 ri1 = ceil(fr->rvdw*scale);
2150 if (fr->vdwtype == evdwSHIFT)
2152 /* Determine the constant energy shift below rvdw_switch.
2153 * Table has a scale factor since we have scaled it down to compensate
2154 * for scaling-up c6/c12 with the derivative factors to save flops in analytical kernels.
2156 fr->enershiftsix = (real)(-1.0/(rc3*rc3)) - 6.0*vdwtab[8*ri0];
2157 fr->enershifttwelve = (real)( 1.0/(rc9*rc3)) - 12.0*vdwtab[8*ri0 + 4];
2159 /* Add the constant part from 0 to rvdw_switch.
2160 * This integration from 0 to rvdw_switch overcounts the number
2161 * of interactions by 1, as it also counts the self interaction.
2162 * We will correct for this later.
2164 eners[0] += 4.0*M_PI*fr->enershiftsix*rc3/3.0;
2165 eners[1] += 4.0*M_PI*fr->enershifttwelve*rc3/3.0;
2167 invscale = 1.0/(scale);
2168 invscale2 = invscale*invscale;
2169 invscale3 = invscale*invscale2;
2171 /* following summation derived from cubic spline definition,
2172 Numerical Recipies in C, second edition, p. 113-116. Exact
2173 for the cubic spline. We first calculate the negative of
2174 the energy from rvdw to rvdw_switch, assuming that g(r)=1,
2175 and then add the more standard, abrupt cutoff correction to
2176 that result, yielding the long-range correction for a
2177 switched function. We perform both the pressure and energy
2178 loops at the same time for simplicity, as the computational
2181 for (i = 0; i < 2; i++)
2183 enersum = 0.0; virsum = 0.0;
2187 /* Since the dispersion table has been scaled down a factor 6.0 and the repulsion
2188 * a factor 12.0 to compensate for the c6/c12 parameters inside nbfp[] being scaled
2189 * up (to save flops in kernels), we need to correct for this.
2198 for (ri = ri0; ri < ri1; ri++)
2202 eb = 2.0*invscale2*r;
2206 pb = 3.0*invscale2*r;
2207 pc = 3.0*invscale*r*r;
2210 /* this "8" is from the packing in the vdwtab array - perhaps should be #define'ed? */
2211 offset = 8*ri + offstart;
2212 y0 = vdwtab[offset];
2213 f = vdwtab[offset+1];
2214 g = vdwtab[offset+2];
2215 h = vdwtab[offset+3];
2217 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);
2218 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);
2221 enersum *= 4.0*M_PI*tabfactor;
2222 virsum *= 4.0*M_PI*tabfactor;
2223 eners[i] -= enersum;
2227 /* now add the correction for rvdw_switch to infinity */
2228 eners[0] += -4.0*M_PI/(3.0*rc3);
2229 eners[1] += 4.0*M_PI/(9.0*rc9);
2230 virs[0] += 8.0*M_PI/rc3;
2231 virs[1] += -16.0*M_PI/(3.0*rc9);
2233 else if ((fr->vdwtype == evdwCUT) || (fr->vdwtype == evdwUSER))
2235 if (fr->vdwtype == evdwUSER && fplog)
2238 "WARNING: using dispersion correction with user tables\n");
2240 rc3 = fr->rvdw*fr->rvdw*fr->rvdw;
2242 /* Contribution beyond the cut-off */
2243 eners[0] += -4.0*M_PI/(3.0*rc3);
2244 eners[1] += 4.0*M_PI/(9.0*rc9);
2245 if (fr->vdw_modifier == eintmodPOTSHIFT)
2247 /* Contribution within the cut-off */
2248 eners[0] += -4.0*M_PI/(3.0*rc3);
2249 eners[1] += 4.0*M_PI/(3.0*rc9);
2251 /* Contribution beyond the cut-off */
2252 virs[0] += 8.0*M_PI/rc3;
2253 virs[1] += -16.0*M_PI/(3.0*rc9);
2258 "Dispersion correction is not implemented for vdw-type = %s",
2259 evdw_names[fr->vdwtype]);
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_large_int_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, const char *confout,
2476 t_inputrec *inputrec,
2477 t_nrnb nrnb[], gmx_wallcycle_t wcycle,
2478 gmx_runtime_t *runtime,
2479 wallclock_gpu_t *gputimes,
2481 gmx_bool bWriteStat)
2484 t_nrnb *nrnb_tot = NULL;
2488 wallcycle_sum(cr, wcycle);
2494 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM,
2495 cr->mpi_comm_mysim);
2503 #if defined(GMX_MPI) && !defined(GMX_THREAD_MPI)
2506 /* reduce nodetime over all MPI processes in the current simulation */
2508 MPI_Allreduce(&runtime->proctime, &sum, 1, MPI_DOUBLE, MPI_SUM,
2509 cr->mpi_comm_mysim);
2510 runtime->proctime = sum;
2516 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
2523 if ((cr->duty & DUTY_PP) && DOMAINDECOMP(cr))
2525 print_dd_statistics(cr, inputrec, fplog);
2537 snew(nrnb_all, cr->nnodes);
2538 nrnb_all[0] = *nrnb;
2539 for (s = 1; s < cr->nnodes; s++)
2541 MPI_Recv(nrnb_all[s].n, eNRNB, MPI_DOUBLE, s, 0,
2542 cr->mpi_comm_mysim, &stat);
2544 pr_load(fplog, cr, nrnb_all);
2549 MPI_Send(nrnb->n, eNRNB, MPI_DOUBLE, MASTERRANK(cr), 0,
2550 cr->mpi_comm_mysim);
2557 wallcycle_print(fplog, cr->nnodes, cr->npmenodes, runtime->realtime,
2560 if (EI_DYNAMICS(inputrec->eI))
2562 delta_t = inputrec->delta_t;
2571 print_perf(fplog, runtime->proctime, runtime->realtime,
2572 runtime->nsteps_done, delta_t, nbfs, mflop);
2576 print_perf(stderr, runtime->proctime, runtime->realtime,
2577 runtime->nsteps_done, delta_t, nbfs, mflop);
2582 extern void initialize_lambdas(FILE *fplog, t_inputrec *ir, int *fep_state, real *lambda, double *lam0)
2584 /* this function works, but could probably use a logic rewrite to keep all the different
2585 types of efep straight. */
2588 t_lambda *fep = ir->fepvals;
2590 if ((ir->efep == efepNO) && (ir->bSimTemp == FALSE))
2592 for (i = 0; i < efptNR; i++)
2604 *fep_state = fep->init_fep_state; /* this might overwrite the checkpoint
2605 if checkpoint is set -- a kludge is in for now
2607 for (i = 0; i < efptNR; i++)
2609 /* overwrite lambda state with init_lambda for now for backwards compatibility */
2610 if (fep->init_lambda >= 0) /* if it's -1, it was never initializd */
2612 lambda[i] = fep->init_lambda;
2615 lam0[i] = lambda[i];
2620 lambda[i] = fep->all_lambda[i][*fep_state];
2623 lam0[i] = lambda[i];
2629 /* need to rescale control temperatures to match current state */
2630 for (i = 0; i < ir->opts.ngtc; i++)
2632 if (ir->opts.ref_t[i] > 0)
2634 ir->opts.ref_t[i] = ir->simtempvals->temperatures[*fep_state];
2640 /* Send to the log the information on the current lambdas */
2643 fprintf(fplog, "Initial vector of lambda components:[ ");
2644 for (i = 0; i < efptNR; i++)
2646 fprintf(fplog, "%10.4f ", lambda[i]);
2648 fprintf(fplog, "]\n");
2654 void init_md(FILE *fplog,
2655 t_commrec *cr, t_inputrec *ir, const output_env_t oenv,
2656 double *t, double *t0,
2657 real *lambda, int *fep_state, double *lam0,
2658 t_nrnb *nrnb, gmx_mtop_t *mtop,
2660 int nfile, const t_filenm fnm[],
2661 gmx_mdoutf_t **outf, t_mdebin **mdebin,
2662 tensor force_vir, tensor shake_vir, rvec mu_tot,
2663 gmx_bool *bSimAnn, t_vcm **vcm, t_state *state, unsigned long Flags)
2668 /* Initial values */
2669 *t = *t0 = ir->init_t;
2672 for (i = 0; i < ir->opts.ngtc; i++)
2674 /* set bSimAnn if any group is being annealed */
2675 if (ir->opts.annealing[i] != eannNO)
2682 update_annealing_target_temp(&(ir->opts), ir->init_t);
2685 /* Initialize lambda variables */
2686 initialize_lambdas(fplog, ir, fep_state, lambda, lam0);
2690 *upd = init_update(fplog, ir);
2696 *vcm = init_vcm(fplog, &mtop->groups, ir);
2699 if (EI_DYNAMICS(ir->eI) && !(Flags & MD_APPENDFILES))
2701 if (ir->etc == etcBERENDSEN)
2703 please_cite(fplog, "Berendsen84a");
2705 if (ir->etc == etcVRESCALE)
2707 please_cite(fplog, "Bussi2007a");
2715 *outf = init_mdoutf(nfile, fnm, Flags, cr, ir, oenv);
2717 *mdebin = init_mdebin((Flags & MD_APPENDFILES) ? NULL : (*outf)->fp_ene,
2718 mtop, ir, (*outf)->fp_dhdl);
2723 please_cite(fplog, "Fritsch12");
2724 please_cite(fplog, "Junghans10");
2726 /* Initiate variables */
2727 clear_mat(force_vir);
2728 clear_mat(shake_vir);