Improved version for builds from modified trees.
[gromacs/qmmm-gamess-us.git] / src / mdlib / calcmu.c
blob994388bbe4da0faa4ec9c833816eea1b7b4cdcf3
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"
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;
56 end = start + homenr;
58 clear_dvec(mu);
59 for(i=start; (i<end); i++)
60 for(m=0; (m<DIM); m++)
61 mu[m] += q[i]*x[i][m];
63 for(m=0; (m<DIM); m++)
64 mu[m] *= ENM2DEBYE;
66 if (nChargePerturbed) {
67 clear_dvec(mu_B);
68 for(i=start; (i<end); i++)
69 for(m=0; (m<DIM); m++)
70 mu_B[m] += qB[i]*x[i][m];
72 for(m=0; (m<DIM); m++)
73 mu_B[m] *= ENM2DEBYE;
74 } else {
75 copy_dvec(mu,mu_B);
79 bool read_mu(FILE *fp,rvec mu,real *vol)
81 /* For backward compatibility */
82 real mmm[4];
84 if (fread(mmm,(size_t)(4*sizeof(real)),1,fp) != 1)
85 return FALSE;
87 copy_rvec(mmm,mu);
88 *vol = mmm[3];
90 return TRUE;