1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * GROwing Monsters And Cloning Shrimps
41 #include<catamount/dclock.h>
47 #ifdef HAVE_SYS_TIME_H
60 #include "chargegroup.h"
83 #include "pull_rotation.h"
84 #include "mpelogging.h"
87 #include "gmx_wallcycle.h"
100 typedef struct gmx_timeprint
{
109 #ifdef HAVE_GETTIMEOFDAY
111 struct timezone tz
= { 0,0 };
114 gettimeofday(&t
,&tz
);
116 seconds
= (double) t
.tv_sec
+ 1e-6*(double)t
.tv_usec
;
122 seconds
= time(NULL
);
129 #define difftime(end,start) ((double)(end)-(double)(start))
131 void print_time(FILE *out
,gmx_runtime_t
*runtime
,gmx_large_int_t step
,
132 t_inputrec
*ir
, t_commrec
*cr
)
143 fprintf(out
,"step %s",gmx_step_str(step
,buf
));
144 if ((step
>= ir
->nstlist
)) {
145 if ((ir
->nstlist
== 0) || ((step
% ir
->nstlist
) == 0)) {
146 /* We have done a full cycle let's update time_per_step */
147 runtime
->last
= gmx_gettime();
148 dt
= difftime(runtime
->last
,runtime
->real
);
149 runtime
->time_per_step
= dt
/(step
- ir
->init_step
+ 1);
151 dt
= (ir
->nsteps
+ ir
->init_step
- step
)*runtime
->time_per_step
;
154 finish
= (time_t) (runtime
->last
+ dt
);
155 sprintf(buf
,"%s",ctime(&finish
));
156 buf
[strlen(buf
)-1]='\0';
157 fprintf(out
,", will finish %s",buf
);
160 fprintf(out
,", remaining runtime: %5d s ",(int)dt
);
174 static double set_proctime(gmx_runtime_t
*runtime
)
180 prev
= runtime
->proc
;
181 runtime
->proc
= dclock();
183 diff
= runtime
->proc
- prev
;
187 prev
= runtime
->proc
;
188 runtime
->proc
= clock();
190 diff
= (double)(runtime
->proc
- prev
)/(double)CLOCKS_PER_SEC
;
194 /* The counter has probably looped, ignore this data */
201 void runtime_start(gmx_runtime_t
*runtime
)
203 runtime
->real
= gmx_gettime();
205 set_proctime(runtime
);
206 runtime
->realtime
= 0;
207 runtime
->proctime
= 0;
209 runtime
->time_per_step
= 0;
212 void runtime_end(gmx_runtime_t
*runtime
)
218 runtime
->proctime
+= set_proctime(runtime
);
219 runtime
->realtime
= now
- runtime
->real
;
223 void runtime_upd_proc(gmx_runtime_t
*runtime
)
225 runtime
->proctime
+= set_proctime(runtime
);
228 void print_date_and_time(FILE *fplog
,int nodeid
,const char *title
,
229 const gmx_runtime_t
*runtime
)
232 char *ts
,time_string
[STRLEN
];
237 tmptime
= (time_t) runtime
->real
;
238 ts
= ctime(&tmptime
);
242 tmptime
= (time_t) gmx_gettime();
243 ts
= ctime(&tmptime
);
245 for(i
=0; ts
[i
]>=' '; i
++)
247 time_string
[i
]=ts
[i
];
252 fprintf(fplog
,"%s on node %d %s\n",title
,nodeid
,time_string
);
256 static void sum_forces(int start
,int end
,rvec f
[],rvec flr
[])
261 pr_rvecs(debug
,0,"fsr",f
+start
,end
-start
);
262 pr_rvecs(debug
,0,"flr",flr
+start
,end
-start
);
264 for(i
=start
; (i
<end
); i
++)
265 rvec_inc(f
[i
],flr
[i
]);
269 * calc_f_el calculates forces due to an electric field.
271 * force is kJ mol^-1 nm^-1 = e * kJ mol^-1 nm^-1 / e
273 * Et[] contains the parameters for the time dependent
274 * part of the field (not yet used).
275 * Ex[] contains the parameters for
276 * the spatial dependent part of the field. You can have cool periodic
277 * fields in principle, but only a constant field is supported
279 * The function should return the energy due to the electric field
280 * (if any) but for now returns 0.
283 * There can be problems with the virial.
284 * Since the field is not self-consistent this is unavoidable.
285 * For neutral molecules the virial is correct within this approximation.
286 * For neutral systems with many charged molecules the error is small.
287 * But for systems with a net charge or a few charged molecules
288 * the error can be significant when the field is high.
289 * Solution: implement a self-consitent electric field into PME.
291 static void calc_f_el(FILE *fp
,int start
,int homenr
,
292 real charge
[],rvec x
[],rvec f
[],
293 t_cosines Ex
[],t_cosines Et
[],double t
)
299 for(m
=0; (m
<DIM
); m
++)
306 Ext
[m
] = cos(Et
[m
].a
[0]*(t
-t0
))*exp(-sqr(t
-t0
)/(2.0*sqr(Et
[m
].a
[2])));
310 Ext
[m
] = cos(Et
[m
].a
[0]*t
);
319 /* Convert the field strength from V/nm to MD-units */
320 Ext
[m
] *= Ex
[m
].a
[0]*FIELDFAC
;
321 for(i
=start
; (i
<start
+homenr
); i
++)
322 f
[i
][m
] += charge
[i
]*Ext
[m
];
331 fprintf(fp
,"%10g %10g %10g %10g #FIELD\n",t
,
332 Ext
[XX
]/FIELDFAC
,Ext
[YY
]/FIELDFAC
,Ext
[ZZ
]/FIELDFAC
);
336 static void calc_virial(FILE *fplog
,int start
,int homenr
,rvec x
[],rvec f
[],
337 tensor vir_part
,t_graph
*graph
,matrix box
,
338 t_nrnb
*nrnb
,const t_forcerec
*fr
,int ePBC
)
343 /* The short-range virial from surrounding boxes */
345 calc_vir(fplog
,SHIFTS
,fr
->shift_vec
,fr
->fshift
,vir_part
,ePBC
==epbcSCREW
,box
);
346 inc_nrnb(nrnb
,eNR_VIRIAL
,SHIFTS
);
348 /* Calculate partial virial, for local atoms only, based on short range.
349 * Total virial is computed in global_stat, called from do_md
351 f_calc_vir(fplog
,start
,start
+homenr
,x
,f
,vir_part
,graph
,box
);
352 inc_nrnb(nrnb
,eNR_VIRIAL
,homenr
);
354 /* Add position restraint contribution */
355 for(i
=0; i
<DIM
; i
++) {
356 vir_part
[i
][i
] += fr
->vir_diag_posres
[i
];
359 /* Add wall contribution */
360 for(i
=0; i
<DIM
; i
++) {
361 vir_part
[i
][ZZ
] += fr
->vir_wall_z
[i
];
365 pr_rvecs(debug
,0,"vir_part",vir_part
,DIM
);
368 static void print_large_forces(FILE *fp
,t_mdatoms
*md
,t_commrec
*cr
,
369 gmx_large_int_t step
,real pforce
,rvec
*x
,rvec
*f
)
373 char buf
[STEPSTRSIZE
];
376 for(i
=md
->start
; i
<md
->start
+md
->homenr
; i
++) {
378 /* We also catch NAN, if the compiler does not optimize this away. */
379 if (fn2
>= pf2
|| fn2
!= fn2
) {
380 fprintf(fp
,"step %s atom %6d x %8.3f %8.3f %8.3f force %12.5e\n",
381 gmx_step_str(step
,buf
),
382 ddglatnr(cr
->dd
,i
),x
[i
][XX
],x
[i
][YY
],x
[i
][ZZ
],sqrt(fn2
));
387 void do_force(FILE *fplog
,t_commrec
*cr
,
388 t_inputrec
*inputrec
,
389 gmx_large_int_t step
,t_nrnb
*nrnb
,gmx_wallcycle_t wcycle
,
392 gmx_groups_t
*groups
,
393 matrix box
,rvec x
[],history_t
*hist
,
397 gmx_enerdata_t
*enerd
,t_fcdata
*fcd
,
398 real lambda
,t_graph
*graph
,
399 t_forcerec
*fr
,gmx_vsite_t
*vsite
,rvec mu_tot
,
400 double t
,FILE *field
,gmx_edsam_t ed
,
407 bool bSepDVDL
,bStateChanged
,bNS
,bFillGrid
,bCalcCGCM
,bBS
;
408 bool bDoLongRange
,bDoForces
,bSepLRF
;
412 float cycles_ppdpme
,cycles_pme
,cycles_seppme
,cycles_force
;
414 start
= mdatoms
->start
;
415 homenr
= mdatoms
->homenr
;
417 bSepDVDL
= (fr
->bSepDVDL
&& do_per_step(step
,inputrec
->nstlog
));
419 clear_mat(vir_force
);
423 pd_cg_range(cr
,&cg0
,&cg1
);
428 if (DOMAINDECOMP(cr
))
430 cg1
= cr
->dd
->ncg_tot
;
442 bStateChanged
= (flags
& GMX_FORCE_STATECHANGED
);
443 bNS
= (flags
& GMX_FORCE_NS
) && (fr
->bAllvsAll
==FALSE
);
444 bFillGrid
= (bNS
&& bStateChanged
);
445 bCalcCGCM
= (bFillGrid
&& !DOMAINDECOMP(cr
));
446 bDoLongRange
= (fr
->bTwinRange
&& bNS
&& (flags
& GMX_FORCE_DOLR
));
447 bDoForces
= (flags
& GMX_FORCE_FORCES
);
448 bSepLRF
= (bDoLongRange
&& bDoForces
&& (flags
& GMX_FORCE_SEPLRF
));
452 update_forcerec(fplog
,fr
,box
);
454 /* Calculate total (local) dipole moment in a temporary common array.
455 * This makes it possible to sum them over nodes faster.
457 calc_mu(start
,homenr
,
458 x
,mdatoms
->chargeA
,mdatoms
->chargeB
,mdatoms
->nChargePerturbed
,
462 if (fr
->ePBC
!= epbcNONE
) {
463 /* Compute shift vectors every step,
464 * because of pressure coupling or box deformation!
466 if ((flags
& GMX_FORCE_DYNAMICBOX
) && bStateChanged
)
467 calc_shifts(box
,fr
->shift_vec
);
470 put_charge_groups_in_box(fplog
,cg0
,cg1
,fr
->ePBC
,box
,
471 &(top
->cgs
),x
,fr
->cg_cm
);
472 inc_nrnb(nrnb
,eNR_CGCM
,homenr
);
473 inc_nrnb(nrnb
,eNR_RESETX
,cg1
-cg0
);
475 else if (EI_ENERGY_MINIMIZATION(inputrec
->eI
) && graph
) {
476 unshift_self(graph
,box
,x
);
479 else if (bCalcCGCM
) {
480 calc_cgcm(fplog
,cg0
,cg1
,&(top
->cgs
),x
,fr
->cg_cm
);
481 inc_nrnb(nrnb
,eNR_CGCM
,homenr
);
486 move_cgcm(fplog
,cr
,fr
->cg_cm
);
489 pr_rvecs(debug
,0,"cgcm",fr
->cg_cm
,top
->cgs
.nr
);
493 if (!(cr
->duty
& DUTY_PME
)) {
494 /* Send particle coordinates to the pme nodes.
495 * Since this is only implemented for domain decomposition
496 * and domain decomposition does not use the graph,
497 * we do not need to worry about shifting.
500 wallcycle_start(wcycle
,ewcPP_PMESENDX
);
501 GMX_MPE_LOG(ev_send_coordinates_start
);
503 bBS
= (inputrec
->nwall
== 2);
506 svmul(inputrec
->wall_ewald_zfac
,boxs
[ZZ
],boxs
[ZZ
]);
509 gmx_pme_send_x(cr
,bBS
? boxs
: box
,x
,mdatoms
->nChargePerturbed
,lambda
,step
);
511 GMX_MPE_LOG(ev_send_coordinates_finish
);
512 wallcycle_stop(wcycle
,ewcPP_PMESENDX
);
516 /* Communicate coordinates and sum dipole if necessary */
519 wallcycle_start(wcycle
,ewcMOVEX
);
520 if (DOMAINDECOMP(cr
))
522 dd_move_x(cr
->dd
,box
,x
);
526 move_x(fplog
,cr
,GMX_LEFT
,GMX_RIGHT
,x
,nrnb
);
528 /* When we don't need the total dipole we sum it in global_stat */
529 if (bStateChanged
&& NEED_MUTOT(*inputrec
))
531 gmx_sumd(2*DIM
,mu
,cr
);
533 wallcycle_stop(wcycle
,ewcMOVEX
);
541 fr
->mu_tot
[i
][j
] = mu
[i
*DIM
+ j
];
545 if (fr
->efep
== efepNO
)
547 copy_rvec(fr
->mu_tot
[0],mu_tot
);
554 (1.0 - lambda
)*fr
->mu_tot
[0][j
] + lambda
*fr
->mu_tot
[1][j
];
559 reset_enerdata(&(inputrec
->opts
),fr
,bNS
,enerd
,MASTER(cr
));
560 clear_rvecs(SHIFTS
,fr
->fshift
);
564 wallcycle_start(wcycle
,ewcNS
);
566 if (graph
&& bStateChanged
)
568 /* Calculate intramolecular shift vectors to make molecules whole */
569 mk_mshift(fplog
,graph
,fr
->ePBC
,box
,x
);
572 /* Reset long range forces if necessary */
575 /* Reset the (long-range) forces if necessary */
576 clear_rvecs(fr
->natoms_force
,bSepLRF
? fr
->f_twin
: f
);
579 /* Do the actual neighbour searching and if twin range electrostatics
580 * also do the calculation of long range forces and energies.
584 groups
,&(inputrec
->opts
),top
,mdatoms
,
585 cr
,nrnb
,lambda
,&dvdl
,&enerd
->grpp
,bFillGrid
,
586 bDoLongRange
,bDoForces
,bSepLRF
? fr
->f_twin
: f
);
589 fprintf(fplog
,sepdvdlformat
,"LR non-bonded",0.0,dvdl
);
591 enerd
->dvdl_lin
+= dvdl
;
593 wallcycle_stop(wcycle
,ewcNS
);
596 if (inputrec
->implicit_solvent
&& bNS
)
598 make_gb_nblist(cr
,inputrec
->gb_algorithm
,inputrec
->rlist
,
599 x
,box
,fr
,&top
->idef
,graph
,fr
->born
);
602 if (DOMAINDECOMP(cr
))
604 if (!(cr
->duty
& DUTY_PME
))
606 wallcycle_start(wcycle
,ewcPPDURINGPME
);
607 dd_force_flop_start(cr
->dd
,nrnb
);
613 /* Enforced rotation has its own cycle counter that starts after the collective
614 * coordinates have been communicated. It is added to ddCyclF */
615 do_rotation(cr
,inputrec
,box
,x
,t
,step
,wcycle
,bNS
);
618 /* Start the force cycle counter.
619 * This counter is stopped in do_forcelow_level.
620 * No parallel communication should occur while this counter is running,
621 * since that will interfere with the dynamic load balancing.
623 wallcycle_start(wcycle
,ewcFORCE
);
624 GMX_MPE_LOG(ev_forcecycles_start
);
628 /* Reset forces for which the virial is calculated separately:
629 * PME/Ewald forces if necessary */
632 if (flags
& GMX_FORCE_VIRIAL
)
634 fr
->f_novirsum
= fr
->f_novirsum_alloc
;
635 GMX_BARRIER(cr
->mpi_comm_mygroup
);
638 clear_rvecs(fr
->f_novirsum_n
,fr
->f_novirsum
);
642 clear_rvecs(homenr
,fr
->f_novirsum
+start
);
644 GMX_BARRIER(cr
->mpi_comm_mygroup
);
648 /* We are not calculating the pressure so we do not need
649 * a separate array for forces that do not contribute
658 /* Add the long range forces to the short range forces */
659 for(i
=0; i
<fr
->natoms_force
; i
++)
661 copy_rvec(fr
->f_twin
[i
],f
[i
]);
664 else if (!(fr
->bTwinRange
&& bNS
))
666 /* Clear the short-range forces */
667 clear_rvecs(fr
->natoms_force
,f
);
670 clear_rvec(fr
->vir_diag_posres
);
672 GMX_BARRIER(cr
->mpi_comm_mygroup
);
674 if (inputrec
->ePull
== epullCONSTRAINT
)
676 clear_pull_forces(inputrec
->pull
);
679 /* update QMMMrec, if necessary */
682 update_QMMMrec(cr
,fr
,x
,mdatoms
,box
,top
);
685 if ((flags
& GMX_FORCE_BONDED
) && top
->idef
.il
[F_POSRES
].nr
> 0)
687 /* Position restraints always require full pbc */
688 set_pbc(&pbc
,inputrec
->ePBC
,box
);
689 v
= posres(top
->idef
.il
[F_POSRES
].nr
,top
->idef
.il
[F_POSRES
].iatoms
,
690 top
->idef
.iparams_posres
,
691 (const rvec
*)x
,fr
->f_novirsum
,fr
->vir_diag_posres
,
692 inputrec
->ePBC
==epbcNONE
? NULL
: &pbc
,lambda
,&dvdl
,
693 fr
->rc_scaling
,fr
->ePBC
,fr
->posres_com
,fr
->posres_comB
);
696 fprintf(fplog
,sepdvdlformat
,
697 interaction_function
[F_POSRES
].longname
,v
,dvdl
);
699 enerd
->term
[F_POSRES
] += v
;
700 /* This linear lambda dependence assumption is only correct
701 * when only k depends on lambda,
702 * not when the reference position depends on lambda.
703 * grompp checks for this.
705 enerd
->dvdl_lin
+= dvdl
;
706 inc_nrnb(nrnb
,eNR_POSRES
,top
->idef
.il
[F_POSRES
].nr
/2);
709 /* Compute the bonded and non-bonded energies and optionally forces */
710 do_force_lowlevel(fplog
,step
,fr
,inputrec
,&(top
->idef
),
711 cr
,nrnb
,wcycle
,mdatoms
,&(inputrec
->opts
),
712 x
,hist
,f
,enerd
,fcd
,mtop
,top
,fr
->born
,
713 &(top
->atomtypes
),bBornRadii
,box
,
714 lambda
,graph
,&(top
->excls
),fr
->mu_tot
,
717 cycles_force
= wallcycle_stop(wcycle
,ewcFORCE
);
718 GMX_BARRIER(cr
->mpi_comm_mygroup
);
722 do_flood(fplog
,cr
,x
,f
,ed
,box
,step
);
725 if (DOMAINDECOMP(cr
))
727 dd_force_flop_stop(cr
->dd
,nrnb
);
730 dd_cycles_add(cr
->dd
,cycles_force
-cycles_pme
,ddCyclF
);
736 if (IR_ELEC_FIELD(*inputrec
))
738 /* Compute forces due to electric field */
739 calc_f_el(MASTER(cr
) ? field
: NULL
,
740 start
,homenr
,mdatoms
->chargeA
,x
,fr
->f_novirsum
,
741 inputrec
->ex
,inputrec
->et
,t
);
744 /* Communicate the forces */
747 wallcycle_start(wcycle
,ewcMOVEF
);
748 if (DOMAINDECOMP(cr
))
750 dd_move_f(cr
->dd
,f
,fr
->fshift
);
751 /* Do we need to communicate the separate force array
752 * for terms that do not contribute to the single sum virial?
753 * Position restraints and electric fields do not introduce
754 * inter-cg forces, only full electrostatics methods do.
755 * When we do not calculate the virial, fr->f_novirsum = f,
756 * so we have already communicated these forces.
758 if (EEL_FULL(fr
->eeltype
) && cr
->dd
->n_intercg_excl
&&
759 (flags
& GMX_FORCE_VIRIAL
))
761 dd_move_f(cr
->dd
,fr
->f_novirsum
,NULL
);
765 /* We should not update the shift forces here,
766 * since f_twin is already included in f.
768 dd_move_f(cr
->dd
,fr
->f_twin
,NULL
);
773 pd_move_f(cr
,f
,nrnb
);
776 pd_move_f(cr
,fr
->f_twin
,nrnb
);
779 wallcycle_stop(wcycle
,ewcMOVEF
);
782 /* If we have NoVirSum forces, but we do not calculate the virial,
783 * we sum fr->f_novirum=f later.
785 if (vsite
&& !(fr
->bF_NoVirSum
&& !(flags
& GMX_FORCE_VIRIAL
)))
787 wallcycle_start(wcycle
,ewcVSITESPREAD
);
788 spread_vsite_f(fplog
,vsite
,x
,f
,fr
->fshift
,nrnb
,
789 &top
->idef
,fr
->ePBC
,fr
->bMolPBC
,graph
,box
,cr
);
790 wallcycle_stop(wcycle
,ewcVSITESPREAD
);
794 wallcycle_start(wcycle
,ewcVSITESPREAD
);
795 spread_vsite_f(fplog
,vsite
,x
,fr
->f_twin
,NULL
,
797 &top
->idef
,fr
->ePBC
,fr
->bMolPBC
,graph
,box
,cr
);
798 wallcycle_stop(wcycle
,ewcVSITESPREAD
);
802 if (flags
& GMX_FORCE_VIRIAL
)
804 /* Calculation of the virial must be done after vsites! */
805 calc_virial(fplog
,mdatoms
->start
,mdatoms
->homenr
,x
,f
,
806 vir_force
,graph
,box
,nrnb
,fr
,inputrec
->ePBC
);
810 if (inputrec
->ePull
== epullUMBRELLA
|| inputrec
->ePull
== epullCONST_F
)
812 /* Calculate the center of mass forces, this requires communication,
813 * which is why pull_potential is called close to other communication.
814 * The virial contribution is calculated directly,
815 * which is why we call pull_potential after calc_virial.
817 set_pbc(&pbc
,inputrec
->ePBC
,box
);
819 enerd
->term
[F_COM_PULL
] =
820 pull_potential(inputrec
->ePull
,inputrec
->pull
,mdatoms
,&pbc
,
821 cr
,t
,lambda
,x
,f
,vir_force
,&dvdl
);
824 fprintf(fplog
,sepdvdlformat
,"Com pull",enerd
->term
[F_COM_PULL
],dvdl
);
826 enerd
->dvdl_lin
+= dvdl
;
829 enerd
->term
[F_COM_PULL
] = 0.0;
831 /* Add the forces from enforced rotation potentials (if any) */
833 enerd
->term
[F_COM_PULL
] += add_rot_forces(inputrec
->rot
, f
, cr
,step
,t
);
836 if (PAR(cr
) && !(cr
->duty
& DUTY_PME
))
838 cycles_ppdpme
= wallcycle_stop(wcycle
,ewcPPDURINGPME
);
839 dd_cycles_add(cr
->dd
,cycles_ppdpme
,ddCyclPPduringPME
);
841 /* In case of node-splitting, the PP nodes receive the long-range
842 * forces, virial and energy from the PME nodes here.
844 wallcycle_start(wcycle
,ewcPP_PMEWAITRECVF
);
846 gmx_pme_receive_f(cr
,fr
->f_novirsum
,fr
->vir_el_recip
,&e
,&dvdl
,
850 fprintf(fplog
,sepdvdlformat
,"PME mesh",e
,dvdl
);
852 enerd
->term
[F_COUL_RECIP
] += e
;
853 enerd
->dvdl_lin
+= dvdl
;
856 dd_cycles_add(cr
->dd
,cycles_seppme
,ddCyclPME
);
858 wallcycle_stop(wcycle
,ewcPP_PMEWAITRECVF
);
861 if (bDoForces
&& fr
->bF_NoVirSum
)
865 /* Spread the mesh force on virtual sites to the other particles...
866 * This is parallellized. MPI communication is performed
867 * if the constructing atoms aren't local.
869 wallcycle_start(wcycle
,ewcVSITESPREAD
);
870 spread_vsite_f(fplog
,vsite
,x
,fr
->f_novirsum
,NULL
,nrnb
,
871 &top
->idef
,fr
->ePBC
,fr
->bMolPBC
,graph
,box
,cr
);
872 wallcycle_stop(wcycle
,ewcVSITESPREAD
);
874 if (flags
& GMX_FORCE_VIRIAL
)
876 /* Now add the forces, this is local */
879 sum_forces(0,fr
->f_novirsum_n
,f
,fr
->f_novirsum
);
883 sum_forces(start
,start
+homenr
,f
,fr
->f_novirsum
);
885 if (EEL_FULL(fr
->eeltype
))
887 /* Add the mesh contribution to the virial */
888 m_add(vir_force
,fr
->vir_el_recip
,vir_force
);
892 pr_rvecs(debug
,0,"vir_force",vir_force
,DIM
);
897 /* Sum the potential energy terms from group contributions */
898 sum_epot(&(inputrec
->opts
),enerd
);
900 if (fr
->print_force
>= 0 && bDoForces
)
902 print_large_forces(stderr
,mdatoms
,cr
,step
,fr
->print_force
,x
,f
);
906 void do_constrain_first(FILE *fplog
,gmx_constr_t constr
,
907 t_inputrec
*ir
,t_mdatoms
*md
,
908 t_state
*state
,rvec
*f
,
909 t_graph
*graph
,t_commrec
*cr
,t_nrnb
*nrnb
,
910 t_forcerec
*fr
, gmx_localtop_t
*top
, tensor shake_vir
)
913 gmx_large_int_t step
;
914 double mass
,tmass
,vcm
[4];
919 snew(savex
,state
->natoms
);
922 end
= md
->homenr
+ start
;
925 fprintf(debug
,"vcm: start=%d, homenr=%d, end=%d\n",
926 start
,md
->homenr
,end
);
927 /* Do a first constrain to reset particles... */
928 step
= ir
->init_step
;
931 char buf
[STEPSTRSIZE
];
932 fprintf(fplog
,"\nConstraining the starting coordinates (step %s)\n",
933 gmx_step_str(step
,buf
));
937 /* constrain the current position */
938 constrain(NULL
,TRUE
,FALSE
,constr
,&(top
->idef
),
939 ir
,NULL
,cr
,step
,0,md
,
940 state
->x
,state
->x
,NULL
,
941 state
->box
,state
->lambda
,&dvdlambda
,
942 NULL
,NULL
,nrnb
,econqCoord
,ir
->epc
==epcMTTK
,state
->veta
,state
->veta
);
945 /* constrain the inital velocity, and save it */
946 /* also may be useful if we need the ekin from the halfstep for velocity verlet */
947 /* might not yet treat veta correctly */
948 constrain(NULL
,TRUE
,FALSE
,constr
,&(top
->idef
),
949 ir
,NULL
,cr
,step
,0,md
,
950 state
->x
,state
->v
,state
->v
,
951 state
->box
,state
->lambda
,&dvdlambda
,
952 NULL
,NULL
,nrnb
,econqVeloc
,ir
->epc
==epcMTTK
,state
->veta
,state
->veta
);
954 /* constrain the inital velocities at t-dt/2 */
955 if (EI_STATE_VELOCITY(ir
->eI
) && ir
->eI
!=eiVV
)
957 for(i
=start
; (i
<end
); i
++)
959 for(m
=0; (m
<DIM
); m
++)
961 /* Reverse the velocity */
962 state
->v
[i
][m
] = -state
->v
[i
][m
];
963 /* Store the position at t-dt in buf */
964 savex
[i
][m
] = state
->x
[i
][m
] + dt
*state
->v
[i
][m
];
967 /* Shake the positions at t=-dt with the positions at t=0
968 * as reference coordinates.
972 char buf
[STEPSTRSIZE
];
973 fprintf(fplog
,"\nConstraining the coordinates at t0-dt (step %s)\n",
974 gmx_step_str(step
,buf
));
977 constrain(NULL
,TRUE
,FALSE
,constr
,&(top
->idef
),
978 ir
,NULL
,cr
,step
,-1,md
,
980 state
->box
,state
->lambda
,&dvdlambda
,
981 state
->v
,NULL
,nrnb
,econqCoord
,ir
->epc
==epcMTTK
,state
->veta
,state
->veta
);
983 for(i
=start
; i
<end
; i
++) {
984 for(m
=0; m
<DIM
; m
++) {
985 /* Re-reverse the velocities */
986 state
->v
[i
][m
] = -state
->v
[i
][m
];
993 for(i
=start
; i
<end
; i
++) {
995 for(m
=0; m
<DIM
; m
++) {
996 vcm
[m
] += state
->v
[i
][m
]*mass
;
1001 if (ir
->nstcomm
!= 0 || debug
) {
1002 /* Compute the global sum of vcm */
1004 fprintf(debug
,"vcm: %8.3f %8.3f %8.3f,"
1005 " total mass = %12.5e\n",vcm
[XX
],vcm
[YY
],vcm
[ZZ
],vcm
[3]);
1009 for(m
=0; (m
<DIM
); m
++)
1012 fprintf(debug
,"vcm: %8.3f %8.3f %8.3f,"
1013 " total mass = %12.5e\n",vcm
[XX
],vcm
[YY
],vcm
[ZZ
],tmass
);
1014 if (ir
->nstcomm
!= 0) {
1015 /* Now we have the velocity of center of mass, let's remove it */
1016 for(i
=start
; (i
<end
); i
++) {
1017 for(m
=0; (m
<DIM
); m
++)
1018 state
->v
[i
][m
] -= vcm
[m
];
1026 void calc_enervirdiff(FILE *fplog
,int eDispCorr
,t_forcerec
*fr
)
1028 double eners
[2],virs
[2],enersum
,virsum
,y0
,f
,g
,h
;
1029 double r0
,r1
,r
,rc3
,rc9
,ea
,eb
,ec
,pa
,pb
,pc
,pd
;
1030 double invscale
,invscale2
,invscale3
;
1031 int ri0
,ri1
,ri
,i
,offstart
,offset
;
1034 fr
->enershiftsix
= 0;
1035 fr
->enershifttwelve
= 0;
1036 fr
->enerdiffsix
= 0;
1037 fr
->enerdifftwelve
= 0;
1039 fr
->virdifftwelve
= 0;
1041 if (eDispCorr
!= edispcNO
) {
1042 for(i
=0; i
<2; i
++) {
1046 if ((fr
->vdwtype
== evdwSWITCH
) || (fr
->vdwtype
== evdwSHIFT
)) {
1047 if (fr
->rvdw_switch
== 0)
1049 "With dispersion correction rvdw-switch can not be zero "
1050 "for vdw-type = %s",evdw_names
[fr
->vdwtype
]);
1052 scale
= fr
->nblists
[0].tab
.scale
;
1053 vdwtab
= fr
->nblists
[0].vdwtab
;
1055 /* Round the cut-offs to exact table values for precision */
1056 ri0
= floor(fr
->rvdw_switch
*scale
);
1057 ri1
= ceil(fr
->rvdw
*scale
);
1063 if (fr
->vdwtype
== evdwSHIFT
) {
1064 /* Determine the constant energy shift below rvdw_switch */
1065 fr
->enershiftsix
= (real
)(-1.0/(rc3
*rc3
)) - vdwtab
[8*ri0
];
1066 fr
->enershifttwelve
= (real
)( 1.0/(rc9
*rc3
)) - vdwtab
[8*ri0
+ 4];
1068 /* Add the constant part from 0 to rvdw_switch.
1069 * This integration from 0 to rvdw_switch overcounts the number
1070 * of interactions by 1, as it also counts the self interaction.
1071 * We will correct for this later.
1073 eners
[0] += 4.0*M_PI
*fr
->enershiftsix
*rc3
/3.0;
1074 eners
[1] += 4.0*M_PI
*fr
->enershifttwelve
*rc3
/3.0;
1076 invscale
= 1.0/(scale
);
1077 invscale2
= invscale
*invscale
;
1078 invscale3
= invscale
*invscale2
;
1080 /* following summation derived from cubic spline definition,
1081 Numerical Recipies in C, second edition, p. 113-116. Exact
1082 for the cubic spline. We first calculate the negative of
1083 the energy from rvdw to rvdw_switch, assuming that g(r)=1,
1084 and then add the more standard, abrupt cutoff correction to
1085 that result, yielding the long-range correction for a
1086 switched function. We perform both the pressure and energy
1087 loops at the same time for simplicity, as the computational
1091 enersum
= 0.0; virsum
= 0.0;
1096 for (ri
=ri0
; ri
<ri1
; ri
++) {
1099 eb
= 2.0*invscale2
*r
;
1103 pb
= 3.0*invscale2
*r
;
1104 pc
= 3.0*invscale
*r
*r
;
1107 /* this "8" is from the packing in the vdwtab array - perhaps
1108 should be #define'ed? */
1109 offset
= 8*ri
+ offstart
;
1110 y0
= vdwtab
[offset
];
1111 f
= vdwtab
[offset
+1];
1112 g
= vdwtab
[offset
+2];
1113 h
= vdwtab
[offset
+3];
1115 enersum
+= y0
*(ea
/3 + eb
/2 + ec
) + f
*(ea
/4 + eb
/3 + ec
/2)+
1116 g
*(ea
/5 + eb
/4 + ec
/3) + h
*(ea
/6 + eb
/5 + ec
/4);
1117 virsum
+= f
*(pa
/4 + pb
/3 + pc
/2 + pd
) +
1118 2*g
*(pa
/5 + pb
/4 + pc
/3 + pd
/2) + 3*h
*(pa
/6 + pb
/5 + pc
/4 + pd
/3);
1121 enersum
*= 4.0*M_PI
;
1123 eners
[i
] -= enersum
;
1127 /* now add the correction for rvdw_switch to infinity */
1128 eners
[0] += -4.0*M_PI
/(3.0*rc3
);
1129 eners
[1] += 4.0*M_PI
/(9.0*rc9
);
1130 virs
[0] += 8.0*M_PI
/rc3
;
1131 virs
[1] += -16.0*M_PI
/(3.0*rc9
);
1133 else if ((fr
->vdwtype
== evdwCUT
) || (fr
->vdwtype
== evdwUSER
)) {
1134 if (fr
->vdwtype
== evdwUSER
&& fplog
)
1136 "WARNING: using dispersion correction with user tables\n");
1137 rc3
= fr
->rvdw
*fr
->rvdw
*fr
->rvdw
;
1139 eners
[0] += -4.0*M_PI
/(3.0*rc3
);
1140 eners
[1] += 4.0*M_PI
/(9.0*rc9
);
1141 virs
[0] += 8.0*M_PI
/rc3
;
1142 virs
[1] += -16.0*M_PI
/(3.0*rc9
);
1145 "Dispersion correction is not implemented for vdw-type = %s",
1146 evdw_names
[fr
->vdwtype
]);
1148 fr
->enerdiffsix
= eners
[0];
1149 fr
->enerdifftwelve
= eners
[1];
1150 /* The 0.5 is due to the Gromacs definition of the virial */
1151 fr
->virdiffsix
= 0.5*virs
[0];
1152 fr
->virdifftwelve
= 0.5*virs
[1];
1156 void calc_dispcorr(FILE *fplog
,t_inputrec
*ir
,t_forcerec
*fr
,
1157 gmx_large_int_t step
,int natoms
,
1158 matrix box
,real lambda
,tensor pres
,tensor virial
,
1159 real
*prescorr
, real
*enercorr
, real
*dvdlcorr
)
1161 bool bCorrAll
,bCorrPres
;
1162 real dvdlambda
,invvol
,dens
,ninter
,avcsix
,avctwelve
,enerdiff
,svir
=0,spres
=0;
1172 if (ir
->eDispCorr
!= edispcNO
) {
1173 bCorrAll
= (ir
->eDispCorr
== edispcAllEner
||
1174 ir
->eDispCorr
== edispcAllEnerPres
);
1175 bCorrPres
= (ir
->eDispCorr
== edispcEnerPres
||
1176 ir
->eDispCorr
== edispcAllEnerPres
);
1178 invvol
= 1/det(box
);
1181 /* Only correct for the interactions with the inserted molecule */
1182 dens
= (natoms
- fr
->n_tpi
)*invvol
;
1187 dens
= natoms
*invvol
;
1188 ninter
= 0.5*natoms
;
1191 if (ir
->efep
== efepNO
)
1193 avcsix
= fr
->avcsix
[0];
1194 avctwelve
= fr
->avctwelve
[0];
1198 avcsix
= (1 - lambda
)*fr
->avcsix
[0] + lambda
*fr
->avcsix
[1];
1199 avctwelve
= (1 - lambda
)*fr
->avctwelve
[0] + lambda
*fr
->avctwelve
[1];
1202 enerdiff
= ninter
*(dens
*fr
->enerdiffsix
- fr
->enershiftsix
);
1203 *enercorr
+= avcsix
*enerdiff
;
1205 if (ir
->efep
!= efepNO
)
1207 dvdlambda
+= (fr
->avcsix
[1] - fr
->avcsix
[0])*enerdiff
;
1211 enerdiff
= ninter
*(dens
*fr
->enerdifftwelve
- fr
->enershifttwelve
);
1212 *enercorr
+= avctwelve
*enerdiff
;
1213 if (fr
->efep
!= efepNO
)
1215 dvdlambda
+= (fr
->avctwelve
[1] - fr
->avctwelve
[0])*enerdiff
;
1221 svir
= ninter
*dens
*avcsix
*fr
->virdiffsix
/3.0;
1222 if (ir
->eDispCorr
== edispcAllEnerPres
)
1224 svir
+= ninter
*dens
*avctwelve
*fr
->virdifftwelve
/3.0;
1226 /* The factor 2 is because of the Gromacs virial definition */
1227 spres
= -2.0*invvol
*svir
*PRESFAC
;
1229 for(m
=0; m
<DIM
; m
++) {
1230 virial
[m
][m
] += svir
;
1231 pres
[m
][m
] += spres
;
1236 /* Can't currently control when it prints, for now, just print when degugging */
1240 fprintf(debug
,"Long Range LJ corr.: <C6> %10.4e, <C12> %10.4e\n",
1246 "Long Range LJ corr.: Epot %10g, Pres: %10g, Vir: %10g\n",
1247 *enercorr
,spres
,svir
);
1251 fprintf(debug
,"Long Range LJ corr.: Epot %10g\n",*enercorr
);
1255 if (fr
->bSepDVDL
&& do_per_step(step
,ir
->nstlog
))
1257 fprintf(fplog
,sepdvdlformat
,"Dispersion correction",
1258 *enercorr
,dvdlambda
);
1260 if (fr
->efep
!= efepNO
)
1262 *dvdlcorr
+= dvdlambda
;
1267 void do_pbc_first(FILE *fplog
,matrix box
,t_forcerec
*fr
,
1268 t_graph
*graph
,rvec x
[])
1271 fprintf(fplog
,"Removing pbc first time\n");
1272 calc_shifts(box
,fr
->shift_vec
);
1274 mk_mshift(fplog
,graph
,fr
->ePBC
,box
,x
);
1276 p_graph(debug
,"do_pbc_first 1",graph
);
1277 shift_self(graph
,box
,x
);
1278 /* By doing an extra mk_mshift the molecules that are broken
1279 * because they were e.g. imported from another software
1280 * will be made whole again. Such are the healing powers
1283 mk_mshift(fplog
,graph
,fr
->ePBC
,box
,x
);
1285 p_graph(debug
,"do_pbc_first 2",graph
);
1288 fprintf(fplog
,"Done rmpbc\n");
1291 static void low_do_pbc_mtop(FILE *fplog
,int ePBC
,matrix box
,
1292 gmx_mtop_t
*mtop
,rvec x
[],
1297 gmx_molblock_t
*molb
;
1299 if (bFirst
&& fplog
)
1300 fprintf(fplog
,"Removing pbc first time\n");
1304 for(mb
=0; mb
<mtop
->nmolblock
; mb
++) {
1305 molb
= &mtop
->molblock
[mb
];
1306 if (molb
->natoms_mol
== 1 ||
1307 (!bFirst
&& mtop
->moltype
[molb
->type
].cgs
.nr
== 1)) {
1308 /* Just one atom or charge group in the molecule, no PBC required */
1309 as
+= molb
->nmol
*molb
->natoms_mol
;
1311 /* Pass NULL iso fplog to avoid graph prints for each molecule type */
1312 mk_graph_ilist(NULL
,mtop
->moltype
[molb
->type
].ilist
,
1313 0,molb
->natoms_mol
,FALSE
,FALSE
,graph
);
1315 for(mol
=0; mol
<molb
->nmol
; mol
++) {
1316 mk_mshift(fplog
,graph
,ePBC
,box
,x
+as
);
1318 shift_self(graph
,box
,x
+as
);
1319 /* The molecule is whole now.
1320 * We don't need the second mk_mshift call as in do_pbc_first,
1321 * since we no longer need this graph.
1324 as
+= molb
->natoms_mol
;
1332 void do_pbc_first_mtop(FILE *fplog
,int ePBC
,matrix box
,
1333 gmx_mtop_t
*mtop
,rvec x
[])
1335 low_do_pbc_mtop(fplog
,ePBC
,box
,mtop
,x
,TRUE
);
1338 void do_pbc_mtop(FILE *fplog
,int ePBC
,matrix box
,
1339 gmx_mtop_t
*mtop
,rvec x
[])
1341 low_do_pbc_mtop(fplog
,ePBC
,box
,mtop
,x
,FALSE
);
1344 void finish_run(FILE *fplog
,t_commrec
*cr
,const char *confout
,
1345 t_inputrec
*inputrec
,
1346 t_nrnb nrnb
[],gmx_wallcycle_t wcycle
,
1347 gmx_runtime_t
*runtime
,
1351 t_nrnb
*nrnb_tot
=NULL
;
1354 double cycles
[ewcNR
];
1356 wallcycle_sum(cr
,wcycle
,cycles
);
1358 if (cr
->nnodes
> 1) {
1362 MPI_Reduce(nrnb
->n
,nrnb_tot
->n
,eNRNB
,MPI_DOUBLE
,MPI_SUM
,
1363 MASTERRANK(cr
),cr
->mpi_comm_mysim
);
1369 if (SIMMASTER(cr
)) {
1370 print_flop(fplog
,nrnb_tot
,&nbfs
,&mflop
);
1371 if (cr
->nnodes
> 1) {
1376 if ((cr
->duty
& DUTY_PP
) && DOMAINDECOMP(cr
)) {
1377 print_dd_statistics(cr
,inputrec
,fplog
);
1380 if (SIMMASTER(cr
)) {
1381 if (PARTDECOMP(cr
)) {
1382 pr_load(fplog
,cr
,nrnb_tot
);
1385 wallcycle_print(fplog
,cr
->nnodes
,cr
->npmenodes
,runtime
->realtime
,
1388 if (EI_DYNAMICS(inputrec
->eI
)) {
1389 delta_t
= inputrec
->delta_t
;
1395 print_perf(fplog
,runtime
->proctime
,runtime
->realtime
,
1396 cr
->nnodes
-cr
->npmenodes
,
1397 runtime
->nsteps_done
,delta_t
,nbfs
,mflop
);
1400 print_perf(stderr
,runtime
->proctime
,runtime
->realtime
,
1401 cr
->nnodes
-cr
->npmenodes
,
1402 runtime
->nsteps_done
,delta_t
,nbfs
,mflop
);
1406 runtime=inputrec->nsteps*inputrec->delta_t;
1408 if (cr->nnodes == 1)
1409 fprintf(stderr,"\n\n");
1410 print_perf(stderr,nodetime,realtime,runtime,&ntot,
1411 cr->nnodes-cr->npmenodes,FALSE);
1413 wallcycle_print(fplog,cr->nnodes,cr->npmenodes,realtime,wcycle,cycles);
1414 print_perf(fplog,nodetime,realtime,runtime,&ntot,cr->nnodes-cr->npmenodes,
1417 pr_load(fplog,cr,nrnb_all);
1424 void init_md(FILE *fplog
,
1425 t_commrec
*cr
,t_inputrec
*ir
,const output_env_t oenv
,
1426 double *t
,double *t0
,
1427 real
*lambda
,double *lam0
,
1428 t_nrnb
*nrnb
,gmx_mtop_t
*mtop
,
1430 int nfile
,const t_filenm fnm
[],
1431 gmx_mdoutf_t
**outf
,t_mdebin
**mdebin
,
1432 tensor force_vir
,tensor shake_vir
,rvec mu_tot
,
1433 bool *bNEMD
,bool *bSimAnn
,t_vcm
**vcm
, t_state
*state
, unsigned long Flags
)
1438 /* Initial values */
1439 *t
= *t0
= ir
->init_t
;
1440 if (ir
->efep
!= efepNO
)
1442 *lam0
= ir
->init_lambda
;
1443 *lambda
= *lam0
+ ir
->init_step
*ir
->delta_lambda
;
1447 *lambda
= *lam0
= 0.0;
1451 for(i
=0;i
<ir
->opts
.ngtc
;i
++)
1453 /* set bSimAnn if any group is being annealed */
1454 if(ir
->opts
.annealing
[i
]!=eannNO
)
1461 update_annealing_target_temp(&(ir
->opts
),ir
->init_t
);
1464 *bNEMD
= (ir
->opts
.ngacc
> 1) || (norm(ir
->opts
.acc
[0]) > 0);
1468 *upd
= init_update(fplog
,ir
);
1473 *vcm
= init_vcm(fplog
,&mtop
->groups
,ir
);
1476 if (EI_DYNAMICS(ir
->eI
) && !(Flags
& MD_APPENDFILES
))
1478 if (ir
->etc
== etcBERENDSEN
)
1480 please_cite(fplog
,"Berendsen84a");
1482 if (ir
->etc
== etcVRESCALE
)
1484 please_cite(fplog
,"Bussi2007a");
1492 *outf
= init_mdoutf(nfile
,fnm
,(Flags
& MD_APPENDFILES
),cr
,ir
,oenv
);
1494 *mdebin
= init_mdebin((Flags
& MD_APPENDFILES
) ? NULL
: (*outf
)->fp_ene
,
1498 /* Initiate variables */
1499 clear_mat(force_vir
);
1500 clear_mat(shake_vir
);