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"
82 #include "pull_rotation.h"
83 #include "gmx_random.h"
84 #include "mpelogging.h"
87 #include "gmx_wallcycle.h"
101 typedef struct gmx_timeprint {
106 /* Portable version of ctime_r implemented in src/gmxlib/string2.c, but we do not want it declared in public installed headers */
108 gmx_ctime_r(const time_t *clock,char *buf, int n);
114 #ifdef HAVE_GETTIMEOFDAY
118 gettimeofday(&t,NULL);
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];
143 #ifndef GMX_THREAD_MPI
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);
180 #ifndef GMX_THREAD_MPI
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;
431 gmx_bool bDoAdressWF;
433 real e,v,dvdlambda[efptNR];
434 real dvdl_dum,lambda_dum;
436 float cycles_ppdpme,cycles_pme,cycles_seppme,cycles_force;
438 start = mdatoms->start;
439 homenr = mdatoms->homenr;
441 bSepDVDL = (fr->bSepDVDL && do_per_step(step,inputrec->nstlog));
443 clear_mat(vir_force);
447 pd_cg_range(cr,&cg0,&cg1);
452 if (DOMAINDECOMP(cr))
454 cg1 = cr->dd->ncg_tot;
466 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
467 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll==FALSE);
468 bFillGrid = (bNS && bStateChanged);
469 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
470 bDoLongRange = (fr->bTwinRange && bNS && (flags & GMX_FORCE_DOLR));
471 bDoForces = (flags & GMX_FORCE_FORCES);
472 bSepLRF = (bDoLongRange && bDoForces && (flags & GMX_FORCE_SEPLRF));
473 /* should probably move this to the forcerec since it doesn't change */
474 bDoAdressWF = ((fr->adress_type!=eAdressOff));
478 update_forcerec(fplog,fr,box);
480 /* Calculate total (local) dipole moment in a temporary common array.
481 * This makes it possible to sum them over nodes faster.
483 calc_mu(start,homenr,
484 x,mdatoms->chargeA,mdatoms->chargeB,mdatoms->nChargePerturbed,
488 if (fr->ePBC != epbcNONE) {
489 /* Compute shift vectors every step,
490 * because of pressure coupling or box deformation!
492 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
493 calc_shifts(box,fr->shift_vec);
496 put_charge_groups_in_box(fplog,cg0,cg1,fr->ePBC,box,
497 &(top->cgs),x,fr->cg_cm);
498 inc_nrnb(nrnb,eNR_CGCM,homenr);
499 inc_nrnb(nrnb,eNR_RESETX,cg1-cg0);
501 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph) {
502 unshift_self(graph,box,x);
505 else if (bCalcCGCM) {
506 calc_cgcm(fplog,cg0,cg1,&(top->cgs),x,fr->cg_cm);
507 inc_nrnb(nrnb,eNR_CGCM,homenr);
512 move_cgcm(fplog,cr,fr->cg_cm);
515 pr_rvecs(debug,0,"cgcm",fr->cg_cm,top->cgs.nr);
519 if (!(cr->duty & DUTY_PME)) {
520 /* Send particle coordinates to the pme nodes.
521 * Since this is only implemented for domain decomposition
522 * and domain decomposition does not use the graph,
523 * we do not need to worry about shifting.
526 wallcycle_start(wcycle,ewcPP_PMESENDX);
527 GMX_MPE_LOG(ev_send_coordinates_start);
529 bBS = (inputrec->nwall == 2);
532 svmul(inputrec->wall_ewald_zfac,boxs[ZZ],boxs[ZZ]);
535 gmx_pme_send_x(cr,bBS ? boxs : box,x,
536 mdatoms->nChargePerturbed,lambda[efptCOUL],
537 ( flags & GMX_FORCE_VIRIAL),step);
539 GMX_MPE_LOG(ev_send_coordinates_finish);
540 wallcycle_stop(wcycle,ewcPP_PMESENDX);
544 /* Communicate coordinates and sum dipole if necessary */
547 wallcycle_start(wcycle,ewcMOVEX);
548 if (DOMAINDECOMP(cr))
550 dd_move_x(cr->dd,box,x);
554 move_x(fplog,cr,GMX_LEFT,GMX_RIGHT,x,nrnb);
556 /* When we don't need the total dipole we sum it in global_stat */
557 if (bStateChanged && NEED_MUTOT(*inputrec))
559 gmx_sumd(2*DIM,mu,cr);
561 wallcycle_stop(wcycle,ewcMOVEX);
566 /* update adress weight beforehand */
569 /* need pbc for adress weight calculation with pbc_dx */
570 set_pbc(&pbc,inputrec->ePBC,box);
571 if(fr->adress_site == eAdressSITEcog)
573 update_adress_weights_cog(top->idef.iparams,top->idef.il,x,fr,mdatoms,
574 inputrec->ePBC==epbcNONE ? NULL : &pbc);
576 else if (fr->adress_site == eAdressSITEcom)
578 update_adress_weights_com(fplog,cg0,cg1,&(top->cgs),x,fr,mdatoms,
579 inputrec->ePBC==epbcNONE ? NULL : &pbc);
581 else if (fr->adress_site == eAdressSITEatomatom){
582 update_adress_weights_atom_per_atom(cg0,cg1,&(top->cgs),x,fr,mdatoms,
583 inputrec->ePBC==epbcNONE ? NULL : &pbc);
587 update_adress_weights_atom(cg0,cg1,&(top->cgs),x,fr,mdatoms,
588 inputrec->ePBC==epbcNONE ? NULL : &pbc);
596 fr->mu_tot[i][j] = mu[i*DIM + j];
600 if (fr->efep == efepNO)
602 copy_rvec(fr->mu_tot[0],mu_tot);
609 (1.0 - lambda[efptCOUL])*fr->mu_tot[0][j] + lambda[efptCOUL]*fr->mu_tot[1][j];
614 reset_enerdata(&(inputrec->opts),fr,bNS,enerd,MASTER(cr));
615 clear_rvecs(SHIFTS,fr->fshift);
619 wallcycle_start(wcycle,ewcNS);
621 if (graph && bStateChanged)
623 /* Calculate intramolecular shift vectors to make molecules whole */
624 mk_mshift(fplog,graph,fr->ePBC,box,x);
627 /* Reset long range forces if necessary */
630 /* Reset the (long-range) forces if necessary */
631 clear_rvecs(fr->natoms_force_constr,bSepLRF ? fr->f_twin : f);
634 /* Do the actual neighbour searching and if twin range electrostatics
635 * also do the calculation of long range forces and energies.
637 for (i=0;i<efptNR;i++) {dvdlambda[i] = 0;}
639 groups,&(inputrec->opts),top,mdatoms,
640 cr,nrnb,lambda,dvdlambda,&enerd->grpp,bFillGrid,
641 bDoLongRange,bDoForces,bSepLRF ? fr->f_twin : f);
644 fprintf(fplog,sepdvdlformat,"LR non-bonded",0.0,dvdlambda);
646 enerd->dvdl_lin[efptVDW] += dvdlambda[efptVDW];
647 enerd->dvdl_lin[efptCOUL] += dvdlambda[efptCOUL];
649 wallcycle_stop(wcycle,ewcNS);
652 if (inputrec->implicit_solvent && bNS)
654 make_gb_nblist(cr,inputrec->gb_algorithm,inputrec->rlist,
655 x,box,fr,&top->idef,graph,fr->born);
658 if (DOMAINDECOMP(cr))
660 if (!(cr->duty & DUTY_PME))
662 wallcycle_start(wcycle,ewcPPDURINGPME);
663 dd_force_flop_start(cr->dd,nrnb);
669 /* Enforced rotation has its own cycle counter that starts after the collective
670 * coordinates have been communicated. It is added to ddCyclF to allow
671 * for proper load-balancing */
672 wallcycle_start(wcycle,ewcROT);
673 do_rotation(cr,inputrec,box,x,t,step,wcycle,bNS);
674 wallcycle_stop(wcycle,ewcROT);
677 /* Start the force cycle counter.
678 * This counter is stopped in do_forcelow_level.
679 * No parallel communication should occur while this counter is running,
680 * since that will interfere with the dynamic load balancing.
682 wallcycle_start(wcycle,ewcFORCE);
686 /* Reset forces for which the virial is calculated separately:
687 * PME/Ewald forces if necessary */
690 if (flags & GMX_FORCE_VIRIAL)
692 fr->f_novirsum = fr->f_novirsum_alloc;
693 GMX_BARRIER(cr->mpi_comm_mygroup);
696 clear_rvecs(fr->f_novirsum_n,fr->f_novirsum);
700 clear_rvecs(homenr,fr->f_novirsum+start);
702 GMX_BARRIER(cr->mpi_comm_mygroup);
706 /* We are not calculating the pressure so we do not need
707 * a separate array for forces that do not contribute
716 /* Add the long range forces to the short range forces */
717 for(i=0; i<fr->natoms_force_constr; i++)
719 copy_rvec(fr->f_twin[i],f[i]);
722 else if (!(fr->bTwinRange && bNS))
724 /* Clear the short-range forces */
725 clear_rvecs(fr->natoms_force_constr,f);
728 clear_rvec(fr->vir_diag_posres);
730 GMX_BARRIER(cr->mpi_comm_mygroup);
732 if (inputrec->ePull == epullCONSTRAINT)
734 clear_pull_forces(inputrec->pull);
737 /* update QMMMrec, if necessary */
740 update_QMMMrec(cr,fr,x,mdatoms,box,top);
743 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
745 /* Position restraints always require full pbc. Check if we already did it for Adress */
746 if(!(bStateChanged && bDoAdressWF))
748 set_pbc(&pbc,inputrec->ePBC,box);
750 v = posres(top->idef.il[F_POSRES].nr,top->idef.il[F_POSRES].iatoms,
751 top->idef.iparams_posres,
752 (const rvec*)x,fr->f_novirsum,fr->vir_diag_posres,
753 inputrec->ePBC==epbcNONE ? NULL : &pbc,lambda[efptRESTRAINT],&(dvdlambda[efptRESTRAINT]),
754 fr->rc_scaling,fr->ePBC,fr->posres_com,fr->posres_comB);
757 fprintf(fplog,sepdvdlformat,
758 interaction_function[F_POSRES].longname,v,dvdlambda);
760 enerd->term[F_POSRES] += v;
761 /* This linear lambda dependence assumption is only correct
762 * when only k depends on lambda,
763 * not when the reference position depends on lambda.
764 * grompp checks for this. (verify this is still the case?)
766 enerd->dvdl_nonlin[efptRESTRAINT] += dvdlambda[efptRESTRAINT]; /* if just the force constant changes, this is linear,
767 but we can't be sure w/o additional checking that is
768 hard to do at this level of code. Otherwise,
769 the dvdl is not differentiable */
770 inc_nrnb(nrnb,eNR_POSRES,top->idef.il[F_POSRES].nr/2);
771 if ((inputrec->fepvals->n_lambda > 0) && (flags & GMX_FORCE_DHDL))
773 for(i=0; i<enerd->n_lambda; i++)
775 lambda_dum = (i==0 ? lambda[efptRESTRAINT] : inputrec->fepvals->all_lambda[efptRESTRAINT][i-1]);
776 v = posres(top->idef.il[F_POSRES].nr,top->idef.il[F_POSRES].iatoms,
777 top->idef.iparams_posres,
778 (const rvec*)x,NULL,NULL,
779 inputrec->ePBC==epbcNONE ? NULL : &pbc,lambda_dum,&dvdl_dum,
780 fr->rc_scaling,fr->ePBC,fr->posres_com,fr->posres_comB);
781 enerd->enerpart_lambda[i] += v;
786 /* Compute the bonded and non-bonded energies and optionally forces */
787 do_force_lowlevel(fplog,step,fr,inputrec,&(top->idef),
788 cr,nrnb,wcycle,mdatoms,&(inputrec->opts),
789 x,hist,f,enerd,fcd,mtop,top,fr->born,
790 &(top->atomtypes),bBornRadii,box,
791 inputrec->fepvals,lambda,graph,&(top->excls),fr->mu_tot,
794 cycles_force = wallcycle_stop(wcycle,ewcFORCE);
795 GMX_BARRIER(cr->mpi_comm_mygroup);
799 do_flood(fplog,cr,x,f,ed,box,step);
802 if (DOMAINDECOMP(cr))
804 dd_force_flop_stop(cr->dd,nrnb);
807 dd_cycles_add(cr->dd,cycles_force-cycles_pme,ddCyclF);
813 if (IR_ELEC_FIELD(*inputrec))
815 /* Compute forces due to electric field */
816 calc_f_el(MASTER(cr) ? field : NULL,
817 start,homenr,mdatoms->chargeA,x,fr->f_novirsum,
818 inputrec->ex,inputrec->et,t);
821 if (bDoAdressWF && fr->adress_icor == eAdressICThermoForce)
823 /* Compute thermodynamic force in hybrid AdResS region */
824 adress_thermo_force(start,homenr,&(top->cgs),x,fr->f_novirsum,fr,mdatoms,
825 inputrec->ePBC==epbcNONE ? NULL : &pbc);
828 /* Communicate the forces */
831 wallcycle_start(wcycle,ewcMOVEF);
832 if (DOMAINDECOMP(cr))
834 dd_move_f(cr->dd,f,fr->fshift);
835 /* Do we need to communicate the separate force array
836 * for terms that do not contribute to the single sum virial?
837 * Position restraints and electric fields do not introduce
838 * inter-cg forces, only full electrostatics methods do.
839 * When we do not calculate the virial, fr->f_novirsum = f,
840 * so we have already communicated these forces.
842 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
843 (flags & GMX_FORCE_VIRIAL))
845 dd_move_f(cr->dd,fr->f_novirsum,NULL);
849 /* We should not update the shift forces here,
850 * since f_twin is already included in f.
852 dd_move_f(cr->dd,fr->f_twin,NULL);
857 pd_move_f(cr,f,nrnb);
860 pd_move_f(cr,fr->f_twin,nrnb);
863 wallcycle_stop(wcycle,ewcMOVEF);
866 /* If we have NoVirSum forces, but we do not calculate the virial,
867 * we sum fr->f_novirum=f later.
869 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
871 wallcycle_start(wcycle,ewcVSITESPREAD);
872 spread_vsite_f(fplog,vsite,x,f,fr->fshift,FALSE,NULL,nrnb,
873 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
874 wallcycle_stop(wcycle,ewcVSITESPREAD);
878 wallcycle_start(wcycle,ewcVSITESPREAD);
879 spread_vsite_f(fplog,vsite,x,fr->f_twin,NULL,FALSE,NULL,
881 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
882 wallcycle_stop(wcycle,ewcVSITESPREAD);
886 if (flags & GMX_FORCE_VIRIAL)
888 /* Calculation of the virial must be done after vsites! */
889 calc_virial(fplog,mdatoms->start,mdatoms->homenr,x,f,
890 vir_force,graph,box,nrnb,fr,inputrec->ePBC);
894 enerd->term[F_COM_PULL] = 0;
895 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
897 /* Calculate the center of mass forces, this requires communication,
898 * which is why pull_potential is called close to other communication.
899 * The virial contribution is calculated directly,
900 * which is why we call pull_potential after calc_virial.
902 set_pbc(&pbc,inputrec->ePBC,box);
903 dvdlambda[efptRESTRAINT] = 0;
904 enerd->term[F_COM_PULL] +=
905 pull_potential(inputrec->ePull,inputrec->pull,mdatoms,&pbc,
906 cr,t,lambda[efptRESTRAINT],x,f,vir_force,&(dvdlambda[efptRESTRAINT]));
909 fprintf(fplog,sepdvdlformat,"Com pull",enerd->term[F_COM_PULL],dvdlambda[efptRESTRAINT]);
911 enerd->dvdl_lin[efptRESTRAINT] += dvdlambda[efptRESTRAINT];
914 /* Add the forces from enforced rotation potentials (if any) */
917 wallcycle_start(wcycle,ewcROTadd);
918 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr,step,t);
919 wallcycle_stop(wcycle,ewcROTadd);
922 if (PAR(cr) && !(cr->duty & DUTY_PME))
924 cycles_ppdpme = wallcycle_stop(wcycle,ewcPPDURINGPME);
925 dd_cycles_add(cr->dd,cycles_ppdpme,ddCyclPPduringPME);
927 /* In case of node-splitting, the PP nodes receive the long-range
928 * forces, virial and energy from the PME nodes here.
930 wallcycle_start(wcycle,ewcPP_PMEWAITRECVF);
931 dvdlambda[efptCOUL] = 0;
932 gmx_pme_receive_f(cr,fr->f_novirsum,fr->vir_el_recip,&e,&dvdlambda[efptCOUL],
936 fprintf(fplog,sepdvdlformat,"PME mesh",e,dvdlambda[efptCOUL]);
938 enerd->term[F_COUL_RECIP] += e;
939 enerd->dvdl_lin[efptCOUL] += dvdlambda[efptCOUL];
942 dd_cycles_add(cr->dd,cycles_seppme,ddCyclPME);
944 wallcycle_stop(wcycle,ewcPP_PMEWAITRECVF);
947 if (bDoForces && fr->bF_NoVirSum)
951 /* Spread the mesh force on virtual sites to the other particles...
952 * This is parallellized. MPI communication is performed
953 * if the constructing atoms aren't local.
955 wallcycle_start(wcycle,ewcVSITESPREAD);
956 spread_vsite_f(fplog,vsite,x,fr->f_novirsum,NULL,
957 (flags & GMX_FORCE_VIRIAL),fr->vir_el_recip,
959 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
960 wallcycle_stop(wcycle,ewcVSITESPREAD);
962 if (flags & GMX_FORCE_VIRIAL)
964 /* Now add the forces, this is local */
967 sum_forces(0,fr->f_novirsum_n,f,fr->f_novirsum);
971 sum_forces(start,start+homenr,f,fr->f_novirsum);
973 if (EEL_FULL(fr->eeltype))
975 /* Add the mesh contribution to the virial */
976 m_add(vir_force,fr->vir_el_recip,vir_force);
980 pr_rvecs(debug,0,"vir_force",vir_force,DIM);
985 /* Sum the potential energy terms from group contributions */
986 sum_epot(&(inputrec->opts),enerd);
988 if (fr->print_force >= 0 && bDoForces)
990 print_large_forces(stderr,mdatoms,cr,step,fr->print_force,x,f);
994 void do_constrain_first(FILE *fplog,gmx_constr_t constr,
995 t_inputrec *ir,t_mdatoms *md,
996 t_state *state,rvec *f,
997 t_graph *graph,t_commrec *cr,t_nrnb *nrnb,
998 t_forcerec *fr, gmx_localtop_t *top, tensor shake_vir)
1001 gmx_large_int_t step;
1002 real dt=ir->delta_t;
1006 snew(savex,state->natoms);
1009 end = md->homenr + start;
1012 fprintf(debug,"vcm: start=%d, homenr=%d, end=%d\n",
1013 start,md->homenr,end);
1014 /* Do a first constrain to reset particles... */
1015 step = ir->init_step;
1018 char buf[STEPSTRSIZE];
1019 fprintf(fplog,"\nConstraining the starting coordinates (step %s)\n",
1020 gmx_step_str(step,buf));
1024 /* constrain the current position */
1025 constrain(NULL,TRUE,FALSE,constr,&(top->idef),
1026 ir,NULL,cr,step,0,md,
1027 state->x,state->x,NULL,
1028 state->box,state->lambda[efptBONDED],&dvdl_dum,
1029 NULL,NULL,nrnb,econqCoord,ir->epc==epcMTTK,state->veta,state->veta);
1032 /* constrain the inital velocity, and save it */
1033 /* also may be useful if we need the ekin from the halfstep for velocity verlet */
1034 /* might not yet treat veta correctly */
1035 constrain(NULL,TRUE,FALSE,constr,&(top->idef),
1036 ir,NULL,cr,step,0,md,
1037 state->x,state->v,state->v,
1038 state->box,state->lambda[efptBONDED],&dvdl_dum,
1039 NULL,NULL,nrnb,econqVeloc,ir->epc==epcMTTK,state->veta,state->veta);
1041 /* constrain the inital velocities at t-dt/2 */
1042 if (EI_STATE_VELOCITY(ir->eI) && ir->eI!=eiVV)
1044 for(i=start; (i<end); i++)
1046 for(m=0; (m<DIM); m++)
1048 /* Reverse the velocity */
1049 state->v[i][m] = -state->v[i][m];
1050 /* Store the position at t-dt in buf */
1051 savex[i][m] = state->x[i][m] + dt*state->v[i][m];
1054 /* Shake the positions at t=-dt with the positions at t=0
1055 * as reference coordinates.
1059 char buf[STEPSTRSIZE];
1060 fprintf(fplog,"\nConstraining the coordinates at t0-dt (step %s)\n",
1061 gmx_step_str(step,buf));
1064 constrain(NULL,TRUE,FALSE,constr,&(top->idef),
1065 ir,NULL,cr,step,-1,md,
1066 state->x,savex,NULL,
1067 state->box,state->lambda[efptBONDED],&dvdl_dum,
1068 state->v,NULL,nrnb,econqCoord,ir->epc==epcMTTK,state->veta,state->veta);
1070 for(i=start; i<end; i++) {
1071 for(m=0; m<DIM; m++) {
1072 /* Re-reverse the velocities */
1073 state->v[i][m] = -state->v[i][m];
1080 void calc_enervirdiff(FILE *fplog,int eDispCorr,t_forcerec *fr)
1082 double eners[2],virs[2],enersum,virsum,y0,f,g,h;
1083 double r0,r1,r,rc3,rc9,ea,eb,ec,pa,pb,pc,pd;
1084 double invscale,invscale2,invscale3;
1085 int ri0,ri1,ri,i,offstart,offset;
1088 fr->enershiftsix = 0;
1089 fr->enershifttwelve = 0;
1090 fr->enerdiffsix = 0;
1091 fr->enerdifftwelve = 0;
1093 fr->virdifftwelve = 0;
1095 if (eDispCorr != edispcNO) {
1096 for(i=0; i<2; i++) {
1100 if ((fr->vdwtype == evdwSWITCH) || (fr->vdwtype == evdwSHIFT)) {
1101 if (fr->rvdw_switch == 0)
1103 "With dispersion correction rvdw-switch can not be zero "
1104 "for vdw-type = %s",evdw_names[fr->vdwtype]);
1106 scale = fr->nblists[0].tab.scale;
1107 vdwtab = fr->nblists[0].vdwtab;
1109 /* Round the cut-offs to exact table values for precision */
1110 ri0 = floor(fr->rvdw_switch*scale);
1111 ri1 = ceil(fr->rvdw*scale);
1117 if (fr->vdwtype == evdwSHIFT) {
1118 /* Determine the constant energy shift below rvdw_switch */
1119 fr->enershiftsix = (real)(-1.0/(rc3*rc3)) - vdwtab[8*ri0];
1120 fr->enershifttwelve = (real)( 1.0/(rc9*rc3)) - vdwtab[8*ri0 + 4];
1122 /* Add the constant part from 0 to rvdw_switch.
1123 * This integration from 0 to rvdw_switch overcounts the number
1124 * of interactions by 1, as it also counts the self interaction.
1125 * We will correct for this later.
1127 eners[0] += 4.0*M_PI*fr->enershiftsix*rc3/3.0;
1128 eners[1] += 4.0*M_PI*fr->enershifttwelve*rc3/3.0;
1130 invscale = 1.0/(scale);
1131 invscale2 = invscale*invscale;
1132 invscale3 = invscale*invscale2;
1134 /* following summation derived from cubic spline definition,
1135 Numerical Recipies in C, second edition, p. 113-116. Exact
1136 for the cubic spline. We first calculate the negative of
1137 the energy from rvdw to rvdw_switch, assuming that g(r)=1,
1138 and then add the more standard, abrupt cutoff correction to
1139 that result, yielding the long-range correction for a
1140 switched function. We perform both the pressure and energy
1141 loops at the same time for simplicity, as the computational
1145 enersum = 0.0; virsum = 0.0;
1150 for (ri=ri0; ri<ri1; ri++) {
1153 eb = 2.0*invscale2*r;
1157 pb = 3.0*invscale2*r;
1158 pc = 3.0*invscale*r*r;
1161 /* this "8" is from the packing in the vdwtab array - perhaps
1162 should be #define'ed? */
1163 offset = 8*ri + offstart;
1164 y0 = vdwtab[offset];
1165 f = vdwtab[offset+1];
1166 g = vdwtab[offset+2];
1167 h = vdwtab[offset+3];
1169 enersum += y0*(ea/3 + eb/2 + ec) + f*(ea/4 + eb/3 + ec/2)+
1170 g*(ea/5 + eb/4 + ec/3) + h*(ea/6 + eb/5 + ec/4);
1171 virsum += f*(pa/4 + pb/3 + pc/2 + pd) +
1172 2*g*(pa/5 + pb/4 + pc/3 + pd/2) + 3*h*(pa/6 + pb/5 + pc/4 + pd/3);
1175 enersum *= 4.0*M_PI;
1177 eners[i] -= enersum;
1181 /* now add the correction for rvdw_switch to infinity */
1182 eners[0] += -4.0*M_PI/(3.0*rc3);
1183 eners[1] += 4.0*M_PI/(9.0*rc9);
1184 virs[0] += 8.0*M_PI/rc3;
1185 virs[1] += -16.0*M_PI/(3.0*rc9);
1187 else if ((fr->vdwtype == evdwCUT) || (fr->vdwtype == evdwUSER)) {
1188 if (fr->vdwtype == evdwUSER && fplog)
1190 "WARNING: using dispersion correction with user tables\n");
1191 rc3 = fr->rvdw*fr->rvdw*fr->rvdw;
1193 eners[0] += -4.0*M_PI/(3.0*rc3);
1194 eners[1] += 4.0*M_PI/(9.0*rc9);
1195 virs[0] += 8.0*M_PI/rc3;
1196 virs[1] += -16.0*M_PI/(3.0*rc9);
1199 "Dispersion correction is not implemented for vdw-type = %s",
1200 evdw_names[fr->vdwtype]);
1202 fr->enerdiffsix = eners[0];
1203 fr->enerdifftwelve = eners[1];
1204 /* The 0.5 is due to the Gromacs definition of the virial */
1205 fr->virdiffsix = 0.5*virs[0];
1206 fr->virdifftwelve = 0.5*virs[1];
1210 void calc_dispcorr(FILE *fplog,t_inputrec *ir,t_forcerec *fr,
1211 gmx_large_int_t step,int natoms,
1212 matrix box,real lambda,tensor pres,tensor virial,
1213 real *prescorr, real *enercorr, real *dvdlcorr)
1215 gmx_bool bCorrAll,bCorrPres;
1216 real dvdlambda,invvol,dens,ninter,avcsix,avctwelve,enerdiff,svir=0,spres=0;
1226 if (ir->eDispCorr != edispcNO) {
1227 bCorrAll = (ir->eDispCorr == edispcAllEner ||
1228 ir->eDispCorr == edispcAllEnerPres);
1229 bCorrPres = (ir->eDispCorr == edispcEnerPres ||
1230 ir->eDispCorr == edispcAllEnerPres);
1232 invvol = 1/det(box);
1235 /* Only correct for the interactions with the inserted molecule */
1236 dens = (natoms - fr->n_tpi)*invvol;
1241 dens = natoms*invvol;
1242 ninter = 0.5*natoms;
1245 if (ir->efep == efepNO)
1247 avcsix = fr->avcsix[0];
1248 avctwelve = fr->avctwelve[0];
1252 avcsix = (1 - lambda)*fr->avcsix[0] + lambda*fr->avcsix[1];
1253 avctwelve = (1 - lambda)*fr->avctwelve[0] + lambda*fr->avctwelve[1];
1256 enerdiff = ninter*(dens*fr->enerdiffsix - fr->enershiftsix);
1257 *enercorr += avcsix*enerdiff;
1259 if (ir->efep != efepNO)
1261 dvdlambda += (fr->avcsix[1] - fr->avcsix[0])*enerdiff;
1265 enerdiff = ninter*(dens*fr->enerdifftwelve - fr->enershifttwelve);
1266 *enercorr += avctwelve*enerdiff;
1267 if (fr->efep != efepNO)
1269 dvdlambda += (fr->avctwelve[1] - fr->avctwelve[0])*enerdiff;
1275 svir = ninter*dens*avcsix*fr->virdiffsix/3.0;
1276 if (ir->eDispCorr == edispcAllEnerPres)
1278 svir += ninter*dens*avctwelve*fr->virdifftwelve/3.0;
1280 /* The factor 2 is because of the Gromacs virial definition */
1281 spres = -2.0*invvol*svir*PRESFAC;
1283 for(m=0; m<DIM; m++) {
1284 virial[m][m] += svir;
1285 pres[m][m] += spres;
1290 /* Can't currently control when it prints, for now, just print when degugging */
1294 fprintf(debug,"Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
1300 "Long Range LJ corr.: Epot %10g, Pres: %10g, Vir: %10g\n",
1301 *enercorr,spres,svir);
1305 fprintf(debug,"Long Range LJ corr.: Epot %10g\n",*enercorr);
1309 if (fr->bSepDVDL && do_per_step(step,ir->nstlog))
1311 fprintf(fplog,sepdvdlformat,"Dispersion correction",
1312 *enercorr,dvdlambda);
1314 if (fr->efep != efepNO)
1316 *dvdlcorr += dvdlambda;
1321 void do_pbc_first(FILE *fplog,matrix box,t_forcerec *fr,
1322 t_graph *graph,rvec x[])
1325 fprintf(fplog,"Removing pbc first time\n");
1326 calc_shifts(box,fr->shift_vec);
1328 mk_mshift(fplog,graph,fr->ePBC,box,x);
1330 p_graph(debug,"do_pbc_first 1",graph);
1331 shift_self(graph,box,x);
1332 /* By doing an extra mk_mshift the molecules that are broken
1333 * because they were e.g. imported from another software
1334 * will be made whole again. Such are the healing powers
1337 mk_mshift(fplog,graph,fr->ePBC,box,x);
1339 p_graph(debug,"do_pbc_first 2",graph);
1342 fprintf(fplog,"Done rmpbc\n");
1345 static void low_do_pbc_mtop(FILE *fplog,int ePBC,matrix box,
1346 gmx_mtop_t *mtop,rvec x[],
1351 gmx_molblock_t *molb;
1353 if (bFirst && fplog)
1354 fprintf(fplog,"Removing pbc first time\n");
1358 for(mb=0; mb<mtop->nmolblock; mb++) {
1359 molb = &mtop->molblock[mb];
1360 if (molb->natoms_mol == 1 ||
1361 (!bFirst && mtop->moltype[molb->type].cgs.nr == 1)) {
1362 /* Just one atom or charge group in the molecule, no PBC required */
1363 as += molb->nmol*molb->natoms_mol;
1365 /* Pass NULL iso fplog to avoid graph prints for each molecule type */
1366 mk_graph_ilist(NULL,mtop->moltype[molb->type].ilist,
1367 0,molb->natoms_mol,FALSE,FALSE,graph);
1369 for(mol=0; mol<molb->nmol; mol++) {
1370 mk_mshift(fplog,graph,ePBC,box,x+as);
1372 shift_self(graph,box,x+as);
1373 /* The molecule is whole now.
1374 * We don't need the second mk_mshift call as in do_pbc_first,
1375 * since we no longer need this graph.
1378 as += molb->natoms_mol;
1386 void do_pbc_first_mtop(FILE *fplog,int ePBC,matrix box,
1387 gmx_mtop_t *mtop,rvec x[])
1389 low_do_pbc_mtop(fplog,ePBC,box,mtop,x,TRUE);
1392 void do_pbc_mtop(FILE *fplog,int ePBC,matrix box,
1393 gmx_mtop_t *mtop,rvec x[])
1395 low_do_pbc_mtop(fplog,ePBC,box,mtop,x,FALSE);
1398 void finish_run(FILE *fplog,t_commrec *cr,const char *confout,
1399 t_inputrec *inputrec,
1400 t_nrnb nrnb[],gmx_wallcycle_t wcycle,
1401 gmx_runtime_t *runtime,
1402 gmx_bool bWriteStat)
1405 t_nrnb *nrnb_tot=NULL;
1408 double cycles[ewcNR];
1410 wallcycle_sum(cr,wcycle,cycles);
1412 if (cr->nnodes > 1) {
1416 MPI_Reduce(nrnb->n,nrnb_tot->n,eNRNB,MPI_DOUBLE,MPI_SUM,
1417 MASTERRANK(cr),cr->mpi_comm_mysim);
1423 if (SIMMASTER(cr)) {
1424 print_flop(fplog,nrnb_tot,&nbfs,&mflop);
1425 if (cr->nnodes > 1) {
1430 if ((cr->duty & DUTY_PP) && DOMAINDECOMP(cr)) {
1431 print_dd_statistics(cr,inputrec,fplog);
1443 snew(nrnb_all,cr->nnodes);
1444 nrnb_all[0] = *nrnb;
1445 for(s=1; s<cr->nnodes; s++)
1447 MPI_Recv(nrnb_all[s].n,eNRNB,MPI_DOUBLE,s,0,
1448 cr->mpi_comm_mysim,&stat);
1450 pr_load(fplog,cr,nrnb_all);
1455 MPI_Send(nrnb->n,eNRNB,MPI_DOUBLE,MASTERRANK(cr),0,
1456 cr->mpi_comm_mysim);
1461 if (SIMMASTER(cr)) {
1462 wallcycle_print(fplog,cr->nnodes,cr->npmenodes,runtime->realtime,
1465 if (EI_DYNAMICS(inputrec->eI)) {
1466 delta_t = inputrec->delta_t;
1472 print_perf(fplog,runtime->proctime,runtime->realtime,
1473 cr->nnodes-cr->npmenodes,
1474 runtime->nsteps_done,delta_t,nbfs,mflop);
1477 print_perf(stderr,runtime->proctime,runtime->realtime,
1478 cr->nnodes-cr->npmenodes,
1479 runtime->nsteps_done,delta_t,nbfs,mflop);
1483 runtime=inputrec->nsteps*inputrec->delta_t;
1485 if (cr->nnodes == 1)
1486 fprintf(stderr,"\n\n");
1487 print_perf(stderr,nodetime,realtime,runtime,&ntot,
1488 cr->nnodes-cr->npmenodes,FALSE);
1490 wallcycle_print(fplog,cr->nnodes,cr->npmenodes,realtime,wcycle,cycles);
1491 print_perf(fplog,nodetime,realtime,runtime,&ntot,cr->nnodes-cr->npmenodes,
1494 pr_load(fplog,cr,nrnb_all);
1501 extern void initialize_lambdas(FILE *fplog,t_inputrec *ir,int *fep_state,real *lambda,double *lam0)
1503 /* this function works, but could probably use a logic rewrite to keep all the different
1504 types of efep straight. */
1507 t_lambda *fep = ir->fepvals;
1509 if ((ir->efep==efepNO) && (ir->bSimTemp == FALSE)) {
1510 for (i=0;i<efptNR;i++) {
1519 *fep_state = fep->init_fep_state; /* this might overwrite the checkpoint
1520 if checkpoint is set -- a kludge is in for now
1522 for (i=0;i<efptNR;i++)
1524 /* overwrite lambda state with init_lambda for now for backwards compatibility */
1525 if (fep->init_lambda>=0) /* if it's -1, it was never initializd */
1527 lambda[i] = fep->init_lambda;
1529 lam0[i] = lambda[i];
1534 lambda[i] = fep->all_lambda[i][*fep_state];
1536 lam0[i] = lambda[i];
1541 /* need to rescale control temperatures to match current state */
1542 for (i=0;i<ir->opts.ngtc;i++) {
1543 if (ir->opts.ref_t[i] > 0) {
1544 ir->opts.ref_t[i] = ir->simtempvals->temperatures[*fep_state];
1550 /* Send to the log the information on the current lambdas */
1553 fprintf(fplog,"Initial vector of lambda components:[ ");
1554 for (i=0;i<efptNR;i++)
1556 fprintf(fplog,"%10.4f ",lambda[i]);
1558 fprintf(fplog,"]\n");
1564 void init_md(FILE *fplog,
1565 t_commrec *cr,t_inputrec *ir,const output_env_t oenv,
1566 double *t,double *t0,
1567 real *lambda, int *fep_state, double *lam0,
1568 t_nrnb *nrnb,gmx_mtop_t *mtop,
1570 int nfile,const t_filenm fnm[],
1571 gmx_mdoutf_t **outf,t_mdebin **mdebin,
1572 tensor force_vir,tensor shake_vir,rvec mu_tot,
1573 gmx_bool *bSimAnn,t_vcm **vcm, t_state *state, unsigned long Flags)
1578 /* Initial values */
1579 *t = *t0 = ir->init_t;
1582 for(i=0;i<ir->opts.ngtc;i++)
1584 /* set bSimAnn if any group is being annealed */
1585 if(ir->opts.annealing[i]!=eannNO)
1592 update_annealing_target_temp(&(ir->opts),ir->init_t);
1595 /* Initialize lambda variables */
1596 initialize_lambdas(fplog,ir,fep_state,lambda,lam0);
1600 *upd = init_update(fplog,ir);
1606 *vcm = init_vcm(fplog,&mtop->groups,ir);
1609 if (EI_DYNAMICS(ir->eI) && !(Flags & MD_APPENDFILES))
1611 if (ir->etc == etcBERENDSEN)
1613 please_cite(fplog,"Berendsen84a");
1615 if (ir->etc == etcVRESCALE)
1617 please_cite(fplog,"Bussi2007a");
1625 *outf = init_mdoutf(nfile,fnm,Flags,cr,ir,oenv);
1627 *mdebin = init_mdebin((Flags & MD_APPENDFILES) ? NULL : (*outf)->fp_ene,
1628 mtop,ir, (*outf)->fp_dhdl);
1633 please_cite(fplog,"Fritsch12");
1634 please_cite(fplog,"Junghans10");
1636 /* Initiate variables */
1637 clear_mat(force_vir);
1638 clear_mat(shake_vir);