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

Public Member Functions

 RK4 ()
 
void initialize () override
 
void local_integrate_step (const typename ParticleHandler< dim >::particle_iterator &begin_particle, const typename ParticleHandler< dim >::particle_iterator &end_particle, const std::vector< Tensor< 1, dim >> &old_velocities, const std::vector< Tensor< 1, dim >> &velocities, const double dt) override
 
bool new_integration_step () override
 
std::array< bool, 3 > required_solution_vectors () const override
 
- Public Member Functions inherited from aspect::Particle::Integrator::Interface< dim >
virtual ~Interface ()=default
 
virtual std::size_t get_data_size () const
 
virtual const void * read_data (const typename ParticleHandler< dim >::particle_iterator &particle, const void *data)
 
virtual void * write_data (const typename ParticleHandler< dim >::particle_iterator &particle, void *data) const
 
virtual void parse_parameters (ParameterHandler &prm)
 
- 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
 
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 Attributes

static const unsigned int n_integrator_properties = 4*dim
 

Private Attributes

unsigned int integrator_substep
 
std::array< unsigned int, 4 > property_indices
 

Additional Inherited Members

- Static Public Member Functions inherited from aspect::Particle::Integrator::Interface< dim >
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::Particle::Integrator::RK4< dim >

Runge Kutta fourth order integrator. This scheme requires storing the original location and intermediate k1, k2, k3 values in the particle properties.

Definition at line 40 of file rk_4.h.

Constructor & Destructor Documentation

§ RK4()

template<int dim>
aspect::Particle::Integrator::RK4< dim >::RK4 ( )

Member Function Documentation

§ initialize()

template<int dim>
void aspect::Particle::Integrator::RK4< dim >::initialize ( )
overridevirtual

Look up where the RK4 data is stored. Done once and cached to avoid repeated lookups.

Reimplemented from aspect::Particle::Integrator::Interface< dim >.

§ local_integrate_step()

template<int dim>
void aspect::Particle::Integrator::RK4< dim >::local_integrate_step ( const typename ParticleHandler< dim >::particle_iterator &  begin_particle,
const typename ParticleHandler< dim >::particle_iterator &  end_particle,
const std::vector< Tensor< 1, dim >> &  old_velocities,
const std::vector< Tensor< 1, dim >> &  velocities,
const double  dt 
)
overridevirtual

Perform an integration step of moving the particles of one cell by the specified timestep dt. This class implements a Runge-Kutta integration scheme that is fourth order accurate in space.

Parameters
[in]begin_particleAn iterator to the first particle to be moved.
[in]end_particleAn iterator to the last particle to be moved.
[in]old_velocitiesThe velocities at t_n, i.e. before the particle movement, for all particles between begin_particle and end_particle at their current position.
[in]velocitiesThe velocities at the particle positions at t_{n+1}, i.e. after the particle movement. Note that this is the velocity at the old positions, but at the new time. It is the responsibility of this function to compute the new location of the particles.
[in]dtThe length of the integration timestep.

Implements aspect::Particle::Integrator::Interface< dim >.

§ new_integration_step()

template<int dim>
bool aspect::Particle::Integrator::RK4< dim >::new_integration_step ( )
overridevirtual

This function is called at the end of every integration step. For the current class it will either increment the step variable and return true to request another integration step, or reset the step variable to 0 and return false if we are already in step 4.

Returns
This function returns true if the integrator requires another integration step. The particle integration will continue to start new integration steps until this function returns false.

Reimplemented from aspect::Particle::Integrator::Interface< dim >.

§ required_solution_vectors()

template<int dim>
std::array<bool, 3> aspect::Particle::Integrator::RK4< dim >::required_solution_vectors ( ) const
overridevirtual

Return a list of boolean values indicating which solution vectors are required for the integration. The first entry indicates if the particle integrator requires the solution vector at the old old time (k-1), the second entry indicates if the particle integrator requires the solution vector at the old time (k), and the third entry indicates if the particle integrator requires the solution vector at the new time (k+1).

The RK4 integrator requires the solution vector at the old time (k) for the first integration step, the solution vector at both the old and new time for the second and third integration steps and the solution vector at the new time (k+1) for the fourth integration step.

Implements aspect::Particle::Integrator::Interface< dim >.

Member Data Documentation

§ n_integrator_properties

template<int dim>
const unsigned int aspect::Particle::Integrator::RK4< dim >::n_integrator_properties = 4*dim
static

We need to tell the property manager how many intermediate properties this integrator requires, so that it can allocate sufficient space for each particle. However, the integrator is not created at the time the property manager is set up and we can not reverse the order of creation, because the integrator needs to know where to store its properties, which requires the property manager to be finished setting up properties. Luckily the number of properties is constant, so we can make it a static property of this class. Therefore, the property manager can access this variable even before any object is constructed.

The Runge-Kutta 4 integrator requires 4 tensors with dim components each.

Definition at line 117 of file rk_4.h.

§ integrator_substep

template<int dim>
unsigned int aspect::Particle::Integrator::RK4< dim >::integrator_substep
private

The current integration step, i.e. for RK4 a number between 0 and 3.

Definition at line 124 of file rk_4.h.

§ property_indices

template<int dim>
std::array<unsigned int,4> aspect::Particle::Integrator::RK4< dim >::property_indices
private

The location of the 4 RK4 data fields stored in the particle properties.

Definition at line 129 of file rk_4.h.


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