![]() |
ASPECT
|
Public Member Functions | |
Interface () | |
~Interface () override=default | |
virtual void | initialize () |
virtual DEAL_II_DEPRECATED void | generate_particles (std::multimap< Particles::internal::LevelInd, Particle< dim >> &particles) |
virtual void | generate_particles (Particles::ParticleHandler< dim > &particle_handler) |
std::pair< Particles::internal::LevelInd, Particle< dim > > | generate_particle (const typename parallel::distributed::Triangulation< dim >::active_cell_iterator &cell, const types::particle_index id) |
virtual void | parse_parameters (ParameterHandler &prm) |
![]() | |
SimulatorAccess () | |
SimulatorAccess (const Simulator< dim > &simulator_object) | |
virtual | ~SimulatorAccess ()=default |
virtual void | initialize_simulator (const Simulator< dim > &simulator_object) |
template<typename PostprocessorType > | |
PostprocessorType * | find_postprocessor () const |
const Introspection< dim > & | introspection () const |
const Simulator< dim > & | get_simulator () const |
const Parameters< dim > & | get_parameters () const |
SimulatorSignals< dim > & | get_signals () const |
MPI_Comm | get_mpi_communicator () const |
TimerOutput & | get_computing_timer () const |
const ConditionalOStream & | get_pcout () const |
double | get_time () const |
double | get_timestep () const |
double | get_old_timestep () const |
unsigned int | get_timestep_number () const |
unsigned int | get_nonlinear_iteration () const |
const parallel::distributed::Triangulation< dim > & | get_triangulation () const |
double | get_volume () const |
const Mapping< dim > & | get_mapping () const |
std::string | get_output_directory () const |
bool | include_adiabatic_heating () const |
bool | include_latent_heat () const |
bool | include_melt_transport () const |
int | get_stokes_velocity_degree () const |
double | get_adiabatic_surface_temperature () const |
double | get_surface_pressure () const |
bool | convert_output_to_years () const |
unsigned int | get_pre_refinement_step () const |
unsigned int | n_compositional_fields () const |
void | get_refinement_criteria (Vector< float > &estimated_error_per_cell) const |
void | get_artificial_viscosity (Vector< float > &viscosity_per_cell, const bool skip_interior_cells=false) const |
void | get_artificial_viscosity_composition (Vector< float > &viscosity_per_cell, const unsigned int compositional_variable) const |
const LinearAlgebra::BlockVector & | get_current_linearization_point () const |
const LinearAlgebra::BlockVector & | get_solution () const |
const LinearAlgebra::BlockVector & | get_old_solution () const |
const LinearAlgebra::BlockVector & | get_old_old_solution () const |
const LinearAlgebra::BlockVector & | get_reaction_vector () const |
const LinearAlgebra::BlockVector & | get_mesh_velocity () const |
const DoFHandler< dim > & | get_dof_handler () const |
const FiniteElement< dim > & | get_fe () const |
const LinearAlgebra::BlockSparseMatrix & | get_system_matrix () const |
const LinearAlgebra::BlockSparseMatrix & | get_system_preconditioner_matrix () const |
const MaterialModel::Interface< dim > & | get_material_model () const |
const GravityModel::Interface< dim > & | get_gravity_model () const |
const InitialTopographyModel::Interface< dim > & | get_initial_topography_model () const |
const GeometryModel::Interface< dim > & | get_geometry_model () const |
const AdiabaticConditions::Interface< dim > & | get_adiabatic_conditions () const |
bool | has_boundary_temperature () const |
DEAL_II_DEPRECATED const BoundaryTemperature::Interface< dim > & | get_boundary_temperature () const |
const BoundaryTemperature::Manager< dim > & | get_boundary_temperature_manager () const |
const BoundaryHeatFlux::Interface< dim > & | get_boundary_heat_flux () const |
bool | has_boundary_composition () const |
DEAL_II_DEPRECATED const BoundaryComposition::Interface< dim > & | get_boundary_composition () const |
const BoundaryComposition::Manager< dim > & | get_boundary_composition_manager () const |
const std::map< types::boundary_id, std::unique_ptr< BoundaryTraction::Interface< dim > > > & | get_boundary_traction () const |
DEAL_II_DEPRECATED const InitialTemperature::Interface< dim > & | get_initial_temperature () const |
std::shared_ptr< const InitialTemperature::Manager< dim > > | get_initial_temperature_manager_pointer () const |
const InitialTemperature::Manager< dim > & | get_initial_temperature_manager () const |
DEAL_II_DEPRECATED const InitialComposition::Interface< dim > & | get_initial_composition () const |
std::shared_ptr< const InitialComposition::Manager< dim > > | get_initial_composition_manager_pointer () const |
const InitialComposition::Manager< dim > & | get_initial_composition_manager () const |
const std::set< types::boundary_id > & | get_fixed_temperature_boundary_indicators () const |
const std::set< types::boundary_id > & | get_fixed_heat_flux_boundary_indicators () const |
const std::set< types::boundary_id > & | get_fixed_composition_boundary_indicators () const |
const std::set< types::boundary_id > & | get_mesh_deformation_boundary_indicators () const |
const BoundaryVelocity::Manager< dim > & | get_boundary_velocity_manager () const |
const HeatingModel::Manager< dim > & | get_heating_model_manager () const |
const MeshRefinement::Manager< dim > & | get_mesh_refinement_manager () const |
const MeltHandler< dim > & | get_melt_handler () const |
const VolumeOfFluidHandler< dim > & | get_volume_of_fluid_handler () const |
const NewtonHandler< dim > & | get_newton_handler () const |
const MeshDeformation::MeshDeformationHandler< dim > & | get_mesh_deformation_handler () const |
const LateralAveraging< dim > & | get_lateral_averaging () const |
const AffineConstraints< double > & | get_current_constraints () const |
bool | simulator_is_past_initialization () const |
double | get_pressure_scaling () const |
bool | pressure_rhs_needs_compatibility_modification () const |
bool | model_has_prescribed_stokes_solution () const |
TableHandler & | get_statistics_object () const |
template<typename PostprocessorType > | |
DEAL_II_DEPRECATED PostprocessorType * | find_postprocessor () const |
const Postprocess::Manager< dim > & | get_postprocess_manager () const |
const Particle::World< dim > & | get_particle_world () const |
Particle::World< dim > & | get_particle_world () |
bool | is_stokes_matrix_free () |
const StokesMatrixFreeHandler< dim > & | get_stokes_matrix_free () const |
RotationProperties< dim > | compute_net_angular_momentum (const bool use_constant_density, const LinearAlgebra::BlockVector &solution, const bool limit_to_top_faces=false) const |
Static Public Member Functions | |
static void | declare_parameters (ParameterHandler &prm) |
![]() | |
static void | get_composition_values_at_q_point (const std::vector< std::vector< double >> &composition_values, const unsigned int q, std::vector< double > &composition_values_at_q_point) |
Protected Member Functions | |
DEAL_II_DEPRECATED std::pair< Particles::internal::LevelInd, Particle< dim > > | generate_particle (const Point< dim > &position, const types::particle_index id) const |
Particles::ParticleIterator< dim > | insert_particle_at_position (const Point< dim > &position, const types::particle_index id, Particles::ParticleHandler< dim > &particle_handler) const |
Protected Attributes | |
std::mt19937 | random_number_generator |
Abstract base class used for classes that generate particles.
aspect::Particle::Generator::Interface< dim >::Interface | ( | ) |
Constructor. Initializes the random number generator.
|
overridedefault |
Destructor. Made virtual so that derived classes can be created and destroyed through pointers to the base class.
|
virtual |
Initialization function. This function is called once at the beginning of the program after parse_parameters is run and after the SimulatorAccess (if applicable) is initialized.
|
virtual |
Generate particles. Every derived class has to decide on the method and number of particles to generate, for example using input parameters declared in their declare_parameters and parse_parameters functions. This function should generate the particles and associate them to their according cells by inserting them into a multimap between cell and particle. This map becomes very large if the particle count per process is large, so we hand it over by reference instead of returning the multimap.
[in,out] | particles | A multimap between cells and their particles. This map will be filled in this function. |
|
virtual |
Generate particles. Every derived class has to decide on the method and number of particles to generate, for example using input parameters declared in their declare_parameters and parse_parameters functions. This function should generate the particles and insert them into particle_handler
by calling its respective functions.
[in,out] | particle_handler | The particle handler into which the generated particles should be inserted. |
Reimplemented in aspect::Particle::Generator::ProbabilityDensityFunction< dim >, aspect::Particle::Generator::UniformRadial< dim >, aspect::Particle::Generator::ReferenceCell< dim >, aspect::Particle::Generator::UniformBox< dim >, aspect::Particle::Generator::QuadraturePoints< dim >, aspect::Particle::Generator::AsciiFile< dim >, and aspect::Particle::Generator::RandomUniform< dim >.
std::pair<Particles::internal::LevelInd,Particle<dim> > aspect::Particle::Generator::Interface< dim >::generate_particle | ( | const typename parallel::distributed::Triangulation< dim >::active_cell_iterator & | cell, |
const types::particle_index | id | ||
) |
Generate one particle in the given cell. This function's main purpose is to provide functionality to fill up cells with too few particles after refinement. Of course it can also be utilized by derived classes to generate the initial particle distribution.
|
static |
Declare the parameters this class takes through input files. The default implementation of this function does not describe any parameters. Consequently, derived classes do not have to overload this function if they do not take any runtime parameters.
|
virtual |
Read the parameters this class declares from the parameter file. The default implementation of this function does not read any parameters. Consequently, derived classes do not have to overload this function if they do not take any runtime parameters.
Reimplemented in aspect::Particle::Generator::ProbabilityDensityFunction< dim >, aspect::Particle::Generator::UniformRadial< dim >, aspect::Particle::Generator::ReferenceCell< dim >, aspect::Particle::Generator::UniformBox< dim >, aspect::Particle::Generator::AsciiFile< dim >, and aspect::Particle::Generator::RandomUniform< dim >.
|
protected |
Generate a particle at the specified position and with the specified id. Many derived classes use this functionality, therefore it is implemented here to avoid duplication. In case the position is not in the local domain this function throws an exception of type ExcParticlePointNotInDomain, which can be caught in the calling plugin.
|
protected |
Generate a particle at the specified position and with the specified id and insert it into the particle_handler
. Many derived classes use this functionality, therefore it is implemented here to avoid duplication. In case the position is not in the local domain this function throws an exception of type ExcParticlePointNotInDomain, which can be caught in the calling plugin. Note that since the cell in which the particle is generated is not known, it has to be found, which is an expensive operation.
position | Position of the particle. |
id | The id of the particle. |
particle_handler | The particle handler into which the particle should be inserted. |
|
protected |
Random number generator. For reproducibility of tests it is initialized in the constructor with a constant.
Definition at line 202 of file interface.h.