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