Add AWH biasing module + tests
[gromacs.git] / src / gromacs / gmxpreprocess / read-conformation.cpp
blob950f751acbb28347abd8c9f0afd5bdca466c2812
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2014,2015,2016,2017, 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.
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.
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.
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.
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.
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.
35 #include "gmxpre.h"
37 #include "read-conformation.h"
39 #include <vector>
41 #include "gromacs/fileio/confio.h"
42 #include "gromacs/topology/atomprop.h"
43 #include "gromacs/topology/atoms.h"
44 #include "gromacs/topology/mtop_util.h"
45 #include "gromacs/topology/topology.h"
46 #include "gromacs/utility/cstringutil.h"
47 #include "gromacs/utility/smalloc.h"
48 #include "gromacs/utility/unique_cptr.h"
50 using gmx::RVec;
52 std::vector<real>
53 makeExclusionDistances(const t_atoms *a, gmx_atomprop_t aps,
54 real defaultDistance, real scaleFactor)
56 std::vector<real> exclusionDistances;
58 exclusionDistances.reserve(a->nr);
59 for (int i = 0; i < a->nr; ++i)
61 real value;
62 if (!gmx_atomprop_query(aps, epropVDW,
63 *(a->resinfo[a->atom[i].resind].name),
64 *(a->atomname[i]), &value))
66 value = defaultDistance;
68 else
70 value *= scaleFactor;
72 exclusionDistances.push_back(value);
74 return exclusionDistances;
77 void readConformation(const char *confin, gmx_mtop_t *top,
78 std::vector<RVec> *x, std::vector<RVec> *v,
79 int *ePBC, matrix box, const char *statusTitle)
81 fprintf(stderr, "Reading %s configuration%s\n", statusTitle,
82 v ? " and velocities" : "");
83 rvec *x_tmp = nullptr, *v_tmp = nullptr;
84 bool dummy;
85 readConfAndTopology(confin, &dummy, top, ePBC, x ? &x_tmp : nullptr, v ? &v_tmp : nullptr, box);
86 const gmx::sfree_guard xguard(x_tmp);
87 const gmx::sfree_guard vguard(v_tmp);
88 if (x && x_tmp)
90 *x = std::vector<RVec>(x_tmp, x_tmp + top->natoms);
92 if (v && v_tmp)
94 *v = std::vector<RVec>(v_tmp, v_tmp + top->natoms);
96 fprintf(stderr, "%s\nContaining %d atoms in %d residues\n",
97 *top->name, top->natoms, gmx_mtop_nres(top));
100 void readConformation(const char *confin, t_topology *top,
101 std::vector<RVec> *x, std::vector<RVec> *v,
102 int *ePBC, matrix box, const char *statusTitle)
104 fprintf(stderr, "Reading %s configuration%s\n", statusTitle,
105 v ? " and velocities" : "");
106 rvec *x_tmp = nullptr, *v_tmp = nullptr;
107 read_tps_conf(confin, top, ePBC, x ? &x_tmp : nullptr, v ? &v_tmp : nullptr, box, FALSE);
108 const gmx::sfree_guard xguard(x_tmp);
109 const gmx::sfree_guard vguard(v_tmp);
110 if (x && x_tmp)
112 *x = std::vector<RVec>(x_tmp, x_tmp + top->atoms.nr);
114 if (v && v_tmp)
116 *v = std::vector<RVec>(v_tmp, v_tmp + top->atoms.nr);
118 fprintf(stderr, "%s\nContaining %d atoms in %d residues\n",
119 *top->name, top->atoms.nr, top->atoms.nres);