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, 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::vector< std::reference_wrapper< ProcessVariable > > > const & getProcessVariables () 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)
 
bool requiresNormalization () const override
 
- Public Member Functions inherited from ProcessLib::SubmeshAssemblySupport
virtual std::vector< 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, 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
std::vector< NumLib::LocalToGlobalIndexMap const * > getDOFTables (int const number_of_processes) const
 
NumLib::ExtrapolatorgetExtrapolator () const
 
NumLib::LocalToGlobalIndexMap const & getSingleComponentDOFTable () const
 
void initializeProcessBoundaryConditionsAndSourceTerms (const NumLib::LocalToGlobalIndexMap &dof_table, const int process_id, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
 
void constructMonolithicProcessDofTable ()
 
void constructDofTableOfSpecifiedProcessStaggeredScheme (const int specified_process_id)
 
virtual std::tuple< NumLib::LocalToGlobalIndexMap *, bool > getDOFTableForExtrapolatorData () const
 
std::vector< GlobalIndexTypegetIndicesOfResiduumWithoutInitialCompensation () const override
 
- Protected Attributes inherited from ProcessLib::Process
MeshLib::Mesh_mesh
 
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
 
std::unique_ptr< NumLib::LocalToGlobalIndexMap_local_to_global_index_map
 
SecondaryVariableCollection _secondary_variables
 
CellAverageData cell_average_data_
 
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 28 of file PhaseFieldProcess.cpp.

39 : Process(std::move(name), mesh, std::move(jacobian_assembler), parameters,
40 integration_order, std::move(process_variables),
41 std::move(secondary_variables), use_monolithic_scheme),
42 _process_data(std::move(process_data))
43{
44 if (use_monolithic_scheme)
45 {
47 "Monolithic scheme is not implemented for the PhaseField process.");
48 }
49
51 mesh, "NodalForces", MeshLib::MeshItemType::Node, DisplacementDim);
52
53 _integration_point_writer.emplace_back(
54 std::make_unique<MeshLib::IntegrationPointWriter>(
55 "sigma_ip",
56 static_cast<int>(mesh.getDimension() == 2 ? 4 : 6) /*n components*/,
57 integration_order, _local_assemblers,
58 &LocalAssemblerInterface::getSigma));
59}
#define OGS_FATAL(...)
Definition Error.h:26
unsigned getDimension() const
Returns the dimension of the mesh (determined by the maximum dimension over all elements).
Definition Mesh.h:88
MeshLib::PropertyVector< double > * _nodal_forces
std::vector< std::unique_ptr< LocalAssemblerInterface > > _local_assemblers
PhaseFieldProcessData< DisplacementDim > _process_data
PhaseFieldLocalAssemblerInterface LocalAssemblerInterface
std::string const name
Definition Process.h:362
std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > _integration_point_writer
Definition Process.h:390
Process(std::string name_, MeshLib::Mesh &mesh, std::unique_ptr< AbstractJacobianAssembler > &&jacobian_assembler, std::vector< std::unique_ptr< ParameterLib::ParameterBase > > const &parameters, unsigned const integration_order, std::vector< std::vector< std::reference_wrapper< ProcessVariable > > > &&process_variables, SecondaryVariableCollection &&secondary_variables, const bool use_monolithic_scheme=true)
Definition Process.cpp:44
PropertyVector< T > * getOrCreateMeshProperty(Mesh &mesh, std::string const &property_name, MeshItemType const item_type, int const number_of_components)

References ProcessLib::Process::_integration_point_writer, ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >::_local_assemblers, ProcessLib::PhaseField::PhaseFieldProcess< DisplacementDim >::_nodal_forces, MeshLib::Mesh::getDimension(), MeshLib::getOrCreateMeshProperty(), ProcessLib::PhaseField::PhaseFieldLocalAssemblerInterface::getSigma(), 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 197 of file PhaseFieldProcess.cpp.

201{
202 DBUG("Assemble PhaseFieldProcess.");
203
204 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
205
206 // For the staggered scheme
207 if (process_id == 1)
208 {
209 DBUG(
210 "Assemble the equations of phase field in "
211 "PhaseFieldProcess for the staggered scheme.");
212 }
213 else
214 {
215 DBUG(
216 "Assemble the equations of deformation in "
217 "PhaseFieldProcess for the staggered scheme.");
218 }
219 dof_tables.emplace_back(_local_to_global_index_map_single_component.get());
220 dof_tables.emplace_back(_local_to_global_index_map.get());
221
222 // Call global assembler for each local assembly item.
225 getActiveElementIDs(), dof_tables, t, dt, x, x_prev, process_id, &M, &K,
226 &b);
227}
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::size_t > const & getActiveElementIDs() const
Definition Process.h:167
VectorMatrixAssembler _global_assembler
Definition Process.h:377
std::unique_ptr< NumLib::LocalToGlobalIndexMap > _local_to_global_index_map
Definition Process.h:368
void assemble(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalMatrix *M, GlobalMatrix *K, GlobalVector *b)
static void executeSelectedMemberDereferenced(Object &object, Method method, Container const &container, std::vector< std::size_t > const &active_container_ids, Args &&... args)

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

