tpsamodel.hpp
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 Copyright 2025 NORCE AS
5
6 This file is part of the Open Porous Media project (OPM).
7
8 OPM is free software: you can redistribute it and/or modify
9 it under the terms of the GNU General Public License as published by
10 the Free Software Foundation, either version 2 of the License, or
11 (at your option) any later version.
12
13 OPM is distributed in the hope that it will be useful,
14 but WITHOUT ANY WARRANTY; without even the implied warranty of
15 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16 GNU General Public License for more details.
17
18 You should have received a copy of the GNU General Public License
19 along with OPM. If not, see <http://www.gnu.org/licenses/>.
20
21 Consult the COPYING file in the top-level source directory of this
22 module for the precise wording of the license and the list of
23 copyright holders.
24*/
25#ifndef TPSA_MODEL_HPP
26#define TPSA_MODEL_HPP
27
28#include <dune/grid/common/gridenums.hh>
29
30#include <opm/grid/utility/ElementChunks.hpp>
31
32#include <opm/material/common/MathToolbox.hpp>
33
37
38#include <array>
39#include <memory>
40
41
42namespace Opm {
43
49template <class TypeTag>
51{
63
64 enum { dimWorld = GridView::dimensionworld };
65 enum { historySize = getPropValue<TypeTag, Properties::SolutionHistorySizeTPSA>() };
66 enum { numEq = getPropValue<TypeTag, Properties::NumEqTPSA>() };
67
68 enum { disp0Idx = Indices::disp0Idx };
69 enum { rot0Idx = Indices::rot0Idx };
70 enum { solidPres0Idx = Indices::solidPres0Idx };
71
72 using MaterialState = MaterialStateTPSA<Evaluation>;
73
74 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
75 using SymTensor = Dune::FieldVector<Scalar, 6>;
76
77public:
82 {
83 protected:
84 SolutionVector blockVector_;
85 public:
92 TpsaBlockVectorWrapper(const std::string&, const std::size_t size)
93 : blockVector_(size)
94 {}
95
100
105 {
106 TpsaBlockVectorWrapper result("dummy", 3);
107 result.blockVector_[0] = 1.0;
108 result.blockVector_[1] = 2.0;
109 result.blockVector_[2] = 3.0;
110
111 return result;
112 }
113
119 SolutionVector& blockVector()
120 { return blockVector_; }
121
127 const SolutionVector& blockVector() const
128 { return blockVector_; }
129
136 bool operator==(const TpsaBlockVectorWrapper& wrapper) const
137 {
138 return std::equal(this->blockVector_.begin(), this->blockVector_.end(),
139 wrapper.blockVector_.begin(), wrapper.blockVector_.end());
140 }
141
147 template<class Serializer>
148 void serializeOp(Serializer& serializer)
149 {
150 serializer(blockVector_);
151 }
152 };
153
154 // ///
155 // Public functions
156 // ///
162 explicit TpsaModel(Simulator& simulator)
163 : linearizer_(std::make_unique<Linearizer>())
164 , newtonMethod_(simulator)
165 , simulator_(simulator)
166 , element_chunks_(simulator.gridView(), Dune::Partitions::all, ThreadManager::maxThreads())
167 {
168 // Initialize equation weights to 1.0
169 eqWeights_.resize(numEq, 1.0);
170
171 // Initialize historic solution vectors
172 // OBS: need at least history size = 2, due to time-derivative of solid-pressure in Flow coupling term
173 const std::size_t numDof = simulator_.model().numGridDof();
174 for (unsigned timeIdx = 0; timeIdx < historySize; ++timeIdx) {
175 solution_[timeIdx] = std::make_unique<TpsaBlockVectorWrapper>("solution", numDof);
176 }
177 }
178
183 {
184 // Initialize the linearizer
185 linearizer_->init(simulator_);
186
187 // Resize material state vector
189 }
190
194 static void registerParameters()
195 {
196 // Newton method parameters
198 }
199
204 {
205 // Update historic solution
206 solution(/*timeIdx=*/1) = solution(/*timeIdx=*/0);
207 }
208
215 {
216 // Syncronize the solution on the ghost and overlap elements
217 using GhostSyncHandle = GridCommHandleGhostSync<PrimaryVariables,
218 SolutionVector,
219 DofMapper,
220 /*commCodim=*/0>;
221
222 auto ghostSync = GhostSyncHandle(solution(/*timeIdx=*/0),
223 simulator_.model().dofMapper());
224
225 simulator_.gridView().communicate(ghostSync,
226 Dune::InteriorBorder_All_Interface,
227 Dune::ForwardCommunication);
228 }
229
237 void updateMaterialState(const unsigned /*timeIdx*/)
238 {
239 // Loop over all elements chuncks and update material state from current solution
240 const auto& elementMapper = simulator_.model().elementMapper();
241#ifdef _OPENMP
242#pragma omp parallel for
243#endif
244 for (const auto& chunk : element_chunks_) {
245 for (const auto& elem : chunk) {
246 const unsigned globalIdx = elementMapper.index(elem);
247 auto& currSol = solution(/*timeIdx=*/0)[globalIdx];
248 setMaterialState_(globalIdx, /*timeIdx=*/0, currSol);
249 }
250 }
251 }
252
253 // ///
254 // Public get and set functions
255 // ///
261 const Linearizer& linearizer() const
262 {
263 return *linearizer_;
264 }
265
271 Linearizer& linearizer()
272 {
273 return *linearizer_;
274 }
275
281 const NewtonMethod& newtonMethod() const
282 {
283 return newtonMethod_;
284 }
285
291 NewtonMethod& newtonMethod()
292 {
293 return newtonMethod_;
294 }
295
302 const SolutionVector& solution(unsigned timeIdx) const
303 {
304 return solution_[timeIdx]->blockVector();
305 }
306
313 SolutionVector& solution(unsigned timeIdx)
314 {
315 return solution_[timeIdx]->blockVector();
316 }
317
322 std::size_t numGridDof() const
323 {
324 return simulator_.model().numGridDof();
325 }
326
331 std::size_t numTotalDof() const
332 {
333 return numGridDof() + numAuxiliaryDof();
334 }
335
342 Scalar dofTotalVolume(unsigned globalIdx) const
343 {
344 return simulator_.model().dofTotalVolume(globalIdx);
345 }
346
354 Scalar eqWeight(unsigned /*dofIdx*/, unsigned eqIdx) const
355 {
356 return eqWeights_[eqIdx];
357 }
358
365 void setEqWeight(unsigned eqIdx, Scalar value)
366 {
367 eqWeights_[eqIdx] = value;
368 }
369
374 std::size_t numAuxiliaryModules() const
375 {
376 return 0;
377 }
378
383 std::size_t numAuxiliaryDof() const
384 {
385 return 0;
386 }
387
397 const MaterialState& materialState(const unsigned globalIdx, unsigned /*timeIdx*/) const
398 {
399 return materialState_[globalIdx];
400 }
401
411 DimVector disp(const unsigned globalIdx, const bool /*with_fracture*/) const
412 {
413 DimVector d;
414 for (std::size_t i = 0; i < 3; ++i) {
415 d[i] = decay<Scalar>(materialState_[globalIdx].displacement(i));
416 }
417 return d;
418 }
419
428 SymTensor delstress(const unsigned /*globalIdx*/) const
429 {
430 SymTensor val;
431 return val;
432 }
433
442 SymTensor fractureStress(const unsigned /*globalIdx*/) const
443 {
444 SymTensor val;
445 return val;
446 }
447
456 SymTensor linstress(const unsigned /*globalIdx*/) const
457 {
458 SymTensor val;
459 return val;
460 }
461
471 SymTensor stress(const unsigned /*globalIdx*/, const bool /*with_fracture*/) const
472 {
473 SymTensor val;
474 return val;
475 }
476
486 SymTensor strain(const unsigned /*globalIdx*/, const bool /*with_fracture*/) const
487 {
488 SymTensor val;
489 return val;
490 }
491
500 Scalar mechPotentialForce(unsigned /*globalIdx*/) const
501 {
502 return Scalar(0.0);
503 }
504
513 Scalar mechPotentialPressForce(unsigned /*globalIdx*/) const
514 {
515 return Scalar(0.0);
516 }
517
526 Scalar mechPotentialTempForce(unsigned /*globalIdx*/) const
527 {
528 return Scalar(0.0);
529 }
530
531protected:
532 // ///
533 // Protected functions
534 // ///
541 {
542 const std::size_t numDof = simulator_.model().numGridDof();
543 materialState_.resize(numDof);
544 }
545
546private:
547 // ///
548 // Private functions
549 // ///
559 void setMaterialState_(const unsigned globalIdx, const unsigned /*timeIdx*/, PrimaryVariables& values)
560 {
561 auto& dofMaterialState = materialState_[globalIdx];
562 for (unsigned dirIdx = 0; dirIdx < 3; ++dirIdx) {
563 dofMaterialState.setDisplacement(dirIdx, values.makeEvaluation(disp0Idx + dirIdx, 0));
564 dofMaterialState.setRotation(dirIdx, values.makeEvaluation(rot0Idx + dirIdx, 0));
565 }
566 dofMaterialState.setSolidPressure(values.makeEvaluation(solidPres0Idx, 0));
567 }
568
569 std::unique_ptr<Linearizer> linearizer_;
570 NewtonMethod newtonMethod_;
571 Simulator& simulator_;
572 ElementChunks<GridView, Dune::Partitions::All> element_chunks_;
573
574 std::array<std::unique_ptr<TpsaBlockVectorWrapper>, historySize> solution_;
575 std::vector<Scalar> eqWeights_;
576 std::vector<MaterialState> materialState_;
577}; // class TpsaModel
578
579} // namespace Opm
580
581
582#endif
Data handle for parallel communication which can be used to set the values values of ghost and overla...
Definition: gridcommhandles.hh:109
static void registerParameters()
Register all run-time parameters for the Newton method.
Definition: newtonmethod.hh:135
Simplifies multi-threaded capabilities.
Definition: threadmanager.hpp:36
Small block vector wrapper class for model solutions.
Definition: tpsamodel.hpp:82
bool operator==(const TpsaBlockVectorWrapper &wrapper) const
Check if incoming block vector is the same as current.
Definition: tpsamodel.hpp:136
static TpsaBlockVectorWrapper serializationTestObject()
Test function for serialization.
Definition: tpsamodel.hpp:104
TpsaBlockVectorWrapper()=default
Default constructor.
TpsaBlockVectorWrapper(const std::string &, const std::size_t size)
Constructor.
Definition: tpsamodel.hpp:92
void serializeOp(Serializer &serializer)
Serializing operation.
Definition: tpsamodel.hpp:148
const SolutionVector & blockVector() const
Get const reference of block vector.
Definition: tpsamodel.hpp:127
SolutionVector & blockVector()
Get reference of block vector.
Definition: tpsamodel.hpp:119
SolutionVector blockVector_
Definition: tpsamodel.hpp:84
TPSA geomechanics model.
Definition: tpsamodel.hpp:51
NewtonMethod & newtonMethod()
Return the Newton method.
Definition: tpsamodel.hpp:291
std::size_t numAuxiliaryDof() const
Return number of auxillary degrees of freedom.
Definition: tpsamodel.hpp:383
Scalar mechPotentialPressForce(unsigned) const
Output potential pressure forces.
Definition: tpsamodel.hpp:513
Scalar mechPotentialForce(unsigned) const
Output potential forces.
Definition: tpsamodel.hpp:500
DimVector disp(const unsigned globalIdx, const bool) const
Output displacement vector.
Definition: tpsamodel.hpp:411
void finishInit()
Initialize TPSA model.
Definition: tpsamodel.hpp:182
void setEqWeight(unsigned eqIdx, Scalar value)
Set weights for equation.
Definition: tpsamodel.hpp:365
SymTensor linstress(const unsigned) const
Output linear stress tensor.
Definition: tpsamodel.hpp:456
std::size_t numGridDof() const
Return number of degrees of freedom in the grid from the Flow model.
Definition: tpsamodel.hpp:322
const Linearizer & linearizer() const
Return the linearizer.
Definition: tpsamodel.hpp:261
const NewtonMethod & newtonMethod() const
Return the Newton method.
Definition: tpsamodel.hpp:281
std::size_t numAuxiliaryModules() const
Return number of auxillary modules.
Definition: tpsamodel.hpp:374
void updateMaterialState(const unsigned)
Update material state for all cells.
Definition: tpsamodel.hpp:237
Scalar mechPotentialTempForce(unsigned) const
Output potential temparature forces.
Definition: tpsamodel.hpp:526
const SolutionVector & solution(unsigned timeIdx) const
Get reference to history solution vector.
Definition: tpsamodel.hpp:302
const MaterialState & materialState(const unsigned globalIdx, unsigned) const
Return current material state.
Definition: tpsamodel.hpp:397
Scalar dofTotalVolume(unsigned globalIdx) const
Return the total grid volume from the Flow model.
Definition: tpsamodel.hpp:342
void resizeMaterialState_()
Resize material state vector.
Definition: tpsamodel.hpp:540
static void registerParameters()
Register runtime parameters.
Definition: tpsamodel.hpp:194
SymTensor delstress(const unsigned) const
Output (del?)stress tensor.
Definition: tpsamodel.hpp:428
Scalar eqWeight(unsigned, unsigned eqIdx) const
Return equation weights.
Definition: tpsamodel.hpp:354
TpsaModel(Simulator &simulator)
Constructor.
Definition: tpsamodel.hpp:162
SymTensor strain(const unsigned, const bool) const
Output strain tensor.
Definition: tpsamodel.hpp:486
SolutionVector & solution(unsigned timeIdx)
Get reference to history solution vector.
Definition: tpsamodel.hpp:313
void prepareTPSA()
Prepare TPSA model for coupled Flow-TPSA scheme.
Definition: tpsamodel.hpp:203
SymTensor fractureStress(const unsigned) const
Output fracture stress tensor.
Definition: tpsamodel.hpp:442
std::size_t numTotalDof() const
Return the total number of degrees of freedom.
Definition: tpsamodel.hpp:331
SymTensor stress(const unsigned, const bool) const
Output stress tensor.
Definition: tpsamodel.hpp:471
Linearizer & linearizer()
Return the linearizer.
Definition: tpsamodel.hpp:271
void syncOverlap()
Sync primary variables in overlapping cells.
Definition: tpsamodel.hpp:214
Definition: fvbaseprimaryvariables.hh:161
Definition: blackoilbioeffectsmodules.hh:45
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:233
The Opm property system, traits with inheritance.