353dc327b0f39a0c4207c7379a080568b583748a
[alexxy/gromacs.git] / src / gromacs / domdec / utility.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2018,2019, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*! \internal \file
36  *
37  * \brief Declares utility functions used in the domain decomposition module.
38  *
39  * \author Berk Hess <hess@kth.se>
40  * \ingroup module_domdec
41  */
42 #ifndef GMX_DOMDEC_DOMDEC_UTILITY_H
43 #define GMX_DOMDEC_DOMDEC_UTILITY_H
44
45 #include "gromacs/gpu_utils/hostallocator.h"
46 #include "gromacs/mdtypes/forcerec.h"
47 #include "gromacs/utility/arrayref.h"
48
49 #include "domdec_internal.h"
50
51 /*! \brief Returns true if the DLB state indicates that the balancer is on. */
52 static inline bool isDlbOn(const gmx_domdec_comm_t* comm)
53 {
54     return (comm->dlbState == DlbState::onCanTurnOff || comm->dlbState == DlbState::onUser);
55 };
56
57 /*! \brief Returns true if the DLB state indicates that the balancer is off/disabled.
58  */
59 static inline bool isDlbDisabled(const DlbState& dlbState)
60 {
61     return (dlbState == DlbState::offUser || dlbState == DlbState::offForever);
62 };
63
64 /*! \brief Returns true if the DLB state indicates that the balancer is off/disabled.
65  */
66 static inline bool isDlbDisabled(const gmx_domdec_comm_t* comm)
67 {
68     return isDlbDisabled(comm->dlbState);
69 };
70
71 /*! \brief Returns the character, x/y/z, corresponding to dimension dim */
72 char dim2char(int dim);
73
74 /*! \brief Sets matrix to convert from Cartesian to lattice coordinates */
75 void make_tric_corr_matrix(int npbcdim, const matrix box, matrix tcm);
76
77 /*! \brief Ensure box obeys the screw restrictions, fatal error if not */
78 void check_screw_box(const matrix box);
79
80 /*! \brief Return the charge group information flags for charge group cg */
81 static inline int ddcginfo(gmx::ArrayRef<const cginfo_mb_t> cginfo_mb, int cg)
82 {
83     size_t index = 0;
84     while (cg >= cginfo_mb[index].cg_end)
85     {
86         index++;
87     }
88     const cginfo_mb_t& cgimb = cginfo_mb[index];
89
90     return cgimb.cginfo[(cg - cgimb.cg_start) % cgimb.cg_mod];
91 };
92
93 /*! \brief Returns the number of MD steps for which load has been recorded */
94 static inline int dd_load_count(const gmx_domdec_comm_t* comm)
95 {
96     return (comm->ddSettings.eFlop ? comm->flop_n : comm->cycl_n[ddCyclF]);
97 }
98
99 /*! \brief Resize the state and f, if !=nullptr, to natoms
100  *
101  * \param[in]  state   Current state
102  * \param[in]  f       The vector of forces to be resized
103  * \param[out] natoms  New number of atoms to accommodate
104  */
105 void dd_resize_state(t_state* state, gmx::PaddedHostVector<gmx::RVec>* f, int natoms);
106
107 /*! \brief Enrsure fr, state and f, if != nullptr, can hold numChargeGroups atoms for the Verlet scheme and charge groups for the group scheme
108  *
109  * \param[in]  fr               Force record
110  * \param[in]  state            Current state
111  * \param[in]  f                The force buffer
112  * \param[out] numChargeGroups  Number of charged groups
113  */
114 void dd_check_alloc_ncg(t_forcerec* fr, t_state* state, gmx::PaddedHostVector<gmx::RVec>* f, int numChargeGroups);
115
116 /*! \brief Returns a domain-to-domain cutoff distance given an atom-to-atom cutoff */
117 static inline real atomToAtomIntoDomainToDomainCutoff(const DDSystemInfo& systemInfo, real cutoff)
118 {
119     if (systemInfo.useUpdateGroups)
120     {
121         cutoff += 2 * systemInfo.maxUpdateGroupRadius;
122     }
123
124     return cutoff;
125 }
126
127 /*! \brief Returns an atom-to-domain cutoff distance given a domain-to-domain cutoff */
128 static inline real domainToDomainIntoAtomToDomainCutoff(const DDSystemInfo& systemInfo, real cutoff)
129 {
130     if (systemInfo.useUpdateGroups)
131     {
132         cutoff -= systemInfo.maxUpdateGroupRadius;
133     }
134
135     return cutoff;
136 }
137
138 #endif