Merge branch 'rotation-4-5' into rotation
[gromacs/adressmacs.git] / src / mdlib / sim_util.c
blobb176ccf6afaa0d008203e1317e52ad68b1e1cd7f
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 "chargegroup.h"
61 #include "vec.h"
62 #include "time.h"
63 #include "nrnb.h"
64 #include "mshift.h"
65 #include "mdrun.h"
66 #include "update.h"
67 #include "physics.h"
68 #include "main.h"
69 #include "mdatoms.h"
70 #include "force.h"
71 #include "bondf.h"
72 #include "pme.h"
73 #include "pppm.h"
74 #include "disre.h"
75 #include "orires.h"
76 #include "network.h"
77 #include "calcmu.h"
78 #include "constr.h"
79 #include "xvgr.h"
80 #include "trnio.h"
81 #include "xtcio.h"
82 #include "copyrite.h"
83 #include "pull_rotation.h"
84 #include "mpelogging.h"
85 #include "domdec.h"
86 #include "partdec.h"
87 #include "gmx_wallcycle.h"
88 #include "genborn.h"
90 #ifdef GMX_LIB_MPI
91 #include <mpi.h>
92 #endif
93 #ifdef GMX_THREADS
94 #include "tmpi.h"
95 #endif
97 #include "qmmm.h"
99 #if 0
100 typedef struct gmx_timeprint {
102 } t_gmx_timeprint;
103 #endif
105 /* Portable version of ctime_r implemented in src/gmxlib/string2.c, but we do not want it declared in public installed headers */
106 char *
107 gmx_ctime_r(const time_t *clock,char *buf, int n);
110 double
111 gmx_gettime()
113 #ifdef HAVE_GETTIMEOFDAY
114 struct timeval t;
115 struct timezone tz = { 0,0 };
116 double seconds;
118 gettimeofday(&t,&tz);
120 seconds = (double) t.tv_sec + 1e-6*(double)t.tv_usec;
122 return seconds;
123 #else
124 double seconds;
126 seconds = time(NULL);
128 return seconds;
129 #endif
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)
138 time_t finish;
139 char timebuf[STRLEN];
140 double dt;
141 char buf[48];
143 #ifndef GMX_THREADS
144 if (!PAR(cr))
145 #endif
147 fprintf(out,"\r");
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;
161 if (ir->nsteps >= 0)
163 if (dt >= 300)
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);
171 else
172 fprintf(out,", remaining runtime: %5d s ",(int)dt);
174 else
176 fprintf(out," performance: %.1f ns/day ",
177 ir->delta_t/1000*24*60*60/runtime->time_per_step);
180 #ifndef GMX_THREADS
181 if (PAR(cr))
183 fprintf(out,"\n");
185 #endif
187 fflush(out);
190 #ifdef NO_CLOCK
191 #define clock() -1
192 #endif
194 static double set_proctime(gmx_runtime_t *runtime)
196 double diff;
197 #ifdef GMX_CRAY_XT3
198 double prev;
200 prev = runtime->proc;
201 runtime->proc = dclock();
203 diff = runtime->proc - prev;
204 #else
205 clock_t prev;
207 prev = runtime->proc;
208 runtime->proc = clock();
210 diff = (double)(runtime->proc - prev)/(double)CLOCKS_PER_SEC;
211 #endif
212 if (diff < 0)
214 /* The counter has probably looped, ignore this data */
215 diff = 0;
218 return diff;
221 void runtime_start(gmx_runtime_t *runtime)
223 runtime->real = gmx_gettime();
224 runtime->proc = 0;
225 set_proctime(runtime);
226 runtime->realtime = 0;
227 runtime->proctime = 0;
228 runtime->last = 0;
229 runtime->time_per_step = 0;
232 void runtime_end(gmx_runtime_t *runtime)
234 double now;
236 now = gmx_gettime();
238 runtime->proctime += set_proctime(runtime);
239 runtime->realtime = now - runtime->real;
240 runtime->real = now;
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)
251 int i;
252 char timebuf[STRLEN];
253 char time_string[STRLEN];
254 time_t tmptime;
256 if (fplog)
258 if (runtime != NULL)
260 tmptime = (time_t) runtime->real;
261 gmx_ctime_r(&tmptime,timebuf,STRLEN);
263 else
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];
272 time_string[i]='\0';
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[])
280 int i;
282 if (gmx_debug_at) {
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
300 * now.
301 * The function should return the energy due to the electric field
302 * (if any) but for now returns 0.
304 * WARNING:
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)
317 rvec Ext;
318 real t0;
319 int i,m;
321 for(m=0; (m<DIM); m++)
323 if (Et[m].n > 0)
325 if (Et[m].n == 3)
327 t0 = Et[m].a[1];
328 Ext[m] = cos(Et[m].a[0]*(t-t0))*exp(-sqr(t-t0)/(2.0*sqr(Et[m].a[2])));
330 else
332 Ext[m] = cos(Et[m].a[0]*t);
335 else
337 Ext[m] = 1.0;
339 if (Ex[m].n > 0)
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];
346 else
348 Ext[m] = 0;
351 if (fp != NULL)
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)
362 int i,j;
363 tensor virtest;
365 /* The short-range virial from surrounding boxes */
366 clear_mat(vir_part);
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];
386 if (debug)
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)
393 int i;
394 real pf2,fn2;
395 char buf[STEPSTRSIZE];
397 pf2 = sqr(pforce);
398 for(i=md->start; i<md->start+md->homenr; i++) {
399 fn2 = norm2(f[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,
412 gmx_localtop_t *top,
413 gmx_mtop_t *mtop,
414 gmx_groups_t *groups,
415 matrix box,rvec x[],history_t *hist,
416 rvec f[],
417 tensor vir_force,
418 t_mdatoms *mdatoms,
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,
423 gmx_bool bBornRadii,
424 int flags)
426 int cg0,cg1,i,j;
427 int start,homenr;
428 double mu[2*DIM];
429 gmx_bool bSepDVDL,bStateChanged,bNS,bFillGrid,bCalcCGCM,bBS;
430 gmx_bool bDoLongRange,bDoForces,bSepLRF;
431 matrix boxs;
432 real e,v,dvdl;
433 t_pbc pbc;
434 float cycles_ppdpme,cycles_pme,cycles_seppme,cycles_force;
436 start = mdatoms->start;
437 homenr = mdatoms->homenr;
439 bSepDVDL = (fr->bSepDVDL && do_per_step(step,inputrec->nstlog));
441 clear_mat(vir_force);
443 if (PARTDECOMP(cr))
445 pd_cg_range(cr,&cg0,&cg1);
447 else
449 cg0 = 0;
450 if (DOMAINDECOMP(cr))
452 cg1 = cr->dd->ncg_tot;
454 else
456 cg1 = top->cgs.nr;
458 if (fr->n_tpi > 0)
460 cg1--;
464 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
465 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll==FALSE);
466 bFillGrid = (bNS && bStateChanged);
467 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
468 bDoLongRange = (fr->bTwinRange && bNS && (flags & GMX_FORCE_DOLR));
469 bDoForces = (flags & GMX_FORCE_FORCES);
470 bSepLRF = (bDoLongRange && bDoForces && (flags & GMX_FORCE_SEPLRF));
472 if (bStateChanged)
474 update_forcerec(fplog,fr,box);
476 /* Calculate total (local) dipole moment in a temporary common array.
477 * This makes it possible to sum them over nodes faster.
479 calc_mu(start,homenr,
480 x,mdatoms->chargeA,mdatoms->chargeB,mdatoms->nChargePerturbed,
481 mu,mu+DIM);
484 if (fr->ePBC != epbcNONE) {
485 /* Compute shift vectors every step,
486 * because of pressure coupling or box deformation!
488 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
489 calc_shifts(box,fr->shift_vec);
491 if (bCalcCGCM) {
492 put_charge_groups_in_box(fplog,cg0,cg1,fr->ePBC,box,
493 &(top->cgs),x,fr->cg_cm);
494 inc_nrnb(nrnb,eNR_CGCM,homenr);
495 inc_nrnb(nrnb,eNR_RESETX,cg1-cg0);
497 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph) {
498 unshift_self(graph,box,x);
501 else if (bCalcCGCM) {
502 calc_cgcm(fplog,cg0,cg1,&(top->cgs),x,fr->cg_cm);
503 inc_nrnb(nrnb,eNR_CGCM,homenr);
506 if (bCalcCGCM) {
507 if (PAR(cr)) {
508 move_cgcm(fplog,cr,fr->cg_cm);
510 if (gmx_debug_at)
511 pr_rvecs(debug,0,"cgcm",fr->cg_cm,top->cgs.nr);
514 #ifdef GMX_MPI
515 if (!(cr->duty & DUTY_PME)) {
516 /* Send particle coordinates to the pme nodes.
517 * Since this is only implemented for domain decomposition
518 * and domain decomposition does not use the graph,
519 * we do not need to worry about shifting.
522 wallcycle_start(wcycle,ewcPP_PMESENDX);
523 GMX_MPE_LOG(ev_send_coordinates_start);
525 bBS = (inputrec->nwall == 2);
526 if (bBS) {
527 copy_mat(box,boxs);
528 svmul(inputrec->wall_ewald_zfac,boxs[ZZ],boxs[ZZ]);
531 gmx_pme_send_x(cr,bBS ? boxs : box,x,
532 mdatoms->nChargePerturbed,lambda,
533 ( flags & GMX_FORCE_VIRIAL),step);
535 GMX_MPE_LOG(ev_send_coordinates_finish);
536 wallcycle_stop(wcycle,ewcPP_PMESENDX);
538 #endif /* GMX_MPI */
540 /* Communicate coordinates and sum dipole if necessary */
541 if (PAR(cr))
543 wallcycle_start(wcycle,ewcMOVEX);
544 if (DOMAINDECOMP(cr))
546 dd_move_x(cr->dd,box,x);
548 else
550 move_x(fplog,cr,GMX_LEFT,GMX_RIGHT,x,nrnb);
552 /* When we don't need the total dipole we sum it in global_stat */
553 if (bStateChanged && NEED_MUTOT(*inputrec))
555 gmx_sumd(2*DIM,mu,cr);
557 wallcycle_stop(wcycle,ewcMOVEX);
559 if (bStateChanged)
561 for(i=0; i<2; i++)
563 for(j=0;j<DIM;j++)
565 fr->mu_tot[i][j] = mu[i*DIM + j];
569 if (fr->efep == efepNO)
571 copy_rvec(fr->mu_tot[0],mu_tot);
573 else
575 for(j=0; j<DIM; j++)
577 mu_tot[j] =
578 (1.0 - lambda)*fr->mu_tot[0][j] + lambda*fr->mu_tot[1][j];
582 /* Reset energies */
583 reset_enerdata(&(inputrec->opts),fr,bNS,enerd,MASTER(cr));
584 clear_rvecs(SHIFTS,fr->fshift);
586 if (bNS)
588 wallcycle_start(wcycle,ewcNS);
590 if (graph && bStateChanged)
592 /* Calculate intramolecular shift vectors to make molecules whole */
593 mk_mshift(fplog,graph,fr->ePBC,box,x);
596 /* Reset long range forces if necessary */
597 if (fr->bTwinRange)
599 /* Reset the (long-range) forces if necessary */
600 clear_rvecs(fr->natoms_force_constr,bSepLRF ? fr->f_twin : f);
603 /* Do the actual neighbour searching and if twin range electrostatics
604 * also do the calculation of long range forces and energies.
606 dvdl = 0;
607 ns(fplog,fr,x,box,
608 groups,&(inputrec->opts),top,mdatoms,
609 cr,nrnb,lambda,&dvdl,&enerd->grpp,bFillGrid,
610 bDoLongRange,bDoForces,bSepLRF ? fr->f_twin : f);
611 if (bSepDVDL)
613 fprintf(fplog,sepdvdlformat,"LR non-bonded",0.0,dvdl);
615 enerd->dvdl_lin += dvdl;
617 wallcycle_stop(wcycle,ewcNS);
620 if (inputrec->implicit_solvent && bNS)
622 make_gb_nblist(cr,inputrec->gb_algorithm,inputrec->rlist,
623 x,box,fr,&top->idef,graph,fr->born);
626 if (DOMAINDECOMP(cr))
628 if (!(cr->duty & DUTY_PME))
630 wallcycle_start(wcycle,ewcPPDURINGPME);
631 dd_force_flop_start(cr->dd,nrnb);
635 if (inputrec->bRot)
637 /* Enforced rotation has its own cycle counter that starts after the collective
638 * coordinates have been communicated. It is added to ddCyclF */
639 do_rotation(cr,inputrec,box,x,t,step,wcycle,bNS);
642 /* Start the force cycle counter.
643 * This counter is stopped in do_forcelow_level.
644 * No parallel communication should occur while this counter is running,
645 * since that will interfere with the dynamic load balancing.
647 wallcycle_start(wcycle,ewcFORCE);
648 GMX_MPE_LOG(ev_forcecycles_start);
650 if (bDoForces)
652 /* Reset forces for which the virial is calculated separately:
653 * PME/Ewald forces if necessary */
654 if (fr->bF_NoVirSum)
656 if (flags & GMX_FORCE_VIRIAL)
658 fr->f_novirsum = fr->f_novirsum_alloc;
659 GMX_BARRIER(cr->mpi_comm_mygroup);
660 if (fr->bDomDec)
662 clear_rvecs(fr->f_novirsum_n,fr->f_novirsum);
664 else
666 clear_rvecs(homenr,fr->f_novirsum+start);
668 GMX_BARRIER(cr->mpi_comm_mygroup);
670 else
672 /* We are not calculating the pressure so we do not need
673 * a separate array for forces that do not contribute
674 * to the pressure.
676 fr->f_novirsum = f;
680 if (bSepLRF)
682 /* Add the long range forces to the short range forces */
683 for(i=0; i<fr->natoms_force_constr; i++)
685 copy_rvec(fr->f_twin[i],f[i]);
688 else if (!(fr->bTwinRange && bNS))
690 /* Clear the short-range forces */
691 clear_rvecs(fr->natoms_force_constr,f);
694 clear_rvec(fr->vir_diag_posres);
696 GMX_BARRIER(cr->mpi_comm_mygroup);
698 if (inputrec->ePull == epullCONSTRAINT)
700 clear_pull_forces(inputrec->pull);
703 /* update QMMMrec, if necessary */
704 if(fr->bQMMM)
706 update_QMMMrec(cr,fr,x,mdatoms,box,top);
709 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
711 /* Position restraints always require full pbc */
712 set_pbc(&pbc,inputrec->ePBC,box);
713 v = posres(top->idef.il[F_POSRES].nr,top->idef.il[F_POSRES].iatoms,
714 top->idef.iparams_posres,
715 (const rvec*)x,fr->f_novirsum,fr->vir_diag_posres,
716 inputrec->ePBC==epbcNONE ? NULL : &pbc,lambda,&dvdl,
717 fr->rc_scaling,fr->ePBC,fr->posres_com,fr->posres_comB);
718 if (bSepDVDL)
720 fprintf(fplog,sepdvdlformat,
721 interaction_function[F_POSRES].longname,v,dvdl);
723 enerd->term[F_POSRES] += v;
724 /* This linear lambda dependence assumption is only correct
725 * when only k depends on lambda,
726 * not when the reference position depends on lambda.
727 * grompp checks for this.
729 enerd->dvdl_lin += dvdl;
730 inc_nrnb(nrnb,eNR_POSRES,top->idef.il[F_POSRES].nr/2);
733 /* Compute the bonded and non-bonded energies and optionally forces */
734 do_force_lowlevel(fplog,step,fr,inputrec,&(top->idef),
735 cr,nrnb,wcycle,mdatoms,&(inputrec->opts),
736 x,hist,f,enerd,fcd,mtop,top,fr->born,
737 &(top->atomtypes),bBornRadii,box,
738 lambda,graph,&(top->excls),fr->mu_tot,
739 flags,&cycles_pme);
741 cycles_force = wallcycle_stop(wcycle,ewcFORCE);
742 GMX_BARRIER(cr->mpi_comm_mygroup);
744 if (ed)
746 do_flood(fplog,cr,x,f,ed,box,step);
749 if (DOMAINDECOMP(cr))
751 dd_force_flop_stop(cr->dd,nrnb);
752 if (wcycle)
754 dd_cycles_add(cr->dd,cycles_force-cycles_pme,ddCyclF);
758 if (bDoForces)
760 if (IR_ELEC_FIELD(*inputrec))
762 /* Compute forces due to electric field */
763 calc_f_el(MASTER(cr) ? field : NULL,
764 start,homenr,mdatoms->chargeA,x,fr->f_novirsum,
765 inputrec->ex,inputrec->et,t);
768 /* Communicate the forces */
769 if (PAR(cr))
771 wallcycle_start(wcycle,ewcMOVEF);
772 if (DOMAINDECOMP(cr))
774 dd_move_f(cr->dd,f,fr->fshift);
775 /* Do we need to communicate the separate force array
776 * for terms that do not contribute to the single sum virial?
777 * Position restraints and electric fields do not introduce
778 * inter-cg forces, only full electrostatics methods do.
779 * When we do not calculate the virial, fr->f_novirsum = f,
780 * so we have already communicated these forces.
782 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
783 (flags & GMX_FORCE_VIRIAL))
785 dd_move_f(cr->dd,fr->f_novirsum,NULL);
787 if (bSepLRF)
789 /* We should not update the shift forces here,
790 * since f_twin is already included in f.
792 dd_move_f(cr->dd,fr->f_twin,NULL);
795 else
797 pd_move_f(cr,f,nrnb);
798 if (bSepLRF)
800 pd_move_f(cr,fr->f_twin,nrnb);
803 wallcycle_stop(wcycle,ewcMOVEF);
806 /* If we have NoVirSum forces, but we do not calculate the virial,
807 * we sum fr->f_novirum=f later.
809 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
811 wallcycle_start(wcycle,ewcVSITESPREAD);
812 spread_vsite_f(fplog,vsite,x,f,fr->fshift,nrnb,
813 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
814 wallcycle_stop(wcycle,ewcVSITESPREAD);
816 if (bSepLRF)
818 wallcycle_start(wcycle,ewcVSITESPREAD);
819 spread_vsite_f(fplog,vsite,x,fr->f_twin,NULL,
820 nrnb,
821 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
822 wallcycle_stop(wcycle,ewcVSITESPREAD);
826 if (flags & GMX_FORCE_VIRIAL)
828 /* Calculation of the virial must be done after vsites! */
829 calc_virial(fplog,mdatoms->start,mdatoms->homenr,x,f,
830 vir_force,graph,box,nrnb,fr,inputrec->ePBC);
834 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
836 /* Calculate the center of mass forces, this requires communication,
837 * which is why pull_potential is called close to other communication.
838 * The virial contribution is calculated directly,
839 * which is why we call pull_potential after calc_virial.
841 set_pbc(&pbc,inputrec->ePBC,box);
842 dvdl = 0;
843 enerd->term[F_COM_PULL] =
844 pull_potential(inputrec->ePull,inputrec->pull,mdatoms,&pbc,
845 cr,t,lambda,x,f,vir_force,&dvdl);
846 if (bSepDVDL)
848 fprintf(fplog,sepdvdlformat,"Com pull",enerd->term[F_COM_PULL],dvdl);
850 enerd->dvdl_lin += dvdl;
852 else
853 enerd->term[F_COM_PULL] = 0.0;
855 /* Add the forces from enforced rotation potentials (if any) */
856 if (inputrec->bRot)
857 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr,step,t);
860 if (PAR(cr) && !(cr->duty & DUTY_PME))
862 cycles_ppdpme = wallcycle_stop(wcycle,ewcPPDURINGPME);
863 dd_cycles_add(cr->dd,cycles_ppdpme,ddCyclPPduringPME);
865 /* In case of node-splitting, the PP nodes receive the long-range
866 * forces, virial and energy from the PME nodes here.
868 wallcycle_start(wcycle,ewcPP_PMEWAITRECVF);
869 dvdl = 0;
870 gmx_pme_receive_f(cr,fr->f_novirsum,fr->vir_el_recip,&e,&dvdl,
871 &cycles_seppme);
872 if (bSepDVDL)
874 fprintf(fplog,sepdvdlformat,"PME mesh",e,dvdl);
876 enerd->term[F_COUL_RECIP] += e;
877 enerd->dvdl_lin += dvdl;
878 if (wcycle)
880 dd_cycles_add(cr->dd,cycles_seppme,ddCyclPME);
882 wallcycle_stop(wcycle,ewcPP_PMEWAITRECVF);
885 if (bDoForces && fr->bF_NoVirSum)
887 if (vsite)
889 /* Spread the mesh force on virtual sites to the other particles...
890 * This is parallellized. MPI communication is performed
891 * if the constructing atoms aren't local.
893 wallcycle_start(wcycle,ewcVSITESPREAD);
894 spread_vsite_f(fplog,vsite,x,fr->f_novirsum,NULL,nrnb,
895 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
896 wallcycle_stop(wcycle,ewcVSITESPREAD);
898 if (flags & GMX_FORCE_VIRIAL)
900 /* Now add the forces, this is local */
901 if (fr->bDomDec)
903 sum_forces(0,fr->f_novirsum_n,f,fr->f_novirsum);
905 else
907 sum_forces(start,start+homenr,f,fr->f_novirsum);
909 if (EEL_FULL(fr->eeltype))
911 /* Add the mesh contribution to the virial */
912 m_add(vir_force,fr->vir_el_recip,vir_force);
914 if (debug)
916 pr_rvecs(debug,0,"vir_force",vir_force,DIM);
921 /* Sum the potential energy terms from group contributions */
922 sum_epot(&(inputrec->opts),enerd);
924 if (fr->print_force >= 0 && bDoForces)
926 print_large_forces(stderr,mdatoms,cr,step,fr->print_force,x,f);
930 void do_constrain_first(FILE *fplog,gmx_constr_t constr,
931 t_inputrec *ir,t_mdatoms *md,
932 t_state *state,rvec *f,
933 t_graph *graph,t_commrec *cr,t_nrnb *nrnb,
934 t_forcerec *fr, gmx_localtop_t *top, tensor shake_vir)
936 int i,m,start,end;
937 gmx_large_int_t step;
938 double mass,tmass,vcm[4];
939 real dt=ir->delta_t;
940 real dvdlambda;
941 rvec *savex;
943 snew(savex,state->natoms);
945 start = md->start;
946 end = md->homenr + start;
948 if (debug)
949 fprintf(debug,"vcm: start=%d, homenr=%d, end=%d\n",
950 start,md->homenr,end);
951 /* Do a first constrain to reset particles... */
952 step = ir->init_step;
953 if (fplog)
955 char buf[STEPSTRSIZE];
956 fprintf(fplog,"\nConstraining the starting coordinates (step %s)\n",
957 gmx_step_str(step,buf));
959 dvdlambda = 0;
961 /* constrain the current position */
962 constrain(NULL,TRUE,FALSE,constr,&(top->idef),
963 ir,NULL,cr,step,0,md,
964 state->x,state->x,NULL,
965 state->box,state->lambda,&dvdlambda,
966 NULL,NULL,nrnb,econqCoord,ir->epc==epcMTTK,state->veta,state->veta);
967 if (EI_VV(ir->eI))
969 /* constrain the inital velocity, and save it */
970 /* also may be useful if we need the ekin from the halfstep for velocity verlet */
971 /* might not yet treat veta correctly */
972 constrain(NULL,TRUE,FALSE,constr,&(top->idef),
973 ir,NULL,cr,step,0,md,
974 state->x,state->v,state->v,
975 state->box,state->lambda,&dvdlambda,
976 NULL,NULL,nrnb,econqVeloc,ir->epc==epcMTTK,state->veta,state->veta);
978 /* constrain the inital velocities at t-dt/2 */
979 if (EI_STATE_VELOCITY(ir->eI) && ir->eI!=eiVV)
981 for(i=start; (i<end); i++)
983 for(m=0; (m<DIM); m++)
985 /* Reverse the velocity */
986 state->v[i][m] = -state->v[i][m];
987 /* Store the position at t-dt in buf */
988 savex[i][m] = state->x[i][m] + dt*state->v[i][m];
991 /* Shake the positions at t=-dt with the positions at t=0
992 * as reference coordinates.
994 if (fplog)
996 char buf[STEPSTRSIZE];
997 fprintf(fplog,"\nConstraining the coordinates at t0-dt (step %s)\n",
998 gmx_step_str(step,buf));
1000 dvdlambda = 0;
1001 constrain(NULL,TRUE,FALSE,constr,&(top->idef),
1002 ir,NULL,cr,step,-1,md,
1003 state->x,savex,NULL,
1004 state->box,state->lambda,&dvdlambda,
1005 state->v,NULL,nrnb,econqCoord,ir->epc==epcMTTK,state->veta,state->veta);
1007 for(i=start; i<end; i++) {
1008 for(m=0; m<DIM; m++) {
1009 /* Re-reverse the velocities */
1010 state->v[i][m] = -state->v[i][m];
1015 for(m=0; (m<4); m++)
1016 vcm[m] = 0;
1017 for(i=start; i<end; i++) {
1018 mass = md->massT[i];
1019 for(m=0; m<DIM; m++) {
1020 vcm[m] += state->v[i][m]*mass;
1022 vcm[3] += mass;
1025 if (ir->nstcomm != 0 || debug) {
1026 /* Compute the global sum of vcm */
1027 if (debug)
1028 fprintf(debug,"vcm: %8.3f %8.3f %8.3f,"
1029 " total mass = %12.5e\n",vcm[XX],vcm[YY],vcm[ZZ],vcm[3]);
1030 if (PAR(cr))
1031 gmx_sumd(4,vcm,cr);
1032 tmass = vcm[3];
1033 for(m=0; (m<DIM); m++)
1034 vcm[m] /= tmass;
1035 if (debug)
1036 fprintf(debug,"vcm: %8.3f %8.3f %8.3f,"
1037 " total mass = %12.5e\n",vcm[XX],vcm[YY],vcm[ZZ],tmass);
1038 if (ir->nstcomm != 0) {
1039 /* Now we have the velocity of center of mass, let's remove it */
1040 for(i=start; (i<end); i++) {
1041 for(m=0; (m<DIM); m++)
1042 state->v[i][m] -= vcm[m];
1047 sfree(savex);
1050 void calc_enervirdiff(FILE *fplog,int eDispCorr,t_forcerec *fr)
1052 double eners[2],virs[2],enersum,virsum,y0,f,g,h;
1053 double r0,r1,r,rc3,rc9,ea,eb,ec,pa,pb,pc,pd;
1054 double invscale,invscale2,invscale3;
1055 int ri0,ri1,ri,i,offstart,offset;
1056 real scale,*vdwtab;
1058 fr->enershiftsix = 0;
1059 fr->enershifttwelve = 0;
1060 fr->enerdiffsix = 0;
1061 fr->enerdifftwelve = 0;
1062 fr->virdiffsix = 0;
1063 fr->virdifftwelve = 0;
1065 if (eDispCorr != edispcNO) {
1066 for(i=0; i<2; i++) {
1067 eners[i] = 0;
1068 virs[i] = 0;
1070 if ((fr->vdwtype == evdwSWITCH) || (fr->vdwtype == evdwSHIFT)) {
1071 if (fr->rvdw_switch == 0)
1072 gmx_fatal(FARGS,
1073 "With dispersion correction rvdw-switch can not be zero "
1074 "for vdw-type = %s",evdw_names[fr->vdwtype]);
1076 scale = fr->nblists[0].tab.scale;
1077 vdwtab = fr->nblists[0].vdwtab;
1079 /* Round the cut-offs to exact table values for precision */
1080 ri0 = floor(fr->rvdw_switch*scale);
1081 ri1 = ceil(fr->rvdw*scale);
1082 r0 = ri0/scale;
1083 r1 = ri1/scale;
1084 rc3 = r0*r0*r0;
1085 rc9 = rc3*rc3*rc3;
1087 if (fr->vdwtype == evdwSHIFT) {
1088 /* Determine the constant energy shift below rvdw_switch */
1089 fr->enershiftsix = (real)(-1.0/(rc3*rc3)) - vdwtab[8*ri0];
1090 fr->enershifttwelve = (real)( 1.0/(rc9*rc3)) - vdwtab[8*ri0 + 4];
1092 /* Add the constant part from 0 to rvdw_switch.
1093 * This integration from 0 to rvdw_switch overcounts the number
1094 * of interactions by 1, as it also counts the self interaction.
1095 * We will correct for this later.
1097 eners[0] += 4.0*M_PI*fr->enershiftsix*rc3/3.0;
1098 eners[1] += 4.0*M_PI*fr->enershifttwelve*rc3/3.0;
1100 invscale = 1.0/(scale);
1101 invscale2 = invscale*invscale;
1102 invscale3 = invscale*invscale2;
1104 /* following summation derived from cubic spline definition,
1105 Numerical Recipies in C, second edition, p. 113-116. Exact
1106 for the cubic spline. We first calculate the negative of
1107 the energy from rvdw to rvdw_switch, assuming that g(r)=1,
1108 and then add the more standard, abrupt cutoff correction to
1109 that result, yielding the long-range correction for a
1110 switched function. We perform both the pressure and energy
1111 loops at the same time for simplicity, as the computational
1112 cost is low. */
1114 for (i=0;i<2;i++) {
1115 enersum = 0.0; virsum = 0.0;
1116 if (i==0)
1117 offstart = 0;
1118 else
1119 offstart = 4;
1120 for (ri=ri0; ri<ri1; ri++) {
1121 r = ri*invscale;
1122 ea = invscale3;
1123 eb = 2.0*invscale2*r;
1124 ec = invscale*r*r;
1126 pa = invscale3;
1127 pb = 3.0*invscale2*r;
1128 pc = 3.0*invscale*r*r;
1129 pd = r*r*r;
1131 /* this "8" is from the packing in the vdwtab array - perhaps
1132 should be #define'ed? */
1133 offset = 8*ri + offstart;
1134 y0 = vdwtab[offset];
1135 f = vdwtab[offset+1];
1136 g = vdwtab[offset+2];
1137 h = vdwtab[offset+3];
1139 enersum += y0*(ea/3 + eb/2 + ec) + f*(ea/4 + eb/3 + ec/2)+
1140 g*(ea/5 + eb/4 + ec/3) + h*(ea/6 + eb/5 + ec/4);
1141 virsum += f*(pa/4 + pb/3 + pc/2 + pd) +
1142 2*g*(pa/5 + pb/4 + pc/3 + pd/2) + 3*h*(pa/6 + pb/5 + pc/4 + pd/3);
1145 enersum *= 4.0*M_PI;
1146 virsum *= 4.0*M_PI;
1147 eners[i] -= enersum;
1148 virs[i] -= virsum;
1151 /* now add the correction for rvdw_switch to infinity */
1152 eners[0] += -4.0*M_PI/(3.0*rc3);
1153 eners[1] += 4.0*M_PI/(9.0*rc9);
1154 virs[0] += 8.0*M_PI/rc3;
1155 virs[1] += -16.0*M_PI/(3.0*rc9);
1157 else if ((fr->vdwtype == evdwCUT) || (fr->vdwtype == evdwUSER)) {
1158 if (fr->vdwtype == evdwUSER && fplog)
1159 fprintf(fplog,
1160 "WARNING: using dispersion correction with user tables\n");
1161 rc3 = fr->rvdw*fr->rvdw*fr->rvdw;
1162 rc9 = rc3*rc3*rc3;
1163 eners[0] += -4.0*M_PI/(3.0*rc3);
1164 eners[1] += 4.0*M_PI/(9.0*rc9);
1165 virs[0] += 8.0*M_PI/rc3;
1166 virs[1] += -16.0*M_PI/(3.0*rc9);
1167 } else {
1168 gmx_fatal(FARGS,
1169 "Dispersion correction is not implemented for vdw-type = %s",
1170 evdw_names[fr->vdwtype]);
1172 fr->enerdiffsix = eners[0];
1173 fr->enerdifftwelve = eners[1];
1174 /* The 0.5 is due to the Gromacs definition of the virial */
1175 fr->virdiffsix = 0.5*virs[0];
1176 fr->virdifftwelve = 0.5*virs[1];
1180 void calc_dispcorr(FILE *fplog,t_inputrec *ir,t_forcerec *fr,
1181 gmx_large_int_t step,int natoms,
1182 matrix box,real lambda,tensor pres,tensor virial,
1183 real *prescorr, real *enercorr, real *dvdlcorr)
1185 gmx_bool bCorrAll,bCorrPres;
1186 real dvdlambda,invvol,dens,ninter,avcsix,avctwelve,enerdiff,svir=0,spres=0;
1187 int m;
1189 *prescorr = 0;
1190 *enercorr = 0;
1191 *dvdlcorr = 0;
1193 clear_mat(virial);
1194 clear_mat(pres);
1196 if (ir->eDispCorr != edispcNO) {
1197 bCorrAll = (ir->eDispCorr == edispcAllEner ||
1198 ir->eDispCorr == edispcAllEnerPres);
1199 bCorrPres = (ir->eDispCorr == edispcEnerPres ||
1200 ir->eDispCorr == edispcAllEnerPres);
1202 invvol = 1/det(box);
1203 if (fr->n_tpi)
1205 /* Only correct for the interactions with the inserted molecule */
1206 dens = (natoms - fr->n_tpi)*invvol;
1207 ninter = fr->n_tpi;
1209 else
1211 dens = natoms*invvol;
1212 ninter = 0.5*natoms;
1215 if (ir->efep == efepNO)
1217 avcsix = fr->avcsix[0];
1218 avctwelve = fr->avctwelve[0];
1220 else
1222 avcsix = (1 - lambda)*fr->avcsix[0] + lambda*fr->avcsix[1];
1223 avctwelve = (1 - lambda)*fr->avctwelve[0] + lambda*fr->avctwelve[1];
1226 enerdiff = ninter*(dens*fr->enerdiffsix - fr->enershiftsix);
1227 *enercorr += avcsix*enerdiff;
1228 dvdlambda = 0.0;
1229 if (ir->efep != efepNO)
1231 dvdlambda += (fr->avcsix[1] - fr->avcsix[0])*enerdiff;
1233 if (bCorrAll)
1235 enerdiff = ninter*(dens*fr->enerdifftwelve - fr->enershifttwelve);
1236 *enercorr += avctwelve*enerdiff;
1237 if (fr->efep != efepNO)
1239 dvdlambda += (fr->avctwelve[1] - fr->avctwelve[0])*enerdiff;
1243 if (bCorrPres)
1245 svir = ninter*dens*avcsix*fr->virdiffsix/3.0;
1246 if (ir->eDispCorr == edispcAllEnerPres)
1248 svir += ninter*dens*avctwelve*fr->virdifftwelve/3.0;
1250 /* The factor 2 is because of the Gromacs virial definition */
1251 spres = -2.0*invvol*svir*PRESFAC;
1253 for(m=0; m<DIM; m++) {
1254 virial[m][m] += svir;
1255 pres[m][m] += spres;
1257 *prescorr += spres;
1260 /* Can't currently control when it prints, for now, just print when degugging */
1261 if (debug)
1263 if (bCorrAll) {
1264 fprintf(debug,"Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
1265 avcsix,avctwelve);
1267 if (bCorrPres)
1269 fprintf(debug,
1270 "Long Range LJ corr.: Epot %10g, Pres: %10g, Vir: %10g\n",
1271 *enercorr,spres,svir);
1273 else
1275 fprintf(debug,"Long Range LJ corr.: Epot %10g\n",*enercorr);
1279 if (fr->bSepDVDL && do_per_step(step,ir->nstlog))
1281 fprintf(fplog,sepdvdlformat,"Dispersion correction",
1282 *enercorr,dvdlambda);
1284 if (fr->efep != efepNO)
1286 *dvdlcorr += dvdlambda;
1291 void do_pbc_first(FILE *fplog,matrix box,t_forcerec *fr,
1292 t_graph *graph,rvec x[])
1294 if (fplog)
1295 fprintf(fplog,"Removing pbc first time\n");
1296 calc_shifts(box,fr->shift_vec);
1297 if (graph) {
1298 mk_mshift(fplog,graph,fr->ePBC,box,x);
1299 if (gmx_debug_at)
1300 p_graph(debug,"do_pbc_first 1",graph);
1301 shift_self(graph,box,x);
1302 /* By doing an extra mk_mshift the molecules that are broken
1303 * because they were e.g. imported from another software
1304 * will be made whole again. Such are the healing powers
1305 * of GROMACS.
1307 mk_mshift(fplog,graph,fr->ePBC,box,x);
1308 if (gmx_debug_at)
1309 p_graph(debug,"do_pbc_first 2",graph);
1311 if (fplog)
1312 fprintf(fplog,"Done rmpbc\n");
1315 static void low_do_pbc_mtop(FILE *fplog,int ePBC,matrix box,
1316 gmx_mtop_t *mtop,rvec x[],
1317 gmx_bool bFirst)
1319 t_graph *graph;
1320 int mb,as,mol;
1321 gmx_molblock_t *molb;
1323 if (bFirst && fplog)
1324 fprintf(fplog,"Removing pbc first time\n");
1326 snew(graph,1);
1327 as = 0;
1328 for(mb=0; mb<mtop->nmolblock; mb++) {
1329 molb = &mtop->molblock[mb];
1330 if (molb->natoms_mol == 1 ||
1331 (!bFirst && mtop->moltype[molb->type].cgs.nr == 1)) {
1332 /* Just one atom or charge group in the molecule, no PBC required */
1333 as += molb->nmol*molb->natoms_mol;
1334 } else {
1335 /* Pass NULL iso fplog to avoid graph prints for each molecule type */
1336 mk_graph_ilist(NULL,mtop->moltype[molb->type].ilist,
1337 0,molb->natoms_mol,FALSE,FALSE,graph);
1339 for(mol=0; mol<molb->nmol; mol++) {
1340 mk_mshift(fplog,graph,ePBC,box,x+as);
1342 shift_self(graph,box,x+as);
1343 /* The molecule is whole now.
1344 * We don't need the second mk_mshift call as in do_pbc_first,
1345 * since we no longer need this graph.
1348 as += molb->natoms_mol;
1350 done_graph(graph);
1353 sfree(graph);
1356 void do_pbc_first_mtop(FILE *fplog,int ePBC,matrix box,
1357 gmx_mtop_t *mtop,rvec x[])
1359 low_do_pbc_mtop(fplog,ePBC,box,mtop,x,TRUE);
1362 void do_pbc_mtop(FILE *fplog,int ePBC,matrix box,
1363 gmx_mtop_t *mtop,rvec x[])
1365 low_do_pbc_mtop(fplog,ePBC,box,mtop,x,FALSE);
1368 void finish_run(FILE *fplog,t_commrec *cr,const char *confout,
1369 t_inputrec *inputrec,
1370 t_nrnb nrnb[],gmx_wallcycle_t wcycle,
1371 gmx_runtime_t *runtime,
1372 gmx_bool bWriteStat)
1374 int i,j;
1375 t_nrnb *nrnb_tot=NULL;
1376 real delta_t;
1377 double nbfs,mflop;
1378 double cycles[ewcNR];
1380 wallcycle_sum(cr,wcycle,cycles);
1382 if (cr->nnodes > 1) {
1383 if (SIMMASTER(cr))
1384 snew(nrnb_tot,1);
1385 #ifdef GMX_MPI
1386 MPI_Reduce(nrnb->n,nrnb_tot->n,eNRNB,MPI_DOUBLE,MPI_SUM,
1387 MASTERRANK(cr),cr->mpi_comm_mysim);
1388 #endif
1389 } else {
1390 nrnb_tot = nrnb;
1393 if (SIMMASTER(cr)) {
1394 print_flop(fplog,nrnb_tot,&nbfs,&mflop);
1395 if (cr->nnodes > 1) {
1396 sfree(nrnb_tot);
1400 if ((cr->duty & DUTY_PP) && DOMAINDECOMP(cr)) {
1401 print_dd_statistics(cr,inputrec,fplog);
1404 #ifdef GMX_MPI
1405 if (PARTDECOMP(cr))
1407 if (MASTER(cr))
1409 t_nrnb *nrnb_all;
1410 int s;
1411 MPI_Status stat;
1413 snew(nrnb_all,cr->nnodes);
1414 nrnb_all[0] = *nrnb;
1415 for(s=1; s<cr->nnodes; s++)
1417 MPI_Recv(nrnb_all[s].n,eNRNB,MPI_DOUBLE,s,0,
1418 cr->mpi_comm_mysim,&stat);
1420 pr_load(fplog,cr,nrnb_all);
1421 sfree(nrnb_all);
1423 else
1425 MPI_Send(nrnb->n,eNRNB,MPI_DOUBLE,MASTERRANK(cr),0,
1426 cr->mpi_comm_mysim);
1429 #endif
1431 if (SIMMASTER(cr)) {
1432 wallcycle_print(fplog,cr->nnodes,cr->npmenodes,runtime->realtime,
1433 wcycle,cycles);
1435 if (EI_DYNAMICS(inputrec->eI)) {
1436 delta_t = inputrec->delta_t;
1437 } else {
1438 delta_t = 0;
1441 if (fplog) {
1442 print_perf(fplog,runtime->proctime,runtime->realtime,
1443 cr->nnodes-cr->npmenodes,
1444 runtime->nsteps_done,delta_t,nbfs,mflop);
1446 if (bWriteStat) {
1447 print_perf(stderr,runtime->proctime,runtime->realtime,
1448 cr->nnodes-cr->npmenodes,
1449 runtime->nsteps_done,delta_t,nbfs,mflop);
1453 runtime=inputrec->nsteps*inputrec->delta_t;
1454 if (bWriteStat) {
1455 if (cr->nnodes == 1)
1456 fprintf(stderr,"\n\n");
1457 print_perf(stderr,nodetime,realtime,runtime,&ntot,
1458 cr->nnodes-cr->npmenodes,FALSE);
1460 wallcycle_print(fplog,cr->nnodes,cr->npmenodes,realtime,wcycle,cycles);
1461 print_perf(fplog,nodetime,realtime,runtime,&ntot,cr->nnodes-cr->npmenodes,
1462 TRUE);
1463 if (PARTDECOMP(cr))
1464 pr_load(fplog,cr,nrnb_all);
1465 if (cr->nnodes > 1)
1466 sfree(nrnb_all);
1471 void init_md(FILE *fplog,
1472 t_commrec *cr,t_inputrec *ir,const output_env_t oenv,
1473 double *t,double *t0,
1474 real *lambda,double *lam0,
1475 t_nrnb *nrnb,gmx_mtop_t *mtop,
1476 gmx_update_t *upd,
1477 int nfile,const t_filenm fnm[],
1478 gmx_mdoutf_t **outf,t_mdebin **mdebin,
1479 tensor force_vir,tensor shake_vir,rvec mu_tot,
1480 gmx_bool *bSimAnn,t_vcm **vcm, t_state *state, unsigned long Flags)
1482 int i,j,n;
1483 real tmpt,mod;
1485 /* Initial values */
1486 *t = *t0 = ir->init_t;
1487 if (ir->efep != efepNO)
1489 *lam0 = ir->init_lambda;
1490 *lambda = *lam0 + ir->init_step*ir->delta_lambda;
1492 else
1494 *lambda = *lam0 = 0.0;
1497 *bSimAnn=FALSE;
1498 for(i=0;i<ir->opts.ngtc;i++)
1500 /* set bSimAnn if any group is being annealed */
1501 if(ir->opts.annealing[i]!=eannNO)
1503 *bSimAnn = TRUE;
1506 if (*bSimAnn)
1508 update_annealing_target_temp(&(ir->opts),ir->init_t);
1511 if (upd)
1513 *upd = init_update(fplog,ir);
1516 if (vcm != NULL)
1518 *vcm = init_vcm(fplog,&mtop->groups,ir);
1521 if (EI_DYNAMICS(ir->eI) && !(Flags & MD_APPENDFILES))
1523 if (ir->etc == etcBERENDSEN)
1525 please_cite(fplog,"Berendsen84a");
1527 if (ir->etc == etcVRESCALE)
1529 please_cite(fplog,"Bussi2007a");
1533 init_nrnb(nrnb);
1535 if (nfile != -1)
1537 *outf = init_mdoutf(nfile,fnm,Flags,cr,ir,oenv);
1539 *mdebin = init_mdebin((Flags & MD_APPENDFILES) ? NULL : (*outf)->fp_ene,
1540 mtop,ir, (*outf)->fp_dhdl);
1543 /* Initiate variables */
1544 clear_mat(force_vir);
1545 clear_mat(shake_vir);
1546 clear_rvec(mu_tot);
1548 debug_gmx();