Update copyright statements and change license to LGPL
[alexxy/gromacs.git] / src / gmxlib / rmpbc.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  * check out http://www.gromacs.org for more information.
7  * Copyright (c) 2012, by the GROMACS development team, led by
8  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9  * others, as listed in the AUTHORS file in the top-level source
10  * directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
41
42 #include "sysstuff.h"
43 #include "typedefs.h"
44 #include "smalloc.h"
45 #include "mshift.h"
46 #include "pbc.h"
47 #include "gstat.h"
48 #include "futil.h"
49 #include "vec.h"        
50
51 typedef struct {
52     int     natoms;
53     t_graph *gr;
54 } rmpbc_graph_t;
55
56 typedef struct gmx_rmpbc {
57     t_idef        *idef;
58     int           natoms_init;
59     int           ePBC;
60     int           ngraph;
61     rmpbc_graph_t *graph;
62 } koeiepoep;
63
64 static t_graph *gmx_rmpbc_get_graph(gmx_rmpbc_t gpbc,int ePBC,int natoms)
65 {
66     int           i;
67     rmpbc_graph_t *gr;
68
69     if (ePBC == epbcNONE
70         || NULL == gpbc
71         || NULL == gpbc->idef
72         || gpbc->idef->ntypes <= 0)
73     {
74         return NULL;
75     }
76
77     gr = NULL;
78     for(i=0; i<gpbc->ngraph; i++)
79     {
80         if (natoms == gpbc->graph[i].natoms)
81         {
82             gr = &gpbc->graph[i];
83         }
84     }
85     if (gr == NULL)
86     {
87         /* We'd like to check with the number of atoms in the topology,
88          * but we don't have that available.
89          * So we check against the number of atoms that gmx_rmpbc_init
90          * was called with.
91          */
92         if (natoms > gpbc->natoms_init)
93         {
94             gmx_fatal(FARGS,"Structure or trajectory file has more atoms (%d) than the topology (%d)",natoms,gpbc->natoms_init);
95         }
96         gpbc->ngraph++;
97         srenew(gpbc->graph,gpbc->ngraph);
98         gr = &gpbc->graph[gpbc->ngraph-1];
99         gr->natoms = natoms;
100         gr->gr     = mk_graph(NULL,gpbc->idef,0,natoms,FALSE,FALSE);
101     }
102
103     return gr->gr;
104 }
105
106 gmx_rmpbc_t gmx_rmpbc_init(t_idef *idef,int ePBC,int natoms,
107                            matrix box)
108 {
109     gmx_rmpbc_t gpbc;
110   
111     snew(gpbc,1);
112
113     gpbc->natoms_init = natoms;
114   
115     /* This sets pbc when we now it,
116      * otherwise we guess it from the instantaneous box in the trajectory.
117      */
118     gpbc->ePBC = ePBC;
119
120     gpbc->idef = idef;
121     if (gpbc->idef->ntypes <= 0)
122     {
123         fprintf(stderr,
124                 "\n"
125                 "WARNING: If there are molecules in the input trajectory file\n"
126                 "         that are broken across periodic boundaries, they\n"
127                 "         cannot be made whole (or treated as whole) without\n"
128                 "         you providing a run input file.\n\n");
129     }
130
131     return gpbc;
132 }
133
134 void gmx_rmpbc_done(gmx_rmpbc_t gpbc)
135 {
136     int i;
137
138     if (NULL != gpbc)
139     {
140         for(i=0; i<gpbc->ngraph; i++)
141         {
142             done_graph(gpbc->graph[i].gr);
143         }
144         if (gpbc->graph != NULL)
145         {
146             sfree(gpbc->graph);
147         }
148     }
149 }
150
151 static int gmx_rmpbc_ePBC(gmx_rmpbc_t gpbc,matrix box)
152 {
153     if (NULL != gpbc && gpbc->ePBC >= 0)
154     {
155         return gpbc->ePBC;
156     }
157     else
158     {
159         return guess_ePBC(box);
160     }
161 }
162
163 void gmx_rmpbc(gmx_rmpbc_t gpbc,int natoms,matrix box,rvec x[])
164 {
165     int     ePBC;
166     t_graph *gr;
167     
168     ePBC = gmx_rmpbc_ePBC(gpbc,box);
169     gr = gmx_rmpbc_get_graph(gpbc,ePBC,natoms);
170     if (gr != NULL)
171     {
172         mk_mshift(stdout,gr,ePBC,box,x);
173         shift_self(gr,box,x);
174     }
175 }
176
177 void gmx_rmpbc_copy(gmx_rmpbc_t gpbc,int natoms,matrix box,rvec x[],rvec x_s[])
178 {
179     int     ePBC;
180     t_graph *gr;
181     int     i;
182
183     ePBC = gmx_rmpbc_ePBC(gpbc,box);
184     gr = gmx_rmpbc_get_graph(gpbc,ePBC,natoms);
185     if (gr != NULL)
186     {
187         mk_mshift(stdout,gr,ePBC,box,x);
188         shift_x(gr,box,x,x_s);
189     }
190     else
191     {
192         for(i=0; i<natoms; i++)
193         {
194             copy_rvec(x[i],x_s[i]);
195         }
196     }
197 }
198
199 void gmx_rmpbc_trxfr(gmx_rmpbc_t gpbc,t_trxframe *fr)
200 {
201     int     ePBC;
202     t_graph *gr;
203
204     if (fr->bX && fr->bBox)
205     {
206         ePBC = gmx_rmpbc_ePBC(gpbc,fr->box);
207         gr = gmx_rmpbc_get_graph(gpbc,ePBC,fr->natoms);
208         if (gr != NULL)
209         {
210             mk_mshift(stdout,gr,ePBC,fr->box,fr->x);
211             shift_self(gr,fr->box,fr->x);
212         }
213     }
214 }
215
216 void rm_gropbc(t_atoms *atoms,rvec x[],matrix box)
217 {
218     real dist;
219     int  n,m,d;
220   
221     /* check periodic boundary */
222     for(n=1;(n<atoms->nr);n++)
223     {
224         for(m=DIM-1; m>=0; m--)
225         {
226             dist = x[n][m]-x[n-1][m];
227             if (fabs(dist) > 0.9*box[m][m])
228             {
229                 if ( dist >  0 )
230                 {
231                     for(d=0; d<=m; d++)
232                     {
233                         x[n][d] -= box[m][d];
234                     }
235                 }
236                 else
237                 {
238                     for(d=0; d<=m; d++)
239                     {
240                         x[n][d] += box[m][d];
241                     }
242                 }
243             }   
244         }
245     }
246 }
247