PRISMS-PF Manual
Loading...
Searching...
No Matches
nucleation_manager.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: © 2025 PRISMS Center at the University of Michigan
2// SPDX-License-Identifier: GNU Lesser General Public Version 2.1
3
4#pragma once
5
6#include <deal.II/base/mpi.h>
7#include <deal.II/fe/fe_values.h>
8
13#include <prismspf/core/types.h>
14
16
18
23
25
26#include <prismspf/config.h>
27
28#include <algorithm>
29#include <list>
30#include <mpi.h>
31#include <random>
32#include <vector>
33
35
36// Note: could make NucleationManager a namespace as everything is static
37
42template <unsigned int dim, unsigned int degree, typename number>
44{
45public:
50 static bool
52 std::vector<Nucleus<dim>> &nuclei);
53
58 static unsigned int
60 const double &delta_t,
61 const double &volume,
62 RNGEngine &rng)
63 {
64 std::poisson_distribution<unsigned int> distribution(rate * delta_t * volume);
65 return distribution(rng);
66 }
67
73 static bool
75 std::vector<Nucleus<dim>> &global_nuclei,
76 const UserInputParameters<dim> &user_inputs,
78
83 static void
85
90 static void
92};
93
94template <unsigned int dim, unsigned int degree, typename number>
95inline bool
97 const SolveContext<dim, degree, number> &solve_context,
98 std::vector<Nucleus<dim>> &nuclei)
99{
100 // Set up references.
101 const UserInputParameters<dim> &user_inputs = solve_context.get_user_inputs();
102 const NucleationParameters &nuc_params = user_inputs.nucleation_parameters;
103 const SimulationTimer &time_info = solve_context.get_simulation_timer();
104 const double delta_t = nuc_params.get_nucleation_period() * time_info.get_timestep();
105 auto &rng = user_inputs.misc_parameters.rng;
106
107 // Set up FEValues
108 unsigned int num_quad_points = SystemWide<dim, degree>::quadrature.size();
109 // Made static because initialization was taking a LOT of time
110 static dealii::FEValues<dim> fe_values(SystemWide<dim, degree>::fe_systems[0],
112 dealii::UpdateFlags::update_values |
113 dealii::UpdateFlags::update_JxW_values);
114 std::list<Nucleus<dim>> new_nuclei_list;
115 // Loop over nucleation rate variables and attempt seeding at each cell
116 for (unsigned int index = 0; index < solve_context.get_field_attributes().size();
117 ++index)
118 {
119 const auto &variable = solve_context.get_field_attributes()[index];
120 if (!variable.is_nucleation_rate_variable)
121 {
122 continue;
123 }
124 std::uniform_int_distribution<unsigned int> nucleating_index_dist(
125 0,
126 variable.nucleating_field_indices.size() - 1);
127 // Perform nucleation logic here
128 // This is where you would check conditions and create nuclei
129 for (const auto &cell : solve_context.get_triangulation_manager()
130 .get_triangulation()
131 .active_cell_iterators())
132 {
133 if (!cell->is_locally_owned())
134 {
135 continue;
136 }
137 std::vector<number> values(num_quad_points, 0.0);
138 // Grab the DoFHandler iterator
139 const auto dof_iterator = cell->as_dof_handler_iterator(
140 solve_context.get_dof_manager().get_field_dof_handler(index));
141
142 // Reinit the cell
143 fe_values.reinit(dof_iterator);
144 // Get the values for a scalar field
145 fe_values.get_function_values(
146 (solve_context.get_solution_indexer().get_solution_vector(index)),
147 values);
148 double nuc_rate = 0.0;
149 double cell_volume = 0.0;
150 for (unsigned int q_point = 0; q_point < num_quad_points; ++q_point)
151 {
152 nuc_rate += values[q_point] * fe_values.get_quadrature().weight(q_point);
153 cell_volume += fe_values.JxW(q_point);
154 }
155 unsigned int num_nuclei_in_cell =
156 calculate_number_of_events(nuc_rate, delta_t, cell_volume, rng);
157 for (unsigned int i = 0; i < num_nuclei_in_cell; ++i)
158 {
159 dealii::Point<dim> nucleus_location_unit_cell;
160 for (unsigned int d = 0; d < dim; ++d)
161 {
162 // Note: if we ever do non-rectangular cells, try a randomly weighted
163 // sum over dealii::GeometryInfo< dim >::unit_cell_vertex
164 static std::uniform_real_distribution<double> uniform_unit_interval(
165 0.0,
166 1.0);
168 }
169 dealii::Point<dim> nucleus_location =
171 .transform_unit_to_real_cell(cell, nucleus_location_unit_cell);
172 double seed_time = time_info.get_time();
173 unsigned int seed_increment = time_info.get_increment();
174 unsigned int nucleating_index =
175 variable.nucleating_field_indices[nucleating_index_dist(rng)];
176
179 seed_time,
180 seed_increment);
181 }
182 }
183 }
184 return gather_exclude_broadcast_nuclei(new_nuclei_list, nuclei, user_inputs, time_info);
185}
186
187template <unsigned int dim, unsigned int degree, typename number>
188inline bool
190 std::list<Nucleus<dim>> &new_nuclei_list,
191 std::vector<Nucleus<dim>> &global_nuclei,
192 const UserInputParameters<dim> &user_inputs,
194{
195 // dont waste time if no nuclei appeared
196 if (!bool(dealii::Utilities::MPI::sum(new_nuclei_list.size(), MPI_COMM_WORLD)))
197 {
198 return false;
199 }
200
201 // Set up refs
202 const NucleationParameters &nuc_params = user_inputs.nucleation_parameters;
203 RNGEngine &rng = user_inputs.misc_parameters.rng;
204
205 // Gather new nuclei to root process
206 std::vector<Nucleus<dim>> new_nuclei(new_nuclei_list.begin(), new_nuclei_list.end());
207 mpi_gather_nuclei(new_nuclei);
208 bool any_nucleation_occurred = false;
209 if (dealii::Utilities::MPI::this_mpi_process(MPI_COMM_WORLD) == 0)
210 {
211 // Remove nuclei within their exclusion distance and add to nuclei list
213 << "[Increment " << time_info.get_increment() << "] : Nucleation\n"
214 << " " << new_nuclei.size() << " nuclei generated before exclusion.\n"
215 << " Excluding nuclei...\n";
216 unsigned int count = 0;
217
218 // remove bias from cell order
219 std::shuffle(new_nuclei.begin(), new_nuclei.end(), rng);
220
221 while (!new_nuclei.empty())
222 {
223 const Nucleus<dim> &nuc = new_nuclei.back();
224 bool valid =
225 std::none_of(global_nuclei.begin(),
226 global_nuclei.end(),
228 {
229 const double distance = prisms::distance<dim, double>(
230 nuc.location,
231 existing_nucleus.location,
232 user_inputs.spatial_discretization.rectangular_mesh);
233 return nuc_params.check_active(existing_nucleus, time_info) &&
234 (distance < nuc_params.get_exclusion_distance() ||
235 (nuc.field_index == existing_nucleus.field_index &&
236 distance <
237 nuc_params.get_same_field_exclusion_distance()));
238 });
239 if (valid)
240 {
241 // Note: Using push_back() in a loop is not good use for
242 // vectors. We also don't want to use reserve() on the upper bound
243 // because that could allocate much more space than needed. I originally
244 // was using a std::list to avoid this issue, but that is unfriendly to
245 // the MPI functions. One solution could be to convert between data
246 // structures as needed, but that also adds overhead. For now, I will
247 // assume that the total number of nuclei is not enough to
248 // cause significant performance issues.
249 global_nuclei.push_back(nuc);
251 << " New nucleus at: " << nuc.location << "\n";
252 ++count;
254 }
255 new_nuclei.pop_back();
256 }
258 << " new nuclei after exclusion.\n"
259 " "
260 << global_nuclei.size() << " total nuclei.\n\n"
261 << std::flush;
262 }
264 mpi_broadcast_nuclei(global_nuclei);
266}
267
268template <unsigned int dim, unsigned int degree, typename number>
269inline void
271 std::vector<Nucleus<dim>> &local_nuclei)
272{
273 // Step 1: Share how many nuclei each rank has
274 int rank = dealii::Utilities::MPI::this_mpi_process(MPI_COMM_WORLD);
275 int num_procs = dealii::Utilities::MPI::n_mpi_processes(MPI_COMM_WORLD);
276 int local_count = local_nuclei.size();
277 std::vector<int> nuclei_counts_per_rank;
278 if (rank == 0)
279 {
281 }
282
284 1,
285 MPI_INT,
287 1,
288 MPI_INT,
289 0,
291
292 // Step 2: Compute displacements and allocate receive buffer on root
293 std::vector<int> recv_displacements;
294 std::vector<Nucleus<dim>> gathered_nuclei;
295 if (rank == 0)
296 {
298 recv_displacements[0] = 0;
299 for (int r = 1; r < num_procs; ++r)
301
304 }
305
306 // Step 3: Gather all nuclei into root's buffer
310 gathered_nuclei.data(),
312 recv_displacements.data(),
314 0,
316 if (rank == 0)
317 {
319 }
320}
321
322template <unsigned int dim, unsigned int degree, typename number>
323inline void
325 std::vector<Nucleus<dim>> &local_nuclei)
326{
327 int rank = dealii::Utilities::MPI::this_mpi_process(MPI_COMM_WORLD);
328 int count = local_nuclei.size();
330
331 if (rank != 0)
332 {
333 local_nuclei.resize(count);
334 }
335
337}
338
static dealii::ConditionalOStream & pout_base()
Generic parallel output stream. Used for essential information in release and debug mode.
Definition conditional_ostreams.cc:43
The class handles the stochastic nucleation in PRISMS-PF.
Definition nucleation_manager.h:44
static bool attempt_nucleation(const SolveContext< dim, degree, number > &solve_context, std::vector< Nucleus< dim > > &nuclei)
Main nucleation function. Iterates over the domain and stochastically adds nuclei to the list.
Definition nucleation_manager.h:96
static void mpi_broadcast_nuclei(std::vector< Nucleus< dim > > &local_nuclei)
Broadcasts nuclei lists from root. Modifies.
Definition nucleation_manager.h:324
static unsigned int calculate_number_of_events(const double &rate, const double &delta_t, const double &volume, RNGEngine &rng)
Samples the poisson distribution to calculate a number of events in a time-volume given a rate.
Definition nucleation_manager.h:59
static bool gather_exclude_broadcast_nuclei(std::list< Nucleus< dim > > &new_nuclei_list, std::vector< Nucleus< dim > > &global_nuclei, const UserInputParameters< dim > &user_inputs, const SimulationTimer &time_info)
Gathers the potential new nuclei from each processor onto the root process, eliminates any nuclei tha...
Definition nucleation_manager.h:189
static void mpi_gather_nuclei(std::vector< Nucleus< dim > > &local_nuclei)
Gathers nuclei lists to root. Modifies.
Definition nucleation_manager.h:270
Definition simulation_timer.h:13
This class provides context for a solver with ptrs to all the relevant dependencies.
Definition solve_context.h:36
const DoFManager< dim, degree > & get_dof_manager() const
Get the dof manager.
Definition solve_context.h:101
const std::vector< FieldAttributes > & get_field_attributes() const
Get the field attributes.
Definition solve_context.h:62
SolutionIndexer< dim, number > & get_solution_indexer() const
Get the solution manager.
Definition solve_context.h:141
const UserInputParameters< dim > & get_user_inputs() const
Get the user-inputs.
Definition solve_context.h:71
const SimulationTimer & get_simulation_timer() const
Get the simulation timer.
Definition solve_context.h:169
const TriangulationManager< dim > & get_triangulation_manager() const
Get the triangulation manager.
Definition solve_context.h:81
This is the main class that handles the construction and solving of user-specified PDEs.
Definition system_wide.h:24
Definition user_input_parameters.h:32
MiscellaneousParameters misc_parameters
Definition user_input_parameters.h:199
NucleationParameters nucleation_parameters
Definition user_input_parameters.h:196
@ Value
Use value of the variable as a criterion for refinement.
Definition grid_refiner_criterion.h:31
std::mt19937 RNGEngine
Definition miscellaneous_parameters.h:19
Definition conditional_ostreams.cc:20
RNGEngine rng
Definition miscellaneous_parameters.h:51
Struct that holds nucleation parameters.
Definition nucleation_parameters.h:27
unsigned int get_nucleation_period() const
Get the nucleation period.
Definition nucleation_parameters.h:60
This class contains mutable utilities for phase field problems.
Definition nucleus.h:23