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

using LocalAssemblerInterface = PhaseFieldLocalAssemblerInterface
 

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, const double t, const double delta_t, 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
 
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 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, 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
 

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 26 of file PhaseFieldProcess.cpp.

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

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 & x_prev,
int const process_id,
GlobalMatrix & M,
GlobalMatrix & K,
GlobalVector & b )
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 185 of file PhaseFieldProcess.cpp.

189{
190 DBUG("Assemble PhaseFieldProcess.");
191
192 std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
193 dof_tables;
194
195 // For the staggered scheme
196 if (process_id == 1)
197 {
198 DBUG(
199 "Assemble the equations of phase field in "
200 "PhaseFieldProcess for the staggered scheme.");
201 }
202 else
203 {
204 DBUG(
205 "Assemble the equations of deformation in "
206 "PhaseFieldProcess for the staggered scheme.");
207 }
208 dof_tables.emplace_back(*_local_to_global_index_map_single_component);
209 dof_tables.emplace_back(*_local_to_global_index_map);
210
211 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
212
213 // Call global assembler for each local assembly item.
216 pv.getActiveElementIDs(), dof_tables, t, dt, x, x_prev, process_id, M,
217 K, b);
218}
void DBUG(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:30
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:155
VectorMatrixAssembler _global_assembler
Definition Process.h:364
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition Process.h:357
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 &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::PhaseField::PhaseFieldProcess< 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 221 of file PhaseFieldProcess.cpp.

225{
226 std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
227 dof_tables;
228
229 // For the staggered scheme
230 if (process_id == 1)
231 {
232 DBUG(
233 "Assemble the Jacobian equations of phase field in "
234 "PhaseFieldProcess for the staggered scheme.");
235 }
236 else
237 {
238 DBUG(
239 "Assemble the Jacobian equations of deformation in "
240 "PhaseFieldProcess for the staggered scheme.");
241 }
242 dof_tables.emplace_back(*_local_to_global_index_map);
243 dof_tables.emplace_back(*_local_to_global_index_map_single_component);
244
245 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
246
247 // Call global assembler for each local assembly item.
250 _local_assemblers, pv.getActiveElementIDs(), dof_tables, t, dt, x,
251 x_prev, process_id, M, K, b, Jac);
252
253 if (process_id == 0)
254 {
256 std::transform(_nodal_forces->begin(), _nodal_forces->end(),
257 _nodal_forces->begin(), [](double val) { return -val; });
258 }
259}
void copyValues(std::vector< double > &u) 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 &x_prev, int const process_id, GlobalMatrix &M, GlobalMatrix &K, GlobalVector &b, GlobalMatrix &Jac)

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

◆ 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 92 of file PhaseFieldProcess.cpp.

93{
94 // For displacement equation.
95 const int mechanics_process_id = 0;
97
98 // TODO move the two data members somewhere else.
99 // for extrapolation of secondary variables of stress or strain
100 std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component{
103 std::make_unique<NumLib::LocalToGlobalIndexMap>(
104 std::move(all_mesh_subsets_single_component),
105 // by location order is needed for output
107
109
110 // For phase field equation.
113}
GlobalSparsityPattern _sparsity_pattern_with_single_component
void constructDofTableOfSpecifiedProcessStaggeredScheme(const int specified_process_id)
Definition Process.cpp:345
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
Definition Process.h:355
MeshLib::Mesh & _mesh
Definition Process.h:354
@ 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 79 of file PhaseFieldProcess.cpp.

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

◆ getMatrixSpecifications()

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

Definition at line 60 of file PhaseFieldProcess.cpp.

62{
63 // For the M process (deformation) in the staggered scheme.
64 if (process_id == 0)
65 {
66 auto const& l = *_local_to_global_index_map;
67 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
68 &l.getGhostIndices(), &this->_sparsity_pattern};
69 }
70
71 // For staggered scheme and phase field process.
73 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
74 &l.getGhostIndices(), &_sparsity_pattern_with_single_component};
75}
GlobalSparsityPattern _sparsity_pattern
Definition Process.h:379

