|
adamantine
|
#include <ThermalOperator.hh>
Public Member Functions | |
| ThermalOperator (MPI_Comm const &communicator, Boundary const &boundary, MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType > &material_properties, std::vector< std::shared_ptr< HeatSource< dim >>> const &heat_sources) | |
| void | reinit (dealii::DoFHandler< dim > const &dof_handler, dealii::AffineConstraints< double > const &affine_constraints, dealii::hp::QCollection< 1 > const &quad) override |
| void | compute_inverse_mass_matrix (dealii::DoFHandler< dim > const &dof_handler, dealii::AffineConstraints< double > const &affine_constraints) override |
| void | clear () override |
| dealii::types::global_dof_index | m () const override |
| dealii::types::global_dof_index | n () const override |
| std::shared_ptr< dealii::LA::distributed::Vector< double, MemorySpaceType > > | get_inverse_mass_matrix () const override |
| dealii::MatrixFree< dim, double > const & | get_matrix_free () const |
| void | vmult (dealii::LA::distributed::Vector< double, MemorySpaceType > &dst, dealii::LA::distributed::Vector< double, MemorySpaceType > const &src) const override |
| void | vmult_add (dealii::LA::distributed::Vector< double, MemorySpaceType > &dst, dealii::LA::distributed::Vector< double, MemorySpaceType > const &src) const override |
| void | initialize_dof_vector (dealii::LA::distributed::Vector< double, MemorySpaceType > &vector) const override |
| void | get_state_from_material_properties () override |
| void | set_state_to_material_properties () override |
| void | set_material_deposition_orientation (std::vector< double > const &deposition_cos, std::vector< double > const &deposition_sin) override |
| void | set_time_and_source_height (double t, double height) override |
Public Member Functions inherited from adamantine::ThermalOperatorBase< dim, MemorySpaceType > | |
| ThermalOperatorBase ()=default | |
| virtual | ~ThermalOperatorBase ()=default |
Private Member Functions | |
| void | update_state_ratios ([[maybe_unused]] unsigned int cell, [[maybe_unused]] unsigned int q, [[maybe_unused]] dealii::VectorizedArray< double > temperature, std::array< dealii::VectorizedArray< double >, MaterialStates::n_material_states > &state_ratios) const |
| void | update_face_state_ratios ([[maybe_unused]] unsigned int face, [[maybe_unused]] unsigned int q, [[maybe_unused]] dealii::VectorizedArray< double > temperature, std::array< dealii::VectorizedArray< double >, MaterialStates::n_material_states > &state_ratios) const |
| dealii::VectorizedArray< double > | get_inv_rho_cp (std::array< dealii::types::material_id, dealii::VectorizedArray< double >::size()> const &material_id, std::array< dealii::VectorizedArray< double >, MaterialStates::n_material_states > const &state_ratios, dealii::VectorizedArray< double > const &temperature, dealii::AlignedVector< dealii::VectorizedArray< double >> const &temperature_powers) const |
| void | cell_local_apply (dealii::MatrixFree< dim, double > const &data, dealii::LA::distributed::Vector< double, MemorySpaceType > &dst, dealii::LA::distributed::Vector< double, MemorySpaceType > const &src, std::pair< unsigned int, unsigned int > const &cell_range) const |
| void | face_local_apply (dealii::MatrixFree< dim, double > const &data, dealii::LA::distributed::Vector< double, MemorySpaceType > &dst, dealii::LA::distributed::Vector< double, MemorySpaceType > const &src, std::pair< unsigned int, unsigned int > const &face_range) const |
| void | cell_local_mass (dealii::MatrixFree< dim, double > const &data, dealii::LA::distributed::Vector< double, MemorySpaceType > &dst, dealii::LA::distributed::Vector< double, MemorySpaceType > const &src, std::pair< unsigned int, unsigned int > const &cell_range) const |
Private Attributes | |
| MPI_Comm const & | _communicator |
| bool | _adiabatic_only_bc = true |
| double | _current_source_height = 0. |
| Boundary | _boundary |
| dealii::MatrixFree< dim, double >::AdditionalData | _matrix_free_data |
| dealii::Table< 2, dealii::VectorizedArray< double > > | _thermal_conductivity |
| MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType > & | _material_properties |
| std::vector< std::shared_ptr< HeatSource< dim > > > | _heat_sources |
| dealii::MatrixFree< dim, double > | _matrix_free |
| dealii::AffineConstraints< double > const * | _affine_constraints |
| std::shared_ptr< dealii::LA::distributed::Vector< double, MemorySpaceType > > | _inverse_mass_matrix |
| std::map< typename dealii::DoFHandler< dim >::cell_iterator, std::pair< unsigned int, unsigned int > > | _cell_it_to_mf_cell_map |
| dealii::Table< 2, dealii::VectorizedArray< double > > | _liquid_ratio |
| dealii::Table< 2, dealii::VectorizedArray< double > > | _powder_ratio |
| dealii::Table< 2, dealii::VectorizedArray< double > > | _face_powder_ratio |
| dealii::Table< 2, std::array< dealii::types::material_id, dealii::VectorizedArray< double >::size()> > | _material_id |
| dealii::Table< 2, std::array< dealii::types::material_id, dealii::VectorizedArray< double >::size()> > | _face_material_id |
| dealii::Table< 2, dealii::VectorizedArray< double > > | _deposition_cos |
| dealii::Table< 2, dealii::VectorizedArray< double > > | _deposition_sin |
This class is the operator associated with the heat equation, i.e., vmult performs
.
Definition at line 26 of file ThermalOperator.hh.
| adamantine::ThermalOperator< dim, n_materials, use_table, p_order, fe_degree, MaterialStates, MemorySpaceType >::ThermalOperator | ( | MPI_Comm const & | communicator, |
| Boundary const & | boundary, | ||
| MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType > & | material_properties, | ||
| std::vector< std::shared_ptr< HeatSource< dim >>> const & | heat_sources | ||
| ) |
Definition at line 29 of file ThermalOperator.templates.hh.
|
private |
Apply the operator on a given set of quadrature points inside each cell.
Definition at line 481 of file ThermalOperator.templates.hh.
|
private |
Apply the mass operator on a given set of quadrature points.
Definition at line 83 of file ThermalOperator.templates.hh.
|
overridevirtual |
Resize the inverse of the mass matrix to zero.
Implements adamantine::ThermalOperatorBase< dim, MemorySpaceType >.
Definition at line 168 of file ThermalOperator.templates.hh.
|
overridevirtual |
Compute the inverse of the mass matrix.
Implements adamantine::ThermalOperatorBase< dim, MemorySpaceType >.
Definition at line 120 of file ThermalOperator.templates.hh.
|
private |
Apply the operator on a given set of quadrature points on each face.
Definition at line 616 of file ThermalOperator.templates.hh.
|
private |
Return the value of
for a given matrix-free cell/face and quadrature point.
Definition at line 418 of file ThermalOperator.templates.hh.
|
inlineoverridevirtual |
Return a shared pointer to the inverse of the mass matrix.
Implements adamantine::ThermalOperatorBase< dim, MemorySpaceType >.
Definition at line 272 of file ThermalOperator.hh.
|
inline |
Return a shared pointer to the underlying MatrixFree object.
Definition at line 281 of file ThermalOperator.hh.
|
overridevirtual |
Implements adamantine::ThermalOperatorBase< dim, MemorySpaceType >.
Definition at line 751 of file ThermalOperator.templates.hh.
|
inlineoverridevirtual |
Implements adamantine::ThermalOperatorBase< dim, MemorySpaceType >.
Definition at line 289 of file ThermalOperator.hh.
|
inlineoverridevirtual |
Return the dimension of the codomain (or range) space. To remember: the matrix is of dimension m×n.
Implements adamantine::ThermalOperatorBase< dim, MemorySpaceType >.
Definition at line 254 of file ThermalOperator.hh.
|
inlineoverridevirtual |
Return the dimension of the domain space. To remember: the matrix is of dimension m×n.
Implements adamantine::ThermalOperatorBase< dim, MemorySpaceType >.
Definition at line 263 of file ThermalOperator.hh.
|
overridevirtual |
Associate the AffineConstraints<double> and the MatrixFree objects to the underlying Triangulation.
Implements adamantine::ThermalOperatorBase< dim, MemorySpaceType >.
Definition at line 58 of file ThermalOperator.templates.hh.
|
overridevirtual |
Set the deposition cosine and sine angles and convert the data from std::vector to dealii::Table<2, dealii::VectorizedArray>
Implements adamantine::ThermalOperatorBase< dim, MemorySpaceType >.
Definition at line 894 of file ThermalOperator.templates.hh.
|
overridevirtual |
Implements adamantine::ThermalOperatorBase< dim, MemorySpaceType >.
Definition at line 884 of file ThermalOperator.templates.hh.
|
inlineoverridevirtual |
Implements adamantine::ThermalOperatorBase< dim, MemorySpaceType >.
Definition at line 300 of file ThermalOperator.hh.
|
private |
Update the ratios of the material state at the face quadrature points.
Definition at line 321 of file ThermalOperator.templates.hh.
|
private |
Update the ratios of the material state.
Definition at line 229 of file ThermalOperator.templates.hh.
|
overridevirtual |
Matrix-vector multiplication. This function applies the operator to the vector src.
| [in] | src | |
| [out] | dst |
Implements adamantine::ThermalOperatorBase< dim, MemorySpaceType >.
Definition at line 178 of file ThermalOperator.templates.hh.
|
overridevirtual |
Matrix-vector multiplication and addition of the result to dst. This function applies the operator to the vector src and add the result to the vector dst.
| [in] | src | |
| [in,out] | dst |
Implements adamantine::ThermalOperatorBase< dim, MemorySpaceType >.
Definition at line 190 of file ThermalOperator.templates.hh.
|
private |
Flag set to true if all the boundary conditions are adiabatic. It is set to false otherwise.
Definition at line 164 of file ThermalOperator.hh.
|
private |
Non-owning pointer to the AffineConstraints from ThermalPhysics.
Definition at line 197 of file ThermalOperator.hh.
|
private |
Boundary ids associated to the domain.
Definition at line 172 of file ThermalOperator.hh.
|
private |
Map between the cell iterator and the position in _inv_rho_cp table.
Definition at line 210 of file ThermalOperator.hh.
|
private |
MPI communicator.
Definition at line 159 of file ThermalOperator.hh.
|
private |
Current height of the heat sources.
Definition at line 168 of file ThermalOperator.hh.
|
private |
Table of the material deposition cosine angles.
Definition at line 243 of file ThermalOperator.hh.
|
private |
Table of the material deposition cosine angles.
Definition at line 247 of file ThermalOperator.hh.
|
mutableprivate |
Table of the material index on faces; mutable so that it can be changed in face_local_apply which is const.
Definition at line 239 of file ThermalOperator.hh.
|
mutableprivate |
Table of the powder fraction on faces; mutable so that it can be changed in face_local_apply which is const.
Definition at line 225 of file ThermalOperator.hh.
|
private |
Vector of heat sources.
Definition at line 189 of file ThermalOperator.hh.
|
private |
The inverse of the mass matrix is computed using an inexact Gauss-Lobatto quadrature. This inexact quadrature makes the mass matrix and therefore also its inverse, a diagonal matrix.
Definition at line 204 of file ThermalOperator.hh.
|
mutableprivate |
Table of the powder fraction inside cells; mutable so that it can be changed in cell_local_apply which is const.
Definition at line 215 of file ThermalOperator.hh.
|
mutableprivate |
Table of the material index inside cells; mutable so that it can be changed in cell_local_apply which is const.
Definition at line 232 of file ThermalOperator.hh.
|
private |
Material properties associated with the domain.
Definition at line 185 of file ThermalOperator.hh.
|
private |
Underlying MatrixFree object.
Definition at line 193 of file ThermalOperator.hh.
|
private |
Data to configure the MatrixFree object.
Definition at line 176 of file ThermalOperator.hh.
|
mutableprivate |
Table of the powder fraction inside cells; mutable so that it can be changed in cell_local_apply which is const.
Definition at line 220 of file ThermalOperator.hh.
|
private |
Table of thermal conductivity coefficient.
Definition at line 180 of file ThermalOperator.hh.