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
45 #ifdef HAVE_SYS_TIME_H
58 #include "chargegroup.h"
80 #include "pull_rotation.h"
81 #include "gmx_random.h"
84 #include "gmx_wallcycle.h"
86 #include "nbnxn_atomdata.h"
87 #include "nbnxn_search.h"
88 #include "nbnxn_kernels/nbnxn_kernel_ref.h"
89 #include "nbnxn_kernels/simd_4xn/nbnxn_kernel_simd_4xn.h"
90 #include "nbnxn_kernels/simd_2xnn/nbnxn_kernel_simd_2xnn.h"
91 #include "nbnxn_kernels/nbnxn_kernel_gpu_ref.h"
93 #include "gromacs/utility/gmxmpi.h"
98 #include "nbnxn_cuda_data_mgmt.h"
99 #include "nbnxn_cuda/nbnxn_cuda.h"
104 #if _POSIX_TIMERS > 0
105 /* Mac and Windows do not support this */
109 clock_gettime(CLOCK_REALTIME, &t);
110 seconds = (double) t.tv_sec + 1e-9*(double)t.tv_nsec;
112 #elif defined HAVE_GETTIMEOFDAY
116 gettimeofday(&t, NULL);
117 seconds = (double) t.tv_sec + 1e-6*(double)t.tv_usec;
123 seconds = time(NULL);
130 gmx_gettime_per_thread()
132 #if _POSIX_THREAD_CPUTIME > 0
136 clock_gettime(CLOCK_THREAD_CPUTIME_ID, &t);
137 seconds = (double) t.tv_sec + 1e-9*(double)t.tv_nsec;
140 return gmx_gettime();
144 void print_time(FILE *out, gmx_runtime_t *runtime, gmx_large_int_t step,
145 t_inputrec *ir, t_commrec gmx_unused *cr)
148 char timebuf[STRLEN];
149 double dt, time_per_step;
152 #ifndef GMX_THREAD_MPI
158 fprintf(out, "step %s", gmx_step_str(step, buf));
159 if ((step >= ir->nstlist))
161 double seconds_since_epoch = gmx_gettime();
162 dt = seconds_since_epoch - runtime->start_time_stamp;
163 time_per_step = dt/(step - ir->init_step + 1);
165 dt = (ir->nsteps + ir->init_step - step) * time_per_step;
171 finish = (time_t) (seconds_since_epoch + dt);
172 gmx_ctime_r(&finish, timebuf, STRLEN);
173 sprintf(buf, "%s", timebuf);
174 buf[strlen(buf)-1] = '\0';
175 fprintf(out, ", will finish %s", buf);
179 fprintf(out, ", remaining runtime: %5d s ", (int)dt);
184 fprintf(out, " performance: %.1f ns/day ",
185 ir->delta_t/1000*24*60*60/time_per_step);
188 #ifndef GMX_THREAD_MPI
198 void runtime_start(gmx_runtime_t *runtime)
200 runtime->start_time_stamp = gmx_gettime();
201 runtime->start_time_stamp_per_thread = gmx_gettime_per_thread();
202 runtime->elapsed_run_time = 0;
205 void runtime_end(gmx_runtime_t *runtime)
207 double now, now_per_thread;
210 now_per_thread = gmx_gettime_per_thread();
212 runtime->elapsed_run_time = now - runtime->start_time_stamp;
213 runtime->elapsed_run_time_per_thread = now_per_thread - runtime->start_time_stamp_per_thread;
216 double runtime_get_elapsed_time(gmx_runtime_t *runtime)
218 return gmx_gettime() - runtime->start_time_stamp;
221 void print_date_and_time(FILE *fplog, int nodeid, const char *title,
222 const gmx_runtime_t *runtime)
225 char timebuf[STRLEN];
226 char time_string[STRLEN];
233 tmptime = (time_t) runtime->start_time_stamp;
234 gmx_ctime_r(&tmptime, timebuf, STRLEN);
238 tmptime = (time_t) gmx_gettime();
239 gmx_ctime_r(&tmptime, timebuf, STRLEN);
241 for (i = 0; timebuf[i] >= ' '; i++)
243 time_string[i] = timebuf[i];
245 time_string[i] = '\0';
247 fprintf(fplog, "%s on node %d %s\n", title, nodeid, time_string);
251 static void sum_forces(int start, int end, rvec f[], rvec flr[])
257 pr_rvecs(debug, 0, "fsr", f+start, end-start);
258 pr_rvecs(debug, 0, "flr", flr+start, end-start);
260 for (i = start; (i < end); i++)
262 rvec_inc(f[i], flr[i]);
267 * calc_f_el calculates forces due to an electric field.
269 * force is kJ mol^-1 nm^-1 = e * kJ mol^-1 nm^-1 / e
271 * Et[] contains the parameters for the time dependent
272 * part of the field (not yet used).
273 * Ex[] contains the parameters for
274 * the spatial dependent part of the field. You can have cool periodic
275 * fields in principle, but only a constant field is supported
277 * The function should return the energy due to the electric field
278 * (if any) but for now returns 0.
281 * There can be problems with the virial.
282 * Since the field is not self-consistent this is unavoidable.
283 * For neutral molecules the virial is correct within this approximation.
284 * For neutral systems with many charged molecules the error is small.
285 * But for systems with a net charge or a few charged molecules
286 * the error can be significant when the field is high.
287 * Solution: implement a self-consitent electric field into PME.
289 static void calc_f_el(FILE *fp, int start, int homenr,
290 real charge[], rvec f[],
291 t_cosines Ex[], t_cosines Et[], double t)
297 for (m = 0; (m < DIM); m++)
304 Ext[m] = cos(Et[m].a[0]*(t-t0))*exp(-sqr(t-t0)/(2.0*sqr(Et[m].a[2])));
308 Ext[m] = cos(Et[m].a[0]*t);
317 /* Convert the field strength from V/nm to MD-units */
318 Ext[m] *= Ex[m].a[0]*FIELDFAC;
319 for (i = start; (i < start+homenr); i++)
321 f[i][m] += charge[i]*Ext[m];
331 fprintf(fp, "%10g %10g %10g %10g #FIELD\n", t,
332 Ext[XX]/FIELDFAC, Ext[YY]/FIELDFAC, Ext[ZZ]/FIELDFAC);
336 static void calc_virial(int start, int homenr, rvec x[], rvec f[],
337 tensor vir_part, t_graph *graph, matrix box,
338 t_nrnb *nrnb, const t_forcerec *fr, int ePBC)
343 /* The short-range virial from surrounding boxes */
345 calc_vir(SHIFTS, fr->shift_vec, fr->fshift, vir_part, ePBC == epbcSCREW, box);
346 inc_nrnb(nrnb, eNR_VIRIAL, SHIFTS);
348 /* Calculate partial virial, for local atoms only, based on short range.
349 * Total virial is computed in global_stat, called from do_md
351 f_calc_vir(start, start+homenr, x, f, vir_part, graph, box);
352 inc_nrnb(nrnb, eNR_VIRIAL, homenr);
354 /* Add position restraint contribution */
355 for (i = 0; i < DIM; i++)
357 vir_part[i][i] += fr->vir_diag_posres[i];
360 /* Add wall contribution */
361 for (i = 0; i < DIM; i++)
363 vir_part[i][ZZ] += fr->vir_wall_z[i];
368 pr_rvecs(debug, 0, "vir_part", vir_part, DIM);
372 static void posres_wrapper(FILE *fplog,
378 matrix box, rvec x[],
379 gmx_enerdata_t *enerd,
387 /* Position restraints always require full pbc */
388 set_pbc(&pbc, ir->ePBC, box);
390 v = posres(top->idef.il[F_POSRES].nr, top->idef.il[F_POSRES].iatoms,
391 top->idef.iparams_posres,
392 (const rvec*)x, fr->f_novirsum, fr->vir_diag_posres,
393 ir->ePBC == epbcNONE ? NULL : &pbc,
394 lambda[efptRESTRAINT], &dvdl,
395 fr->rc_scaling, fr->ePBC, fr->posres_com, fr->posres_comB);
398 gmx_print_sepdvdl(fplog, interaction_function[F_POSRES].longname, v, dvdl);
400 enerd->term[F_POSRES] += v;
401 /* If just the force constant changes, the FEP term is linear,
402 * but if k changes, it is not.
404 enerd->dvdl_nonlin[efptRESTRAINT] += dvdl;
405 inc_nrnb(nrnb, eNR_POSRES, top->idef.il[F_POSRES].nr/2);
407 if ((ir->fepvals->n_lambda > 0) && (flags & GMX_FORCE_DHDL))
409 for (i = 0; i < enerd->n_lambda; i++)
411 real dvdl_dum, lambda_dum;
413 lambda_dum = (i == 0 ? lambda[efptRESTRAINT] : ir->fepvals->all_lambda[efptRESTRAINT][i-1]);
414 v = posres(top->idef.il[F_POSRES].nr, top->idef.il[F_POSRES].iatoms,
415 top->idef.iparams_posres,
416 (const rvec*)x, NULL, NULL,
417 ir->ePBC == epbcNONE ? NULL : &pbc, lambda_dum, &dvdl,
418 fr->rc_scaling, fr->ePBC, fr->posres_com, fr->posres_comB);
419 enerd->enerpart_lambda[i] += v;
424 static void pull_potential_wrapper(FILE *fplog,
428 matrix box, rvec x[],
432 gmx_enerdata_t *enerd,
439 /* Calculate the center of mass forces, this requires communication,
440 * which is why pull_potential is called close to other communication.
441 * The virial contribution is calculated directly,
442 * which is why we call pull_potential after calc_virial.
444 set_pbc(&pbc, ir->ePBC, box);
446 enerd->term[F_COM_PULL] +=
447 pull_potential(ir->ePull, ir->pull, mdatoms, &pbc,
448 cr, t, lambda[efptRESTRAINT], x, f, vir_force, &dvdl);
451 gmx_print_sepdvdl(fplog, "Com pull", enerd->term[F_COM_PULL], dvdl);
453 enerd->dvdl_lin[efptRESTRAINT] += dvdl;
456 static void pme_receive_force_ener(FILE *fplog,
459 gmx_wallcycle_t wcycle,
460 gmx_enerdata_t *enerd,
464 float cycles_ppdpme, cycles_seppme;
466 cycles_ppdpme = wallcycle_stop(wcycle, ewcPPDURINGPME);
467 dd_cycles_add(cr->dd, cycles_ppdpme, ddCyclPPduringPME);
469 /* In case of node-splitting, the PP nodes receive the long-range
470 * forces, virial and energy from the PME nodes here.
472 wallcycle_start(wcycle, ewcPP_PMEWAITRECVF);
474 gmx_pme_receive_f(cr, fr->f_novirsum, fr->vir_el_recip, &e, &dvdl,
478 gmx_print_sepdvdl(fplog, "PME mesh", e, dvdl);
480 enerd->term[F_COUL_RECIP] += e;
481 enerd->dvdl_lin[efptCOUL] += dvdl;
484 dd_cycles_add(cr->dd, cycles_seppme, ddCyclPME);
486 wallcycle_stop(wcycle, ewcPP_PMEWAITRECVF);
489 static void print_large_forces(FILE *fp, t_mdatoms *md, t_commrec *cr,
490 gmx_large_int_t step, real pforce, rvec *x, rvec *f)
494 char buf[STEPSTRSIZE];
497 for (i = md->start; i < md->start+md->homenr; i++)
500 /* We also catch NAN, if the compiler does not optimize this away. */
501 if (fn2 >= pf2 || fn2 != fn2)
503 fprintf(fp, "step %s atom %6d x %8.3f %8.3f %8.3f force %12.5e\n",
504 gmx_step_str(step, buf),
505 ddglatnr(cr->dd, i), x[i][XX], x[i][YY], x[i][ZZ], sqrt(fn2));
510 static void post_process_forces(t_commrec *cr,
511 gmx_large_int_t step,
512 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
514 matrix box, rvec x[],
519 t_forcerec *fr, gmx_vsite_t *vsite,
526 /* Spread the mesh force on virtual sites to the other particles...
527 * This is parallellized. MPI communication is performed
528 * if the constructing atoms aren't local.
530 wallcycle_start(wcycle, ewcVSITESPREAD);
531 spread_vsite_f(vsite, x, fr->f_novirsum, NULL,
532 (flags & GMX_FORCE_VIRIAL), fr->vir_el_recip,
534 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
535 wallcycle_stop(wcycle, ewcVSITESPREAD);
537 if (flags & GMX_FORCE_VIRIAL)
539 /* Now add the forces, this is local */
542 sum_forces(0, fr->f_novirsum_n, f, fr->f_novirsum);
546 sum_forces(mdatoms->start, mdatoms->start+mdatoms->homenr,
549 if (EEL_FULL(fr->eeltype))
551 /* Add the mesh contribution to the virial */
552 m_add(vir_force, fr->vir_el_recip, vir_force);
556 pr_rvecs(debug, 0, "vir_force", vir_force, DIM);
561 if (fr->print_force >= 0)
563 print_large_forces(stderr, mdatoms, cr, step, fr->print_force, x, f);
567 static void do_nb_verlet(t_forcerec *fr,
568 interaction_const_t *ic,
569 gmx_enerdata_t *enerd,
570 int flags, int ilocality,
574 int nnbl, kernel_type, enr_nbnxn_kernel_ljc, enr_nbnxn_kernel_lj;
576 nonbonded_verlet_group_t *nbvg;
579 if (!(flags & GMX_FORCE_NONBONDED))
581 /* skip non-bonded calculation */
585 nbvg = &fr->nbv->grp[ilocality];
587 /* CUDA kernel launch overhead is already timed separately */
588 if (fr->cutoff_scheme != ecutsVERLET)
590 gmx_incons("Invalid cut-off scheme passed!");
593 bCUDA = (nbvg->kernel_type == nbnxnk8x8x8_CUDA);
597 wallcycle_sub_start(wcycle, ewcsNONBONDED);
599 switch (nbvg->kernel_type)
601 case nbnxnk4x4_PlainC:
602 nbnxn_kernel_ref(&nbvg->nbl_lists,
608 enerd->grpp.ener[egCOULSR],
610 enerd->grpp.ener[egBHAMSR] :
611 enerd->grpp.ener[egLJSR]);
614 case nbnxnk4xN_SIMD_4xN:
615 nbnxn_kernel_simd_4xn(&nbvg->nbl_lists,
622 enerd->grpp.ener[egCOULSR],
624 enerd->grpp.ener[egBHAMSR] :
625 enerd->grpp.ener[egLJSR]);
627 case nbnxnk4xN_SIMD_2xNN:
628 nbnxn_kernel_simd_2xnn(&nbvg->nbl_lists,
635 enerd->grpp.ener[egCOULSR],
637 enerd->grpp.ener[egBHAMSR] :
638 enerd->grpp.ener[egLJSR]);
641 case nbnxnk8x8x8_CUDA:
642 nbnxn_cuda_launch_kernel(fr->nbv->cu_nbv, nbvg->nbat, flags, ilocality);
645 case nbnxnk8x8x8_PlainC:
646 nbnxn_kernel_gpu_ref(nbvg->nbl_lists.nbl[0],
651 nbvg->nbat->out[0].f,
653 enerd->grpp.ener[egCOULSR],
655 enerd->grpp.ener[egBHAMSR] :
656 enerd->grpp.ener[egLJSR]);
660 gmx_incons("Invalid nonbonded kernel type passed!");
665 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
668 if (EEL_RF(ic->eeltype) || ic->eeltype == eelCUT)
670 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_RF;
672 else if ((!bCUDA && nbvg->ewald_excl == ewaldexclAnalytical) ||
673 (bCUDA && nbnxn_cuda_is_kernel_ewald_analytical(fr->nbv->cu_nbv)))
675 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_EWALD;
679 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_TAB;
681 enr_nbnxn_kernel_lj = eNR_NBNXN_LJ;
682 if (flags & GMX_FORCE_ENERGY)
684 /* In eNR_??? the nbnxn F+E kernels are always the F kernel + 1 */
685 enr_nbnxn_kernel_ljc += 1;
686 enr_nbnxn_kernel_lj += 1;
689 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc,
690 nbvg->nbl_lists.natpair_ljq);
691 inc_nrnb(nrnb, enr_nbnxn_kernel_lj,
692 nbvg->nbl_lists.natpair_lj);
693 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc-eNR_NBNXN_LJ_RF+eNR_NBNXN_RF,
694 nbvg->nbl_lists.natpair_q);
697 void do_force_cutsVERLET(FILE *fplog, t_commrec *cr,
698 t_inputrec *inputrec,
699 gmx_large_int_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
701 gmx_groups_t gmx_unused *groups,
702 matrix box, rvec x[], history_t *hist,
706 gmx_enerdata_t *enerd, t_fcdata *fcd,
707 real *lambda, t_graph *graph,
708 t_forcerec *fr, interaction_const_t *ic,
709 gmx_vsite_t *vsite, rvec mu_tot,
710 double t, FILE *field, gmx_edsam_t ed,
718 gmx_bool bSepDVDL, bStateChanged, bNS, bFillGrid, bCalcCGCM, bBS;
719 gmx_bool bDoLongRange, bDoForces, bSepLRF, bUseGPU, bUseOrEmulGPU;
720 gmx_bool bDiffKernels = FALSE;
722 rvec vzero, box_diag;
724 float cycles_pme, cycles_force;
725 nonbonded_verlet_t *nbv;
729 nb_kernel_type = fr->nbv->grp[0].kernel_type;
731 start = mdatoms->start;
732 homenr = mdatoms->homenr;
734 bSepDVDL = (fr->bSepDVDL && do_per_step(step, inputrec->nstlog));
736 clear_mat(vir_force);
739 if (DOMAINDECOMP(cr))
741 cg1 = cr->dd->ncg_tot;
752 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
753 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
754 bFillGrid = (bNS && bStateChanged);
755 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
756 bDoLongRange = (fr->bTwinRange && bNS && (flags & GMX_FORCE_DO_LR));
757 bDoForces = (flags & GMX_FORCE_FORCES);
758 bSepLRF = (bDoLongRange && bDoForces && (flags & GMX_FORCE_SEPLRF));
759 bUseGPU = fr->nbv->bUseGPU;
760 bUseOrEmulGPU = bUseGPU || (nbv->grp[0].kernel_type == nbnxnk8x8x8_PlainC);
764 update_forcerec(fr, box);
766 if (NEED_MUTOT(*inputrec))
768 /* Calculate total (local) dipole moment in a temporary common array.
769 * This makes it possible to sum them over nodes faster.
771 calc_mu(start, homenr,
772 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
777 if (fr->ePBC != epbcNONE)
779 /* Compute shift vectors every step,
780 * because of pressure coupling or box deformation!
782 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
784 calc_shifts(box, fr->shift_vec);
789 put_atoms_in_box_omp(fr->ePBC, box, homenr, x);
790 inc_nrnb(nrnb, eNR_SHIFTX, homenr);
792 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
794 unshift_self(graph, box, x);
798 nbnxn_atomdata_copy_shiftvec(flags & GMX_FORCE_DYNAMICBOX,
799 fr->shift_vec, nbv->grp[0].nbat);
802 if (!(cr->duty & DUTY_PME))
804 /* Send particle coordinates to the pme nodes.
805 * Since this is only implemented for domain decomposition
806 * and domain decomposition does not use the graph,
807 * we do not need to worry about shifting.
810 wallcycle_start(wcycle, ewcPP_PMESENDX);
812 bBS = (inputrec->nwall == 2);
816 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
819 gmx_pme_send_x(cr, bBS ? boxs : box, x,
820 mdatoms->nChargePerturbed, lambda[efptCOUL],
821 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)), step);
823 wallcycle_stop(wcycle, ewcPP_PMESENDX);
827 /* do gridding for pair search */
830 if (graph && bStateChanged)
832 /* Calculate intramolecular shift vectors to make molecules whole */
833 mk_mshift(fplog, graph, fr->ePBC, box, x);
837 box_diag[XX] = box[XX][XX];
838 box_diag[YY] = box[YY][YY];
839 box_diag[ZZ] = box[ZZ][ZZ];
841 wallcycle_start(wcycle, ewcNS);
844 wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
845 nbnxn_put_on_grid(nbv->nbs, fr->ePBC, box,
847 0, mdatoms->homenr, -1, fr->cginfo, x,
849 nbv->grp[eintLocal].kernel_type,
850 nbv->grp[eintLocal].nbat);
851 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
855 wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
856 nbnxn_put_on_grid_nonlocal(nbv->nbs, domdec_zones(cr->dd),
858 nbv->grp[eintNonlocal].kernel_type,
859 nbv->grp[eintNonlocal].nbat);
860 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
863 if (nbv->ngrp == 1 ||
864 nbv->grp[eintNonlocal].nbat == nbv->grp[eintLocal].nbat)
866 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatAll,
867 nbv->nbs, mdatoms, fr->cginfo);
871 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatLocal,
872 nbv->nbs, mdatoms, fr->cginfo);
873 nbnxn_atomdata_set(nbv->grp[eintNonlocal].nbat, eatAll,
874 nbv->nbs, mdatoms, fr->cginfo);
876 wallcycle_stop(wcycle, ewcNS);
879 /* initialize the GPU atom data and copy shift vector */
884 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
885 nbnxn_cuda_init_atomdata(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
886 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
889 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
890 nbnxn_cuda_upload_shiftvec(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
891 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
894 /* do local pair search */
897 wallcycle_start_nocount(wcycle, ewcNS);
898 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
899 nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintLocal].nbat,
902 nbv->min_ci_balanced,
903 &nbv->grp[eintLocal].nbl_lists,
905 nbv->grp[eintLocal].kernel_type,
907 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
911 /* initialize local pair-list on the GPU */
912 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
913 nbv->grp[eintLocal].nbl_lists.nbl[0],
916 wallcycle_stop(wcycle, ewcNS);
920 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
921 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
922 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, FALSE, x,
923 nbv->grp[eintLocal].nbat);
924 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
925 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
930 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
931 /* launch local nonbonded F on GPU */
932 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFNo,
934 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
937 /* Communicate coordinates and sum dipole if necessary +
938 do non-local pair search */
939 if (DOMAINDECOMP(cr))
941 bDiffKernels = (nbv->grp[eintNonlocal].kernel_type !=
942 nbv->grp[eintLocal].kernel_type);
946 /* With GPU+CPU non-bonded calculations we need to copy
947 * the local coordinates to the non-local nbat struct
948 * (in CPU format) as the non-local kernel call also
949 * calculates the local - non-local interactions.
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, TRUE, x,
954 nbv->grp[eintNonlocal].nbat);
955 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
956 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
961 wallcycle_start_nocount(wcycle, ewcNS);
962 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
966 nbnxn_grid_add_simple(nbv->nbs, nbv->grp[eintNonlocal].nbat);
969 nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintNonlocal].nbat,
972 nbv->min_ci_balanced,
973 &nbv->grp[eintNonlocal].nbl_lists,
975 nbv->grp[eintNonlocal].kernel_type,
978 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
980 if (nbv->grp[eintNonlocal].kernel_type == nbnxnk8x8x8_CUDA)
982 /* initialize non-local pair-list on the GPU */
983 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
984 nbv->grp[eintNonlocal].nbl_lists.nbl[0],
987 wallcycle_stop(wcycle, ewcNS);
991 wallcycle_start(wcycle, ewcMOVEX);
992 dd_move_x(cr->dd, box, x);
994 /* When we don't need the total dipole we sum it in global_stat */
995 if (bStateChanged && NEED_MUTOT(*inputrec))
997 gmx_sumd(2*DIM, mu, cr);
999 wallcycle_stop(wcycle, ewcMOVEX);
1001 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1002 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
1003 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatNonlocal, FALSE, x,
1004 nbv->grp[eintNonlocal].nbat);
1005 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
1006 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1009 if (bUseGPU && !bDiffKernels)
1011 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
1012 /* launch non-local nonbonded F on GPU */
1013 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFNo,
1015 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1021 /* launch D2H copy-back F */
1022 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1023 if (DOMAINDECOMP(cr) && !bDiffKernels)
1025 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintNonlocal].nbat,
1026 flags, eatNonlocal);
1028 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintLocal].nbat,
1030 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1033 if (bStateChanged && NEED_MUTOT(*inputrec))
1037 gmx_sumd(2*DIM, mu, cr);
1040 for (i = 0; i < 2; i++)
1042 for (j = 0; j < DIM; j++)
1044 fr->mu_tot[i][j] = mu[i*DIM + j];
1048 if (fr->efep == efepNO)
1050 copy_rvec(fr->mu_tot[0], mu_tot);
1054 for (j = 0; j < DIM; j++)
1057 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] +
1058 lambda[efptCOUL]*fr->mu_tot[1][j];
1062 /* Reset energies */
1063 reset_enerdata(fr, bNS, enerd, MASTER(cr));
1064 clear_rvecs(SHIFTS, fr->fshift);
1066 if (DOMAINDECOMP(cr))
1068 if (!(cr->duty & DUTY_PME))
1070 wallcycle_start(wcycle, ewcPPDURINGPME);
1071 dd_force_flop_start(cr->dd, nrnb);
1077 /* Enforced rotation has its own cycle counter that starts after the collective
1078 * coordinates have been communicated. It is added to ddCyclF to allow
1079 * for proper load-balancing */
1080 wallcycle_start(wcycle, ewcROT);
1081 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1082 wallcycle_stop(wcycle, ewcROT);
1085 /* Start the force cycle counter.
1086 * This counter is stopped in do_forcelow_level.
1087 * No parallel communication should occur while this counter is running,
1088 * since that will interfere with the dynamic load balancing.
1090 wallcycle_start(wcycle, ewcFORCE);
1093 /* Reset forces for which the virial is calculated separately:
1094 * PME/Ewald forces if necessary */
1095 if (fr->bF_NoVirSum)
1097 if (flags & GMX_FORCE_VIRIAL)
1099 fr->f_novirsum = fr->f_novirsum_alloc;
1102 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1106 clear_rvecs(homenr, fr->f_novirsum+start);
1111 /* We are not calculating the pressure so we do not need
1112 * a separate array for forces that do not contribute
1119 /* Clear the short- and long-range forces */
1120 clear_rvecs(fr->natoms_force_constr, f);
1121 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1123 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1126 clear_rvec(fr->vir_diag_posres);
1129 if (inputrec->ePull == epullCONSTRAINT)
1131 clear_pull_forces(inputrec->pull);
1134 /* We calculate the non-bonded forces, when done on the CPU, here.
1135 * We do this before calling do_force_lowlevel, as in there bondeds
1136 * forces are calculated before PME, which does communication.
1137 * With this order, non-bonded and bonded force calculation imbalance
1138 * can be balanced out by the domain decomposition load balancing.
1143 /* Maybe we should move this into do_force_lowlevel */
1144 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFYes,
1148 if (!bUseOrEmulGPU || bDiffKernels)
1152 if (DOMAINDECOMP(cr))
1154 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal,
1155 bDiffKernels ? enbvClearFYes : enbvClearFNo,
1165 aloc = eintNonlocal;
1168 /* Add all the non-bonded force to the normal force array.
1169 * This can be split into a local a non-local part when overlapping
1170 * communication with calculation with domain decomposition.
1172 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1173 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1174 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1175 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatAll, nbv->grp[aloc].nbat, f);
1176 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1177 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1178 wallcycle_start_nocount(wcycle, ewcFORCE);
1180 /* if there are multiple fshift output buffers reduce them */
1181 if ((flags & GMX_FORCE_VIRIAL) &&
1182 nbv->grp[aloc].nbl_lists.nnbl > 1)
1184 nbnxn_atomdata_add_nbat_fshift_to_fshift(nbv->grp[aloc].nbat,
1189 /* update QMMMrec, if necessary */
1192 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1195 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1197 posres_wrapper(fplog, flags, bSepDVDL, inputrec, nrnb, top, box, x,
1201 /* Compute the bonded and non-bonded energies and optionally forces */
1202 do_force_lowlevel(fplog, step, fr, inputrec, &(top->idef),
1203 cr, nrnb, wcycle, mdatoms,
1204 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, top, fr->born,
1205 &(top->atomtypes), bBornRadii, box,
1206 inputrec->fepvals, lambda, graph, &(top->excls), fr->mu_tot,
1207 flags, &cycles_pme);
1211 if (do_per_step(step, inputrec->nstcalclr))
1213 /* Add the long range forces to the short range forces */
1214 for (i = 0; i < fr->natoms_force_constr; i++)
1216 rvec_add(fr->f_twin[i], f[i], f[i]);
1221 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1225 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1228 if (bUseOrEmulGPU && !bDiffKernels)
1230 /* wait for non-local forces (or calculate in emulation mode) */
1231 if (DOMAINDECOMP(cr))
1235 wallcycle_start(wcycle, ewcWAIT_GPU_NB_NL);
1236 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1237 nbv->grp[eintNonlocal].nbat,
1239 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1241 cycles_force += wallcycle_stop(wcycle, ewcWAIT_GPU_NB_NL);
1245 wallcycle_start_nocount(wcycle, ewcFORCE);
1246 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFYes,
1248 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1250 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1251 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1252 /* skip the reduction if there was no non-local work to do */
1253 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1255 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatNonlocal,
1256 nbv->grp[eintNonlocal].nbat, f);
1258 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1259 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1265 /* Communicate the forces */
1268 wallcycle_start(wcycle, ewcMOVEF);
1269 if (DOMAINDECOMP(cr))
1271 dd_move_f(cr->dd, f, fr->fshift);
1272 /* Do we need to communicate the separate force array
1273 * for terms that do not contribute to the single sum virial?
1274 * Position restraints and electric fields do not introduce
1275 * inter-cg forces, only full electrostatics methods do.
1276 * When we do not calculate the virial, fr->f_novirsum = f,
1277 * so we have already communicated these forces.
1279 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1280 (flags & GMX_FORCE_VIRIAL))
1282 dd_move_f(cr->dd, fr->f_novirsum, NULL);
1286 /* We should not update the shift forces here,
1287 * since f_twin is already included in f.
1289 dd_move_f(cr->dd, fr->f_twin, NULL);
1292 wallcycle_stop(wcycle, ewcMOVEF);
1298 /* wait for local forces (or calculate in emulation mode) */
1301 wallcycle_start(wcycle, ewcWAIT_GPU_NB_L);
1302 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1303 nbv->grp[eintLocal].nbat,
1305 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1307 wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
1309 /* now clear the GPU outputs while we finish the step on the CPU */
1311 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1312 nbnxn_cuda_clear_outputs(nbv->cu_nbv, flags);
1313 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1317 wallcycle_start_nocount(wcycle, ewcFORCE);
1318 do_nb_verlet(fr, ic, enerd, flags, eintLocal,
1319 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
1321 wallcycle_stop(wcycle, ewcFORCE);
1323 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1324 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1325 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1327 /* skip the reduction if there was no non-local work to do */
1328 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatLocal,
1329 nbv->grp[eintLocal].nbat, f);
1331 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1332 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1335 if (DOMAINDECOMP(cr))
1337 dd_force_flop_stop(cr->dd, nrnb);
1340 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1346 if (IR_ELEC_FIELD(*inputrec))
1348 /* Compute forces due to electric field */
1349 calc_f_el(MASTER(cr) ? field : NULL,
1350 start, homenr, mdatoms->chargeA, fr->f_novirsum,
1351 inputrec->ex, inputrec->et, t);
1354 /* If we have NoVirSum forces, but we do not calculate the virial,
1355 * we sum fr->f_novirum=f later.
1357 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1359 wallcycle_start(wcycle, ewcVSITESPREAD);
1360 spread_vsite_f(vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1361 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1362 wallcycle_stop(wcycle, ewcVSITESPREAD);
1366 wallcycle_start(wcycle, ewcVSITESPREAD);
1367 spread_vsite_f(vsite, x, fr->f_twin, NULL, FALSE, NULL,
1369 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1370 wallcycle_stop(wcycle, ewcVSITESPREAD);
1374 if (flags & GMX_FORCE_VIRIAL)
1376 /* Calculation of the virial must be done after vsites! */
1377 calc_virial(mdatoms->start, mdatoms->homenr, x, f,
1378 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1382 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1384 pull_potential_wrapper(fplog, bSepDVDL, cr, inputrec, box, x,
1385 f, vir_force, mdatoms, enerd, lambda, t);
1388 /* Add the forces from enforced rotation potentials (if any) */
1391 wallcycle_start(wcycle, ewcROTadd);
1392 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
1393 wallcycle_stop(wcycle, ewcROTadd);
1396 if (PAR(cr) && !(cr->duty & DUTY_PME))
1398 /* In case of node-splitting, the PP nodes receive the long-range
1399 * forces, virial and energy from the PME nodes here.
1401 pme_receive_force_ener(fplog, bSepDVDL, cr, wcycle, enerd, fr);
1406 post_process_forces(cr, step, nrnb, wcycle,
1407 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
1411 /* Sum the potential energy terms from group contributions */
1412 sum_epot(&(enerd->grpp), enerd->term);
1415 void do_force_cutsGROUP(FILE *fplog, t_commrec *cr,
1416 t_inputrec *inputrec,
1417 gmx_large_int_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1418 gmx_localtop_t *top,
1419 gmx_groups_t *groups,
1420 matrix box, rvec x[], history_t *hist,
1424 gmx_enerdata_t *enerd, t_fcdata *fcd,
1425 real *lambda, t_graph *graph,
1426 t_forcerec *fr, gmx_vsite_t *vsite, rvec mu_tot,
1427 double t, FILE *field, gmx_edsam_t ed,
1428 gmx_bool bBornRadii,
1434 gmx_bool bSepDVDL, bStateChanged, bNS, bFillGrid, bCalcCGCM, bBS;
1435 gmx_bool bDoLongRangeNS, bDoForces, bDoPotential, bSepLRF;
1436 gmx_bool bDoAdressWF;
1438 rvec vzero, box_diag;
1439 real e, v, dvdlambda[efptNR];
1441 float cycles_pme, cycles_force;
1443 start = mdatoms->start;
1444 homenr = mdatoms->homenr;
1446 bSepDVDL = (fr->bSepDVDL && do_per_step(step, inputrec->nstlog));
1448 clear_mat(vir_force);
1452 pd_cg_range(cr, &cg0, &cg1);
1457 if (DOMAINDECOMP(cr))
1459 cg1 = cr->dd->ncg_tot;
1471 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
1472 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
1473 /* Should we update the long-range neighborlists at this step? */
1474 bDoLongRangeNS = fr->bTwinRange && bNS;
1475 /* Should we perform the long-range nonbonded evaluation inside the neighborsearching? */
1476 bFillGrid = (bNS && bStateChanged);
1477 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
1478 bDoForces = (flags & GMX_FORCE_FORCES);
1479 bDoPotential = (flags & GMX_FORCE_ENERGY);
1480 bSepLRF = ((inputrec->nstcalclr > 1) && bDoForces &&
1481 (flags & GMX_FORCE_SEPLRF) && (flags & GMX_FORCE_DO_LR));
1483 /* should probably move this to the forcerec since it doesn't change */
1484 bDoAdressWF = ((fr->adress_type != eAdressOff));
1488 update_forcerec(fr, box);
1490 if (NEED_MUTOT(*inputrec))
1492 /* Calculate total (local) dipole moment in a temporary common array.
1493 * This makes it possible to sum them over nodes faster.
1495 calc_mu(start, homenr,
1496 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
1501 if (fr->ePBC != epbcNONE)
1503 /* Compute shift vectors every step,
1504 * because of pressure coupling or box deformation!
1506 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
1508 calc_shifts(box, fr->shift_vec);
1513 put_charge_groups_in_box(fplog, cg0, cg1, fr->ePBC, box,
1514 &(top->cgs), x, fr->cg_cm);
1515 inc_nrnb(nrnb, eNR_CGCM, homenr);
1516 inc_nrnb(nrnb, eNR_RESETX, cg1-cg0);
1518 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
1520 unshift_self(graph, box, x);
1525 calc_cgcm(fplog, cg0, cg1, &(top->cgs), x, fr->cg_cm);
1526 inc_nrnb(nrnb, eNR_CGCM, homenr);
1533 move_cgcm(fplog, cr, fr->cg_cm);
1537 pr_rvecs(debug, 0, "cgcm", fr->cg_cm, top->cgs.nr);
1542 if (!(cr->duty & DUTY_PME))
1544 /* Send particle coordinates to the pme nodes.
1545 * Since this is only implemented for domain decomposition
1546 * and domain decomposition does not use the graph,
1547 * we do not need to worry about shifting.
1550 wallcycle_start(wcycle, ewcPP_PMESENDX);
1552 bBS = (inputrec->nwall == 2);
1555 copy_mat(box, boxs);
1556 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
1559 gmx_pme_send_x(cr, bBS ? boxs : box, x,
1560 mdatoms->nChargePerturbed, lambda[efptCOUL],
1561 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)), step);
1563 wallcycle_stop(wcycle, ewcPP_PMESENDX);
1565 #endif /* GMX_MPI */
1567 /* Communicate coordinates and sum dipole if necessary */
1570 wallcycle_start(wcycle, ewcMOVEX);
1571 if (DOMAINDECOMP(cr))
1573 dd_move_x(cr->dd, box, x);
1577 move_x(cr, x, nrnb);
1579 wallcycle_stop(wcycle, ewcMOVEX);
1582 /* update adress weight beforehand */
1583 if (bStateChanged && bDoAdressWF)
1585 /* need pbc for adress weight calculation with pbc_dx */
1586 set_pbc(&pbc, inputrec->ePBC, box);
1587 if (fr->adress_site == eAdressSITEcog)
1589 update_adress_weights_cog(top->idef.iparams, top->idef.il, x, fr, mdatoms,
1590 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1592 else if (fr->adress_site == eAdressSITEcom)
1594 update_adress_weights_com(fplog, cg0, cg1, &(top->cgs), x, fr, mdatoms,
1595 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1597 else if (fr->adress_site == eAdressSITEatomatom)
1599 update_adress_weights_atom_per_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1600 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1604 update_adress_weights_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1605 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1609 if (NEED_MUTOT(*inputrec))
1616 gmx_sumd(2*DIM, mu, cr);
1618 for (i = 0; i < 2; i++)
1620 for (j = 0; j < DIM; j++)
1622 fr->mu_tot[i][j] = mu[i*DIM + j];
1626 if (fr->efep == efepNO)
1628 copy_rvec(fr->mu_tot[0], mu_tot);
1632 for (j = 0; j < DIM; j++)
1635 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] + lambda[efptCOUL]*fr->mu_tot[1][j];
1640 /* Reset energies */
1641 reset_enerdata(fr, bNS, enerd, MASTER(cr));
1642 clear_rvecs(SHIFTS, fr->fshift);
1646 wallcycle_start(wcycle, ewcNS);
1648 if (graph && bStateChanged)
1650 /* Calculate intramolecular shift vectors to make molecules whole */
1651 mk_mshift(fplog, graph, fr->ePBC, box, x);
1654 /* Do the actual neighbour searching */
1656 groups, top, mdatoms,
1657 cr, nrnb, bFillGrid,
1660 wallcycle_stop(wcycle, ewcNS);
1663 if (inputrec->implicit_solvent && bNS)
1665 make_gb_nblist(cr, inputrec->gb_algorithm,
1666 x, box, fr, &top->idef, graph, fr->born);
1669 if (DOMAINDECOMP(cr))
1671 if (!(cr->duty & DUTY_PME))
1673 wallcycle_start(wcycle, ewcPPDURINGPME);
1674 dd_force_flop_start(cr->dd, nrnb);
1680 /* Enforced rotation has its own cycle counter that starts after the collective
1681 * coordinates have been communicated. It is added to ddCyclF to allow
1682 * for proper load-balancing */
1683 wallcycle_start(wcycle, ewcROT);
1684 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1685 wallcycle_stop(wcycle, ewcROT);
1688 /* Start the force cycle counter.
1689 * This counter is stopped in do_forcelow_level.
1690 * No parallel communication should occur while this counter is running,
1691 * since that will interfere with the dynamic load balancing.
1693 wallcycle_start(wcycle, ewcFORCE);
1697 /* Reset forces for which the virial is calculated separately:
1698 * PME/Ewald forces if necessary */
1699 if (fr->bF_NoVirSum)
1701 if (flags & GMX_FORCE_VIRIAL)
1703 fr->f_novirsum = fr->f_novirsum_alloc;
1706 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1710 clear_rvecs(homenr, fr->f_novirsum+start);
1715 /* We are not calculating the pressure so we do not need
1716 * a separate array for forces that do not contribute
1723 /* Clear the short- and long-range forces */
1724 clear_rvecs(fr->natoms_force_constr, f);
1725 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1727 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1730 clear_rvec(fr->vir_diag_posres);
1732 if (inputrec->ePull == epullCONSTRAINT)
1734 clear_pull_forces(inputrec->pull);
1737 /* update QMMMrec, if necessary */
1740 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1743 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1745 posres_wrapper(fplog, flags, bSepDVDL, inputrec, nrnb, top, box, x,
1749 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_FBPOSRES].nr > 0)
1751 /* Flat-bottomed position restraints always require full pbc */
1752 if (!(bStateChanged && bDoAdressWF))
1754 set_pbc(&pbc, inputrec->ePBC, box);
1756 v = fbposres(top->idef.il[F_FBPOSRES].nr, top->idef.il[F_FBPOSRES].iatoms,
1757 top->idef.iparams_fbposres,
1758 (const rvec*)x, fr->f_novirsum, fr->vir_diag_posres,
1759 inputrec->ePBC == epbcNONE ? NULL : &pbc,
1760 fr->rc_scaling, fr->ePBC, fr->posres_com);
1761 enerd->term[F_FBPOSRES] += v;
1762 inc_nrnb(nrnb, eNR_FBPOSRES, top->idef.il[F_FBPOSRES].nr/2);
1765 /* Compute the bonded and non-bonded energies and optionally forces */
1766 do_force_lowlevel(fplog, step, fr, inputrec, &(top->idef),
1767 cr, nrnb, wcycle, mdatoms,
1768 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, top, fr->born,
1769 &(top->atomtypes), bBornRadii, box,
1770 inputrec->fepvals, lambda,
1771 graph, &(top->excls), fr->mu_tot,
1777 if (do_per_step(step, inputrec->nstcalclr))
1779 /* Add the long range forces to the short range forces */
1780 for (i = 0; i < fr->natoms_force_constr; i++)
1782 rvec_add(fr->f_twin[i], f[i], f[i]);
1787 cycles_force = wallcycle_stop(wcycle, ewcFORCE);
1791 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1794 if (DOMAINDECOMP(cr))
1796 dd_force_flop_stop(cr->dd, nrnb);
1799 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1805 if (IR_ELEC_FIELD(*inputrec))
1807 /* Compute forces due to electric field */
1808 calc_f_el(MASTER(cr) ? field : NULL,
1809 start, homenr, mdatoms->chargeA, fr->f_novirsum,
1810 inputrec->ex, inputrec->et, t);
1813 if (bDoAdressWF && fr->adress_icor == eAdressICThermoForce)
1815 /* Compute thermodynamic force in hybrid AdResS region */
1816 adress_thermo_force(start, homenr, &(top->cgs), x, fr->f_novirsum, fr, mdatoms,
1817 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1820 /* Communicate the forces */
1823 wallcycle_start(wcycle, ewcMOVEF);
1824 if (DOMAINDECOMP(cr))
1826 dd_move_f(cr->dd, f, fr->fshift);
1827 /* Do we need to communicate the separate force array
1828 * for terms that do not contribute to the single sum virial?
1829 * Position restraints and electric fields do not introduce
1830 * inter-cg forces, only full electrostatics methods do.
1831 * When we do not calculate the virial, fr->f_novirsum = f,
1832 * so we have already communicated these forces.
1834 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1835 (flags & GMX_FORCE_VIRIAL))
1837 dd_move_f(cr->dd, fr->f_novirsum, NULL);
1841 /* We should not update the shift forces here,
1842 * since f_twin is already included in f.
1844 dd_move_f(cr->dd, fr->f_twin, NULL);
1849 pd_move_f(cr, f, nrnb);
1852 pd_move_f(cr, fr->f_twin, nrnb);
1855 wallcycle_stop(wcycle, ewcMOVEF);
1858 /* If we have NoVirSum forces, but we do not calculate the virial,
1859 * we sum fr->f_novirum=f later.
1861 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1863 wallcycle_start(wcycle, ewcVSITESPREAD);
1864 spread_vsite_f(vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1865 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1866 wallcycle_stop(wcycle, ewcVSITESPREAD);
1870 wallcycle_start(wcycle, ewcVSITESPREAD);
1871 spread_vsite_f(vsite, x, fr->f_twin, NULL, FALSE, NULL,
1873 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1874 wallcycle_stop(wcycle, ewcVSITESPREAD);
1878 if (flags & GMX_FORCE_VIRIAL)
1880 /* Calculation of the virial must be done after vsites! */
1881 calc_virial(mdatoms->start, mdatoms->homenr, x, f,
1882 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1886 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1888 pull_potential_wrapper(fplog, bSepDVDL, cr, inputrec, box, x,
1889 f, vir_force, mdatoms, enerd, lambda, t);
1892 /* Add the forces from enforced rotation potentials (if any) */
1895 wallcycle_start(wcycle, ewcROTadd);
1896 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
1897 wallcycle_stop(wcycle, ewcROTadd);
1900 if (PAR(cr) && !(cr->duty & DUTY_PME))
1902 /* In case of node-splitting, the PP nodes receive the long-range
1903 * forces, virial and energy from the PME nodes here.
1905 pme_receive_force_ener(fplog, bSepDVDL, cr, wcycle, enerd, fr);
1910 post_process_forces(cr, step, nrnb, wcycle,
1911 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
1915 /* Sum the potential energy terms from group contributions */
1916 sum_epot(&(enerd->grpp), enerd->term);
1919 void do_force(FILE *fplog, t_commrec *cr,
1920 t_inputrec *inputrec,
1921 gmx_large_int_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1922 gmx_localtop_t *top,
1923 gmx_groups_t *groups,
1924 matrix box, rvec x[], history_t *hist,
1928 gmx_enerdata_t *enerd, t_fcdata *fcd,
1929 real *lambda, t_graph *graph,
1931 gmx_vsite_t *vsite, rvec mu_tot,
1932 double t, FILE *field, gmx_edsam_t ed,
1933 gmx_bool bBornRadii,
1936 /* modify force flag if not doing nonbonded */
1937 if (!fr->bNonbonded)
1939 flags &= ~GMX_FORCE_NONBONDED;
1942 switch (inputrec->cutoff_scheme)
1945 do_force_cutsVERLET(fplog, cr, inputrec,
1961 do_force_cutsGROUP(fplog, cr, inputrec,
1976 gmx_incons("Invalid cut-off scheme passed!");
1981 void do_constrain_first(FILE *fplog, gmx_constr_t constr,
1982 t_inputrec *ir, t_mdatoms *md,
1983 t_state *state, t_commrec *cr, t_nrnb *nrnb,
1984 t_forcerec *fr, gmx_localtop_t *top)
1986 int i, m, start, end;
1987 gmx_large_int_t step;
1988 real dt = ir->delta_t;
1992 snew(savex, state->natoms);
1995 end = md->homenr + start;
1999 fprintf(debug, "vcm: start=%d, homenr=%d, end=%d\n",
2000 start, md->homenr, end);
2002 /* Do a first constrain to reset particles... */
2003 step = ir->init_step;
2006 char buf[STEPSTRSIZE];
2007 fprintf(fplog, "\nConstraining the starting coordinates (step %s)\n",
2008 gmx_step_str(step, buf));
2012 /* constrain the current position */
2013 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2014 ir, NULL, cr, step, 0, md,
2015 state->x, state->x, NULL,
2016 fr->bMolPBC, state->box,
2017 state->lambda[efptBONDED], &dvdl_dum,
2018 NULL, NULL, nrnb, econqCoord,
2019 ir->epc == epcMTTK, state->veta, state->veta);
2022 /* constrain the inital velocity, and save it */
2023 /* also may be useful if we need the ekin from the halfstep for velocity verlet */
2024 /* might not yet treat veta correctly */
2025 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2026 ir, NULL, cr, step, 0, md,
2027 state->x, state->v, state->v,
2028 fr->bMolPBC, state->box,
2029 state->lambda[efptBONDED], &dvdl_dum,
2030 NULL, NULL, nrnb, econqVeloc,
2031 ir->epc == epcMTTK, state->veta, state->veta);
2033 /* constrain the inital velocities at t-dt/2 */
2034 if (EI_STATE_VELOCITY(ir->eI) && ir->eI != eiVV)
2036 for (i = start; (i < end); i++)
2038 for (m = 0; (m < DIM); m++)
2040 /* Reverse the velocity */
2041 state->v[i][m] = -state->v[i][m];
2042 /* Store the position at t-dt in buf */
2043 savex[i][m] = state->x[i][m] + dt*state->v[i][m];
2046 /* Shake the positions at t=-dt with the positions at t=0
2047 * as reference coordinates.
2051 char buf[STEPSTRSIZE];
2052 fprintf(fplog, "\nConstraining the coordinates at t0-dt (step %s)\n",
2053 gmx_step_str(step, buf));
2056 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2057 ir, NULL, cr, step, -1, md,
2058 state->x, savex, NULL,
2059 fr->bMolPBC, state->box,
2060 state->lambda[efptBONDED], &dvdl_dum,
2061 state->v, NULL, nrnb, econqCoord,
2062 ir->epc == epcMTTK, state->veta, state->veta);
2064 for (i = start; i < end; i++)
2066 for (m = 0; m < DIM; m++)
2068 /* Re-reverse the velocities */
2069 state->v[i][m] = -state->v[i][m];
2076 void calc_enervirdiff(FILE *fplog, int eDispCorr, t_forcerec *fr)
2078 double eners[2], virs[2], enersum, virsum, y0, f, g, h;
2079 double r0, r1, r, rc3, rc9, ea, eb, ec, pa, pb, pc, pd;
2080 double invscale, invscale2, invscale3;
2081 int ri0, ri1, ri, i, offstart, offset;
2082 real scale, *vdwtab, tabfactor, tmp;
2084 fr->enershiftsix = 0;
2085 fr->enershifttwelve = 0;
2086 fr->enerdiffsix = 0;
2087 fr->enerdifftwelve = 0;
2089 fr->virdifftwelve = 0;
2091 if (eDispCorr != edispcNO)
2093 for (i = 0; i < 2; i++)
2098 if ((fr->vdwtype == evdwSWITCH) || (fr->vdwtype == evdwSHIFT))
2100 if (fr->rvdw_switch == 0)
2103 "With dispersion correction rvdw-switch can not be zero "
2104 "for vdw-type = %s", evdw_names[fr->vdwtype]);
2107 scale = fr->nblists[0].table_elec_vdw.scale;
2108 vdwtab = fr->nblists[0].table_vdw.data;
2110 /* Round the cut-offs to exact table values for precision */
2111 ri0 = floor(fr->rvdw_switch*scale);
2112 ri1 = ceil(fr->rvdw*scale);
2118 if (fr->vdwtype == evdwSHIFT)
2120 /* Determine the constant energy shift below rvdw_switch.
2121 * Table has a scale factor since we have scaled it down to compensate
2122 * for scaling-up c6/c12 with the derivative factors to save flops in analytical kernels.
2124 fr->enershiftsix = (real)(-1.0/(rc3*rc3)) - 6.0*vdwtab[8*ri0];
2125 fr->enershifttwelve = (real)( 1.0/(rc9*rc3)) - 12.0*vdwtab[8*ri0 + 4];
2127 /* Add the constant part from 0 to rvdw_switch.
2128 * This integration from 0 to rvdw_switch overcounts the number
2129 * of interactions by 1, as it also counts the self interaction.
2130 * We will correct for this later.
2132 eners[0] += 4.0*M_PI*fr->enershiftsix*rc3/3.0;
2133 eners[1] += 4.0*M_PI*fr->enershifttwelve*rc3/3.0;
2135 invscale = 1.0/(scale);
2136 invscale2 = invscale*invscale;
2137 invscale3 = invscale*invscale2;
2139 /* following summation derived from cubic spline definition,
2140 Numerical Recipies in C, second edition, p. 113-116. Exact
2141 for the cubic spline. We first calculate the negative of
2142 the energy from rvdw to rvdw_switch, assuming that g(r)=1,
2143 and then add the more standard, abrupt cutoff correction to
2144 that result, yielding the long-range correction for a
2145 switched function. We perform both the pressure and energy
2146 loops at the same time for simplicity, as the computational
2149 for (i = 0; i < 2; i++)
2151 enersum = 0.0; virsum = 0.0;
2155 /* Since the dispersion table has been scaled down a factor 6.0 and the repulsion
2156 * a factor 12.0 to compensate for the c6/c12 parameters inside nbfp[] being scaled
2157 * up (to save flops in kernels), we need to correct for this.
2166 for (ri = ri0; ri < ri1; ri++)
2170 eb = 2.0*invscale2*r;
2174 pb = 3.0*invscale2*r;
2175 pc = 3.0*invscale*r*r;
2178 /* this "8" is from the packing in the vdwtab array - perhaps should be #define'ed? */
2179 offset = 8*ri + offstart;
2180 y0 = vdwtab[offset];
2181 f = vdwtab[offset+1];
2182 g = vdwtab[offset+2];
2183 h = vdwtab[offset+3];
2185 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);
2186 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);
2189 enersum *= 4.0*M_PI*tabfactor;
2190 virsum *= 4.0*M_PI*tabfactor;
2191 eners[i] -= enersum;
2195 /* now add the correction for rvdw_switch to infinity */
2196 eners[0] += -4.0*M_PI/(3.0*rc3);
2197 eners[1] += 4.0*M_PI/(9.0*rc9);
2198 virs[0] += 8.0*M_PI/rc3;
2199 virs[1] += -16.0*M_PI/(3.0*rc9);
2201 else if ((fr->vdwtype == evdwCUT) || (fr->vdwtype == evdwUSER))
2203 if (fr->vdwtype == evdwUSER && fplog)
2206 "WARNING: using dispersion correction with user tables\n");
2208 rc3 = fr->rvdw*fr->rvdw*fr->rvdw;
2210 /* Contribution beyond the cut-off */
2211 eners[0] += -4.0*M_PI/(3.0*rc3);
2212 eners[1] += 4.0*M_PI/(9.0*rc9);
2213 if (fr->vdw_modifier == eintmodPOTSHIFT)
2215 /* Contribution within the cut-off */
2216 eners[0] += -4.0*M_PI/(3.0*rc3);
2217 eners[1] += 4.0*M_PI/(3.0*rc9);
2219 /* Contribution beyond the cut-off */
2220 virs[0] += 8.0*M_PI/rc3;
2221 virs[1] += -16.0*M_PI/(3.0*rc9);
2226 "Dispersion correction is not implemented for vdw-type = %s",
2227 evdw_names[fr->vdwtype]);
2229 fr->enerdiffsix = eners[0];
2230 fr->enerdifftwelve = eners[1];
2231 /* The 0.5 is due to the Gromacs definition of the virial */
2232 fr->virdiffsix = 0.5*virs[0];
2233 fr->virdifftwelve = 0.5*virs[1];
2237 void calc_dispcorr(FILE *fplog, t_inputrec *ir, t_forcerec *fr,
2238 gmx_large_int_t step, int natoms,
2239 matrix box, real lambda, tensor pres, tensor virial,
2240 real *prescorr, real *enercorr, real *dvdlcorr)
2242 gmx_bool bCorrAll, bCorrPres;
2243 real dvdlambda, invvol, dens, ninter, avcsix, avctwelve, enerdiff, svir = 0, spres = 0;
2253 if (ir->eDispCorr != edispcNO)
2255 bCorrAll = (ir->eDispCorr == edispcAllEner ||
2256 ir->eDispCorr == edispcAllEnerPres);
2257 bCorrPres = (ir->eDispCorr == edispcEnerPres ||
2258 ir->eDispCorr == edispcAllEnerPres);
2260 invvol = 1/det(box);
2263 /* Only correct for the interactions with the inserted molecule */
2264 dens = (natoms - fr->n_tpi)*invvol;
2269 dens = natoms*invvol;
2270 ninter = 0.5*natoms;
2273 if (ir->efep == efepNO)
2275 avcsix = fr->avcsix[0];
2276 avctwelve = fr->avctwelve[0];
2280 avcsix = (1 - lambda)*fr->avcsix[0] + lambda*fr->avcsix[1];
2281 avctwelve = (1 - lambda)*fr->avctwelve[0] + lambda*fr->avctwelve[1];
2284 enerdiff = ninter*(dens*fr->enerdiffsix - fr->enershiftsix);
2285 *enercorr += avcsix*enerdiff;
2287 if (ir->efep != efepNO)
2289 dvdlambda += (fr->avcsix[1] - fr->avcsix[0])*enerdiff;
2293 enerdiff = ninter*(dens*fr->enerdifftwelve - fr->enershifttwelve);
2294 *enercorr += avctwelve*enerdiff;
2295 if (fr->efep != efepNO)
2297 dvdlambda += (fr->avctwelve[1] - fr->avctwelve[0])*enerdiff;
2303 svir = ninter*dens*avcsix*fr->virdiffsix/3.0;
2304 if (ir->eDispCorr == edispcAllEnerPres)
2306 svir += ninter*dens*avctwelve*fr->virdifftwelve/3.0;
2308 /* The factor 2 is because of the Gromacs virial definition */
2309 spres = -2.0*invvol*svir*PRESFAC;
2311 for (m = 0; m < DIM; m++)
2313 virial[m][m] += svir;
2314 pres[m][m] += spres;
2319 /* Can't currently control when it prints, for now, just print when degugging */
2324 fprintf(debug, "Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
2330 "Long Range LJ corr.: Epot %10g, Pres: %10g, Vir: %10g\n",
2331 *enercorr, spres, svir);
2335 fprintf(debug, "Long Range LJ corr.: Epot %10g\n", *enercorr);
2339 if (fr->bSepDVDL && do_per_step(step, ir->nstlog))
2341 gmx_print_sepdvdl(fplog, "Dispersion correction", *enercorr, dvdlambda);
2343 if (fr->efep != efepNO)
2345 *dvdlcorr += dvdlambda;
2350 void do_pbc_first(FILE *fplog, matrix box, t_forcerec *fr,
2351 t_graph *graph, rvec x[])
2355 fprintf(fplog, "Removing pbc first time\n");
2357 calc_shifts(box, fr->shift_vec);
2360 mk_mshift(fplog, graph, fr->ePBC, box, x);
2363 p_graph(debug, "do_pbc_first 1", graph);
2365 shift_self(graph, box, x);
2366 /* By doing an extra mk_mshift the molecules that are broken
2367 * because they were e.g. imported from another software
2368 * will be made whole again. Such are the healing powers
2371 mk_mshift(fplog, graph, fr->ePBC, box, x);
2374 p_graph(debug, "do_pbc_first 2", graph);
2379 fprintf(fplog, "Done rmpbc\n");
2383 static void low_do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2384 gmx_mtop_t *mtop, rvec x[],
2389 gmx_molblock_t *molb;
2391 if (bFirst && fplog)
2393 fprintf(fplog, "Removing pbc first time\n");
2398 for (mb = 0; mb < mtop->nmolblock; mb++)
2400 molb = &mtop->molblock[mb];
2401 if (molb->natoms_mol == 1 ||
2402 (!bFirst && mtop->moltype[molb->type].cgs.nr == 1))
2404 /* Just one atom or charge group in the molecule, no PBC required */
2405 as += molb->nmol*molb->natoms_mol;
2409 /* Pass NULL iso fplog to avoid graph prints for each molecule type */
2410 mk_graph_ilist(NULL, mtop->moltype[molb->type].ilist,
2411 0, molb->natoms_mol, FALSE, FALSE, graph);
2413 for (mol = 0; mol < molb->nmol; mol++)
2415 mk_mshift(fplog, graph, ePBC, box, x+as);
2417 shift_self(graph, box, x+as);
2418 /* The molecule is whole now.
2419 * We don't need the second mk_mshift call as in do_pbc_first,
2420 * since we no longer need this graph.
2423 as += molb->natoms_mol;
2431 void do_pbc_first_mtop(FILE *fplog, int ePBC, matrix box,
2432 gmx_mtop_t *mtop, rvec x[])
2434 low_do_pbc_mtop(fplog, ePBC, box, mtop, x, TRUE);
2437 void do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2438 gmx_mtop_t *mtop, rvec x[])
2440 low_do_pbc_mtop(fplog, ePBC, box, mtop, x, FALSE);
2443 void finish_run(FILE *fplog, t_commrec *cr,
2444 t_inputrec *inputrec,
2445 t_nrnb nrnb[], gmx_wallcycle_t wcycle,
2446 gmx_runtime_t *runtime,
2447 wallclock_gpu_t *gputimes,
2448 gmx_bool bWriteStat)
2451 t_nrnb *nrnb_tot = NULL;
2454 double elapsed_run_time_over_all_ranks = 0;
2455 double elapsed_run_time_per_thread_over_all_ranks = 0;
2456 wallcycle_sum(cr, wcycle);
2462 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM,
2463 cr->mpi_comm_mysim);
2474 /* reduce elapsed_run_time over all MPI ranks in the current simulation */
2475 MPI_Allreduce(&runtime->elapsed_run_time,
2476 &elapsed_run_time_over_all_ranks,
2477 1, MPI_DOUBLE, MPI_SUM,
2478 cr->mpi_comm_mysim);
2479 elapsed_run_time_over_all_ranks /= cr->nnodes;
2480 /* reduce elapsed_run_time_per_thread over all MPI ranks in the current simulation */
2481 MPI_Allreduce(&runtime->elapsed_run_time_per_thread,
2482 &elapsed_run_time_per_thread_over_all_ranks,
2483 1, MPI_DOUBLE, MPI_SUM,
2484 cr->mpi_comm_mysim);
2489 elapsed_run_time_over_all_ranks = runtime->elapsed_run_time;
2490 elapsed_run_time_per_thread_over_all_ranks = runtime->elapsed_run_time_per_thread;
2495 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
2502 if ((cr->duty & DUTY_PP) && DOMAINDECOMP(cr))
2504 print_dd_statistics(cr, inputrec, fplog);
2516 snew(nrnb_all, cr->nnodes);
2517 nrnb_all[0] = *nrnb;
2518 for (s = 1; s < cr->nnodes; s++)
2520 MPI_Recv(nrnb_all[s].n, eNRNB, MPI_DOUBLE, s, 0,
2521 cr->mpi_comm_mysim, &stat);
2523 pr_load(fplog, cr, nrnb_all);
2528 MPI_Send(nrnb->n, eNRNB, MPI_DOUBLE, MASTERRANK(cr), 0,
2529 cr->mpi_comm_mysim);
2536 wallcycle_print(fplog, cr->nnodes, cr->npmenodes, runtime->elapsed_run_time,
2539 if (EI_DYNAMICS(inputrec->eI))
2541 delta_t = inputrec->delta_t;
2550 print_perf(fplog, elapsed_run_time_per_thread_over_all_ranks,
2551 elapsed_run_time_over_all_ranks,
2552 runtime->nsteps_done, delta_t, nbfs, mflop);
2556 print_perf(stderr, elapsed_run_time_per_thread_over_all_ranks,
2557 elapsed_run_time_over_all_ranks,
2558 runtime->nsteps_done, delta_t, nbfs, mflop);
2563 extern void initialize_lambdas(FILE *fplog, t_inputrec *ir, int *fep_state, real *lambda, double *lam0)
2565 /* this function works, but could probably use a logic rewrite to keep all the different
2566 types of efep straight. */
2569 t_lambda *fep = ir->fepvals;
2571 if ((ir->efep == efepNO) && (ir->bSimTemp == FALSE))
2573 for (i = 0; i < efptNR; i++)
2585 *fep_state = fep->init_fep_state; /* this might overwrite the checkpoint
2586 if checkpoint is set -- a kludge is in for now
2588 for (i = 0; i < efptNR; i++)
2590 /* overwrite lambda state with init_lambda for now for backwards compatibility */
2591 if (fep->init_lambda >= 0) /* if it's -1, it was never initializd */
2593 lambda[i] = fep->init_lambda;
2596 lam0[i] = lambda[i];
2601 lambda[i] = fep->all_lambda[i][*fep_state];
2604 lam0[i] = lambda[i];
2610 /* need to rescale control temperatures to match current state */
2611 for (i = 0; i < ir->opts.ngtc; i++)
2613 if (ir->opts.ref_t[i] > 0)
2615 ir->opts.ref_t[i] = ir->simtempvals->temperatures[*fep_state];
2621 /* Send to the log the information on the current lambdas */
2624 fprintf(fplog, "Initial vector of lambda components:[ ");
2625 for (i = 0; i < efptNR; i++)
2627 fprintf(fplog, "%10.4f ", lambda[i]);
2629 fprintf(fplog, "]\n");
2635 void init_md(FILE *fplog,
2636 t_commrec *cr, t_inputrec *ir, const output_env_t oenv,
2637 double *t, double *t0,
2638 real *lambda, int *fep_state, double *lam0,
2639 t_nrnb *nrnb, gmx_mtop_t *mtop,
2641 int nfile, const t_filenm fnm[],
2642 gmx_mdoutf_t **outf, t_mdebin **mdebin,
2643 tensor force_vir, tensor shake_vir, rvec mu_tot,
2644 gmx_bool *bSimAnn, t_vcm **vcm, unsigned long Flags)
2649 /* Initial values */
2650 *t = *t0 = ir->init_t;
2653 for (i = 0; i < ir->opts.ngtc; i++)
2655 /* set bSimAnn if any group is being annealed */
2656 if (ir->opts.annealing[i] != eannNO)
2663 update_annealing_target_temp(&(ir->opts), ir->init_t);
2666 /* Initialize lambda variables */
2667 initialize_lambdas(fplog, ir, fep_state, lambda, lam0);
2671 *upd = init_update(ir);
2677 *vcm = init_vcm(fplog, &mtop->groups, ir);
2680 if (EI_DYNAMICS(ir->eI) && !(Flags & MD_APPENDFILES))
2682 if (ir->etc == etcBERENDSEN)
2684 please_cite(fplog, "Berendsen84a");
2686 if (ir->etc == etcVRESCALE)
2688 please_cite(fplog, "Bussi2007a");
2696 *outf = init_mdoutf(nfile, fnm, Flags, cr, ir, oenv);
2698 *mdebin = init_mdebin((Flags & MD_APPENDFILES) ? NULL : (*outf)->fp_ene,
2699 mtop, ir, (*outf)->fp_dhdl);
2704 please_cite(fplog, "Fritsch12");
2705 please_cite(fplog, "Junghans10");
2707 /* Initiate variables */
2708 clear_mat(force_vir);
2709 clear_mat(shake_vir);