2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
37 /* This file is completely threadsafe - keep it that way! */
40 #include "gromacs/legacyheaders/typedefs.h"
44 #include "gromacs/legacyheaders/macros.h"
45 #include "gromacs/math/vec.h"
46 #include "gromacs/pbcutil/pbc.h"
47 #include "gromacs/random/random.h"
48 #include "gromacs/utility/smalloc.h"
50 /* The source code in this file should be thread-safe.
51 Please keep it that way. */
53 int gmx_int64_to_int(gmx_int64_t step, const char *warn)
59 if (warn != NULL && (step < INT_MIN || step > INT_MAX))
61 fprintf(stderr, "\nWARNING during %s:\n", warn);
62 fprintf(stderr, "step value ");
63 fprintf(stderr, "%"GMX_PRId64, step);
64 fprintf(stderr, " does not fit in int, converted to %d\n\n", i);
70 void init_inputrec(t_inputrec *ir)
72 memset(ir, 0, (size_t)sizeof(*ir));
74 snew(ir->expandedvals, 1);
75 snew(ir->simtempvals, 1);
78 static void done_pull_group(t_pull_group *pgrp)
87 static void done_pull_params(pull_params_t *pull)
91 for (i = 0; i < pull->ngroup+1; i++)
93 done_pull_group(pull->group);
100 void done_inputrec(t_inputrec *ir)
104 for (m = 0; (m < DIM); m++)
112 sfree(ir->ex[m].phi);
120 sfree(ir->et[m].phi);
124 sfree(ir->opts.nrdf);
125 sfree(ir->opts.ref_t);
126 sfree(ir->opts.annealing);
127 sfree(ir->opts.anneal_npoints);
128 sfree(ir->opts.anneal_time);
129 sfree(ir->opts.anneal_temp);
130 sfree(ir->opts.tau_t);
132 sfree(ir->opts.nFreeze);
133 sfree(ir->opts.QMmethod);
134 sfree(ir->opts.QMbasis);
135 sfree(ir->opts.QMcharge);
136 sfree(ir->opts.QMmult);
138 sfree(ir->opts.CASorbitals);
139 sfree(ir->opts.CASelectrons);
140 sfree(ir->opts.SAon);
141 sfree(ir->opts.SAoff);
142 sfree(ir->opts.SAsteps);
143 sfree(ir->opts.bOPT);
148 done_pull_params(ir->pull);
153 static void zero_history(history_t *hist)
155 hist->disre_initf = 0;
156 hist->ndisrepairs = 0;
157 hist->disre_rm3tav = NULL;
158 hist->orire_initf = 0;
159 hist->norire_Dtav = 0;
160 hist->orire_Dtav = NULL;
163 static void zero_ekinstate(ekinstate_t *eks)
168 eks->ekinh_old = NULL;
169 eks->ekinscalef_nhc = NULL;
170 eks->ekinscaleh_nhc = NULL;
171 eks->vscale_nhc = NULL;
176 static void init_swapstate(swapstate_t *swapstate)
180 swapstate->eSwapCoords = 0;
181 swapstate->nAverage = 0;
183 /* Ion/water position swapping */
184 for (ic = 0; ic < eCompNR; ic++)
186 for (ii = 0; ii < eIonNR; ii++)
188 swapstate->nat_req[ic][ii] = 0;
189 swapstate->nat_req_p[ic][ii] = NULL;
190 swapstate->inflow_netto[ic][ii] = 0;
191 swapstate->inflow_netto_p[ic][ii] = NULL;
192 swapstate->nat_past[ic][ii] = NULL;
193 swapstate->nat_past_p[ic][ii] = NULL;
194 swapstate->fluxfromAtoB[ic][ii] = 0;
195 swapstate->fluxfromAtoB_p[ic][ii] = NULL;
198 swapstate->fluxleak = NULL;
199 swapstate->nions = 0;
200 swapstate->comp_from = NULL;
201 swapstate->channel_label = NULL;
202 swapstate->bFromCpt = 0;
203 swapstate->nat[eChan0] = 0;
204 swapstate->nat[eChan1] = 0;
205 swapstate->xc_old_whole[eChan0] = NULL;
206 swapstate->xc_old_whole[eChan1] = NULL;
207 swapstate->xc_old_whole_p[eChan0] = NULL;
208 swapstate->xc_old_whole_p[eChan1] = NULL;
211 void init_energyhistory(energyhistory_t * enerhist)
215 enerhist->ener_ave = NULL;
216 enerhist->ener_sum = NULL;
217 enerhist->ener_sum_sim = NULL;
218 enerhist->dht = NULL;
220 enerhist->nsteps = 0;
222 enerhist->nsteps_sim = 0;
223 enerhist->nsum_sim = 0;
225 enerhist->dht = NULL;
228 static void done_delta_h_history(delta_h_history_t *dht)
232 for (i = 0; i < dht->nndh; i++)
240 void done_energyhistory(energyhistory_t * enerhist)
242 sfree(enerhist->ener_ave);
243 sfree(enerhist->ener_sum);
244 sfree(enerhist->ener_sum_sim);
246 if (enerhist->dht != NULL)
248 done_delta_h_history(enerhist->dht);
249 sfree(enerhist->dht);
253 void init_gtc_state(t_state *state, int ngtc, int nnhpres, int nhchainlength)
258 state->nnhpres = nnhpres;
259 state->nhchainlength = nhchainlength;
262 snew(state->nosehoover_xi, state->nhchainlength*state->ngtc);
263 snew(state->nosehoover_vxi, state->nhchainlength*state->ngtc);
264 snew(state->therm_integral, state->ngtc);
265 for (i = 0; i < state->ngtc; i++)
267 for (j = 0; j < state->nhchainlength; j++)
269 state->nosehoover_xi[i*state->nhchainlength + j] = 0.0;
270 state->nosehoover_vxi[i*state->nhchainlength + j] = 0.0;
273 for (i = 0; i < state->ngtc; i++)
275 state->therm_integral[i] = 0.0;
280 state->nosehoover_xi = NULL;
281 state->nosehoover_vxi = NULL;
282 state->therm_integral = NULL;
285 if (state->nnhpres > 0)
287 snew(state->nhpres_xi, state->nhchainlength*nnhpres);
288 snew(state->nhpres_vxi, state->nhchainlength*nnhpres);
289 for (i = 0; i < nnhpres; i++)
291 for (j = 0; j < state->nhchainlength; j++)
293 state->nhpres_xi[i*nhchainlength + j] = 0.0;
294 state->nhpres_vxi[i*nhchainlength + j] = 0.0;
300 state->nhpres_xi = NULL;
301 state->nhpres_vxi = NULL;
306 void init_state(t_state *state, int natoms, int ngtc, int nnhpres, int nhchainlength, int nlambda)
310 state->natoms = natoms;
312 state->fep_state = 0;
314 snew(state->lambda, efptNR);
315 for (i = 0; i < efptNR; i++)
317 state->lambda[i] = 0;
320 clear_mat(state->box);
321 clear_mat(state->box_rel);
322 clear_mat(state->boxv);
323 clear_mat(state->pres_prev);
324 clear_mat(state->svir_prev);
325 clear_mat(state->fvir_prev);
326 init_gtc_state(state, ngtc, nnhpres, nhchainlength);
327 state->nalloc = state->natoms;
328 if (state->nalloc > 0)
330 snew(state->x, state->nalloc);
331 snew(state->v, state->nalloc);
340 zero_history(&state->hist);
341 zero_ekinstate(&state->ekinstate);
342 init_energyhistory(&state->enerhist);
343 init_df_history(&state->dfhist, nlambda);
344 init_swapstate(&state->swapstate);
345 state->ddp_count = 0;
346 state->ddp_count_cg_gl = 0;
348 state->cg_gl_nalloc = 0;
351 void done_state(t_state *state)
374 state->cg_gl_nalloc = 0;
377 sfree(state->lambda);
381 sfree(state->nosehoover_xi);
382 sfree(state->nosehoover_vxi);
383 sfree(state->therm_integral);
387 t_state *serial_init_local_state(t_state *state_global)
390 t_state *state_local;
392 snew(state_local, 1);
394 /* Copy all the contents */
395 *state_local = *state_global;
396 snew(state_local->lambda, efptNR);
397 /* local storage for lambda */
398 for (i = 0; i < efptNR; i++)
400 state_local->lambda[i] = state_global->lambda[i];
406 static void do_box_rel(t_inputrec *ir, matrix box_rel, matrix b, gmx_bool bInit)
410 for (d = YY; d <= ZZ; d++)
412 for (d2 = XX; d2 <= (ir->epct == epctSEMIISOTROPIC ? YY : ZZ); d2++)
414 /* We need to check if this box component is deformed
415 * or if deformation of another component might cause
416 * changes in this component due to box corrections.
418 if (ir->deform[d][d2] == 0 &&
419 !(d == ZZ && d2 == XX && ir->deform[d][YY] != 0 &&
420 (b[YY][d2] != 0 || ir->deform[YY][d2] != 0)))
424 box_rel[d][d2] = b[d][d2]/b[XX][XX];
428 b[d][d2] = b[XX][XX]*box_rel[d][d2];
435 void set_box_rel(t_inputrec *ir, t_state *state)
437 /* Make sure the box obeys the restrictions before we fix the ratios */
438 correct_box(NULL, 0, state->box, NULL);
440 clear_mat(state->box_rel);
442 if (PRESERVE_SHAPE(*ir))
444 do_box_rel(ir, state->box_rel, state->box, TRUE);
448 void preserve_box_shape(t_inputrec *ir, matrix box_rel, matrix b)
450 if (PRESERVE_SHAPE(*ir))
452 do_box_rel(ir, box_rel, b, FALSE);
456 real max_cutoff(real cutoff1, real cutoff2)
458 if (cutoff1 == 0 || cutoff2 == 0)
464 return max(cutoff1, cutoff2);
468 void init_df_history(df_history_t *dfhist, int nlambda)
472 dfhist->nlambda = nlambda;
474 dfhist->wl_delta = 0;
478 snew(dfhist->sum_weights, dfhist->nlambda);
479 snew(dfhist->sum_dg, dfhist->nlambda);
480 snew(dfhist->sum_minvar, dfhist->nlambda);
481 snew(dfhist->sum_variance, dfhist->nlambda);
482 snew(dfhist->n_at_lam, dfhist->nlambda);
483 snew(dfhist->wl_histo, dfhist->nlambda);
485 /* allocate transition matrices here */
486 snew(dfhist->Tij, dfhist->nlambda);
487 snew(dfhist->Tij_empirical, dfhist->nlambda);
489 /* allocate accumulators for various transition matrix
490 free energy methods here */
491 snew(dfhist->accum_p, dfhist->nlambda);
492 snew(dfhist->accum_m, dfhist->nlambda);
493 snew(dfhist->accum_p2, dfhist->nlambda);
494 snew(dfhist->accum_m2, dfhist->nlambda);
496 for (i = 0; i < dfhist->nlambda; i++)
498 snew(dfhist->Tij[i], dfhist->nlambda);
499 snew(dfhist->Tij_empirical[i], dfhist->nlambda);
500 snew((dfhist->accum_p)[i], dfhist->nlambda);
501 snew((dfhist->accum_m)[i], dfhist->nlambda);
502 snew((dfhist->accum_p2)[i], dfhist->nlambda);
503 snew((dfhist->accum_m2)[i], dfhist->nlambda);
508 extern void copy_df_history(df_history_t *df_dest, df_history_t *df_source)
512 /* Currently, there should not be any difference in nlambda between the two,
513 but this is included for completeness for potential later functionality */
514 df_dest->nlambda = df_source->nlambda;
515 df_dest->bEquil = df_source->bEquil;
516 df_dest->wl_delta = df_source->wl_delta;
518 for (i = 0; i < df_dest->nlambda; i++)
520 df_dest->sum_weights[i] = df_source->sum_weights[i];
521 df_dest->sum_dg[i] = df_source->sum_dg[i];
522 df_dest->sum_minvar[i] = df_source->sum_minvar[i];
523 df_dest->sum_variance[i] = df_source->sum_variance[i];
524 df_dest->n_at_lam[i] = df_source->n_at_lam[i];
525 df_dest->wl_histo[i] = df_source->wl_histo[i];
528 for (i = 0; i < df_dest->nlambda; i++)
530 for (j = 0; j < df_dest->nlambda; j++)
532 df_dest->accum_p[i][j] = df_source->accum_p[i][j];
533 df_dest->accum_m[i][j] = df_source->accum_m[i][j];
534 df_dest->accum_p2[i][j] = df_source->accum_p2[i][j];
535 df_dest->accum_m2[i][j] = df_source->accum_m2[i][j];
536 df_dest->Tij[i][j] = df_source->Tij[i][j];
537 df_dest->Tij_empirical[i][j] = df_source->Tij_empirical[i][j];
542 void done_df_history(df_history_t *dfhist)
546 if (dfhist->nlambda > 0)
548 sfree(dfhist->n_at_lam);
549 sfree(dfhist->wl_histo);
550 sfree(dfhist->sum_weights);
551 sfree(dfhist->sum_dg);
552 sfree(dfhist->sum_minvar);
553 sfree(dfhist->sum_variance);
555 for (i = 0; i < dfhist->nlambda; i++)
557 sfree(dfhist->Tij[i]);
558 sfree(dfhist->Tij_empirical[i]);
559 sfree(dfhist->accum_p[i]);
560 sfree(dfhist->accum_m[i]);
561 sfree(dfhist->accum_p2[i]);
562 sfree(dfhist->accum_m2[i]);
567 dfhist->wl_delta = 0;