Add nblib virials function with test
[alexxy/gromacs.git] / api / nblib / virials.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 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 /*! \inpublicapi \file
36  * \brief
37  * Implements functionality to compute virials from force data
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  */
44
45 #include "gromacs/mdlib/calcvir.h"
46 #include "gromacs/utility/arrayref.h"
47 #include "nblib/box.h"
48 #include "nblib/exception.h"
49 #include "nblib/virials.h"
50
51 namespace nblib
52 {
53
54 void computeVirialTensor(gmx::ArrayRef<const Vec3> coordinates,
55                          gmx::ArrayRef<const Vec3> forces,
56                          gmx::ArrayRef<const Vec3> shiftVectors,
57                          gmx::ArrayRef<const Vec3> shiftForces,
58                          const Box&                box,
59                          gmx::ArrayRef<real>       virialOutput)
60 {
61     if (virialOutput.size() != 9)
62     {
63         throw InputException("Virial array size incorrect, should be 9");
64     }
65
66     // set virial output array to zero
67     std::fill(virialOutput.begin(), virialOutput.end(), 0.0);
68     // use legacy tensor format
69     rvec* virial = reinterpret_cast<rvec*>(virialOutput.data());
70
71     // compute the virials from surrounding boxes
72     const rvec* fshift    = as_rvec_array(shiftForces.data());
73     const rvec* shift_vec = as_rvec_array(shiftVectors.data());
74     calc_vir(numShiftVectors, shift_vec, fshift, virial, false, box.legacyMatrix());
75
76     // calculate partial virial, for local atoms only, based on short range
77     // NOTE: GROMACS uses a variable called mdatoms.homenr which is basically number of atoms in current processor
78     //       Used numAtoms from the coordinate array in its place
79     auto        numAtoms = coordinates.size();
80     const rvec* f        = as_rvec_array(forces.data());
81     const rvec* x        = as_rvec_array(coordinates.data());
82     calc_vir(numAtoms, x, f, virial, false, box.legacyMatrix());
83 }
84
85 } // namespace nblib