1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * GROwing Monsters And Cloning Shrimps
41 #include <catamount/dclock.h>
47 #ifdef HAVE_SYS_TIME_H
60 #include "chargegroup.h"
83 #include "pull_rotation.h"
84 #include "gmx_random.h"
87 #include "gmx_wallcycle.h"
89 #include "nbnxn_atomdata.h"
90 #include "nbnxn_search.h"
91 #include "nbnxn_kernels/nbnxn_kernel_ref.h"
92 #include "nbnxn_kernels/nbnxn_kernel_simd_4xn.h"
93 #include "nbnxn_kernels/nbnxn_kernel_simd_2xnn.h"
94 #include "nbnxn_kernels/nbnxn_kernel_gpu_ref.h"
96 #include "gromacs/utility/gmxmpi.h"
101 #include "nbnxn_cuda_data_mgmt.h"
102 #include "nbnxn_cuda/nbnxn_cuda.h"
107 #ifdef HAVE_GETTIMEOFDAY
111 gettimeofday(&t, NULL);
113 seconds = (double) t.tv_sec + 1e-6*(double)t.tv_usec;
119 seconds = time(NULL);
126 #define difftime(end, start) ((double)(end)-(double)(start))
128 void print_time(FILE *out, gmx_runtime_t *runtime, gmx_large_int_t step,
129 t_inputrec *ir, t_commrec *cr)
132 char timebuf[STRLEN];
136 #ifndef GMX_THREAD_MPI
142 fprintf(out, "step %s", gmx_step_str(step, buf));
143 if ((step >= ir->nstlist))
145 runtime->last = gmx_gettime();
146 dt = difftime(runtime->last, runtime->real);
147 runtime->time_per_step = dt/(step - ir->init_step + 1);
149 dt = (ir->nsteps + ir->init_step - step)*runtime->time_per_step;
155 finish = (time_t) (runtime->last + dt);
156 gmx_ctime_r(&finish, timebuf, STRLEN);
157 sprintf(buf, "%s", timebuf);
158 buf[strlen(buf)-1] = '\0';
159 fprintf(out, ", will finish %s", buf);
163 fprintf(out, ", remaining runtime: %5d s ", (int)dt);
168 fprintf(out, " performance: %.1f ns/day ",
169 ir->delta_t/1000*24*60*60/runtime->time_per_step);
172 #ifndef GMX_THREAD_MPI
186 static double set_proctime(gmx_runtime_t *runtime)
192 prev = runtime->proc;
193 runtime->proc = dclock();
195 diff = runtime->proc - prev;
199 prev = runtime->proc;
200 runtime->proc = clock();
202 diff = (double)(runtime->proc - prev)/(double)CLOCKS_PER_SEC;
206 /* The counter has probably looped, ignore this data */
213 void runtime_start(gmx_runtime_t *runtime)
215 runtime->real = gmx_gettime();
217 set_proctime(runtime);
218 runtime->realtime = 0;
219 runtime->proctime = 0;
221 runtime->time_per_step = 0;
224 void runtime_end(gmx_runtime_t *runtime)
230 runtime->proctime += set_proctime(runtime);
231 runtime->realtime = now - runtime->real;
235 void runtime_upd_proc(gmx_runtime_t *runtime)
237 runtime->proctime += set_proctime(runtime);
240 void print_date_and_time(FILE *fplog, int nodeid, const char *title,
241 const gmx_runtime_t *runtime)
244 char timebuf[STRLEN];
245 char time_string[STRLEN];
252 tmptime = (time_t) runtime->real;
253 gmx_ctime_r(&tmptime, timebuf, STRLEN);
257 tmptime = (time_t) gmx_gettime();
258 gmx_ctime_r(&tmptime, timebuf, STRLEN);
260 for (i = 0; timebuf[i] >= ' '; i++)
262 time_string[i] = timebuf[i];
264 time_string[i] = '\0';
266 fprintf(fplog, "%s on node %d %s\n", title, nodeid, time_string);
270 static void sum_forces(int start, int end, rvec f[], rvec flr[])
276 pr_rvecs(debug, 0, "fsr", f+start, end-start);
277 pr_rvecs(debug, 0, "flr", flr+start, end-start);
279 for (i = start; (i < end); i++)
281 rvec_inc(f[i], flr[i]);
286 * calc_f_el calculates forces due to an electric field.
288 * force is kJ mol^-1 nm^-1 = e * kJ mol^-1 nm^-1 / e
290 * Et[] contains the parameters for the time dependent
291 * part of the field (not yet used).
292 * Ex[] contains the parameters for
293 * the spatial dependent part of the field. You can have cool periodic
294 * fields in principle, but only a constant field is supported
296 * The function should return the energy due to the electric field
297 * (if any) but for now returns 0.
300 * There can be problems with the virial.
301 * Since the field is not self-consistent this is unavoidable.
302 * For neutral molecules the virial is correct within this approximation.
303 * For neutral systems with many charged molecules the error is small.
304 * But for systems with a net charge or a few charged molecules
305 * the error can be significant when the field is high.
306 * Solution: implement a self-consitent electric field into PME.
308 static void calc_f_el(FILE *fp, int start, int homenr,
309 real charge[], rvec x[], rvec f[],
310 t_cosines Ex[], t_cosines Et[], double t)
316 for (m = 0; (m < DIM); m++)
323 Ext[m] = cos(Et[m].a[0]*(t-t0))*exp(-sqr(t-t0)/(2.0*sqr(Et[m].a[2])));
327 Ext[m] = cos(Et[m].a[0]*t);
336 /* Convert the field strength from V/nm to MD-units */
337 Ext[m] *= Ex[m].a[0]*FIELDFAC;
338 for (i = start; (i < start+homenr); i++)
340 f[i][m] += charge[i]*Ext[m];
350 fprintf(fp, "%10g %10g %10g %10g #FIELD\n", t,
351 Ext[XX]/FIELDFAC, Ext[YY]/FIELDFAC, Ext[ZZ]/FIELDFAC);
355 static void calc_virial(FILE *fplog, int start, int homenr, rvec x[], rvec f[],
356 tensor vir_part, t_graph *graph, matrix box,
357 t_nrnb *nrnb, const t_forcerec *fr, int ePBC)
362 /* The short-range virial from surrounding boxes */
364 calc_vir(fplog, SHIFTS, fr->shift_vec, fr->fshift, vir_part, ePBC == epbcSCREW, box);
365 inc_nrnb(nrnb, eNR_VIRIAL, SHIFTS);
367 /* Calculate partial virial, for local atoms only, based on short range.
368 * Total virial is computed in global_stat, called from do_md
370 f_calc_vir(fplog, start, start+homenr, x, f, vir_part, graph, box);
371 inc_nrnb(nrnb, eNR_VIRIAL, homenr);
373 /* Add position restraint contribution */
374 for (i = 0; i < DIM; i++)
376 vir_part[i][i] += fr->vir_diag_posres[i];
379 /* Add wall contribution */
380 for (i = 0; i < DIM; i++)
382 vir_part[i][ZZ] += fr->vir_wall_z[i];
387 pr_rvecs(debug, 0, "vir_part", vir_part, DIM);
391 static void posres_wrapper(FILE *fplog,
397 matrix box, rvec x[],
399 gmx_enerdata_t *enerd,
407 /* Position restraints always require full pbc */
408 set_pbc(&pbc, ir->ePBC, box);
410 v = posres(top->idef.il[F_POSRES].nr, top->idef.il[F_POSRES].iatoms,
411 top->idef.iparams_posres,
412 (const rvec*)x, fr->f_novirsum, fr->vir_diag_posres,
413 ir->ePBC == epbcNONE ? NULL : &pbc,
414 lambda[efptRESTRAINT], &dvdl,
415 fr->rc_scaling, fr->ePBC, fr->posres_com, fr->posres_comB);
418 gmx_print_sepdvdl(fplog, interaction_function[F_POSRES].longname, v, dvdl);
420 enerd->term[F_POSRES] += v;
421 /* If just the force constant changes, the FEP term is linear,
422 * but if k changes, it is not.
424 enerd->dvdl_nonlin[efptRESTRAINT] += dvdl;
425 inc_nrnb(nrnb, eNR_POSRES, top->idef.il[F_POSRES].nr/2);
427 if ((ir->fepvals->n_lambda > 0) && (flags & GMX_FORCE_DHDL))
429 for (i = 0; i < enerd->n_lambda; i++)
431 real dvdl_dum, lambda_dum;
433 lambda_dum = (i == 0 ? lambda[efptRESTRAINT] : ir->fepvals->all_lambda[efptRESTRAINT][i-1]);
434 v = posres(top->idef.il[F_POSRES].nr, top->idef.il[F_POSRES].iatoms,
435 top->idef.iparams_posres,
436 (const rvec*)x, NULL, NULL,
437 ir->ePBC == epbcNONE ? NULL : &pbc, lambda_dum, &dvdl,
438 fr->rc_scaling, fr->ePBC, fr->posres_com, fr->posres_comB);
439 enerd->enerpart_lambda[i] += v;
444 static void pull_potential_wrapper(FILE *fplog,
448 matrix box, rvec x[],
452 gmx_enerdata_t *enerd,
459 /* Calculate the center of mass forces, this requires communication,
460 * which is why pull_potential is called close to other communication.
461 * The virial contribution is calculated directly,
462 * which is why we call pull_potential after calc_virial.
464 set_pbc(&pbc, ir->ePBC, box);
466 enerd->term[F_COM_PULL] +=
467 pull_potential(ir->ePull, ir->pull, mdatoms, &pbc,
468 cr, t, lambda[efptRESTRAINT], x, f, vir_force, &dvdl);
471 gmx_print_sepdvdl(fplog, "Com pull", enerd->term[F_COM_PULL], dvdl);
473 enerd->dvdl_lin[efptRESTRAINT] += dvdl;
476 static void pme_receive_force_ener(FILE *fplog,
479 gmx_wallcycle_t wcycle,
480 gmx_enerdata_t *enerd,
484 float cycles_ppdpme, cycles_seppme;
486 cycles_ppdpme = wallcycle_stop(wcycle, ewcPPDURINGPME);
487 dd_cycles_add(cr->dd, cycles_ppdpme, ddCyclPPduringPME);
489 /* In case of node-splitting, the PP nodes receive the long-range
490 * forces, virial and energy from the PME nodes here.
492 wallcycle_start(wcycle, ewcPP_PMEWAITRECVF);
494 gmx_pme_receive_f(cr, fr->f_novirsum, fr->vir_el_recip, &e, &dvdl,
498 gmx_print_sepdvdl(fplog, "PME mesh", e, dvdl);
500 enerd->term[F_COUL_RECIP] += e;
501 enerd->dvdl_lin[efptCOUL] += dvdl;
504 dd_cycles_add(cr->dd, cycles_seppme, ddCyclPME);
506 wallcycle_stop(wcycle, ewcPP_PMEWAITRECVF);
509 static void print_large_forces(FILE *fp, t_mdatoms *md, t_commrec *cr,
510 gmx_large_int_t step, real pforce, rvec *x, rvec *f)
514 char buf[STEPSTRSIZE];
517 for (i = md->start; i < md->start+md->homenr; i++)
520 /* We also catch NAN, if the compiler does not optimize this away. */
521 if (fn2 >= pf2 || fn2 != fn2)
523 fprintf(fp, "step %s atom %6d x %8.3f %8.3f %8.3f force %12.5e\n",
524 gmx_step_str(step, buf),
525 ddglatnr(cr->dd, i), x[i][XX], x[i][YY], x[i][ZZ], sqrt(fn2));
530 static void post_process_forces(FILE *fplog,
532 gmx_large_int_t step,
533 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
535 matrix box, rvec x[],
540 t_forcerec *fr, gmx_vsite_t *vsite,
547 /* Spread the mesh force on virtual sites to the other particles...
548 * This is parallellized. MPI communication is performed
549 * if the constructing atoms aren't local.
551 wallcycle_start(wcycle, ewcVSITESPREAD);
552 spread_vsite_f(fplog, vsite, x, fr->f_novirsum, NULL,
553 (flags & GMX_FORCE_VIRIAL), fr->vir_el_recip,
555 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
556 wallcycle_stop(wcycle, ewcVSITESPREAD);
558 if (flags & GMX_FORCE_VIRIAL)
560 /* Now add the forces, this is local */
563 sum_forces(0, fr->f_novirsum_n, f, fr->f_novirsum);
567 sum_forces(mdatoms->start, mdatoms->start+mdatoms->homenr,
570 if (EEL_FULL(fr->eeltype))
572 /* Add the mesh contribution to the virial */
573 m_add(vir_force, fr->vir_el_recip, vir_force);
577 pr_rvecs(debug, 0, "vir_force", vir_force, DIM);
582 if (fr->print_force >= 0)
584 print_large_forces(stderr, mdatoms, cr, step, fr->print_force, x, f);
588 static void do_nb_verlet(t_forcerec *fr,
589 interaction_const_t *ic,
590 gmx_enerdata_t *enerd,
591 int flags, int ilocality,
594 gmx_wallcycle_t wcycle)
596 int nnbl, kernel_type, enr_nbnxn_kernel_ljc, enr_nbnxn_kernel_lj;
598 nonbonded_verlet_group_t *nbvg;
600 if (!(flags & GMX_FORCE_NONBONDED))
602 /* skip non-bonded calculation */
606 nbvg = &fr->nbv->grp[ilocality];
608 /* CUDA kernel launch overhead is already timed separately */
609 if (fr->cutoff_scheme != ecutsVERLET)
611 gmx_incons("Invalid cut-off scheme passed!");
614 if (nbvg->kernel_type != nbnxnk8x8x8_CUDA)
616 wallcycle_sub_start(wcycle, ewcsNONBONDED);
618 switch (nbvg->kernel_type)
620 case nbnxnk4x4_PlainC:
621 nbnxn_kernel_ref(&nbvg->nbl_lists,
627 enerd->grpp.ener[egCOULSR],
629 enerd->grpp.ener[egBHAMSR] :
630 enerd->grpp.ener[egLJSR]);
633 case nbnxnk4xN_SIMD_4xN:
634 nbnxn_kernel_simd_4xn(&nbvg->nbl_lists,
641 enerd->grpp.ener[egCOULSR],
643 enerd->grpp.ener[egBHAMSR] :
644 enerd->grpp.ener[egLJSR]);
646 case nbnxnk4xN_SIMD_2xNN:
647 nbnxn_kernel_simd_2xnn(&nbvg->nbl_lists,
654 enerd->grpp.ener[egCOULSR],
656 enerd->grpp.ener[egBHAMSR] :
657 enerd->grpp.ener[egLJSR]);
660 case nbnxnk8x8x8_CUDA:
661 nbnxn_cuda_launch_kernel(fr->nbv->cu_nbv, nbvg->nbat, flags, ilocality);
664 case nbnxnk8x8x8_PlainC:
665 nbnxn_kernel_gpu_ref(nbvg->nbl_lists.nbl[0],
670 nbvg->nbat->out[0].f,
672 enerd->grpp.ener[egCOULSR],
674 enerd->grpp.ener[egBHAMSR] :
675 enerd->grpp.ener[egLJSR]);
679 gmx_incons("Invalid nonbonded kernel type passed!");
682 if (nbvg->kernel_type != nbnxnk8x8x8_CUDA)
684 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
687 if (EEL_RF(ic->eeltype) || ic->eeltype == eelCUT)
689 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_RF;
691 else if (nbvg->ewald_excl == ewaldexclTable)
693 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_TAB;
697 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_EWALD;
699 enr_nbnxn_kernel_lj = eNR_NBNXN_LJ;
700 if (flags & GMX_FORCE_ENERGY)
702 /* In eNR_??? the nbnxn F+E kernels are always the F kernel + 1 */
703 enr_nbnxn_kernel_ljc += 1;
704 enr_nbnxn_kernel_lj += 1;
707 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc,
708 nbvg->nbl_lists.natpair_ljq);
709 inc_nrnb(nrnb, enr_nbnxn_kernel_lj,
710 nbvg->nbl_lists.natpair_lj);
711 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc-eNR_NBNXN_LJ_RF+eNR_NBNXN_RF,
712 nbvg->nbl_lists.natpair_q);
715 void do_force_cutsVERLET(FILE *fplog, t_commrec *cr,
716 t_inputrec *inputrec,
717 gmx_large_int_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
720 gmx_groups_t *groups,
721 matrix box, rvec x[], history_t *hist,
725 gmx_enerdata_t *enerd, t_fcdata *fcd,
726 real *lambda, t_graph *graph,
727 t_forcerec *fr, interaction_const_t *ic,
728 gmx_vsite_t *vsite, rvec mu_tot,
729 double t, FILE *field, gmx_edsam_t ed,
737 gmx_bool bSepDVDL, bStateChanged, bNS, bFillGrid, bCalcCGCM, bBS;
738 gmx_bool bDoLongRange, bDoForces, bSepLRF, bUseGPU, bUseOrEmulGPU;
739 gmx_bool bDiffKernels = FALSE;
741 rvec vzero, box_diag;
743 float cycles_pme, cycles_force;
744 nonbonded_verlet_t *nbv;
748 nb_kernel_type = fr->nbv->grp[0].kernel_type;
750 start = mdatoms->start;
751 homenr = mdatoms->homenr;
753 bSepDVDL = (fr->bSepDVDL && do_per_step(step, inputrec->nstlog));
755 clear_mat(vir_force);
758 if (DOMAINDECOMP(cr))
760 cg1 = cr->dd->ncg_tot;
771 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
772 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
773 bFillGrid = (bNS && bStateChanged);
774 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
775 bDoLongRange = (fr->bTwinRange && bNS && (flags & GMX_FORCE_DO_LR));
776 bDoForces = (flags & GMX_FORCE_FORCES);
777 bSepLRF = (bDoLongRange && bDoForces && (flags & GMX_FORCE_SEPLRF));
778 bUseGPU = fr->nbv->bUseGPU;
779 bUseOrEmulGPU = bUseGPU || (nbv->grp[0].kernel_type == nbnxnk8x8x8_PlainC);
783 update_forcerec(fplog, fr, box);
785 if (NEED_MUTOT(*inputrec))
787 /* Calculate total (local) dipole moment in a temporary common array.
788 * This makes it possible to sum them over nodes faster.
790 calc_mu(start, homenr,
791 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
796 if (fr->ePBC != epbcNONE)
798 /* Compute shift vectors every step,
799 * because of pressure coupling or box deformation!
801 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
803 calc_shifts(box, fr->shift_vec);
808 put_atoms_in_box_omp(fr->ePBC, box, homenr, x);
809 inc_nrnb(nrnb, eNR_SHIFTX, homenr);
811 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
813 unshift_self(graph, box, x);
817 nbnxn_atomdata_copy_shiftvec(flags & GMX_FORCE_DYNAMICBOX,
818 fr->shift_vec, nbv->grp[0].nbat);
821 if (!(cr->duty & DUTY_PME))
823 /* Send particle coordinates to the pme nodes.
824 * Since this is only implemented for domain decomposition
825 * and domain decomposition does not use the graph,
826 * we do not need to worry about shifting.
829 wallcycle_start(wcycle, ewcPP_PMESENDX);
831 bBS = (inputrec->nwall == 2);
835 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
838 gmx_pme_send_x(cr, bBS ? boxs : box, x,
839 mdatoms->nChargePerturbed, lambda[efptCOUL],
840 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)), step);
842 wallcycle_stop(wcycle, ewcPP_PMESENDX);
846 /* do gridding for pair search */
849 if (graph && bStateChanged)
851 /* Calculate intramolecular shift vectors to make molecules whole */
852 mk_mshift(fplog, graph, fr->ePBC, box, x);
856 box_diag[XX] = box[XX][XX];
857 box_diag[YY] = box[YY][YY];
858 box_diag[ZZ] = box[ZZ][ZZ];
860 wallcycle_start(wcycle, ewcNS);
863 wallcycle_sub_start(wcycle, ewcsNBS_GRID_LOCAL);
864 nbnxn_put_on_grid(nbv->nbs, fr->ePBC, box,
866 0, mdatoms->homenr, -1, fr->cginfo, x,
868 nbv->grp[eintLocal].kernel_type,
869 nbv->grp[eintLocal].nbat);
870 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_LOCAL);
874 wallcycle_sub_start(wcycle, ewcsNBS_GRID_NONLOCAL);
875 nbnxn_put_on_grid_nonlocal(nbv->nbs, domdec_zones(cr->dd),
877 nbv->grp[eintNonlocal].kernel_type,
878 nbv->grp[eintNonlocal].nbat);
879 wallcycle_sub_stop(wcycle, ewcsNBS_GRID_NONLOCAL);
882 if (nbv->ngrp == 1 ||
883 nbv->grp[eintNonlocal].nbat == nbv->grp[eintLocal].nbat)
885 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatAll,
886 nbv->nbs, mdatoms, fr->cginfo);
890 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat, eatLocal,
891 nbv->nbs, mdatoms, fr->cginfo);
892 nbnxn_atomdata_set(nbv->grp[eintNonlocal].nbat, eatAll,
893 nbv->nbs, mdatoms, fr->cginfo);
895 wallcycle_stop(wcycle, ewcNS);
898 /* initialize the GPU atom data and copy shift vector */
903 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
904 nbnxn_cuda_init_atomdata(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
905 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
908 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
909 nbnxn_cuda_upload_shiftvec(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
910 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
913 /* do local pair search */
916 wallcycle_start_nocount(wcycle, ewcNS);
917 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_LOCAL);
918 nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintLocal].nbat,
921 nbv->min_ci_balanced,
922 &nbv->grp[eintLocal].nbl_lists,
924 nbv->grp[eintLocal].kernel_type,
926 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_LOCAL);
930 /* initialize local pair-list on the GPU */
931 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
932 nbv->grp[eintLocal].nbl_lists.nbl[0],
935 wallcycle_stop(wcycle, ewcNS);
939 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
940 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
941 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, FALSE, x,
942 nbv->grp[eintLocal].nbat);
943 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
944 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
949 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
950 /* launch local nonbonded F on GPU */
951 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFNo,
953 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
956 /* Communicate coordinates and sum dipole if necessary +
957 do non-local pair search */
958 if (DOMAINDECOMP(cr))
960 bDiffKernels = (nbv->grp[eintNonlocal].kernel_type !=
961 nbv->grp[eintLocal].kernel_type);
965 /* With GPU+CPU non-bonded calculations we need to copy
966 * the local coordinates to the non-local nbat struct
967 * (in CPU format) as the non-local kernel call also
968 * calculates the local - non-local interactions.
970 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
971 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
972 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatLocal, TRUE, x,
973 nbv->grp[eintNonlocal].nbat);
974 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
975 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
980 wallcycle_start_nocount(wcycle, ewcNS);
981 wallcycle_sub_start(wcycle, ewcsNBS_SEARCH_NONLOCAL);
985 nbnxn_grid_add_simple(nbv->nbs, nbv->grp[eintNonlocal].nbat);
988 nbnxn_make_pairlist(nbv->nbs, nbv->grp[eintNonlocal].nbat,
991 nbv->min_ci_balanced,
992 &nbv->grp[eintNonlocal].nbl_lists,
994 nbv->grp[eintNonlocal].kernel_type,
997 wallcycle_sub_stop(wcycle, ewcsNBS_SEARCH_NONLOCAL);
999 if (nbv->grp[eintNonlocal].kernel_type == nbnxnk8x8x8_CUDA)
1001 /* initialize non-local pair-list on the GPU */
1002 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
1003 nbv->grp[eintNonlocal].nbl_lists.nbl[0],
1006 wallcycle_stop(wcycle, ewcNS);
1010 wallcycle_start(wcycle, ewcMOVEX);
1011 dd_move_x(cr->dd, box, x);
1013 /* When we don't need the total dipole we sum it in global_stat */
1014 if (bStateChanged && NEED_MUTOT(*inputrec))
1016 gmx_sumd(2*DIM, mu, cr);
1018 wallcycle_stop(wcycle, ewcMOVEX);
1020 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1021 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
1022 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs, eatNonlocal, FALSE, x,
1023 nbv->grp[eintNonlocal].nbat);
1024 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
1025 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1028 if (bUseGPU && !bDiffKernels)
1030 wallcycle_start(wcycle, ewcLAUNCH_GPU_NB);
1031 /* launch non-local nonbonded F on GPU */
1032 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFNo,
1034 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1040 /* launch D2H copy-back F */
1041 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1042 if (DOMAINDECOMP(cr) && !bDiffKernels)
1044 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintNonlocal].nbat,
1045 flags, eatNonlocal);
1047 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintLocal].nbat,
1049 cycles_force += wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1052 if (bStateChanged && NEED_MUTOT(*inputrec))
1056 gmx_sumd(2*DIM, mu, cr);
1059 for (i = 0; i < 2; i++)
1061 for (j = 0; j < DIM; j++)
1063 fr->mu_tot[i][j] = mu[i*DIM + j];
1067 if (fr->efep == efepNO)
1069 copy_rvec(fr->mu_tot[0], mu_tot);
1073 for (j = 0; j < DIM; j++)
1076 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] +
1077 lambda[efptCOUL]*fr->mu_tot[1][j];
1081 /* Reset energies */
1082 reset_enerdata(&(inputrec->opts), fr, bNS, enerd, MASTER(cr));
1083 clear_rvecs(SHIFTS, fr->fshift);
1085 if (DOMAINDECOMP(cr))
1087 if (!(cr->duty & DUTY_PME))
1089 wallcycle_start(wcycle, ewcPPDURINGPME);
1090 dd_force_flop_start(cr->dd, nrnb);
1096 /* Enforced rotation has its own cycle counter that starts after the collective
1097 * coordinates have been communicated. It is added to ddCyclF to allow
1098 * for proper load-balancing */
1099 wallcycle_start(wcycle, ewcROT);
1100 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1101 wallcycle_stop(wcycle, ewcROT);
1104 /* Start the force cycle counter.
1105 * This counter is stopped in do_forcelow_level.
1106 * No parallel communication should occur while this counter is running,
1107 * since that will interfere with the dynamic load balancing.
1109 wallcycle_start(wcycle, ewcFORCE);
1112 /* Reset forces for which the virial is calculated separately:
1113 * PME/Ewald forces if necessary */
1114 if (fr->bF_NoVirSum)
1116 if (flags & GMX_FORCE_VIRIAL)
1118 fr->f_novirsum = fr->f_novirsum_alloc;
1121 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1125 clear_rvecs(homenr, fr->f_novirsum+start);
1130 /* We are not calculating the pressure so we do not need
1131 * a separate array for forces that do not contribute
1138 /* Clear the short- and long-range forces */
1139 clear_rvecs(fr->natoms_force_constr, f);
1140 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1142 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1145 clear_rvec(fr->vir_diag_posres);
1148 if (inputrec->ePull == epullCONSTRAINT)
1150 clear_pull_forces(inputrec->pull);
1153 /* We calculate the non-bonded forces, when done on the CPU, here.
1154 * We do this before calling do_force_lowlevel, as in there bondeds
1155 * forces are calculated before PME, which does communication.
1156 * With this order, non-bonded and bonded force calculation imbalance
1157 * can be balanced out by the domain decomposition load balancing.
1162 /* Maybe we should move this into do_force_lowlevel */
1163 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFYes,
1167 if (!bUseOrEmulGPU || bDiffKernels)
1171 if (DOMAINDECOMP(cr))
1173 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal,
1174 bDiffKernels ? enbvClearFYes : enbvClearFNo,
1184 aloc = eintNonlocal;
1187 /* Add all the non-bonded force to the normal force array.
1188 * This can be split into a local a non-local part when overlapping
1189 * communication with calculation with domain decomposition.
1191 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1192 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1193 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1194 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatAll, nbv->grp[aloc].nbat, f);
1195 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1196 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1197 wallcycle_start_nocount(wcycle, ewcFORCE);
1199 /* if there are multiple fshift output buffers reduce them */
1200 if ((flags & GMX_FORCE_VIRIAL) &&
1201 nbv->grp[aloc].nbl_lists.nnbl > 1)
1203 nbnxn_atomdata_add_nbat_fshift_to_fshift(nbv->grp[aloc].nbat,
1208 /* update QMMMrec, if necessary */
1211 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1214 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1216 posres_wrapper(fplog, flags, bSepDVDL, inputrec, nrnb, top, box, x,
1217 f, enerd, lambda, fr);
1220 /* Compute the bonded and non-bonded energies and optionally forces */
1221 do_force_lowlevel(fplog, step, fr, inputrec, &(top->idef),
1222 cr, nrnb, wcycle, mdatoms, &(inputrec->opts),
1223 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, mtop, top, fr->born,
1224 &(top->atomtypes), bBornRadii, box,
1225 inputrec->fepvals, lambda, graph, &(top->excls), fr->mu_tot,
1226 flags, &cycles_pme);
1230 if (do_per_step(step, inputrec->nstcalclr))
1232 /* Add the long range forces to the short range forces */
1233 for (i = 0; i < fr->natoms_force_constr; i++)
1235 rvec_add(fr->f_twin[i], f[i], f[i]);
1240 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1244 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1247 if (bUseOrEmulGPU && !bDiffKernels)
1249 /* wait for non-local forces (or calculate in emulation mode) */
1250 if (DOMAINDECOMP(cr))
1254 wallcycle_start(wcycle, ewcWAIT_GPU_NB_NL);
1255 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1256 nbv->grp[eintNonlocal].nbat,
1258 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1260 cycles_force += wallcycle_stop(wcycle, ewcWAIT_GPU_NB_NL);
1264 wallcycle_start_nocount(wcycle, ewcFORCE);
1265 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFYes,
1267 cycles_force += wallcycle_stop(wcycle, ewcFORCE);
1269 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1270 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1271 /* skip the reduction if there was no non-local work to do */
1272 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1274 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatNonlocal,
1275 nbv->grp[eintNonlocal].nbat, f);
1277 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1278 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1284 /* Communicate the forces */
1287 wallcycle_start(wcycle, ewcMOVEF);
1288 if (DOMAINDECOMP(cr))
1290 dd_move_f(cr->dd, f, fr->fshift);
1291 /* Do we need to communicate the separate force array
1292 * for terms that do not contribute to the single sum virial?
1293 * Position restraints and electric fields do not introduce
1294 * inter-cg forces, only full electrostatics methods do.
1295 * When we do not calculate the virial, fr->f_novirsum = f,
1296 * so we have already communicated these forces.
1298 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1299 (flags & GMX_FORCE_VIRIAL))
1301 dd_move_f(cr->dd, fr->f_novirsum, NULL);
1305 /* We should not update the shift forces here,
1306 * since f_twin is already included in f.
1308 dd_move_f(cr->dd, fr->f_twin, NULL);
1311 wallcycle_stop(wcycle, ewcMOVEF);
1317 /* wait for local forces (or calculate in emulation mode) */
1320 wallcycle_start(wcycle, ewcWAIT_GPU_NB_L);
1321 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1322 nbv->grp[eintLocal].nbat,
1324 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1326 wallcycle_stop(wcycle, ewcWAIT_GPU_NB_L);
1328 /* now clear the GPU outputs while we finish the step on the CPU */
1330 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1331 nbnxn_cuda_clear_outputs(nbv->cu_nbv, flags);
1332 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
1336 wallcycle_start_nocount(wcycle, ewcFORCE);
1337 do_nb_verlet(fr, ic, enerd, flags, eintLocal,
1338 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
1340 wallcycle_stop(wcycle, ewcFORCE);
1342 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1343 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1344 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1346 /* skip the reduction if there was no non-local work to do */
1347 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs, eatLocal,
1348 nbv->grp[eintLocal].nbat, f);
1350 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1351 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1354 if (DOMAINDECOMP(cr))
1356 dd_force_flop_stop(cr->dd, nrnb);
1359 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1365 if (IR_ELEC_FIELD(*inputrec))
1367 /* Compute forces due to electric field */
1368 calc_f_el(MASTER(cr) ? field : NULL,
1369 start, homenr, mdatoms->chargeA, x, fr->f_novirsum,
1370 inputrec->ex, inputrec->et, t);
1373 /* If we have NoVirSum forces, but we do not calculate the virial,
1374 * we sum fr->f_novirum=f later.
1376 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1378 wallcycle_start(wcycle, ewcVSITESPREAD);
1379 spread_vsite_f(fplog, vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1380 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1381 wallcycle_stop(wcycle, ewcVSITESPREAD);
1385 wallcycle_start(wcycle, ewcVSITESPREAD);
1386 spread_vsite_f(fplog, vsite, x, fr->f_twin, NULL, FALSE, NULL,
1388 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1389 wallcycle_stop(wcycle, ewcVSITESPREAD);
1393 if (flags & GMX_FORCE_VIRIAL)
1395 /* Calculation of the virial must be done after vsites! */
1396 calc_virial(fplog, mdatoms->start, mdatoms->homenr, x, f,
1397 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1401 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1403 pull_potential_wrapper(fplog, bSepDVDL, cr, inputrec, box, x,
1404 f, vir_force, mdatoms, enerd, lambda, t);
1407 /* Add the forces from enforced rotation potentials (if any) */
1410 wallcycle_start(wcycle, ewcROTadd);
1411 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
1412 wallcycle_stop(wcycle, ewcROTadd);
1415 if (PAR(cr) && !(cr->duty & DUTY_PME))
1417 /* In case of node-splitting, the PP nodes receive the long-range
1418 * forces, virial and energy from the PME nodes here.
1420 pme_receive_force_ener(fplog, bSepDVDL, cr, wcycle, enerd, fr);
1425 post_process_forces(fplog, cr, step, nrnb, wcycle,
1426 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
1430 /* Sum the potential energy terms from group contributions */
1431 sum_epot(&(inputrec->opts), &(enerd->grpp), enerd->term);
1434 void do_force_cutsGROUP(FILE *fplog, t_commrec *cr,
1435 t_inputrec *inputrec,
1436 gmx_large_int_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1437 gmx_localtop_t *top,
1439 gmx_groups_t *groups,
1440 matrix box, rvec x[], history_t *hist,
1444 gmx_enerdata_t *enerd, t_fcdata *fcd,
1445 real *lambda, t_graph *graph,
1446 t_forcerec *fr, gmx_vsite_t *vsite, rvec mu_tot,
1447 double t, FILE *field, gmx_edsam_t ed,
1448 gmx_bool bBornRadii,
1454 gmx_bool bSepDVDL, bStateChanged, bNS, bFillGrid, bCalcCGCM, bBS;
1455 gmx_bool bDoLongRangeNS, bDoForces, bDoPotential, bSepLRF;
1456 gmx_bool bDoAdressWF;
1458 rvec vzero, box_diag;
1459 real e, v, dvdlambda[efptNR];
1461 float cycles_pme, cycles_force;
1463 start = mdatoms->start;
1464 homenr = mdatoms->homenr;
1466 bSepDVDL = (fr->bSepDVDL && do_per_step(step, inputrec->nstlog));
1468 clear_mat(vir_force);
1472 pd_cg_range(cr, &cg0, &cg1);
1477 if (DOMAINDECOMP(cr))
1479 cg1 = cr->dd->ncg_tot;
1491 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
1492 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll == FALSE);
1493 /* Should we update the long-range neighborlists at this step? */
1494 bDoLongRangeNS = fr->bTwinRange && bNS;
1495 /* Should we perform the long-range nonbonded evaluation inside the neighborsearching? */
1496 bFillGrid = (bNS && bStateChanged);
1497 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
1498 bDoForces = (flags & GMX_FORCE_FORCES);
1499 bDoPotential = (flags & GMX_FORCE_ENERGY);
1500 bSepLRF = ((inputrec->nstcalclr > 1) && bDoForces &&
1501 (flags & GMX_FORCE_SEPLRF) && (flags & GMX_FORCE_DO_LR));
1503 /* should probably move this to the forcerec since it doesn't change */
1504 bDoAdressWF = ((fr->adress_type != eAdressOff));
1508 update_forcerec(fplog, fr, box);
1510 if (NEED_MUTOT(*inputrec))
1512 /* Calculate total (local) dipole moment in a temporary common array.
1513 * This makes it possible to sum them over nodes faster.
1515 calc_mu(start, homenr,
1516 x, mdatoms->chargeA, mdatoms->chargeB, mdatoms->nChargePerturbed,
1521 if (fr->ePBC != epbcNONE)
1523 /* Compute shift vectors every step,
1524 * because of pressure coupling or box deformation!
1526 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
1528 calc_shifts(box, fr->shift_vec);
1533 put_charge_groups_in_box(fplog, cg0, cg1, fr->ePBC, box,
1534 &(top->cgs), x, fr->cg_cm);
1535 inc_nrnb(nrnb, eNR_CGCM, homenr);
1536 inc_nrnb(nrnb, eNR_RESETX, cg1-cg0);
1538 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph)
1540 unshift_self(graph, box, x);
1545 calc_cgcm(fplog, cg0, cg1, &(top->cgs), x, fr->cg_cm);
1546 inc_nrnb(nrnb, eNR_CGCM, homenr);
1553 move_cgcm(fplog, cr, fr->cg_cm);
1557 pr_rvecs(debug, 0, "cgcm", fr->cg_cm, top->cgs.nr);
1562 if (!(cr->duty & DUTY_PME))
1564 /* Send particle coordinates to the pme nodes.
1565 * Since this is only implemented for domain decomposition
1566 * and domain decomposition does not use the graph,
1567 * we do not need to worry about shifting.
1570 wallcycle_start(wcycle, ewcPP_PMESENDX);
1572 bBS = (inputrec->nwall == 2);
1575 copy_mat(box, boxs);
1576 svmul(inputrec->wall_ewald_zfac, boxs[ZZ], boxs[ZZ]);
1579 gmx_pme_send_x(cr, bBS ? boxs : box, x,
1580 mdatoms->nChargePerturbed, lambda[efptCOUL],
1581 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)), step);
1583 wallcycle_stop(wcycle, ewcPP_PMESENDX);
1585 #endif /* GMX_MPI */
1587 /* Communicate coordinates and sum dipole if necessary */
1590 wallcycle_start(wcycle, ewcMOVEX);
1591 if (DOMAINDECOMP(cr))
1593 dd_move_x(cr->dd, box, x);
1597 move_x(fplog, cr, GMX_LEFT, GMX_RIGHT, x, nrnb);
1599 wallcycle_stop(wcycle, ewcMOVEX);
1602 /* update adress weight beforehand */
1603 if (bStateChanged && bDoAdressWF)
1605 /* need pbc for adress weight calculation with pbc_dx */
1606 set_pbc(&pbc, inputrec->ePBC, box);
1607 if (fr->adress_site == eAdressSITEcog)
1609 update_adress_weights_cog(top->idef.iparams, top->idef.il, x, fr, mdatoms,
1610 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1612 else if (fr->adress_site == eAdressSITEcom)
1614 update_adress_weights_com(fplog, cg0, cg1, &(top->cgs), x, fr, mdatoms,
1615 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1617 else if (fr->adress_site == eAdressSITEatomatom)
1619 update_adress_weights_atom_per_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1620 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1624 update_adress_weights_atom(cg0, cg1, &(top->cgs), x, fr, mdatoms,
1625 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1629 if (NEED_MUTOT(*inputrec))
1636 gmx_sumd(2*DIM, mu, cr);
1638 for (i = 0; i < 2; i++)
1640 for (j = 0; j < DIM; j++)
1642 fr->mu_tot[i][j] = mu[i*DIM + j];
1646 if (fr->efep == efepNO)
1648 copy_rvec(fr->mu_tot[0], mu_tot);
1652 for (j = 0; j < DIM; j++)
1655 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] + lambda[efptCOUL]*fr->mu_tot[1][j];
1660 /* Reset energies */
1661 reset_enerdata(&(inputrec->opts), fr, bNS, enerd, MASTER(cr));
1662 clear_rvecs(SHIFTS, fr->fshift);
1666 wallcycle_start(wcycle, ewcNS);
1668 if (graph && bStateChanged)
1670 /* Calculate intramolecular shift vectors to make molecules whole */
1671 mk_mshift(fplog, graph, fr->ePBC, box, x);
1674 /* Do the actual neighbour searching and if twin range electrostatics
1675 * also do the calculation of long range forces and energies.
1677 for (i = 0; i < efptNR; i++)
1681 ns(fplog, fr, x, box,
1682 groups, &(inputrec->opts), top, mdatoms,
1683 cr, nrnb, lambda, dvdlambda, &enerd->grpp, bFillGrid,
1686 wallcycle_stop(wcycle, ewcNS);
1689 if (inputrec->implicit_solvent && bNS)
1691 make_gb_nblist(cr, inputrec->gb_algorithm, inputrec->rlist,
1692 x, box, fr, &top->idef, graph, fr->born);
1695 if (DOMAINDECOMP(cr))
1697 if (!(cr->duty & DUTY_PME))
1699 wallcycle_start(wcycle, ewcPPDURINGPME);
1700 dd_force_flop_start(cr->dd, nrnb);
1706 /* Enforced rotation has its own cycle counter that starts after the collective
1707 * coordinates have been communicated. It is added to ddCyclF to allow
1708 * for proper load-balancing */
1709 wallcycle_start(wcycle, ewcROT);
1710 do_rotation(cr, inputrec, box, x, t, step, wcycle, bNS);
1711 wallcycle_stop(wcycle, ewcROT);
1714 /* Start the force cycle counter.
1715 * This counter is stopped in do_forcelow_level.
1716 * No parallel communication should occur while this counter is running,
1717 * since that will interfere with the dynamic load balancing.
1719 wallcycle_start(wcycle, ewcFORCE);
1723 /* Reset forces for which the virial is calculated separately:
1724 * PME/Ewald forces if necessary */
1725 if (fr->bF_NoVirSum)
1727 if (flags & GMX_FORCE_VIRIAL)
1729 fr->f_novirsum = fr->f_novirsum_alloc;
1732 clear_rvecs(fr->f_novirsum_n, fr->f_novirsum);
1736 clear_rvecs(homenr, fr->f_novirsum+start);
1741 /* We are not calculating the pressure so we do not need
1742 * a separate array for forces that do not contribute
1749 /* Clear the short- and long-range forces */
1750 clear_rvecs(fr->natoms_force_constr, f);
1751 if (bSepLRF && do_per_step(step, inputrec->nstcalclr))
1753 clear_rvecs(fr->natoms_force_constr, fr->f_twin);
1756 clear_rvec(fr->vir_diag_posres);
1758 if (inputrec->ePull == epullCONSTRAINT)
1760 clear_pull_forces(inputrec->pull);
1763 /* update QMMMrec, if necessary */
1766 update_QMMMrec(cr, fr, x, mdatoms, box, top);
1769 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1771 posres_wrapper(fplog, flags, bSepDVDL, inputrec, nrnb, top, box, x,
1772 f, enerd, lambda, fr);
1775 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_FBPOSRES].nr > 0)
1777 /* Flat-bottomed position restraints always require full pbc */
1778 if (!(bStateChanged && bDoAdressWF))
1780 set_pbc(&pbc, inputrec->ePBC, box);
1782 v = fbposres(top->idef.il[F_FBPOSRES].nr, top->idef.il[F_FBPOSRES].iatoms,
1783 top->idef.iparams_fbposres,
1784 (const rvec*)x, fr->f_novirsum, fr->vir_diag_posres,
1785 inputrec->ePBC == epbcNONE ? NULL : &pbc,
1786 fr->rc_scaling, fr->ePBC, fr->posres_com);
1787 enerd->term[F_FBPOSRES] += v;
1788 inc_nrnb(nrnb, eNR_FBPOSRES, top->idef.il[F_FBPOSRES].nr/2);
1791 /* Compute the bonded and non-bonded energies and optionally forces */
1792 do_force_lowlevel(fplog, step, fr, inputrec, &(top->idef),
1793 cr, nrnb, wcycle, mdatoms, &(inputrec->opts),
1794 x, hist, f, bSepLRF ? fr->f_twin : f, enerd, fcd, mtop, top, fr->born,
1795 &(top->atomtypes), bBornRadii, box,
1796 inputrec->fepvals, lambda,
1797 graph, &(top->excls), fr->mu_tot,
1803 if (do_per_step(step, inputrec->nstcalclr))
1805 /* Add the long range forces to the short range forces */
1806 for (i = 0; i < fr->natoms_force_constr; i++)
1808 rvec_add(fr->f_twin[i], f[i], f[i]);
1813 cycles_force = wallcycle_stop(wcycle, ewcFORCE);
1817 do_flood(cr, inputrec, x, f, ed, box, step, bNS);
1820 if (DOMAINDECOMP(cr))
1822 dd_force_flop_stop(cr->dd, nrnb);
1825 dd_cycles_add(cr->dd, cycles_force-cycles_pme, ddCyclF);
1831 if (IR_ELEC_FIELD(*inputrec))
1833 /* Compute forces due to electric field */
1834 calc_f_el(MASTER(cr) ? field : NULL,
1835 start, homenr, mdatoms->chargeA, x, fr->f_novirsum,
1836 inputrec->ex, inputrec->et, t);
1839 if (bDoAdressWF && fr->adress_icor == eAdressICThermoForce)
1841 /* Compute thermodynamic force in hybrid AdResS region */
1842 adress_thermo_force(start, homenr, &(top->cgs), x, fr->f_novirsum, fr, mdatoms,
1843 inputrec->ePBC == epbcNONE ? NULL : &pbc);
1846 /* Communicate the forces */
1849 wallcycle_start(wcycle, ewcMOVEF);
1850 if (DOMAINDECOMP(cr))
1852 dd_move_f(cr->dd, f, fr->fshift);
1853 /* Do we need to communicate the separate force array
1854 * for terms that do not contribute to the single sum virial?
1855 * Position restraints and electric fields do not introduce
1856 * inter-cg forces, only full electrostatics methods do.
1857 * When we do not calculate the virial, fr->f_novirsum = f,
1858 * so we have already communicated these forces.
1860 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1861 (flags & GMX_FORCE_VIRIAL))
1863 dd_move_f(cr->dd, fr->f_novirsum, NULL);
1867 /* We should not update the shift forces here,
1868 * since f_twin is already included in f.
1870 dd_move_f(cr->dd, fr->f_twin, NULL);
1875 pd_move_f(cr, f, nrnb);
1878 pd_move_f(cr, fr->f_twin, nrnb);
1881 wallcycle_stop(wcycle, ewcMOVEF);
1884 /* If we have NoVirSum forces, but we do not calculate the virial,
1885 * we sum fr->f_novirum=f later.
1887 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1889 wallcycle_start(wcycle, ewcVSITESPREAD);
1890 spread_vsite_f(fplog, vsite, x, f, fr->fshift, FALSE, NULL, nrnb,
1891 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1892 wallcycle_stop(wcycle, ewcVSITESPREAD);
1896 wallcycle_start(wcycle, ewcVSITESPREAD);
1897 spread_vsite_f(fplog, vsite, x, fr->f_twin, NULL, FALSE, NULL,
1899 &top->idef, fr->ePBC, fr->bMolPBC, graph, box, cr);
1900 wallcycle_stop(wcycle, ewcVSITESPREAD);
1904 if (flags & GMX_FORCE_VIRIAL)
1906 /* Calculation of the virial must be done after vsites! */
1907 calc_virial(fplog, mdatoms->start, mdatoms->homenr, x, f,
1908 vir_force, graph, box, nrnb, fr, inputrec->ePBC);
1912 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1914 pull_potential_wrapper(fplog, bSepDVDL, cr, inputrec, box, x,
1915 f, vir_force, mdatoms, enerd, lambda, t);
1918 /* Add the forces from enforced rotation potentials (if any) */
1921 wallcycle_start(wcycle, ewcROTadd);
1922 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
1923 wallcycle_stop(wcycle, ewcROTadd);
1926 if (PAR(cr) && !(cr->duty & DUTY_PME))
1928 /* In case of node-splitting, the PP nodes receive the long-range
1929 * forces, virial and energy from the PME nodes here.
1931 pme_receive_force_ener(fplog, bSepDVDL, cr, wcycle, enerd, fr);
1936 post_process_forces(fplog, cr, step, nrnb, wcycle,
1937 top, box, x, f, vir_force, mdatoms, graph, fr, vsite,
1941 /* Sum the potential energy terms from group contributions */
1942 sum_epot(&(inputrec->opts), &(enerd->grpp), enerd->term);
1945 void do_force(FILE *fplog, t_commrec *cr,
1946 t_inputrec *inputrec,
1947 gmx_large_int_t step, t_nrnb *nrnb, gmx_wallcycle_t wcycle,
1948 gmx_localtop_t *top,
1950 gmx_groups_t *groups,
1951 matrix box, rvec x[], history_t *hist,
1955 gmx_enerdata_t *enerd, t_fcdata *fcd,
1956 real *lambda, t_graph *graph,
1958 gmx_vsite_t *vsite, rvec mu_tot,
1959 double t, FILE *field, gmx_edsam_t ed,
1960 gmx_bool bBornRadii,
1963 /* modify force flag if not doing nonbonded */
1964 if (!fr->bNonbonded)
1966 flags &= ~GMX_FORCE_NONBONDED;
1969 switch (inputrec->cutoff_scheme)
1972 do_force_cutsVERLET(fplog, cr, inputrec,
1988 do_force_cutsGROUP(fplog, cr, inputrec,
2003 gmx_incons("Invalid cut-off scheme passed!");
2008 void do_constrain_first(FILE *fplog, gmx_constr_t constr,
2009 t_inputrec *ir, t_mdatoms *md,
2010 t_state *state, rvec *f,
2011 t_graph *graph, t_commrec *cr, t_nrnb *nrnb,
2012 t_forcerec *fr, gmx_localtop_t *top, tensor shake_vir)
2014 int i, m, start, end;
2015 gmx_large_int_t step;
2016 real dt = ir->delta_t;
2020 snew(savex, state->natoms);
2023 end = md->homenr + start;
2027 fprintf(debug, "vcm: start=%d, homenr=%d, end=%d\n",
2028 start, md->homenr, end);
2030 /* Do a first constrain to reset particles... */
2031 step = ir->init_step;
2034 char buf[STEPSTRSIZE];
2035 fprintf(fplog, "\nConstraining the starting coordinates (step %s)\n",
2036 gmx_step_str(step, buf));
2040 /* constrain the current position */
2041 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2042 ir, NULL, cr, step, 0, md,
2043 state->x, state->x, NULL,
2044 fr->bMolPBC, state->box,
2045 state->lambda[efptBONDED], &dvdl_dum,
2046 NULL, NULL, nrnb, econqCoord,
2047 ir->epc == epcMTTK, state->veta, state->veta);
2050 /* constrain the inital velocity, and save it */
2051 /* also may be useful if we need the ekin from the halfstep for velocity verlet */
2052 /* might not yet treat veta correctly */
2053 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2054 ir, NULL, cr, step, 0, md,
2055 state->x, state->v, state->v,
2056 fr->bMolPBC, state->box,
2057 state->lambda[efptBONDED], &dvdl_dum,
2058 NULL, NULL, nrnb, econqVeloc,
2059 ir->epc == epcMTTK, state->veta, state->veta);
2061 /* constrain the inital velocities at t-dt/2 */
2062 if (EI_STATE_VELOCITY(ir->eI) && ir->eI != eiVV)
2064 for (i = start; (i < end); i++)
2066 for (m = 0; (m < DIM); m++)
2068 /* Reverse the velocity */
2069 state->v[i][m] = -state->v[i][m];
2070 /* Store the position at t-dt in buf */
2071 savex[i][m] = state->x[i][m] + dt*state->v[i][m];
2074 /* Shake the positions at t=-dt with the positions at t=0
2075 * as reference coordinates.
2079 char buf[STEPSTRSIZE];
2080 fprintf(fplog, "\nConstraining the coordinates at t0-dt (step %s)\n",
2081 gmx_step_str(step, buf));
2084 constrain(NULL, TRUE, FALSE, constr, &(top->idef),
2085 ir, NULL, cr, step, -1, md,
2086 state->x, savex, NULL,
2087 fr->bMolPBC, state->box,
2088 state->lambda[efptBONDED], &dvdl_dum,
2089 state->v, NULL, nrnb, econqCoord,
2090 ir->epc == epcMTTK, state->veta, state->veta);
2092 for (i = start; i < end; i++)
2094 for (m = 0; m < DIM; m++)
2096 /* Re-reverse the velocities */
2097 state->v[i][m] = -state->v[i][m];
2104 void calc_enervirdiff(FILE *fplog, int eDispCorr, t_forcerec *fr)
2106 double eners[2], virs[2], enersum, virsum, y0, f, g, h;
2107 double r0, r1, r, rc3, rc9, ea, eb, ec, pa, pb, pc, pd;
2108 double invscale, invscale2, invscale3;
2109 int ri0, ri1, ri, i, offstart, offset;
2110 real scale, *vdwtab, tabfactor, tmp;
2112 fr->enershiftsix = 0;
2113 fr->enershifttwelve = 0;
2114 fr->enerdiffsix = 0;
2115 fr->enerdifftwelve = 0;
2117 fr->virdifftwelve = 0;
2119 if (eDispCorr != edispcNO)
2121 for (i = 0; i < 2; i++)
2126 if ((fr->vdwtype == evdwSWITCH) || (fr->vdwtype == evdwSHIFT))
2128 if (fr->rvdw_switch == 0)
2131 "With dispersion correction rvdw-switch can not be zero "
2132 "for vdw-type = %s", evdw_names[fr->vdwtype]);
2135 scale = fr->nblists[0].table_elec_vdw.scale;
2136 vdwtab = fr->nblists[0].table_vdw.data;
2138 /* Round the cut-offs to exact table values for precision */
2139 ri0 = floor(fr->rvdw_switch*scale);
2140 ri1 = ceil(fr->rvdw*scale);
2146 if (fr->vdwtype == evdwSHIFT)
2148 /* Determine the constant energy shift below rvdw_switch.
2149 * Table has a scale factor since we have scaled it down to compensate
2150 * for scaling-up c6/c12 with the derivative factors to save flops in analytical kernels.
2152 fr->enershiftsix = (real)(-1.0/(rc3*rc3)) - 6.0*vdwtab[8*ri0];
2153 fr->enershifttwelve = (real)( 1.0/(rc9*rc3)) - 12.0*vdwtab[8*ri0 + 4];
2155 /* Add the constant part from 0 to rvdw_switch.
2156 * This integration from 0 to rvdw_switch overcounts the number
2157 * of interactions by 1, as it also counts the self interaction.
2158 * We will correct for this later.
2160 eners[0] += 4.0*M_PI*fr->enershiftsix*rc3/3.0;
2161 eners[1] += 4.0*M_PI*fr->enershifttwelve*rc3/3.0;
2163 invscale = 1.0/(scale);
2164 invscale2 = invscale*invscale;
2165 invscale3 = invscale*invscale2;
2167 /* following summation derived from cubic spline definition,
2168 Numerical Recipies in C, second edition, p. 113-116. Exact
2169 for the cubic spline. We first calculate the negative of
2170 the energy from rvdw to rvdw_switch, assuming that g(r)=1,
2171 and then add the more standard, abrupt cutoff correction to
2172 that result, yielding the long-range correction for a
2173 switched function. We perform both the pressure and energy
2174 loops at the same time for simplicity, as the computational
2177 for (i = 0; i < 2; i++)
2179 enersum = 0.0; virsum = 0.0;
2183 /* Since the dispersion table has been scaled down a factor 6.0 and the repulsion
2184 * a factor 12.0 to compensate for the c6/c12 parameters inside nbfp[] being scaled
2185 * up (to save flops in kernels), we need to correct for this.
2194 for (ri = ri0; ri < ri1; ri++)
2198 eb = 2.0*invscale2*r;
2202 pb = 3.0*invscale2*r;
2203 pc = 3.0*invscale*r*r;
2206 /* this "8" is from the packing in the vdwtab array - perhaps should be #define'ed? */
2207 offset = 8*ri + offstart;
2208 y0 = vdwtab[offset];
2209 f = vdwtab[offset+1];
2210 g = vdwtab[offset+2];
2211 h = vdwtab[offset+3];
2213 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);
2214 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);
2217 enersum *= 4.0*M_PI*tabfactor;
2218 virsum *= 4.0*M_PI*tabfactor;
2219 eners[i] -= enersum;
2223 /* now add the correction for rvdw_switch to infinity */
2224 eners[0] += -4.0*M_PI/(3.0*rc3);
2225 eners[1] += 4.0*M_PI/(9.0*rc9);
2226 virs[0] += 8.0*M_PI/rc3;
2227 virs[1] += -16.0*M_PI/(3.0*rc9);
2229 else if ((fr->vdwtype == evdwCUT) || (fr->vdwtype == evdwUSER))
2231 if (fr->vdwtype == evdwUSER && fplog)
2234 "WARNING: using dispersion correction with user tables\n");
2236 rc3 = fr->rvdw*fr->rvdw*fr->rvdw;
2238 /* Contribution beyond the cut-off */
2239 eners[0] += -4.0*M_PI/(3.0*rc3);
2240 eners[1] += 4.0*M_PI/(9.0*rc9);
2241 if (fr->vdw_modifier == eintmodPOTSHIFT)
2243 /* Contribution within the cut-off */
2244 eners[0] += -4.0*M_PI/(3.0*rc3);
2245 eners[1] += 4.0*M_PI/(3.0*rc9);
2247 /* Contribution beyond the cut-off */
2248 virs[0] += 8.0*M_PI/rc3;
2249 virs[1] += -16.0*M_PI/(3.0*rc9);
2254 "Dispersion correction is not implemented for vdw-type = %s",
2255 evdw_names[fr->vdwtype]);
2257 fr->enerdiffsix = eners[0];
2258 fr->enerdifftwelve = eners[1];
2259 /* The 0.5 is due to the Gromacs definition of the virial */
2260 fr->virdiffsix = 0.5*virs[0];
2261 fr->virdifftwelve = 0.5*virs[1];
2265 void calc_dispcorr(FILE *fplog, t_inputrec *ir, t_forcerec *fr,
2266 gmx_large_int_t step, int natoms,
2267 matrix box, real lambda, tensor pres, tensor virial,
2268 real *prescorr, real *enercorr, real *dvdlcorr)
2270 gmx_bool bCorrAll, bCorrPres;
2271 real dvdlambda, invvol, dens, ninter, avcsix, avctwelve, enerdiff, svir = 0, spres = 0;
2281 if (ir->eDispCorr != edispcNO)
2283 bCorrAll = (ir->eDispCorr == edispcAllEner ||
2284 ir->eDispCorr == edispcAllEnerPres);
2285 bCorrPres = (ir->eDispCorr == edispcEnerPres ||
2286 ir->eDispCorr == edispcAllEnerPres);
2288 invvol = 1/det(box);
2291 /* Only correct for the interactions with the inserted molecule */
2292 dens = (natoms - fr->n_tpi)*invvol;
2297 dens = natoms*invvol;
2298 ninter = 0.5*natoms;
2301 if (ir->efep == efepNO)
2303 avcsix = fr->avcsix[0];
2304 avctwelve = fr->avctwelve[0];
2308 avcsix = (1 - lambda)*fr->avcsix[0] + lambda*fr->avcsix[1];
2309 avctwelve = (1 - lambda)*fr->avctwelve[0] + lambda*fr->avctwelve[1];
2312 enerdiff = ninter*(dens*fr->enerdiffsix - fr->enershiftsix);
2313 *enercorr += avcsix*enerdiff;
2315 if (ir->efep != efepNO)
2317 dvdlambda += (fr->avcsix[1] - fr->avcsix[0])*enerdiff;
2321 enerdiff = ninter*(dens*fr->enerdifftwelve - fr->enershifttwelve);
2322 *enercorr += avctwelve*enerdiff;
2323 if (fr->efep != efepNO)
2325 dvdlambda += (fr->avctwelve[1] - fr->avctwelve[0])*enerdiff;
2331 svir = ninter*dens*avcsix*fr->virdiffsix/3.0;
2332 if (ir->eDispCorr == edispcAllEnerPres)
2334 svir += ninter*dens*avctwelve*fr->virdifftwelve/3.0;
2336 /* The factor 2 is because of the Gromacs virial definition */
2337 spres = -2.0*invvol*svir*PRESFAC;
2339 for (m = 0; m < DIM; m++)
2341 virial[m][m] += svir;
2342 pres[m][m] += spres;
2347 /* Can't currently control when it prints, for now, just print when degugging */
2352 fprintf(debug, "Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
2358 "Long Range LJ corr.: Epot %10g, Pres: %10g, Vir: %10g\n",
2359 *enercorr, spres, svir);
2363 fprintf(debug, "Long Range LJ corr.: Epot %10g\n", *enercorr);
2367 if (fr->bSepDVDL && do_per_step(step, ir->nstlog))
2369 gmx_print_sepdvdl(fplog, "Dispersion correction", *enercorr, dvdlambda);
2371 if (fr->efep != efepNO)
2373 *dvdlcorr += dvdlambda;
2378 void do_pbc_first(FILE *fplog, matrix box, t_forcerec *fr,
2379 t_graph *graph, rvec x[])
2383 fprintf(fplog, "Removing pbc first time\n");
2385 calc_shifts(box, fr->shift_vec);
2388 mk_mshift(fplog, graph, fr->ePBC, box, x);
2391 p_graph(debug, "do_pbc_first 1", graph);
2393 shift_self(graph, box, x);
2394 /* By doing an extra mk_mshift the molecules that are broken
2395 * because they were e.g. imported from another software
2396 * will be made whole again. Such are the healing powers
2399 mk_mshift(fplog, graph, fr->ePBC, box, x);
2402 p_graph(debug, "do_pbc_first 2", graph);
2407 fprintf(fplog, "Done rmpbc\n");
2411 static void low_do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2412 gmx_mtop_t *mtop, rvec x[],
2417 gmx_molblock_t *molb;
2419 if (bFirst && fplog)
2421 fprintf(fplog, "Removing pbc first time\n");
2426 for (mb = 0; mb < mtop->nmolblock; mb++)
2428 molb = &mtop->molblock[mb];
2429 if (molb->natoms_mol == 1 ||
2430 (!bFirst && mtop->moltype[molb->type].cgs.nr == 1))
2432 /* Just one atom or charge group in the molecule, no PBC required */
2433 as += molb->nmol*molb->natoms_mol;
2437 /* Pass NULL iso fplog to avoid graph prints for each molecule type */
2438 mk_graph_ilist(NULL, mtop->moltype[molb->type].ilist,
2439 0, molb->natoms_mol, FALSE, FALSE, graph);
2441 for (mol = 0; mol < molb->nmol; mol++)
2443 mk_mshift(fplog, graph, ePBC, box, x+as);
2445 shift_self(graph, box, x+as);
2446 /* The molecule is whole now.
2447 * We don't need the second mk_mshift call as in do_pbc_first,
2448 * since we no longer need this graph.
2451 as += molb->natoms_mol;
2459 void do_pbc_first_mtop(FILE *fplog, int ePBC, matrix box,
2460 gmx_mtop_t *mtop, rvec x[])
2462 low_do_pbc_mtop(fplog, ePBC, box, mtop, x, TRUE);
2465 void do_pbc_mtop(FILE *fplog, int ePBC, matrix box,
2466 gmx_mtop_t *mtop, rvec x[])
2468 low_do_pbc_mtop(fplog, ePBC, box, mtop, x, FALSE);
2471 void finish_run(FILE *fplog, t_commrec *cr, const char *confout,
2472 t_inputrec *inputrec,
2473 t_nrnb nrnb[], gmx_wallcycle_t wcycle,
2474 gmx_runtime_t *runtime,
2475 wallclock_gpu_t *gputimes,
2477 gmx_bool bWriteStat)
2480 t_nrnb *nrnb_tot = NULL;
2484 wallcycle_sum(cr, wcycle);
2490 MPI_Allreduce(nrnb->n, nrnb_tot->n, eNRNB, MPI_DOUBLE, MPI_SUM,
2491 cr->mpi_comm_mysim);
2499 #if defined(GMX_MPI) && !defined(GMX_THREAD_MPI)
2502 /* reduce nodetime over all MPI processes in the current simulation */
2504 MPI_Allreduce(&runtime->proctime, &sum, 1, MPI_DOUBLE, MPI_SUM,
2505 cr->mpi_comm_mysim);
2506 runtime->proctime = sum;
2512 print_flop(fplog, nrnb_tot, &nbfs, &mflop);
2519 if ((cr->duty & DUTY_PP) && DOMAINDECOMP(cr))
2521 print_dd_statistics(cr, inputrec, fplog);
2533 snew(nrnb_all, cr->nnodes);
2534 nrnb_all[0] = *nrnb;
2535 for (s = 1; s < cr->nnodes; s++)
2537 MPI_Recv(nrnb_all[s].n, eNRNB, MPI_DOUBLE, s, 0,
2538 cr->mpi_comm_mysim, &stat);
2540 pr_load(fplog, cr, nrnb_all);
2545 MPI_Send(nrnb->n, eNRNB, MPI_DOUBLE, MASTERRANK(cr), 0,
2546 cr->mpi_comm_mysim);
2553 wallcycle_print(fplog, cr->nnodes, cr->npmenodes, runtime->realtime,
2556 if (EI_DYNAMICS(inputrec->eI))
2558 delta_t = inputrec->delta_t;
2567 print_perf(fplog, runtime->proctime, runtime->realtime,
2568 runtime->nsteps_done, delta_t, nbfs, mflop);
2572 print_perf(stderr, runtime->proctime, runtime->realtime,
2573 runtime->nsteps_done, delta_t, nbfs, mflop);
2578 extern void initialize_lambdas(FILE *fplog, t_inputrec *ir, int *fep_state, real *lambda, double *lam0)
2580 /* this function works, but could probably use a logic rewrite to keep all the different
2581 types of efep straight. */
2584 t_lambda *fep = ir->fepvals;
2586 if ((ir->efep == efepNO) && (ir->bSimTemp == FALSE))
2588 for (i = 0; i < efptNR; i++)
2600 *fep_state = fep->init_fep_state; /* this might overwrite the checkpoint
2601 if checkpoint is set -- a kludge is in for now
2603 for (i = 0; i < efptNR; i++)
2605 /* overwrite lambda state with init_lambda for now for backwards compatibility */
2606 if (fep->init_lambda >= 0) /* if it's -1, it was never initializd */
2608 lambda[i] = fep->init_lambda;
2611 lam0[i] = lambda[i];
2616 lambda[i] = fep->all_lambda[i][*fep_state];
2619 lam0[i] = lambda[i];
2625 /* need to rescale control temperatures to match current state */
2626 for (i = 0; i < ir->opts.ngtc; i++)
2628 if (ir->opts.ref_t[i] > 0)
2630 ir->opts.ref_t[i] = ir->simtempvals->temperatures[*fep_state];
2636 /* Send to the log the information on the current lambdas */
2639 fprintf(fplog, "Initial vector of lambda components:[ ");
2640 for (i = 0; i < efptNR; i++)
2642 fprintf(fplog, "%10.4f ", lambda[i]);
2644 fprintf(fplog, "]\n");
2650 void init_md(FILE *fplog,
2651 t_commrec *cr, t_inputrec *ir, const output_env_t oenv,
2652 double *t, double *t0,
2653 real *lambda, int *fep_state, double *lam0,
2654 t_nrnb *nrnb, gmx_mtop_t *mtop,
2656 int nfile, const t_filenm fnm[],
2657 gmx_mdoutf_t **outf, t_mdebin **mdebin,
2658 tensor force_vir, tensor shake_vir, rvec mu_tot,
2659 gmx_bool *bSimAnn, t_vcm **vcm, t_state *state, unsigned long Flags)
2664 /* Initial values */
2665 *t = *t0 = ir->init_t;
2668 for (i = 0; i < ir->opts.ngtc; i++)
2670 /* set bSimAnn if any group is being annealed */
2671 if (ir->opts.annealing[i] != eannNO)
2678 update_annealing_target_temp(&(ir->opts), ir->init_t);
2681 /* Initialize lambda variables */
2682 initialize_lambdas(fplog, ir, fep_state, lambda, lam0);
2686 *upd = init_update(fplog, ir);
2692 *vcm = init_vcm(fplog, &mtop->groups, ir);
2695 if (EI_DYNAMICS(ir->eI) && !(Flags & MD_APPENDFILES))
2697 if (ir->etc == etcBERENDSEN)
2699 please_cite(fplog, "Berendsen84a");
2701 if (ir->etc == etcVRESCALE)
2703 please_cite(fplog, "Bussi2007a");
2711 *outf = init_mdoutf(nfile, fnm, Flags, cr, ir, oenv);
2713 *mdebin = init_mdebin((Flags & MD_APPENDFILES) ? NULL : (*outf)->fp_ene,
2714 mtop, ir, (*outf)->fp_dhdl);
2719 please_cite(fplog, "Fritsch12");
2720 please_cite(fplog, "Junghans10");
2722 /* Initiate variables */
2723 clear_mat(force_vir);
2724 clear_mat(shake_vir);