opm-simulators
Loading...
Searching...
No Matches
FlowProblemTPSA.hpp
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 FLOW_PROBLEM_TPSA_HPP
26#define FLOW_PROBLEM_TPSA_HPP
27
28#include <dune/common/fvector.hh>
29
30#include <dune/grid/common/rangegenerators.hh>
31
32#include <opm/common/OpmLog/OpmLog.hpp>
33
34#include <opm/input/eclipse/EclipseState/Grid/FaceDir.hpp>
35
36#include <opm/material/common/MathToolbox.hpp>
37#include <opm/material/materialstates/MaterialStateTPSA.hpp>
38
39#include <opm/models/io/vtktpsamodule.hpp>
40#include <opm/models/tpsa/tpsabaseproperties.hpp>
41#include <opm/models/tpsa/tpsamodel.hpp>
43
44#include <opm/simulators/flow/FacePropertiesTPSA.hpp>
46
47#include <cmath>
48#include <memory>
49#include <stdexcept>
50#include <string>
51#include <utility>
52
53#include <fmt/format.h>
54
55
56namespace Opm {
57
61template <class TypeTag>
62class FlowProblemTPSA : public FlowProblemBlackoil<TypeTag>
63{
64public:
65 using ParentType = FlowProblemBlackoil<TypeTag>;
66
78
79 enum { dimWorld = GridView::dimensionworld };
83 enum { numPhases = FluidSystem::numPhases };
84
85 enum { contiRotEqIdx = Indices::contiRotEqIdx };
86 enum { contiSolidPresEqIdx = Indices::contiSolidPresEqIdx };
87 enum { solidPres0Idx = Indices::solidPres0Idx };
88
89 using CartesianIndexMapper = Dune::CartesianIndexMapper<Grid>;
90 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
92 using InitialMaterialState = MaterialStateTPSA<Scalar>;
93 using Toolbox = MathToolbox<Evaluation>;
94
95 // ///
96 // Public functions
97 // ///
103 FlowProblemTPSA(Simulator& simulator)
104 : ParentType(simulator)
105 , faceProps_(simulator.vanguard().eclState(),
106 simulator.vanguard().gridView(),
107 simulator.vanguard().cartesianIndexMapper(),
108 simulator.vanguard().grid(),
109 simulator.vanguard().cellCentroids())
110 , geoMechModel_(simulator)
111 {
112 if constexpr(enableMech) {
113 // Add VTK TPSA to output module
114 this->model().addOutputModule(std::make_unique<VtkTpsaModule<TypeTag>>(simulator));
115
116 // Sanity check
117 const auto& mechSolver = simulator.vanguard().eclState().runspec().mechSolver();
118 if (!mechSolver.tpsa()) {
119 std::string msg = "Simulator with Tpsa-geomechanics enabled compile time, but deck does not contain "
120 "TPSA keyword!";
121 OpmLog::error(msg);
122 throw std::runtime_error(msg);
123 }
124 }
125 else {
126 // Sanity check
127 const auto& mechSolver = simulator.vanguard().eclState().runspec().mechSolver();
128 if (mechSolver.tpsa()) {
129 std::string msg = "TPSA keyword in deck, but Tpsa-geomechanics disabled compile-time!";
130 OpmLog::error(msg);
131 throw std::runtime_error(msg);
132 }
133 }
134 }
135
139 static void registerParameters()
140 {
141 // Register parameters for parent class
143
144 // Geomech model parameters
145 GeomechModel::registerParameters();
146
147 // VTK output parameters
149 }
150
155 {
156 // FlowProblemBlackoil::finishInit()
158
159 // Read initial conditions and set material state
161
162 // Calculate face properties
163 faceProps_.finishInit();
164
165 // Set equation weights
167
168 // Initialize the TPSA model
169 geoMechModel_.finishInit();
170 }
171
178 {
179 // Set up initial solution for the Flow model
181
182 // Initialize soultions as zero
183 auto& uCur = geoMechModel_.solution(/*timeIdx=*/0);
184 uCur = Scalar(0.0);
185
186 // Loop through grid and set initial solution from material state
187 ElementContext elemCtx(this->simulator());
188 for (const auto& elem : elements(this->gridView())) {
189 // Ignore everything which is not in the interior if the current process' piece of the grid
190 if (elem.partitionType() != Dune::InteriorEntity) {
191 continue;
192 }
193
194 // Loop over sub control volumes and set initial solutions
195 elemCtx.updateStencil(elem);
196 for (unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(/*timeIdx=*/0); ++dofIdx) {
197 const unsigned globalIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
198 auto& priVars = uCur[globalIdx];
199 priVars.assignNaive(initialMaterialState_[globalIdx]);
200 priVars.checkDefined();
201 }
202 }
203
204 // synchronize the ghost/overlapping DOFs (if necessary)
205 geoMechModel_.syncOverlap();
206
207 // Set history solutions to the initial solution.
208 for (unsigned timeIdx = 1; timeIdx < historySize; ++timeIdx) {
209 geoMechModel_.solution(timeIdx) = geoMechModel_.solution(/*timeIdx=*/0);
210 }
211
212 // Set material state
213 geoMechModel_.updateMaterialState(/*timeIdx=*/0);
214 }
215
220 {
221 // Average shear modulus over domain
222 Scalar avgSmodulus = 0.0;
223 const auto& gridView = this->gridView();
224 ElementContext elemCtx(this->simulator());
225 for(const auto& elem: elements(gridView, Dune::Partitions::interior)) {
226 elemCtx.updatePrimaryStencil(elem);
227 int elemIdx = elemCtx.globalSpaceIndex(/*spaceIdx=*/0, /*timeIdx=*/0);
228 avgSmodulus += this->shearModulus(elemIdx);
229 }
230 std::size_t numDof = this->model().numGridDof();
231 const auto& comm = this->simulator().vanguard().grid().comm();
232 avgSmodulus = comm.sum(avgSmodulus);
233 Scalar numTotalDof = comm.sum(numDof);
234 avgSmodulus /= numTotalDof;
235 avgSmodulus = std::sqrt(avgSmodulus);
236
237 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
238 if (eqIdx < contiRotEqIdx) {
239 geoMechModel_.setEqWeight(eqIdx, 1 / avgSmodulus);
240 }
241 else {
242 geoMechModel_.setEqWeight(eqIdx, avgSmodulus);
243 }
244 }
245 }
246
250 void beginTimeStep() override
251 {
252 // Call parent class beginTimeStep()
254
255 // Update mechanics boundary conditions.
256 // NOTE: Flow boundary conditions should be updated in ParentType::beginTimeStep()
257 if (this->nonTrivialBoundaryConditions()) {
258 geoMechModel_.linearizer().updateBoundaryConditionData();
259 }
260 }
261
274 std::pair<BCMECHType, Dune::FieldVector<Evaluation, 3>>
275 mechBoundaryCondition(const unsigned int globalSpaceIdx, const int directionId)
276 {
277 // Default boundary conditions if BCCON/BCPROP not defined
278 if (!this->nonTrivialBoundaryConditions()) {
279 return { BCMECHType::NONE, Dune::FieldVector<Evaluation, 3>{0.0, 0.0, 0.0} };
280 }
281
282 // Default for BCPROP index = 0 or no BCPROP defined at current episode
283 FaceDir::DirEnum dir = FaceDir::FromIntersectionIndex(directionId);
284 const auto& schedule = this->simulator().vanguard().schedule();
285 if (this->bcindex_(dir)[globalSpaceIdx] == 0 || schedule[this->episodeIndex()].bcprop.size() == 0) {
286 return { BCMECHType::NONE, Dune::FieldVector<Evaluation, 3>{0.0, 0.0, 0.0} };
287 }
288
289 // Get current BC
290 const auto& bc = schedule[this->episodeIndex()].bcprop[this->bcindex_(dir)[globalSpaceIdx]];
291 if (bc.bcmechtype == BCMECHType::FREE) {
292 return { BCMECHType::FREE, Dune::FieldVector<Evaluation, 3>{0.0, 0.0, 0.0} };
293 }
294 else {
295 return { bc.bcmechtype, Dune::FieldVector<Evaluation, 3>{0.0, 0.0, 0.0} };
296 }
297 }
298
308 void tpsaSource(Dune::FieldVector<Evaluation, numEq>& sourceTerm,
309 unsigned globalSpaceIdx,
310 unsigned timeIdx)
311 {
312 sourceTerm = 0.0;
313
314 // Coupling term Flow -> TPSA
315 const auto biot = this->biotCoeff(globalSpaceIdx);
316 const auto lameParam = this->lame(globalSpaceIdx);
317
318 const auto& iq = this->model().intensiveQuantities(globalSpaceIdx, timeIdx);
319 const auto& fs = iq.fluidState();
320 const auto pres = decay<Scalar>(fs.pressure(this->refPressurePhaseIdx_()));
321 const auto initPres = this->initialFluidState(globalSpaceIdx).pressure(this->refPressurePhaseIdx_());
322
323 auto sourceFromFlow = -biot / lameParam * (pres - initPres);
324 sourceTerm[contiSolidPresEqIdx] += sourceFromFlow;
325 }
326
336 Scalar rockMechPoroChange(unsigned elementIdx, unsigned timeIdx) const
337 {
338 // TODO: get timeIdx=1 solid pressure from a cached materialState (or intensiveQuantities) if/when implemented
339 assert (timeIdx <= historySize);
340 const auto solidPres = (timeIdx == 0) ?
341 decay<Scalar>( geoMechModel_.materialState(elementIdx, /*timeIdx=*/timeIdx).solidPressure()) :
342 geoMechModel_.solution(/*timeIdx=*/timeIdx)[elementIdx][solidPres0Idx];
343 const auto biot = this->biotCoeff(elementIdx);
344 const auto lameParam = this->lame(elementIdx);
345
346 return biot / lameParam * solidPres;
347 }
348
349 // ///
350 // Public get functions
351 // ///
359 Scalar weightAverage(unsigned globalElemIdxIn, unsigned globalElemIdxOut)
360 {
361 return faceProps_.weightAverage(globalElemIdxIn, globalElemIdxOut);
362 }
363
371 Scalar weightAverageBoundary(unsigned globalElemIdxIn, unsigned boundaryFaceIdx) const
372 {
373 return faceProps_.weightAverageBoundary(globalElemIdxIn, boundaryFaceIdx);
374 }
375
383 Scalar weightProduct(unsigned globalElemIdxIn, unsigned globalElemIdxOut) const
384 {
385 return faceProps_.weightProduct(globalElemIdxIn, globalElemIdxOut);
386 }
387
395 Scalar normalDistance(unsigned globalElemIdxIn, unsigned globalElemIdxOut) const
396 {
397 return faceProps_.normalDistance(globalElemIdxIn, globalElemIdxOut);
398 }
399
407 Scalar normalDistanceBoundary(unsigned globalElemIdxIn, unsigned boundaryFaceIdx) const
408 {
409 return faceProps_.normalDistanceBoundary(globalElemIdxIn, boundaryFaceIdx);
410 }
411
419 DimVector cellFaceNormal(unsigned globalElemIdxIn, unsigned globalElemIdxOut)
420 {
421 return faceProps_.cellFaceNormal(globalElemIdxIn, globalElemIdxOut);
422 }
423
431 const DimVector& cellFaceNormalBoundary(unsigned globalElemIdxIn, unsigned boundaryFaceIdx) const
432 {
433 return faceProps_.cellFaceNormalBoundary(globalElemIdxIn, boundaryFaceIdx);
434 }
435
442 Scalar shearModulus(unsigned globalElemIdx) const
443 {
444 return faceProps_.shearModulus(globalElemIdx);
445 }
446
452 bool laggedScheme() const
453 {
454 const auto& mechSolver = this->simulator().vanguard().eclState().runspec().mechSolver();
455 return mechSolver.laggedScheme();
456 }
457
463 bool fixedStressScheme() const
464 {
465 const auto& mechSolver = this->simulator().vanguard().eclState().runspec().mechSolver();
466 return mechSolver.fixedStressScheme();
467 }
468
474 const GeomechModel& geoMechModel() const
475 {
476 return geoMechModel_;
477 }
478
484 GeomechModel& geoMechModel()
485 {
486 return geoMechModel_;
487 }
488
494 std::pair<int, int> fixedStressParameters() const
495 {
496 const auto& mechSolver = this->simulator().vanguard().eclState().runspec().mechSolver();
497 return std::make_pair(mechSolver.fixedStressMinIter(), mechSolver.fixedStressMaxIter());
498 }
499
500protected:
501 // ///
502 // Protected functions
503 // ///
508 {
509 // ///
510 // OBS: No equilibration keywords (e.g., STREQUIL) implemented yet!
511 // ///
512
513 // Set all initial material state variables to zero
514 std::size_t numDof = this->model().numGridDof();
515 initialMaterialState_.resize(numDof);
516 for (std::size_t dofIdx = 0; dofIdx < numDof; ++dofIdx) {
517 auto& dofMaterialState = initialMaterialState_[dofIdx];
518 for (unsigned dirIdx = 0; dirIdx < 3; ++dirIdx) {
519 dofMaterialState.setDisplacement(dirIdx, 0.0);
520 dofMaterialState.setRotation(dirIdx, 0.0);
521 }
522 dofMaterialState.setSolidPressure(0.0);
523 }
524 }
525
526private:
527 FaceProperties faceProps_;
528 GeomechModel geoMechModel_;
529
530 std::vector<Scalar> biotcoeff_;
531 std::vector<InitialMaterialState> initialMaterialState_;
532}; // class FlowProblemTPSA
533
534} // namespace Opm
535
536#endif
This problem simulates an input file given in the data format used by the commercial ECLiPSE simulato...
Definition CollectDataOnIORank.hpp:49
Cell face properties needed in TPSA equation calculations.
Definition FacePropertiesTPSA.hpp:48
FlowProblemBlackoil(Simulator &simulator)
Definition FlowProblemBlackoil.hpp:187
void finishInit()
Called by the Opm::Simulator in order to initialize the problem.
Definition FlowProblemBlackoil.hpp:303
void initialSolutionApplied() override
Callback used by the model to indicate that the initial solution has been determined for all degrees ...
Definition FlowProblemBlackoil.hpp:614
void beginTimeStep() override
Called by the simulator before each time integration.
Definition FlowProblemBlackoil.hpp:294
static void registerParameters()
Registers all available parameters for the problem and the model.
Definition FlowProblemBlackoil.hpp:173
Scalar shearModulus(unsigned globalElemIdx) const
Direct access to shear modulus in an element.
Definition FlowProblemTPSA.hpp:442
Scalar weightAverageBoundary(unsigned globalElemIdxIn, unsigned boundaryFaceIdx) const
Direct access to normal distance at the boundary.
Definition FlowProblemTPSA.hpp:371
static void registerParameters()
Register runtime parameters.
Definition FlowProblemTPSA.hpp:139
const DimVector & cellFaceNormalBoundary(unsigned globalElemIdxIn, unsigned boundaryFaceIdx) const
Direct access to face normal at the boundary.
Definition FlowProblemTPSA.hpp:431
std::pair< BCMECHType, Dune::FieldVector< Evaluation, 3 > > mechBoundaryCondition(const unsigned int globalSpaceIdx, const int directionId)
Organize mechanics boundary conditions.
Definition FlowProblemTPSA.hpp:275
void initialSolutionApplied() override
Set initial solution for the problem.
Definition FlowProblemTPSA.hpp:177
std::pair< int, int > fixedStressParameters() const
Get fixed-stress iteration parameters.
Definition FlowProblemTPSA.hpp:494
void beginTimeStep() override
Called by the simulator before each time integration.
Definition FlowProblemTPSA.hpp:250
DimVector cellFaceNormal(unsigned globalElemIdxIn, unsigned globalElemIdxOut)
Direct access to face normal between two elements.
Definition FlowProblemTPSA.hpp:419
bool fixedStressScheme() const
Flow-TPSA fixed-stress coupling scheme activated?
Definition FlowProblemTPSA.hpp:463
Scalar normalDistance(unsigned globalElemIdxIn, unsigned globalElemIdxOut) const
Direct access to normal distance between two elements.
Definition FlowProblemTPSA.hpp:395
void computeAndSetEqWeights_()
Compute weights to rescale the TPSA equations.
Definition FlowProblemTPSA.hpp:219
const GeomechModel & geoMechModel() const
Get TPSA model.
Definition FlowProblemTPSA.hpp:474
Scalar weightAverage(unsigned globalElemIdxIn, unsigned globalElemIdxOut)
Direct access to average (half-)weight at interface between two elements.
Definition FlowProblemTPSA.hpp:359
void tpsaSource(Dune::FieldVector< Evaluation, numEq > &sourceTerm, unsigned globalSpaceIdx, unsigned timeIdx)
Set mechanics source term, in particular coupling terms.
Definition FlowProblemTPSA.hpp:308
GeomechModel & geoMechModel()
Get TPSA model.
Definition FlowProblemTPSA.hpp:484
FlowProblemTPSA(Simulator &simulator)
Constructor.
Definition FlowProblemTPSA.hpp:103
Scalar rockMechPoroChange(unsigned elementIdx, unsigned timeIdx) const
Pore volume change due to geomechanics.
Definition FlowProblemTPSA.hpp:336
void finishInit()
Initialize the problem.
Definition FlowProblemTPSA.hpp:154
bool laggedScheme() const
Flow-TPSA lagged coupling scheme activated?
Definition FlowProblemTPSA.hpp:452
Scalar normalDistanceBoundary(unsigned globalElemIdxIn, unsigned boundaryFaceIdx) const
Direct access to normal distance at the boundary.
Definition FlowProblemTPSA.hpp:407
Scalar weightProduct(unsigned globalElemIdxIn, unsigned globalElemIdxOut) const
Direct access to product of weights at interface between two elements.
Definition FlowProblemTPSA.hpp:383
void readInitalConditionsTPSA_()
Read initial conditions and generate material state for TPSA model.
Definition FlowProblemTPSA.hpp:507
Scalar lame(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition FlowProblem.hpp:735
Scalar biotCoeff(const Context &context, unsigned spaceIdx, unsigned timeIdx) const
Definition FlowProblem.hpp:745
VTK output module for TPSA quantities.
Definition vtktpsamodule.hpp:47
static void registerParameters()
Register runtime parameters.
Definition vtktpsamodule.hpp:81
This file contains a set of helper functions used by VFPProd / VFPInj.
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
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240
This file provides the infrastructure to retrieve run-time parameters.