votca 2024.1-dev
Loading...
Searching...
No Matches
csg_stat_imc.h
Go to the documentation of this file.
1/*
2 * Copyright 2009-2021 The VOTCA Development Team (http://www.votca.org)
3 *
4 * Licensed under the Apache License, Version 2.0 (the "License");
5 * you may not use this file except in compliance with the License.
6 * You may obtain a copy of the License at
7 *
8 * http://www.apache.org/licenses/LICENSE-2.0
9 *
10 * Unless required by applicable law or agreed to in writing, software
11 * distributed under the License is distributed on an "AS IS" BASIS,
12 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
13 * See the License for the specific language governing permissions and
14 * limitations under the License.
15 *
16 */
17
18#ifndef VOTCA_CSG_CSG_STAT_IMC_H
19#define VOTCA_CSG_CSG_STAT_IMC_H
20
21// VOTCA includes
22#include <votca/tools/average.h>
25
26// Local VOTCA includes
28
29namespace votca {
30namespace csg {
39class Imc {
40 public:
41 void Initialize(void);
42
44 void LoadOptions(const std::string &file);
45
47 void BeginEvaluate(Topology *top, Topology *top_atom);
48
50 void EndEvaluate();
51
52 void BlockLength(votca::Index length) { block_length_ = length; }
53 void DoImc(bool do_imc) { do_imc_ = do_imc; }
54 void IncludeIntra(bool include_intra) { include_intra_ = include_intra; }
55 void Extension(std::string ext) { extension_ = ext; }
56
57 protected:
59
60 using group_matrix = Eigen::MatrixXd;
61 using pair_matrix = Eigen::Block<group_matrix>;
62
76
77 // a pair of interactions which are correlated
86
88 struct group_t {
89 std::vector<interaction_t *> interactions_;
91 std::vector<pair_t> pairs_;
92 };
93
96 // length of the block to write out and averages are clear after every write
98 // calculate the inverse monte carlos parameters (cross correlations)
99 bool do_imc_ = false;
100 // include the intramolecular neighbors
101 bool include_intra_ = false;
102
103 // file extension for the distributions
104 std::string extension_;
105
106 // number of frames we processed
109
111 std::vector<tools::Property *> bonded_;
113 std::vector<tools::Property *> nonbonded_;
114
116 std::map<std::string, std::unique_ptr<interaction_t> > interactions_;
118 std::map<std::string, std::unique_ptr<group_t> > groups_;
119
121 interaction_t *AddInteraction(tools::Property *p, bool is_bonded);
122
124 group_t *getGroup(const std::string &name);
125
127 void InitializeGroups();
128
129 void WriteDist(const std::string &suffix = "");
130 void WriteIMCData(const std::string &suffix = "");
131 void WriteIMCBlock(const std::string &suffix);
132
133 void CalcDeltaS(interaction_t *interaction,
134 Eigen::VectorBlock<Eigen::VectorXd> &dS);
135
136 void ClearAverages();
137
139 public:
140 std::vector<tools::HistogramNew> current_hists_;
141 std::vector<tools::HistogramNew> current_hists_force_;
143 double cur_vol_;
144
146 void EvalConfiguration(Topology *top, Topology *top_atom) override;
148 void DoNonbonded(Topology *top);
150 void DoBonded(Topology *top);
151 };
153 void DoCorrelations(Imc::Worker *worker);
154
156
157 public:
158 std::unique_ptr<CsgApplication::Worker> ForkWorker();
159 void MergeWorker(CsgApplication::Worker *worker_);
160};
161
163 votca::Index offset_i, votca::Index offset_j,
164 const pair_matrix &corr)
165 : i1_(i1), i2_(i2), offset_i_(offset_i), offset_j_(offset_j), corr_(corr) {}
166
167} // namespace csg
168} // namespace votca
169
170#endif // VOTCA_CSG_CSG_STAT_IMC_H
Worker, derived from Thread, does the work.
void EvalConfiguration(Topology *top, Topology *top_atom) override
evaluate current conformation
std::vector< tools::HistogramNew > current_hists_
void DoNonbonded(Topology *top)
process non-bonded interactions for given frame
std::vector< tools::HistogramNew > current_hists_force_
void DoBonded(Topology *top)
process bonded interactions for given frame
class to calculate distribution functions and cross correlations for inverse monte carlo
Eigen::MatrixXd group_matrix
void BlockLength(votca::Index length)
interaction_t * AddInteraction(tools::Property *p, bool is_bonded)
create a new interaction entry based on given options
votca::Index block_length_
std::unique_ptr< CsgApplication::Worker > ForkWorker()
void WriteIMCData(const std::string &suffix="")
void EndEvaluate()
end coarse graining a trajectory
tools::Average< double > avg_vol_
tools::Property options_
the options parsed from cg definition file
void Extension(std::string ext)
void DoImc(bool do_imc)
void InitializeGroups()
initializes the group structs after interactions were added
std::map< std::string, std::unique_ptr< group_t > > groups_
map group-name to group
void IncludeIntra(bool include_intra)
std::vector< tools::Property * > nonbonded_
list of non-bonded interactions
std::string extension_
std::map< std::string, std::unique_ptr< interaction_t > > interactions_
map interaction-name to interaction
votca::Index nblock_
void MergeWorker(CsgApplication::Worker *worker_)
void DoCorrelations(Imc::Worker *worker)
update the correlations after interations were processed
void Initialize(void)
bool processed_some_frames_
votca::Index nframes_
std::vector< tools::Property * > bonded_
list of bonded interactions
void WriteDist(const std::string &suffix="")
void BeginEvaluate(Topology *top, Topology *top_atom)
begin coarse graining a trajectory
void CalcDeltaS(interaction_t *interaction, Eigen::VectorBlock< Eigen::VectorXd > &dS)
void WriteIMCBlock(const std::string &suffix)
void LoadOptions(const std::string &file)
load cg definitions file
Eigen::Block< group_matrix > pair_matrix
group_t * getGroup(const std::string &name)
get group by name, creates one if it doesn't exist
topology of the whole system
Definition topology.h:81
class to generate histograms
class to manage program options with xml serialization functionality
Definition property.h:55
base class for all analysis tools
Definition basebead.h:33
Eigen::Index Index
Definition types.h:26
struct to store collected information for groups (e.g. crosscorrelations)
std::vector< pair_t > pairs_
std::vector< interaction_t * > interactions_
struct to store collected information for interactions
tools::HistogramNew average_
tools::HistogramNew average_force_
pair_t(interaction_t *i1, interaction_t *i2, votca::Index offset_i, votca::Index offset_j, const pair_matrix &corr)
interaction_t * i1_
interaction_t * i2_