OGS
ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim > Class Template Referencefinal

Detailed Description

template<int DisplacementDim>
class ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim >

A class to simulate thermo-mechanical fracturing process using phase-field approach in solids described by.

\[ \mathrm{div} \left[ \left(d^2 + k \right) \boldsymbol{\sigma}_0^+ + \boldsymbol{\sigma}_0^- \right] + \varrho \boldsymbol{b} = \boldsymbol{0} \]

\[ 2d \psi^+(\boldsymbol{\epsilon}_\mathrm{el}) - \frac{1 - d}{2 \varepsilon} g_\mathrm{c} - 2 \varepsilon g_\mathrm{c} \mathrm{div}(\mathrm{grad} d) = 0 \]

\[ (\varrho c_\mathrm{p})_\mathrm{eff} \dfrac{\partial \vartheta}{\partial t} - \mathrm{div} \left(\boldsymbol{\kappa}_\mathrm{eff} \mathrm{grad}\,\vartheta \right) = 0 \]

where

\begin{eqnarray*} &d:& \mbox{order parameter,}\\ &\varrho:& \mbox{density,}\\ &g_\mathrm{c}:& \mbox{fracture energy,}\\ &\varepsilon:& \mbox{length scale}\\ &c_\mathrm{p}:& \mbox{specific heat capacity in constant pressure}\\ &\kappa_\mathrm{eff}:& \mbox{effectiev thermal conductivity}\\ \end{eqnarray*}

Detailed model description can refer [22].

Definition at line 53 of file ThermoMechanicalPhaseFieldProcess.h.

#include <ThermoMechanicalPhaseFieldProcess.h>

Inheritance diagram for ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim >:
[legend]
Collaboration diagram for ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim >:
[legend]

Public Member Functions

 ThermoMechanicalPhaseFieldProcess (std::string name, MeshLib::Mesh &mesh, std::unique_ptr< ProcessLib::AbstractJacobianAssembler > &&jacobian_assembler, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, unsigned const integration_order, std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > &&process_variables, ThermoMechanicalPhaseFieldProcessData< DisplacementDim > &&process_data, SecondaryVariableCollection &&secondary_variables, int const mechanics_related_process_id, int const phase_field_process_id, int const heat_conduction_process_id)
 
MathLib::MatrixSpecifications getMatrixSpecifications (const int process_id) const override
 
NumLib::LocalToGlobalIndexMap const & getDOFTable (const int process_id) const override
 
ODESystem interface
bool isLinear () const override
 
- Public Member Functions inherited from ProcessLib::Process
 Process (std::string name_, MeshLib::Mesh &mesh, std::unique_ptr< AbstractJacobianAssembler > &&jacobian_assembler, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, unsigned const integration_order, std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > &&process_variables, SecondaryVariableCollection &&secondary_variables, const bool use_monolithic_scheme=true)
 
void preTimestep (std::vector< GlobalVector * > const &x, const double t, const double delta_t, const int process_id)
 Preprocessing before starting assembly for new timestep.
 
void postTimestep (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, const double delta_t, int const process_id)
 Postprocessing after a complete timestep.
 
void postNonLinearSolver (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, double const dt, int const process_id)
 
void preIteration (const unsigned iter, GlobalVector const &x) final
 
void computeSecondaryVariable (double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_prev, int const process_id)
 compute secondary variables for the coupled equations or for output.
 
NumLib::IterationResult postIteration (GlobalVector const &x) final
 
