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

Public Member Functions

void execute (internal::Assembly::Scratch::ScratchBase< dim > &scratch_base, internal::Assembly::CopyData::CopyDataBase< dim > &data_base) const override
 
std::vector< double > compute_residual (internal::Assembly::Scratch::ScratchBase< dim > &scratch_base) const override
 
virtual std::vector< double > advection_prefactors (internal::Assembly::Scratch::ScratchBase< dim > &scratch_base) const override
 
virtual std::vector< double > diffusion_prefactors (internal::Assembly::Scratch::ScratchBase< dim > &scratch_base) const override
 
- Public Member Functions inherited from aspect::Assemblers::Interface< dim >
virtual void create_additional_material_model_outputs (MaterialModel::MaterialModelOutputs< dim > &) const
 
- Public Member Functions inherited from aspect::Plugins::InterfaceBase
virtual ~InterfaceBase ()=default
 
virtual void initialize ()
 
virtual void update ()
 
virtual void parse_parameters (ParameterHandler &prm)
 
- Public Member Functions inherited from aspect::Assemblers::AdvectionStabilizationInterface< dim >
virtual ~AdvectionStabilizationInterface ()
 
- Public Member Functions inherited from aspect::SimulatorAccess< dim >
 SimulatorAccess ()
 
 SimulatorAccess (const Simulator< dim > &simulator_object)
 
virtual ~SimulatorAccess ()=default
 
virtual void initialize_simulator (const Simulator< dim > &simulator_object)
 
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
 
unsigned int n_particle_managers () const
 
const Particle::Manager< dim > & get_particle_manager (const unsigned int particle_manager_index) const
 
Particle::Manager< dim > & get_particle_manager (const unsigned int particle_manager_index)
 
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
 

Additional Inherited Members

- Static Public Member Functions inherited from aspect::Plugins::InterfaceBase
static void declare_parameters (ParameterHandler &prm)
 
- 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)
 

Detailed Description

template<int dim>
class aspect::Assemblers::DiffusionSystem< dim >

This class assembles the terms for the matrix and right-hand-side equation for the current cell in case we only want to solve the diffusion equation.

Definition at line 70 of file advection.h.

Member Function Documentation

§ execute()

template<int dim>
void aspect::Assemblers::DiffusionSystem< dim >::execute ( internal::Assembly::Scratch::ScratchBase< dim > &  scratch,
internal::Assembly::CopyData::CopyDataBase< dim > &  data 
) const
overridevirtual

Execute this assembler object. This function performs the primary work of an assembler. More precisely, it uses information for the current cell that is stored in scratch (like the material properties on this cell and the position of quadrature points) and computes the matrix and right hand side contributions for a set of terms for the given cell. These contributions are stored in data. Note, that the data in scratch and data is shared between all active assemblers so that each assembler should only add contributions to data, not overwrite entries in the matrix. After all assemblers have finished, the final content of data is distributed into the global matrix and right hand side vector.

Implements aspect::Assemblers::Interface< dim >.

§ compute_residual()

template<int dim>
std::vector<double> aspect::Assemblers::DiffusionSystem< dim >::compute_residual ( internal::Assembly::Scratch::ScratchBase< dim > &  ) const
overridevirtual

A required function for objects that implement the assembly of terms in an equation that requires the computation of residuals (in particular the advection equation in ASPECT). Just like the assemblers itself, the residual that we use to compute the necessary entropy viscosity depend on the equation (i.e. which terms are actually included in the equation). Thus different objects compute different residuals (i.e. the residual for a melt advection equation looks different from the residual for a passive compositional field). For assemblers for the Stokes system, an implementation of this function is not necessary.

Reimplemented from aspect::Assemblers::Interface< dim >.

§ advection_prefactors()

template<int dim>
virtual std::vector<double> aspect::Assemblers::DiffusionSystem< dim >::advection_prefactors ( internal::Assembly::Scratch::ScratchBase< dim > &  scratch_base) const
overridevirtual

This function returns a representative prefactor for the advection term of the equation for each quadrature point of the current cell. In the non-dimensional case this is simply 1.0, but for other quantities like temperature it is computed using physical units (like density and specific heat capacity). This information is useful for algorithms that depend on the magnitude of individual terms, like stabilization methods.

Reimplemented from aspect::Assemblers::AdvectionStabilizationInterface< dim >.

§ diffusion_prefactors()

template<int dim>
virtual std::vector<double> aspect::Assemblers::DiffusionSystem< dim >::diffusion_prefactors ( internal::Assembly::Scratch::ScratchBase< dim > &  scratch_base) const
overridevirtual

This function returns a representative conductivity for the diffusion part of the equation for each quadrature point of the current cell. For the pure advection case this factor is 0.0, but for other quantities like temperature it is computed using physical units (like thermal conductivity). This information is useful for algorithms depending on the magnitude of individual terms, like stabilization methods.

Reimplemented from aspect::Assemblers::AdvectionStabilizationInterface< dim >.


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