ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
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 if (name == "particle_neighbor_pids") {
139 auto &system = get_system();
140 system.on_observable_calc();
141 std::unordered_map<int, std::vector<int>> dict;
142 context()->parallel_try_catch([&]() {
143 auto neighbor_pids = get_neighbor_pids(system);
144 Utils::Mpi::gather_buffer(neighbor_pids, context()->get_comm());
145 std::ranges::for_each(neighbor_pids, [&dict](auto const &nbhood) {
146 dict[nbhood.pid] = nbhood.neighbor_pids;
147 });
148 });
150 }
151#ifdef ESPRESSO_DPD
152 if (name == "dpd_stress") {
153 auto const result = dpd_stress(context()->get_comm());
154 return result.as_vector();
155 }
156#endif // ESPRESSO_DPD
157 if (name == "min_dist") {
158 auto const p_types1 = get_value<std::vector<int>>(parameters, "p_types1");
159 auto const p_types2 = get_value<std::vector<int>>(parameters, "p_types2");
160 for (auto const p_type : p_types1) {
161 context()->parallel_try_catch([&]() { check_particle_type(p_type); });
162 }
163 for (auto const p_type : p_types2) {
164 context()->parallel_try_catch([&]() { check_particle_type(p_type); });
165 }
167 }
168 if (name == "center_of_mass") {
169 auto const p_type = get_value<int>(parameters, "p_type");
170 context()->parallel_try_catch([&]() { check_particle_type(p_type); });
171 auto const local = center_of_mass(get_system(), p_type);
172 return mpi_reduce_sum(context()->get_comm(), local).as_vector();
173 }
174 if (name == "angular_momentum") {
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 = angular_momentum(get_system(), p_type);
178 return mpi_reduce_sum(context()->get_comm(), local).as_vector();
179 }
180 if (name == "nbhood") {
181 auto const pos = get_value<Utils::Vector3d>(parameters, "pos");
182 auto const radius = get_value<double>(parameters, "r_catch");
183 auto const result = nbhood(get_system(), pos, radius);
184 return result;
185 }
186 if (name == "calc_re") {
187 auto const &system = get_system();
188 auto const chain_start = get_value<int>(parameters, "chain_start");
189 auto const chain_length = get_value<int>(parameters, "chain_length");
190 auto const n_chains = get_value<int>(parameters, "number_of_chains");
192 auto const result = calc_re(system, chain_start, chain_length, n_chains);
193 return std::vector<double>(result.begin(), result.end());
194 }
195 if (name == "calc_rg") {
196 auto const &system = get_system();
197 auto const chain_start = get_value<int>(parameters, "chain_start");
198 auto const chain_length = get_value<int>(parameters, "chain_length");
199 auto const n_chains = get_value<int>(parameters, "number_of_chains");
202 context()->parallel_try_catch([&]() {
203 auto const result = calc_rg(system, chain_start, chain_length, n_chains);
204 output = Variant{std::vector<double>(result.begin(), result.end())};
205 });
206 return output;
207 }
208 if (name == "calc_rh") {
209 auto const &system = get_system();
210 auto const chain_start = get_value<int>(parameters, "chain_start");
211 auto const chain_length = get_value<int>(parameters, "chain_length");
212 auto const n_chains = get_value<int>(parameters, "number_of_chains");
214 auto const result = calc_rh(system, chain_start, chain_length, n_chains);
215 return std::vector<double>(result.begin(), result.end());
216 }
217 if (name == "gyration_tensor") {
218 auto const p_types = get_value<std::vector<int>>(parameters, "p_types");
219 for (auto const p_type : p_types) {
220 context()->parallel_try_catch([&]() { check_particle_type(p_type); });
221 }
222 auto const mat = gyration_tensor(get_system(), p_types);
223 return std::vector<double>(mat.begin(), mat.end());
224 }
225 if (name == "moment_of_inertia_matrix") {
226 auto const p_type = get_value<int>(parameters, "p_type");
227 context()->parallel_try_catch([&]() { check_particle_type(p_type); });
228 auto const local = moment_of_inertia_matrix(get_system(), p_type);
229 return mpi_reduce_sum(context()->get_comm(), local).as_vector();
230 }
231 if (name == "structure_factor") {
232 auto const order = get_value<int>(parameters, "sf_order");
233 auto const p_types = get_value<std::vector<int>>(parameters, "sf_types");
235 if (order < 1)
236 throw std::domain_error("order has to be a strictly positive number");
237 });
238 for (auto const p_type : p_types) {
239 context()->parallel_try_catch([&]() { check_particle_type(p_type); });
240 }
241 auto const result = structure_factor(get_system(), p_types, order);
242 return make_vector_of_variants(result);
243 }
244 if (name == "distribution") {
245 auto const &box_l = get_system().box_geo->length();
246 auto const r_max_limit =
247 0.5 * std::min(std::min(box_l[0], box_l[1]), box_l[2]);
248 auto const r_min = get_value_or<double>(parameters, "r_min", 0.);
249 auto const r_max = get_value_or<double>(parameters, "r_max", r_max_limit);
250 auto const r_bins = get_value_or<int>(parameters, "r_bins", 100);
251 auto const log_flag = get_value_or<bool>(parameters, "log_flag", false);
252 auto const int_flag = get_value_or<bool>(parameters, "int_flag", false);
253 context()->parallel_try_catch([=]() {
254 if (log_flag and r_min <= 0.) {
255 throw std::domain_error("Parameter 'r_min' must be > 0");
256 }
257 if (r_min < 0.) {
258 throw std::domain_error("Parameter 'r_min' must be >= 0");
259 }
260 if (r_min >= r_max) {
261 throw std::domain_error("Parameter 'r_max' must be > 'r_min'");
262 }
263 if (r_max > r_max_limit) {
264 throw std::domain_error("Parameter 'r_max' must be <= box_l / 2");
265 }
266 if (r_bins <= 0) {
267 throw std::domain_error("Parameter 'r_bins' must be >= 1");
268 }
269 });
270 auto const p_types1 =
271 get_value<std::vector<int>>(parameters, "type_list_a");
272 auto const p_types2 =
273 get_value<std::vector<int>>(parameters, "type_list_b");
274 for (auto const p_type : p_types1) {
275 context()->parallel_try_catch([&]() { check_particle_type(p_type); });
276 }
277 for (auto const p_type : p_types2) {
278 context()->parallel_try_catch([&]() { check_particle_type(p_type); });
279 }
283 }
284 if (name == "calculate_energy") {
285 return m_obs_stat->do_call_method("calculate_energy", {});
286 }
287 if (name == "calculate_scalar_pressure") {
288 return m_obs_stat->do_call_method("calculate_scalar_pressure", {});
289 }
290 if (name == "calculate_pressure_tensor") {
291 return m_obs_stat->do_call_method("calculate_pressure_tensor", {});
292 }
293#ifdef ESPRESSO_NPT
294 if (name == "get_instantaneous_pressure") {
295 return get_system().npt_inst_pressure->p_inst[0];
296 }
297 if (name == "get_instantaneous_pressure_virial") {
298 return get_system().npt_inst_pressure->p_inst[1];
299 }
300#endif // ESPRESSO_NPT
301 return {};
302}
303
304} // namespace Analysis
305} // 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.
Describes a cell structure / cell system.
Particle * get_local_particle(int id)
Get a local particle by id.
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
Context * context() const
Responsible context.
std::string_view name() const
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:184
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:133
auto make_unordered_map_of_variants(std::unordered_map< K, V > const &v)
Definition Variant.hpp:144
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:148
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.
Recursive variant implementation.
Definition Variant.hpp:84