SYCL: Avoid using no_init read accessor in rocFFT
[alexxy/gromacs.git] / src / gromacs / pbcutil / pbc.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_PBCUTIL_PBC_H
39 #define GMX_PBCUTIL_PBC_H
40
41 #include <stdio.h>
42
43 #include <string>
44
45 #include "gromacs/math/vectypes.h"
46 #include "gromacs/mdtypes/md_enums.h"
47 #include "gromacs/utility/basedefinitions.h"
48 #include "gromacs/utility/enumerationhelpers.h"
49 #include "gromacs/utility/real.h"
50
51 struct gmx_domdec_t;
52 struct gmx_mtop_t;
53
54 namespace gmx
55 {
56 template<typename>
57 class ArrayRef;
58 } // namespace gmx
59
60 //! Names for all values in PBC types enumeration
61 extern const gmx::EnumerationArray<PbcType, std::string> c_pbcTypeNames;
62
63 /* Maximum number of combinations of single triclinic box vectors
64  * required to shift atoms that are within a brick of the size of
65  * the diagonal of the box to within the maximum cut-off distance.
66  */
67 #define MAX_NTRICVEC 12
68
69 /*! \brief Structure containing info on periodic boundary conditions */
70 typedef struct t_pbc
71 {
72     //! The PBC type
73     PbcType pbcType;
74     //! Number of dimensions in which PBC is exerted
75     int ndim_ePBC;
76     /*! \brief Determines how to compute distance vectors.
77      *
78      *  Indicator of how to compute distance vectors, depending
79      *  on PBC type (depends on pbcType and dimensions with(out) DD)
80      *  and the box angles.
81      */
82     int pbcTypeDX;
83     /*! \brief Used for selecting which dimensions to use in PBC.
84      *
85      *  In case of 1-D PBC this indicates which dimension is used,
86      *  in case of 2-D PBC this indicates the opposite
87      */
88     int dim;
89     //! The simulation box
90     matrix box;
91     //! The lengths of the diagonal of the full box
92     rvec fbox_diag;
93     //! Halve of the above
94     rvec hbox_diag;
95     //! Negative of the above
96     rvec mhbox_diag;
97     //! Maximum allowed cutoff squared for the box and PBC used
98     real max_cutoff2;
99     /*! \brief Number of triclinic shift vectors.
100      *
101      *  Number of triclinic shift vectors depends on the skewedness
102      *  of the box, that is mostly on the angles. For triclinic boxes
103      *  we first take the closest image along each Cartesian dimension
104      *  independently. When the resulting distance^2 is larger than
105      *  max_cutoff2, up to ntric_vec triclinic shift vectors need to
106      *  be tried. Because of the restrictions imposed on the unit-cell
107      *  by GROMACS, ntric_vec <= MAX_NTRICVEC = 12.
108      */
109     int ntric_vec;
110     //! The triclinic shift vectors in grid cells. Internal use only.
111     ivec tric_shift[MAX_NTRICVEC];
112     //!  The triclinic shift vectors in length units
113     rvec tric_vec[MAX_NTRICVEC];
114 } t_pbc;
115
116 #define TRICLINIC(box) ((box)[YY][XX] != 0 || (box)[ZZ][XX] != 0 || (box)[ZZ][YY] != 0)
117
118 #define NTRICIMG 14
119 #define NCUCVERT 24
120 #define NCUCEDGE 36
121
122 enum
123 {
124     ecenterTRIC, /* 0.5*(a+b+c)                  */
125     ecenterRECT, /* (0.5*a[x],0.5*b[y],0.5*c[z]) */
126     ecenterZERO, /* (0,0,0)                      */
127     ecenterDEF = ecenterTRIC
128 };
129
130 /*! \brief Returns the number of dimensions that use pbc
131  *
132  * \param[in] pbcType The periodic boundary condition type
133  * \return the number of dimensions that use pbc, starting at X
134  */
135 int numPbcDimensions(PbcType pbcType);
136
137 /*! \brief Dump the contents of the pbc structure to the file
138  *
139  * \param[in] fp  The file pointer to write to
140  * \param[in] pbc The periodic boundary condition information structure
141  */
142 void dump_pbc(FILE* fp, t_pbc* pbc);
143
144 /*! \brief Check the box for consistency
145  *
146  * When \p pbcType=PbcTypes::Unset, the type of pbc is guessed from the box matrix.
147  *
148  * \param[in] pbcType The pbc identifier
149  * \param[in] box     The box matrix
150  * \return NULL if the box is supported by Gromacs.
151  *         Otherwise returns a string with the problem.
152  */
153 const char* check_box(PbcType pbcType, const matrix box);
154
155 /*! \brief Creates box matrix from edge lengths and angles.
156  *
157  * \param[in,out] box        The box matrix
158  * \param[in] vec            The edge lengths
159  * \param[in] angleInDegrees The angles
160  */
161 void matrix_convert(matrix box, const rvec vec, const rvec angleInDegrees);
162
163 /*! \brief Compute the maximum cutoff for the box
164
165  * Returns the square of the maximum cut-off allowed for the box,
166  * taking into account that the grid neighborsearch code and pbc_dx
167  * only check combinations of single box-vector shifts.
168  *
169  * \param[in] pbcType The pbc identifier
170  * \param[in] box  The box matrix
171  * \return the maximum cut-off.
172  */
173 real max_cutoff2(PbcType pbcType, const matrix box);
174
175 /*! \brief Guess PBC type
176  *
177  * Guesses the type of periodic boundary conditions using the box
178  *
179  * \param[in] box  The box matrix
180  * \return The pbc type identifier
181  */
182 PbcType guessPbcType(const matrix box);
183
184 /*! \brief Corrects the box if necessary
185  *
186  * Checks for un-allowed box angles and corrects the box.
187  *
188  * \param[in] fplog File for debug output
189  * \param[in] step  The MD step number
190  * \param[in] box   The simulation cell
191  * \return TRUE when the box was corrected.
192  */
193 bool correct_box(FILE* fplog, int step, tensor box);
194
195 /*! \brief Initiate the periodic boundary condition algorithms.
196  *
197  * pbc_dx will not use pbc and return the normal difference vector
198  * when one or more of the diagonal elements of box are zero.
199  * When \p pbcType=PbcType::Unset, the type of pbc is guessed from the box matrix.
200  *
201  * \param[in,out] pbc The pbc information structure
202  * \param[in] pbcType The PBC identifier
203  * \param[in] box     The box tensor
204  */
205 void set_pbc(t_pbc* pbc, PbcType pbcType, const matrix box);
206
207 /*! \brief Initiate the periodic boundary condition algorithms.
208  *
209  * As set_pbc, but additionally sets that correct distances can
210  * be obtained using (combinations of) single box-vector shifts.
211  * Should be used with pbc_dx_aiuc.
212  * If domdecCells!=NULL pbc is not used for directions
213  * with dd->nc[i]==1 with bSingleDir==TRUE or
214  * with dd->nc[i]<=2 with bSingleDir==FALSE.
215  * Note that when no PBC is required only pbc->pbcType is set,
216  * the rest of the struct will be invalid.
217  *
218  * \param[in,out] pbc     The pbc information structure
219  * \param[in] pbcType     The PBC identifier
220  * \param[in] domdecCells 3D integer vector describing the number of DD cells
221  *                        or nullptr if not using DD.
222  * \param[in] bSingleDir  TRUE if DD communicates only in one direction along dimensions
223  * \param[in] box         The box tensor
224  * \return the pbc structure when pbc operations are required, NULL otherwise.
225  */
226 t_pbc* set_pbc_dd(t_pbc* pbc, PbcType pbcType, const ivec domdecCells, bool bSingleDir, const matrix box);
227
228 /*! \brief Compute distance with PBC
229  *
230  * Calculate the correct distance vector from x2 to x1 and put it in dx.
231  * set_pbc must be called before ever calling this routine.
232  *
233  * Note that for triclinic boxes that do not obey the GROMACS unit-cell
234  * restrictions, pbc_dx and pbc_dx_aiuc will not correct for PBC.
235  * \param[in,out] pbc The pbc information structure
236  * \param[in]    x1  Coordinates for particle 1
237  * \param[in]    x2  Coordinates for particle 2
238  * \param[out]   dx  Distance vector
239  */
240 void pbc_dx(const t_pbc* pbc, const rvec x1, const rvec x2, rvec dx);
241
242 /*! \brief Compute distance vector for simple PBC types
243  *
244  * Calculate the correct distance vector from x2 to x1 and put it in dx,
245  * This function can only be used when all atoms are in the rectangular
246  * or triclinic unit-cell.
247  * set_pbc_dd or set_pbc must be called before ever calling this routine.
248  * \param[in,out] pbc The pbc information structure
249  * \param[in]    x1  Coordinates for particle 1
250  * \param[in]    x2  Coordinates for particle 2
251  * \param[out]   dx  Distance vector
252  * \return the ishift required to shift x1 at closest distance to x2;
253  * i.e. if 0<=ishift<c_numShiftVectors then x1 - x2 + shift_vec[ishift] = dx
254  * (see calc_shifts below on how to obtain shift_vec)
255  */
256 int pbc_dx_aiuc(const t_pbc* pbc, const rvec x1, const rvec x2, rvec dx);
257
258 /*! \brief Compute distance with PBC
259  *
260  * As pbc_dx, but for double precision vectors.
261  * set_pbc must be called before ever calling this routine.
262  * \param[in,out] pbc The pbc information structure
263  * \param[in]    x1  Coordinates for particle 1
264  * \param[in]    x2  Coordinates for particle 2
265  * \param[out]   dx  Distance vector
266  */
267 void pbc_dx_d(const t_pbc* pbc, const dvec x1, const dvec x2, dvec dx);
268
269 /*! \brief Computes shift vectors
270  *
271  * This routine calculates ths shift vectors necessary to use the
272  * neighbor searching routine.
273  * \param[in]  box       The simulation box
274  * \param[out] shift_vec The shifting vectors
275  */
276 void calc_shifts(const matrix box, gmx::ArrayRef<gmx::RVec> shift_vec);
277
278 /*! \brief Calculates the center of the box.
279  *
280  * See the description for the enum ecenter above.
281  * \param[in]  ecenter    Description of center type
282  * \param[in]  box        The simulation box
283  * \param[out] box_center The center of the box
284  */
285 void calc_box_center(int ecenter, const matrix box, rvec box_center);
286
287 /*! \brief Calculates the NTRICIMG box images
288  *
289  * \param[in]  box The simulation box
290  * \param[out] img The triclinic box images
291  */
292 void calc_triclinic_images(const matrix box, rvec img[]);
293
294 /*! \brief Calculates the NCUCVERT vertices of a compact unitcell
295  *
296  * \param[in]  ecenter The center type
297  * \param[in]  box     The simulation box
298  * \param[out] vert    The vertices
299  */
300 void calc_compact_unitcell_vertices(int ecenter, const matrix box, rvec vert[]);
301
302 /*! \brief Compute unitcell edges
303  *
304  * \return an array of unitcell edges of length NCUCEDGE*2,
305  * this is an index in vert[], which is calculated by calc_unitcell_vertices.
306  * The index consists of NCUCEDGE pairs of vertex indices.
307  * The index does not change, so it needs to be retrieved only once.
308  */
309 int* compact_unitcell_edges();
310
311 /*! \brief Put atoms inside the simulations box
312  *
313  * These routines puts ONE or ALL atoms in the box, not caring
314  * about charge groups!
315  * Also works for triclinic cells.
316  *
317  * \param[in]     pbcType The pbc type
318  * \param[in]     box     The simulation box
319  * \param[in,out] x       The coordinates of the atoms
320  */
321 void put_atoms_in_box(PbcType pbcType, const matrix box, gmx::ArrayRef<gmx::RVec> x);
322
323 /*! \brief Parallellizes put_atoms_in_box()
324  *
325  * This wrapper function around put_atoms_in_box() with the ugly manual
326  * workload splitting is needed to avoid silently introducing multithreading
327  * in tools.
328  *
329  * \param[in]     pbcType    The pbc type
330  * \param[in]     box        The simulation box
331  * \param[in,out] x          The coordinates of the atoms
332  * \param[in]     nth        number of threads to be used in the given module
333  */
334 void put_atoms_in_box_omp(PbcType pbcType, const matrix box, gmx::ArrayRef<gmx::RVec> x, gmx_unused int nth);
335
336 /*! \brief Put atoms inside triclinic box
337  *
338  * This puts ALL atoms in the triclinic unit cell, centered around the
339  * box center as calculated by calc_box_center.
340  * \param[in]    ecenter The pbc center type
341  * \param[in]    box     The simulation box
342  * \param[in,out] x       The coordinates of the atoms
343  */
344 void put_atoms_in_triclinic_unitcell(int ecenter, const matrix box, gmx::ArrayRef<gmx::RVec> x);
345
346 /*! \brief Put atoms inside the unitcell
347  *
348  * This puts ALL atoms at the closest distance for the center of the box
349  * as calculated by calc_box_center.
350  * When \p pbcType=PbcType::Unset, the type of pbc is guessed from the box matrix.
351  *
352  * \param[in]    pbcType The pbc type
353  * \param[in]    ecenter The pbc center type
354  * \param[in]    box     The simulation box
355  * \param[in,out] x      The coordinates of the atoms
356  */
357 void put_atoms_in_compact_unitcell(PbcType pbcType, int ecenter, const matrix box, gmx::ArrayRef<gmx::RVec> x);
358
359 /*! \brief Make all molecules whole by shifting positions
360  *
361  * \param[in]     fplog     Log file
362  * \param[in]     pbcType   The PBC type
363  * \param[in]     box       The simulation box
364  * \param[in]     mtop      System topology definition
365  * \param[in,out] x         The coordinates of the atoms
366  */
367 void do_pbc_first_mtop(FILE* fplog, PbcType pbcType, const matrix box, const gmx_mtop_t* mtop, rvec x[]);
368
369 /*! \brief Make molecules consisting of multiple charge groups whole by shifting positions
370  *
371  * \param[in]     pbcType   The PBC type
372  * \param[in]     box       The simulation box
373  * \param[in]     mtop      System topology definition
374  * \param[in,out] x         The coordinates of the atoms
375  */
376 void do_pbc_mtop(PbcType pbcType, const matrix box, const gmx_mtop_t* mtop, rvec x[]);
377
378 #endif