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