Tests that verify serializer behaviours
[gromacs.git] / src / gromacs / utility / inmemoryserializer.cpp
blob2dec0caec9ae0345c3609152cf09e265d63be4d0
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2016,2017,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.
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
36 * \brief
37 * Defines gmx::ISerializer implementation for in-memory serialization.
39 * \author Teemu Murtola <teemu.murtola@gmail.com>
40 * \ingroup module_utility
42 #include "gmxpre.h"
44 #include "inmemoryserializer.h"
46 #include "config.h"
48 #include <algorithm>
49 #include <vector>
51 namespace gmx
54 namespace
57 template<typename T>
58 class CharBuffer
60 public:
61 static constexpr size_t ValueSize = sizeof(T);
63 explicit CharBuffer(T value) { u.v = value; }
64 explicit CharBuffer(const char buffer[]) { std::copy(buffer, buffer + ValueSize, u.c); }
66 T value() const { return u.v; }
68 void appendTo(std::vector<char>* buffer)
70 buffer->insert(buffer->end(), u.c, u.c + ValueSize);
73 private:
74 union {
75 char c[ValueSize];
76 T v;
77 } u;
80 //! Return \c value with the byte order swapped.
81 template<typename T>
82 T swapEndian(const T& value)
84 union {
85 T value_;
86 std::array<char, sizeof(T)> valueAsCharArray_;
87 } endianessSwappedValue;
89 endianessSwappedValue.value_ = value;
90 int hiByte = sizeof(T) - 1;
91 for (int loByte = 0; hiByte > loByte; loByte++, hiByte--)
93 std::swap(endianessSwappedValue.valueAsCharArray_[loByte],
94 endianessSwappedValue.valueAsCharArray_[hiByte]);
97 return endianessSwappedValue.value_;
100 /*! \brief Change the host-dependent endian settings to either Swap or DoNotSwap.
102 * \param endianSwapBehavior input swap behavior, might depend on host.
104 * \return Host-independent setting, either Swap or DoNotSwap. */
105 EndianSwapBehavior setEndianSwapBehaviorFromHost(EndianSwapBehavior endianSwapBehavior)
107 if (endianSwapBehavior == EndianSwapBehavior::SwapIfHostIsBigEndian)
109 return GMX_INTEGER_BIG_ENDIAN ? EndianSwapBehavior::Swap : EndianSwapBehavior::DoNotSwap;
111 else if (endianSwapBehavior == EndianSwapBehavior::SwapIfHostIsLittleEndian)
113 return GMX_INTEGER_BIG_ENDIAN ? EndianSwapBehavior::DoNotSwap : EndianSwapBehavior::Swap;
115 else
117 return endianSwapBehavior;
121 } // namespace
123 /********************************************************************
124 * InMemorySerializer
127 class InMemorySerializer::Impl
129 public:
130 Impl(EndianSwapBehavior endianSwapBehavior) :
131 endianSwapBehavior_(setEndianSwapBehaviorFromHost(endianSwapBehavior))
135 template<typename T>
136 void doValue(T value)
138 if (endianSwapBehavior_ == EndianSwapBehavior::Swap)
140 CharBuffer<T>(swapEndian(value)).appendTo(&buffer_);
142 else
144 CharBuffer<T>(value).appendTo(&buffer_);
147 void doString(const std::string& value)
149 doValue<uint64_t>(value.size());
150 buffer_.insert(buffer_.end(), value.begin(), value.end());
152 void doOpaque(const char* data, std::size_t size)
154 buffer_.insert(buffer_.end(), data, data + size);
157 std::vector<char> buffer_;
158 EndianSwapBehavior endianSwapBehavior_;
161 InMemorySerializer::InMemorySerializer(EndianSwapBehavior endianSwapBehavior) :
162 impl_(new Impl(endianSwapBehavior))
166 InMemorySerializer::~InMemorySerializer() = default;
168 std::vector<char> InMemorySerializer::finishAndGetBuffer()
170 return std::move(impl_->buffer_);
173 void InMemorySerializer::doBool(bool* value)
175 impl_->doValue(*value);
178 void InMemorySerializer::doUChar(unsigned char* value)
180 impl_->doValue(*value);
183 void InMemorySerializer::doChar(char* value)
185 impl_->doValue(*value);
188 void InMemorySerializer::doUShort(unsigned short* value)
190 impl_->doValue(*value);
193 void InMemorySerializer::doInt(int* value)
195 impl_->doValue(*value);
198 void InMemorySerializer::doInt32(int32_t* value)
200 impl_->doValue(*value);
203 void InMemorySerializer::doInt64(int64_t* value)
205 impl_->doValue(*value);
208 void InMemorySerializer::doFloat(float* value)
210 impl_->doValue(*value);
213 void InMemorySerializer::doDouble(double* value)
215 impl_->doValue(*value);
218 void InMemorySerializer::doReal(real* value)
220 impl_->doValue(*value);
223 void InMemorySerializer::doRvec(rvec* value)
225 for (int d = 0; d < DIM; d++)
227 doReal(&(*value)[d]);
231 void InMemorySerializer::doIvec(ivec* value)
233 for (int d = 0; d < DIM; d++)
235 doInt(&(*value)[d]);
239 void InMemorySerializer::doString(std::string* value)
241 impl_->doString(*value);
244 void InMemorySerializer::doOpaque(char* data, std::size_t size)
246 impl_->doOpaque(data, size);
249 /********************************************************************
250 * InMemoryDeserializer
253 class InMemoryDeserializer::Impl
255 public:
256 explicit Impl(ArrayRef<const char> buffer, bool sourceIsDouble, EndianSwapBehavior endianSwapBehavior) :
257 buffer_(buffer),
258 sourceIsDouble_(sourceIsDouble),
259 pos_(0),
260 endianSwapBehavior_(setEndianSwapBehaviorFromHost(endianSwapBehavior))
264 template<typename T>
265 void doValue(T* value)
267 if (endianSwapBehavior_ == EndianSwapBehavior::Swap)
269 *value = swapEndian(CharBuffer<T>(&buffer_[pos_]).value());
271 else
273 *value = CharBuffer<T>(&buffer_[pos_]).value();
275 pos_ += CharBuffer<T>::ValueSize;
277 void doString(std::string* value)
279 uint64_t size;
280 doValue<uint64_t>(&size);
281 *value = std::string(&buffer_[pos_], size);
282 pos_ += size;
284 void doOpaque(char* data, std::size_t size)
286 std::copy(&buffer_[pos_], &buffer_[pos_ + size], data);
287 pos_ += size;
290 ArrayRef<const char> buffer_;
291 bool sourceIsDouble_;
292 size_t pos_;
293 EndianSwapBehavior endianSwapBehavior_;
296 InMemoryDeserializer::InMemoryDeserializer(ArrayRef<const char> buffer,
297 bool sourceIsDouble,
298 EndianSwapBehavior endianSwapBehavior) :
299 impl_(new Impl(buffer, sourceIsDouble, endianSwapBehavior))
303 InMemoryDeserializer::~InMemoryDeserializer() = default;
305 bool InMemoryDeserializer::sourceIsDouble() const
307 return impl_->sourceIsDouble_;
310 void InMemoryDeserializer::doBool(bool* value)
312 impl_->doValue(value);
315 void InMemoryDeserializer::doUChar(unsigned char* value)
317 impl_->doValue(value);
320 void InMemoryDeserializer::doChar(char* value)
322 impl_->doValue(value);
325 void InMemoryDeserializer::doUShort(unsigned short* value)
327 impl_->doValue(value);
330 void InMemoryDeserializer::doInt(int* value)
332 impl_->doValue(value);
335 void InMemoryDeserializer::doInt32(int32_t* value)
337 impl_->doValue(value);
340 void InMemoryDeserializer::doInt64(int64_t* value)
342 impl_->doValue(value);
345 void InMemoryDeserializer::doFloat(float* value)
347 impl_->doValue(value);
350 void InMemoryDeserializer::doDouble(double* value)
352 impl_->doValue(value);
355 void InMemoryDeserializer::doReal(real* value)
357 if (sourceIsDouble())
359 double temp = 0.0;
360 doDouble(&temp);
361 *value = temp;
363 else
365 float temp = 0.0;
366 doFloat(&temp);
367 *value = temp;
371 void InMemoryDeserializer::doRvec(rvec* value)
373 for (int d = 0; d < DIM; d++)
375 doReal(&(*value)[d]);
379 void InMemoryDeserializer::doIvec(ivec* value)
381 for (int d = 0; d < DIM; d++)
383 doInt(&(*value)[d]);
387 void InMemoryDeserializer::doString(std::string* value)
389 impl_->doString(value);
392 void InMemoryDeserializer::doOpaque(char* data, std::size_t size)
394 impl_->doOpaque(data, size);
397 } // namespace gmx