votca 2024.2-dev
Loading...
Searching...
No Matches
rdf_calculator.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_RDF_CALCULATOR_H
19#define VOTCA_CSG_RDF_CALCULATOR_H
20
21// Standard includes
22#include <cmath>
23#include <memory>
24
25// Third party includes
26#include <boost/numeric/ublas/io.hpp>
27#include <boost/numeric/ublas/matrix.hpp>
28#include <boost/numeric/ublas/matrix_proxy.hpp>
29#include <boost/numeric/ublas/symmetric.hpp>
30#include <boost/numeric/ublas/vector.hpp>
31#include <boost/numeric/ublas/vector_proxy.hpp>
32
33// VOTCA includes
34#include <votca/tools/average.h>
37
38// Local VOTCA includes
40
41namespace votca {
42namespace csg {
43using namespace votca::tools;
44
54 public:
57
58 void Initialize(void);
59
61 void LoadOptions(const std::string &file);
62
64 void BeginEvaluate(Topology *top, Topology *top_atom);
65
67 void EndEvaluate();
68
69 void WriteEvery(Index write_every) { write_every_ = write_every; }
70 void DoBlocks(bool do_blocks) { do_blocks_ = do_blocks; }
71 void DoVolumeCorrection(bool do_vol_corr) { do_vol_corr_ = do_vol_corr; }
72 void SetSubvolRadius(double r) { subvol_rad_ = r; }
73 double AnalyticVolumeCorrection(double t) {
74
75 std::cout << "DBG " << t << " "
76 << 1.0 / 24.0 *
77 (16.0 * t * t - 12.0 * t * t * t + t * t * t * t * t)
78 << std::endl;
79 return 1.0 / 24.0 * (16.0 * t * t - 12.0 * t * t * t + t * t * t * t * t);
80 }
81
82 protected:
84
85 using group_matrix = Eigen::MatrixXd;
86 using pair_matrix = Eigen::Block<group_matrix>;
87
99
100 // a pair of interactions which are correlated
109
111 struct group_t {
112 std::list<interaction_t *> interactions_;
114 std::vector<pair_t> pairs_;
115 };
116
119 // we want to write out every so many frames
121 // we want do do block averaging -> clear averagings every write out
123
124 // number of frames we processed
128 Eigen::Vector3d boxc_; // center of box
130
132 std::vector<Property *> bonded_;
134 std::vector<Property *> nonbonded_;
135
137 std::map<std::string, std::unique_ptr<interaction_t>> interactions_;
139 std::map<std::string, std::unique_ptr<group_t>> groups_;
140
143
145 group_t *getGroup(const std::string &name);
146
147 void WriteDist(const std::string &suffix = "");
148
149 void ClearAverages();
150
152 public:
153 std::vector<HistogramNew> current_hists_;
155 double cur_vol_;
156 double cur_beadlist_1_count_; // need to normalize to avg density for
157 // subvol
159
161 void EvalConfiguration(Topology *top, Topology *top_atom) override;
163 void DoNonbonded(Topology *top);
165 void DoBonded(Topology *top);
166 };
169
171
172 public:
173 std::unique_ptr<CsgApplication::Worker> ForkWorker();
174 void MergeWorker(CsgApplication::Worker *worker_);
175};
176
179 Index offset_i, Index offset_j,
180 const pair_matrix &corr)
181 : i1_(i1), i2_(i2), offset_i_(offset_i), offset_j_(offset_j), corr_(corr) {}
182
183} // namespace csg
184} // namespace votca
185
186#endif // VOTCA_CSG_RDF_CALCULATOR_H
Worker, derived from Thread, does the work.
void DoNonbonded(Topology *top)
process non-bonded interactions for given frame
void DoBonded(Topology *top)
process bonded interactions for given frame
void EvalConfiguration(Topology *top, Topology *top_atom) override
evaluate current conformation
std::vector< HistogramNew > current_hists_
class to calculate distribution functions and cross correlations for inverse monte carlo
Eigen::Block< group_matrix > pair_matrix
std::unique_ptr< CsgApplication::Worker > ForkWorker()
double AnalyticVolumeCorrection(double t)
std::vector< Property * > nonbonded_
list of non-bonded interactions
Average< double > avg_vol_
interaction_t * AddInteraction(Property *p)
create a new interaction entry based on given options
void SetSubvolRadius(double r)
void DoVolumeCorrection(bool do_vol_corr)
void DoBlocks(bool do_blocks)
std::vector< Property * > bonded_
list of bonded interactions
Eigen::MatrixXd group_matrix
void BeginEvaluate(Topology *top, Topology *top_atom)
begin coarse graining a trajectory
void EndEvaluate()
end coarse graining a trajectory
void WriteEvery(Index write_every)
std::map< std::string, std::unique_ptr< interaction_t > > interactions_
std::map ineteractionm-name to interaction
void WriteDist(const std::string &suffix="")
Property options_
the options parsed from cg definition file
std::map< std::string, std::unique_ptr< group_t > > groups_
std::map group-name to group
void DoCorrelations(RDFCalculator::Worker *worker)
update the correlations after interations were processed
void LoadOptions(const std::string &file)
load cg definitions file
void MergeWorker(CsgApplication::Worker *worker_)
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::list< interaction_t * > interactions_
struct to store collected information for interactions
pair_t(interaction_t *i1, interaction_t *i2, Index offset_i, Index offset_j, const pair_matrix &corr)