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
,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
;
114 md
->nr
= mtop
->natoms
;
119 if (md
->nr
> md
->nalloc
) {
120 md
->nalloc
= over_alloc_dd(md
->nr
);
122 if (md
->nMassPerturbed
) {
123 srenew(md
->massA
,md
->nalloc
);
124 srenew(md
->massB
,md
->nalloc
);
126 srenew(md
->massT
,md
->nalloc
);
127 srenew(md
->invmass
,md
->nalloc
);
128 srenew(md
->chargeA
,md
->nalloc
);
129 if (md
->nPerturbed
) {
130 srenew(md
->chargeB
,md
->nalloc
);
132 srenew(md
->typeA
,md
->nalloc
);
133 if (md
->nPerturbed
) {
134 srenew(md
->typeB
,md
->nalloc
);
136 srenew(md
->ptype
,md
->nalloc
);
137 if (opts
->ngtc
> 1) {
138 srenew(md
->cTC
,md
->nalloc
);
139 /* We always copy cTC with domain decomposition */
141 srenew(md
->cENER
,md
->nalloc
);
143 srenew(md
->cACC
,md
->nalloc
);
146 opts
->nFreeze
[0][XX
] || opts
->nFreeze
[0][YY
] || opts
->nFreeze
[0][ZZ
]))
147 srenew(md
->cFREEZE
,md
->nalloc
);
149 srenew(md
->cVCM
,md
->nalloc
);
151 srenew(md
->cORF
,md
->nalloc
);
153 srenew(md
->bPerturbed
,md
->nalloc
);
155 /* Note that these user t_mdatoms array pointers are NULL
156 * when there is only one group present.
157 * Therefore, when adding code, the user should use something like:
158 * gprnrU1 = (md->cU1==NULL ? 0 : md->cU1[localatindex])
160 if (mtop
->groups
.grpnr
[egcUser1
] != NULL
)
161 srenew(md
->cU1
,md
->nalloc
);
162 if (mtop
->groups
.grpnr
[egcUser2
] != NULL
)
163 srenew(md
->cU2
,md
->nalloc
);
166 srenew(md
->bQM
,md
->nalloc
);
169 for(i
=0; (i
<md
->nr
); i
++) {
172 gmx_mtop_atomnr_to_atom(mtop
,ag
,&atom
);
180 ae
= as
+ molblock
[molb
].nmol
*molblock
[molb
].natoms_mol
;
182 atoms_mol
= &mtop
->moltype
[molblock
[molb
].type
].atoms
;
183 atom
= &atoms_mol
->atom
[(ag
- as
) % atoms_mol
->nr
];
187 md
->cFREEZE
[i
] = ggrpnr(groups
,egcFREEZE
,ag
);
189 if (EI_ENERGY_MINIMIZATION(ir
->eI
)) {
192 } else if (ir
->eI
== eiBD
) {
193 /* Make the mass proportional to the friction coefficient for BD.
194 * This is necessary for the constraint algorithms.
197 mA
= ir
->bd_fric
*ir
->delta_t
;
198 mB
= ir
->bd_fric
*ir
->delta_t
;
200 fac
= ir
->delta_t
/opts
->tau_t
[md
->cTC
? groups
->grpnr
[egcTC
][ag
] : 0];
208 if (md
->nMassPerturbed
) {
215 } else if (md
->cFREEZE
) {
217 if (opts
->nFreeze
[g
][XX
] && opts
->nFreeze
[g
][YY
] && opts
->nFreeze
[g
][ZZ
])
218 /* Set the mass of completely frozen particles to ALMOST_ZERO iso 0
219 * to avoid div by zero in lincs or shake.
220 * Note that constraints can still move a partially frozen particle.
222 md
->invmass
[i
] = ALMOST_ZERO
;
224 md
->invmass
[i
] = 1.0/mA
;
226 md
->invmass
[i
] = 1.0/mA
;
228 md
->chargeA
[i
] = atom
->q
;
229 md
->typeA
[i
] = atom
->type
;
230 if (md
->nPerturbed
) {
231 md
->chargeB
[i
] = atom
->qB
;
232 md
->typeB
[i
] = atom
->typeB
;
233 md
->bPerturbed
[i
] = PERTURBED(*atom
);
235 md
->ptype
[i
] = atom
->ptype
;
237 md
->cTC
[i
] = groups
->grpnr
[egcTC
][ag
];
239 (groups
->grpnr
[egcENER
] ? groups
->grpnr
[egcENER
][ag
] : 0);
241 md
->cACC
[i
] = groups
->grpnr
[egcACC
][ag
];
243 md
->cVCM
[i
] = groups
->grpnr
[egcVCM
][ag
];
245 md
->cORF
[i
] = groups
->grpnr
[egcORFIT
][ag
];
248 md
->cU1
[i
] = groups
->grpnr
[egcUser1
][ag
];
250 md
->cU2
[i
] = groups
->grpnr
[egcUser2
][ag
];
253 if (groups
->grpnr
[egcQMMM
] == 0 ||
254 groups
->grpnr
[egcQMMM
][ag
] < groups
->grps
[egcQMMM
].nr
-1) {
267 void update_mdatoms(t_mdatoms
*md
,real lambda
)
274 if (md
->nMassPerturbed
) {
275 for(al
=0; (al
<end
); al
++) {
276 if (md
->bPerturbed
[al
]) {
277 md
->massT
[al
] = L1
*md
->massA
[al
]+ lambda
*md
->massB
[al
];
278 if (md
->invmass
[al
] > 1.1*ALMOST_ZERO
)
279 md
->invmass
[al
] = 1.0/md
->massT
[al
];
282 md
->tmass
= L1
*md
->tmassA
+ lambda
*md
->tmassB
;
284 md
->tmass
= md
->tmassA
;