927db7ac8910ff2b14f24e1c5872c932d1164f26
[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, 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 "gromacs/math/vectypes.h"
42 #include "gromacs/topology/idef.h"
43 #include "gromacs/utility/basedefinitions.h"
44 #include "gromacs/utility/real.h"
45
46 typedef real rvec5[5];
47
48 /* Distance restraining stuff */
49 typedef struct t_disresdata
50 {
51     int      dr_weighting; /* Weighting of pairs in one restraint              */
52     gmx_bool dr_bMixed;    /* Use sqrt of the instantaneous times              *
53                             * the time averaged violation                      */
54     real dr_fc;            /* Force constant for disres,                       *
55                             * which is multiplied by a (possibly)              *
56                             * different factor for each restraint              */
57     real  dr_tau;          /* Time constant for disres                    */
58     real  ETerm;           /* multiplication factor for time averaging         */
59     real  ETerm1;          /* 1 - ETerm1                                       */
60     real  exp_min_t_tau;   /* Factor for slowly switching on the force         */
61     int   nres;            /* The number of distance restraints                */
62     int   npair;           /* The number of distance restraint pairs           */
63     int   type_min;        /* The minimum iparam type index for restraints     */
64     real  sumviol;         /* The sum of violations                            */
65     real* rt;              /* The instantaneous distance (npair)               */
66     real* rm3tav;          /* The time averaged distance (npair)               */
67     real* Rtl_6;           /* The instantaneous r^-6 (nres)                    */
68     real* Rt_6;            /* The instantaneous ensemble averaged r^-6 (nres)  */
69     real* Rtav_6;          /* The time and ensemble averaged r^-6 (nres)       */
70     int   nsystems;        /* The number of systems for ensemble averaging     */
71
72     /* TODO: Implement a proper solution for parallel disre indexing */
73     const t_iatom* forceatomsStart; /* Pointer to the start of the disre forceatoms */
74 } t_disresdata;
75
76 /* All coefficients for the matrix equation for the orientation tensor */
77 struct OriresMatEq
78 {
79     real rhs[5];    /* The right hand side of the matrix equation */
80     real mat[5][5]; /* The matrix                                 */
81 };
82
83 /* Orientation restraining stuff */
84 typedef struct t_oriresdata
85 {
86     real         fc;            /* Force constant for the restraints                  */
87     real         edt;           /* Multiplication factor for time averaging           */
88     real         edt_1;         /* 1 - edt                                            */
89     real         exp_min_t_tau; /* Factor for slowly switching on the force         */
90     int          nr;            /* The number of orientation restraints               */
91     int          nex;           /* The number of experiments                          */
92     int          typeMin;       /* The minimum iparam type index for restraints       */
93     int          nref;          /* The number of atoms for the fit                    */
94     real*        mref;          /* The masses of the reference atoms                  */
95     rvec*        xref;          /* The reference coordinates for the fit (nref)       */
96     rvec*        xtmp;          /* Temporary array for fitting (nref)                 */
97     matrix       R;             /* Rotation matrix to rotate to the reference coor.   */
98     tensor*      S;             /* Array of order tensors for each experiment (nexp)  */
99     rvec5*       Dinsl;         /* The order matrix D for all restraints (nr x 5)     */
100     rvec5*       Dins;          /* The ensemble averaged D (nr x 5)                   */
101     rvec5*       Dtav;          /* The time and ensemble averaged D (nr x 5)          */
102     real*        oinsl;         /* The calculated instantaneous orientations          */
103     real*        oins;          /* The calculated emsemble averaged orientations      */
104     real*        otav;          /* The calculated time and ensemble averaged orient.  */
105     real         rmsdev;        /* The weighted (using kfac) RMS deviation            */
106     OriresMatEq* tmpEq;         /* An temporary array of matrix + rhs                 */
107     real*        eig;           /* Eigenvalues/vectors, for output only (nex x 12)    */
108
109     /* variables for diagonalization with diagonalize_orires_tensors()*/
110     double** M;
111     double*  eig_diag;
112     double** v;
113 } t_oriresdata;
114
115 typedef struct bondedtable_t
116 {
117     int   n;     /* n+1 is the number of points */
118     real  scale; /* distance between two points */
119     real* data;  /* the actual table data, per point there are 4 numbers */
120 } bondedtable_t;
121
122 /*
123  * Data struct used in the force calculation routines
124  * for storing the tables for bonded interactions and
125  * for storing information which is needed in following steps
126  * (for instance for time averaging in distance retraints)
127  * or for storing output, since force routines only return the potential.
128  */
129 typedef struct t_fcdata
130 {
131     bondedtable_t* bondtab;
132     bondedtable_t* angletab;
133     bondedtable_t* dihtab;
134
135     t_disresdata disres;
136     t_oriresdata orires;
137 } t_fcdata;
138
139 #endif