added Verlet scheme and NxN non-bonded functionality
[gromacs.git] / src / mdlib / calcmu.c
blobc3de19f57125d26546557e43b613e243792d7d55
1 /*
2 *
3 * This source code is part of
4 *
5 * G R O M A C S
6 *
7 * GROningen MAchine for Chemical Simulations
8 *
9 * VERSION 3.2.0
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
32 * And Hey:
33 * GROwing Monsters And Cloning Shrimps
35 /* This file is completely threadsafe - keep it that way! */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
40 #include <stdio.h>
41 #include <stdlib.h>
42 #include "typedefs.h"
43 #include "network.h"
44 #include "vec.h"
45 #include "gmx_fatal.h"
46 #include "physics.h"
47 #include "main.h"
48 #include "calcmu.h"
49 #include "gmx_omp_nthreads.h"
51 void calc_mu(int start,int homenr,rvec x[],real q[],real qB[],
52 int nChargePerturbed,
53 dvec mu,dvec mu_B)
55 int i,end,m;
56 double mu_x, mu_y, mu_z;
58 end = start + homenr;
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++)
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 mu_x += qB[i]*x[i][XX];
86 mu_y += qB[i]*x[i][YY];
87 mu_z += qB[i]*x[i][ZZ];
89 mu_B[XX] = mu_x * ENM2DEBYE;
90 mu_B[YY] = mu_y * ENM2DEBYE;
91 mu_B[ZZ] = mu_z * ENM2DEBYE;
93 else
95 copy_dvec(mu,mu_B);
99 gmx_bool read_mu(FILE *fp,rvec mu,real *vol)
101 /* For backward compatibility */
102 real mmm[4];
104 if (fread(mmm,(size_t)(4*sizeof(real)),1,fp) != 1)
106 return FALSE;
109 copy_rvec(mmm,mu);
110 *vol = mmm[3];
112 return TRUE;