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,2018,2019, 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.
39 #include "gromacs/fileio/xdrf.h"
40 #include "gromacs/math/vectypes.h"
41 #include "gromacs/utility/fatalerror.h"
42 #include "gromacs/utility/smalloc.h"
44 int xdr_real(XDR *xdrs, real *r)
51 ret = xdr_float(xdrs, &f);
56 return xdr_float(xdrs, static_cast<float *>(r));
60 int xdr3drcoord(XDR *xdrs, real *fp, int *size, real *precision)
70 gmx_fatal(FARGS, "Don't know what to malloc for ffp, isize = %d", isize);
75 for (i = 0; (i < isize); i++)
80 ret = xdr3dfcoord(xdrs, ffp, size, &fprec);
83 for (i = 0; (i < isize); i++)
91 return xdr3dfcoord(xdrs, reinterpret_cast<float *>(fp), size, reinterpret_cast<float *>(precision));
95 int xdr_int32(XDR *xdrs, int32_t *i)
97 // Note that this implementation assumes that an int is at least
98 // 32 bits, which is not strictly required by the language, but
99 // good enough in practice on 32- or 64-bit systems. GROMACS
100 // requires 64-bit systems.
101 static_assert(sizeof(int) >= 4, "XDR handling assumes that an int32_t can be stored in an int");
102 int temporary = static_cast<int>(*i);
103 int ret = xdr_int(xdrs, &temporary);
104 *i = static_cast<int32_t>(temporary);
109 int xdr_int64(XDR *xdrs, int64_t *i)
111 // Note that this implementation assumes that an int is at least
112 // 32 bits, which is not strictly required by the language, but
113 // good enough in practice on 32- or 64-bit systems. GROMACS
114 // requires 64-bit systems.
115 static_assert(2*sizeof(int) >= 8, "XDR handling assumes that an int64_t can be stored in two ints");
119 static const int64_t two_p32_m1 = 0xFFFFFFFF;
120 int64_t imaj64, imin64;
122 imaj64 = ((*i)>>32) & two_p32_m1;
123 imin64 = (*i) & two_p32_m1;
124 imaj = static_cast<int>(imaj64);
125 imin = static_cast<int>(imin64);
126 ret = xdr_int(xdrs, &imaj);
127 ret |= xdr_int(xdrs, &imin);
129 *i = ((static_cast<int64_t>(imaj) << 32) | (static_cast<int64_t>(imin) & two_p32_m1));