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 [16].

Definition at line 54 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. More...
 
void postTimestep (std::vector< GlobalVector * > const &x, const double t, const double delta_t, int const process_id)
 Postprocessing after a complete timestep. More...
 
void postNonLinearSolver (GlobalVector const &x, GlobalVector const &xdot, 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_dot, int const process_id)
 compute secondary variables for the coupled equations or for output. More...
 
NumLib::IterationResult postIteration (GlobalVector const &x) final
 
void initialize ()
 
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 setCoupledSolutionsForStaggeredScheme (CoupledSolutionsForStaggeredScheme *const coupled_solutions)
 
void updateDeactivatedSubdomains (double const time, const int process_id)
 
bool isMonolithicSchemeUsed () const
 
virtual void setCoupledTermForTheStaggeredSchemeToLocalAssemblers (int 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 &xdot, 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 &xdot, const double dxdot_dx, const double dx_dx, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac) final
 
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
 
SecondaryVariableCollection const & getSecondaryVariables () const
 
std::vector< std::unique_ptr< IntegrationPointWriter > > const * getIntegrationPointWriter (MeshLib::Mesh const &mesh) 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)
 

Private Member Functions

void constructDofTable () override
 
void initializeBoundaryConditions () override
 
void initializeConcreteProcess (NumLib::LocalToGlobalIndexMap const &dof_table, MeshLib::Mesh const &mesh, unsigned const integration_order) override
 Process specific initialization called by initialize(). More...
 
void assembleConcreteProcess (const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, 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 &xdot, const double dxdot_dx, const double dx_dx, 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, double const t, double const dt, int const process_id) override
 
void postNonLinearSolverConcreteProcess (GlobalVector const &x, GlobalVector const &xdot, 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. More...
 
int const _mechanics_related_process_id
 ID of the processes that contains mechanical process. More...
 
int const _phase_field_process_id
 ID of phase field process. More...
 
int const _heat_conduction_process_id
 ID of heat conduction process. More...
 

Additional Inherited Members

- Public Types inherited from ProcessLib::Process
using NonlinearSolver = NumLib::NonlinearSolverBase
 
using TimeDiscretization = NumLib::TimeDiscretization
 
- Public Attributes inherited from ProcessLib::Process
std::string const name
 
- Protected Member Functions inherited from ProcessLib::Process
NumLib::ExtrapolatorgetExtrapolator () const
 
NumLib::LocalToGlobalIndexMap const & getSingleComponentDOFTable () const
 
void initializeProcessBoundaryConditionsAndSourceTerms (const NumLib::LocalToGlobalIndexMap &dof_table, const int process_id)
 
void constructMonolithicProcessDofTable ()
 
void constructDofTableOfSpecifiedProcessStaggeredScheme (const int specified_prosess_id)
 
virtual std::tuple< NumLib::LocalToGlobalIndexMap *, bool > getDOFTableForExtrapolatorData () const
 
- 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
 
VectorMatrixAssembler _global_assembler
 
const bool _use_monolithic_scheme
 
CoupledSolutionsForStaggeredScheme_coupled_solutions
 
unsigned const _integration_order
 
std::vector< std::unique_ptr< 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:323
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:22
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 &  xdot,
int const  process_id,
GlobalMatrix M,
GlobalMatrix K,
GlobalVector b 
)
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 190 of file ThermoMechanicalPhaseFieldProcess.cpp.

196 {
197  DBUG("Assemble the equations for ThermoMechanicalPhaseFieldProcess.");
198 
199  std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
200  dof_table = {std::ref(*_local_to_global_index_map)};
201  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
202 
203  // Call global assembler for each local assembly item.
206  pv.getActiveElementIDs(), dof_table, t, dt, x, xdot, process_id, M, K,
207  b);
208 }
void DBUG(char const *fmt, Args const &... args)
Definition: Logging.h:27
std::vector< std::size_t > const & getActiveElementIDs() const
std::vector< std::reference_wrapper< ProcessVariable > > const & getProcessVariables(const int process_id) const
Definition: Process.h:145
VectorMatrixAssembler _global_assembler
Definition: Process.h:333
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition: Process.h:329
std::vector< std::unique_ptr< ThermoMechanicalPhaseFieldLocalAssemblerInterface > > _local_assemblers
void assemble(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< std::reference_wrapper< NumLib::LocalToGlobalIndexMap >> const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, 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 &  xdot,
const double  dxdot_dx,
const double  dx_dx,
int const  process_id,
GlobalMatrix M,
GlobalMatrix K,
GlobalVector b,
GlobalMatrix Jac 
)
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 211 of file ThermoMechanicalPhaseFieldProcess.cpp.

