Consistently change densityfitting spreading weight to charge
[gromacs.git] / src / gromacs / applied_forces / densityfittingamplitudelookup.cpp
blob04c56616adc317bc298fc853c6305d0ccecc31fe
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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.
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 * Implements amplitude lookup for density fitting
39 * \author Christian Blau <blau@kth.se>
40 * \ingroup module_applied_forces
42 #include "gmxpre.h"
44 #include "densityfittingamplitudelookup.h"
46 #include <algorithm>
47 #include <vector>
49 #include "gromacs/mdtypes/mdatom.h"
51 namespace gmx
53 /********************************************************************
54 * DensityFittingAmplitudeLookup::Impl
57 class DensityFittingAmplitudeLookupImpl
59 public:
60 DensityFittingAmplitudeLookupImpl() = default;
61 DensityFittingAmplitudeLookupImpl(const DensityFittingAmplitudeLookupImpl &) = default;
62 virtual ~DensityFittingAmplitudeLookupImpl() = default;
64 virtual const std::vector<real> &operator()(const t_mdatoms &atoms,
65 ArrayRef<const int> localIndex) = 0;
66 virtual std::unique_ptr<DensityFittingAmplitudeLookupImpl> clone() = 0;
69 namespace
71 class UnitAmplitudes final : public DensityFittingAmplitudeLookupImpl
73 public:
74 UnitAmplitudes() = default;
75 UnitAmplitudes(const UnitAmplitudes &) = default;
76 ~UnitAmplitudes() override = default;
77 std::unique_ptr<DensityFittingAmplitudeLookupImpl> clone() override;
78 const std::vector<real> &operator()(const t_mdatoms &atoms,
79 ArrayRef<const int> localIndex) override;
81 private:
82 std::vector<real> amplitude_;
85 std::unique_ptr<DensityFittingAmplitudeLookupImpl> UnitAmplitudes::clone()
87 return std::make_unique<UnitAmplitudes>(*this);
90 const std::vector<real> &UnitAmplitudes::operator()(const t_mdatoms & /*atoms*/,
91 ArrayRef<const int> localIndex)
93 if (amplitude_.size() != localIndex.size())
95 amplitude_.resize(localIndex.size(), 1.);
98 return amplitude_;
101 class ChargesAsAmplitudes final : public DensityFittingAmplitudeLookupImpl
103 public:
104 ChargesAsAmplitudes() = default;
105 ChargesAsAmplitudes(const ChargesAsAmplitudes &) = default;
106 ~ChargesAsAmplitudes() override = default;
107 std::unique_ptr<DensityFittingAmplitudeLookupImpl> clone() override;
108 const std::vector<real> &operator()(const t_mdatoms &atoms,
109 ArrayRef<const int> localIndex) override;
111 private:
112 std::vector<real> amplitude_;
115 std::unique_ptr<DensityFittingAmplitudeLookupImpl> ChargesAsAmplitudes::clone()
117 return std::make_unique<ChargesAsAmplitudes>(*this);
120 const std::vector<real> &ChargesAsAmplitudes::operator()(const t_mdatoms &atoms,
121 ArrayRef<const int> localIndex)
123 if (amplitude_.size() != localIndex.size())
125 amplitude_.resize(localIndex.size());
128 std::transform(std::begin(localIndex), std::end(localIndex), std::begin(amplitude_),
129 [&atoms](gmx::index index){return atoms.chargeA[index]; }
131 return amplitude_;
134 class MassesAsAmplitudes final : public DensityFittingAmplitudeLookupImpl
136 public:
137 MassesAsAmplitudes() = default;
138 MassesAsAmplitudes(const MassesAsAmplitudes &) = default;
139 ~MassesAsAmplitudes() override = default;
140 std::unique_ptr<DensityFittingAmplitudeLookupImpl> clone() override;
141 const std::vector<real> &operator()(const t_mdatoms &atoms,
142 ArrayRef<const int> localIndex) override;
144 private:
145 std::vector<real> amplitude_;
148 std::unique_ptr<DensityFittingAmplitudeLookupImpl> MassesAsAmplitudes::clone()
150 return std::make_unique<MassesAsAmplitudes>(*this);
153 const std::vector<real> &MassesAsAmplitudes::operator()(const t_mdatoms &atoms,
154 ArrayRef<const int> localIndex)
156 if (amplitude_.size() != localIndex.size())
158 amplitude_.resize(localIndex.size());
161 std::transform(std::begin(localIndex), std::end(localIndex), std::begin(amplitude_),
162 [&atoms](gmx::index index) { return atoms.massT[index]; });
163 return amplitude_;
166 } // namespace
168 /********************************************************************
169 * DensityFittingAmplitudeLookup
173 DensityFittingAmplitudeLookup::DensityFittingAmplitudeLookup(const DensityFittingAmplitudeMethod &method)
175 switch (method)
177 case DensityFittingAmplitudeMethod::Unity:
178 impl_ = std::make_unique<UnitAmplitudes>();
179 break;
180 case DensityFittingAmplitudeMethod::Mass:
181 impl_ = std::make_unique<MassesAsAmplitudes>();
182 break;
183 case DensityFittingAmplitudeMethod::Charge:
184 impl_ = std::make_unique<ChargesAsAmplitudes>();
185 break;
186 default:
187 break;
191 const std::vector<real> &DensityFittingAmplitudeLookup::operator()(const t_mdatoms &atoms,
192 ArrayRef<const int> localIndex)
194 return (*impl_)(atoms, localIndex);
197 DensityFittingAmplitudeLookup::~DensityFittingAmplitudeLookup() = default;
199 DensityFittingAmplitudeLookup::DensityFittingAmplitudeLookup(const DensityFittingAmplitudeLookup &other)
200 : impl_(other.impl_->clone())
204 DensityFittingAmplitudeLookup &DensityFittingAmplitudeLookup::operator=(const DensityFittingAmplitudeLookup &other)
206 impl_ = other.impl_->clone();
207 return *this;
210 DensityFittingAmplitudeLookup::DensityFittingAmplitudeLookup(DensityFittingAmplitudeLookup &&) noexcept = default;
212 DensityFittingAmplitudeLookup &DensityFittingAmplitudeLookup::operator=(DensityFittingAmplitudeLookup &&) noexcept = default;
214 } // namespace gmx