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>
150diffusivity(
unsigned elemIdx1,
unsigned elemIdx2)
const
152 if (diffusivity_.empty())
158template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
162 if (dispersivity_.empty())
168template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
171 const std::function<
unsigned int(
unsigned int)>& map,
const bool applyNncMultregT)
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());
180 unsigned numElements = elemMapper.size();
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();
186 const bool disableNNC = eclState_.getSimulationConfig().useNONNC();
189 extractPermeability_(map);
192 extractPermeability_();
202 if (num_threads == 1) {
203 trans_.reserve(numElements*3*1.05);
206 transBoundary_.clear();
209 if ( enableEnergy_ && !onlyTrans) {
210 thermalHalfTrans_.clear();
211 if (num_threads == 1) {
212 thermalHalfTrans_.reserve(numElements*6*1.05);
215 thermalHalfTransBoundary_.clear();
219 if (updateDiffusivity && !onlyTrans) {
220 diffusivity_.clear();
221 if (num_threads == 1) {
222 diffusivity_.reserve(numElements*3*1.05);
228 if (updateDispersivity && !onlyTrans) {
229 dispersivity_.clear();
230 if (num_threads == 1) {
231 dispersivity_.reserve(numElements*3*1.05);
233 extractDispersion_();
239 bool useSmallestMultiplier;
240 bool pinchOption4ALL;
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;
252 if (global && comm.size() > 1) {
253 comm.broadcast(&useSmallestMultiplier, 1, 0);
254 comm.broadcast(&pinchOption4ALL, 1, 0);
255 comm.broadcast(&pinchActive, 1, 0);
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);
265 auto harmonicMean = [](
const Scalar x1,
const Scalar x2)
267 return (std::abs(x1) < 1e-30 || std::abs(x2) < 1e-30)
269 : 1.0 / (1.0 / x1 + 1.0 / x2);
272 auto faceIdToDir = [](
int insideFaceIdx)
274 switch (insideFaceIdx) {
277 return FaceDir::XPlus;
280 return FaceDir::YPlus;
284 return FaceDir::ZPlus;
286 throw std::logic_error(
"Could not determine a face direction");
290 auto halfDiff = [](
const DimVector& faceAreaNormal,
295 return computeHalfDiffusivity_(faceAreaNormal,
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);
313 const auto& nnc_input = eclState_.getInputNNC().input();
316#pragma omp parallel for
318 for (
const auto& chunk : ElementChunks(gridView_, Dune::Partitions::all, num_threads)) {
319 for (
const auto& elem : chunk) {
324 inside.
elemIdx = elemMapper.index(elem);
328 template getFieldPropCartesianIdx<Grid>(inside.
elemIdx);
330 auto computeHalf = [
this, &faceAreaNormal, &inside, &outside]
331 (
const auto& halfComputer,
332 const auto& prop1,
const auto& prop2) -> std::array<Scalar,2>
335 halfComputer(faceAreaNormal,
339 halfComputer(faceAreaNormal,
346 auto computeHalfMean = [&inside, &outside, &computeHalf, &ntg, &harmonicMean]
347 (
const auto& halfComputer,
const auto& prop)
349 auto onesided = computeHalf(halfComputer, prop[inside.
elemIdx], prop[outside.
elemIdx]);
350 applyNtg_(onesided[0], inside, ntg);
351 applyNtg_(onesided[1], outside, ntg);
354 return harmonicMean(onesided[0], onesided[1]);
357 unsigned boundaryIsIdx = 0;
358 for (
const auto& intersection : intersections(gridView_, elem)) {
360 if (intersection.boundary()) {
362 const auto& geometry = intersection.geometry();
365 faceAreaNormal = intersection.centerUnitOuterNormal();
366 faceAreaNormal *= geometry.volume();
368 Scalar transBoundaryIs =
369 computeHalfTrans_(faceAreaNormal,
370 intersection.indexInInside(),
372 permeability_[inside.
elemIdx]);
377 applyMultipliers_(transBoundaryIs, intersection.indexInInside(), inside.
cartElemIdx, transMult);
378 transBoundary.insert_or_assign(std::make_pair(inside.
elemIdx, boundaryIsIdx), transBoundaryIs);
382 if (enableEnergy_ && !onlyTrans) {
383 Scalar transBoundaryEnergyIs =
384 computeHalfDiffusivity_(faceAreaNormal,
387 thermalHalfTransBoundary.insert_or_assign(std::make_pair(inside.
elemIdx, boundaryIsIdx),
388 transBoundaryEnergyIs);
395 if (!intersection.neighbor()) {
402 const auto& outsideElem = intersection.outside();
403 outside.
elemIdx = elemMapper.index(outsideElem);
408 template getFieldPropCartesianIdx<Grid>(outside.
elemIdx);
422 inside.
faceIdx = intersection.indexInInside();
423 outside.
faceIdx = intersection.indexInOutside();
430 if (enableEnergy_ && !onlyTrans) {
435 if (updateDiffusivity && !onlyTrans) {
438 if (updateDispersivity && !onlyTrans) {
444 typename std::is_same<Grid, Dune::CpGrid>::type isCpGrid;
445 computeFaceProperties(intersection,
451 Scalar trans = computeHalfMean(computeHalfTrans_, permeability_);
457 auto find_layer = [&cartDims](std::size_t cell) {
459 auto k = cell / cartDims[1];
467 if (std::abs(kup -kdown) > 1) {
473 if (useSmallestMultiplier) {
477 applyAllZMultipliers_(trans, inside, outside, transMult, cartDims);
485 bool foundInputNNC =
false;
486 if (! nnc_input.empty()) {
488 auto it = std::lower_bound(nnc_input.begin(), nnc_input.end(),
489 NNCdata { inside.cartElemIdx, outside.cartElemIdx, 0.0 });
490 foundInputNNC = it != nnc_input.end() && it->cell1 == inside.
cartElemIdx && it->cell2 == outside.
cartElemIdx;
492 if (! foundInputNNC) {
494 trans *= transMult.getRegionMultiplier(inside.
cartElemIdx,
502 if (enableEnergy_ && !onlyTrans) {
503 const auto half = computeHalf(halfDiff, 1.0, 1.0);
512 if (updateDiffusivity && !onlyTrans) {
514 computeHalfMean(halfDiff, porosity_));
518 if (updateDispersivity && !onlyTrans) {
520 computeHalfMean(halfDiff, dispersion_));
525 centroids_cache_.clear();
528#pragma omp parallel sections
538 transBoundary.finalize();
542 thermalHalfTransBoundary.finalize();
546 thermalHalfTrans.finalize();
550 diffusivity.finalize();
554 dispersivity.finalize();
558 this->updateFromEclState_(global);
561 std::unordered_map<std::size_t,int> globalToLocal;
564 for (
const auto& elem : elements(grid_.leafGridView())) {
565 int elemIdx = elemMapper.index(elem);
566 int cartElemIdx = cartMapper_.cartesianIndex(elemIdx);
567 globalToLocal[cartElemIdx] = elemIdx;
576 this->applyPinchNncToGridTrans_(globalToLocal, applyNncMultregT);
577 this->applyNncToGridTrans_(globalToLocal);
578 this->applyEditNncToGridTrans_(globalToLocal);
579 this->applyEditNncrToGridTrans_(globalToLocal);
580 if (applyNncMultregT) {
581 this->applyNncMultreg_(globalToLocal);
583 warnEditNNC_ =
false;
588 this->removeNonCartesianTransmissibilities_(disableNNC);
591template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
595 unsigned numElem = gridView_.size(0);
596 permeability_.resize(numElem);
602 const auto& fp = eclState_.fieldProps();
603 if (fp.has_double(
"PERMX")) {
604 const std::vector<double>& permxData = this-> lookUpData_.assignFieldPropsDoubleOnLeaf(fp,
"PERMX");
606 std::vector<double> permyData;
607 if (fp.has_double(
"PERMY"))
608 permyData = this-> lookUpData_.assignFieldPropsDoubleOnLeaf(fp,
"PERMY");
610 permyData = permxData;
612 std::vector<double> permzData;
613 if (fp.has_double(
"PERMZ"))
614 permzData = this-> lookUpData_.assignFieldPropsDoubleOnLeaf(fp,
"PERMZ");
616 permzData = permxData;
618 for (std::size_t dofIdx = 0; dofIdx < numElem; ++ dofIdx) {
619 permeability_[dofIdx] = 0.0;
620 permeability_[dofIdx][0][0] = permxData[dofIdx];
621 permeability_[dofIdx][1][1] = permyData[dofIdx];
622 permeability_[dofIdx][2][2] = permzData[dofIdx];
629 throw std::logic_error(
"Can't read the intrinsic permeability from the ecl state. "
630 "(The PERM{X,Y,Z} keywords are missing)");
633template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
637 unsigned numElem = gridView_.size(0);
638 permeability_.resize(numElem);
644 const auto& fp = eclState_.fieldProps();
645 if (fp.has_double(
"PERMX")) {
646 const std::vector<double>& permxData =
647 this->lookUpData_.assignFieldPropsDoubleOnLeaf(fp,
"PERMX");
649 std::vector<double> permyData;
650 if (fp.has_double(
"PERMY")){
651 permyData = this->lookUpData_.assignFieldPropsDoubleOnLeaf(fp,
"PERMY");
654 permyData = permxData;
657 std::vector<double> permzData;
658 if (fp.has_double(
"PERMZ")) {
659 permzData = this->lookUpData_.assignFieldPropsDoubleOnLeaf(fp,
"PERMZ");
662 permzData = permxData;
665 for (std::size_t dofIdx = 0; dofIdx < numElem; ++ dofIdx) {
666 permeability_[dofIdx] = 0.0;
667 std::size_t inputDofIdx = map(dofIdx);
668 permeability_[dofIdx][0][0] = permxData[inputDofIdx];
669 permeability_[dofIdx][1][1] = permyData[inputDofIdx];
670 permeability_[dofIdx][2][2] = permzData[inputDofIdx];
676 throw std::logic_error(
"Can't read the intrinsic permeability from the ecl state. "
677 "(The PERM{X,Y,Z} keywords are missing)");
681template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
689 const auto& fp = eclState_.fieldProps();
690 if (fp.has_double(
"PORO")) {
691 if constexpr (std::is_same_v<Scalar,double>) {
692 porosity_ = this->lookUpData_.assignFieldPropsDoubleOnLeaf(fp,
"PORO");
695 const auto por = this->lookUpData_.assignFieldPropsDoubleOnLeaf(fp,
"PORO");
696 porosity_.resize(por.size());
697 std::ranges::copy(por, porosity_.begin());
701 throw std::logic_error(
"Can't read the porosity from the ecl state. "
702 "(The PORO keywords are missing)");
706template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
710 if (!enableDispersivity_) {
711 throw std::runtime_error(
"Dispersion disabled at compile time, but the deck "
712 "contains the DISPERC keyword.");
714 const auto& fp = eclState_.fieldProps();
715 if constexpr (std::is_same_v<Scalar,double>) {
716 dispersion_ = this->lookUpData_.assignFieldPropsDoubleOnLeaf(fp,
"DISPERC");
719 const auto disp = this->lookUpData_.assignFieldPropsDoubleOnLeaf(fp,
"DISPERC");
720 dispersion_.resize(disp.size());
721 std::ranges::copy(disp, dispersion_.begin());
725template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
729 const auto& cartDims = cartMapper_.cartesianDimensions();
730 for (
auto&& trans: trans_) {
732 if (removeAll || trans.second < transmissibilityThreshold_) {
733 const auto&
id = trans.first;
735 int gc1 = std::min(cartMapper_.cartesianIndex(elements.first), cartMapper_.cartesianIndex(elements.second));
736 int gc2 = std::max(cartMapper_.cartesianIndex(elements.first), cartMapper_.cartesianIndex(elements.second));
741 if (gc2 - gc1 == 1 || gc2 - gc1 == cartDims[0] || gc2 - gc1 == cartDims[0]*cartDims[1] || gc2 - gc1 == 0) {
750template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
755 const TransMult& transMult,
756 const std::array<int, dimWorld>& cartDims)
758 if (grid_.maxLevel() > 0) {
759 OPM_THROW(std::invalid_argument,
"MULTZ not support with LGRS, yet.");
765 unsigned lastCartElemIdx;
770 lastCartElemIdx = outside.
cartElemIdx - cartDims[0]*cartDims[1];
773 Scalar mult = transMult.getMultiplier(lastCartElemIdx , FaceDir::ZPlus) *
774 transMult.getMultiplier(outside.
cartElemIdx , FaceDir::ZMinus);
778 for (
auto cartElemIdx = inside.
cartElemIdx; cartElemIdx < lastCartElemIdx;) {
779 auto multiplier = transMult.getMultiplier(cartElemIdx, FaceDir::ZPlus);
780 cartElemIdx += cartDims[0]*cartDims[1];
781 multiplier *= transMult.getMultiplier(cartElemIdx, FaceDir::ZMinus);
782 mult = std::min(mult,
static_cast<Scalar
>(multiplier));
793template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
797 const FieldPropsManager* fp =
798 (global) ? &(eclState_.fieldProps()) :
799 &(eclState_.globalFieldProps());
801 std::array<bool,3> is_tran {fp->tran_active(
"TRANX"),
802 fp->tran_active(
"TRANY"),
803 fp->tran_active(
"TRANZ")};
805 if (!(is_tran[0] || is_tran[1] || is_tran[2])) {
810 std::array<std::string, 3> keywords {
"TRANX",
"TRANY",
"TRANZ"};
811 std::array<std::vector<double>,3> trans = createTransmissibilityArrays_(is_tran);
812 auto key = keywords.begin();
813 auto perform = is_tran.begin();
815 for (
auto it = trans.begin(); it != trans.end(); ++it, ++key, ++perform) {
817 if (grid_.maxLevel() > 0) {
818 OPM_THROW(std::invalid_argument,
"Calculations on TRANX/TRANY/TRANZ arrays are not support with LGRS, yet.");
820 fp->apply_tran(*key, *it);
824 resetTransmissibilityFromArrays_(is_tran, trans);
827template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
828std::array<std::vector<double>,3>
832 const auto& cartDims = cartMapper_.cartesianDimensions();
833 ElementMapper elemMapper(gridView_, Dune::mcmgElementLayout());
835 auto numElem = gridView_.size(0);
836 std::array<std::vector<double>,3> trans = {
837 std::vector<double>(is_tran[0] ? numElem : 0, 0),
838 std::vector<double>(is_tran[1] ? numElem : 0, 0),
839 std::vector<double>(is_tran[2] ? numElem : 0, 0)
843 for (
const auto& elem : elements(gridView_)) {
844 for (
const auto& intersection : intersections(gridView_, elem)) {
846 if (!intersection.neighbor()) {
861 unsigned c1 = elemMapper.index(intersection.inside());
862 unsigned c2 = elemMapper.index(intersection.outside());
863 int gc1 = cartMapper_.cartesianIndex(c1);
864 int gc2 = cartMapper_.cartesianIndex(c2);
865 if (std::tie(gc1, c1) > std::tie(gc2, c2)) {
879 if ((gc2 - gc1 == 1 || (gc2 == gc1 && (intersection.indexInInside() == 0 || intersection.indexInInside() == 1)))
884 trans[0][c1] = trans_[isID];
887 else if ((gc2 - gc1 == cartDims[0] || (gc2 == gc1 && (intersection.indexInInside() == 2 || intersection.indexInInside() == 3)))
892 trans[1][c1] = trans_[isID];
895 else if (gc2 - gc1 == cartDims[0]*cartDims[1] ||
896 (gc2 == gc1 && (intersection.indexInInside() == 4 || intersection.indexInInside() == 5)))
900 trans[2][c1] = trans_[isID];
910template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
913 const std::array<std::vector<double>,3>& trans)
915 const auto& cartDims = cartMapper_.cartesianDimensions();
916 ElementMapper elemMapper(gridView_, Dune::mcmgElementLayout());
919 for (
const auto& elem : elements(gridView_)) {
920 for (
const auto& intersection : intersections(gridView_, elem)) {
921 if (!intersection.neighbor()) {
936 unsigned c1 = elemMapper.index(intersection.inside());
937 unsigned c2 = elemMapper.index(intersection.outside());
938 int gc1 = cartMapper_.cartesianIndex(c1);
939 int gc2 = cartMapper_.cartesianIndex(c2);
940 if (std::tie(gc1, c1) > std::tie(gc2, c2)) {
954 if ((gc2 - gc1 == 1 || (gc2 == gc1 && (intersection.indexInInside() == 0 || intersection.indexInInside() == 1)))
959 trans_[isID] = trans[0][c1];
962 else if ((gc2 - gc1 == cartDims[0] || (gc2 == gc1 && (intersection.indexInInside() == 2|| intersection.indexInInside() == 3)))
967 trans_[isID] = trans[1][c1];
970 else if (gc2 - gc1 == cartDims[0]*cartDims[1] ||
971 (gc2 == gc1 && (intersection.indexInInside() == 4 || intersection.indexInInside() == 5)))
975 trans_[isID] = trans[2][c1];
984template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
985template<
class Intersection>
991 std::false_type)
const
994 const auto& geometry = intersection.geometry();
997 faceAreaNormal = intersection.centerUnitOuterNormal();
998 faceAreaNormal *= geometry.volume();
1001template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
1002template<
class Intersection>
1008 std::true_type)
const
1010 int faceIdx = intersection.id();
1012 if (grid_.maxLevel() == 0) {
1015 faceAreaNormal = grid_.faceAreaNormalEcl(faceIdx);
1018 if ((intersection.inside().level() != intersection.outside().level())) {
1025 const auto& parentIntersection =
1026 grid_.getParentIntersectionFromLgrBoundaryFace(intersection);
1027 const auto& parentIntersectionGeometry = parentIntersection.geometry();
1031 inside.
faceCenter = (intersection.inside().level() == 0)
1032 ? parentIntersectionGeometry.center()
1033 : grid_.faceCenterEcl(inside.
elemIdx, inside.
faceIdx, intersection);
1034 outside.
faceCenter = (intersection.outside().level() == 0)
1035 ? parentIntersectionGeometry.center()
1036 : grid_.faceCenterEcl(outside.
elemIdx, outside.
faceIdx, intersection);
1044 faceAreaNormal = intersection.centerUnitOuterNormal();
1045 faceAreaNormal *= intersection.geometry().volume();
1048 assert(intersection.inside().level() == intersection.outside().level());
1054 if (intersection.inside().level() > 0) {
1055 faceAreaNormal = intersection.centerUnitOuterNormal();
1056 faceAreaNormal *= intersection.geometry().volume();
1059 faceAreaNormal = grid_.faceAreaNormalEcl(faceIdx);
1065template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
1069 const bool applyNncMultregT)
1071 const auto& pinchNnc = eclState_.getPinchNNC();
1072 const auto& transMult = this->eclState_.getTransMult();
1074 for (
const auto& nncEntry : pinchNnc) {
1075 auto c1 = nncEntry.cell1;
1076 auto c2 = nncEntry.cell2;
1077 auto lowIt = cartesianToCompressed.find(c1);
1078 auto highIt = cartesianToCompressed.find(c2);
1079 int low = (lowIt == cartesianToCompressed.end())? -1 : lowIt->second;
1080 int high = (highIt == cartesianToCompressed.end())? -1 : highIt->second;
1083 std::swap(low, high);
1086 if (low == -1 && high == -1) {
1091 if (low == -1 || high == -1) {
1100 if (candidate != trans_.end()) {
1104 candidate->second = nncEntry.trans;
1105 if (applyNncMultregT) {
1106 const auto mult = transMult.getRegionMultiplierNNC(c1, c2);
1107 candidate->second *= mult;
1114template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
1120 const auto& nnc_input = eclState_.getInputNNC().input();
1122 for (
const auto& nncEntry : nnc_input) {
1123 auto c1 = nncEntry.cell1;
1124 auto c2 = nncEntry.cell2;
1125 auto lowIt = cartesianToCompressed.find(c1);
1126 auto highIt = cartesianToCompressed.find(c2);
1127 int low = (lowIt == cartesianToCompressed.end())? -1 : lowIt->second;
1128 int high = (highIt == cartesianToCompressed.end())? -1 : highIt->second;
1131 std::swap(low, high);
1134 if (low == -1 && high == -1) {
1139 if (low == -1 || high == -1) {
1141 std::ostringstream sstr;
1142 sstr <<
"NNC between active and inactive cells ("
1143 << low <<
" -> " << high <<
") with globalcell is (" << c1 <<
"->" << c2 <<
")";
1144 OpmLog::warning(sstr.str());
1148 if (
auto candidate = trans_.find(
details::isId(low, high)); candidate != trans_.end()) {
1152 candidate->second += nncEntry.trans;
1182template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
1186 const auto& input = eclState_.getInputNNC();
1187 applyEditNncToGridTransHelper_(globalToLocal,
"EDITNNC",
1189 [&input](
const NNCdata& nnc){
1190 return input.edit_location(nnc);},
1192 [](Scalar& trans,
const Scalar& rhs){ trans *= rhs;});
1195template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
1199 const auto& input = eclState_.getInputNNC();
1200 applyEditNncToGridTransHelper_(globalToLocal,
"EDITNNCR",
1202 [&input](
const NNCdata& nnc){
1203 return input.editr_location(nnc);},
1205 [](Scalar& trans,
const Scalar& rhs){ trans = rhs;});
1208template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
1211 const std::string& keyword,
1212 const std::vector<NNCdata>& nncs,
1213 const std::function<KeywordLocation(
const NNCdata&)>& getLocation,
1214 const std::function<
void(Scalar&,
const Scalar&)>& apply)
1219 const auto& cartDims = cartMapper_.cartesianDimensions();
1221 auto format_ijk = [&cartDims](std::size_t cell) -> std::string
1223 auto i = cell % cartDims[0]; cell /= cartDims[0];
1224 auto j = cell % cartDims[1];
1225 auto k = cell / cartDims[1];
1227 return fmt::format(
"({},{},{})", i + 1,j + 1,k + 1);
1230 auto print_warning = [&format_ijk, &getLocation, &keyword] (
const NNCdata& nnc)
1232 const auto& location = getLocation( nnc );
1233 auto warning = fmt::format(
"Problem with {} keyword\n"
1235 "No NNC defined for connection {} -> {}", keyword, location.filename,
1236 location.lineno, format_ijk(nnc.cell1), format_ijk(nnc.cell2));
1237 OpmLog::warning(keyword, warning);
1243 auto nnc = nncs.begin();
1244 auto end = nncs.end();
1245 std::size_t warning_count = 0;
1246 while (nnc != end) {
1247 auto c1 = nnc->cell1;
1248 auto c2 = nnc->cell2;
1249 auto lowIt = globalToLocal.find(c1);
1250 auto highIt = globalToLocal.find(c2);
1252 if (lowIt == globalToLocal.end() || highIt == globalToLocal.end()) {
1254 if (lowIt != highIt && warnEditNNC_) {
1255 print_warning(*nnc);
1262 auto low = lowIt->second, high = highIt->second;
1265 std::swap(low, high);
1269 if (candidate == trans_.end() && warnEditNNC_) {
1270 print_warning(*nnc);
1276 while (nnc != end && c1 == nnc->cell1 && c2 == nnc->cell2) {
1277 apply(candidate->second, nnc->trans);
1283 if (warning_count > 0) {
1284 auto warning = fmt::format(
"Problems with {} keyword\n"
1285 "A total of {} connections not defined in grid", keyword, warning_count);
1286 OpmLog::warning(warning);
1290template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
1293applyNncMultreg_(
const std::unordered_map<std::size_t,int>& cartesianToCompressed)
1295 const auto& inputNNC = this->eclState_.getInputNNC();
1296 const auto& transMult = this->eclState_.getTransMult();
1298 auto compressedIdx = [&cartesianToCompressed](
const std::size_t globIdx)
1300 auto ixPos = cartesianToCompressed.find(globIdx);
1301 return (ixPos == cartesianToCompressed.end()) ? -1 : ixPos->second;
1315 for (
const auto& nncList : {&NNC::input, &NNC::editr}) {
1316 for (
const auto& nncEntry : (inputNNC.*nncList)()) {
1317 const auto c1 = nncEntry.cell1;
1318 const auto c2 = nncEntry.cell2;
1320 auto low = compressedIdx(c1);
1321 auto high = compressedIdx(c2);
1323 if ((low == -1) || (high == -1)) {
1328 std::swap(low, high);
1331 auto candidate = this->trans_.find(
details::isId(low, high));
1332 if (candidate != this->trans_.end()) {
1333 candidate->second *= transMult.getRegionMultiplierNNC(c1, c2);
1339template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
1347 assert(faceIdx >= 0);
1348 unsigned dimIdx = faceIdx / 2;
1349 assert(dimIdx < dimWorld);
1350 Scalar halfTrans = perm[dimIdx][dimIdx];
1351 halfTrans *= std::abs(Dune::dot(areaNormal, distance));
1352 halfTrans /= distance.two_norm2();
1357template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
1364 Scalar halfDiff = poro;
1365 halfDiff *= std::abs(Dune::dot(areaNormal, distance));
1366 halfDiff /= distance.two_norm2();
1371template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
1375 const unsigned& cellIdx)
const
1377 const auto& cellCenter = centroids_cache_.empty() ? centroids_(cellIdx)
1378 : centroids_cache_[cellIdx];
1380 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
1381 x[dimIdx] -= cellCenter[dimIdx];
1387template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
1391 unsigned cartElemIdx,
1392 const TransMult& transMult)
const
1397 trans *= transMult.getMultiplier(cartElemIdx,
1398 FaceDir::FromIntersectionIndex(faceIdx));
1401template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
1405 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:1389
Scalar diffusivity(unsigned elemIdx1, unsigned elemIdx2) const
Return the diffusivity for the intersection between two elements.
Definition: Transmissibility_impl.hpp:150
void applyEditNncToGridTrans_(const std::unordered_map< std::size_t, int > &globalToLocal)
Multiplies the grid transmissibilities according to EDITNNC.
Definition: Transmissibility_impl.hpp:1184
void computeFaceProperties(const Intersection &intersection, FaceInfo &inside, FaceInfo &outside, DimVector &faceAreaNormal, std::false_type) const
Definition: Transmissibility_impl.hpp:987
Dune::FieldMatrix< Scalar, dimWorld, dimWorld > DimMatrix
Definition: Transmissibility.hpp:59
void update(bool global, TransUpdateQuantities update_quantities=TransUpdateQuantities::All, const std::function< unsigned int(unsigned int)> &map={}, bool applyNncMultRegT=false)
Definition: Transmissibility_impl.hpp:170
DimVector distanceVector_(const DimVector &faceCenter, const unsigned &cellIdx) const
Definition: Transmissibility_impl.hpp:1374
void applyNncMultreg_(const std::unordered_map< std::size_t, int > &globalToLocal)
Definition: Transmissibility_impl.hpp:1293
void applyEditNncrToGridTrans_(const std::unordered_map< std::size_t, int > &globalToLocal)
Resets the grid transmissibilities according to EDITNNCR.
Definition: Transmissibility_impl.hpp:1197
void extractDispersion_()
Definition: Transmissibility_impl.hpp:708
Scalar thermalHalfTrans(unsigned insideElemIdx, unsigned outsideElemIdx) const
Return the thermal "half transmissibility" for the intersection between two elements.
Definition: Transmissibility_impl.hpp:136
Scalar transmissibilityBoundary(unsigned elemIdx, unsigned boundaryFaceIdx) const
Return the transmissibility for a given boundary segment.
Definition: Transmissibility_impl.hpp:129
std::array< std::vector< double >, 3 > createTransmissibilityArrays_(const std::array< bool, 3 > &is_tran)
Creates TRANS{XYZ} arrays for modification by FieldProps data.
Definition: Transmissibility_impl.hpp:830
void applyEditNncToGridTransHelper_(const std::unordered_map< std::size_t, int > &globalToLocal, const std::string &keyword, const std::vector< NNCdata > &nncs, const std::function< KeywordLocation(const NNCdata &)> &getLocation, const std::function< void(Scalar &, const Scalar &)> &apply)
Definition: Transmissibility_impl.hpp:1210
void updateFromEclState_(bool global)
Definition: Transmissibility_impl.hpp:795
void extractPermeability_()
Definition: Transmissibility_impl.hpp:593
void applyAllZMultipliers_(Scalar &trans, const FaceInfo &inside, const FaceInfo &outside, const TransMult &transMult, const std::array< int, dimWorld > &cartDims)
Apply the Multipliers for the case PINCH(4)==TOPBOT.
Definition: Transmissibility_impl.hpp:752
TransUpdateQuantities
Compute all transmissibilities.
Definition: Transmissibility.hpp:157
void resetTransmissibilityFromArrays_(const std::array< bool, 3 > &is_tran, const std::array< std::vector< double >, 3 > &trans)
overwrites calculated transmissibilities
Definition: Transmissibility_impl.hpp:912
Dune::FieldVector< Scalar, dimWorld > DimVector
Definition: Transmissibility.hpp:60
static Scalar computeHalfDiffusivity_(const DimVector &areaNormal, const DimVector &distance, const Scalar poro)
Definition: Transmissibility_impl.hpp:1360
Scalar transmissibilityThreshold_
Definition: Transmissibility.hpp:289
Scalar transmissibility(unsigned elemIdx1, unsigned elemIdx2) const
Return the transmissibility for the intersection between two elements.
Definition: Transmissibility_impl.hpp:122
static void applyNtg_(Scalar &trans, const FaceInfo &face, const std::vector< double > &ntg)
Definition: Transmissibility_impl.hpp:1403
const EclipseState & eclState_
Definition: Transmissibility.hpp:283
Scalar dispersivity(unsigned elemIdx1, unsigned elemIdx2) const
Return the dispersivity for the intersection between two elements.
Definition: Transmissibility_impl.hpp:160
Transmissibility(const EclipseState &eclState, const GridView &gridView, const CartesianIndexMapper &cartMapper, const Grid &grid, std::function< std::array< double, dimWorld >(int)> centroids, bool enableEnergy, bool enableDiffusivity, bool enableDispersivity)
Definition: Transmissibility_impl.hpp:97
void applyNncToGridTrans_(const std::unordered_map< std::size_t, int > &cartesianToCompressed)
Definition: Transmissibility_impl.hpp:1117
void applyPinchNncToGridTrans_(const std::unordered_map< std::size_t, int > &cartesianToCompressed, bool applyNncMultregT)
Applies the previous calculate transmissibilities to the NNCs created via PINCH.
Definition: Transmissibility_impl.hpp:1068
void extractPorosity_()
Definition: Transmissibility_impl.hpp:683
void removeNonCartesianTransmissibilities_(bool removeAll)
Definition: Transmissibility_impl.hpp:727
static Scalar computeHalfTrans_(const DimVector &areaNormal, int faceIdx, const DimVector &distance, const DimMatrix &perm)
Definition: Transmissibility_impl.hpp:1342
Scalar thermalHalfTransBoundary(unsigned insideElemIdx, unsigned boundaryFaceIdx) const
Definition: Transmissibility_impl.hpp:143
constexpr unsigned elemIdxShift
Definition: FacePropertiesTPSA_impl.hpp:52
std::pair< std::uint32_t, std::uint32_t > isIdReverse(const std::uint64_t &id)
Definition: Transmissibility_impl.hpp:78
std::uint64_t directionalIsId(std::uint32_t elemIdx1, std::uint32_t elemIdx2)
Definition: Transmissibility_impl.hpp:89
std::uint64_t isId(std::uint32_t elemIdx1, std::uint32_t elemIdx2)
Definition: Transmissibility_impl.hpp:70
Definition: blackoilbioeffectsmodules.hh:45
Definition: Transmissibility.hpp:167
unsigned cartElemIdx
Definition: Transmissibility.hpp:171
DimVector faceCenter
Definition: Transmissibility.hpp:168
unsigned elemIdx
Definition: Transmissibility.hpp:170
int faceIdx
Definition: Transmissibility.hpp:169