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 * GROningen Mixture of Alchemy and Childrens' Stories
36 /* This file is completely threadsafe - keep it that way! */
49 #include "thread_mpi.h"
52 /* The source code in this file should be thread-safe.
53 Please keep it that way. */
57 static gmx_bool bOverAllocDD = FALSE;
59 static tMPI_Thread_mutex_t over_alloc_mutex = TMPI_THREAD_MUTEX_INITIALIZER;
63 void set_over_alloc_dd(gmx_bool set)
66 tMPI_Thread_mutex_lock(&over_alloc_mutex);
67 /* we just make sure that we don't set this at the same time.
68 We don't worry too much about reading this rarely-set variable */
72 tMPI_Thread_mutex_unlock(&over_alloc_mutex);
76 int over_alloc_dd(int n)
80 return OVER_ALLOC_FAC*n + 100;
88 int gmx_large_int_to_int(gmx_large_int_t step, const char *warn)
94 if (warn != NULL && (step < INT_MIN || step > INT_MAX))
96 fprintf(stderr, "\nWARNING during %s:\n", warn);
97 fprintf(stderr, "step value ");
98 fprintf(stderr, gmx_large_int_pfmt, step);
99 fprintf(stderr, " does not fit in int, converted to %d\n\n", i);
105 char *gmx_step_str(gmx_large_int_t i, char *buf)
107 sprintf(buf, gmx_large_int_pfmt, i);
112 void init_block(t_block *block)
117 block->nalloc_index = 1;
118 snew(block->index, block->nalloc_index);
122 void init_blocka(t_blocka *block)
128 block->nalloc_index = 1;
129 snew(block->index, block->nalloc_index);
135 void init_atom(t_atoms *at)
145 at->atomtypeB = NULL;
149 void init_atomtypes(t_atomtypes *at)
154 at->atomnumber = NULL;
155 at->gb_radius = NULL;
159 void init_groups(gmx_groups_t *groups)
163 groups->ngrpname = 0;
164 groups->grpname = NULL;
165 for (g = 0; (g < egcNR); g++)
167 groups->grps[g].nm_ind = NULL;
168 groups->ngrpnr[g] = 0;
169 groups->grpnr[g] = NULL;
174 void init_mtop(gmx_mtop_t *mtop)
178 mtop->moltype = NULL;
180 mtop->molblock = NULL;
181 mtop->maxres_renum = 0;
183 init_groups(&mtop->groups);
184 init_block(&mtop->mols);
185 open_symtab(&mtop->symtab);
188 void init_top (t_topology *top)
193 init_atom (&(top->atoms));
194 init_atomtypes(&(top->atomtypes));
195 init_block(&top->cgs);
196 init_block(&top->mols);
197 init_blocka(&top->excls);
198 open_symtab(&top->symtab);
201 void init_inputrec(t_inputrec *ir)
203 memset(ir, 0, (size_t)sizeof(*ir));
204 snew(ir->fepvals, 1);
205 snew(ir->expandedvals, 1);
206 snew(ir->simtempvals, 1);
209 void stupid_fill_block(t_block *grp, int natom, gmx_bool bOneIndexGroup)
215 grp->nalloc_index = 2;
216 snew(grp->index, grp->nalloc_index);
218 grp->index[1] = natom;
223 grp->nalloc_index = natom+1;
224 snew(grp->index, grp->nalloc_index);
225 snew(grp->index, natom+1);
226 for (i = 0; (i <= natom); i++)
234 void stupid_fill_blocka(t_blocka *grp, int natom)
238 grp->nalloc_a = natom;
239 snew(grp->a, grp->nalloc_a);
240 for (i = 0; (i < natom); i++)
246 grp->nalloc_index = natom + 1;
247 snew(grp->index, grp->nalloc_index);
248 for (i = 0; (i <= natom); i++)
255 void copy_blocka(const t_blocka *src, t_blocka *dest)
260 dest->nalloc_index = dest->nr + 1;
261 snew(dest->index, dest->nalloc_index);
262 for (i = 0; i < dest->nr+1; i++)
264 dest->index[i] = src->index[i];
266 dest->nra = src->nra;
267 dest->nalloc_a = dest->nra + 1;
268 snew(dest->a, dest->nalloc_a);
269 for (i = 0; i < dest->nra+1; i++)
271 dest->a[i] = src->a[i];
275 void done_block(t_block *block)
279 block->nalloc_index = 0;
282 void done_blocka(t_blocka *block)
291 block->nalloc_index = 0;
295 void done_atom (t_atoms *at)
303 sfree(at->atomtypeB);
310 void done_atomtypes(t_atomtypes *atype)
313 sfree(atype->radius);
315 sfree(atype->surftens);
316 sfree(atype->atomnumber);
317 sfree(atype->gb_radius);
321 void done_moltype(gmx_moltype_t *molt)
325 done_atom(&molt->atoms);
326 done_block(&molt->cgs);
327 done_blocka(&molt->excls);
329 for (f = 0; f < F_NRE; f++)
331 sfree(molt->ilist[f].iatoms);
332 molt->ilist[f].nalloc = 0;
336 void done_molblock(gmx_molblock_t *molb)
338 if (molb->nposres_xA > 0)
340 molb->nposres_xA = 0;
341 free(molb->posres_xA);
343 if (molb->nposres_xB > 0)
345 molb->nposres_xB = 0;
346 free(molb->posres_xB);
350 void done_mtop(gmx_mtop_t *mtop, gmx_bool bDoneSymtab)
356 done_symtab(&mtop->symtab);
359 sfree(mtop->ffparams.functype);
360 sfree(mtop->ffparams.iparams);
362 for (i = 0; i < mtop->nmoltype; i++)
364 done_moltype(&mtop->moltype[i]);
366 sfree(mtop->moltype);
367 for (i = 0; i < mtop->nmolblock; i++)
369 done_molblock(&mtop->molblock[i]);
371 sfree(mtop->molblock);
372 done_block(&mtop->mols);
375 void done_top(t_topology *top)
379 sfree(top->idef.functype);
380 sfree(top->idef.iparams);
381 for (f = 0; f < F_NRE; ++f)
383 sfree(top->idef.il[f].iatoms);
384 top->idef.il[f].iatoms = NULL;
385 top->idef.il[f].nalloc = 0;
388 done_atom (&(top->atoms));
391 done_atomtypes(&(top->atomtypes));
393 done_symtab(&(top->symtab));
394 done_block(&(top->cgs));
395 done_block(&(top->mols));
396 done_blocka(&(top->excls));
399 static void done_pullgrp(t_pullgrp *pgrp)
402 sfree(pgrp->ind_loc);
404 sfree(pgrp->weight_loc);
407 static void done_pull(t_pull *pull)
411 for (i = 0; i < pull->ngrp+1; i++)
413 done_pullgrp(pull->grp);
414 done_pullgrp(pull->dyna);
418 void done_inputrec(t_inputrec *ir)
422 for (m = 0; (m < DIM); m++)
430 sfree(ir->ex[m].phi);
438 sfree(ir->et[m].phi);
442 sfree(ir->opts.nrdf);
443 sfree(ir->opts.ref_t);
444 sfree(ir->opts.annealing);
445 sfree(ir->opts.anneal_npoints);
446 sfree(ir->opts.anneal_time);
447 sfree(ir->opts.anneal_temp);
448 sfree(ir->opts.tau_t);
450 sfree(ir->opts.nFreeze);
451 sfree(ir->opts.QMmethod);
452 sfree(ir->opts.QMbasis);
453 sfree(ir->opts.QMcharge);
454 sfree(ir->opts.QMmult);
456 sfree(ir->opts.CASorbitals);
457 sfree(ir->opts.CASelectrons);
458 sfree(ir->opts.SAon);
459 sfree(ir->opts.SAoff);
460 sfree(ir->opts.SAsteps);
461 sfree(ir->opts.bOPT);
471 static void zero_ekinstate(ekinstate_t *eks)
476 eks->ekinh_old = NULL;
477 eks->ekinscalef_nhc = NULL;
478 eks->ekinscaleh_nhc = NULL;
479 eks->vscale_nhc = NULL;
484 void init_energyhistory(energyhistory_t * enerhist)
488 enerhist->ener_ave = NULL;
489 enerhist->ener_sum = NULL;
490 enerhist->ener_sum_sim = NULL;
491 enerhist->dht = NULL;
493 enerhist->nsteps = 0;
495 enerhist->nsteps_sim = 0;
496 enerhist->nsum_sim = 0;
498 enerhist->dht = NULL;
501 static void done_delta_h_history(delta_h_history_t *dht)
505 for (i = 0; i < dht->nndh; i++)
513 void done_energyhistory(energyhistory_t * enerhist)
515 sfree(enerhist->ener_ave);
516 sfree(enerhist->ener_sum);
517 sfree(enerhist->ener_sum_sim);
519 if (enerhist->dht != NULL)
521 done_delta_h_history(enerhist->dht);
522 sfree(enerhist->dht);
526 void init_gtc_state(t_state *state, int ngtc, int nnhpres, int nhchainlength)
531 state->nnhpres = nnhpres;
532 state->nhchainlength = nhchainlength;
535 snew(state->nosehoover_xi, state->nhchainlength*state->ngtc);
536 snew(state->nosehoover_vxi, state->nhchainlength*state->ngtc);
537 snew(state->therm_integral, state->ngtc);
538 for (i = 0; i < state->ngtc; i++)
540 for (j = 0; j < state->nhchainlength; j++)
542 state->nosehoover_xi[i*state->nhchainlength + j] = 0.0;
543 state->nosehoover_vxi[i*state->nhchainlength + j] = 0.0;
546 for (i = 0; i < state->ngtc; i++)
548 state->therm_integral[i] = 0.0;
553 state->nosehoover_xi = NULL;
554 state->nosehoover_vxi = NULL;
555 state->therm_integral = NULL;
558 if (state->nnhpres > 0)
560 snew(state->nhpres_xi, state->nhchainlength*nnhpres);
561 snew(state->nhpres_vxi, state->nhchainlength*nnhpres);
562 for (i = 0; i < nnhpres; i++)
564 for (j = 0; j < state->nhchainlength; j++)
566 state->nhpres_xi[i*nhchainlength + j] = 0.0;
567 state->nhpres_vxi[i*nhchainlength + j] = 0.0;
573 state->nhpres_xi = NULL;
574 state->nhpres_vxi = NULL;
579 void init_state(t_state *state, int natoms, int ngtc, int nnhpres, int nhchainlength, int nlambda)
583 state->natoms = natoms;
587 snew(state->lambda, efptNR);
588 for (i = 0; i < efptNR; i++)
590 state->lambda[i] = 0;
593 clear_mat(state->box);
594 clear_mat(state->box_rel);
595 clear_mat(state->boxv);
596 clear_mat(state->pres_prev);
597 clear_mat(state->svir_prev);
598 clear_mat(state->fvir_prev);
599 init_gtc_state(state, ngtc, nnhpres, nhchainlength);
600 state->nalloc = state->natoms;
601 if (state->nalloc > 0)
603 snew(state->x, state->nalloc);
604 snew(state->v, state->nalloc);
614 zero_ekinstate(&state->ekinstate);
616 init_energyhistory(&state->enerhist);
618 init_df_history(&state->dfhist, nlambda, 0);
620 state->ddp_count = 0;
621 state->ddp_count_cg_gl = 0;
623 state->cg_gl_nalloc = 0;
626 void done_state(t_state *state)
628 if (state->nosehoover_xi)
630 sfree(state->nosehoover_xi);
653 state->cg_gl_nalloc = 0;
656 sfree(state->lambda);
660 sfree(state->nosehoover_xi);
661 sfree(state->nosehoover_vxi);
662 sfree(state->therm_integral);
666 static void do_box_rel(t_inputrec *ir, matrix box_rel, matrix b, gmx_bool bInit)
670 for (d = YY; d <= ZZ; d++)
672 for (d2 = XX; d2 <= (ir->epct == epctSEMIISOTROPIC ? YY : ZZ); d2++)
674 /* We need to check if this box component is deformed
675 * or if deformation of another component might cause
676 * changes in this component due to box corrections.
678 if (ir->deform[d][d2] == 0 &&
679 !(d == ZZ && d2 == XX && ir->deform[d][YY] != 0 &&
680 (b[YY][d2] != 0 || ir->deform[YY][d2] != 0)))
684 box_rel[d][d2] = b[d][d2]/b[XX][XX];
688 b[d][d2] = b[XX][XX]*box_rel[d][d2];
695 void set_box_rel(t_inputrec *ir, t_state *state)
697 /* Make sure the box obeys the restrictions before we fix the ratios */
698 correct_box(NULL, 0, state->box, NULL);
700 clear_mat(state->box_rel);
702 if (PRESERVE_SHAPE(*ir))
704 do_box_rel(ir, state->box_rel, state->box, TRUE);
708 void preserve_box_shape(t_inputrec *ir, matrix box_rel, matrix b)
710 if (PRESERVE_SHAPE(*ir))
712 do_box_rel(ir, box_rel, b, FALSE);
716 void add_t_atoms(t_atoms *atoms, int natom_extra, int nres_extra)
722 srenew(atoms->atomname, atoms->nr+natom_extra);
723 srenew(atoms->atom, atoms->nr+natom_extra);
724 if (NULL != atoms->pdbinfo)
726 srenew(atoms->pdbinfo, atoms->nr+natom_extra);
728 if (NULL != atoms->atomtype)
730 srenew(atoms->atomtype, atoms->nr+natom_extra);
732 if (NULL != atoms->atomtypeB)
734 srenew(atoms->atomtypeB, atoms->nr+natom_extra);
736 for (i = atoms->nr; (i < atoms->nr+natom_extra); i++)
738 atoms->atomname[i] = NULL;
739 memset(&atoms->atom[i], 0, sizeof(atoms->atom[i]));
740 if (NULL != atoms->pdbinfo)
742 memset(&atoms->pdbinfo[i], 0, sizeof(atoms->pdbinfo[i]));
744 if (NULL != atoms->atomtype)
746 atoms->atomtype[i] = NULL;
748 if (NULL != atoms->atomtypeB)
750 atoms->atomtypeB[i] = NULL;
753 atoms->nr += natom_extra;
757 srenew(atoms->resinfo, atoms->nres+nres_extra);
758 for (i = atoms->nres; (i < atoms->nres+nres_extra); i++)
760 memset(&atoms->resinfo[i], 0, sizeof(atoms->resinfo[i]));
762 atoms->nres += nres_extra;
766 void init_t_atoms(t_atoms *atoms, int natoms, gmx_bool bPdbinfo)
770 snew(atoms->atomname, natoms);
771 atoms->atomtype = NULL;
772 atoms->atomtypeB = NULL;
773 snew(atoms->resinfo, natoms);
774 snew(atoms->atom, natoms);
777 snew(atoms->pdbinfo, natoms);
781 atoms->pdbinfo = NULL;
785 t_atoms *copy_t_atoms(t_atoms *src)
791 init_t_atoms(dst, src->nr, (NULL != src->pdbinfo));
793 if (NULL != src->atomname)
795 snew(dst->atomname, src->nr);
797 if (NULL != src->atomtype)
799 snew(dst->atomtype, src->nr);
801 if (NULL != src->atomtypeB)
803 snew(dst->atomtypeB, src->nr);
805 for (i = 0; (i < src->nr); i++)
807 dst->atom[i] = src->atom[i];
808 if (NULL != src->pdbinfo)
810 dst->pdbinfo[i] = src->pdbinfo[i];
812 if (NULL != src->atomname)
814 dst->atomname[i] = src->atomname[i];
816 if (NULL != src->atomtype)
818 dst->atomtype[i] = src->atomtype[i];
820 if (NULL != src->atomtypeB)
822 dst->atomtypeB[i] = src->atomtypeB[i];
825 dst->nres = src->nres;
826 for (i = 0; (i < src->nres); i++)
828 dst->resinfo[i] = src->resinfo[i];
833 void t_atoms_set_resinfo(t_atoms *atoms, int atom_ind, t_symtab *symtab,
834 const char *resname, int resnr, unsigned char ic,
835 int chainnum, char chainid)
839 ri = &atoms->resinfo[atoms->atom[atom_ind].resind];
840 ri->name = put_symtab(symtab, resname);
844 ri->chainnum = chainnum;
845 ri->chainid = chainid;
848 void free_t_atoms(t_atoms *atoms, gmx_bool bFreeNames)
854 for (i = 0; i < atoms->nr; i++)
856 sfree(*atoms->atomname[i]);
857 *atoms->atomname[i] = NULL;
859 for (i = 0; i < atoms->nres; i++)
861 sfree(*atoms->resinfo[i].name);
862 *atoms->resinfo[i].name = NULL;
865 sfree(atoms->atomname);
866 /* Do we need to free atomtype and atomtypeB as well ? */
867 sfree(atoms->resinfo);
871 sfree(atoms->pdbinfo);
875 atoms->atomname = NULL;
876 atoms->resinfo = NULL;
878 atoms->pdbinfo = NULL;
881 real max_cutoff(real cutoff1, real cutoff2)
883 if (cutoff1 == 0 || cutoff2 == 0)
889 return max(cutoff1, cutoff2);
893 extern void init_df_history(df_history_t *dfhist, int nlambda, real wl_delta)
898 dfhist->nlambda = nlambda;
899 dfhist->wl_delta = wl_delta;
900 snew(dfhist->sum_weights, dfhist->nlambda);
901 snew(dfhist->sum_dg, dfhist->nlambda);
902 snew(dfhist->sum_minvar, dfhist->nlambda);
903 snew(dfhist->sum_variance, dfhist->nlambda);
904 snew(dfhist->n_at_lam, dfhist->nlambda);
905 snew(dfhist->wl_histo, dfhist->nlambda);
907 /* allocate transition matrices here */
908 snew(dfhist->Tij, dfhist->nlambda);
909 snew(dfhist->Tij_empirical, dfhist->nlambda);
911 for (i = 0; i < dfhist->nlambda; i++)
913 snew(dfhist->Tij[i], dfhist->nlambda);
914 snew(dfhist->Tij_empirical[i], dfhist->nlambda);
917 snew(dfhist->accum_p, dfhist->nlambda);
918 snew(dfhist->accum_m, dfhist->nlambda);
919 snew(dfhist->accum_p2, dfhist->nlambda);
920 snew(dfhist->accum_m2, dfhist->nlambda);
922 for (i = 0; i < dfhist->nlambda; i++)
924 snew((dfhist->accum_p)[i], dfhist->nlambda);
925 snew((dfhist->accum_m)[i], dfhist->nlambda);
926 snew((dfhist->accum_p2)[i], dfhist->nlambda);
927 snew((dfhist->accum_m2)[i], dfhist->nlambda);
931 extern void copy_df_history(df_history_t *df_dest, df_history_t *df_source)
935 init_df_history(df_dest, df_source->nlambda, df_source->wl_delta);
936 df_dest->nlambda = df_source->nlambda;
937 df_dest->bEquil = df_source->bEquil;
938 for (i = 0; i < df_dest->nlambda; i++)
940 df_dest->sum_weights[i] = df_source->sum_weights[i];
941 df_dest->sum_dg[i] = df_source->sum_dg[i];
942 df_dest->sum_minvar[i] = df_source->sum_minvar[i];
943 df_dest->sum_variance[i] = df_source->sum_variance[i];
944 df_dest->n_at_lam[i] = df_source->n_at_lam[i];
945 df_dest->wl_histo[i] = df_source->wl_histo[i];
946 df_dest->accum_p[i] = df_source->accum_p[i];
947 df_dest->accum_m[i] = df_source->accum_m[i];
948 df_dest->accum_p2[i] = df_source->accum_p2[i];
949 df_dest->accum_m2[i] = df_source->accum_m2[i];
952 for (i = 0; i < df_dest->nlambda; i++)
954 for (j = 0; j < df_dest->nlambda; j++)
956 df_dest->Tij[i][j] = df_source->Tij[i][j];
957 df_dest->Tij_empirical[i][j] = df_source->Tij_empirical[i][j];