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