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

Public Member Functions

 AsciiDataBoundary ()
 
virtual void initialize (const std::set< types::boundary_id > &boundary_ids, const unsigned int components)
 
void update ()
 
double get_data_component (const types::boundary_id boundary_indicator, const Point< dim > &position, const unsigned int component) const
 
double get_maximum_component_value (const types::boundary_id boundary_indicator, const unsigned int component) const
 
Tensor< 1, dim-1 > vector_gradient (const types::boundary_id boundary_indicator, const Point< dim > &p, const unsigned int component) const
 
void parse_parameters (ParameterHandler &prm, const std::string &subsection_name="Ascii data model", const bool parse_time_dependent_parameters=true)
 
- Public Member Functions inherited from aspect::Utilities::AsciiDataBase< dim >
 AsciiDataBase ()
 
void parse_parameters (ParameterHandler &prm, const std::string &subsection_name="Ascii data model")
 
- 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 Member Functions

static void declare_parameters (ParameterHandler &prm, const std::string &default_directory, const std::string &default_filename, const std::string &subsection_name="Ascii data model", const bool declare_time_dependent_parameters=true)
 
- Static Public Member Functions inherited from aspect::Utilities::AsciiDataBase< dim >
static void declare_parameters (ParameterHandler &prm, const std::string &default_directory, const std::string &default_filename, const std::string &subsection_name="Ascii data model")
 
- 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)
 

Protected Member Functions

void update_data (const types::boundary_id boundary_id, const bool reload_both_files)
 
void end_time_dependence ()
 
std::string create_filename (const int filenumber, const types::boundary_id boundary_id) const
 

Protected Attributes

int current_file_number
 
int first_data_file_number
 
bool decreasing_file_order
 
double data_file_time_step
 
double time_weight
 
bool time_dependent
 
std::map< types::boundary_id, std::unique_ptr< aspect::Utilities::StructuredDataLookup< dim-1 > > > lookups
 
std::map< types::boundary_id, std::unique_ptr< aspect::Utilities::StructuredDataLookup< dim-1 > > > old_lookups
 

Additional Inherited Members

- Public Attributes inherited from aspect::Utilities::AsciiDataBase< dim >
std::string data_directory
 
std::string data_file_name
 
double scale_factor
 

Detailed Description

template<int dim>
class aspect::Utilities::AsciiDataBoundary< dim >

A base class that implements boundary conditions determined from a AsciiData input file.

Definition at line 363 of file structured_data.h.

Constructor & Destructor Documentation

§ AsciiDataBoundary()

Constructor

Member Function Documentation

§ initialize()

template<int dim>
virtual void aspect::Utilities::AsciiDataBoundary< dim >::initialize ( const std::set< types::boundary_id > &  boundary_ids,
const unsigned int  components 
)
virtual

Initialization function. This function is called once at the beginning of the program. Checks preconditions.

§ update()

template<int dim>
void aspect::Utilities::AsciiDataBoundary< dim >::update ( )

A function that is called at the beginning of each time step. For the current plugin, this function loads the next data files if necessary and outputs a warning if the end of the set of data files is reached.

§ get_data_component()

template<int dim>
double aspect::Utilities::AsciiDataBoundary< dim >::get_data_component ( const types::boundary_id  boundary_indicator,
const Point< dim > &  position,
const unsigned int  component 
) const

Returns the data component at the given position.

§ get_maximum_component_value()

template<int dim>
double aspect::Utilities::AsciiDataBoundary< dim >::get_maximum_component_value ( const types::boundary_id  boundary_indicator,
const unsigned int  component 
) const

Returns the maximum value of the given data component.

§ vector_gradient()

template<int dim>
Tensor<1,dim-1> aspect::Utilities::AsciiDataBoundary< dim >::vector_gradient ( const types::boundary_id  boundary_indicator,
const Point< dim > &  p,
const unsigned int  component 
) const

Return the gradients of the parameters from the parameter file.

§ declare_parameters()

template<int dim>
static void aspect::Utilities::AsciiDataBoundary< dim >::declare_parameters ( ParameterHandler &  prm,
const std::string &  default_directory,
const std::string &  default_filename,
const std::string &  subsection_name = "Ascii data model",
const bool  declare_time_dependent_parameters = true 
)
static

Declare the parameters all derived classes take from input files.

