ASPECT
Public Member Functions | Static Public Member Functions | Public Attributes | Private Attributes | List of all members
aspect::MeltHandler< dim > Class Template Reference
Inheritance diagram for aspect::MeltHandler< dim >:
Inheritance graph
[legend]

Public Member Functions

 MeltHandler (ParameterHandler &prm)
 
void edit_finite_element_variables (const Parameters< dim > &parameters, std::vector< VariableDeclaration< dim >> &variables)
 
void set_assemblers (Assemblers::Manager< dim > &assemblers) const
 
void initialize () const
 
void initialize_simulator (const Simulator< dim > &simulator_object) override
 
void compute_melt_variables (LinearAlgebra::BlockSparseMatrix &system_matrix, LinearAlgebra::BlockVector &solution, LinearAlgebra::BlockVector &system_rhs)
 
bool is_porosity (const typename Simulator< dim >::AdvectionField &advection_field) const
 
void apply_free_surface_stabilization_with_melt (const double free_surface_theta, const typename DoFHandler< dim >::active_cell_iterator &cell, internal::Assembly::Scratch::StokesSystem< dim > &scratch, internal::Assembly::CopyData::StokesSystem< dim > &data) const
 
void add_current_constraints (AffineConstraints< double > &constraints)
 
bool is_melt_cell (const typename DoFHandler< dim >::active_cell_iterator &cell) const
 
double limited_darcy_coefficient (const double K_D, const bool is_melt_cell) const
 
const BoundaryFluidPressure::Interface< dim > & get_boundary_fluid_pressure () const
 
- Public Member Functions inherited from aspect::SimulatorAccess< dim >
 SimulatorAccess ()
 
 SimulatorAccess (const Simulator< dim > &simulator_object)
 
virtual ~SimulatorAccess ()=default
 
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
 
const TimeStepping::Manager< dim > & get_timestepping_manager () 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
 
double get_end_time () 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::BlockVectorget_current_linearization_point () const
 
const LinearAlgebra::BlockVectorget_solution () const
 
const LinearAlgebra::BlockVectorget_old_solution () const
 
const LinearAlgebra::BlockVectorget_old_old_solution () const
 
const LinearAlgebra::BlockVectorget_reaction_vector () const
 
const LinearAlgebra::BlockVectorget_mesh_velocity () const
 
const DoFHandler< dim > & get_dof_handler () const
 
const FiniteElement< dim > & get_fe () const
 
const LinearAlgebra::BlockSparseMatrixget_system_matrix () const
 
const LinearAlgebra::BlockSparseMatrixget_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
 
const BoundaryTemperature::Manager< dim > & get_boundary_temperature_manager () const
 
const BoundaryHeatFlux::Interface< dim > & get_boundary_heat_flux () const
 
bool has_boundary_composition () const
 
const BoundaryComposition::Manager< dim > & get_boundary_composition_manager () const
 
const BoundaryTraction::Manager< dim > & get_boundary_traction_manager () const
 
std::shared_ptr< const InitialTemperature::Manager< dim > > get_initial_temperature_manager_pointer () const
 
const InitialTemperature::Manager< dim > & get_initial_temperature_manager () 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
 
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 create_material_model_outputs (MaterialModel::MaterialModelOutputs< dim > &output)
 
- Static Public Member Functions inherited from aspect::SimulatorAccess< dim >
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)
 

Public Attributes

Melt::Parameters< dim > melt_parameters
 

Private Attributes

const std::unique_ptr< aspect::BoundaryFluidPressure::Interface< dim > > boundary_fluid_pressure
 
std::vector< bool > is_melt_cell_vector
 
AffineConstraints< double > current_constraints
 

Detailed Description

template<int dim>
class aspect::MeltHandler< dim >

Class that contains all runtime parameters and other helper functions related to melt transport. A global instance can be retrieved with SimulatorAccess<dim>::get_melt_handler(), but keep in mind that it only exists if parameters.include_melt_transport is true.

Definition at line 431 of file melt.h.

Constructor & Destructor Documentation

§ MeltHandler()

template<int dim>
aspect::MeltHandler< dim >::MeltHandler ( ParameterHandler &  prm)

Member Function Documentation

§ create_material_model_outputs()

template<int dim>
static void aspect::MeltHandler< dim >::create_material_model_outputs ( MaterialModel::MaterialModelOutputs< dim > &  output)
static

Create an additional material model output object that contains the additional output variables needed in simulation with melt transport, and attaches a pointer to it to the corresponding vector in the MaterialModel::MaterialModelOutputs structure.

§ edit_finite_element_variables()

template<int dim>
void aspect::MeltHandler< dim >::edit_finite_element_variables ( const Parameters< dim > &  parameters,
std::vector< VariableDeclaration< dim >> &  variables 
)

Add the additional variables we need in simulations with melt migration to the list of variables, which will be used later to set up the introspection object.

§ set_assemblers()

template<int dim>
void aspect::MeltHandler< dim >::set_assemblers ( Assemblers::Manager< dim > &  assemblers) const

