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"
105 typedef struct gmx_timeprint {
110 /* Portable version of ctime_r implemented in src/gmxlib/string2.c, but we do not want it declared in public installed headers */
112 gmx_ctime_r(const time_t *clock, char *buf, int n);
118 #ifdef HAVE_GETTIMEOFDAY
122 gettimeofday(&t, NULL);
124 seconds = (double) t.tv_sec + 1e-6*(double)t.tv_usec;
130 seconds = time(NULL);
137 #define difftime(end, start) ((double)(end)-(double)(start))
139 void print_time(FILE *out, gmx_runtime_t *runtime, gmx_large_int_t step,
140 t_inputrec *ir, t_commrec *cr)
143 char timebuf[STRLEN];
147 #ifndef GMX_THREAD_MPI
153 fprintf(out, "step %s", gmx_step_str(step, buf));
154 if ((step >= ir->nstlist))
156 runtime->last = gmx_gettime();
157 dt = difftime(runtime->last, runtime->real);
158 runtime->time_per_step = dt/(step - ir->init_step + 1);
160 dt = (ir->nsteps + ir->init_step - step)*runtime->time_per_step;
166 finish = (time_t) (runtime->last + dt);
167 gmx_ctime_r(&finish, timebuf, STRLEN);
168 sprintf(buf, "%s", timebuf);
169 buf[strlen(buf)-1] = '\0';
170 fprintf(out, ", will finish %s", buf);
174 fprintf(out, ", remaining runtime: %5d s ", (int)dt);
179 fprintf(out, " performance: %.1f ns/day ",
180 ir->delta_t/1000*24*60*60/runtime->time_per_step);
183 #ifndef GMX_THREAD_MPI
197 static double set_proctime(gmx_runtime_t *runtime)
203 prev = runtime->proc;
204 runtime->proc = dclock();
206 diff = runtime->proc - prev;
210 prev = runtime->proc;
211 runtime->proc = clock();
213 diff = (double)(runtime->proc - prev)/(double)CLOCKS_PER_SEC;
217 /* The counter has probably looped, ignore this data */
224 void runtime_start(gmx_runtime_t *runtime)
226 runtime->real = gmx_gettime();
228 set_proctime(runtime);
229 runtime->realtime = 0;
230 runtime->proctime = 0;
232 runtime->time_per_step = 0;
235 void runtime_end(gmx_runtime_t *runtime)
241 runtime->proctime += set_proctime(runtime);
242 runtime->realtime = now - runtime->real;
246 void runtime_upd_proc(gmx_runtime_t *runtime)
248 runtime->proctime += set_proctime(runtime);
251 void print_date_and_time(FILE *fplog, int nodeid, const char *title,
252 const gmx_runtime_t *runtime)
255 char timebuf[STRLEN];
256 char time_string[STRLEN];
263 tmptime = (time_t) runtime->real;
264 gmx_ctime_r(&tmptime, timebuf, STRLEN);
268 tmptime = (time_t) gmx_gettime();
269 gmx_ctime_r(&tmptime, timebuf, STRLEN);
271 for (i = 0; timebuf[i] >= ' '; i++)
273 time_string[i] = timebuf[i];
275 time_string[i] = '\0';
277 fprintf(fplog, "%s on node %d %s\n", title, nodeid, time_string);
281 static void sum_forces(int start, int end, rvec f[], rvec flr[])
287 pr_rvecs(debug, 0, "fsr", f+start, end-start);
288 pr_rvecs(debug, 0, "flr", flr+start, end-start);
290 for (i = start; (i < end); i++)
292 rvec_inc(f[i], flr[i]);
297 * calc_f_el calculates forces due to an electric field.
299 * force is kJ mol^-1 nm^-1 = e * kJ mol^-1 nm^-1 / e
301 * Et[] contains the parameters for the time dependent
302 * part of the field (not yet used).
303 * Ex[] contains the parameters for
304 * the spatial dependent part of the field. You can have cool periodic
305 * fields in principle, but only a constant field is supported
307 * The function should return the energy due to the electric field
308 * (if any) but for now returns 0.
311 * There can be problems with the virial.
312 * Since the field is not self-consistent this is unavoidable.
313 * For neutral molecules the virial is correct within this approximation.
314 * For neutral systems with many charged molecules the error is small.
315 * But for systems with a net charge or a few charged molecules
316 * the error can be significant when the field is high.
317 * Solution: implement a self-consitent electric field into PME.
319 static void calc_f_el(FILE *fp, int start, int homenr,
320 real charge[], rvec x[], rvec f[],
321 t_cosines Ex[], t_cosines Et[], double t)
327 for (m = 0; (m < DIM); m++)
334 Ext[m] = cos(Et[m].a[0]*(t-t0))*exp(-sqr(t-t0)/(2.0*sqr(Et[m].a[2])));
338 Ext[m] = cos(Et[m].a[0]*t);
347 /* Convert the field strength from V/nm to MD-units */
348 Ext[m] *= Ex[m].a[0]*FIELDFAC;
349 for (i = start; (i < start+homenr); i++)
351 f[i][m] += charge[i]*Ext[m];
361 fprintf(fp, "%10g %10g %10g %10g #FIELD\n", t,
362 Ext[XX]/FIELDFAC, Ext[YY]/FIELDFAC, Ext[ZZ]/FIELDFAC);
366 static void calc_virial(FILE *fplog, int start, int homenr, rvec x[], rvec f[],
367 tensor vir_part, t_graph *graph, matrix box,
368 t_nrnb *nrnb, const t_forcerec *fr, int ePBC)
373 /* The short-range virial from surrounding boxes */
375 calc_vir(fplog, SHIFTS, fr->shift_vec, fr->fshift, vir_part, ePBC == epbcSCREW, box);
376 inc_nrnb(nrnb, eNR_VIRIAL, SHIFTS);
378 /* Calculate partial virial, for local atoms only, based on short range.
379 * Total virial is computed in global_stat, called from do_md
381 f_calc_vir(fplog, start, start+homenr, x, f, vir_part, graph, box);
382 inc_nrnb(nrnb, eNR_VIRIAL, homenr);
384 /* Add position restraint contribution */
385 for (i = 0; i < DIM; i++)
387 vir_part[i][i] += fr->vir_diag_posres[i];
390 /* Add wall contribution */
391 for (i = 0; i < DIM; i++)
393 vir_part[i][ZZ] += fr->vir_wall_z[i];
398 pr_rvecs(debug, 0, "vir_part", vir_part, DIM);
402 static void posres_wrapper(FILE *fplog,
408 matrix box, rvec x[],
410 gmx_enerdata_t *enerd,
418 /* Position restraints always require full pbc */
419 set_pbc(&pbc, ir->ePBC, box);
421 v = posres(top->idef.il[F_POSRES].nr, top->idef.il[F_POSRES].iatoms,
422 top->idef.iparams_posres,
423 (const rvec*)x, fr->f_novirsum, fr->vir_diag_posres,
424 ir->ePBC == epbcNONE ? NULL : &pbc,
425 lambda[efptRESTRAINT], &dvdl,
426 fr->rc_scaling, fr->ePBC, fr->posres_com, fr->posres_comB);
429 gmx_print_sepdvdl(fplog, interaction_function[F_POSRES].longname, v, dvdl);
431 enerd->term[F_POSRES] += v;
432 /* If just the force constant changes, the FEP term is linear,
433 * but if k changes, it is not.
435 enerd->dvdl_nonlin[efptRESTRAINT] += dvdl;
436 inc_nrnb(nrnb, eNR_POSRES, top->idef.il[F_POSRES].nr/2);
438 if ((ir->fepvals->n_lambda > 0) && (flags & GMX_FORCE_DHDL))
440 for (i = 0; i < enerd->n_lambda; i++)
442 real dvdl_dum, lambda_dum;
444 lambda_dum = (i == 0 ? lambda[efptRESTRAINT] : ir->fepvals->all_lambda[efptRESTRAINT][i-1]);
445 v = posres(top->idef.il[F_POSRES].nr, top->idef.il[F_POSRES].iatoms,
446 top->idef.iparams_posres,
447 (const rvec*)x, NULL, NULL,
448 ir->ePBC == epbcNONE ? NULL : &pbc, lambda_dum, &dvdl,
449 fr->rc_scaling, fr->ePBC, fr->posres_com, fr->posres_comB);
450 enerd->enerpart_lambda[i] += v;
455 static void pull_potential_wrapper(FILE *fplog,
459 matrix box, rvec x[],
463 gmx_enerdata_t *enerd,
470 /* Calculate the center of mass forces, this requires communication,
471 * which is why pull_potential is called close to other communication.
472 * The virial contribution is calculated directly,
473 * which is why we call pull_potential after calc_virial.
475 set_pbc(&pbc, ir->ePBC, box);
477 enerd->term[F_COM_PULL] +=
478 pull_potential(ir->ePull, ir->pull, mdatoms, &pbc,
479 cr, t, lambda[efptRESTRAINT], x, f, vir_force, &dvdl);
482 gmx_print_sepdvdl(fplog, "Com pull", enerd->term[F_COM_PULL], dvdl);
484 enerd->dvdl_lin[efptRESTRAINT] += dvdl;
487 static void pme_receive_force_ener(FILE *fplog,
490 gmx_wallcycle_t wcycle,
491 gmx_enerdata_t *enerd,
495 float cycles_ppdpme, cycles_seppme;
497 cycles_ppdpme = wallcycle_stop(wcycle, ewcPPDURINGPME);
498 dd_cycles_add(cr->dd, cycles_ppdpme, ddCyclPPduringPME);
500 /* In case of node-splitting, the PP nodes receive the long-range
501 * forces, virial and energy from the PME nodes here.
503 wallcycle_start(wcycle, ewcPP_PMEWAITRECVF);
505 gmx_pme_receive_f(cr, fr->f_novirsum, fr->vir_el_recip, &e, &dvdl,
509 gmx_print_sepdvdl(fplog, "PME mesh", e, dvdl);
511 enerd->term[F_COUL_RECIP] += e;
512 enerd->dvdl_lin[efptCOUL] += dvdl;
515 dd_cycles_add(cr->dd, cycles_seppme, ddCyclPME);
517 wallcycle_stop(wcycle, ewcPP_PMEWAITRECVF);
520 static void print_large_forces(FILE *fp, t_mdatoms *md, t_commrec *cr,
521 gmx_large_int_t step, real pforce, rvec *x, rvec *f)
525 char buf[STEPSTRSIZE];
528 for (i = md->start; i < md->start+md->homenr; i++)
531 /* We also catch NAN, if the compiler does not optimize this away. */
532 if (fn2 >= pf2 || fn2 != fn2)
534 fprintf(fp, "step %s atom %6d x %8.3f %8.3f %8.3f force %12.5e\n",
535 gmx_step_str(step, buf),
536 ddglatnr(cr->dd, i), x[i][XX], x[i][YY], x[i][ZZ], sqrt(fn2));
541 static void post_process_forces(FILE *fplog,
543 gmx_large_int_t step,
544 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
546 matrix box, rvec x[],
551 t_forcerec *fr, gmx_vsite_t *vsite,
558 /* Spread the mesh force on virtual sites to the other particles...
559 * This is parallellized. MPI communication is performed
560 * if the constructing atoms aren't local.
562 wallcycle_start(wcycle, ewcVSITESPREAD);
563 spread_vsite_f(fplog, vsite, x, fr->f_novirsum, NULL,
564 (flags & GMX_FORCE_VIRIAL), fr->vir_el_recip,
566 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
567 wallcycle_stop(wcycle, ewcVSITESPREAD);
569 if (flags & GMX_FORCE_VIRIAL)
571 /* Now add the forces, this is local */
574 sum_forces(0, fr->f_novirsum_n, f, fr->f_novirsum);
578 sum_forces(mdatoms->start, mdatoms->start+mdatoms->homenr,
581 if (EEL_FULL(fr->eeltype))
583 /* Add the mesh contribution to the virial */
584 m_add(vir_force, fr->vir_el_recip, vir_force);
588 pr_rvecs(debug, 0, "vir_force", vir_force, DIM);
593 if (fr->print_force >= 0)
595 print_large_forces(stderr, mdatoms, cr, step, fr->print_force, x, f);
599 static void do_nb_verlet(t_forcerec *fr,
600 interaction_const_t *ic,
601 gmx_enerdata_t *enerd,
602 int flags, int ilocality,
605 gmx_wallcycle_t wcycle)
607 int nnbl, kernel_type, enr_nbnxn_kernel_ljc, enr_nbnxn_kernel_lj;
609 nonbonded_verlet_group_t *nbvg;
611 if (!(flags & GMX_FORCE_NONBONDED))
613 /* skip non-bonded calculation */
617 nbvg = &fr->nbv->grp[ilocality];
619 /* CUDA kernel launch overhead is already timed separately */
620 if (fr->cutoff_scheme != ecutsVERLET)
622 gmx_incons("Invalid cut-off scheme passed!");
625 if (nbvg->kernel_type != nbnxnk8x8x8_CUDA)
627 wallcycle_sub_start(wcycle, ewcsNONBONDED);
629 switch (nbvg->kernel_type)
631 case nbnxnk4x4_PlainC:
632 nbnxn_kernel_ref(&nbvg->nbl_lists,
638 enerd->grpp.ener[egCOULSR],
640 enerd->grpp.ener[egBHAMSR] :
641 enerd->grpp.ener[egLJSR]);
644 case nbnxnk4xN_SIMD_4xN:
645 nbnxn_kernel_simd_4xn(&nbvg->nbl_lists,
652 enerd->grpp.ener[egCOULSR],
654 enerd->grpp.ener[egBHAMSR] :
655 enerd->grpp.ener[egLJSR]);
657 case nbnxnk4xN_SIMD_2xNN:
658 nbnxn_kernel_simd_2xnn(&nbvg->nbl_lists,
665 enerd->grpp.ener[egCOULSR],
667 enerd->grpp.ener[egBHAMSR] :
668 enerd->grpp.ener[egLJSR]);
671 case nbnxnk8x8x8_CUDA:
672 nbnxn_cuda_launch_kernel(fr->nbv->cu_nbv, nbvg->nbat, flags, ilocality);
675 case nbnxnk8x8x8_PlainC:
676 nbnxn_kernel_gpu_ref(nbvg->nbl_lists.nbl[0],
681 nbvg->nbat->out[0].f,
683 enerd->grpp.ener[egCOULSR],
685 enerd->grpp.ener[egBHAMSR] :
686 enerd->grpp.ener[egLJSR]);
690 gmx_incons("Invalid nonbonded kernel type passed!");
693 if (nbvg->kernel_type != nbnxnk8x8x8_CUDA)
695 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
698 if (EEL_RF(ic->eeltype) || ic->eeltype == eelCUT)
700 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_RF;
702 else if (nbvg->ewald_excl == ewaldexclTable)
704 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_TAB;
708 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_EWALD;
710 enr_nbnxn_kernel_lj = eNR_NBNXN_LJ;
711 if (flags & GMX_FORCE_ENERGY)
713 /* In eNR_??? the nbnxn F+E kernels are always the F kernel + 1 */
714 enr_nbnxn_kernel_ljc += 1;
715 enr_nbnxn_kernel_lj += 1;
718 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc,
719 nbvg->nbl_lists.natpair_ljq);
720 inc_nrnb(nrnb, enr_nbnxn_kernel_lj,
721 nbvg->nbl_lists.natpair_lj);
722 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc-eNR_NBNXN_LJ_RF+eNR_NBNXN_RF,
723 nbvg->nbl_lists.natpair_q);
726 void do_force_cutsVERLET(FILE *fplog, t_commrec *cr,
727 t_inputrec *inputrec,
728 gmx_large_int_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
731 gmx_groups_t *groups,
732 matrix box, rvec x[], history_t *hist,
736 gmx_enerdata_t *enerd, t_fcdata *fcd,
737 real *lambda, t_graph *graph,
738 t_forcerec *fr, interaction_const_t *ic,
739 gmx_vsite_t *vsite, rvec mu_tot,
740 double t, FILE *field, gmx_edsam_t ed,
748 gmx_bool bSepDVDL, bStateChanged, bNS, bFillGrid, bCalcCGCM, bBS;
749 gmx_bool bDoLongRange, bDoForces, bSepLRF, bUseGPU, bUseOrEmulGPU;
750 gmx_bool bDiffKernels = FALSE;
752 rvec vzero, box_diag;
754 float cycles_pme, cycles_force;
755 nonbonded_verlet_t *nbv;
759 nb_kernel_type = fr->nbv->grp[0].kernel_type;
761 start = mdatoms->start;
762 homenr = mdatoms->homenr;
764 bSepDVDL = (fr->bSepDVDL && do_per_step(step, inputrec->nstlog));
766 clear_mat(vir_force);
769 if (DOMAINDECOMP(cr))
771 cg1 = cr->dd->ncg_tot;
782 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
783 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
784 bFillGrid = (bNS && bStateChanged);
785 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
786 bDoLongRange = (fr->bTwinRange && bNS && (flags & GMX_FORCE_DO_LR));
787 bDoForces = (flags & GMX_FORCE_FORCES);
788 bSepLRF = (bDoLongRange && bDoForces && (flags & GMX_FORCE_SEPLRF));
789 bUseGPU = fr->nbv->bUseGPU;
790 bUseOrEmulGPU = bUseGPU || (nbv->grp[0].kernel_type == nbnxnk8x8x8_PlainC);
794 update_forcerec(fplog, fr, box);
796 if (NEED_MUTOT(*inputrec))
798 /* Calculate total (local) dipole moment in a temporary common array.
799 * This makes it possible to sum them over nodes faster.
801 calc_mu(start, homenr,
802 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
807 if (fr->ePBC != epbcNONE)
809 /* Compute shift vectors every step,
810 * because of pressure coupling or box deformation!
812 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
814 calc_shifts(box, fr->shift_vec);
819 put_atoms_in_box_omp(fr->ePBC, box, homenr, x);
820 inc_nrnb(nrnb, eNR_SHIFTX, homenr);
822 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
824 unshift_self(graph, box, x);
828 nbnxn_atomdata_copy_shiftvec(flags & GMX_FORCE_DYNAMICBOX,
829 fr->shift_vec, nbv->grp[0].nbat);
832 if (!(cr->duty & DUTY_PME))
834 /* Send particle coordinates to the pme nodes.
835 * Since this is only implemented for domain decomposition
836 * and domain decomposition does not use the graph,
837 * we do not need to worry about shifting.
840 wallcycle_start(wcycle, ewcPP_PMESENDX);
842 bBS = (inputrec->nwall == 2);
846 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
849 gmx_pme_send_x(cr, bBS ? boxs : box, x,
850 mdatoms->nChargePerturbed, lambda[efptCOUL],
851 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)), step);
853 wallcycle_stop(wcycle, ewcPP_PMESENDX);
857 /* do gridding for pair search */
860 if (graph && bStateChanged)
862 /* Calculate intramolecular shift vectors to make molecules whole */
863 mk_mshift(fplog, graph, fr->ePBC, box, x);
867 box_diag[XX] = box[XX][XX];
868 box_diag[YY] = box[YY][YY];
869 box_diag[ZZ] = box[ZZ][ZZ];
871 wallcycle_start(wcycle, ewcNS);
874 wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
875 nbnxn_put_on_grid(nbv->nbs, fr->ePBC, box,
877 0, mdatoms->homenr, -1, fr->cginfo, x,
879 nbv->grp[eintLocal].kernel_type,
880 nbv->grp[eintLocal].nbat);
881 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
885 wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
886 nbnxn_put_on_grid_nonlocal(nbv->nbs, domdec_zones(cr->dd),
888 nbv->grp[eintNonlocal].kernel_type,
889 nbv->grp[eintNonlocal].nbat);
890 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
893 if (nbv->ngrp == 1 ||
894 nbv->grp[eintNonlocal].nbat == nbv->grp[eintLocal].nbat)
896 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatAll,
897 nbv->nbs, mdatoms, fr->cginfo);
901 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatLocal,
902 nbv->nbs, mdatoms, fr->cginfo);
903 nbnxn_atomdata_set(nbv->grp[eintNonlocal].nbat, eatAll,
904 nbv->nbs, mdatoms, fr->cginfo);
906 wallcycle_stop(wcycle, ewcNS);
909 /* initialize the GPU atom data and copy shift vector */
914 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
915 nbnxn_cuda_init_atomdata(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
916 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
919 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
920 nbnxn_cuda_upload_shiftvec(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
921 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
924 /* do local pair search */
927 wallcycle_start_nocount(wcycle, ewcNS);
928 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
929 nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintLocal].nbat,
932 nbv->min_ci_balanced,
933 &nbv->grp[eintLocal].nbl_lists,
935 nbv->grp[eintLocal].kernel_type,
937 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
941 /* initialize local pair-list on the GPU */
942 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
943 nbv->grp[eintLocal].nbl_lists.nbl[0],
946 wallcycle_stop(wcycle, ewcNS);
950 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
951 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
952 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, FALSE, x,
953 nbv->grp[eintLocal].nbat);
954 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
955 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
960 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
961 /* launch local nonbonded F on GPU */
962 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFNo,
964 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
967 /* Communicate coordinates and sum dipole if necessary +
968 do non-local pair search */
969 if (DOMAINDECOMP(cr))
971 bDiffKernels = (nbv->grp[eintNonlocal].kernel_type !=
972 nbv->grp[eintLocal].kernel_type);
976 /* With GPU+CPU non-bonded calculations we need to copy
977 * the local coordinates to the non-local nbat struct
978 * (in CPU format) as the non-local kernel call also
979 * calculates the local - non-local interactions.
981 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
982 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
983 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, TRUE, x,
984 nbv->grp[eintNonlocal].nbat);
985 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
986 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
991 wallcycle_start_nocount(wcycle, ewcNS);
992 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
996 nbnxn_grid_add_simple(nbv->nbs, nbv->grp[eintNonlocal].nbat);
999 nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintNonlocal].nbat,
1002 nbv->min_ci_balanced,
1003 &nbv->grp[eintNonlocal].nbl_lists,
1005 nbv->grp[eintNonlocal].kernel_type,
1008 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
1010 if (nbv->grp[eintNonlocal].kernel_type == nbnxnk8x8x8_CUDA)
1012 /* initialize non-local pair-list on the GPU */
1013 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
1014 nbv->grp[eintNonlocal].nbl_lists.nbl[0],
1017 wallcycle_stop(wcycle, ewcNS);
1021 wallcycle_start(wcycle, ewcMOVEX);
1022 dd_move_x(cr->dd, box, x);
1024 /* When we don't need the total dipole we sum it in global_stat */
1025 if (bStateChanged && NEED_MUTOT(*inputrec))
1027 gmx_sumd(2*DIM, mu, cr);
1029 wallcycle_stop(wcycle, ewcMOVEX);
1031 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1032 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
1033 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatNonlocal, FALSE, x,
1034 nbv->grp[eintNonlocal].nbat);
1035 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
1036 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1039 if (bUseGPU && !bDiffKernels)
1041 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
1042 /* launch non-local nonbonded F on GPU */
1043 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFNo,
1045 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1051 /* launch D2H copy-back F */
1052 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1053 if (DOMAINDECOMP(cr) && !bDiffKernels)
1055 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintNonlocal].nbat,
1056 flags, eatNonlocal);
1058 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintLocal].nbat,
1060 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1063 if (bStateChanged && NEED_MUTOT(*inputrec))
1067 gmx_sumd(2*DIM, mu, cr);
1070 for (i = 0; i < 2; i++)
1072 for (j = 0; j < DIM; j++)
1074 fr->mu_tot[i][j] = mu[i*DIM + j];
1078 if (fr->efep == efepNO)
1080 copy_rvec(fr->mu_tot[0], mu_tot);
1084 for (j = 0; j < DIM; j++)
1087 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] +
1088 lambda[efptCOUL]*fr->mu_tot[1][j];
1092 /* Reset energies */
1093 reset_enerdata(&(inputrec->opts), fr, bNS, enerd, MASTER(cr));
1094 clear_rvecs(SHIFTS, fr->fshift);
1096 if (DOMAINDECOMP(cr))
1098 if (!(cr->duty & DUTY_PME))
1100 wallcycle_start(wcycle, ewcPPDURINGPME);
1101 dd_force_flop_start(cr->dd, nrnb);
1107 /* Enforced rotation has its own cycle counter that starts after the collective
1108 * coordinates have been communicated. It is added to ddCyclF to allow
1109 * for proper load-balancing */
1110 wallcycle_start(wcycle, ewcROT);
1111 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1112 wallcycle_stop(wcycle, ewcROT);
1115 /* Start the force cycle counter.
1116 * This counter is stopped in do_forcelow_level.
1117 * No parallel communication should occur while this counter is running,
1118 * since that will interfere with the dynamic load balancing.
1120 wallcycle_start(wcycle, ewcFORCE);
1123 /* Reset forces for which the virial is calculated separately:
1124 * PME/Ewald forces if necessary */
1125 if (fr->bF_NoVirSum)
1127 if (flags & GMX_FORCE_VIRIAL)
1129 fr->f_novirsum = fr->f_novirsum_alloc;
1132 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1136 clear_rvecs(homenr, fr->f_novirsum+start);
1141 /* We are not calculating the pressure so we do not need
1142 * a separate array for forces that do not contribute
1149 /* Clear the short- and long-range forces */
1150 clear_rvecs(fr->natoms_force_constr, f);
1151 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1153 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1156 clear_rvec(fr->vir_diag_posres);
1159 if (inputrec->ePull == epullCONSTRAINT)
1161 clear_pull_forces(inputrec->pull);
1164 /* We calculate the non-bonded forces, when done on the CPU, here.
1165 * We do this before calling do_force_lowlevel, as in there bondeds
1166 * forces are calculated before PME, which does communication.
1167 * With this order, non-bonded and bonded force calculation imbalance
1168 * can be balanced out by the domain decomposition load balancing.
1173 /* Maybe we should move this into do_force_lowlevel */
1174 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFYes,
1178 if (!bUseOrEmulGPU || bDiffKernels)
1182 if (DOMAINDECOMP(cr))
1184 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal,
1185 bDiffKernels ? enbvClearFYes : enbvClearFNo,
1195 aloc = eintNonlocal;
1198 /* Add all the non-bonded force to the normal force array.
1199 * This can be split into a local a non-local part when overlapping
1200 * communication with calculation with domain decomposition.
1202 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1203 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1204 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1205 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatAll, nbv->grp[aloc].nbat, f);
1206 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1207 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1208 wallcycle_start_nocount(wcycle, ewcFORCE);
1210 /* if there are multiple fshift output buffers reduce them */
1211 if ((flags & GMX_FORCE_VIRIAL) &&
1212 nbv->grp[aloc].nbl_lists.nnbl > 1)
1214 nbnxn_atomdata_add_nbat_fshift_to_fshift(nbv->grp[aloc].nbat,
1219 /* update QMMMrec, if necessary */
1222 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1225 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1227 posres_wrapper(fplog, flags, bSepDVDL, inputrec, nrnb, top, box, x,
1228 f, enerd, lambda, fr);
1231 /* Compute the bonded and non-bonded energies and optionally forces */
1232 do_force_lowlevel(fplog, step, fr, inputrec, &(top->idef),
1233 cr, nrnb, wcycle, mdatoms, &(inputrec->opts),
1234 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, mtop, top, fr->born,
1235 &(top->atomtypes), bBornRadii, box,
1236 inputrec->fepvals, lambda, graph, &(top->excls), fr->mu_tot,
1237 flags, &cycles_pme);
1241 if (do_per_step(step, inputrec->nstcalclr))
1243 /* Add the long range forces to the short range forces */
1244 for (i = 0; i < fr->natoms_force_constr; i++)
1246 rvec_add(fr->f_twin[i], f[i], f[i]);
1251 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1255 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1258 if (bUseOrEmulGPU && !bDiffKernels)
1260 /* wait for non-local forces (or calculate in emulation mode) */
1261 if (DOMAINDECOMP(cr))
1265 wallcycle_start(wcycle, ewcWAIT_GPU_NB_NL);
1266 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1267 nbv->grp[eintNonlocal].nbat,
1269 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1271 cycles_force += wallcycle_stop(wcycle, ewcWAIT_GPU_NB_NL);
1275 wallcycle_start_nocount(wcycle, ewcFORCE);
1276 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFYes,
1278 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1280 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1281 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1282 /* skip the reduction if there was no non-local work to do */
1283 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1285 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatNonlocal,
1286 nbv->grp[eintNonlocal].nbat, f);
1288 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1289 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1295 /* Communicate the forces */
1298 wallcycle_start(wcycle, ewcMOVEF);
1299 if (DOMAINDECOMP(cr))
1301 dd_move_f(cr->dd, f, fr->fshift);
1302 /* Do we need to communicate the separate force array
1303 * for terms that do not contribute to the single sum virial?
1304 * Position restraints and electric fields do not introduce
1305 * inter-cg forces, only full electrostatics methods do.
1306 * When we do not calculate the virial, fr->f_novirsum = f,
1307 * so we have already communicated these forces.
1309 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1310 (flags & GMX_FORCE_VIRIAL))
1312 dd_move_f(cr->dd, fr->f_novirsum, NULL);
1316 /* We should not update the shift forces here,
1317 * since f_twin is already included in f.
1319 dd_move_f(cr->dd, fr->f_twin, NULL);
1322 wallcycle_stop(wcycle, ewcMOVEF);
1328 /* wait for local forces (or calculate in emulation mode) */
1331 wallcycle_start(wcycle, ewcWAIT_GPU_NB_L);
1332 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1333 nbv->grp[eintLocal].nbat,
1335 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1337 wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
1339 /* now clear the GPU outputs while we finish the step on the CPU */
1341 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1342 nbnxn_cuda_clear_outputs(nbv->cu_nbv, flags);
1343 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1347 wallcycle_start_nocount(wcycle, ewcFORCE);
1348 do_nb_verlet(fr, ic, enerd, flags, eintLocal,
1349 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
1351 wallcycle_stop(wcycle, ewcFORCE);
1353 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1354 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1355 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1357 /* skip the reduction if there was no non-local work to do */
1358 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatLocal,
1359 nbv->grp[eintLocal].nbat, f);
1361 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1362 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1365 if (DOMAINDECOMP(cr))
1367 dd_force_flop_stop(cr->dd, nrnb);
1370 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1376 if (IR_ELEC_FIELD(*inputrec))
1378 /* Compute forces due to electric field */
1379 calc_f_el(MASTER(cr) ? field : NULL,
1380 start, homenr, mdatoms->chargeA, x, fr->f_novirsum,
1381 inputrec->ex, inputrec->et, t);
1384 /* If we have NoVirSum forces, but we do not calculate the virial,
1385 * we sum fr->f_novirum=f later.
1387 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1389 wallcycle_start(wcycle, ewcVSITESPREAD);
1390 spread_vsite_f(fplog, vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1391 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1392 wallcycle_stop(wcycle, ewcVSITESPREAD);
1396 wallcycle_start(wcycle, ewcVSITESPREAD);
1397 spread_vsite_f(fplog, vsite, x, fr->f_twin, NULL, FALSE, NULL,
1399 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1400 wallcycle_stop(wcycle, ewcVSITESPREAD);
1404 if (flags & GMX_FORCE_VIRIAL)
1406 /* Calculation of the virial must be done after vsites! */
1407 calc_virial(fplog, mdatoms->start, mdatoms->homenr, x, f,
1408 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1412 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1414 pull_potential_wrapper(fplog, bSepDVDL, cr, inputrec, box, x,
1415 f, vir_force, mdatoms, enerd, lambda, t);
1418 /* Add the forces from enforced rotation potentials (if any) */
1421 wallcycle_start(wcycle, ewcROTadd);
1422 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
1423 wallcycle_stop(wcycle, ewcROTadd);
1426 if (PAR(cr) && !(cr->duty & DUTY_PME))
1428 /* In case of node-splitting, the PP nodes receive the long-range
1429 * forces, virial and energy from the PME nodes here.
1431 pme_receive_force_ener(fplog, bSepDVDL, cr, wcycle, enerd, fr);
1436 post_process_forces(fplog, cr, step, nrnb, wcycle,
1437 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
1441 /* Sum the potential energy terms from group contributions */
1442 sum_epot(&(inputrec->opts), &(enerd->grpp), enerd->term);
1445 void do_force_cutsGROUP(FILE *fplog, t_commrec *cr,
1446 t_inputrec *inputrec,
1447 gmx_large_int_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1448 gmx_localtop_t *top,
1450 gmx_groups_t *groups,
1451 matrix box, rvec x[], history_t *hist,
1455 gmx_enerdata_t *enerd, t_fcdata *fcd,
1456 real *lambda, t_graph *graph,
1457 t_forcerec *fr, gmx_vsite_t *vsite, rvec mu_tot,
1458 double t, FILE *field, gmx_edsam_t ed,
1459 gmx_bool bBornRadii,
1465 gmx_bool bSepDVDL, bStateChanged, bNS, bFillGrid, bCalcCGCM, bBS;
1466 gmx_bool bDoLongRangeNS, bDoForces, bDoPotential, bSepLRF;
1467 gmx_bool bDoAdressWF;
1469 rvec vzero, box_diag;
1470 real e, v, dvdlambda[efptNR];
1472 float cycles_pme, cycles_force;
1474 start = mdatoms->start;
1475 homenr = mdatoms->homenr;
1477 bSepDVDL = (fr->bSepDVDL && do_per_step(step, inputrec->nstlog));
1479 clear_mat(vir_force);
1483 pd_cg_range(cr, &cg0, &cg1);
1488 if (DOMAINDECOMP(cr))
1490 cg1 = cr->dd->ncg_tot;
1502 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
1503 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
1504 /* Should we update the long-range neighborlists at this step? */
1505 bDoLongRangeNS = fr->bTwinRange && bNS;
1506 /* Should we perform the long-range nonbonded evaluation inside the neighborsearching? */
1507 bFillGrid = (bNS && bStateChanged);
1508 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
1509 bDoForces = (flags & GMX_FORCE_FORCES);
1510 bDoPotential = (flags & GMX_FORCE_ENERGY);
1511 bSepLRF = ((inputrec->nstcalclr > 1) && bDoForces &&
1512 (flags & GMX_FORCE_SEPLRF) && (flags & GMX_FORCE_DO_LR));
1514 /* should probably move this to the forcerec since it doesn't change */
1515 bDoAdressWF = ((fr->adress_type != eAdressOff));
1519 update_forcerec(fplog, fr, box);
1521 if (NEED_MUTOT(*inputrec))
1523 /* Calculate total (local) dipole moment in a temporary common array.
1524 * This makes it possible to sum them over nodes faster.
1526 calc_mu(start, homenr,
1527 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
1532 if (fr->ePBC != epbcNONE)
1534 /* Compute shift vectors every step,
1535 * because of pressure coupling or box deformation!
1537 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
1539 calc_shifts(box, fr->shift_vec);
1544 put_charge_groups_in_box(fplog, cg0, cg1, fr->ePBC, box,
1545 &(top->cgs), x, fr->cg_cm);
1546 inc_nrnb(nrnb, eNR_CGCM, homenr);
1547 inc_nrnb(nrnb, eNR_RESETX, cg1-cg0);
1549 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
1551 unshift_self(graph, box, x);
1556 calc_cgcm(fplog, cg0, cg1, &(top->cgs), x, fr->cg_cm);
1557 inc_nrnb(nrnb, eNR_CGCM, homenr);
1564 move_cgcm(fplog, cr, fr->cg_cm);
1568 pr_rvecs(debug, 0, "cgcm", fr->cg_cm, top->cgs.nr);
1573 if (!(cr->duty & DUTY_PME))
1575 /* Send particle coordinates to the pme nodes.
1576 * Since this is only implemented for domain decomposition
1577 * and domain decomposition does not use the graph,
1578 * we do not need to worry about shifting.
1581 wallcycle_start(wcycle, ewcPP_PMESENDX);
1583 bBS = (inputrec->nwall == 2);
1586 copy_mat(box, boxs);
1587 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
1590 gmx_pme_send_x(cr, bBS ? boxs : box, x,
1591 mdatoms->nChargePerturbed, lambda[efptCOUL],
1592 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)), step);
1594 wallcycle_stop(wcycle, ewcPP_PMESENDX);
1596 #endif /* GMX_MPI */
1598 /* Communicate coordinates and sum dipole if necessary */
1601 wallcycle_start(wcycle, ewcMOVEX);
1602 if (DOMAINDECOMP(cr))
1604 dd_move_x(cr->dd, box, x);
1608 move_x(fplog, cr, GMX_LEFT, GMX_RIGHT, x, nrnb);
1610 wallcycle_stop(wcycle, ewcMOVEX);
1613 /* update adress weight beforehand */
1614 if (bStateChanged && bDoAdressWF)
1616 /* need pbc for adress weight calculation with pbc_dx */
1617 set_pbc(&pbc, inputrec->ePBC, box);
1618 if (fr->adress_site == eAdressSITEcog)
1620 update_adress_weights_cog(top->idef.iparams, top->idef.il, x, fr, mdatoms,
1621 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1623 else if (fr->adress_site == eAdressSITEcom)
1625 update_adress_weights_com(fplog, cg0, cg1, &(top->cgs), x, fr, mdatoms,
1626 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1628 else if (fr->adress_site == eAdressSITEatomatom)
1630 update_adress_weights_atom_per_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1631 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1635 update_adress_weights_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1636 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1640 if (NEED_MUTOT(*inputrec))
1647 gmx_sumd(2*DIM, mu, cr);
1649 for (i = 0; i < 2; i++)
1651 for (j = 0; j < DIM; j++)
1653 fr->mu_tot[i][j] = mu[i*DIM + j];
1657 if (fr->efep == efepNO)
1659 copy_rvec(fr->mu_tot[0], mu_tot);
1663 for (j = 0; j < DIM; j++)
1666 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] + lambda[efptCOUL]*fr->mu_tot[1][j];
1671 /* Reset energies */
1672 reset_enerdata(&(inputrec->opts), fr, bNS, enerd, MASTER(cr));
1673 clear_rvecs(SHIFTS, fr->fshift);
1677 wallcycle_start(wcycle, ewcNS);
1679 if (graph && bStateChanged)
1681 /* Calculate intramolecular shift vectors to make molecules whole */
1682 mk_mshift(fplog, graph, fr->ePBC, box, x);
1685 /* Do the actual neighbour searching and if twin range electrostatics
1686 * also do the calculation of long range forces and energies.
1688 for (i = 0; i < efptNR; i++)
1692 ns(fplog, fr, x, box,
1693 groups, &(inputrec->opts), top, mdatoms,
1694 cr, nrnb, lambda, dvdlambda, &enerd->grpp, bFillGrid,
1697 wallcycle_stop(wcycle, ewcNS);
1700 if (inputrec->implicit_solvent && bNS)
1702 make_gb_nblist(cr, inputrec->gb_algorithm, inputrec->rlist,
1703 x, box, fr, &top->idef, graph, fr->born);
1706 if (DOMAINDECOMP(cr))
1708 if (!(cr->duty & DUTY_PME))
1710 wallcycle_start(wcycle, ewcPPDURINGPME);
1711 dd_force_flop_start(cr->dd, nrnb);
1717 /* Enforced rotation has its own cycle counter that starts after the collective
1718 * coordinates have been communicated. It is added to ddCyclF to allow
1719 * for proper load-balancing */
1720 wallcycle_start(wcycle, ewcROT);
1721 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1722 wallcycle_stop(wcycle, ewcROT);
1725 /* Start the force cycle counter.
1726 * This counter is stopped in do_forcelow_level.
1727 * No parallel communication should occur while this counter is running,
1728 * since that will interfere with the dynamic load balancing.
1730 wallcycle_start(wcycle, ewcFORCE);
1734 /* Reset forces for which the virial is calculated separately:
1735 * PME/Ewald forces if necessary */
1736 if (fr->bF_NoVirSum)
1738 if (flags & GMX_FORCE_VIRIAL)
1740 fr->f_novirsum = fr->f_novirsum_alloc;
1743 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1747 clear_rvecs(homenr, fr->f_novirsum+start);
1752 /* We are not calculating the pressure so we do not need
1753 * a separate array for forces that do not contribute
1760 /* Clear the short- and long-range forces */
1761 clear_rvecs(fr->natoms_force_constr, f);
1762 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1764 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1767 clear_rvec(fr->vir_diag_posres);
1769 if (inputrec->ePull == epullCONSTRAINT)
1771 clear_pull_forces(inputrec->pull);
1774 /* update QMMMrec, if necessary */
1777 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1780 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1782 posres_wrapper(fplog, flags, bSepDVDL, inputrec, nrnb, top, box, x,
1783 f, enerd, lambda, fr);
1786 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_FBPOSRES].nr > 0)
1788 /* Flat-bottomed position restraints always require full pbc */
1789 if (!(bStateChanged && bDoAdressWF))
1791 set_pbc(&pbc, inputrec->ePBC, box);
1793 v = fbposres(top->idef.il[F_FBPOSRES].nr, top->idef.il[F_FBPOSRES].iatoms,
1794 top->idef.iparams_fbposres,
1795 (const rvec*)x, fr->f_novirsum, fr->vir_diag_posres,
1796 inputrec->ePBC == epbcNONE ? NULL : &pbc,
1797 fr->rc_scaling, fr->ePBC, fr->posres_com);
1798 enerd->term[F_FBPOSRES] += v;
1799 inc_nrnb(nrnb, eNR_FBPOSRES, top->idef.il[F_FBPOSRES].nr/2);
1802 /* Compute the bonded and non-bonded energies and optionally forces */
1803 do_force_lowlevel(fplog, step, fr, inputrec, &(top->idef),
1804 cr, nrnb, wcycle, mdatoms, &(inputrec->opts),
1805 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, mtop, top, fr->born,
1806 &(top->atomtypes), bBornRadii, box,
1807 inputrec->fepvals, lambda,
1808 graph, &(top->excls), fr->mu_tot,
1814 if (do_per_step(step, inputrec->nstcalclr))
1816 /* Add the long range forces to the short range forces */
1817 for (i = 0; i < fr->natoms_force_constr; i++)
1819 rvec_add(fr->f_twin[i], f[i], f[i]);
1824 cycles_force = wallcycle_stop(wcycle, ewcFORCE);
1828 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1831 if (DOMAINDECOMP(cr))
1833 dd_force_flop_stop(cr->dd, nrnb);
1836 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1842 if (IR_ELEC_FIELD(*inputrec))
1844 /* Compute forces due to electric field */
1845 calc_f_el(MASTER(cr) ? field : NULL,
1846 start, homenr, mdatoms->chargeA, x, fr->f_novirsum,
1847 inputrec->ex, inputrec->et, t);
1850 if (bDoAdressWF && fr->adress_icor == eAdressICThermoForce)
1852 /* Compute thermodynamic force in hybrid AdResS region */
1853 adress_thermo_force(start, homenr, &(top->cgs), x, fr->f_novirsum, fr, mdatoms,
1854 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1857 /* Communicate the forces */
1860 wallcycle_start(wcycle, ewcMOVEF);
1861 if (DOMAINDECOMP(cr))
1863 dd_move_f(cr->dd, f, fr->fshift);
1864 /* Do we need to communicate the separate force array
1865 * for terms that do not contribute to the single sum virial?
1866 * Position restraints and electric fields do not introduce
1867 * inter-cg forces, only full electrostatics methods do.
1868 * When we do not calculate the virial, fr->f_novirsum = f,
1869 * so we have already communicated these forces.
1871 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1872 (flags & GMX_FORCE_VIRIAL))
1874 dd_move_f(cr->dd, fr->f_novirsum, NULL);
1878 /* We should not update the shift forces here,
1879 * since f_twin is already included in f.
1881 dd_move_f(cr->dd, fr->f_twin, NULL);
1886 pd_move_f(cr, f, nrnb);
1889 pd_move_f(cr, fr->f_twin, nrnb);
1892 wallcycle_stop(wcycle, ewcMOVEF);
1895 /* If we have NoVirSum forces, but we do not calculate the virial,
1896 * we sum fr->f_novirum=f later.
1898 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1900 wallcycle_start(wcycle, ewcVSITESPREAD);
1901 spread_vsite_f(fplog, vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1902 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1903 wallcycle_stop(wcycle, ewcVSITESPREAD);
1907 wallcycle_start(wcycle, ewcVSITESPREAD);
1908 spread_vsite_f(fplog, vsite, x, fr->f_twin, NULL, FALSE, NULL,
1910 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1911 wallcycle_stop(wcycle, ewcVSITESPREAD);
1915 if (flags & GMX_FORCE_VIRIAL)
1917 /* Calculation of the virial must be done after vsites! */
1918 calc_virial(fplog, mdatoms->start, mdatoms->homenr, x, f,
1919 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1923 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1925 pull_potential_wrapper(fplog, bSepDVDL, cr, inputrec, box, x,
1926 f, vir_force, mdatoms, enerd, lambda, t);
1929 /* Add the forces from enforced rotation potentials (if any) */
1932 wallcycle_start(wcycle, ewcROTadd);
1933 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
1934 wallcycle_stop(wcycle, ewcROTadd);
1937 if (PAR(cr) && !(cr->duty & DUTY_PME))
1939 /* In case of node-splitting, the PP nodes receive the long-range
1940 * forces, virial and energy from the PME nodes here.
1942 pme_receive_force_ener(fplog, bSepDVDL, cr, wcycle, enerd, fr);
1947 post_process_forces(fplog, cr, step, nrnb, wcycle,
1948 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
1952 /* Sum the potential energy terms from group contributions */
1953 sum_epot(&(inputrec->opts), &(enerd->grpp), enerd->term);
1956 void do_force(FILE *fplog, t_commrec *cr,
1957 t_inputrec *inputrec,
1958 gmx_large_int_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1959 gmx_localtop_t *top,
1961 gmx_groups_t *groups,
1962 matrix box, rvec x[], history_t *hist,
1966 gmx_enerdata_t *enerd, t_fcdata *fcd,
1967 real *lambda, t_graph *graph,
1969 gmx_vsite_t *vsite, rvec mu_tot,
1970 double t, FILE *field, gmx_edsam_t ed,
1971 gmx_bool bBornRadii,
1974 /* modify force flag if not doing nonbonded */
1975 if (!fr->bNonbonded)
1977 flags &= ~GMX_FORCE_NONBONDED;
1980 switch (inputrec->cutoff_scheme)
1983 do_force_cutsVERLET(fplog, cr, inputrec,
1999 do_force_cutsGROUP(fplog, cr, inputrec,
2014 gmx_incons("Invalid cut-off scheme passed!");
2019 void do_constrain_first(FILE *fplog, gmx_constr_t constr,
2020 t_inputrec *ir, t_mdatoms *md,
2021 t_state *state, rvec *f,
2022 t_graph *graph, t_commrec *cr, t_nrnb *nrnb,
2023 t_forcerec *fr, gmx_localtop_t *top, tensor shake_vir)
2025 int i, m, start, end;
2026 gmx_large_int_t step;
2027 real dt = ir->delta_t;
2031 snew(savex, state->natoms);
2034 end = md->homenr + start;
2038 fprintf(debug, "vcm: start=%d, homenr=%d, end=%d\n",
2039 start, md->homenr, end);
2041 /* Do a first constrain to reset particles... */
2042 step = ir->init_step;
2045 char buf[STEPSTRSIZE];
2046 fprintf(fplog, "\nConstraining the starting coordinates (step %s)\n",
2047 gmx_step_str(step, buf));
2051 /* constrain the current position */
2052 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2053 ir, NULL, cr, step, 0, md,
2054 state->x, state->x, NULL,
2055 fr->bMolPBC, state->box,
2056 state->lambda[efptBONDED], &dvdl_dum,
2057 NULL, NULL, nrnb, econqCoord,
2058 ir->epc == epcMTTK, state->veta, state->veta);
2061 /* constrain the inital velocity, and save it */
2062 /* also may be useful if we need the ekin from the halfstep for velocity verlet */
2063 /* might not yet treat veta correctly */
2064 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2065 ir, NULL, cr, step, 0, md,
2066 state->x, state->v, state->v,
2067 fr->bMolPBC, state->box,
2068 state->lambda[efptBONDED], &dvdl_dum,
2069 NULL, NULL, nrnb, econqVeloc,
2070 ir->epc == epcMTTK, state->veta, state->veta);
2072 /* constrain the inital velocities at t-dt/2 */
2073 if (EI_STATE_VELOCITY(ir->eI) && ir->eI != eiVV)
2075 for (i = start; (i < end); i++)
2077 for (m = 0; (m < DIM); m++)
2079 /* Reverse the velocity */
2080 state->v[i][m] = -state->v[i][m];
2081 /* Store the position at t-dt in buf */
2082 savex[i][m] = state->x[i][m] + dt*state->v[i][m];
2085 /* Shake the positions at t=-dt with the positions at t=0
2086 * as reference coordinates.
2090 char buf[STEPSTRSIZE];
2091 fprintf(fplog, "\nConstraining the coordinates at t0-dt (step %s)\n",
2092 gmx_step_str(step, buf));
2095 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2096 ir, NULL, cr, step, -1, md,
2097 state->x, savex, NULL,
2098 fr->bMolPBC, state->box,
2099 state->lambda[efptBONDED], &dvdl_dum,
2100 state->v, NULL, nrnb, econqCoord,
2101 ir->epc == epcMTTK, state->veta, state->veta);
2103 for (i = start; i < end; i++)
2105 for (m = 0; m < DIM; m++)
2107 /* Re-reverse the velocities */
2108 state->v[i][m] = -state->v[i][m];
2115 void calc_enervirdiff(FILE *fplog, int eDispCorr, t_forcerec *fr)
2117 double eners[2], virs[2], enersum, virsum, y0, f, g, h;
2118 double r0, r1, r, rc3, rc9, ea, eb, ec, pa, pb, pc, pd;
2119 double invscale, invscale2, invscale3;
2120 int ri0, ri1, ri, i, offstart, offset;
2121 real scale, *vdwtab, tabfactor, tmp;
2123 fr->enershiftsix = 0;
2124 fr->enershifttwelve = 0;
2125 fr->enerdiffsix = 0;
2126 fr->enerdifftwelve = 0;
2128 fr->virdifftwelve = 0;
2130 if (eDispCorr != edispcNO)
2132 for (i = 0; i < 2; i++)
2137 if ((fr->vdwtype == evdwSWITCH) || (fr->vdwtype == evdwSHIFT))
2139 if (fr->rvdw_switch == 0)
2142 "With dispersion correction rvdw-switch can not be zero "
2143 "for vdw-type = %s", evdw_names[fr->vdwtype]);
2146 scale = fr->nblists[0].table_elec_vdw.scale;
2147 vdwtab = fr->nblists[0].table_vdw.data;
2149 /* Round the cut-offs to exact table values for precision */
2150 ri0 = floor(fr->rvdw_switch*scale);
2151 ri1 = ceil(fr->rvdw*scale);
2157 if (fr->vdwtype == evdwSHIFT)
2159 /* Determine the constant energy shift below rvdw_switch.
2160 * Table has a scale factor since we have scaled it down to compensate
2161 * for scaling-up c6/c12 with the derivative factors to save flops in analytical kernels.
2163 fr->enershiftsix = (real)(-1.0/(rc3*rc3)) - 6.0*vdwtab[8*ri0];
2164 fr->enershifttwelve = (real)( 1.0/(rc9*rc3)) - 12.0*vdwtab[8*ri0 + 4];
2166 /* Add the constant part from 0 to rvdw_switch.
2167 * This integration from 0 to rvdw_switch overcounts the number
2168 * of interactions by 1, as it also counts the self interaction.
2169 * We will correct for this later.
2171 eners[0] += 4.0*M_PI*fr->enershiftsix*rc3/3.0;
2172 eners[1] += 4.0*M_PI*fr->enershifttwelve*rc3/3.0;
2174 invscale = 1.0/(scale);
2175 invscale2 = invscale*invscale;
2176 invscale3 = invscale*invscale2;
2178 /* following summation derived from cubic spline definition,
2179 Numerical Recipies in C, second edition, p. 113-116. Exact
2180 for the cubic spline. We first calculate the negative of
2181 the energy from rvdw to rvdw_switch, assuming that g(r)=1,
2182 and then add the more standard, abrupt cutoff correction to
2183 that result, yielding the long-range correction for a
2184 switched function. We perform both the pressure and energy
2185 loops at the same time for simplicity, as the computational
2188 for (i = 0; i < 2; i++)
2190 enersum = 0.0; virsum = 0.0;
2194 /* Since the dispersion table has been scaled down a factor 6.0 and the repulsion
2195 * a factor 12.0 to compensate for the c6/c12 parameters inside nbfp[] being scaled
2196 * up (to save flops in kernels), we need to correct for this.
2205 for (ri = ri0; ri < ri1; ri++)
2209 eb = 2.0*invscale2*r;
2213 pb = 3.0*invscale2*r;
2214 pc = 3.0*invscale*r*r;
2217 /* this "8" is from the packing in the vdwtab array - perhaps should be #define'ed? */
2218 offset = 8*ri + offstart;
2219 y0 = vdwtab[offset];
2220 f = vdwtab[offset+1];
2221 g = vdwtab[offset+2];
2222 h = vdwtab[offset+3];
2224 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);
2225 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);
2228 enersum *= 4.0*M_PI*tabfactor;
2229 virsum *= 4.0*M_PI*tabfactor;
2230 eners[i] -= enersum;
2234 /* now add the correction for rvdw_switch to infinity */
2235 eners[0] += -4.0*M_PI/(3.0*rc3);
2236 eners[1] += 4.0*M_PI/(9.0*rc9);
2237 virs[0] += 8.0*M_PI/rc3;
2238 virs[1] += -16.0*M_PI/(3.0*rc9);
2240 else if ((fr->vdwtype == evdwCUT) || (fr->vdwtype == evdwUSER))
2242 if (fr->vdwtype == evdwUSER && fplog)
2245 "WARNING: using dispersion correction with user tables\n");
2247 rc3 = fr->rvdw*fr->rvdw*fr->rvdw;
2249 /* Contribution beyond the cut-off */
2250 eners[0] += -4.0*M_PI/(3.0*rc3);
2251 eners[1] += 4.0*M_PI/(9.0*rc9);
2252 if (fr->vdw_modifier == eintmodPOTSHIFT)
2254 /* Contribution within the cut-off */
2255 eners[0] += -4.0*M_PI/(3.0*rc3);
2256 eners[1] += 4.0*M_PI/(3.0*rc9);
2258 /* Contribution beyond the cut-off */
2259 virs[0] += 8.0*M_PI/rc3;
2260 virs[1] += -16.0*M_PI/(3.0*rc9);
2265 "Dispersion correction is not implemented for vdw-type = %s",
2266 evdw_names[fr->vdwtype]);
2268 fr->enerdiffsix = eners[0];
2269 fr->enerdifftwelve = eners[1];
2270 /* The 0.5 is due to the Gromacs definition of the virial */
2271 fr->virdiffsix = 0.5*virs[0];
2272 fr->virdifftwelve = 0.5*virs[1];
2276 void calc_dispcorr(FILE *fplog, t_inputrec *ir, t_forcerec *fr,
2277 gmx_large_int_t step, int natoms,
2278 matrix box, real lambda, tensor pres, tensor virial,
2279 real *prescorr, real *enercorr, real *dvdlcorr)
2281 gmx_bool bCorrAll, bCorrPres;
2282 real dvdlambda, invvol, dens, ninter, avcsix, avctwelve, enerdiff, svir = 0, spres = 0;
2292 if (ir->eDispCorr != edispcNO)
2294 bCorrAll = (ir->eDispCorr == edispcAllEner ||
2295 ir->eDispCorr == edispcAllEnerPres);
2296 bCorrPres = (ir->eDispCorr == edispcEnerPres ||
2297 ir->eDispCorr == edispcAllEnerPres);
2299 invvol = 1/det(box);
2302 /* Only correct for the interactions with the inserted molecule */
2303 dens = (natoms - fr->n_tpi)*invvol;
2308 dens = natoms*invvol;
2309 ninter = 0.5*natoms;
2312 if (ir->efep == efepNO)
2314 avcsix = fr->avcsix[0];
2315 avctwelve = fr->avctwelve[0];
2319 avcsix = (1 - lambda)*fr->avcsix[0] + lambda*fr->avcsix[1];
2320 avctwelve = (1 - lambda)*fr->avctwelve[0] + lambda*fr->avctwelve[1];
2323 enerdiff = ninter*(dens*fr->enerdiffsix - fr->enershiftsix);
2324 *enercorr += avcsix*enerdiff;
2326 if (ir->efep != efepNO)
2328 dvdlambda += (fr->avcsix[1] - fr->avcsix[0])*enerdiff;
2332 enerdiff = ninter*(dens*fr->enerdifftwelve - fr->enershifttwelve);
2333 *enercorr += avctwelve*enerdiff;
2334 if (fr->efep != efepNO)
2336 dvdlambda += (fr->avctwelve[1] - fr->avctwelve[0])*enerdiff;
2342 svir = ninter*dens*avcsix*fr->virdiffsix/3.0;
2343 if (ir->eDispCorr == edispcAllEnerPres)
2345 svir += ninter*dens*avctwelve*fr->virdifftwelve/3.0;
2347 /* The factor 2 is because of the Gromacs virial definition */
2348 spres = -2.0*invvol*svir*PRESFAC;
2350 for (m = 0; m < DIM; m++)
2352 virial[m][m] += svir;
2353 pres[m][m] += spres;
2358 /* Can't currently control when it prints, for now, just print when degugging */
2363 fprintf(debug, "Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
2369 "Long Range LJ corr.: Epot %10g, Pres: %10g, Vir: %10g\n",
2370 *enercorr, spres, svir);
2374 fprintf(debug, "Long Range LJ corr.: Epot %10g\n", *enercorr);
2378 if (fr->bSepDVDL && do_per_step(step, ir->nstlog))
2380 gmx_print_sepdvdl(fplog, "Dispersion correction", *enercorr, dvdlambda);
2382 if (fr->efep != efepNO)
2384 *dvdlcorr += dvdlambda;
2389 void do_pbc_first(FILE *fplog, matrix box, t_forcerec *fr,
2390 t_graph *graph, rvec x[])
2394 fprintf(fplog, "Removing pbc first time\n");
2396 calc_shifts(box, fr->shift_vec);
2399 mk_mshift(fplog, graph, fr->ePBC, box, x);
2402 p_graph(debug, "do_pbc_first 1", graph);
2404 shift_self(graph, box, x);
2405 /* By doing an extra mk_mshift the molecules that are broken
2406 * because they were e.g. imported from another software
2407 * will be made whole again. Such are the healing powers
2410 mk_mshift(fplog, graph, fr->ePBC, box, x);
2413 p_graph(debug, "do_pbc_first 2", graph);
2418 fprintf(fplog, "Done rmpbc\n");
2422 static void low_do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2423 gmx_mtop_t *mtop, rvec x[],
2428 gmx_molblock_t *molb;
2430 if (bFirst && fplog)
2432 fprintf(fplog, "Removing pbc first time\n");
2437 for (mb = 0; mb < mtop->nmolblock; mb++)
2439 molb = &mtop->molblock[mb];
2440 if (molb->natoms_mol == 1 ||
2441 (!bFirst && mtop->moltype[molb->type].cgs.nr == 1))
2443 /* Just one atom or charge group in the molecule, no PBC required */
2444 as += molb->nmol*molb->natoms_mol;
2448 /* Pass NULL iso fplog to avoid graph prints for each molecule type */
2449 mk_graph_ilist(NULL, mtop->moltype[molb->type].ilist,
2450 0, molb->natoms_mol, FALSE, FALSE, graph);
2452 for (mol = 0; mol < molb->nmol; mol++)
2454 mk_mshift(fplog, graph, ePBC, box, x+as);
2456 shift_self(graph, box, x+as);
2457 /* The molecule is whole now.
2458 * We don't need the second mk_mshift call as in do_pbc_first,
2459 * since we no longer need this graph.
2462 as += molb->natoms_mol;
2470 void do_pbc_first_mtop(FILE *fplog, int ePBC, matrix box,
2471 gmx_mtop_t *mtop, rvec x[])
2473 low_do_pbc_mtop(fplog, ePBC, box, mtop, x, TRUE);
2476 void do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2477 gmx_mtop_t *mtop, rvec x[])
2479 low_do_pbc_mtop(fplog, ePBC, box, mtop, x, FALSE);
2482 void finish_run(FILE *fplog, t_commrec *cr, const char *confout,
2483 t_inputrec *inputrec,
2484 t_nrnb nrnb[], gmx_wallcycle_t wcycle,
2485 gmx_runtime_t *runtime,
2486 wallclock_gpu_t *gputimes,
2488 gmx_bool bWriteStat)
2491 t_nrnb *nrnb_tot = NULL;
2495 wallcycle_sum(cr, wcycle);
2501 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM,
2502 cr->mpi_comm_mysim);
2510 #if defined(GMX_MPI) && !defined(GMX_THREAD_MPI)
2513 /* reduce nodetime over all MPI processes in the current simulation */
2515 MPI_Allreduce(&runtime->proctime, &sum, 1, MPI_DOUBLE, MPI_SUM,
2516 cr->mpi_comm_mysim);
2517 runtime->proctime = sum;
2523 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
2530 if ((cr->duty & DUTY_PP) && DOMAINDECOMP(cr))
2532 print_dd_statistics(cr, inputrec, fplog);
2544 snew(nrnb_all, cr->nnodes);
2545 nrnb_all[0] = *nrnb;
2546 for (s = 1; s < cr->nnodes; s++)
2548 MPI_Recv(nrnb_all[s].n, eNRNB, MPI_DOUBLE, s, 0,
2549 cr->mpi_comm_mysim, &stat);
2551 pr_load(fplog, cr, nrnb_all);
2556 MPI_Send(nrnb->n, eNRNB, MPI_DOUBLE, MASTERRANK(cr), 0,
2557 cr->mpi_comm_mysim);
2564 wallcycle_print(fplog, cr->nnodes, cr->npmenodes, runtime->realtime,
2567 if (EI_DYNAMICS(inputrec->eI))
2569 delta_t = inputrec->delta_t;
2578 print_perf(fplog, runtime->proctime, runtime->realtime,
2579 runtime->nsteps_done, delta_t, nbfs, mflop);
2583 print_perf(stderr, runtime->proctime, runtime->realtime,
2584 runtime->nsteps_done, delta_t, nbfs, mflop);
2589 extern void initialize_lambdas(FILE *fplog, t_inputrec *ir, int *fep_state, real *lambda, double *lam0)
2591 /* this function works, but could probably use a logic rewrite to keep all the different
2592 types of efep straight. */
2595 t_lambda *fep = ir->fepvals;
2597 if ((ir->efep == efepNO) && (ir->bSimTemp == FALSE))
2599 for (i = 0; i < efptNR; i++)
2611 *fep_state = fep->init_fep_state; /* this might overwrite the checkpoint
2612 if checkpoint is set -- a kludge is in for now
2614 for (i = 0; i < efptNR; i++)
2616 /* overwrite lambda state with init_lambda for now for backwards compatibility */
2617 if (fep->init_lambda >= 0) /* if it's -1, it was never initializd */
2619 lambda[i] = fep->init_lambda;
2622 lam0[i] = lambda[i];
2627 lambda[i] = fep->all_lambda[i][*fep_state];
2630 lam0[i] = lambda[i];
2636 /* need to rescale control temperatures to match current state */
2637 for (i = 0; i < ir->opts.ngtc; i++)
2639 if (ir->opts.ref_t[i] > 0)
2641 ir->opts.ref_t[i] = ir->simtempvals->temperatures[*fep_state];
2647 /* Send to the log the information on the current lambdas */
2650 fprintf(fplog, "Initial vector of lambda components:[ ");
2651 for (i = 0; i < efptNR; i++)
2653 fprintf(fplog, "%10.4f ", lambda[i]);
2655 fprintf(fplog, "]\n");
2661 void init_md(FILE *fplog,
2662 t_commrec *cr, t_inputrec *ir, const output_env_t oenv,
2663 double *t, double *t0,
2664 real *lambda, int *fep_state, double *lam0,
2665 t_nrnb *nrnb, gmx_mtop_t *mtop,
2667 int nfile, const t_filenm fnm[],
2668 gmx_mdoutf_t **outf, t_mdebin **mdebin,
2669 tensor force_vir, tensor shake_vir, rvec mu_tot,
2670 gmx_bool *bSimAnn, t_vcm **vcm, t_state *state, unsigned long Flags)
2675 /* Initial values */
2676 *t = *t0 = ir->init_t;
2679 for (i = 0; i < ir->opts.ngtc; i++)
2681 /* set bSimAnn if any group is being annealed */
2682 if (ir->opts.annealing[i] != eannNO)
2689 update_annealing_target_temp(&(ir->opts), ir->init_t);
2692 /* Initialize lambda variables */
2693 initialize_lambdas(fplog, ir, fep_state, lambda, lam0);
2697 *upd = init_update(fplog, ir);
2703 *vcm = init_vcm(fplog, &mtop->groups, ir);
2706 if (EI_DYNAMICS(ir->eI) && !(Flags & MD_APPENDFILES))
2708 if (ir->etc == etcBERENDSEN)
2710 please_cite(fplog, "Berendsen84a");
2712 if (ir->etc == etcVRESCALE)
2714 please_cite(fplog, "Bussi2007a");
2722 *outf = init_mdoutf(nfile, fnm, Flags, cr, ir, oenv);
2724 *mdebin = init_mdebin((Flags & MD_APPENDFILES) ? NULL : (*outf)->fp_ene,
2725 mtop, ir, (*outf)->fp_dhdl);
2730 please_cite(fplog, "Fritsch12");
2731 please_cite(fplog, "Junghans10");
2733 /* Initiate variables */
2734 clear_mat(force_vir);
2735 clear_mat(shake_vir);