4fe1eaa10450d39d36c945e49b81005a9bc873db
[alexxy/gromacs.git] / src / gromacs / mdlib / calcmu.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2012,2014,2015,2017,2018,2019, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 /* This file is completely threadsafe - keep it that way! */
38 #include "gmxpre.h"
39
40 #include "calcmu.h"
41
42 #include <cstdio>
43 #include <cstdlib>
44
45 #include "gromacs/gmxlib/network.h"
46 #include "gromacs/math/units.h"
47 #include "gromacs/math/vec.h"
48 #include "gromacs/mdlib/gmx_omp_nthreads.h"
49
50 void calc_mu(int                      start,
51              int                      homenr,
52              gmx::ArrayRef<gmx::RVec> x,
53              const real               q[],
54              const real               qB[],
55              int                      nChargePerturbed,
56              dvec                     mu,
57              dvec                     mu_B)
58 {
59     int    end, m;
60     double mu_x, mu_y, mu_z;
61
62     end = start + homenr;
63
64     mu_x = mu_y = mu_z = 0.0;
65 #pragma omp parallel for reduction(+: mu_x, mu_y, mu_z) schedule(static) \
66     num_threads(gmx_omp_nthreads_get(emntDefault))
67     for (int i = start; i < end; i++)
68     {
69         // Trivial OpenMP region that cannot throw
70         mu_x += q[i] * x[i][XX];
71         mu_y += q[i] * x[i][YY];
72         mu_z += q[i] * x[i][ZZ];
73     }
74     mu[XX] = mu_x;
75     mu[YY] = mu_y;
76     mu[ZZ] = mu_z;
77
78     for (m = 0; (m < DIM); m++)
79     {
80         mu[m] *= ENM2DEBYE;
81     }
82
83     if (nChargePerturbed)
84     {
85         mu_x = mu_y = mu_z = 0.0;
86 #pragma omp parallel for reduction(+: mu_x, mu_y, mu_z) schedule(static) \
87         num_threads(gmx_omp_nthreads_get(emntDefault))
88         for (int i = start; i < end; i++)
89         {
90             // Trivial OpenMP region that cannot throw
91             mu_x += qB[i] * x[i][XX];
92             mu_y += qB[i] * x[i][YY];
93             mu_z += qB[i] * x[i][ZZ];
94         }
95         mu_B[XX] = mu_x * ENM2DEBYE;
96         mu_B[YY] = mu_y * ENM2DEBYE;
97         mu_B[ZZ] = mu_z * ENM2DEBYE;
98     }
99     else
100     {
101         copy_dvec(mu, mu_B);
102     }
103 }
104
105 gmx_bool read_mu(FILE* fp, rvec mu, real* vol)
106 {
107     /* For backward compatibility */
108     real mmm[4];
109
110     if (fread(mmm, static_cast<size_t>(4 * sizeof(real)), 1, fp) != 1)
111     {
112         return FALSE;
113     }
114
115     copy_rvec(mmm, mu);
116     *vol = mmm[3];
117
118     return TRUE;
119 }