Revised wording in pdb2gmx.c, hopefully clearer now.
[gromacs/rigid-bodies.git] / src / gmxlib / rmpbc.c
blob86b713968cc2cf39373cb1daec6e46ea3cafb260
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
68 || NULL == gpbc
69 || NULL == gpbc->idef
70 || gpbc->idef->ntypes <= 0)
72 return NULL;
75 gr = NULL;
76 for(i=0; i<gpbc->ngraph; i++)
78 if (natoms == gpbc->graph[i].natoms)
80 gr = &gpbc->graph[i];
83 if (gr == NULL)
85 /* We'd like to check with the number of atoms in the topology,
86 * but we don't have that available.
87 * So we check against the number of atoms that gmx_rmpbc_init
88 * was called with.
90 if (natoms > gpbc->natoms_init)
92 gmx_fatal(FARGS,"Structure or trajectory file has more atoms (%d) than the topology (%d)",natoms,gpbc->natoms_init);
94 gpbc->ngraph++;
95 srenew(gpbc->graph,gpbc->ngraph);
96 gr = &gpbc->graph[gpbc->ngraph-1];
97 gr->natoms = natoms;
98 gr->gr = mk_graph(NULL,gpbc->idef,0,natoms,FALSE,FALSE);
101 return gr->gr;
104 gmx_rmpbc_t gmx_rmpbc_init(t_idef *idef,int ePBC,int natoms,
105 matrix box)
107 gmx_rmpbc_t gpbc;
109 snew(gpbc,1);
111 gpbc->natoms_init = natoms;
113 /* This sets pbc when we now it,
114 * otherwise we guess it from the instantaneous box in the trajectory.
116 gpbc->ePBC = ePBC;
118 gpbc->idef = idef;
119 if (gpbc->idef->ntypes <= 0)
121 fprintf(stderr,
122 "\n"
123 "WARNING: if there are broken molecules in the trajectory file,\n"
124 " they can not be made whole without a run input file\n\n");
127 return gpbc;
130 void gmx_rmpbc_done(gmx_rmpbc_t gpbc)
132 int i;
134 if (NULL != gpbc)
136 for(i=0; i<gpbc->ngraph; i++)
138 done_graph(gpbc->graph[i].gr);
140 if (gpbc->graph != NULL)
142 sfree(gpbc->graph);
147 static int gmx_rmpbc_ePBC(gmx_rmpbc_t gpbc,matrix box)
149 if (NULL != gpbc && gpbc->ePBC >= 0)
151 return gpbc->ePBC;
153 else
155 return guess_ePBC(box);
159 void gmx_rmpbc(gmx_rmpbc_t gpbc,int natoms,matrix box,rvec x[])
161 int ePBC;
162 t_graph *gr;
164 ePBC = gmx_rmpbc_ePBC(gpbc,box);
165 gr = gmx_rmpbc_get_graph(gpbc,ePBC,natoms);
166 if (gr != NULL)
168 mk_mshift(stdout,gr,ePBC,box,x);
169 shift_x(gr,box,x,x);
173 void gmx_rmpbc_copy(gmx_rmpbc_t gpbc,int natoms,matrix box,rvec x[],rvec x_s[])
175 int ePBC;
176 t_graph *gr;
177 int i;
179 ePBC = gmx_rmpbc_ePBC(gpbc,box);
180 gr = gmx_rmpbc_get_graph(gpbc,ePBC,natoms);
181 if (gr != NULL)
183 mk_mshift(stdout,gr,ePBC,box,x);
184 shift_x(gr,box,x,x_s);
186 else
188 for(i=0; i<natoms; i++)
190 copy_rvec(x[i],x_s[i]);
195 void gmx_rmpbc_trxfr(gmx_rmpbc_t gpbc,t_trxframe *fr)
197 int ePBC;
198 t_graph *gr;
200 if (fr->bX && fr->bBox)
202 ePBC = gmx_rmpbc_ePBC(gpbc,fr->box);
203 gr = gmx_rmpbc_get_graph(gpbc,ePBC,fr->natoms);
204 if (gr != NULL)
206 mk_mshift(stdout,gr,ePBC,fr->box,fr->x);
207 shift_x(gr,fr->box,fr->x,fr->x);
212 void rm_gropbc(t_atoms *atoms,rvec x[],matrix box)
214 real dist;
215 int n,m,d;
217 /* check periodic boundary */
218 for(n=1;(n<atoms->nr);n++)
220 for(m=DIM-1; m>=0; m--)
222 dist = x[n][m]-x[n-1][m];
223 if (fabs(dist) > 0.9*box[m][m])
225 if ( dist > 0 )
227 for(d=0; d<=m; d++)
229 x[n][d] -= box[m][d];
232 else
234 for(d=0; d<=m; d++)
236 x[n][d] += box[m][d];