◆ 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,
GlobalVector & b,
GlobalMatrix & Jac )
overrideprivatevirtual

Implements ProcessLib::Process.

Definition at line 230 of file PhaseFieldProcess.cpp.

234{
235 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
236
237 // For the staggered scheme
238 if (process_id == 1)
239 {
240 DBUG(
241 "Assemble the Jacobian equations of phase field in "
242 "PhaseFieldProcess for the staggered scheme.");
243 }
244 else
245 {
246 DBUG(
247 "Assemble the Jacobian equations of deformation in "
248 "PhaseFieldProcess for the staggered scheme.");
249 }
250 dof_tables.emplace_back(_local_to_global_index_map.get());
251 dof_tables.emplace_back(_local_to_global_index_map_single_component.get());
252
253 // Call global assembler for each local assembly item.
256 _local_assemblers, getActiveElementIDs(), dof_tables, t, dt, x, x_prev,
257 process_id, &b, &Jac);
258
259 if (process_id == 0)
260 {
262 std::transform(_nodal_forces->begin(), _nodal_forces->end(),
263 _nodal_forces->begin(), [](double val) { return -val; });
264 }
265}
void copyValues(std::vector< double > &u) const
void assembleWithJacobian(std::size_t const mesh_item_id, LocalAssemblerInterface &local_assembler, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, const double t, double const dt, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, int const process_id, GlobalVector *b, GlobalMatrix *Jac)

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

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

102{
103 // For displacement equation.
104 const int mechanics_process_id = 0;
106
107 // TODO move the two data members somewhere else.
108 // for extrapolation of secondary variables of stress or strain
109 std::vector<MeshLib::MeshSubset> all_mesh_subsets_single_component{
112 std::make_unique<NumLib::LocalToGlobalIndexMap>(
113 std::move(all_mesh_subsets_single_component),
114 // by location order is needed for output
116
118
119 // For phase field equation.
122}
GlobalSparsityPattern _sparsity_pattern_with_single_component
void constructDofTableOfSpecifiedProcessStaggeredScheme(const int specified_process_id)
Definition Process.cpp:354
std::unique_ptr< MeshLib::MeshSubset const > _mesh_subset_all_nodes
Definition Process.h:366
MeshLib::Mesh & _mesh
Definition Process.h:365
@ 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 88 of file PhaseFieldProcess.cpp.

89{
90 // For the M process (deformation) in the staggered scheme.
91 if (process_id == 0)
92 {
94 }
95
96 // For the equation of phasefield
98}

◆ getMatrixSpecifications()

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

Definition at line 69 of file PhaseFieldProcess.cpp.

71{
72 // For the M process (deformation) in the staggered scheme.
73 if (process_id == 0)
74 {
75 auto const& l = *_local_to_global_index_map;
76 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
77 &l.getGhostIndices(), &this->_sparsity_pattern};
78 }
79
80 // For staggered scheme and phase field process.
82 return {l.dofSizeWithoutGhosts(), l.dofSizeWithoutGhosts(),
83 &l.getGhostIndices(), &_sparsity_pattern_with_single_component};
84}
GlobalSparsityPattern _sparsity_pattern
Definition Process.h:392

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

183{
184 // Staggered scheme:
185 // for the equations of deformation.
186 const int mechanical_process_id = 0;
188 *_local_to_global_index_map, mechanical_process_id, media);
189 // for the phase field
190 const int phasefield_process_id = 1;
192 *_local_to_global_index_map_single_component, phasefield_process_id,
193 media);
194}
void initializeProcessBoundaryConditionsAndSourceTerms(const NumLib::LocalToGlobalIndexMap &dof_table, const int process_id, std::map< int, std::shared_ptr< MaterialPropertyLib::Medium > > const &media)
Definition Process.cpp:90

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

129{
131 DisplacementDim, PhaseFieldLocalAssembler>(
132 mesh.getElements(), dof_table, _local_assemblers,
133 NumLib::IntegrationOrder{integration_order}, mesh.isAxiallySymmetric(),
135
137 "sigma",
139 DisplacementDim>::RowsAtCompileTime,
142
144 "epsilon",
146 DisplacementDim>::RowsAtCompileTime,
149
151 "sigma_tensile",
153 DisplacementDim>::RowsAtCompileTime,
156
158 "sigma_compressive",
160 DisplacementDim>::RowsAtCompileTime,
163
165 "eps_tensile",
167 DisplacementDim>::RowsAtCompileTime,
170
173
174 // Initialize local assemblers after all variables have been set.
178}
virtual void initialize(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)
SecondaryVariableCollection _secondary_variables
Definition Process.h:370
NumLib::Extrapolator & getExtrapolator() const
Definition Process.h:208
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)
void setIPDataInitialConditions(std::vector< std::unique_ptr< MeshLib::IntegrationPointWriter > > const &_integration_point_writer, MeshLib::Properties const &mesh_properties, LocalAssemblersVector &local_assemblers)
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(), MeshLib::Mesh::getProperties(), ProcessLib::LocalAssemblerInterface::initialize(), MeshLib::Mesh::isAxiallySymmetric(), ProcessLib::makeExtrapolator(), and ProcessLib::setIPDataInitialConditions().

