elasticitylocalresidualtpsa.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 ELASTICITY_LOCAL_RESIDUAL_TPSA_HPP
26#define ELASTICITY_LOCAL_RESIDUAL_TPSA_HPP
27
28#include <dune/common/fvector.hh>
29
30#include <opm/input/eclipse/Schedule/BCProp.hpp>
31
32#include <opm/material/common/MathToolbox.hpp>
33#include <opm/material/materialstates/MaterialStateTPSA.hpp>
34
38
40
41#include <stdexcept>
42
43
44namespace Opm {
45
67template <class TypeTag>
69{
74
75 using MaterialState = MaterialStateTPSA<Evaluation>;
76 using Toolbox = MathToolbox<Evaluation>;
77
78 enum { conti0EqIdx = Indices::conti0EqIdx };
79 enum { contiRotEqIdx = Indices::contiRotEqIdx };
80 enum { contiSolidPresEqIdx = Indices::contiSolidPresEqIdx };
81 enum { numEq = getPropValue<TypeTag, Properties::NumEqTPSA>() };
82
83public:
95 template <class LhsEval>
96 static void computeVolumeTerm(Dune::FieldVector<LhsEval, numEq>& volTerm,
97 const MaterialState& materialState,
98 const Problem& problem,
99 const unsigned globalIndex)
100 {
101 // Reset volume terms
102 volTerm = 0.0;
103
104 // Rotation equations (one per direction)
105 const Scalar sModulus = problem.shearModulus(globalIndex);
106 for (unsigned dirIdx = 0; dirIdx < 3; ++dirIdx) {
107 volTerm[contiRotEqIdx + dirIdx] +=
108 Toolbox::template decay<LhsEval>(materialState.rotation(dirIdx))
109 / sModulus;
110 }
111
112 // Solid pressure equation
113 const Scalar lame = problem.lame(globalIndex);
114 volTerm[contiSolidPresEqIdx] +=
115 Toolbox::template decay<LhsEval>(materialState.solidPressure())
116 / lame;
117 }
118
132 static void computeFaceTerm(Dune::FieldVector<Evaluation, numEq>& faceTerm,
133 const MaterialState& materialStateIn,
134 const MaterialState& materialStateEx,
135 Problem& problem,
136 const unsigned globalIndexIn,
137 const unsigned globalIndexEx)
138 {
139 // Reset face terms
140 faceTerm = 0.0;
141
142 // Extract some face properties
143 const Scalar weightAvgIn = problem.weightAverage(globalIndexIn, globalIndexEx);
144 const Scalar weightAvgEx = problem.weightAverage(globalIndexEx, globalIndexIn);
145 const Scalar weightProd = problem.weightProduct(globalIndexIn, globalIndexEx);
146 const Scalar normDist = problem.normalDistance(globalIndexIn, globalIndexEx);
147 const auto& faceNormal = problem.cellFaceNormal(globalIndexIn, globalIndexEx);
148
149 // Effective shear modulus
150 const Scalar sModulusIn = problem.shearModulus(globalIndexIn);
151 const Scalar sModulusEx = problem.shearModulus(globalIndexEx);
152 const Scalar eff_sModulus = weightAvgIn * sModulusIn + weightAvgEx * sModulusEx;
153
154 // Distance ratio
155 const Scalar distRatio = 0.5 * weightProd / normDist;
156
157 // Solid pressures
158 const Evaluation& solidPIn = materialStateIn.solidPressure();
159 const Scalar solidPEx = decay<Scalar>(materialStateEx.solidPressure());
160
161 // ///
162 // Solid pressure equation (direction-independent equation)
163 // ///
164 faceTerm[contiSolidPresEqIdx] +=
165 distRatio * eff_sModulus * (solidPIn - solidPEx);
166
167 // ///
168 // Displacement, rotation and solid pressure (directional-dependent) equations
169 // ///
170 // Lambda function for computing modulo 3 of possibly negative integers to get indices in cross product.
171 // E.g. if i = x-dir(=0), we want y-dir(=1) and z-dir(=2), hence -1 mod 3 must equal 2 and not -1
172 auto modNeg = [](int i) { return ((i % 3) + 3) % 3; };
173
174 // Loop over x-, y- and z-dir (corresponding to dirIdx = 0, 1, 2)
175 for (int dirIdx = 0; dirIdx < 3; ++dirIdx) {
176 // Direction indices in cross-product
177 unsigned dirIdxNeg = modNeg(dirIdx - 1);
178 unsigned dirIdxPos = modNeg(dirIdx + 1);
179
180 // Displacement equation
181 const Scalar faceNormalDir = faceNormal[dirIdx];
182 const Scalar faceNormalNeg = faceNormal[dirIdxNeg];
183 const Scalar faceNormalPos = faceNormal[dirIdxPos];
184
185 const Evaluation& dispIn = materialStateIn.displacement(dirIdx);
186 const Scalar dispEx = decay<Scalar>(materialStateEx.displacement(dirIdx));
187
188 const Evaluation& rotInNeg = materialStateIn.rotation(dirIdxNeg);
189 const Evaluation& rotInPos = materialStateIn.rotation(dirIdxPos);
190 const Scalar rotExNeg = decay<Scalar>(materialStateEx.rotation(dirIdxNeg));
191 const Scalar rotExPos = decay<Scalar>(materialStateEx.rotation(dirIdxPos));
192
193 faceTerm[conti0EqIdx + dirIdx] +=
194 2.0 * (eff_sModulus / normDist) * (dispIn - dispEx)
195 - weightAvgIn * (faceNormalNeg * rotInPos - faceNormalPos * rotInNeg)
196 - weightAvgEx * (faceNormalNeg * rotExPos - faceNormalPos * rotExNeg)
197 - faceNormalDir * (weightAvgIn * solidPIn + weightAvgEx * solidPEx);
198
199 // Rotation equation
200 const Evaluation& dispInNeg = materialStateIn.displacement(dirIdxNeg);
201 const Evaluation& dispInPos = materialStateIn.displacement(dirIdxPos);
202 const Scalar dispExNeg = decay<Scalar>(materialStateEx.displacement(dirIdxNeg));
203 const Scalar dispExPos = decay<Scalar>(materialStateEx.displacement(dirIdxPos));
204
205 faceTerm[contiRotEqIdx + dirIdx] +=
206 - weightAvgEx * (faceNormalNeg * dispInPos - faceNormalPos * dispInNeg)
207 - weightAvgIn * (faceNormalNeg * dispExPos - faceNormalPos * dispExNeg);
208
209 // Solid pressure equation
210 faceTerm[contiSolidPresEqIdx] +=
211 - faceNormalDir * (weightAvgEx * dispIn + weightAvgIn * dispEx);
212 }
213 }
214
224 template <class BoundaryConditionData>
225 static void computeBoundaryTerm(Dune::FieldVector<Evaluation, numEq>& bndryTerm,
226 const MaterialState& materialState,
227 const BoundaryConditionData& bdyInfo,
228 Problem& problem,
229 unsigned globalIndex)
230 {
231 // Switch between possible boundary conditions
232 switch (bdyInfo.type) {
233 // OBS: NONE is interpreted as FIXED with zero displacement
234 case BCMECHType::NONE:
235 computeBoundaryTermFixed(bndryTerm,
236 materialState,
237 bdyInfo,
238 problem,
239 globalIndex);
240 break;
241 case BCMECHType::FIXED:
242 throw std::runtime_error("BCTYPE FIXED has not been implemented in TPSA");
243 case BCMECHType::FREE:
244 computeBoundaryTermFree(bndryTerm,
245 materialState,
246 bdyInfo,
247 problem,
248 globalIndex);
249 break;
250 default:
251 throw std::logic_error("Unknown boundary condition type " +
252 std::to_string(static_cast<int>(bdyInfo.type)) +
253 " in computeBoundaryFlux()." );
254 }
255 }
256
269 template <class BoundaryConditionData>
270 static void computeBoundaryTermFixed(Dune::FieldVector<Evaluation, numEq>& bndryTerm,
271 const MaterialState& materialState,
272 const BoundaryConditionData& bdyInfo,
273 Problem& problem,
274 unsigned globalIndex)
275 {
276 // !!!!
277 // Only BCMECHType::NONE, where we have zero displacement on boundary face, have been implemented!
278 // !!!!
279
280 // Reset bondary term
281 bndryTerm = 0.0;
282
283 // Extract some face properties
284 const unsigned bfIdx = bdyInfo.boundaryFaceIndex;
285 const Scalar weightAvg = problem.weightAverageBoundary(globalIndex, bfIdx);
286 const Scalar normDist = problem.normalDistanceBoundary(globalIndex, bfIdx);
287 const auto& faceNormal = problem.cellFaceNormalBoundary(globalIndex, bfIdx);
288
289 // Effective shear modulus (= cell shear modulus)
290 const Scalar eff_sModulus = problem.shearModulus(globalIndex);
291
292 // Solid pressure
293 const Evaluation& solidP = materialState.solidPressure();
294
295 // ///
296 // Displacement equation
297 // ///
298 // Lambda function for computing modulo 3 of possibly negative integers to get indices in cross product.
299 // E.g. if i = x-dir(=0), we want y-dir(=1) and z-dir(=2), hence -1 mod 3 must equal 2 and not -1
300 auto modNeg = [](int i) { return ((i % 3) + 3) % 3; };
301
302 // Loop over x-, y- and z-dir (corresponding to dirIdx = 0, 1, 2)
303 for (int dirIdx = 0; dirIdx < 3; ++dirIdx) {
304 // Direction indices in cross-product
305 unsigned dirIdxNeg = modNeg(dirIdx - 1);
306 unsigned dirIdxPos = modNeg(dirIdx + 1);
307
308 // Displacement equation
309 const Scalar faceNormalDir = faceNormal[dirIdx];
310 const Scalar faceNormalNeg = faceNormal[dirIdxNeg];
311 const Scalar faceNormalPos = faceNormal[dirIdxPos];
312
313 const Evaluation& disp = materialState.displacement(dirIdx);
314
315 const Evaluation& rotNeg = materialState.rotation(dirIdxNeg);
316 const Evaluation& rotPos = materialState.rotation(dirIdxPos);
317
318 bndryTerm[conti0EqIdx + dirIdx] +=
319 2.0 * (eff_sModulus / normDist) * disp
320 - weightAvg * (faceNormalNeg * rotPos - faceNormalPos * rotNeg)
321 - faceNormalDir * weightAvg * solidP;
322 }
323 }
324
337 template <class BoundaryConditionData>
338 static void computeBoundaryTermFree(Dune::FieldVector<Evaluation, numEq>& bndryTerm,
339 const MaterialState& materialState,
340 const BoundaryConditionData& bdyInfo,
341 Problem& problem,
342 unsigned globalIndex)
343 {
344 // Reset bondary term
345 bndryTerm = 0.0;
346
347 // Face properties
348 const unsigned bfIdx = bdyInfo.boundaryFaceIndex;
349 const Scalar weightAvg = 1.0;
350 const Scalar normDist = problem.normalDistanceBoundary(globalIndex, bfIdx);
351 const auto& faceNormal = problem.cellFaceNormalBoundary(globalIndex, bfIdx);
352
353 const Scalar sModulus = problem.shearModulus(globalIndex);
354
355 // ///
356 // Solid pressure equation (direction-independent equation)
357 // ///
358 const Evaluation& solidP = materialState.solidPressure();
359 bndryTerm[contiSolidPresEqIdx] +=
360 0.5 * (normDist / sModulus) * solidP;
361
362 // ///
363 // Rotation and solid pressure (directional-dependent) equations
364 // ///
365 // Lambda function for computing modulo 3 of possibly negative integers to get indices in cross product.
366 // E.g. if i = x-dir(=0), we want y-dir(=1) and z-dir(=2), hence -1 mod 3 must equal 2 and not -1
367 auto modNeg = [](int i) { return ((i % 3) + 3) % 3; };
368
369 // Pre-compute dot product for rotation equation
370 Evaluation dotProd = 0;
371 for (int dirIdx = 0; dirIdx < 3; ++dirIdx) {
372 dotProd += faceNormal[dirIdx] * materialState.rotation(dirIdx);
373 }
374
375 // Loop over x-, y- and z-dir (corresponding to dirIdx = 0, 1, 2)
376 for (int dirIdx = 0; dirIdx < 3; ++dirIdx) {
377 // Direction indices in cross-product
378 unsigned dirIdxNeg = modNeg(dirIdx - 1);
379 unsigned dirIdxPos = modNeg(dirIdx + 1);
380
381 // Rotation equation
382 const Scalar faceNormalDir = faceNormal[dirIdx];
383 const Scalar faceNormalNeg = faceNormal[dirIdxNeg];
384 const Scalar faceNormalPos = faceNormal[dirIdxPos];
385
386 const Evaluation& dispNeg = materialState.displacement(dirIdxNeg);
387 const Evaluation& dispPos = materialState.displacement(dirIdxPos);
388
389 const Evaluation& rot = materialState.rotation(dirIdx);
390
391 bndryTerm[contiRotEqIdx + dirIdx] +=
392 - weightAvg * (faceNormalNeg * dispPos - faceNormalPos * dispNeg)
393 + 0.5 * (normDist / sModulus) * (dotProd * faceNormalDir - rot);
394
395 // Solid pressure (directional-dependent) equation
396 const Evaluation& disp = materialState.displacement(dirIdx);
397
398 bndryTerm[contiSolidPresEqIdx] +=
399 - faceNormalDir * weightAvg * disp;
400 }
401 }
402
411 static void computeSourceTerm(Dune::FieldVector<Evaluation, numEq>& sourceTerm,
412 Problem& problem,
413 unsigned globalSpaceIdex,
414 unsigned timeIdx)
415 {
416 // Reset source terms
417 sourceTerm = 0.0;
418
419 // Get source term from problem
420 // NOTE: separate source function than Flow source(...)!
421 problem.tpsaSource(sourceTerm, globalSpaceIdex, timeIdx);
422 }
423}; // class ElasticityLocalResidual
424
425} // namespace Opm
426
427#endif
Defines a type tags and some fundamental properties all models.
Calculation of (linear) elasticity model terms for the residual.
Definition: elasticitylocalresidualtpsa.hpp:69
static void computeFaceTerm(Dune::FieldVector< Evaluation, numEq > &faceTerm, const MaterialState &materialStateIn, const MaterialState &materialStateEx, Problem &problem, const unsigned globalIndexIn, const unsigned globalIndexEx)
Calculate terms across cell faces in TPSA formulation.
Definition: elasticitylocalresidualtpsa.hpp:132
static void computeBoundaryTermFixed(Dune::FieldVector< Evaluation, numEq > &bndryTerm, const MaterialState &materialState, const BoundaryConditionData &bdyInfo, Problem &problem, unsigned globalIndex)
Calculate fixed displacement boundary condition in TPSA formulation.
Definition: elasticitylocalresidualtpsa.hpp:270
static void computeBoundaryTermFree(Dune::FieldVector< Evaluation, numEq > &bndryTerm, const MaterialState &materialState, const BoundaryConditionData &bdyInfo, Problem &problem, unsigned globalIndex)
Calculate free (or zero traction) boundary condition in TPSA formulation.
Definition: elasticitylocalresidualtpsa.hpp:338
static void computeVolumeTerm(Dune::FieldVector< LhsEval, numEq > &volTerm, const MaterialState &materialState, const Problem &problem, const unsigned globalIndex)
Calculate volume terms in TPSA formulation.
Definition: elasticitylocalresidualtpsa.hpp:96
static void computeSourceTerm(Dune::FieldVector< Evaluation, numEq > &sourceTerm, Problem &problem, unsigned globalSpaceIdex, unsigned timeIdx)
Calculate source term in TPSA formulation.
Definition: elasticitylocalresidualtpsa.hpp:411
static void computeBoundaryTerm(Dune::FieldVector< Evaluation, numEq > &bndryTerm, const MaterialState &materialState, const BoundaryConditionData &bdyInfo, Problem &problem, unsigned globalIndex)
Calculate boundary conditions in TPSA formulation given by BCCON/BCPROP.
Definition: elasticitylocalresidualtpsa.hpp:225
Declare the properties used by the infrastructure code of the finite volume discretizations.
@ NONE
Definition: DeferredLogger.hpp:46
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
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)