votca 2024.2-dev
Loading...
Searching...
No Matches
rdf_calculator.cc
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// Standard includes
19#include <fstream>
20#include <iomanip>
21#include <iostream>
22
23// Third party includes
24#include <boost/lexical_cast.hpp>
25
26// VOTCA includes
29
30// Local VOTCA includes
31#include <votca/csg/beadlist.h>
32#include <votca/csg/imcio.h>
34
35// Local private includes
36#include "rdf_calculator.h"
37
38namespace votca {
39namespace csg {
40
42 : write_every_(0),
43 do_blocks_(false),
44 nblock_(0),
45 subvol_rad_(0),
46 do_vol_corr_(false),
47 processed_some_frames_(false) {}
48
50
51// begin the coarse graining process
52// here the data structures are prepared to handle all the data
54 // do some output
55 std::cout << "begin to calculate distribution functions\n";
56 std::cout << "# of bonded interactions: " << bonded_.size() << std::endl;
57 std::cout << "# of non-bonded interactions: " << nonbonded_.size()
58 << std::endl;
59
60 if (bonded_.size() + nonbonded_.size() == 0) {
61 throw std::runtime_error(
62 "No interactions defined in options xml-file - nothing to be done");
63 }
64
65 // initialize non-bonded structures
66 for (Property *prop : nonbonded_) {
68 i->is_bonded_ = false;
69 }
70}
71
73 Eigen::Matrix3d box = top->getBox();
74 boxc_ = box.rowwise().sum() / 2.0;
75
76 std::cout << "Using center of box: " << boxc_ << std::endl;
77 // we didn't process any frames so far
78 nframes_ = 0;
79 nblock_ = 0;
81
82 // initialize non-bonded structures
83 for (Property *prop : nonbonded_) {
84 std::string name = prop->get("name").value();
85
86 interaction_t &i = *interactions_[name];
87
88 // count total species for ideal densities
89
90 BeadList allbeads1, allbeads2;
91 allbeads1.Generate(*top, prop->get("type1").value());
92 allbeads2.Generate(*top, prop->get("type2").value());
93
94 if (allbeads1.size() == 0) {
95 throw std::runtime_error("Topology does not have beads of type \"" +
96 prop->get("type1").value() +
97 "\"\n"
98 "This was specified in type1 of interaction \"" +
99 name + "\"");
100 }
101 if (allbeads2.size() == 0) {
102 throw std::runtime_error("Topology does not have beads of type \"" +
103 prop->get("type2").value() +
104 "\"\n"
105 "This was specified in type2 of interaction \"" +
106 name + "\"");
107 }
108 // calculate normalization factor for rdf
109
110 if (do_vol_corr_) {
111 std::cout << "Volume correction on" << std::endl;
112 i.norm_ = 1. / (4.0 * votca::tools::conv::Pi * i.step_);
113 } else {
114 std::cout << "Volume correction off" << std::endl;
115 i.norm_ = 1. / (4. * votca::tools::conv::Pi * i.step_);
116 }
117 }
118}
119// create an entry for interactions
121 std::string name = p->get("name").value();
122 std::string group;
123
124 group = "none";
125
126 auto inter = std::make_unique<interaction_t>();
127 interaction_t *i = inter.get();
128 inter->index_ = interactions_.size();
129 interactions_[name] = std::move(inter);
130 getGroup(group)->interactions_.push_back(i);
131
132 i->step_ = p->get("step").as<double>();
133 i->min_ = p->get("min").as<double>();
134 i->max_ = p->get("max").as<double>();
135 i->norm_ = 1.0;
136 i->p_ = p;
137
138 // initialize the current and average histogram
139 Index n = (Index)((i->max_ - i->min_) / i->step_ + 1.000000001);
140
141 i->average_.Initialize(i->min_, i->max_ + i->step_, n);
142
143 return i;
144}
145
146// end of trajectory, post processing data
148 if (nframes_ > 0) {
149 if (!do_blocks_) {
150 WriteDist();
151 }
152 }
153 // clear interactions and groups
154 interactions_.clear();
155 groups_.clear();
157 throw std::runtime_error(
158 "no frames were processed. Please check your input");
159 }
160}
161
162// load options from xml file
163void RDFCalculator::LoadOptions(const std::string &file) {
164 options_.LoadFromXML(file);
165 bonded_ = options_.Select("cg.bonded");
166 nonbonded_ = options_.Select("cg.non-bonded");
167}
168
169// evaluate current conformation
173 // process non-bonded interactions
174 DoNonbonded(top);
175 // process bonded interactions
176 DoBonded(top);
177}
178
180
181 nframes_ = 0;
182 for (auto &interaction_ : interactions_) {
183 interaction_.second->average_.Clear();
184 }
185
186 for (auto &group_ : groups_) {
187 group_.second->corr_.setZero();
188 }
189}
190
191class IMCNBSearchHandler {
192 public:
193 IMCNBSearchHandler(HistogramNew *hist, double subvol_rad,
194 Eigen::Vector3d boxc, bool do_vol_corr)
195 : hist_(hist),
196 subvol_rad_(subvol_rad),
197 boxc_(boxc),
198 do_vol_corr_(do_vol_corr) {}
199
202 Eigen::Vector3d boxc_; // center of box
204
205 bool FoundPair(Bead *b1, Bead *, const Eigen::Vector3d &, const double dist) {
206
207 if (do_vol_corr_) {
208 double dr = (b1->Pos() - boxc_).norm();
209 if (dist + dr > subvol_rad_) {
210 // 2.0 is because everything is normalized to 4 PI
211 hist_->Process(dist, 2.0 / SurfaceRatio(dist, dr));
212 } else {
213 hist_->Process(dist);
214 }
215
216 } else {
217 hist_->Process(dist);
218 }
219 return false;
220 }
221
222 double SurfaceRatio(double dist, double r) {
223 // r: distance of particle from ex center
224 // dist: distance between particles
225 return (1.0 + (subvol_rad_ * subvol_rad_ - r * r - dist * dist) /
226 (2.0 * dist * r));
227 }
228};
229
230// process non-bonded interactions for current frame
232 for (Property *prop : rdfcalculator_->nonbonded_) {
233 std::string name = prop->get("name").value();
234
235 interaction_t &i = *rdfcalculator_->interactions_[name];
236
237 // generate the bead lists
238 BeadList beads1, beads2;
239
240 beads1.GenerateInSphericalSubvolume(*top, prop->get("type1").value(),
241 rdfcalculator_->boxc_,
242 rdfcalculator_->subvol_rad_);
243 beads2.GenerateInSphericalSubvolume(*top, prop->get("type2").value(),
244 rdfcalculator_->boxc_,
245 rdfcalculator_->subvol_rad_);
246
247 cur_beadlist_1_count_ = (double)beads1.size();
248 cur_beadlist_2_count_ = (double)beads2.size();
249
250 // same types, so put factor 1/2 because of already counted interactions
251 if (prop->get("type1").value() == prop->get("type2").value()) {
252 cur_beadlist_2_count_ /= 2.0;
253 }
254
255 // generate the neighbour list
256 std::unique_ptr<NBList> nb;
257
258 bool gridsearch = true;
259
260 if (rdfcalculator_->options_.exists("cg.nbsearch")) {
261 if (rdfcalculator_->options_.get("cg.nbsearch").as<std::string>() ==
262 "grid") {
263 gridsearch = true;
264 } else if (rdfcalculator_->options_.get("cg.nbsearch")
265 .as<std::string>() == "simple") {
266 gridsearch = false;
267 } else {
268 throw std::runtime_error("cg.nbsearch invalid, can be grid or simple");
269 }
270 }
271 if (gridsearch) {
272 nb = std::make_unique<NBListGrid>();
273 } else {
274 nb = std::make_unique<NBList>();
275 }
276
277 nb->setCutoff(i.max_ + i.step_);
278
279 // clear the current histogram
280 current_hists_[i.index_].Clear();
281
282 IMCNBSearchHandler h(&(current_hists_[i.index_]),
283 rdfcalculator_->subvol_rad_, rdfcalculator_->boxc_,
284 rdfcalculator_->do_vol_corr_);
285 nb->SetMatchFunction(&h, &IMCNBSearchHandler::FoundPair);
286
287 // is it same types or different types?
288 if (prop->get("type1").value() == prop->get("type2").value()) {
289 nb->Generate(beads1);
290 } else {
291 nb->Generate(beads1, beads2);
292 }
293
294 // store particle number in subvolume for each interaction
295 i.avg_beadlist_1_count_.Process(cur_beadlist_1_count_);
296 i.avg_beadlist_2_count_.Process(cur_beadlist_2_count_);
297 }
298}
299
300// process non-bonded interactions for current frame
302 for (Property *prop : rdfcalculator_->bonded_) {
303 std::string name = prop->get("name").value();
304
305 interaction_t &i = *rdfcalculator_->interactions_[name];
306
307 // clear the current histogram
308 current_hists_[i.index_].Clear();
309
310 // now fill with new data
311 std::vector<Interaction *> vec = top->InteractionsInGroup(name);
312
313 for (auto ic : vec) {
314 double v = ic->EvaluateVar(*top);
315 current_hists_[i.index_].Process(v);
316 }
317 }
318}
319
320// returns a group, creates it if doesn't exist
322 std::map<std::string, std::unique_ptr<group_t>>::iterator iter;
323 iter = groups_.find(name);
324 if (iter == groups_.end()) {
325 return (groups_[name] = std::make_unique<group_t>()).get();
326 }
327 return (*iter).second.get();
328}
329
330// write the distribution function
331void RDFCalculator::WriteDist(const std::string &suffix) {
332
333 // for all interactions
334 for (auto &interaction_ : interactions_) {
335 // calculate the rdf
336 Table &t = interaction_.second->average_.data();
337 Table dist(t);
338
339 interaction_.second->norm_ /=
340 (interaction_.second->avg_beadlist_1_count_.getAvg() *
341 interaction_.second->avg_beadlist_2_count_.getAvg());
342 dist.y() = avg_vol_.getAvg() * interaction_.second->norm_ *
343 dist.y().cwiseQuotient(dist.x().cwiseAbs2());
344
345 dist.Save((interaction_.first) + suffix + ".dist.new");
346 std::cout << "written " << (interaction_.first) + suffix + ".dist.new\n";
347
348 std::cout << "Avg. number of particles in subvol for "
349 << (interaction_.first) << std::endl;
350 std::cout << "beadlist 1: "
351 << interaction_.second->avg_beadlist_1_count_.getAvg()
352 << std::endl;
353 std::cout << "beadlist 2: "
354 << interaction_.second->avg_beadlist_2_count_.getAvg()
355 << std::endl;
356 }
357
358 std::cout << "Volume used for normalization: " << avg_vol_.getAvg()
359 << std::endl;
360}
361
362std::unique_ptr<CsgApplication::Worker> RDFCalculator::ForkWorker() {
363 auto worker = std::make_unique<RDFCalculator::Worker>();
364
365 worker->current_hists_.resize(interactions_.size());
366 worker->rdfcalculator_ = this;
367
368 for (auto &interaction_ : interactions_) {
369 interaction_t *i = interaction_.second.get();
370 worker->current_hists_[i->index_].Initialize(
372 }
373 return worker;
374}
375
378 RDFCalculator::Worker *worker =
379 dynamic_cast<RDFCalculator::Worker *>(worker_);
380 // update the average
381
382 ++nframes_;
383
384 avg_vol_.Process(worker->cur_vol_);
385
386 for (auto &interaction_ : interactions_) {
387 interaction_t *i = interaction_.second.get();
388 i->average_.data().y() =
389 (((double)nframes_ - 1.0) * i->average_.data().y() +
390 worker->current_hists_[i->index_].data().y()) /
391 (double)nframes_;
392 }
393
394 if (write_every_ != 0) {
395 if ((nframes_ % write_every_) == 0) {
396 nblock_++;
397 std::string suffix =
398 std::string("_") + boost::lexical_cast<std::string>(nblock_);
399 WriteDist(suffix);
400 if (do_blocks_) {
402 }
403 }
404 }
405}
406
407} // namespace csg
408} // namespace votca
virtual Eigen::Vector3d & Pos()
Definition basebead.h:128
Index Generate(Topology &top, const std::string &select)
Select all beads of type "select".
Definition beadlist.cc:30
Index size() const
Definition beadlist.h:52
Index GenerateInSphericalSubvolume(Topology &top, const std::string &select, Eigen::Vector3d ref, double radius)
Select all beads of type "select" withn a radius "radius" of reference vector "ref".
Definition beadlist.cc:57
information about a bead
Definition bead.h:50
Worker, derived from Thread, does the work.
bool FoundPair(Bead *b1, Bead *, const Eigen::Vector3d &, const double dist)
double SurfaceRatio(double dist, double r)
IMCNBSearchHandler(HistogramNew *hist, double subvol_rad, Eigen::Vector3d boxc, bool do_vol_corr)
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_
std::unique_ptr< CsgApplication::Worker > ForkWorker()
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
std::vector< Property * > bonded_
list of bonded interactions
void BeginEvaluate(Topology *top, Topology *top_atom)
begin coarse graining a trajectory
void EndEvaluate()
end coarse graining a trajectory
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 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
std::vector< Interaction * > InteractionsInGroup(const std::string &group)
Definition topology.cc:201
const Eigen::Matrix3d & getBox() const
Definition topology.h:298
void Process(const T &value)
Definition average.h:48
const T & getAvg() const
Definition average.h:83
class to generate histograms
double getMax() const
get the upper bound of the histogram intervaö
Index getNBins() const
Get number of grid points.
double getMin() const
get the lower bound of the histogram intervaö
Table & data()
get access to content of histogram
void Initialize(double min, double max, Index nbins)
Initialize the HistogramNew.
void Process(const double &v, double scale=1.0)
process a data point
class to manage program options with xml serialization functionality
Definition property.h:55
Property & get(const std::string &key)
get existing property
Definition property.cc:79
std::string & value()
reference to value of property
Definition property.h:153
T as() const
return value as type
Definition property.h:283
std::vector< Property * > Select(const std::string &filter)
select property based on a filter
Definition property.cc:185
void LoadFromXML(std::string filename)
Definition property.cc:238
class to store tables like rdfs, tabulated potentials, etc
Definition table.h:36
double & x(Index i)
Definition table.h:44
double & y(Index i)
Definition table.h:45
void Save(std::string filename) const
Definition table.cc:59
const double Pi
Definition constants.h:36
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