◆ isLinear()

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

Definition at line 62 of file PhaseFieldProcess.cpp.

63{
64 return false;
65}

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

409{
410 return process_id == 1;
411}

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

337{
338 _process_data.crack_volume = 0.0;
339
340 if (isPhaseFieldProcess(process_id))
341 {
342 if (_process_data.propagating_pressurized_crack)
343 {
344 auto& u = *x[0];
345 MathLib::LinAlg::scale(const_cast<GlobalVector&>(u),
346 1 / _process_data.pressure);
347 }
348 return;
349 }
350
351 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
352
353 dof_tables.emplace_back(_local_to_global_index_map.get());
354 dof_tables.emplace_back(_local_to_global_index_map_single_component.get());
355
356 DBUG("PostNonLinearSolver crack volume computation.");
357
360 getActiveElementIDs(), dof_tables, x, t, _process_data.crack_volume);
361
362#ifdef USE_PETSC
364 _process_data.crack_volume, MPI_SUM, BaseLib::MPI::Mpi{});
365#endif
366
367 INFO("Integral of crack: {:g}", _process_data.crack_volume);
368
369 if (_process_data.propagating_pressurized_crack)
370 {
371 _process_data.pressure_old = _process_data.pressure;
372 _process_data.pressure =
373 _process_data.injected_volume / _process_data.crack_volume;
374 _process_data.pressure_error =
375 std::abs(_process_data.pressure_old - _process_data.pressure) /
376 _process_data.pressure;
377 INFO("Internal pressure: {:g} and Pressure error: {:.4e}",
378 _process_data.pressure, _process_data.pressure_error);
379
380 auto& u = *x[0];
381 MathLib::LinAlg::scale(const_cast<GlobalVector&>(u),
382 _process_data.pressure);
383 }
384}
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
static T allreduce(T const &value, MPI_Op const &mpi_op, Mpi const &mpi)
Definition MPI.h:128
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< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, double &crack_volume)=0

References BaseLib::MPI::allreduce(), DBUG(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), 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 285 of file PhaseFieldProcess.cpp.

289{
290 if (isPhaseFieldProcess(process_id))
291 {
292 DBUG("PostTimestep PhaseFieldProcess.");
293
294 _process_data.elastic_energy = 0.0;
295 _process_data.surface_energy = 0.0;
296 _process_data.pressure_work = 0.0;
297
298 std::vector<NumLib::LocalToGlobalIndexMap const*> dof_tables;
299
300 dof_tables.emplace_back(_local_to_global_index_map.get());
301 dof_tables.emplace_back(
303
306 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 BaseLib::MPI::Mpi mpi{};
312 _process_data.elastic_energy =
313 BaseLib::MPI::allreduce(_process_data.elastic_energy, MPI_SUM, mpi);
314 _process_data.surface_energy =
315 BaseLib::MPI::allreduce(_process_data.surface_energy, MPI_SUM, mpi);
316 _process_data.pressure_work =
317 BaseLib::MPI::allreduce(_process_data.pressure_work, MPI_SUM, mpi);
318#endif
319
320 INFO(
321 "Elastic energy: {:g} Surface energy: {:g} Pressure work: {:g} at "
322 "time: {:g} ",
323 _process_data.elastic_energy, _process_data.surface_energy,
324 _process_data.pressure_work, t);
325 if (_process_data.propagating_pressurized_crack)
326 {
327 INFO("Pressure: {:g} at time: {:g} ", _process_data.pressure, t);
328 }
329 }
330}
virtual void computeEnergy(std::size_t mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, double &elastic_energy, double &surface_energy, double &pressure_work)=0

References BaseLib::MPI::allreduce(), DBUG(), NumLib::SerialExecutor::executeSelectedMemberOnDereferenced(), 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 268 of file PhaseFieldProcess.cpp.

271{
272 DBUG("PreTimestep PhaseFieldProcess {:d}.", process_id);
273
274 _process_data.injected_volume = t;
275
278
281 getActiveElementIDs(), getDOFTable(process_id), *x[process_id], t, dt);
282}
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(), 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 387 of file PhaseFieldProcess.cpp.

389{
390 lower.setZero();
393
394 GlobalIndexType const x_begin = _x_previous_timestep->getRangeBegin();
395 GlobalIndexType const x_end = _x_previous_timestep->getRangeEnd();
396
397 for (GlobalIndexType i = x_begin; i < x_end; i++)
398 {
399 if ((*_x_previous_timestep)[i] > _process_data.irreversible_threshold)
400 {
401 upper.set(i, 1.0);
402 }
403 }
404}
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

◆ _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 98 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 93 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 104 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 106 of file PhaseFieldProcess.h.


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