2 * This file is part of the GROMACS molecular simulation package.
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.
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.
37 * Implements gmx::analysismodules::Freevolume.
39 * \author Titov Anatoly <Wapuk-cobaka@yandex.ru>
40 * \ingroup module_trajectoryanalysis
49 #include <gromacs/trajectoryanalysis.h>
50 #include <gromacs/pbcutil/pbc.h>
51 #include <gromacs/pbcutil/rmpbc.h>
52 #include <gromacs/utility/smalloc.h>
53 #include <gromacs/math/vectypes.h>
54 #include <gromacs/math/vec.h>
55 #include <gromacs/math/do_fit.h>
57 #define MAX_NTRICVEC 12
69 void make_graph(int mgwi_natoms, rvec *mgwi_x, std::vector< std::vector< node > > &mgwi_graph)
71 mgwi_graph.resize(mgwi_natoms);
72 for (int i = 0; i < mgwi_natoms; i++) {
73 mgwi_graph[i].resize(mgwi_natoms);
75 for (int i = 0; i < mgwi_natoms; i++) {
76 for (int j = 0; j < mgwi_natoms; j++) {
77 rvec_sub(mgwi_x[i], mgwi_x[j], mgwi_graph[i][j].r);
78 mgwi_graph[i][j].n = 0;
83 void update_graph(std::vector< std::vector< node > > &ugwi_graph, rvec *ugwi_x, long double ugwi_epsi) {
85 int ugwi_for = ugwi_graph.size();
86 for (int i = 0; i < ugwi_for; i++) {
87 for (int j = i; j < ugwi_for; j++) {
88 rvec_sub(ugwi_x[i], ugwi_x[j], ugwi_temp);
89 rvec_dec(ugwi_temp, ugwi_graph[i][j].r.as_vec());
90 if (norm(ugwi_temp) <= ugwi_epsi) {
103 void check_domains(long double cd_delta, int cd_frames, std::vector< std::vector< std::vector< node > > > &cd_graph) {
104 int cd_for1 = cd_graph.size(), cd_for2 = cd_graph[1].size();
105 for (int k = 0; k < cd_for1; k++) {
106 for (int i = 0; i < cd_for2; i++) {
107 for (int j = 0; j < cd_for2; j++) {
108 if (cd_graph[k][i][j].n >= cd_frames * cd_delta) {
109 cd_graph[k][i][j].yep = true;
112 cd_graph[k][i][j].yep = false;
119 void find_domain_sizes(std::vector< std::vector< std::vector< node > > > fds_graph, std::vector< std::vector< int > > &fds_domsizes) {
120 fds_domsizes.resize(fds_graph.size());
121 int fds_for1 = fds_graph.size(), fds_for2 = fds_graph[1].size();
122 for (int i = 0; i < fds_for1; i++) {
123 fds_domsizes[i].resize(fds_for2, 0);
124 for (int j = 0; j < fds_for2; j++) {
125 for (int k = 0; k < fds_for2; k++) {
126 if (fds_graph[i][j][k].yep) {
127 fds_domsizes[i][j]++;
134 void get_maxsized_domain(std::vector< int > &gmd_max_d, std::vector< std::vector< std::vector< node > > > gmd_graph, std::vector< std::vector< int > > gmd_domsizes) {
135 int gmd_number1 = 0, gmd_number2 = 0;
136 int gmd_for1 = gmd_domsizes.size(), gmd_for2 = gmd_domsizes[0].size();
137 for (int i = 0; i < gmd_for1; i++) {
138 for (int j = 0; j < gmd_for2; j++) {
139 if (gmd_domsizes[i][j] >= gmd_domsizes[gmd_number1][gmd_number2]) {
146 int gmd_for3 = gmd_graph[gmd_number1][gmd_number2].size();
147 for (int i = 0; i < gmd_for3; i++) {
148 if (gmd_graph[gmd_number1][gmd_number2][i].yep) {
149 gmd_max_d.push_back(i);
154 void get_min_domain(std::vector< int > &gmd_min_d, std::vector< std::vector< std::vector< node > > > gmd_graph, std::vector< std::vector< int > > gmd_domsizes, int gmd_min_dom_size) {
155 int gmd_number1 = 0, gmd_number2 = 0;
156 int gmd_for1 = gmd_domsizes.size(), gmd_for2 = gmd_domsizes[0].size();
157 for (int i = 0; i < gmd_for1; i++) {
158 for (int j = 0; j < gmd_for2; j++) {
159 if (gmd_domsizes[gmd_number1][gmd_number2] < gmd_min_dom_size && gmd_domsizes[i][j] >= gmd_min_dom_size) {
163 if (gmd_domsizes[i][j] <= gmd_domsizes[gmd_number1][gmd_number2] && gmd_domsizes[i][j] >= gmd_min_dom_size) {
170 int gmd_for3 = gmd_graph[gmd_number1][gmd_number2].size();
171 for (int i = 0; i < gmd_for3; i++) {
172 if (gmd_graph[gmd_number1][gmd_number2][i].yep) {
173 gmd_min_d.push_back(i);
178 void delete_domain_from_graph(std::vector< std::vector< std::vector< node > > > &ddf_graph, std::vector< int > ddf_domain) {
179 int ddfg_for1 = ddf_domain.size(), ddfg_for2 = ddf_graph.size(), ddfg_for3 = ddf_graph[1].size();
180 for (int i = 0; i < ddfg_for1; i++) {
181 for (int k = 0; k < ddfg_for2; k++) {
182 for (int j = 0; j < ddfg_for3; j++) {
183 if (ddf_graph[k][ddf_domain[i]][j].yep) {
184 ddf_graph[k][ddf_domain[i]][j].yep = false;
186 if (ddf_graph[k][j][ddf_domain[i]].yep) {
187 ddf_graph[k][j][ddf_domain[i]].yep = false;
194 bool check_domsizes(std::vector< std::vector< int > > cd_domsizes, int cd_domain_min_size) {
195 int cd_for1 = cd_domsizes.size(), cd_for2 = cd_domsizes[0].size();
196 for (int i = 0; i < cd_for1; i++) {
197 for (int j = 0; j < cd_for2; j++) {
198 if (cd_domsizes[i][j] >= cd_domain_min_size) {
206 void print_domains(std::vector< std::vector< int > > pd_domains, std::vector< int > index, std::string fnNdx_) {
208 fpNdx_ = std::fopen(fnNdx_.c_str(), "w+");
210 for (int i = 0; i < pd_domains.size(); i++) {
211 std::fprintf(fpNdx_, "[domain_%d]\n", i + 1);
213 for (int j = 0; j < pd_domains[i].size(); j++) {
215 if (write_count > 20) {
217 std::fprintf(fpNdx_, "\n");
219 std::fprintf(fpNdx_, "%5d ", index[pd_domains[i][j]] + 1);
221 std::fprintf(fpNdx_,"\n\n");
223 std::fprintf(fpNdx_,"\n");
229 * Class used to compute free volume in a simulations box.
231 * Inherits TrajectoryAnalysisModule and all functions from there.
232 * Does not implement any new functionality.
234 * \ingroup module_trajectoryanalysis
236 class Domains : public TrajectoryAnalysisModule
243 //! Set the options and setting
244 virtual void initOptions(IOptionsContainer *options,
245 TrajectoryAnalysisSettings *settings);
247 //! First routine called by the analysis framework
248 // virtual void initAnalysis(const t_trxframe &fr, t_pbc *pbc);
249 virtual void initAnalysis(const TrajectoryAnalysisSettings &settings,
250 const TopologyInformation &top);
252 virtual void initAfterFirstFrame(const TrajectoryAnalysisSettings &settings,
253 const t_trxframe &fr);
255 //! Call for each frame of the trajectory
256 // virtual void analyzeFrame(const t_trxframe &fr, t_pbc *pbc);
257 virtual void analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
258 TrajectoryAnalysisModuleData *pdata);
260 //! Last routine called by the analysis framework
261 // virtual void finishAnalysis(t_pbc *pbc);
262 virtual void finishAnalysis(int nframes);
264 //! Routine to write output, that is additional over the built-in
265 virtual void writeOutput();
271 std::vector< std::vector< std::vector< node > > > graph;
273 std::vector< std::vector< int > > domains;
274 std::vector< std::vector< int > > domsizes;
276 std::vector< int > index;
277 std::vector< int > numbers;
278 std::vector< std::vector < RVec > > trajectory;
281 int domain_min_size = 5; // selectable
285 double delta = 0.90; //0.95 // selectable
286 double epsi = 0.15; //0.30 // selectable
290 int DomainSearchingAlgorythm = 0; // selectable
291 // Copy and assign disallowed by base.
294 Domains::Domains(): TrajectoryAnalysisModule()
303 Domains::initOptions(IOptionsContainer *options,
304 TrajectoryAnalysisSettings *settings)
306 static const char *const desc[] = {
307 "[THISMODULE] to be done"
309 // Add the descriptive text (program help text) to the options
310 settings->setHelpText(desc);
311 // Add option for selecting a subset of atoms
312 options->addOption(SelectionOption("select")
313 .store(&selec).required()
314 .description("Atoms that are considered as part of the excluded volume"));
315 // Add option for output file name
316 options->addOption(FileNameOption("on").filetype(eftIndex).outputFile()
317 .store(&fnNdx_).defaultBasename("domains")
318 .description("Index file from the domains"));
319 // Add option for domain min size constant
320 options->addOption(gmx::IntegerOption("dms")
321 .store(&domain_min_size)
322 .description("minimum domain size"));
323 // Add option for Domain's Searching Algorythm
324 options->addOption(gmx::IntegerOption("DSA")
325 .store(&DomainSearchingAlgorythm)
326 .description("Domain's Searching Algorythm: 0 == default (from bigger to smaller) | 1 == (from smaller to bigger)"));
327 // Add option for epsi constant
328 options->addOption(DoubleOption("epsilon")
330 .description("thermal vibrations' constant"));
331 // Add option for delta constant
332 options->addOption(DoubleOption("delta")
334 .description("domain membership probability"));
335 // Control input settings
336 settings->setFlags(TrajectoryAnalysisSettings::efNoUserPBC);
337 settings->setPBC(true);
341 Domains::initAnalysis(const TrajectoryAnalysisSettings &settings,
342 const TopologyInformation &top)
344 domains_ePBC = top.ePBC();
348 Domains::initAfterFirstFrame(const TrajectoryAnalysisSettings &settings,
349 const t_trxframe &fr)
352 t_pbc *ppbc = settings.hasPBC() ? &pbc : NULL;
354 copy_mat(fr.box, boxx);
356 set_pbc(ppbc, domains_ePBC, boxx);
358 ConstArrayRef< int > atomind = selec.atomIndices();
360 for (ConstArrayRef<int>::iterator ai = atomind.begin(); (ai < atomind.end()); ai++) {
361 index.push_back(*ai);
363 trajectory.resize(2);
364 trajectory[0].resize(selec.atomCount());
366 for (int i = 0; i < selec.atomCount(); i++) {
367 trajectory[0][i] = fr.x[index[i]];
370 bone = index.size() - domain_min_size + 1;
373 for (int i = 0; i < bone; i++) {
374 snew(w_rls[i], index.size());
375 for (int j = 0; j < index.size(); j++) {
376 if (j >= i && j <= i + domain_min_size - 1) {
383 snew(etalon, index.size());
384 for (int j = 0; j < index.size(); j++) {
385 copy_rvec(trajectory[0][j].as_vec(), etalon[j]);
387 reset_x(index.size(), NULL, index.size(), NULL, etalon, w_rls[i]);
388 make_graph(index.size(), etalon, graph[i]);
391 trajectory[1].resize(index.size());
395 Domains::analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
396 TrajectoryAnalysisModuleData *pdata)
398 for (int i = 0; i < index.size(); i++) {
399 trajectory[1][i] = fr.x[index[i]];
405 #pragma omp for schedule(dynamic)
406 for (int j = 0; j < bone; j++) {
408 snew(etalon, index.size());
409 for (int k = 0; k < index.size(); k++) {
410 copy_rvec(trajectory[0][k].as_vec(), etalon[k]);
412 snew(traj, index.size());
413 for (int k = 0; k < index.size(); k++) {
414 copy_rvec(trajectory[1][k].as_vec(), traj[k]);
416 reset_x(index.size(), NULL, index.size(), NULL, etalon, w_rls[j]);
417 reset_x(index.size(), NULL, index.size(), NULL, traj, w_rls[j]);
418 do_fit(index.size(), w_rls[j], etalon, traj);
419 update_graph(graph[j], traj, epsi);
424 std::cout << "frame: " << frames << " analyzed\n";
427 //domains -s '/home/toluk/Develop/samples/reca_rd/reca_rd.mono.tpr' -f '/home/toluk/Develop/samples/reca_rd/reca_rd.mono.xtc' -select 'name CA'
428 //domains -s '/home/toluk/Develop/samples/banana_phone/pgk.md.non-sol.tpr' -f '/home/toluk/Develop/samples/banana_phone/pgk.md.non-sol.10th.xtc' -select 'name CA'
429 //domains -s '/home/toluk/Data/reca_rd/reca_rd.mono.tpr' -f '/home/toluk/Data/reca_rd/reca_rd.mono.xtc' - select 'name CA'
432 Domains::finishAnalysis(int nframes)
436 std::cout << "final cheking\n";
437 check_domains(delta, frames, graph);
439 std::cout << "finding domains' sizes\n";
440 find_domain_sizes(graph, domsizes);
442 std::cout << "finding domains\n";
443 std::vector< int > a;
445 while (check_domsizes(domsizes, domain_min_size)) {
446 domains.push_back(a);
447 if (DomainSearchingAlgorythm == 0) {
448 get_maxsized_domain(domains.back(), graph, domsizes);
449 } else if (DomainSearchingAlgorythm == 1) {
450 get_min_domain(domains.back(), graph, domsizes, domain_min_size);
452 std::cout << "new domain: " << domains.back().size() << " atoms\n";
453 delete_domain_from_graph(graph, domains.back());
455 find_domain_sizes(graph, domsizes);
457 for (int i = 0; i < bone; i++) {
464 Domains::writeOutput()
466 std::cout << "making output file\n";
467 print_domains(domains, index, fnNdx_); // see function for details | numbers from index
468 std::cout << "\n END \n";
473 * The main function for the analysis template.
476 main(int argc, char *argv[])
478 return gmx::TrajectoryAnalysisCommandLineRunner::runAsMain<Domains>(argc, argv);