Transmissibility_impl.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 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
23#ifndef OPM_TRANSMISSIBILITY_IMPL_HPP
24#define OPM_TRANSMISSIBILITY_IMPL_HPP
25
26#ifndef OPM_TRANSMISSIBILITY_HPP
27#include <config.h>
29#endif
30
31#include <dune/common/version.hh>
32#include <dune/grid/common/mcmgmapper.hh>
33
34#include <opm/common/OpmLog/KeywordLocation.hpp>
35#include <opm/common/utility/ThreadSafeMapBuilder.hpp>
36
37#include <opm/grid/CpGrid.hpp>
38#include <opm/grid/utility/ElementChunks.hpp>
39
40#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
41#include <opm/input/eclipse/EclipseState/Grid/FaceDir.hpp>
42#include <opm/input/eclipse/EclipseState/Grid/FieldPropsManager.hpp>
43#include <opm/input/eclipse/EclipseState/Grid/TransMult.hpp>
44#include <opm/input/eclipse/Units/Units.hpp>
45
47
48#include <algorithm>
49#include <array>
50#include <cassert>
51#include <cmath>
52#include <cstddef>
53#include <cstdint>
54#include <functional>
55#include <initializer_list>
56#include <sstream>
57#include <stdexcept>
58#include <type_traits>
59#include <utility>
60#include <vector>
61
62#include <fmt/format.h>
63
64namespace Opm {
65
66namespace details {
67
68 constexpr unsigned elemIdxShift = 32; // bits
69
70 std::uint64_t isId(std::uint32_t elemIdx1, std::uint32_t elemIdx2)
71 {
72 const std::uint32_t elemAIdx = std::min(elemIdx1, elemIdx2);
73 const std::uint64_t elemBIdx = std::max(elemIdx1, elemIdx2);
74
75 return (elemBIdx << elemIdxShift) + elemAIdx;
76 }
77
78 std::pair<std::uint32_t, std::uint32_t> isIdReverse(const std::uint64_t& id)
79 {
80 // Assigning an unsigned integer to a narrower type discards the most significant bits.
81 // See "The C programming language", section A.6.2.
82 // NOTE that the ordering of element A and B may have changed
83 const std::uint32_t elemAIdx = static_cast<uint32_t>(id);
84 const std::uint32_t elemBIdx = (id - elemAIdx) >> elemIdxShift;
85
86 return std::make_pair(elemAIdx, elemBIdx);
87 }
88
89 std::uint64_t directionalIsId(std::uint32_t elemIdx1, std::uint32_t elemIdx2)
90 {
91 return (std::uint64_t(elemIdx1) << elemIdxShift) + elemIdx2;
92 }
93}
94
95template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
97Transmissibility(const EclipseState& eclState,
98 const GridView& gridView,
99 const CartesianIndexMapper& cartMapper,
100 const Grid& grid,
101 std::function<std::array<double,dimWorld>(int)> centroids,
102 bool enableEnergy,
103 bool enableDiffusivity,
104 bool enableDispersivity)
105 : eclState_(eclState)
106 , gridView_(gridView)
107 , cartMapper_(cartMapper)
108 , grid_(grid)
109 , centroids_(centroids)
110 , enableEnergy_(enableEnergy)
111 , enableDiffusivity_(enableDiffusivity)
112 , enableDispersivity_(enableDispersivity)
113 , lookUpData_(gridView)
114 , lookUpCartesianData_(gridView, cartMapper)
115{
116 const UnitSystem& unitSystem = eclState_.getDeckUnitSystem();
117 transmissibilityThreshold_ = unitSystem.parse("Transmissibility").getSIScaling() * 1e-6;
118}
119
120template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
122transmissibility(unsigned elemIdx1, unsigned elemIdx2) const
123{
124 return trans_.at(details::isId(elemIdx1, elemIdx2));
125}
126
127template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
129transmissibilityBoundary(unsigned elemIdx, unsigned boundaryFaceIdx) const
130{
131 return transBoundary_.at(std::make_pair(elemIdx, boundaryFaceIdx));
132}
133
134template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
136thermalHalfTrans(unsigned insideElemIdx, unsigned outsideElemIdx) const
137{
138 return thermalHalfTrans_.at(details::directionalIsId(insideElemIdx, outsideElemIdx));
139}
140
141template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
143thermalHalfTransBoundary(unsigned insideElemIdx, unsigned boundaryFaceIdx) const
144{
145 return thermalHalfTransBoundary_.at(std::make_pair(insideElemIdx, boundaryFaceIdx));
146}
147
148template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
149const std::map<std::pair<unsigned, unsigned>, Scalar>& Transmissibility<Grid,GridView,ElementMapper,CartesianIndexMapper,Scalar>::
151{
152 return thermalHalfTransBoundary_;
153}
154
155template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
157diffusivity(unsigned elemIdx1, unsigned elemIdx2) const
158{
159 if (diffusivity_.empty())
160 return 0.0;
161
162 return diffusivity_.at(details::isId(elemIdx1, elemIdx2));
163}
164
165template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
167dispersivity(unsigned elemIdx1, unsigned elemIdx2) const
168{
169 if (dispersivity_.empty())
170 return 0.0;
171
172 return dispersivity_.at(details::isId(elemIdx1, elemIdx2));
173}
174
175template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
177update(bool global, const TransUpdateQuantities update_quantities,
178 const std::function<unsigned int(unsigned int)>& map, const bool applyNncMultregT)
179{
180 // whether only update the permeability related transmissibility
181 const bool onlyTrans = (update_quantities == TransUpdateQuantities::Trans);
182 const auto& cartDims = cartMapper_.cartesianDimensions();
183 const auto& transMult = eclState_.getTransMult();
184 const auto& comm = gridView_.comm();
185 ElementMapper elemMapper(gridView_, Dune::mcmgElementLayout());
186
187 unsigned numElements = elemMapper.size();
188 // get the ntg values, the ntg values are modified for the cells merged with minpv
189 const std::vector<double>& ntg = this->lookUpData_.assignFieldPropsDoubleOnLeaf(eclState_.fieldProps(), "NTG");
190 const bool updateDiffusivity = eclState_.getSimulationConfig().isDiffusive();
191 const bool updateDispersivity = eclState_.getSimulationConfig().rock_config().dispersion();
192
193 const bool disableNNC = eclState_.getSimulationConfig().useNONNC();
194
195 if (map) {
196 extractPermeability_(map);
197 }
198 else {
199 extractPermeability_();
200 }
201
202 const int num_threads = ThreadManager::maxThreads();
203
204 // reserving some space in the hashmap upfront saves quite a bit of time because
205 // resizes are costly for hashmaps and there would be quite a few of them if we
206 // would not have a rough idea of how large the final map will be (the rough idea
207 // is a conforming Cartesian grid).
208 trans_.clear();
209 if (num_threads == 1) {
210 trans_.reserve(numElements*3*1.05);
211 }
212
213 transBoundary_.clear();
214
215 // if energy is enabled, let's do the same for the "thermal half transmissibilities"
216 if ( enableEnergy_ && !onlyTrans) {
217 thermalHalfTrans_.clear();
218 if (num_threads == 1) {
219 thermalHalfTrans_.reserve(numElements*6*1.05);
220 }
221
222 thermalHalfTransBoundary_.clear();
223 }
224
225 // if diffusion is enabled, let's do the same for the "diffusivity"
226 if (updateDiffusivity && !onlyTrans) {
227 diffusivity_.clear();
228 if (num_threads == 1) {
229 diffusivity_.reserve(numElements*3*1.05);
230 }
231 extractPorosity_();
232 }
233
234 // if dispersion is enabled, let's do the same for the "dispersivity"
235 if (updateDispersivity && !onlyTrans) {
236 dispersivity_.clear();
237 if (num_threads == 1) {
238 dispersivity_.reserve(numElements*3*1.05);
239 }
240 extractDispersion_();
241 }
242
243 // The MULTZ needs special case if the option is ALL
244 // Then the smallest multiplier is applied.
245 // Default is to apply the top and bottom multiplier
246 bool useSmallestMultiplier;
247 bool pinchOption4ALL;
248 bool pinchActive;
249 if (comm.rank() == 0) {
250 const auto& eclGrid = eclState_.getInputGrid();
251 pinchActive = eclGrid.isPinchActive();
252 auto pinchTransCalcMode = eclGrid.getPinchOption();
253 useSmallestMultiplier = eclGrid.getMultzOption() == PinchMode::ALL;
254 pinchOption4ALL = (pinchTransCalcMode == PinchMode::ALL);
255 if (pinchOption4ALL) {
256 useSmallestMultiplier = false;
257 }
258 }
259 if (global && comm.size() > 1) {
260 comm.broadcast(&useSmallestMultiplier, 1, 0);
261 comm.broadcast(&pinchOption4ALL, 1, 0);
262 comm.broadcast(&pinchActive, 1, 0);
263 }
264
265 // fill the centroids cache to avoid repeated calculations in loops below
266 centroids_cache_.resize(gridView_.size(0));
267 for (const auto& elem : elements(gridView_)) {
268 const unsigned elemIdx = elemMapper.index(elem);
269 centroids_cache_[elemIdx] = centroids_(elemIdx);
270 }
271
272 auto harmonicMean = [](const Scalar x1, const Scalar x2)
273 {
274 return (std::abs(x1) < 1e-30 || std::abs(x2) < 1e-30)
275 ? 0.0
276 : 1.0 / (1.0 / x1 + 1.0 / x2);
277 };
278
279 auto faceIdToDir = [](int insideFaceIdx)
280 {
281 switch (insideFaceIdx) {
282 case 0:
283 case 1:
284 return FaceDir::XPlus;
285 case 2:
286 case 3:
287 return FaceDir::YPlus;
288 break;
289 case 4:
290 case 5:
291 return FaceDir::ZPlus;
292 default:
293 throw std::logic_error("Could not determine a face direction");
294 }
295 };
296
297 auto halfDiff = [](const DimVector& faceAreaNormal,
298 const unsigned,
299 const DimVector& distVector,
300 const Scalar prop)
301 {
302 return computeHalfDiffusivity_(faceAreaNormal,
303 distVector,
304 prop);
305 };
306
307 ThreadSafeMapBuilder transBoundary(transBoundary_, num_threads,
308 MapBuilderInsertionMode::Insert_Or_Assign);
309 ThreadSafeMapBuilder transMap(trans_, num_threads,
310 MapBuilderInsertionMode::Insert_Or_Assign);
311 ThreadSafeMapBuilder thermalHalfTransBoundary(thermalHalfTransBoundary_, num_threads,
312 MapBuilderInsertionMode::Insert_Or_Assign);
313 ThreadSafeMapBuilder thermalHalfTrans(thermalHalfTrans_, num_threads,
314 MapBuilderInsertionMode::Insert_Or_Assign);
315 ThreadSafeMapBuilder diffusivity(diffusivity_, num_threads,
316 MapBuilderInsertionMode::Insert_Or_Assign);
317 ThreadSafeMapBuilder dispersivity(dispersivity_, num_threads,
318 MapBuilderInsertionMode::Insert_Or_Assign);
319
320 const auto& nnc_input = eclState_.getInputNNC().input();
321
322#ifdef _OPENMP
323#pragma omp parallel for
324#endif
325 for (const auto& chunk : ElementChunks(gridView_, Dune::Partitions::all, num_threads)) {
326 for (const auto& elem : chunk) {
327 FaceInfo inside;
328 FaceInfo outside;
329 DimVector faceAreaNormal;
330
331 inside.elemIdx = elemMapper.index(elem);
332 // Get the Cartesian index of the origin cells (parent or equivalent cell on level zero),
333 // for CpGrid with LGRs. For general grids and no LGRs, get the usual Cartesian Index.
334 inside.cartElemIdx = this->lookUpCartesianData_.
335 template getFieldPropCartesianIdx<Grid>(inside.elemIdx);
336
337 auto computeHalf = [this, &faceAreaNormal, &inside, &outside]
338 (const auto& halfComputer,
339 const auto& prop1, const auto& prop2) -> std::array<Scalar,2>
340 {
341 return {
342 halfComputer(faceAreaNormal,
343 inside.faceIdx,
344 distanceVector_(inside.faceCenter, inside.elemIdx),
345 prop1),
346 halfComputer(faceAreaNormal,
347 outside.faceIdx,
348 distanceVector_(outside.faceCenter, outside.elemIdx),
349 prop2)
350 };
351 };
352
353 auto computeHalfMean = [&inside, &outside, &computeHalf, &ntg, &harmonicMean]
354 (const auto& halfComputer, const auto& prop)
355 {
356 auto onesided = computeHalf(halfComputer, prop[inside.elemIdx], prop[outside.elemIdx]);
357 applyNtg_(onesided[0], inside, ntg);
358 applyNtg_(onesided[1], outside, ntg);
359
360 //TODO Add support for multipliers
361 return harmonicMean(onesided[0], onesided[1]);
362 };
363
364 unsigned boundaryIsIdx = 0;
365 for (const auto& intersection : intersections(gridView_, elem)) {
366 // deal with grid boundaries
367 if (intersection.boundary()) {
368 // compute the transmissibilty for the boundary intersection
369 const auto& geometry = intersection.geometry();
370 inside.faceCenter = geometry.center();
371
372 faceAreaNormal = intersection.centerUnitOuterNormal();
373 faceAreaNormal *= geometry.volume();
374
375 Scalar transBoundaryIs =
376 computeHalfTrans_(faceAreaNormal,
377 intersection.indexInInside(),
378 distanceVector_(inside.faceCenter, inside.elemIdx),
379 permeability_[inside.elemIdx]);
380
381 // normally there would be two half-transmissibilities that would be
382 // averaged. on the grid boundary there only is the half
383 // transmissibility of the interior element.
384 applyMultipliers_(transBoundaryIs, intersection.indexInInside(), inside.cartElemIdx, transMult);
385 transBoundary.insert_or_assign(std::make_pair(inside.elemIdx, boundaryIsIdx), transBoundaryIs);
386
387 // for boundary intersections we also need to compute the thermal
388 // half transmissibilities
389 if (enableEnergy_ && !onlyTrans) {
390 Scalar transBoundaryEnergyIs =
391 computeHalfDiffusivity_(faceAreaNormal,
392 distanceVector_(inside.faceCenter, inside.elemIdx),
393 1.0);
394 thermalHalfTransBoundary.insert_or_assign(std::make_pair(inside.elemIdx, boundaryIsIdx),
395 transBoundaryEnergyIs);
396 }
397
398 ++boundaryIsIdx;
399 continue;
400 }
401
402 if (!intersection.neighbor()) {
403 // elements can be on process boundaries, i.e. they are not on the
404 // domain boundary yet they don't have neighbors.
405 ++boundaryIsIdx;
406 continue;
407 }
408
409 const auto& outsideElem = intersection.outside();
410 outside.elemIdx = elemMapper.index(outsideElem);
411
412 // Get the Cartesian index of the origin cells (parent or equivalent cell on level zero),
413 // for CpGrid with LGRs. For general grids and no LGRs, get the usual Cartesian Index.
414 outside.cartElemIdx = this->lookUpCartesianData_.
415 template getFieldPropCartesianIdx<Grid>(outside.elemIdx);
416
417 // we only need to calculate a face's transmissibility
418 // once...
419 // In a parallel run inside.cartElemIdx > outside.cartElemIdx does not imply inside.elemIdx > outside.elemIdx for
420 // ghost cells and we need to use the cartesian index as this will be used when applying Z multipliers
421 // To cover the case where both cells are part of an LGR and as a consequence might have
422 // the same cartesian index, we tie their Cartesian indices and the ones on the leaf grid view.
423 if (std::tie(inside.cartElemIdx, inside.elemIdx) > std::tie(outside.cartElemIdx, outside.elemIdx)) {
424 continue;
425 }
426
427 // local indices of the faces of the inside and
428 // outside elements which contain the intersection
429 inside.faceIdx = intersection.indexInInside();
430 outside.faceIdx = intersection.indexInOutside();
431
432 if (inside.faceIdx == -1) {
433 // NNC. Set zero transmissibility, as it will be
434 // *added to* by applyNncToGridTrans_() later.
435 assert(outside.faceIdx == -1);
436 transMap.insert_or_assign(details::isId(inside.elemIdx, outside.elemIdx), 0.0);
437 if (enableEnergy_ && !onlyTrans) {
438 thermalHalfTrans.insert_or_assign(details::directionalIsId(inside.elemIdx, outside.elemIdx), 0.0);
439 thermalHalfTrans.insert_or_assign(details::directionalIsId(outside.elemIdx, inside.elemIdx), 0.0);
440 }
441
442 if (updateDiffusivity && !onlyTrans) {
443 diffusivity.insert_or_assign(details::isId(inside.elemIdx, outside.elemIdx), 0.0);
444 }
445 if (updateDispersivity && !onlyTrans) {
446 dispersivity.insert_or_assign(details::isId(inside.elemIdx, outside.elemIdx), 0.0);
447 }
448 continue;
449 }
450
451 typename std::is_same<Grid, Dune::CpGrid>::type isCpGrid;
452 computeFaceProperties(intersection,
453 inside,
454 outside,
455 faceAreaNormal,
456 isCpGrid);
457
458 Scalar trans = computeHalfMean(computeHalfTrans_, permeability_);
459
460 // apply the full face transmissibility multipliers
461 // for the inside ...
462 if (!pinchActive) {
463 if (inside.faceIdx > 3) { // top or bottom
464 auto find_layer = [&cartDims](std::size_t cell) {
465 cell /= cartDims[0];
466 auto k = cell / cartDims[1];
467 return k;
468 };
469 int kup = find_layer(inside.cartElemIdx);
470 int kdown = find_layer(outside.cartElemIdx);
471 // When a grid is a CpGrid with LGRs, insideCartElemIdx coincides with outsideCartElemIdx
472 // for cells on the leaf with the same parent cell on level zero.
473 assert((kup != kdown) || (inside.cartElemIdx == outside.cartElemIdx));
474 if (std::abs(kup -kdown) > 1) {
475 trans = 0.0;
476 }
477 }
478 }
479
480 if (useSmallestMultiplier) {
481 // PINCH(4) == TOPBOT is assumed here as we set useSmallestMultipliers
482 // to false if PINCH(4) == ALL holds
483 // In contrast to the name this will also apply
484 applyAllZMultipliers_(trans, inside, outside, transMult, cartDims);
485 }
486 else {
487 applyMultipliers_(trans, inside.faceIdx, inside.cartElemIdx, transMult);
488 // ... and outside elements
489 applyMultipliers_(trans, outside.faceIdx, outside.cartElemIdx, transMult);
490 }
491
492 bool foundInputNNC = false;
493 if (! nnc_input.empty()) {
494 // Skip region multipliers for overlapping input NNCs (they are handled later)
495 auto it = std::lower_bound(nnc_input.begin(), nnc_input.end(),
496 NNCdata { inside.cartElemIdx, outside.cartElemIdx, 0.0 });
497 foundInputNNC = it != nnc_input.end() && it->cell1 == inside.cartElemIdx && it->cell2 == outside.cartElemIdx;
498 }
499 if (! foundInputNNC) {
500 // apply the region multipliers (cf. the MULTREGT keyword)
501 trans *= transMult.getRegionMultiplier(inside.cartElemIdx,
502 outside.cartElemIdx,
503 faceIdToDir(inside.faceIdx));
504 }
505
506 transMap.insert_or_assign(details::isId(inside.elemIdx, outside.elemIdx), trans);
507
508 // update the "thermal half transmissibility" for the intersection
509 if (enableEnergy_ && !onlyTrans) {
510 const auto half = computeHalf(halfDiff, 1.0, 1.0);
511 // TODO Add support for multipliers
512 thermalHalfTrans.insert_or_assign(details::directionalIsId(inside.elemIdx, outside.elemIdx),
513 half[0]);
514 thermalHalfTrans.insert_or_assign(details::directionalIsId(outside.elemIdx, inside.elemIdx),
515 half[1]);
516 }
517
518 // update the "diffusive half transmissibility" for the intersection
519 if (updateDiffusivity && !onlyTrans) {
520 diffusivity.insert_or_assign(details::isId(inside.elemIdx, outside.elemIdx),
521 computeHalfMean(halfDiff, porosity_));
522 }
523
524 // update the "dispersivity half transmissibility" for the intersection
525 if (updateDispersivity && !onlyTrans) {
526 dispersivity.insert_or_assign(details::isId(inside.elemIdx, outside.elemIdx),
527 computeHalfMean(halfDiff, dispersion_));
528 }
529 }
530 }
531 }
532 centroids_cache_.clear();
533
534#ifdef _OPENMP
535#pragma omp parallel sections
536#endif
537 {
538#ifdef _OPENMP
539#pragma omp section
540#endif
541 transMap.finalize();
542#ifdef _OPENMP
543#pragma omp section
544#endif
545 transBoundary.finalize();
546#ifdef _OPENMP
547#pragma omp section
548#endif
549 thermalHalfTransBoundary.finalize();
550#ifdef _OPENMP
551#pragma omp section
552#endif
553 thermalHalfTrans.finalize();
554#ifdef _OPENMP
555#pragma omp section
556#endif
557 diffusivity.finalize();
558#ifdef _OPENMP
559#pragma omp section
560#endif
561 dispersivity.finalize();
562 }
563
564 // Potentially overwrite and/or modify transmissibilities based on input from deck
565 this->updateFromEclState_(global);
566
567 // Create mapping from global to local index
568 std::unordered_map<std::size_t,int> globalToLocal;
569
570 // Loop over all elements (global grid) and store Cartesian index
571 for (const auto& elem : elements(grid_.leafGridView())) {
572 int elemIdx = elemMapper.index(elem);
573 int cartElemIdx = cartMapper_.cartesianIndex(elemIdx);
574 globalToLocal[cartElemIdx] = elemIdx;
575 }
576
577 if (!disableNNC) {
578 // For EDITNNC and EDITNNCR we warn only once
579 // If transmissibility is used for load balancing this will be done
580 // when computing the gobal transmissibilities and all warnings will
581 // be seen in a parallel. Unfortunately, when we do not use transmissibilities
582 // we will only see warnings for the partition of process 0 and also false positives.
583 this->applyPinchNncToGridTrans_(globalToLocal, applyNncMultregT);
584 this->applyNncToGridTrans_(globalToLocal);
585 this->applyEditNncToGridTrans_(globalToLocal);
586 this->applyEditNncrToGridTrans_(globalToLocal);
587 if (applyNncMultregT) {
588 this->applyNncMultreg_(globalToLocal);
589 }
590 warnEditNNC_ = false;
591 }
592
593 // If disableNNC == true, remove all non-neighbouring transmissibilities.
594 // If disableNNC == false, remove very small non-neighbouring transmissibilities.
595 this->removeNonCartesianTransmissibilities_(disableNNC);
596}
597
598template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
601{
602 unsigned numElem = gridView_.size(/*codim=*/0);
603 permeability_.resize(numElem);
604
605 // read the intrinsic permeabilities from the eclState. Note that all arrays
606 // provided by eclState are one-per-cell of "uncompressed" grid, whereas the
607 // simulation grid might remove a few elements. (e.g. because it is distributed
608 // over several processes.)
609 const auto& fp = eclState_.fieldProps();
610 if (fp.has_double("PERMX")) {
611 const std::vector<double>& permxData = this-> lookUpData_.assignFieldPropsDoubleOnLeaf(fp, "PERMX");
612
613 std::vector<double> permyData;
614 if (fp.has_double("PERMY"))
615 permyData = this-> lookUpData_.assignFieldPropsDoubleOnLeaf(fp,"PERMY");
616 else
617 permyData = permxData;
618
619 std::vector<double> permzData;
620 if (fp.has_double("PERMZ"))
621 permzData = this-> lookUpData_.assignFieldPropsDoubleOnLeaf(fp,"PERMZ");
622 else
623 permzData = permxData;
624
625 for (std::size_t dofIdx = 0; dofIdx < numElem; ++ dofIdx) {
626 permeability_[dofIdx] = 0.0;
627 permeability_[dofIdx][0][0] = permxData[dofIdx];
628 permeability_[dofIdx][1][1] = permyData[dofIdx];
629 permeability_[dofIdx][2][2] = permzData[dofIdx];
630 }
631
632 // for now we don't care about non-diagonal entries
633
634 }
635 else
636 throw std::logic_error("Can't read the intrinsic permeability from the ecl state. "
637 "(The PERM{X,Y,Z} keywords are missing)");
638}
639
640template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
642extractPermeability_(const std::function<unsigned int(unsigned int)>& map)
643{
644 unsigned numElem = gridView_.size(/*codim=*/0);
645 permeability_.resize(numElem);
646
647 // read the intrinsic permeabilities from the eclState. Note that all arrays
648 // provided by eclState are one-per-cell of "uncompressed" grid, whereas the
649 // simulation grid might remove a few elements. (e.g. because it is distributed
650 // over several processes.)
651 const auto& fp = eclState_.fieldProps();
652 if (fp.has_double("PERMX")) {
653 const std::vector<double>& permxData =
654 this->lookUpData_.assignFieldPropsDoubleOnLeaf(fp,"PERMX");
655
656 std::vector<double> permyData;
657 if (fp.has_double("PERMY")){
658 permyData = this->lookUpData_.assignFieldPropsDoubleOnLeaf(fp,"PERMY");
659 }
660 else {
661 permyData = permxData;
662 }
663
664 std::vector<double> permzData;
665 if (fp.has_double("PERMZ")) {
666 permzData = this->lookUpData_.assignFieldPropsDoubleOnLeaf(fp,"PERMZ");
667 }
668 else {
669 permzData = permxData;
670 }
671
672 for (std::size_t dofIdx = 0; dofIdx < numElem; ++ dofIdx) {
673 permeability_[dofIdx] = 0.0;
674 std::size_t inputDofIdx = map(dofIdx);
675 permeability_[dofIdx][0][0] = permxData[inputDofIdx];
676 permeability_[dofIdx][1][1] = permyData[inputDofIdx];
677 permeability_[dofIdx][2][2] = permzData[inputDofIdx];
678 }
679
680 // for now we don't care about non-diagonal entries
681 }
682 else {
683 throw std::logic_error("Can't read the intrinsic permeability from the ecl state. "
684 "(The PERM{X,Y,Z} keywords are missing)");
685 }
686}
687
688template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
691{
692 // read the intrinsic porosity from the eclState. Note that all arrays
693 // provided by eclState are one-per-cell of "uncompressed" grid, whereas the
694 // simulation grid might remove a few elements. (e.g. because it is distributed
695 // over several processes.)
696 const auto& fp = eclState_.fieldProps();
697 if (fp.has_double("PORO")) {
698 if constexpr (std::is_same_v<Scalar,double>) {
699 porosity_ = this->lookUpData_.assignFieldPropsDoubleOnLeaf(fp,"PORO");
700 }
701 else {
702 const auto por = this->lookUpData_.assignFieldPropsDoubleOnLeaf(fp,"PORO");
703 porosity_.resize(por.size());
704 std::ranges::copy(por, porosity_.begin());
705 }
706 }
707 else {
708 throw std::logic_error("Can't read the porosity from the ecl state. "
709 "(The PORO keywords are missing)");
710 }
711}
712
713template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
716{
717 if (!enableDispersivity_) {
718 throw std::runtime_error("Dispersion disabled at compile time, but the deck "
719 "contains the DISPERC keyword.");
720 }
721 const auto& fp = eclState_.fieldProps();
722 if constexpr (std::is_same_v<Scalar,double>) {
723 dispersion_ = this->lookUpData_.assignFieldPropsDoubleOnLeaf(fp,"DISPERC");
724 }
725 else {
726 const auto disp = this->lookUpData_.assignFieldPropsDoubleOnLeaf(fp,"DISPERC");
727 dispersion_.resize(disp.size());
728 std::ranges::copy(disp, dispersion_.begin());
729 }
730}
731
732template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
735{
736 const auto& cartDims = cartMapper_.cartesianDimensions();
737 for (auto&& trans: trans_) {
738 //either remove all NNC transmissibilities or those less than the threshold (by default 1e-6 in the deck's unit system)
739 if (removeAll || trans.second < transmissibilityThreshold_) {
740 const auto& id = trans.first;
741 const auto& elements = details::isIdReverse(id);
742 int gc1 = std::min(cartMapper_.cartesianIndex(elements.first), cartMapper_.cartesianIndex(elements.second));
743 int gc2 = std::max(cartMapper_.cartesianIndex(elements.first), cartMapper_.cartesianIndex(elements.second));
744
745 // only adjust the NNCs
746 // When LGRs, all neighbors in the LGR are cartesian neighbours on the level grid representing the LGR.
747 // When elements on the leaf grid view have the same parent cell, gc1 and gc2 coincide.
748 if (gc2 - gc1 == 1 || gc2 - gc1 == cartDims[0] || gc2 - gc1 == cartDims[0]*cartDims[1] || gc2 - gc1 == 0) {
749 continue;
750 }
751
752 trans.second = 0.0;
753 }
754 }
755}
756
757template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
759applyAllZMultipliers_(Scalar& trans,
760 const FaceInfo& inside,
761 const FaceInfo& outside,
762 const TransMult& transMult,
763 const std::array<int, dimWorld>& cartDims)
764{
765 if (grid_.maxLevel() > 0) {
766 OPM_THROW(std::invalid_argument, "MULTZ not support with LGRS, yet.");
767 }
768 if (inside.faceIdx > 3) { // top or or bottom
769 assert(inside.faceIdx == 5); // as insideCartElemIdx < outsideCartElemIdx holds for the Z column
770 // For CpGrid with LGRs, insideCartElemIdx == outsideCartElemIdx when cells on the leaf have the same parent cell on level zero.
771 assert(outside.cartElemIdx >= inside.cartElemIdx);
772 unsigned lastCartElemIdx;
773 if (outside.cartElemIdx == inside.cartElemIdx) {
774 lastCartElemIdx = outside.cartElemIdx;
775 }
776 else {
777 lastCartElemIdx = outside.cartElemIdx - cartDims[0]*cartDims[1];
778 }
779 // Last multiplier using (Z+)*(Z-)
780 Scalar mult = transMult.getMultiplier(lastCartElemIdx , FaceDir::ZPlus) *
781 transMult.getMultiplier(outside.cartElemIdx , FaceDir::ZMinus);
782
783 // pick the smallest multiplier using (Z+)*(Z-) while looking down
784 // the pillar until reaching the other end of the connection
785 for (auto cartElemIdx = inside.cartElemIdx; cartElemIdx < lastCartElemIdx;) {
786 auto multiplier = transMult.getMultiplier(cartElemIdx, FaceDir::ZPlus);
787 cartElemIdx += cartDims[0]*cartDims[1];
788 multiplier *= transMult.getMultiplier(cartElemIdx, FaceDir::ZMinus);
789 mult = std::min(mult, static_cast<Scalar>(multiplier));
790 }
791
792 trans *= mult;
793 }
794 else {
795 applyMultipliers_(trans, inside.faceIdx, inside.cartElemIdx, transMult);
796 applyMultipliers_(trans, outside.faceIdx, outside.cartElemIdx, transMult);
797 }
798}
799
800template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
802updateFromEclState_(bool global)
803{
804 const FieldPropsManager* fp =
805 (global) ? &(eclState_.fieldProps()) :
806 &(eclState_.globalFieldProps());
807
808 std::array<bool,3> is_tran {fp->tran_active("TRANX"),
809 fp->tran_active("TRANY"),
810 fp->tran_active("TRANZ")};
811
812 if (!(is_tran[0] || is_tran[1] || is_tran[2])) {
813 // Skip unneeded expensive traversals
814 return;
815 }
816
817 std::array<std::string, 3> keywords {"TRANX", "TRANY", "TRANZ"};
818 std::array<std::vector<double>,3> trans = createTransmissibilityArrays_(is_tran);
819 auto key = keywords.begin();
820 auto perform = is_tran.begin();
821
822 for (auto it = trans.begin(); it != trans.end(); ++it, ++key, ++perform) {
823 if (*perform) {
824 if (grid_.maxLevel() > 0) {
825 OPM_THROW(std::invalid_argument, "Calculations on TRANX/TRANY/TRANZ arrays are not support with LGRS, yet.");
826 }
827 fp->apply_tran(*key, *it);
828 }
829 }
830
831 resetTransmissibilityFromArrays_(is_tran, trans);
832}
833
834template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
835std::array<std::vector<double>,3>
837createTransmissibilityArrays_(const std::array<bool,3>& is_tran)
838{
839 const auto& cartDims = cartMapper_.cartesianDimensions();
840 ElementMapper elemMapper(gridView_, Dune::mcmgElementLayout());
841
842 auto numElem = gridView_.size(/*codim=*/0);
843 std::array<std::vector<double>,3> trans = {
844 std::vector<double>(is_tran[0] ? numElem : 0, 0),
845 std::vector<double>(is_tran[1] ? numElem : 0, 0),
846 std::vector<double>(is_tran[2] ? numElem : 0, 0)
847 };
848
849 // compute the transmissibilities for all intersections
850 for (const auto& elem : elements(gridView_)) {
851 for (const auto& intersection : intersections(gridView_, elem)) {
852 // store intersection, this might be costly
853 if (!intersection.neighbor()) {
854 continue; // intersection is on the domain boundary
855 }
856
857 // In the EclState TRANX[c1] is transmissibility in X+
858 // direction. we only store transmissibilities in the +
859 // direction. Same for Y and Z. Ordering of compressed (c1,c2) and cartesian index
860 // (gc1, gc2) is coherent (c1 < c2 <=> gc1 < gc2) only in a serial run.
861 // In a parallel run this only holds in the interior as elements in the
862 // ghost overlap region might be ordered after the others. Hence we need
863 // to use the cartesian index to select the compressed index where to store
864 // the transmissibility value.
865 // c1 < c2 <=> gc1 < gc2 is no longer true (even in serial) when the grid is a
866 // CpGrid with LGRs. When cells c1 and c2 have the same parent
867 // cell on level zero, then gc1 == gc2.
868 unsigned c1 = elemMapper.index(intersection.inside());
869 unsigned c2 = elemMapper.index(intersection.outside());
870 int gc1 = cartMapper_.cartesianIndex(c1);
871 int gc2 = cartMapper_.cartesianIndex(c2);
872 if (std::tie(gc1, c1) > std::tie(gc2, c2)) {
873 // we only need to handle each connection once, thank you.
874 // We do this when gc1 is smaller than the other to find the
875 // correct place to store in parallel when ghost/overlap elements
876 // are ordered last
877 continue;
878 }
879
880 auto isID = details::isId(c1, c2);
881
882 // For CpGrid with LGRs, when leaf grid view cells with indices c1 and c2
883 // have the same parent cell on level zero, then gc2 - gc1 == 0. In that case,
884 // 'intersection.indexInSIde()' needed to be checked to determine the direction, i.e.
885 // add in the if/else-if 'gc2 == gc1 && intersection.indexInInside() == ... '
886 if ((gc2 - gc1 == 1 || (gc2 == gc1 && (intersection.indexInInside() == 0 || intersection.indexInInside() == 1)))
887 && cartDims[0] > 1)
888 {
889 if (is_tran[0]) {
890 // set simulator internal transmissibilities to values from inputTranx
891 trans[0][c1] = trans_[isID];
892 }
893 }
894 else if ((gc2 - gc1 == cartDims[0] || (gc2 == gc1 && (intersection.indexInInside() == 2 || intersection.indexInInside() == 3)))
895 && cartDims[1] > 1)
896 {
897 if (is_tran[1]) {
898 // set simulator internal transmissibilities to values from inputTrany
899 trans[1][c1] = trans_[isID];
900 }
901 }
902 else if (gc2 - gc1 == cartDims[0]*cartDims[1] ||
903 (gc2 == gc1 && (intersection.indexInInside() == 4 || intersection.indexInInside() == 5)))
904 {
905 if (is_tran[2]) {
906 // set simulator internal transmissibilities to values from inputTranz
907 trans[2][c1] = trans_[isID];
908 }
909 }
910 // else.. We don't support modification of NNC at the moment.
911 }
912 }
913
914 return trans;
915}
916
917template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
919resetTransmissibilityFromArrays_(const std::array<bool,3>& is_tran,
920 const std::array<std::vector<double>,3>& trans)
921{
922 const auto& cartDims = cartMapper_.cartesianDimensions();
923 ElementMapper elemMapper(gridView_, Dune::mcmgElementLayout());
924
925 // compute the transmissibilities for all intersections
926 for (const auto& elem : elements(gridView_)) {
927 for (const auto& intersection : intersections(gridView_, elem)) {
928 if (!intersection.neighbor()) {
929 continue; // intersection is on the domain boundary
930 }
931
932 // In the EclState TRANX[c1] is transmissibility in X+
933 // direction. we only store transmissibilities in the +
934 // direction. Same for Y and Z. Ordering of compressed (c1,c2) and cartesian index
935 // (gc1, gc2) is coherent (c1 < c2 <=> gc1 < gc2) only in a serial run.
936 // In a parallel run this only holds in the interior as elements in the
937 // ghost overlap region might be ordered after the others. Hence we need
938 // to use the cartesian index to select the compressed index where to store
939 // the transmissibility value.
940 // c1 < c2 <=> gc1 < gc2 is no longer true (even in serial) when the grid is a
941 // CpGrid with LGRs. When cells c1 and c2 have the same parent
942 // cell on level zero, then gc1 == gc2.
943 unsigned c1 = elemMapper.index(intersection.inside());
944 unsigned c2 = elemMapper.index(intersection.outside());
945 int gc1 = cartMapper_.cartesianIndex(c1);
946 int gc2 = cartMapper_.cartesianIndex(c2);
947 if (std::tie(gc1, c1) > std::tie(gc2, c2)) {
948 // we only need to handle each connection once, thank you.
949 // We do this when gc1 is smaller than the other to find the
950 // correct place to read in parallel when ghost/overlap elements
951 // are ordered last
952 continue;
953 }
954
955 auto isID = details::isId(c1, c2);
956
957 // For CpGrid with LGRs, when leaf grid view cells with indices c1 and c2
958 // have the same parent cell on level zero, then gc2 - gc1 == 0. In that case,
959 // 'intersection.indexInSIde()' needed to be checked to determine the direction, i.e.
960 // add in the if/else-if 'gc2 == gc1 && intersection.indexInInside() == ... '
961 if ((gc2 - gc1 == 1 || (gc2 == gc1 && (intersection.indexInInside() == 0 || intersection.indexInInside() == 1)))
962 && cartDims[0] > 1)
963 {
964 if (is_tran[0]) {
965 // set simulator internal transmissibilities to values from inputTranx
966 trans_[isID] = trans[0][c1];
967 }
968 }
969 else if ((gc2 - gc1 == cartDims[0] || (gc2 == gc1 && (intersection.indexInInside() == 2|| intersection.indexInInside() == 3)))
970 && cartDims[1] > 1)
971 {
972 if (is_tran[1]) {
973 // set simulator internal transmissibilities to values from inputTrany
974 trans_[isID] = trans[1][c1];
975 }
976 }
977 else if (gc2 - gc1 == cartDims[0]*cartDims[1] ||
978 (gc2 == gc1 && (intersection.indexInInside() == 4 || intersection.indexInInside() == 5)))
979 {
980 if (is_tran[2]) {
981 // set simulator internal transmissibilities to values from inputTranz
982 trans_[isID] = trans[2][c1];
983 }
984 }
985
986 // else.. We don't support modification of NNC at the moment.
987 }
988 }
989}
990
991template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
992template<class Intersection>
994computeFaceProperties(const Intersection& intersection,
995 FaceInfo& inside,
996 FaceInfo& outside,
997 DimVector& faceAreaNormal,
998 /*isCpGrid=*/std::false_type) const
999{
1000 // default implementation for DUNE grids
1001 const auto& geometry = intersection.geometry();
1002 outside.faceCenter = inside.faceCenter = geometry.center();
1003
1004 faceAreaNormal = intersection.centerUnitOuterNormal();
1005 faceAreaNormal *= geometry.volume();
1006}
1007
1008template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
1009template<class Intersection>
1011computeFaceProperties(const Intersection& intersection,
1012 FaceInfo& inside,
1013 FaceInfo& outside,
1014 DimVector& faceAreaNormal,
1015 /*isCpGrid=*/std::true_type) const
1016{
1017 int faceIdx = intersection.id();
1018
1019 if (grid_.maxLevel() == 0) {
1020 inside.faceCenter = grid_.faceCenterEcl(inside.elemIdx, inside.faceIdx, intersection);
1021 outside.faceCenter = grid_.faceCenterEcl(outside.elemIdx, outside.faceIdx, intersection);
1022 faceAreaNormal = grid_.faceAreaNormalEcl(faceIdx);
1023 }
1024 else {
1025 if ((intersection.inside().level() != intersection.outside().level())) {
1026 // For CpGrid with LGRs, intersection laying on the boundary of an LGR, having two neighboring cells:
1027 // one coarse neighboring cell and one refined neighboring cell, we get the corresponding parent
1028 // intersection (from level 0), and use the center of the parent intersection for the coarse
1029 // neighboring cell.
1030
1031 // Get parent intersection and its geometry
1032 const auto& parentIntersection =
1033 grid_.getParentIntersectionFromLgrBoundaryFace(intersection);
1034 const auto& parentIntersectionGeometry = parentIntersection.geometry();
1035
1036 // For the coarse neighboring cell, take the center of the parent intersection.
1037 // For the refined neighboring cell, take the 'usual' center.
1038 inside.faceCenter = (intersection.inside().level() == 0)
1039 ? parentIntersectionGeometry.center()
1040 : grid_.faceCenterEcl(inside.elemIdx, inside.faceIdx, intersection);
1041 outside.faceCenter = (intersection.outside().level() == 0)
1042 ? parentIntersectionGeometry.center()
1043 : grid_.faceCenterEcl(outside.elemIdx, outside.faceIdx, intersection);
1044
1045 // For some computations, it seems to be benefitial to replace the actual area of the refined face, by
1046 // the area of its parent face.
1047 // faceAreaNormal = parentIntersection.centerUnitOuterNormal();
1048 // faceAreaNormal *= parentIntersectionGeometry.volume();
1049
1051 faceAreaNormal = intersection.centerUnitOuterNormal();
1052 faceAreaNormal *= intersection.geometry().volume();
1053 }
1054 else {
1055 assert(intersection.inside().level() == intersection.outside().level());
1056
1057 inside.faceCenter = grid_.faceCenterEcl(inside.elemIdx, inside.faceIdx, intersection);
1058 outside.faceCenter = grid_.faceCenterEcl(outside.elemIdx, outside.faceIdx, intersection);
1059
1060 // When the CpGrid has LGRs, we compute the face area normal differently.
1061 if (intersection.inside().level() > 0) { // remove intersection.inside().level() > 0
1062 faceAreaNormal = intersection.centerUnitOuterNormal();
1063 faceAreaNormal *= intersection.geometry().volume();
1064 }
1065 else {
1066 faceAreaNormal = grid_.faceAreaNormalEcl(faceIdx);
1067 }
1068 }
1069 }
1070}
1071
1072template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
1073void
1075applyPinchNncToGridTrans_(const std::unordered_map<std::size_t,int>& cartesianToCompressed,
1076 const bool applyNncMultregT)
1077{
1078 const auto& pinchNnc = eclState_.getPinchNNC();
1079 const auto& transMult = this->eclState_.getTransMult();
1080
1081 for (const auto& nncEntry : pinchNnc) {
1082 auto c1 = nncEntry.cell1;
1083 auto c2 = nncEntry.cell2;
1084 auto lowIt = cartesianToCompressed.find(c1);
1085 auto highIt = cartesianToCompressed.find(c2);
1086 int low = (lowIt == cartesianToCompressed.end())? -1 : lowIt->second;
1087 int high = (highIt == cartesianToCompressed.end())? -1 : highIt->second;
1088
1089 if (low > high) {
1090 std::swap(low, high);
1091 }
1092
1093 if (low == -1 && high == -1) {
1094 // Silently discard as it is not between active cells
1095 continue;
1096 }
1097
1098 if (low == -1 || high == -1) {
1099 // We can end up here if one of the cells is overlap/ghost, because those
1100 // are lacking connections to other cells in the ghost/overlap.
1101 // Hence discard the NNC if it is between active cell and inactive cell
1102 continue;
1103 }
1104
1105 {
1106 auto candidate = trans_.find(details::isId(low, high));
1107 if (candidate != trans_.end()) {
1108 // the correctly calculated transmissibility is stored in
1109 // the NNC. Overwrite previous value with it.
1110 // taking the region multiplier into account.
1111 candidate->second = nncEntry.trans;
1112 if (applyNncMultregT) {
1113 const auto mult = transMult.getRegionMultiplierNNC(c1, c2);
1114 candidate->second *= mult;
1115 }
1116 }
1117 }
1118 }
1119}
1120
1121template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
1122void
1124applyNncToGridTrans_(const std::unordered_map<std::size_t,int>& cartesianToCompressed)
1125{
1126 // First scale NNCs with EDITNNC.
1127 const auto& nnc_input = eclState_.getInputNNC().input();
1128
1129 for (const auto& nncEntry : nnc_input) {
1130 auto c1 = nncEntry.cell1;
1131 auto c2 = nncEntry.cell2;
1132 auto lowIt = cartesianToCompressed.find(c1);
1133 auto highIt = cartesianToCompressed.find(c2);
1134 int low = (lowIt == cartesianToCompressed.end())? -1 : lowIt->second;
1135 int high = (highIt == cartesianToCompressed.end())? -1 : highIt->second;
1136
1137 if (low > high) {
1138 std::swap(low, high);
1139 }
1140
1141 if (low == -1 && high == -1) {
1142 // Silently discard as it is not between active cells
1143 continue;
1144 }
1145
1146 if (low == -1 || high == -1) {
1147 // Discard the NNC if it is between active cell and inactive cell
1148 std::ostringstream sstr;
1149 sstr << "NNC between active and inactive cells ("
1150 << low << " -> " << high << ") with globalcell is (" << c1 << "->" << c2 <<")";
1151 OpmLog::warning(sstr.str());
1152 continue;
1153 }
1154
1155 if (auto candidate = trans_.find(details::isId(low, high)); candidate != trans_.end()) {
1156 // NNC is represented by the grid and might be a neighboring connection
1157 // In this case the transmissibilty is added to the value already
1158 // set or computed.
1159 candidate->second += nncEntry.trans;
1160 }
1161 // if (enableEnergy_) {
1162 // auto candidate = thermalHalfTrans_.find(details::directionalIsId(low, high));
1163 // if (candidate != trans_.end()) {
1164 // // NNC is represented by the grid and might be a neighboring connection
1165 // // In this case the transmissibilty is added to the value already
1166 // // set or computed.
1167 // candidate->second += nncEntry.transEnergy1;
1168 // }
1169 // auto candidate = thermalHalfTrans_.find(details::directionalIsId(high, low));
1170 // if (candidate != trans_.end()) {
1171 // // NNC is represented by the grid and might be a neighboring connection
1172 // // In this case the transmissibilty is added to the value already
1173 // // set or computed.
1174 // candidate->second += nncEntry.transEnergy2;
1175 // }
1176 // }
1177 // if (enableDiffusivity_) {
1178 // auto candidate = diffusivity_.find(details::isId(low, high));
1179 // if (candidate != trans_.end()) {
1180 // // NNC is represented by the grid and might be a neighboring connection
1181 // // In this case the transmissibilty is added to the value already
1182 // // set or computed.
1183 // candidate->second += nncEntry.transDiffusion;
1184 // }
1185 // }
1186 }
1187}
1188
1189template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
1191applyEditNncToGridTrans_(const std::unordered_map<std::size_t,int>& globalToLocal)
1192{
1193 const auto& input = eclState_.getInputNNC();
1194 applyEditNncToGridTransHelper_(globalToLocal, "EDITNNC",
1195 input.edit(),
1196 [&input](const NNCdata& nnc){
1197 return input.edit_location(nnc);},
1198 // Multiply transmissibility with EDITNNC value
1199 [](Scalar& trans, const Scalar& rhs){ trans *= rhs;});
1200}
1201
1202template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
1204applyEditNncrToGridTrans_(const std::unordered_map<std::size_t,int>& globalToLocal)
1205{
1206 const auto& input = eclState_.getInputNNC();
1207 applyEditNncToGridTransHelper_(globalToLocal, "EDITNNCR",
1208 input.editr(),
1209 [&input](const NNCdata& nnc){
1210 return input.editr_location(nnc);},
1211 // Replace Transmissibility with EDITNNCR value
1212 [](Scalar& trans, const Scalar& rhs){ trans = rhs;});
1213}
1214
1215template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
1217applyEditNncToGridTransHelper_(const std::unordered_map<std::size_t,int>& globalToLocal,
1218 const std::string& keyword,
1219 const std::vector<NNCdata>& nncs,
1220 const std::function<KeywordLocation(const NNCdata&)>& getLocation,
1221 const std::function<void(Scalar&, const Scalar&)>& apply)
1222{
1223 if (nncs.empty()) {
1224 return;
1225 }
1226 const auto& cartDims = cartMapper_.cartesianDimensions();
1227
1228 auto format_ijk = [&cartDims](std::size_t cell) -> std::string
1229 {
1230 auto i = cell % cartDims[0]; cell /= cartDims[0];
1231 auto j = cell % cartDims[1];
1232 auto k = cell / cartDims[1];
1233
1234 return fmt::format("({},{},{})", i + 1,j + 1,k + 1);
1235 };
1236
1237 auto print_warning = [&format_ijk, &getLocation, &keyword] (const NNCdata& nnc)
1238 {
1239 const auto& location = getLocation( nnc );
1240 auto warning = fmt::format("Problem with {} keyword\n"
1241 "In {} line {} \n"
1242 "No NNC defined for connection {} -> {}", keyword, location.filename,
1243 location.lineno, format_ijk(nnc.cell1), format_ijk(nnc.cell2));
1244 OpmLog::warning(keyword, warning);
1245 };
1246
1247 // editNnc is supposed to only reference non-neighboring connections and not
1248 // neighboring connections. Use all entries for scaling if there is an NNC.
1249 // variable nnc incremented in loop body.
1250 auto nnc = nncs.begin();
1251 auto end = nncs.end();
1252 std::size_t warning_count = 0;
1253 while (nnc != end) {
1254 auto c1 = nnc->cell1;
1255 auto c2 = nnc->cell2;
1256 auto lowIt = globalToLocal.find(c1);
1257 auto highIt = globalToLocal.find(c2);
1258
1259 if (lowIt == globalToLocal.end() || highIt == globalToLocal.end()) {
1260 // Prevent warnings for NNCs stored on other processes in parallel (both cells inactive)
1261 if (lowIt != highIt && warnEditNNC_) {
1262 print_warning(*nnc);
1263 warning_count++;
1264 }
1265 ++nnc;
1266 continue;
1267 }
1268
1269 auto low = lowIt->second, high = highIt->second;
1270
1271 if (low > high) {
1272 std::swap(low, high);
1273 }
1274
1275 auto candidate = trans_.find(details::isId(low, high));
1276 if (candidate == trans_.end() && warnEditNNC_) {
1277 print_warning(*nnc);
1278 ++nnc;
1279 warning_count++;
1280 }
1281 else {
1282 // NNC exists
1283 while (nnc != end && c1 == nnc->cell1 && c2 == nnc->cell2) {
1284 apply(candidate->second, nnc->trans);
1285 ++nnc;
1286 }
1287 }
1288 }
1289
1290 if (warning_count > 0) {
1291 auto warning = fmt::format("Problems with {} keyword\n"
1292 "A total of {} connections not defined in grid", keyword, warning_count);
1293 OpmLog::warning(warning);
1294 }
1295}
1296
1297template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
1298void
1300applyNncMultreg_(const std::unordered_map<std::size_t,int>& cartesianToCompressed)
1301{
1302 const auto& inputNNC = this->eclState_.getInputNNC();
1303 const auto& transMult = this->eclState_.getTransMult();
1304
1305 auto compressedIdx = [&cartesianToCompressed](const std::size_t globIdx)
1306 {
1307 auto ixPos = cartesianToCompressed.find(globIdx);
1308 return (ixPos == cartesianToCompressed.end()) ? -1 : ixPos->second;
1309 };
1310
1311 // Apply region-based transmissibility multipliers (i.e., the MULTREGT
1312 // keyword) to those transmissibilities that are directly assigned from
1313 // the input.
1314 //
1315 // * NNC::input() covers the NNC keyword and any numerical aquifers
1316 // * NNC::editr() covers the EDITNNCR keyword
1317 //
1318 // Note: We do not apply MULTREGT to the entries in NNC::edit() since
1319 // those act as regular multipliers and have already been fully
1320 // accounted for in the multiplier part of the main loop of update() and
1321 // the applyEditNncToGridTrans_() member function.
1322 for (const auto& nncList : {&NNC::input, &NNC::editr}) {
1323 for (const auto& nncEntry : (inputNNC.*nncList)()) {
1324 const auto c1 = nncEntry.cell1;
1325 const auto c2 = nncEntry.cell2;
1326
1327 auto low = compressedIdx(c1);
1328 auto high = compressedIdx(c2);
1329
1330 if ((low == -1) || (high == -1)) {
1331 continue;
1332 }
1333
1334 if (low > high) {
1335 std::swap(low, high);
1336 }
1337
1338 auto candidate = this->trans_.find(details::isId(low, high));
1339 if (candidate != this->trans_.end()) {
1340 candidate->second *= transMult.getRegionMultiplierNNC(c1, c2);
1341 }
1342 }
1343 }
1344}
1345
1346template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
1347Scalar
1349computeHalfTrans_(const DimVector& areaNormal,
1350 int faceIdx, // in the reference element that contains the intersection
1351 const DimVector& distance,
1352 const DimMatrix& perm)
1353{
1354 assert(faceIdx >= 0);
1355 unsigned dimIdx = faceIdx / 2;
1356 assert(dimIdx < dimWorld);
1357 Scalar halfTrans = perm[dimIdx][dimIdx];
1358 halfTrans *= std::abs(Dune::dot(areaNormal, distance));
1359 halfTrans /= distance.two_norm2();
1360
1361 return halfTrans;
1362}
1363
1364template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
1365Scalar
1367computeHalfDiffusivity_(const DimVector& areaNormal,
1368 const DimVector& distance,
1369 const Scalar poro)
1370{
1371 Scalar halfDiff = poro;
1372 halfDiff *= std::abs(Dune::dot(areaNormal, distance));
1373 halfDiff /= distance.two_norm2();
1374
1375 return halfDiff;
1376}
1377
1378template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
1381distanceVector_(const DimVector& faceCenter,
1382 const unsigned& cellIdx) const
1383{
1384 const auto& cellCenter = centroids_cache_.empty() ? centroids_(cellIdx)
1385 : centroids_cache_[cellIdx];
1386 DimVector x = faceCenter;
1387 for (unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
1388 x[dimIdx] -= cellCenter[dimIdx];
1389 }
1390
1391 return x;
1392}
1393
1394template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
1396applyMultipliers_(Scalar& trans,
1397 unsigned faceIdx,
1398 unsigned cartElemIdx,
1399 const TransMult& transMult) const
1400{
1401 // apply multiplier for the transmissibility of the face. (the
1402 // face index is the index of the reference-element face which
1403 // contains the intersection of interest.)
1404 trans *= transMult.getMultiplier(cartElemIdx,
1405 FaceDir::FromIntersectionIndex(faceIdx));
1406}
1407
1408template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
1410applyNtg_(Scalar& trans,
1411 const FaceInfo& face,
1412 const std::vector<double>& ntg)
1413{
1414 // apply multiplier for the transmissibility of the face. (the
1415 // face index is the index of the reference-element face which
1416 // contains the intersection of interest.)
1417 // NTG does not apply to top and bottom faces
1418 if (face.faceIdx >= 0 && face.faceIdx <= 3) {
1419 trans *= ntg[face.elemIdx];
1420 }
1421}
1422
1423} // namespace Opm
1424
1425#endif // OPM_TRANSMISSIBILITY_IMPL_HPP
static unsigned maxThreads()
Return the maximum number of threads of the current process.
Definition: threadmanager.hpp:66
void applyMultipliers_(Scalar &trans, unsigned faceIdx, unsigned cartElemIdx, const TransMult &transMult) const
Definition: Transmissibility_impl.hpp:1396
Scalar diffusivity(unsigned elemIdx1, unsigned elemIdx2) const
Return the diffusivity for the intersection between two elements.
Definition: Transmissibility_impl.hpp:157
void applyEditNncToGridTrans_(const std::unordered_map< std::size_t, int > &globalToLocal)
Multiplies the grid transmissibilities according to EDITNNC.
Definition: Transmissibility_impl.hpp:1191
void computeFaceProperties(const Intersection &intersection, FaceInfo &inside, FaceInfo &outside, DimVector &faceAreaNormal, std::false_type) const
Definition: Transmissibility_impl.hpp:994
Dune::FieldMatrix< Scalar, dimWorld, dimWorld > DimMatrix
Definition: Transmissibility.hpp:59
void update(bool global, TransUpdateQuantities update_quantities=TransUpdateQuantities::All, const std::function< unsigned int(unsigned int)> &map={}, bool applyNncMultRegT=false)
Definition: Transmissibility_impl.hpp:177
DimVector distanceVector_(const DimVector &faceCenter, const unsigned &cellIdx) const
Definition: Transmissibility_impl.hpp:1381
void applyNncMultreg_(const std::unordered_map< std::size_t, int > &globalToLocal)
Definition: Transmissibility_impl.hpp:1300
void applyEditNncrToGridTrans_(const std::unordered_map< std::size_t, int > &globalToLocal)
Resets the grid transmissibilities according to EDITNNCR.
Definition: Transmissibility_impl.hpp:1204
void extractDispersion_()
Definition: Transmissibility_impl.hpp:715
const std::map< std::pair< unsigned, unsigned >, Scalar > & getThermalHalfTransBoundary() const
Definition: Transmissibility_impl.hpp:150
Scalar thermalHalfTrans(unsigned insideElemIdx, unsigned outsideElemIdx) const
Return the thermal "half transmissibility" for the intersection between two elements.
Definition: Transmissibility_impl.hpp:136
Scalar transmissibilityBoundary(unsigned elemIdx, unsigned boundaryFaceIdx) const
Return the transmissibility for a given boundary segment.
Definition: Transmissibility_impl.hpp:129
std::array< std::vector< double >, 3 > createTransmissibilityArrays_(const std::array< bool, 3 > &is_tran)
Creates TRANS{XYZ} arrays for modification by FieldProps data.
Definition: Transmissibility_impl.hpp:837
void applyEditNncToGridTransHelper_(const std::unordered_map< std::size_t, int > &globalToLocal, const std::string &keyword, const std::vector< NNCdata > &nncs, const std::function< KeywordLocation(const NNCdata &)> &getLocation, const std::function< void(Scalar &, const Scalar &)> &apply)
Definition: Transmissibility_impl.hpp:1217
void updateFromEclState_(bool global)
Definition: Transmissibility_impl.hpp:802
void extractPermeability_()
Definition: Transmissibility_impl.hpp:600
void applyAllZMultipliers_(Scalar &trans, const FaceInfo &inside, const FaceInfo &outside, const TransMult &transMult, const std::array< int, dimWorld > &cartDims)
Apply the Multipliers for the case PINCH(4)==TOPBOT.
Definition: Transmissibility_impl.hpp:759
TransUpdateQuantities
Compute all transmissibilities.
Definition: Transmissibility.hpp:159
void resetTransmissibilityFromArrays_(const std::array< bool, 3 > &is_tran, const std::array< std::vector< double >, 3 > &trans)
overwrites calculated transmissibilities
Definition: Transmissibility_impl.hpp:919
Dune::FieldVector< Scalar, dimWorld > DimVector
Definition: Transmissibility.hpp:60
static Scalar computeHalfDiffusivity_(const DimVector &areaNormal, const DimVector &distance, const Scalar poro)
Definition: Transmissibility_impl.hpp:1367
Scalar transmissibilityThreshold_
Definition: Transmissibility.hpp:291
Scalar transmissibility(unsigned elemIdx1, unsigned elemIdx2) const
Return the transmissibility for the intersection between two elements.
Definition: Transmissibility_impl.hpp:122
static void applyNtg_(Scalar &trans, const FaceInfo &face, const std::vector< double > &ntg)
Definition: Transmissibility_impl.hpp:1410
const EclipseState & eclState_
Definition: Transmissibility.hpp:285
Scalar dispersivity(unsigned elemIdx1, unsigned elemIdx2) const
Return the dispersivity for the intersection between two elements.
Definition: Transmissibility_impl.hpp:167
Transmissibility(const EclipseState &eclState, const GridView &gridView, const CartesianIndexMapper &cartMapper, const Grid &grid, std::function< std::array< double, dimWorld >(int)> centroids, bool enableEnergy, bool enableDiffusivity, bool enableDispersivity)
Definition: Transmissibility_impl.hpp:97
void applyNncToGridTrans_(const std::unordered_map< std::size_t, int > &cartesianToCompressed)
Definition: Transmissibility_impl.hpp:1124
void applyPinchNncToGridTrans_(const std::unordered_map< std::size_t, int > &cartesianToCompressed, bool applyNncMultregT)
Applies the previous calculate transmissibilities to the NNCs created via PINCH.
Definition: Transmissibility_impl.hpp:1075
void extractPorosity_()
Definition: Transmissibility_impl.hpp:690
void removeNonCartesianTransmissibilities_(bool removeAll)
Definition: Transmissibility_impl.hpp:734
static Scalar computeHalfTrans_(const DimVector &areaNormal, int faceIdx, const DimVector &distance, const DimMatrix &perm)
Definition: Transmissibility_impl.hpp:1349
Scalar thermalHalfTransBoundary(unsigned insideElemIdx, unsigned boundaryFaceIdx) const
Definition: Transmissibility_impl.hpp:143
constexpr unsigned elemIdxShift
Definition: FacePropertiesTPSA_impl.hpp:52
std::pair< std::uint32_t, std::uint32_t > isIdReverse(const std::uint64_t &id)
Definition: Transmissibility_impl.hpp:78
std::uint64_t directionalIsId(std::uint32_t elemIdx1, std::uint32_t elemIdx2)
Definition: Transmissibility_impl.hpp:89
std::uint64_t isId(std::uint32_t elemIdx1, std::uint32_t elemIdx2)
Definition: Transmissibility_impl.hpp:70
Definition: blackoilbioeffectsmodules.hh:45
Definition: Transmissibility.hpp:169
unsigned cartElemIdx
Definition: Transmissibility.hpp:173
DimVector faceCenter
Definition: Transmissibility.hpp:170
unsigned elemIdx
Definition: Transmissibility.hpp:172
int faceIdx
Definition: Transmissibility.hpp:171