Parameters
prmThe parameter handler in which the parameters are declared.
default_directoryThe default value for the data directory parameter.
default_filenameThe default value for the filename parameter.
subsection_nameThe name of the parameter file subsection all parameters will be declared in. The function will enter this subsection, declare all parameters, then leave this subsection, to return prm in the same subsection it was in before.
declare_time_dependent_parametersWhether to declare the parameter that are only needed for time dependent AsciiDataBoundary objects. If the caller already knows time dependence is not supported for the current application, disabling this parameter avoids introducing these parameters to the parameter handler.

§ parse_parameters()

template<int dim>
void aspect::Utilities::AsciiDataBoundary< dim >::parse_parameters ( ParameterHandler &  prm,
const std::string &  subsection_name = "Ascii data model",
const bool  parse_time_dependent_parameters = true 
)

Read the parameters from the parameter file.

Parameters
prmThe parameter handler from which the parameters are parsed.
subsection_nameThe name of the parameter file subsection all parameters will be parsed from. The function will enter this subsection, parse all parameters, then leave this subsection, to return prm in the same subsection it was in before.
parse_time_dependent_parametersWhether to parse the parameter that are only needed for time dependent AsciiDataBoundary objects. This parameter always needs to be set to the same value that was handed over to declare_parameters().

§ update_data()

template<int dim>
void aspect::Utilities::AsciiDataBoundary< dim >::update_data ( const types::boundary_id  boundary_id,
const bool  reload_both_files 
)
protected

Handles the update of the data in lookup.

§ end_time_dependence()

template<int dim>
void aspect::Utilities::AsciiDataBoundary< dim >::end_time_dependence ( )
protected

Handles settings and user notification in case the time-dependent part of the boundary condition is over.

§ create_filename()

template<int dim>
std::string aspect::Utilities::AsciiDataBoundary< dim >::create_filename ( const int  filenumber,
const types::boundary_id  boundary_id 
) const
protected

Create a filename out of the name template.

Member Data Documentation

§ current_file_number

template<int dim>
int aspect::Utilities::AsciiDataBoundary< dim >::current_file_number
protected

A variable that stores the currently used data file of a series. It gets updated if necessary by update().

Definition at line 459 of file structured_data.h.

§ first_data_file_number

template<int dim>
int aspect::Utilities::AsciiDataBoundary< dim >::first_data_file_number
protected

Number of the first data file to be loaded when the model time is larger than 'First data file model time'.

Definition at line 465 of file structured_data.h.

§ decreasing_file_order

template<int dim>
bool aspect::Utilities::AsciiDataBoundary< dim >::decreasing_file_order
protected

In some cases the boundary files are not numbered in increasing but in decreasing order (e.g. 'Ma BP'). If this flag is set to 'True' the plugin will first load the file with the number 'First data file number' and decrease the file number during the model run.

Definition at line 473 of file structured_data.h.

§ data_file_time_step

template<int dim>
double aspect::Utilities::AsciiDataBoundary< dim >::data_file_time_step
protected

Time in model units (depends on other model inputs) between two data files.

Definition at line 479 of file structured_data.h.

§ time_weight

template<int dim>
double aspect::Utilities::AsciiDataBoundary< dim >::time_weight
protected

Weight between data file n and n+1 while the current time is between the two values t(n) and t(n+1).

Definition at line 485 of file structured_data.h.

§ time_dependent

template<int dim>
bool aspect::Utilities::AsciiDataBoundary< dim >::time_dependent
protected

State whether we have time_dependent boundary conditions. Switched off after finding no more data files to suppress attempts to read in new files.

Definition at line 492 of file structured_data.h.

§ lookups

template<int dim>
std::map<types::boundary_id, std::unique_ptr<aspect::Utilities::StructuredDataLookup<dim-1> > > aspect::Utilities::AsciiDataBoundary< dim >::lookups
protected

Map between the boundary id and an object that reads and processes data we get from text files.

Definition at line 499 of file structured_data.h.

§ old_lookups

template<int dim>
std::map<types::boundary_id, std::unique_ptr<aspect::Utilities::StructuredDataLookup<dim-1> > > aspect::Utilities::AsciiDataBoundary< dim >::old_lookups
protected

Map between the boundary id and the old data objects.

Definition at line 505 of file structured_data.h.


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