Move main.*, splitter.*, gmx_omp_nthreads.* to mdlib
[gromacs.git] / src / gromacs / mdlib / calcmu.cpp
blobc29e37005262a5c1ca1b0a8413865e358d9a319c
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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, 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.
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.
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.
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.
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.
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.
37 /* This file is completely threadsafe - keep it that way! */
38 #include "gmxpre.h"
40 #include "calcmu.h"
42 #include <stdio.h>
43 #include <stdlib.h>
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"
50 void calc_mu(int start, int homenr, rvec x[], real q[], real qB[],
51 int nChargePerturbed,
52 dvec mu, dvec mu_B)
54 int i, end, m;
55 double mu_x, mu_y, mu_z;
57 end = start + homenr;
59 mu_x = mu_y = mu_z = 0.0;
60 #pragma omp parallel for reduction(+: mu_x, mu_y, mu_z) schedule(static) \
61 num_threads(gmx_omp_nthreads_get(emntDefault))
62 for (i = start; i < end; i++)
64 // Trivial OpenMP region that cannot throw
65 mu_x += q[i]*x[i][XX];
66 mu_y += q[i]*x[i][YY];
67 mu_z += q[i]*x[i][ZZ];
69 mu[XX] = mu_x;
70 mu[YY] = mu_y;
71 mu[ZZ] = mu_z;
73 for (m = 0; (m < DIM); m++)
75 mu[m] *= ENM2DEBYE;
78 if (nChargePerturbed)
80 mu_x = mu_y = mu_z = 0.0;
81 #pragma omp parallel for reduction(+: mu_x, mu_y, mu_z) schedule(static) \
82 num_threads(gmx_omp_nthreads_get(emntDefault))
83 for (i = start; i < end; i++)
85 // Trivial OpenMP region that cannot throw
86 mu_x += qB[i]*x[i][XX];
87 mu_y += qB[i]*x[i][YY];
88 mu_z += qB[i]*x[i][ZZ];
90 mu_B[XX] = mu_x * ENM2DEBYE;
91 mu_B[YY] = mu_y * ENM2DEBYE;
92 mu_B[ZZ] = mu_z * ENM2DEBYE;
94 else
96 copy_dvec(mu, mu_B);
100 gmx_bool read_mu(FILE *fp, rvec mu, real *vol)
102 /* For backward compatibility */
103 real mmm[4];
105 if (fread(mmm, (size_t)(4*sizeof(real)), 1, fp) != 1)
107 return FALSE;
110 copy_rvec(mmm, mu);
111 *vol = mmm[3];
113 return TRUE;