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

Public Member Functions

 RK2 ()
 
virtual 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)
 
virtual bool new_integration_step ()
 
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
 
- Public Member Functions inherited from aspect::Particle::Integrator::Interface< dim >
virtual ~Interface ()
 
virtual void parse_parameters (ParameterHandler &prm)
 
- Public Member Functions inherited from aspect::SimulatorAccess< dim >
 SimulatorAccess ()
 
 SimulatorAccess (const Simulator< dim > &simulator_object)
 
virtual ~SimulatorAccess ()
 
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
 
TimerOutputget_computing_timer () const
 
const ConditionalOStreamget_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::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
 
void compute_material_model_input_values (const LinearAlgebra::BlockVector &input_solution, const FEValuesBase< dim, dim > &input_finite_element_values, const typename DoFHandler< dim >::active_cell_iterator &cell, const bool compute_strainrate, MaterialModel::MaterialModelInputs< dim > &material_model_inputs) 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
 
const InitialTemperature::Manager< dim > & get_initial_temperature_manager () const
 
DEAL_II_DEPRECATED const InitialComposition::Interface< dim > & get_initial_composition () 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_free_surface_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 WorldBuilder::World & get_world_builder () const
 
const FreeSurfaceHandler< dim > & get_free_surface_handler () const
 
const LateralAveraging< dim > & get_lateral_averaging () const
 
const ConstraintMatrix & get_current_constraints () const
 
bool simulator_is_initialized () const
 
double get_pressure_scaling () const
 
bool pressure_rhs_needs_compatibility_modification () const
 
bool model_has_prescribed_stokes_solution () const
 
TableHandlerget_statistics_object () const
 
template<typename PostprocessorType >
DEAL_II_DEPRECATED PostprocessorType * find_postprocessor () const
 
const Postprocess::Manager< dim > & get_postprocess_manager () const
 

Private Attributes

unsigned int integrator_substep
 
std::map< types::particle_index, Point< dim > > loc0
 

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::RK2< dim >

Runge Kutta second order integrator. This scheme requires storing the original location, and the read/write_data functions reflect this.

Definition at line 42 of file rk_2.h.

Constructor & Destructor Documentation

§ RK2()

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

Member Function Documentation

§ local_integrate_step()

template<int dim>
virtual void aspect::Particle::Integrator::RK2< 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 
)
virtual

Perform an integration step of moving the particles of one cell by the specified timestep dt. This function implements a second-order accurate Runge-Kutta integration scheme.

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>
virtual bool aspect::Particle::Integrator::RK2< dim >::new_integration_step ( )
virtual

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 2.

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 >.

§ get_data_size()

template<int dim>
virtual std::size_t aspect::Particle::Integrator::RK2< dim >::get_data_size ( ) const
virtual

Return data length of the integration related data required for communication in terms of number of bytes. When data about particles is transported from one processor to another, or stored on disk for snapshots, integrators get the chance to store whatever information they need with each particle. This function returns how many pieces of additional information a concrete integrator class needs to store for each particle.

Returns
The number of bytes required to store the relevant integrator data for one particle.

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

§ read_data()

template<int dim>
virtual const void* aspect::Particle::Integrator::RK2< dim >::read_data ( const typename ParticleHandler< dim >::particle_iterator &  particle,
const void *  data 
)
virtual

Read integration related data for a particle specified by particle_id from the data array. This function is called after transferring a particle to the local domain during an integration step.

Parameters
[in]dataA pointer into the data array. The pointer marks the position where this function starts reading.
[in]particleAn iterator pointing to the particle to read the data for.
Returns
The updated position of the pointer into the data array. The return value is data advanced by get_data_size() bytes.

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

§ write_data()

template<int dim>
virtual void* aspect::Particle::Integrator::RK2< dim >::write_data ( const typename ParticleHandler< dim >::particle_iterator &  particle,
void *  data 
) const
virtual

Write integration related data to a vector for a particle specified by particle_id. This function is called in cases where particles leave the local domain during an integration step to transfer this data to another process.

Parameters
[in]particleAn iterator pointing to the particle to write the data for.
[in]dataA pointer into the array of integrator data.
Returns
The updated position of the pointer into the data array. The return value is data advanced by get_data_size() bytes.

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

Member Data Documentation

§ integrator_substep

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

The current integration step, i.e for RK2 a number that is either 0 or 1.

Definition at line 119 of file rk_2.h.

§ loc0

template<int dim>
std::map<types::particle_index, Point<dim> > aspect::Particle::Integrator::RK2< dim >::loc0
private

The particle location before the first integration step. This is used in the second step and transferred to another process if the particle leaves the domain during the first step.

Definition at line 126 of file rk_2.h.


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