Apply clang-format-11
[alexxy/gromacs.git] / src / gromacs / applied_forces / densityfitting / densityfittingamplitudelookup.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2019,2020,2021, 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/utility/basedefinitions.h"
50 #include "gromacs/utility/arrayref.h"
51
52 namespace gmx
53 {
54 /********************************************************************
55  * DensityFittingAmplitudeLookup::Impl
56  */
57
58 class DensityFittingAmplitudeLookupImpl
59 {
60 public:
61     DensityFittingAmplitudeLookupImpl()                                         = default;
62     DensityFittingAmplitudeLookupImpl(const DensityFittingAmplitudeLookupImpl&) = default;
63     virtual ~DensityFittingAmplitudeLookupImpl()                                = default;
64
65     virtual const std::vector<real>& operator()(ArrayRef<const real> /*chargeA*/,
66                                                 ArrayRef<const real> /*massT*/,
67                                                 ArrayRef<const int> localIndex) = 0;
68     virtual std::unique_ptr<DensityFittingAmplitudeLookupImpl> clone()          = 0;
69 };
70
71 namespace
72 {
73 class UnitAmplitudes final : public DensityFittingAmplitudeLookupImpl
74 {
75 public:
76     UnitAmplitudes()                      = default;
77     UnitAmplitudes(const UnitAmplitudes&) = default;
78     ~UnitAmplitudes() override            = default;
79     std::unique_ptr<DensityFittingAmplitudeLookupImpl> clone() override;
80     const std::vector<real>&                           operator()(ArrayRef<const real> /*chargeA*/,
81                                         ArrayRef<const real> /*massT*/,
82                                         ArrayRef<const int> localIndex) override;
83
84 private:
85     std::vector<real> amplitude_;
86 };
87
88 std::unique_ptr<DensityFittingAmplitudeLookupImpl> UnitAmplitudes::clone()
89 {
90     return std::make_unique<UnitAmplitudes>(*this);
91 };
92
93 const std::vector<real>& UnitAmplitudes::operator()(ArrayRef<const real> /*chargeA*/,
94                                                     ArrayRef<const real> /*massT*/,
95                                                     ArrayRef<const int> localIndex)
96 {
97     if (amplitude_.size() != localIndex.size())
98     {
99         amplitude_.resize(localIndex.size(), 1.);
100     }
101
102     return amplitude_;
103 }
104
105 class ChargesAsAmplitudes final : public DensityFittingAmplitudeLookupImpl
106 {
107 public:
108     ChargesAsAmplitudes()                           = default;
109     ChargesAsAmplitudes(const ChargesAsAmplitudes&) = default;
110     ~ChargesAsAmplitudes() override                 = default;
111     std::unique_ptr<DensityFittingAmplitudeLookupImpl> clone() override;
112     const std::vector<real>&                           operator()(ArrayRef<const real> chargeA,
113                                         ArrayRef<const real> /*massT*/,
114                                         ArrayRef<const int> localIndex) override;
115
116 private:
117     std::vector<real> amplitude_;
118 };
119
120 std::unique_ptr<DensityFittingAmplitudeLookupImpl> ChargesAsAmplitudes::clone()
121 {
122     return std::make_unique<ChargesAsAmplitudes>(*this);
123 };
124
125 const std::vector<real>& ChargesAsAmplitudes::operator()(ArrayRef<const real> chargeA,
126                                                          ArrayRef<const real> /*massT*/,
127                                                          ArrayRef<const int> localIndex)
128 {
129     if (amplitude_.size() != localIndex.size())
130     {
131         amplitude_.resize(localIndex.size());
132     }
133
134     std::transform(std::begin(localIndex),
135                    std::end(localIndex),
136                    std::begin(amplitude_),
137                    [&chargeA](gmx::index index) { return chargeA[index]; });
138     return amplitude_;
139 }
140
141 class MassesAsAmplitudes final : public DensityFittingAmplitudeLookupImpl
142 {
143 public:
144     MassesAsAmplitudes()                          = default;
145     MassesAsAmplitudes(const MassesAsAmplitudes&) = default;
146     ~MassesAsAmplitudes() override                = default;
147     std::unique_ptr<DensityFittingAmplitudeLookupImpl> clone() override;
148     const std::vector<real>&                           operator()(ArrayRef<const real> /*chargeA*/,
149                                         ArrayRef<const real> massT,
150                                         ArrayRef<const int>  localIndex) override;
151
152 private:
153     std::vector<real> amplitude_;
154 };
155
156 std::unique_ptr<DensityFittingAmplitudeLookupImpl> MassesAsAmplitudes::clone()
157 {
158     return std::make_unique<MassesAsAmplitudes>(*this);
159 };
160
161 const std::vector<real>& MassesAsAmplitudes::operator()(ArrayRef<const real> /*chargeA*/,
162                                                         ArrayRef<const real> massT,
163                                                         ArrayRef<const int>  localIndex)
164 {
165     if (amplitude_.size() != localIndex.size())
166     {
167         amplitude_.resize(localIndex.size());
168     }
169
170     std::transform(std::begin(localIndex),
171                    std::end(localIndex),
172                    std::begin(amplitude_),
173                    [&massT](gmx::index index) { return massT[index]; });
174     return amplitude_;
175 }
176
177 } // namespace
178
179 /********************************************************************
180  * DensityFittingAmplitudeLookup
181  */
182
183
184 DensityFittingAmplitudeLookup::DensityFittingAmplitudeLookup(const DensityFittingAmplitudeMethod& method)
185 {
186     switch (method)
187     {
188         case DensityFittingAmplitudeMethod::Unity:
189             impl_ = std::make_unique<UnitAmplitudes>();
190             break;
191         case DensityFittingAmplitudeMethod::Mass:
192             impl_ = std::make_unique<MassesAsAmplitudes>();
193             break;
194         case DensityFittingAmplitudeMethod::Charge:
195             impl_ = std::make_unique<ChargesAsAmplitudes>();
196             break;
197         default: break;
198     }
199 }
200
201 const std::vector<real>& DensityFittingAmplitudeLookup::operator()(ArrayRef<const real> chargeA,
202                                                                    ArrayRef<const real> massT,
203                                                                    ArrayRef<const int>  localIndex)
204 {
205     return (*impl_)(chargeA, massT, localIndex);
206 }
207
208 DensityFittingAmplitudeLookup::~DensityFittingAmplitudeLookup() = default;
209
210 DensityFittingAmplitudeLookup::DensityFittingAmplitudeLookup(const DensityFittingAmplitudeLookup& other) :
211     impl_(other.impl_->clone())
212 {
213 }
214
215 DensityFittingAmplitudeLookup& DensityFittingAmplitudeLookup::operator=(const DensityFittingAmplitudeLookup& other)
216 {
217     impl_ = other.impl_->clone();
218     return *this;
219 }
220
221 DensityFittingAmplitudeLookup::DensityFittingAmplitudeLookup(DensityFittingAmplitudeLookup&&) noexcept = default;
222
223 DensityFittingAmplitudeLookup&
224 DensityFittingAmplitudeLookup::operator=(DensityFittingAmplitudeLookup&&) noexcept = default;
225
226 } // namespace gmx