adamantine
Public Member Functions | Static Public Member Functions | Static Public Attributes | Private Member Functions | Private Attributes | List of all members
adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType > Class Template Reference

#include <MaterialProperty.hh>

Public Member Functions

 MaterialProperty (MPI_Comm const &communicator, dealii::parallel::distributed::Triangulation< dim > const &tria, boost::property_tree::ptree const &database)
 
bool properties_use_table () const
 
double get_cell_value (typename dealii::Triangulation< dim >::active_cell_iterator const &cell, StateProperty prop) const
 
double get_cell_value (typename dealii::Triangulation< dim >::active_cell_iterator const &cell, Property prop) const
 
double get (dealii::types::material_id material_id, Property prop) const
 
double get_mechanical_property (typename dealii::Triangulation< dim >::active_cell_iterator const &cell, StateProperty prop) const
 
Kokkos::View< typename internal::Data2D< n_materials, g_n_properties >::type, typename MemorySpaceType::kokkos_space > get_properties ()
 
Kokkos::View< typename internal::Data5D< n_materials, g_n_thermal_state_properties, MaterialStates::n_material_states, table_size, 2 >::type, typename MemorySpaceType::kokkos_space > get_state_property_tables ()
 
Kokkos::View< typename internal::Data4D< n_materials, g_n_thermal_state_properties, MaterialStates::n_material_states, p_order+1 >::type, typename MemorySpaceType::kokkos_space > get_state_property_polynomials ()
 
void reinit_dofs ()
 
void update (dealii::DoFHandler< dim > const &temperature_dof_handler, dealii::LA::distributed::Vector< double, MemorySpaceType > const &temperature)
 
void update_boundary_material_properties (dealii::DoFHandler< dim > const &temperature_dof_handler, dealii::LA::distributed::Vector< double, MemorySpaceType > const &temperature)
 
