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