void initialize (std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
 
void setInitialConditions (std::vector< GlobalVector * > &process_solutions, std::vector< GlobalVector * > const &process_solutions_prev, double const t, int const process_id)
 
MathLib::MatrixSpecifications getMatrixSpecifications (const int process_id) const override
 
void updateDeactivatedSubdomains (double const time, const int process_id)
 
virtual bool isMonolithicSchemeUsed () const
 
virtual void extrapolateIntegrationPointValuesToNodes (const double, std::vector< GlobalVector * > const &, std::vector< GlobalVector * > &)
 
void preAssemble (const double t, double const dt, GlobalVector const &x) final
 
void assemble (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b) final
 
void assembleWithJacobian (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac) final
 
void preOutput (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id)
 
std::vector< NumLib::IndexValueVector< GlobalIndexType > > const * getKnownSolutions (double const t, GlobalVector const &x, int const process_id) const final
 
MeshLib::MeshgetMesh () const
 
std::vector< std::reference_wrapper< ProcessVariable > > const & getProcessVariables (const int process_id) const
 
std::vector< std::size_t > const & getActiveElementIDs () const
 
SecondaryVariableCollection const & getSecondaryVariables () const
 
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > const & getIntegrationPointWriters () const
 
virtual Eigen::Vector3d getFlux (std::size_t, MathLib::Point3d const &, double const, std::vector< GlobalVector * > const &) const
 
virtual void solveReactionEquation (std::vector< GlobalVector * > &, std::vector< GlobalVector * > const &, double const, double const, NumLib::EquationSystem &, int const)
 
- Public Member Functions inherited from ProcessLib::SubmeshAssemblySupport
virtual std::vector< std::string > initializeAssemblyOnSubmeshes (std::vector< std::reference_wrapper< MeshLib::Mesh > > const &meshes)
 
virtual ~SubmeshAssemblySupport ()=default
 

Private Member Functions

void constructDofTable () override
 
void initializeBoundaryConditions (std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media) override
 
void initializeConcreteProcess (NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override
 Process specific initialization called by initialize().
 
void assembleConcreteProcess (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b) override
 
void assembleWithJacobianConcreteProcess (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac) override
 
void preTimestepConcreteProcess (std::vector< GlobalVector * > const &x, double const t, double const dt, const int process_id) override
 
void postTimestepConcreteProcess (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id) override
 
void postNonLinearSolverConcreteProcess (std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, const double t, double const dt, int const process_id) override
 
NumLib::LocalToGlobalIndexMapgetDOFTableByProcessID (const int process_id) const
 

Private Attributes

ThermoMechanicalPhaseFieldProcessData< DisplacementDim > _process_data
 
std::vector< std::unique_ptr< ThermoMechanicalPhaseFieldLocalAssemblerInterface > > _local_assemblers
 
std::unique_ptr< NumLib::LocalToGlobalIndexMap_local_to_global_index_map_single_component
 
GlobalSparsityPattern _sparsity_pattern_with_single_component
 Sparsity pattern for the phase field equation, and it is initialized.
 
int const _mechanics_related_process_id
 ID of the processes that contains mechanical process.
 
int const _phase_field_process_id
 ID of phase field process.
 
int const _heat_conduction_process_id
 ID of heat conduction process.
 

Additional Inherited Members

- Public Attributes inherited from ProcessLib::Process
std::string const name
 
- Static Public Attributes inherited from ProcessLib::Process
static PROCESSLIB_EXPORT const std::string constant_one_parameter_name = "constant_one"
 
- Protected Member Functions inherited from ProcessLib::Process
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables (int const number_of_processes) const
 
NumLib::ExtrapolatorgetExtrapolator () const
 
NumLib::LocalToGlobalIndexMap const & getSingleComponentDOFTable () const
 
void initializeProcessBoundaryConditionsAndSourceTerms (const NumLib::LocalToGlobalIndexMap &dof_table, const int process_id, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
 
void constructMonolithicProcessDofTable ()
 
void constructDofTableOfSpecifiedProcessStaggeredScheme (const int specified_process_id)
 
virtual std::tuple< NumLib::LocalToGlobalIndexMap *, bool > getDOFTableForExtrapolatorData () const
 
std::vector< GlobalIndexTypegetIndicesOfResiduumWithoutInitialCompensation () const override
 
- Protected Attributes inherited from ProcessLib::Process
MeshLib::Mesh_mesh
 
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
 
std::unique_ptr< NumLib::LocalToGlobalIndexMap_local_to_global_index_map
 
SecondaryVariableCollection _secondary_variables
 
std::unique_ptr< ProcessLib::AbstractJacobianAssembler_jacobian_assembler
 
VectorMatrixAssembler _global_assembler
 
const bool _use_monolithic_scheme
 
unsigned const _integration_order
 
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > _integration_point_writer
 
GlobalSparsityPattern _sparsity_pattern
 
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > _process_variables
 
std::vector< BoundaryConditionCollection_boundary_conditions
 

Constructor & Destructor Documentation

◆ ThermoMechanicalPhaseFieldProcess()

template<int DisplacementDim>
ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim >::ThermoMechanicalPhaseFieldProcess ( std::string name,
MeshLib::Mesh & mesh,
std::unique_ptr< ProcessLib::AbstractJacobianAssembler > && jacobian_assembler,
std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const & parameters,
unsigned const integration_order,
std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > && process_variables,
ThermoMechanicalPhaseFieldProcessData< DisplacementDim > && process_data,
SecondaryVariableCollection && secondary_variables,
int const mechanics_related_process_id,
int const phase_field_process_id,
int const heat_conduction_process_id )

Definition at line 26 of file ThermoMechanicalPhaseFieldProcess.cpp.

42 : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
43 integration_order, std::move(process_variables),
44 std::move(secondary_variables), false),
45 _process_data(std::move(process_data)),
46 _mechanics_related_process_id(mechanics_related_process_id),
47 _phase_field_process_id(phase_field_process_id),
48 _heat_conduction_process_id(heat_conduction_process_id)
49{
50}
std::string const name
Definition Process.h:354
Process(std::string name_, MeshLib::Mesh &mesh, std::unique_ptr< AbstractJacobianAssembler > &&jacobian_assembler, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, unsigned const integration_order, std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > &&process_variables, SecondaryVariableCollection &&secondary_variables, const bool use_monolithic_scheme=true)
Definition Process.cpp:44
int const _mechanics_related_process_id
ID of the processes that contains mechanical process.

Member Function Documentation

◆ assembleConcreteProcess()

template<int DisplacementDim>
void ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim >::assembleConcreteProcess ( const double t,
double const dt,
std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
int const process_id,
GlobalMatrix & M,
GlobalMatrix & K,
GlobalVector & b )
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 192 of file ThermoMechanicalPhaseFieldProcess.cpp.

198{
199 DBUG("Assemble the equations for ThermoMechanicalPhaseFieldProcess.");
200
201 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_table = {
203 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
204
205 // Call global assembler for each local assembly item.
208 pv.getActiveElementIDs(), dof_table, t, dt, x, x_prev, process_id, M, K,
209 b);
210}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
std::vector< std::size_t > const & getActiveElementIDs() const
std::vector< std::reference_wrapper< ProcessVariable > > const & getProcessVariables(const int process_id) const
Definition Process.h:155
VectorMatrixAssembler _global_assembler
Definition Process.h:367
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition Process.h:360
std::vector< std::unique_ptr< ThermoMechanicalPhaseFieldLocalAssemblerInterface > > _local_assemblers
void assemble(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b)
static void executeSelectedMemberDereferenced(Object &object, Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)

References ProcessLib::VectorMatrixAssembler::assemble(), DBUG(), NumLib::SerialExecutor::executeSelectedMemberDereferenced(), and ProcessLib::ProcessVariable::getActiveElementIDs().

◆ assembleWithJacobianConcreteProcess()

template<int DisplacementDim>
void ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim >::assembleWithJacobianConcreteProcess ( const double t,
double const dt,
std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
int const process_id,
GlobalMatrix & M,
GlobalMatrix & K,
GlobalVector & b,
GlobalMatrix & Jac )
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 213 of file ThermoMechanicalPhaseFieldProcess.cpp.

218{
219 // For the staggered scheme
220 if (process_id == _mechanics_related_process_id)
221 {
222 DBUG(
223 "Assemble the Jacobian equations of "
224 "temperature-deformation in ThermoMechanicalPhaseFieldProcess for "
225 "the staggered scheme.");
226 }
227
228 if (process_id == _phase_field_process_id)
229 {
230 DBUG(
231 "Assemble the Jacobian equations ofphase field in "
232 "ThermoMechanicalPhaseFieldProcess for the staggered scheme.");
233 }
234 else
235 {
236 DBUG(
237 "Assemble the Jacobian equations of heat conduction in "
238 "ThermoMechanicalPhaseFieldProcess for the staggered scheme.");
239 }
240
241 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
242 dof_tables.emplace_back(
244 dof_tables.emplace_back(
246 dof_tables.emplace_back(&getDOFTableByProcessID(_phase_field_process_id));
247
248 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
249
252 _local_assemblers, pv.getActiveElementIDs(), dof_tables, t, dt, x,
253 x_prev, process_id, M, K, b, Jac);
254}
NumLib::LocalToGlobalIndexMap & getDOFTableByProcessID(const int process_id) const
void assembleWithJacobian(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac)

References ProcessLib::VectorMatrixAssembler::assembleWithJacobian(), DBUG(), NumLib::SerialExecutor::executeSelectedMemberDereferenced(), and ProcessLib::ProcessVariable::getActiveElementIDs().

◆ constructDofTable()

template<int DisplacementDim>
void ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim >::constructDofTable ( )
overrideprivatevirtual

This function is for general cases, in which all equations of the coupled processes have the same number of unknowns. For the general cases with the staggered scheme, all equations of the coupled processes share one DOF table hold by _local_to_global_index_map. Other cases can be considered by overloading this member function in the derived class.

Reimplemented from ProcessLib::Process.

Definition at line 105 of file ThermoMechanicalPhaseFieldProcess.cpp.

106{
107 // For displacement equation.
110
111 // TODO move the two data members somewhere else.
112 // for extrapolation of secondary variables of stress or strain
113 std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component{
116 std::make_unique<NumLib::LocalToGlobalIndexMap>(
117 std::move(all_mesh_subsets_single_component),
118 // by location order is needed for output
120
122
123 // For phase field equation or the heat conduction.
126}
void constructDofTableOfSpecifiedProcessStaggeredScheme(const int specified_process_id)
Definition Process.cpp:355
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
Definition Process.h:358
MeshLib::Mesh & _mesh
Definition Process.h:357
GlobalSparsityPattern _sparsity_pattern_with_single_component
Sparsity pattern for the phase field equation, and it is initialized.
@ BY_LOCATION
Ordering data by spatial location.
GlobalSparsityPattern computeSparsityPattern(LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh)
Computes a sparsity pattern for the given inputs.

References NumLib::BY_LOCATION, and NumLib::computeSparsityPattern().

◆ getDOFTable()

template<int DisplacementDim>
NumLib::LocalToGlobalIndexMap const & ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim >::getDOFTable ( const int process_id) const
overridevirtual

Reimplemented from ProcessLib::Process.

Definition at line 78 of file ThermoMechanicalPhaseFieldProcess.cpp.

80{
81 if (process_id == _mechanics_related_process_id)
82 {
84 }
85
86 // For the equation of phasefield or heat conduction.
88}

◆ getDOFTableByProcessID()

template<int DisplacementDim>
NumLib::LocalToGlobalIndexMap & ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim >::getDOFTableByProcessID ( const int process_id) const
private

Definition at line 92 of file ThermoMechanicalPhaseFieldProcess.cpp.

94{
95 if (process_id == _mechanics_related_process_id)
96 {
98 }
99
100 // For the equation of phasefield or heat conduction.
102}

◆ getMatrixSpecifications()

template<int DisplacementDim>
MathLib::MatrixSpecifications ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim >::getMatrixSpecifications ( const int process_id) const
override

Definition at line 60 of file ThermoMechanicalPhaseFieldProcess.cpp.

62{
63 if (process_id == _mechanics_related_process_id)
64 {
65 auto const& l = *_local_to_global_index_map;
66 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
67 &l.getGhostIndices(), &this->_sparsity_pattern};
68 }
69
70 // For staggered scheme and phase field process or heat conduction.
72 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
73 &l.getGhostIndices(), &_sparsity_pattern_with_single_component};
74}
GlobalSparsityPattern _sparsity_pattern
Definition Process.h:382

◆ initializeBoundaryConditions()

template<int DisplacementDim>
void ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim >::initializeBoundaryConditions ( std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const & media)
overrideprivatevirtual

Member function to initialize the boundary conditions for all coupled processes. It is called by initialize().

Reimplemented from ProcessLib::Process.

Definition at line 171 of file ThermoMechanicalPhaseFieldProcess.cpp.

175{
176 // Staggered scheme:
177 // for the equations of temperature-deformation.
181 // for the phase field
185 // for heat conduction
189}
void initializeProcessBoundaryConditionsAndSourceTerms(const NumLib::LocalToGlobalIndexMap &dof_table, const int process_id, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
Definition Process.cpp:89

◆ initializeConcreteProcess()

template<int DisplacementDim>
void ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim >::initializeConcreteProcess ( NumLib::LocalToGlobalIndexMap const & dof_table,
MeshLib::Mesh const & mesh,
unsigned const integration_order )
overrideprivatevirtual

Process specific initialization called by initialize().

Implements ProcessLib::Process.

Definition at line 129 of file ThermoMechanicalPhaseFieldProcess.cpp.

133{
135 DisplacementDim, ThermoMechanicalPhaseFieldLocalAssembler>(
136 mesh.getElements(), dof_table, _local_assemblers,
137 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
140
142 "sigma",
145 DisplacementDim>::RowsAtCompileTime,
148
150 "epsilon",
152 DisplacementDim>::RowsAtCompileTime,
154 &ThermoMechanicalPhaseFieldLocalAssemblerInterface::
155 getIntPtEpsilon));
156
158 "heat_flux",
159 makeExtrapolator(mesh.getDimension(), getExtrapolator(),
161 &ThermoMechanicalPhaseFieldLocalAssemblerInterface::
162 getIntPtHeatFlux));
163
164 // Initialize local assemblers after all variables have been set.
168}
virtual void initialize(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
SecondaryVariableCollection _secondary_variables
Definition Process.h:362
NumLib::Extrapolator & getExtrapolator() const
Definition Process.h:199
void addSecondaryVariable(std::string const &internal_name, SecondaryVariableFunctions &&fcts)
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
void createLocalAssemblers(std::vector< MeshLib::Element * > const &mesh_elements, NumLib::LocalToGlobalIndexMap const &dof_table, std::vector< std::unique_ptr< LocalAssemblerInterface > > &local_assemblers, ProviderOrOrder const &provider_or_order, ExtraCtorArgs &&... extra_ctor_args)
SecondaryVariableFunctions makeExtrapolator(const unsigned num_components, NumLib::Extrapolator &extrapolator, LocalAssemblerCollection const &local_assemblers, typename NumLib::ExtrapolatableLocalAssemblerCollection< LocalAssemblerCollection >::IntegrationPointValuesMethod integration_point_values_method)
static void executeMemberOnDereferenced(Method method, Container const &container, Args &&... args)
virtual std::vector< double > const & getIntPtSigma(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0

References ProcessLib::SmallDeformation::createLocalAssemblers(), NumLib::SerialExecutor::executeMemberOnDereferenced(), MeshLib::Mesh::getDimension(), MeshLib::Mesh::getElements(), ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssemblerInterface::getIntPtHeatFlux(), ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldLocalAssemblerInterface::getIntPtSigma(), ProcessLib::LocalAssemblerInterface::initialize(), MeshLib::Mesh::isAxiallySymmetric(), and ProcessLib::makeExtrapolator().

◆ isLinear()

template<int DisplacementDim>
bool ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim >::isLinear ( ) const
override

Definition at line 53 of file ThermoMechanicalPhaseFieldProcess.cpp.

54{
55 return false;
56}

◆ postNonLinearSolverConcreteProcess()

template<int DisplacementDim>
void ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim >::postNonLinearSolverConcreteProcess ( std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
const double t,
double const dt,
int const process_id )
overrideprivatevirtual

Reimplemented from ProcessLib::Process.

Definition at line 301 of file ThermoMechanicalPhaseFieldProcess.cpp.

306{
307 if (process_id != _mechanics_related_process_id)
308 {
309 return;
310 }
311
312 DBUG("PostNonLinearSolver ThermoMechanicalPhaseFieldProcess.");
313 // Calculate strain, stress or other internal variables of mechanics.
314 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
315
318 pv.getActiveElementIDs(), getDOFTables(x.size()), x, x_prev, t, dt,
319 process_id);
320}
void postNonLinearSolver(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables(int const number_of_processes) const
Definition Process.cpp:386
static void executeSelectedMemberOnDereferenced(Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)

References DBUG(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::ProcessVariable::getActiveElementIDs(), and ProcessLib::LocalAssemblerInterface::postNonLinearSolver().

◆ postTimestepConcreteProcess()

template<int DisplacementDim>
void ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim >::postTimestepConcreteProcess ( std::vector< GlobalVector * > const & x,
std::vector< GlobalVector * > const & x_prev,
double const t,
double const dt,
int const process_id )
overrideprivatevirtual

Reimplemented from ProcessLib::Process.

Definition at line 279 of file ThermoMechanicalPhaseFieldProcess.cpp.

285{
286 if (process_id != 0)
287 {
288 return;
289 }
290
291 DBUG("PostTimestep ThermoMechanicalPhaseFieldProcess.");
292
293 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
297 x_prev, t, dt, process_id);
298}
virtual void postTimestep(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id)

References DBUG(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::ProcessVariable::getActiveElementIDs(), and ProcessLib::LocalAssemblerInterface::postTimestep().

◆ preTimestepConcreteProcess()

template<int DisplacementDim>
void ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim >::preTimestepConcreteProcess ( std::vector< GlobalVector * > const & x,
double const t,
double const dt,
const int process_id )
overrideprivatevirtual

Reimplemented from ProcessLib::Process.

Definition at line 257 of file ThermoMechanicalPhaseFieldProcess.cpp.

262{
263 DBUG("PreTimestep ThermoMechanicalPhaseFieldProcess.");
264
265 if (process_id != _mechanics_related_process_id)
266 {
267 return;
268 }
269
270 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
271
275 *x[process_id], t, dt);
276}
virtual void preTimestep(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, double const delta_t)
NumLib::LocalToGlobalIndexMap const & getDOFTable(const int process_id) const override

References DBUG(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::ProcessVariable::getActiveElementIDs(), and ProcessLib::LocalAssemblerInterface::preTimestep().

Member Data Documentation

◆ _heat_conduction_process_id

template<int DisplacementDim>
int const ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim >::_heat_conduction_process_id
private

ID of heat conduction process.

Definition at line 146 of file ThermoMechanicalPhaseFieldProcess.h.

◆ _local_assemblers

template<int DisplacementDim>
std::vector< std::unique_ptr<ThermoMechanicalPhaseFieldLocalAssemblerInterface> > ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim >::_local_assemblers
private

Definition at line 130 of file ThermoMechanicalPhaseFieldProcess.h.

◆ _local_to_global_index_map_single_component

template<int DisplacementDim>
std::unique_ptr<NumLib::LocalToGlobalIndexMap> ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim >::_local_to_global_index_map_single_component
private

Definition at line 133 of file ThermoMechanicalPhaseFieldProcess.h.

◆ _mechanics_related_process_id

template<int DisplacementDim>
int const ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim >::_mechanics_related_process_id
private

ID of the processes that contains mechanical process.

Definition at line 140 of file ThermoMechanicalPhaseFieldProcess.h.

◆ _phase_field_process_id

template<int DisplacementDim>
int const ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim >::_phase_field_process_id
private

ID of phase field process.

Definition at line 143 of file ThermoMechanicalPhaseFieldProcess.h.

◆ _process_data

template<int DisplacementDim>
ThermoMechanicalPhaseFieldProcessData<DisplacementDim> ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim >::_process_data
private

Definition at line 126 of file ThermoMechanicalPhaseFieldProcess.h.

◆ _sparsity_pattern_with_single_component

template<int DisplacementDim>
GlobalSparsityPattern ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim >::_sparsity_pattern_with_single_component
private

Sparsity pattern for the phase field equation, and it is initialized.

Definition at line 137 of file ThermoMechanicalPhaseFieldProcess.h.


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