Sort all includes in src/gromacs
[alexxy/gromacs.git] / src / gromacs / selection / tests / toputils.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2013,2014, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*! \internal \file
36  * \brief
37  * Implements test helper routines from toputils.h.
38  *
39  * \author Teemu Murtola <teemu.murtola@gmail.com>
40  * \ingroup module_selection
41  */
42 #include "gmxpre.h"
43
44 #include "toputils.h"
45
46 #include <cstring>
47
48 #include "gromacs/fileio/tpxio.h"
49 #include "gromacs/fileio/trx.h"
50 #include "gromacs/fileio/trxio.h"
51 #include "gromacs/math/vec.h"
52 #include "gromacs/topology/atoms.h"
53 #include "gromacs/topology/topology.h"
54 #include "gromacs/utility/gmxassert.h"
55 #include "gromacs/utility/smalloc.h"
56
57 #include "testutils/testfilemanager.h"
58
59 namespace gmx
60 {
61 namespace test
62 {
63
64 TopologyManager::TopologyManager()
65     : top_(NULL), frame_(NULL)
66 {
67 }
68
69 TopologyManager::~TopologyManager()
70 {
71     if (top_ != NULL)
72     {
73         free_t_atoms(&top_->atoms, TRUE);
74         done_top(top_);
75         sfree(top_);
76     }
77
78     if (frame_ != NULL)
79     {
80         sfree(frame_->x);
81         sfree(frame_->v);
82         sfree(frame_->f);
83         sfree(frame_);
84     }
85 }
86
87 void TopologyManager::requestFrame()
88 {
89     GMX_RELEASE_ASSERT(top_ == NULL,
90                        "Frame must be requested before initializing topology");
91     if (frame_ == NULL)
92     {
93         snew(frame_, 1);
94     }
95 }
96
97 void TopologyManager::requestVelocities()
98 {
99     GMX_RELEASE_ASSERT(frame_ != NULL,
100                        "Velocities requested before requesting a frame");
101     frame_->bV = TRUE;
102     if (frame_->natoms > 0)
103     {
104         snew(frame_->v, frame_->natoms);
105     }
106 }
107
108 void TopologyManager::requestForces()
109 {
110     GMX_RELEASE_ASSERT(frame_ != NULL,
111                        "Forces requested before requesting a frame");
112     frame_->bF = TRUE;
113     if (frame_->natoms > 0)
114     {
115         snew(frame_->f, frame_->natoms);
116     }
117 }
118
119 void TopologyManager::loadTopology(const char *filename)
120 {
121     char    title[STRLEN];
122     int     ePBC;
123     rvec   *xtop = NULL;
124     matrix  box;
125
126     GMX_RELEASE_ASSERT(top_ == NULL, "Topology initialized more than once");
127     snew(top_, 1);
128     read_tps_conf(gmx::test::TestFileManager::getInputFilePath(filename).c_str(),
129                   title, top_, &ePBC, frame_ != NULL ? &xtop : NULL,
130                   NULL, box, FALSE);
131
132     if (frame_ != NULL)
133     {
134         frame_->flags  = TRX_NEED_X;
135         frame_->natoms = top_->atoms.nr;
136         frame_->bX     = TRUE;
137         snew(frame_->x, frame_->natoms);
138         std::memcpy(frame_->x, xtop, sizeof(*frame_->x) * frame_->natoms);
139         frame_->bBox   = TRUE;
140         copy_mat(box, frame_->box);
141     }
142
143     sfree(xtop);
144 }
145
146 void TopologyManager::initAtoms(int count)
147 {
148     GMX_RELEASE_ASSERT(top_ == NULL, "Topology initialized more than once");
149     snew(top_, 1);
150     init_t_atoms(&top_->atoms, count, FALSE);
151     for (int i = 0; i < count; ++i)
152     {
153         top_->atoms.atom[i].m = (i % 3 == 0 ? 2.0 : 1.0);
154     }
155     if (frame_ != NULL)
156     {
157         frame_->flags  = TRX_NEED_X;
158         frame_->natoms = count;
159         frame_->bX     = TRUE;
160         snew(frame_->x, count);
161         if (frame_->bV)
162         {
163             snew(frame_->v, count);
164         }
165         if (frame_->bF)
166         {
167             snew(frame_->f, count);
168         }
169     }
170 }
171
172 void TopologyManager::initAtomTypes(int count, const char *const types[])
173 {
174     GMX_RELEASE_ASSERT(top_ != NULL, "Topology not initialized");
175     atomtypes_.reserve(count);
176     for (int i = 0; i < count; ++i)
177     {
178         atomtypes_.push_back(gmx_strdup(types[i]));
179     }
180     snew(top_->atoms.atomtype, top_->atoms.nr);
181     for (int i = 0, j = 0; i < top_->atoms.nr; ++i, ++j)
182     {
183         if (j == count)
184         {
185             j = 0;
186         }
187         top_->atoms.atomtype[i] = &atomtypes_[j];
188     }
189 }
190
191 void TopologyManager::initUniformResidues(int residueSize)
192 {
193     GMX_RELEASE_ASSERT(top_ != NULL, "Topology not initialized");
194     int residueIndex = -1;
195     for (int i = 0; i < top_->atoms.nr; ++i)
196     {
197         if (i % residueSize == 0)
198         {
199             ++residueIndex;
200         }
201         top_->atoms.atom[i].resind = residueIndex;
202     }
203 }
204
205 void TopologyManager::initUniformMolecules(int moleculeSize)
206 {
207     GMX_RELEASE_ASSERT(top_ != NULL, "Topology not initialized");
208     int index = 0;
209     top_->mols.nalloc_index = (top_->atoms.nr + moleculeSize - 1) / moleculeSize + 1;
210     snew(top_->mols.index, top_->mols.nalloc_index);
211     top_->mols.nr = 0;
212     while (index < top_->atoms.nr)
213     {
214         top_->mols.index[top_->mols.nr] = index;
215         ++top_->mols.nr;
216         index += moleculeSize;
217     }
218     top_->mols.index[top_->mols.nr] = top_->atoms.nr;
219 }
220
221 } // namespace test
222 } // namespace gmx