◆ initializeBoundaryConditions()

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

171{
172 // Staggered scheme:
173 // for the equations of deformation.
174 const int mechanical_process_id = 0;
176 *_local_to_global_index_map, mechanical_process_id, media);
177 // for the phase field
178 const int phasefield_process_id = 1;
180 *_local_to_global_index_map_single_component, phasefield_process_id,
181 media);
182}
void initializeProcessBoundaryConditionsAndSourceTerms(const NumLib::LocalToGlobalIndexMap &dof_table, const int process_id, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
Definition Process.cpp:73

◆ 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 116 of file PhaseFieldProcess.cpp.

120{
122 DisplacementDim, PhaseFieldLocalAssembler>(
123 mesh.getElements(), dof_table, _local_assemblers,
124 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
126
128 "sigma",
130 DisplacementDim>::RowsAtCompileTime,
133
135 "epsilon",
137 DisplacementDim>::RowsAtCompileTime,
140
142 "sigma_tensile",
144 DisplacementDim>::RowsAtCompileTime,
147
149 "sigma_compressive",
151 DisplacementDim>::RowsAtCompileTime,
154
156 "eps_tensile",
158 DisplacementDim>::RowsAtCompileTime,
161
162 // Initialize local assemblers after all variables have been set.
166}
virtual void initialize(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
SecondaryVariableCollection _secondary_variables
Definition Process.h:359
NumLib::Extrapolator & getExtrapolator() const
Definition Process.h:196
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 & getIntPtSigmaTensile(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 & getIntPtEpsilonTensile(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 & 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
virtual std::vector< double > const & getIntPtSigmaCompressive(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 ( ) const
override

Definition at line 53 of file PhaseFieldProcess.cpp.

54{
55 return false;
56}

◆ 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 412 of file PhaseFieldProcess.cpp.

414{
415 return process_id == 1;
416}

◆ postNonLinearSolverConcreteProcess()

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

339{
340 _process_data.crack_volume = 0.0;
341
342 if (isPhaseFieldProcess(process_id))
343 {
344 if (_process_data.propagating_pressurized_crack)
345 {
346 auto& u = *x[0];
347 MathLib::LinAlg::scale(const_cast<GlobalVector&>(u),
348 1 / _process_data.pressure);
349 }
350 return;
351 }
352
353 std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
354 dof_tables;
355
356 dof_tables.emplace_back(*_local_to_global_index_map);
357 dof_tables.emplace_back(*_local_to_global_index_map_single_component);
358
359 DBUG("PostNonLinearSolver crack volume computation.");
360
361 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
364 pv.getActiveElementIDs(), dof_tables, x, t, _process_data.crack_volume);
365
366#ifdef USE_PETSC
367 double const crack_volume = _process_data.crack_volume;
368 MPI_Allreduce(&crack_volume, &_process_data.crack_volume, 1, MPI_DOUBLE,
369 MPI_SUM, PETSC_COMM_WORLD);
370#endif
371
372 INFO("Integral of crack: {:g}", _process_data.crack_volume);
373
374 if (_process_data.propagating_pressurized_crack)
375 {
376 _process_data.pressure_old = _process_data.pressure;
377 _process_data.pressure =
378 _process_data.injected_volume / _process_data.crack_volume;
379 _process_data.pressure_error =
380 std::abs(_process_data.pressure_old - _process_data.pressure) /
381 _process_data.pressure;
382 INFO("Internal pressure: {:g} and Pressure error: {:.4e}",
383 _process_data.pressure, _process_data.pressure_error);
384
385 auto& u = *x[0];
386 MathLib::LinAlg::scale(const_cast<GlobalVector&>(u),
387 _process_data.pressure);
388 }
389}
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:35
Global vector based on Eigen vector.
Definition EigenVector.h:25
bool isPhaseFieldProcess(int const process_id) const
void scale(PETScVector &x, PetscScalar const a)
Definition LinAlg.cpp:44
static void executeSelectedMemberOnDereferenced(Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)
virtual void computeCrackIntegral(std::size_t mesh_item_id, std::vector< std::reference_wrapper< NumLib::LocalToGlobalIndexMap > > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, double &crack_volume)=0

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

◆ postTimestepConcreteProcess()

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

Reimplemented from ProcessLib::Process.

Definition at line 282 of file PhaseFieldProcess.cpp.

286{
287 if (isPhaseFieldProcess(process_id))
288 {
289 DBUG("PostTimestep PhaseFieldProcess.");
290
291 _process_data.elastic_energy = 0.0;
292 _process_data.surface_energy = 0.0;
293 _process_data.pressure_work = 0.0;
294
295 std::vector<std::reference_wrapper<NumLib::LocalToGlobalIndexMap>>
296 dof_tables;
297
298 dof_tables.emplace_back(*_local_to_global_index_map);
299 dof_tables.emplace_back(*_local_to_global_index_map_single_component);
300
302 getProcessVariables(process_id)[0];
303
306 pv.getActiveElementIDs(), dof_tables, x, t,
307 _process_data.elastic_energy, _process_data.surface_energy,
308 _process_data.pressure_work);
309
310#ifdef USE_PETSC
311 double const elastic_energy = _process_data.elastic_energy;
312 MPI_Allreduce(&elastic_energy, &_process_data.elastic_energy, 1,
313 MPI_DOUBLE, MPI_SUM, PETSC_COMM_WORLD);
314 double const surface_energy = _process_data.surface_energy;
315 MPI_Allreduce(&surface_energy, &_process_data.surface_energy, 1,
316 MPI_DOUBLE, MPI_SUM, PETSC_COMM_WORLD);
317 double const pressure_work = _process_data.pressure_work;
318 MPI_Allreduce(&pressure_work, &_process_data.pressure_work, 1,
319 MPI_DOUBLE, MPI_SUM, PETSC_COMM_WORLD);
320#endif
321
322 INFO(
323 "Elastic energy: {:g} Surface energy: {:g} Pressure work: {:g} at "
324 "time: {:g} ",
325 _process_data.elastic_energy, _process_data.surface_energy,
326 _process_data.pressure_work, t);
327 if (_process_data.propagating_pressurized_crack)
328 {
329 INFO("Pressure: {:g} at time: {:g} ", _process_data.pressure, t);
330 }
331 }
332}
virtual void computeEnergy(std::size_t mesh_item_id, std::vector< std::reference_wrapper< NumLib::LocalToGlobalIndexMap > > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, double &elastic_energy, double &surface_energy, double &pressure_work)=0

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

◆ 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 262 of file PhaseFieldProcess.cpp.

265{
266 DBUG("PreTimestep PhaseFieldProcess {:d}.", process_id);
267
268 _process_data.injected_volume = t;
269
272
273 ProcessLib::ProcessVariable const& pv = getProcessVariables(process_id)[0];
274
277 pv.getActiveElementIDs(), getDOFTable(process_id), *x[process_id], t,
278 dt);
279}
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(), and ProcessLib::LocalAssemblerInterface::preTimestep().

◆ updateConstraints()

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

Definition at line 392 of file PhaseFieldProcess.cpp.

394{
395 lower.setZero();
398
399 GlobalIndexType const x_begin = _x_previous_timestep->getRangeBegin();
400 GlobalIndexType const x_end = _x_previous_timestep->getRangeEnd();
401
402 for (GlobalIndexType i = x_begin; i < x_end; i++)
403 {
404 if ((*_x_previous_timestep)[i] > _process_data.irreversible_threshold)
405 {
406 upper.set(i, 1.0);
407 }
408 }
409}
GlobalMatrix::IndexType GlobalIndexType
void set(IndexType rowId, double v)
set entry
Definition EigenVector.h:73
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 96 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 99 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 94 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 105 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 107 of file PhaseFieldProcess.h.


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