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>
94template <
unsigned int dim,
unsigned int degree,
typename number>
112 dealii::UpdateFlags::update_values |
113 dealii::UpdateFlags::update_JxW_values);
120 if (!
variable.is_nucleation_rate_variable)
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);
150 for (
unsigned int q_point = 0; q_point < num_quad_points; ++q_point)
160 for (
unsigned int d = 0;
d <
dim; ++
d)
173 unsigned int seed_increment =
time_info.get_increment();
187template <
unsigned int dim,
unsigned int degree,
typename number>
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;
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()));
251 <<
" New nucleus at: " <<
nuc.location <<
"\n";
258 <<
" new nuclei after exclusion.\n"
268template <
unsigned int dim,
unsigned int degree,
typename number>
322template <
unsigned int dim,
unsigned int degree,
typename number>
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
@ 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