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"
86 #include "gmx_wallcycle.h"
99 typedef struct gmx_timeprint {
104 /* Portable version of ctime_r implemented in src/gmxlib/string2.c, but we do not want it declared in public installed headers */
106 gmx_ctime_r(const time_t *clock,char *buf, int n);
112 #ifdef HAVE_GETTIMEOFDAY
116 gettimeofday(&t,NULL);
118 seconds = (double) t.tv_sec + 1e-6*(double)t.tv_usec;
124 seconds = time(NULL);
131 #define difftime(end,start) ((double)(end)-(double)(start))
133 void print_time(FILE *out,gmx_runtime_t *runtime,gmx_large_int_t step,
134 t_inputrec *ir, t_commrec *cr)
137 char timebuf[STRLEN];
147 fprintf(out,"step %s",gmx_step_str(step,buf));
148 if ((step >= ir->nstlist))
150 if ((ir->nstlist == 0) || ((step % ir->nstlist) == 0))
152 /* We have done a full cycle let's update time_per_step */
153 runtime->last = gmx_gettime();
154 dt = difftime(runtime->last,runtime->real);
155 runtime->time_per_step = dt/(step - ir->init_step + 1);
157 dt = (ir->nsteps + ir->init_step - step)*runtime->time_per_step;
163 finish = (time_t) (runtime->last + dt);
164 gmx_ctime_r(&finish,timebuf,STRLEN);
165 sprintf(buf,"%s",timebuf);
166 buf[strlen(buf)-1]='\0';
167 fprintf(out,", will finish %s",buf);
170 fprintf(out,", remaining runtime: %5d s ",(int)dt);
174 fprintf(out," performance: %.1f ns/day ",
175 ir->delta_t/1000*24*60*60/runtime->time_per_step);
192 static double set_proctime(gmx_runtime_t *runtime)
198 prev = runtime->proc;
199 runtime->proc = dclock();
201 diff = runtime->proc - prev;
205 prev = runtime->proc;
206 runtime->proc = clock();
208 diff = (double)(runtime->proc - prev)/(double)CLOCKS_PER_SEC;
212 /* The counter has probably looped, ignore this data */
219 void runtime_start(gmx_runtime_t *runtime)
221 runtime->real = gmx_gettime();
223 set_proctime(runtime);
224 runtime->realtime = 0;
225 runtime->proctime = 0;
227 runtime->time_per_step = 0;
230 void runtime_end(gmx_runtime_t *runtime)
236 runtime->proctime += set_proctime(runtime);
237 runtime->realtime = now - runtime->real;
241 void runtime_upd_proc(gmx_runtime_t *runtime)
243 runtime->proctime += set_proctime(runtime);
246 void print_date_and_time(FILE *fplog,int nodeid,const char *title,
247 const gmx_runtime_t *runtime)
250 char timebuf[STRLEN];
251 char time_string[STRLEN];
258 tmptime = (time_t) runtime->real;
259 gmx_ctime_r(&tmptime,timebuf,STRLEN);
263 tmptime = (time_t) gmx_gettime();
264 gmx_ctime_r(&tmptime,timebuf,STRLEN);
266 for(i=0; timebuf[i]>=' '; i++)
268 time_string[i]=timebuf[i];
272 fprintf(fplog,"%s on node %d %s\n",title,nodeid,time_string);
276 static void sum_forces(int start,int end,rvec f[],rvec flr[])
281 pr_rvecs(debug,0,"fsr",f+start,end-start);
282 pr_rvecs(debug,0,"flr",flr+start,end-start);
284 for(i=start; (i<end); i++)
285 rvec_inc(f[i],flr[i]);
289 * calc_f_el calculates forces due to an electric field.
291 * force is kJ mol^-1 nm^-1 = e * kJ mol^-1 nm^-1 / e
293 * Et[] contains the parameters for the time dependent
294 * part of the field (not yet used).
295 * Ex[] contains the parameters for
296 * the spatial dependent part of the field. You can have cool periodic
297 * fields in principle, but only a constant field is supported
299 * The function should return the energy due to the electric field
300 * (if any) but for now returns 0.
303 * There can be problems with the virial.
304 * Since the field is not self-consistent this is unavoidable.
305 * For neutral molecules the virial is correct within this approximation.
306 * For neutral systems with many charged molecules the error is small.
307 * But for systems with a net charge or a few charged molecules
308 * the error can be significant when the field is high.
309 * Solution: implement a self-consitent electric field into PME.
311 static void calc_f_el(FILE *fp,int start,int homenr,
312 real charge[],rvec x[],rvec f[],
313 t_cosines Ex[],t_cosines Et[],double t)
319 for(m=0; (m<DIM); m++)
326 Ext[m] = cos(Et[m].a[0]*(t-t0))*exp(-sqr(t-t0)/(2.0*sqr(Et[m].a[2])));
330 Ext[m] = cos(Et[m].a[0]*t);
339 /* Convert the field strength from V/nm to MD-units */
340 Ext[m] *= Ex[m].a[0]*FIELDFAC;
341 for(i=start; (i<start+homenr); i++)
342 f[i][m] += charge[i]*Ext[m];
351 fprintf(fp,"%10g %10g %10g %10g #FIELD\n",t,
352 Ext[XX]/FIELDFAC,Ext[YY]/FIELDFAC,Ext[ZZ]/FIELDFAC);
356 static void calc_virial(FILE *fplog,int start,int homenr,rvec x[],rvec f[],
357 tensor vir_part,t_graph *graph,matrix box,
358 t_nrnb *nrnb,const t_forcerec *fr,int ePBC)
363 /* The short-range virial from surrounding boxes */
365 calc_vir(fplog,SHIFTS,fr->shift_vec,fr->fshift,vir_part,ePBC==epbcSCREW,box);
366 inc_nrnb(nrnb,eNR_VIRIAL,SHIFTS);
368 /* Calculate partial virial, for local atoms only, based on short range.
369 * Total virial is computed in global_stat, called from do_md
371 f_calc_vir(fplog,start,start+homenr,x,f,vir_part,graph,box);
372 inc_nrnb(nrnb,eNR_VIRIAL,homenr);
374 /* Add position restraint contribution */
375 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++) {
381 vir_part[i][ZZ] += fr->vir_wall_z[i];
385 pr_rvecs(debug,0,"vir_part",vir_part,DIM);
388 static void print_large_forces(FILE *fp,t_mdatoms *md,t_commrec *cr,
389 gmx_large_int_t step,real pforce,rvec *x,rvec *f)
393 char buf[STEPSTRSIZE];
396 for(i=md->start; i<md->start+md->homenr; i++) {
398 /* We also catch NAN, if the compiler does not optimize this away. */
399 if (fn2 >= pf2 || fn2 != fn2) {
400 fprintf(fp,"step %s atom %6d x %8.3f %8.3f %8.3f force %12.5e\n",
401 gmx_step_str(step,buf),
402 ddglatnr(cr->dd,i),x[i][XX],x[i][YY],x[i][ZZ],sqrt(fn2));
407 void do_force(FILE *fplog,t_commrec *cr,
408 t_inputrec *inputrec,
409 gmx_large_int_t step,t_nrnb *nrnb,gmx_wallcycle_t wcycle,
412 gmx_groups_t *groups,
413 matrix box,rvec x[],history_t *hist,
417 gmx_enerdata_t *enerd,t_fcdata *fcd,
418 real lambda,t_graph *graph,
419 t_forcerec *fr,gmx_vsite_t *vsite,rvec mu_tot,
420 double t,FILE *field,gmx_edsam_t ed,
427 gmx_bool bSepDVDL,bStateChanged,bNS,bFillGrid,bCalcCGCM,bBS;
428 gmx_bool bDoLongRange,bDoForces,bSepLRF;
432 float cycles_ppdpme,cycles_pme,cycles_seppme,cycles_force;
434 start = mdatoms->start;
435 homenr = mdatoms->homenr;
437 bSepDVDL = (fr->bSepDVDL && do_per_step(step,inputrec->nstlog));
439 clear_mat(vir_force);
443 pd_cg_range(cr,&cg0,&cg1);
448 if (DOMAINDECOMP(cr))
450 cg1 = cr->dd->ncg_tot;
462 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
463 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll==FALSE);
464 bFillGrid = (bNS && bStateChanged);
465 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
466 bDoLongRange = (fr->bTwinRange && bNS && (flags & GMX_FORCE_DOLR));
467 bDoForces = (flags & GMX_FORCE_FORCES);
468 bSepLRF = (bDoLongRange && bDoForces && (flags & GMX_FORCE_SEPLRF));
472 update_forcerec(fplog,fr,box);
474 /* Calculate total (local) dipole moment in a temporary common array.
475 * This makes it possible to sum them over nodes faster.
477 calc_mu(start,homenr,
478 x,mdatoms->chargeA,mdatoms->chargeB,mdatoms->nChargePerturbed,
482 if (fr->ePBC != epbcNONE) {
483 /* Compute shift vectors every step,
484 * because of pressure coupling or box deformation!
486 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
487 calc_shifts(box,fr->shift_vec);
490 put_charge_groups_in_box(fplog,cg0,cg1,fr->ePBC,box,
491 &(top->cgs),x,fr->cg_cm);
492 inc_nrnb(nrnb,eNR_CGCM,homenr);
493 inc_nrnb(nrnb,eNR_RESETX,cg1-cg0);
495 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph) {
496 unshift_self(graph,box,x);
499 else if (bCalcCGCM) {
500 calc_cgcm(fplog,cg0,cg1,&(top->cgs),x,fr->cg_cm);
501 inc_nrnb(nrnb,eNR_CGCM,homenr);
506 move_cgcm(fplog,cr,fr->cg_cm);
509 pr_rvecs(debug,0,"cgcm",fr->cg_cm,top->cgs.nr);
513 if (!(cr->duty & DUTY_PME)) {
514 /* Send particle coordinates to the pme nodes.
515 * Since this is only implemented for domain decomposition
516 * and domain decomposition does not use the graph,
517 * we do not need to worry about shifting.
520 wallcycle_start(wcycle,ewcPP_PMESENDX);
522 bBS = (inputrec->nwall == 2);
525 svmul(inputrec->wall_ewald_zfac,boxs[ZZ],boxs[ZZ]);
528 gmx_pme_send_x(cr,bBS ? boxs : box,x,
529 mdatoms->nChargePerturbed,lambda,
530 ( flags & GMX_FORCE_VIRIAL),step);
532 wallcycle_stop(wcycle,ewcPP_PMESENDX);
536 /* Communicate coordinates and sum dipole if necessary */
539 wallcycle_start(wcycle,ewcMOVEX);
540 if (DOMAINDECOMP(cr))
542 dd_move_x(cr->dd,box,x);
546 move_x(fplog,cr,GMX_LEFT,GMX_RIGHT,x,nrnb);
548 /* When we don't need the total dipole we sum it in global_stat */
549 if (bStateChanged && NEED_MUTOT(*inputrec))
551 gmx_sumd(2*DIM,mu,cr);
553 wallcycle_stop(wcycle,ewcMOVEX);
561 fr->mu_tot[i][j] = mu[i*DIM + j];
565 if (fr->efep == efepNO)
567 copy_rvec(fr->mu_tot[0],mu_tot);
574 (1.0 - lambda)*fr->mu_tot[0][j] + lambda*fr->mu_tot[1][j];
579 reset_enerdata(&(inputrec->opts),fr,bNS,enerd,MASTER(cr));
580 clear_rvecs(SHIFTS,fr->fshift);
584 wallcycle_start(wcycle,ewcNS);
586 if (graph && bStateChanged)
588 /* Calculate intramolecular shift vectors to make molecules whole */
589 mk_mshift(fplog,graph,fr->ePBC,box,x);
592 /* Reset long range forces if necessary */
595 /* Reset the (long-range) forces if necessary */
596 clear_rvecs(fr->natoms_force_constr,bSepLRF ? fr->f_twin : f);
599 /* Do the actual neighbour searching and if twin range electrostatics
600 * also do the calculation of long range forces and energies.
604 groups,&(inputrec->opts),top,mdatoms,
605 cr,nrnb,lambda,&dvdl,&enerd->grpp,bFillGrid,
606 bDoLongRange,bDoForces,bSepLRF ? fr->f_twin : f);
609 fprintf(fplog,sepdvdlformat,"LR non-bonded",0.0,dvdl);
611 enerd->dvdl_lin += dvdl;
613 wallcycle_stop(wcycle,ewcNS);
616 if (inputrec->implicit_solvent && bNS)
618 make_gb_nblist(cr,inputrec->gb_algorithm,inputrec->rlist,
619 x,box,fr,&top->idef,graph,fr->born);
622 if (DOMAINDECOMP(cr))
624 if (!(cr->duty & DUTY_PME))
626 wallcycle_start(wcycle,ewcPPDURINGPME);
627 dd_force_flop_start(cr->dd,nrnb);
633 /* Enforced rotation has its own cycle counter that starts after the collective
634 * coordinates have been communicated. It is added to ddCyclF */
635 do_rotation(cr,inputrec,box,x,t,step,wcycle,bNS);
638 /* Start the force cycle counter.
639 * This counter is stopped in do_forcelow_level.
640 * No parallel communication should occur while this counter is running,
641 * since that will interfere with the dynamic load balancing.
643 wallcycle_start(wcycle,ewcFORCE);
647 /* Reset forces for which the virial is calculated separately:
648 * PME/Ewald forces if necessary */
651 if (flags & GMX_FORCE_VIRIAL)
653 fr->f_novirsum = fr->f_novirsum_alloc;
656 clear_rvecs(fr->f_novirsum_n,fr->f_novirsum);
660 clear_rvecs(homenr,fr->f_novirsum+start);
665 /* We are not calculating the pressure so we do not need
666 * a separate array for forces that do not contribute
675 /* Add the long range forces to the short range forces */
676 for(i=0; i<fr->natoms_force_constr; i++)
678 copy_rvec(fr->f_twin[i],f[i]);
681 else if (!(fr->bTwinRange && bNS))
683 /* Clear the short-range forces */
684 clear_rvecs(fr->natoms_force_constr,f);
687 clear_rvec(fr->vir_diag_posres);
689 if (inputrec->ePull == epullCONSTRAINT)
691 clear_pull_forces(inputrec->pull);
694 /* update QMMMrec, if necessary */
697 update_QMMMrec(cr,fr,x,mdatoms,box,top);
700 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
702 /* Position restraints always require full pbc */
703 set_pbc(&pbc,inputrec->ePBC,box);
704 v = posres(top->idef.il[F_POSRES].nr,top->idef.il[F_POSRES].iatoms,
705 top->idef.iparams_posres,
706 (const rvec*)x,fr->f_novirsum,fr->vir_diag_posres,
707 inputrec->ePBC==epbcNONE ? NULL : &pbc,lambda,&dvdl,
708 fr->rc_scaling,fr->ePBC,fr->posres_com,fr->posres_comB);
711 fprintf(fplog,sepdvdlformat,
712 interaction_function[F_POSRES].longname,v,dvdl);
714 enerd->term[F_POSRES] += v;
715 /* This linear lambda dependence assumption is only correct
716 * when only k depends on lambda,
717 * not when the reference position depends on lambda.
718 * grompp checks for this.
720 enerd->dvdl_lin += dvdl;
721 inc_nrnb(nrnb,eNR_POSRES,top->idef.il[F_POSRES].nr/2);
724 /* Compute the bonded and non-bonded energies and optionally forces */
725 do_force_lowlevel(fplog,step,fr,inputrec,&(top->idef),
726 cr,nrnb,wcycle,mdatoms,&(inputrec->opts),
727 x,hist,f,enerd,fcd,mtop,top,fr->born,
728 &(top->atomtypes),bBornRadii,box,
729 lambda,graph,&(top->excls),fr->mu_tot,
732 cycles_force = wallcycle_stop(wcycle,ewcFORCE);
736 do_flood(fplog,cr,x,f,ed,box,step);
739 if (DOMAINDECOMP(cr))
741 dd_force_flop_stop(cr->dd,nrnb);
744 dd_cycles_add(cr->dd,cycles_force-cycles_pme,ddCyclF);
750 if (IR_ELEC_FIELD(*inputrec))
752 /* Compute forces due to electric field */
753 calc_f_el(MASTER(cr) ? field : NULL,
754 start,homenr,mdatoms->chargeA,x,fr->f_novirsum,
755 inputrec->ex,inputrec->et,t);
758 /* Communicate the forces */
761 wallcycle_start(wcycle,ewcMOVEF);
762 if (DOMAINDECOMP(cr))
764 dd_move_f(cr->dd,f,fr->fshift);
765 /* Do we need to communicate the separate force array
766 * for terms that do not contribute to the single sum virial?
767 * Position restraints and electric fields do not introduce
768 * inter-cg forces, only full electrostatics methods do.
769 * When we do not calculate the virial, fr->f_novirsum = f,
770 * so we have already communicated these forces.
772 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
773 (flags & GMX_FORCE_VIRIAL))
775 dd_move_f(cr->dd,fr->f_novirsum,NULL);
779 /* We should not update the shift forces here,
780 * since f_twin is already included in f.
782 dd_move_f(cr->dd,fr->f_twin,NULL);
787 pd_move_f(cr,f,nrnb);
790 pd_move_f(cr,fr->f_twin,nrnb);
793 wallcycle_stop(wcycle,ewcMOVEF);
796 /* If we have NoVirSum forces, but we do not calculate the virial,
797 * we sum fr->f_novirum=f later.
799 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
801 wallcycle_start(wcycle,ewcVSITESPREAD);
802 spread_vsite_f(fplog,vsite,x,f,fr->fshift,nrnb,
803 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
804 wallcycle_stop(wcycle,ewcVSITESPREAD);
808 wallcycle_start(wcycle,ewcVSITESPREAD);
809 spread_vsite_f(fplog,vsite,x,fr->f_twin,NULL,
811 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
812 wallcycle_stop(wcycle,ewcVSITESPREAD);
816 if (flags & GMX_FORCE_VIRIAL)
818 /* Calculation of the virial must be done after vsites! */
819 calc_virial(fplog,mdatoms->start,mdatoms->homenr,x,f,
820 vir_force,graph,box,nrnb,fr,inputrec->ePBC);
824 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
826 /* Calculate the center of mass forces, this requires communication,
827 * which is why pull_potential is called close to other communication.
828 * The virial contribution is calculated directly,
829 * which is why we call pull_potential after calc_virial.
831 set_pbc(&pbc,inputrec->ePBC,box);
833 enerd->term[F_COM_PULL] =
834 pull_potential(inputrec->ePull,inputrec->pull,mdatoms,&pbc,
835 cr,t,lambda,x,f,vir_force,&dvdl);
838 fprintf(fplog,sepdvdlformat,"Com pull",enerd->term[F_COM_PULL],dvdl);
840 enerd->dvdl_lin += dvdl;
843 enerd->term[F_COM_PULL] = 0.0;
845 /* Add the forces from enforced rotation potentials (if any) */
847 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr, step, t);
850 if (PAR(cr) && !(cr->duty & DUTY_PME))
852 cycles_ppdpme = wallcycle_stop(wcycle,ewcPPDURINGPME);
853 dd_cycles_add(cr->dd,cycles_ppdpme,ddCyclPPduringPME);
855 /* In case of node-splitting, the PP nodes receive the long-range
856 * forces, virial and energy from the PME nodes here.
858 wallcycle_start(wcycle,ewcPP_PMEWAITRECVF);
860 gmx_pme_receive_f(cr,fr->f_novirsum,fr->vir_el_recip,&e,&dvdl,
864 fprintf(fplog,sepdvdlformat,"PME mesh",e,dvdl);
866 enerd->term[F_COUL_RECIP] += e;
867 enerd->dvdl_lin += dvdl;
870 dd_cycles_add(cr->dd,cycles_seppme,ddCyclPME);
872 wallcycle_stop(wcycle,ewcPP_PMEWAITRECVF);
875 if (bDoForces && fr->bF_NoVirSum)
879 /* Spread the mesh force on virtual sites to the other particles...
880 * This is parallellized. MPI communication is performed
881 * if the constructing atoms aren't local.
883 wallcycle_start(wcycle,ewcVSITESPREAD);
884 spread_vsite_f(fplog,vsite,x,fr->f_novirsum,NULL,nrnb,
885 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
886 wallcycle_stop(wcycle,ewcVSITESPREAD);
888 if (flags & GMX_FORCE_VIRIAL)
890 /* Now add the forces, this is local */
893 sum_forces(0,fr->f_novirsum_n,f,fr->f_novirsum);
897 sum_forces(start,start+homenr,f,fr->f_novirsum);
899 if (EEL_FULL(fr->eeltype))
901 /* Add the mesh contribution to the virial */
902 m_add(vir_force,fr->vir_el_recip,vir_force);
906 pr_rvecs(debug,0,"vir_force",vir_force,DIM);
911 /* Sum the potential energy terms from group contributions */
912 sum_epot(&(inputrec->opts),enerd);
914 if (fr->print_force >= 0 && bDoForces)
916 print_large_forces(stderr,mdatoms,cr,step,fr->print_force,x,f);
920 void do_constrain_first(FILE *fplog,gmx_constr_t constr,
921 t_inputrec *ir,t_mdatoms *md,
922 t_state *state,rvec *f,
923 t_graph *graph,t_commrec *cr,t_nrnb *nrnb,
924 t_forcerec *fr, gmx_localtop_t *top, tensor shake_vir)
927 gmx_large_int_t step;
928 double mass,tmass,vcm[4];
933 snew(savex,state->natoms);
936 end = md->homenr + start;
939 fprintf(debug,"vcm: start=%d, homenr=%d, end=%d\n",
940 start,md->homenr,end);
941 /* Do a first constrain to reset particles... */
942 step = ir->init_step;
945 char buf[STEPSTRSIZE];
946 fprintf(fplog,"\nConstraining the starting coordinates (step %s)\n",
947 gmx_step_str(step,buf));
951 /* constrain the current position */
952 constrain(NULL,TRUE,FALSE,constr,&(top->idef),
953 ir,NULL,cr,step,0,md,
954 state->x,state->x,NULL,
955 state->box,state->lambda,&dvdlambda,
956 NULL,NULL,nrnb,econqCoord,ir->epc==epcMTTK,state->veta,state->veta);
959 /* constrain the inital velocity, and save it */
960 /* also may be useful if we need the ekin from the halfstep for velocity verlet */
961 /* might not yet treat veta correctly */
962 constrain(NULL,TRUE,FALSE,constr,&(top->idef),
963 ir,NULL,cr,step,0,md,
964 state->x,state->v,state->v,
965 state->box,state->lambda,&dvdlambda,
966 NULL,NULL,nrnb,econqVeloc,ir->epc==epcMTTK,state->veta,state->veta);
968 /* constrain the inital velocities at t-dt/2 */
969 if (EI_STATE_VELOCITY(ir->eI) && ir->eI!=eiVV)
971 for(i=start; (i<end); i++)
973 for(m=0; (m<DIM); m++)
975 /* Reverse the velocity */
976 state->v[i][m] = -state->v[i][m];
977 /* Store the position at t-dt in buf */
978 savex[i][m] = state->x[i][m] + dt*state->v[i][m];
981 /* Shake the positions at t=-dt with the positions at t=0
982 * as reference coordinates.
986 char buf[STEPSTRSIZE];
987 fprintf(fplog,"\nConstraining the coordinates at t0-dt (step %s)\n",
988 gmx_step_str(step,buf));
991 constrain(NULL,TRUE,FALSE,constr,&(top->idef),
992 ir,NULL,cr,step,-1,md,
994 state->box,state->lambda,&dvdlambda,
995 state->v,NULL,nrnb,econqCoord,ir->epc==epcMTTK,state->veta,state->veta);
997 for(i=start; i<end; i++) {
998 for(m=0; m<DIM; m++) {
999 /* Re-reverse the velocities */
1000 state->v[i][m] = -state->v[i][m];
1005 for(m=0; (m<4); m++)
1007 for(i=start; i<end; i++) {
1008 mass = md->massT[i];
1009 for(m=0; m<DIM; m++) {
1010 vcm[m] += state->v[i][m]*mass;
1015 if (ir->nstcomm != 0 || debug) {
1016 /* Compute the global sum of vcm */
1018 fprintf(debug,"vcm: %8.3f %8.3f %8.3f,"
1019 " total mass = %12.5e\n",vcm[XX],vcm[YY],vcm[ZZ],vcm[3]);
1023 for(m=0; (m<DIM); m++)
1026 fprintf(debug,"vcm: %8.3f %8.3f %8.3f,"
1027 " total mass = %12.5e\n",vcm[XX],vcm[YY],vcm[ZZ],tmass);
1028 if (ir->nstcomm != 0) {
1029 /* Now we have the velocity of center of mass, let's remove it */
1030 for(i=start; (i<end); i++) {
1031 for(m=0; (m<DIM); m++)
1032 state->v[i][m] -= vcm[m];
1040 void calc_enervirdiff(FILE *fplog,int eDispCorr,t_forcerec *fr)
1042 double eners[2],virs[2],enersum,virsum,y0,f,g,h;
1043 double r0,r1,r,rc3,rc9,ea,eb,ec,pa,pb,pc,pd;
1044 double invscale,invscale2,invscale3;
1045 int ri0,ri1,ri,i,offstart,offset;
1048 fr->enershiftsix = 0;
1049 fr->enershifttwelve = 0;
1050 fr->enerdiffsix = 0;
1051 fr->enerdifftwelve = 0;
1053 fr->virdifftwelve = 0;
1055 if (eDispCorr != edispcNO) {
1056 for(i=0; i<2; i++) {
1060 if ((fr->vdwtype == evdwSWITCH) || (fr->vdwtype == evdwSHIFT)) {
1061 if (fr->rvdw_switch == 0)
1063 "With dispersion correction rvdw-switch can not be zero "
1064 "for vdw-type = %s",evdw_names[fr->vdwtype]);
1066 scale = fr->nblists[0].tab.scale;
1067 vdwtab = fr->nblists[0].vdwtab;
1069 /* Round the cut-offs to exact table values for precision */
1070 ri0 = floor(fr->rvdw_switch*scale);
1071 ri1 = ceil(fr->rvdw*scale);
1077 if (fr->vdwtype == evdwSHIFT) {
1078 /* Determine the constant energy shift below rvdw_switch */
1079 fr->enershiftsix = (real)(-1.0/(rc3*rc3)) - vdwtab[8*ri0];
1080 fr->enershifttwelve = (real)( 1.0/(rc9*rc3)) - vdwtab[8*ri0 + 4];
1082 /* Add the constant part from 0 to rvdw_switch.
1083 * This integration from 0 to rvdw_switch overcounts the number
1084 * of interactions by 1, as it also counts the self interaction.
1085 * We will correct for this later.
1087 eners[0] += 4.0*M_PI*fr->enershiftsix*rc3/3.0;
1088 eners[1] += 4.0*M_PI*fr->enershifttwelve*rc3/3.0;
1090 invscale = 1.0/(scale);
1091 invscale2 = invscale*invscale;
1092 invscale3 = invscale*invscale2;
1094 /* following summation derived from cubic spline definition,
1095 Numerical Recipies in C, second edition, p. 113-116. Exact
1096 for the cubic spline. We first calculate the negative of
1097 the energy from rvdw to rvdw_switch, assuming that g(r)=1,
1098 and then add the more standard, abrupt cutoff correction to
1099 that result, yielding the long-range correction for a
1100 switched function. We perform both the pressure and energy
1101 loops at the same time for simplicity, as the computational
1105 enersum = 0.0; virsum = 0.0;
1110 for (ri=ri0; ri<ri1; ri++) {
1113 eb = 2.0*invscale2*r;
1117 pb = 3.0*invscale2*r;
1118 pc = 3.0*invscale*r*r;
1121 /* this "8" is from the packing in the vdwtab array - perhaps
1122 should be #define'ed? */
1123 offset = 8*ri + offstart;
1124 y0 = vdwtab[offset];
1125 f = vdwtab[offset+1];
1126 g = vdwtab[offset+2];
1127 h = vdwtab[offset+3];
1129 enersum += y0*(ea/3 + eb/2 + ec) + f*(ea/4 + eb/3 + ec/2)+
1130 g*(ea/5 + eb/4 + ec/3) + h*(ea/6 + eb/5 + ec/4);
1131 virsum += f*(pa/4 + pb/3 + pc/2 + pd) +
1132 2*g*(pa/5 + pb/4 + pc/3 + pd/2) + 3*h*(pa/6 + pb/5 + pc/4 + pd/3);
1135 enersum *= 4.0*M_PI;
1137 eners[i] -= enersum;
1141 /* now add the correction for rvdw_switch to infinity */
1142 eners[0] += -4.0*M_PI/(3.0*rc3);
1143 eners[1] += 4.0*M_PI/(9.0*rc9);
1144 virs[0] += 8.0*M_PI/rc3;
1145 virs[1] += -16.0*M_PI/(3.0*rc9);
1147 else if ((fr->vdwtype == evdwCUT) || (fr->vdwtype == evdwUSER)) {
1148 if (fr->vdwtype == evdwUSER && fplog)
1150 "WARNING: using dispersion correction with user tables\n");
1151 rc3 = fr->rvdw*fr->rvdw*fr->rvdw;
1153 eners[0] += -4.0*M_PI/(3.0*rc3);
1154 eners[1] += 4.0*M_PI/(9.0*rc9);
1155 virs[0] += 8.0*M_PI/rc3;
1156 virs[1] += -16.0*M_PI/(3.0*rc9);
1159 "Dispersion correction is not implemented for vdw-type = %s",
1160 evdw_names[fr->vdwtype]);
1162 fr->enerdiffsix = eners[0];
1163 fr->enerdifftwelve = eners[1];
1164 /* The 0.5 is due to the Gromacs definition of the virial */
1165 fr->virdiffsix = 0.5*virs[0];
1166 fr->virdifftwelve = 0.5*virs[1];
1170 void calc_dispcorr(FILE *fplog,t_inputrec *ir,t_forcerec *fr,
1171 gmx_large_int_t step,int natoms,
1172 matrix box,real lambda,tensor pres,tensor virial,
1173 real *prescorr, real *enercorr, real *dvdlcorr)
1175 gmx_bool bCorrAll,bCorrPres;
1176 real dvdlambda,invvol,dens,ninter,avcsix,avctwelve,enerdiff,svir=0,spres=0;
1186 if (ir->eDispCorr != edispcNO) {
1187 bCorrAll = (ir->eDispCorr == edispcAllEner ||
1188 ir->eDispCorr == edispcAllEnerPres);
1189 bCorrPres = (ir->eDispCorr == edispcEnerPres ||
1190 ir->eDispCorr == edispcAllEnerPres);
1192 invvol = 1/det(box);
1195 /* Only correct for the interactions with the inserted molecule */
1196 dens = (natoms - fr->n_tpi)*invvol;
1201 dens = natoms*invvol;
1202 ninter = 0.5*natoms;
1205 if (ir->efep == efepNO)
1207 avcsix = fr->avcsix[0];
1208 avctwelve = fr->avctwelve[0];
1212 avcsix = (1 - lambda)*fr->avcsix[0] + lambda*fr->avcsix[1];
1213 avctwelve = (1 - lambda)*fr->avctwelve[0] + lambda*fr->avctwelve[1];
1216 enerdiff = ninter*(dens*fr->enerdiffsix - fr->enershiftsix);
1217 *enercorr += avcsix*enerdiff;
1219 if (ir->efep != efepNO)
1221 dvdlambda += (fr->avcsix[1] - fr->avcsix[0])*enerdiff;
1225 enerdiff = ninter*(dens*fr->enerdifftwelve - fr->enershifttwelve);
1226 *enercorr += avctwelve*enerdiff;
1227 if (fr->efep != efepNO)
1229 dvdlambda += (fr->avctwelve[1] - fr->avctwelve[0])*enerdiff;
1235 svir = ninter*dens*avcsix*fr->virdiffsix/3.0;
1236 if (ir->eDispCorr == edispcAllEnerPres)
1238 svir += ninter*dens*avctwelve*fr->virdifftwelve/3.0;
1240 /* The factor 2 is because of the Gromacs virial definition */
1241 spres = -2.0*invvol*svir*PRESFAC;
1243 for(m=0; m<DIM; m++) {
1244 virial[m][m] += svir;
1245 pres[m][m] += spres;
1250 /* Can't currently control when it prints, for now, just print when degugging */
1254 fprintf(debug,"Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
1260 "Long Range LJ corr.: Epot %10g, Pres: %10g, Vir: %10g\n",
1261 *enercorr,spres,svir);
1265 fprintf(debug,"Long Range LJ corr.: Epot %10g\n",*enercorr);
1269 if (fr->bSepDVDL && do_per_step(step,ir->nstlog))
1271 fprintf(fplog,sepdvdlformat,"Dispersion correction",
1272 *enercorr,dvdlambda);
1274 if (fr->efep != efepNO)
1276 *dvdlcorr += dvdlambda;
1281 void do_pbc_first(FILE *fplog,matrix box,t_forcerec *fr,
1282 t_graph *graph,rvec x[])
1285 fprintf(fplog,"Removing pbc first time\n");
1286 calc_shifts(box,fr->shift_vec);
1288 mk_mshift(fplog,graph,fr->ePBC,box,x);
1290 p_graph(debug,"do_pbc_first 1",graph);
1291 shift_self(graph,box,x);
1292 /* By doing an extra mk_mshift the molecules that are broken
1293 * because they were e.g. imported from another software
1294 * will be made whole again. Such are the healing powers
1297 mk_mshift(fplog,graph,fr->ePBC,box,x);
1299 p_graph(debug,"do_pbc_first 2",graph);
1302 fprintf(fplog,"Done rmpbc\n");
1305 static void low_do_pbc_mtop(FILE *fplog,int ePBC,matrix box,
1306 gmx_mtop_t *mtop,rvec x[],
1311 gmx_molblock_t *molb;
1313 if (bFirst && fplog)
1314 fprintf(fplog,"Removing pbc first time\n");
1318 for(mb=0; mb<mtop->nmolblock; mb++) {
1319 molb = &mtop->molblock[mb];
1320 if (molb->natoms_mol == 1 ||
1321 (!bFirst && mtop->moltype[molb->type].cgs.nr == 1)) {
1322 /* Just one atom or charge group in the molecule, no PBC required */
1323 as += molb->nmol*molb->natoms_mol;
1325 /* Pass NULL iso fplog to avoid graph prints for each molecule type */
1326 mk_graph_ilist(NULL,mtop->moltype[molb->type].ilist,
1327 0,molb->natoms_mol,FALSE,FALSE,graph);
1329 for(mol=0; mol<molb->nmol; mol++) {
1330 mk_mshift(fplog,graph,ePBC,box,x+as);
1332 shift_self(graph,box,x+as);
1333 /* The molecule is whole now.
1334 * We don't need the second mk_mshift call as in do_pbc_first,
1335 * since we no longer need this graph.
1338 as += molb->natoms_mol;
1346 void do_pbc_first_mtop(FILE *fplog,int ePBC,matrix box,
1347 gmx_mtop_t *mtop,rvec x[])
1349 low_do_pbc_mtop(fplog,ePBC,box,mtop,x,TRUE);
1352 void do_pbc_mtop(FILE *fplog,int ePBC,matrix box,
1353 gmx_mtop_t *mtop,rvec x[])
1355 low_do_pbc_mtop(fplog,ePBC,box,mtop,x,FALSE);
1358 void finish_run(FILE *fplog,t_commrec *cr,const char *confout,
1359 t_inputrec *inputrec,
1360 t_nrnb nrnb[],gmx_wallcycle_t wcycle,
1361 gmx_runtime_t *runtime,
1362 gmx_bool bWriteStat)
1365 t_nrnb *nrnb_tot=NULL;
1368 double cycles[ewcNR];
1370 wallcycle_sum(cr,wcycle,cycles);
1372 if (cr->nnodes > 1) {
1376 MPI_Reduce(nrnb->n,nrnb_tot->n,eNRNB,MPI_DOUBLE,MPI_SUM,
1377 MASTERRANK(cr),cr->mpi_comm_mysim);
1383 if (SIMMASTER(cr)) {
1384 print_flop(fplog,nrnb_tot,&nbfs,&mflop);
1385 if (cr->nnodes > 1) {
1390 if ((cr->duty & DUTY_PP) && DOMAINDECOMP(cr)) {
1391 print_dd_statistics(cr,inputrec,fplog);
1403 snew(nrnb_all,cr->nnodes);
1404 nrnb_all[0] = *nrnb;
1405 for(s=1; s<cr->nnodes; s++)
1407 MPI_Recv(nrnb_all[s].n,eNRNB,MPI_DOUBLE,s,0,
1408 cr->mpi_comm_mysim,&stat);
1410 pr_load(fplog,cr,nrnb_all);
1415 MPI_Send(nrnb->n,eNRNB,MPI_DOUBLE,MASTERRANK(cr),0,
1416 cr->mpi_comm_mysim);
1421 if (SIMMASTER(cr)) {
1422 wallcycle_print(fplog,cr->nnodes,cr->npmenodes,runtime->realtime,
1425 if (EI_DYNAMICS(inputrec->eI)) {
1426 delta_t = inputrec->delta_t;
1432 print_perf(fplog,runtime->proctime,runtime->realtime,
1433 cr->nnodes-cr->npmenodes,
1434 runtime->nsteps_done,delta_t,nbfs,mflop);
1437 print_perf(stderr,runtime->proctime,runtime->realtime,
1438 cr->nnodes-cr->npmenodes,
1439 runtime->nsteps_done,delta_t,nbfs,mflop);
1443 runtime=inputrec->nsteps*inputrec->delta_t;
1445 if (cr->nnodes == 1)
1446 fprintf(stderr,"\n\n");
1447 print_perf(stderr,nodetime,realtime,runtime,&ntot,
1448 cr->nnodes-cr->npmenodes,FALSE);
1450 wallcycle_print(fplog,cr->nnodes,cr->npmenodes,realtime,wcycle,cycles);
1451 print_perf(fplog,nodetime,realtime,runtime,&ntot,cr->nnodes-cr->npmenodes,
1454 pr_load(fplog,cr,nrnb_all);
1461 void init_md(FILE *fplog,
1462 t_commrec *cr,t_inputrec *ir,const output_env_t oenv,
1463 double *t,double *t0,
1464 real *lambda,double *lam0,
1465 t_nrnb *nrnb,gmx_mtop_t *mtop,
1467 int nfile,const t_filenm fnm[],
1468 gmx_mdoutf_t **outf,t_mdebin **mdebin,
1469 tensor force_vir,tensor shake_vir,rvec mu_tot,
1470 gmx_bool *bSimAnn,t_vcm **vcm, t_state *state, unsigned long Flags)
1475 /* Initial values */
1476 *t = *t0 = ir->init_t;
1477 if (ir->efep != efepNO)
1479 *lam0 = ir->init_lambda;
1480 *lambda = *lam0 + ir->init_step*ir->delta_lambda;
1484 *lambda = *lam0 = 0.0;
1488 for(i=0;i<ir->opts.ngtc;i++)
1490 /* set bSimAnn if any group is being annealed */
1491 if(ir->opts.annealing[i]!=eannNO)
1498 update_annealing_target_temp(&(ir->opts),ir->init_t);
1503 *upd = init_update(fplog,ir);
1508 *vcm = init_vcm(fplog,&mtop->groups,ir);
1511 if (EI_DYNAMICS(ir->eI) && !(Flags & MD_APPENDFILES))
1513 if (ir->etc == etcBERENDSEN)
1515 please_cite(fplog,"Berendsen84a");
1517 if (ir->etc == etcVRESCALE)
1519 please_cite(fplog,"Bussi2007a");
1527 *outf = init_mdoutf(nfile,fnm,Flags,cr,ir,oenv);
1529 *mdebin = init_mdebin((Flags & MD_APPENDFILES) ? NULL : (*outf)->fp_ene,
1530 mtop,ir, (*outf)->fp_dhdl);
1533 /* Initiate variables */
1534 clear_mat(force_vir);
1535 clear_mat(shake_vir);