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, 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! */
46 #include "thread_mpi/threads.h"
49 #include "gromacs/math/vec.h"
50 #include "gromacs/pbcutil/pbc.h"
51 #include "gromacs/random/random.h"
52 #include "gromacs/topology/atoms.h"
53 #include "gromacs/topology/block.h"
54 #include "gromacs/topology/symtab.h"
55 #include "gromacs/utility/smalloc.h"
57 /* The source code in this file should be thread-safe.
58 Please keep it that way. */
62 static gmx_bool bOverAllocDD = FALSE;
63 static tMPI_Thread_mutex_t over_alloc_mutex = TMPI_THREAD_MUTEX_INITIALIZER;
66 void set_over_alloc_dd(gmx_bool set)
68 tMPI_Thread_mutex_lock(&over_alloc_mutex);
69 /* we just make sure that we don't set this at the same time.
70 We don't worry too much about reading this rarely-set variable */
72 tMPI_Thread_mutex_unlock(&over_alloc_mutex);
75 int over_alloc_dd(int n)
79 return OVER_ALLOC_FAC*n + 100;
87 int gmx_int64_to_int(gmx_int64_t step, const char *warn)
93 if (warn != NULL && (step < INT_MIN || step > INT_MAX))
95 fprintf(stderr, "\nWARNING during %s:\n", warn);
96 fprintf(stderr, "step value ");
97 fprintf(stderr, "%"GMX_PRId64, step);
98 fprintf(stderr, " does not fit in int, converted to %d\n\n", i);
104 char *gmx_step_str(gmx_int64_t i, char *buf)
106 sprintf(buf, "%"GMX_PRId64, i);
111 void init_groups(gmx_groups_t *groups)
115 groups->ngrpname = 0;
116 groups->grpname = NULL;
117 for (g = 0; (g < egcNR); g++)
119 groups->grps[g].nm_ind = NULL;
120 groups->ngrpnr[g] = 0;
121 groups->grpnr[g] = NULL;
126 void init_mtop(gmx_mtop_t *mtop)
130 mtop->moltype = NULL;
132 mtop->molblock = NULL;
133 mtop->maxres_renum = 0;
135 init_groups(&mtop->groups);
136 init_block(&mtop->mols);
137 open_symtab(&mtop->symtab);
140 void init_top(t_topology *top)
145 init_atom (&(top->atoms));
146 init_atomtypes(&(top->atomtypes));
147 init_block(&top->cgs);
148 init_block(&top->mols);
149 init_blocka(&top->excls);
150 open_symtab(&top->symtab);
153 void init_inputrec(t_inputrec *ir)
155 memset(ir, 0, (size_t)sizeof(*ir));
156 snew(ir->fepvals, 1);
157 snew(ir->expandedvals, 1);
158 snew(ir->simtempvals, 1);
161 void done_moltype(gmx_moltype_t *molt)
165 done_atom(&molt->atoms);
166 done_block(&molt->cgs);
167 done_blocka(&molt->excls);
169 for (f = 0; f < F_NRE; f++)
171 sfree(molt->ilist[f].iatoms);
172 molt->ilist[f].nalloc = 0;
176 void done_molblock(gmx_molblock_t *molb)
178 if (molb->nposres_xA > 0)
180 molb->nposres_xA = 0;
181 sfree(molb->posres_xA);
183 if (molb->nposres_xB > 0)
185 molb->nposres_xB = 0;
186 sfree(molb->posres_xB);
190 void done_mtop(gmx_mtop_t *mtop, gmx_bool bDoneSymtab)
196 done_symtab(&mtop->symtab);
199 sfree(mtop->ffparams.functype);
200 sfree(mtop->ffparams.iparams);
202 for (i = 0; i < mtop->nmoltype; i++)
204 done_moltype(&mtop->moltype[i]);
206 sfree(mtop->moltype);
207 for (i = 0; i < mtop->nmolblock; i++)
209 done_molblock(&mtop->molblock[i]);
211 sfree(mtop->molblock);
212 done_block(&mtop->mols);
215 void done_top(t_topology *top)
219 sfree(top->idef.functype);
220 sfree(top->idef.iparams);
221 for (f = 0; f < F_NRE; ++f)
223 sfree(top->idef.il[f].iatoms);
224 top->idef.il[f].iatoms = NULL;
225 top->idef.il[f].nalloc = 0;
228 done_atom (&(top->atoms));
231 done_atomtypes(&(top->atomtypes));
233 done_symtab(&(top->symtab));
234 done_block(&(top->cgs));
235 done_block(&(top->mols));
236 done_blocka(&(top->excls));
239 static void done_pull_group(t_pull_group *pgrp)
244 sfree(pgrp->ind_loc);
246 sfree(pgrp->weight_loc);
250 static void done_pull(t_pull *pull)
254 for (i = 0; i < pull->ngroup+1; i++)
256 done_pull_group(pull->group);
257 done_pull_group(pull->dyna);
261 void done_inputrec(t_inputrec *ir)
265 for (m = 0; (m < DIM); m++)
273 sfree(ir->ex[m].phi);
281 sfree(ir->et[m].phi);
285 sfree(ir->opts.nrdf);
286 sfree(ir->opts.ref_t);
287 sfree(ir->opts.annealing);
288 sfree(ir->opts.anneal_npoints);
289 sfree(ir->opts.anneal_time);
290 sfree(ir->opts.anneal_temp);
291 sfree(ir->opts.tau_t);
293 sfree(ir->opts.nFreeze);
294 sfree(ir->opts.QMmethod);
295 sfree(ir->opts.QMbasis);
296 sfree(ir->opts.QMcharge);
297 sfree(ir->opts.QMmult);
299 sfree(ir->opts.CASorbitals);
300 sfree(ir->opts.CASelectrons);
301 sfree(ir->opts.SAon);
302 sfree(ir->opts.SAoff);
303 sfree(ir->opts.SAsteps);
304 sfree(ir->opts.bOPT);
314 static void zero_history(history_t *hist)
316 hist->disre_initf = 0;
317 hist->ndisrepairs = 0;
318 hist->disre_rm3tav = NULL;
319 hist->orire_initf = 0;
320 hist->norire_Dtav = 0;
321 hist->orire_Dtav = NULL;
324 static void zero_ekinstate(ekinstate_t *eks)
329 eks->ekinh_old = NULL;
330 eks->ekinscalef_nhc = NULL;
331 eks->ekinscaleh_nhc = NULL;
332 eks->vscale_nhc = NULL;
337 static void init_swapstate(swapstate_t *swapstate)
341 swapstate->eSwapCoords = 0;
342 swapstate->nAverage = 0;
344 /* Ion/water position swapping */
345 for (ic = 0; ic < eCompNR; ic++)
347 for (ii = 0; ii < eIonNR; ii++)
349 swapstate->nat_req[ic][ii] = 0;
350 swapstate->nat_req_p[ic][ii] = NULL;
351 swapstate->inflow_netto[ic][ii] = 0;
352 swapstate->inflow_netto_p[ic][ii] = NULL;
353 swapstate->nat_past[ic][ii] = NULL;
354 swapstate->nat_past_p[ic][ii] = NULL;
355 swapstate->fluxfromAtoB[ic][ii] = 0;
356 swapstate->fluxfromAtoB_p[ic][ii] = NULL;
359 swapstate->fluxleak = NULL;
360 swapstate->nions = 0;
361 swapstate->comp_from = NULL;
362 swapstate->channel_label = NULL;
363 swapstate->bFromCpt = 0;
364 swapstate->nat[eChan0] = 0;
365 swapstate->nat[eChan1] = 0;
366 swapstate->xc_old_whole[eChan0] = NULL;
367 swapstate->xc_old_whole[eChan1] = NULL;
368 swapstate->xc_old_whole_p[eChan0] = NULL;
369 swapstate->xc_old_whole_p[eChan1] = NULL;
372 void init_energyhistory(energyhistory_t * enerhist)
376 enerhist->ener_ave = NULL;
377 enerhist->ener_sum = NULL;
378 enerhist->ener_sum_sim = NULL;
379 enerhist->dht = NULL;
381 enerhist->nsteps = 0;
383 enerhist->nsteps_sim = 0;
384 enerhist->nsum_sim = 0;
386 enerhist->dht = NULL;
389 static void done_delta_h_history(delta_h_history_t *dht)
393 for (i = 0; i < dht->nndh; i++)
401 void done_energyhistory(energyhistory_t * enerhist)
403 sfree(enerhist->ener_ave);
404 sfree(enerhist->ener_sum);
405 sfree(enerhist->ener_sum_sim);
407 if (enerhist->dht != NULL)
409 done_delta_h_history(enerhist->dht);
410 sfree(enerhist->dht);
414 void init_gtc_state(t_state *state, int ngtc, int nnhpres, int nhchainlength)
419 state->nnhpres = nnhpres;
420 state->nhchainlength = nhchainlength;
423 snew(state->nosehoover_xi, state->nhchainlength*state->ngtc);
424 snew(state->nosehoover_vxi, state->nhchainlength*state->ngtc);
425 snew(state->therm_integral, state->ngtc);
426 for (i = 0; i < state->ngtc; i++)
428 for (j = 0; j < state->nhchainlength; j++)
430 state->nosehoover_xi[i*state->nhchainlength + j] = 0.0;
431 state->nosehoover_vxi[i*state->nhchainlength + j] = 0.0;
434 for (i = 0; i < state->ngtc; i++)
436 state->therm_integral[i] = 0.0;
441 state->nosehoover_xi = NULL;
442 state->nosehoover_vxi = NULL;
443 state->therm_integral = NULL;
446 if (state->nnhpres > 0)
448 snew(state->nhpres_xi, state->nhchainlength*nnhpres);
449 snew(state->nhpres_vxi, state->nhchainlength*nnhpres);
450 for (i = 0; i < nnhpres; i++)
452 for (j = 0; j < state->nhchainlength; j++)
454 state->nhpres_xi[i*nhchainlength + j] = 0.0;
455 state->nhpres_vxi[i*nhchainlength + j] = 0.0;
461 state->nhpres_xi = NULL;
462 state->nhpres_vxi = NULL;
467 void init_state(t_state *state, int natoms, int ngtc, int nnhpres, int nhchainlength, int nlambda)
471 state->natoms = natoms;
474 snew(state->lambda, efptNR);
475 for (i = 0; i < efptNR; i++)
477 state->lambda[i] = 0;
480 clear_mat(state->box);
481 clear_mat(state->box_rel);
482 clear_mat(state->boxv);
483 clear_mat(state->pres_prev);
484 clear_mat(state->svir_prev);
485 clear_mat(state->fvir_prev);
486 init_gtc_state(state, ngtc, nnhpres, nhchainlength);
487 state->nalloc = state->natoms;
488 if (state->nalloc > 0)
490 snew(state->x, state->nalloc);
491 snew(state->v, state->nalloc);
500 zero_history(&state->hist);
501 zero_ekinstate(&state->ekinstate);
502 init_energyhistory(&state->enerhist);
503 init_df_history(&state->dfhist, nlambda);
504 init_swapstate(&state->swapstate);
505 state->ddp_count = 0;
506 state->ddp_count_cg_gl = 0;
508 state->cg_gl_nalloc = 0;
511 void done_state(t_state *state)
534 state->cg_gl_nalloc = 0;
537 sfree(state->lambda);
541 sfree(state->nosehoover_xi);
542 sfree(state->nosehoover_vxi);
543 sfree(state->therm_integral);
547 t_state *serial_init_local_state(t_state *state_global)
550 t_state *state_local;
552 snew(state_local, 1);
554 /* Copy all the contents */
555 *state_local = *state_global;
556 snew(state_local->lambda, efptNR);
557 /* local storage for lambda */
558 for (i = 0; i < efptNR; i++)
560 state_local->lambda[i] = state_global->lambda[i];
566 static void do_box_rel(t_inputrec *ir, matrix box_rel, matrix b, gmx_bool bInit)
570 for (d = YY; d <= ZZ; d++)
572 for (d2 = XX; d2 <= (ir->epct == epctSEMIISOTROPIC ? YY : ZZ); d2++)
574 /* We need to check if this box component is deformed
575 * or if deformation of another component might cause
576 * changes in this component due to box corrections.
578 if (ir->deform[d][d2] == 0 &&
579 !(d == ZZ && d2 == XX && ir->deform[d][YY] != 0 &&
580 (b[YY][d2] != 0 || ir->deform[YY][d2] != 0)))
584 box_rel[d][d2] = b[d][d2]/b[XX][XX];
588 b[d][d2] = b[XX][XX]*box_rel[d][d2];
595 void set_box_rel(t_inputrec *ir, t_state *state)
597 /* Make sure the box obeys the restrictions before we fix the ratios */
598 correct_box(NULL, 0, state->box, NULL);
600 clear_mat(state->box_rel);
602 if (PRESERVE_SHAPE(*ir))
604 do_box_rel(ir, state->box_rel, state->box, TRUE);
608 void preserve_box_shape(t_inputrec *ir, matrix box_rel, matrix b)
610 if (PRESERVE_SHAPE(*ir))
612 do_box_rel(ir, box_rel, b, FALSE);
616 real max_cutoff(real cutoff1, real cutoff2)
618 if (cutoff1 == 0 || cutoff2 == 0)
624 return max(cutoff1, cutoff2);
628 void init_df_history(df_history_t *dfhist, int nlambda)
632 dfhist->nlambda = nlambda;
634 dfhist->wl_delta = 0;
638 snew(dfhist->sum_weights, dfhist->nlambda);
639 snew(dfhist->sum_dg, dfhist->nlambda);
640 snew(dfhist->sum_minvar, dfhist->nlambda);
641 snew(dfhist->sum_variance, dfhist->nlambda);
642 snew(dfhist->n_at_lam, dfhist->nlambda);
643 snew(dfhist->wl_histo, dfhist->nlambda);
645 /* allocate transition matrices here */
646 snew(dfhist->Tij, dfhist->nlambda);
647 snew(dfhist->Tij_empirical, dfhist->nlambda);
649 /* allocate accumulators for various transition matrix
650 free energy methods here */
651 snew(dfhist->accum_p, dfhist->nlambda);
652 snew(dfhist->accum_m, dfhist->nlambda);
653 snew(dfhist->accum_p2, dfhist->nlambda);
654 snew(dfhist->accum_m2, dfhist->nlambda);
656 for (i = 0; i < dfhist->nlambda; i++)
658 snew(dfhist->Tij[i], dfhist->nlambda);
659 snew(dfhist->Tij_empirical[i], dfhist->nlambda);
660 snew((dfhist->accum_p)[i], dfhist->nlambda);
661 snew((dfhist->accum_m)[i], dfhist->nlambda);
662 snew((dfhist->accum_p2)[i], dfhist->nlambda);
663 snew((dfhist->accum_m2)[i], dfhist->nlambda);
668 extern void copy_df_history(df_history_t *df_dest, df_history_t *df_source)
672 /* Currently, there should not be any difference in nlambda between the two,
673 but this is included for completeness for potential later functionality */
674 df_dest->nlambda = df_source->nlambda;
675 df_dest->bEquil = df_source->bEquil;
676 df_dest->wl_delta = df_source->wl_delta;
678 for (i = 0; i < df_dest->nlambda; i++)
680 df_dest->sum_weights[i] = df_source->sum_weights[i];
681 df_dest->sum_dg[i] = df_source->sum_dg[i];
682 df_dest->sum_minvar[i] = df_source->sum_minvar[i];
683 df_dest->sum_variance[i] = df_source->sum_variance[i];
684 df_dest->n_at_lam[i] = df_source->n_at_lam[i];
685 df_dest->wl_histo[i] = df_source->wl_histo[i];
688 for (i = 0; i < df_dest->nlambda; i++)
690 for (j = 0; j < df_dest->nlambda; j++)
692 df_dest->accum_p[i][j] = df_source->accum_p[i][j];
693 df_dest->accum_m[i][j] = df_source->accum_m[i][j];
694 df_dest->accum_p2[i][j] = df_source->accum_p2[i][j];
695 df_dest->accum_m2[i][j] = df_source->accum_m2[i][j];
696 df_dest->Tij[i][j] = df_source->Tij[i][j];
697 df_dest->Tij_empirical[i][j] = df_source->Tij_empirical[i][j];
702 void done_df_history(df_history_t *dfhist)
706 if (dfhist->nlambda > 0)
708 sfree(dfhist->n_at_lam);
709 sfree(dfhist->wl_histo);
710 sfree(dfhist->sum_weights);
711 sfree(dfhist->sum_dg);
712 sfree(dfhist->sum_minvar);
713 sfree(dfhist->sum_variance);
715 for (i = 0; i < dfhist->nlambda; i++)
717 sfree(dfhist->Tij[i]);
718 sfree(dfhist->Tij_empirical[i]);
719 sfree(dfhist->accum_p[i]);
720 sfree(dfhist->accum_m[i]);
721 sfree(dfhist->accum_p2[i]);
722 sfree(dfhist->accum_m2[i]);
727 dfhist->wl_delta = 0;