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
51 #include "visibility.h"
61 #include "chargegroup.h"
84 #include "pull_rotation.h"
85 #include "gmx_random.h"
86 #include "mpelogging.h"
89 #include "gmx_wallcycle.h"
91 #include "nbnxn_atomdata.h"
92 #include "nbnxn_search.h"
93 #include "nbnxn_kernels/nbnxn_kernel_ref.h"
94 #include "nbnxn_kernels/nbnxn_kernel_x86_simd128.h"
95 #include "nbnxn_kernels/nbnxn_kernel_x86_simd256.h"
96 #include "nbnxn_kernels/nbnxn_kernel_gpu_ref.h"
101 #ifdef GMX_THREAD_MPI
108 #include "nbnxn_cuda_data_mgmt.h"
109 #include "nbnxn_cuda/nbnxn_cuda.h"
112 typedef struct gmx_timeprint {
117 /* Portable version of ctime_r implemented in src/gmxlib/string2.c, but we do not want it declared in public installed headers */
120 gmx_ctime_r(const time_t *clock,char *buf, int n);
126 #ifdef HAVE_GETTIMEOFDAY
130 gettimeofday(&t,NULL);
132 seconds = (double) t.tv_sec + 1e-6*(double)t.tv_usec;
138 seconds = time(NULL);
145 #define difftime(end,start) ((double)(end)-(double)(start))
147 void print_time(FILE *out,gmx_runtime_t *runtime,gmx_large_int_t step,
148 t_inputrec *ir, t_commrec *cr)
151 char timebuf[STRLEN];
155 #ifndef GMX_THREAD_MPI
161 fprintf(out,"step %s",gmx_step_str(step,buf));
162 if ((step >= ir->nstlist))
164 runtime->last = gmx_gettime();
165 dt = difftime(runtime->last,runtime->real);
166 runtime->time_per_step = dt/(step - ir->init_step + 1);
168 dt = (ir->nsteps + ir->init_step - step)*runtime->time_per_step;
174 finish = (time_t) (runtime->last + dt);
175 gmx_ctime_r(&finish,timebuf,STRLEN);
176 sprintf(buf,"%s",timebuf);
177 buf[strlen(buf)-1]='\0';
178 fprintf(out,", will finish %s",buf);
181 fprintf(out,", remaining runtime: %5d s ",(int)dt);
185 fprintf(out," performance: %.1f ns/day ",
186 ir->delta_t/1000*24*60*60/runtime->time_per_step);
189 #ifndef GMX_THREAD_MPI
203 static double set_proctime(gmx_runtime_t *runtime)
209 prev = runtime->proc;
210 runtime->proc = dclock();
212 diff = runtime->proc - prev;
216 prev = runtime->proc;
217 runtime->proc = clock();
219 diff = (double)(runtime->proc - prev)/(double)CLOCKS_PER_SEC;
223 /* The counter has probably looped, ignore this data */
230 void runtime_start(gmx_runtime_t *runtime)
232 runtime->real = gmx_gettime();
234 set_proctime(runtime);
235 runtime->realtime = 0;
236 runtime->proctime = 0;
238 runtime->time_per_step = 0;
241 void runtime_end(gmx_runtime_t *runtime)
247 runtime->proctime += set_proctime(runtime);
248 runtime->realtime = now - runtime->real;
252 void runtime_upd_proc(gmx_runtime_t *runtime)
254 runtime->proctime += set_proctime(runtime);
257 void print_date_and_time(FILE *fplog,int nodeid,const char *title,
258 const gmx_runtime_t *runtime)
261 char timebuf[STRLEN];
262 char time_string[STRLEN];
269 tmptime = (time_t) runtime->real;
270 gmx_ctime_r(&tmptime,timebuf,STRLEN);
274 tmptime = (time_t) gmx_gettime();
275 gmx_ctime_r(&tmptime,timebuf,STRLEN);
277 for(i=0; timebuf[i]>=' '; i++)
279 time_string[i]=timebuf[i];
283 fprintf(fplog,"%s on node %d %s\n",title,nodeid,time_string);
287 static void sum_forces(int start,int end,rvec f[],rvec flr[])
292 pr_rvecs(debug,0,"fsr",f+start,end-start);
293 pr_rvecs(debug,0,"flr",flr+start,end-start);
295 for(i=start; (i<end); i++)
296 rvec_inc(f[i],flr[i]);
300 * calc_f_el calculates forces due to an electric field.
302 * force is kJ mol^-1 nm^-1 = e * kJ mol^-1 nm^-1 / e
304 * Et[] contains the parameters for the time dependent
305 * part of the field (not yet used).
306 * Ex[] contains the parameters for
307 * the spatial dependent part of the field. You can have cool periodic
308 * fields in principle, but only a constant field is supported
310 * The function should return the energy due to the electric field
311 * (if any) but for now returns 0.
314 * There can be problems with the virial.
315 * Since the field is not self-consistent this is unavoidable.
316 * For neutral molecules the virial is correct within this approximation.
317 * For neutral systems with many charged molecules the error is small.
318 * But for systems with a net charge or a few charged molecules
319 * the error can be significant when the field is high.
320 * Solution: implement a self-consitent electric field into PME.
322 static void calc_f_el(FILE *fp,int start,int homenr,
323 real charge[],rvec x[],rvec f[],
324 t_cosines Ex[],t_cosines Et[],double t)
330 for(m=0; (m<DIM); m++)
337 Ext[m] = cos(Et[m].a[0]*(t-t0))*exp(-sqr(t-t0)/(2.0*sqr(Et[m].a[2])));
341 Ext[m] = cos(Et[m].a[0]*t);
350 /* Convert the field strength from V/nm to MD-units */
351 Ext[m] *= Ex[m].a[0]*FIELDFAC;
352 for(i=start; (i<start+homenr); i++)
353 f[i][m] += charge[i]*Ext[m];
362 fprintf(fp,"%10g %10g %10g %10g #FIELD\n",t,
363 Ext[XX]/FIELDFAC,Ext[YY]/FIELDFAC,Ext[ZZ]/FIELDFAC);
367 static void calc_virial(FILE *fplog,int start,int homenr,rvec x[],rvec f[],
368 tensor vir_part,t_graph *graph,matrix box,
369 t_nrnb *nrnb,const t_forcerec *fr,int ePBC)
374 /* The short-range virial from surrounding boxes */
376 calc_vir(fplog,SHIFTS,fr->shift_vec,fr->fshift,vir_part,ePBC==epbcSCREW,box);
377 inc_nrnb(nrnb,eNR_VIRIAL,SHIFTS);
379 /* Calculate partial virial, for local atoms only, based on short range.
380 * Total virial is computed in global_stat, called from do_md
382 f_calc_vir(fplog,start,start+homenr,x,f,vir_part,graph,box);
383 inc_nrnb(nrnb,eNR_VIRIAL,homenr);
385 /* Add position restraint contribution */
386 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++) {
392 vir_part[i][ZZ] += fr->vir_wall_z[i];
396 pr_rvecs(debug,0,"vir_part",vir_part,DIM);
399 static void posres_wrapper(FILE *fplog,
407 gmx_enerdata_t *enerd,
415 /* Position restraints always require full pbc */
416 set_pbc(&pbc,ir->ePBC,box);
418 v = posres(top->idef.il[F_POSRES].nr,top->idef.il[F_POSRES].iatoms,
419 top->idef.iparams_posres,
420 (const rvec*)x,fr->f_novirsum,fr->vir_diag_posres,
421 ir->ePBC==epbcNONE ? NULL : &pbc,
422 lambda[efptRESTRAINT],&dvdl,
423 fr->rc_scaling,fr->ePBC,fr->posres_com,fr->posres_comB);
426 fprintf(fplog,sepdvdlformat,
427 interaction_function[F_POSRES].longname,v,dvdl);
429 enerd->term[F_POSRES] += v;
430 /* If just the force constant changes, the FEP term is linear,
431 * but if k changes, it is not.
433 enerd->dvdl_nonlin[efptRESTRAINT] += dvdl;
434 inc_nrnb(nrnb,eNR_POSRES,top->idef.il[F_POSRES].nr/2);
436 if ((ir->fepvals->n_lambda > 0) && (flags & GMX_FORCE_DHDL))
438 for(i=0; i<enerd->n_lambda; i++)
440 real dvdl_dum,lambda_dum;
442 lambda_dum = (i==0 ? lambda[efptRESTRAINT] : ir->fepvals->all_lambda[efptRESTRAINT][i-1]);
443 v = posres(top->idef.il[F_POSRES].nr,top->idef.il[F_POSRES].iatoms,
444 top->idef.iparams_posres,
445 (const rvec*)x,NULL,NULL,
446 ir->ePBC==epbcNONE ? NULL : &pbc,lambda_dum,&dvdl,
447 fr->rc_scaling,fr->ePBC,fr->posres_com,fr->posres_comB);
448 enerd->enerpart_lambda[i] += v;
453 static void pull_potential_wrapper(FILE *fplog,
461 gmx_enerdata_t *enerd,
468 /* Calculate the center of mass forces, this requires communication,
469 * which is why pull_potential is called close to other communication.
470 * The virial contribution is calculated directly,
471 * which is why we call pull_potential after calc_virial.
473 set_pbc(&pbc,ir->ePBC,box);
475 enerd->term[F_COM_PULL] +=
476 pull_potential(ir->ePull,ir->pull,mdatoms,&pbc,
477 cr,t,lambda[efptRESTRAINT],x,f,vir_force,&dvdl);
480 fprintf(fplog,sepdvdlformat,"Com pull",enerd->term[F_COM_PULL],dvdl);
482 enerd->dvdl_lin[efptRESTRAINT] += dvdl;
485 static void pme_receive_force_ener(FILE *fplog,
488 gmx_wallcycle_t wcycle,
489 gmx_enerdata_t *enerd,
493 float cycles_ppdpme,cycles_seppme;
495 cycles_ppdpme = wallcycle_stop(wcycle,ewcPPDURINGPME);
496 dd_cycles_add(cr->dd,cycles_ppdpme,ddCyclPPduringPME);
498 /* In case of node-splitting, the PP nodes receive the long-range
499 * forces, virial and energy from the PME nodes here.
501 wallcycle_start(wcycle,ewcPP_PMEWAITRECVF);
503 gmx_pme_receive_f(cr,fr->f_novirsum,fr->vir_el_recip,&e,&dvdl,
507 fprintf(fplog,sepdvdlformat,"PME mesh",e,dvdl);
509 enerd->term[F_COUL_RECIP] += e;
510 enerd->dvdl_lin[efptCOUL] += dvdl;
513 dd_cycles_add(cr->dd,cycles_seppme,ddCyclPME);
515 wallcycle_stop(wcycle,ewcPP_PMEWAITRECVF);
518 static void print_large_forces(FILE *fp,t_mdatoms *md,t_commrec *cr,
519 gmx_large_int_t step,real pforce,rvec *x,rvec *f)
523 char buf[STEPSTRSIZE];
526 for(i=md->start; i<md->start+md->homenr; i++) {
528 /* We also catch NAN, if the compiler does not optimize this away. */
529 if (fn2 >= pf2 || fn2 != fn2) {
530 fprintf(fp,"step %s atom %6d x %8.3f %8.3f %8.3f force %12.5e\n",
531 gmx_step_str(step,buf),
532 ddglatnr(cr->dd,i),x[i][XX],x[i][YY],x[i][ZZ],sqrt(fn2));
537 static void post_process_forces(FILE *fplog,
539 gmx_large_int_t step,
540 t_nrnb *nrnb,gmx_wallcycle_t wcycle,
547 t_forcerec *fr,gmx_vsite_t *vsite,
554 /* Spread the mesh force on virtual sites to the other particles...
555 * This is parallellized. MPI communication is performed
556 * if the constructing atoms aren't local.
558 wallcycle_start(wcycle,ewcVSITESPREAD);
559 spread_vsite_f(fplog,vsite,x,fr->f_novirsum,NULL,
560 (flags & GMX_FORCE_VIRIAL),fr->vir_el_recip,
562 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
563 wallcycle_stop(wcycle,ewcVSITESPREAD);
565 if (flags & GMX_FORCE_VIRIAL)
567 /* Now add the forces, this is local */
570 sum_forces(0,fr->f_novirsum_n,f,fr->f_novirsum);
574 sum_forces(mdatoms->start,mdatoms->start+mdatoms->homenr,
577 if (EEL_FULL(fr->eeltype))
579 /* Add the mesh contribution to the virial */
580 m_add(vir_force,fr->vir_el_recip,vir_force);
584 pr_rvecs(debug,0,"vir_force",vir_force,DIM);
589 if (fr->print_force >= 0)
591 print_large_forces(stderr,mdatoms,cr,step,fr->print_force,x,f);
595 static void do_nb_verlet(t_forcerec *fr,
596 interaction_const_t *ic,
597 gmx_enerdata_t *enerd,
598 int flags, int ilocality,
601 gmx_wallcycle_t wcycle)
603 int nnbl, kernel_type, enr_nbnxn_kernel_ljc, enr_nbnxn_kernel_lj;
605 nonbonded_verlet_group_t *nbvg;
607 if (!(flags & GMX_FORCE_NONBONDED))
609 /* skip non-bonded calculation */
613 nbvg = &fr->nbv->grp[ilocality];
615 /* CUDA kernel launch overhead is already timed separately */
616 if (fr->cutoff_scheme != ecutsVERLET)
618 gmx_incons("Invalid cut-off scheme passed!");
621 if (nbvg->kernel_type != nbk8x8x8_CUDA)
623 wallcycle_sub_start(wcycle, ewcsNONBONDED);
625 switch (nbvg->kernel_type)
628 nbnxn_kernel_ref(&nbvg->nbl_lists,
634 enerd->grpp.ener[egCOULSR],
636 enerd->grpp.ener[egBHAMSR] :
637 enerd->grpp.ener[egLJSR]);
640 case nbk4xN_X86_SIMD128:
641 nbnxn_kernel_x86_simd128(&nbvg->nbl_lists,
648 enerd->grpp.ener[egCOULSR],
650 enerd->grpp.ener[egBHAMSR] :
651 enerd->grpp.ener[egLJSR]);
653 case nbk4xN_X86_SIMD256:
654 nbnxn_kernel_x86_simd256(&nbvg->nbl_lists,
661 enerd->grpp.ener[egCOULSR],
663 enerd->grpp.ener[egBHAMSR] :
664 enerd->grpp.ener[egLJSR]);
668 nbnxn_cuda_launch_kernel(fr->nbv->cu_nbv, nbvg->nbat, flags, ilocality);
671 case nbk8x8x8_PlainC:
672 nbnxn_kernel_gpu_ref(nbvg->nbl_lists.nbl[0],
677 nbvg->nbat->out[0].f,
679 enerd->grpp.ener[egCOULSR],
681 enerd->grpp.ener[egBHAMSR] :
682 enerd->grpp.ener[egLJSR]);
686 gmx_incons("Invalid nonbonded kernel type passed!");
689 if (nbvg->kernel_type != nbk8x8x8_CUDA)
691 wallcycle_sub_stop(wcycle, ewcsNONBONDED);
694 if (EEL_RF(ic->eeltype) || ic->eeltype == eelCUT)
696 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_RF;
698 else if (nbvg->ewald_excl == ewaldexclTable)
700 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_TAB;
704 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_EWALD;
706 enr_nbnxn_kernel_lj = eNR_NBNXN_LJ;
707 if (flags & GMX_FORCE_ENERGY)
709 /* In eNR_??? the nbnxn F+E kernels are always the F kernel + 1 */
710 enr_nbnxn_kernel_ljc += 1;
711 enr_nbnxn_kernel_lj += 1;
714 inc_nrnb(nrnb,enr_nbnxn_kernel_ljc,
715 nbvg->nbl_lists.natpair_ljq);
716 inc_nrnb(nrnb,enr_nbnxn_kernel_lj,
717 nbvg->nbl_lists.natpair_lj);
718 inc_nrnb(nrnb,enr_nbnxn_kernel_ljc-eNR_NBNXN_LJ_RF+eNR_NBNXN_RF,
719 nbvg->nbl_lists.natpair_q);
722 void do_force_cutsVERLET(FILE *fplog,t_commrec *cr,
723 t_inputrec *inputrec,
724 gmx_large_int_t step,t_nrnb *nrnb,gmx_wallcycle_t wcycle,
727 gmx_groups_t *groups,
728 matrix box,rvec x[],history_t *hist,
732 gmx_enerdata_t *enerd,t_fcdata *fcd,
733 real *lambda,t_graph *graph,
734 t_forcerec *fr, interaction_const_t *ic,
735 gmx_vsite_t *vsite,rvec mu_tot,
736 double t,FILE *field,gmx_edsam_t ed,
744 gmx_bool bSepDVDL,bStateChanged,bNS,bFillGrid,bCalcCGCM,bBS;
745 gmx_bool bDoLongRange,bDoForces,bSepLRF,bUseGPU,bUseOrEmulGPU;
746 gmx_bool bDiffKernels=FALSE;
750 float cycles_pme,cycles_force;
751 nonbonded_verlet_t *nbv;
755 nb_kernel_type = fr->nbv->grp[0].kernel_type;
757 start = mdatoms->start;
758 homenr = mdatoms->homenr;
760 bSepDVDL = (fr->bSepDVDL && do_per_step(step,inputrec->nstlog));
762 clear_mat(vir_force);
765 if (DOMAINDECOMP(cr))
767 cg1 = cr->dd->ncg_tot;
778 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
779 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll==FALSE);
780 bFillGrid = (bNS && bStateChanged);
781 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
782 bDoLongRange = (fr->bTwinRange && bNS && (flags & GMX_FORCE_DO_LR));
783 bDoForces = (flags & GMX_FORCE_FORCES);
784 bSepLRF = (bDoLongRange && bDoForces && (flags & GMX_FORCE_SEPLRF));
785 bUseGPU = fr->nbv->bUseGPU;
786 bUseOrEmulGPU = bUseGPU || (nbv->grp[0].kernel_type == nbk8x8x8_PlainC);
790 update_forcerec(fplog,fr,box);
792 if (NEED_MUTOT(*inputrec))
794 /* Calculate total (local) dipole moment in a temporary common array.
795 * This makes it possible to sum them over nodes faster.
797 calc_mu(start,homenr,
798 x,mdatoms->chargeA,mdatoms->chargeB,mdatoms->nChargePerturbed,
803 if (fr->ePBC != epbcNONE) {
804 /* Compute shift vectors every step,
805 * because of pressure coupling or box deformation!
807 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
808 calc_shifts(box,fr->shift_vec);
811 put_atoms_in_box_omp(fr->ePBC,box,homenr,x);
812 inc_nrnb(nrnb,eNR_SHIFTX,homenr);
814 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph) {
815 unshift_self(graph,box,x);
819 nbnxn_atomdata_copy_shiftvec(flags & GMX_FORCE_DYNAMICBOX,
820 fr->shift_vec,nbv->grp[0].nbat);
823 if (!(cr->duty & DUTY_PME)) {
824 /* Send particle coordinates to the pme nodes.
825 * Since this is only implemented for domain decomposition
826 * and domain decomposition does not use the graph,
827 * we do not need to worry about shifting.
830 wallcycle_start(wcycle,ewcPP_PMESENDX);
831 GMX_MPE_LOG(ev_send_coordinates_start);
833 bBS = (inputrec->nwall == 2);
836 svmul(inputrec->wall_ewald_zfac,boxs[ZZ],boxs[ZZ]);
839 gmx_pme_send_x(cr,bBS ? boxs : box,x,
840 mdatoms->nChargePerturbed,lambda[efptCOUL],
841 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),step);
843 GMX_MPE_LOG(ev_send_coordinates_finish);
844 wallcycle_stop(wcycle,ewcPP_PMESENDX);
848 /* do gridding for pair search */
851 if (graph && bStateChanged)
853 /* Calculate intramolecular shift vectors to make molecules whole */
854 mk_mshift(fplog,graph,fr->ePBC,box,x);
858 box_diag[XX] = box[XX][XX];
859 box_diag[YY] = box[YY][YY];
860 box_diag[ZZ] = box[ZZ][ZZ];
862 wallcycle_start(wcycle,ewcNS);
865 wallcycle_sub_start(wcycle,ewcsNBS_GRID_LOCAL);
866 nbnxn_put_on_grid(nbv->nbs,fr->ePBC,box,
868 0,mdatoms->homenr,-1,fr->cginfo,x,
870 nbv->grp[eintLocal].kernel_type,
871 nbv->grp[eintLocal].nbat);
872 wallcycle_sub_stop(wcycle,ewcsNBS_GRID_LOCAL);
876 wallcycle_sub_start(wcycle,ewcsNBS_GRID_NONLOCAL);
877 nbnxn_put_on_grid_nonlocal(nbv->nbs,domdec_zones(cr->dd),
879 nbv->grp[eintNonlocal].kernel_type,
880 nbv->grp[eintNonlocal].nbat);
881 wallcycle_sub_stop(wcycle,ewcsNBS_GRID_NONLOCAL);
884 if (nbv->ngrp == 1 ||
885 nbv->grp[eintNonlocal].nbat == nbv->grp[eintLocal].nbat)
887 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat,eatAll,
888 nbv->nbs,mdatoms,fr->cginfo);
892 nbnxn_atomdata_set(nbv->grp[eintLocal].nbat,eatLocal,
893 nbv->nbs,mdatoms,fr->cginfo);
894 nbnxn_atomdata_set(nbv->grp[eintNonlocal].nbat,eatAll,
895 nbv->nbs,mdatoms,fr->cginfo);
897 wallcycle_stop(wcycle, ewcNS);
900 /* initialize the GPU atom data and copy shift vector */
905 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
906 nbnxn_cuda_init_atomdata(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
907 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
910 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
911 nbnxn_cuda_upload_shiftvec(nbv->cu_nbv, nbv->grp[eintLocal].nbat);
912 wallcycle_stop(wcycle, ewcLAUNCH_GPU_NB);
915 /* do local pair search */
918 wallcycle_start_nocount(wcycle,ewcNS);
919 wallcycle_sub_start(wcycle,ewcsNBS_SEARCH_LOCAL);
920 nbnxn_make_pairlist(nbv->nbs,nbv->grp[eintLocal].nbat,
923 nbv->min_ci_balanced,
924 &nbv->grp[eintLocal].nbl_lists,
926 nbv->grp[eintLocal].kernel_type,
928 wallcycle_sub_stop(wcycle,ewcsNBS_SEARCH_LOCAL);
932 /* initialize local pair-list on the GPU */
933 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
934 nbv->grp[eintLocal].nbl_lists.nbl[0],
937 wallcycle_stop(wcycle, ewcNS);
941 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
942 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
943 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs,eatLocal,FALSE,x,
944 nbv->grp[eintLocal].nbat);
945 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
946 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
951 wallcycle_start(wcycle,ewcLAUNCH_GPU_NB);
952 /* launch local nonbonded F on GPU */
953 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFNo,
955 wallcycle_stop(wcycle,ewcLAUNCH_GPU_NB);
958 /* Communicate coordinates and sum dipole if necessary +
959 do non-local pair search */
960 if (DOMAINDECOMP(cr))
962 bDiffKernels = (nbv->grp[eintNonlocal].kernel_type !=
963 nbv->grp[eintLocal].kernel_type);
967 /* With GPU+CPU non-bonded calculations we need to copy
968 * the local coordinates to the non-local nbat struct
969 * (in CPU format) as the non-local kernel call also
970 * calculates the local - non-local interactions.
972 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
973 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
974 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs,eatLocal,TRUE,x,
975 nbv->grp[eintNonlocal].nbat);
976 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
977 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
982 wallcycle_start_nocount(wcycle,ewcNS);
983 wallcycle_sub_start(wcycle,ewcsNBS_SEARCH_NONLOCAL);
987 nbnxn_grid_add_simple(nbv->nbs,nbv->grp[eintNonlocal].nbat);
990 nbnxn_make_pairlist(nbv->nbs,nbv->grp[eintNonlocal].nbat,
993 nbv->min_ci_balanced,
994 &nbv->grp[eintNonlocal].nbl_lists,
996 nbv->grp[eintNonlocal].kernel_type,
999 wallcycle_sub_stop(wcycle,ewcsNBS_SEARCH_NONLOCAL);
1001 if (nbv->grp[eintNonlocal].kernel_type == nbk8x8x8_CUDA)
1003 /* initialize non-local pair-list on the GPU */
1004 nbnxn_cuda_init_pairlist(nbv->cu_nbv,
1005 nbv->grp[eintNonlocal].nbl_lists.nbl[0],
1008 wallcycle_stop(wcycle,ewcNS);
1012 wallcycle_start(wcycle,ewcMOVEX);
1013 dd_move_x(cr->dd,box,x);
1015 /* When we don't need the total dipole we sum it in global_stat */
1016 if (bStateChanged && NEED_MUTOT(*inputrec))
1018 gmx_sumd(2*DIM,mu,cr);
1020 wallcycle_stop(wcycle,ewcMOVEX);
1022 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1023 wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
1024 nbnxn_atomdata_copy_x_to_nbat_x(nbv->nbs,eatNonlocal,FALSE,x,
1025 nbv->grp[eintNonlocal].nbat);
1026 wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
1027 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1030 if (bUseGPU && !bDiffKernels)
1032 wallcycle_start(wcycle,ewcLAUNCH_GPU_NB);
1033 /* launch non-local nonbonded F on GPU */
1034 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFNo,
1036 cycles_force += wallcycle_stop(wcycle,ewcLAUNCH_GPU_NB);
1042 /* launch D2H copy-back F */
1043 wallcycle_start_nocount(wcycle, ewcLAUNCH_GPU_NB);
1044 if (DOMAINDECOMP(cr) && !bDiffKernels)
1046 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintNonlocal].nbat,
1047 flags, eatNonlocal);
1049 nbnxn_cuda_launch_cpyback(nbv->cu_nbv, nbv->grp[eintLocal].nbat,
1051 cycles_force += wallcycle_stop(wcycle,ewcLAUNCH_GPU_NB);
1054 if (bStateChanged && NEED_MUTOT(*inputrec))
1058 gmx_sumd(2*DIM,mu,cr);
1065 fr->mu_tot[i][j] = mu[i*DIM + j];
1069 if (fr->efep == efepNO)
1071 copy_rvec(fr->mu_tot[0],mu_tot);
1075 for(j=0; j<DIM; j++)
1078 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] +
1079 lambda[efptCOUL]*fr->mu_tot[1][j];
1083 /* Reset energies */
1084 reset_enerdata(&(inputrec->opts),fr,bNS,enerd,MASTER(cr));
1085 clear_rvecs(SHIFTS,fr->fshift);
1087 if (DOMAINDECOMP(cr))
1089 if (!(cr->duty & DUTY_PME))
1091 wallcycle_start(wcycle,ewcPPDURINGPME);
1092 dd_force_flop_start(cr->dd,nrnb);
1096 /* Start the force cycle counter.
1097 * This counter is stopped in do_forcelow_level.
1098 * No parallel communication should occur while this counter is running,
1099 * since that will interfere with the dynamic load balancing.
1101 wallcycle_start(wcycle,ewcFORCE);
1104 /* Reset forces for which the virial is calculated separately:
1105 * PME/Ewald forces if necessary */
1106 if (fr->bF_NoVirSum)
1108 if (flags & GMX_FORCE_VIRIAL)
1110 fr->f_novirsum = fr->f_novirsum_alloc;
1111 GMX_BARRIER(cr->mpi_comm_mygroup);
1114 clear_rvecs(fr->f_novirsum_n,fr->f_novirsum);
1118 clear_rvecs(homenr,fr->f_novirsum+start);
1120 GMX_BARRIER(cr->mpi_comm_mygroup);
1124 /* We are not calculating the pressure so we do not need
1125 * a separate array for forces that do not contribute
1132 /* Clear the short- and long-range forces */
1133 clear_rvecs(fr->natoms_force_constr,f);
1134 if(bSepLRF && do_per_step(step,inputrec->nstcalclr))
1136 clear_rvecs(fr->natoms_force_constr,fr->f_twin);
1139 clear_rvec(fr->vir_diag_posres);
1141 GMX_BARRIER(cr->mpi_comm_mygroup);
1143 if (inputrec->ePull == epullCONSTRAINT)
1145 clear_pull_forces(inputrec->pull);
1148 /* update QMMMrec, if necessary */
1151 update_QMMMrec(cr,fr,x,mdatoms,box,top);
1154 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1156 posres_wrapper(fplog,flags,bSepDVDL,inputrec,nrnb,top,box,x,
1160 /* Compute the bonded and non-bonded energies and optionally forces */
1161 do_force_lowlevel(fplog,step,fr,inputrec,&(top->idef),
1162 cr,nrnb,wcycle,mdatoms,&(inputrec->opts),
1163 x,hist,f, bSepLRF ? fr->f_twin : f,enerd,fcd,mtop,top,fr->born,
1164 &(top->atomtypes),bBornRadii,box,
1165 inputrec->fepvals,lambda,graph,&(top->excls),fr->mu_tot,
1166 flags, &cycles_pme);
1170 if (do_per_step(step,inputrec->nstcalclr))
1172 /* Add the long range forces to the short range forces */
1173 for(i=0; i<fr->natoms_force_constr; i++)
1175 rvec_add(fr->f_twin[i],f[i],f[i]);
1182 /* Maybe we should move this into do_force_lowlevel */
1183 do_nb_verlet(fr, ic, enerd, flags, eintLocal, enbvClearFYes,
1188 if (!bUseOrEmulGPU || bDiffKernels)
1192 if (DOMAINDECOMP(cr))
1194 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal,
1195 bDiffKernels ? enbvClearFYes : enbvClearFNo,
1205 aloc = eintNonlocal;
1208 /* Add all the non-bonded force to the normal force array.
1209 * This can be split into a local a non-local part when overlapping
1210 * communication with calculation with domain decomposition.
1212 cycles_force += wallcycle_stop(wcycle,ewcFORCE);
1213 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1214 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1215 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs,eatAll,nbv->grp[aloc].nbat,f);
1216 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1217 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1218 wallcycle_start_nocount(wcycle,ewcFORCE);
1220 /* if there are multiple fshift output buffers reduce them */
1221 if ((flags & GMX_FORCE_VIRIAL) &&
1222 nbv->grp[aloc].nbl_lists.nnbl > 1)
1224 nbnxn_atomdata_add_nbat_fshift_to_fshift(nbv->grp[aloc].nbat,
1229 cycles_force += wallcycle_stop(wcycle,ewcFORCE);
1230 GMX_BARRIER(cr->mpi_comm_mygroup);
1234 do_flood(fplog,cr,x,f,ed,box,step,bNS);
1237 if (bUseOrEmulGPU && !bDiffKernels)
1239 /* wait for non-local forces (or calculate in emulation mode) */
1240 if (DOMAINDECOMP(cr))
1244 wallcycle_start(wcycle,ewcWAIT_GPU_NB_NL);
1245 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1246 nbv->grp[eintNonlocal].nbat,
1248 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1250 cycles_force += wallcycle_stop(wcycle,ewcWAIT_GPU_NB_NL);
1254 wallcycle_start_nocount(wcycle,ewcFORCE);
1255 do_nb_verlet(fr, ic, enerd, flags, eintNonlocal, enbvClearFYes,
1257 cycles_force += wallcycle_stop(wcycle,ewcFORCE);
1259 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1260 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1261 /* skip the reduction if there was no non-local work to do */
1262 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1264 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs,eatNonlocal,
1265 nbv->grp[eintNonlocal].nbat,f);
1267 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1268 cycles_force += wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1274 /* Communicate the forces */
1277 wallcycle_start(wcycle,ewcMOVEF);
1278 if (DOMAINDECOMP(cr))
1280 dd_move_f(cr->dd,f,fr->fshift);
1281 /* Do we need to communicate the separate force array
1282 * for terms that do not contribute to the single sum virial?
1283 * Position restraints and electric fields do not introduce
1284 * inter-cg forces, only full electrostatics methods do.
1285 * When we do not calculate the virial, fr->f_novirsum = f,
1286 * so we have already communicated these forces.
1288 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1289 (flags & GMX_FORCE_VIRIAL))
1291 dd_move_f(cr->dd,fr->f_novirsum,NULL);
1295 /* We should not update the shift forces here,
1296 * since f_twin is already included in f.
1298 dd_move_f(cr->dd,fr->f_twin,NULL);
1301 wallcycle_stop(wcycle,ewcMOVEF);
1307 /* wait for local forces (or calculate in emulation mode) */
1310 wallcycle_start(wcycle,ewcWAIT_GPU_NB_L);
1311 nbnxn_cuda_wait_gpu(nbv->cu_nbv,
1312 nbv->grp[eintLocal].nbat,
1314 enerd->grpp.ener[egLJSR], enerd->grpp.ener[egCOULSR],
1316 wallcycle_stop(wcycle,ewcWAIT_GPU_NB_L);
1318 /* now clear the GPU outputs while we finish the step on the CPU */
1319 nbnxn_cuda_clear_outputs(nbv->cu_nbv, flags);
1323 wallcycle_start_nocount(wcycle,ewcFORCE);
1324 do_nb_verlet(fr, ic, enerd, flags, eintLocal,
1325 DOMAINDECOMP(cr) ? enbvClearFNo : enbvClearFYes,
1327 wallcycle_stop(wcycle,ewcFORCE);
1329 wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
1330 wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
1331 if (nbv->grp[eintLocal].nbl_lists.nbl[0]->nsci > 0)
1333 /* skip the reduction if there was no non-local work to do */
1334 nbnxn_atomdata_add_nbat_f_to_f(nbv->nbs,eatLocal,
1335 nbv->grp[eintLocal].nbat,f);
1337 wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
1338 wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
1341 if (DOMAINDECOMP(cr))
1343 dd_force_flop_stop(cr->dd,nrnb);
1346 dd_cycles_add(cr->dd,cycles_force-cycles_pme,ddCyclF);
1352 if (IR_ELEC_FIELD(*inputrec))
1354 /* Compute forces due to electric field */
1355 calc_f_el(MASTER(cr) ? field : NULL,
1356 start,homenr,mdatoms->chargeA,x,fr->f_novirsum,
1357 inputrec->ex,inputrec->et,t);
1360 /* If we have NoVirSum forces, but we do not calculate the virial,
1361 * we sum fr->f_novirum=f later.
1363 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1365 wallcycle_start(wcycle,ewcVSITESPREAD);
1366 spread_vsite_f(fplog,vsite,x,f,fr->fshift,FALSE,NULL,nrnb,
1367 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
1368 wallcycle_stop(wcycle,ewcVSITESPREAD);
1372 wallcycle_start(wcycle,ewcVSITESPREAD);
1373 spread_vsite_f(fplog,vsite,x,fr->f_twin,NULL,FALSE,NULL,
1375 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
1376 wallcycle_stop(wcycle,ewcVSITESPREAD);
1380 if (flags & GMX_FORCE_VIRIAL)
1382 /* Calculation of the virial must be done after vsites! */
1383 calc_virial(fplog,mdatoms->start,mdatoms->homenr,x,f,
1384 vir_force,graph,box,nrnb,fr,inputrec->ePBC);
1388 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1390 pull_potential_wrapper(fplog,bSepDVDL,cr,inputrec,box,x,
1391 f,vir_force,mdatoms,enerd,lambda,t);
1394 if (PAR(cr) && !(cr->duty & DUTY_PME))
1396 /* In case of node-splitting, the PP nodes receive the long-range
1397 * forces, virial and energy from the PME nodes here.
1399 pme_receive_force_ener(fplog,bSepDVDL,cr,wcycle,enerd,fr);
1404 post_process_forces(fplog,cr,step,nrnb,wcycle,
1405 top,box,x,f,vir_force,mdatoms,graph,fr,vsite,
1409 /* Sum the potential energy terms from group contributions */
1410 sum_epot(&(inputrec->opts),&(enerd->grpp),enerd->term);
1413 void do_force_cutsGROUP(FILE *fplog,t_commrec *cr,
1414 t_inputrec *inputrec,
1415 gmx_large_int_t step,t_nrnb *nrnb,gmx_wallcycle_t wcycle,
1416 gmx_localtop_t *top,
1418 gmx_groups_t *groups,
1419 matrix box,rvec x[],history_t *hist,
1423 gmx_enerdata_t *enerd,t_fcdata *fcd,
1424 real *lambda,t_graph *graph,
1425 t_forcerec *fr,gmx_vsite_t *vsite,rvec mu_tot,
1426 double t,FILE *field,gmx_edsam_t ed,
1427 gmx_bool bBornRadii,
1433 gmx_bool bSepDVDL,bStateChanged,bNS,bFillGrid,bCalcCGCM,bBS;
1434 gmx_bool bDoLongRangeNS,bDoForces,bDoPotential,bSepLRF;
1435 gmx_bool bDoAdressWF;
1437 rvec vzero,box_diag;
1438 real e,v,dvdlambda[efptNR];
1440 float cycles_pme,cycles_force;
1442 start = mdatoms->start;
1443 homenr = mdatoms->homenr;
1445 bSepDVDL = (fr->bSepDVDL && do_per_step(step,inputrec->nstlog));
1447 clear_mat(vir_force);
1451 pd_cg_range(cr,&cg0,&cg1);
1456 if (DOMAINDECOMP(cr))
1458 cg1 = cr->dd->ncg_tot;
1470 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
1471 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll==FALSE);
1472 /* Should we update the long-range neighborlists at this step? */
1473 bDoLongRangeNS = fr->bTwinRange && bNS;
1474 /* Should we perform the long-range nonbonded evaluation inside the neighborsearching? */
1475 bFillGrid = (bNS && bStateChanged);
1476 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
1477 bDoForces = (flags & GMX_FORCE_FORCES);
1478 bDoPotential = (flags & GMX_FORCE_ENERGY);
1479 bSepLRF = ((inputrec->nstcalclr>1) && bDoForces &&
1480 (flags & GMX_FORCE_SEPLRF) && (flags & GMX_FORCE_DO_LR));
1482 /* should probably move this to the forcerec since it doesn't change */
1483 bDoAdressWF = ((fr->adress_type!=eAdressOff));
1487 update_forcerec(fplog,fr,box);
1489 if (NEED_MUTOT(*inputrec))
1491 /* Calculate total (local) dipole moment in a temporary common array.
1492 * This makes it possible to sum them over nodes faster.
1494 calc_mu(start,homenr,
1495 x,mdatoms->chargeA,mdatoms->chargeB,mdatoms->nChargePerturbed,
1500 if (fr->ePBC != epbcNONE) {
1501 /* Compute shift vectors every step,
1502 * because of pressure coupling or box deformation!
1504 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
1505 calc_shifts(box,fr->shift_vec);
1508 put_charge_groups_in_box(fplog,cg0,cg1,fr->ePBC,box,
1509 &(top->cgs),x,fr->cg_cm);
1510 inc_nrnb(nrnb,eNR_CGCM,homenr);
1511 inc_nrnb(nrnb,eNR_RESETX,cg1-cg0);
1513 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph) {
1514 unshift_self(graph,box,x);
1517 else if (bCalcCGCM) {
1518 calc_cgcm(fplog,cg0,cg1,&(top->cgs),x,fr->cg_cm);
1519 inc_nrnb(nrnb,eNR_CGCM,homenr);
1524 move_cgcm(fplog,cr,fr->cg_cm);
1527 pr_rvecs(debug,0,"cgcm",fr->cg_cm,top->cgs.nr);
1531 if (!(cr->duty & DUTY_PME)) {
1532 /* Send particle coordinates to the pme nodes.
1533 * Since this is only implemented for domain decomposition
1534 * and domain decomposition does not use the graph,
1535 * we do not need to worry about shifting.
1538 wallcycle_start(wcycle,ewcPP_PMESENDX);
1539 GMX_MPE_LOG(ev_send_coordinates_start);
1541 bBS = (inputrec->nwall == 2);
1544 svmul(inputrec->wall_ewald_zfac,boxs[ZZ],boxs[ZZ]);
1547 gmx_pme_send_x(cr,bBS ? boxs : box,x,
1548 mdatoms->nChargePerturbed,lambda[efptCOUL],
1549 (flags & (GMX_FORCE_VIRIAL | GMX_FORCE_ENERGY)),step);
1551 GMX_MPE_LOG(ev_send_coordinates_finish);
1552 wallcycle_stop(wcycle,ewcPP_PMESENDX);
1554 #endif /* GMX_MPI */
1556 /* Communicate coordinates and sum dipole if necessary */
1559 wallcycle_start(wcycle,ewcMOVEX);
1560 if (DOMAINDECOMP(cr))
1562 dd_move_x(cr->dd,box,x);
1566 move_x(fplog,cr,GMX_LEFT,GMX_RIGHT,x,nrnb);
1568 wallcycle_stop(wcycle,ewcMOVEX);
1571 /* update adress weight beforehand */
1572 if(bStateChanged && bDoAdressWF)
1574 /* need pbc for adress weight calculation with pbc_dx */
1575 set_pbc(&pbc,inputrec->ePBC,box);
1576 if(fr->adress_site == eAdressSITEcog)
1578 update_adress_weights_cog(top->idef.iparams,top->idef.il,x,fr,mdatoms,
1579 inputrec->ePBC==epbcNONE ? NULL : &pbc);
1581 else if (fr->adress_site == eAdressSITEcom)
1583 update_adress_weights_com(fplog,cg0,cg1,&(top->cgs),x,fr,mdatoms,
1584 inputrec->ePBC==epbcNONE ? NULL : &pbc);
1586 else if (fr->adress_site == eAdressSITEatomatom){
1587 update_adress_weights_atom_per_atom(cg0,cg1,&(top->cgs),x,fr,mdatoms,
1588 inputrec->ePBC==epbcNONE ? NULL : &pbc);
1592 update_adress_weights_atom(cg0,cg1,&(top->cgs),x,fr,mdatoms,
1593 inputrec->ePBC==epbcNONE ? NULL : &pbc);
1597 if (NEED_MUTOT(*inputrec))
1604 gmx_sumd(2*DIM,mu,cr);
1610 fr->mu_tot[i][j] = mu[i*DIM + j];
1614 if (fr->efep == efepNO)
1616 copy_rvec(fr->mu_tot[0],mu_tot);
1620 for(j=0; j<DIM; j++)
1623 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] + lambda[efptCOUL]*fr->mu_tot[1][j];
1628 /* Reset energies */
1629 reset_enerdata(&(inputrec->opts),fr,bNS,enerd,MASTER(cr));
1630 clear_rvecs(SHIFTS,fr->fshift);
1634 wallcycle_start(wcycle,ewcNS);
1636 if (graph && bStateChanged)
1638 /* Calculate intramolecular shift vectors to make molecules whole */
1639 mk_mshift(fplog,graph,fr->ePBC,box,x);
1642 /* Do the actual neighbour searching and if twin range electrostatics
1643 * also do the calculation of long range forces and energies.
1645 for (i=0;i<efptNR;i++) {dvdlambda[i] = 0;}
1647 groups,&(inputrec->opts),top,mdatoms,
1648 cr,nrnb,lambda,dvdlambda,&enerd->grpp,bFillGrid,
1652 fprintf(fplog,sepdvdlformat,"LR non-bonded",0.0,dvdlambda);
1654 enerd->dvdl_lin[efptVDW] += dvdlambda[efptVDW];
1655 enerd->dvdl_lin[efptCOUL] += dvdlambda[efptCOUL];
1657 wallcycle_stop(wcycle,ewcNS);
1660 if (inputrec->implicit_solvent && bNS)
1662 make_gb_nblist(cr,inputrec->gb_algorithm,inputrec->rlist,
1663 x,box,fr,&top->idef,graph,fr->born);
1666 if (DOMAINDECOMP(cr))
1668 if (!(cr->duty & DUTY_PME))
1670 wallcycle_start(wcycle,ewcPPDURINGPME);
1671 dd_force_flop_start(cr->dd,nrnb);
1677 /* Enforced rotation has its own cycle counter that starts after the collective
1678 * coordinates have been communicated. It is added to ddCyclF to allow
1679 * for proper load-balancing */
1680 wallcycle_start(wcycle,ewcROT);
1681 do_rotation(cr,inputrec,box,x,t,step,wcycle,bNS);
1682 wallcycle_stop(wcycle,ewcROT);
1685 /* Start the force cycle counter.
1686 * This counter is stopped in do_forcelow_level.
1687 * No parallel communication should occur while this counter is running,
1688 * since that will interfere with the dynamic load balancing.
1690 wallcycle_start(wcycle,ewcFORCE);
1694 /* Reset forces for which the virial is calculated separately:
1695 * PME/Ewald forces if necessary */
1696 if (fr->bF_NoVirSum)
1698 if (flags & GMX_FORCE_VIRIAL)
1700 fr->f_novirsum = fr->f_novirsum_alloc;
1701 GMX_BARRIER(cr->mpi_comm_mygroup);
1704 clear_rvecs(fr->f_novirsum_n,fr->f_novirsum);
1708 clear_rvecs(homenr,fr->f_novirsum+start);
1710 GMX_BARRIER(cr->mpi_comm_mygroup);
1714 /* We are not calculating the pressure so we do not need
1715 * a separate array for forces that do not contribute
1722 /* Clear the short- and long-range forces */
1723 clear_rvecs(fr->natoms_force_constr,f);
1724 if(bSepLRF && do_per_step(step,inputrec->nstcalclr))
1726 clear_rvecs(fr->natoms_force_constr,fr->f_twin);
1729 clear_rvec(fr->vir_diag_posres);
1731 GMX_BARRIER(cr->mpi_comm_mygroup);
1733 if (inputrec->ePull == epullCONSTRAINT)
1735 clear_pull_forces(inputrec->pull);
1738 /* update QMMMrec, if necessary */
1741 update_QMMMrec(cr,fr,x,mdatoms,box,top);
1744 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
1746 posres_wrapper(fplog,flags,bSepDVDL,inputrec,nrnb,top,box,x,
1750 /* Compute the bonded and non-bonded energies and optionally forces */
1751 do_force_lowlevel(fplog,step,fr,inputrec,&(top->idef),
1752 cr,nrnb,wcycle,mdatoms,&(inputrec->opts),
1753 x,hist,f, bSepLRF ? fr->f_twin : f,enerd,fcd,mtop,top,fr->born,
1754 &(top->atomtypes),bBornRadii,box,
1755 inputrec->fepvals,lambda,
1756 graph,&(top->excls),fr->mu_tot,
1762 if (do_per_step(step,inputrec->nstcalclr))
1764 /* Add the long range forces to the short range forces */
1765 for(i=0; i<fr->natoms_force_constr; i++)
1767 rvec_add(fr->f_twin[i],f[i],f[i]);
1772 cycles_force = wallcycle_stop(wcycle,ewcFORCE);
1773 GMX_BARRIER(cr->mpi_comm_mygroup);
1777 do_flood(fplog,cr,x,f,ed,box,step,bNS);
1780 if (DOMAINDECOMP(cr))
1782 dd_force_flop_stop(cr->dd,nrnb);
1785 dd_cycles_add(cr->dd,cycles_force-cycles_pme,ddCyclF);
1791 if (IR_ELEC_FIELD(*inputrec))
1793 /* Compute forces due to electric field */
1794 calc_f_el(MASTER(cr) ? field : NULL,
1795 start,homenr,mdatoms->chargeA,x,fr->f_novirsum,
1796 inputrec->ex,inputrec->et,t);
1799 if (bDoAdressWF && fr->adress_icor == eAdressICThermoForce)
1801 /* Compute thermodynamic force in hybrid AdResS region */
1802 adress_thermo_force(start,homenr,&(top->cgs),x,fr->f_novirsum,fr,mdatoms,
1803 inputrec->ePBC==epbcNONE ? NULL : &pbc);
1806 /* Communicate the forces */
1809 wallcycle_start(wcycle,ewcMOVEF);
1810 if (DOMAINDECOMP(cr))
1812 dd_move_f(cr->dd,f,fr->fshift);
1813 /* Do we need to communicate the separate force array
1814 * for terms that do not contribute to the single sum virial?
1815 * Position restraints and electric fields do not introduce
1816 * inter-cg forces, only full electrostatics methods do.
1817 * When we do not calculate the virial, fr->f_novirsum = f,
1818 * so we have already communicated these forces.
1820 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
1821 (flags & GMX_FORCE_VIRIAL))
1823 dd_move_f(cr->dd,fr->f_novirsum,NULL);
1827 /* We should not update the shift forces here,
1828 * since f_twin is already included in f.
1830 dd_move_f(cr->dd,fr->f_twin,NULL);
1835 pd_move_f(cr,f,nrnb);
1838 pd_move_f(cr,fr->f_twin,nrnb);
1841 wallcycle_stop(wcycle,ewcMOVEF);
1844 /* If we have NoVirSum forces, but we do not calculate the virial,
1845 * we sum fr->f_novirum=f later.
1847 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
1849 wallcycle_start(wcycle,ewcVSITESPREAD);
1850 spread_vsite_f(fplog,vsite,x,f,fr->fshift,FALSE,NULL,nrnb,
1851 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
1852 wallcycle_stop(wcycle,ewcVSITESPREAD);
1856 wallcycle_start(wcycle,ewcVSITESPREAD);
1857 spread_vsite_f(fplog,vsite,x,fr->f_twin,NULL,FALSE,NULL,
1859 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
1860 wallcycle_stop(wcycle,ewcVSITESPREAD);
1864 if (flags & GMX_FORCE_VIRIAL)
1866 /* Calculation of the virial must be done after vsites! */
1867 calc_virial(fplog,mdatoms->start,mdatoms->homenr,x,f,
1868 vir_force,graph,box,nrnb,fr,inputrec->ePBC);
1872 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
1874 pull_potential_wrapper(fplog,bSepDVDL,cr,inputrec,box,x,
1875 f,vir_force,mdatoms,enerd,lambda,t);
1878 /* Add the forces from enforced rotation potentials (if any) */
1881 wallcycle_start(wcycle,ewcROTadd);
1882 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr,step,t);
1883 wallcycle_stop(wcycle,ewcROTadd);
1886 if (PAR(cr) && !(cr->duty & DUTY_PME))
1888 /* In case of node-splitting, the PP nodes receive the long-range
1889 * forces, virial and energy from the PME nodes here.
1891 pme_receive_force_ener(fplog,bSepDVDL,cr,wcycle,enerd,fr);
1896 post_process_forces(fplog,cr,step,nrnb,wcycle,
1897 top,box,x,f,vir_force,mdatoms,graph,fr,vsite,
1901 /* Sum the potential energy terms from group contributions */
1902 sum_epot(&(inputrec->opts),&(enerd->grpp),enerd->term);
1905 void do_force(FILE *fplog,t_commrec *cr,
1906 t_inputrec *inputrec,
1907 gmx_large_int_t step,t_nrnb *nrnb,gmx_wallcycle_t wcycle,
1908 gmx_localtop_t *top,
1910 gmx_groups_t *groups,
1911 matrix box,rvec x[],history_t *hist,
1915 gmx_enerdata_t *enerd,t_fcdata *fcd,
1916 real *lambda,t_graph *graph,
1918 gmx_vsite_t *vsite,rvec mu_tot,
1919 double t,FILE *field,gmx_edsam_t ed,
1920 gmx_bool bBornRadii,
1923 /* modify force flag if not doing nonbonded */
1924 if (!fr->bNonbonded)
1926 flags &= ~GMX_FORCE_NONBONDED;
1929 switch (inputrec->cutoff_scheme)
1932 do_force_cutsVERLET(fplog, cr, inputrec,
1948 do_force_cutsGROUP(fplog, cr, inputrec,
1963 gmx_incons("Invalid cut-off scheme passed!");
1968 void do_constrain_first(FILE *fplog,gmx_constr_t constr,
1969 t_inputrec *ir,t_mdatoms *md,
1970 t_state *state,rvec *f,
1971 t_graph *graph,t_commrec *cr,t_nrnb *nrnb,
1972 t_forcerec *fr, gmx_localtop_t *top, tensor shake_vir)
1975 gmx_large_int_t step;
1976 real dt=ir->delta_t;
1980 snew(savex,state->natoms);
1983 end = md->homenr + start;
1986 fprintf(debug,"vcm: start=%d, homenr=%d, end=%d\n",
1987 start,md->homenr,end);
1988 /* Do a first constrain to reset particles... */
1989 step = ir->init_step;
1992 char buf[STEPSTRSIZE];
1993 fprintf(fplog,"\nConstraining the starting coordinates (step %s)\n",
1994 gmx_step_str(step,buf));
1998 /* constrain the current position */
1999 constrain(NULL,TRUE,FALSE,constr,&(top->idef),
2000 ir,NULL,cr,step,0,md,
2001 state->x,state->x,NULL,
2002 fr->bMolPBC,state->box,
2003 state->lambda[efptBONDED],&dvdl_dum,
2004 NULL,NULL,nrnb,econqCoord,
2005 ir->epc==epcMTTK,state->veta,state->veta);
2008 /* constrain the inital velocity, and save it */
2009 /* also may be useful if we need the ekin from the halfstep for velocity verlet */
2010 /* might not yet treat veta correctly */
2011 constrain(NULL,TRUE,FALSE,constr,&(top->idef),
2012 ir,NULL,cr,step,0,md,
2013 state->x,state->v,state->v,
2014 fr->bMolPBC,state->box,
2015 state->lambda[efptBONDED],&dvdl_dum,
2016 NULL,NULL,nrnb,econqVeloc,
2017 ir->epc==epcMTTK,state->veta,state->veta);
2019 /* constrain the inital velocities at t-dt/2 */
2020 if (EI_STATE_VELOCITY(ir->eI) && ir->eI!=eiVV)
2022 for(i=start; (i<end); i++)
2024 for(m=0; (m<DIM); m++)
2026 /* Reverse the velocity */
2027 state->v[i][m] = -state->v[i][m];
2028 /* Store the position at t-dt in buf */
2029 savex[i][m] = state->x[i][m] + dt*state->v[i][m];
2032 /* Shake the positions at t=-dt with the positions at t=0
2033 * as reference coordinates.
2037 char buf[STEPSTRSIZE];
2038 fprintf(fplog,"\nConstraining the coordinates at t0-dt (step %s)\n",
2039 gmx_step_str(step,buf));
2042 constrain(NULL,TRUE,FALSE,constr,&(top->idef),
2043 ir,NULL,cr,step,-1,md,
2044 state->x,savex,NULL,
2045 fr->bMolPBC,state->box,
2046 state->lambda[efptBONDED],&dvdl_dum,
2047 state->v,NULL,nrnb,econqCoord,
2048 ir->epc==epcMTTK,state->veta,state->veta);
2050 for(i=start; i<end; i++) {
2051 for(m=0; m<DIM; m++) {
2052 /* Re-reverse the velocities */
2053 state->v[i][m] = -state->v[i][m];
2060 void calc_enervirdiff(FILE *fplog,int eDispCorr,t_forcerec *fr)
2062 double eners[2],virs[2],enersum,virsum,y0,f,g,h;
2063 double r0,r1,r,rc3,rc9,ea,eb,ec,pa,pb,pc,pd;
2064 double invscale,invscale2,invscale3;
2065 int ri0,ri1,ri,i,offstart,offset;
2066 real scale,*vdwtab,tabfactor,tmp;
2068 fr->enershiftsix = 0;
2069 fr->enershifttwelve = 0;
2070 fr->enerdiffsix = 0;
2071 fr->enerdifftwelve = 0;
2073 fr->virdifftwelve = 0;
2075 if (eDispCorr != edispcNO) {
2076 for(i=0; i<2; i++) {
2080 if ((fr->vdwtype == evdwSWITCH) || (fr->vdwtype == evdwSHIFT)) {
2081 if (fr->rvdw_switch == 0)
2083 "With dispersion correction rvdw-switch can not be zero "
2084 "for vdw-type = %s",evdw_names[fr->vdwtype]);
2086 scale = fr->nblists[0].table_elec_vdw.scale;
2087 vdwtab = fr->nblists[0].table_vdw.data;
2089 /* Round the cut-offs to exact table values for precision */
2090 ri0 = floor(fr->rvdw_switch*scale);
2091 ri1 = ceil(fr->rvdw*scale);
2097 if (fr->vdwtype == evdwSHIFT)
2099 /* Determine the constant energy shift below rvdw_switch.
2100 * Table has a scale factor since we have scaled it down to compensate
2101 * for scaling-up c6/c12 with the derivative factors to save flops in analytical kernels.
2103 fr->enershiftsix = (real)(-1.0/(rc3*rc3)) - 6.0*vdwtab[8*ri0];
2104 fr->enershifttwelve = (real)( 1.0/(rc9*rc3)) - 12.0*vdwtab[8*ri0 + 4];
2106 /* Add the constant part from 0 to rvdw_switch.
2107 * This integration from 0 to rvdw_switch overcounts the number
2108 * of interactions by 1, as it also counts the self interaction.
2109 * We will correct for this later.
2111 eners[0] += 4.0*M_PI*fr->enershiftsix*rc3/3.0;
2112 eners[1] += 4.0*M_PI*fr->enershifttwelve*rc3/3.0;
2114 invscale = 1.0/(scale);
2115 invscale2 = invscale*invscale;
2116 invscale3 = invscale*invscale2;
2118 /* following summation derived from cubic spline definition,
2119 Numerical Recipies in C, second edition, p. 113-116. Exact
2120 for the cubic spline. We first calculate the negative of
2121 the energy from rvdw to rvdw_switch, assuming that g(r)=1,
2122 and then add the more standard, abrupt cutoff correction to
2123 that result, yielding the long-range correction for a
2124 switched function. We perform both the pressure and energy
2125 loops at the same time for simplicity, as the computational
2129 enersum = 0.0; virsum = 0.0;
2133 /* Since the dispersion table has been scaled down a factor 6.0 and the repulsion
2134 * a factor 12.0 to compensate for the c6/c12 parameters inside nbfp[] being scaled
2135 * up (to save flops in kernels), we need to correct for this.
2144 for (ri=ri0; ri<ri1; ri++) {
2147 eb = 2.0*invscale2*r;
2151 pb = 3.0*invscale2*r;
2152 pc = 3.0*invscale*r*r;
2155 /* this "8" is from the packing in the vdwtab array - perhaps should be #define'ed? */
2156 offset = 8*ri + offstart;
2157 y0 = vdwtab[offset];
2158 f = vdwtab[offset+1];
2159 g = vdwtab[offset+2];
2160 h = vdwtab[offset+3];
2162 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);
2163 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);
2166 enersum *= 4.0*M_PI*tabfactor;
2167 virsum *= 4.0*M_PI*tabfactor;
2168 eners[i] -= enersum;
2172 /* now add the correction for rvdw_switch to infinity */
2173 eners[0] += -4.0*M_PI/(3.0*rc3);
2174 eners[1] += 4.0*M_PI/(9.0*rc9);
2175 virs[0] += 8.0*M_PI/rc3;
2176 virs[1] += -16.0*M_PI/(3.0*rc9);
2178 else if ((fr->vdwtype == evdwCUT) || (fr->vdwtype == evdwUSER)) {
2179 if (fr->vdwtype == evdwUSER && fplog)
2181 "WARNING: using dispersion correction with user tables\n");
2182 rc3 = fr->rvdw*fr->rvdw*fr->rvdw;
2184 /* Contribution beyond the cut-off */
2185 eners[0] += -4.0*M_PI/(3.0*rc3);
2186 eners[1] += 4.0*M_PI/(9.0*rc9);
2187 if (fr->vdw_modifier==eintmodPOTSHIFT) {
2188 /* Contribution within the cut-off */
2189 eners[0] += -4.0*M_PI/(3.0*rc3);
2190 eners[1] += 4.0*M_PI/(3.0*rc9);
2192 /* Contribution beyond the cut-off */
2193 virs[0] += 8.0*M_PI/rc3;
2194 virs[1] += -16.0*M_PI/(3.0*rc9);
2197 "Dispersion correction is not implemented for vdw-type = %s",
2198 evdw_names[fr->vdwtype]);
2200 fr->enerdiffsix = eners[0];
2201 fr->enerdifftwelve = eners[1];
2202 /* The 0.5 is due to the Gromacs definition of the virial */
2203 fr->virdiffsix = 0.5*virs[0];
2204 fr->virdifftwelve = 0.5*virs[1];
2208 void calc_dispcorr(FILE *fplog,t_inputrec *ir,t_forcerec *fr,
2209 gmx_large_int_t step,int natoms,
2210 matrix box,real lambda,tensor pres,tensor virial,
2211 real *prescorr, real *enercorr, real *dvdlcorr)
2213 gmx_bool bCorrAll,bCorrPres;
2214 real dvdlambda,invvol,dens,ninter,avcsix,avctwelve,enerdiff,svir=0,spres=0;
2224 if (ir->eDispCorr != edispcNO) {
2225 bCorrAll = (ir->eDispCorr == edispcAllEner ||
2226 ir->eDispCorr == edispcAllEnerPres);
2227 bCorrPres = (ir->eDispCorr == edispcEnerPres ||
2228 ir->eDispCorr == edispcAllEnerPres);
2230 invvol = 1/det(box);
2233 /* Only correct for the interactions with the inserted molecule */
2234 dens = (natoms - fr->n_tpi)*invvol;
2239 dens = natoms*invvol;
2240 ninter = 0.5*natoms;
2243 if (ir->efep == efepNO)
2245 avcsix = fr->avcsix[0];
2246 avctwelve = fr->avctwelve[0];
2250 avcsix = (1 - lambda)*fr->avcsix[0] + lambda*fr->avcsix[1];
2251 avctwelve = (1 - lambda)*fr->avctwelve[0] + lambda*fr->avctwelve[1];
2254 enerdiff = ninter*(dens*fr->enerdiffsix - fr->enershiftsix);
2255 *enercorr += avcsix*enerdiff;
2257 if (ir->efep != efepNO)
2259 dvdlambda += (fr->avcsix[1] - fr->avcsix[0])*enerdiff;
2263 enerdiff = ninter*(dens*fr->enerdifftwelve - fr->enershifttwelve);
2264 *enercorr += avctwelve*enerdiff;
2265 if (fr->efep != efepNO)
2267 dvdlambda += (fr->avctwelve[1] - fr->avctwelve[0])*enerdiff;
2273 svir = ninter*dens*avcsix*fr->virdiffsix/3.0;
2274 if (ir->eDispCorr == edispcAllEnerPres)
2276 svir += ninter*dens*avctwelve*fr->virdifftwelve/3.0;
2278 /* The factor 2 is because of the Gromacs virial definition */
2279 spres = -2.0*invvol*svir*PRESFAC;
2281 for(m=0; m<DIM; m++) {
2282 virial[m][m] += svir;
2283 pres[m][m] += spres;
2288 /* Can't currently control when it prints, for now, just print when degugging */
2292 fprintf(debug,"Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
2298 "Long Range LJ corr.: Epot %10g, Pres: %10g, Vir: %10g\n",
2299 *enercorr,spres,svir);
2303 fprintf(debug,"Long Range LJ corr.: Epot %10g\n",*enercorr);
2307 if (fr->bSepDVDL && do_per_step(step,ir->nstlog))
2309 fprintf(fplog,sepdvdlformat,"Dispersion correction",
2310 *enercorr,dvdlambda);
2312 if (fr->efep != efepNO)
2314 *dvdlcorr += dvdlambda;
2319 void do_pbc_first(FILE *fplog,matrix box,t_forcerec *fr,
2320 t_graph *graph,rvec x[])
2323 fprintf(fplog,"Removing pbc first time\n");
2324 calc_shifts(box,fr->shift_vec);
2326 mk_mshift(fplog,graph,fr->ePBC,box,x);
2328 p_graph(debug,"do_pbc_first 1",graph);
2329 shift_self(graph,box,x);
2330 /* By doing an extra mk_mshift the molecules that are broken
2331 * because they were e.g. imported from another software
2332 * will be made whole again. Such are the healing powers
2335 mk_mshift(fplog,graph,fr->ePBC,box,x);
2337 p_graph(debug,"do_pbc_first 2",graph);
2340 fprintf(fplog,"Done rmpbc\n");
2343 static void low_do_pbc_mtop(FILE *fplog,int ePBC,matrix box,
2344 gmx_mtop_t *mtop,rvec x[],
2349 gmx_molblock_t *molb;
2351 if (bFirst && fplog)
2352 fprintf(fplog,"Removing pbc first time\n");
2356 for(mb=0; mb<mtop->nmolblock; mb++) {
2357 molb = &mtop->molblock[mb];
2358 if (molb->natoms_mol == 1 ||
2359 (!bFirst && mtop->moltype[molb->type].cgs.nr == 1)) {
2360 /* Just one atom or charge group in the molecule, no PBC required */
2361 as += molb->nmol*molb->natoms_mol;
2363 /* Pass NULL iso fplog to avoid graph prints for each molecule type */
2364 mk_graph_ilist(NULL,mtop->moltype[molb->type].ilist,
2365 0,molb->natoms_mol,FALSE,FALSE,graph);
2367 for(mol=0; mol<molb->nmol; mol++) {
2368 mk_mshift(fplog,graph,ePBC,box,x+as);
2370 shift_self(graph,box,x+as);
2371 /* The molecule is whole now.
2372 * We don't need the second mk_mshift call as in do_pbc_first,
2373 * since we no longer need this graph.
2376 as += molb->natoms_mol;
2384 void do_pbc_first_mtop(FILE *fplog,int ePBC,matrix box,
2385 gmx_mtop_t *mtop,rvec x[])
2387 low_do_pbc_mtop(fplog,ePBC,box,mtop,x,TRUE);
2390 void do_pbc_mtop(FILE *fplog,int ePBC,matrix box,
2391 gmx_mtop_t *mtop,rvec x[])
2393 low_do_pbc_mtop(fplog,ePBC,box,mtop,x,FALSE);
2396 void finish_run(FILE *fplog,t_commrec *cr,const char *confout,
2397 t_inputrec *inputrec,
2398 t_nrnb nrnb[],gmx_wallcycle_t wcycle,
2399 gmx_runtime_t *runtime,
2400 wallclock_gpu_t *gputimes,
2402 gmx_bool bWriteStat)
2405 t_nrnb *nrnb_tot=NULL;
2409 wallcycle_sum(cr,wcycle);
2415 MPI_Allreduce(nrnb->n,nrnb_tot->n,eNRNB,MPI_DOUBLE,MPI_SUM,
2416 cr->mpi_comm_mysim);
2424 #if defined(GMX_MPI) && !defined(GMX_THREAD_MPI)
2427 /* reduce nodetime over all MPI processes in the current simulation */
2429 MPI_Allreduce(&runtime->proctime,&sum,1,MPI_DOUBLE,MPI_SUM,
2430 cr->mpi_comm_mysim);
2431 runtime->proctime = sum;
2437 print_flop(fplog,nrnb_tot,&nbfs,&mflop);
2444 if ((cr->duty & DUTY_PP) && DOMAINDECOMP(cr))
2446 print_dd_statistics(cr,inputrec,fplog);
2458 snew(nrnb_all,cr->nnodes);
2459 nrnb_all[0] = *nrnb;
2460 for(s=1; s<cr->nnodes; s++)
2462 MPI_Recv(nrnb_all[s].n,eNRNB,MPI_DOUBLE,s,0,
2463 cr->mpi_comm_mysim,&stat);
2465 pr_load(fplog,cr,nrnb_all);
2470 MPI_Send(nrnb->n,eNRNB,MPI_DOUBLE,MASTERRANK(cr),0,
2471 cr->mpi_comm_mysim);
2478 wallcycle_print(fplog,cr->nnodes,cr->npmenodes,runtime->realtime,
2481 if (EI_DYNAMICS(inputrec->eI))
2483 delta_t = inputrec->delta_t;
2492 print_perf(fplog,runtime->proctime,runtime->realtime,
2493 cr->nnodes-cr->npmenodes,
2494 runtime->nsteps_done,delta_t,nbfs,mflop,
2499 print_perf(stderr,runtime->proctime,runtime->realtime,
2500 cr->nnodes-cr->npmenodes,
2501 runtime->nsteps_done,delta_t,nbfs,mflop,
2507 extern void initialize_lambdas(FILE *fplog,t_inputrec *ir,int *fep_state,real *lambda,double *lam0)
2509 /* this function works, but could probably use a logic rewrite to keep all the different
2510 types of efep straight. */
2513 t_lambda *fep = ir->fepvals;
2515 if ((ir->efep==efepNO) && (ir->bSimTemp == FALSE)) {
2516 for (i=0;i<efptNR;i++) {
2525 *fep_state = fep->init_fep_state; /* this might overwrite the checkpoint
2526 if checkpoint is set -- a kludge is in for now
2528 for (i=0;i<efptNR;i++)
2530 /* overwrite lambda state with init_lambda for now for backwards compatibility */
2531 if (fep->init_lambda>=0) /* if it's -1, it was never initializd */
2533 lambda[i] = fep->init_lambda;
2535 lam0[i] = lambda[i];
2540 lambda[i] = fep->all_lambda[i][*fep_state];
2542 lam0[i] = lambda[i];
2547 /* need to rescale control temperatures to match current state */
2548 for (i=0;i<ir->opts.ngtc;i++) {
2549 if (ir->opts.ref_t[i] > 0) {
2550 ir->opts.ref_t[i] = ir->simtempvals->temperatures[*fep_state];
2556 /* Send to the log the information on the current lambdas */
2559 fprintf(fplog,"Initial vector of lambda components:[ ");
2560 for (i=0;i<efptNR;i++)
2562 fprintf(fplog,"%10.4f ",lambda[i]);
2564 fprintf(fplog,"]\n");
2570 void init_md(FILE *fplog,
2571 t_commrec *cr,t_inputrec *ir,const output_env_t oenv,
2572 double *t,double *t0,
2573 real *lambda, int *fep_state, double *lam0,
2574 t_nrnb *nrnb,gmx_mtop_t *mtop,
2576 int nfile,const t_filenm fnm[],
2577 gmx_mdoutf_t **outf,t_mdebin **mdebin,
2578 tensor force_vir,tensor shake_vir,rvec mu_tot,
2579 gmx_bool *bSimAnn,t_vcm **vcm, t_state *state, unsigned long Flags)
2584 /* Initial values */
2585 *t = *t0 = ir->init_t;
2588 for(i=0;i<ir->opts.ngtc;i++)
2590 /* set bSimAnn if any group is being annealed */
2591 if(ir->opts.annealing[i]!=eannNO)
2598 update_annealing_target_temp(&(ir->opts),ir->init_t);
2601 /* Initialize lambda variables */
2602 initialize_lambdas(fplog,ir,fep_state,lambda,lam0);
2606 *upd = init_update(fplog,ir);
2612 *vcm = init_vcm(fplog,&mtop->groups,ir);
2615 if (EI_DYNAMICS(ir->eI) && !(Flags & MD_APPENDFILES))
2617 if (ir->etc == etcBERENDSEN)
2619 please_cite(fplog,"Berendsen84a");
2621 if (ir->etc == etcVRESCALE)
2623 please_cite(fplog,"Bussi2007a");
2631 *outf = init_mdoutf(nfile,fnm,Flags,cr,ir,oenv);
2633 *mdebin = init_mdebin((Flags & MD_APPENDFILES) ? NULL : (*outf)->fp_ene,
2634 mtop,ir, (*outf)->fp_dhdl);
2639 please_cite(fplog,"Fritsch12");
2640 please_cite(fplog,"Junghans10");
2642 /* Initiate variables */
2643 clear_mat(force_vir);
2644 clear_mat(shake_vir);