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);
314#pragma omp parallel for
316 for (
const auto& chunk : ElementChunks(gridView_, Dune::Partitions::all, num_threads)) {
317 for (
const auto& elem : chunk) {
322 inside.
elemIdx = elemMapper.index(elem);
326 template getFieldPropCartesianIdx<Grid>(inside.
elemIdx);
328 auto computeHalf = [
this, &faceAreaNormal, &inside, &outside]
329 (
const auto& halfComputer,
330 const auto& prop1,
const auto& prop2) -> std::array<Scalar,2>
333 halfComputer(faceAreaNormal,
337 halfComputer(faceAreaNormal,
344 auto computeHalfMean = [&inside, &outside, &computeHalf, &ntg, &harmonicMean]
345 (
const auto& halfComputer,
const auto& prop)
347 auto onesided = computeHalf(halfComputer, prop[inside.
elemIdx], prop[outside.
elemIdx]);
348 applyNtg_(onesided[0], inside, ntg);
349 applyNtg_(onesided[1], outside, ntg);
352 return harmonicMean(onesided[0], onesided[1]);
355 unsigned boundaryIsIdx = 0;
356 for (
const auto& intersection : intersections(gridView_, elem)) {
358 if (intersection.boundary()) {
360 const auto& geometry = intersection.geometry();
363 faceAreaNormal = intersection.centerUnitOuterNormal();
364 faceAreaNormal *= geometry.volume();
366 Scalar transBoundaryIs =
367 computeHalfTrans_(faceAreaNormal,
368 intersection.indexInInside(),
370 permeability_[inside.
elemIdx]);
375 applyMultipliers_(transBoundaryIs, intersection.indexInInside(), inside.
cartElemIdx, transMult);
376 transBoundary.insert_or_assign(std::make_pair(inside.
elemIdx, boundaryIsIdx), transBoundaryIs);
380 if (enableEnergy_ && !onlyTrans) {
381 Scalar transBoundaryEnergyIs =
382 computeHalfDiffusivity_(faceAreaNormal,
385 thermalHalfTransBoundary.insert_or_assign(std::make_pair(inside.
elemIdx, boundaryIsIdx),
386 transBoundaryEnergyIs);
393 if (!intersection.neighbor()) {
400 const auto& outsideElem = intersection.outside();
401 outside.
elemIdx = elemMapper.index(outsideElem);
406 template getFieldPropCartesianIdx<Grid>(outside.
elemIdx);
420 inside.
faceIdx = intersection.indexInInside();
421 outside.
faceIdx = intersection.indexInOutside();
428 if (enableEnergy_ && !onlyTrans) {
433 if (updateDiffusivity && !onlyTrans) {
436 if (updateDispersivity && !onlyTrans) {
442 typename std::is_same<Grid, Dune::CpGrid>::type isCpGrid;
443 computeFaceProperties(intersection,
449 Scalar trans = computeHalfMean(computeHalfTrans_, permeability_);
455 auto find_layer = [&cartDims](std::size_t cell) {
457 auto k = cell / cartDims[1];
465 if (std::abs(kup -kdown) > 1) {
471 if (useSmallestMultiplier) {
475 applyAllZMultipliers_(trans, inside, outside, transMult, cartDims);
484 trans *= transMult.getRegionMultiplier(inside.
cartElemIdx,
491 if (enableEnergy_ && !onlyTrans) {
492 const auto half = computeHalf(halfDiff, 1.0, 1.0);
501 if (updateDiffusivity && !onlyTrans) {
503 computeHalfMean(halfDiff, porosity_));
507 if (updateDispersivity && !onlyTrans) {
509 computeHalfMean(halfDiff, dispersion_));
514 centroids_cache_.clear();
517#pragma omp parallel sections
527 transBoundary.finalize();
531 thermalHalfTransBoundary.finalize();
535 thermalHalfTrans.finalize();
539 diffusivity.finalize();
543 dispersivity.finalize();
547 this->updateFromEclState_(global);
550 std::unordered_map<std::size_t,int> globalToLocal;
553 for (
const auto& elem : elements(grid_.leafGridView())) {
554 int elemIdx = elemMapper.index(elem);
555 int cartElemIdx = cartMapper_.cartesianIndex(elemIdx);
556 globalToLocal[cartElemIdx] = elemIdx;
565 this->applyEditNncToGridTrans_(globalToLocal);
566 this->applyPinchNncToGridTrans_(globalToLocal);
567 this->applyNncToGridTrans_(globalToLocal);
568 this->applyEditNncrToGridTrans_(globalToLocal);
569 if (applyNncMultregT) {
570 this->applyNncMultreg_(globalToLocal);
572 warnEditNNC_ =
false;
577 this->removeNonCartesianTransmissibilities_(disableNNC);
580template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
584 unsigned numElem = gridView_.size(0);
585 permeability_.resize(numElem);
591 const auto& fp = eclState_.fieldProps();
592 if (fp.has_double(
"PERMX")) {
593 const std::vector<double>& permxData = this-> lookUpData_.assignFieldPropsDoubleOnLeaf(fp,
"PERMX");
595 std::vector<double> permyData;
596 if (fp.has_double(
"PERMY"))
597 permyData = this-> lookUpData_.assignFieldPropsDoubleOnLeaf(fp,
"PERMY");
599 permyData = permxData;
601 std::vector<double> permzData;
602 if (fp.has_double(
"PERMZ"))
603 permzData = this-> lookUpData_.assignFieldPropsDoubleOnLeaf(fp,
"PERMZ");
605 permzData = permxData;
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];
618 throw std::logic_error(
"Can't read the intrinsic permeability from the ecl state. "
619 "(The PERM{X,Y,Z} keywords are missing)");
622template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
626 unsigned numElem = gridView_.size(0);
627 permeability_.resize(numElem);
633 const auto& fp = eclState_.fieldProps();
634 if (fp.has_double(
"PERMX")) {
635 const std::vector<double>& permxData =
636 this->lookUpData_.assignFieldPropsDoubleOnLeaf(fp,
"PERMX");
638 std::vector<double> permyData;
639 if (fp.has_double(
"PERMY")){
640 permyData = this->lookUpData_.assignFieldPropsDoubleOnLeaf(fp,
"PERMY");
643 permyData = permxData;
646 std::vector<double> permzData;
647 if (fp.has_double(
"PERMZ")) {
648 permzData = this->lookUpData_.assignFieldPropsDoubleOnLeaf(fp,
"PERMZ");
651 permzData = permxData;
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];
665 throw std::logic_error(
"Can't read the intrinsic permeability from the ecl state. "
666 "(The PERM{X,Y,Z} keywords are missing)");
670template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
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");
684 const auto por = this->lookUpData_.assignFieldPropsDoubleOnLeaf(fp,
"PORO");
685 porosity_.resize(por.size());
686 std::copy(por.begin(), por.end(), porosity_.begin());
690 throw std::logic_error(
"Can't read the porosity from the ecl state. "
691 "(The PORO keywords are missing)");
695template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
699 if (!enableDispersivity_) {
700 throw std::runtime_error(
"Dispersion disabled at compile time, but the deck "
701 "contains the DISPERC keyword.");
703 const auto& fp = eclState_.fieldProps();
704 if constexpr (std::is_same_v<Scalar,double>) {
705 dispersion_ = this->lookUpData_.assignFieldPropsDoubleOnLeaf(fp,
"DISPERC");
708 const auto disp = this->lookUpData_.assignFieldPropsDoubleOnLeaf(fp,
"DISPERC");
709 dispersion_.resize(disp.size());
710 std::copy(disp.begin(), disp.end(), dispersion_.begin());
714template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
718 const auto& cartDims = cartMapper_.cartesianDimensions();
719 for (
auto&& trans: trans_) {
721 if (removeAll || trans.second < transmissibilityThreshold_) {
722 const auto&
id = trans.first;
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));
730 if (gc2 - gc1 == 1 || gc2 - gc1 == cartDims[0] || gc2 - gc1 == cartDims[0]*cartDims[1] || gc2 - gc1 == 0) {
739template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
744 const TransMult& transMult,
745 const std::array<int, dimWorld>& cartDims)
747 if (grid_.maxLevel() > 0) {
748 OPM_THROW(std::invalid_argument,
"MULTZ not support with LGRS, yet.");
754 unsigned lastCartElemIdx;
759 lastCartElemIdx = outside.
cartElemIdx - cartDims[0]*cartDims[1];
762 Scalar mult = transMult.getMultiplier(lastCartElemIdx , FaceDir::ZPlus) *
763 transMult.getMultiplier(outside.
cartElemIdx , FaceDir::ZMinus);
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));
782template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
786 const FieldPropsManager* fp =
787 (global) ? &(eclState_.fieldProps()) :
788 &(eclState_.globalFieldProps());
790 std::array<bool,3> is_tran {fp->tran_active(
"TRANX"),
791 fp->tran_active(
"TRANY"),
792 fp->tran_active(
"TRANZ")};
794 if (!(is_tran[0] || is_tran[1] || is_tran[2])) {
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();
804 for (
auto it = trans.begin(); it != trans.end(); ++it, ++key, ++perform) {
806 if (grid_.maxLevel() > 0) {
807 OPM_THROW(std::invalid_argument,
"Calculations on TRANX/TRANY/TRANZ arrays are not support with LGRS, yet.");
809 fp->apply_tran(*key, *it);
813 resetTransmissibilityFromArrays_(is_tran, trans);
816template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
817std::array<std::vector<double>,3>
821 const auto& cartDims = cartMapper_.cartesianDimensions();
822 ElementMapper elemMapper(gridView_, Dune::mcmgElementLayout());
824 auto numElem = gridView_.size(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)
832 for (
const auto& elem : elements(gridView_)) {
833 for (
const auto& intersection : intersections(gridView_, elem)) {
835 if (!intersection.neighbor()) {
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)) {
868 if ((gc2 - gc1 == 1 || (gc2 == gc1 && (intersection.indexInInside() == 0 || intersection.indexInInside() == 1)))
873 trans[0][c1] = trans_[isID];
876 else if ((gc2 - gc1 == cartDims[0] || (gc2 == gc1 && (intersection.indexInInside() == 2 || intersection.indexInInside() == 3)))
881 trans[1][c1] = trans_[isID];
884 else if (gc2 - gc1 == cartDims[0]*cartDims[1] ||
885 (gc2 == gc1 && (intersection.indexInInside() == 4 || intersection.indexInInside() == 5)))
889 trans[2][c1] = trans_[isID];
899template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
902 const std::array<std::vector<double>,3>& trans)
904 const auto& cartDims = cartMapper_.cartesianDimensions();
905 ElementMapper elemMapper(gridView_, Dune::mcmgElementLayout());
908 for (
const auto& elem : elements(gridView_)) {
909 for (
const auto& intersection : intersections(gridView_, elem)) {
910 if (!intersection.neighbor()) {
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)) {
943 if ((gc2 - gc1 == 1 || (gc2 == gc1 && (intersection.indexInInside() == 0 || intersection.indexInInside() == 1)))
948 trans_[isID] = trans[0][c1];
951 else if ((gc2 - gc1 == cartDims[0] || (gc2 == gc1 && (intersection.indexInInside() == 2|| intersection.indexInInside() == 3)))
956 trans_[isID] = trans[1][c1];
959 else if (gc2 - gc1 == cartDims[0]*cartDims[1] ||
960 (gc2 == gc1 && (intersection.indexInInside() == 4 || intersection.indexInInside() == 5)))
964 trans_[isID] = trans[2][c1];
973template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
974template<
class Intersection>
980 std::false_type)
const
983 const auto& geometry = intersection.geometry();
986 faceAreaNormal = intersection.centerUnitOuterNormal();
987 faceAreaNormal *= geometry.volume();
990template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
991template<
class Intersection>
997 std::true_type)
const
999 int faceIdx = intersection.id();
1001 if (grid_.maxLevel() == 0) {
1004 faceAreaNormal = grid_.faceAreaNormalEcl(faceIdx);
1007 if ((intersection.inside().level() != intersection.outside().level())) {
1014 const auto& parentIntersection =
1015 grid_.getParentIntersectionFromLgrBoundaryFace(intersection);
1016 const auto& parentIntersectionGeometry = parentIntersection.geometry();
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);
1033 faceAreaNormal = intersection.centerUnitOuterNormal();
1034 faceAreaNormal *= intersection.geometry().volume();
1037 assert(intersection.inside().level() == intersection.outside().level());
1043 if (intersection.inside().level() > 0) {
1044 faceAreaNormal = intersection.centerUnitOuterNormal();
1045 faceAreaNormal *= intersection.geometry().volume();
1048 faceAreaNormal = grid_.faceAreaNormalEcl(faceIdx);
1054template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
1060 const auto& nnc_input = eclState_.getPinchNNC();
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;
1071 std::swap(low, high);
1074 if (low == -1 && high == -1) {
1079 if (low == -1 || high == -1) {
1088 if (candidate != trans_.end()) {
1091 candidate->second = nncEntry.trans;
1097template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
1103 const auto& nnc_input = eclState_.getInputNNC().input();
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;
1114 std::swap(low, high);
1117 if (low == -1 && high == -1) {
1122 if (low == -1 || high == -1) {
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());
1131 if (
auto candidate = trans_.find(
details::isId(low, high)); candidate != trans_.end()) {
1135 candidate->second += nncEntry.trans;
1165template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
1169 const auto& input = eclState_.getInputNNC();
1170 applyEditNncToGridTransHelper_(globalToLocal,
"EDITNNC",
1172 [&input](
const NNCdata& nnc){
1173 return input.edit_location(nnc);},
1175 [](Scalar& trans,
const Scalar& rhs){ trans *= rhs;});
1178template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
1182 const auto& input = eclState_.getInputNNC();
1183 applyEditNncToGridTransHelper_(globalToLocal,
"EDITNNCR",
1185 [&input](
const NNCdata& nnc){
1186 return input.editr_location(nnc);},
1188 [](Scalar& trans,
const Scalar& rhs){ trans = rhs;});
1191template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
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)
1202 const auto& cartDims = cartMapper_.cartesianDimensions();
1204 auto format_ijk = [&cartDims](std::size_t cell) -> std::string
1206 auto i = cell % cartDims[0]; cell /= cartDims[0];
1207 auto j = cell % cartDims[1];
1208 auto k = cell / cartDims[1];
1210 return fmt::format(
"({},{},{})", i + 1,j + 1,k + 1);
1213 auto print_warning = [&format_ijk, &getLocation, &keyword] (
const NNCdata& nnc)
1215 const auto& location = getLocation( nnc );
1216 auto warning = fmt::format(
"Problem with {} keyword\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);
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);
1235 if (lowIt == globalToLocal.end() || highIt == globalToLocal.end()) {
1237 if (lowIt != highIt && warnEditNNC_) {
1238 print_warning(*nnc);
1245 auto low = lowIt->second, high = highIt->second;
1248 std::swap(low, high);
1252 if (candidate == trans_.end() && warnEditNNC_) {
1253 print_warning(*nnc);
1259 while (nnc != end && c1 == nnc->cell1 && c2 == nnc->cell2) {
1260 apply(candidate->second, nnc->trans);
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);
1273template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
1276applyNncMultreg_(
const std::unordered_map<std::size_t,int>& cartesianToCompressed)
1278 const auto& inputNNC = this->eclState_.getInputNNC();
1279 const auto& transMult = this->eclState_.getTransMult();
1281 auto compressedIdx = [&cartesianToCompressed](
const std::size_t globIdx)
1283 auto ixPos = cartesianToCompressed.find(globIdx);
1284 return (ixPos == cartesianToCompressed.end()) ? -1 : ixPos->second;
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;
1303 auto low = compressedIdx(c1);
1304 auto high = compressedIdx(c2);
1306 if ((low == -1) || (high == -1)) {
1311 std::swap(low, high);
1314 auto candidate = this->trans_.find(
details::isId(low, high));
1315 if (candidate != this->trans_.end()) {
1316 candidate->second *= transMult.getRegionMultiplierNNC(c1, c2);
1322template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
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();
1340template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
1347 Scalar halfDiff = poro;
1348 halfDiff *= std::abs(Dune::dot(areaNormal, distance));
1349 halfDiff /= distance.two_norm2();
1354template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
1358 const unsigned& cellIdx)
const
1360 const auto& cellCenter = centroids_cache_.empty() ? centroids_(cellIdx)
1361 : centroids_cache_[cellIdx];
1363 for (
unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx) {
1364 x[dimIdx] -= cellCenter[dimIdx];
1370template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
1374 unsigned cartElemIdx,
1375 const TransMult& transMult)
const
1380 trans *= transMult.getMultiplier(cartElemIdx,
1381 FaceDir::FromIntersectionIndex(faceIdx));
1384template<
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper,
class Scalar>
1388 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: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