2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2019,2020, 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/mdtypes/mdatom.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()(const t_mdatoms& atoms, ArrayRef<const int> localIndex) = 0;
66 virtual std::unique_ptr<DensityFittingAmplitudeLookupImpl> clone() = 0;
71 class UnitAmplitudes final : public DensityFittingAmplitudeLookupImpl
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, ArrayRef<const int> localIndex) override;
81 std::vector<real> amplitude_;
84 std::unique_ptr<DensityFittingAmplitudeLookupImpl> UnitAmplitudes::clone()
86 return std::make_unique<UnitAmplitudes>(*this);
89 const std::vector<real>& UnitAmplitudes::operator()(const t_mdatoms& /*atoms*/, ArrayRef<const int> localIndex)
91 if (amplitude_.size() != localIndex.size())
93 amplitude_.resize(localIndex.size(), 1.);
99 class ChargesAsAmplitudes final : public DensityFittingAmplitudeLookupImpl
102 ChargesAsAmplitudes() = default;
103 ChargesAsAmplitudes(const ChargesAsAmplitudes&) = default;
104 ~ChargesAsAmplitudes() override = default;
105 std::unique_ptr<DensityFittingAmplitudeLookupImpl> clone() override;
106 const std::vector<real>& operator()(const t_mdatoms& atoms, ArrayRef<const int> localIndex) override;
109 std::vector<real> amplitude_;
112 std::unique_ptr<DensityFittingAmplitudeLookupImpl> ChargesAsAmplitudes::clone()
114 return std::make_unique<ChargesAsAmplitudes>(*this);
117 const std::vector<real>& ChargesAsAmplitudes::operator()(const t_mdatoms& atoms, ArrayRef<const int> localIndex)
119 if (amplitude_.size() != localIndex.size())
121 amplitude_.resize(localIndex.size());
124 std::transform(std::begin(localIndex), std::end(localIndex), std::begin(amplitude_),
125 [&atoms](gmx::index index) { return atoms.chargeA[index]; });
129 class MassesAsAmplitudes final : public DensityFittingAmplitudeLookupImpl
132 MassesAsAmplitudes() = default;
133 MassesAsAmplitudes(const MassesAsAmplitudes&) = default;
134 ~MassesAsAmplitudes() override = default;
135 std::unique_ptr<DensityFittingAmplitudeLookupImpl> clone() override;
136 const std::vector<real>& operator()(const t_mdatoms& atoms, ArrayRef<const int> localIndex) override;
139 std::vector<real> amplitude_;
142 std::unique_ptr<DensityFittingAmplitudeLookupImpl> MassesAsAmplitudes::clone()
144 return std::make_unique<MassesAsAmplitudes>(*this);
147 const std::vector<real>& MassesAsAmplitudes::operator()(const t_mdatoms& atoms, ArrayRef<const int> localIndex)
149 if (amplitude_.size() != localIndex.size())
151 amplitude_.resize(localIndex.size());
154 std::transform(std::begin(localIndex), std::end(localIndex), std::begin(amplitude_),
155 [&atoms](gmx::index index) { return atoms.massT[index]; });
161 /********************************************************************
162 * DensityFittingAmplitudeLookup
166 DensityFittingAmplitudeLookup::DensityFittingAmplitudeLookup(const DensityFittingAmplitudeMethod& method)
170 case DensityFittingAmplitudeMethod::Unity:
171 impl_ = std::make_unique<UnitAmplitudes>();
173 case DensityFittingAmplitudeMethod::Mass:
174 impl_ = std::make_unique<MassesAsAmplitudes>();
176 case DensityFittingAmplitudeMethod::Charge:
177 impl_ = std::make_unique<ChargesAsAmplitudes>();
183 const std::vector<real>& DensityFittingAmplitudeLookup::operator()(const t_mdatoms& atoms,
184 ArrayRef<const int> localIndex)
186 return (*impl_)(atoms, localIndex);
189 DensityFittingAmplitudeLookup::~DensityFittingAmplitudeLookup() = default;
191 DensityFittingAmplitudeLookup::DensityFittingAmplitudeLookup(const DensityFittingAmplitudeLookup& other) :
192 impl_(other.impl_->clone())
196 DensityFittingAmplitudeLookup& DensityFittingAmplitudeLookup::operator=(const DensityFittingAmplitudeLookup& other)
198 impl_ = other.impl_->clone();
202 DensityFittingAmplitudeLookup::DensityFittingAmplitudeLookup(DensityFittingAmplitudeLookup&&) noexcept = default;
204 DensityFittingAmplitudeLookup& DensityFittingAmplitudeLookup::
205 operator=(DensityFittingAmplitudeLookup&&) noexcept = default;