217 {
218  std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
219  dof_tables;
220  // For the staggered scheme
221  if (process_id == _mechanics_related_process_id)
222  {
223  DBUG(
224  "Assemble the Jacobian equations of "
225  "temperature-deformation in ThermoMechanicalPhaseFieldProcess for "
226  "the staggered scheme.");
227  }
228 
229  if (process_id == _phase_field_process_id)
230  {
231  DBUG(
232  "Assemble the Jacobian equations ofphase field in "
233  "ThermoMechanicalPhaseFieldProcess for the staggered scheme.");
234  }
235  else
236  {
237  DBUG(
238  "Assemble the Jacobian equations of heat conduction in "
239  "ThermoMechanicalPhaseFieldProcess for the staggered scheme.");
240  }
241  dof_tables.emplace_back(
243  dof_tables.emplace_back(
245  dof_tables.emplace_back(getDOFTableByProcessID(_phase_field_process_id));
246 
247  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
248 
251  _local_assemblers, pv.getActiveElementIDs(), dof_tables, t, dt, x, xdot,
252  dxdot_dx, dx_dx, process_id, M, K, b, Jac);
253 }
NumLib::LocalToGlobalIndexMap & getDOFTableByProcessID(const int process_id) const
void assembleWithJacobian(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< std::reference_wrapper< NumLib::LocalToGlobalIndexMap >> const &dof_tables, const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &xdot, const double dxdot_dx, const double dx_dx, 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_prosess_id)
Definition: Process.cpp:303
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
Definition: Process.h:327
MeshLib::Mesh & _mesh
Definition: Process.h:326
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:352

◆ initializeBoundaryConditions()

template<int DisplacementDim>
void ProcessLib::ThermoMechanicalPhaseField::ThermoMechanicalPhaseFieldProcess< DisplacementDim >::initializeBoundaryConditions
overrideprivatevirtual

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

Reimplemented from ProcessLib::Process.

Definition at line 172 of file ThermoMechanicalPhaseFieldProcess.cpp.

173 {
174  // Staggered scheme:
175  // for the equations of temperature-deformation.
179  // for the phase field
183  // for heat conduction
187 }
void initializeProcessBoundaryConditionsAndSourceTerms(const NumLib::LocalToGlobalIndexMap &dof_table, const int process_id)
Definition: Process.cpp:67

◆ 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  mesh.isAxiallySymmetric(), integration_order, _process_data,
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)
NumLib::Extrapolator & getExtrapolator() const
Definition: Process.h:185
SecondaryVariableCollection _secondary_variables
Definition: Process.h:331
void addSecondaryVariable(std::string const &internal_name, SecondaryVariableFunctions &&fcts)
Eigen::Matrix< double, kelvin_vector_dimensions(DisplacementDim), 1, Eigen::ColMajor > KelvinVectorType
Definition: KelvinVector.h:48
void createLocalAssemblers(std::vector< MeshLib::Element * > const &mesh_elements, NumLib::LocalToGlobalIndexMap const &dof_table, std::vector< std::unique_ptr< LocalAssemblerInterface >> &local_assemblers, 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::getIntPtSigma(), ProcessLib::LocalAssemblerInterface::initialize(), MeshLib::Mesh::isAxiallySymmetric(), and ProcessLib::makeExtrapolator().

◆ isLinear()

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

Definition at line 53 of file ThermoMechanicalPhaseFieldProcess.cpp.

54 {
55  return false;
56 }

◆ postNonLinearSolverConcreteProcess()

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

Reimplemented from ProcessLib::Process.

Definition at line 305 of file ThermoMechanicalPhaseFieldProcess.cpp.

309 {
310  if (process_id != _mechanics_related_process_id)
311  {
312  return;
313  }
314 
315  DBUG("PostNonLinearSolver ThermoMechanicalPhaseFieldProcess.");
316  // Calculate strain, stress or other internal variables of mechanics.
317  const bool use_monolithic_scheme = false;
318  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
319 
322  pv.getActiveElementIDs(), getDOFTable(process_id), x, xdot, t, dt,
323  use_monolithic_scheme, process_id);
324 }
void postNonLinearSolver(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, GlobalVector const &xdot, double const t, double const dt, bool const use_monolithic_scheme, int const process_id)
NumLib::LocalToGlobalIndexMap const & getDOFTable(const int process_id) const override
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,
double const  t,
double const  dt,
int const  process_id 
)
overrideprivatevirtual

Reimplemented from ProcessLib::Process.

Definition at line 278 of file ThermoMechanicalPhaseFieldProcess.cpp.

283 {
284  if (process_id != 0)
285  {
286  return;
287  }
288 
289  DBUG("PostTimestep ThermoMechanicalPhaseFieldProcess.");
290  std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
291  auto const n_processes = x.size();
292  dof_tables.reserve(n_processes);
293  for (std::size_t process_id = 0; process_id < n_processes; ++process_id)
294  {
295  dof_tables.push_back(&getDOFTable(process_id));
296  }
297 
298  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
301  _local_assemblers, pv.getActiveElementIDs(), dof_tables, x, t, dt);
302 }
virtual void postTimestep(std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, double const dt)

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 256 of file ThermoMechanicalPhaseFieldProcess.cpp.

261 {
262  DBUG("PreTimestep ThermoMechanicalPhaseFieldProcess.");
263 
264  if (process_id != _mechanics_related_process_id)
265  {
266  return;
267  }
268 
269  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
270 
274  *x[process_id], t, dt);
275 }
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)

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 144 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 128 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 131 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 138 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 141 of file ThermoMechanicalPhaseFieldProcess.h.

◆ _process_data

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

Definition at line 124 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 135 of file ThermoMechanicalPhaseFieldProcess.h.


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