Merge branch 'master' into rotation
[gromacs/adressmacs.git] / src / mdlib / sim_util.c
blobea53069dec5ddad8e89b9b71197869b8b11774ab
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"
82 #include "pull_rotation.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_THREADS
93 #include "tmpi.h"
94 #endif
96 #include "qmmm.h"
98 #if 0
99 typedef struct gmx_timeprint {
101 } t_gmx_timeprint;
102 #endif
105 double
106 gmx_gettime()
108 #ifdef HAVE_GETTIMEOFDAY
109 struct timeval t;
110 struct timezone tz = { 0,0 };
111 double seconds;
113 gettimeofday(&t,&tz);
115 seconds = (double) t.tv_sec + 1e-6*(double)t.tv_usec;
117 return seconds;
118 #else
119 double seconds;
121 seconds = time(NULL);
123 return seconds;
124 #endif
128 #define difftime(end,start) ((double)(end)-(double)(start))
130 void print_time(FILE *out,gmx_runtime_t *runtime,gmx_large_int_t step,
131 t_inputrec *ir, t_commrec *cr)
133 time_t finish;
135 double dt;
136 char buf[48];
138 #ifndef GMX_THREADS
139 if (!PAR(cr))
140 #endif
141 fprintf(out,"\r");
142 fprintf(out,"step %s",gmx_step_str(step,buf));
143 if ((step >= ir->nstlist)) {
144 if ((ir->nstlist == 0) || ((step % ir->nstlist) == 0)) {
145 /* We have done a full cycle let's update time_per_step */
146 runtime->last = gmx_gettime();
147 dt = difftime(runtime->last,runtime->real);
148 runtime->time_per_step = dt/(step - ir->init_step + 1);
150 dt = (ir->nsteps + ir->init_step - step)*runtime->time_per_step;
152 if (dt >= 300) {
153 finish = (time_t) (runtime->last + dt);
154 sprintf(buf,"%s",ctime(&finish));
155 buf[strlen(buf)-1]='\0';
156 fprintf(out,", will finish %s",buf);
158 else
159 fprintf(out,", remaining runtime: %5d s ",(int)dt);
161 #ifndef GMX_THREADS
162 if (PAR(cr))
163 fprintf(out,"\n");
164 #endif
166 fflush(out);
169 #ifdef NO_CLOCK
170 #define clock() -1
171 #endif
173 static double set_proctime(gmx_runtime_t *runtime)
175 double diff;
176 #ifdef GMX_CRAY_XT3
177 double prev;
179 prev = runtime->proc;
180 runtime->proc = dclock();
182 diff = runtime->proc - prev;
183 #else
184 clock_t prev;
186 prev = runtime->proc;
187 runtime->proc = clock();
189 diff = (double)(runtime->proc - prev)/(double)CLOCKS_PER_SEC;
190 #endif
191 if (diff < 0)
193 /* The counter has probably looped, ignore this data */
194 diff = 0;
197 return diff;
200 void runtime_start(gmx_runtime_t *runtime)
202 runtime->real = gmx_gettime();
203 runtime->proc = 0;
204 set_proctime(runtime);
205 runtime->realtime = 0;
206 runtime->proctime = 0;
207 runtime->last = 0;
208 runtime->time_per_step = 0;
211 void runtime_end(gmx_runtime_t *runtime)
213 double now;
215 now = gmx_gettime();
217 runtime->proctime += set_proctime(runtime);
218 runtime->realtime = now - runtime->real;
219 runtime->real = now;
222 void runtime_upd_proc(gmx_runtime_t *runtime)
224 runtime->proctime += set_proctime(runtime);
227 void print_date_and_time(FILE *fplog,int nodeid,const char *title,
228 const gmx_runtime_t *runtime)
230 int i;
231 char *ts,time_string[STRLEN];
232 time_t tmptime;
234 if (runtime != NULL)
236 tmptime = (time_t) runtime->real;
237 ts = ctime(&tmptime);
239 else
241 tmptime = (time_t) gmx_gettime();
242 ts = ctime(&tmptime);
244 for(i=0; ts[i]>=' '; i++)
246 time_string[i]=ts[i];
248 time_string[i]='\0';
249 if (fplog)
251 fprintf(fplog,"%s on node %d %s\n",title,nodeid,time_string);
255 static void sum_forces(int start,int end,rvec f[],rvec flr[])
257 int i;
259 if (gmx_debug_at) {
260 pr_rvecs(debug,0,"fsr",f+start,end-start);
261 pr_rvecs(debug,0,"flr",flr+start,end-start);
263 for(i=start; (i<end); i++)
264 rvec_inc(f[i],flr[i]);
268 * calc_f_el calculates forces due to an electric field.
270 * force is kJ mol^-1 nm^-1 = e * kJ mol^-1 nm^-1 / e
272 * Et[] contains the parameters for the time dependent
273 * part of the field (not yet used).
274 * Ex[] contains the parameters for
275 * the spatial dependent part of the field. You can have cool periodic
276 * fields in principle, but only a constant field is supported
277 * now.
278 * The function should return the energy due to the electric field
279 * (if any) but for now returns 0.
281 * WARNING:
282 * There can be problems with the virial.
283 * Since the field is not self-consistent this is unavoidable.
284 * For neutral molecules the virial is correct within this approximation.
285 * For neutral systems with many charged molecules the error is small.
286 * But for systems with a net charge or a few charged molecules
287 * the error can be significant when the field is high.
288 * Solution: implement a self-consitent electric field into PME.
290 static void calc_f_el(FILE *fp,int start,int homenr,
291 real charge[],rvec x[],rvec f[],
292 t_cosines Ex[],t_cosines Et[],double t)
294 rvec Ext;
295 real t0;
296 int i,m;
298 for(m=0; (m<DIM); m++)
300 if (Et[m].n > 0)
302 if (Et[m].n == 3)
304 t0 = Et[m].a[1];
305 Ext[m] = cos(Et[m].a[0]*(t-t0))*exp(-sqr(t-t0)/(2.0*sqr(Et[m].a[2])));
307 else
309 Ext[m] = cos(Et[m].a[0]*t);
312 else
314 Ext[m] = 1.0;
316 if (Ex[m].n > 0)
318 /* Convert the field strength from V/nm to MD-units */
319 Ext[m] *= Ex[m].a[0]*FIELDFAC;
320 for(i=start; (i<start+homenr); i++)
321 f[i][m] += charge[i]*Ext[m];
323 else
325 Ext[m] = 0;
328 if (fp != NULL)
330 fprintf(fp,"%10g %10g %10g %10g #FIELD\n",t,
331 Ext[XX]/FIELDFAC,Ext[YY]/FIELDFAC,Ext[ZZ]/FIELDFAC);
335 static void calc_virial(FILE *fplog,int start,int homenr,rvec x[],rvec f[],
336 tensor vir_part,t_graph *graph,matrix box,
337 t_nrnb *nrnb,const t_forcerec *fr,int ePBC)
339 int i,j;
340 tensor virtest;
342 /* The short-range virial from surrounding boxes */
343 clear_mat(vir_part);
344 calc_vir(fplog,SHIFTS,fr->shift_vec,fr->fshift,vir_part,ePBC==epbcSCREW,box);
345 inc_nrnb(nrnb,eNR_VIRIAL,SHIFTS);
347 /* Calculate partial virial, for local atoms only, based on short range.
348 * Total virial is computed in global_stat, called from do_md
350 f_calc_vir(fplog,start,start+homenr,x,f,vir_part,graph,box);
351 inc_nrnb(nrnb,eNR_VIRIAL,homenr);
353 /* Add position restraint contribution */
354 for(i=0; i<DIM; i++) {
355 vir_part[i][i] += fr->vir_diag_posres[i];
358 /* Add wall contribution */
359 for(i=0; i<DIM; i++) {
360 vir_part[i][ZZ] += fr->vir_wall_z[i];
363 if (debug)
364 pr_rvecs(debug,0,"vir_part",vir_part,DIM);
367 static void print_large_forces(FILE *fp,t_mdatoms *md,t_commrec *cr,
368 gmx_large_int_t step,real pforce,rvec *x,rvec *f)
370 int i;
371 real pf2,fn2;
372 char buf[22];
374 pf2 = sqr(pforce);
375 for(i=md->start; i<md->start+md->homenr; i++) {
376 fn2 = norm2(f[i]);
377 if (fn2 >= pf2) {
378 fprintf(fp,"step %s atom %6d x %8.3f %8.3f %8.3f force %12.5e\n",
379 gmx_step_str(step,buf),
380 ddglatnr(cr->dd,i),x[i][XX],x[i][YY],x[i][ZZ],sqrt(fn2));
385 void do_force(FILE *fplog,t_commrec *cr,
386 t_inputrec *inputrec,
387 gmx_large_int_t step,t_nrnb *nrnb,gmx_wallcycle_t wcycle,
388 gmx_localtop_t *top,
389 gmx_mtop_t *mtop,
390 gmx_groups_t *groups,
391 matrix box,rvec x[],history_t *hist,
392 rvec f[],
393 tensor vir_force,
394 t_mdatoms *mdatoms,
395 gmx_enerdata_t *enerd,t_fcdata *fcd,
396 real lambda,t_graph *graph,
397 t_forcerec *fr,gmx_vsite_t *vsite,rvec mu_tot,
398 double t,FILE *field,gmx_edsam_t ed,
399 bool bBornRadii,
400 int flags)
402 int cg0,cg1,i,j;
403 int start,homenr;
404 double mu[2*DIM];
405 bool bSepDVDL,bStateChanged,bNS,bFillGrid,bCalcCGCM,bBS;
406 bool bDoLongRange,bDoForces,bSepLRF;
407 matrix boxs;
408 real e,v,dvdl;
409 t_pbc pbc;
410 float cycles_ppdpme,cycles_pme,cycles_seppme,cycles_force;
412 start = mdatoms->start;
413 homenr = mdatoms->homenr;
415 bSepDVDL = (fr->bSepDVDL && do_per_step(step,inputrec->nstlog));
417 clear_mat(vir_force);
419 if (PARTDECOMP(cr))
421 pd_cg_range(cr,&cg0,&cg1);
423 else
425 cg0 = 0;
426 if (DOMAINDECOMP(cr))
428 cg1 = cr->dd->ncg_tot;
430 else
432 cg1 = top->cgs.nr;
434 if (fr->n_tpi > 0)
436 cg1--;
440 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
441 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll==FALSE);
442 bFillGrid = (bNS && bStateChanged);
443 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
444 bDoLongRange = (fr->bTwinRange && bNS && (flags & GMX_FORCE_DOLR));
445 bDoForces = (flags & GMX_FORCE_FORCES);
446 bSepLRF = (bDoLongRange && bDoForces && (flags & GMX_FORCE_SEPLRF));
448 if (bStateChanged)
450 update_forcerec(fplog,fr,box);
452 /* Calculate total (local) dipole moment in a temporary common array.
453 * This makes it possible to sum them over nodes faster.
455 calc_mu(start,homenr,
456 x,mdatoms->chargeA,mdatoms->chargeB,mdatoms->nChargePerturbed,
457 mu,mu+DIM);
460 if (fr->ePBC != epbcNONE) {
461 /* Compute shift vectors every step,
462 * because of pressure coupling or box deformation!
464 if ((flags & GMX_FORCE_DYNAMICBOX) && bStateChanged)
465 calc_shifts(box,fr->shift_vec);
467 if (bCalcCGCM) {
468 put_charge_groups_in_box(fplog,cg0,cg1,fr->ePBC,box,
469 &(top->cgs),x,fr->cg_cm);
470 inc_nrnb(nrnb,eNR_CGCM,homenr);
471 inc_nrnb(nrnb,eNR_RESETX,cg1-cg0);
473 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph) {
474 unshift_self(graph,box,x);
477 else if (bCalcCGCM) {
478 calc_cgcm(fplog,cg0,cg1,&(top->cgs),x,fr->cg_cm);
479 inc_nrnb(nrnb,eNR_CGCM,homenr);
482 if (bCalcCGCM) {
483 if (PAR(cr)) {
484 move_cgcm(fplog,cr,fr->cg_cm);
486 if (gmx_debug_at)
487 pr_rvecs(debug,0,"cgcm",fr->cg_cm,top->cgs.nr);
490 #ifdef GMX_MPI
491 if (!(cr->duty & DUTY_PME)) {
492 /* Send particle coordinates to the pme nodes.
493 * Since this is only implemented for domain decomposition
494 * and domain decomposition does not use the graph,
495 * we do not need to worry about shifting.
498 wallcycle_start(wcycle,ewcPP_PMESENDX);
499 GMX_MPE_LOG(ev_send_coordinates_start);
501 bBS = (inputrec->nwall == 2);
502 if (bBS) {
503 copy_mat(box,boxs);
504 svmul(inputrec->wall_ewald_zfac,boxs[ZZ],boxs[ZZ]);
507 gmx_pme_send_x(cr,bBS ? boxs : box,x,mdatoms->nChargePerturbed,lambda,step);
509 GMX_MPE_LOG(ev_send_coordinates_finish);
510 wallcycle_stop(wcycle,ewcPP_PMESENDX);
512 #endif /* GMX_MPI */
514 /* Communicate coordinates and sum dipole if necessary */
515 if (PAR(cr))
517 wallcycle_start(wcycle,ewcMOVEX);
518 if (DOMAINDECOMP(cr))
520 dd_move_x(cr->dd,box,x);
522 else
524 move_x(fplog,cr,GMX_LEFT,GMX_RIGHT,x,nrnb);
526 /* When we don't need the total dipole we sum it in global_stat */
527 if (bStateChanged && NEED_MUTOT(*inputrec))
529 gmx_sumd(2*DIM,mu,cr);
531 wallcycle_stop(wcycle,ewcMOVEX);
533 if (bStateChanged)
535 for(i=0; i<2; i++)
537 for(j=0;j<DIM;j++)
539 fr->mu_tot[i][j] = mu[i*DIM + j];
543 if (fr->efep == efepNO)
545 copy_rvec(fr->mu_tot[0],mu_tot);
547 else
549 for(j=0; j<DIM; j++)
551 mu_tot[j] =
552 (1.0 - lambda)*fr->mu_tot[0][j] + lambda*fr->mu_tot[1][j];
556 /* Reset energies */
557 reset_enerdata(&(inputrec->opts),fr,bNS,enerd,MASTER(cr));
558 clear_rvecs(SHIFTS,fr->fshift);
560 if (bNS)
562 wallcycle_start(wcycle,ewcNS);
564 if (graph && bStateChanged)
566 /* Calculate intramolecular shift vectors to make molecules whole */
567 mk_mshift(fplog,graph,fr->ePBC,box,x);
570 /* Reset long range forces if necessary */
571 if (fr->bTwinRange)
573 /* Reset the (long-range) forces if necessary */
574 clear_rvecs(fr->natoms_force,bSepLRF ? fr->f_twin : f);
577 /* Do the actual neighbour searching and if twin range electrostatics
578 * also do the calculation of long range forces and energies.
580 dvdl = 0;
581 ns(fplog,fr,x,box,
582 groups,&(inputrec->opts),top,mdatoms,
583 cr,nrnb,lambda,&dvdl,&enerd->grpp,bFillGrid,
584 bDoLongRange,bDoForces,bSepLRF ? fr->f_twin : f);
585 if (bSepDVDL)
587 fprintf(fplog,sepdvdlformat,"LR non-bonded",0.0,dvdl);
589 enerd->dvdl_lin += dvdl;
591 wallcycle_stop(wcycle,ewcNS);
594 if (inputrec->implicit_solvent && bNS)
596 make_gb_nblist(cr,mtop->natoms,inputrec->gb_algorithm,inputrec->rlist,
597 x,box,fr,&top->idef,graph,fr->born);
600 if (DOMAINDECOMP(cr))
602 if (!(cr->duty & DUTY_PME))
604 wallcycle_start(wcycle,ewcPPDURINGPME);
605 dd_force_flop_start(cr->dd,nrnb);
609 if (inputrec->bRot)
611 /* Enforced rotation has its own cycle counter that starts after the collective
612 * coordinates have been communicated. It is added to ddCyclF */
613 do_rotation(cr,inputrec,box,x,t,step,wcycle,bNS);
616 /* Start the force cycle counter.
617 * This counter is stopped in do_forcelow_level.
618 * No parallel communication should occur while this counter is running,
619 * since that will interfere with the dynamic load balancing.
621 wallcycle_start(wcycle,ewcFORCE);
622 GMX_MPE_LOG(ev_forcecycles_start);
624 if (bDoForces)
626 /* Reset forces for which the virial is calculated separately:
627 * PME/Ewald forces if necessary */
628 if (fr->bF_NoVirSum)
630 if (flags & GMX_FORCE_VIRIAL)
632 fr->f_novirsum = fr->f_novirsum_alloc;
633 GMX_BARRIER(cr->mpi_comm_mygroup);
634 if (fr->bDomDec)
636 clear_rvecs(fr->f_novirsum_n,fr->f_novirsum);
638 else
640 clear_rvecs(homenr,fr->f_novirsum+start);
642 GMX_BARRIER(cr->mpi_comm_mygroup);
644 else
646 /* We are not calculating the pressure so we do not need
647 * a separate array for forces that do not contribute
648 * to the pressure.
650 fr->f_novirsum = f;
654 if (bSepLRF)
656 /* Add the long range forces to the short range forces */
657 for(i=0; i<fr->natoms_force; i++)
659 copy_rvec(fr->f_twin[i],f[i]);
662 else if (!(fr->bTwinRange && bNS))
664 /* Clear the short-range forces */
665 clear_rvecs(fr->natoms_force,f);
668 clear_rvec(fr->vir_diag_posres);
670 GMX_BARRIER(cr->mpi_comm_mygroup);
672 if (inputrec->ePull == epullCONSTRAINT)
674 clear_pull_forces(inputrec->pull);
677 /* update QMMMrec, if necessary */
678 if(fr->bQMMM)
680 update_QMMMrec(cr,fr,x,mdatoms,box,top);
683 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
685 /* Position restraints always require full pbc */
686 set_pbc(&pbc,inputrec->ePBC,box);
687 v = posres(top->idef.il[F_POSRES].nr,top->idef.il[F_POSRES].iatoms,
688 top->idef.iparams_posres,
689 (const rvec*)x,fr->f_novirsum,fr->vir_diag_posres,
690 inputrec->ePBC==epbcNONE ? NULL : &pbc,lambda,&dvdl,
691 fr->rc_scaling,fr->ePBC,fr->posres_com,fr->posres_comB);
692 if (bSepDVDL)
694 fprintf(fplog,sepdvdlformat,
695 interaction_function[F_POSRES].longname,v,dvdl);
697 enerd->term[F_POSRES] += v;
698 /* This linear lambda dependence assumption is only correct
699 * when only k depends on lambda,
700 * not when the reference position depends on lambda.
701 * grompp checks for this.
703 enerd->dvdl_lin += dvdl;
704 inc_nrnb(nrnb,eNR_POSRES,top->idef.il[F_POSRES].nr/2);
707 /* Compute the bonded and non-bonded energies and optionally forces */
708 do_force_lowlevel(fplog,step,fr,inputrec,&(top->idef),
709 cr,nrnb,wcycle,mdatoms,&(inputrec->opts),
710 x,hist,f,enerd,fcd,mtop,top,fr->born,
711 &(top->atomtypes),bBornRadii,box,
712 lambda,graph,&(top->excls),fr->mu_tot,
713 flags,&cycles_pme);
715 cycles_force = wallcycle_stop(wcycle,ewcFORCE);
716 GMX_BARRIER(cr->mpi_comm_mygroup);
718 if (ed)
720 do_flood(fplog,cr,x,f,ed,box,step);
723 if (DOMAINDECOMP(cr))
725 dd_force_flop_stop(cr->dd,nrnb);
726 if (wcycle)
728 dd_cycles_add(cr->dd,cycles_force-cycles_pme,ddCyclF);
732 if (bDoForces)
734 if (IR_ELEC_FIELD(*inputrec))
736 /* Compute forces due to electric field */
737 calc_f_el(MASTER(cr) ? field : NULL,
738 start,homenr,mdatoms->chargeA,x,fr->f_novirsum,
739 inputrec->ex,inputrec->et,t);
742 /* Communicate the forces */
743 if (PAR(cr))
745 wallcycle_start(wcycle,ewcMOVEF);
746 if (DOMAINDECOMP(cr))
748 dd_move_f(cr->dd,f,fr->fshift);
749 /* Do we need to communicate the separate force array
750 * for terms that do not contribute to the single sum virial?
751 * Position restraints and electric fields do not introduce
752 * inter-cg forces, only full electrostatics methods do.
753 * When we do not calculate the virial, fr->f_novirsum = f,
754 * so we have already communicated these forces.
756 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
757 (flags & GMX_FORCE_VIRIAL))
759 dd_move_f(cr->dd,fr->f_novirsum,NULL);
761 if (bSepLRF)
763 /* We should not update the shift forces here,
764 * since f_twin is already included in f.
766 dd_move_f(cr->dd,fr->f_twin,NULL);
769 else
771 pd_move_f(cr,f,nrnb);
772 if (bSepLRF)
774 pd_move_f(cr,fr->f_twin,nrnb);
777 wallcycle_stop(wcycle,ewcMOVEF);
780 /* If we have NoVirSum forces, but we do not calculate the virial,
781 * we sum fr->f_novirum=f later.
783 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
785 wallcycle_start(wcycle,ewcVSITESPREAD);
786 spread_vsite_f(fplog,vsite,x,f,fr->fshift,nrnb,
787 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
788 wallcycle_stop(wcycle,ewcVSITESPREAD);
790 if (bSepLRF)
792 wallcycle_start(wcycle,ewcVSITESPREAD);
793 spread_vsite_f(fplog,vsite,x,fr->f_twin,NULL,
794 nrnb,
795 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
796 wallcycle_stop(wcycle,ewcVSITESPREAD);
800 if (flags & GMX_FORCE_VIRIAL)
802 /* Calculation of the virial must be done after vsites! */
803 calc_virial(fplog,mdatoms->start,mdatoms->homenr,x,f,
804 vir_force,graph,box,nrnb,fr,inputrec->ePBC);
808 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
810 /* Calculate the center of mass forces, this requires communication,
811 * which is why pull_potential is called close to other communication.
812 * The virial contribution is calculated directly,
813 * which is why we call pull_potential after calc_virial.
815 set_pbc(&pbc,inputrec->ePBC,box);
816 dvdl = 0;
817 enerd->term[F_COM_PULL] =
818 pull_potential(inputrec->ePull,inputrec->pull,mdatoms,&pbc,
819 cr,t,lambda,x,f,vir_force,&dvdl);
820 if (bSepDVDL)
822 fprintf(fplog,sepdvdlformat,"Com pull",enerd->term[F_COM_PULL],dvdl);
824 enerd->dvdl_lin += dvdl;
826 else
827 enerd->term[F_COM_PULL] = 0.0;
829 /* Add the forces from enforced rotation potentials (if any) */
830 if (inputrec->bRot)
831 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr,step,t);
834 if (PAR(cr) && !(cr->duty & DUTY_PME))
836 cycles_ppdpme = wallcycle_stop(wcycle,ewcPPDURINGPME);
837 dd_cycles_add(cr->dd,cycles_ppdpme,ddCyclPPduringPME);
839 /* In case of node-splitting, the PP nodes receive the long-range
840 * forces, virial and energy from the PME nodes here.
842 wallcycle_start(wcycle,ewcPP_PMEWAITRECVF);
843 dvdl = 0;
844 gmx_pme_receive_f(cr,fr->f_novirsum,fr->vir_el_recip,&e,&dvdl,
845 &cycles_seppme);
846 if (bSepDVDL)
848 fprintf(fplog,sepdvdlformat,"PME mesh",e,dvdl);
850 enerd->term[F_COUL_RECIP] += e;
851 enerd->dvdl_lin += dvdl;
852 if (wcycle)
854 dd_cycles_add(cr->dd,cycles_seppme,ddCyclPME);
856 wallcycle_stop(wcycle,ewcPP_PMEWAITRECVF);
859 if (bDoForces && fr->bF_NoVirSum)
861 if (vsite)
863 /* Spread the mesh force on virtual sites to the other particles...
864 * This is parallellized. MPI communication is performed
865 * if the constructing atoms aren't local.
867 wallcycle_start(wcycle,ewcVSITESPREAD);
868 spread_vsite_f(fplog,vsite,x,fr->f_novirsum,NULL,nrnb,
869 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
870 wallcycle_stop(wcycle,ewcVSITESPREAD);
872 if (flags & GMX_FORCE_VIRIAL)
874 /* Now add the forces, this is local */
875 if (fr->bDomDec)
877 sum_forces(0,fr->f_novirsum_n,f,fr->f_novirsum);
879 else
881 sum_forces(start,start+homenr,f,fr->f_novirsum);
883 if (EEL_FULL(fr->eeltype))
885 /* Add the mesh contribution to the virial */
886 m_add(vir_force,fr->vir_el_recip,vir_force);
888 if (debug)
890 pr_rvecs(debug,0,"vir_force",vir_force,DIM);
895 /* Sum the potential energy terms from group contributions */
896 sum_epot(&(inputrec->opts),enerd);
898 if (fr->print_force >= 0 && bDoForces)
900 print_large_forces(stderr,mdatoms,cr,step,fr->print_force,x,f);
904 void do_constrain_first(FILE *fplog,gmx_constr_t constr,
905 t_inputrec *ir,t_mdatoms *md,
906 t_state *state,rvec *f,
907 t_graph *graph,t_commrec *cr,t_nrnb *nrnb,
908 t_forcerec *fr, gmx_localtop_t *top, tensor shake_vir)
910 int i,m,start,end;
911 gmx_large_int_t step;
912 double mass,tmass,vcm[4];
913 real dt=ir->delta_t;
914 real dvdlambda;
915 rvec *savex;
917 snew(savex,state->natoms);
919 start = md->start;
920 end = md->homenr + start;
922 if (debug)
923 fprintf(debug,"vcm: start=%d, homenr=%d, end=%d\n",
924 start,md->homenr,end);
925 /* Do a first constrain to reset particles... */
926 step = ir->init_step;
927 if (fplog)
928 fprintf(fplog,"\nConstraining the starting coordinates (step %d)\n",step);
929 dvdlambda = 0;
931 /* constrain the current position */
932 constrain(NULL,TRUE,FALSE,constr,&(top->idef),
933 ir,NULL,cr,step,0,md,
934 state->x,state->x,NULL,
935 state->box,state->lambda,&dvdlambda,
936 NULL,NULL,nrnb,econqCoord,ir->epc==epcMTTK,state->veta,state->veta);
937 if (EI_VV(ir->eI))
939 /* constrain the inital velocity, and save it */
940 /* also may be useful if we need the ekin from the halfstep for velocity verlet */
941 /* might not yet treat veta correctly */
942 constrain(NULL,TRUE,FALSE,constr,&(top->idef),
943 ir,NULL,cr,step,0,md,
944 state->x,state->v,state->v,
945 state->box,state->lambda,&dvdlambda,
946 NULL,NULL,nrnb,econqVeloc,ir->epc==epcMTTK,state->veta,state->veta);
948 /* constrain the inital velocities at t-dt/2 */
949 if (EI_STATE_VELOCITY(ir->eI) && ir->eI!=eiVV)
951 for(i=start; (i<end); i++)
953 for(m=0; (m<DIM); m++)
955 /* Reverse the velocity */
956 state->v[i][m] = -state->v[i][m];
957 /* Store the position at t-dt in buf */
958 savex[i][m] = state->x[i][m] + dt*state->v[i][m];
961 /* Shake the positions at t=-dt with the positions at t=0
962 * as reference coordinates.
964 if (fplog)
966 fprintf(fplog,"\nConstraining the coordinates at t0-dt (step %d)\n",
967 step);
969 dvdlambda = 0;
970 constrain(NULL,TRUE,FALSE,constr,&(top->idef),
971 ir,NULL,cr,step,-1,md,
972 state->x,savex,NULL,
973 state->box,state->lambda,&dvdlambda,
974 state->v,NULL,nrnb,econqCoord,ir->epc==epcMTTK,state->veta,state->veta);
976 for(i=start; i<end; i++) {
977 for(m=0; m<DIM; m++) {
978 /* Re-reverse the velocities */
979 state->v[i][m] = -state->v[i][m];
984 for(m=0; (m<4); m++)
985 vcm[m] = 0;
986 for(i=start; i<end; i++) {
987 mass = md->massT[i];
988 for(m=0; m<DIM; m++) {
989 vcm[m] += state->v[i][m]*mass;
991 vcm[3] += mass;
994 if (ir->nstcomm != 0 || debug) {
995 /* Compute the global sum of vcm */
996 if (debug)
997 fprintf(debug,"vcm: %8.3f %8.3f %8.3f,"
998 " total mass = %12.5e\n",vcm[XX],vcm[YY],vcm[ZZ],vcm[3]);
999 if (PAR(cr))
1000 gmx_sumd(4,vcm,cr);
1001 tmass = vcm[3];
1002 for(m=0; (m<DIM); m++)
1003 vcm[m] /= tmass;
1004 if (debug)
1005 fprintf(debug,"vcm: %8.3f %8.3f %8.3f,"
1006 " total mass = %12.5e\n",vcm[XX],vcm[YY],vcm[ZZ],tmass);
1007 if (ir->nstcomm != 0) {
1008 /* Now we have the velocity of center of mass, let's remove it */
1009 for(i=start; (i<end); i++) {
1010 for(m=0; (m<DIM); m++)
1011 state->v[i][m] -= vcm[m];
1016 sfree(savex);
1019 void calc_enervirdiff(FILE *fplog,int eDispCorr,t_forcerec *fr)
1021 double eners[2],virs[2],enersum,virsum,y0,f,g,h;
1022 double r0,r1,r,rc3,rc9,ea,eb,ec,pa,pb,pc,pd;
1023 double invscale,invscale2,invscale3;
1024 int ri0,ri1,ri,i,offstart,offset;
1025 real scale,*vdwtab;
1027 fr->enershiftsix = 0;
1028 fr->enershifttwelve = 0;
1029 fr->enerdiffsix = 0;
1030 fr->enerdifftwelve = 0;
1031 fr->virdiffsix = 0;
1032 fr->virdifftwelve = 0;
1034 if (eDispCorr != edispcNO) {
1035 for(i=0; i<2; i++) {
1036 eners[i] = 0;
1037 virs[i] = 0;
1039 if ((fr->vdwtype == evdwSWITCH) || (fr->vdwtype == evdwSHIFT)) {
1040 if (fr->rvdw_switch == 0)
1041 gmx_fatal(FARGS,
1042 "With dispersion correction rvdw-switch can not be zero "
1043 "for vdw-type = %s",evdw_names[fr->vdwtype]);
1045 scale = fr->nblists[0].tab.scale;
1046 vdwtab = fr->nblists[0].vdwtab;
1048 /* Round the cut-offs to exact table values for precision */
1049 ri0 = floor(fr->rvdw_switch*scale);
1050 ri1 = ceil(fr->rvdw*scale);
1051 r0 = ri0/scale;
1052 r1 = ri1/scale;
1053 rc3 = r0*r0*r0;
1054 rc9 = rc3*rc3*rc3;
1056 if (fr->vdwtype == evdwSHIFT) {
1057 /* Determine the constant energy shift below rvdw_switch */
1058 fr->enershiftsix = (real)(-1.0/(rc3*rc3)) - vdwtab[8*ri0];
1059 fr->enershifttwelve = (real)( 1.0/(rc9*rc3)) - vdwtab[8*ri0 + 4];
1061 /* Add the constant part from 0 to rvdw_switch.
1062 * This integration from 0 to rvdw_switch overcounts the number
1063 * of interactions by 1, as it also counts the self interaction.
1064 * We will correct for this later.
1066 eners[0] += 4.0*M_PI*fr->enershiftsix*rc3/3.0;
1067 eners[1] += 4.0*M_PI*fr->enershifttwelve*rc3/3.0;
1069 invscale = 1.0/(scale);
1070 invscale2 = invscale*invscale;
1071 invscale3 = invscale*invscale2;
1073 /* following summation derived from cubic spline definition,
1074 Numerical Recipies in C, second edition, p. 113-116. Exact
1075 for the cubic spline. We first calculate the negative of
1076 the energy from rvdw to rvdw_switch, assuming that g(r)=1,
1077 and then add the more standard, abrupt cutoff correction to
1078 that result, yielding the long-range correction for a
1079 switched function. We perform both the pressure and energy
1080 loops at the same time for simplicity, as the computational
1081 cost is low. */
1083 for (i=0;i<2;i++) {
1084 enersum = 0.0; virsum = 0.0;
1085 if (i==0)
1086 offstart = 0;
1087 else
1088 offstart = 4;
1089 for (ri=ri0; ri<ri1; ri++) {
1090 r = ri*invscale;
1091 ea = invscale3;
1092 eb = 2.0*invscale2*r;
1093 ec = invscale*r*r;
1095 pa = invscale3;
1096 pb = 3.0*invscale2*r;
1097 pc = 3.0*invscale*r*r;
1098 pd = r*r*r;
1100 /* this "8" is from the packing in the vdwtab array - perhaps
1101 should be #define'ed? */
1102 offset = 8*ri + offstart;
1103 y0 = vdwtab[offset];
1104 f = vdwtab[offset+1];
1105 g = vdwtab[offset+2];
1106 h = vdwtab[offset+3];
1108 enersum += y0*(ea/3 + eb/2 + ec) + f*(ea/4 + eb/3 + ec/2)+
1109 g*(ea/5 + eb/4 + ec/3) + h*(ea/6 + eb/5 + ec/4);
1110 virsum += f*(pa/4 + pb/3 + pc/2 + pd) +
1111 2*g*(pa/5 + pb/4 + pc/3 + pd/2) + 3*h*(pa/6 + pb/5 + pc/4 + pd/3);
1114 enersum *= 4.0*M_PI;
1115 virsum *= 4.0*M_PI;
1116 eners[i] -= enersum;
1117 virs[i] -= virsum;
1120 /* now add the correction for rvdw_switch to infinity */
1121 eners[0] += -4.0*M_PI/(3.0*rc3);
1122 eners[1] += 4.0*M_PI/(9.0*rc9);
1123 virs[0] += 8.0*M_PI/rc3;
1124 virs[1] += -16.0*M_PI/(3.0*rc9);
1126 else if ((fr->vdwtype == evdwCUT) || (fr->vdwtype == evdwUSER)) {
1127 if (fr->vdwtype == evdwUSER && fplog)
1128 fprintf(fplog,
1129 "WARNING: using dispersion correction with user tables\n");
1130 rc3 = fr->rvdw*fr->rvdw*fr->rvdw;
1131 rc9 = rc3*rc3*rc3;
1132 eners[0] += -4.0*M_PI/(3.0*rc3);
1133 eners[1] += 4.0*M_PI/(9.0*rc9);
1134 virs[0] += 8.0*M_PI/rc3;
1135 virs[1] += -16.0*M_PI/(3.0*rc9);
1136 } else {
1137 gmx_fatal(FARGS,
1138 "Dispersion correction is not implemented for vdw-type = %s",
1139 evdw_names[fr->vdwtype]);
1141 fr->enerdiffsix = eners[0];
1142 fr->enerdifftwelve = eners[1];
1143 /* The 0.5 is due to the Gromacs definition of the virial */
1144 fr->virdiffsix = 0.5*virs[0];
1145 fr->virdifftwelve = 0.5*virs[1];
1149 void calc_dispcorr(FILE *fplog,t_inputrec *ir,t_forcerec *fr,
1150 gmx_large_int_t step,int natoms,
1151 matrix box,real lambda,tensor pres,tensor virial,
1152 real *prescorr, real *enercorr, real *dvdlcorr)
1154 bool bCorrAll,bCorrPres;
1155 real dvdlambda,invvol,dens,ninter,avcsix,avctwelve,enerdiff,svir=0,spres=0;
1156 int m;
1158 *prescorr = 0;
1159 *enercorr = 0;
1160 *dvdlcorr = 0;
1162 clear_mat(virial);
1163 clear_mat(pres);
1165 if (ir->eDispCorr != edispcNO) {
1166 bCorrAll = (ir->eDispCorr == edispcAllEner ||
1167 ir->eDispCorr == edispcAllEnerPres);
1168 bCorrPres = (ir->eDispCorr == edispcEnerPres ||
1169 ir->eDispCorr == edispcAllEnerPres);
1171 invvol = 1/det(box);
1172 if (fr->n_tpi)
1174 /* Only correct for the interactions with the inserted molecule */
1175 dens = (natoms - fr->n_tpi)*invvol;
1176 ninter = fr->n_tpi;
1178 else
1180 dens = natoms*invvol;
1181 ninter = 0.5*natoms;
1184 if (ir->efep == efepNO)
1186 avcsix = fr->avcsix[0];
1187 avctwelve = fr->avctwelve[0];
1189 else
1191 avcsix = (1 - lambda)*fr->avcsix[0] + lambda*fr->avcsix[1];
1192 avctwelve = (1 - lambda)*fr->avctwelve[0] + lambda*fr->avctwelve[1];
1195 enerdiff = ninter*(dens*fr->enerdiffsix - fr->enershiftsix);
1196 *enercorr += avcsix*enerdiff;
1197 dvdlambda = 0.0;
1198 if (ir->efep != efepNO)
1200 dvdlambda += (fr->avcsix[1] - fr->avcsix[0])*enerdiff;
1202 if (bCorrAll)
1204 enerdiff = ninter*(dens*fr->enerdifftwelve - fr->enershifttwelve);
1205 *enercorr += avctwelve*enerdiff;
1206 if (fr->efep != efepNO)
1208 dvdlambda += (fr->avctwelve[1] - fr->avctwelve[0])*enerdiff;
1212 if (bCorrPres)
1214 svir = ninter*dens*avcsix*fr->virdiffsix/3.0;
1215 if (ir->eDispCorr == edispcAllEnerPres)
1217 svir += ninter*dens*avctwelve*fr->virdifftwelve/3.0;
1219 /* The factor 2 is because of the Gromacs virial definition */
1220 spres = -2.0*invvol*svir*PRESFAC;
1222 for(m=0; m<DIM; m++) {
1223 virial[m][m] += svir;
1224 pres[m][m] += spres;
1226 *prescorr += spres;
1229 /* Can't currently control when it prints, for now, just print when degugging */
1230 if (debug)
1232 if (bCorrAll) {
1233 fprintf(debug,"Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
1234 avcsix,avctwelve);
1236 if (bCorrPres)
1238 fprintf(debug,
1239 "Long Range LJ corr.: Epot %10g, Pres: %10g, Vir: %10g\n",
1240 *enercorr,spres,svir);
1242 else
1244 fprintf(debug,"Long Range LJ corr.: Epot %10g\n",*enercorr);
1248 if (fr->bSepDVDL && do_per_step(step,ir->nstlog))
1249 fprintf(fplog,sepdvdlformat,"Dispersion correction",
1250 *enercorr,dvdlambda);
1252 if (fr->efep != efepNO)
1254 *dvdlcorr += dvdlambda;
1259 void do_pbc_first(FILE *fplog,matrix box,t_forcerec *fr,
1260 t_graph *graph,rvec x[])
1262 if (fplog)
1263 fprintf(fplog,"Removing pbc first time\n");
1264 calc_shifts(box,fr->shift_vec);
1265 if (graph) {
1266 mk_mshift(fplog,graph,fr->ePBC,box,x);
1267 if (gmx_debug_at)
1268 p_graph(debug,"do_pbc_first 1",graph);
1269 shift_self(graph,box,x);
1270 /* By doing an extra mk_mshift the molecules that are broken
1271 * because they were e.g. imported from another software
1272 * will be made whole again. Such are the healing powers
1273 * of GROMACS.
1275 mk_mshift(fplog,graph,fr->ePBC,box,x);
1276 if (gmx_debug_at)
1277 p_graph(debug,"do_pbc_first 2",graph);
1279 if (fplog)
1280 fprintf(fplog,"Done rmpbc\n");
1283 static void low_do_pbc_mtop(FILE *fplog,int ePBC,matrix box,
1284 gmx_mtop_t *mtop,rvec x[],
1285 bool bFirst)
1287 t_graph *graph;
1288 int mb,as,mol;
1289 gmx_molblock_t *molb;
1291 if (bFirst && fplog)
1292 fprintf(fplog,"Removing pbc first time\n");
1294 snew(graph,1);
1295 as = 0;
1296 for(mb=0; mb<mtop->nmolblock; mb++) {
1297 molb = &mtop->molblock[mb];
1298 if (molb->natoms_mol == 1 ||
1299 (!bFirst && mtop->moltype[molb->type].cgs.nr == 1)) {
1300 /* Just one atom or charge group in the molecule, no PBC required */
1301 as += molb->nmol*molb->natoms_mol;
1302 } else {
1303 /* Pass NULL iso fplog to avoid graph prints for each molecule type */
1304 mk_graph_ilist(NULL,mtop->moltype[molb->type].ilist,
1305 0,molb->natoms_mol,FALSE,FALSE,graph);
1307 for(mol=0; mol<molb->nmol; mol++) {
1308 mk_mshift(fplog,graph,ePBC,box,x+as);
1310 shift_self(graph,box,x+as);
1311 /* The molecule is whole now.
1312 * We don't need the second mk_mshift call as in do_pbc_first,
1313 * since we no longer need this graph.
1316 as += molb->natoms_mol;
1318 done_graph(graph);
1321 sfree(graph);
1324 void do_pbc_first_mtop(FILE *fplog,int ePBC,matrix box,
1325 gmx_mtop_t *mtop,rvec x[])
1327 low_do_pbc_mtop(fplog,ePBC,box,mtop,x,TRUE);
1330 void do_pbc_mtop(FILE *fplog,int ePBC,matrix box,
1331 gmx_mtop_t *mtop,rvec x[])
1333 low_do_pbc_mtop(fplog,ePBC,box,mtop,x,FALSE);
1336 void finish_run(FILE *fplog,t_commrec *cr,const char *confout,
1337 t_inputrec *inputrec,
1338 t_nrnb nrnb[],gmx_wallcycle_t wcycle,
1339 gmx_runtime_t *runtime,
1340 bool bWriteStat)
1342 int i,j;
1343 t_nrnb *nrnb_tot=NULL;
1344 real delta_t;
1345 double nbfs,mflop;
1346 double cycles[ewcNR];
1348 wallcycle_sum(cr,wcycle,cycles);
1350 if (cr->nnodes > 1) {
1351 if (SIMMASTER(cr))
1352 snew(nrnb_tot,1);
1353 #ifdef GMX_MPI
1354 MPI_Reduce(nrnb->n,nrnb_tot->n,eNRNB,MPI_DOUBLE,MPI_SUM,
1355 MASTERRANK(cr),cr->mpi_comm_mysim);
1356 #endif
1357 } else {
1358 nrnb_tot = nrnb;
1361 if (SIMMASTER(cr)) {
1362 print_flop(fplog,nrnb_tot,&nbfs,&mflop);
1363 if (cr->nnodes > 1) {
1364 sfree(nrnb_tot);
1368 if ((cr->duty & DUTY_PP) && DOMAINDECOMP(cr)) {
1369 print_dd_statistics(cr,inputrec,fplog);
1372 if (SIMMASTER(cr)) {
1373 if (PARTDECOMP(cr)) {
1374 pr_load(fplog,cr,nrnb_tot);
1377 wallcycle_print(fplog,cr->nnodes,cr->npmenodes,runtime->realtime,
1378 wcycle,cycles);
1380 if (EI_DYNAMICS(inputrec->eI)) {
1381 delta_t = inputrec->delta_t;
1382 } else {
1383 delta_t = 0;
1386 if (fplog) {
1387 print_perf(fplog,runtime->proctime,runtime->realtime,
1388 cr->nnodes-cr->npmenodes,
1389 runtime->nsteps_done,delta_t,nbfs,mflop);
1391 if (bWriteStat) {
1392 print_perf(stderr,runtime->proctime,runtime->realtime,
1393 cr->nnodes-cr->npmenodes,
1394 runtime->nsteps_done,delta_t,nbfs,mflop);
1398 runtime=inputrec->nsteps*inputrec->delta_t;
1399 if (bWriteStat) {
1400 if (cr->nnodes == 1)
1401 fprintf(stderr,"\n\n");
1402 print_perf(stderr,nodetime,realtime,runtime,&ntot,
1403 cr->nnodes-cr->npmenodes,FALSE);
1405 wallcycle_print(fplog,cr->nnodes,cr->npmenodes,realtime,wcycle,cycles);
1406 print_perf(fplog,nodetime,realtime,runtime,&ntot,cr->nnodes-cr->npmenodes,
1407 TRUE);
1408 if (PARTDECOMP(cr))
1409 pr_load(fplog,cr,nrnb_all);
1410 if (cr->nnodes > 1)
1411 sfree(nrnb_all);
1416 void init_md(FILE *fplog,
1417 t_commrec *cr,t_inputrec *ir,const output_env_t oenv,
1418 double *t,double *t0,
1419 real *lambda,double *lam0,
1420 t_nrnb *nrnb,gmx_mtop_t *mtop,
1421 gmx_update_t *upd,
1422 int nfile,const t_filenm fnm[],
1423 int *fp_trn,int *fp_xtc,ener_file_t *fp_ene,const char **fn_cpt,
1424 FILE **fp_dhdl,FILE **fp_field,
1425 t_mdebin **mdebin,
1426 tensor force_vir,tensor shake_vir,rvec mu_tot,
1427 bool *bNEMD,bool *bSimAnn,t_vcm **vcm, t_state *state, unsigned long Flags)
1429 int i,j,n;
1430 real tmpt,mod;
1431 char filemode[3];
1433 sprintf(filemode, (Flags & MD_APPENDFILES) ? "a+" : "w+");
1435 /* Initial values */
1436 *t = *t0 = ir->init_t;
1437 if (ir->efep != efepNO)
1439 *lam0 = ir->init_lambda;
1440 *lambda = *lam0 + ir->init_step*ir->delta_lambda;
1442 else
1444 *lambda = *lam0 = 0.0;
1447 *bSimAnn=FALSE;
1448 for(i=0;i<ir->opts.ngtc;i++)
1450 /* set bSimAnn if any group is being annealed */
1451 if(ir->opts.annealing[i]!=eannNO)
1453 *bSimAnn = TRUE;
1456 if (*bSimAnn)
1458 update_annealing_target_temp(&(ir->opts),ir->init_t);
1461 *bNEMD = (ir->opts.ngacc > 1) || (norm(ir->opts.acc[0]) > 0);
1463 if (upd)
1465 *upd = init_update(fplog,ir);
1468 if (vcm != NULL)
1470 *vcm = init_vcm(fplog,&mtop->groups,ir);
1473 if (EI_DYNAMICS(ir->eI) && !(Flags & MD_APPENDFILES))
1475 if (ir->etc == etcBERENDSEN)
1477 please_cite(fplog,"Berendsen84a");
1479 if (ir->etc == etcVRESCALE)
1481 please_cite(fplog,"Bussi2007a");
1485 init_nrnb(nrnb);
1487 if (nfile != -1)
1489 *fp_trn = -1;
1490 *fp_ene = NULL;
1491 *fp_xtc = -1;
1493 if (MASTER(cr))
1495 *fp_trn = open_trn(ftp2fn(efTRN,nfile,fnm), filemode);
1496 if (ir->nstxtcout > 0)
1498 *fp_xtc = open_xtc(ftp2fn(efXTC,nfile,fnm), filemode);
1500 *fp_ene = open_enx(ftp2fn(efEDR,nfile,fnm), filemode);
1501 *fn_cpt = opt2fn("-cpo",nfile,fnm);
1503 if ((fp_dhdl != NULL) && ir->efep != efepNO && ir->nstdhdl > 0)
1505 if(Flags & MD_APPENDFILES)
1507 *fp_dhdl= gmx_fio_fopen(opt2fn("-dhdl",nfile,fnm),filemode);
1509 else
1511 *fp_dhdl = open_dhdl(opt2fn("-dhdl",nfile,fnm),ir,oenv);
1515 if ((fp_field != NULL) &&
1516 (ir->ex[XX].n || ir->ex[YY].n ||ir->ex[ZZ].n))
1518 if(Flags & MD_APPENDFILES)
1520 *fp_dhdl=gmx_fio_fopen(opt2fn("-field",nfile,fnm),filemode);
1522 else
1524 *fp_field = xvgropen(opt2fn("-field",nfile,fnm),
1525 "Applied electric field","Time (ps)",
1526 "E (V/nm)",oenv);
1530 *mdebin = init_mdebin( (Flags & MD_APPENDFILES) ? NULL : *fp_ene,
1531 mtop,ir);
1534 /* Initiate variables */
1535 clear_mat(force_vir);
1536 clear_mat(shake_vir);
1537 clear_rvec(mu_tot);
1539 debug_gmx();