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
47 #include <gromacs/trajectoryanalysis.h>
48 #include <gromacs/pbcutil/pbc.h>
49 #include <gromacs/utility/smalloc.h>
50 #include <gromacs/math/do_fit.h>
52 #define MAX_NTRICVEC 12
64 void make_graph(int mgwi_natoms, rvec *mgwi_x, std::vector< std::vector< node > > &mgwi_graph)
66 mgwi_graph.resize(mgwi_natoms);
67 for (int i = 0; i < mgwi_natoms; i++) {
68 mgwi_graph[i].resize(mgwi_natoms);
70 for (int i = 0; i < mgwi_natoms; i++) {
71 for (int j = 0; j < mgwi_natoms; j++) {
72 rvec_sub(mgwi_x[i], mgwi_x[j], mgwi_graph[i][j].r);
73 mgwi_graph[i][j].n = 0;
78 void update_graph(std::vector< std::vector< node > > &ugwi_graph, rvec *ugwi_x, long double ugwi_epsi) {
80 int ugwi_for = ugwi_graph.size();
81 for (int i = 0; i < ugwi_for; i++) {
82 for (int j = i; j < ugwi_for; j++) {
83 rvec_sub(ugwi_x[i], ugwi_x[j], ugwi_temp);
84 rvec_dec(ugwi_temp, ugwi_graph[i][j].r.as_vec());
85 if (norm(ugwi_temp) <= ugwi_epsi) {
98 void check_domains(long double cd_delta, int cd_frames, std::vector< std::vector< std::vector< node > > > &cd_graph) {
99 int cd_for1 = cd_graph.size(), cd_for2 = cd_graph[1].size();
100 for (int k = 0; k < cd_for1; k++) {
101 for (int i = 0; i < cd_for2; i++) {
102 for (int j = 0; j < cd_for2; j++) {
103 if (cd_graph[k][i][j].n >= cd_frames * cd_delta) {
104 cd_graph[k][i][j].yep = true;
107 cd_graph[k][i][j].yep = false;
114 void find_domain_sizes(std::vector< std::vector< std::vector< node > > > fds_graph, std::vector< std::vector< int > > &fds_domsizes) {
115 fds_domsizes.resize(fds_graph.size());
116 int fds_for1 = fds_graph.size(), fds_for2 = fds_graph[1].size();
117 for (int i = 0; i < fds_for1; i++) {
118 fds_domsizes[i].resize(fds_for2, 0);
119 for (int j = 0; j < fds_for2; j++) {
120 for (int k = 0; k < fds_for2; k++) {
121 if (fds_graph[i][j][k].yep) {
122 fds_domsizes[i][j]++;
129 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) {
130 int gmd_number1 = 0, gmd_number2 = 0;
131 int gmd_for1 = gmd_domsizes.size(), gmd_for2 = gmd_domsizes[0].size();
132 for (int i = 0; i < gmd_for1; i++) {
133 for (int j = 0; j < gmd_for2; j++) {
134 if (gmd_domsizes[i][j] >= gmd_domsizes[gmd_number1][gmd_number2]) {
141 int gmd_for3 = gmd_graph[gmd_number1][gmd_number2].size();
142 for (int i = 0; i < gmd_for3; i++) {
143 if (gmd_graph[gmd_number1][gmd_number2][i].yep) {
144 gmd_max_d.push_back(i);
149 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) {
150 int gmd_number1 = 0, gmd_number2 = 0;
151 int gmd_for1 = gmd_domsizes.size(), gmd_for2 = gmd_domsizes[0].size();
152 for (int i = 0; i < gmd_for1; i++) {
153 for (int j = 0; j < gmd_for2; j++) {
154 if (gmd_domsizes[gmd_number1][gmd_number2] < gmd_min_dom_size && gmd_domsizes[i][j] >= gmd_min_dom_size) {
158 if (gmd_domsizes[i][j] <= gmd_domsizes[gmd_number1][gmd_number2] && gmd_domsizes[i][j] >= gmd_min_dom_size) {
165 int gmd_for3 = gmd_graph[gmd_number1][gmd_number2].size();
166 for (int i = 0; i < gmd_for3; i++) {
167 if (gmd_graph[gmd_number1][gmd_number2][i].yep) {
168 gmd_min_d.push_back(i);
173 void delete_domain_from_graph(std::vector< std::vector< std::vector< node > > > &ddf_graph, std::vector< int > ddf_domain) {
174 int ddfg_for1 = ddf_domain.size(), ddfg_for2 = ddf_graph.size(), ddfg_for3 = ddf_graph[1].size();
175 for (int i = 0; i < ddfg_for1; i++) {
176 for (int k = 0; k < ddfg_for2; k++) {
177 for (int j = 0; j < ddfg_for3; j++) {
178 if (ddf_graph[k][ddf_domain[i]][j].yep) {
179 ddf_graph[k][ddf_domain[i]][j].yep = false;
181 if (ddf_graph[k][j][ddf_domain[i]].yep) {
182 ddf_graph[k][j][ddf_domain[i]].yep = false;
189 bool check_domsizes(std::vector< std::vector< int > > cd_domsizes, int cd_domain_min_size) {
190 int cd_for1 = cd_domsizes.size(), cd_for2 = cd_domsizes[0].size();
191 for (int i = 0; i < cd_for1; i++) {
192 for (int j = 0; j < cd_for2; j++) {
193 if (cd_domsizes[i][j] >= cd_domain_min_size) {
201 void print_domains(std::vector< std::vector< int > > pd_domains, std::vector< int > index, std::string fnNdx_, int dms, double epsi, double delta) {
203 fpNdx_ = std::fopen(fnNdx_.c_str(), "w+");
205 for (int i = 0; i < pd_domains.size(); i++) {
206 std::fprintf(fpNdx_, "[domain_№_%3d_%3d_%4.3f_%4.3f]\n", i + 1, dms, epsi, delta);
208 for (int j = 0; j < pd_domains[i].size(); j++) {
210 if (write_count > 20) {
212 std::fprintf(fpNdx_, "\n");
214 std::fprintf(fpNdx_, "%5d ", index[pd_domains[i][j]] + 1);
216 std::fprintf(fpNdx_,"\n\n");
218 std::fprintf(fpNdx_,"\n");
223 * \ingroup module_trajectoryanalysis
225 class Domains : public TrajectoryAnalysisModule
232 //! Set the options and setting
233 virtual void initOptions(IOptionsContainer *options,
234 TrajectoryAnalysisSettings *settings);
236 //! First routine called by the analysis framework
237 // virtual void initAnalysis(const t_trxframe &fr, t_pbc *pbc);
238 virtual void initAnalysis(const TrajectoryAnalysisSettings &settings,
239 const TopologyInformation &top);
241 virtual void initAfterFirstFrame(const TrajectoryAnalysisSettings &settings,
242 const t_trxframe &fr);
244 //! Call for each frame of the trajectory
245 // virtual void analyzeFrame(const t_trxframe &fr, t_pbc *pbc);
246 virtual void analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
247 TrajectoryAnalysisModuleData *pdata);
249 //! Last routine called by the analysis framework
250 // virtual void finishAnalysis(t_pbc *pbc);
251 virtual void finishAnalysis(int nframes);
253 //! Routine to write output, that is additional over the built-in
254 virtual void writeOutput();
260 std::vector< std::vector< std::vector< node > > > graph;
262 std::vector< std::vector< int > > domains;
263 std::vector< std::vector< int > > domsizes;
265 std::vector< int > index;
266 std::vector< int > numbers;
267 std::vector< std::vector < RVec > > trajectory;
270 int domain_min_size = 5; // selectable
274 double delta = 0.90; // selectable
275 double epsi = 0.15; // selectable
278 int DomainSearchingAlgorythm = 0; // selectable
279 // Copy and assign disallowed by base.
282 Domains::Domains(): TrajectoryAnalysisModule()
291 Domains::initOptions(IOptionsContainer *options,
292 TrajectoryAnalysisSettings *settings)
294 static const char *const desc[] = {
295 "[THISMODULE] to be done"
297 // Add the descriptive text (program help text) to the options
298 settings->setHelpText(desc);
299 // Add option for selecting a subset of atoms
300 options->addOption(SelectionOption("select")
301 .store(&selec).required()
302 .description("Atoms that are considered as part of the excluded volume"));
303 // Add option for output file name
304 options->addOption(FileNameOption("on").filetype(eftIndex).outputFile()
305 .store(&fnNdx_).defaultBasename("domains")
306 .description("Index file from the domains"));
307 // Add option for domain min size constant
308 options->addOption(gmx::IntegerOption("dms")
309 .store(&domain_min_size)
310 .description("minimum domain size"));
311 // Add option for Domain's Searching Algorythm
312 options->addOption(gmx::IntegerOption("DSA")
313 .store(&DomainSearchingAlgorythm)
314 .description("Domain's Searching Algorythm: 0 == default (from bigger to smaller) | 1 == (from smaller to bigger)"));
315 // Add option for epsi constant
316 options->addOption(DoubleOption("epsilon")
318 .description("thermal vibrations' constant"));
319 // Add option for delta constant
320 options->addOption(DoubleOption("delta")
322 .description("domain membership probability"));
323 // Control input settings
324 settings->setFlags(TrajectoryAnalysisSettings::efNoUserPBC);
325 settings->setPBC(true);
329 Domains::initAnalysis(const TrajectoryAnalysisSettings &settings,
330 const TopologyInformation &top)
332 domains_ePBC = top.ePBC();
336 Domains::initAfterFirstFrame(const TrajectoryAnalysisSettings &settings,
337 const t_trxframe &fr)
340 t_pbc *ppbc = settings.hasPBC() ? &pbc : NULL;
342 copy_mat(fr.box, boxx);
344 set_pbc(ppbc, domains_ePBC, boxx);
346 ConstArrayRef< int > atomind = selec.atomIndices();
348 for (ConstArrayRef<int>::iterator ai = atomind.begin(); (ai < atomind.end()); ai++) {
349 index.push_back(*ai);
351 trajectory.resize(2);
352 trajectory[0].resize(selec.atomCount());
354 for (int i = 0; i < selec.atomCount(); i++) {
355 trajectory[0][i] = fr.x[index[i]];
358 bone = index.size() - domain_min_size + 1;
361 for (int i = 0; i < bone; i++) {
362 snew(w_rls[i], index.size());
363 for (int j = 0; j < index.size(); j++) {
364 if (j >= i && j <= i + domain_min_size - 1) {
371 snew(etalon, index.size());
372 for (int j = 0; j < index.size(); j++) {
373 copy_rvec(trajectory[0][j].as_vec(), etalon[j]);
375 reset_x(index.size(), NULL, index.size(), NULL, etalon, w_rls[i]);
376 make_graph(index.size(), etalon, graph[i]);
379 trajectory[1].resize(index.size());
383 Domains::analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
384 TrajectoryAnalysisModuleData *pdata)
386 for (int i = 0; i < index.size(); i++) {
387 trajectory[1][i] = fr.x[index[i]];
393 #pragma omp for schedule(dynamic)
394 for (int j = 0; j < bone; j++) {
396 snew(etalon, index.size());
397 for (int k = 0; k < index.size(); k++) {
398 copy_rvec(trajectory[0][k].as_vec(), etalon[k]);
400 snew(traj, index.size());
401 for (int k = 0; k < index.size(); k++) {
402 copy_rvec(trajectory[1][k].as_vec(), traj[k]);
404 reset_x(index.size(), NULL, index.size(), NULL, etalon, w_rls[j]);
405 reset_x(index.size(), NULL, index.size(), NULL, traj, w_rls[j]);
406 do_fit(index.size(), w_rls[j], etalon, traj);
407 update_graph(graph[j], traj, epsi);
412 std::cout << "frame: " << frames << " analyzed\n";
416 Domains::finishAnalysis(int nframes)
420 std::cout << "final cheking\n";
421 check_domains(delta, frames, graph);
423 std::cout << "finding domains' sizes\n";
424 find_domain_sizes(graph, domsizes);
426 std::cout << "finding domains\n";
427 std::vector< int > a;
429 while (check_domsizes(domsizes, domain_min_size)) {
430 domains.push_back(a);
431 if (DomainSearchingAlgorythm == 0) {
432 get_maxsized_domain(domains.back(), graph, domsizes);
433 } else if (DomainSearchingAlgorythm == 1) {
434 get_min_domain(domains.back(), graph, domsizes, domain_min_size);
436 std::cout << "new domain: " << domains.back().size() << " atoms\n";
437 delete_domain_from_graph(graph, domains.back());
439 find_domain_sizes(graph, domsizes);
441 for (int i = 0; i < bone; i++) {
448 Domains::writeOutput()
450 std::cout << "making output file\n";
451 print_domains(domains, index, fnNdx_, domain_min_size, epsi, delta); // see function for details | numbers from index
452 std::cout << "\n END \n";
456 * The main function for the analysis template.
459 main(int argc, char *argv[])
461 return gmx::TrajectoryAnalysisCommandLineRunner::runAsMain<Domains>(argc, argv);