e77586e1572e52ade5e422de61daeb6507deb48f
[alexxy/gromacs.git] / src / gromacs / mdtypes / state.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014,2015,2016, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 /* This file is completely threadsafe - keep it that way! */
38 #include "gmxpre.h"
39
40 #include "state.h"
41
42 #include <cstring>
43
44 #include <algorithm>
45
46 #include "gromacs/math/vec.h"
47 #include "gromacs/math/veccompare.h"
48 #include "gromacs/mdtypes/df_history.h"
49 #include "gromacs/mdtypes/inputrec.h"
50 #include "gromacs/mdtypes/md_enums.h"
51 #include "gromacs/utility/compare.h"
52 #include "gromacs/utility/gmxassert.h"
53 #include "gromacs/utility/smalloc.h"
54
55 /* The source code in this file should be thread-safe.
56       Please keep it that way. */
57
58 static void zero_history(history_t *hist)
59 {
60     hist->disre_initf  = 0;
61     hist->ndisrepairs  = 0;
62     hist->disre_rm3tav = NULL;
63     hist->orire_initf  = 0;
64     hist->norire_Dtav  = 0;
65     hist->orire_Dtav   = NULL;
66 }
67
68 static void zero_ekinstate(ekinstate_t *eks)
69 {
70     eks->ekin_n         = 0;
71     eks->ekinh          = NULL;
72     eks->ekinf          = NULL;
73     eks->ekinh_old      = NULL;
74     eks->ekinscalef_nhc.resize(0);
75     eks->ekinscaleh_nhc.resize(0);
76     eks->vscale_nhc.resize(0);
77     eks->dekindl        = 0;
78     eks->mvcos          = 0;
79 }
80
81 static void init_swapstate(swapstate_t *swapstate)
82 {
83     /* Ion/water position swapping */
84     swapstate->eSwapCoords            = 0;
85     swapstate->nIonTypes              = 0;
86     swapstate->nAverage               = 0;
87     swapstate->fluxleak               = 0;
88     swapstate->fluxleak_p             = NULL;
89     swapstate->bFromCpt               = 0;
90     swapstate->nat[eChan0]            = 0;
91     swapstate->nat[eChan1]            = 0;
92     swapstate->xc_old_whole[eChan0]   = NULL;
93     swapstate->xc_old_whole[eChan1]   = NULL;
94     swapstate->xc_old_whole_p[eChan0] = NULL;
95     swapstate->xc_old_whole_p[eChan1] = NULL;
96     swapstate->ionType                = NULL;
97 }
98
99 void init_gtc_state(t_state *state, int ngtc, int nnhpres, int nhchainlength)
100 {
101     state->ngtc          = ngtc;
102     state->nnhpres       = nnhpres;
103     state->nhchainlength = nhchainlength;
104     state->nosehoover_xi.resize(state->nhchainlength*state->ngtc, 0);
105     state->nosehoover_vxi.resize(state->nhchainlength*state->ngtc, 0);
106     state->therm_integral.resize(state->ngtc, 0);
107     state->nhpres_xi.resize(state->nhchainlength*nnhpres, 0);
108     state->nhpres_vxi.resize(state->nhchainlength*nnhpres, 0);
109 }
110
111
112 void init_state(t_state *state, int natoms, int ngtc, int nnhpres, int nhchainlength, int dfhistNumLambda)
113 {
114     state->natoms    = natoms;
115     state->flags     = 0;
116     state->fep_state = 0;
117     state->lambda.resize(efptNR, 0);
118     state->veta   = 0;
119     clear_mat(state->box);
120     clear_mat(state->box_rel);
121     clear_mat(state->boxv);
122     clear_mat(state->pres_prev);
123     clear_mat(state->svir_prev);
124     clear_mat(state->fvir_prev);
125     init_gtc_state(state, ngtc, nnhpres, nhchainlength);
126     if (state->natoms > 0)
127     {
128         /* We need to allocate one element extra, since we might use
129          * (unaligned) 4-wide SIMD loads to access rvec entries.
130          */
131         state->x.resize(state->natoms + 1);
132         state->v.resize(state->natoms + 1);
133     }
134     else
135     {
136         state->x.resize(0);
137         state->v.resize(0);
138     }
139     state->cg_p.resize(0);
140     zero_history(&state->hist);
141     zero_ekinstate(&state->ekinstate);
142     if (dfhistNumLambda > 0)
143     {
144         snew(state->dfhist, 1);
145         init_df_history(state->dfhist, dfhistNumLambda);
146     }
147     else
148     {
149         state->dfhist = NULL;
150     }
151     state->swapstate       = NULL;
152     state->edsamstate      = NULL;
153     state->ddp_count       = 0;
154     state->ddp_count_cg_gl = 0;
155     state->cg_gl.resize(0);
156 }
157
158 void comp_state(const t_state *st1, const t_state *st2,
159                 gmx_bool bRMSD, real ftol, real abstol)
160 {
161     int i, j, nc;
162
163     fprintf(stdout, "comparing flags\n");
164     cmp_int(stdout, "flags", -1, st1->flags, st2->flags);
165     fprintf(stdout, "comparing box\n");
166     cmp_rvecs(stdout, "box", DIM, st1->box, st2->box, FALSE, ftol, abstol);
167     fprintf(stdout, "comparing box_rel\n");
168     cmp_rvecs(stdout, "box_rel", DIM, st1->box_rel, st2->box_rel, FALSE, ftol, abstol);
169     fprintf(stdout, "comparing boxv\n");
170     cmp_rvecs(stdout, "boxv", DIM, st1->boxv, st2->boxv, FALSE, ftol, abstol);
171     if (st1->flags & (1<<estSVIR_PREV))
172     {
173         fprintf(stdout, "comparing shake vir_prev\n");
174         cmp_rvecs(stdout, "svir_prev", DIM, st1->svir_prev, st2->svir_prev, FALSE, ftol, abstol);
175     }
176     if (st1->flags & (1<<estFVIR_PREV))
177     {
178         fprintf(stdout, "comparing force vir_prev\n");
179         cmp_rvecs(stdout, "fvir_prev", DIM, st1->fvir_prev, st2->fvir_prev, FALSE, ftol, abstol);
180     }
181     if (st1->flags & (1<<estPRES_PREV))
182     {
183         fprintf(stdout, "comparing prev_pres\n");
184         cmp_rvecs(stdout, "pres_prev", DIM, st1->pres_prev, st2->pres_prev, FALSE, ftol, abstol);
185     }
186     cmp_int(stdout, "ngtc", -1, st1->ngtc, st2->ngtc);
187     cmp_int(stdout, "nhchainlength", -1, st1->nhchainlength, st2->nhchainlength);
188     if (st1->ngtc == st2->ngtc && st1->nhchainlength == st2->nhchainlength)
189     {
190         for (i = 0; i < st1->ngtc; i++)
191         {
192             nc = i*st1->nhchainlength;
193             for (j = 0; j < nc; j++)
194             {
195                 cmp_real(stdout, "nosehoover_xi",
196                          i, st1->nosehoover_xi[nc+j], st2->nosehoover_xi[nc+j], ftol, abstol);
197             }
198         }
199     }
200     cmp_int(stdout, "nnhpres", -1, st1->nnhpres, st2->nnhpres);
201     if (st1->nnhpres == st2->nnhpres && st1->nhchainlength == st2->nhchainlength)
202     {
203         for (i = 0; i < st1->nnhpres; i++)
204         {
205             nc = i*st1->nhchainlength;
206             for (j = 0; j < nc; j++)
207             {
208                 cmp_real(stdout, "nosehoover_xi",
209                          i, st1->nhpres_xi[nc+j], st2->nhpres_xi[nc+j], ftol, abstol);
210             }
211         }
212     }
213
214     cmp_int(stdout, "natoms", -1, st1->natoms, st2->natoms);
215     if (st1->natoms == st2->natoms)
216     {
217         if ((st1->flags & (1<<estX)) && (st2->flags & (1<<estX)))
218         {
219             fprintf(stdout, "comparing x\n");
220             cmp_rvecs(stdout, "x", st1->natoms, as_rvec_array(st1->x.data()), as_rvec_array(st2->x.data()), bRMSD, ftol, abstol);
221         }
222         if ((st1->flags & (1<<estV)) && (st2->flags & (1<<estV)))
223         {
224             fprintf(stdout, "comparing v\n");
225             cmp_rvecs(stdout, "v", st1->natoms, as_rvec_array(st1->v.data()), as_rvec_array(st2->v.data()), bRMSD, ftol, abstol);
226         }
227     }
228 }
229
230 rvec *getRvecArrayFromPaddedRVecVector(const PaddedRVecVector *v,
231                                        unsigned int            n)
232 {
233     GMX_ASSERT(v->size() >= n, "We can't copy more elements than the vector size");
234
235     rvec *dest;
236
237     snew(dest, n);
238
239     const rvec *vPtr = as_rvec_array(v->data());
240     for (unsigned int i = 0; i < n; i++)
241     {
242         copy_rvec(vPtr[i], dest[i]);
243     }
244
245     return dest;
246 }