5ded958e743d7edd5bef5af185ee3c53fca0d0dc
[alexxy/gromacs.git] / api / nblib / listed_forces / calculator.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 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 a bonded force calculator
38  *
39  * \author Victor Holanda <victor.holanda@cscs.ch>
40  * \author Joe Jordan <ejjordan@kth.se>
41  * \author Prashanth Kanduri <kanduri@cscs.ch>
42  * \author Sebastian Keller <keller@cscs.ch>
43  * \author Artem Zhmurov <zhmurov@gmail.com>
44  */
45 #include "nblib/box.h"
46 #include "nblib/exception.h"
47 #include "nblib/pbc.hpp"
48 #include "nblib/listed_forces/calculator.h"
49 #include "nblib/listed_forces/dataflow.hpp"
50 #include "nblib/listed_forces/helpers.hpp"
51
52 namespace nblib
53 {
54
55 ListedForceCalculator::~ListedForceCalculator() = default;
56
57 ListedForceCalculator::ListedForceCalculator(const ListedInteractionData& interactions,
58                                              size_t                       bufferSize,
59                                              int                          nthr,
60                                              const Box&                   box) :
61     numThreads(nthr),
62     masterForceBuffer_(bufferSize, Vec3{ 0, 0, 0 }),
63     pbcHolder_(std::make_unique<PbcHolder>(PbcType::Xyz, box))
64 {
65     // split up the work
66     threadedInteractions_ = splitListedWork(interactions, bufferSize, numThreads);
67
68     // set up the reduction buffers
69     int threadRange = bufferSize / numThreads;
70     for (int i = 0; i < numThreads; ++i)
71     {
72         int rangeStart = i * threadRange;
73         int rangeEnd   = rangeStart + threadRange;
74         // the last range is a bit bigger due to integer division
75         if (i == numThreads - 1)
76         {
77             rangeEnd = bufferSize;
78         }
79
80         threadedForceBuffers_.push_back(std::make_unique<ForceBuffer<Vec3>>(
81                 masterForceBuffer_.data(), rangeStart, rangeEnd));
82     }
83 }
84
85 void ListedForceCalculator::computeForcesAndEnergies(gmx::ArrayRef<const Vec3> x, bool usePbc)
86 {
87     energyBuffer_.fill(0);
88     std::vector<std::array<real, std::tuple_size<ListedInteractionData>::value>> energiesPerThread(numThreads);
89
90 #pragma omp parallel for num_threads(numThreads) schedule(static)
91     for (int thread = 0; thread < numThreads; ++thread)
92     {
93         if (usePbc)
94         {
95             energiesPerThread[thread] = reduceListedForces(
96                     threadedInteractions_[thread], x, threadedForceBuffers_[thread].get(), *pbcHolder_);
97         }
98         else
99         {
100             energiesPerThread[thread] = reduceListedForces(
101                     threadedInteractions_[thread], x, threadedForceBuffers_[thread].get(), NoPbc{});
102         }
103     }
104
105     // reduce energies
106     for (int thread = 0; thread < numThreads; ++thread)
107     {
108         for (int type = 0; type < int(energyBuffer_.size()); ++type)
109         {
110             energyBuffer_[type] += energiesPerThread[thread][type];
111         }
112     }
113
114     // reduce forces
115 #pragma omp parallel for num_threads(numThreads) schedule(static)
116     for (int thread = 0; thread < numThreads; ++thread)
117     {
118         auto& thisBuffer = *threadedForceBuffers_[thread];
119         // access outliers from other threads
120         for (int otherThread = 0; otherThread < numThreads; ++otherThread)
121         {
122             auto& otherBuffer = *threadedForceBuffers_[otherThread];
123             for (const auto& outlier : otherBuffer)
124             {
125                 int index = outlier.first;
126                 // check whether this outlier falls within the range of <thread>
127                 if (thisBuffer.inRange(index))
128                 {
129                     auto force = outlier.second;
130                     masterForceBuffer_[index] += force;
131                 }
132             }
133         }
134     }
135 }
136
137 void ListedForceCalculator::compute(gmx::ArrayRef<const Vec3> coordinates, gmx::ArrayRef<Vec3> forces, bool usePbc)
138 {
139     if (coordinates.size() != forces.size())
140     {
141         throw InputException("Coordinates array and force buffer size mismatch");
142     }
143
144     // check if the force buffers have the same size
145     if (masterForceBuffer_.size() != forces.size())
146     {
147         throw InputException("Input force buffer size mismatch with listed forces buffer");
148     }
149
150     // compute forces and fill in local buffers
151     computeForcesAndEnergies(coordinates, usePbc);
152
153     // add forces to output force buffers
154     for (int pIndex = 0; pIndex < int(forces.size()); pIndex++)
155     {
156         forces[pIndex] += masterForceBuffer_[pIndex];
157     }
158 }
159 void ListedForceCalculator::compute(gmx::ArrayRef<const Vec3> coordinates,
160                                     gmx::ArrayRef<Vec3>       forces,
161                                     EnergyType&               energies,
162                                     bool                      usePbc)
163 {
164     compute(coordinates, forces, usePbc);
165
166     energies = energyBuffer_;
167 }
168
169 } // namespace nblib