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 /*! \libinternal \file
36  * \brief
37  * Declares classes to calculate density fitting forces.
38  *
39  * \author Christian Blau <blau@kth.se>
40  * \inlibraryapi
41  * \ingroup module_math
42  */
43 
44 #include "gromacs/math/gausstransform.h"
45 #include "gromacs/math/vectypes.h"
46 #include "gromacs/mdspan/extensions.h"
47 #include "gromacs/mdspan/mdspan.h"
48 #include "gromacs/utility/classhelpers.h"
49 
50 namespace gmx
51 {
52 
53 
54 /*! \libinternal
55  * \brief
56  * Manages evaluation of density-fitting forces for particles that were spread
57  * with a kernel.
58  *
59  * The density fitting force is the inner product of the derivative of the
60  * density-density goodness-of-fit function at all density values and the
61  * gradient of the spreading kernel used to simulate the density.
62  *
63  * For a discrete Gauss transform spreading kernel and density derivative \f$D\f$
64  * it is the Gauss transform of a single coordinate, weighted by a distance vector
65  * and density derivative at each lattice point
66  * \f[
67  *  \sum_v D_v \frac{\vec{x_v} - \vec{x}}{\sigma^2}
68  *  \frac{1}{\sigma^3 \sqrt(2^3*\pi^3)} * \exp(-\frac{(x_v-x)^2}{2 \sigma^2})
69  * \f]
70  *
71  */
72 class DensityFittingForce
73 {
74 public:
75     /*! \brief Construct density fitting force evaluator.
76      * \param[in] kernelShapeParameters the global parameters of the density spreading kernel
77      */
78     DensityFittingForce(const GaussianSpreadKernelParameters::Shape& kernelShapeParameters);
79 
80     ~DensityFittingForce();
81     //! Copy constructor
82     DensityFittingForce(const DensityFittingForce& other);
83     //! Copy assignment
84     DensityFittingForce& operator=(const DensityFittingForce& other);
85     //! Move constructor
86     DensityFittingForce(DensityFittingForce&& other) noexcept;
87     //! Move assignment
88     DensityFittingForce& operator=(DensityFittingForce&& other) noexcept;
89     /*! \brief
90      * Evaluate density-fitting force between a reference density and the
91      * density generated by a kernel.
92      *
93      * Implements the equation in the class description.
94      *
95      * \param[in] localParameters   local parameters of the spreading kernel
96      * \param[in] densityDerivative the spatial derivative of the similarity
97      *                              measure between the reference density and
98      *                              a density that was generated using the
99      *                              spreading kernel
100      *
101      * \returns the force that increases the measure of the goodness of fit used to calculate the density derivative
102      */
103     RVec evaluateForce(const GaussianSpreadKernelParameters::PositionAndAmplitude& localParameters,
104                        basic_mdspan<const float, dynamicExtents3D> densityDerivative);
105 
106 private:
107     class Impl;
108     PrivateImplPointer<Impl> impl_;
109 };
110 
111 } // namespace gmx
112