26fb504dcdf318749f505c37aa4b78f48280198e
[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,2015, 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/cstringutil.h"
55 #include "gromacs/utility/gmxassert.h"
56 #include "gromacs/utility/smalloc.h"
57
58 #include "testutils/testfilemanager.h"
59
60 namespace gmx
61 {
62 namespace test
63 {
64
65 TopologyManager::TopologyManager()
66     : top_(NULL), frame_(NULL)
67 {
68 }
69
70 TopologyManager::~TopologyManager()
71 {
72     if (top_ != NULL)
73     {
74         free_t_atoms(&top_->atoms, TRUE);
75         done_top(top_);
76         sfree(top_);
77     }
78
79     if (frame_ != NULL)
80     {
81         sfree(frame_->x);
82         sfree(frame_->v);
83         sfree(frame_->f);
84         sfree(frame_);
85     }
86 }
87
88 void TopologyManager::requestFrame()
89 {
90     GMX_RELEASE_ASSERT(top_ == NULL,
91                        "Frame must be requested before initializing topology");
92     if (frame_ == NULL)
93     {
94         snew(frame_, 1);
95     }
96 }
97
98 void TopologyManager::requestVelocities()
99 {
100     GMX_RELEASE_ASSERT(frame_ != NULL,
101                        "Velocities requested before requesting a frame");
102     frame_->bV = TRUE;
103     if (frame_->natoms > 0)
104     {
105         snew(frame_->v, frame_->natoms);
106     }
107 }
108
109 void TopologyManager::requestForces()
110 {
111     GMX_RELEASE_ASSERT(frame_ != NULL,
112                        "Forces requested before requesting a frame");
113     frame_->bF = TRUE;
114     if (frame_->natoms > 0)
115     {
116         snew(frame_->f, frame_->natoms);
117     }
118 }
119
120 void TopologyManager::loadTopology(const char *filename)
121 {
122     char    title[STRLEN];
123     int     ePBC;
124     rvec   *xtop = NULL;
125     matrix  box;
126
127     GMX_RELEASE_ASSERT(top_ == NULL, "Topology initialized more than once");
128     snew(top_, 1);
129     read_tps_conf(gmx::test::TestFileManager::getInputFilePath(filename).c_str(),
130                   title, top_, &ePBC, frame_ != NULL ? &xtop : NULL,
131                   NULL, box, FALSE);
132
133     if (frame_ != NULL)
134     {
135         frame_->flags  = TRX_NEED_X;
136         frame_->natoms = top_->atoms.nr;
137         frame_->bX     = TRUE;
138         snew(frame_->x, frame_->natoms);
139         std::memcpy(frame_->x, xtop, sizeof(*frame_->x) * frame_->natoms);
140         frame_->bBox   = TRUE;
141         copy_mat(box, frame_->box);
142     }
143
144     sfree(xtop);
145 }
146
147 void TopologyManager::initAtoms(int count)
148 {
149     GMX_RELEASE_ASSERT(top_ == NULL, "Topology initialized more than once");
150     snew(top_, 1);
151     init_t_atoms(&top_->atoms, count, FALSE);
152     for (int i = 0; i < count; ++i)
153     {
154         top_->atoms.atom[i].m = (i % 3 == 0 ? 2.0 : 1.0);
155     }
156     if (frame_ != NULL)
157     {
158         frame_->flags  = TRX_NEED_X;
159         frame_->natoms = count;
160         frame_->bX     = TRUE;
161         snew(frame_->x, count);
162         if (frame_->bV)
163         {
164             snew(frame_->v, count);
165         }
166         if (frame_->bF)
167         {
168             snew(frame_->f, count);
169         }
170     }
171 }
172
173 void TopologyManager::initAtomTypes(int count, const char *const types[])
174 {
175     GMX_RELEASE_ASSERT(top_ != NULL, "Topology not initialized");
176     atomtypes_.reserve(count);
177     for (int i = 0; i < count; ++i)
178     {
179         atomtypes_.push_back(gmx_strdup(types[i]));
180     }
181     snew(top_->atoms.atomtype, top_->atoms.nr);
182     for (int i = 0, j = 0; i < top_->atoms.nr; ++i, ++j)
183     {
184         if (j == count)
185         {
186             j = 0;
187         }
188         top_->atoms.atomtype[i] = &atomtypes_[j];
189     }
190 }
191
192 void TopologyManager::initUniformResidues(int residueSize)
193 {
194     GMX_RELEASE_ASSERT(top_ != NULL, "Topology not initialized");
195     int residueIndex = -1;
196     for (int i = 0; i < top_->atoms.nr; ++i)
197     {
198         if (i % residueSize == 0)
199         {
200             ++residueIndex;
201         }
202         top_->atoms.atom[i].resind = residueIndex;
203     }
204 }
205
206 void TopologyManager::initUniformMolecules(int moleculeSize)
207 {
208     GMX_RELEASE_ASSERT(top_ != NULL, "Topology not initialized");
209     int index = 0;
210     top_->mols.nalloc_index = (top_->atoms.nr + moleculeSize - 1) / moleculeSize + 1;
211     snew(top_->mols.index, top_->mols.nalloc_index);
212     top_->mols.nr = 0;
213     while (index < top_->atoms.nr)
214     {
215         top_->mols.index[top_->mols.nr] = index;
216         ++top_->mols.nr;
217         index += moleculeSize;
218     }
219     top_->mols.index[top_->mols.nr] = top_->atoms.nr;
220 }
221
222 } // namespace test
223 } // namespace gmx