template<bool use_table, StateProperty state_property>
dealii::VectorizedArray< double > compute_material_property (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
 
template<bool use_table>
KOKKOS_FUNCTION double compute_material_property (StateProperty state_property, dealii::types::material_id const material_id, double const *state_ratios, double temperature) const
 
Kokkos::View< double **, typename MemorySpaceType::kokkos_space > get_state () const
 
double get_state_ratio (typename dealii::Triangulation< dim >::active_cell_iterator const &cell, typename MaterialStates::State material_state) const
 
void set_initial_state ()
 
void set_state (dealii::Table< 2, dealii::VectorizedArray< double >> const &liquid_ratio, dealii::Table< 2, dealii::VectorizedArray< double >> const &powder_ratio, std::map< typename dealii::DoFHandler< dim >::cell_iterator, std::pair< unsigned int, unsigned int >> &cell_it_to_mf_cell_map, dealii::DoFHandler< dim > const &dof_handler)
 
void set_state_device (Kokkos::View< double *, typename MemorySpaceType::kokkos_space > liquid_ratio, Kokkos::View< double *, typename MemorySpaceType::kokkos_space > powder_ratio, std::map< typename dealii::DoFHandler< dim >::cell_iterator, std::vector< unsigned int >> const &_cell_it_to_mf_pos, dealii::DoFHandler< dim > const &dof_handler)
 
void set_cell_state (std::vector< std::array< double, MaterialStates::n_material_states >> const &cell_state)
 
dealii::DoFHandler< dim > const & get_dof_handler () const
 
std::unordered_map< dealii::types::global_dof_index, unsigned int > get_dofs_map () const
 

Static Public Member Functions

static KOKKOS_FUNCTION double compute_property_from_table (Kokkos::View< typename internal::Data5D< n_materials, g_n_thermal_state_properties, MaterialStates::n_material_states, table_size, 2 >::type, typename MemorySpaceType::kokkos_space > state_property_tables, unsigned int const material_id, unsigned int const material_state, unsigned int const property, double const temperature)
 

Static Public Attributes

static unsigned constexpr int table_size = 12
 

Private Member Functions

void fill_properties (boost::property_tree::ptree const &database)
 
dealii::types::global_dof_index get_dof_index (typename dealii::Triangulation< dim >::active_cell_iterator const &cell) const
 
dealii::LA::distributed::Vector< double, MemorySpaceType > compute_average_temperature (dealii::DoFHandler< dim > const &temperature_dof_handler, dealii::LA::distributed::Vector< double, MemorySpaceType > const &temperature) const
 

Private Attributes

MPI_Comm _communicator
 
bool _use_table
 
Kokkos::View< typename internal::Data5D< n_materials, g_n_thermal_state_properties, MaterialStates::n_material_states, table_size, 2 >::type, typename MemorySpaceType::kokkos_space > _state_property_tables
 
Kokkos::View< typename internal::Data4D< n_materials, g_n_thermal_state_properties, MaterialStates::n_material_states, p_order+1 >::type, typename MemorySpaceType::kokkos_space > _state_property_polynomials
 
Kokkos::View< typename internal::Data2D< n_materials, g_n_properties >::type, typename MemorySpaceType::kokkos_space > _properties
 
Kokkos::View< double **, typename MemorySpaceType::kokkos_space > _state
 
Kokkos::View< double **, typename MemorySpaceType::kokkos_space > _property_values
 
Kokkos::View< typename internal::Data4D< n_materials, g_n_mechanical_state_properties, table_size, 2 >::type, dealii::MemorySpace::Host::kokkos_space > _mechanical_properties_tables_host
 
Kokkos::View< typename internal::Data3D< n_materials, g_n_mechanical_state_properties, p_order+1 >::type, dealii::MemorySpace::Host::kokkos_space > _mechanical_properties_polynomials_host
 
Kokkos::View< typename internal::Data2D< n_materials, g_n_mechanical_state_properties >::type, dealii::MemorySpace::Host::kokkos_space > _mechanical_properties_host
 
dealii::FE_DGQ< dim > _fe
 
dealii::DoFHandler< dim > _mp_dof_handler
 
std::unordered_map< dealii::types::global_dof_index, unsigned int > _dofs_map
 

Detailed Description

template<int dim, int n_materials, int p_order, typename MaterialStates, typename MemorySpaceType>
class adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >

This class stores the material properties for all the materials

Definition at line 126 of file MaterialProperty.hh.

Constructor & Destructor Documentation

◆ MaterialProperty()

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::MaterialProperty ( MPI_Comm const &  communicator,
dealii::parallel::distributed::Triangulation< dim > const &  tria,
boost::property_tree::ptree const &  database 
)

Constructor.

Definition at line 126 of file MaterialProperty.templates.hh.

Member Function Documentation

◆ compute_average_temperature()

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
dealii::LA::distributed::Vector< double, MemorySpaceType > adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::compute_average_temperature ( dealii::DoFHandler< dim > const &  temperature_dof_handler,
dealii::LA::distributed::Vector< double, MemorySpaceType > const &  temperature 
) const
private

Compute the average of the temperature on every cell.

Definition at line 1113 of file MaterialProperty.templates.hh.

◆ compute_material_property() [1/2]

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
template<bool use_table>
KOKKOS_FUNCTION double adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::compute_material_property ( StateProperty  state_property,
dealii::types::material_id const  material_id,
double const *  state_ratios,
double  temperature 
) const

Compute a material property at a quadrature point for a mix of states.

Definition at line 806 of file MaterialProperty.hh.

◆ compute_material_property() [2/2]

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
template<bool use_table, StateProperty state_property>
dealii::VectorizedArray< double > adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::compute_material_property ( 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

Compute a material property at a quadrature point for a mix of states.

Definition at line 533 of file MaterialProperty.hh.

◆ compute_property_from_table()

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
KOKKOS_FUNCTION double adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::compute_property_from_table ( Kokkos::View< typename internal::Data5D< n_materials, g_n_thermal_state_properties, MaterialStates::n_material_states, table_size, 2 >::type, typename MemorySpaceType::kokkos_space >  state_property_tables,
unsigned int const  material_id,
unsigned int const  material_state,
unsigned int const  property,
double const  temperature 
)
static

Compute a property from a table given the temperature.

Definition at line 1149 of file MaterialProperty.templates.hh.

◆ fill_properties()

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
void adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::fill_properties ( boost::property_tree::ptree const &  database)
private

Fill the _properties map.

Definition at line 785 of file MaterialProperty.templates.hh.

◆ get()

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
double adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::get ( dealii::types::material_id  material_id,
Property  prop 
) const
inline

Return the value of a given Property for a given material id.

Definition at line 440 of file MaterialProperty.hh.

◆ get_cell_value() [1/2]

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
double adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::get_cell_value ( typename dealii::Triangulation< dim >::active_cell_iterator const &  cell,
Property  prop 
) const

Return the value of the given Property for a given cell.

Definition at line 165 of file MaterialProperty.templates.hh.

◆ get_cell_value() [2/2]

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
double adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::get_cell_value ( typename dealii::Triangulation< dim >::active_cell_iterator const &  cell,
StateProperty  prop 
) const

Return the value of the given StateProperty for a given cell.

Definition at line 149 of file MaterialProperty.templates.hh.

◆ get_dof_handler()

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
dealii::DoFHandler< dim > const & adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::get_dof_handler
inline

Return the underlying the DoFHandler.

Definition at line 524 of file MaterialProperty.hh.

◆ get_dof_index()

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
dealii::types::global_dof_index adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::get_dof_index ( typename dealii::Triangulation< dim >::active_cell_iterator const &  cell) const
inlineprivate

Return the index of the dof associated to the cell.

Definition at line 505 of file MaterialProperty.hh.

◆ get_dofs_map()

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
std::unordered_map<dealii::types::global_dof_index, unsigned int> adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::get_dofs_map ( ) const
inline

Return the mapping between the degrees of freedom and the local index of the cells.

Definition at line 312 of file MaterialProperty.hh.

◆ get_mechanical_property()

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
double adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::get_mechanical_property ( typename dealii::Triangulation< dim >::active_cell_iterator const &  cell,
StateProperty  prop 
) const

Return the values of the given mechanical StateProperty for a given cell.

Definition at line 181 of file MaterialProperty.templates.hh.

◆ get_properties()

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
Kokkos::View< typename internal::Data2D< n_materials, g_n_properties >::type, typename MemorySpaceType::kokkos_space > adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::get_properties
inline

Return the properties of the material that are independent of the state of the material.

Definition at line 452 of file MaterialProperty.hh.

◆ get_state()

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
Kokkos::View< double **, typename MemorySpaceType::kokkos_space > adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::get_state
inline

Get the array of material state vectors. The order of the different state vectors is given by the MaterialState enum. Each entry in the vector correspond to a cell in the mesh and has a value between 0 and 1. The sum of the states for a given cell is equal to 1.

Definition at line 497 of file MaterialProperty.hh.

◆ get_state_property_polynomials()

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
Kokkos::View< typename internal::Data4D< n_materials, g_n_thermal_state_properties, MaterialStates::n_material_states, p_order+1 >::type, typename MemorySpaceType::kokkos_space > adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::get_state_property_polynomials
inline

Return the properties of the material that are dependent of the state of the material and which have beese set using polynomials.

Definition at line 488 of file MaterialProperty.hh.

◆ get_state_property_tables()

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
Kokkos::View< typename internal::Data5D< n_materials, g_n_thermal_state_properties, MaterialStates::n_material_states, MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::table_size, 2 >::type, typename MemorySpaceType::kokkos_space > adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::get_state_property_tables
inline

Return the properties of the material that are dependent of the state of the material and which have been set using tables.

Definition at line 476 of file MaterialProperty.hh.

◆ get_state_ratio()

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
double adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::get_state_ratio ( typename dealii::Triangulation< dim >::active_cell_iterator const &  cell,
typename MaterialStates::State  material_state 
) const

Get the ratio of a given MaterialState for a given cell. The sum of the states for a given cell is equal to 1.

Definition at line 196 of file MaterialProperty.templates.hh.

◆ properties_use_table()

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
bool adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::properties_use_table
inline

Return true if the material properties are given in table format. Return false if they are given in polynomial format.

Definition at line 460 of file MaterialProperty.hh.

◆ reinit_dofs()

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
void adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::reinit_dofs

Reinitialize the DoFHandler associated with MaterialProperty and resize the state vectors.

Definition at line 212 of file MaterialProperty.templates.hh.

◆ set_cell_state()

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
void adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::set_cell_state ( std::vector< std::array< double, MaterialStates::n_material_states >> const &  cell_state)

Set the ratio of the material states at the cell level.

Definition at line 724 of file MaterialProperty.templates.hh.

◆ set_initial_state()

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
void adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::set_initial_state

Set the values in _state from the values of the user index of the Triangulation.

Definition at line 744 of file MaterialProperty.templates.hh.

◆ set_state()

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
void adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::set_state ( dealii::Table< 2, dealii::VectorizedArray< double >> const &  liquid_ratio,
dealii::Table< 2, dealii::VectorizedArray< double >> const &  powder_ratio,
std::map< typename dealii::DoFHandler< dim >::cell_iterator, std::pair< unsigned int, unsigned int >> &  cell_it_to_mf_cell_map,
dealii::DoFHandler< dim > const &  dof_handler 
)

Set the ratio of the material states from ThermalOperator.

Definition at line 526 of file MaterialProperty.templates.hh.

◆ set_state_device()

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
void adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::set_state_device ( Kokkos::View< double *, typename MemorySpaceType::kokkos_space >  liquid_ratio,
Kokkos::View< double *, typename MemorySpaceType::kokkos_space >  powder_ratio,
std::map< typename dealii::DoFHandler< dim >::cell_iterator, std::vector< unsigned int >> const &  _cell_it_to_mf_pos,
dealii::DoFHandler< dim > const &  dof_handler 
)

Set the ratio of the material states from ThermalOperatorDevice.

Definition at line 612 of file MaterialProperty.templates.hh.

◆ update()

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
void adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::update ( dealii::DoFHandler< dim > const &  temperature_dof_handler,
dealii::LA::distributed::Vector< double, MemorySpaceType > const &  temperature 
)

Update the material state, i.e, the ratio of liquid, powder, and solid and the material properties given the field of temperature.

Definition at line 242 of file MaterialProperty.templates.hh.

◆ update_boundary_material_properties()

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
void adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::update_boundary_material_properties ( dealii::DoFHandler< dim > const &  temperature_dof_handler,
dealii::LA::distributed::Vector< double, MemorySpaceType > const &  temperature 
)

Update the material properties necessary to compute the radiative and convective boundary conditions given the field of temperature.

Definition at line 435 of file MaterialProperty.templates.hh.

Member Data Documentation

◆ _communicator

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
MPI_Comm adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::_communicator
private

MPI communicator.

Definition at line 354 of file MaterialProperty.hh.

◆ _dofs_map

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
std::unordered_map<dealii::types::global_dof_index, unsigned int> adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::_dofs_map
private

Mapping between the degrees of freedom and the local index of the cells.

Definition at line 433 of file MaterialProperty.hh.

◆ _fe

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
dealii::FE_DGQ<dim> adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::_fe
private

Discontinuous piecewise constant finite element.

Definition at line 425 of file MaterialProperty.hh.

◆ _mechanical_properties_host

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
Kokkos::View<typename internal::Data2D<n_materials, g_n_mechanical_state_properties>::type, dealii::MemorySpace::Host::kokkos_space> adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::_mechanical_properties_host
private

Temperature independent mechanical properties.

Definition at line 421 of file MaterialProperty.hh.

◆ _mechanical_properties_polynomials_host

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
Kokkos::View< typename internal::Data3D<n_materials, g_n_mechanical_state_properties, p_order + 1>::type, dealii::MemorySpace::Host::kokkos_space> adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::_mechanical_properties_polynomials_host
private

Mechanical properties which have been set using polynomials.

Definition at line 414 of file MaterialProperty.hh.

◆ _mechanical_properties_tables_host

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
Kokkos::View< typename internal::Data4D<n_materials, g_n_mechanical_state_properties, table_size, 2>::type, dealii::MemorySpace::Host::kokkos_space> adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::_mechanical_properties_tables_host
private

Mechanical properties which have been set using tables.

Definition at line 406 of file MaterialProperty.hh.

◆ _mp_dof_handler

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
dealii::DoFHandler<dim> adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::_mp_dof_handler
private

DoFHandler associated to the _state array.

Definition at line 429 of file MaterialProperty.hh.

◆ _properties

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
Kokkos::View<typename internal::Data2D<n_materials, g_n_properties>::type, typename MemorySpaceType::kokkos_space> adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::_properties
private

Properties of the material that are independent of the state of the material.

Definition at line 383 of file MaterialProperty.hh.

◆ _property_values

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
Kokkos::View<double **, typename MemorySpaceType::kokkos_space> adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::_property_values
private

Thermal properties of the material that are dependent of the state of the material.

Definition at line 395 of file MaterialProperty.hh.

◆ _state

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
Kokkos::View<double **, typename MemorySpaceType::kokkos_space> adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::_state
private

Ratio of each in MaterarialState in each cell.

Definition at line 389 of file MaterialProperty.hh.

◆ _state_property_polynomials

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
Kokkos::View<typename internal::Data4D< n_materials, g_n_thermal_state_properties, MaterialStates::n_material_states, p_order + 1>::type, typename MemorySpaceType::kokkos_space> adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::_state_property_polynomials
private

Thermal material properties which have been set using polynomials.

Definition at line 376 of file MaterialProperty.hh.

◆ _state_property_tables

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
Kokkos::View<typename internal::Data5D< n_materials, g_n_thermal_state_properties, MaterialStates::n_material_states, table_size, 2>::type, typename MemorySpaceType::kokkos_space> adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::_state_property_tables
private

Thermal material properties which have been set using tables.

Definition at line 367 of file MaterialProperty.hh.

◆ _use_table

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
bool adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::_use_table
private

If the flag is true the material properties are saved under a table. Otherwise the material properties are saved as polynomials.

Definition at line 359 of file MaterialProperty.hh.

◆ table_size

template<int dim, int n_materials, int p_order, typename MaterialStates , typename MemorySpaceType >
unsigned constexpr int adamantine::MaterialProperty< dim, n_materials, p_order, MaterialStates, MemorySpaceType >::table_size = 12
staticconstexpr

Size of the table, i.e. number of temperature/property pairs, used to describe the material properties.

Definition at line 133 of file MaterialProperty.hh.