Determine, based on the run-time parameters of the current simulation, which functions need to be called in order to assemble linear systems, matrices, and right hand side vectors.

§ initialize()

template<int dim>
void aspect::MeltHandler< dim >::initialize ( ) const

Initialize function. This is mainly to check that the melt transport parameters chosen in the input file are consistent with the rest of the options. We can not do this in the parse_parameters function, as we do not have simulator access at that point.

§ initialize_simulator()

template<int dim>
void aspect::MeltHandler< dim >::initialize_simulator ( const Simulator< dim > &  simulator_object)
overridevirtual

Setup SimulatorAccess for the plugins related to melt transport.

Reimplemented from aspect::SimulatorAccess< dim >.

§ compute_melt_variables()

template<int dim>
void aspect::MeltHandler< dim >::compute_melt_variables ( LinearAlgebra::BlockSparseMatrix system_matrix,
LinearAlgebra::BlockVector solution,
LinearAlgebra::BlockVector system_rhs 
)

Compute fluid velocity and solid pressure in this ghosted solution vector. The fluid velocity is computed by solving a mass matrix problem, and the solid pressure is computed algebraically.

Parameters
system_matrixThe system matrix with an already set up sparsity pattern that will be used by this function to compute the melt variables.
solutionThe existing solution vector that contains the values for porosity, compaction pressure, fluid pressure and solid velocity obtained by solving the Stokes and advection system, and that will be updated with the computed values for fluid velocity and solid pressure.
system_rhsThe right-hand side vector that will be used by this function to compute the melt variables.

§ is_porosity()

template<int dim>
bool aspect::MeltHandler< dim >::is_porosity ( const typename Simulator< dim >::AdvectionField &  advection_field) const

Return whether this object refers to the porosity field.

§ apply_free_surface_stabilization_with_melt()

template<int dim>
void aspect::MeltHandler< dim >::apply_free_surface_stabilization_with_melt ( const double  free_surface_theta,
const typename DoFHandler< dim >::active_cell_iterator &  cell,
internal::Assembly::Scratch::StokesSystem< dim > &  scratch,
internal::Assembly::CopyData::StokesSystem< dim > &  data 
) const

Apply free surface stabilization to a cell of the system matrix when melt transport is used in the computation. Called during assembly of the system matrix.

§ add_current_constraints()

template<int dim>
void aspect::MeltHandler< dim >::add_current_constraints ( AffineConstraints< double > &  constraints)

Constrain the compaction pressure to zero in all cells that are not "melt cells" (cells where the porosity is above a given threshold). This reverts the system of equations we solve back to the Stokes system without melt transport for these cells.

§ is_melt_cell()

template<int dim>
bool aspect::MeltHandler< dim >::is_melt_cell ( const typename DoFHandler< dim >::active_cell_iterator &  cell) const

Returns the entry of the private variable is_melt_cell_vector for the cell given in the input, describing if we have melt transport in this cell or not.

§ limited_darcy_coefficient()

template<int dim>
double aspect::MeltHandler< dim >::limited_darcy_coefficient ( const double  K_D,
const bool  is_melt_cell 
) const

Given the Darcy coefficient K_D as computed by the material model, limit the coefficient to a minimum value (computed as the K_D variation threshold given in the input file times the reference Darcy coefficient) in melt cells and return this value. If is_melt_cell is false, return zero.

§ get_boundary_fluid_pressure()

template<int dim>
const BoundaryFluidPressure::Interface<dim>& aspect::MeltHandler< dim >::get_boundary_fluid_pressure ( ) const

Return a pointer to the boundary fluid pressure.

Member Data Documentation

§ melt_parameters

template<int dim>
Melt::Parameters<dim> aspect::MeltHandler< dim >::melt_parameters

The object that stores the run-time parameters that control the how the melt transport equations are solved.

Definition at line 540 of file melt.h.

§ boundary_fluid_pressure

template<int dim>
const std::unique_ptr<aspect::BoundaryFluidPressure::Interface<dim> > aspect::MeltHandler< dim >::boundary_fluid_pressure
private

Store a pointer to the fluid pressure boundary plugin, so that the initialization can be done together with the other objects related to melt transport.

Definition at line 548 of file melt.h.

§ is_melt_cell_vector

template<int dim>
std::vector<bool> aspect::MeltHandler< dim >::is_melt_cell_vector
private

is_melt_cell_vector[cell->active_cell_index()] says whether we want to solve the melt transport equations (as opposed to the Stokes equations without melt) in this cell or not. The value is set to true or false based on the porosity in the cell (only cells where the porosity is above a threshold are considered melt cells).

Definition at line 557 of file melt.h.

§ current_constraints

template<int dim>
AffineConstraints<double> aspect::MeltHandler< dim >::current_constraints
private

Constraint object. We need to save the current constraints at the start of every time step so that we can add the melt constraints, which depend on the solution of the porosity field, later after we have computed this solution.

Definition at line 565 of file melt.h.


The documentation for this class was generated from the following file: