ASPECT
|
Classes | |
struct | BaseElements |
struct | BlockIndices |
struct | ComponentIndices |
struct | ComponentMasks |
struct | Extractors |
struct | FaceQuadratures |
struct | IndexSets |
struct | PolynomialDegree |
struct | Quadratures |
Public Member Functions | |
Introspection (const std::vector< VariableDeclaration< dim >> &variables, const Parameters< dim > ¶meters) | |
~Introspection () | |
const std::vector< unsigned int > & | get_composition_base_element_indices () const |
const std::vector< unsigned int > & | get_compositional_field_indices_with_base_element (const unsigned int base_element_index) const |
unsigned int | compositional_index_for_name (const std::string &name) const |
std::string | name_for_compositional_index (const unsigned int index) const |
const std::vector< std::string > & | get_composition_names () const |
const std::vector< CompositionalFieldDescription > & | get_composition_descriptions () const |
const std::vector< std::string > & | chemical_composition_field_names () const |
const std::vector< unsigned int > & | chemical_composition_field_indices () const |
unsigned int | n_chemical_composition_fields () const |
bool | composition_type_exists (const CompositionalFieldDescription::Type &type) const |
unsigned int | find_composition_type (const CompositionalFieldDescription::Type &type) const |
bool | compositional_name_exists (const std::string &name) const |
const std::vector< unsigned int > & | get_indices_for_fields_of_type (const CompositionalFieldDescription::Type &type) const |
const std::vector< std::string > & | get_names_for_fields_of_type (const CompositionalFieldDescription::Type &type) const |
unsigned int | get_number_of_fields_of_type (const CompositionalFieldDescription::Type &type) const |
bool | is_stokes_component (const unsigned int component_index) const |
bool | is_composition_component (const unsigned int component_index) const |
Public Member Functions inherited from aspect::FEVariableCollection< dim > | |
FEVariableCollection () | |
FEVariableCollection (const std::vector< VariableDeclaration< dim >> &variable_definitions) | |
void | initialize (const std::vector< VariableDeclaration< dim >> &variable_definitions) |
const FEVariable< dim > & | variable (const std::string &name) const |
std::vector< const FEVariable< dim > * > | variables_with_name (const std::string &name) const |
bool | variable_exists (const std::string &name) const |
const std::vector< FEVariable< dim > > & | get_variables () const |
unsigned int | n_components () const |
unsigned int | n_blocks () const |
const std::vector< const FiniteElement< dim > * > & | get_fes () const |
const std::vector< unsigned int > & | get_multiplicities () const |
const std::vector< unsigned int > & | get_components_to_blocks () const |
Public Attributes | |
Things that are independent of the current mesh | |
const unsigned int | n_components |
const unsigned int | n_compositional_fields |
const bool | use_discontinuous_temperature_discretization |
const std::vector< bool > | use_discontinuous_composition_discretization |
const ComponentIndices | component_indices |
const unsigned int | n_blocks |
const BlockIndices | block_indices |
const Extractors | extractors |
const BaseElements | base_elements |
const PolynomialDegree | polynomial_degree |
const Quadratures | quadratures |
const FaceQuadratures | face_quadratures |
const ComponentMasks | component_masks |
Things that depend on the current mesh | |
std::vector< types::global_dof_index > | system_dofs_per_block |
IndexSets | index_sets |
Parameters< dim >::AdvectionFieldMethod::Kind | temperature_method |
std::vector< typename Parameters< dim >::AdvectionFieldMethod::Kind > | compositional_field_methods |
Private Attributes | |
std::vector< std::string > | composition_names |
std::vector< CompositionalFieldDescription > | composition_descriptions |
std::vector< std::vector< std::string > > | composition_names_for_type |
std::vector< std::vector< unsigned int > > | composition_indices_for_type |
std::vector< unsigned int > | composition_base_element_indices |
std::map< unsigned int, std::vector< unsigned int > > | compositional_field_indices_with_base_element |
Additional Inherited Members | |
Protected Attributes inherited from aspect::FEVariableCollection< dim > | |
std::vector< FEVariable< dim > > | variables |
unsigned int | n_components_ |
unsigned int | n_blocks_ |
std::vector< const FiniteElement< dim > * > | fes |
std::vector< unsigned int > | multiplicities |
std::vector< unsigned int > | components_to_blocks |
The introspection class provides information about the simulation as a whole. In particular, it provides symbolic names for things like the velocity, pressure and other variables, along with their corresponding vector and scalar FEValues extractors, component masks, etc.
The purpose of this class is primarily to provide these symbolic names so that we do not have to use implicit knowledge about the ordering of variables (e.g., earlier versions of ASPECT had many places where we built a scalar FEValues extractor at component 'dim' since that is where we knew that the pressure lies in the finite element; this kind of implicit knowledge is no longer necessary with the Introspection class). The Introspection class is used both by the Simulator class itself, but is also exported to plugins via the SimulatorAccess class.
Definition at line 128 of file introspection.h.
aspect::Introspection< dim >::Introspection | ( | const std::vector< VariableDeclaration< dim >> & | variables, |
const Parameters< dim > & | parameters | ||
) |
Constructor.
aspect::Introspection< dim >::~Introspection | ( | ) |
Destructor.
const std::vector<unsigned int>& aspect::Introspection< dim >::get_composition_base_element_indices | ( | ) | const |
Return a vector that contains the base element indices of the deal.II FiniteElement that are used for compositional fields. Note that compositional fields can share the same base element, so this vector can (and usually will) be smaller than the number of compositional fields. The function get_compositional_field_indices_with_base_element() can be used to translate from base element index to all compositional field indices using the specified base element. If several fields are the finite element type (same degree and both continuous or both discontinuous), they share base elements. If you have no compositional fields, the vector returned has length 0. If all compositional fields have the same finite element space, the length is 1.
const std::vector<unsigned int>& aspect::Introspection< dim >::get_compositional_field_indices_with_base_element | ( | const unsigned int | base_element_index | ) | const |
Return a vector with all compositional field indices that belong to a given base element index as returned by get_composition_base_element_indices(). The indices returned are therefore between 0 and n_compositional_fields-1. If you have a single compositional field, this function returns {0} when passing in base_element_index=0.
unsigned int aspect::Introspection< dim >::compositional_index_for_name | ( | const std::string & | name | ) | const |
A function that gets the name of a compositional field as an input parameter and returns its index. If the name is not found, an exception is thrown.
name | The name of compositional field (as specified in the input file) |
std::string aspect::Introspection< dim >::name_for_compositional_index | ( | const unsigned int | index | ) | const |
A function that gets the index of a compositional field as an input parameter and returns its name.
index | The index of compositional field |
const std::vector<std::string>& aspect::Introspection< dim >::get_composition_names | ( | ) | const |
A function that returns the full list of compositional field names.
const std::vector<CompositionalFieldDescription>& aspect::Introspection< dim >::get_composition_descriptions | ( | ) | const |
A function that returns the full vector of compositional field descriptions.
const std::vector<std::string>& aspect::Introspection< dim >::chemical_composition_field_names | ( | ) | const |
A function that returns the names of compositional fields that correspond to chemical compositions.
This function is shorthand for get_names_for_fields_of_type(CompositionalFieldDescription::chemical_composition).
const std::vector<unsigned int>& aspect::Introspection< dim >::chemical_composition_field_indices | ( | ) | const |
A function that returns the indices of compositional fields that correspond to chemical compositions.
This function is shorthand for get_indices_for_fields_of_type(CompositionalFieldDescription::chemical_composition).
unsigned int aspect::Introspection< dim >::n_chemical_composition_fields | ( | ) | const |
A function that returns the number of compositional fields that correspond to chemical compositions.
This function is shorthand for get_number_of_fields_of_type(CompositionalFieldDescription::chemical_composition).
bool aspect::Introspection< dim >::composition_type_exists | ( | const CompositionalFieldDescription::Type & | type | ) | const |
A function that gets the type of a compositional field as an input parameter and returns if any compositional field of that type is used in this simulation.
type | The type of compositional field (as specified in the input file) |
unsigned int aspect::Introspection< dim >::find_composition_type | ( | const CompositionalFieldDescription::Type & | type | ) | const |
A function that gets the type of a compositional field as an input parameter and returns the index of the first compositional field of this type used in this simulation. If no such field is found, the function returns the number of compositional fields.
type | The type of compositional field (as specified in the input file) |
bool aspect::Introspection< dim >::compositional_name_exists | ( | const std::string & | name | ) | const |
A function that gets the name of a compositional field as an input parameter and returns if the compositional field is used in this simulation.
name | The name of compositional field (as specified in the input file) |
const std::vector<unsigned int>& aspect::Introspection< dim >::get_indices_for_fields_of_type | ( | const CompositionalFieldDescription::Type & | type | ) | const |
Get the indices of the compositional fields which are of a particular type (chemical composition, porosity, etc.).
const std::vector<std::string>& aspect::Introspection< dim >::get_names_for_fields_of_type | ( | const CompositionalFieldDescription::Type & | type | ) | const |
Get the names of the compositional fields which are of a particular type (chemical composition, porosity, etc.).
unsigned int aspect::Introspection< dim >::get_number_of_fields_of_type | ( | const CompositionalFieldDescription::Type & | type | ) | const |
Get the number of compositional fields which are of a particular type (chemical composition, porosity, etc.).
bool aspect::Introspection< dim >::is_stokes_component | ( | const unsigned int | component_index | ) | const |
A function that gets a component index as an input parameter and returns if the component is one of the stokes system (i.e. if it is the pressure or one of the velocity components).
component_index | The component index to check. |
bool aspect::Introspection< dim >::is_composition_component | ( | const unsigned int | component_index | ) | const |
A function that gets a component index as an input parameter and returns if the component is one of the compositional fields.
component_index | The component index to check. |
const unsigned int aspect::Introspection< dim >::n_components |
The number of vector components used by the finite element description of this problem. It equals \(d+2+n_c\) where \(d\) is the dimension and equals the number of velocity components, and \(n_c\) is the number of advected (compositional) fields. The remaining components are the scalar pressure and temperature fields.
Definition at line 156 of file introspection.h.
const unsigned int aspect::Introspection< dim >::n_compositional_fields |
The number of compositional fields.
Definition at line 161 of file introspection.h.
const bool aspect::Introspection< dim >::use_discontinuous_temperature_discretization |
A variable that holds whether the temperature field should use a discontinuous discretization.
Definition at line 167 of file introspection.h.
const std::vector<bool> aspect::Introspection< dim >::use_discontinuous_composition_discretization |
A variable that holds whether the composition field(s) should use a discontinuous discretization.
Definition at line 173 of file introspection.h.
const ComponentIndices aspect::Introspection< dim >::component_indices |
A variable that enumerates the vector components of the finite element that correspond to each of the variables in this problem.
Definition at line 190 of file introspection.h.
const unsigned int aspect::Introspection< dim >::n_blocks |
The number of vector blocks. This equals \(3+n_c\) where, in comparison to the n_components field, the velocity components form a single block.
Definition at line 197 of file introspection.h.
const BlockIndices aspect::Introspection< dim >::block_indices |
A variable that enumerates the vector blocks of the finite element that correspond to each of the variables in this problem.
Definition at line 224 of file introspection.h.
const Extractors aspect::Introspection< dim >::extractors |
A variable that contains extractors for every block of the finite element used in the overall description.
Definition at line 244 of file introspection.h.
const BaseElements aspect::Introspection< dim >::base_elements |
A variable that enumerates the base elements of the finite element that correspond to each of the variables in this problem.
Definition at line 267 of file introspection.h.
const PolynomialDegree aspect::Introspection< dim >::polynomial_degree |
A variable that enumerates the polynomial degree of the finite element that correspond to each of the variables in this problem.
Definition at line 289 of file introspection.h.
const Quadratures aspect::Introspection< dim >::quadratures |
A variable that enumerates the polynomial degree of the finite element that correspond to each of the variables in this problem.
Definition at line 324 of file introspection.h.
const FaceQuadratures aspect::Introspection< dim >::face_quadratures |
A variable that enumerates the polynomial degree of the finite element that correspond to each of the variables in this problem.
Definition at line 347 of file introspection.h.
const ComponentMasks aspect::Introspection< dim >::component_masks |
A variable that contains component masks for each of the variables in this problem. Component masks are a deal.II concept, see the deal.II glossary.
Definition at line 369 of file introspection.h.
std::vector<types::global_dof_index> aspect::Introspection< dim >::system_dofs_per_block |
A variable that describes how many of the degrees of freedom on the current mesh belong to each of the n_blocks blocks of the finite element.
Definition at line 384 of file introspection.h.
IndexSets aspect::Introspection< dim >::index_sets |
A variable that contains index sets describing which of the globally enumerated degrees of freedom are owned by the current processor in a parallel computation.
Definition at line 453 of file introspection.h.
Parameters<dim>::AdvectionFieldMethod::Kind aspect::Introspection< dim >::temperature_method |
A variable that contains the field method for the temperature field and is used to determine how to solve it when solving a timestep.
Definition at line 459 of file introspection.h.
std::vector<typename Parameters<dim>::AdvectionFieldMethod::Kind> aspect::Introspection< dim >::compositional_field_methods |
A vector that contains a field method for every compositional field and is used to determine how to solve a particular field when solving a timestep.
Definition at line 466 of file introspection.h.
|
private |
A vector that stores the names of the compositional fields that will be used in the simulation.
Definition at line 640 of file introspection.h.
|
private |
A vector that stores descriptions of each compositional field, including its type (i.e. whether the compositional field corresponds to chemical composition, porosity etc.).
Definition at line 647 of file introspection.h.
|
private |
A vector of vectors of composition names that stores the names of the compositional fields corresponding to each field type given in CompositionalFieldDescription.
Definition at line 654 of file introspection.h.
|
private |
A vector of vectors of composition indices that stores the indices of the compositional fields corresponding to each field type given in CompositionalFieldDescription.
Definition at line 661 of file introspection.h.
|
private |
List of base element indices used by compositional fields. Cached result returned by get_composition_base_element_indices().
Definition at line 667 of file introspection.h.
|
private |
Map base_element_index to list of compositional field indices that use that base element. Cached result returned by get_compositional_field_indices_with_base_element();
Definition at line 674 of file introspection.h.