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
45 #include "mtop_util.h"
47 #define ALMOST_ZERO 1e-30
49 t_mdatoms *init_mdatoms(FILE *fp,gmx_mtop_t *mtop,gmx_bool bFreeEnergy)
55 gmx_mtop_atomloop_all_t aloop;
60 md->nenergrp = mtop->groups.grps[egcENER].nr;
65 aloop = gmx_mtop_atomloop_all_init(mtop);
66 while(gmx_mtop_atomloop_all_next(aloop,&a,&atom)) {
67 if (ggrpnr(&mtop->groups,egcVCM,a) > 0)
70 if (bFreeEnergy && PERTURBED(*atom)) {
72 if (atom->mB != atom->m)
74 if (atom->qB != atom->q)
75 md->nChargePerturbed++;
85 if (bFreeEnergy && fp)
87 "There are %d atoms and %d charges for free energy perturbation\n",
88 md->nPerturbed,md->nChargePerturbed);
90 md->bOrires = gmx_mtop_ftype_count(mtop,F_ORIRES);
95 void atoms2md(gmx_mtop_t *mtop,t_inputrec *ir,
96 int nindex,int *index,
101 int i,g,ag,as,ae,molb;
105 gmx_groups_t *groups;
106 gmx_molblock_t *molblock;
110 groups = &mtop->groups;
112 molblock = mtop->molblock;
114 /* Index==NULL indicates particle decomposition,
115 * unless we have an empty DD node, so also check for homenr and start.
116 * This should be signaled properly with an extra parameter or nindex==-1.
118 if (index == NULL && (homenr > 0 || start > 0)) {
119 md->nr = mtop->natoms;
124 if (md->nr > md->nalloc) {
125 md->nalloc = over_alloc_dd(md->nr);
127 if (md->nMassPerturbed) {
128 srenew(md->massA,md->nalloc);
129 srenew(md->massB,md->nalloc);
131 srenew(md->massT,md->nalloc);
132 srenew(md->invmass,md->nalloc);
133 srenew(md->chargeA,md->nalloc);
134 if (md->nPerturbed) {
135 srenew(md->chargeB,md->nalloc);
137 srenew(md->typeA,md->nalloc);
138 if (md->nPerturbed) {
139 srenew(md->typeB,md->nalloc);
141 srenew(md->ptype,md->nalloc);
142 if (opts->ngtc > 1) {
143 srenew(md->cTC,md->nalloc);
144 /* We always copy cTC with domain decomposition */
146 srenew(md->cENER,md->nalloc);
148 srenew(md->cACC,md->nalloc);
151 opts->nFreeze[0][XX] || opts->nFreeze[0][YY] || opts->nFreeze[0][ZZ]))
152 srenew(md->cFREEZE,md->nalloc);
154 srenew(md->cVCM,md->nalloc);
156 srenew(md->cORF,md->nalloc);
158 srenew(md->bPerturbed,md->nalloc);
160 /* Note that these user t_mdatoms array pointers are NULL
161 * when there is only one group present.
162 * Therefore, when adding code, the user should use something like:
163 * gprnrU1 = (md->cU1==NULL ? 0 : md->cU1[localatindex])
165 if (mtop->groups.grpnr[egcUser1] != NULL)
166 srenew(md->cU1,md->nalloc);
167 if (mtop->groups.grpnr[egcUser2] != NULL)
168 srenew(md->cU2,md->nalloc);
171 srenew(md->bQM,md->nalloc);
173 srenew(md->wf,md->nalloc);
174 srenew(md->tf_table_index,md->nalloc);
180 for(i=0; (i<md->nr); i++) {
183 gmx_mtop_atomnr_to_atom(mtop,ag,&atom);
191 ae = as + molblock[molb].nmol*molblock[molb].natoms_mol;
193 atoms_mol = &mtop->moltype[molblock[molb].type].atoms;
194 atom = &atoms_mol->atom[(ag - as) % atoms_mol->nr];
198 md->cFREEZE[i] = ggrpnr(groups,egcFREEZE,ag);
200 if (EI_ENERGY_MINIMIZATION(ir->eI))
202 /* Displacement is proportional to F, masses used for constraints */
206 else if (ir->eI == eiBD)
208 /* With BD the physical masses are irrelevant.
209 * To keep the code simple we use most of the normal MD code path
210 * for BD. Thus for constraining the masses should be proportional
211 * to the friction coefficient. We set the absolute value such that
212 * m/2<(dx/dt)^2> = m/2*2kT/fric*dt = kT/2 => m=fric*dt/2
213 * Then if we set the (meaningless) velocity to v=dx/dt, we get the
214 * correct kinetic energy and temperature using the usual code path.
215 * Thus with BD v*dt will give the displacement and the reported
216 * temperature can signal bad integration (too large time step).
220 mA = 0.5*ir->bd_fric*ir->delta_t;
221 mB = 0.5*ir->bd_fric*ir->delta_t;
225 /* The friction coefficient is mass/tau_t */
226 fac = ir->delta_t/opts->tau_t[md->cTC ? groups->grpnr[egcTC][ag] : 0];
227 mA = 0.5*atom->m*fac;
228 mB = 0.5*atom->mB*fac;
236 if (md->nMassPerturbed) {
243 } else if (md->cFREEZE) {
245 if (opts->nFreeze[g][XX] && opts->nFreeze[g][YY] && opts->nFreeze[g][ZZ])
246 /* Set the mass of completely frozen particles to ALMOST_ZERO iso 0
247 * to avoid div by zero in lincs or shake.
248 * Note that constraints can still move a partially frozen particle.
250 md->invmass[i] = ALMOST_ZERO;
252 md->invmass[i] = 1.0/mA;
254 md->invmass[i] = 1.0/mA;
256 md->chargeA[i] = atom->q;
257 md->typeA[i] = atom->type;
258 if (md->nPerturbed) {
259 md->chargeB[i] = atom->qB;
260 md->typeB[i] = atom->typeB;
261 md->bPerturbed[i] = PERTURBED(*atom);
263 md->ptype[i] = atom->ptype;
265 md->cTC[i] = groups->grpnr[egcTC][ag];
267 (groups->grpnr[egcENER] ? groups->grpnr[egcENER][ag] : 0);
269 md->cACC[i] = groups->grpnr[egcACC][ag];
271 md->cVCM[i] = groups->grpnr[egcVCM][ag];
273 md->cORF[i] = groups->grpnr[egcORFIT][ag];
276 md->cU1[i] = groups->grpnr[egcUser1][ag];
278 md->cU2[i] = groups->grpnr[egcUser2][ag];
281 if (groups->grpnr[egcQMMM] == 0 ||
282 groups->grpnr[egcQMMM][ag] < groups->grps[egcQMMM].nr-1) {
288 /* Initialize AdResS weighting functions to adressw */
291 /* if no tf table groups specified, use default table */
292 md->tf_table_index[i] = DEFAULT_TF_TABLE;
293 if (ir->adress->n_tf_grps > 0){
294 /* if tf table groups specified, tf is only applied to thoose energy groups*/
295 md->tf_table_index[i] = NO_TF_TABLE;
296 /* check wether atom is in one of the relevant energy groups and assign a table index */
297 for (g=0; g<ir->adress->n_tf_grps; g++){
298 if (md->cENER[i] == ir->adress->tf_table_index[g]){
299 md->tf_table_index[i] = g;
311 void update_mdatoms(t_mdatoms *md,real lambda)
318 if (md->nMassPerturbed) {
319 for(al=0; (al<end); al++) {
320 if (md->bPerturbed[al]) {
321 md->massT[al] = L1*md->massA[al]+ lambda*md->massB[al];
322 if (md->invmass[al] > 1.1*ALMOST_ZERO)
323 md->invmass[al] = 1.0/md->massT[al];
326 md->tmass = L1*md->tmassA + lambda*md->tmassB;
328 md->tmass = md->tmassA;