Sort all includes in src/gromacs
[alexxy/gromacs.git] / src / gromacs / mdlib / calcmu.c
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, 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 "gromacs/legacyheaders/calcmu.h"
41
42 #include <stdio.h>
43 #include <stdlib.h>
44
45 #include "gromacs/legacyheaders/gmx_omp_nthreads.h"
46 #include "gromacs/legacyheaders/network.h"
47 #include "gromacs/legacyheaders/typedefs.h"
48 #include "gromacs/math/units.h"
49 #include "gromacs/math/vec.h"
50
51 void calc_mu(int start, int homenr, rvec x[], real q[], real qB[],
52              int nChargePerturbed,
53              dvec mu, dvec mu_B)
54 {
55     int    i, end, m;
56     double mu_x, mu_y, mu_z;
57
58     end   = start + homenr;
59
60     mu_x = mu_y = mu_z = 0.0;
61 #pragma omp parallel for reduction(+: mu_x, mu_y, mu_z) schedule(static) \
62     num_threads(gmx_omp_nthreads_get(emntDefault))
63     for (i = start; i < end; i++)
64     {
65         mu_x += q[i]*x[i][XX];
66         mu_y += q[i]*x[i][YY];
67         mu_z += q[i]*x[i][ZZ];
68     }
69     mu[XX] = mu_x;
70     mu[YY] = mu_y;
71     mu[ZZ] = mu_z;
72
73     for (m = 0; (m < DIM); m++)
74     {
75         mu[m] *= ENM2DEBYE;
76     }
77
78     if (nChargePerturbed)
79     {
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++)
84         {
85             mu_x += qB[i]*x[i][XX];
86             mu_y += qB[i]*x[i][YY];
87             mu_z += qB[i]*x[i][ZZ];
88         }
89         mu_B[XX] = mu_x * ENM2DEBYE;
90         mu_B[YY] = mu_y * ENM2DEBYE;
91         mu_B[ZZ] = mu_z * ENM2DEBYE;
92     }
93     else
94     {
95         copy_dvec(mu, mu_B);
96     }
97 }
98
99 gmx_bool read_mu(FILE *fp, rvec mu, real *vol)
100 {
101     /* For backward compatibility */
102     real mmm[4];
103
104     if (fread(mmm, (size_t)(4*sizeof(real)), 1, fp) != 1)
105     {
106         return FALSE;
107     }
108
109     copy_rvec(mmm, mu);
110     *vol = mmm[3];
111
112     return TRUE;
113 }