3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Gromacs Runs On Most of All Computer Systems
54 #include "types/membedt.h"
55 #include "types/globsig.h"
59 #include "thread_mpi/threads.h"
66 #define MD_POLARISE (1<<2)
67 #define MD_IONIZE (1<<3)
68 #define MD_RERUN (1<<4)
69 #define MD_RERUN_VSITE (1<<5)
70 #define MD_FFSCAN (1<<6)
71 #define MD_SEPPOT (1<<7)
72 #define MD_PARTDEC (1<<9)
73 #define MD_DDBONDCHECK (1<<10)
74 #define MD_DDBONDCOMM (1<<11)
75 #define MD_CONFOUT (1<<12)
76 #define MD_REPRODUCIBLE (1<<13)
77 #define MD_READ_RNG (1<<14)
78 #define MD_APPENDFILES (1<<15)
79 #define MD_APPENDFILESSET (1<<21)
80 #define MD_KEEPANDNUMCPT (1<<16)
81 #define MD_READ_EKIN (1<<17)
82 #define MD_STARTFROMCPT (1<<18)
83 #define MD_RESETCOUNTERSHALFWAY (1<<19)
84 #define MD_TUNEPME (1<<20)
85 #define MD_TESTVERLET (1<<22)
87 /* The options for the domain decomposition MPI task ordering */
89 ddnoSEL, ddnoINTERLEAVE, ddnoPP_PME, ddnoCARTESIAN, ddnoNR
92 /* The options for the thread affinity setting, default: auto */
94 threadaffSEL, threadaffAUTO, threadaffON, threadaffOFF, threadaffNR
98 int nthreads_tot; /* Total number of threads requested (TMPI) */
99 int nthreads_tmpi; /* Number of TMPI threads requested */
100 int nthreads_omp; /* Number of OpenMP threads requested */
101 int nthreads_omp_pme; /* As nthreads_omp, but for PME only nodes */
102 int thread_affinity; /* Thread affinity switch, see enum above */
103 int core_pinning_stride; /* Logical core pinning stride */
104 int core_pinning_offset; /* Logical core pinning offset */
105 char *gpu_id; /* GPU id's to use, each specified as chars */
108 /* Variables for temporary use with the deform option,
109 * used in runner.c and md.c.
110 * (These variables should be stored in the tpx file.)
112 extern gmx_large_int_t deform_init_init_step_tpx;
113 extern matrix deform_init_box_tpx;
114 #ifdef GMX_THREAD_MPI
115 extern tMPI_Thread_mutex_t deform_init_box_mutex;
117 /* The minimum number of atoms per tMPI thread. With fewer atoms than this,
118 * the number of threads will get lowered.
120 #define MIN_ATOMS_PER_MPI_THREAD 90
121 #define MIN_ATOMS_PER_GPU 900
125 typedef double gmx_integrator_t (FILE *log, t_commrec *cr,
126 int nfile, const t_filenm fnm[],
127 const output_env_t oenv, gmx_bool bVerbose,
128 gmx_bool bCompact, int nstglobalcomm,
129 gmx_vsite_t *vsite, gmx_constr_t constr,
131 t_inputrec *inputrec,
132 gmx_mtop_t *mtop, t_fcdata *fcd,
135 t_nrnb *nrnb, gmx_wallcycle_t wcycle,
138 int repl_ex_nst, int repl_ex_nex, int repl_ex_seed,
140 real cpt_period, real max_hours,
141 const char *deviceOptions,
143 gmx_runtime_t *runtime);
145 /* ROUTINES from md.c */
147 gmx_integrator_t do_md;
150 /* ROUTINES from minimize.c */
152 gmx_integrator_t do_steep;
153 /* Do steepest descents EM */
155 gmx_integrator_t do_cg;
156 /* Do conjugate gradient EM */
158 gmx_integrator_t do_lbfgs;
159 /* Do conjugate gradient L-BFGS */
161 gmx_integrator_t do_nm;
162 /* Do normal mode analysis */
164 /* ROUTINES from tpi.c */
166 gmx_integrator_t do_tpi;
167 /* Do test particle insertion */
169 void init_npt_masses(t_inputrec *ir, t_state *state, t_extmass *MassQ, gmx_bool bInit);
171 int ExpandedEnsembleDynamics(FILE *log, t_inputrec *ir, gmx_enerdata_t *enerd,
172 t_state *state, t_extmass *MassQ, df_history_t *dfhist,
173 gmx_large_int_t step, gmx_rng_t mcrng,
174 rvec *v, t_mdatoms *mdatoms);
176 void PrintFreeEnergyInfoToFile(FILE *outfile, t_lambda *fep, t_expanded *expand, t_simtemp *simtemp, df_history_t *dfhist,
177 int nlam, int frequency, gmx_large_int_t step);
179 void get_mc_state(gmx_rng_t rng, t_state *state);
181 void set_mc_state(gmx_rng_t rng, t_state *state);
183 /* check the version */
184 void check_ir_old_tpx_versions(t_commrec *cr, FILE *fplog,
185 t_inputrec *ir, gmx_mtop_t *mtop);
187 /* Allocate and initialize node-local state entries. */
188 void set_state_entries(t_state *state, const t_inputrec *ir, int nnodes);
190 /* Broadcast the data for a simulation, and allocate node-specific settings
191 such as rng generators. */
192 void init_parallel(FILE *log, t_commrec *cr, t_inputrec *inputrec,
195 int mdrunner(gmx_hw_opt_t *hw_opt,
196 FILE *fplog, t_commrec *cr, int nfile,
197 const t_filenm fnm[], const output_env_t oenv, gmx_bool bVerbose,
198 gmx_bool bCompact, int nstglobalcomm, ivec ddxyz, int dd_node_order,
199 real rdd, real rconstr, const char *dddlb_opt, real dlb_scale,
200 const char *ddcsx, const char *ddcsy, const char *ddcsz,
201 const char *nbpu_opt,
202 int nsteps_cmdline, int nstepout, int resetstep,
203 int nmultisim, int repl_ex_nst, int repl_ex_nex,
204 int repl_ex_seed, real pforce, real cpt_period, real max_hours,
205 const char *deviceOptions, unsigned long Flags);
206 /* Driver routine, that calls the different methods */
212 #endif /* _mdrun_h */