Merge branch 'master' of git://git.gromacs.org/gromacs
[gromacs/adressmacs.git] / src / gmxlib / rmpbc.c
blob0f7d005fb354bb90e23156e6a174cf438bd84f8e
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
3 *
4 * This source code is part of
5 *
6 * G R O M A C S
7 *
8 * GROningen MAchine for Chemical Simulations
9 *
10 * VERSION 3.2.0
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
33 * And Hey:
34 * GROningen Mixture of Alchemy and Childrens' Stories
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
40 #include "sysstuff.h"
41 #include "typedefs.h"
42 #include "smalloc.h"
43 #include "mshift.h"
44 #include "pbc.h"
45 #include "gstat.h"
46 #include "futil.h"
47 #include "vec.h"
49 typedef struct {
50 int natoms;
51 t_graph *gr;
52 } rmpbc_graph_t;
54 typedef struct gmx_rmpbc {
55 t_idef *idef;
56 int natoms_init;
57 int ePBC;
58 int ngraph;
59 rmpbc_graph_t *graph;
60 } koeiepoep;
62 static t_graph *gmx_rmpbc_get_graph(gmx_rmpbc_t gpbc,int ePBC,int natoms)
64 int i;
65 rmpbc_graph_t *gr;
67 if (ePBC == epbcNONE || gpbc->idef->ntypes <= 0)
69 return NULL;
72 gr = NULL;
73 for(i=0; i<gpbc->ngraph; i++)
75 if (natoms == gpbc->graph[i].natoms)
77 gr = &gpbc->graph[i];
80 if (gr == NULL)
82 /* We'd like to check with the number of atoms in the topology,
83 * but we don't have that available.
84 * So we check against the number of atoms that gmx_rmpbc_init
85 * was called with.
87 if (natoms > gpbc->natoms_init)
89 gmx_fatal(FARGS,"Structure or trajectory file has more atoms (%d) than the topology (%d)",natoms,gpbc->natoms_init);
91 gpbc->ngraph++;
92 srenew(gpbc->graph,gpbc->ngraph);
93 gr = &gpbc->graph[gpbc->ngraph-1];
94 gr->natoms = natoms;
95 gr->gr = mk_graph(NULL,gpbc->idef,0,natoms,FALSE,FALSE);
98 return gr->gr;
101 gmx_rmpbc_t gmx_rmpbc_init(t_idef *idef,int ePBC,int natoms,
102 matrix box)
104 gmx_rmpbc_t gpbc;
106 snew(gpbc,1);
108 gpbc->natoms_init = natoms;
110 /* This sets pbc when we now it,
111 * otherwise we guess it from the instantaneous box in the trajectory.
113 gpbc->ePBC = ePBC;
115 gpbc->idef = idef;
116 if (gpbc->idef->ntypes <= 0)
118 fprintf(stderr,
119 "\n"
120 "WARNING: if there are broken molecules in the trajectory file,\n"
121 " they can not be made whole without a run input file\n\n");
124 return gpbc;
127 void gmx_rmpbc_done(gmx_rmpbc_t gpbc)
129 int i;
131 for(i=0; i<gpbc->ngraph; i++)
133 done_graph(gpbc->graph[i].gr);
135 if (gpbc->graph != NULL)
137 sfree(gpbc->graph);
141 static int gmx_rmpbc_ePBC(gmx_rmpbc_t gpbc,matrix box)
143 if (gpbc->ePBC >= 0)
145 return gpbc->ePBC;
147 else
149 return guess_ePBC(box);
153 void gmx_rmpbc(gmx_rmpbc_t gpbc,int natoms,matrix box,rvec x[])
155 int ePBC;
156 t_graph *gr;
158 ePBC = gmx_rmpbc_ePBC(gpbc,box);
159 gr = gmx_rmpbc_get_graph(gpbc,ePBC,natoms);
160 if (gr != NULL)
162 mk_mshift(stdout,gr,ePBC,box,x);
163 shift_x(gr,box,x,x);
167 void gmx_rmpbc_copy(gmx_rmpbc_t gpbc,int natoms,matrix box,rvec x[],rvec x_s[])
169 int ePBC;
170 t_graph *gr;
171 int i;
173 ePBC = gmx_rmpbc_ePBC(gpbc,box);
174 gr = gmx_rmpbc_get_graph(gpbc,ePBC,natoms);
175 if (gr != NULL)
177 mk_mshift(stdout,gr,ePBC,box,x);
178 shift_x(gr,box,x,x_s);
180 else
182 for(i=0; i<natoms; i++)
184 copy_rvec(x[i],x_s[i]);
189 void gmx_rmpbc_trxfr(gmx_rmpbc_t gpbc,t_trxframe *fr)
191 int ePBC;
192 t_graph *gr;
194 if (fr->bX && fr->bBox)
196 ePBC = gmx_rmpbc_ePBC(gpbc,fr->box);
197 gr = gmx_rmpbc_get_graph(gpbc,ePBC,fr->natoms);
198 if (gr != NULL)
200 mk_mshift(stdout,gr,ePBC,fr->box,fr->x);
201 shift_x(gr,fr->box,fr->x,fr->x);
206 void rm_gropbc(t_atoms *atoms,rvec x[],matrix box)
208 real dist;
209 int n,m,d;
211 /* check periodic boundary */
212 for(n=1;(n<atoms->nr);n++)
214 for(m=DIM-1; m>=0; m--)
216 dist = x[n][m]-x[n-1][m];
217 if (fabs(dist) > 0.9*box[m][m])
219 if ( dist > 0 )
221 for(d=0; d<=m; d++)
223 x[n][d] -= box[m][d];
226 else
228 for(d=0; d<=m; d++)
230 x[n][d] += box[m][d];