Introduce ObservablesHistory container
[gromacs.git] / src / gromacs / essentialdynamics / edsam.h
blobbd3a6b3e0cdb9333bfbbe6caf1c8fd02b18d3c21
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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, 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.
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.
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.
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.
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.
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.
37 /*! \libinternal \file
39 * \brief
40 * Declares functions to calculate both essential dynamics constraints
41 * as well as flooding potentials and forces.
43 * \authors Bert de Groot <bgroot@gwdg.de>, Oliver Lange <oliver.lange@tum.de>,
44 * Carsten Kutzner <ckutzne@gwdg.de>
46 * \inlibraryapi
48 #ifndef GMX_ESSENTIALDYNAMICS_EDSAM_H
49 #define GMX_ESSENTIALDYNAMICS_EDSAM_H
51 #include "gromacs/math/vectypes.h"
52 #include "gromacs/utility/basedefinitions.h"
55 /*! \brief Abstract type for essential dynamics
57 * The main type is defined only in edsam.c
59 typedef struct gmx_edsam *gmx_edsam_t;
61 struct edsamhistory_t;
62 struct gmx_domdec_t;
63 struct gmx_mtop_t;
64 struct gmx_output_env_t;
65 struct t_commrec;
66 struct t_filenm;
67 struct t_inputrec;
69 /*! \brief Applies essential dynamics constrains as defined in the .edi input file.
71 * \param ir MD input parameter record.
72 * \param step Number of the time step.
73 * \param cr Data needed for MPI communication.
74 * \param xs The local positions on this processor.
75 * \param v The local velocities.
76 * \param box The simulation box.
77 * \param ed The essential dynamics data.
79 void do_edsam(const t_inputrec *ir, gmx_int64_t step,
80 t_commrec *cr, rvec xs[], rvec v[], matrix box, gmx_edsam_t ed);
83 /*! \brief Reads in the .edi file containing the essential dynamics (ED) and flooding data.
85 * This function opens the ED input and output files, reads in all datasets it finds
86 * in the input file, and cross-checks whether the .edi file information is consistent
87 * with the ED data found in the checkpoint file (if present).
88 * gmx make_edi can be used to create an .edi input file.
90 * \param natoms Number of atoms of the whole MD system.
91 * \param oh Observables history, contains ED observables history
92 * \param nfile Number of entries (files) in the fnm structure.
93 * \param fnm The filenames struct; it contains also the names of the
94 * essential dynamics and flooding in + output files.
95 * \param Flags Flags passed over from main, used to determine
96 * whether we are appending.
97 * \param oenv Needed to open the output xvgr file.
98 * \param cr Data needed for MPI communication.
99 * \returns Pointer to the initialized essential dynamics / flooding data.
101 gmx_edsam_t ed_open(int natoms, struct ObservablesHistory *oh, int nfile, const t_filenm fnm[],
102 unsigned long Flags, const gmx_output_env_t *oenv, t_commrec *cr);
104 /*! \brief Initializes the essential dynamics and flooding module.
106 * \param mtop Molecular topology.
107 * \param ir MD input parameter record.
108 * \param cr Data needed for MPI communication.
109 * \param ed The essential dynamics data.
110 * \param x Positions of the whole MD system.
111 * \param box The simulation box.
112 * \param EDstate ED data stored in the checkpoint file.
114 void init_edsam(const gmx_mtop_t *mtop, const t_inputrec *ir, t_commrec *cr,
115 gmx_edsam_t ed, rvec x[], matrix box, edsamhistory_t *EDstate);
118 /*! \brief Make a selection of the home atoms for the ED groups.
120 * Should be called at every domain decomposition.
122 * \param dd Domain decomposition data.
123 * \param ed Essential dynamics and flooding data.
125 void dd_make_local_ed_indices(gmx_domdec_t *dd, gmx_edsam_t ed);
128 /*! \brief Evaluate the flooding potential(s) and forces as requested in the .edi input file.
130 * \param cr Data needed for MPI communication.
131 * \param ir MD input parameter record.
132 * \param x Positions on the local processor.
133 * \param force Forcefield forces to which the flooding forces are added.
134 * \param ed The essential dynamics data.
135 * \param box The simulation box.
136 * \param step Number of the time step.
137 * \param bNS Are we in a neighbor searching step?
139 void do_flood(t_commrec *cr, const t_inputrec *ir, rvec x[], rvec force[], gmx_edsam_t ed,
140 matrix box, gmx_int64_t step, gmx_bool bNS);
142 /*! \brief Clean up
144 * \param ed The essential dynamics data
146 void done_ed(gmx_edsam_t *ed);
148 #endif