6#include <deal.II/base/mpi.h>
7#include <deal.II/fe/fe_values.h>
26#include <prismspf/config.h>
42template <
unsigned int dim,
unsigned int degree,
typename number>
60 const double &delta_t,
64 std::poisson_distribution<unsigned int> distribution(rate * delta_t * volume);
65 return distribution(rng);
94template <
unsigned int dim,
unsigned int degree,
typename number>
112 dealii::UpdateFlags::update_values |
113 dealii::UpdateFlags::update_JxW_values);
114 std::list<Nucleus<dim>> new_nuclei_list;
120 if (!variable.is_nucleation_rate_variable)
124 std::uniform_int_distribution<unsigned int> nucleating_index_dist(
126 variable.nucleating_field_indices.size() - 1);
131 .active_cell_iterators())
133 if (!cell->is_locally_owned())
137 std::vector<number> values(num_quad_points, 0.0);
139 const auto dof_iterator = cell->as_dof_handler_iterator(
143 fe_values.reinit(dof_iterator);
145 fe_values.get_function_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)
152 nuc_rate += values[q_point] * fe_values.get_quadrature().weight(q_point);
153 cell_volume += fe_values.JxW(q_point);
155 unsigned int num_nuclei_in_cell =
157 for (
unsigned int i = 0; i < num_nuclei_in_cell; ++i)
159 dealii::Point<dim> nucleus_location_unit_cell;
160 for (
unsigned int d = 0; d < dim; ++d)
164 static std::uniform_real_distribution<double> uniform_unit_interval(
167 nucleus_location_unit_cell[d] = uniform_unit_interval(rng);
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();
174 unsigned int nucleating_index =
175 variable.nucleating_field_indices[nucleating_index_dist(rng)];
177 new_nuclei_list.emplace_back(nucleating_index,
187template <
unsigned int dim,
unsigned int degree,
typename number>
196 if (!
bool(dealii::Utilities::MPI::sum(new_nuclei_list.size(), MPI_COMM_WORLD)))
206 std::vector<Nucleus<dim>> new_nuclei(new_nuclei_list.begin(), new_nuclei_list.end());
208 bool any_nucleation_occurred =
false;
209 if (dealii::Utilities::MPI::this_mpi_process(MPI_COMM_WORLD) == 0)
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;
219 std::shuffle(new_nuclei.begin(), new_nuclei.end(), rng);
221 while (!new_nuclei.empty())
225 std::none_of(global_nuclei.begin(),
229 const double distance = prisms::distance<dim, double>(
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 &&
237 nuc_params.get_same_field_exclusion_distance()));
249 global_nuclei.push_back(nuc);
251 <<
" New nucleus at: " << nuc.
location <<
"\n";
253 any_nucleation_occurred =
true;
255 new_nuclei.pop_back();
258 <<
" new nuclei after exclusion.\n"
260 << global_nuclei.size() <<
" total nuclei.\n\n"
263 MPI_Bcast(&any_nucleation_occurred, 1, MPI_CXX_BOOL, 0, MPI_COMM_WORLD);
265 return any_nucleation_occurred;
268template <
unsigned int dim,
unsigned int degree,
typename number>
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;
280 nuclei_counts_per_rank.resize(num_procs);
283 MPI_Gather(&local_count,
286 nuclei_counts_per_rank.data(),
293 std::vector<int> recv_displacements;
294 std::vector<Nucleus<dim>> gathered_nuclei;
297 recv_displacements.resize(num_procs);
298 recv_displacements[0] = 0;
299 for (
int r = 1; r < num_procs; ++r)
300 recv_displacements[r] = recv_displacements[r - 1] + nuclei_counts_per_rank[r - 1];
302 int total_count = recv_displacements.back() + nuclei_counts_per_rank.back();
303 gathered_nuclei.resize(total_count);
307 MPI_Gatherv(local_nuclei.data(),
310 gathered_nuclei.data(),
311 nuclei_counts_per_rank.data(),
312 recv_displacements.data(),
318 local_nuclei = gathered_nuclei;
322template <
unsigned int dim,
unsigned int degree,
typename number>
327 int rank = dealii::Utilities::MPI::this_mpi_process(MPI_COMM_WORLD);
328 int count = local_nuclei.size();
329 MPI_Bcast(&count, 1, MPI_INT, 0, MPI_COMM_WORLD);
333 local_nuclei.resize(count);
339PRISMS_PF_END_NAMESPACE
static dealii::ConditionalOStream & pout_base()
Generic parallel output stream. Used for essential information in release and debug mode.
Definition conditional_ostreams.cc:44
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
unsigned int get_increment() const
Definition simulation_timer.h:22
double get_timestep() const
Definition simulation_timer.h:34
double get_time() const
Definition simulation_timer.h:28
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
static const std::array< const dealii::FESystem< dim >, 2 > fe_systems
Scalar and Vector FE systems.
Definition system_wide.h:29
static const dealii::QGaussLobatto< dim > quadrature
Quadrature rule.
Definition system_wide.h:41
static const dealii::MappingQ1< dim > mapping
Mappings to and from reference cell.
Definition system_wide.h:36
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
static MPI_Datatype mpi_datatype()
Definition nucleus.h:77
dealii::Point< dim > location
Definition nucleus.h:44