2 * This file is part of the GROMACS molecular simulation package.
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.
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.
37 * Implements amplitude lookup for density fitting
39 * \author Christian Blau <blau@kth.se>
40 * \ingroup module_applied_forces
44 #include "densityfittingamplitudelookup.h"
49 #include "gromacs/utility/basedefinitions.h"
50 #include "gromacs/utility/arrayref.h"
54 /********************************************************************
55 * DensityFittingAmplitudeLookup::Impl
58 class DensityFittingAmplitudeLookupImpl
61 DensityFittingAmplitudeLookupImpl() = default;
62 DensityFittingAmplitudeLookupImpl(const DensityFittingAmplitudeLookupImpl&) = default;
63 virtual ~DensityFittingAmplitudeLookupImpl() = default;
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;
73 class UnitAmplitudes final : public DensityFittingAmplitudeLookupImpl
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;
85 std::vector<real> amplitude_;
88 std::unique_ptr<DensityFittingAmplitudeLookupImpl> UnitAmplitudes::clone()
90 return std::make_unique<UnitAmplitudes>(*this);
93 const std::vector<real>& UnitAmplitudes::operator()(ArrayRef<const real> /*chargeA*/,
94 ArrayRef<const real> /*massT*/,
95 ArrayRef<const int> localIndex)
97 if (amplitude_.size() != localIndex.size())
99 amplitude_.resize(localIndex.size(), 1.);
105 class ChargesAsAmplitudes final : public DensityFittingAmplitudeLookupImpl
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;
117 std::vector<real> amplitude_;
120 std::unique_ptr<DensityFittingAmplitudeLookupImpl> ChargesAsAmplitudes::clone()
122 return std::make_unique<ChargesAsAmplitudes>(*this);
125 const std::vector<real>& ChargesAsAmplitudes::operator()(ArrayRef<const real> chargeA,
126 ArrayRef<const real> /*massT*/,
127 ArrayRef<const int> localIndex)
129 if (amplitude_.size() != localIndex.size())
131 amplitude_.resize(localIndex.size());
134 std::transform(std::begin(localIndex),
135 std::end(localIndex),
136 std::begin(amplitude_),
137 [&chargeA](gmx::index index) { return chargeA[index]; });
141 class MassesAsAmplitudes final : public DensityFittingAmplitudeLookupImpl
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;
153 std::vector<real> amplitude_;
156 std::unique_ptr<DensityFittingAmplitudeLookupImpl> MassesAsAmplitudes::clone()
158 return std::make_unique<MassesAsAmplitudes>(*this);
161 const std::vector<real>& MassesAsAmplitudes::operator()(ArrayRef<const real> /*chargeA*/,
162 ArrayRef<const real> massT,
163 ArrayRef<const int> localIndex)
165 if (amplitude_.size() != localIndex.size())
167 amplitude_.resize(localIndex.size());
170 std::transform(std::begin(localIndex),
171 std::end(localIndex),
172 std::begin(amplitude_),
173 [&massT](gmx::index index) { return massT[index]; });
179 /********************************************************************
180 * DensityFittingAmplitudeLookup
184 DensityFittingAmplitudeLookup::DensityFittingAmplitudeLookup(const DensityFittingAmplitudeMethod& method)
188 case DensityFittingAmplitudeMethod::Unity:
189 impl_ = std::make_unique<UnitAmplitudes>();
191 case DensityFittingAmplitudeMethod::Mass:
192 impl_ = std::make_unique<MassesAsAmplitudes>();
194 case DensityFittingAmplitudeMethod::Charge:
195 impl_ = std::make_unique<ChargesAsAmplitudes>();
201 const std::vector<real>& DensityFittingAmplitudeLookup::operator()(ArrayRef<const real> chargeA,
202 ArrayRef<const real> massT,
203 ArrayRef<const int> localIndex)
205 return (*impl_)(chargeA, massT, localIndex);
208 DensityFittingAmplitudeLookup::~DensityFittingAmplitudeLookup() = default;
210 DensityFittingAmplitudeLookup::DensityFittingAmplitudeLookup(const DensityFittingAmplitudeLookup& other) :
211 impl_(other.impl_->clone())
215 DensityFittingAmplitudeLookup& DensityFittingAmplitudeLookup::operator=(const DensityFittingAmplitudeLookup& other)
217 impl_ = other.impl_->clone();
221 DensityFittingAmplitudeLookup::DensityFittingAmplitudeLookup(DensityFittingAmplitudeLookup&&) noexcept = default;
223 DensityFittingAmplitudeLookup&
224 DensityFittingAmplitudeLookup::operator=(DensityFittingAmplitudeLookup&&) noexcept = default;