OGS
ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim > Class Template Referencefinal

Detailed Description

template<int DisplacementDim>
class ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >

Definition at line 22 of file PhaseFieldProcess.h.

#include <PhaseFieldProcess.h>

Inheritance diagram for ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >:
[legend]
Collaboration diagram for ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >:
[legend]

Public Member Functions

 PhaseFieldProcess (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, PhaseFieldProcessData< DisplacementDim > &&process_data, SecondaryVariableCollection &&secondary_variables, bool const use_monolithic_scheme)
 
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, 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 & 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)
 

Private Types

using LocalAssemblerInterface = PhaseFieldLocalAssemblerInterface
 

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, 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, const double t, const double delta_t, int const process_id) override
 
void postNonLinearSolverConcreteProcess (GlobalVector const &x, GlobalVector const &xdot, const double t, double const dt, int const process_id) override
 
void updateConstraints (GlobalVector &lower, GlobalVector &upper, int const process_id) override
 
bool isPhaseFieldProcess (int const process_id) const
 

Private Attributes

PhaseFieldProcessData< DisplacementDim > _process_data
 
std::vector< std::unique_ptr< LocalAssemblerInterface > > _local_assemblers
 
std::unique_ptr< NumLib::LocalToGlobalIndexMap_local_to_global_index_map_single_component
 
MeshLib::PropertyVector< double > * _nodal_forces = nullptr
 
GlobalSparsityPattern _sparsity_pattern_with_single_component
 
std::unique_ptr< GlobalVector_x_previous_timestep
 

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
 
- 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
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
 
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
 
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
 

Member Typedef Documentation

◆ LocalAssemblerInterface

template<int DisplacementDim>
using ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >::LocalAssemblerInterface = PhaseFieldLocalAssemblerInterface
private

Definition at line 51 of file PhaseFieldProcess.h.

Constructor & Destructor Documentation

◆ PhaseFieldProcess()

template<int DisplacementDim>
ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >::PhaseFieldProcess ( 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,
PhaseFieldProcessData< DisplacementDim > &&  process_data,
SecondaryVariableCollection &&  secondary_variables,
bool const  use_monolithic_scheme 
)

Definition at line 25 of file PhaseFieldProcess.cpp.

36  : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
37  integration_order, std::move(process_variables),
38  std::move(secondary_variables), use_monolithic_scheme),
39  _process_data(std::move(process_data))
40 {
41  if (use_monolithic_scheme)
42  {
43  OGS_FATAL(
44  "Monolithic scheme is not implemented for the PhaseField process.");
45  }
46 
47  _nodal_forces = MeshLib::getOrCreateMeshProperty<double>(
48  mesh, "NodalForces", MeshLib::MeshItemType::Node, DisplacementDim);
49 }
#define OGS_FATAL(...)
Definition: Error.h:26
MeshLib::PropertyVector< double > * _nodal_forces
PhaseFieldProcessData< DisplacementDim > _process_data
std::string const name
Definition: Process.h:327
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:24

References ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >::_nodal_forces, MeshLib::Node, and OGS_FATAL.

Member Function Documentation

◆ assembleConcreteProcess()

template<int DisplacementDim>
void ProcessLib::PhaseField::PhaseFieldProcess< 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 161 of file PhaseFieldProcess.cpp.

165 {
166  DBUG("Assemble PhaseFieldProcess.");
167 
168  std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
169  dof_tables;
170 
171  // For the staggered scheme
172  if (process_id == 1)
173  {
174  DBUG(
175  "Assemble the equations of phase field in "
176  "PhaseFieldProcess for the staggered scheme.");
177  }
178  else
179  {
180  DBUG(
181  "Assemble the equations of deformation in "
182  "PhaseFieldProcess for the staggered scheme.");
183  }
184  dof_tables.emplace_back(*_local_to_global_index_map_single_component);
185  dof_tables.emplace_back(*_local_to_global_index_map);
186 
187  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
188 
189  // Call global assembler for each local assembly item.
192  pv.getActiveElementIDs(), dof_tables, t, dt, x, xdot, process_id, M, K,
193  b);
194 }
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition: Logging.h:29
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map_single_component
std::vector< std::unique_ptr< LocalAssemblerInterface > > _local_assemblers
std::vector< std::size_t > const & getActiveElementIDs() const
std::vector< std::reference_wrapper< ProcessVariable > > const & getProcessVariables(const int process_id) const
Definition: Process.h:148
VectorMatrixAssembler _global_assembler
Definition: Process.h:337
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition: Process.h:333
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 const double t
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(), ProcessLib::ProcessVariable::getActiveElementIDs(), and MathLib::t.

