Apply clang-format to source tree
[alexxy/gromacs.git] / src / gromacs / applied_forces / densityfittingamplitudelookup.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 /*! \internal \file
36  * \brief
37  * Implements amplitude lookup for density fitting
38  *
39  * \author Christian Blau <blau@kth.se>
40  * \ingroup module_applied_forces
41  */
42 #include "gmxpre.h"
43
44 #include "densityfittingamplitudelookup.h"
45
46 #include <algorithm>
47 #include <vector>
48
49 #include "gromacs/mdtypes/mdatom.h"
50
51 namespace gmx
52 {
53 /********************************************************************
54  * DensityFittingAmplitudeLookup::Impl
55  */
56
57 class DensityFittingAmplitudeLookupImpl
58 {
59 public:
60     DensityFittingAmplitudeLookupImpl()                                         = default;
61     DensityFittingAmplitudeLookupImpl(const DensityFittingAmplitudeLookupImpl&) = default;
62     virtual ~DensityFittingAmplitudeLookupImpl()                                = default;
63
64     virtual const std::vector<real>& operator()(const t_mdatoms& atoms, ArrayRef<const int> localIndex) = 0;
65     virtual std::unique_ptr<DensityFittingAmplitudeLookupImpl> clone() = 0;
66 };
67
68 namespace
69 {
70 class UnitAmplitudes final : public DensityFittingAmplitudeLookupImpl
71 {
72 public:
73     UnitAmplitudes()                      = default;
74     UnitAmplitudes(const UnitAmplitudes&) = default;
75     ~UnitAmplitudes() override            = default;
76     std::unique_ptr<DensityFittingAmplitudeLookupImpl> clone() override;
77     const std::vector<real>& operator()(const t_mdatoms& atoms, ArrayRef<const int> localIndex) override;
78
79 private:
80     std::vector<real> amplitude_;
81 };
82
83 std::unique_ptr<DensityFittingAmplitudeLookupImpl> UnitAmplitudes::clone()
84 {
85     return std::make_unique<UnitAmplitudes>(*this);
86 };
87
88 const std::vector<real>& UnitAmplitudes::operator()(const t_mdatoms& /*atoms*/, ArrayRef<const int> localIndex)
89 {
90     if (amplitude_.size() != localIndex.size())
91     {
92         amplitude_.resize(localIndex.size(), 1.);
93     }
94
95     return amplitude_;
96 }
97
98 class ChargesAsAmplitudes final : public DensityFittingAmplitudeLookupImpl
99 {
100 public:
101     ChargesAsAmplitudes()                           = default;
102     ChargesAsAmplitudes(const ChargesAsAmplitudes&) = default;
103     ~ChargesAsAmplitudes() override                 = default;
104     std::unique_ptr<DensityFittingAmplitudeLookupImpl> clone() override;
105     const std::vector<real>& operator()(const t_mdatoms& atoms, ArrayRef<const int> localIndex) override;
106
107 private:
108     std::vector<real> amplitude_;
109 };
110
111 std::unique_ptr<DensityFittingAmplitudeLookupImpl> ChargesAsAmplitudes::clone()
112 {
113     return std::make_unique<ChargesAsAmplitudes>(*this);
114 };
115
116 const std::vector<real>& ChargesAsAmplitudes::operator()(const t_mdatoms& atoms, ArrayRef<const int> localIndex)
117 {
118     if (amplitude_.size() != localIndex.size())
119     {
120         amplitude_.resize(localIndex.size());
121     }
122
123     std::transform(std::begin(localIndex), std::end(localIndex), std::begin(amplitude_),
124                    [&atoms](gmx::index index) { return atoms.chargeA[index]; });
125     return amplitude_;
126 }
127
128 class MassesAsAmplitudes final : public DensityFittingAmplitudeLookupImpl
129 {
130 public:
131     MassesAsAmplitudes()                          = default;
132     MassesAsAmplitudes(const MassesAsAmplitudes&) = default;
133     ~MassesAsAmplitudes() override                = default;
134     std::unique_ptr<DensityFittingAmplitudeLookupImpl> clone() override;
135     const std::vector<real>& operator()(const t_mdatoms& atoms, ArrayRef<const int> localIndex) override;
136
137 private:
138     std::vector<real> amplitude_;
139 };
140
141 std::unique_ptr<DensityFittingAmplitudeLookupImpl> MassesAsAmplitudes::clone()
142 {
143     return std::make_unique<MassesAsAmplitudes>(*this);
144 };
145
146 const std::vector<real>& MassesAsAmplitudes::operator()(const t_mdatoms& atoms, ArrayRef<const int> localIndex)
147 {
148     if (amplitude_.size() != localIndex.size())
149     {
150         amplitude_.resize(localIndex.size());
151     }
152
153     std::transform(std::begin(localIndex), std::end(localIndex), std::begin(amplitude_),
154                    [&atoms](gmx::index index) { return atoms.massT[index]; });
155     return amplitude_;
156 }
157
158 } // namespace
159
160 /********************************************************************
161  * DensityFittingAmplitudeLookup
162  */
163
164
165 DensityFittingAmplitudeLookup::DensityFittingAmplitudeLookup(const DensityFittingAmplitudeMethod& method)
166 {
167     switch (method)
168     {
169         case DensityFittingAmplitudeMethod::Unity:
170             impl_ = std::make_unique<UnitAmplitudes>();
171             break;
172         case DensityFittingAmplitudeMethod::Mass:
173             impl_ = std::make_unique<MassesAsAmplitudes>();
174             break;
175         case DensityFittingAmplitudeMethod::Charge:
176             impl_ = std::make_unique<ChargesAsAmplitudes>();
177             break;
178         default: break;
179     }
180 }
181
182 const std::vector<real>& DensityFittingAmplitudeLookup::operator()(const t_mdatoms&    atoms,
183                                                                    ArrayRef<const int> localIndex)
184 {
185     return (*impl_)(atoms, localIndex);
186 }
187
188 DensityFittingAmplitudeLookup::~DensityFittingAmplitudeLookup() = default;
189
190 DensityFittingAmplitudeLookup::DensityFittingAmplitudeLookup(const DensityFittingAmplitudeLookup& other) :
191     impl_(other.impl_->clone())
192 {
193 }
194
195 DensityFittingAmplitudeLookup& DensityFittingAmplitudeLookup::operator=(const DensityFittingAmplitudeLookup& other)
196 {
197     impl_ = other.impl_->clone();
198     return *this;
199 }
200
201 DensityFittingAmplitudeLookup::DensityFittingAmplitudeLookup(DensityFittingAmplitudeLookup&&) noexcept = default;
202
203 DensityFittingAmplitudeLookup& DensityFittingAmplitudeLookup::
204                                operator=(DensityFittingAmplitudeLookup&&) noexcept = default;
205
206 } // namespace gmx