GEOS
PhaseFieldDamageFEM.hpp
Go to the documentation of this file.
1 /*
2  * ------------------------------------------------------------------------------------------------------------
3  * SPDX-License-Identifier: LGPL-2.1-only
4  *
5  * Copyright (c) 2016-2024 Lawrence Livermore National Security LLC
6  * Copyright (c) 2018-2024 Total, S.A
7  * Copyright (c) 2018-2024 The Board of Trustees of the Leland Stanford Junior University
8  * Copyright (c) 2023-2024 Chevron
9  * Copyright (c) 2019- GEOS/GEOSX Contributors
10  * All rights reserved
11  *
12  * See top level LICENSE, COPYRIGHT, CONTRIBUTORS, NOTICE, and ACKNOWLEDGEMENTS files for details.
13  * ------------------------------------------------------------------------------------------------------------
14  */
15 
20 #ifndef GEOS_PHYSICSSOLVERS_SIMPLEPDE_PHASEFIELDDAMAGE_HPP_
21 #define GEOS_PHYSICSSOLVERS_SIMPLEPDE_PHASEFIELDDAMAGE_HPP_
22 
27 
28 struct stabledt
29 {
30  double m_maxdt;
31 };
32 
33 namespace geos
34 {
35 namespace dataRepository
36 {
37 class Group;
38 }
39 class FieldSpecificationBase;
40 class FiniteElementBase;
41 class DomainPartition;
42 
44 {
45 public:
46  PhaseFieldDamageFEM( const string & name, Group * const parent );
47 
48  virtual ~PhaseFieldDamageFEM() override;
49 
50  static string catalogName()
51  {
52  return "PhaseFieldDamageFEM";
53  }
57  string getCatalogName() const override { return catalogName(); }
58 
59  static string coupledSolverAttributePrefix() { return "damage"; }
60 
61  virtual void registerDataOnMesh( Group & meshBodies ) override final;
62 
71  virtual real64 solverStep( real64 const & time_n,
72  real64 const & dt,
73  integer const cycleNumber,
74  DomainPartition & domain ) override;
75 
76  virtual real64 explicitStep( real64 const & time_n,
77  real64 const & dt,
78  integer const cycleNumber,
79  DomainPartition & domain ) override;
80 
81  virtual void setupDofs( DomainPartition const & domain,
82  DofManager & dofManager ) const override;
83 
84  virtual void assembleSystem( real64 const time, real64 const dt,
85  DomainPartition & domain,
86  DofManager const & dofManager,
88  arrayView1d< real64 > const & localRhs ) override;
89 
90  virtual void applyBoundaryConditions( real64 const time, real64 const dt,
91  DomainPartition & domain,
92  DofManager const & dofManager,
94  arrayView1d< real64 > const & localRhs ) override;
95 
96  virtual real64 calculateResidualNorm( real64 const & time_n,
97  real64 const & dt,
98  DomainPartition const & domain,
99  DofManager const & dofManager,
100  arrayView1d< real64 const > const & localRhs ) override;
101 
102  virtual void applySystemSolution( DofManager const & dofManager,
103  arrayView1d< real64 const > const & localSolution,
104  real64 const scalingFactor,
105  real64 const dt,
106  DomainPartition & domain ) override;
107 
108  virtual void updateState( DomainPartition & domain ) override final;
109 
110  virtual void
112  real64 const &,
113  DomainPartition & ) override {}
114 
115  virtual void implicitStepComplete( real64 const & time,
116  real64 const & dt,
117  DomainPartition & domain ) override;
118 
119  virtual void
121 
124  void applyDirichletBCImplicit( real64 const time,
125  DofManager const & dofManager,
126  DomainPartition & domain,
127  CRSMatrixView< real64, globalIndex const > const & localMatrix,
128  arrayView1d< real64 > const & localRhs );
129 
130  void applyIrreversibilityConstraint( DofManager const & dofManager,
131  DomainPartition & domain,
132  CRSMatrixView< real64, globalIndex const > const & localMatrix,
133  arrayView1d< real64 > const & localRhs );
134 
135  virtual void saveSequentialIterationState( DomainPartition & domain ) override;
136 
137  enum class TimeIntegrationOption
138  {
139  SteadyState,
140  ImplicitTransient,
141  ExplicitTransient
142  };
143 
144  enum class LocalDissipation
145  {
146  Linear,
147  Quadratic,
148  };
149 
151  {
152  static constexpr char const * coeffNameString() { return "coeffField"; }
153  static constexpr char const * localDissipationOptionString() { return "localDissipation"; }
154  static constexpr char const * irreversibilityFlagString() { return "irreversibilityFlag"; }
155  static constexpr char const * damageUpperBoundString() { return "damageUpperBound"; }
156  static constexpr char const * solidModelNamesString() { return "solidMaterialNames"; }
157 
158  dataRepository::ViewKey timeIntegrationOption = { "timeIntegrationOption" };
159  dataRepository::ViewKey fieldVarName = { "fieldName" };
160  } PhaseFieldDamageFEMViewKeys;
161 
162  inline ParallelVector const * getSolution() const
163  {
164  return &m_solution;
165  }
166 
167  inline globalIndex getSize() const
168  {
169  return m_matrix.numGlobalRows();
170  }
171 
172  string const & getFieldName() const
173  {
174  return m_fieldName;
175  }
176 
177 protected:
178  virtual void postInputInitialization() override final;
179 
180 private:
181  string m_fieldName;
182  stabledt m_stabledt;
183  TimeIntegrationOption m_timeIntegrationOption;
184  LocalDissipation m_localDissipationOption;
185  integer m_irreversibilityFlag;
186  real64 m_damageUpperBound;
187 
188  array1d< real64 > m_coeff;
189 
191 };
192 
195  "Linear",
196  "Quadratic" );
197 ENUM_STRINGS( PhaseFieldDamageFEM::TimeIntegrationOption,
198  "SteadyState",
199  "ImplicitTransient",
200  "ExplicitTransient" );
201 
202 } /* namespace geos */
203 
204 #endif /* GEOS_PHYSICSSOLVERS_SIMPLEPDE_PHASEFIELDDAMAGE_HPP_ */
The DoFManager is responsible for allocating global dofs, constructing sparsity patterns,...
Definition: DofManager.hpp:44
Partition of the decomposed physical domain. It also manages the connexion information to its neighbo...
virtual void postInputInitialization() override final
string getCatalogName() const override
virtual void registerDataOnMesh(Group &meshBodies) override final
Register wrappers that contain data on the mesh objects.
virtual void saveSequentialIterationState(DomainPartition &domain) override
Save the state of the solver for sequential iteration.
Base class for all physics solvers.
ParallelMatrix m_matrix
System matrix.
ParallelVector m_solution
System solution vector.
virtual void resetStateToBeginningOfStep(DomainPartition &) override
reset state of physics back to the beginning of the step.
virtual void updateState(DomainPartition &domain) override final
Recompute all dependent quantities from primary variables (including constitutive models)
virtual void applyBoundaryConditions(real64 const time, real64 const dt, DomainPartition &domain, DofManager const &dofManager, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs) override
apply boundary condition to system
virtual void applySystemSolution(DofManager const &dofManager, arrayView1d< real64 const > const &localSolution, real64 const scalingFactor, real64 const dt, DomainPartition &domain) override
Function to apply the solution vector to the state.
virtual real64 solverStep(real64 const &time_n, real64 const &dt, integer const cycleNumber, DomainPartition &domain) override
entry function to perform a solver step
virtual void setupDofs(DomainPartition const &domain, DofManager &dofManager) const override
Populate degree-of-freedom manager with fields relevant to this solver.
virtual void assembleSystem(real64 const time, real64 const dt, DomainPartition &domain, DofManager const &dofManager, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs) override
function to assemble the linear system matrix and rhs
virtual void implicitStepComplete(real64 const &time, real64 const &dt, DomainPartition &domain) override
perform cleanup for implicit timestep
virtual real64 explicitStep(real64 const &time_n, real64 const &dt, integer const cycleNumber, DomainPartition &domain) override
Entry function for an explicit time integration step.
virtual real64 calculateResidualNorm(real64 const &time_n, real64 const &dt, DomainPartition const &domain, DofManager const &dofManager, arrayView1d< real64 const > const &localRhs) override
calculate the norm of the global system residual
virtual void implicitStepSetup(real64 const &, real64 const &, DomainPartition &) override
function to perform setup for implicit timestep
Group::wrapperMap::KeyIndex ViewKey
Type alias for KeyIndexT type used for wrapper lookups.
Definition: Group.hpp:1662
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
Definition: DataTypes.hpp:180
LvArray::CRSMatrixView< T, COL_INDEX, localIndex const, LvArray::ChaiBuffer > CRSMatrixView
Alias for CRS Matrix View.
Definition: DataTypes.hpp:310
GEOS_GLOBALINDEX_TYPE globalIndex
Global index type (for indexing objects across MPI partitions).
Definition: DataTypes.hpp:88
ENUM_STRINGS(LinearSolverParameters::SolverType, "direct", "cg", "gmres", "fgmres", "bicgstab", "preconditioner")
Declare strings associated with enumeration values.
double real64
64-bit floating point type.
Definition: DataTypes.hpp:99
std::int32_t integer
Signed integer type.
Definition: DataTypes.hpp:82
Array< T, 1 > array1d
Alias for 1D array.
Definition: DataTypes.hpp:176
LAInterface::ParallelVector ParallelVector
Alias for ParallelVector.
Structure to hold scoped key names.