Extend serialization to gmx_int64_t
[gromacs.git] / src / gromacs / utility / inmemoryserializer.cpp
bloba741b6760773068b43d2287ab939378e0c4043f8
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2016,2017, 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 #include "gmxpre.h"
37 #include "inmemoryserializer.h"
39 #include <algorithm>
40 #include <vector>
42 #include "gromacs/utility/gmxassert.h"
44 namespace gmx
47 namespace
50 //! Returns offset to add to `pos` to get it aligned at `alignment` bytes.
51 size_t alignedOffset(size_t pos, size_t alignment)
53 return (alignment - pos % alignment) % alignment;
56 } // namespace
58 /********************************************************************
59 * InMemorySerializer
62 class InMemorySerializer::Impl
64 public:
65 template <typename T>
66 void doValue(T value)
68 // Here, we assume that the vector memory buffer start is aligned,
69 // similar to what malloc() guarantees.
70 const size_t size = buffer_.size();
71 const size_t pos = size + alignedOffset(size, alignof(T));
72 buffer_.resize(pos + sizeof(T));
73 *reinterpret_cast<T *>(&buffer_[pos]) = value;
75 void doString(const std::string &value)
77 doValue<size_t>(value.size());
78 const size_t pos = buffer_.size();
79 buffer_.resize(pos + value.size());
80 std::copy(value.begin(), value.end(), &buffer_[pos]);
83 std::vector<char> buffer_;
86 InMemorySerializer::InMemorySerializer()
87 : impl_(new Impl)
91 InMemorySerializer::~InMemorySerializer()
95 std::vector<char> InMemorySerializer::finishAndGetBuffer()
97 return std::move(impl_->buffer_);
100 void InMemorySerializer::doUChar(unsigned char *value)
102 impl_->doValue(*value);
105 void InMemorySerializer::doInt(int *value)
107 impl_->doValue(*value);
110 void InMemorySerializer::doInt64(gmx_int64_t *value)
112 impl_->doValue(*value);
115 void InMemorySerializer::doFloat(float *value)
117 impl_->doValue(*value);
120 void InMemorySerializer::doDouble(double *value)
122 impl_->doValue(*value);
125 void InMemorySerializer::doString(std::string *value)
127 impl_->doString(*value);
130 /********************************************************************
131 * InMemoryDeserializer
134 class InMemoryDeserializer::Impl
136 public:
137 explicit Impl(const std::vector<char> &buffer)
138 : buffer_(buffer), pos_(0)
142 template <typename T>
143 void doValue(T *value)
145 pos_ += alignedOffset(pos_, alignof(T));
146 *value = *reinterpret_cast<const T *>(&buffer_[pos_]);
147 pos_ += sizeof(T);
149 void doString(std::string *value)
151 size_t size;
152 doValue<size_t>(&size);
153 *value = std::string(&buffer_[pos_], &buffer_[pos_ + size]);
154 pos_ += size;
157 const std::vector<char> &buffer_;
158 size_t pos_;
161 InMemoryDeserializer::InMemoryDeserializer(const std::vector<char> &buffer)
162 : impl_(new Impl(buffer))
166 InMemoryDeserializer::~InMemoryDeserializer()
170 void InMemoryDeserializer::doUChar(unsigned char *value)
172 impl_->doValue(value);
175 void InMemoryDeserializer::doInt(int *value)
177 impl_->doValue(value);
180 void InMemoryDeserializer::doInt64(gmx_int64_t *value)
182 impl_->doValue(value);
185 void InMemoryDeserializer::doFloat(float *value)
187 impl_->doValue(value);
190 void InMemoryDeserializer::doDouble(double *value)
192 impl_->doValue(value);
195 void InMemoryDeserializer::doString(std::string *value)
197 impl_->doString(value);
200 } // namespace gmx