Counter resetting also on the PME only nodes
[gromacs/AngularHB.git] / src / mdlib / sim_util.c
blobcd088061d174f5bc8e5b1ab0ce798084fab2d337
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
3 *
4 * This source code is part of
5 *
6 * G R O M A C S
7 *
8 * GROningen MAchine for Chemical Simulations
9 *
10 * VERSION 3.2.0
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
33 * And Hey:
34 * GROwing Monsters And Cloning Shrimps
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
40 #ifdef GMX_CRAY_XT3
41 #include<catamount/dclock.h>
42 #endif
45 #include <stdio.h>
46 #include <time.h>
47 #ifdef HAVE_SYS_TIME_H
48 #include <sys/time.h>
49 #endif
50 #include <math.h>
51 #include "typedefs.h"
52 #include "string2.h"
53 #include "gmxfio.h"
54 #include "smalloc.h"
55 #include "names.h"
56 #include "confio.h"
57 #include "mvdata.h"
58 #include "txtdump.h"
59 #include "pbc.h"
60 #include "vec.h"
61 #include "time.h"
62 #include "nrnb.h"
63 #include "mshift.h"
64 #include "mdrun.h"
65 #include "update.h"
66 #include "physics.h"
67 #include "main.h"
68 #include "mdatoms.h"
69 #include "force.h"
70 #include "bondf.h"
71 #include "pme.h"
72 #include "pppm.h"
73 #include "disre.h"
74 #include "orires.h"
75 #include "network.h"
76 #include "calcmu.h"
77 #include "constr.h"
78 #include "xvgr.h"
79 #include "trnio.h"
80 #include "xtcio.h"
81 #include "copyrite.h"
83 #include "mpelogging.h"
84 #include "domdec.h"
85 #include "partdec.h"
86 #include "gmx_wallcycle.h"
87 #include "genborn.h"
89 #ifdef GMX_LIB_MPI
90 #include <mpi.h>
91 #endif
92 #ifdef GMX_THREAD_MPI
93 #include "thread_mpi.h"
94 #endif
97 #include "qmmm.h"
99 typedef struct gmx_timeprint {
101 } t_gmx_timeprint;
104 double
105 gmx_gettime()
107 #ifdef HAVE_GETTIMEOFDAY
108 struct timeval t;
109 struct timezone tz = { 0,0 };
110 double seconds;
112 gettimeofday(&t,&tz);
114 seconds = (double) t.tv_sec + 1e-6*(double)t.tv_usec;
116 return seconds;
117 #else
118 double seconds;
120 seconds = time(NULL);
122 return seconds;
123 #endif
127 #define difftime(end,start) ((double)(end)-(double)(start))
129 void print_time(FILE *out,gmx_runtime_t *runtime,gmx_step_t step,t_inputrec *ir)
131 time_t finish;
133 double dt;
134 char buf[48];
136 if (!gmx_parallel_env)
137 fprintf(out,"\r");
138 fprintf(out,"step %s",gmx_step_str(step,buf));
139 if ((step >= ir->nstlist)) {
140 if ((ir->nstlist == 0) || ((step % ir->nstlist) == 0)) {
141 /* We have done a full cycle let's update time_per_step */
142 runtime->last = gmx_gettime();
143 dt = difftime(runtime->last,runtime->real);
144 runtime->time_per_step = dt/(step - ir->init_step + 1);
146 dt = (ir->nsteps + ir->init_step - step)*runtime->time_per_step;
148 if (dt >= 300) {
149 finish = (time_t) (runtime->last + dt);
150 sprintf(buf,"%s",ctime(&finish));
151 buf[strlen(buf)-1]='\0';
152 fprintf(out,", will finish %s",buf);
154 else
155 fprintf(out,", remaining runtime: %5d s ",(int)dt);
157 if (gmx_parallel_env)
158 fprintf(out,"\n");
160 fflush(out);
163 #ifdef NO_CLOCK
164 #define clock() -1
165 #endif
167 static double set_proctime(gmx_runtime_t *runtime)
169 double diff;
170 #ifdef GMX_CRAY_XT3
171 double prev;
173 prev = runtime->proc;
174 runtime->proc = dclock();
176 diff = runtime->proc - prev;
177 #else
178 clock_t prev;
180 prev = runtime->proc;
181 runtime->proc = clock();
183 diff = (double)(runtime->proc - prev)/(double)CLOCKS_PER_SEC;
184 #endif
185 if (diff < 0)
187 /* The counter has probably looped, ignore this data */
188 diff = 0;
191 return diff;
194 void runtime_start(gmx_runtime_t *runtime)
196 runtime->real = gmx_gettime();
197 set_proctime(runtime);
198 runtime->realtime = 0;
199 runtime->proctime = 0;
200 runtime->last = 0;
201 runtime->time_per_step = 0;
204 void runtime_end(gmx_runtime_t *runtime)
206 double now;
208 now = gmx_gettime();
210 runtime->proctime += set_proctime(runtime);
211 runtime->realtime = now - runtime->real;
212 runtime->real = now;
215 void runtime_upd_proc(gmx_runtime_t *runtime)
217 runtime->proctime += set_proctime(runtime);
220 void print_date_and_time(FILE *fplog,int nodeid,const char *title,
221 const gmx_runtime_t *runtime)
223 int i;
224 char *ts,time_string[STRLEN];
225 time_t tmptime;
227 if (runtime != NULL)
229 tmptime = (time_t) runtime->real;
230 ts = ctime(&tmptime);
232 else
234 tmptime = (time_t) gmx_gettime();
235 ts = ctime(&tmptime);
237 for(i=0; ts[i]>=' '; i++)
239 time_string[i]=ts[i];
241 time_string[i]='\0';
242 if (fplog)
244 fprintf(fplog,"%s on node %d %s\n",title,nodeid,time_string);
248 static void sum_forces(int start,int end,rvec f[],rvec flr[])
250 int i;
252 if (gmx_debug_at) {
253 pr_rvecs(debug,0,"fsr",f+start,end-start);
254 pr_rvecs(debug,0,"flr",flr+start,end-start);
256 for(i=start; (i<end); i++)
257 rvec_inc(f[i],flr[i]);
261 * calc_f_el calculates forces due to an electric field.
263 * force is kJ mol^-1 nm^-1 = e * kJ mol^-1 nm^-1 / e
265 * Et[] contains the parameters for the time dependent
266 * part of the field (not yet used).
267 * Ex[] contains the parameters for
268 * the spatial dependent part of the field. You can have cool periodic
269 * fields in principle, but only a constant field is supported
270 * now.
271 * The function should return the energy due to the electric field
272 * (if any) but for now returns 0.
274 * WARNING:
275 * There can be problems with the virial.
276 * Since the field is not self-consistent this is unavoidable.
277 * For neutral molecules the virial is correct within this approximation.
278 * For neutral systems with many charged molecules the error is small.
279 * But for systems with a net charge or a few charged molecules
280 * the error can be significant when the field is high.
281 * Solution: implement a self-consitent electric field into PME.
283 static void calc_f_el(FILE *fp,int start,int homenr,
284 real charge[],rvec x[],rvec f[],
285 t_cosines Ex[],t_cosines Et[],double t)
287 rvec Ext;
288 real t0;
289 int i,m;
291 for(m=0; (m<DIM); m++) {
292 if (Et[m].n) {
293 if (Et[m].n == 3) {
294 t0 = Et[m].a[1];
295 Ext[m] = cos(Et[m].a[0]*(t-t0))*exp(-sqr(t-t0)/(2.0*sqr(Et[m].a[2])));
297 else
298 Ext[m] = cos(Et[m].a[0]*t);
300 else
301 Ext[m] = 1.0;
302 if (Ex[m].n) {
303 /* Convert the field strength from V/nm to MD-units */
304 Ext[m] *= Ex[m].a[0]*FIELDFAC;
305 for(i=start; (i<start+homenr); i++)
306 f[i][m] += charge[i]*Ext[m];
308 else
309 Ext[m] = 0;
311 if (fp)
312 fprintf(fp,"%10g %10g %10g %10g #FIELD\n",t,
313 Ext[XX]/FIELDFAC,Ext[YY]/FIELDFAC,Ext[ZZ]/FIELDFAC);
316 static void calc_virial(FILE *fplog,int start,int homenr,rvec x[],rvec f[],
317 tensor vir_part,t_graph *graph,matrix box,
318 t_nrnb *nrnb,const t_forcerec *fr,int ePBC)
320 int i,j;
321 tensor virtest;
323 /* The short-range virial from surrounding boxes */
324 clear_mat(vir_part);
325 calc_vir(fplog,SHIFTS,fr->shift_vec,fr->fshift,vir_part,ePBC==epbcSCREW,box);
326 inc_nrnb(nrnb,eNR_VIRIAL,SHIFTS);
328 /* Calculate partial virial, for local atoms only, based on short range.
329 * Total virial is computed in global_stat, called from do_md
331 f_calc_vir(fplog,start,start+homenr,x,f,vir_part,graph,box);
332 inc_nrnb(nrnb,eNR_VIRIAL,homenr);
334 /* Add position restraint contribution */
335 for(i=0; i<DIM; i++) {
336 vir_part[i][i] += fr->vir_diag_posres[i];
339 /* Add wall contribution */
340 for(i=0; i<DIM; i++) {
341 vir_part[i][ZZ] += fr->vir_wall_z[i];
344 if (debug)
345 pr_rvecs(debug,0,"vir_part",vir_part,DIM);
348 static void print_large_forces(FILE *fp,t_mdatoms *md,t_commrec *cr,
349 gmx_step_t step,real pforce,rvec *x,rvec *f)
351 int i;
352 real pf2,fn2;
353 char buf[22];
355 pf2 = sqr(pforce);
356 for(i=md->start; i<md->start+md->homenr; i++) {
357 fn2 = norm2(f[i]);
358 if (fn2 >= pf2) {
359 fprintf(fp,"step %s atom %6d x %8.3f %8.3f %8.3f force %12.5e\n",
360 gmx_step_str(step,buf),
361 ddglatnr(cr->dd,i),x[i][XX],x[i][YY],x[i][ZZ],sqrt(fn2));
366 void do_force(FILE *fplog,t_commrec *cr,
367 t_inputrec *inputrec,
368 gmx_step_t step,t_nrnb *nrnb,gmx_wallcycle_t wcycle,
369 gmx_localtop_t *top,
370 gmx_mtop_t *mtop,
371 gmx_groups_t *groups,
372 matrix box,rvec x[],history_t *hist,
373 rvec f[],
374 tensor vir_force,
375 t_mdatoms *mdatoms,
376 gmx_enerdata_t *enerd,t_fcdata *fcd,
377 real lambda,t_graph *graph,
378 t_forcerec *fr,gmx_vsite_t *vsite,rvec mu_tot,
379 double t,FILE *field,gmx_edsam_t ed,
380 bool bBornRadii,
381 int flags)
383 int cg0,cg1,i,j;
384 int start,homenr;
385 double mu[2*DIM];
386 bool bSepDVDL,bStateChanged,bNS,bFillGrid,bCalcCGCM,bBS,bDoForces;
387 matrix boxs;
388 real e,v,dvdl;
389 t_pbc pbc;
390 float cycles_ppdpme,cycles_pme,cycles_seppme,cycles_force;
392 start = mdatoms->start;
393 homenr = mdatoms->homenr;
395 bSepDVDL = (fr->bSepDVDL && do_per_step(step,inputrec->nstlog));
397 clear_mat(vir_force);
399 if (PARTDECOMP(cr)) {
400 pd_cg_range(cr,&cg0,&cg1);
401 } else {
402 cg0 = 0;
403 if (DOMAINDECOMP(cr))
404 cg1 = cr->dd->ncg_tot;
405 else
406 cg1 = top->cgs.nr;
407 if (fr->n_tpi > 0)
408 cg1--;
411 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
412 bNS = (flags & GMX_FORCE_NS);
413 bFillGrid = (bNS && bStateChanged);
414 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
415 bDoForces = (flags & GMX_FORCE_FORCES);
417 if (bStateChanged)
419 update_forcerec(fplog,fr,box);
421 /* Calculate total (local) dipole moment in a temporary common array.
422 * This makes it possible to sum them over nodes faster.
424 calc_mu(start,homenr,
425 x,mdatoms->chargeA,mdatoms->chargeB,mdatoms->nChargePerturbed,
426 mu,mu+DIM);
429 if (fr->ePBC != epbcNONE) {
430 /* Compute shift vectors every step,
431 * because of pressure coupling or box deformation!
433 if (DYNAMIC_BOX(*inputrec) && bStateChanged)
434 calc_shifts(box,fr->shift_vec);
436 if (bCalcCGCM) {
437 put_charge_groups_in_box(fplog,cg0,cg1,fr->ePBC,box,
438 &(top->cgs),x,fr->cg_cm);
439 inc_nrnb(nrnb,eNR_CGCM,homenr);
440 inc_nrnb(nrnb,eNR_RESETX,cg1-cg0);
442 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph) {
443 unshift_self(graph,box,x);
446 else if (bCalcCGCM) {
447 calc_cgcm(fplog,cg0,cg1,&(top->cgs),x,fr->cg_cm);
448 inc_nrnb(nrnb,eNR_CGCM,homenr);
451 if (bCalcCGCM) {
452 if (PAR(cr)) {
453 move_cgcm(fplog,cr,fr->cg_cm);
455 if (gmx_debug_at)
456 pr_rvecs(debug,0,"cgcm",fr->cg_cm,top->cgs.nr);
459 #ifdef GMX_MPI
460 if (!(cr->duty & DUTY_PME)) {
461 /* Send particle coordinates to the pme nodes.
462 * Since this is only implemented for domain decomposition
463 * and domain decomposition does not use the graph,
464 * we do not need to worry about shifting.
467 wallcycle_start(wcycle,ewcPP_PMESENDX);
468 GMX_MPE_LOG(ev_send_coordinates_start);
470 bBS = (inputrec->nwall == 2);
471 if (bBS) {
472 copy_mat(box,boxs);
473 svmul(inputrec->wall_ewald_zfac,boxs[ZZ],boxs[ZZ]);
476 gmx_pme_send_x(cr,bBS ? boxs : box,x,mdatoms->nChargePerturbed,lambda,step);
478 GMX_MPE_LOG(ev_send_coordinates_finish);
479 wallcycle_stop(wcycle,ewcPP_PMESENDX);
481 #endif /* GMX_MPI */
483 /* Communicate coordinates and sum dipole if necessary */
484 if (PAR(cr))
486 wallcycle_start(wcycle,ewcMOVEX);
487 if (DOMAINDECOMP(cr))
489 dd_move_x(cr->dd,box,x);
491 else
493 move_x(fplog,cr,GMX_LEFT,GMX_RIGHT,x,nrnb);
495 /* When we don't need the total dipole we sum it in global_stat */
496 if (bStateChanged && NEED_MUTOT(*inputrec))
498 gmx_sumd(2*DIM,mu,cr);
500 wallcycle_stop(wcycle,ewcMOVEX);
502 if (bStateChanged)
504 for(i=0; i<2; i++)
506 for(j=0;j<DIM;j++)
508 fr->mu_tot[i][j] = mu[i*DIM + j];
512 if (fr->efep == efepNO)
514 copy_rvec(fr->mu_tot[0],mu_tot);
516 else
518 for(j=0; j<DIM; j++)
520 mu_tot[j] =
521 (1.0 - lambda)*fr->mu_tot[0][j] + lambda*fr->mu_tot[1][j];
525 /* Reset energies */
526 reset_enerdata(&(inputrec->opts),fr,bNS,enerd,MASTER(cr));
527 if (bNS) {
528 wallcycle_start(wcycle,ewcNS);
530 if (graph && bStateChanged)
531 /* Calculate intramolecular shift vectors to make molecules whole */
532 mk_mshift(fplog,graph,fr->ePBC,box,x);
534 /* Reset long range forces if necessary */
535 if (fr->bTwinRange) {
536 clear_rvecs(fr->natoms_force,fr->f_twin);
537 clear_rvecs(SHIFTS,fr->fshift_twin);
539 /* Do the actual neighbour searching and if twin range electrostatics
540 * also do the calculation of long range forces and energies.
542 dvdl = 0;
543 ns(fplog,fr,x,f,box,groups,&(inputrec->opts),top,mdatoms,
544 cr,nrnb,lambda,&dvdl,&enerd->grpp,bFillGrid,bDoForces);
545 if (bSepDVDL)
546 fprintf(fplog,sepdvdlformat,"LR non-bonded",0,dvdl);
547 enerd->dvdl_lr = dvdl;
549 wallcycle_stop(wcycle,ewcNS);
552 if(inputrec->implicit_solvent && bNS)
554 make_gb_nblist(cr,mtop->natoms,inputrec->gb_algorithm, inputrec->rlist,x,fr,&top->idef,fr->born);
557 if (DOMAINDECOMP(cr)) {
558 if (!(cr->duty & DUTY_PME)) {
559 wallcycle_start(wcycle,ewcPPDURINGPME);
560 dd_force_flop_start(cr->dd,nrnb);
564 /* Start the force cycle counter.
565 * This counter is stopped in do_forcelow_level.
566 * No parallel communication should occur while this counter is running,
567 * since that will interfere with the dynamic load balancing.
569 wallcycle_start(wcycle,ewcFORCE);
571 if (bDoForces) {
572 /* Reset PME/Ewald forces if necessary */
573 if (fr->bF_NoVirSum)
575 if (flags & GMX_FORCE_VIRIAL) {
576 fr->f_novirsum = fr->f_novirsum_alloc;
577 GMX_BARRIER(cr->mpi_comm_mygroup);
578 if (fr->bDomDec)
579 clear_rvecs(fr->f_novirsum_n,fr->f_novirsum);
580 else
581 clear_rvecs(homenr,fr->f_novirsum+start);
582 GMX_BARRIER(cr->mpi_comm_mygroup);
583 } else {
584 /* We are not calculating the pressure so we do not need
585 * a separate array for forces that do not contribute to the pressure.
587 fr->f_novirsum = f;
590 /* Copy long range forces into normal buffers */
591 if (fr->bTwinRange) {
592 for(i=0; i<fr->natoms_force; i++)
593 copy_rvec(fr->f_twin[i],f[i]);
594 for(i=0; i<SHIFTS; i++)
595 copy_rvec(fr->fshift_twin[i],fr->fshift[i]);
597 else {
598 if (DOMAINDECOMP(cr))
599 clear_rvecs(cr->dd->nat_tot,f);
600 else
601 clear_rvecs(mdatoms->nr,f);
602 clear_rvecs(SHIFTS,fr->fshift);
604 clear_rvec(fr->vir_diag_posres);
605 GMX_BARRIER(cr->mpi_comm_mygroup);
607 if (inputrec->ePull == epullCONSTRAINT)
608 clear_pull_forces(inputrec->pull);
610 /* update QMMMrec, if necessary */
611 if(fr->bQMMM)
612 update_QMMMrec(cr,fr,x,mdatoms,box,top);
614 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0) {
615 /* Position restraints always require full pbc */
616 set_pbc(&pbc,inputrec->ePBC,box);
617 v = posres(top->idef.il[F_POSRES].nr,top->idef.il[F_POSRES].iatoms,
618 top->idef.iparams_posres,
619 (const rvec*)x,fr->f_novirsum,fr->vir_diag_posres,
620 inputrec->ePBC==epbcNONE ? NULL : &pbc,lambda,&dvdl,
621 fr->rc_scaling,fr->ePBC,fr->posres_com,fr->posres_comB);
622 if (bSepDVDL) {
623 fprintf(fplog,sepdvdlformat,
624 interaction_function[F_POSRES].longname,v,dvdl);
626 enerd->term[F_POSRES] += v;
627 /* This linear lambda dependence assumption is only correct
628 * when only k depends on lambda,
629 * not when the reference position depends on lambda.
630 * grompp checks for this.
632 enerd->dvdl_lin += dvdl;
633 inc_nrnb(nrnb,eNR_POSRES,top->idef.il[F_POSRES].nr/2);
635 /* Compute the bonded and non-bonded forces */
636 do_force_lowlevel(fplog,step,fr,inputrec,&(top->idef),
637 cr,nrnb,wcycle,mdatoms,&(inputrec->opts),
638 x,hist,f,enerd,fcd,mtop,top,fr->born,
639 &(top->atomtypes),bBornRadii,box,
640 lambda,graph,&(top->excls),fr->mu_tot,
641 flags,&cycles_pme);
643 cycles_force = wallcycle_stop(wcycle,ewcFORCE);
644 GMX_BARRIER(cr->mpi_comm_mygroup);
646 if (ed) {
647 do_flood(fplog,cr,x,f,ed,box,step);
650 if (DOMAINDECOMP(cr)) {
651 dd_force_flop_stop(cr->dd,nrnb);
652 if (wcycle)
653 dd_cycles_add(cr->dd,cycles_force-cycles_pme,ddCyclF);
656 if (bDoForces)
658 /* Compute forces due to electric field */
659 calc_f_el(MASTER(cr) ? field : NULL,
660 start,homenr,mdatoms->chargeA,x,f,
661 inputrec->ex,inputrec->et,t);
663 /* When using PME/Ewald we compute the long range virial there.
664 * otherwise we do it based on long range forces from twin range
665 * cut-off based calculation (or not at all).
668 /* Communicate the forces */
669 if (PAR(cr))
671 wallcycle_start(wcycle,ewcMOVEF);
672 if (DOMAINDECOMP(cr))
674 dd_move_f(cr->dd,f,fr->fshift);
675 /* Position restraint do not introduce inter-cg forces.
676 * When we do not calculate the virial, fr->f_novirsum = f.
678 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
679 (flags & GMX_FORCE_VIRIAL))
681 dd_move_f(cr->dd,fr->f_novirsum,NULL);
684 else
686 pd_move_f(cr,f,nrnb);
688 wallcycle_stop(wcycle,ewcMOVEF);
692 if (bDoForces) {
693 /* If we have NoVirSum forces, but we do not calculate the virial,
694 * we sum fr->f_novirum=f later.
696 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL))) {
697 wallcycle_start(wcycle,ewcVSITESPREAD);
698 spread_vsite_f(fplog,vsite,x,f,fr->fshift,nrnb,
699 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
700 wallcycle_stop(wcycle,ewcVSITESPREAD);
703 if (flags & GMX_FORCE_VIRIAL) {
704 /* Calculation of the virial must be done after vsites! */
705 calc_virial(fplog,mdatoms->start,mdatoms->homenr,x,f,
706 vir_force,graph,box,nrnb,fr,inputrec->ePBC);
710 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F) {
711 /* Calculate the center of mass forces, this requires communication,
712 * which is why pull_potential is called close to other communication.
713 * The virial contribution is calculated directly,
714 * which is why we call pull_potential after calc_virial.
716 set_pbc(&pbc,inputrec->ePBC,box);
717 dvdl = 0;
718 enerd->term[F_COM_PULL] =
719 pull_potential(inputrec->ePull,inputrec->pull,mdatoms,&pbc,
720 cr,t,lambda,x,f,vir_force,&dvdl);
721 if (bSepDVDL)
722 fprintf(fplog,sepdvdlformat,"Com pull",enerd->term[F_COM_PULL],dvdl);
723 enerd->dvdl_lin += dvdl;
726 if (PAR(cr) && !(cr->duty & DUTY_PME))
728 cycles_ppdpme = wallcycle_stop(wcycle,ewcPPDURINGPME);
729 dd_cycles_add(cr->dd,cycles_ppdpme,ddCyclPPduringPME);
731 /* In case of node-splitting, the PP nodes receive the long-range
732 * forces, virial and energy from the PME nodes here.
734 wallcycle_start(wcycle,ewcPP_PMEWAITRECVF);
735 dvdl = 0;
736 gmx_pme_receive_f(cr,fr->f_novirsum,fr->vir_el_recip,&e,&dvdl,
737 &cycles_seppme);
738 if (bSepDVDL)
740 fprintf(fplog,sepdvdlformat,"PME mesh",e,dvdl);
742 enerd->term[F_COUL_RECIP] += e;
743 enerd->dvdl_lin += dvdl;
744 if (wcycle)
746 dd_cycles_add(cr->dd,cycles_seppme,ddCyclPME);
748 wallcycle_stop(wcycle,ewcPP_PMEWAITRECVF);
751 if (bDoForces && fr->bF_NoVirSum) {
752 if (vsite) {
753 /* Spread the mesh force on virtual sites to the other particles...
754 * This is parallellized. MPI communication is performed
755 * if the constructing atoms aren't local.
757 wallcycle_start(wcycle,ewcVSITESPREAD);
758 spread_vsite_f(fplog,vsite,x,fr->f_novirsum,NULL,nrnb,
759 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
760 wallcycle_stop(wcycle,ewcVSITESPREAD);
762 if (flags & GMX_FORCE_VIRIAL) {
763 /* Now add the forces, this is local */
764 if (fr->bDomDec) {
765 sum_forces(0,fr->f_novirsum_n,f,fr->f_novirsum);
766 } else {
767 sum_forces(start,start+homenr,f,fr->f_novirsum);
769 if (EEL_FULL(fr->eeltype)) {
770 /* Add the mesh contribution to the virial */
771 m_add(vir_force,fr->vir_el_recip,vir_force);
773 if (debug)
774 pr_rvecs(debug,0,"vir_force",vir_force,DIM);
778 /* Sum the potential energy terms from group contributions */
779 sum_epot(&(inputrec->opts),enerd);
781 if (fr->print_force >= 0 && bDoForces)
782 print_large_forces(stderr,mdatoms,cr,step,fr->print_force,x,f);
785 void do_constrain_first(FILE *fplog,gmx_constr_t constr,
786 t_inputrec *inputrec,t_mdatoms *md,
787 t_state *state,
788 t_graph *graph,t_commrec *cr,t_nrnb *nrnb,
789 t_forcerec *fr,t_idef *idef)
791 int i,m,start,end;
792 gmx_step_t step;
793 double mass,tmass,vcm[4];
794 real dt=inputrec->delta_t;
795 real dvdlambda;
796 rvec *xcon;
797 char buf[22];
799 start = md->start;
800 end = md->homenr + start;
801 if (debug)
803 fprintf(debug,"vcm: start=%d, homenr=%d, end=%d\n",
804 start,md->homenr,end);
806 snew(xcon,state->nalloc);
808 /* Do a first constraining to reset particles... */
809 step = inputrec->init_step;
810 if (fplog)
812 fprintf(fplog,"\nConstraining the starting coordinates (step %s)\n",
813 gmx_step_str(step,buf));
815 dvdlambda = 0;
816 constrain(NULL,TRUE,FALSE,constr,idef,
817 inputrec,cr,step,0,md,
818 state->x,state->x,NULL,
819 state->box,state->lambda,&dvdlambda,
820 NULL,NULL,nrnb,econqCoord);
822 if (EI_STATE_VELOCITY(inputrec->eI)) {
823 for(i=start; (i<end); i++) {
824 for(m=0; (m<DIM); m++) {
825 /* Reverse the velocity */
826 state->v[i][m] = -state->v[i][m];
827 /* Store the position at t-dt in xcon */
828 xcon[i][m] = state->x[i][m] + dt*state->v[i][m];
832 /* Constrain the positions at t=-dt with the positions at t=0
833 * as reference coordinates.
835 if (fplog)
837 fprintf(fplog,"\nConstraining the coordinates at t0-dt (step %s)\n",
838 gmx_step_str(step,buf));
840 dvdlambda = 0;
841 constrain(NULL,TRUE,FALSE,constr,idef,
842 inputrec,cr,step,-1,md,
843 state->x,xcon,NULL,
844 state->box,state->lambda,&dvdlambda,
845 state->v,NULL,nrnb,econqCoord);
847 sfree(xcon);
849 for(m=0; (m<4); m++)
851 vcm[m] = 0;
853 for(i=start; i<end; i++)
855 mass = md->massT[i];
856 for(m=0; m<DIM; m++)
858 /* Re-reverse the velocities */
859 state->v[i][m] = -state->v[i][m];
860 vcm[m] += state->v[i][m]*mass;
862 vcm[3] += mass;
865 if (inputrec->nstcomm != 0 || debug)
867 /* Compute the global sum of vcm */
868 if (debug)
870 fprintf(debug,
871 "vcm: %8.3f %8.3f %8.3f,"
872 " total mass = %12.5e\n",
873 vcm[XX],vcm[YY],vcm[ZZ],vcm[3]);
875 if (PAR(cr))
877 gmx_sumd(4,vcm,cr);
879 tmass = vcm[3];
880 for(m=0; (m<DIM); m++)
882 vcm[m] /= tmass;
884 if (debug)
886 fprintf(debug,"vcm: %8.3f %8.3f %8.3f,"
887 " total mass = %12.5e\n",vcm[XX],vcm[YY],vcm[ZZ],tmass);
889 if (inputrec->nstcomm != 0)
891 /* Now we have the velocity of center of mass,
892 * let's remove it.
894 for(i=start; (i<end); i++)
896 for(m=0; (m<DIM); m++)
898 state->v[i][m] -= vcm[m];
906 void calc_enervirdiff(FILE *fplog,int eDispCorr,t_forcerec *fr)
908 double eners[2],virs[2],enersum,virsum,y0,f,g,h;
909 double r0,r1,r,rc3,rc9,ea,eb,ec,pa,pb,pc,pd;
910 double invscale,invscale2,invscale3;
911 int ri0,ri1,ri,i,offstart,offset;
912 real scale,*vdwtab;
914 fr->enershiftsix = 0;
915 fr->enershifttwelve = 0;
916 fr->enerdiffsix = 0;
917 fr->enerdifftwelve = 0;
918 fr->virdiffsix = 0;
919 fr->virdifftwelve = 0;
921 if (eDispCorr != edispcNO) {
922 for(i=0; i<2; i++) {
923 eners[i] = 0;
924 virs[i] = 0;
926 if ((fr->vdwtype == evdwSWITCH) || (fr->vdwtype == evdwSHIFT)) {
927 if (fr->rvdw_switch == 0)
928 gmx_fatal(FARGS,
929 "With dispersion correction rvdw-switch can not be zero "
930 "for vdw-type = %s",evdw_names[fr->vdwtype]);
932 scale = fr->nblists[0].tab.scale;
933 vdwtab = fr->nblists[0].vdwtab;
935 /* Round the cut-offs to exact table values for precision */
936 ri0 = floor(fr->rvdw_switch*scale);
937 ri1 = ceil(fr->rvdw*scale);
938 r0 = ri0/scale;
939 r1 = ri1/scale;
940 rc3 = r0*r0*r0;
941 rc9 = rc3*rc3*rc3;
943 if (fr->vdwtype == evdwSHIFT) {
944 /* Determine the constant energy shift below rvdw_switch */
945 fr->enershiftsix = (real)(-1.0/(rc3*rc3)) - vdwtab[8*ri0];
946 fr->enershifttwelve = (real)( 1.0/(rc9*rc3)) - vdwtab[8*ri0 + 4];
948 /* Add the constant part from 0 to rvdw_switch.
949 * This integration from 0 to rvdw_switch overcounts the number
950 * of interactions by 1, as it also counts the self interaction.
951 * We will correct for this later.
953 eners[0] += 4.0*M_PI*fr->enershiftsix*rc3/3.0;
954 eners[1] += 4.0*M_PI*fr->enershifttwelve*rc3/3.0;
956 invscale = 1.0/(scale);
957 invscale2 = invscale*invscale;
958 invscale3 = invscale*invscale2;
960 /* following summation derived from cubic spline definition,
961 Numerical Recipies in C, second edition, p. 113-116. Exact
962 for the cubic spline. We first calculate the negative of
963 the energy from rvdw to rvdw_switch, assuming that g(r)=1,
964 and then add the more standard, abrupt cutoff correction to
965 that result, yielding the long-range correction for a
966 switched function. We perform both the pressure and energy
967 loops at the same time for simplicity, as the computational
968 cost is low. */
970 for (i=0;i<2;i++) {
971 enersum = 0.0; virsum = 0.0;
972 if (i==0)
973 offstart = 0;
974 else
975 offstart = 4;
976 for (ri=ri0; ri<ri1; ri++) {
977 r = ri*invscale;
978 ea = invscale3;
979 eb = 2.0*invscale2*r;
980 ec = invscale*r*r;
982 pa = invscale3;
983 pb = 3.0*invscale2*r;
984 pc = 3.0*invscale*r*r;
985 pd = r*r*r;
987 /* this "8" is from the packing in the vdwtab array - perhaps
988 should be #define'ed? */
989 offset = 8*ri + offstart;
990 y0 = vdwtab[offset];
991 f = vdwtab[offset+1];
992 g = vdwtab[offset+2];
993 h = vdwtab[offset+3];
995 enersum += y0*(ea/3 + eb/2 + ec) + f*(ea/4 + eb/3 + ec/2)+
996 g*(ea/5 + eb/4 + ec/3) + h*(ea/6 + eb/5 + ec/4);
997 virsum += f*(pa/4 + pb/3 + pc/2 + pd) +
998 2*g*(pa/5 + pb/4 + pc/3 + pd/2) + 3*h*(pa/6 + pb/5 + pc/4 + pd/3);
1001 enersum *= 4.0*M_PI;
1002 virsum *= 4.0*M_PI;
1003 eners[i] -= enersum;
1004 virs[i] -= virsum;
1007 /* now add the correction for rvdw_switch to infinity */
1008 eners[0] += -4.0*M_PI/(3.0*rc3);
1009 eners[1] += 4.0*M_PI/(9.0*rc9);
1010 virs[0] += 8.0*M_PI/rc3;
1011 virs[1] += -16.0*M_PI/(3.0*rc9);
1013 else if ((fr->vdwtype == evdwCUT) || (fr->vdwtype == evdwUSER)) {
1014 if (fr->vdwtype == evdwUSER && fplog)
1015 fprintf(fplog,
1016 "WARNING: using dispersion correction with user tables\n");
1017 rc3 = fr->rvdw*fr->rvdw*fr->rvdw;
1018 rc9 = rc3*rc3*rc3;
1019 eners[0] += -4.0*M_PI/(3.0*rc3);
1020 eners[1] += 4.0*M_PI/(9.0*rc9);
1021 virs[0] += 8.0*M_PI/rc3;
1022 virs[1] += -16.0*M_PI/(3.0*rc9);
1023 } else {
1024 gmx_fatal(FARGS,
1025 "Dispersion correction is not implemented for vdw-type = %s",
1026 evdw_names[fr->vdwtype]);
1028 fr->enerdiffsix = eners[0];
1029 fr->enerdifftwelve = eners[1];
1030 /* The 0.5 is due to the Gromacs definition of the virial */
1031 fr->virdiffsix = 0.5*virs[0];
1032 fr->virdifftwelve = 0.5*virs[1];
1036 void calc_dispcorr(FILE *fplog,t_inputrec *ir,t_forcerec *fr,
1037 gmx_step_t step,int natoms,matrix box,real lambda,
1038 tensor pres,tensor virial,gmx_enerdata_t *enerd)
1040 bool bCorrAll,bCorrPres;
1041 real dvdlambda,invvol,dens,ninter,avcsix,avctwelve,enerdiff,svir=0,spres=0;
1042 int m;
1044 enerd->term[F_DISPCORR] = 0.0;
1045 enerd->term[F_PDISPCORR] = 0.0;
1047 if (ir->eDispCorr != edispcNO)
1049 bCorrAll = (ir->eDispCorr == edispcAllEner ||
1050 ir->eDispCorr == edispcAllEnerPres);
1051 bCorrPres = (ir->eDispCorr == edispcEnerPres ||
1052 ir->eDispCorr == edispcAllEnerPres);
1054 invvol = 1/det(box);
1055 if (fr->n_tpi)
1057 /* Only correct for the interactions with the inserted molecule */
1058 dens = (natoms - fr->n_tpi)*invvol;
1059 ninter = fr->n_tpi;
1061 else
1063 dens = natoms*invvol;
1064 ninter = 0.5*natoms;
1067 if (ir->efep == efepNO)
1069 avcsix = fr->avcsix[0];
1070 avctwelve = fr->avctwelve[0];
1072 else
1074 avcsix = (1 - lambda)*fr->avcsix[0] + lambda*fr->avcsix[1];
1075 avctwelve = (1 - lambda)*fr->avctwelve[0] + lambda*fr->avctwelve[1];
1078 enerdiff = ninter*(dens*fr->enerdiffsix - fr->enershiftsix);
1079 enerd->term[F_DISPCORR] += avcsix*enerdiff;
1080 dvdlambda = 0.0;
1081 if (ir->efep != efepNO)
1083 dvdlambda += (fr->avcsix[1] - fr->avcsix[0])*enerdiff;
1086 if (bCorrAll)
1088 enerdiff = ninter*(dens*fr->enerdifftwelve - fr->enershifttwelve);
1089 enerd->term[F_DISPCORR] += avctwelve*enerdiff;
1090 if (fr->efep != efepNO)
1092 dvdlambda += (fr->avctwelve[1] - fr->avctwelve[0])*enerdiff;
1096 if (bCorrPres)
1098 svir = ninter*dens*avcsix*fr->virdiffsix/3.0;
1099 if (ir->eDispCorr == edispcAllEnerPres)
1101 svir += ninter*dens*avctwelve*fr->virdifftwelve/3.0;
1103 /* The factor 2 is because of the Gromacs virial definition */
1104 spres = -2.0*invvol*svir*PRESFAC;
1106 for(m=0; m<DIM; m++) {
1107 virial[m][m] += svir;
1108 pres[m][m] += spres;
1110 enerd->term[F_PDISPCORR] = spres;
1111 enerd->term[F_PRES] += spres;
1114 if (fr->bSepDVDL && do_per_step(step,ir->nstlog))
1116 fprintf(fplog,sepdvdlformat,"Dispersion correction",
1117 enerd->term[F_DISPCORR],dvdlambda);
1120 enerd->term[F_EPOT] += enerd->term[F_DISPCORR];
1121 if (fr->efep != efepNO)
1123 enerd->dvdl_lin += dvdlambda;
1129 void do_pbc_first(FILE *fplog,matrix box,t_forcerec *fr,
1130 t_graph *graph,rvec x[])
1132 if (fplog)
1133 fprintf(fplog,"Removing pbc first time\n");
1134 calc_shifts(box,fr->shift_vec);
1135 if (graph) {
1136 mk_mshift(fplog,graph,fr->ePBC,box,x);
1137 if (gmx_debug_at)
1138 p_graph(debug,"do_pbc_first 1",graph);
1139 shift_self(graph,box,x);
1140 /* By doing an extra mk_mshift the molecules that are broken
1141 * because they were e.g. imported from another software
1142 * will be made whole again. Such are the healing powers
1143 * of GROMACS.
1145 mk_mshift(fplog,graph,fr->ePBC,box,x);
1146 if (gmx_debug_at)
1147 p_graph(debug,"do_pbc_first 2",graph);
1149 if (fplog)
1150 fprintf(fplog,"Done rmpbc\n");
1153 static void low_do_pbc_mtop(FILE *fplog,int ePBC,matrix box,
1154 gmx_mtop_t *mtop,rvec x[],
1155 bool bFirst)
1157 t_graph *graph;
1158 int mb,as,mol;
1159 gmx_molblock_t *molb;
1161 if (bFirst && fplog)
1162 fprintf(fplog,"Removing pbc first time\n");
1164 snew(graph,1);
1165 as = 0;
1166 for(mb=0; mb<mtop->nmolblock; mb++) {
1167 molb = &mtop->molblock[mb];
1168 if (molb->natoms_mol == 1 ||
1169 (!bFirst && mtop->moltype[molb->type].cgs.nr == 1)) {
1170 /* Just one atom or charge group in the molecule, no PBC required */
1171 as += molb->nmol*molb->natoms_mol;
1172 } else {
1173 /* Pass NULL iso fplog to avoid graph prints for each molecule type */
1174 mk_graph_ilist(NULL,mtop->moltype[molb->type].ilist,
1175 0,molb->natoms_mol,FALSE,FALSE,graph);
1177 for(mol=0; mol<molb->nmol; mol++) {
1178 mk_mshift(fplog,graph,ePBC,box,x+as);
1180 shift_self(graph,box,x+as);
1181 /* The molecule is whole now.
1182 * We don't need the second mk_mshift call as in do_pbc_first,
1183 * since we no longer need this graph.
1186 as += molb->natoms_mol;
1188 done_graph(graph);
1191 sfree(graph);
1194 void do_pbc_first_mtop(FILE *fplog,int ePBC,matrix box,
1195 gmx_mtop_t *mtop,rvec x[])
1197 low_do_pbc_mtop(fplog,ePBC,box,mtop,x,TRUE);
1200 void do_pbc_mtop(FILE *fplog,int ePBC,matrix box,
1201 gmx_mtop_t *mtop,rvec x[])
1203 low_do_pbc_mtop(fplog,ePBC,box,mtop,x,FALSE);
1206 void finish_run(FILE *fplog,t_commrec *cr,char *confout,
1207 t_inputrec *inputrec,
1208 t_nrnb nrnb[],gmx_wallcycle_t wcycle,
1209 gmx_runtime_t *runtime,
1210 bool bWriteStat)
1212 int i,j;
1213 t_nrnb *nrnb_tot=NULL;
1214 real delta_t;
1215 double nbfs,mflop;
1216 double cycles[ewcNR];
1218 wallcycle_sum(cr,wcycle,cycles);
1220 if (cr->nnodes > 1) {
1221 if (SIMMASTER(cr))
1222 snew(nrnb_tot,1);
1223 #ifdef GMX_MPI
1224 MPI_Reduce(nrnb->n,nrnb_tot->n,eNRNB,MPI_DOUBLE,MPI_SUM,
1225 MASTERRANK(cr),cr->mpi_comm_mysim);
1226 #endif
1227 } else {
1228 nrnb_tot = nrnb;
1231 if (SIMMASTER(cr)) {
1232 print_flop(fplog,nrnb_tot,&nbfs,&mflop);
1233 if (cr->nnodes > 1) {
1234 sfree(nrnb_tot);
1238 if ((cr->duty & DUTY_PP) && DOMAINDECOMP(cr)) {
1239 print_dd_statistics(cr,inputrec,fplog);
1242 if (SIMMASTER(cr)) {
1243 if (PARTDECOMP(cr)) {
1244 pr_load(fplog,cr,nrnb_tot);
1247 wallcycle_print(fplog,cr->nnodes,cr->npmenodes,runtime->realtime,
1248 wcycle,cycles);
1250 if (EI_DYNAMICS(inputrec->eI)) {
1251 delta_t = inputrec->delta_t;
1252 } else {
1253 delta_t = 0;
1256 if (fplog) {
1257 print_perf(fplog,runtime->proctime,runtime->realtime,
1258 cr->nnodes-cr->npmenodes,
1259 runtime->nsteps_done,delta_t,nbfs,mflop);
1261 if (bWriteStat) {
1262 print_perf(stderr,runtime->proctime,runtime->realtime,
1263 cr->nnodes-cr->npmenodes,
1264 runtime->nsteps_done,delta_t,nbfs,mflop);
1268 runtime=inputrec->nsteps*inputrec->delta_t;
1269 if (bWriteStat) {
1270 if (cr->nnodes == 1)
1271 fprintf(stderr,"\n\n");
1272 print_perf(stderr,nodetime,realtime,runtime,&ntot,
1273 cr->nnodes-cr->npmenodes,FALSE);
1275 wallcycle_print(fplog,cr->nnodes,cr->npmenodes,realtime,wcycle,cycles);
1276 print_perf(fplog,nodetime,realtime,runtime,&ntot,cr->nnodes-cr->npmenodes,
1277 TRUE);
1278 if (PARTDECOMP(cr))
1279 pr_load(fplog,cr,nrnb_all);
1280 if (cr->nnodes > 1)
1281 sfree(nrnb_all);
1286 void init_md(FILE *fplog,
1287 t_commrec *cr,t_inputrec *ir,
1288 double *t,double *t0,
1289 real *lambda,double *lam0,
1290 t_nrnb *nrnb,gmx_mtop_t *mtop,
1291 gmx_update_t *upd,
1292 int nfile,t_filenm fnm[],
1293 int *fp_trn,int *fp_xtc,int *fp_ene,char **fn_cpt,
1294 FILE **fp_dhdl,FILE **fp_field,
1295 t_mdebin **mdebin,
1296 tensor force_vir,tensor shake_vir,rvec mu_tot,
1297 bool *bNEMD,bool *bSimAnn,t_vcm **vcm, unsigned long Flags)
1299 int i,j,n;
1300 real tmpt,mod;
1301 char filemode[2];
1303 sprintf(filemode, (Flags & MD_APPENDFILES) ? "a" : "w");
1305 /* Initial values */
1306 *t = *t0 = ir->init_t;
1307 if (ir->efep != efepNO)
1309 *lam0 = ir->init_lambda;
1310 *lambda = *lam0 + ir->init_step*ir->delta_lambda;
1312 else
1314 *lambda = *lam0 = 0.0;
1317 *bSimAnn=FALSE;
1318 for(i=0;i<ir->opts.ngtc;i++)
1320 /* set bSimAnn if any group is being annealed */
1321 if(ir->opts.annealing[i]!=eannNO)
1323 *bSimAnn = TRUE;
1326 if (*bSimAnn)
1328 update_annealing_target_temp(&(ir->opts),ir->init_t);
1331 *bNEMD = (ir->opts.ngacc > 1) || (norm(ir->opts.acc[0]) > 0);
1333 if (upd)
1335 *upd = init_update(fplog,ir);
1338 if (vcm != NULL)
1340 *vcm = init_vcm(fplog,&mtop->groups,ir);
1343 if (EI_DYNAMICS(ir->eI) && !(Flags & MD_APPENDFILES))
1345 if (ir->etc == etcBERENDSEN)
1347 please_cite(fplog,"Berendsen84a");
1349 if (ir->etc == etcVRESCALE)
1351 please_cite(fplog,"Bussi2007a");
1355 init_nrnb(nrnb);
1357 if (nfile != -1)
1359 *fp_trn = -1;
1360 *fp_ene = -1;
1361 *fp_xtc = -1;
1363 if (MASTER(cr))
1365 *fp_trn = open_trn(ftp2fn(efTRN,nfile,fnm), filemode);
1366 if (ir->nstxtcout > 0)
1368 *fp_xtc = open_xtc(ftp2fn(efXTC,nfile,fnm), filemode);
1370 *fp_ene = open_enx(ftp2fn(efEDR,nfile,fnm), filemode);
1371 *fn_cpt = opt2fn("-cpo",nfile,fnm);
1373 if ((fp_dhdl != NULL) && ir->efep != efepNO && ir->nstdhdl > 0)
1375 if(Flags & MD_APPENDFILES)
1377 *fp_dhdl= gmx_fio_fopen(opt2fn("-dhdl",nfile,fnm),filemode);
1379 else
1381 *fp_dhdl = open_dhdl(opt2fn("-dhdl",nfile,fnm),ir);
1385 if ((fp_field != NULL) &&
1386 (ir->ex[XX].n || ir->ex[YY].n ||ir->ex[ZZ].n))
1388 if(Flags & MD_APPENDFILES)
1390 *fp_dhdl=gmx_fio_fopen(opt2fn("-field",nfile,fnm),filemode);
1392 else
1394 *fp_field = xvgropen(opt2fn("-field",nfile,fnm),
1395 "Applied electric field","Time (ps)",
1396 "E (V/nm)");
1400 *mdebin = init_mdebin( (Flags & MD_APPENDFILES) ? -1 : *fp_ene,mtop,ir);
1403 /* Initiate variables */
1404 clear_mat(force_vir);
1405 clear_mat(shake_vir);
1406 clear_rvec(mu_tot);
1408 debug_gmx();