ASPECT
|
Public Member Functions | |
SolutionEvaluator (const SimulatorAccess< dim > &simulator, const UpdateFlags update_flags) | |
void | reinit (const typename DoFHandler< dim >::active_cell_iterator &cell, const ArrayView< Point< dim >> &positions) |
void | evaluate (const ArrayView< double > &solution_values, const std::vector< EvaluationFlags::EvaluationFlags > &evaluation_flags) |
void | get_solution (const unsigned int evaluation_point, const ArrayView< double > &solution, const std::vector< EvaluationFlags::EvaluationFlags > &evaluation_flags) const |
void | get_gradients (const unsigned int evaluation_point, const ArrayView< Tensor< 1, dim >> &gradients, const std::vector< EvaluationFlags::EvaluationFlags > &evaluation_flags) const |
FEPointEvaluation< dim, dim > & | get_velocity_or_fluid_velocity_evaluator (const bool use_fluid_velocity) |
NonMatching::MappingInfo< dim > & | get_mapping_info () |
unsigned int | n_components () const |
Private Attributes | |
NonMatching::MappingInfo< dim > | mapping_info |
FEPointEvaluation< dim, dim > | velocity |
std::unique_ptr< FEPointEvaluation< 1, dim > > | pressure |
FEPointEvaluation< 1, dim > | temperature |
std::vector< std::unique_ptr< internal::DynamicFEPointEvaluation< dim > > > | compositions |
std::unique_ptr< FEPointEvaluation< dim, dim > > | fluid_velocity |
std::unique_ptr< FEPointEvaluation< 1, dim > > | compaction_pressure |
std::unique_ptr< FEPointEvaluation< 1, dim > > | fluid_pressure |
std::array< unsigned int, 3 > | melt_component_indices |
const SimulatorAccess< dim > & | simulator_access |
This class evaluates the solution vector at arbitrary positions inside a cell. This base class only provides the interface for SolutionEvaluatorImplementation. See there for more details.
Definition at line 141 of file solution_evaluator.h.
aspect::SolutionEvaluator< dim >::SolutionEvaluator | ( | const SimulatorAccess< dim > & | simulator, |
const UpdateFlags | update_flags | ||
) |
Constructor. Create the member variables given a simulator access simulator
and a set of update flags update_flags
. The update_flags
control if only the solution or also the gradients should be evaluated.
void aspect::SolutionEvaluator< dim >::reinit | ( | const typename DoFHandler< dim >::active_cell_iterator & | cell, |
const ArrayView< Point< dim >> & | positions | ||
) |
Reinitialize all variables to prepare for evaluation for the given cell
and at the given positions
in the reference coordinate system of that cell.
void aspect::SolutionEvaluator< dim >::evaluate | ( | const ArrayView< double > & | solution_values, |
const std::vector< EvaluationFlags::EvaluationFlags > & | evaluation_flags | ||
) |
Evaluate all variables in the cell and at the positions controlled by a previous call to reinit().
solution_values
contains the values of the degrees of freedom. evaluation_flags
controls if nothing, the solution values, and/or the gradients should be computed. The size of solution_flags
has to be equal to the number of components as returned by n_components().
void aspect::SolutionEvaluator< dim >::get_solution | ( | const unsigned int | evaluation_point, |
const ArrayView< double > & | solution, | ||
const std::vector< EvaluationFlags::EvaluationFlags > & | evaluation_flags | ||
) | const |
Fill solution
with all solution components at the given evaluation_point
. Note that this function only works after a successful call to reinit() and evaluate() because this function only returns the results of the computation that happened in those functions.
evaluation_point | The index of the evaluation point in the positions array. |
solution | The array to fill with the solution values. This array has to be of size n_components(). |
evaluation_flags | The flags that indicate which values should be copied into the solution array. This vector has to be of size n_components(). |
void aspect::SolutionEvaluator< dim >::get_gradients | ( | const unsigned int | evaluation_point, |
const ArrayView< Tensor< 1, dim >> & | gradients, | ||
const std::vector< EvaluationFlags::EvaluationFlags > & | evaluation_flags | ||
) | const |
Fill gradients
with all solution gradients at the given evaluation_point
. Note that this function only works after a successful call to reinit() and evaluate() because this function only returns the results of the computation that happened in those functions.
evaluation_point | The index of the evaluation point in the positions array. |
gradients | The array to fill with the solution gradients. This array has to be of size n_components(). |
evaluation_flags | The flags that indicate which gradients should be copied into the solution array. This vector has to be of size n_components(). |
FEPointEvaluation<dim, dim>& aspect::SolutionEvaluator< dim >::get_velocity_or_fluid_velocity_evaluator | ( | const bool | use_fluid_velocity | ) |
Return the evaluator for velocity or fluid velocity. This is the only information necessary for advecting particles.
NonMatching::MappingInfo<dim>& aspect::SolutionEvaluator< dim >::get_mapping_info | ( | ) |
Return the cached mapping information.
unsigned int aspect::SolutionEvaluator< dim >::n_components | ( | ) | const |
Return the number of components in the solution vector.
|
private |
MappingInfo object for the FEPointEvaluation objects
Definition at line 228 of file solution_evaluator.h.
|
private |
FEPointEvaluation objects for all common components of ASPECT's finite element solution. These objects are used inside of the member functions of this class.
Definition at line 235 of file solution_evaluator.h.
|
private |
Definition at line 236 of file solution_evaluator.h.
|
private |
Definition at line 237 of file solution_evaluator.h.
|
private |
We group compositions by type (FiniteElement) and evaluate them together.
Definition at line 243 of file solution_evaluator.h.
|
private |
Pointers to FEPointEvaluation objects for all melt components of ASPECT's finite element solution, which only point to valid objects in case we use melt transport. Other documentation like for the objects directly above.
Definition at line 251 of file solution_evaluator.h.
|
private |
Definition at line 252 of file solution_evaluator.h.
|
private |
Definition at line 253 of file solution_evaluator.h.
|
private |
The component indices for the three melt formulation variables fluid velocity, compaction pressure, and fluid pressure (in this order). They are cached to avoid repeated expensive lookups.
Definition at line 261 of file solution_evaluator.h.
|
private |
Reference to the active simulator access object. Provides access to the general simulation variables.
Definition at line 267 of file solution_evaluator.h.