ac6e3808d2876929d896bb5b3bbae72902977c87
[alexxy/gromacs-fitng.git] / src / fitng.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 gmx::analysismodules::Freevolume.
38  *
39  * \author Titov Anatoly <Wapuk-cobaka@yandex.ru>
40  * \ingroup module_trajectoryanalysis
41  */
42
43 #include <iostream>
44 #include <chrono>
45 #include <omp.h>
46 #include <thread>
47 #include <string>
48
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>
56
57 #define MAX_NTRICVEC 12
58
59 using namespace gmx;
60
61 using gmx::RVec;
62
63 struct node {
64     short int n;
65     RVec r;
66     bool yep;
67 };
68
69 void make_graph(int mgwi_natoms, rvec *mgwi_x, std::vector< std::vector< node > > &mgwi_graph)
70 {
71     mgwi_graph.resize(mgwi_natoms);
72     for (int i = 0; i < mgwi_natoms; i++) {
73         mgwi_graph[i].resize(mgwi_natoms);
74     }
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;
79         }
80     }
81 }
82
83 void update_graph(std::vector< std::vector< node > > &ugwi_graph, rvec *ugwi_x, long double ugwi_epsi) {
84     rvec ugwi_temp;
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) {
91                 if (i == j) {
92                     ugwi_graph[i][j].n++;
93                 }
94                 else {
95                     ugwi_graph[i][j].n++;
96                     ugwi_graph[j][i].n++;
97                 }
98             }
99         }
100     }
101 }
102
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;
110                 }
111                 else {
112                     cd_graph[k][i][j].yep = false;
113                 }
114             }
115         }
116     }
117 }
118
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]++;
128                 }
129             }
130         }
131     }
132 }
133
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(), gmd_for3 = gmd_graph[1][1].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]) {
140                 gmd_number1 = i;
141                 gmd_number2 = j;
142             }
143         }
144     }
145     gmd_max_d.resize(0);
146     for (int i = 0; i < gmd_for3; i++) {
147         if (gmd_graph[gmd_number1][gmd_number2][i].yep) {
148             gmd_max_d.push_back(i);
149         }
150     }
151 }
152
153 void delete_domain_from_graph(std::vector< std::vector< std::vector< node > > > &ddf_graph, std::vector< int > ddf_domain) {
154     int ddfg_for1 = ddf_domain.size(), ddfg_for2 = ddf_graph.size(), ddfg_for3 = ddf_graph[1].size();
155     for (int i = 0; i < ddfg_for1; i++) {
156         for (int k = 0; k < ddfg_for2; k++) {
157             for (int j = 0; j < ddfg_for3; j++) {
158                 if (ddf_graph[k][ddf_domain[i]][j].yep) {
159                     ddf_graph[k][ddf_domain[i]][j].yep = false;
160                 }
161                 if (ddf_graph[k][j][ddf_domain[i]].yep) {
162                     ddf_graph[k][j][ddf_domain[i]].yep = false;
163                 }
164             }
165         }
166     }
167 }
168
169 bool check_domsizes(std::vector< std::vector< int > > cd_domsizes, int cd_domain_min_size) {
170     int cd_for1 = cd_domsizes.size(), cd_for2 = cd_domsizes[0].size();
171     for (int i = 0; i < cd_for1; i++) {
172         for (int j = 0; j < cd_for2; j++) {
173             if (cd_domsizes[i][j] >= cd_domain_min_size) {
174                 return true;
175             }
176         }
177     }
178     return false;
179 }
180
181 void print_domains(std::vector< std::vector< int > > pd_domains, std::vector< int > index, std::string fnNdx_) {
182     FILE                *fpNdx_;
183     fpNdx_ = std::fopen(fnNdx_.c_str(), "w+");
184     int write_count;
185     for (int i = 0; i < pd_domains.size(); i++) {
186         std::fprintf(fpNdx_, "[domain_%d]\n", i + 1);
187         write_count = 0;
188         for (int j = 0; j < pd_domains[i].size(); j++) {
189             write_count++;
190             if (write_count > 20) {
191                 write_count -= 20;
192                 std::fprintf(fpNdx_, "\n");
193             }
194             std::fprintf(fpNdx_, "%5d ", index[pd_domains[i][j]] + 1);
195         }
196         std::fprintf(fpNdx_,"\n\n");
197     }
198     std::fprintf(fpNdx_,"\n");
199     std::fclose(fpNdx_);
200 }
201
202
203 /*! \brief
204  * Class used to compute free volume in a simulations box.
205  *
206  * Inherits TrajectoryAnalysisModule and all functions from there.
207  * Does not implement any new functionality.
208  *
209  * \ingroup module_trajectoryanalysis
210  */
211 class Domains : public TrajectoryAnalysisModule
212 {
213     public:
214
215         Domains();
216         virtual ~Domains();
217
218         //! Set the options and setting
219         virtual void initOptions(IOptionsContainer          *options,
220                                  TrajectoryAnalysisSettings *settings);
221
222         //! First routine called by the analysis framework
223         // virtual void initAnalysis(const t_trxframe &fr, t_pbc *pbc);
224         virtual void initAnalysis(const TrajectoryAnalysisSettings &settings,
225                                   const TopologyInformation        &top);
226
227         virtual void initAfterFirstFrame(const TrajectoryAnalysisSettings   &settings,
228                                          const t_trxframe                   &fr);
229
230         //! Call for each frame of the trajectory
231         // virtual void analyzeFrame(const t_trxframe &fr, t_pbc *pbc);
232         virtual void analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
233                                   TrajectoryAnalysisModuleData *pdata);
234
235         //! Last routine called by the analysis framework
236         // virtual void finishAnalysis(t_pbc *pbc);
237         virtual void finishAnalysis(int nframes);
238
239         //! Routine to write output, that is additional over the built-in
240         virtual void writeOutput();
241
242     private:
243
244         std::string                                                 fnNdx_;
245
246         std::vector< std::vector< std::vector< node > > >           graph;
247
248         std::vector< std::vector< int > >                           domains;
249         std::vector< std::vector< int > >                           domsizes;
250
251         std::vector< int >                                          index;
252         std::vector< int >                                          numbers;
253         std::vector< std::vector < RVec > >                         trajectory;
254         Selection                                                   selec;
255         int                                                         frames              = 0;
256         int                                                         domain_min_size     = 5; // should be selectable
257
258         real                                                        **w_rls;
259         int                                                         bone;
260         double                                                      delta               = 0.90; //0.95 // should be selectable
261         double                                                      epsi                = 0.15; //0.3 колебания внутри домена // should be selectable
262
263         int                                                         domains_ePBC;
264         // Copy and assign disallowed by base.
265 };
266
267 Domains::Domains(): TrajectoryAnalysisModule()
268 {
269 }
270
271 Domains::~Domains()
272 {
273 }
274
275 void
276 Domains::initOptions(IOptionsContainer          *options,
277                      TrajectoryAnalysisSettings *settings)
278 {
279     static const char *const desc[] = {
280         "[THISMODULE] to be done"
281     };
282     // Add the descriptive text (program help text) to the options
283     settings->setHelpText(desc);
284     // Add option for selecting a subset of atoms
285     options->addOption(SelectionOption("select")
286                            .store(&selec).required()
287                            .description("Atoms that are considered as part of the excluded volume"));
288     // Add option for output file name
289     options->addOption(FileNameOption("on").filetype(eftIndex).outputFile()
290                             .store(&fnNdx_).defaultBasename("domains")
291                             .description("Index file from the domains"));
292     // Add option for etalon_frame constant
293     options->addOption(gmx::IntegerOption("dms")
294                             .store(&domain_min_size)
295                             .description("minimum domain size"));
296     // Add option for epsi constant
297     options->addOption(DoubleOption("epsilon")
298                             .store(&epsi)
299                             .description("thermal vibrations' constant"));
300     // Add option for delta constant
301     options->addOption(DoubleOption("delta")
302                             .store(&delta)
303                             .description("domain membership probability"));
304     // Control input settings
305     settings->setFlags(TrajectoryAnalysisSettings::efNoUserPBC);
306     settings->setPBC(true);
307 }
308
309 void
310 Domains::initAnalysis(const TrajectoryAnalysisSettings &settings,
311                       const TopologyInformation        &top)
312 {
313     domains_ePBC = top.ePBC();
314 }
315
316 void
317 Domains::initAfterFirstFrame(const TrajectoryAnalysisSettings       &settings,
318                              const t_trxframe                       &fr)
319 {
320     t_pbc  pbc;
321     t_pbc *ppbc = settings.hasPBC() ? &pbc : NULL;
322     matrix boxx;
323     copy_mat(fr.box, boxx);
324     if (ppbc != NULL) {
325         set_pbc(ppbc, domains_ePBC, boxx);
326     }
327     ConstArrayRef< int > atomind  = selec.atomIndices();
328     index.resize(0);
329     for (ConstArrayRef<int>::iterator ai = atomind.begin(); (ai < atomind.end()); ai++) {
330         index.push_back(*ai);
331     }
332     trajectory.resize(2);
333     trajectory[0].resize(selec.atomCount());
334
335     for (int i = 0; i < selec.atomCount(); i++) {
336         trajectory[0][i] = fr.x[index[i]];
337     }
338
339     bone = index.size() - domain_min_size + 1;
340     graph.resize(bone);
341     snew(w_rls, bone);
342     for (int i = 0; i < bone; i++) {
343         snew(w_rls[i], index.size());
344         for (int j = 0; j < index.size(); j++) {
345             if (j >= i && j <= i + domain_min_size - 1) {
346                 w_rls[i][j] = 1;
347             } else {
348                 w_rls[i][j] = 0;
349             }
350         }
351         rvec *etalon;
352         snew(etalon, index.size());
353         for (int j = 0; j < index.size(); j++) {
354             copy_rvec(trajectory[0][j].as_vec(), etalon[j]);
355         }
356         reset_x(index.size(), NULL, index.size(), NULL, etalon, w_rls[i]);
357         make_graph(index.size(), etalon, graph[i]);
358         sfree(etalon);
359     }
360     trajectory[1].resize(index.size());
361 }
362
363 void
364 Domains::analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
365                       TrajectoryAnalysisModuleData *pdata)
366 {
367     for (int i = 0; i < index.size(); i++) {
368         trajectory[1][i] = fr.x[index[i]];
369     }
370     frames++;
371
372     #pragma omp parallel
373     {
374         #pragma omp for schedule(dynamic)
375         for (int j = 0; j < bone; j++) {
376             rvec *etalon, *traj;
377             snew(etalon, index.size());
378             for (int k = 0; k < index.size(); k++) {
379                 copy_rvec(trajectory[0][k].as_vec(), etalon[k]);
380             }
381             snew(traj, index.size());
382             for (int k = 0; k < index.size(); k++) {
383                 copy_rvec(trajectory[1][k].as_vec(), traj[k]);
384             }
385             reset_x(index.size(), NULL, index.size(), NULL, etalon, w_rls[j]);
386             reset_x(index.size(), NULL, index.size(), NULL, traj, w_rls[j]);
387             do_fit(index.size(), w_rls[j], etalon, traj);
388             update_graph(graph[j], traj, epsi);
389             sfree(etalon);
390             sfree(traj);
391         }
392     }
393     std::cout << "frame: " << frames << " analyzed\n";
394 }
395
396 //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'
397 //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'
398 //domains -s '/home/toluk/Data/reca_rd/reca_rd.mono.tpr' -f '/home/toluk/Data/reca_rd/reca_rd.mono.xtc' - select 'name CA'
399
400 void
401 Domains::finishAnalysis(int nframes)
402 {
403     frames -= 1;
404
405     std::cout << "final cheking\n";
406     check_domains(delta, frames, graph);
407
408     std::cout << "finding domains' sizes\n";
409     find_domain_sizes(graph, domsizes);
410
411     std::cout << "finding domains\n";
412     std::vector< int > a;
413     a.resize(0);
414     while (check_domsizes(domsizes, domain_min_size)) {
415         domains.push_back(a);
416         get_maxsized_domain(domains.back(), graph, domsizes);
417         delete_domain_from_graph(graph, domains.back());
418         domsizes.resize(0);
419         find_domain_sizes(graph, domsizes);
420     }
421     for (int i = 0; i < bone; i++) {
422         sfree(w_rls[i]);
423     }
424     sfree(w_rls);
425 }
426
427 void
428 Domains::writeOutput()
429 {
430     std::cout << "making output file\n";
431     print_domains(domains, index, fnNdx_); // see function for details | numbers from index
432     std::cout << "\n END \n";
433 }
434
435
436 /*! \brief
437  * The main function for the analysis template.
438  */
439 int
440 main(int argc, char *argv[])
441 {
442     return gmx::TrajectoryAnalysisCommandLineRunner::runAsMain<Domains>(argc, argv);
443 }