23#ifndef OPM_TRANSMISSIBILITY_IMPL_HPP
24#define OPM_TRANSMISSIBILITY_IMPL_HPP
26#ifndef OPM_TRANSMISSIBILITY_HPP
31#include <dune/common/version.hh>
32#include <dune/grid/common/mcmgmapper.hh>
34#include <opm/common/OpmLog/KeywordLocation.hpp>
35#include <opm/common/utility/ThreadSafeMapBuilder.hpp>
37#include <opm/grid/CpGrid.hpp>
38#include <opm/grid/utility/ElementChunks.hpp>
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>
55#include <initializer_list>
62#include <fmt/format.h>
70 std::uint64_t
isId(std::uint32_t elemIdx1, std::uint32_t elemIdx2)
72 const std::uint32_t elemAIdx = std::min(elemIdx1, elemIdx2);
73 const std::uint64_t elemBIdx = std::max(elemIdx1, elemIdx2);
78 std::pair<std::uint32_t, std::uint32_t>
isIdReverse(
const std::uint64_t&
id)
83 const std::uint32_t elemAIdx =
static_cast<uint32_t
>(id);
84 const std::uint32_t elemBIdx = (
id - elemAIdx) >>
elemIdxShift;
86 return std::make_pair(elemAIdx, elemBIdx);
91 return (std::uint64_t(elemIdx1) <<
elemIdxShift) + elemIdx2;
95template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
98 const GridView& gridView,
99 const CartesianIndexMapper& cartMapper,
101 std::function<std::array<double,dimWorld>(
int)> centroids,
103 bool enableDiffusivity,
104 bool enableDispersivity)
105 : eclState_(eclState)
106 , gridView_(gridView)
107 , cartMapper_(cartMapper)
109 , centroids_(centroids)
110 , enableEnergy_(enableEnergy)
111 , enableDiffusivity_(enableDiffusivity)
112 , enableDispersivity_(enableDispersivity)
113 , lookUpData_(gridView)
114 , lookUpCartesianData_(gridView, cartMapper)
116 const UnitSystem& unitSystem =
eclState_.getDeckUnitSystem();
120template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
127template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
131 return transBoundary_.at(std::make_pair(elemIdx, boundaryFaceIdx));
134template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
141template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
145 return thermalHalfTransBoundary_.at(std::make_pair(insideElemIdx, boundaryFaceIdx));
148template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
152 return thermalHalfTransBoundary_;
155template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
157diffusivity(
unsigned elemIdx1,
unsigned elemIdx2)
const
159 if (diffusivity_.empty())
165template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
169 if (dispersivity_.empty())
175template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
178 const std::function<
unsigned int(
unsigned int)>& map,
const bool applyNncMultregT)
181 const bool onlyTrans = (update_quantities == TransUpdateQuantities::Trans);
182 const auto& cartDims = cartMapper_.cartesianDimensions();
183 const auto& transMult = eclState_.getTransMult();
184 const auto& comm = gridView_.comm();
185 ElementMapper elemMapper(gridView_, Dune::mcmgElementLayout());
187 unsigned numElements = elemMapper.size();
189 const std::vector<double>& ntg = this->lookUpData_.assignFieldPropsDoubleOnLeaf(eclState_.fieldProps(),
"NTG");
190 const bool updateDiffusivity = eclState_.getSimulationConfig().isDiffusive();
191 const bool updateDispersivity = eclState_.getSimulationConfig().rock_config().dispersion();
193 const bool disableNNC = eclState_.getSimulationConfig().useNONNC();
196 extractPermeability_(map);
199 extractPermeability_();
209 if (num_threads == 1) {
210 trans_.reserve(numElements*3*1.05);
213 transBoundary_.clear();
216 if ( enableEnergy_ && !onlyTrans) {
217 thermalHalfTrans_.clear();
218 if (num_threads == 1) {
219 thermalHalfTrans_.reserve(numElements*6*1.05);
222 thermalHalfTransBoundary_.clear();
226 if (updateDiffusivity && !onlyTrans) {
227 diffusivity_.clear();
228 if (num_threads == 1) {
229 diffusivity_.reserve(numElements*3*1.05);
235 if (updateDispersivity && !onlyTrans) {
236 dispersivity_.clear();
237 if (num_threads == 1) {
238 dispersivity_.reserve(numElements*3*1.05);
240 extractDispersion_();
246 bool useSmallestMultiplier;
247 bool pinchOption4ALL;
249 if (comm.rank() == 0) {
250 const auto& eclGrid = eclState_.getInputGrid();
251 pinchActive = eclGrid.isPinchActive();
252 auto pinchTransCalcMode = eclGrid.getPinchOption();
253 useSmallestMultiplier = eclGrid.getMultzOption() == PinchMode::ALL;
254 pinchOption4ALL = (pinchTransCalcMode == PinchMode::ALL);
255 if (pinchOption4ALL) {
256 useSmallestMultiplier =
false;
259 if (global && comm.size() > 1) {
260 comm.broadcast(&useSmallestMultiplier, 1, 0);
261 comm.broadcast(&pinchOption4ALL, 1, 0);
262 comm.broadcast(&pinchActive, 1, 0);
266 centroids_cache_.resize(gridView_.size(0));
267 for (
const auto& elem : elements(gridView_)) {
268 const unsigned elemIdx = elemMapper.index(elem);
269 centroids_cache_[elemIdx] = centroids_(elemIdx);
272 auto harmonicMean = [](
const Scalar x1,
const Scalar x2)
274 return (std::abs(x1) < 1e-30 || std::abs(x2) < 1e-30)
276 : 1.0 / (1.0 / x1 + 1.0 / x2);
279 auto faceIdToDir = [](
int insideFaceIdx)
281 switch (insideFaceIdx) {
284 return FaceDir::XPlus;
287 return FaceDir::YPlus;
291 return FaceDir::ZPlus;
293 throw std::logic_error(
"Could not determine a face direction");
297 auto halfDiff = [](
const DimVector& faceAreaNormal,
302 return computeHalfDiffusivity_(faceAreaNormal,
307 ThreadSafeMapBuilder transBoundary(transBoundary_, num_threads,
308 MapBuilderInsertionMode::Insert_Or_Assign);
309 ThreadSafeMapBuilder transMap(trans_, num_threads,
310 MapBuilderInsertionMode::Insert_Or_Assign);
311 ThreadSafeMapBuilder thermalHalfTransBoundary(thermalHalfTransBoundary_, num_threads,
312 MapBuilderInsertionMode::Insert_Or_Assign);
313 ThreadSafeMapBuilder thermalHalfTrans(thermalHalfTrans_, num_threads,
314 MapBuilderInsertionMode::Insert_Or_Assign);
315 ThreadSafeMapBuilder diffusivity(diffusivity_, num_threads,
316 MapBuilderInsertionMode::Insert_Or_Assign);
317 ThreadSafeMapBuilder dispersivity(dispersivity_, num_threads,
318 MapBuilderInsertionMode::Insert_Or_Assign);
320 const auto& nnc_input = eclState_.getInputNNC().input();
323#pragma omp parallel for
325 for (
const auto& chunk : ElementChunks(gridView_, Dune::Partitions::all, num_threads)) {
326 for (
const auto& elem : chunk) {
331 inside.
elemIdx = elemMapper.index(elem);
335 template getFieldPropCartesianIdx<Grid>(inside.
elemIdx);
337 auto computeHalf = [
this, &faceAreaNormal, &inside, &outside]
338 (
const auto& halfComputer,
339 const auto& prop1,
const auto& prop2) -> std::array<Scalar,2>
342 halfComputer(faceAreaNormal,
346 halfComputer(faceAreaNormal,
353 auto computeHalfMean = [&inside, &outside, &computeHalf, &ntg, &harmonicMean]
354 (
const auto& halfComputer,
const auto& prop)
356 auto onesided = computeHalf(halfComputer, prop[inside.
elemIdx], prop[outside.
elemIdx]);
357 applyNtg_(onesided[0], inside, ntg);
358 applyNtg_(onesided[1], outside, ntg);
361 return harmonicMean(onesided[0], onesided[1]);
364 unsigned boundaryIsIdx = 0;
365 for (
const auto& intersection : intersections(gridView_, elem)) {
367 if (intersection.boundary()) {
369 const auto& geometry = intersection.geometry();
372 faceAreaNormal = intersection.centerUnitOuterNormal();
373 faceAreaNormal *= geometry.volume();
375 Scalar transBoundaryIs =
376 computeHalfTrans_(faceAreaNormal,
377 intersection.indexInInside(),
379 permeability_[inside.
elemIdx]);
384 applyMultipliers_(transBoundaryIs, intersection.indexInInside(), inside.
cartElemIdx, transMult);
385 transBoundary.insert_or_assign(std::make_pair(inside.
elemIdx, boundaryIsIdx), transBoundaryIs);
389 if (enableEnergy_ && !onlyTrans) {
390 Scalar transBoundaryEnergyIs =
391 computeHalfDiffusivity_(faceAreaNormal,
394 thermalHalfTransBoundary.insert_or_assign(std::make_pair(inside.
elemIdx, boundaryIsIdx),
395 transBoundaryEnergyIs);
402 if (!intersection.neighbor()) {
409 const auto& outsideElem = intersection.outside();
410 outside.
elemIdx = elemMapper.index(outsideElem);
415 template getFieldPropCartesianIdx<Grid>(outside.
elemIdx);
429 inside.
faceIdx = intersection.indexInInside();
430 outside.
faceIdx = intersection.indexInOutside();
437 if (enableEnergy_ && !onlyTrans) {
442 if (updateDiffusivity && !onlyTrans) {
445 if (updateDispersivity && !onlyTrans) {
451 typename std::is_same<Grid, Dune::CpGrid>::type isCpGrid;
452 computeFaceProperties(intersection,
458 Scalar trans = computeHalfMean(computeHalfTrans_, permeability_);
464 auto find_layer = [&cartDims](std::size_t cell) {
466 auto k = cell / cartDims[1];
474 if (std::abs(kup -kdown) > 1) {
480 if (useSmallestMultiplier) {
484 applyAllZMultipliers_(trans, inside, outside, transMult, cartDims);
492 bool foundInputNNC =
false;
493 if (! nnc_input.empty()) {
495 auto it = std::lower_bound(nnc_input.begin(), nnc_input.end(),
496 NNCdata { inside.cartElemIdx, outside.cartElemIdx, 0.0 });
497 foundInputNNC = it != nnc_input.end() && it->cell1 == inside.
cartElemIdx && it->cell2 == outside.
cartElemIdx;
499 if (! foundInputNNC) {
501 trans *= transMult.getRegionMultiplier(inside.
cartElemIdx,
509 if (enableEnergy_ && !onlyTrans) {
510 const auto half = computeHalf(halfDiff, 1.0, 1.0);
519 if (updateDiffusivity && !onlyTrans) {
521 computeHalfMean(halfDiff, porosity_));
525 if (updateDispersivity && !onlyTrans) {
527 computeHalfMean(halfDiff, dispersion_));
532 centroids_cache_.clear();
535#pragma omp parallel sections
545 transBoundary.finalize();
549 thermalHalfTransBoundary.finalize();
553 thermalHalfTrans.finalize();
557 diffusivity.finalize();
561 dispersivity.finalize();
565 this->updateFromEclState_(global);
568 std::unordered_map<std::size_t,int> globalToLocal;
571 for (
const auto& elem : elements(grid_.leafGridView())) {
572 int elemIdx = elemMapper.index(elem);
573 int cartElemIdx = cartMapper_.cartesianIndex(elemIdx);
574 globalToLocal[cartElemIdx] = elemIdx;
583 this->applyPinchNncToGridTrans_(globalToLocal, applyNncMultregT);
584 this->applyNncToGridTrans_(globalToLocal);
585 this->applyEditNncToGridTrans_(globalToLocal);
586 this->applyEditNncrToGridTrans_(globalToLocal);
587 if (applyNncMultregT) {
588 this->applyNncMultreg_(globalToLocal);
590 warnEditNNC_ =
false;
595 this->removeNonCartesianTransmissibilities_(disableNNC);
598template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
602 unsigned numElem = gridView_.size(0);
603 permeability_.resize(numElem);
609 const auto& fp = eclState_.fieldProps();
610 if (fp.has_double(
"PERMX")) {
611 const std::vector<double>& permxData = this-> lookUpData_.assignFieldPropsDoubleOnLeaf(fp,
"PERMX");
613 std::vector<double> permyData;
614 if (fp.has_double(
"PERMY"))
615 permyData = this-> lookUpData_.assignFieldPropsDoubleOnLeaf(fp,
"PERMY");
617 permyData = permxData;
619 std::vector<double> permzData;
620 if (fp.has_double(
"PERMZ"))
621 permzData = this-> lookUpData_.assignFieldPropsDoubleOnLeaf(fp,
"PERMZ");
623 permzData = permxData;
625 for (std::size_t dofIdx = 0; dofIdx < numElem; ++ dofIdx) {
626 permeability_[dofIdx] = 0.0;
627 permeability_[dofIdx][0][0] = permxData[dofIdx];
628 permeability_[dofIdx][1][1] = permyData[dofIdx];
629 permeability_[dofIdx][2][2] = permzData[dofIdx];
636 throw std::logic_error(
"Can't read the intrinsic permeability from the ecl state. "
637 "(The PERM{X,Y,Z} keywords are missing)");
640template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
644 unsigned numElem = gridView_.size(0);
645 permeability_.resize(numElem);
651 const auto& fp = eclState_.fieldProps();
652 if (fp.has_double(
"PERMX")) {
653 const std::vector<double>& permxData =
654 this->lookUpData_.assignFieldPropsDoubleOnLeaf(fp,
"PERMX");
656 std::vector<double> permyData;
657 if (fp.has_double(
"PERMY")){
658 permyData = this->lookUpData_.assignFieldPropsDoubleOnLeaf(fp,
"PERMY");
661 permyData = permxData;
664 std::vector<double> permzData;
665 if (fp.has_double(
"PERMZ")) {
666 permzData = this->lookUpData_.assignFieldPropsDoubleOnLeaf(fp,
"PERMZ");
669 permzData = permxData;
672 for (std::size_t dofIdx = 0; dofIdx < numElem; ++ dofIdx) {
673 permeability_[dofIdx] = 0.0;
674 std::size_t inputDofIdx = map(dofIdx);
675 permeability_[dofIdx][0][0] = permxData[inputDofIdx];
676 permeability_[dofIdx][1][1] = permyData[inputDofIdx];
677 permeability_[dofIdx][2][2] = permzData[inputDofIdx];
683 throw std::logic_error(
"Can't read the intrinsic permeability from the ecl state. "
684 "(The PERM{X,Y,Z} keywords are missing)");
688template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
696 const auto& fp = eclState_.fieldProps();
697 if (fp.has_double(
"PORO")) {
698 if constexpr (std::is_same_v<Scalar,double>) {
699 porosity_ = this->lookUpData_.assignFieldPropsDoubleOnLeaf(fp,
"PORO");
702 const auto por = this->lookUpData_.assignFieldPropsDoubleOnLeaf(fp,
"PORO");
703 porosity_.resize(por.size());
704 std::ranges::copy(por, porosity_.begin());
708 throw std::logic_error(
"Can't read the porosity from the ecl state. "
709 "(The PORO keywords are missing)");
713template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
717 if (!enableDispersivity_) {
718 throw std::runtime_error(
"Dispersion disabled at compile time, but the deck "
719 "contains the DISPERC keyword.");
721 const auto& fp = eclState_.fieldProps();
722 if constexpr (std::is_same_v<Scalar,double>) {
723 dispersion_ = this->lookUpData_.assignFieldPropsDoubleOnLeaf(fp,
"DISPERC");
726 const auto disp = this->lookUpData_.assignFieldPropsDoubleOnLeaf(fp,
"DISPERC");
727 dispersion_.resize(disp.size());
728 std::ranges::copy(disp, dispersion_.begin());
732template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
736 const auto& cartDims = cartMapper_.cartesianDimensions();
737 for (
auto&& trans: trans_) {
739 if (removeAll || trans.second < transmissibilityThreshold_) {
740 const auto&
id = trans.first;
742 int gc1 = std::min(cartMapper_.cartesianIndex(elements.first), cartMapper_.cartesianIndex(elements.second));
743 int gc2 = std::max(cartMapper_.cartesianIndex(elements.first), cartMapper_.cartesianIndex(elements.second));
748 if (gc2 - gc1 == 1 || gc2 - gc1 == cartDims[0] || gc2 - gc1 == cartDims[0]*cartDims[1] || gc2 - gc1 == 0) {
757template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
762 const TransMult& transMult,
763 const std::array<int, dimWorld>& cartDims)
765 if (grid_.maxLevel() > 0) {
766 OPM_THROW(std::invalid_argument,
"MULTZ not support with LGRS, yet.");
772 unsigned lastCartElemIdx;
777 lastCartElemIdx = outside.
cartElemIdx - cartDims[0]*cartDims[1];
780 Scalar mult = transMult.getMultiplier(lastCartElemIdx , FaceDir::ZPlus) *
781 transMult.getMultiplier(outside.
cartElemIdx , FaceDir::ZMinus);
785 for (
auto cartElemIdx = inside.
cartElemIdx; cartElemIdx < lastCartElemIdx;) {
786 auto multiplier = transMult.getMultiplier(cartElemIdx, FaceDir::ZPlus);
787 cartElemIdx += cartDims[0]*cartDims[1];
788 multiplier *= transMult.getMultiplier(cartElemIdx, FaceDir::ZMinus);
789 mult = std::min(mult,
static_cast<Scalar
>(multiplier));
800template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
804 const FieldPropsManager* fp =
805 (global) ? &(eclState_.fieldProps()) :
806 &(eclState_.globalFieldProps());
808 std::array<bool,3> is_tran {fp->tran_active(
"TRANX"),
809 fp->tran_active(
"TRANY"),
810 fp->tran_active(
"TRANZ")};
812 if (!(is_tran[0] || is_tran[1] || is_tran[2])) {
817 std::array<std::string, 3> keywords {
"TRANX",
"TRANY",
"TRANZ"};
818 std::array<std::vector<double>,3> trans = createTransmissibilityArrays_(is_tran);
819 auto key = keywords.begin();
820 auto perform = is_tran.begin();
822 for (
auto it = trans.begin(); it != trans.end(); ++it, ++key, ++perform) {
824 if (grid_.maxLevel() > 0) {
825 OPM_THROW(std::invalid_argument,
"Calculations on TRANX/TRANY/TRANZ arrays are not support with LGRS, yet.");
827 fp->apply_tran(*key, *it);
831 resetTransmissibilityFromArrays_(is_tran, trans);
834template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
835std::array<std::vector<double>,3>
839 const auto& cartDims = cartMapper_.cartesianDimensions();
840 ElementMapper elemMapper(gridView_, Dune::mcmgElementLayout());
842 auto numElem = gridView_.size(0);
843 std::array<std::vector<double>,3> trans = {
844 std::vector<double>(is_tran[0] ? numElem : 0, 0),
845 std::vector<double>(is_tran[1] ? numElem : 0, 0),
846 std::vector<double>(is_tran[2] ? numElem : 0, 0)
850 for (
const auto& elem : elements(gridView_)) {
851 for (
const auto& intersection : intersections(gridView_, elem)) {
853 if (!intersection.neighbor()) {
868 unsigned c1 = elemMapper.index(intersection.inside());
869 unsigned c2 = elemMapper.index(intersection.outside());
870 int gc1 = cartMapper_.cartesianIndex(c1);
871 int gc2 = cartMapper_.cartesianIndex(c2);
872 if (std::tie(gc1, c1) > std::tie(gc2, c2)) {
886 if ((gc2 - gc1 == 1 || (gc2 == gc1 && (intersection.indexInInside() == 0 || intersection.indexInInside() == 1)))
891 trans[0][c1] = trans_[isID];
894 else if ((gc2 - gc1 == cartDims[0] || (gc2 == gc1 && (intersection.indexInInside() == 2 || intersection.indexInInside() == 3)))
899 trans[1][c1] = trans_[isID];
902 else if (gc2 - gc1 == cartDims[0]*cartDims[1] ||
903 (gc2 == gc1 && (intersection.indexInInside() == 4 || intersection.indexInInside() == 5)))
907 trans[2][c1] = trans_[isID];
917template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
920 const std::array<std::vector<double>,3>& trans)
922 const auto& cartDims = cartMapper_.cartesianDimensions();
923 ElementMapper elemMapper(gridView_, Dune::mcmgElementLayout());
926 for (
const auto& elem : elements(gridView_)) {
927 for (
const auto& intersection : intersections(gridView_, elem)) {
928 if (!intersection.neighbor()) {
943 unsigned c1 = elemMapper.index(intersection.inside());
944 unsigned c2 = elemMapper.index(intersection.outside());
945 int gc1 = cartMapper_.cartesianIndex(c1);
946 int gc2 = cartMapper_.cartesianIndex(c2);
947 if (std::tie(gc1, c1) > std::tie(gc2, c2)) {
961 if ((gc2 - gc1 == 1 || (gc2 == gc1 && (intersection.indexInInside() == 0 || intersection.indexInInside() == 1)))
966 trans_[isID] = trans[0][c1];
969 else if ((gc2 - gc1 == cartDims[0] || (gc2 == gc1 && (intersection.indexInInside() == 2|| intersection.indexInInside() == 3)))
974 trans_[isID] = trans[1][c1];
977 else if (gc2 - gc1 == cartDims[0]*cartDims[1] ||
978 (gc2 == gc1 && (intersection.indexInInside() == 4 || intersection.indexInInside() == 5)))
982 trans_[isID] = trans[2][c1];
991template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
992template<
class Intersection>
998 std::false_type)
const
1001 const auto& geometry = intersection.geometry();
1004 faceAreaNormal = intersection.centerUnitOuterNormal();
1005 faceAreaNormal *= geometry.volume();
1008template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
1009template<
class Intersection>
1015 std::true_type)
const
1017 int faceIdx = intersection.id();
1019 if (grid_.maxLevel() == 0) {
1022 faceAreaNormal = grid_.faceAreaNormalEcl(faceIdx);
1025 if ((intersection.inside().level() != intersection.outside().level())) {
1032 const auto& parentIntersection =
1033 grid_.getParentIntersectionFromLgrBoundaryFace(intersection);
1034 const auto& parentIntersectionGeometry = parentIntersection.geometry();
1038 inside.
faceCenter = (intersection.inside().level() == 0)
1039 ? parentIntersectionGeometry.center()
1040 : grid_.faceCenterEcl(inside.
elemIdx, inside.
faceIdx, intersection);
1041 outside.
faceCenter = (intersection.outside().level() == 0)
1042 ? parentIntersectionGeometry.center()
1043 : grid_.faceCenterEcl(outside.
elemIdx, outside.
faceIdx, intersection);
1051 faceAreaNormal = intersection.centerUnitOuterNormal();
1052 faceAreaNormal *= intersection.geometry().volume();
1055 assert(intersection.inside().level() == intersection.outside().level());
1061 if (intersection.inside().level() > 0) {
1062 faceAreaNormal = intersection.centerUnitOuterNormal();
1063 faceAreaNormal *= intersection.geometry().volume();
1066 faceAreaNormal = grid_.faceAreaNormalEcl(faceIdx);
1072template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
1076 const bool applyNncMultregT)
1078 const auto& pinchNnc = eclState_.getPinchNNC();
1079 const auto& transMult = this->eclState_.getTransMult();
1081 for (
const auto& nncEntry : pinchNnc) {
1082 auto c1 = nncEntry.cell1;
1083 auto c2 = nncEntry.cell2;
1084 auto lowIt = cartesianToCompressed.find(c1);
1085 auto highIt = cartesianToCompressed.find(c2);
1086 int low = (lowIt == cartesianToCompressed.end())? -1 : lowIt->second;
1087 int high = (highIt == cartesianToCompressed.end())? -1 : highIt->second;
1090 std::swap(low, high);
1093 if (low == -1 && high == -1) {
1098 if (low == -1 || high == -1) {
1107 if (candidate != trans_.end()) {
1111 candidate->second = nncEntry.trans;
1112 if (applyNncMultregT) {
1113 const auto mult = transMult.getRegionMultiplierNNC(c1, c2);
1114 candidate->second *= mult;
1121template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
1127 const auto& nnc_input = eclState_.getInputNNC().input();
1129 for (
const auto& nncEntry : nnc_input) {
1130 auto c1 = nncEntry.cell1;
1131 auto c2 = nncEntry.cell2;
1132 auto lowIt = cartesianToCompressed.find(c1);
1133 auto highIt = cartesianToCompressed.find(c2);
1134 int low = (lowIt == cartesianToCompressed.end())? -1 : lowIt->second;
1135 int high = (highIt == cartesianToCompressed.end())? -1 : highIt->second;
1138 std::swap(low, high);
1141 if (low == -1 && high == -1) {
1146 if (low == -1 || high == -1) {
1148 std::ostringstream sstr;
1149 sstr <<
"NNC between active and inactive cells ("
1150 << low <<
" -> " << high <<
") with globalcell is (" << c1 <<
"->" << c2 <<
")";
1151 OpmLog::warning(sstr.str());
1155 if (
auto candidate = trans_.find(
details::isId(low, high)); candidate != trans_.end()) {
1159 candidate->second += nncEntry.trans;
1189template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
1193 const auto& input = eclState_.getInputNNC();
1194 applyEditNncToGridTransHelper_(globalToLocal,
"EDITNNC",
1196 [&input](
const NNCdata& nnc){
1197 return input.edit_location(nnc);},
1199 [](Scalar& trans,
const Scalar& rhs){ trans *= rhs;});
1202template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
1206 const auto& input = eclState_.getInputNNC();
1207 applyEditNncToGridTransHelper_(globalToLocal,
"EDITNNCR",
1209 [&input](
const NNCdata& nnc){
1210 return input.editr_location(nnc);},
1212 [](Scalar& trans,
const Scalar& rhs){ trans = rhs;});
1215template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
1218 const std::string& keyword,
1219 const std::vector<NNCdata>& nncs,
1220 const std::function<KeywordLocation(
const NNCdata&)>& getLocation,
1221 const std::function<
void(Scalar&,
const Scalar&)>& apply)
1226 const auto& cartDims = cartMapper_.cartesianDimensions();
1228 auto format_ijk = [&cartDims](std::size_t cell) -> std::string
1230 auto i = cell % cartDims[0]; cell /= cartDims[0];
1231 auto j = cell % cartDims[1];
1232 auto k = cell / cartDims[1];
1234 return fmt::format(
"({},{},{})", i + 1,j + 1,k + 1);
1237 auto print_warning = [&format_ijk, &getLocation, &keyword] (
const NNCdata& nnc)
1239 const auto& location = getLocation( nnc );
1240 auto warning = fmt::format(
"Problem with {} keyword\n"
1242 "No NNC defined for connection {} -> {}", keyword, location.filename,
1243 location.lineno, format_ijk(nnc.cell1), format_ijk(nnc.cell2));
1244 OpmLog::warning(keyword, warning);
1250 auto nnc = nncs.begin();
1251 auto end = nncs.end();
1252 std::size_t warning_count = 0;
1253 while (nnc != end) {
1254 auto c1 = nnc->cell1;
1255 auto c2 = nnc->cell2;
1256 auto lowIt = globalToLocal.find(c1);
1257 auto highIt = globalToLocal.find(c2);
1259 if (lowIt == globalToLocal.end() || highIt == globalToLocal.end()) {
1261 if (lowIt != highIt && warnEditNNC_) {
1262 print_warning(*nnc);
1269 auto low = lowIt->second, high = highIt->second;
1272 std::swap(low, high);
1276 if (candidate == trans_.end() && warnEditNNC_) {
1277 print_warning(*nnc);
1283 while (nnc != end && c1 == nnc->cell1 && c2 == nnc->cell2) {
1284 apply(candidate->second, nnc->trans);
1290 if (warning_count > 0) {
1291 auto warning = fmt::format(
"Problems with {} keyword\n"
1292 "A total of {} connections not defined in grid", keyword, warning_count);
1293 OpmLog::warning(warning);
1297template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
1300applyNncMultreg_(
const std::unordered_map<std::size_t,int>& cartesianToCompressed)
1302 const auto& inputNNC = this->eclState_.getInputNNC();
1303 const auto& transMult = this->eclState_.getTransMult();
1305 auto compressedIdx = [&cartesianToCompressed](
const std::size_t globIdx)
1307 auto ixPos = cartesianToCompressed.find(globIdx);
1308 return (ixPos == cartesianToCompressed.end()) ? -1 : ixPos->second;
1322 for (
const auto& nncList : {&NNC::input, &NNC::editr}) {
1323 for (
const auto& nncEntry : (inputNNC.*nncList)()) {
1324 const auto c1 = nncEntry.cell1;
1325 const auto c2 = nncEntry.cell2;
1327 auto low = compressedIdx(c1);
1328 auto high = compressedIdx(c2);
1330 if ((low == -1) || (high == -1)) {
1335 std::swap(low, high);
1338 auto candidate = this->trans_.find(
details::isId(low, high));
1339 if (candidate != this->trans_.end()) {
1340 candidate->second *= transMult.getRegionMultiplierNNC(c1, c2);
1346template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
1354 assert(faceIdx >= 0);
1355 unsigned dimIdx = faceIdx / 2;
1356 assert(dimIdx < dimWorld);
1357 Scalar halfTrans = perm[dimIdx][dimIdx];
1358 halfTrans *= std::abs(Dune::dot(areaNormal, distance));
1359 halfTrans /= distance.two_norm2();
1364template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
1371 Scalar halfDiff = poro;
1372 halfDiff *= std::abs(Dune::dot(areaNormal, distance));
1373 halfDiff /= distance.two_norm2();
1378template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
1382 const unsigned& cellIdx)
const
1384 const auto& cellCenter = centroids_cache_.empty() ? centroids_(cellIdx)
1385 : centroids_cache_[cellIdx];
1387 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
1388 x[dimIdx] -= cellCenter[dimIdx];
1394template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
1398 unsigned cartElemIdx,
1399 const TransMult& transMult)
const
1404 trans *= transMult.getMultiplier(cartElemIdx,
1405 FaceDir::FromIntersectionIndex(faceIdx));
1408template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
1412 const std::vector<double>& ntg)
static unsigned maxThreads()
Return the maximum number of threads of the current process.
Definition: threadmanager.hpp:66
void applyMultipliers_(Scalar &trans, unsigned faceIdx, unsigned cartElemIdx, const TransMult &transMult) const
Definition: Transmissibility_impl.hpp:1396
Scalar diffusivity(unsigned elemIdx1, unsigned elemIdx2) const
Return the diffusivity for the intersection between two elements.
Definition: Transmissibility_impl.hpp:157
void applyEditNncToGridTrans_(const std::unordered_map< std::size_t, int > &globalToLocal)
Multiplies the grid transmissibilities according to EDITNNC.
Definition: Transmissibility_impl.hpp:1191
void computeFaceProperties(const Intersection &intersection, FaceInfo &inside, FaceInfo &outside, DimVector &faceAreaNormal, std::false_type) const
Definition: Transmissibility_impl.hpp:994
Dune::FieldMatrix< Scalar, dimWorld, dimWorld > DimMatrix
Definition: Transmissibility.hpp:59
void update(bool global, TransUpdateQuantities update_quantities=TransUpdateQuantities::All, const std::function< unsigned int(unsigned int)> &map={}, bool applyNncMultRegT=false)
Definition: Transmissibility_impl.hpp:177
DimVector distanceVector_(const DimVector &faceCenter, const unsigned &cellIdx) const
Definition: Transmissibility_impl.hpp:1381
void applyNncMultreg_(const std::unordered_map< std::size_t, int > &globalToLocal)
Definition: Transmissibility_impl.hpp:1300
void applyEditNncrToGridTrans_(const std::unordered_map< std::size_t, int > &globalToLocal)
Resets the grid transmissibilities according to EDITNNCR.
Definition: Transmissibility_impl.hpp:1204
void extractDispersion_()
Definition: Transmissibility_impl.hpp:715
const std::map< std::pair< unsigned, unsigned >, Scalar > & getThermalHalfTransBoundary() const
Definition: Transmissibility_impl.hpp:150
Scalar thermalHalfTrans(unsigned insideElemIdx, unsigned outsideElemIdx) const
Return the thermal "half transmissibility" for the intersection between two elements.
Definition: Transmissibility_impl.hpp:136
Scalar transmissibilityBoundary(unsigned elemIdx, unsigned boundaryFaceIdx) const
Return the transmissibility for a given boundary segment.
Definition: Transmissibility_impl.hpp:129
std::array< std::vector< double >, 3 > createTransmissibilityArrays_(const std::array< bool, 3 > &is_tran)
Creates TRANS{XYZ} arrays for modification by FieldProps data.
Definition: Transmissibility_impl.hpp:837
void applyEditNncToGridTransHelper_(const std::unordered_map< std::size_t, int > &globalToLocal, const std::string &keyword, const std::vector< NNCdata > &nncs, const std::function< KeywordLocation(const NNCdata &)> &getLocation, const std::function< void(Scalar &, const Scalar &)> &apply)
Definition: Transmissibility_impl.hpp:1217
void updateFromEclState_(bool global)
Definition: Transmissibility_impl.hpp:802
void extractPermeability_()
Definition: Transmissibility_impl.hpp:600
void applyAllZMultipliers_(Scalar &trans, const FaceInfo &inside, const FaceInfo &outside, const TransMult &transMult, const std::array< int, dimWorld > &cartDims)
Apply the Multipliers for the case PINCH(4)==TOPBOT.
Definition: Transmissibility_impl.hpp:759
TransUpdateQuantities
Compute all transmissibilities.
Definition: Transmissibility.hpp:159
void resetTransmissibilityFromArrays_(const std::array< bool, 3 > &is_tran, const std::array< std::vector< double >, 3 > &trans)
overwrites calculated transmissibilities
Definition: Transmissibility_impl.hpp:919
Dune::FieldVector< Scalar, dimWorld > DimVector
Definition: Transmissibility.hpp:60
static Scalar computeHalfDiffusivity_(const DimVector &areaNormal, const DimVector &distance, const Scalar poro)
Definition: Transmissibility_impl.hpp:1367
Scalar transmissibilityThreshold_
Definition: Transmissibility.hpp:291
Scalar transmissibility(unsigned elemIdx1, unsigned elemIdx2) const
Return the transmissibility for the intersection between two elements.
Definition: Transmissibility_impl.hpp:122
static void applyNtg_(Scalar &trans, const FaceInfo &face, const std::vector< double > &ntg)
Definition: Transmissibility_impl.hpp:1410
const EclipseState & eclState_
Definition: Transmissibility.hpp:285
Scalar dispersivity(unsigned elemIdx1, unsigned elemIdx2) const
Return the dispersivity for the intersection between two elements.
Definition: Transmissibility_impl.hpp:167
Transmissibility(const EclipseState &eclState, const GridView &gridView, const CartesianIndexMapper &cartMapper, const Grid &grid, std::function< std::array< double, dimWorld >(int)> centroids, bool enableEnergy, bool enableDiffusivity, bool enableDispersivity)
Definition: Transmissibility_impl.hpp:97
void applyNncToGridTrans_(const std::unordered_map< std::size_t, int > &cartesianToCompressed)
Definition: Transmissibility_impl.hpp:1124
void applyPinchNncToGridTrans_(const std::unordered_map< std::size_t, int > &cartesianToCompressed, bool applyNncMultregT)
Applies the previous calculate transmissibilities to the NNCs created via PINCH.
Definition: Transmissibility_impl.hpp:1075
void extractPorosity_()
Definition: Transmissibility_impl.hpp:690
void removeNonCartesianTransmissibilities_(bool removeAll)
Definition: Transmissibility_impl.hpp:734
static Scalar computeHalfTrans_(const DimVector &areaNormal, int faceIdx, const DimVector &distance, const DimMatrix &perm)
Definition: Transmissibility_impl.hpp:1349
Scalar thermalHalfTransBoundary(unsigned insideElemIdx, unsigned boundaryFaceIdx) const
Definition: Transmissibility_impl.hpp:143
constexpr unsigned elemIdxShift
Definition: FacePropertiesTPSA_impl.hpp:52
std::pair< std::uint32_t, std::uint32_t > isIdReverse(const std::uint64_t &id)
Definition: Transmissibility_impl.hpp:78
std::uint64_t directionalIsId(std::uint32_t elemIdx1, std::uint32_t elemIdx2)
Definition: Transmissibility_impl.hpp:89
std::uint64_t isId(std::uint32_t elemIdx1, std::uint32_t elemIdx2)
Definition: Transmissibility_impl.hpp:70
Definition: blackoilbioeffectsmodules.hh:45
Definition: Transmissibility.hpp:169
unsigned cartElemIdx
Definition: Transmissibility.hpp:173
DimVector faceCenter
Definition: Transmissibility.hpp:170
unsigned elemIdx
Definition: Transmissibility.hpp:172
int faceIdx
Definition: Transmissibility.hpp:171