Loading [MathJax]/extensions/TeX/AMSmath.js
ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages Concepts
Analysis.cpp
Go to the documentation of this file.
1/*
2 * Copyright (C) 2013-2022 The ESPResSo project
3 *
4 * This file is part of ESPResSo.
5 *
6 * ESPResSo is free software: you can redistribute it and/or modify
7 * it under the terms of the GNU General Public License as published by
8 * the Free Software Foundation, either version 3 of the License, or
9 * (at your option) any later version.
10 *
11 * ESPResSo is distributed in the hope that it will be useful,
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 * GNU General Public License for more details.
15 *
16 * You should have received a copy of the GNU General Public License
17 * along with this program. If not, see <http://www.gnu.org/licenses/>.
18 */
19
20#include "Analysis.hpp"
21#include "ObservableStat.hpp"
22
23#include "core/BoxGeometry.hpp"
27#include "core/cells.hpp"
29#include "core/dpd.hpp"
31#include "core/npt.hpp"
32
34
35#include <utils/Vector.hpp>
36#include <utils/contains.hpp>
39
40#include <boost/mpi/collectives/all_reduce.hpp>
41
42#include <algorithm>
43#include <cmath>
44#include <functional>
45#include <memory>
46#include <sstream>
47#include <stdexcept>
48#include <string>
49#include <vector>
50
51namespace ScriptInterface {
52namespace Analysis {
53
54/** @brief Check if a contiguous range of particle ids exists. */
55static void check_topology(CellStructure const &cell_structure, int chain_start,
56 int chain_length, int n_chains) {
57 try {
58 if (n_chains <= 0) {
59 throw std::domain_error("Chain analysis needs at least 1 chain");
60 }
61 if (chain_length <= 0) {
62 throw std::domain_error("Chain analysis needs at least 1 bead per chain");
63 }
64
65 auto n_particles_local = 0;
66 for (int i = 0; i < n_chains; ++i) {
67 for (int j = 0; j < chain_length; ++j) {
68 auto const pid = chain_start + i * chain_length + j;
69 auto ptr = cell_structure.get_local_particle(pid);
70 if (ptr != nullptr and not ptr->is_ghost()) {
72 }
73 }
74 }
75 auto const n_particles_total = boost::mpi::all_reduce(
76 ::comm_cart, n_particles_local, std::plus<int>{});
77
79 for (int i = 0; i < chain_length * n_chains; ++i) {
80 auto const pid = chain_start + i;
81 auto ptr = cell_structure.get_local_particle(pid);
82 int local_count = 0;
83 if (ptr != nullptr and not ptr->is_ghost()) {
84 local_count = 1;
85 }
86 auto const total_count =
87 boost::mpi::all_reduce(::comm_cart, local_count, std::plus<int>{});
88 if (total_count == 0) {
89 std::stringstream error_msg;
90 error_msg << "Particle with id " << pid << " does not exist; "
91 << "cannot perform analysis on the range chain_start="
92 << chain_start << ", number_of_chains=" << n_chains
93 << ", chain_length=" << chain_length << ". "
94 << "Please provide a contiguous range of particle ids.";
95 throw std::runtime_error(error_msg.str());
96 }
97 }
98 }
99 } catch (...) {
100 if (::comm_cart.rank() == 0) {
101 throw;
102 }
103 throw Exception("");
104 }
105}
106
107void Analysis::check_particle_type(int p_type) const {
108 auto const &nonbonded_ias = get_system().nonbonded_ias;
109 if (p_type < 0 or p_type > nonbonded_ias->get_max_seen_particle_type()) {
110 std::stringstream error_msg;
111 error_msg << "Particle type " << p_type << " does not exist";
112 throw std::invalid_argument(error_msg.str());
113 }
114}
115
116Variant Analysis::do_call_method(std::string const &name,
117 VariantMap const &parameters) {
118 if (name == "linear_momentum") {
119 auto const local = calc_linear_momentum(
120 get_system(), get_value_or<bool>(parameters, "include_particles", true),
121 get_value_or<bool>(parameters, "include_lbfluid", true));
122 return mpi_reduce_sum(context()->get_comm(), local).as_vector();
123 }
124 if (name == "particle_energy") {
125 auto &system = get_system();
126 auto const pid = get_value<int>(parameters, "pid");
127 auto const local = system.particle_short_range_energy_contribution(pid);
128 return mpi_reduce_sum(context()->get_comm(), local);
129 }
130 if (name == "particle_bond_energy") {
131 auto &system = get_system();
132 auto const pid = get_value<int>(parameters, "pid");
133 auto const bond_id = get_value<int>(parameters, "bond_id");
134 auto const partners = get_value<std::vector<int>>(parameters, "partners");
135 auto const local = system.particle_bond_energy(pid, bond_id, partners);
136 return Utils::Mpi::reduce_optional(context()->get_comm(), local);
137 }
138#ifdef DIPOLE_FIELD_TRACKING
139 if (name == "calc_long_range_fields") {
140 get_system().calculate_long_range_fields();
141 return {};
142 }
143#endif
144 if (name == "particle_neighbor_pids") {
145 auto &system = get_system();
146 system.on_observable_calc();
147 std::unordered_map<int, std::vector<int>> dict;
148 context()->parallel_try_catch([&]() {
149 auto neighbor_pids = get_neighbor_pids(system);
150 Utils::Mpi::gather_buffer(neighbor_pids, context()->get_comm());
151 std::ranges::for_each(neighbor_pids, [&dict](auto const &nbhood) {
152 dict[nbhood.pid] = nbhood.neighbor_pids;
153 });
154 });
156 }
157#ifdef DPD
158 if (name == "dpd_stress") {
159 auto const result = dpd_stress(context()->get_comm());
160 return result.as_vector();
161 }
162#endif // DPD
163 if (name == "min_dist") {
164 auto const p_types1 = get_value<std::vector<int>>(parameters, "p_types1");
165 auto const p_types2 = get_value<std::vector<int>>(parameters, "p_types2");
166 for (auto const p_type : p_types1) {
167 context()->parallel_try_catch([&]() { check_particle_type(p_type); });
168 }
169 for (auto const p_type : p_types2) {
170 context()->parallel_try_catch([&]() { check_particle_type(p_type); });
171 }
173 }
174 if (name == "center_of_mass") {
175 auto const p_type = get_value<int>(parameters, "p_type");
176 context()->parallel_try_catch([&]() { check_particle_type(p_type); });
177 auto const local = center_of_mass(get_system(), p_type);
178 return mpi_reduce_sum(context()->get_comm(), local).as_vector();
179 }
180 if (name == "angular_momentum") {
181 auto const p_type = get_value<int>(parameters, "p_type");
182 context()->parallel_try_catch([&]() { check_particle_type(p_type); });
183 auto const local = angular_momentum(get_system(), p_type);
184 return mpi_reduce_sum(context()->get_comm(), local).as_vector();
185 }
186 if (name == "nbhood") {
187 auto const pos = get_value<Utils::Vector3d>(parameters, "pos");
188 auto const radius = get_value<double>(parameters, "r_catch");
189 auto const result = nbhood(get_system(), pos, radius);
190 return result;
191 }
192 if (name == "calc_re") {
193 auto const &system = get_system();
194 auto const chain_start = get_value<int>(parameters, "chain_start");
195 auto const chain_length = get_value<int>(parameters, "chain_length");
196 auto const n_chains = get_value<int>(parameters, "number_of_chains");
198 auto const result = calc_re(system, chain_start, chain_length, n_chains);
199 return std::vector<double>(result.begin(), result.end());
200 }
201 if (name == "calc_rg") {
202 auto const &system = get_system();
203 auto const chain_start = get_value<int>(parameters, "chain_start");
204 auto const chain_length = get_value<int>(parameters, "chain_length");
205 auto const n_chains = get_value<int>(parameters, "number_of_chains");
208 context()->parallel_try_catch([&]() {
209 auto const result = calc_rg(system, chain_start, chain_length, n_chains);
210 output = Variant{std::vector<double>(result.begin(), result.end())};
211 });
212 return output;
213 }
214 if (name == "calc_rh") {
215 auto const &system = get_system();
216 auto const chain_start = get_value<int>(parameters, "chain_start");
217 auto const chain_length = get_value<int>(parameters, "chain_length");
218 auto const n_chains = get_value<int>(parameters, "number_of_chains");
220 auto const result = calc_rh(system, chain_start, chain_length, n_chains);
221 return std::vector<double>(result.begin(), result.end());
222 }
223 if (name == "gyration_tensor") {
224 auto const p_types = get_value<std::vector<int>>(parameters, "p_types");
225 for (auto const p_type : p_types) {
226 context()->parallel_try_catch([&]() { check_particle_type(p_type); });
227 }
228 auto const mat = gyration_tensor(get_system(), p_types);
229 return std::vector<double>(mat.begin(), mat.end());
230 }
231 if (name == "moment_of_inertia_matrix") {
232 auto const p_type = get_value<int>(parameters, "p_type");
233 context()->parallel_try_catch([&]() { check_particle_type(p_type); });
234 auto const local = moment_of_inertia_matrix(get_system(), p_type);
235 return mpi_reduce_sum(context()->get_comm(), local).as_vector();
236 }
237 if (name == "structure_factor") {
238 auto const order = get_value<int>(parameters, "sf_order");
239 auto const p_types = get_value<std::vector<int>>(parameters, "sf_types");
241 if (order < 1)
242 throw std::domain_error("order has to be a strictly positive number");
243 });
244 for (auto const p_type : p_types) {
245 context()->parallel_try_catch([&]() { check_particle_type(p_type); });
246 }
247 auto const result = structure_factor(get_system(), p_types, order);
248 return make_vector_of_variants(result);
249 }
250 if (name == "distribution") {
251 auto const &box_l = get_system().box_geo->length();
252 auto const r_max_limit =
253 0.5 * std::min(std::min(box_l[0], box_l[1]), box_l[2]);
254 auto const r_min = get_value_or<double>(parameters, "r_min", 0.);
255 auto const r_max = get_value_or<double>(parameters, "r_max", r_max_limit);
256 auto const r_bins = get_value_or<int>(parameters, "r_bins", 100);
257 auto const log_flag = get_value_or<bool>(parameters, "log_flag", false);
258 auto const int_flag = get_value_or<bool>(parameters, "int_flag", false);
259 context()->parallel_try_catch([=]() {
260 if (log_flag and r_min <= 0.) {
261 throw std::domain_error("Parameter 'r_min' must be > 0");
262 }
263 if (r_min < 0.) {
264 throw std::domain_error("Parameter 'r_min' must be >= 0");
265 }
266 if (r_min >= r_max) {
267 throw std::domain_error("Parameter 'r_max' must be > 'r_min'");
268 }
269 if (r_max > r_max_limit) {
270 throw std::domain_error("Parameter 'r_max' must be <= box_l / 2");
271 }
272 if (r_bins <= 0) {
273 throw std::domain_error("Parameter 'r_bins' must be >= 1");
274 }
275 });
276 auto const p_types1 =
277 get_value<std::vector<int>>(parameters, "type_list_a");
278 auto const p_types2 =
279 get_value<std::vector<int>>(parameters, "type_list_b");
280 for (auto const p_type : p_types1) {
281 context()->parallel_try_catch([&]() { check_particle_type(p_type); });
282 }
283 for (auto const p_type : p_types2) {
284 context()->parallel_try_catch([&]() { check_particle_type(p_type); });
285 }
289 }
290 if (name == "calculate_energy") {
291 return m_obs_stat->do_call_method("calculate_energy", {});
292 }
293 if (name == "calculate_scalar_pressure") {
294 return m_obs_stat->do_call_method("calculate_scalar_pressure", {});
295 }
296 if (name == "calculate_pressure_tensor") {
297 return m_obs_stat->do_call_method("calculate_pressure_tensor", {});
298 }
299#ifdef NPT
300 if (name == "get_instantaneous_pressure") {
301 return get_system().npt_inst_pressure->p_inst[0];
302 }
303 if (name == "get_instantaneous_pressure_virial") {
304 return get_system().npt_inst_pressure->p_inst[1];
305 }
306#endif // NPT
307 return {};
308}
309
310} // namespace Analysis
311} // namespace ScriptInterface
Vector implementation and trait types for boost qvm interoperability.
std::vector< NeighborPIDs > get_neighbor_pids(System::System const &system)
Returns pairs of particle ids and neighbor particle id lists.
Definition cells.cpp:188
This file contains everything related to the global cell structure / cell system.
Variant do_call_method(std::string const &name, VariantMap const &parameters) override
Definition Analysis.cpp:116
virtual void parallel_try_catch(std::function< void()> const &cb) const =0
boost::string_ref name() const
Context * context() const
Responsible context.
boost::mpi::communicator comm_cart
The communicator.
This file contains the asynchronous MPI communication.
Utils::Vector9d dpd_stress(boost::mpi::communicator const &comm)
Viscous stress tensor of the DPD interaction.
Definition dpd.cpp:179
Routines to use DPD as thermostat or pair force .
static void check_topology(CellStructure const &cell_structure, int chain_start, int chain_length, int n_chains)
Check if a contiguous range of particle ids exists.
Definition Analysis.cpp:55
T get_value(Variant const &v)
Extract value of specific type T from a Variant.
std::unordered_map< std::string, Variant > VariantMap
Definition Variant.hpp:69
auto make_unordered_map_of_variants(std::unordered_map< K, V > const &v)
Definition Variant.hpp:80
T mpi_reduce_sum(boost::mpi::communicator const &comm, T const &result)
Reduce object by sum on the head node.
auto make_vector_of_variants(std::vector< T > const &v)
Definition Variant.hpp:88
boost::make_recursive_variant< None, bool, int, std::size_t, double, std::string, ObjectRef, Utils::Vector3b, Utils::Vector3i, Utils::Vector2d, Utils::Vector3d, Utils::Vector4d, std::vector< int >, std::vector< double >, std::vector< boost::recursive_variant_ >, std::unordered_map< int, boost::recursive_variant_ >, std::unordered_map< std::string, boost::recursive_variant_ > >::type Variant
Possible types for parameters.
Definition Variant.hpp:67
void gather_buffer(std::vector< T, Allocator > &buffer, boost::mpi::communicator const &comm, int root=0)
Gather buffer with different size on each node.
T reduce_optional(boost::mpi::communicator const &comm, std::optional< T > const &result)
Reduce an optional on the head node.
Various procedures concerning interactions between particles.
Exports for the NpT code.
Utils::Vector3d center_of_mass(System::System const &system, int p_type)
Calculate the center of mass of particles of a certain type.
Utils::Vector3d angular_momentum(System::System const &system, int p_type)
Calculate the angular momentum of particles of a certain type.
std::vector< int > nbhood(System::System const &system, Utils::Vector3d const &pos, double dist)
Find all particles within a given radius dist around a position pos.
Utils::Vector9d moment_of_inertia_matrix(System::System const &system, int p_type)
Calculate the moment of inertia of particles of a certain type.
Utils::Vector9d gyration_tensor(System::System const &system, std::vector< int > const &p_types)
Calculate the gyration tensor of particles of certain types.
Utils::Vector3d calc_linear_momentum(System::System const &system, bool include_particles, bool include_lbfluid)
Calculate total momentum of the system (particles & LB fluid).
std::vector< std::vector< double > > structure_factor(System::System const &system, std::vector< int > const &p_types, int order)
Calculate the spherically averaged structure factor.
double mindist(System::System const &system, std::vector< int > const &set1, std::vector< int > const &set2)
Calculate the minimal distance of two particles with types in set1 and set2, respectively.
std::vector< std::vector< double > > calc_part_distribution(System::System const &system, std::vector< int > const &p1_types, std::vector< int > const &p2_types, double r_min, double r_max, int r_bins, bool log_flag, bool int_flag)
Calculate the distribution of particles around others.
Statistical tools to analyze simulations.
std::array< double, 4 > calc_rg(System::System const &system, int chain_start, int chain_length, int n_chains)
Calculate the radius of gyration.
std::array< double, 2 > calc_rh(System::System const &system, int chain_start, int chain_length, int n_chains)
Calculate the hydrodynamic radius (ref.
std::array< double, 4 > calc_re(System::System const &system, int chain_start, int chain_length, int n_chains)
Calculate the end-to-end-distance.
This file contains the code for statistics on chains.
Describes a cell structure / cell system.
Particle * get_local_particle(int id)
Get a local particle by id.