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"
84 #include "mpelogging.h"
87 #include "gmx_wallcycle.h"
100 typedef struct gmx_timeprint {
105 /* Portable version of ctime_r implemented in src/gmxlib/string2.c, but we do not want it declared in public installed headers */
107 gmx_ctime_r(const time_t *clock,char *buf, int n);
113 #ifdef HAVE_GETTIMEOFDAY
115 struct timezone tz = { 0,0 };
118 gettimeofday(&t,&tz);
120 seconds = (double) t.tv_sec + 1e-6*(double)t.tv_usec;
126 seconds = time(NULL);
133 #define difftime(end,start) ((double)(end)-(double)(start))
135 void print_time(FILE *out,gmx_runtime_t *runtime,gmx_large_int_t step,
136 t_inputrec *ir, t_commrec *cr)
139 char timebuf[STRLEN];
149 fprintf(out,"step %s",gmx_step_str(step,buf));
150 if ((step >= ir->nstlist))
152 if ((ir->nstlist == 0) || ((step % ir->nstlist) == 0))
154 /* We have done a full cycle let's update time_per_step */
155 runtime->last = gmx_gettime();
156 dt = difftime(runtime->last,runtime->real);
157 runtime->time_per_step = dt/(step - ir->init_step + 1);
159 dt = (ir->nsteps + ir->init_step - step)*runtime->time_per_step;
165 finish = (time_t) (runtime->last + dt);
166 gmx_ctime_r(&finish,timebuf,STRLEN);
167 sprintf(buf,"%s",timebuf);
168 buf[strlen(buf)-1]='\0';
169 fprintf(out,", will finish %s",buf);
172 fprintf(out,", remaining runtime: %5d s ",(int)dt);
176 fprintf(out," performance: %.1f ns/day ",
177 ir->delta_t/1000*24*60*60/runtime->time_per_step);
194 static double set_proctime(gmx_runtime_t *runtime)
200 prev = runtime->proc;
201 runtime->proc = dclock();
203 diff = runtime->proc - prev;
207 prev = runtime->proc;
208 runtime->proc = clock();
210 diff = (double)(runtime->proc - prev)/(double)CLOCKS_PER_SEC;
214 /* The counter has probably looped, ignore this data */
221 void runtime_start(gmx_runtime_t *runtime)
223 runtime->real = gmx_gettime();
225 set_proctime(runtime);
226 runtime->realtime = 0;
227 runtime->proctime = 0;
229 runtime->time_per_step = 0;
232 void runtime_end(gmx_runtime_t *runtime)
238 runtime->proctime += set_proctime(runtime);
239 runtime->realtime = now - runtime->real;
243 void runtime_upd_proc(gmx_runtime_t *runtime)
245 runtime->proctime += set_proctime(runtime);
248 void print_date_and_time(FILE *fplog,int nodeid,const char *title,
249 const gmx_runtime_t *runtime)
252 char timebuf[STRLEN];
253 char time_string[STRLEN];
260 tmptime = (time_t) runtime->real;
261 gmx_ctime_r(&tmptime,timebuf,STRLEN);
265 tmptime = (time_t) gmx_gettime();
266 gmx_ctime_r(&tmptime,timebuf,STRLEN);
268 for(i=0; timebuf[i]>=' '; i++)
270 time_string[i]=timebuf[i];
274 fprintf(fplog,"%s on node %d %s\n",title,nodeid,time_string);
278 static void sum_forces(int start,int end,rvec f[],rvec flr[])
283 pr_rvecs(debug,0,"fsr",f+start,end-start);
284 pr_rvecs(debug,0,"flr",flr+start,end-start);
286 for(i=start; (i<end); i++)
287 rvec_inc(f[i],flr[i]);
291 * calc_f_el calculates forces due to an electric field.
293 * force is kJ mol^-1 nm^-1 = e * kJ mol^-1 nm^-1 / e
295 * Et[] contains the parameters for the time dependent
296 * part of the field (not yet used).
297 * Ex[] contains the parameters for
298 * the spatial dependent part of the field. You can have cool periodic
299 * fields in principle, but only a constant field is supported
301 * The function should return the energy due to the electric field
302 * (if any) but for now returns 0.
305 * There can be problems with the virial.
306 * Since the field is not self-consistent this is unavoidable.
307 * For neutral molecules the virial is correct within this approximation.
308 * For neutral systems with many charged molecules the error is small.
309 * But for systems with a net charge or a few charged molecules
310 * the error can be significant when the field is high.
311 * Solution: implement a self-consitent electric field into PME.
313 static void calc_f_el(FILE *fp,int start,int homenr,
314 real charge[],rvec x[],rvec f[],
315 t_cosines Ex[],t_cosines Et[],double t)
321 for(m=0; (m<DIM); m++)
328 Ext[m] = cos(Et[m].a[0]*(t-t0))*exp(-sqr(t-t0)/(2.0*sqr(Et[m].a[2])));
332 Ext[m] = cos(Et[m].a[0]*t);
341 /* Convert the field strength from V/nm to MD-units */
342 Ext[m] *= Ex[m].a[0]*FIELDFAC;
343 for(i=start; (i<start+homenr); i++)
344 f[i][m] += charge[i]*Ext[m];
353 fprintf(fp,"%10g %10g %10g %10g #FIELD\n",t,
354 Ext[XX]/FIELDFAC,Ext[YY]/FIELDFAC,Ext[ZZ]/FIELDFAC);
358 static void calc_virial(FILE *fplog,int start,int homenr,rvec x[],rvec f[],
359 tensor vir_part,t_graph *graph,matrix box,
360 t_nrnb *nrnb,const t_forcerec *fr,int ePBC)
365 /* The short-range virial from surrounding boxes */
367 calc_vir(fplog,SHIFTS,fr->shift_vec,fr->fshift,vir_part,ePBC==epbcSCREW,box);
368 inc_nrnb(nrnb,eNR_VIRIAL,SHIFTS);
370 /* Calculate partial virial, for local atoms only, based on short range.
371 * Total virial is computed in global_stat, called from do_md
373 f_calc_vir(fplog,start,start+homenr,x,f,vir_part,graph,box);
374 inc_nrnb(nrnb,eNR_VIRIAL,homenr);
376 /* Add position restraint contribution */
377 for(i=0; i<DIM; i++) {
378 vir_part[i][i] += fr->vir_diag_posres[i];
381 /* Add wall contribution */
382 for(i=0; i<DIM; i++) {
383 vir_part[i][ZZ] += fr->vir_wall_z[i];
387 pr_rvecs(debug,0,"vir_part",vir_part,DIM);
390 static void print_large_forces(FILE *fp,t_mdatoms *md,t_commrec *cr,
391 gmx_large_int_t step,real pforce,rvec *x,rvec *f)
395 char buf[STEPSTRSIZE];
398 for(i=md->start; i<md->start+md->homenr; i++) {
400 /* We also catch NAN, if the compiler does not optimize this away. */
401 if (fn2 >= pf2 || fn2 != fn2) {
402 fprintf(fp,"step %s atom %6d x %8.3f %8.3f %8.3f force %12.5e\n",
403 gmx_step_str(step,buf),
404 ddglatnr(cr->dd,i),x[i][XX],x[i][YY],x[i][ZZ],sqrt(fn2));
409 void do_force(FILE *fplog,t_commrec *cr,
410 t_inputrec *inputrec,
411 gmx_large_int_t step,t_nrnb *nrnb,gmx_wallcycle_t wcycle,
414 gmx_groups_t *groups,
415 matrix box,rvec x[],history_t *hist,
419 gmx_enerdata_t *enerd,t_fcdata *fcd,
420 real lambda,t_graph *graph,
421 t_forcerec *fr,gmx_vsite_t *vsite,rvec mu_tot,
422 double t,FILE *field,gmx_edsam_t ed,
429 gmx_bool bSepDVDL,bStateChanged,bNS,bFillGrid,bCalcCGCM,bBS;
430 gmx_bool bDoLongRange,bDoForces,bSepLRF;
434 float cycles_ppdpme,cycles_pme,cycles_seppme,cycles_force;
436 start = mdatoms->start;
437 homenr = mdatoms->homenr;
439 bSepDVDL = (fr->bSepDVDL && do_per_step(step,inputrec->nstlog));
441 clear_mat(vir_force);
445 pd_cg_range(cr,&cg0,&cg1);
450 if (DOMAINDECOMP(cr))
452 cg1 = cr->dd->ncg_tot;
464 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
465 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll==FALSE);
466 bFillGrid = (bNS && bStateChanged);
467 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
468 bDoLongRange = (fr->bTwinRange && bNS && (flags & GMX_FORCE_DOLR));
469 bDoForces = (flags & GMX_FORCE_FORCES);
470 bSepLRF = (bDoLongRange && bDoForces && (flags & GMX_FORCE_SEPLRF));
474 update_forcerec(fplog,fr,box);
476 /* Calculate total (local) dipole moment in a temporary common array.
477 * This makes it possible to sum them over nodes faster.
479 calc_mu(start,homenr,
480 x,mdatoms->chargeA,mdatoms->chargeB,mdatoms->nChargePerturbed,
484 if (fr->ePBC != epbcNONE) {
485 /* Compute shift vectors every step,
486 * because of pressure coupling or box deformation!
488 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
489 calc_shifts(box,fr->shift_vec);
492 put_charge_groups_in_box(fplog,cg0,cg1,fr->ePBC,box,
493 &(top->cgs),x,fr->cg_cm);
494 inc_nrnb(nrnb,eNR_CGCM,homenr);
495 inc_nrnb(nrnb,eNR_RESETX,cg1-cg0);
497 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph) {
498 unshift_self(graph,box,x);
501 else if (bCalcCGCM) {
502 calc_cgcm(fplog,cg0,cg1,&(top->cgs),x,fr->cg_cm);
503 inc_nrnb(nrnb,eNR_CGCM,homenr);
508 move_cgcm(fplog,cr,fr->cg_cm);
511 pr_rvecs(debug,0,"cgcm",fr->cg_cm,top->cgs.nr);
515 if (!(cr->duty & DUTY_PME)) {
516 /* Send particle coordinates to the pme nodes.
517 * Since this is only implemented for domain decomposition
518 * and domain decomposition does not use the graph,
519 * we do not need to worry about shifting.
522 wallcycle_start(wcycle,ewcPP_PMESENDX);
523 GMX_MPE_LOG(ev_send_coordinates_start);
525 bBS = (inputrec->nwall == 2);
528 svmul(inputrec->wall_ewald_zfac,boxs[ZZ],boxs[ZZ]);
531 gmx_pme_send_x(cr,bBS ? boxs : box,x,
532 mdatoms->nChargePerturbed,lambda,
533 /* FIX ME after 4.5 */
534 /* we are using gmx_bool of type char */
535 ( flags & GMX_FORCE_VIRIAL) != 0,step);
537 GMX_MPE_LOG(ev_send_coordinates_finish);
538 wallcycle_stop(wcycle,ewcPP_PMESENDX);
542 /* Communicate coordinates and sum dipole if necessary */
545 wallcycle_start(wcycle,ewcMOVEX);
546 if (DOMAINDECOMP(cr))
548 dd_move_x(cr->dd,box,x);
552 move_x(fplog,cr,GMX_LEFT,GMX_RIGHT,x,nrnb);
554 /* When we don't need the total dipole we sum it in global_stat */
555 if (bStateChanged && NEED_MUTOT(*inputrec))
557 gmx_sumd(2*DIM,mu,cr);
559 wallcycle_stop(wcycle,ewcMOVEX);
567 fr->mu_tot[i][j] = mu[i*DIM + j];
571 if (fr->efep == efepNO)
573 copy_rvec(fr->mu_tot[0],mu_tot);
580 (1.0 - lambda)*fr->mu_tot[0][j] + lambda*fr->mu_tot[1][j];
585 reset_enerdata(&(inputrec->opts),fr,bNS,enerd,MASTER(cr));
586 clear_rvecs(SHIFTS,fr->fshift);
590 wallcycle_start(wcycle,ewcNS);
592 if (graph && bStateChanged)
594 /* Calculate intramolecular shift vectors to make molecules whole */
595 mk_mshift(fplog,graph,fr->ePBC,box,x);
598 /* Reset long range forces if necessary */
601 /* Reset the (long-range) forces if necessary */
602 clear_rvecs(fr->natoms_force_constr,bSepLRF ? fr->f_twin : f);
605 /* Do the actual neighbour searching and if twin range electrostatics
606 * also do the calculation of long range forces and energies.
610 groups,&(inputrec->opts),top,mdatoms,
611 cr,nrnb,lambda,&dvdl,&enerd->grpp,bFillGrid,
612 bDoLongRange,bDoForces,bSepLRF ? fr->f_twin : f);
615 fprintf(fplog,sepdvdlformat,"LR non-bonded",0.0,dvdl);
617 enerd->dvdl_lin += dvdl;
619 wallcycle_stop(wcycle,ewcNS);
622 if (inputrec->implicit_solvent && bNS)
624 make_gb_nblist(cr,inputrec->gb_algorithm,inputrec->rlist,
625 x,box,fr,&top->idef,graph,fr->born);
628 if (DOMAINDECOMP(cr))
630 if (!(cr->duty & DUTY_PME))
632 wallcycle_start(wcycle,ewcPPDURINGPME);
633 dd_force_flop_start(cr->dd,nrnb);
637 /* Start the force cycle counter.
638 * This counter is stopped in do_forcelow_level.
639 * No parallel communication should occur while this counter is running,
640 * since that will interfere with the dynamic load balancing.
642 wallcycle_start(wcycle,ewcFORCE);
646 /* Reset forces for which the virial is calculated separately:
647 * PME/Ewald forces if necessary */
650 if (flags & GMX_FORCE_VIRIAL)
652 fr->f_novirsum = fr->f_novirsum_alloc;
653 GMX_BARRIER(cr->mpi_comm_mygroup);
656 clear_rvecs(fr->f_novirsum_n,fr->f_novirsum);
660 clear_rvecs(homenr,fr->f_novirsum+start);
662 GMX_BARRIER(cr->mpi_comm_mygroup);
666 /* We are not calculating the pressure so we do not need
667 * a separate array for forces that do not contribute
676 /* Add the long range forces to the short range forces */
677 for(i=0; i<fr->natoms_force_constr; i++)
679 copy_rvec(fr->f_twin[i],f[i]);
682 else if (!(fr->bTwinRange && bNS))
684 /* Clear the short-range forces */
685 clear_rvecs(fr->natoms_force_constr,f);
688 clear_rvec(fr->vir_diag_posres);
690 GMX_BARRIER(cr->mpi_comm_mygroup);
692 if (inputrec->ePull == epullCONSTRAINT)
694 clear_pull_forces(inputrec->pull);
697 /* update QMMMrec, if necessary */
700 update_QMMMrec(cr,fr,x,mdatoms,box,top);
703 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
705 /* Position restraints always require full pbc */
706 set_pbc(&pbc,inputrec->ePBC,box);
707 v = posres(top->idef.il[F_POSRES].nr,top->idef.il[F_POSRES].iatoms,
708 top->idef.iparams_posres,
709 (const rvec*)x,fr->f_novirsum,fr->vir_diag_posres,
710 inputrec->ePBC==epbcNONE ? NULL : &pbc,lambda,&dvdl,
711 fr->rc_scaling,fr->ePBC,fr->posres_com,fr->posres_comB);
714 fprintf(fplog,sepdvdlformat,
715 interaction_function[F_POSRES].longname,v,dvdl);
717 enerd->term[F_POSRES] += v;
718 /* This linear lambda dependence assumption is only correct
719 * when only k depends on lambda,
720 * not when the reference position depends on lambda.
721 * grompp checks for this.
723 enerd->dvdl_lin += dvdl;
724 inc_nrnb(nrnb,eNR_POSRES,top->idef.il[F_POSRES].nr/2);
727 /* Compute the bonded and non-bonded energies and optionally forces */
728 do_force_lowlevel(fplog,step,fr,inputrec,&(top->idef),
729 cr,nrnb,wcycle,mdatoms,&(inputrec->opts),
730 x,hist,f,enerd,fcd,mtop,top,fr->born,
731 &(top->atomtypes),bBornRadii,box,
732 lambda,graph,&(top->excls),fr->mu_tot,
735 cycles_force = wallcycle_stop(wcycle,ewcFORCE);
736 GMX_BARRIER(cr->mpi_comm_mygroup);
740 do_flood(fplog,cr,x,f,ed,box,step);
743 if (DOMAINDECOMP(cr))
745 dd_force_flop_stop(cr->dd,nrnb);
748 dd_cycles_add(cr->dd,cycles_force-cycles_pme,ddCyclF);
754 if (IR_ELEC_FIELD(*inputrec))
756 /* Compute forces due to electric field */
757 calc_f_el(MASTER(cr) ? field : NULL,
758 start,homenr,mdatoms->chargeA,x,fr->f_novirsum,
759 inputrec->ex,inputrec->et,t);
762 /* Communicate the forces */
765 wallcycle_start(wcycle,ewcMOVEF);
766 if (DOMAINDECOMP(cr))
768 dd_move_f(cr->dd,f,fr->fshift);
769 /* Do we need to communicate the separate force array
770 * for terms that do not contribute to the single sum virial?
771 * Position restraints and electric fields do not introduce
772 * inter-cg forces, only full electrostatics methods do.
773 * When we do not calculate the virial, fr->f_novirsum = f,
774 * so we have already communicated these forces.
776 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
777 (flags & GMX_FORCE_VIRIAL))
779 dd_move_f(cr->dd,fr->f_novirsum,NULL);
783 /* We should not update the shift forces here,
784 * since f_twin is already included in f.
786 dd_move_f(cr->dd,fr->f_twin,NULL);
791 pd_move_f(cr,f,nrnb);
794 pd_move_f(cr,fr->f_twin,nrnb);
797 wallcycle_stop(wcycle,ewcMOVEF);
800 /* If we have NoVirSum forces, but we do not calculate the virial,
801 * we sum fr->f_novirum=f later.
803 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
805 wallcycle_start(wcycle,ewcVSITESPREAD);
806 spread_vsite_f(fplog,vsite,x,f,fr->fshift,nrnb,
807 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
808 wallcycle_stop(wcycle,ewcVSITESPREAD);
812 wallcycle_start(wcycle,ewcVSITESPREAD);
813 spread_vsite_f(fplog,vsite,x,fr->f_twin,NULL,
815 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
816 wallcycle_stop(wcycle,ewcVSITESPREAD);
820 if (flags & GMX_FORCE_VIRIAL)
822 /* Calculation of the virial must be done after vsites! */
823 calc_virial(fplog,mdatoms->start,mdatoms->homenr,x,f,
824 vir_force,graph,box,nrnb,fr,inputrec->ePBC);
828 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
830 /* Calculate the center of mass forces, this requires communication,
831 * which is why pull_potential is called close to other communication.
832 * The virial contribution is calculated directly,
833 * which is why we call pull_potential after calc_virial.
835 set_pbc(&pbc,inputrec->ePBC,box);
837 enerd->term[F_COM_PULL] =
838 pull_potential(inputrec->ePull,inputrec->pull,mdatoms,&pbc,
839 cr,t,lambda,x,f,vir_force,&dvdl);
842 fprintf(fplog,sepdvdlformat,"Com pull",enerd->term[F_COM_PULL],dvdl);
844 enerd->dvdl_lin += dvdl;
847 if (PAR(cr) && !(cr->duty & DUTY_PME))
849 cycles_ppdpme = wallcycle_stop(wcycle,ewcPPDURINGPME);
850 dd_cycles_add(cr->dd,cycles_ppdpme,ddCyclPPduringPME);
852 /* In case of node-splitting, the PP nodes receive the long-range
853 * forces, virial and energy from the PME nodes here.
855 wallcycle_start(wcycle,ewcPP_PMEWAITRECVF);
857 gmx_pme_receive_f(cr,fr->f_novirsum,fr->vir_el_recip,&e,&dvdl,
861 fprintf(fplog,sepdvdlformat,"PME mesh",e,dvdl);
863 enerd->term[F_COUL_RECIP] += e;
864 enerd->dvdl_lin += dvdl;
867 dd_cycles_add(cr->dd,cycles_seppme,ddCyclPME);
869 wallcycle_stop(wcycle,ewcPP_PMEWAITRECVF);
872 if (bDoForces && fr->bF_NoVirSum)
876 /* Spread the mesh force on virtual sites to the other particles...
877 * This is parallellized. MPI communication is performed
878 * if the constructing atoms aren't local.
880 wallcycle_start(wcycle,ewcVSITESPREAD);
881 spread_vsite_f(fplog,vsite,x,fr->f_novirsum,NULL,nrnb,
882 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
883 wallcycle_stop(wcycle,ewcVSITESPREAD);
885 if (flags & GMX_FORCE_VIRIAL)
887 /* Now add the forces, this is local */
890 sum_forces(0,fr->f_novirsum_n,f,fr->f_novirsum);
894 sum_forces(start,start+homenr,f,fr->f_novirsum);
896 if (EEL_FULL(fr->eeltype))
898 /* Add the mesh contribution to the virial */
899 m_add(vir_force,fr->vir_el_recip,vir_force);
903 pr_rvecs(debug,0,"vir_force",vir_force,DIM);
908 /* Sum the potential energy terms from group contributions */
909 sum_epot(&(inputrec->opts),enerd);
911 if (fr->print_force >= 0 && bDoForces)
913 print_large_forces(stderr,mdatoms,cr,step,fr->print_force,x,f);
917 void do_constrain_first(FILE *fplog,gmx_constr_t constr,
918 t_inputrec *ir,t_mdatoms *md,
919 t_state *state,rvec *f,
920 t_graph *graph,t_commrec *cr,t_nrnb *nrnb,
921 t_forcerec *fr, gmx_localtop_t *top, tensor shake_vir)
924 gmx_large_int_t step;
925 double mass,tmass,vcm[4];
930 snew(savex,state->natoms);
933 end = md->homenr + start;
936 fprintf(debug,"vcm: start=%d, homenr=%d, end=%d\n",
937 start,md->homenr,end);
938 /* Do a first constrain to reset particles... */
939 step = ir->init_step;
942 char buf[STEPSTRSIZE];
943 fprintf(fplog,"\nConstraining the starting coordinates (step %s)\n",
944 gmx_step_str(step,buf));
948 /* constrain the current position */
949 constrain(NULL,TRUE,FALSE,constr,&(top->idef),
950 ir,NULL,cr,step,0,md,
951 state->x,state->x,NULL,
952 state->box,state->lambda,&dvdlambda,
953 NULL,NULL,nrnb,econqCoord,ir->epc==epcMTTK,state->veta,state->veta);
956 /* constrain the inital velocity, and save it */
957 /* also may be useful if we need the ekin from the halfstep for velocity verlet */
958 /* might not yet treat veta correctly */
959 constrain(NULL,TRUE,FALSE,constr,&(top->idef),
960 ir,NULL,cr,step,0,md,
961 state->x,state->v,state->v,
962 state->box,state->lambda,&dvdlambda,
963 NULL,NULL,nrnb,econqVeloc,ir->epc==epcMTTK,state->veta,state->veta);
965 /* constrain the inital velocities at t-dt/2 */
966 if (EI_STATE_VELOCITY(ir->eI) && ir->eI!=eiVV)
968 for(i=start; (i<end); i++)
970 for(m=0; (m<DIM); m++)
972 /* Reverse the velocity */
973 state->v[i][m] = -state->v[i][m];
974 /* Store the position at t-dt in buf */
975 savex[i][m] = state->x[i][m] + dt*state->v[i][m];
978 /* Shake the positions at t=-dt with the positions at t=0
979 * as reference coordinates.
983 char buf[STEPSTRSIZE];
984 fprintf(fplog,"\nConstraining the coordinates at t0-dt (step %s)\n",
985 gmx_step_str(step,buf));
988 constrain(NULL,TRUE,FALSE,constr,&(top->idef),
989 ir,NULL,cr,step,-1,md,
991 state->box,state->lambda,&dvdlambda,
992 state->v,NULL,nrnb,econqCoord,ir->epc==epcMTTK,state->veta,state->veta);
994 for(i=start; i<end; i++) {
995 for(m=0; m<DIM; m++) {
996 /* Re-reverse the velocities */
997 state->v[i][m] = -state->v[i][m];
1002 for(m=0; (m<4); m++)
1004 for(i=start; i<end; i++) {
1005 mass = md->massT[i];
1006 for(m=0; m<DIM; m++) {
1007 vcm[m] += state->v[i][m]*mass;
1012 if (ir->nstcomm != 0 || debug) {
1013 /* Compute the global sum of vcm */
1015 fprintf(debug,"vcm: %8.3f %8.3f %8.3f,"
1016 " total mass = %12.5e\n",vcm[XX],vcm[YY],vcm[ZZ],vcm[3]);
1020 for(m=0; (m<DIM); m++)
1023 fprintf(debug,"vcm: %8.3f %8.3f %8.3f,"
1024 " total mass = %12.5e\n",vcm[XX],vcm[YY],vcm[ZZ],tmass);
1025 if (ir->nstcomm != 0) {
1026 /* Now we have the velocity of center of mass, let's remove it */
1027 for(i=start; (i<end); i++) {
1028 for(m=0; (m<DIM); m++)
1029 state->v[i][m] -= vcm[m];
1037 void calc_enervirdiff(FILE *fplog,int eDispCorr,t_forcerec *fr)
1039 double eners[2],virs[2],enersum,virsum,y0,f,g,h;
1040 double r0,r1,r,rc3,rc9,ea,eb,ec,pa,pb,pc,pd;
1041 double invscale,invscale2,invscale3;
1042 int ri0,ri1,ri,i,offstart,offset;
1045 fr->enershiftsix = 0;
1046 fr->enershifttwelve = 0;
1047 fr->enerdiffsix = 0;
1048 fr->enerdifftwelve = 0;
1050 fr->virdifftwelve = 0;
1052 if (eDispCorr != edispcNO) {
1053 for(i=0; i<2; i++) {
1057 if ((fr->vdwtype == evdwSWITCH) || (fr->vdwtype == evdwSHIFT)) {
1058 if (fr->rvdw_switch == 0)
1060 "With dispersion correction rvdw-switch can not be zero "
1061 "for vdw-type = %s",evdw_names[fr->vdwtype]);
1063 scale = fr->nblists[0].tab.scale;
1064 vdwtab = fr->nblists[0].vdwtab;
1066 /* Round the cut-offs to exact table values for precision */
1067 ri0 = floor(fr->rvdw_switch*scale);
1068 ri1 = ceil(fr->rvdw*scale);
1074 if (fr->vdwtype == evdwSHIFT) {
1075 /* Determine the constant energy shift below rvdw_switch */
1076 fr->enershiftsix = (real)(-1.0/(rc3*rc3)) - vdwtab[8*ri0];
1077 fr->enershifttwelve = (real)( 1.0/(rc9*rc3)) - vdwtab[8*ri0 + 4];
1079 /* Add the constant part from 0 to rvdw_switch.
1080 * This integration from 0 to rvdw_switch overcounts the number
1081 * of interactions by 1, as it also counts the self interaction.
1082 * We will correct for this later.
1084 eners[0] += 4.0*M_PI*fr->enershiftsix*rc3/3.0;
1085 eners[1] += 4.0*M_PI*fr->enershifttwelve*rc3/3.0;
1087 invscale = 1.0/(scale);
1088 invscale2 = invscale*invscale;
1089 invscale3 = invscale*invscale2;
1091 /* following summation derived from cubic spline definition,
1092 Numerical Recipies in C, second edition, p. 113-116. Exact
1093 for the cubic spline. We first calculate the negative of
1094 the energy from rvdw to rvdw_switch, assuming that g(r)=1,
1095 and then add the more standard, abrupt cutoff correction to
1096 that result, yielding the long-range correction for a
1097 switched function. We perform both the pressure and energy
1098 loops at the same time for simplicity, as the computational
1102 enersum = 0.0; virsum = 0.0;
1107 for (ri=ri0; ri<ri1; ri++) {
1110 eb = 2.0*invscale2*r;
1114 pb = 3.0*invscale2*r;
1115 pc = 3.0*invscale*r*r;
1118 /* this "8" is from the packing in the vdwtab array - perhaps
1119 should be #define'ed? */
1120 offset = 8*ri + offstart;
1121 y0 = vdwtab[offset];
1122 f = vdwtab[offset+1];
1123 g = vdwtab[offset+2];
1124 h = vdwtab[offset+3];
1126 enersum += y0*(ea/3 + eb/2 + ec) + f*(ea/4 + eb/3 + ec/2)+
1127 g*(ea/5 + eb/4 + ec/3) + h*(ea/6 + eb/5 + ec/4);
1128 virsum += f*(pa/4 + pb/3 + pc/2 + pd) +
1129 2*g*(pa/5 + pb/4 + pc/3 + pd/2) + 3*h*(pa/6 + pb/5 + pc/4 + pd/3);
1132 enersum *= 4.0*M_PI;
1134 eners[i] -= enersum;
1138 /* now add the correction for rvdw_switch to infinity */
1139 eners[0] += -4.0*M_PI/(3.0*rc3);
1140 eners[1] += 4.0*M_PI/(9.0*rc9);
1141 virs[0] += 8.0*M_PI/rc3;
1142 virs[1] += -16.0*M_PI/(3.0*rc9);
1144 else if ((fr->vdwtype == evdwCUT) || (fr->vdwtype == evdwUSER)) {
1145 if (fr->vdwtype == evdwUSER && fplog)
1147 "WARNING: using dispersion correction with user tables\n");
1148 rc3 = fr->rvdw*fr->rvdw*fr->rvdw;
1150 eners[0] += -4.0*M_PI/(3.0*rc3);
1151 eners[1] += 4.0*M_PI/(9.0*rc9);
1152 virs[0] += 8.0*M_PI/rc3;
1153 virs[1] += -16.0*M_PI/(3.0*rc9);
1156 "Dispersion correction is not implemented for vdw-type = %s",
1157 evdw_names[fr->vdwtype]);
1159 fr->enerdiffsix = eners[0];
1160 fr->enerdifftwelve = eners[1];
1161 /* The 0.5 is due to the Gromacs definition of the virial */
1162 fr->virdiffsix = 0.5*virs[0];
1163 fr->virdifftwelve = 0.5*virs[1];
1167 void calc_dispcorr(FILE *fplog,t_inputrec *ir,t_forcerec *fr,
1168 gmx_large_int_t step,int natoms,
1169 matrix box,real lambda,tensor pres,tensor virial,
1170 real *prescorr, real *enercorr, real *dvdlcorr)
1172 gmx_bool bCorrAll,bCorrPres;
1173 real dvdlambda,invvol,dens,ninter,avcsix,avctwelve,enerdiff,svir=0,spres=0;
1183 if (ir->eDispCorr != edispcNO) {
1184 bCorrAll = (ir->eDispCorr == edispcAllEner ||
1185 ir->eDispCorr == edispcAllEnerPres);
1186 bCorrPres = (ir->eDispCorr == edispcEnerPres ||
1187 ir->eDispCorr == edispcAllEnerPres);
1189 invvol = 1/det(box);
1192 /* Only correct for the interactions with the inserted molecule */
1193 dens = (natoms - fr->n_tpi)*invvol;
1198 dens = natoms*invvol;
1199 ninter = 0.5*natoms;
1202 if (ir->efep == efepNO)
1204 avcsix = fr->avcsix[0];
1205 avctwelve = fr->avctwelve[0];
1209 avcsix = (1 - lambda)*fr->avcsix[0] + lambda*fr->avcsix[1];
1210 avctwelve = (1 - lambda)*fr->avctwelve[0] + lambda*fr->avctwelve[1];
1213 enerdiff = ninter*(dens*fr->enerdiffsix - fr->enershiftsix);
1214 *enercorr += avcsix*enerdiff;
1216 if (ir->efep != efepNO)
1218 dvdlambda += (fr->avcsix[1] - fr->avcsix[0])*enerdiff;
1222 enerdiff = ninter*(dens*fr->enerdifftwelve - fr->enershifttwelve);
1223 *enercorr += avctwelve*enerdiff;
1224 if (fr->efep != efepNO)
1226 dvdlambda += (fr->avctwelve[1] - fr->avctwelve[0])*enerdiff;
1232 svir = ninter*dens*avcsix*fr->virdiffsix/3.0;
1233 if (ir->eDispCorr == edispcAllEnerPres)
1235 svir += ninter*dens*avctwelve*fr->virdifftwelve/3.0;
1237 /* The factor 2 is because of the Gromacs virial definition */
1238 spres = -2.0*invvol*svir*PRESFAC;
1240 for(m=0; m<DIM; m++) {
1241 virial[m][m] += svir;
1242 pres[m][m] += spres;
1247 /* Can't currently control when it prints, for now, just print when degugging */
1251 fprintf(debug,"Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
1257 "Long Range LJ corr.: Epot %10g, Pres: %10g, Vir: %10g\n",
1258 *enercorr,spres,svir);
1262 fprintf(debug,"Long Range LJ corr.: Epot %10g\n",*enercorr);
1266 if (fr->bSepDVDL && do_per_step(step,ir->nstlog))
1268 fprintf(fplog,sepdvdlformat,"Dispersion correction",
1269 *enercorr,dvdlambda);
1271 if (fr->efep != efepNO)
1273 *dvdlcorr += dvdlambda;
1278 void do_pbc_first(FILE *fplog,matrix box,t_forcerec *fr,
1279 t_graph *graph,rvec x[])
1282 fprintf(fplog,"Removing pbc first time\n");
1283 calc_shifts(box,fr->shift_vec);
1285 mk_mshift(fplog,graph,fr->ePBC,box,x);
1287 p_graph(debug,"do_pbc_first 1",graph);
1288 shift_self(graph,box,x);
1289 /* By doing an extra mk_mshift the molecules that are broken
1290 * because they were e.g. imported from another software
1291 * will be made whole again. Such are the healing powers
1294 mk_mshift(fplog,graph,fr->ePBC,box,x);
1296 p_graph(debug,"do_pbc_first 2",graph);
1299 fprintf(fplog,"Done rmpbc\n");
1302 static void low_do_pbc_mtop(FILE *fplog,int ePBC,matrix box,
1303 gmx_mtop_t *mtop,rvec x[],
1308 gmx_molblock_t *molb;
1310 if (bFirst && fplog)
1311 fprintf(fplog,"Removing pbc first time\n");
1315 for(mb=0; mb<mtop->nmolblock; mb++) {
1316 molb = &mtop->molblock[mb];
1317 if (molb->natoms_mol == 1 ||
1318 (!bFirst && mtop->moltype[molb->type].cgs.nr == 1)) {
1319 /* Just one atom or charge group in the molecule, no PBC required */
1320 as += molb->nmol*molb->natoms_mol;
1322 /* Pass NULL iso fplog to avoid graph prints for each molecule type */
1323 mk_graph_ilist(NULL,mtop->moltype[molb->type].ilist,
1324 0,molb->natoms_mol,FALSE,FALSE,graph);
1326 for(mol=0; mol<molb->nmol; mol++) {
1327 mk_mshift(fplog,graph,ePBC,box,x+as);
1329 shift_self(graph,box,x+as);
1330 /* The molecule is whole now.
1331 * We don't need the second mk_mshift call as in do_pbc_first,
1332 * since we no longer need this graph.
1335 as += molb->natoms_mol;
1343 void do_pbc_first_mtop(FILE *fplog,int ePBC,matrix box,
1344 gmx_mtop_t *mtop,rvec x[])
1346 low_do_pbc_mtop(fplog,ePBC,box,mtop,x,TRUE);
1349 void do_pbc_mtop(FILE *fplog,int ePBC,matrix box,
1350 gmx_mtop_t *mtop,rvec x[])
1352 low_do_pbc_mtop(fplog,ePBC,box,mtop,x,FALSE);
1355 void finish_run(FILE *fplog,t_commrec *cr,const char *confout,
1356 t_inputrec *inputrec,
1357 t_nrnb nrnb[],gmx_wallcycle_t wcycle,
1358 gmx_runtime_t *runtime,
1359 gmx_bool bWriteStat)
1362 t_nrnb *nrnb_tot=NULL;
1365 double cycles[ewcNR];
1367 wallcycle_sum(cr,wcycle,cycles);
1369 if (cr->nnodes > 1) {
1373 MPI_Reduce(nrnb->n,nrnb_tot->n,eNRNB,MPI_DOUBLE,MPI_SUM,
1374 MASTERRANK(cr),cr->mpi_comm_mysim);
1380 if (SIMMASTER(cr)) {
1381 print_flop(fplog,nrnb_tot,&nbfs,&mflop);
1382 if (cr->nnodes > 1) {
1387 if ((cr->duty & DUTY_PP) && DOMAINDECOMP(cr)) {
1388 print_dd_statistics(cr,inputrec,fplog);
1400 snew(nrnb_all,cr->nnodes);
1401 nrnb_all[0] = *nrnb;
1402 for(s=1; s<cr->nnodes; s++)
1404 MPI_Recv(nrnb_all[s].n,eNRNB,MPI_DOUBLE,s,0,
1405 cr->mpi_comm_mysim,&stat);
1407 pr_load(fplog,cr,nrnb_all);
1412 MPI_Send(nrnb->n,eNRNB,MPI_DOUBLE,MASTERRANK(cr),0,
1413 cr->mpi_comm_mysim);
1418 if (SIMMASTER(cr)) {
1419 wallcycle_print(fplog,cr->nnodes,cr->npmenodes,runtime->realtime,
1422 if (EI_DYNAMICS(inputrec->eI)) {
1423 delta_t = inputrec->delta_t;
1429 print_perf(fplog,runtime->proctime,runtime->realtime,
1430 cr->nnodes-cr->npmenodes,
1431 runtime->nsteps_done,delta_t,nbfs,mflop);
1434 print_perf(stderr,runtime->proctime,runtime->realtime,
1435 cr->nnodes-cr->npmenodes,
1436 runtime->nsteps_done,delta_t,nbfs,mflop);
1440 runtime=inputrec->nsteps*inputrec->delta_t;
1442 if (cr->nnodes == 1)
1443 fprintf(stderr,"\n\n");
1444 print_perf(stderr,nodetime,realtime,runtime,&ntot,
1445 cr->nnodes-cr->npmenodes,FALSE);
1447 wallcycle_print(fplog,cr->nnodes,cr->npmenodes,realtime,wcycle,cycles);
1448 print_perf(fplog,nodetime,realtime,runtime,&ntot,cr->nnodes-cr->npmenodes,
1451 pr_load(fplog,cr,nrnb_all);
1458 void init_md(FILE *fplog,
1459 t_commrec *cr,t_inputrec *ir,const output_env_t oenv,
1460 double *t,double *t0,
1461 real *lambda,double *lam0,
1462 t_nrnb *nrnb,gmx_mtop_t *mtop,
1464 int nfile,const t_filenm fnm[],
1465 gmx_mdoutf_t **outf,t_mdebin **mdebin,
1466 tensor force_vir,tensor shake_vir,rvec mu_tot,
1467 gmx_bool *bSimAnn,t_vcm **vcm, t_state *state, unsigned long Flags)
1472 /* Initial values */
1473 *t = *t0 = ir->init_t;
1474 if (ir->efep != efepNO)
1476 *lam0 = ir->init_lambda;
1477 *lambda = *lam0 + ir->init_step*ir->delta_lambda;
1481 *lambda = *lam0 = 0.0;
1485 for(i=0;i<ir->opts.ngtc;i++)
1487 /* set bSimAnn if any group is being annealed */
1488 if(ir->opts.annealing[i]!=eannNO)
1495 update_annealing_target_temp(&(ir->opts),ir->init_t);
1500 *upd = init_update(fplog,ir);
1505 *vcm = init_vcm(fplog,&mtop->groups,ir);
1508 if (EI_DYNAMICS(ir->eI) && !(Flags & MD_APPENDFILES))
1510 if (ir->etc == etcBERENDSEN)
1512 please_cite(fplog,"Berendsen84a");
1514 if (ir->etc == etcVRESCALE)
1516 please_cite(fplog,"Bussi2007a");
1524 *outf = init_mdoutf(nfile,fnm,Flags,cr,ir,oenv);
1526 *mdebin = init_mdebin((Flags & MD_APPENDFILES) ? NULL : (*outf)->fp_ene,
1527 mtop,ir, (*outf)->fp_dhdl);
1530 /* Initiate variables */
1531 clear_mat(force_vir);
1532 clear_mat(shake_vir);