Add replacements for pbc enumerations
[gromacs.git] / src / gromacs / domdec / utility.h
blobe4d2519f876018f84b70570046c79d2eea98dc9f
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2018, 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.
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.
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.
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.
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.
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.
35 /*! \internal \file
37 * \brief Declares utility functions used in the domain decomposition module.
39 * \author Berk Hess <hess@kth.se>
40 * \ingroup module_domdec
42 #ifndef GMX_DOMDEC_DOMDEC_UTILITY_H
43 #define GMX_DOMDEC_DOMDEC_UTILITY_H
45 #include "gromacs/math/paddedvector.h"
46 #include "gromacs/mdtypes/forcerec.h"
48 #include "domdec_internal.h"
50 /*! \brief Returns true if the DLB state indicates that the balancer is on. */
51 static inline bool isDlbOn(const gmx_domdec_comm_t *comm)
53 return (comm->dlbState == DlbState::onCanTurnOff ||
54 comm->dlbState == DlbState::onUser);
57 /*! \brief Returns true if the DLB state indicates that the balancer is off/disabled.
59 static inline bool isDlbDisabled(const gmx_domdec_comm_t *comm)
61 return (comm->dlbState == DlbState::offUser ||
62 comm->dlbState == DlbState::offForever);
65 /*! \brief Returns the character, x/y/z, corresponding to dimension dim */
66 char dim2char(int dim);
68 /*! \brief Sets matrix to convert from Cartesian to lattice coordinates */
69 void make_tric_corr_matrix(int npbcdim, const matrix box, matrix tcm);
71 /*! \brief Ensure box obeys the screw restrictions, fatal error if not */
72 void check_screw_box(const matrix box);
74 /*! \brief Return the charge group information flags for charge group cg */
75 static inline int ddcginfo(const cginfo_mb_t *cginfo_mb,
76 int cg)
78 while (cg >= cginfo_mb->cg_end)
80 cginfo_mb++;
83 return cginfo_mb->cginfo[(cg - cginfo_mb->cg_start) % cginfo_mb->cg_mod];
86 /*! \brief Returns the number of MD steps for which load has been recorded */
87 static inline int dd_load_count(const gmx_domdec_comm_t *comm)
89 return (comm->eFlop ? comm->flop_n : comm->cycl_n[ddCyclF]);
92 /*! \brief Resize the state and f, if !=nullptr, to natoms */
93 void dd_resize_state(t_state *state,
94 PaddedVector<gmx::RVec> *f,
95 int natoms);
97 /*! \brief Enrsure fr, state and f, if != nullptr, can hold numChargeGroups atoms for the Verlet scheme and charge groups for the group scheme */
98 void dd_check_alloc_ncg(t_forcerec *fr,
99 t_state *state,
100 PaddedVector<gmx::RVec> *f,
101 int numChargeGroups);
103 /*! \brief Returns a domain-to-domain cutoff distance given an atom-to-atom cutoff */
104 static inline real
105 atomToAtomIntoDomainToDomainCutoff(const gmx_domdec_comm_t &comm,
106 real cutoff)
108 if (comm.useUpdateGroups)
110 GMX_ASSERT(comm.updateGroupsCog, "updateGroupsCog should be initialized here");
112 cutoff += 2*comm.updateGroupsCog->maxUpdateGroupRadius();
115 return cutoff;
118 /*! \brief Returns an atom-to-domain cutoff distance given a domain-to-domain cutoff */
119 static inline real
120 domainToDomainIntoAtomToDomainCutoff(const gmx_domdec_comm_t &comm,
121 real cutoff)
123 if (comm.useUpdateGroups)
125 GMX_ASSERT(comm.updateGroupsCog, "updateGroupsCog should be initialized here");
127 cutoff -= comm.updateGroupsCog->maxUpdateGroupRadius();
130 return cutoff;
133 #endif