Merge branch 'master' into rotation
[gromacs/adressmacs.git] / src / mdlib / sim_util.c
blob3b343969f68890520688b2f83235182a6994ae62
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
97 #include "qmmm.h"
99 #if 0
100 typedef struct gmx_timeprint {
102 } t_gmx_timeprint;
103 #endif
106 double
107 gmx_gettime()
109 #ifdef HAVE_GETTIMEOFDAY
110 struct timeval t;
111 struct timezone tz = { 0,0 };
112 double seconds;
114 gettimeofday(&t,&tz);
116 seconds = (double) t.tv_sec + 1e-6*(double)t.tv_usec;
118 return seconds;
119 #else
120 double seconds;
122 seconds = time(NULL);
124 return seconds;
125 #endif
129 #define difftime(end,start) ((double)(end)-(double)(start))
131 void print_time(FILE *out,gmx_runtime_t *runtime,gmx_large_int_t step,t_inputrec *ir)
133 time_t finish;
135 double dt;
136 char buf[48];
138 if (!gmx_parallel_env())
139 fprintf(out,"\r");
140 fprintf(out,"step %s",gmx_step_str(step,buf));
141 if ((step >= ir->nstlist)) {
142 if ((ir->nstlist == 0) || ((step % ir->nstlist) == 0)) {
143 /* We have done a full cycle let's update time_per_step */
144 runtime->last = gmx_gettime();
145 dt = difftime(runtime->last,runtime->real);
146 runtime->time_per_step = dt/(step - ir->init_step + 1);
148 dt = (ir->nsteps + ir->init_step - step)*runtime->time_per_step;
150 if (dt >= 300) {
151 finish = (time_t) (runtime->last + dt);
152 sprintf(buf,"%s",ctime(&finish));
153 buf[strlen(buf)-1]='\0';
154 fprintf(out,", will finish %s",buf);
156 else
157 fprintf(out,", remaining runtime: %5d s ",(int)dt);
159 if (gmx_parallel_env())
160 fprintf(out,"\n");
162 fflush(out);
165 #ifdef NO_CLOCK
166 #define clock() -1
167 #endif
169 static double set_proctime(gmx_runtime_t *runtime)
171 double diff;
172 #ifdef GMX_CRAY_XT3
173 double prev;
175 prev = runtime->proc;
176 runtime->proc = dclock();
178 diff = runtime->proc - prev;
179 #else
180 clock_t prev;
182 prev = runtime->proc;
183 runtime->proc = clock();
185 diff = (double)(runtime->proc - prev)/(double)CLOCKS_PER_SEC;
186 #endif
187 if (diff < 0)
189 /* The counter has probably looped, ignore this data */
190 diff = 0;
193 return diff;
196 void runtime_start(gmx_runtime_t *runtime)
198 runtime->real = gmx_gettime();
199 runtime->proc = 0;
200 set_proctime(runtime);
201 runtime->realtime = 0;
202 runtime->proctime = 0;
203 runtime->last = 0;
204 runtime->time_per_step = 0;
207 void runtime_end(gmx_runtime_t *runtime)
209 double now;
211 now = gmx_gettime();
213 runtime->proctime += set_proctime(runtime);
214 runtime->realtime = now - runtime->real;
215 runtime->real = now;
218 void runtime_upd_proc(gmx_runtime_t *runtime)
220 runtime->proctime += set_proctime(runtime);
223 void print_date_and_time(FILE *fplog,int nodeid,const char *title,
224 const gmx_runtime_t *runtime)
226 int i;
227 char *ts,time_string[STRLEN];
228 time_t tmptime;
230 if (runtime != NULL)
232 tmptime = (time_t) runtime->real;
233 ts = ctime(&tmptime);
235 else
237 tmptime = (time_t) gmx_gettime();
238 ts = ctime(&tmptime);
240 for(i=0; ts[i]>=' '; i++)
242 time_string[i]=ts[i];
244 time_string[i]='\0';
245 if (fplog)
247 fprintf(fplog,"%s on node %d %s\n",title,nodeid,time_string);
251 static void sum_forces(int start,int end,rvec f[],rvec flr[])
253 int i;
255 if (gmx_debug_at) {
256 pr_rvecs(debug,0,"fsr",f+start,end-start);
257 pr_rvecs(debug,0,"flr",flr+start,end-start);
259 for(i=start; (i<end); i++)
260 rvec_inc(f[i],flr[i]);
264 * calc_f_el calculates forces due to an electric field.
266 * force is kJ mol^-1 nm^-1 = e * kJ mol^-1 nm^-1 / e
268 * Et[] contains the parameters for the time dependent
269 * part of the field (not yet used).
270 * Ex[] contains the parameters for
271 * the spatial dependent part of the field. You can have cool periodic
272 * fields in principle, but only a constant field is supported
273 * now.
274 * The function should return the energy due to the electric field
275 * (if any) but for now returns 0.
277 * WARNING:
278 * There can be problems with the virial.
279 * Since the field is not self-consistent this is unavoidable.
280 * For neutral molecules the virial is correct within this approximation.
281 * For neutral systems with many charged molecules the error is small.
282 * But for systems with a net charge or a few charged molecules
283 * the error can be significant when the field is high.
284 * Solution: implement a self-consitent electric field into PME.
286 static void calc_f_el(FILE *fp,int start,int homenr,
287 real charge[],rvec x[],rvec f[],
288 t_cosines Ex[],t_cosines Et[],double t)
290 rvec Ext;
291 real t0;
292 int i,m;
294 for(m=0; (m<DIM); m++)
296 if (Et[m].n > 0)
298 if (Et[m].n == 3)
300 t0 = Et[m].a[1];
301 Ext[m] = cos(Et[m].a[0]*(t-t0))*exp(-sqr(t-t0)/(2.0*sqr(Et[m].a[2])));
303 else
305 Ext[m] = cos(Et[m].a[0]*t);
308 else
310 Ext[m] = 1.0;
312 if (Ex[m].n > 0)
314 /* Convert the field strength from V/nm to MD-units */
315 Ext[m] *= Ex[m].a[0]*FIELDFAC;
316 for(i=start; (i<start+homenr); i++)
317 f[i][m] += charge[i]*Ext[m];
319 else
321 Ext[m] = 0;
324 if (fp != NULL)
326 fprintf(fp,"%10g %10g %10g %10g #FIELD\n",t,
327 Ext[XX]/FIELDFAC,Ext[YY]/FIELDFAC,Ext[ZZ]/FIELDFAC);
331 static void calc_virial(FILE *fplog,int start,int homenr,rvec x[],rvec f[],
332 tensor vir_part,t_graph *graph,matrix box,
333 t_nrnb *nrnb,const t_forcerec *fr,int ePBC)
335 int i,j;
336 tensor virtest;
338 /* The short-range virial from surrounding boxes */
339 clear_mat(vir_part);
340 calc_vir(fplog,SHIFTS,fr->shift_vec,fr->fshift,vir_part,ePBC==epbcSCREW,box);
341 inc_nrnb(nrnb,eNR_VIRIAL,SHIFTS);
343 /* Calculate partial virial, for local atoms only, based on short range.
344 * Total virial is computed in global_stat, called from do_md
346 f_calc_vir(fplog,start,start+homenr,x,f,vir_part,graph,box);
347 inc_nrnb(nrnb,eNR_VIRIAL,homenr);
349 /* Add position restraint contribution */
350 for(i=0; i<DIM; i++) {
351 vir_part[i][i] += fr->vir_diag_posres[i];
354 /* Add wall contribution */
355 for(i=0; i<DIM; i++) {
356 vir_part[i][ZZ] += fr->vir_wall_z[i];
359 if (debug)
360 pr_rvecs(debug,0,"vir_part",vir_part,DIM);
363 static void print_large_forces(FILE *fp,t_mdatoms *md,t_commrec *cr,
364 gmx_large_int_t step,real pforce,rvec *x,rvec *f)
366 int i;
367 real pf2,fn2;
368 char buf[22];
370 pf2 = sqr(pforce);
371 for(i=md->start; i<md->start+md->homenr; i++) {
372 fn2 = norm2(f[i]);
373 if (fn2 >= pf2) {
374 fprintf(fp,"step %s atom %6d x %8.3f %8.3f %8.3f force %12.5e\n",
375 gmx_step_str(step,buf),
376 ddglatnr(cr->dd,i),x[i][XX],x[i][YY],x[i][ZZ],sqrt(fn2));
381 void do_force(FILE *fplog,t_commrec *cr,
382 t_inputrec *inputrec,
383 gmx_large_int_t step,t_nrnb *nrnb,gmx_wallcycle_t wcycle,
384 gmx_localtop_t *top,
385 gmx_mtop_t *mtop,
386 gmx_groups_t *groups,
387 matrix box,rvec x[],history_t *hist,
388 rvec f[],
389 tensor vir_force,
390 t_mdatoms *mdatoms,
391 gmx_enerdata_t *enerd,t_fcdata *fcd,
392 real lambda,t_graph *graph,
393 t_forcerec *fr,gmx_vsite_t *vsite,rvec mu_tot,
394 double t,FILE *field,gmx_edsam_t ed,
395 bool bBornRadii,
396 int flags)
398 int cg0,cg1,i,j;
399 int start,homenr;
400 double mu[2*DIM];
401 bool bSepDVDL,bStateChanged,bNS,bFillGrid,bCalcCGCM,bBS;
402 bool bDoLongRange,bDoForces,bSepLRF;
403 matrix boxs;
404 real e,v,dvdl;
405 t_pbc pbc;
406 float cycles_ppdpme,cycles_pme,cycles_seppme,cycles_force;
408 start = mdatoms->start;
409 homenr = mdatoms->homenr;
411 bSepDVDL = (fr->bSepDVDL && do_per_step(step,inputrec->nstlog));
413 clear_mat(vir_force);
415 if (PARTDECOMP(cr))
417 pd_cg_range(cr,&cg0,&cg1);
419 else
421 cg0 = 0;
422 if (DOMAINDECOMP(cr))
424 cg1 = cr->dd->ncg_tot;
426 else
428 cg1 = top->cgs.nr;
430 if (fr->n_tpi > 0)
432 cg1--;
436 bStateChanged = (flags & GMX_FORCE_STATECHANGED);
437 bNS = (flags & GMX_FORCE_NS) && (fr->bAllvsAll==FALSE);
438 bFillGrid = (bNS && bStateChanged);
439 bCalcCGCM = (bFillGrid && !DOMAINDECOMP(cr));
440 bDoLongRange = (fr->bTwinRange && bNS && (flags & GMX_FORCE_DOLR));
441 bDoForces = (flags & GMX_FORCE_FORCES);
442 bSepLRF = (bDoLongRange && bDoForces && (flags & GMX_FORCE_SEPLRF));
444 if (bStateChanged)
446 update_forcerec(fplog,fr,box);
448 /* Calculate total (local) dipole moment in a temporary common array.
449 * This makes it possible to sum them over nodes faster.
451 calc_mu(start,homenr,
452 x,mdatoms->chargeA,mdatoms->chargeB,mdatoms->nChargePerturbed,
453 mu,mu+DIM);
456 if (fr->ePBC != epbcNONE) {
457 /* Compute shift vectors every step,
458 * because of pressure coupling or box deformation!
460 if (DYNAMIC_BOX(*inputrec) && bStateChanged)
461 calc_shifts(box,fr->shift_vec);
463 if (bCalcCGCM) {
464 put_charge_groups_in_box(fplog,cg0,cg1,fr->ePBC,box,
465 &(top->cgs),x,fr->cg_cm);
466 inc_nrnb(nrnb,eNR_CGCM,homenr);
467 inc_nrnb(nrnb,eNR_RESETX,cg1-cg0);
469 else if (EI_ENERGY_MINIMIZATION(inputrec->eI) && graph) {
470 unshift_self(graph,box,x);
473 else if (bCalcCGCM) {
474 calc_cgcm(fplog,cg0,cg1,&(top->cgs),x,fr->cg_cm);
475 inc_nrnb(nrnb,eNR_CGCM,homenr);
478 if (bCalcCGCM) {
479 if (PAR(cr)) {
480 move_cgcm(fplog,cr,fr->cg_cm);
482 if (gmx_debug_at)
483 pr_rvecs(debug,0,"cgcm",fr->cg_cm,top->cgs.nr);
486 #ifdef GMX_MPI
487 if (!(cr->duty & DUTY_PME)) {
488 /* Send particle coordinates to the pme nodes.
489 * Since this is only implemented for domain decomposition
490 * and domain decomposition does not use the graph,
491 * we do not need to worry about shifting.
494 wallcycle_start(wcycle,ewcPP_PMESENDX);
495 GMX_MPE_LOG(ev_send_coordinates_start);
497 bBS = (inputrec->nwall == 2);
498 if (bBS) {
499 copy_mat(box,boxs);
500 svmul(inputrec->wall_ewald_zfac,boxs[ZZ],boxs[ZZ]);
503 gmx_pme_send_x(cr,bBS ? boxs : box,x,mdatoms->nChargePerturbed,lambda,step);
505 GMX_MPE_LOG(ev_send_coordinates_finish);
506 wallcycle_stop(wcycle,ewcPP_PMESENDX);
508 #endif /* GMX_MPI */
510 /* Communicate coordinates and sum dipole if necessary */
511 if (PAR(cr))
513 wallcycle_start(wcycle,ewcMOVEX);
514 if (DOMAINDECOMP(cr))
516 dd_move_x(cr->dd,box,x);
518 else
520 move_x(fplog,cr,GMX_LEFT,GMX_RIGHT,x,nrnb);
522 /* When we don't need the total dipole we sum it in global_stat */
523 if (bStateChanged && NEED_MUTOT(*inputrec))
525 gmx_sumd(2*DIM,mu,cr);
527 wallcycle_stop(wcycle,ewcMOVEX);
529 if (bStateChanged)
531 for(i=0; i<2; i++)
533 for(j=0;j<DIM;j++)
535 fr->mu_tot[i][j] = mu[i*DIM + j];
539 if (fr->efep == efepNO)
541 copy_rvec(fr->mu_tot[0],mu_tot);
543 else
545 for(j=0; j<DIM; j++)
547 mu_tot[j] =
548 (1.0 - lambda)*fr->mu_tot[0][j] + lambda*fr->mu_tot[1][j];
552 /* Reset energies */
553 reset_enerdata(&(inputrec->opts),fr,bNS,enerd,MASTER(cr));
554 clear_rvecs(SHIFTS,fr->fshift);
556 if (bNS)
558 wallcycle_start(wcycle,ewcNS);
560 if (graph && bStateChanged)
562 /* Calculate intramolecular shift vectors to make molecules whole */
563 mk_mshift(fplog,graph,fr->ePBC,box,x);
566 /* Reset long range forces if necessary */
567 if (fr->bTwinRange)
569 /* Reset the (long-range) forces if necessary */
570 clear_rvecs(fr->natoms_force,bSepLRF ? fr->f_twin : f);
573 /* Do the actual neighbour searching and if twin range electrostatics
574 * also do the calculation of long range forces and energies.
576 dvdl = 0;
577 ns(fplog,fr,x,box,
578 groups,&(inputrec->opts),top,mdatoms,
579 cr,nrnb,lambda,&dvdl,&enerd->grpp,bFillGrid,
580 bDoLongRange,bDoForces,bSepLRF ? fr->f_twin : f);
581 if (bSepDVDL)
583 fprintf(fplog,sepdvdlformat,"LR non-bonded",0.0,dvdl);
585 enerd->dvdl_lin += dvdl;
587 wallcycle_stop(wcycle,ewcNS);
590 if (inputrec->implicit_solvent && bNS)
592 make_gb_nblist(cr,mtop->natoms,inputrec->gb_algorithm,inputrec->rlist,
593 x,box,fr,&top->idef,graph,fr->born);
596 if (DOMAINDECOMP(cr))
598 if (!(cr->duty & DUTY_PME))
600 wallcycle_start(wcycle,ewcPPDURINGPME);
601 dd_force_flop_start(cr->dd,nrnb);
605 if (inputrec->bRot)
607 /* Enforced rotation has its own cycle counter that starts after the collective
608 * coordinates have been communicated. It is added to ddCyclF */
609 do_rotation(cr,inputrec,box,x,t,step,wcycle,bNS);
612 /* Start the force cycle counter.
613 * This counter is stopped in do_forcelow_level.
614 * No parallel communication should occur while this counter is running,
615 * since that will interfere with the dynamic load balancing.
617 wallcycle_start(wcycle,ewcFORCE);
618 GMX_MPE_LOG(ev_forcecycles_start);
620 if (bDoForces)
622 /* Reset forces for which the virial is calculated separately:
623 * PME/Ewald forces if necessary */
624 if (fr->bF_NoVirSum)
626 if (flags & GMX_FORCE_VIRIAL)
628 fr->f_novirsum = fr->f_novirsum_alloc;
629 GMX_BARRIER(cr->mpi_comm_mygroup);
630 if (fr->bDomDec)
632 clear_rvecs(fr->f_novirsum_n,fr->f_novirsum);
634 else
636 clear_rvecs(homenr,fr->f_novirsum+start);
638 GMX_BARRIER(cr->mpi_comm_mygroup);
640 else
642 /* We are not calculating the pressure so we do not need
643 * a separate array for forces that do not contribute
644 * to the pressure.
646 fr->f_novirsum = f;
650 if (bSepLRF)
652 /* Add the long range forces to the short range forces */
653 for(i=0; i<fr->natoms_force; i++)
655 copy_rvec(fr->f_twin[i],f[i]);
658 else if (!(fr->bTwinRange && bNS))
660 /* Clear the short-range forces */
661 clear_rvecs(fr->natoms_force,f);
664 clear_rvec(fr->vir_diag_posres);
666 GMX_BARRIER(cr->mpi_comm_mygroup);
668 if (inputrec->ePull == epullCONSTRAINT)
670 clear_pull_forces(inputrec->pull);
673 /* update QMMMrec, if necessary */
674 if(fr->bQMMM)
676 update_QMMMrec(cr,fr,x,mdatoms,box,top);
679 if ((flags & GMX_FORCE_BONDED) && top->idef.il[F_POSRES].nr > 0)
681 /* Position restraints always require full pbc */
682 set_pbc(&pbc,inputrec->ePBC,box);
683 v = posres(top->idef.il[F_POSRES].nr,top->idef.il[F_POSRES].iatoms,
684 top->idef.iparams_posres,
685 (const rvec*)x,fr->f_novirsum,fr->vir_diag_posres,
686 inputrec->ePBC==epbcNONE ? NULL : &pbc,lambda,&dvdl,
687 fr->rc_scaling,fr->ePBC,fr->posres_com,fr->posres_comB);
688 if (bSepDVDL)
690 fprintf(fplog,sepdvdlformat,
691 interaction_function[F_POSRES].longname,v,dvdl);
693 enerd->term[F_POSRES] += v;
694 /* This linear lambda dependence assumption is only correct
695 * when only k depends on lambda,
696 * not when the reference position depends on lambda.
697 * grompp checks for this.
699 enerd->dvdl_lin += dvdl;
700 inc_nrnb(nrnb,eNR_POSRES,top->idef.il[F_POSRES].nr/2);
703 /* Compute the bonded and non-bonded energies and optionally forces */
704 do_force_lowlevel(fplog,step,fr,inputrec,&(top->idef),
705 cr,nrnb,wcycle,mdatoms,&(inputrec->opts),
706 x,hist,f,enerd,fcd,mtop,top,fr->born,
707 &(top->atomtypes),bBornRadii,box,
708 lambda,graph,&(top->excls),fr->mu_tot,
709 flags,&cycles_pme);
711 cycles_force = wallcycle_stop(wcycle,ewcFORCE);
712 GMX_BARRIER(cr->mpi_comm_mygroup);
714 if (ed)
716 do_flood(fplog,cr,x,f,ed,box,step);
719 if (DOMAINDECOMP(cr))
721 dd_force_flop_stop(cr->dd,nrnb);
722 if (wcycle)
724 dd_cycles_add(cr->dd,cycles_force-cycles_pme,ddCyclF);
728 if (bDoForces)
730 if (IR_ELEC_FIELD(*inputrec))
732 /* Compute forces due to electric field */
733 calc_f_el(MASTER(cr) ? field : NULL,
734 start,homenr,mdatoms->chargeA,x,fr->f_novirsum,
735 inputrec->ex,inputrec->et,t);
738 /* Communicate the forces */
739 if (PAR(cr))
741 wallcycle_start(wcycle,ewcMOVEF);
742 if (DOMAINDECOMP(cr))
744 dd_move_f(cr->dd,f,fr->fshift);
745 /* Do we need to communicate the separate force array
746 * for terms that do not contribute to the single sum virial?
747 * Position restraints and electric fields do not introduce
748 * inter-cg forces, only full electrostatics methods do.
749 * When we do not calculate the virial, fr->f_novirsum = f,
750 * so we have already communicated these forces.
752 if (EEL_FULL(fr->eeltype) && cr->dd->n_intercg_excl &&
753 (flags & GMX_FORCE_VIRIAL))
755 dd_move_f(cr->dd,fr->f_novirsum,NULL);
757 if (bSepLRF)
759 /* We should not update the shift forces here,
760 * since f_twin is already included in f.
762 dd_move_f(cr->dd,fr->f_twin,NULL);
765 else
767 pd_move_f(cr,f,nrnb);
768 if (bSepLRF)
770 pd_move_f(cr,fr->f_twin,nrnb);
773 wallcycle_stop(wcycle,ewcMOVEF);
776 /* If we have NoVirSum forces, but we do not calculate the virial,
777 * we sum fr->f_novirum=f later.
779 if (vsite && !(fr->bF_NoVirSum && !(flags & GMX_FORCE_VIRIAL)))
781 wallcycle_start(wcycle,ewcVSITESPREAD);
782 spread_vsite_f(fplog,vsite,x,f,fr->fshift,nrnb,
783 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
784 wallcycle_stop(wcycle,ewcVSITESPREAD);
786 if (bSepLRF)
788 wallcycle_start(wcycle,ewcVSITESPREAD);
789 spread_vsite_f(fplog,vsite,x,fr->f_twin,NULL,
790 nrnb,
791 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
792 wallcycle_stop(wcycle,ewcVSITESPREAD);
796 if (flags & GMX_FORCE_VIRIAL)
798 /* Calculation of the virial must be done after vsites! */
799 calc_virial(fplog,mdatoms->start,mdatoms->homenr,x,f,
800 vir_force,graph,box,nrnb,fr,inputrec->ePBC);
804 if (inputrec->ePull == epullUMBRELLA || inputrec->ePull == epullCONST_F)
806 /* Calculate the center of mass forces, this requires communication,
807 * which is why pull_potential is called close to other communication.
808 * The virial contribution is calculated directly,
809 * which is why we call pull_potential after calc_virial.
811 set_pbc(&pbc,inputrec->ePBC,box);
812 dvdl = 0;
813 enerd->term[F_COM_PULL] =
814 pull_potential(inputrec->ePull,inputrec->pull,mdatoms,&pbc,
815 cr,t,lambda,x,f,vir_force,&dvdl);
816 if (bSepDVDL)
818 fprintf(fplog,sepdvdlformat,"Com pull",enerd->term[F_COM_PULL],dvdl);
820 enerd->dvdl_lin += dvdl;
822 else
823 enerd->term[F_COM_PULL] = 0.0;
825 /* Add the forces from enforced rotation potentials (if any) */
826 if (inputrec->bRot)
827 enerd->term[F_COM_PULL] += add_rot_forces(inputrec->rot, f, cr,step,t);
830 if (PAR(cr) && !(cr->duty & DUTY_PME))
832 cycles_ppdpme = wallcycle_stop(wcycle,ewcPPDURINGPME);
833 dd_cycles_add(cr->dd,cycles_ppdpme,ddCyclPPduringPME);
835 /* In case of node-splitting, the PP nodes receive the long-range
836 * forces, virial and energy from the PME nodes here.
838 wallcycle_start(wcycle,ewcPP_PMEWAITRECVF);
839 dvdl = 0;
840 gmx_pme_receive_f(cr,fr->f_novirsum,fr->vir_el_recip,&e,&dvdl,
841 &cycles_seppme);
842 if (bSepDVDL)
844 fprintf(fplog,sepdvdlformat,"PME mesh",e,dvdl);
846 enerd->term[F_COUL_RECIP] += e;
847 enerd->dvdl_lin += dvdl;
848 if (wcycle)
850 dd_cycles_add(cr->dd,cycles_seppme,ddCyclPME);
852 wallcycle_stop(wcycle,ewcPP_PMEWAITRECVF);
855 if (bDoForces && fr->bF_NoVirSum)
857 if (vsite)
859 /* Spread the mesh force on virtual sites to the other particles...
860 * This is parallellized. MPI communication is performed
861 * if the constructing atoms aren't local.
863 wallcycle_start(wcycle,ewcVSITESPREAD);
864 spread_vsite_f(fplog,vsite,x,fr->f_novirsum,NULL,nrnb,
865 &top->idef,fr->ePBC,fr->bMolPBC,graph,box,cr);
866 wallcycle_stop(wcycle,ewcVSITESPREAD);
868 if (flags & GMX_FORCE_VIRIAL)
870 /* Now add the forces, this is local */
871 if (fr->bDomDec)
873 sum_forces(0,fr->f_novirsum_n,f,fr->f_novirsum);
875 else
877 sum_forces(start,start+homenr,f,fr->f_novirsum);
879 if (EEL_FULL(fr->eeltype))
881 /* Add the mesh contribution to the virial */
882 m_add(vir_force,fr->vir_el_recip,vir_force);
884 if (debug)
886 pr_rvecs(debug,0,"vir_force",vir_force,DIM);
891 /* Sum the potential energy terms from group contributions */
892 sum_epot(&(inputrec->opts),enerd);
894 if (fr->print_force >= 0 && bDoForces)
896 print_large_forces(stderr,mdatoms,cr,step,fr->print_force,x,f);
900 void do_constrain_first(FILE *fplog,gmx_constr_t constr,
901 t_inputrec *inputrec,t_mdatoms *md,
902 t_state *state,
903 t_graph *graph,t_commrec *cr,t_nrnb *nrnb,
904 t_forcerec *fr,t_idef *idef)
906 int i,m,start,end;
907 gmx_large_int_t step;
908 double mass,tmass,vcm[4];
909 real dt=inputrec->delta_t;
910 real dvdlambda;
911 rvec *xcon;
912 char buf[22];
914 start = md->start;
915 end = md->homenr + start;
916 if (debug)
918 fprintf(debug,"vcm: start=%d, homenr=%d, end=%d\n",
919 start,md->homenr,end);
921 snew(xcon,state->nalloc);
923 /* Do a first constraining to reset particles... */
924 step = inputrec->init_step;
925 if (fplog)
927 fprintf(fplog,"\nConstraining the starting coordinates (step %s)\n",
928 gmx_step_str(step,buf));
930 dvdlambda = 0;
931 constrain(NULL,TRUE,FALSE,constr,idef,
932 inputrec,cr,step,0,md,
933 state->x,state->x,NULL,
934 state->box,state->lambda,&dvdlambda,
935 NULL,NULL,nrnb,econqCoord);
937 if (EI_STATE_VELOCITY(inputrec->eI)) {
938 for(i=start; (i<end); i++) {
939 for(m=0; (m<DIM); m++) {
940 /* Reverse the velocity */
941 state->v[i][m] = -state->v[i][m];
942 /* Store the position at t-dt in xcon */
943 xcon[i][m] = state->x[i][m] + dt*state->v[i][m];
947 /* Constrain the positions at t=-dt with the positions at t=0
948 * as reference coordinates.
950 if (fplog)
952 fprintf(fplog,"\nConstraining the coordinates at t0-dt (step %s)\n",
953 gmx_step_str(step,buf));
955 dvdlambda = 0;
956 constrain(NULL,TRUE,FALSE,constr,idef,
957 inputrec,cr,step,-1,md,
958 state->x,xcon,NULL,
959 state->box,state->lambda,&dvdlambda,
960 state->v,NULL,nrnb,econqCoord);
962 sfree(xcon);
964 for(m=0; (m<4); m++)
966 vcm[m] = 0;
968 for(i=start; i<end; i++)
970 mass = md->massT[i];
971 for(m=0; m<DIM; m++)
973 /* Re-reverse the velocities */
974 state->v[i][m] = -state->v[i][m];
975 vcm[m] += state->v[i][m]*mass;
977 vcm[3] += mass;
980 if (inputrec->nstcomm != 0 || debug)
982 /* Compute the global sum of vcm */
983 if (debug)
985 fprintf(debug,
986 "vcm: %8.3f %8.3f %8.3f,"
987 " total mass = %12.5e\n",
988 vcm[XX],vcm[YY],vcm[ZZ],vcm[3]);
990 if (PAR(cr))
992 gmx_sumd(4,vcm,cr);
994 tmass = vcm[3];
995 for(m=0; (m<DIM); m++)
997 vcm[m] /= tmass;
999 if (debug)
1001 fprintf(debug,"vcm: %8.3f %8.3f %8.3f,"
1002 " total mass = %12.5e\n",vcm[XX],vcm[YY],vcm[ZZ],tmass);
1004 if (inputrec->nstcomm != 0)
1006 /* Now we have the velocity of center of mass,
1007 * let's remove it.
1009 for(i=start; (i<end); i++)
1011 for(m=0; (m<DIM); m++)
1013 state->v[i][m] -= vcm[m];
1021 void calc_enervirdiff(FILE *fplog,int eDispCorr,t_forcerec *fr)
1023 double eners[2],virs[2],enersum,virsum,y0,f,g,h;
1024 double r0,r1,r,rc3,rc9,ea,eb,ec,pa,pb,pc,pd;
1025 double invscale,invscale2,invscale3;
1026 int ri0,ri1,ri,i,offstart,offset;
1027 real scale,*vdwtab;
1029 fr->enershiftsix = 0;
1030 fr->enershifttwelve = 0;
1031 fr->enerdiffsix = 0;
1032 fr->enerdifftwelve = 0;
1033 fr->virdiffsix = 0;
1034 fr->virdifftwelve = 0;
1036 if (eDispCorr != edispcNO) {
1037 for(i=0; i<2; i++) {
1038 eners[i] = 0;
1039 virs[i] = 0;
1041 if ((fr->vdwtype == evdwSWITCH) || (fr->vdwtype == evdwSHIFT)) {
1042 if (fr->rvdw_switch == 0)
1043 gmx_fatal(FARGS,
1044 "With dispersion correction rvdw-switch can not be zero "
1045 "for vdw-type = %s",evdw_names[fr->vdwtype]);
1047 scale = fr->nblists[0].tab.scale;
1048 vdwtab = fr->nblists[0].vdwtab;
1050 /* Round the cut-offs to exact table values for precision */
1051 ri0 = floor(fr->rvdw_switch*scale);
1052 ri1 = ceil(fr->rvdw*scale);
1053 r0 = ri0/scale;
1054 r1 = ri1/scale;
1055 rc3 = r0*r0*r0;
1056 rc9 = rc3*rc3*rc3;
1058 if (fr->vdwtype == evdwSHIFT) {
1059 /* Determine the constant energy shift below rvdw_switch */
1060 fr->enershiftsix = (real)(-1.0/(rc3*rc3)) - vdwtab[8*ri0];
1061 fr->enershifttwelve = (real)( 1.0/(rc9*rc3)) - vdwtab[8*ri0 + 4];
1063 /* Add the constant part from 0 to rvdw_switch.
1064 * This integration from 0 to rvdw_switch overcounts the number
1065 * of interactions by 1, as it also counts the self interaction.
1066 * We will correct for this later.
1068 eners[0] += 4.0*M_PI*fr->enershiftsix*rc3/3.0;
1069 eners[1] += 4.0*M_PI*fr->enershifttwelve*rc3/3.0;
1071 invscale = 1.0/(scale);
1072 invscale2 = invscale*invscale;
1073 invscale3 = invscale*invscale2;
1075 /* following summation derived from cubic spline definition,
1076 Numerical Recipies in C, second edition, p. 113-116. Exact
1077 for the cubic spline. We first calculate the negative of
1078 the energy from rvdw to rvdw_switch, assuming that g(r)=1,
1079 and then add the more standard, abrupt cutoff correction to
1080 that result, yielding the long-range correction for a
1081 switched function. We perform both the pressure and energy
1082 loops at the same time for simplicity, as the computational
1083 cost is low. */
1085 for (i=0;i<2;i++) {
1086 enersum = 0.0; virsum = 0.0;
1087 if (i==0)
1088 offstart = 0;
1089 else
1090 offstart = 4;
1091 for (ri=ri0; ri<ri1; ri++) {
1092 r = ri*invscale;
1093 ea = invscale3;
1094 eb = 2.0*invscale2*r;
1095 ec = invscale*r*r;
1097 pa = invscale3;
1098 pb = 3.0*invscale2*r;
1099 pc = 3.0*invscale*r*r;
1100 pd = r*r*r;
1102 /* this "8" is from the packing in the vdwtab array - perhaps
1103 should be #define'ed? */
1104 offset = 8*ri + offstart;
1105 y0 = vdwtab[offset];
1106 f = vdwtab[offset+1];
1107 g = vdwtab[offset+2];
1108 h = vdwtab[offset+3];
1110 enersum += y0*(ea/3 + eb/2 + ec) + f*(ea/4 + eb/3 + ec/2)+
1111 g*(ea/5 + eb/4 + ec/3) + h*(ea/6 + eb/5 + ec/4);
1112 virsum += f*(pa/4 + pb/3 + pc/2 + pd) +
1113 2*g*(pa/5 + pb/4 + pc/3 + pd/2) + 3*h*(pa/6 + pb/5 + pc/4 + pd/3);
1116 enersum *= 4.0*M_PI;
1117 virsum *= 4.0*M_PI;
1118 eners[i] -= enersum;
1119 virs[i] -= virsum;
1122 /* now add the correction for rvdw_switch to infinity */
1123 eners[0] += -4.0*M_PI/(3.0*rc3);
1124 eners[1] += 4.0*M_PI/(9.0*rc9);
1125 virs[0] += 8.0*M_PI/rc3;
1126 virs[1] += -16.0*M_PI/(3.0*rc9);
1128 else if ((fr->vdwtype == evdwCUT) || (fr->vdwtype == evdwUSER)) {
1129 if (fr->vdwtype == evdwUSER && fplog)
1130 fprintf(fplog,
1131 "WARNING: using dispersion correction with user tables\n");
1132 rc3 = fr->rvdw*fr->rvdw*fr->rvdw;
1133 rc9 = rc3*rc3*rc3;
1134 eners[0] += -4.0*M_PI/(3.0*rc3);
1135 eners[1] += 4.0*M_PI/(9.0*rc9);
1136 virs[0] += 8.0*M_PI/rc3;
1137 virs[1] += -16.0*M_PI/(3.0*rc9);
1138 } else {
1139 gmx_fatal(FARGS,
1140 "Dispersion correction is not implemented for vdw-type = %s",
1141 evdw_names[fr->vdwtype]);
1143 fr->enerdiffsix = eners[0];
1144 fr->enerdifftwelve = eners[1];
1145 /* The 0.5 is due to the Gromacs definition of the virial */
1146 fr->virdiffsix = 0.5*virs[0];
1147 fr->virdifftwelve = 0.5*virs[1];
1151 void calc_dispcorr(FILE *fplog,t_inputrec *ir,t_forcerec *fr,
1152 gmx_large_int_t step,int natoms,matrix box,real lambda,
1153 tensor pres,tensor virial,gmx_enerdata_t *enerd)
1155 bool bCorrAll,bCorrPres;
1156 real dvdlambda,invvol,dens,ninter,avcsix,avctwelve,enerdiff,svir=0,spres=0;
1157 int m;
1159 enerd->term[F_DISPCORR] = 0.0;
1160 enerd->term[F_PDISPCORR] = 0.0;
1162 if (ir->eDispCorr != edispcNO)
1164 bCorrAll = (ir->eDispCorr == edispcAllEner ||
1165 ir->eDispCorr == edispcAllEnerPres);
1166 bCorrPres = (ir->eDispCorr == edispcEnerPres ||
1167 ir->eDispCorr == edispcAllEnerPres);
1169 invvol = 1/det(box);
1170 if (fr->n_tpi)
1172 /* Only correct for the interactions with the inserted molecule */
1173 dens = (natoms - fr->n_tpi)*invvol;
1174 ninter = fr->n_tpi;
1176 else
1178 dens = natoms*invvol;
1179 ninter = 0.5*natoms;
1182 if (ir->efep == efepNO)
1184 avcsix = fr->avcsix[0];
1185 avctwelve = fr->avctwelve[0];
1187 else
1189 avcsix = (1 - lambda)*fr->avcsix[0] + lambda*fr->avcsix[1];
1190 avctwelve = (1 - lambda)*fr->avctwelve[0] + lambda*fr->avctwelve[1];
1193 enerdiff = ninter*(dens*fr->enerdiffsix - fr->enershiftsix);
1194 enerd->term[F_DISPCORR] += avcsix*enerdiff;
1195 dvdlambda = 0.0;
1196 if (ir->efep != efepNO)
1198 dvdlambda += (fr->avcsix[1] - fr->avcsix[0])*enerdiff;
1201 if (bCorrAll)
1203 enerdiff = ninter*(dens*fr->enerdifftwelve - fr->enershifttwelve);
1204 enerd->term[F_DISPCORR] += avctwelve*enerdiff;
1205 if (fr->efep != efepNO)
1207 dvdlambda += (fr->avctwelve[1] - fr->avctwelve[0])*enerdiff;
1211 if (bCorrPres)
1213 svir = ninter*dens*avcsix*fr->virdiffsix/3.0;
1214 if (ir->eDispCorr == edispcAllEnerPres)
1216 svir += ninter*dens*avctwelve*fr->virdifftwelve/3.0;
1218 /* The factor 2 is because of the Gromacs virial definition */
1219 spres = -2.0*invvol*svir*PRESFAC;
1221 for(m=0; m<DIM; m++) {
1222 virial[m][m] += svir;
1223 pres[m][m] += spres;
1225 enerd->term[F_PDISPCORR] = spres;
1226 enerd->term[F_PRES] += spres;
1229 if (fr->bSepDVDL && do_per_step(step,ir->nstlog))
1231 fprintf(fplog,sepdvdlformat,"Dispersion correction",
1232 enerd->term[F_DISPCORR],dvdlambda);
1235 enerd->term[F_EPOT] += enerd->term[F_DISPCORR];
1236 if (fr->efep != efepNO)
1238 enerd->dvdl_lin += dvdlambda;
1244 void do_pbc_first(FILE *fplog,matrix box,t_forcerec *fr,
1245 t_graph *graph,rvec x[])
1247 if (fplog)
1248 fprintf(fplog,"Removing pbc first time\n");
1249 calc_shifts(box,fr->shift_vec);
1250 if (graph) {
1251 mk_mshift(fplog,graph,fr->ePBC,box,x);
1252 if (gmx_debug_at)
1253 p_graph(debug,"do_pbc_first 1",graph);
1254 shift_self(graph,box,x);
1255 /* By doing an extra mk_mshift the molecules that are broken
1256 * because they were e.g. imported from another software
1257 * will be made whole again. Such are the healing powers
1258 * of GROMACS.
1260 mk_mshift(fplog,graph,fr->ePBC,box,x);
1261 if (gmx_debug_at)
1262 p_graph(debug,"do_pbc_first 2",graph);
1264 if (fplog)
1265 fprintf(fplog,"Done rmpbc\n");
1268 static void low_do_pbc_mtop(FILE *fplog,int ePBC,matrix box,
1269 gmx_mtop_t *mtop,rvec x[],
1270 bool bFirst)
1272 t_graph *graph;
1273 int mb,as,mol;
1274 gmx_molblock_t *molb;
1276 if (bFirst && fplog)
1277 fprintf(fplog,"Removing pbc first time\n");
1279 snew(graph,1);
1280 as = 0;
1281 for(mb=0; mb<mtop->nmolblock; mb++) {
1282 molb = &mtop->molblock[mb];
1283 if (molb->natoms_mol == 1 ||
1284 (!bFirst && mtop->moltype[molb->type].cgs.nr == 1)) {
1285 /* Just one atom or charge group in the molecule, no PBC required */
1286 as += molb->nmol*molb->natoms_mol;
1287 } else {
1288 /* Pass NULL iso fplog to avoid graph prints for each molecule type */
1289 mk_graph_ilist(NULL,mtop->moltype[molb->type].ilist,
1290 0,molb->natoms_mol,FALSE,FALSE,graph);
1292 for(mol=0; mol<molb->nmol; mol++) {
1293 mk_mshift(fplog,graph,ePBC,box,x+as);
1295 shift_self(graph,box,x+as);
1296 /* The molecule is whole now.
1297 * We don't need the second mk_mshift call as in do_pbc_first,
1298 * since we no longer need this graph.
1301 as += molb->natoms_mol;
1303 done_graph(graph);
1306 sfree(graph);
1309 void do_pbc_first_mtop(FILE *fplog,int ePBC,matrix box,
1310 gmx_mtop_t *mtop,rvec x[])
1312 low_do_pbc_mtop(fplog,ePBC,box,mtop,x,TRUE);
1315 void do_pbc_mtop(FILE *fplog,int ePBC,matrix box,
1316 gmx_mtop_t *mtop,rvec x[])
1318 low_do_pbc_mtop(fplog,ePBC,box,mtop,x,FALSE);
1321 void finish_run(FILE *fplog,t_commrec *cr,const char *confout,
1322 t_inputrec *inputrec,
1323 t_nrnb nrnb[],gmx_wallcycle_t wcycle,
1324 gmx_runtime_t *runtime,
1325 bool bWriteStat)
1327 int i,j;
1328 t_nrnb *nrnb_tot=NULL;
1329 real delta_t;
1330 double nbfs,mflop;
1331 double cycles[ewcNR];
1333 wallcycle_sum(cr,wcycle,cycles);
1335 if (cr->nnodes > 1) {
1336 if (SIMMASTER(cr))
1337 snew(nrnb_tot,1);
1338 #ifdef GMX_MPI
1339 MPI_Reduce(nrnb->n,nrnb_tot->n,eNRNB,MPI_DOUBLE,MPI_SUM,
1340 MASTERRANK(cr),cr->mpi_comm_mysim);
1341 #endif
1342 } else {
1343 nrnb_tot = nrnb;
1346 if (SIMMASTER(cr)) {
1347 print_flop(fplog,nrnb_tot,&nbfs,&mflop);
1348 if (cr->nnodes > 1) {
1349 sfree(nrnb_tot);
1353 if ((cr->duty & DUTY_PP) && DOMAINDECOMP(cr)) {
1354 print_dd_statistics(cr,inputrec,fplog);
1357 if (SIMMASTER(cr)) {
1358 if (PARTDECOMP(cr)) {
1359 pr_load(fplog,cr,nrnb_tot);
1362 wallcycle_print(fplog,cr->nnodes,cr->npmenodes,runtime->realtime,
1363 wcycle,cycles);
1365 if (EI_DYNAMICS(inputrec->eI)) {
1366 delta_t = inputrec->delta_t;
1367 } else {
1368 delta_t = 0;
1371 if (fplog) {
1372 print_perf(fplog,runtime->proctime,runtime->realtime,
1373 cr->nnodes-cr->npmenodes,
1374 runtime->nsteps_done,delta_t,nbfs,mflop);
1376 if (bWriteStat) {
1377 print_perf(stderr,runtime->proctime,runtime->realtime,
1378 cr->nnodes-cr->npmenodes,
1379 runtime->nsteps_done,delta_t,nbfs,mflop);
1383 runtime=inputrec->nsteps*inputrec->delta_t;
1384 if (bWriteStat) {
1385 if (cr->nnodes == 1)
1386 fprintf(stderr,"\n\n");
1387 print_perf(stderr,nodetime,realtime,runtime,&ntot,
1388 cr->nnodes-cr->npmenodes,FALSE);
1390 wallcycle_print(fplog,cr->nnodes,cr->npmenodes,realtime,wcycle,cycles);
1391 print_perf(fplog,nodetime,realtime,runtime,&ntot,cr->nnodes-cr->npmenodes,
1392 TRUE);
1393 if (PARTDECOMP(cr))
1394 pr_load(fplog,cr,nrnb_all);
1395 if (cr->nnodes > 1)
1396 sfree(nrnb_all);
1401 void init_md(FILE *fplog,
1402 t_commrec *cr,t_inputrec *ir,const output_env_t oenv,
1403 double *t,double *t0,
1404 real *lambda,double *lam0,
1405 t_nrnb *nrnb,gmx_mtop_t *mtop,
1406 gmx_update_t *upd,
1407 int nfile,const t_filenm fnm[],
1408 int *fp_trn,int *fp_xtc,ener_file_t *fp_ene,const char **fn_cpt,
1409 FILE **fp_dhdl,FILE **fp_field,
1410 t_mdebin **mdebin,
1411 tensor force_vir,tensor shake_vir,rvec mu_tot,
1412 bool *bNEMD,bool *bSimAnn,t_vcm **vcm, unsigned long Flags)
1414 int i,j,n;
1415 real tmpt,mod;
1416 char filemode[2];
1418 sprintf(filemode, (Flags & MD_APPENDFILES) ? "a" : "w");
1420 /* Initial values */
1421 *t = *t0 = ir->init_t;
1422 if (ir->efep != efepNO)
1424 *lam0 = ir->init_lambda;
1425 *lambda = *lam0 + ir->init_step*ir->delta_lambda;
1427 else
1429 *lambda = *lam0 = 0.0;
1432 *bSimAnn=FALSE;
1433 for(i=0;i<ir->opts.ngtc;i++)
1435 /* set bSimAnn if any group is being annealed */
1436 if(ir->opts.annealing[i]!=eannNO)
1438 *bSimAnn = TRUE;
1441 if (*bSimAnn)
1443 update_annealing_target_temp(&(ir->opts),ir->init_t);
1446 *bNEMD = (ir->opts.ngacc > 1) || (norm(ir->opts.acc[0]) > 0);
1448 if (upd)
1450 *upd = init_update(fplog,ir);
1453 if (vcm != NULL)
1455 *vcm = init_vcm(fplog,&mtop->groups,ir);
1458 if (EI_DYNAMICS(ir->eI) && !(Flags & MD_APPENDFILES))
1460 if (ir->etc == etcBERENDSEN)
1462 please_cite(fplog,"Berendsen84a");
1464 if (ir->etc == etcVRESCALE)
1466 please_cite(fplog,"Bussi2007a");
1470 init_nrnb(nrnb);
1472 if (nfile != -1)
1474 *fp_trn = -1;
1475 *fp_ene = NULL;
1476 *fp_xtc = -1;
1478 if (MASTER(cr))
1480 *fp_trn = open_trn(ftp2fn(efTRN,nfile,fnm), filemode);
1481 if (ir->nstxtcout > 0)
1483 *fp_xtc = open_xtc(ftp2fn(efXTC,nfile,fnm), filemode);
1485 *fp_ene = open_enx(ftp2fn(efEDR,nfile,fnm), filemode);
1486 *fn_cpt = opt2fn("-cpo",nfile,fnm);
1488 if ((fp_dhdl != NULL) && ir->efep != efepNO && ir->nstdhdl > 0)
1490 if(Flags & MD_APPENDFILES)
1492 *fp_dhdl= gmx_fio_fopen(opt2fn("-dhdl",nfile,fnm),filemode);
1494 else
1496 *fp_dhdl = open_dhdl(opt2fn("-dhdl",nfile,fnm),ir,oenv);
1500 if ((fp_field != NULL) &&
1501 (ir->ex[XX].n || ir->ex[YY].n ||ir->ex[ZZ].n))
1503 if(Flags & MD_APPENDFILES)
1505 *fp_dhdl=gmx_fio_fopen(opt2fn("-field",nfile,fnm),filemode);
1507 else
1509 *fp_field = xvgropen(opt2fn("-field",nfile,fnm),
1510 "Applied electric field","Time (ps)",
1511 "E (V/nm)",oenv);
1515 *mdebin = init_mdebin( (Flags & MD_APPENDFILES) ? NULL : *fp_ene,
1516 mtop,ir);
1519 /* Initiate variables */
1520 clear_mat(force_vir);
1521 clear_mat(shake_vir);
1522 clear_rvec(mu_tot);
1524 debug_gmx();