◆ assembleWithJacobianConcreteProcess()

template<int DisplacementDim>
void ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >::assembleWithJacobianConcreteProcess ( 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,
GlobalMatrix Jac 
)
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 197 of file PhaseFieldProcess.cpp.

201 {
202  std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
203  dof_tables;
204 
205  // For the staggered scheme
206  if (process_id == 1)
207  {
208  DBUG(
209  "Assemble the Jacobian equations of phase field in "
210  "PhaseFieldProcess for the staggered scheme.");
211  }
212  else
213  {
214  DBUG(
215  "Assemble the Jacobian equations of deformation in "
216  "PhaseFieldProcess for the staggered scheme.");
217  }
218  dof_tables.emplace_back(*_local_to_global_index_map);
219  dof_tables.emplace_back(*_local_to_global_index_map_single_component);
220 
221  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
222 
223  // Call global assembler for each local assembly item.
226  _local_assemblers, pv.getActiveElementIDs(), dof_tables, t, dt, x, xdot,
227  process_id, M, K, b, Jac);
228 
229  if (process_id == 0)
230  {
232  std::transform(_nodal_forces->begin(), _nodal_forces->end(),
233  _nodal_forces->begin(), [](double val) { return -val; });
234  }
235 }
void copyValues(std::vector< double > &u) const
Definition: EigenVector.h:106
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, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac)

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

◆ constructDofTable()

template<int DisplacementDim>
void ProcessLib::PhaseField::PhaseFieldProcess< 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 91 of file PhaseFieldProcess.cpp.

92 {
93  // For displacement equation.
94  const int mechanics_process_id = 0;
96 
97  // TODO move the two data members somewhere else.
98  // for extrapolation of secondary variables of stress or strain
99  std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component{
102  std::make_unique<NumLib::LocalToGlobalIndexMap>(
103  std::move(all_mesh_subsets_single_component),
104  // by location order is needed for output
106 
108 
109  // For phase field equation.
112 }
GlobalSparsityPattern _sparsity_pattern_with_single_component
void constructDofTableOfSpecifiedProcessStaggeredScheme(const int specified_prosess_id)
Definition: Process.cpp:307
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
Definition: Process.h:331
MeshLib::Mesh & _mesh
Definition: Process.h:330
@ 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::PhaseField::PhaseFieldProcess< DisplacementDim >::getDOFTable ( const int  process_id) const
overridevirtual

Reimplemented from ProcessLib::Process.

Definition at line 78 of file PhaseFieldProcess.cpp.

79 {
80  // For the M process (deformation) in the staggered scheme.
81  if (process_id == 0)
82  {
84  }
85 
86  // For the equation of phasefield
88 }

◆ getMatrixSpecifications()

template<int DisplacementDim>
MathLib::MatrixSpecifications ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >::getMatrixSpecifications ( const int  process_id) const
override

Definition at line 59 of file PhaseFieldProcess.cpp.

61 {
62  // For the M process (deformation) in the staggered scheme.
63  if (process_id == 0)
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.
72  return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
73  &l.getGhostIndices(), &_sparsity_pattern_with_single_component};
74 }
GlobalSparsityPattern _sparsity_pattern
Definition: Process.h:356

◆ initializeBoundaryConditions()

template<int DisplacementDim>
void ProcessLib::PhaseField::PhaseFieldProcess< 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 147 of file PhaseFieldProcess.cpp.

148 {
149  // Staggered scheme:
150  // for the equations of deformation.
151  const int mechanical_process_id = 0;
153  *_local_to_global_index_map, mechanical_process_id);
154  // for the phase field
155  const int phasefield_process_id = 1;
157  *_local_to_global_index_map_single_component, phasefield_process_id);
158 }
void initializeProcessBoundaryConditionsAndSourceTerms(const NumLib::LocalToGlobalIndexMap &dof_table, const int process_id)
Definition: Process.cpp:69

◆ initializeConcreteProcess()

template<int DisplacementDim>
void ProcessLib::PhaseField::PhaseFieldProcess< 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 115 of file PhaseFieldProcess.cpp.

