Allow using COM of previous step as PBC reference
[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,2017,2018, 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/paddedvector.h"
47 #include "gromacs/math/vec.h"
48 #include "gromacs/math/veccompare.h"
49 #include "gromacs/mdtypes/awh-history.h"
50 #include "gromacs/mdtypes/df_history.h"
51 #include "gromacs/mdtypes/inputrec.h"
52 #include "gromacs/mdtypes/md_enums.h"
53 #include "gromacs/mdtypes/pull-params.h"
54 #include "gromacs/mdtypes/swaphistory.h"
55 #include "gromacs/pbcutil/boxutilities.h"
56 #include "gromacs/pbcutil/pbc.h"
57 #include "gromacs/utility/compare.h"
58 #include "gromacs/utility/gmxassert.h"
59 #include "gromacs/utility/smalloc.h"
60
61 /* The source code in this file should be thread-safe.
62       Please keep it that way. */
63
64 history_t::history_t() : disre_initf(0),
65                          ndisrepairs(0),
66                          disre_rm3tav(nullptr),
67                          orire_initf(0),
68                          norire_Dtav(0),
69                          orire_Dtav(nullptr)
70 {
71 };
72
73 ekinstate_t::ekinstate_t() : bUpToDate(FALSE),
74                              ekin_n(0),
75                              ekinh(nullptr),
76                              ekinf(nullptr),
77                              ekinh_old(nullptr),
78                              ekin_total(),
79
80                              dekindl(0),
81                              mvcos(0)
82 {
83     clear_mat(ekin_total);
84 };
85
86 void init_gtc_state(t_state *state, int ngtc, int nnhpres, int nhchainlength)
87 {
88     state->ngtc          = ngtc;
89     state->nnhpres       = nnhpres;
90     state->nhchainlength = nhchainlength;
91     state->nosehoover_xi.resize(state->nhchainlength*state->ngtc, 0);
92     state->nosehoover_vxi.resize(state->nhchainlength*state->ngtc, 0);
93     state->therm_integral.resize(state->ngtc, 0);
94     state->baros_integral = 0.0;
95     state->nhpres_xi.resize(state->nhchainlength*nnhpres, 0);
96     state->nhpres_vxi.resize(state->nhchainlength*nnhpres, 0);
97 }
98
99
100 /* Checkpoint code relies on this function having no effect if
101    state->natoms is > 0 and passed as natoms. */
102 void state_change_natoms(t_state *state, int natoms)
103 {
104     state->natoms = natoms;
105
106     /* We need padding, since we might use SIMD access */
107     const size_t paddedSize = gmx::paddedRVecVectorSize(state->natoms);
108
109     if (state->flags & (1 << estX))
110     {
111         state->x.resize(paddedSize);
112     }
113     if (state->flags & (1 << estV))
114     {
115         state->v.resize(paddedSize);
116     }
117     if (state->flags & (1 << estCGP))
118     {
119         state->cg_p.resize(paddedSize);
120     }
121 }
122
123 void init_dfhist_state(t_state *state, int dfhistNumLambda)
124 {
125     if (dfhistNumLambda > 0)
126     {
127         snew(state->dfhist, 1);
128         init_df_history(state->dfhist, dfhistNumLambda);
129     }
130     else
131     {
132         state->dfhist = nullptr;
133     }
134 }
135
136 void comp_state(const t_state *st1, const t_state *st2,
137                 gmx_bool bRMSD, real ftol, real abstol)
138 {
139     int i, j, nc;
140
141     fprintf(stdout, "comparing flags\n");
142     cmp_int(stdout, "flags", -1, st1->flags, st2->flags);
143     fprintf(stdout, "comparing box\n");
144     cmp_rvecs(stdout, "box", DIM, st1->box, st2->box, FALSE, ftol, abstol);
145     fprintf(stdout, "comparing box_rel\n");
146     cmp_rvecs(stdout, "box_rel", DIM, st1->box_rel, st2->box_rel, FALSE, ftol, abstol);
147     fprintf(stdout, "comparing boxv\n");
148     cmp_rvecs(stdout, "boxv", DIM, st1->boxv, st2->boxv, FALSE, ftol, abstol);
149     if (st1->flags & (1<<estSVIR_PREV))
150     {
151         fprintf(stdout, "comparing shake vir_prev\n");
152         cmp_rvecs(stdout, "svir_prev", DIM, st1->svir_prev, st2->svir_prev, FALSE, ftol, abstol);
153     }
154     if (st1->flags & (1<<estFVIR_PREV))
155     {
156         fprintf(stdout, "comparing force vir_prev\n");
157         cmp_rvecs(stdout, "fvir_prev", DIM, st1->fvir_prev, st2->fvir_prev, FALSE, ftol, abstol);
158     }
159     if (st1->flags & (1<<estPRES_PREV))
160     {
161         fprintf(stdout, "comparing prev_pres\n");
162         cmp_rvecs(stdout, "pres_prev", DIM, st1->pres_prev, st2->pres_prev, FALSE, ftol, abstol);
163     }
164     cmp_int(stdout, "ngtc", -1, st1->ngtc, st2->ngtc);
165     cmp_int(stdout, "nhchainlength", -1, st1->nhchainlength, st2->nhchainlength);
166     if (st1->ngtc == st2->ngtc && st1->nhchainlength == st2->nhchainlength)
167     {
168         for (i = 0; i < st1->ngtc; i++)
169         {
170             nc = i*st1->nhchainlength;
171             for (j = 0; j < nc; j++)
172             {
173                 cmp_real(stdout, "nosehoover_xi",
174                          i, st1->nosehoover_xi[nc+j], st2->nosehoover_xi[nc+j], ftol, abstol);
175             }
176         }
177     }
178     cmp_int(stdout, "nnhpres", -1, st1->nnhpres, st2->nnhpres);
179     if (st1->nnhpres == st2->nnhpres && st1->nhchainlength == st2->nhchainlength)
180     {
181         for (i = 0; i < st1->nnhpres; i++)
182         {
183             nc = i*st1->nhchainlength;
184             for (j = 0; j < nc; j++)
185             {
186                 cmp_real(stdout, "nosehoover_xi",
187                          i, st1->nhpres_xi[nc+j], st2->nhpres_xi[nc+j], ftol, abstol);
188             }
189         }
190     }
191
192     cmp_int(stdout, "natoms", -1, st1->natoms, st2->natoms);
193     if (st1->natoms == st2->natoms)
194     {
195         if ((st1->flags & (1<<estX)) && (st2->flags & (1<<estX)))
196         {
197             fprintf(stdout, "comparing x\n");
198             cmp_rvecs(stdout, "x", st1->natoms, as_rvec_array(st1->x.data()), as_rvec_array(st2->x.data()), bRMSD, ftol, abstol);
199         }
200         if ((st1->flags & (1<<estV)) && (st2->flags & (1<<estV)))
201         {
202             fprintf(stdout, "comparing v\n");
203             cmp_rvecs(stdout, "v", st1->natoms, as_rvec_array(st1->v.data()), as_rvec_array(st2->v.data()), bRMSD, ftol, abstol);
204         }
205     }
206 }
207
208 rvec *makeRvecArray(gmx::ArrayRef<const gmx::RVec> v,
209                     gmx::index                     n)
210 {
211     GMX_ASSERT(v.size() >= n, "We can't copy more elements than the vector size");
212
213     rvec *dest;
214
215     snew(dest, n);
216
217     const rvec *vPtr = as_rvec_array(v.data());
218     for (gmx::index i = 0; i < n; i++)
219     {
220         copy_rvec(vPtr[i], dest[i]);
221     }
222
223     return dest;
224 }
225
226 t_state::t_state() : natoms(0),
227                      ngtc(0),
228                      nnhpres(0),
229                      nhchainlength(0),
230                      flags(0),
231                      fep_state(0),
232                      lambda(),
233
234                      baros_integral(0),
235                      veta(0),
236                      vol0(0),
237
238                      ekinstate(),
239                      hist(),
240                      dfhist(nullptr),
241                      awhHistory(nullptr),
242                      ddp_count(0),
243                      ddp_count_cg_gl(0)
244
245 {
246     // It would be nicer to initialize these with {} or {{0}} in the
247     // above initialization list, but uncrustify doesn't understand
248     // that.
249     // TODO Fix this if we switch to clang-format some time.
250     lambda = {{ 0 }};
251     clear_mat(box);
252     clear_mat(box_rel);
253     clear_mat(boxv);
254     clear_mat(pres_prev);
255     clear_mat(svir_prev);
256     clear_mat(fvir_prev);
257 }
258
259 void set_box_rel(const t_inputrec *ir, t_state *state)
260 {
261     /* Make sure the box obeys the restrictions before we fix the ratios */
262     correct_box(nullptr, 0, state->box, nullptr);
263
264     clear_mat(state->box_rel);
265
266     if (inputrecPreserveShape(ir))
267     {
268         const int ndim = ir->epct == epctSEMIISOTROPIC ? 2 : 3;
269         do_box_rel(ndim, ir->deform, state->box_rel, state->box, true);
270     }
271 }
272
273 void preserve_box_shape(const t_inputrec *ir, matrix box_rel, matrix box)
274 {
275     if (inputrecPreserveShape(ir))
276     {
277         const int ndim = ir->epct == epctSEMIISOTROPIC ? 2 : 3;
278         do_box_rel(ndim, ir->deform, box_rel, box, false);
279     }
280 }