1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * GROwing Monsters And Cloning Shrimps
41 #include<catamount/dclock.h>
47 #ifdef HAVE_SYS_TIME_H
60 #include "chargegroup.h"
83 #include "pull_rotation.h"
84 #include "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
117 gettimeofday(&t,NULL);
119 seconds = (double) t.tv_sec + 1e-6*(double)t.tv_usec;
125 seconds = time(NULL);
132 #define difftime(end,start) ((double)(end)-(double)(start))
134 void print_time(FILE *out,gmx_runtime_t *runtime,gmx_large_int_t step,
135 t_inputrec *ir, t_commrec *cr)
138 char timebuf[STRLEN];
148 fprintf(out,"step %s",gmx_step_str(step,buf));
149 if ((step >= ir->nstlist))
151 if ((ir->nstlist == 0) || ((step % ir->nstlist) == 0))
153 /* We have done a full cycle let's update time_per_step */
154 runtime->last = gmx_gettime();
155 dt = difftime(runtime->last,runtime->real);
156 runtime->time_per_step = dt/(step - ir->init_step + 1);
158 dt = (ir->nsteps + ir->init_step - step)*runtime->time_per_step;
164 finish = (time_t) (runtime->last + dt);
165 gmx_ctime_r(&finish,timebuf,STRLEN);
166 sprintf(buf,"%s",timebuf);
167 buf[strlen(buf)-1]='\0';
168 fprintf(out,", will finish %s",buf);
171 fprintf(out,", remaining runtime: %5d s ",(int)dt);
175 fprintf(out," performance: %.1f ns/day ",
176 ir->delta_t/1000*24*60*60/runtime->time_per_step);
193 static double set_proctime(gmx_runtime_t *runtime)
199 prev = runtime->proc;
200 runtime->proc = dclock();
202 diff = runtime->proc - prev;
206 prev = runtime->proc;
207 runtime->proc = clock();
209 diff = (double)(runtime->proc - prev)/(double)CLOCKS_PER_SEC;
213 /* The counter has probably looped, ignore this data */
220 void runtime_start(gmx_runtime_t *runtime)
222 runtime->real = gmx_gettime();
224 set_proctime(runtime);
225 runtime->realtime = 0;
226 runtime->proctime = 0;
228 runtime->time_per_step = 0;
231 void runtime_end(gmx_runtime_t *runtime)
237 runtime->proctime += set_proctime(runtime);
238 runtime->realtime = now - runtime->real;
242 void runtime_upd_proc(gmx_runtime_t *runtime)
244 runtime->proctime += set_proctime(runtime);
247 void print_date_and_time(FILE *fplog,int nodeid,const char *title,
248 const gmx_runtime_t *runtime)
251 char timebuf[STRLEN];
252 char time_string[STRLEN];
259 tmptime = (time_t) runtime->real;
260 gmx_ctime_r(&tmptime,timebuf,STRLEN);
264 tmptime = (time_t) gmx_gettime();
265 gmx_ctime_r(&tmptime,timebuf,STRLEN);
267 for(i=0; timebuf[i]>=' '; i++)
269 time_string[i]=timebuf[i];
273 fprintf(fplog,"%s on node %d %s\n",title,nodeid,time_string);
277 static void sum_forces(int start,int end,rvec f[],rvec flr[])
282 pr_rvecs(debug,0,"fsr",f+start,end-start);
283 pr_rvecs(debug,0,"flr",flr+start,end-start);
285 for(i=start; (i<end); i++)
286 rvec_inc(f[i],flr[i]);
290 * calc_f_el calculates forces due to an electric field.
292 * force is kJ mol^-1 nm^-1 = e * kJ mol^-1 nm^-1 / e
294 * Et[] contains the parameters for the time dependent
295 * part of the field (not yet used).
296 * Ex[] contains the parameters for
297 * the spatial dependent part of the field. You can have cool periodic
298 * fields in principle, but only a constant field is supported
300 * The function should return the energy due to the electric field
301 * (if any) but for now returns 0.
304 * There can be problems with the virial.
305 * Since the field is not self-consistent this is unavoidable.
306 * For neutral molecules the virial is correct within this approximation.
307 * For neutral systems with many charged molecules the error is small.
308 * But for systems with a net charge or a few charged molecules
309 * the error can be significant when the field is high.
310 * Solution: implement a self-consitent electric field into PME.
312 static void calc_f_el(FILE *fp,int start,int homenr,
313 real charge[],rvec x[],rvec f[],
314 t_cosines Ex[],t_cosines Et[],double t)
320 for(m=0; (m<DIM); m++)
327 Ext[m] = cos(Et[m].a[0]*(t-t0))*exp(-sqr(t-t0)/(2.0*sqr(Et[m].a[2])));
331 Ext[m] = cos(Et[m].a[0]*t);
340 /* Convert the field strength from V/nm to MD-units */
341 Ext[m] *= Ex[m].a[0]*FIELDFAC;
342 for(i=start; (i<start+homenr); i++)
343 f[i][m] += charge[i]*Ext[m];
352 fprintf(fp,"%10g %10g %10g %10g #FIELD\n",t,
353 Ext[XX]/FIELDFAC,Ext[YY]/FIELDFAC,Ext[ZZ]/FIELDFAC);
357 static void calc_virial(FILE *fplog,int start,int homenr,rvec x[],rvec f[],
358 tensor vir_part,t_graph *graph,matrix box,
359 t_nrnb *nrnb,const t_forcerec *fr,int ePBC)
364 /* The short-range virial from surrounding boxes */
366 calc_vir(fplog,SHIFTS,fr->shift_vec,fr->fshift,vir_part,ePBC==epbcSCREW,box);
367 inc_nrnb(nrnb,eNR_VIRIAL,SHIFTS);
369 /* Calculate partial virial, for local atoms only, based on short range.
370 * Total virial is computed in global_stat, called from do_md
372 f_calc_vir(fplog,start,start+homenr,x,f,vir_part,graph,box);
373 inc_nrnb(nrnb,eNR_VIRIAL,homenr);
375 /* Add position restraint contribution */
376 for(i=0; i<DIM; i++) {
377 vir_part[i][i] += fr->vir_diag_posres[i];
380 /* Add wall contribution */
381 for(i=0; i<DIM; i++) {
382 vir_part[i][ZZ] += fr->vir_wall_z[i];
386 pr_rvecs(debug,0,"vir_part",vir_part,DIM);
389 static void print_large_forces(FILE *fp,t_mdatoms *md,t_commrec *cr,
390 gmx_large_int_t step,real pforce,rvec *x,rvec *f)
394 char buf[STEPSTRSIZE];
397 for(i=md->start; i<md->start+md->homenr; i++) {
399 /* We also catch NAN, if the compiler does not optimize this away. */
400 if (fn2 >= pf2 || fn2 != fn2) {
401 fprintf(fp,"step %s atom %6d x %8.3f %8.3f %8.3f force %12.5e\n",
402 gmx_step_str(step,buf),
403 ddglatnr(cr->dd,i),x[i][XX],x[i][YY],x[i][ZZ],sqrt(fn2));
408 void do_force(FILE *fplog,t_commrec *cr,
409 t_inputrec *inputrec,
410 gmx_large_int_t step,t_nrnb *nrnb,gmx_wallcycle_t wcycle,
413 gmx_groups_t *groups,
414 matrix box,rvec x[],history_t *hist,
418 gmx_enerdata_t *enerd,t_fcdata *fcd,
419 real lambda,t_graph *graph,
420 t_forcerec *fr,gmx_vsite_t *vsite,rvec mu_tot,
421 double t,FILE *field,gmx_edsam_t ed,
428 gmx_bool bSepDVDL,bStateChanged,bNS,bFillGrid,bCalcCGCM,bBS;
429 gmx_bool bDoLongRange,bDoForces,bSepLRF;
433 float cycles_ppdpme,cycles_pme,cycles_seppme,cycles_force;
435 start = mdatoms->start;
436 homenr = mdatoms->homenr;
438 bSepDVDL = (fr->bSepDVDL && do_per_step(step,inputrec->nstlog));
440 clear_mat(vir_force);
444 pd_cg_range(cr,&cg0,&cg1);
449 if (DOMAINDECOMP(cr))
451 cg1 = cr->dd->ncg_tot;
463 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
464 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll==FALSE);
465 bFillGrid = (bNS && bStateChanged);
466 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
467 bDoLongRange = (fr->bTwinRange && bNS && (flags & GMX_FORCE_DOLR));
468 bDoForces = (flags & GMX_FORCE_FORCES);
469 bSepLRF = (bDoLongRange && bDoForces && (flags & GMX_FORCE_SEPLRF));
473 update_forcerec(fplog,fr,box);
475 /* Calculate total (local) dipole moment in a temporary common array.
476 * This makes it possible to sum them over nodes faster.
478 calc_mu(start,homenr,
479 x,mdatoms->chargeA,mdatoms->chargeB,mdatoms->nChargePerturbed,
483 if (fr->ePBC != epbcNONE) {
484 /* Compute shift vectors every step,
485 * because of pressure coupling or box deformation!
487 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
488 calc_shifts(box,fr->shift_vec);
491 put_charge_groups_in_box(fplog,cg0,cg1,fr->ePBC,box,
492 &(top->cgs),x,fr->cg_cm);
493 inc_nrnb(nrnb,eNR_CGCM,homenr);
494 inc_nrnb(nrnb,eNR_RESETX,cg1-cg0);
496 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph) {
497 unshift_self(graph,box,x);
500 else if (bCalcCGCM) {
501 calc_cgcm(fplog,cg0,cg1,&(top->cgs),x,fr->cg_cm);
502 inc_nrnb(nrnb,eNR_CGCM,homenr);
507 move_cgcm(fplog,cr,fr->cg_cm);
510 pr_rvecs(debug,0,"cgcm",fr->cg_cm,top->cgs.nr);
514 if (!(cr->duty & DUTY_PME)) {
515 /* Send particle coordinates to the pme nodes.
516 * Since this is only implemented for domain decomposition
517 * and domain decomposition does not use the graph,
518 * we do not need to worry about shifting.
521 wallcycle_start(wcycle,ewcPP_PMESENDX);
522 GMX_MPE_LOG(ev_send_coordinates_start);
524 bBS = (inputrec->nwall == 2);
527 svmul(inputrec->wall_ewald_zfac,boxs[ZZ],boxs[ZZ]);
530 gmx_pme_send_x(cr,bBS ? boxs : box,x,
531 mdatoms->nChargePerturbed,lambda,
532 ( flags & GMX_FORCE_VIRIAL),step);
534 GMX_MPE_LOG(ev_send_coordinates_finish);
535 wallcycle_stop(wcycle,ewcPP_PMESENDX);
539 /* Communicate coordinates and sum dipole if necessary */
542 wallcycle_start(wcycle,ewcMOVEX);
543 if (DOMAINDECOMP(cr))
545 dd_move_x(cr->dd,box,x);
549 move_x(fplog,cr,GMX_LEFT,GMX_RIGHT,x,nrnb);
551 /* When we don't need the total dipole we sum it in global_stat */
552 if (bStateChanged && NEED_MUTOT(*inputrec))
554 gmx_sumd(2*DIM,mu,cr);
556 wallcycle_stop(wcycle,ewcMOVEX);
564 fr->mu_tot[i][j] = mu[i*DIM + j];
568 if (fr->efep == efepNO)
570 copy_rvec(fr->mu_tot[0],mu_tot);
577 (1.0 - lambda)*fr->mu_tot[0][j] + lambda*fr->mu_tot[1][j];
582 reset_enerdata(&(inputrec->opts),fr,bNS,enerd,MASTER(cr));
583 clear_rvecs(SHIFTS,fr->fshift);
587 wallcycle_start(wcycle,ewcNS);
589 if (graph && bStateChanged)
591 /* Calculate intramolecular shift vectors to make molecules whole */
592 mk_mshift(fplog,graph,fr->ePBC,box,x);
595 /* Reset long range forces if necessary */
598 /* Reset the (long-range) forces if necessary */
599 clear_rvecs(fr->natoms_force_constr,bSepLRF ? fr->f_twin : f);
602 /* Do the actual neighbour searching and if twin range electrostatics
603 * also do the calculation of long range forces and energies.
607 groups,&(inputrec->opts),top,mdatoms,
608 cr,nrnb,lambda,&dvdl,&enerd->grpp,bFillGrid,
609 bDoLongRange,bDoForces,bSepLRF ? fr->f_twin : f);
612 fprintf(fplog,sepdvdlformat,"LR non-bonded",0.0,dvdl);
614 enerd->dvdl_lin += dvdl;
616 wallcycle_stop(wcycle,ewcNS);
619 if (inputrec->implicit_solvent && bNS)
621 make_gb_nblist(cr,inputrec->gb_algorithm,inputrec->rlist,
622 x,box,fr,&top->idef,graph,fr->born);
625 if (DOMAINDECOMP(cr))
627 if (!(cr->duty & DUTY_PME))
629 wallcycle_start(wcycle,ewcPPDURINGPME);
630 dd_force_flop_start(cr->dd,nrnb);
636 /* Enforced rotation has its own cycle counter that starts after the collective
637 * coordinates have been communicated. It is added to ddCyclF to allow
638 * for proper load-balancing */
639 wallcycle_start(wcycle,ewcROT);
640 do_rotation(cr,inputrec,box,x,t,step,wcycle,bNS);
641 wallcycle_stop(wcycle,ewcROT);
644 /* Start the force cycle counter.
645 * This counter is stopped in do_forcelow_level.
646 * No parallel communication should occur while this counter is running,
647 * since that will interfere with the dynamic load balancing.
649 wallcycle_start(wcycle,ewcFORCE);
653 /* Reset forces for which the virial is calculated separately:
654 * PME/Ewald forces if necessary */
657 if (flags & GMX_FORCE_VIRIAL)
659 fr->f_novirsum = fr->f_novirsum_alloc;
660 GMX_BARRIER(cr->mpi_comm_mygroup);
663 clear_rvecs(fr->f_novirsum_n,fr->f_novirsum);
667 clear_rvecs(homenr,fr->f_novirsum+start);
669 GMX_BARRIER(cr->mpi_comm_mygroup);
673 /* We are not calculating the pressure so we do not need
674 * a separate array for forces that do not contribute
683 /* Add the long range forces to the short range forces */
684 for(i=0; i<fr->natoms_force_constr; i++)
686 copy_rvec(fr->f_twin[i],f[i]);
689 else if (!(fr->bTwinRange && bNS))
691 /* Clear the short-range forces */
692 clear_rvecs(fr->natoms_force_constr,f);
695 clear_rvec(fr->vir_diag_posres);
697 GMX_BARRIER(cr->mpi_comm_mygroup);
699 if (inputrec->ePull == epullCONSTRAINT)
701 clear_pull_forces(inputrec->pull);
704 /* update QMMMrec, if necessary */
707 update_QMMMrec(cr,fr,x,mdatoms,box,top);
710 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
712 /* Position restraints always require full pbc */
713 set_pbc(&pbc,inputrec->ePBC,box);
714 v = posres(top->idef.il[F_POSRES].nr,top->idef.il[F_POSRES].iatoms,
715 top->idef.iparams_posres,
716 (const rvec*)x,fr->f_novirsum,fr->vir_diag_posres,
717 inputrec->ePBC==epbcNONE ? NULL : &pbc,lambda,&dvdl,
718 fr->rc_scaling,fr->ePBC,fr->posres_com,fr->posres_comB);
721 fprintf(fplog,sepdvdlformat,
722 interaction_function[F_POSRES].longname,v,dvdl);
724 enerd->term[F_POSRES] += v;
725 /* This linear lambda dependence assumption is only correct
726 * when only k depends on lambda,
727 * not when the reference position depends on lambda.
728 * grompp checks for this.
730 enerd->dvdl_lin += dvdl;
731 inc_nrnb(nrnb,eNR_POSRES,top->idef.il[F_POSRES].nr/2);
734 /* Compute the bonded and non-bonded energies and optionally forces */
735 do_force_lowlevel(fplog,step,fr,inputrec,&(top->idef),
736 cr,nrnb,wcycle,mdatoms,&(inputrec->opts),
737 x,hist,f,enerd,fcd,mtop,top,fr->born,
738 &(top->atomtypes),bBornRadii,box,
739 lambda,graph,&(top->excls),fr->mu_tot,
742 cycles_force = wallcycle_stop(wcycle,ewcFORCE);
743 GMX_BARRIER(cr->mpi_comm_mygroup);
747 do_flood(fplog,cr,x,f,ed,box,step);
750 if (DOMAINDECOMP(cr))
752 dd_force_flop_stop(cr->dd,nrnb);
755 dd_cycles_add(cr->dd,cycles_force-cycles_pme,ddCyclF);
761 if (IR_ELEC_FIELD(*inputrec))
763 /* Compute forces due to electric field */
764 calc_f_el(MASTER(cr) ? field : NULL,
765 start,homenr,mdatoms->chargeA,x,fr->f_novirsum,
766 inputrec->ex,inputrec->et,t);
769 /* Communicate the forces */
772 wallcycle_start(wcycle,ewcMOVEF);
773 if (DOMAINDECOMP(cr))
775 dd_move_f(cr->dd,f,fr->fshift);
776 /* Do we need to communicate the separate force array
777 * for terms that do not contribute to the single sum virial?
778 * Position restraints and electric fields do not introduce
779 * inter-cg forces, only full electrostatics methods do.
780 * When we do not calculate the virial, fr->f_novirsum = f,
781 * so we have already communicated these forces.
783 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
784 (flags & GMX_FORCE_VIRIAL))
786 dd_move_f(cr->dd,fr->f_novirsum,NULL);
790 /* We should not update the shift forces here,
791 * since f_twin is already included in f.
793 dd_move_f(cr->dd,fr->f_twin,NULL);
798 pd_move_f(cr,f,nrnb);
801 pd_move_f(cr,fr->f_twin,nrnb);
804 wallcycle_stop(wcycle,ewcMOVEF);
807 /* If we have NoVirSum forces, but we do not calculate the virial,
808 * we sum fr->f_novirum=f later.
810 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
812 wallcycle_start(wcycle,ewcVSITESPREAD);
813 spread_vsite_f(fplog,vsite,x,f,fr->fshift,nrnb,
814 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
815 wallcycle_stop(wcycle,ewcVSITESPREAD);
819 wallcycle_start(wcycle,ewcVSITESPREAD);
820 spread_vsite_f(fplog,vsite,x,fr->f_twin,NULL,
822 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
823 wallcycle_stop(wcycle,ewcVSITESPREAD);
827 if (flags & GMX_FORCE_VIRIAL)
829 /* Calculation of the virial must be done after vsites! */
830 calc_virial(fplog,mdatoms->start,mdatoms->homenr,x,f,
831 vir_force,graph,box,nrnb,fr,inputrec->ePBC);
835 enerd->term[F_COM_PULL] = 0;
836 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
838 /* Calculate the center of mass forces, this requires communication,
839 * which is why pull_potential is called close to other communication.
840 * The virial contribution is calculated directly,
841 * which is why we call pull_potential after calc_virial.
843 set_pbc(&pbc,inputrec->ePBC,box);
845 enerd->term[F_COM_PULL] +=
846 pull_potential(inputrec->ePull,inputrec->pull,mdatoms,&pbc,
847 cr,t,lambda,x,f,vir_force,&dvdl);
850 fprintf(fplog,sepdvdlformat,"Com pull",enerd->term[F_COM_PULL],dvdl);
852 enerd->dvdl_lin += dvdl;
855 /* Add the forces from enforced rotation potentials (if any) */
858 wallcycle_start(wcycle,ewcROTadd);
859 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr,step,t);
860 wallcycle_stop(wcycle,ewcROTadd);
863 if (PAR(cr) && !(cr->duty & DUTY_PME))
865 cycles_ppdpme = wallcycle_stop(wcycle,ewcPPDURINGPME);
866 dd_cycles_add(cr->dd,cycles_ppdpme,ddCyclPPduringPME);
868 /* In case of node-splitting, the PP nodes receive the long-range
869 * forces, virial and energy from the PME nodes here.
871 wallcycle_start(wcycle,ewcPP_PMEWAITRECVF);
873 gmx_pme_receive_f(cr,fr->f_novirsum,fr->vir_el_recip,&e,&dvdl,
877 fprintf(fplog,sepdvdlformat,"PME mesh",e,dvdl);
879 enerd->term[F_COUL_RECIP] += e;
880 enerd->dvdl_lin += dvdl;
883 dd_cycles_add(cr->dd,cycles_seppme,ddCyclPME);
885 wallcycle_stop(wcycle,ewcPP_PMEWAITRECVF);
888 if (bDoForces && fr->bF_NoVirSum)
892 /* Spread the mesh force on virtual sites to the other particles...
893 * This is parallellized. MPI communication is performed
894 * if the constructing atoms aren't local.
896 wallcycle_start(wcycle,ewcVSITESPREAD);
897 spread_vsite_f(fplog,vsite,x,fr->f_novirsum,NULL,nrnb,
898 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
899 wallcycle_stop(wcycle,ewcVSITESPREAD);
901 if (flags & GMX_FORCE_VIRIAL)
903 /* Now add the forces, this is local */
906 sum_forces(0,fr->f_novirsum_n,f,fr->f_novirsum);
910 sum_forces(start,start+homenr,f,fr->f_novirsum);
912 if (EEL_FULL(fr->eeltype))
914 /* Add the mesh contribution to the virial */
915 m_add(vir_force,fr->vir_el_recip,vir_force);
919 pr_rvecs(debug,0,"vir_force",vir_force,DIM);
924 /* Sum the potential energy terms from group contributions */
925 sum_epot(&(inputrec->opts),enerd);
927 if (fr->print_force >= 0 && bDoForces)
929 print_large_forces(stderr,mdatoms,cr,step,fr->print_force,x,f);
933 void do_constrain_first(FILE *fplog,gmx_constr_t constr,
934 t_inputrec *ir,t_mdatoms *md,
935 t_state *state,rvec *f,
936 t_graph *graph,t_commrec *cr,t_nrnb *nrnb,
937 t_forcerec *fr, gmx_localtop_t *top, tensor shake_vir)
940 gmx_large_int_t step;
941 double mass,tmass,vcm[4];
946 snew(savex,state->natoms);
949 end = md->homenr + start;
952 fprintf(debug,"vcm: start=%d, homenr=%d, end=%d\n",
953 start,md->homenr,end);
954 /* Do a first constrain to reset particles... */
955 step = ir->init_step;
958 char buf[STEPSTRSIZE];
959 fprintf(fplog,"\nConstraining the starting coordinates (step %s)\n",
960 gmx_step_str(step,buf));
964 /* constrain the current position */
965 constrain(NULL,TRUE,FALSE,constr,&(top->idef),
966 ir,NULL,cr,step,0,md,
967 state->x,state->x,NULL,
968 state->box,state->lambda,&dvdlambda,
969 NULL,NULL,nrnb,econqCoord,ir->epc==epcMTTK,state->veta,state->veta);
972 /* constrain the inital velocity, and save it */
973 /* also may be useful if we need the ekin from the halfstep for velocity verlet */
974 /* might not yet treat veta correctly */
975 constrain(NULL,TRUE,FALSE,constr,&(top->idef),
976 ir,NULL,cr,step,0,md,
977 state->x,state->v,state->v,
978 state->box,state->lambda,&dvdlambda,
979 NULL,NULL,nrnb,econqVeloc,ir->epc==epcMTTK,state->veta,state->veta);
981 /* constrain the inital velocities at t-dt/2 */
982 if (EI_STATE_VELOCITY(ir->eI) && ir->eI!=eiVV)
984 for(i=start; (i<end); i++)
986 for(m=0; (m<DIM); m++)
988 /* Reverse the velocity */
989 state->v[i][m] = -state->v[i][m];
990 /* Store the position at t-dt in buf */
991 savex[i][m] = state->x[i][m] + dt*state->v[i][m];
994 /* Shake the positions at t=-dt with the positions at t=0
995 * as reference coordinates.
999 char buf[STEPSTRSIZE];
1000 fprintf(fplog,"\nConstraining the coordinates at t0-dt (step %s)\n",
1001 gmx_step_str(step,buf));
1004 constrain(NULL,TRUE,FALSE,constr,&(top->idef),
1005 ir,NULL,cr,step,-1,md,
1006 state->x,savex,NULL,
1007 state->box,state->lambda,&dvdlambda,
1008 state->v,NULL,nrnb,econqCoord,ir->epc==epcMTTK,state->veta,state->veta);
1010 for(i=start; i<end; i++) {
1011 for(m=0; m<DIM; m++) {
1012 /* Re-reverse the velocities */
1013 state->v[i][m] = -state->v[i][m];
1018 for(m=0; (m<4); m++)
1020 for(i=start; i<end; i++) {
1021 mass = md->massT[i];
1022 for(m=0; m<DIM; m++) {
1023 vcm[m] += state->v[i][m]*mass;
1028 if (ir->nstcomm != 0 || debug) {
1029 /* Compute the global sum of vcm */
1031 fprintf(debug,"vcm: %8.3f %8.3f %8.3f,"
1032 " total mass = %12.5e\n",vcm[XX],vcm[YY],vcm[ZZ],vcm[3]);
1036 for(m=0; (m<DIM); m++)
1039 fprintf(debug,"vcm: %8.3f %8.3f %8.3f,"
1040 " total mass = %12.5e\n",vcm[XX],vcm[YY],vcm[ZZ],tmass);
1041 if (ir->nstcomm != 0) {
1042 /* Now we have the velocity of center of mass, let's remove it */
1043 for(i=start; (i<end); i++) {
1044 for(m=0; (m<DIM); m++)
1045 state->v[i][m] -= vcm[m];
1053 void calc_enervirdiff(FILE *fplog,int eDispCorr,t_forcerec *fr)
1055 double eners[2],virs[2],enersum,virsum,y0,f,g,h;
1056 double r0,r1,r,rc3,rc9,ea,eb,ec,pa,pb,pc,pd;
1057 double invscale,invscale2,invscale3;
1058 int ri0,ri1,ri,i,offstart,offset;
1061 fr->enershiftsix = 0;
1062 fr->enershifttwelve = 0;
1063 fr->enerdiffsix = 0;
1064 fr->enerdifftwelve = 0;
1066 fr->virdifftwelve = 0;
1068 if (eDispCorr != edispcNO) {
1069 for(i=0; i<2; i++) {
1073 if ((fr->vdwtype == evdwSWITCH) || (fr->vdwtype == evdwSHIFT)) {
1074 if (fr->rvdw_switch == 0)
1076 "With dispersion correction rvdw-switch can not be zero "
1077 "for vdw-type = %s",evdw_names[fr->vdwtype]);
1079 scale = fr->nblists[0].tab.scale;
1080 vdwtab = fr->nblists[0].vdwtab;
1082 /* Round the cut-offs to exact table values for precision */
1083 ri0 = floor(fr->rvdw_switch*scale);
1084 ri1 = ceil(fr->rvdw*scale);
1090 if (fr->vdwtype == evdwSHIFT) {
1091 /* Determine the constant energy shift below rvdw_switch */
1092 fr->enershiftsix = (real)(-1.0/(rc3*rc3)) - vdwtab[8*ri0];
1093 fr->enershifttwelve = (real)( 1.0/(rc9*rc3)) - vdwtab[8*ri0 + 4];
1095 /* Add the constant part from 0 to rvdw_switch.
1096 * This integration from 0 to rvdw_switch overcounts the number
1097 * of interactions by 1, as it also counts the self interaction.
1098 * We will correct for this later.
1100 eners[0] += 4.0*M_PI*fr->enershiftsix*rc3/3.0;
1101 eners[1] += 4.0*M_PI*fr->enershifttwelve*rc3/3.0;
1103 invscale = 1.0/(scale);
1104 invscale2 = invscale*invscale;
1105 invscale3 = invscale*invscale2;
1107 /* following summation derived from cubic spline definition,
1108 Numerical Recipies in C, second edition, p. 113-116. Exact
1109 for the cubic spline. We first calculate the negative of
1110 the energy from rvdw to rvdw_switch, assuming that g(r)=1,
1111 and then add the more standard, abrupt cutoff correction to
1112 that result, yielding the long-range correction for a
1113 switched function. We perform both the pressure and energy
1114 loops at the same time for simplicity, as the computational
1118 enersum = 0.0; virsum = 0.0;
1123 for (ri=ri0; ri<ri1; ri++) {
1126 eb = 2.0*invscale2*r;
1130 pb = 3.0*invscale2*r;
1131 pc = 3.0*invscale*r*r;
1134 /* this "8" is from the packing in the vdwtab array - perhaps
1135 should be #define'ed? */
1136 offset = 8*ri + offstart;
1137 y0 = vdwtab[offset];
1138 f = vdwtab[offset+1];
1139 g = vdwtab[offset+2];
1140 h = vdwtab[offset+3];
1142 enersum += y0*(ea/3 + eb/2 + ec) + f*(ea/4 + eb/3 + ec/2)+
1143 g*(ea/5 + eb/4 + ec/3) + h*(ea/6 + eb/5 + ec/4);
1144 virsum += f*(pa/4 + pb/3 + pc/2 + pd) +
1145 2*g*(pa/5 + pb/4 + pc/3 + pd/2) + 3*h*(pa/6 + pb/5 + pc/4 + pd/3);
1148 enersum *= 4.0*M_PI;
1150 eners[i] -= enersum;
1154 /* now add the correction for rvdw_switch to infinity */
1155 eners[0] += -4.0*M_PI/(3.0*rc3);
1156 eners[1] += 4.0*M_PI/(9.0*rc9);
1157 virs[0] += 8.0*M_PI/rc3;
1158 virs[1] += -16.0*M_PI/(3.0*rc9);
1160 else if ((fr->vdwtype == evdwCUT) || (fr->vdwtype == evdwUSER)) {
1161 if (fr->vdwtype == evdwUSER && fplog)
1163 "WARNING: using dispersion correction with user tables\n");
1164 rc3 = fr->rvdw*fr->rvdw*fr->rvdw;
1166 eners[0] += -4.0*M_PI/(3.0*rc3);
1167 eners[1] += 4.0*M_PI/(9.0*rc9);
1168 virs[0] += 8.0*M_PI/rc3;
1169 virs[1] += -16.0*M_PI/(3.0*rc9);
1172 "Dispersion correction is not implemented for vdw-type = %s",
1173 evdw_names[fr->vdwtype]);
1175 fr->enerdiffsix = eners[0];
1176 fr->enerdifftwelve = eners[1];
1177 /* The 0.5 is due to the Gromacs definition of the virial */
1178 fr->virdiffsix = 0.5*virs[0];
1179 fr->virdifftwelve = 0.5*virs[1];
1183 void calc_dispcorr(FILE *fplog,t_inputrec *ir,t_forcerec *fr,
1184 gmx_large_int_t step,int natoms,
1185 matrix box,real lambda,tensor pres,tensor virial,
1186 real *prescorr, real *enercorr, real *dvdlcorr)
1188 gmx_bool bCorrAll,bCorrPres;
1189 real dvdlambda,invvol,dens,ninter,avcsix,avctwelve,enerdiff,svir=0,spres=0;
1199 if (ir->eDispCorr != edispcNO) {
1200 bCorrAll = (ir->eDispCorr == edispcAllEner ||
1201 ir->eDispCorr == edispcAllEnerPres);
1202 bCorrPres = (ir->eDispCorr == edispcEnerPres ||
1203 ir->eDispCorr == edispcAllEnerPres);
1205 invvol = 1/det(box);
1208 /* Only correct for the interactions with the inserted molecule */
1209 dens = (natoms - fr->n_tpi)*invvol;
1214 dens = natoms*invvol;
1215 ninter = 0.5*natoms;
1218 if (ir->efep == efepNO)
1220 avcsix = fr->avcsix[0];
1221 avctwelve = fr->avctwelve[0];
1225 avcsix = (1 - lambda)*fr->avcsix[0] + lambda*fr->avcsix[1];
1226 avctwelve = (1 - lambda)*fr->avctwelve[0] + lambda*fr->avctwelve[1];
1229 enerdiff = ninter*(dens*fr->enerdiffsix - fr->enershiftsix);
1230 *enercorr += avcsix*enerdiff;
1232 if (ir->efep != efepNO)
1234 dvdlambda += (fr->avcsix[1] - fr->avcsix[0])*enerdiff;
1238 enerdiff = ninter*(dens*fr->enerdifftwelve - fr->enershifttwelve);
1239 *enercorr += avctwelve*enerdiff;
1240 if (fr->efep != efepNO)
1242 dvdlambda += (fr->avctwelve[1] - fr->avctwelve[0])*enerdiff;
1248 svir = ninter*dens*avcsix*fr->virdiffsix/3.0;
1249 if (ir->eDispCorr == edispcAllEnerPres)
1251 svir += ninter*dens*avctwelve*fr->virdifftwelve/3.0;
1253 /* The factor 2 is because of the Gromacs virial definition */
1254 spres = -2.0*invvol*svir*PRESFAC;
1256 for(m=0; m<DIM; m++) {
1257 virial[m][m] += svir;
1258 pres[m][m] += spres;
1263 /* Can't currently control when it prints, for now, just print when degugging */
1267 fprintf(debug,"Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
1273 "Long Range LJ corr.: Epot %10g, Pres: %10g, Vir: %10g\n",
1274 *enercorr,spres,svir);
1278 fprintf(debug,"Long Range LJ corr.: Epot %10g\n",*enercorr);
1282 if (fr->bSepDVDL && do_per_step(step,ir->nstlog))
1284 fprintf(fplog,sepdvdlformat,"Dispersion correction",
1285 *enercorr,dvdlambda);
1287 if (fr->efep != efepNO)
1289 *dvdlcorr += dvdlambda;
1294 void do_pbc_first(FILE *fplog,matrix box,t_forcerec *fr,
1295 t_graph *graph,rvec x[])
1298 fprintf(fplog,"Removing pbc first time\n");
1299 calc_shifts(box,fr->shift_vec);
1301 mk_mshift(fplog,graph,fr->ePBC,box,x);
1303 p_graph(debug,"do_pbc_first 1",graph);
1304 shift_self(graph,box,x);
1305 /* By doing an extra mk_mshift the molecules that are broken
1306 * because they were e.g. imported from another software
1307 * will be made whole again. Such are the healing powers
1310 mk_mshift(fplog,graph,fr->ePBC,box,x);
1312 p_graph(debug,"do_pbc_first 2",graph);
1315 fprintf(fplog,"Done rmpbc\n");
1318 static void low_do_pbc_mtop(FILE *fplog,int ePBC,matrix box,
1319 gmx_mtop_t *mtop,rvec x[],
1324 gmx_molblock_t *molb;
1326 if (bFirst && fplog)
1327 fprintf(fplog,"Removing pbc first time\n");
1331 for(mb=0; mb<mtop->nmolblock; mb++) {
1332 molb = &mtop->molblock[mb];
1333 if (molb->natoms_mol == 1 ||
1334 (!bFirst && mtop->moltype[molb->type].cgs.nr == 1)) {
1335 /* Just one atom or charge group in the molecule, no PBC required */
1336 as += molb->nmol*molb->natoms_mol;
1338 /* Pass NULL iso fplog to avoid graph prints for each molecule type */
1339 mk_graph_ilist(NULL,mtop->moltype[molb->type].ilist,
1340 0,molb->natoms_mol,FALSE,FALSE,graph);
1342 for(mol=0; mol<molb->nmol; mol++) {
1343 mk_mshift(fplog,graph,ePBC,box,x+as);
1345 shift_self(graph,box,x+as);
1346 /* The molecule is whole now.
1347 * We don't need the second mk_mshift call as in do_pbc_first,
1348 * since we no longer need this graph.
1351 as += molb->natoms_mol;
1359 void do_pbc_first_mtop(FILE *fplog,int ePBC,matrix box,
1360 gmx_mtop_t *mtop,rvec x[])
1362 low_do_pbc_mtop(fplog,ePBC,box,mtop,x,TRUE);
1365 void do_pbc_mtop(FILE *fplog,int ePBC,matrix box,
1366 gmx_mtop_t *mtop,rvec x[])
1368 low_do_pbc_mtop(fplog,ePBC,box,mtop,x,FALSE);
1371 void finish_run(FILE *fplog,t_commrec *cr,const char *confout,
1372 t_inputrec *inputrec,
1373 t_nrnb nrnb[],gmx_wallcycle_t wcycle,
1374 gmx_runtime_t *runtime,
1375 gmx_bool bWriteStat)
1378 t_nrnb *nrnb_tot=NULL;
1381 double cycles[ewcNR];
1383 wallcycle_sum(cr,wcycle,cycles);
1385 if (cr->nnodes > 1) {
1389 MPI_Reduce(nrnb->n,nrnb_tot->n,eNRNB,MPI_DOUBLE,MPI_SUM,
1390 MASTERRANK(cr),cr->mpi_comm_mysim);
1396 if (SIMMASTER(cr)) {
1397 print_flop(fplog,nrnb_tot,&nbfs,&mflop);
1398 if (cr->nnodes > 1) {
1403 if ((cr->duty & DUTY_PP) && DOMAINDECOMP(cr)) {
1404 print_dd_statistics(cr,inputrec,fplog);
1416 snew(nrnb_all,cr->nnodes);
1417 nrnb_all[0] = *nrnb;
1418 for(s=1; s<cr->nnodes; s++)
1420 MPI_Recv(nrnb_all[s].n,eNRNB,MPI_DOUBLE,s,0,
1421 cr->mpi_comm_mysim,&stat);
1423 pr_load(fplog,cr,nrnb_all);
1428 MPI_Send(nrnb->n,eNRNB,MPI_DOUBLE,MASTERRANK(cr),0,
1429 cr->mpi_comm_mysim);
1434 if (SIMMASTER(cr)) {
1435 wallcycle_print(fplog,cr->nnodes,cr->npmenodes,runtime->realtime,
1438 if (EI_DYNAMICS(inputrec->eI)) {
1439 delta_t = inputrec->delta_t;
1445 print_perf(fplog,runtime->proctime,runtime->realtime,
1446 cr->nnodes-cr->npmenodes,
1447 runtime->nsteps_done,delta_t,nbfs,mflop);
1450 print_perf(stderr,runtime->proctime,runtime->realtime,
1451 cr->nnodes-cr->npmenodes,
1452 runtime->nsteps_done,delta_t,nbfs,mflop);
1456 runtime=inputrec->nsteps*inputrec->delta_t;
1458 if (cr->nnodes == 1)
1459 fprintf(stderr,"\n\n");
1460 print_perf(stderr,nodetime,realtime,runtime,&ntot,
1461 cr->nnodes-cr->npmenodes,FALSE);
1463 wallcycle_print(fplog,cr->nnodes,cr->npmenodes,realtime,wcycle,cycles);
1464 print_perf(fplog,nodetime,realtime,runtime,&ntot,cr->nnodes-cr->npmenodes,
1467 pr_load(fplog,cr,nrnb_all);
1474 void init_md(FILE *fplog,
1475 t_commrec *cr,t_inputrec *ir,const output_env_t oenv,
1476 double *t,double *t0,
1477 real *lambda,double *lam0,
1478 t_nrnb *nrnb,gmx_mtop_t *mtop,
1480 int nfile,const t_filenm fnm[],
1481 gmx_mdoutf_t **outf,t_mdebin **mdebin,
1482 tensor force_vir,tensor shake_vir,rvec mu_tot,
1483 gmx_bool *bSimAnn,t_vcm **vcm, t_state *state, unsigned long Flags)
1488 /* Initial values */
1489 *t = *t0 = ir->init_t;
1490 if (ir->efep != efepNO)
1492 *lam0 = ir->init_lambda;
1493 *lambda = *lam0 + ir->init_step*ir->delta_lambda;
1497 *lambda = *lam0 = 0.0;
1501 for(i=0;i<ir->opts.ngtc;i++)
1503 /* set bSimAnn if any group is being annealed */
1504 if(ir->opts.annealing[i]!=eannNO)
1511 update_annealing_target_temp(&(ir->opts),ir->init_t);
1516 *upd = init_update(fplog,ir);
1521 *vcm = init_vcm(fplog,&mtop->groups,ir);
1524 if (EI_DYNAMICS(ir->eI) && !(Flags & MD_APPENDFILES))
1526 if (ir->etc == etcBERENDSEN)
1528 please_cite(fplog,"Berendsen84a");
1530 if (ir->etc == etcVRESCALE)
1532 please_cite(fplog,"Bussi2007a");
1540 *outf = init_mdoutf(nfile,fnm,Flags,cr,ir,oenv);
1542 *mdebin = init_mdebin((Flags & MD_APPENDFILES) ? NULL : (*outf)->fp_ene,
1543 mtop,ir, (*outf)->fp_dhdl);
1546 /* Initiate variables */
1547 clear_mat(force_vir);
1548 clear_mat(shake_vir);