8e708c737d7be5e818ba77b8ab25f461375ab014
[alexxy/gromacs.git] / src / gromacs / mdtypes / group.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2019,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 classes from group.h.
38  *
39  * \author Kevin Boyd <kevin.boyd@uconn.edu>
40  * \ingroup module_mdtypes
41  */
42 #include "gmxpre.h"
43
44 #include "group.h"
45
46 #include "gromacs/utility/exceptions.h"
47
48 gmx_ekindata_t::gmx_ekindata_t(int numTempCoupleGroups, real cos_accel, int numThreads) :
49     ngtc(numTempCoupleGroups),
50     nthreads_(numThreads)
51 {
52     tcstat.resize(ngtc);
53     /* Set Berendsen tcoupl lambda's to 1,
54      * so runs without Berendsen coupling are not affected.
55      */
56     for (int i = 0; i < ngtc; i++)
57     {
58         tcstat[i].lambda         = 1.0;
59         tcstat[i].vscale_nhc     = 1.0;
60         tcstat[i].ekinscaleh_nhc = 1.0;
61         tcstat[i].ekinscalef_nhc = 1.0;
62     }
63
64     snew(ekin_work_alloc, nthreads_);
65     snew(ekin_work, nthreads_);
66     snew(dekindl_work, nthreads_);
67
68 #pragma omp parallel for num_threads(nthreads_) schedule(static)
69     for (int thread = 0; thread < nthreads_; thread++)
70     {
71         try
72         {
73             constexpr int EKIN_WORK_BUFFER_SIZE = 2;
74             /* Allocate 2 extra elements on both sides, so in single
75              * precision we have
76              * EKIN_WORK_BUFFER_SIZE*DIM*DIM*sizeof(real) = 72/144 bytes
77              * buffer on both sides to avoid cache pollution.
78              */
79             snew(ekin_work_alloc[thread], ngtc + 2 * EKIN_WORK_BUFFER_SIZE);
80             ekin_work[thread] = ekin_work_alloc[thread] + EKIN_WORK_BUFFER_SIZE;
81             /* Nasty hack so we can have the per-thread accumulation
82              * variable for dekindl in the same thread-local cache lines
83              * as the per-thread accumulation tensors for ekin[fh],
84              * because they are accumulated in the same loop. */
85             dekindl_work[thread] = &(ekin_work[thread][ngtc][0][0]);
86         }
87         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
88     }
89
90     cosacc.cos_accel = cos_accel;
91 }
92
93 gmx_ekindata_t::~gmx_ekindata_t()
94 {
95     for (int i = 0; i < nthreads_; i++)
96     {
97         sfree(ekin_work_alloc[i]);
98     }
99     sfree(ekin_work_alloc);
100     sfree(ekin_work);
101     sfree(dekindl_work);
102 }