119 {
121  DisplacementDim, PhaseFieldLocalAssembler>(
122  mesh.getElements(), dof_table, _local_assemblers,
123  NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
124  _process_data);
125 
127  "sigma",
129  DisplacementDim>::RowsAtCompileTime,
132 
134  "epsilon",
136  DisplacementDim>::RowsAtCompileTime,
139 
140  // Initialize local assemblers after all variables have been set.
144 }
virtual void initialize(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
NumLib::Extrapolator & getExtrapolator() const
Definition: Process.h:184
SecondaryVariableCollection _secondary_variables
Definition: Process.h:335
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, NumLib::IntegrationOrder const integration_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 & getIntPtEpsilon(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
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::getElements(), ProcessLib::LocalAssemblerInterface::initialize(), MeshLib::Mesh::isAxiallySymmetric(), and ProcessLib::makeExtrapolator().

◆ isLinear()

template<int DisplacementDim>
bool ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >::isLinear
override

Definition at line 52 of file PhaseFieldProcess.cpp.

53 {
54  return false;
55 }

◆ isPhaseFieldProcess()

template<int DisplacementDim>
bool ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >::isPhaseFieldProcess ( int const  process_id) const
private

Check whether the process represented by process_id is/has mechanical process. In the present implementation, the mechanical process has process_id == 0 in the staggered scheme.

Definition at line 364 of file PhaseFieldProcess.cpp.

366 {
367  return process_id == 1;
368 }

◆ postNonLinearSolverConcreteProcess()

template<int DisplacementDim>
void ProcessLib::PhaseField::PhaseFieldProcess< 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 292 of file PhaseFieldProcess.cpp.

295 {
296  _process_data.crack_volume = 0.0;
297 
298  if (!isPhaseFieldProcess(process_id))
299  {
300  std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
301  dof_tables;
302 
303  dof_tables.emplace_back(*_local_to_global_index_map);
304  dof_tables.emplace_back(*_local_to_global_index_map_single_component);
305 
306  DBUG("PostNonLinearSolver crack volume computation.");
307 
308  ProcessLib::ProcessVariable const& pv =
309  getProcessVariables(process_id)[0];
312  pv.getActiveElementIDs(), dof_tables, x, t,
313  _process_data.crack_volume, _coupled_solutions);
314 
315  INFO("Integral of crack: {:g}", _process_data.crack_volume);
316 
317  if (_process_data.hydro_crack)
318  {
319  _process_data.pressure_old = _process_data.pressure;
320  _process_data.pressure =
321  _process_data.injected_volume / _process_data.crack_volume;
322  _process_data.pressure_error =
323  std::abs(_process_data.pressure_old - _process_data.pressure) /
324  _process_data.pressure;
325  INFO("Internal pressure: {:g} and Pressure error: {:.4e}",
326  _process_data.pressure, _process_data.pressure_error);
327  auto& u = *_coupled_solutions->coupled_xs[0];
328  MathLib::LinAlg::scale(const_cast<GlobalVector&>(u),
329  _process_data.pressure);
330  }
331  }
332  else
333  {
334  if (_process_data.hydro_crack)
335  {
336  auto& u = *_coupled_solutions->coupled_xs[0];
337  MathLib::LinAlg::scale(const_cast<GlobalVector&>(u),
338  1 / _process_data.pressure);
339  }
340  }
341 }
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition: Logging.h:34
Global vector based on Eigen vector.
Definition: EigenVector.h:28
bool isPhaseFieldProcess(int const process_id) const
CoupledSolutionsForStaggeredScheme * _coupled_solutions
Definition: Process.h:343
void scale(PETScVector &x, PetscScalar const a)
Definition: LinAlg.cpp:44
static const double u
static void executeSelectedMemberOnDereferenced(Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)
std::vector< GlobalVector * > const & coupled_xs
References to the current solutions of the coupled processes.
virtual void computeCrackIntegral(std::size_t mesh_item_id, std::vector< std::reference_wrapper< NumLib::LocalToGlobalIndexMap >> const &dof_tables, GlobalVector const &x, double const t, double &crack_volume, CoupledSolutionsForStaggeredScheme const *const cpl_xs)=0

References DBUG(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::ProcessVariable::getActiveElementIDs(), INFO(), MathLib::LinAlg::scale(), MathLib::t, and MathLib::u.

◆ postTimestepConcreteProcess()

template<int DisplacementDim>
void ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >::postTimestepConcreteProcess ( std::vector< GlobalVector * > const &  x,
const double  t,
const double  delta_t,
int const  process_id 
)
overrideprivatevirtual

Reimplemented from ProcessLib::Process.

Definition at line 258 of file PhaseFieldProcess.cpp.

261 {
262  if (isPhaseFieldProcess(process_id))
263  {
264  DBUG("PostTimestep PhaseFieldProcess.");
265 
266  _process_data.elastic_energy = 0.0;
267  _process_data.surface_energy = 0.0;
268  _process_data.pressure_work = 0.0;
269 
270  std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
271  dof_tables;
272 
273  dof_tables.emplace_back(*_local_to_global_index_map);
274  dof_tables.emplace_back(*_local_to_global_index_map_single_component);
275 
276  ProcessLib::ProcessVariable const& pv =
277  getProcessVariables(process_id)[0];
278 
281  pv.getActiveElementIDs(), dof_tables, *x[process_id], t,
282  _process_data.elastic_energy, _process_data.surface_energy,
283  _process_data.pressure_work, _coupled_solutions);
284 
285  INFO("Elastic energy: {:g} Surface energy: {:g} Pressure work: {:g} ",
286  _process_data.elastic_energy, _process_data.surface_energy,
287  _process_data.pressure_work);
288  }
289 }
virtual void computeEnergy(std::size_t mesh_item_id, std::vector< std::reference_wrapper< NumLib::LocalToGlobalIndexMap >> const &dof_tables, GlobalVector const &x, double const t, double &elastic_energy, double &surface_energy, double &pressure_work, CoupledSolutionsForStaggeredScheme const *const cpl_xs)=0

References DBUG(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), ProcessLib::ProcessVariable::getActiveElementIDs(), INFO(), and MathLib::t.

◆ preTimestepConcreteProcess()

template<int DisplacementDim>
void ProcessLib::PhaseField::PhaseFieldProcess< 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 238 of file PhaseFieldProcess.cpp.

241 {
242  DBUG("PreTimestep PhaseFieldProcess {:d}.", process_id);
243 
244  _process_data.injected_volume = t;
245 
248 
249  ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
250 
253  pv.getActiveElementIDs(), getDOFTable(process_id), *x[process_id], t,
254  dt);
255 }
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
std::unique_ptr< GlobalVector > _x_previous_timestep

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

