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 * GROwing Monsters And Cloning Shrimps
44 #include "mtop_util.h"
46 #define ALMOST_ZERO 1e-30
48 t_mdatoms *init_mdatoms(FILE *fp,gmx_mtop_t *mtop,gmx_bool bFreeEnergy)
54 gmx_mtop_atomloop_all_t aloop;
59 md->nenergrp = mtop->groups.grps[egcENER].nr;
64 aloop = gmx_mtop_atomloop_all_init(mtop);
65 while(gmx_mtop_atomloop_all_next(aloop,&a,&atom)) {
66 if (ggrpnr(&mtop->groups,egcVCM,a) > 0)
69 if (bFreeEnergy && PERTURBED(*atom)) {
71 if (atom->mB != atom->m)
73 if (atom->qB != atom->q)
74 md->nChargePerturbed++;
84 if (bFreeEnergy && fp)
86 "There are %d atoms and %d charges for free energy perturbation\n",
87 md->nPerturbed,md->nChargePerturbed);
89 md->bOrires = gmx_mtop_ftype_count(mtop,F_ORIRES);
94 void atoms2md(gmx_mtop_t *mtop,t_inputrec *ir,
95 int nindex,int *index,
100 int i,g,ag,as,ae,molb;
104 gmx_groups_t *groups;
105 gmx_molblock_t *molblock;
109 groups = &mtop->groups;
111 molblock = mtop->molblock;
113 /* Index==NULL indicates particle decomposition,
114 * unless we have an empty DD node, so also check for homenr and start.
115 * This should be signaled properly with an extra parameter or nindex==-1.
117 if (index == NULL && (homenr > 0 || start > 0)) {
118 md->nr = mtop->natoms;
123 if (md->nr > md->nalloc) {
124 md->nalloc = over_alloc_dd(md->nr);
126 if (md->nMassPerturbed) {
127 srenew(md->massA,md->nalloc);
128 srenew(md->massB,md->nalloc);
130 srenew(md->massT,md->nalloc);
131 srenew(md->invmass,md->nalloc);
132 srenew(md->chargeA,md->nalloc);
133 if (md->nPerturbed) {
134 srenew(md->chargeB,md->nalloc);
136 srenew(md->typeA,md->nalloc);
137 if (md->nPerturbed) {
138 srenew(md->typeB,md->nalloc);
140 srenew(md->ptype,md->nalloc);
141 if (opts->ngtc > 1) {
142 srenew(md->cTC,md->nalloc);
143 /* We always copy cTC with domain decomposition */
145 srenew(md->cENER,md->nalloc);
147 srenew(md->cACC,md->nalloc);
150 opts->nFreeze[0][XX] || opts->nFreeze[0][YY] || opts->nFreeze[0][ZZ]))
151 srenew(md->cFREEZE,md->nalloc);
153 srenew(md->cVCM,md->nalloc);
155 srenew(md->cORF,md->nalloc);
157 srenew(md->bPerturbed,md->nalloc);
159 /* Note that these user t_mdatoms array pointers are NULL
160 * when there is only one group present.
161 * Therefore, when adding code, the user should use something like:
162 * gprnrU1 = (md->cU1==NULL ? 0 : md->cU1[localatindex])
164 if (mtop->groups.grpnr[egcUser1] != NULL)
165 srenew(md->cU1,md->nalloc);
166 if (mtop->groups.grpnr[egcUser2] != NULL)
167 srenew(md->cU2,md->nalloc);
170 srenew(md->bQM,md->nalloc);
173 for(i=0; (i<md->nr); i++) {
176 gmx_mtop_atomnr_to_atom(mtop,ag,&atom);
184 ae = as + molblock[molb].nmol*molblock[molb].natoms_mol;
186 atoms_mol = &mtop->moltype[molblock[molb].type].atoms;
187 atom = &atoms_mol->atom[(ag - as) % atoms_mol->nr];
191 md->cFREEZE[i] = ggrpnr(groups,egcFREEZE,ag);
193 if (EI_ENERGY_MINIMIZATION(ir->eI)) {
196 } else if (ir->eI == eiBD) {
197 /* Make the mass proportional to the friction coefficient for BD.
198 * This is necessary for the constraint algorithms.
201 mA = ir->bd_fric*ir->delta_t;
202 mB = ir->bd_fric*ir->delta_t;
204 fac = ir->delta_t/opts->tau_t[md->cTC ? groups->grpnr[egcTC][ag] : 0];
212 if (md->nMassPerturbed) {
219 } else if (md->cFREEZE) {
221 if (opts->nFreeze[g][XX] && opts->nFreeze[g][YY] && opts->nFreeze[g][ZZ])
222 /* Set the mass of completely frozen particles to ALMOST_ZERO iso 0
223 * to avoid div by zero in lincs or shake.
224 * Note that constraints can still move a partially frozen particle.
226 md->invmass[i] = ALMOST_ZERO;
228 md->invmass[i] = 1.0/mA;
230 md->invmass[i] = 1.0/mA;
232 md->chargeA[i] = atom->q;
233 md->typeA[i] = atom->type;
234 if (md->nPerturbed) {
235 md->chargeB[i] = atom->qB;
236 md->typeB[i] = atom->typeB;
237 md->bPerturbed[i] = PERTURBED(*atom);
239 md->ptype[i] = atom->ptype;
241 md->cTC[i] = groups->grpnr[egcTC][ag];
243 (groups->grpnr[egcENER] ? groups->grpnr[egcENER][ag] : 0);
245 md->cACC[i] = groups->grpnr[egcACC][ag];
247 md->cVCM[i] = groups->grpnr[egcVCM][ag];
249 md->cORF[i] = groups->grpnr[egcORFIT][ag];
252 md->cU1[i] = groups->grpnr[egcUser1][ag];
254 md->cU2[i] = groups->grpnr[egcUser2][ag];
257 if (groups->grpnr[egcQMMM] == 0 ||
258 groups->grpnr[egcQMMM][ag] < groups->grps[egcQMMM].nr-1) {
271 void update_mdatoms(t_mdatoms *md,real lambda)
278 if (md->nMassPerturbed) {
279 for(al=0; (al<end); al++) {
280 if (md->bPerturbed[al]) {
281 md->massT[al] = L1*md->massA[al]+ lambda*md->massB[al];
282 if (md->invmass[al] > 1.1*ALMOST_ZERO)
283 md->invmass[al] = 1.0/md->massT[al];
286 md->tmass = L1*md->tmassA + lambda*md->tmassB;
288 md->tmass = md->tmassA;