0dc0c23d789e30411e9aec8651015d93aca22e60
[alexxy/gromacs.git] / src / gromacs / mdtypes / fcdata.h
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) 2012,2014,2015,2016,2017 by the GROMACS development team.
7  * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #ifndef GMX_MDTYPES_FCDATA_H
39 #define GMX_MDTYPES_FCDATA_H
40
41 #include <memory>
42 #include <vector>
43
44 #include "gromacs/math/vectypes.h"
45 #include "gromacs/topology/idef.h"
46 #include "gromacs/utility/arrayref.h"
47 #include "gromacs/utility/classhelpers.h"
48 #include "gromacs/utility/real.h"
49
50 enum class DistanceRestraintWeighting : int;
51 struct gmx_mtop_t;
52 struct gmx_multisim_t;
53 struct t_commrec;
54 struct t_inputrec;
55 class t_state;
56
57 typedef real rvec5[5];
58
59 /* Distance restraining stuff */
60 typedef struct t_disresdata
61 {
62     DistanceRestraintWeighting dr_weighting; /* Weighting of pairs in one restraint              */
63     bool                       dr_bMixed;    /* Use sqrt of the instantaneous times              *
64                                               * the time averaged violation                      */
65     real dr_fc;                              /* Force constant for disres,                       *
66                                               * which is multiplied by a (possibly)              *
67                                               * different factor for each restraint              */
68     real  dr_tau;                            /* Time constant for disres                          */
69     real  ETerm;                             /* multiplication factor for time averaging         */
70     real  ETerm1;                            /* 1 - ETerm1                                       */
71     real  exp_min_t_tau;                     /* Factor for slowly switching on the force         */
72     int   nres;                              /* The number of distance restraints                */
73     int   npair;                             /* The number of distance restraint pairs           */
74     int   type_min;                          /* The minimum iparam type index for restraints     */
75     real  sumviol;                           /* The sum of violations                            */
76     real* rt;                                /* The instantaneous distance (npair)               */
77     real* rm3tav;                            /* The time averaged distance (npair)               */
78     real* Rtl_6;                             /* The instantaneous r^-6 (nres)                    */
79     real* Rt_6;                              /* The instantaneous ensemble averaged r^-6 (nres)  */
80     real* Rtav_6;                            /* The time and ensemble averaged r^-6 (nres)       */
81     int   nsystems;                          /* The number of systems for ensemble averaging     */
82
83     /* TODO: Implement a proper solution for parallel disre indexing */
84     const t_iatom* forceatomsStart; /* Pointer to the start of the disre forceatoms */
85 } t_disresdata;
86
87 /* All coefficients for the matrix equation for the orientation tensor */
88 struct OriresMatEq
89 {
90     real rhs[5];    /* The right hand side of the matrix equation */
91     real mat[5][5]; /* The matrix                                 */
92 };
93
94 //! \brief Orientation restraining stuff
95 struct t_oriresdata
96 {
97     /*! \brief Constructor
98      *
99      * \param[in] fplog  Log file, can be nullptr
100      * \param[in] mtop   The global topology
101      * \param[in] ir     The input record
102      * \param[in] cr     The commrec, can be nullptr when not running in parallel
103      * \param[in] ms     The multisim communicator, pass nullptr to avoid ensemble averaging
104      * \param[in,out] globalState  The global state, orientation restraint entires are added
105      *
106      * \throws InvalidInputError when there is domain decomposition, fewer than 5 restraints,
107      *         periodic molecules or more than 1 molecule for a moleculetype with restraints.
108      */
109     t_oriresdata(FILE*                 fplog,
110                  const gmx_mtop_t&     mtop,
111                  const t_inputrec&     ir,
112                  const t_commrec*      cr,
113                  const gmx_multisim_t* ms,
114                  t_state*              globalState);
115
116     //! Destructor
117     ~t_oriresdata();
118
119     //! Force constant for the restraints
120     real fc;
121     //! Multiplication factor for time averaging
122     real edt;
123     //! 1 - edt
124     real edt_1;
125     //! Factor for slowly switching on the force
126     real exp_min_t_tau;
127     //! The number of orientation restraints
128     const int numRestraints;
129     //! The number of experiments
130     int numExperiments;
131     //! The minimum iparam type index for restraints
132     int typeMin;
133     //! The number of atoms for the fit
134     int numReferenceAtoms;
135     //! The masses of the reference atoms
136     std::vector<real> mref;
137     //! The reference coordinates for the fit
138     std::vector<gmx::RVec> xref;
139     //! Temporary array for fitting
140     std::vector<gmx::RVec> xtmp;
141     //! Rotation matrix to rotate to the reference coordinates
142     matrix rotationMatrix;
143     //! Array of order tensors, one for each experiment
144     tensor* orderTensors = nullptr;
145     //! The order tensor D for all restraints
146     rvec5* DTensors = nullptr;
147     //! The ensemble averaged D for all restraints
148     rvec5* DTensorsEnsembleAv = nullptr;
149     //! The time and ensemble averaged D restraints
150     rvec5* DTensorsTimeAndEnsembleAv = nullptr;
151     //! The calculated instantaneous orientations
152     std::vector<real> orientations;
153     //! The calculated emsemble averaged orientations
154     gmx::ArrayRef<real> orientationsEnsembleAv;
155     //! Buffer for the calculated emsemble averaged orientations, only used with ensemble averaging
156     std::vector<real> orientationsEnsembleAvBuffer;
157     //! The calculated time and ensemble averaged orientations
158     gmx::ArrayRef<real> orientationsTimeAndEnsembleAv;
159     //! The weighted (using kfac) RMS deviation
160     std::vector<real> orientationsTimeAndEnsembleAvBuffer;
161     //! Buffer for the weighted (using kfac) RMS deviation, only used with time averaging
162     real rmsdev;
163     //! An temporary array of matrix + rhs
164     std::vector<OriresMatEq> tmpEq;
165     //! The number of eigenvalues + eigenvectors per experiment
166     static constexpr int c_numEigenRealsPerExperiment = 12;
167     //! Eigenvalues/vectors, for output only (numExperiments x 12)
168     std::vector<real> eigenOutput;
169
170     // variables for diagonalization with diagonalize_orires_tensors()
171     //! Tensor to diagonalize
172     std::array<gmx::DVec, DIM> M;
173     //! Eigenvalues
174     std::array<double, DIM> eig_diag;
175     //! Eigenvectors
176     std::array<gmx::DVec, DIM> v;
177
178     // Default copy and assign would be incorrect and manual versions are not yet implemented.
179     GMX_DISALLOW_COPY_AND_ASSIGN(t_oriresdata);
180 };
181
182 /* Cubic spline table for tabulated bonded interactions */
183 struct bondedtable_t
184 {
185     int               n;     /* n+1 is the number of points */
186     real              scale; /* distance between two points */
187     std::vector<real> data;  /* the actual table data, per point there are 4 numbers */
188 };
189
190 /*
191  * Data struct used in the force calculation routines
192  * for storing the tables for bonded interactions and
193  * for storing information which is needed in following steps
194  * (for instance for time averaging in distance retraints)
195  * or for storing output, since force routines only return the potential.
196  */
197 struct t_fcdata
198 {
199     std::vector<bondedtable_t> bondtab;
200     std::vector<bondedtable_t> angletab;
201     std::vector<bondedtable_t> dihtab;
202
203     // TODO: Convert to C++ and unique_ptr (currently this data is not freed)
204     t_disresdata*                 disres = nullptr;
205     std::unique_ptr<t_oriresdata> orires;
206 };
207
208 #endif