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"
46 #include "gmx_omp_nthreads.h"
48 #define ALMOST_ZERO 1e-30
50 t_mdatoms *init_mdatoms(FILE *fp,gmx_mtop_t *mtop,gmx_bool bFreeEnergy)
56 gmx_mtop_atomloop_all_t aloop;
61 md->nenergrp = mtop->groups.grps[egcENER].nr;
66 aloop = gmx_mtop_atomloop_all_init(mtop);
67 while(gmx_mtop_atomloop_all_next(aloop,&a,&atom)) {
68 if (ggrpnr(&mtop->groups,egcVCM,a) > 0)
71 if (bFreeEnergy && PERTURBED(*atom)) {
73 if (atom->mB != atom->m)
75 if (atom->qB != atom->q)
76 md->nChargePerturbed++;
86 if (bFreeEnergy && fp)
88 "There are %d atoms and %d charges for free energy perturbation\n",
89 md->nPerturbed,md->nChargePerturbed);
91 md->bOrires = gmx_mtop_ftype_count(mtop,F_ORIRES);
96 void atoms2md(gmx_mtop_t *mtop,t_inputrec *ir,
97 int nindex,int *index,
101 gmx_mtop_atomlookup_t alook;
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);
172 srenew(md->wf,md->nalloc);
173 srenew(md->tf_table_index,md->nalloc);
177 alook = gmx_mtop_atomlookup_init(mtop);
179 #pragma omp parallel for num_threads(gmx_omp_nthreads_get(emntDefault)) schedule(static)
180 for(i=0; i<md->nr; i++) {
190 gmx_mtop_atomnr_to_atom(alook,ag,&atom);
193 md->cFREEZE[i] = ggrpnr(groups,egcFREEZE,ag);
195 if (EI_ENERGY_MINIMIZATION(ir->eI))
197 /* Displacement is proportional to F, masses used for constraints */
201 else if (ir->eI == eiBD)
203 /* With BD the physical masses are irrelevant.
204 * To keep the code simple we use most of the normal MD code path
205 * for BD. Thus for constraining the masses should be proportional
206 * to the friction coefficient. We set the absolute value such that
207 * m/2<(dx/dt)^2> = m/2*2kT/fric*dt = kT/2 => m=fric*dt/2
208 * Then if we set the (meaningless) velocity to v=dx/dt, we get the
209 * correct kinetic energy and temperature using the usual code path.
210 * Thus with BD v*dt will give the displacement and the reported
211 * temperature can signal bad integration (too large time step).
215 mA = 0.5*ir->bd_fric*ir->delta_t;
216 mB = 0.5*ir->bd_fric*ir->delta_t;
220 /* The friction coefficient is mass/tau_t */
221 fac = ir->delta_t/opts->tau_t[md->cTC ? groups->grpnr[egcTC][ag] : 0];
222 mA = 0.5*atom->m*fac;
223 mB = 0.5*atom->mB*fac;
231 if (md->nMassPerturbed) {
238 } else if (md->cFREEZE) {
240 if (opts->nFreeze[g][XX] && opts->nFreeze[g][YY] && opts->nFreeze[g][ZZ])
241 /* Set the mass of completely frozen particles to ALMOST_ZERO iso 0
242 * to avoid div by zero in lincs or shake.
243 * Note that constraints can still move a partially frozen particle.
245 md->invmass[i] = ALMOST_ZERO;
247 md->invmass[i] = 1.0/mA;
249 md->invmass[i] = 1.0/mA;
251 md->chargeA[i] = atom->q;
252 md->typeA[i] = atom->type;
253 if (md->nPerturbed) {
254 md->chargeB[i] = atom->qB;
255 md->typeB[i] = atom->typeB;
256 md->bPerturbed[i] = PERTURBED(*atom);
258 md->ptype[i] = atom->ptype;
260 md->cTC[i] = groups->grpnr[egcTC][ag];
262 (groups->grpnr[egcENER] ? groups->grpnr[egcENER][ag] : 0);
264 md->cACC[i] = groups->grpnr[egcACC][ag];
266 md->cVCM[i] = groups->grpnr[egcVCM][ag];
268 md->cORF[i] = groups->grpnr[egcORFIT][ag];
271 md->cU1[i] = groups->grpnr[egcUser1][ag];
273 md->cU2[i] = groups->grpnr[egcUser2][ag];
276 if (groups->grpnr[egcQMMM] == 0 ||
277 groups->grpnr[egcQMMM][ag] < groups->grps[egcQMMM].nr-1) {
283 /* Initialize AdResS weighting functions to adressw */
286 /* if no tf table groups specified, use default table */
287 md->tf_table_index[i] = DEFAULT_TF_TABLE;
288 if (ir->adress->n_tf_grps > 0){
289 /* if tf table groups specified, tf is only applied to thoose energy groups*/
290 md->tf_table_index[i] = NO_TF_TABLE;
291 /* check wether atom is in one of the relevant energy groups and assign a table index */
292 for (g=0; g<ir->adress->n_tf_grps; g++){
293 if (md->cENER[i] == ir->adress->tf_table_index[g]){
294 md->tf_table_index[i] = g;
301 gmx_mtop_atomlookup_destroy(alook);
308 void update_mdatoms(t_mdatoms *md,real lambda)
315 if (md->nMassPerturbed) {
316 for(al=0; (al<end); al++) {
317 if (md->bPerturbed[al]) {
318 md->massT[al] = L1*md->massA[al]+ lambda*md->massB[al];
319 if (md->invmass[al] > 1.1*ALMOST_ZERO)
320 md->invmass[al] = 1.0/md->massT[al];
323 md->tmass = L1*md->tmassA + lambda*md->tmassB;
325 md->tmass = md->tmassA;