◆ updateConstraints()

template<int DisplacementDim>
void ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >::updateConstraints ( GlobalVector lower,
GlobalVector upper,
int const  process_id 
)
overrideprivate

Definition at line 344 of file PhaseFieldProcess.cpp.

346 {
347  lower.setZero();
350 
351  GlobalIndexType const x_begin = _x_previous_timestep->getRangeBegin();
352  GlobalIndexType const x_end = _x_previous_timestep->getRangeEnd();
353 
354  for (GlobalIndexType i = x_begin; i < x_end; i++)
355  {
356  if ((*_x_previous_timestep)[i] > _process_data.irreversible_threshold)
357  {
358  upper.set(i, 1.0);
359  }
360  }
361 }
GlobalMatrix::IndexType GlobalIndexType
void set(IndexType rowId, double v)
set entry
Definition: EigenVector.h:76
void copy(PETScVector const &x, PETScVector &y)
Definition: LinAlg.cpp:37
void setLocalAccessibleVector(PETScVector const &x)
Definition: LinAlg.cpp:27

References MathLib::LinAlg::copy(), MathLib::EigenVector::set(), MathLib::LinAlg::setLocalAccessibleVector(), and MathLib::EigenVector::setZero().

Member Data Documentation

◆ _local_assemblers

template<int DisplacementDim>
std::vector<std::unique_ptr<LocalAssemblerInterface> > ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >::_local_assemblers
private

Definition at line 93 of file PhaseFieldProcess.h.

◆ _local_to_global_index_map_single_component

template<int DisplacementDim>
std::unique_ptr<NumLib::LocalToGlobalIndexMap> ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >::_local_to_global_index_map_single_component
private

Definition at line 96 of file PhaseFieldProcess.h.

◆ _nodal_forces

template<int DisplacementDim>
MeshLib::PropertyVector<double>* ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >::_nodal_forces = nullptr
private

◆ _process_data

template<int DisplacementDim>
PhaseFieldProcessData<DisplacementDim> ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >::_process_data
private

Definition at line 91 of file PhaseFieldProcess.h.

◆ _sparsity_pattern_with_single_component

template<int DisplacementDim>
GlobalSparsityPattern ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >::_sparsity_pattern_with_single_component
private

Sparsity pattern for the phase field equation, and it is initialized only if the staggered scheme is used.

Definition at line 102 of file PhaseFieldProcess.h.

◆ _x_previous_timestep

template<int DisplacementDim>
std::unique_ptr<GlobalVector> ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >::_x_previous_timestep
private

Definition at line 104 of file PhaseFieldProcess.h.


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