Transmissibility_impl.hpp
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
23#ifndef OPM_TRANSMISSIBILITY_IMPL_HPP
24#define OPM_TRANSMISSIBILITY_IMPL_HPP
25
26#include <dune/common/version.hh>
27#include <dune/grid/common/mcmgmapper.hh>
28
29#include <opm/grid/CpGrid.hpp>
30
31#include <opm/common/OpmLog/KeywordLocation.hpp>
32#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
33#include <opm/input/eclipse/EclipseState/Grid/FaceDir.hpp>
34#include <opm/input/eclipse/EclipseState/Grid/FieldPropsManager.hpp>
35#include <opm/input/eclipse/EclipseState/Grid/TransMult.hpp>
36#include <opm/input/eclipse/Units/Units.hpp>
37
39
40#include <fmt/format.h>
41
42#include <algorithm>
43#include <array>
44#include <cassert>
45#include <cmath>
46#include <cstddef>
47#include <cstdint>
48#include <functional>
49#include <initializer_list>
50#include <sstream>
51#include <stdexcept>
52#include <type_traits>
53#include <utility>
54#include <vector>
55
56namespace Opm {
57
58namespace details {
59
60 constexpr unsigned elemIdxShift = 32; // bits
61
62 std::uint64_t isId(std::uint32_t elemIdx1, std::uint32_t elemIdx2)
63 {
64 std::uint32_t elemAIdx = std::min(elemIdx1, elemIdx2);
65 std::uint64_t elemBIdx = std::max(elemIdx1, elemIdx2);
66
67 return (elemBIdx<<elemIdxShift) + elemAIdx;
68 }
69
70 std::pair<std::uint32_t, std::uint32_t> isIdReverse(const std::uint64_t& id)
71 {
72 // Assigning an unsigned integer to a narrower type discards the most significant bits.
73 // See "The C programming language", section A.6.2.
74 // NOTE that the ordering of element A and B may have changed
75 std::uint32_t elemAIdx = id;
76 std::uint32_t elemBIdx = (id - elemAIdx) >> elemIdxShift;
77
78 return std::make_pair(elemAIdx, elemBIdx);
79 }
80
81 std::uint64_t directionalIsId(std::uint32_t elemIdx1, std::uint32_t elemIdx2)
82 {
83 return (std::uint64_t(elemIdx1)<<elemIdxShift) + elemIdx2;
84 }
85}
86
87template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
89Transmissibility(const EclipseState& eclState,
90 const GridView& gridView,
91 const CartesianIndexMapper& cartMapper,
92 const Grid& grid,
93 std::function<std::array<double,dimWorld>(int)> centroids,
94 bool enableEnergy,
95 bool enableDiffusivity,
96 bool enableDispersivity)
97 : eclState_(eclState)
98 , gridView_(gridView)
99 , cartMapper_(cartMapper)
100 , grid_(grid)
101 , centroids_(centroids)
102 , enableEnergy_(enableEnergy)
103 , enableDiffusivity_(enableDiffusivity)
104 , enableDispersivity_(enableDispersivity)
105 , lookUpData_(gridView)
106 , lookUpCartesianData_(gridView, cartMapper)
107{
108 const UnitSystem& unitSystem = eclState_.getDeckUnitSystem();
109 transmissibilityThreshold_ = unitSystem.parse("Transmissibility").getSIScaling() * 1e-6;
110}
111
112template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
114transmissibility(unsigned elemIdx1, unsigned elemIdx2) const
115{
116 return trans_.at(details::isId(elemIdx1, elemIdx2));
117}
118
119template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
121transmissibilityBoundary(unsigned elemIdx, unsigned boundaryFaceIdx) const
122{
123 return transBoundary_.at(std::make_pair(elemIdx, boundaryFaceIdx));
124}
125
126template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
128thermalHalfTrans(unsigned insideElemIdx, unsigned outsideElemIdx) const
129{
130 return thermalHalfTrans_.at(details::directionalIsId(insideElemIdx, outsideElemIdx));
131}
132
133template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
135thermalHalfTransBoundary(unsigned insideElemIdx, unsigned boundaryFaceIdx) const
136{
137 return thermalHalfTransBoundary_.at(std::make_pair(insideElemIdx, boundaryFaceIdx));
138}
139
140template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
142diffusivity(unsigned elemIdx1, unsigned elemIdx2) const
143{
144 if (diffusivity_.empty())
145 return 0.0;
146
147 return diffusivity_.at(details::isId(elemIdx1, elemIdx2));
148}
149
150template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
152dispersivity(unsigned elemIdx1, unsigned elemIdx2) const
153{
154 if (dispersivity_.empty())
155 return 0.0;
156
157 return dispersivity_.at(details::isId(elemIdx1, elemIdx2));
158}
159
160template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
162update(bool global, const std::function<unsigned int(unsigned int)>& map, const bool applyNncMultregT)
163{
164 const auto& cartDims = cartMapper_.cartesianDimensions();
165 const auto& transMult = eclState_.getTransMult();
166 const auto& comm = gridView_.comm();
167 ElementMapper elemMapper(gridView_, Dune::mcmgElementLayout());
168
169 unsigned numElements = elemMapper.size();
170 // get the ntg values, the ntg values are modified for the cells merged with minpv
171 const std::vector<double>& ntg = this->lookUpData_.assignFieldPropsDoubleOnLeaf(eclState_.fieldProps(), "NTG");
172 const bool updateDiffusivity = eclState_.getSimulationConfig().isDiffusive();
173 const bool updateDispersivity = eclState_.getSimulationConfig().rock_config().dispersion();
174
175 const bool disableNNC = eclState_.getSimulationConfig().useNONNC();
176
177 if (map)
178 extractPermeability_(map);
179 else
180 extractPermeability_();
181
182 // calculate the axis specific centroids of all elements
183 std::array<std::vector<DimVector>, dimWorld> axisCentroids;
184
185 for (unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx)
186 axisCentroids[dimIdx].resize(numElements);
187
188 for (const auto& elem : elements(gridView_)) {
189 unsigned elemIdx = elemMapper.index(elem);
190
191 // compute the axis specific "centroids" used for the transmissibilities. for
192 // consistency with the flow simulator, we use the element centers as
193 // computed by opm-parser's Opm::EclipseGrid class for all axes.
194 std::array<double, dimWorld> centroid = centroids_(elemIdx);
195
196 for (unsigned axisIdx = 0; axisIdx < dimWorld; ++axisIdx)
197 for (unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx)
198 axisCentroids[axisIdx][elemIdx][dimIdx] = centroid[dimIdx];
199 }
200
201 // reserving some space in the hashmap upfront saves quite a bit of time because
202 // resizes are costly for hashmaps and there would be quite a few of them if we
203 // would not have a rough idea of how large the final map will be (the rough idea
204 // is a conforming Cartesian grid).
205 trans_.clear();
206 trans_.reserve(numElements*3*1.05);
207
208 transBoundary_.clear();
209
210 // if energy is enabled, let's do the same for the "thermal half transmissibilities"
211 if (enableEnergy_) {
212 thermalHalfTrans_.clear();
213 thermalHalfTrans_.reserve(numElements*6*1.05);
214
215 thermalHalfTransBoundary_.clear();
216 }
217
218 // if diffusion is enabled, let's do the same for the "diffusivity"
219 if (updateDiffusivity) {
220 diffusivity_.clear();
221 diffusivity_.reserve(numElements*3*1.05);
222 extractPorosity_();
223 }
224
225 // if dispersion is enabled, let's do the same for the "dispersivity"
226 if (updateDispersivity) {
227 dispersivity_.clear();
228 dispersivity_.reserve(numElements*3*1.05);
229 extractDispersion_();
230 }
231
232 // The MULTZ needs special case if the option is ALL
233 // Then the smallest multiplier is applied.
234 // Default is to apply the top and bottom multiplier
235 bool useSmallestMultiplier;
236 bool pinchActive;
237 if (comm.rank() == 0) {
238 const auto& eclGrid = eclState_.getInputGrid();
239 pinchActive = eclGrid.isPinchActive();
240 useSmallestMultiplier = eclGrid.getMultzOption() == PinchMode::ALL;
241 }
242 if (global && comm.size() > 1) {
243 comm.broadcast(&useSmallestMultiplier, 1, 0);
244 comm.broadcast(&pinchActive, 1, 0);
245 }
246
247 // compute the transmissibilities for all intersections
248 for (const auto& elem : elements(gridView_)) {
249 unsigned elemIdx = elemMapper.index(elem);
250
251 auto isIt = gridView_.ibegin(elem);
252 const auto& isEndIt = gridView_.iend(elem);
253 unsigned boundaryIsIdx = 0;
254 for (; isIt != isEndIt; ++ isIt) {
255 // store intersection, this might be costly
256 const auto& intersection = *isIt;
257
258 // deal with grid boundaries
259 if (intersection.boundary()) {
260 // compute the transmissibilty for the boundary intersection
261 const auto& geometry = intersection.geometry();
262 const auto& faceCenterInside = geometry.center();
263
264 auto faceAreaNormal = intersection.centerUnitOuterNormal();
265 faceAreaNormal *= geometry.volume();
266
267 Scalar transBoundaryIs;
268 computeHalfTrans_(transBoundaryIs,
269 faceAreaNormal,
270 intersection.indexInInside(),
271 distanceVector_(faceCenterInside,
272 intersection.indexInInside(),
273 elemIdx,
274 axisCentroids),
275 permeability_[elemIdx]);
276
277 // normally there would be two half-transmissibilities that would be
278 // averaged. on the grid boundary there only is the half
279 // transmissibility of the interior element.
280 unsigned insideCartElemIdx = cartMapper_.cartesianIndex(elemIdx);
281 applyMultipliers_(transBoundaryIs, intersection.indexInInside(), insideCartElemIdx, transMult);
282 transBoundary_[std::make_pair(elemIdx, boundaryIsIdx)] = transBoundaryIs;
283
284 // for boundary intersections we also need to compute the thermal
285 // half transmissibilities
286 if (enableEnergy_) {
287 Scalar transBoundaryEnergyIs;
288 computeHalfDiffusivity_(transBoundaryEnergyIs,
289 faceAreaNormal,
290 distanceVector_(faceCenterInside,
291 intersection.indexInInside(),
292 elemIdx,
293 axisCentroids),
294 1.0);
295 thermalHalfTransBoundary_[std::make_pair(elemIdx, boundaryIsIdx)] =
296 transBoundaryEnergyIs;
297 }
298
299 ++ boundaryIsIdx;
300 continue;
301 }
302
303 if (!intersection.neighbor()) {
304 // elements can be on process boundaries, i.e. they are not on the
305 // domain boundary yet they don't have neighbors.
306 ++ boundaryIsIdx;
307 continue;
308 }
309
310 const auto& outsideElem = intersection.outside();
311 unsigned outsideElemIdx = elemMapper.index(outsideElem);
312
313 // Get the Cartesian indices of the origen cells (parent or equivalent cell on level zero), for CpGrid with LGRs.
314 // For genral grids and no LGRs, get the usual Cartesian Index.
315 unsigned insideCartElemIdx = this-> lookUpCartesianData_.template getFieldPropCartesianIdx<Grid>(elemIdx);
316 unsigned outsideCartElemIdx = this-> lookUpCartesianData_.template getFieldPropCartesianIdx<Grid>(outsideElemIdx);
317
318 // we only need to calculate a face's transmissibility
319 // once...
320 // In a parallel run insideCartElemIdx>outsideCartElemIdx does not imply elemIdx>outsideElemIdx for
321 // ghost cells and we need to use the cartesian index as this will be used when applying Z multipliers
322 // To cover the case where both cells are part of an LGR and as a consequence might have
323 // the same cartesian index, we tie their Cartesian indices and the ones on the leaf grid view.
324 if (std::tie(insideCartElemIdx, elemIdx) > std::tie(outsideCartElemIdx, outsideElemIdx))
325 continue;
326
327 // local indices of the faces of the inside and
328 // outside elements which contain the intersection
329 int insideFaceIdx = intersection.indexInInside();
330 int outsideFaceIdx = intersection.indexInOutside();
331
332 if (insideFaceIdx == -1) {
333 // NNC. Set zero transmissibility, as it will be
334 // *added to* by applyNncToGridTrans_() later.
335 assert(outsideFaceIdx == -1);
336 trans_[details::isId(elemIdx, outsideElemIdx)] = 0.0;
337 if (enableEnergy_){
338 thermalHalfTrans_[details::directionalIsId(elemIdx, outsideElemIdx)] = 0.0;
339 thermalHalfTrans_[details::directionalIsId(outsideElemIdx, elemIdx)] = 0.0;
340 }
341
342 if (updateDiffusivity) {
343 diffusivity_[details::isId(elemIdx, outsideElemIdx)] = 0.0;
344 }
345 if (updateDispersivity) {
346 dispersivity_[details::isId(elemIdx, outsideElemIdx)] = 0.0;
347 }
348 continue;
349 }
350
351 DimVector faceCenterInside;
352 DimVector faceCenterOutside;
353 DimVector faceAreaNormal;
354
355 typename std::is_same<Grid, Dune::CpGrid>::type isCpGrid;
356 computeFaceProperties(intersection,
357 elemIdx,
358 insideFaceIdx,
359 outsideElemIdx,
360 outsideFaceIdx,
361 faceCenterInside,
362 faceCenterOutside,
363 faceAreaNormal,
364 isCpGrid);
365
366 Scalar halfTrans1;
367 Scalar halfTrans2;
368
369 computeHalfTrans_(halfTrans1,
370 faceAreaNormal,
371 insideFaceIdx,
372 distanceVector_(faceCenterInside,
373 intersection.indexInInside(),
374 elemIdx,
375 axisCentroids),
376 permeability_[elemIdx]);
377 computeHalfTrans_(halfTrans2,
378 faceAreaNormal,
379 outsideFaceIdx,
380 distanceVector_(faceCenterOutside,
381 intersection.indexInOutside(),
382 outsideElemIdx,
383 axisCentroids),
384 permeability_[outsideElemIdx]);
385
386 applyNtg_(halfTrans1, insideFaceIdx, elemIdx, ntg);
387 applyNtg_(halfTrans2, outsideFaceIdx, outsideElemIdx, ntg);
388
389 // convert half transmissibilities to full face
390 // transmissibilities using the harmonic mean
391 Scalar trans;
392 if (std::abs(halfTrans1) < 1e-30 || std::abs(halfTrans2) < 1e-30)
393 // avoid division by zero
394 trans = 0.0;
395 else
396 trans = 1.0 / (1.0/halfTrans1 + 1.0/halfTrans2);
397
398 // apply the full face transmissibility multipliers
399 // for the inside ...
400 if(!pinchActive){
401 if (insideFaceIdx > 3){// top or bottom
402 auto find_layer = [&cartDims](std::size_t cell){
403 cell /= cartDims[0];
404 auto k = cell / cartDims[1];
405 return k;
406 };
407 int kup = find_layer(insideCartElemIdx);
408 int kdown=find_layer(outsideCartElemIdx);
409 // When a grid is a CpGrid with LGRs, insideCartElemIdx coincides with outsideCartElemIdx
410 // for cells on the leaf with the same parent cell on level zero.
411 assert((kup != kdown) || (insideCartElemIdx == outsideCartElemIdx));
412 if(std::abs(kup -kdown) > 1){
413 trans = 0.0;
414 }
415 }
416 }
417
418 if (useSmallestMultiplier)
419 {
420 // Currently PINCH(4) is never queries and hence PINCH(4) == TOPBOT is assumed
421 // and in this branch PINCH(5) == ALL holds
422 applyAllZMultipliers_(trans, insideFaceIdx, outsideFaceIdx, insideCartElemIdx,
423 outsideCartElemIdx, transMult, cartDims,
424 /* pinchTop= */ false);
425 }
426 else
427 {
428 applyMultipliers_(trans, insideFaceIdx, insideCartElemIdx, transMult);
429 // ... and outside elements
430 applyMultipliers_(trans, outsideFaceIdx, outsideCartElemIdx, transMult);
431 }
432
433 // apply the region multipliers (cf. the MULTREGT keyword)
434 FaceDir::DirEnum faceDir;
435 switch (insideFaceIdx) {
436 case 0:
437 case 1:
438 faceDir = FaceDir::XPlus;
439 break;
440
441 case 2:
442 case 3:
443 faceDir = FaceDir::YPlus;
444 break;
445
446 case 4:
447 case 5:
448 faceDir = FaceDir::ZPlus;
449 break;
450
451 default:
452 throw std::logic_error("Could not determine a face direction");
453 }
454
455 trans *= transMult.getRegionMultiplier(insideCartElemIdx,
456 outsideCartElemIdx,
457 faceDir);
458
459 trans_[details::isId(elemIdx, outsideElemIdx)] = trans;
460
461 // update the "thermal half transmissibility" for the intersection
462 if (enableEnergy_) {
463
464 Scalar halfDiffusivity1;
465 Scalar halfDiffusivity2;
466
467 computeHalfDiffusivity_(halfDiffusivity1,
468 faceAreaNormal,
469 distanceVector_(faceCenterInside,
470 intersection.indexInInside(),
471 elemIdx,
472 axisCentroids),
473 1.0);
474 computeHalfDiffusivity_(halfDiffusivity2,
475 faceAreaNormal,
476 distanceVector_(faceCenterOutside,
477 intersection.indexInOutside(),
478 outsideElemIdx,
479 axisCentroids),
480 1.0);
481 //TODO Add support for multipliers
482 thermalHalfTrans_[details::directionalIsId(elemIdx, outsideElemIdx)] = halfDiffusivity1;
483 thermalHalfTrans_[details::directionalIsId(outsideElemIdx, elemIdx)] = halfDiffusivity2;
484 }
485
486 // update the "diffusive half transmissibility" for the intersection
487 if (updateDiffusivity) {
488
489 Scalar halfDiffusivity1;
490 Scalar halfDiffusivity2;
491
492 computeHalfDiffusivity_(halfDiffusivity1,
493 faceAreaNormal,
494 distanceVector_(faceCenterInside,
495 intersection.indexInInside(),
496 elemIdx,
497 axisCentroids),
498 porosity_[elemIdx]);
499 computeHalfDiffusivity_(halfDiffusivity2,
500 faceAreaNormal,
501 distanceVector_(faceCenterOutside,
502 intersection.indexInOutside(),
503 outsideElemIdx,
504 axisCentroids),
505 porosity_[outsideElemIdx]);
506
507 applyNtg_(halfDiffusivity1, insideFaceIdx, elemIdx, ntg);
508 applyNtg_(halfDiffusivity2, outsideFaceIdx, outsideElemIdx, ntg);
509
510 //TODO Add support for multipliers
511 Scalar diffusivity;
512 if (std::abs(halfDiffusivity1) < 1e-30 || std::abs(halfDiffusivity2) < 1e-30)
513 // avoid division by zero
514 diffusivity = 0.0;
515 else
516 diffusivity = 1.0 / (1.0/halfDiffusivity1 + 1.0/halfDiffusivity2);
517
518
519 diffusivity_[details::isId(elemIdx, outsideElemIdx)] = diffusivity;
520 }
521
522 // update the "dispersivity half transmissibility" for the intersection
523 if (updateDispersivity) {
524
525 Scalar halfDispersivity1;
526 Scalar halfDispersivity2;
527
528 computeHalfDiffusivity_(halfDispersivity1,
529 faceAreaNormal,
530 distanceVector_(faceCenterInside,
531 intersection.indexInInside(),
532 elemIdx,
533 axisCentroids),
534 dispersion_[elemIdx]);
535 computeHalfDiffusivity_(halfDispersivity2,
536 faceAreaNormal,
537 distanceVector_(faceCenterOutside,
538 intersection.indexInOutside(),
539 outsideElemIdx,
540 axisCentroids),
541 dispersion_[outsideElemIdx]);
542
543 applyNtg_(halfDispersivity1, insideFaceIdx, elemIdx, ntg);
544 applyNtg_(halfDispersivity2, outsideFaceIdx, outsideElemIdx, ntg);
545
546 //TODO Add support for multipliers
547 Scalar dispersivity;
548 if (std::abs(halfDispersivity1) < 1e-30 || std::abs(halfDispersivity2) < 1e-30)
549 // avoid division by zero
550 dispersivity = 0.0;
551 else
552 dispersivity = 1.0 / (1.0/halfDispersivity1 + 1.0/halfDispersivity2);
553
554
555 dispersivity_[details::isId(elemIdx, outsideElemIdx)] = dispersivity;
556 }
557 }
558 }
559
560 // Potentially overwrite and/or modify transmissibilities based on input from deck
561 this->updateFromEclState_(global);
562
563 // Create mapping from global to local index
564 std::unordered_map<std::size_t,int> globalToLocal;
565
566 // Loop over all elements (global grid) and store Cartesian index
567 for (const auto& elem : elements(grid_.leafGridView())) {
568 int elemIdx = elemMapper.index(elem);
569 int cartElemIdx = cartMapper_.cartesianIndex(elemIdx);
570 globalToLocal[cartElemIdx] = elemIdx;
571 }
572
573 if (!disableNNC) {
574 this->applyEditNncToGridTrans_(globalToLocal);
575 this->applyNncToGridTrans_(globalToLocal);
576 this->applyEditNncrToGridTrans_(globalToLocal);
577 if (applyNncMultregT) {
578 this->applyNncMultreg_(globalToLocal);
579 }
580 }
581
582 // If disableNNC == true, remove all non-neighbouring transmissibilities.
583 // If disableNNC == false, remove very small non-neighbouring transmissibilities.
584 this->removeNonCartesianTransmissibilities_(disableNNC);
585}
586
587template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
590{
591 unsigned numElem = gridView_.size(/*codim=*/0);
592 permeability_.resize(numElem);
593
594 // read the intrinsic permeabilities from the eclState. Note that all arrays
595 // provided by eclState are one-per-cell of "uncompressed" grid, whereas the
596 // simulation grid might remove a few elements. (e.g. because it is distributed
597 // over several processes.)
598 const auto& fp = eclState_.fieldProps();
599 if (fp.has_double("PERMX")) {
600 const std::vector<double>& permxData = this-> lookUpData_.assignFieldPropsDoubleOnLeaf(fp, "PERMX");
601
602 std::vector<double> permyData;
603 if (fp.has_double("PERMY"))
604 permyData = this-> lookUpData_.assignFieldPropsDoubleOnLeaf(fp,"PERMY");
605 else
606 permyData = permxData;
607
608 std::vector<double> permzData;
609 if (fp.has_double("PERMZ"))
610 permzData = this-> lookUpData_.assignFieldPropsDoubleOnLeaf(fp,"PERMZ");
611 else
612 permzData = permxData;
613
614 for (std::size_t dofIdx = 0; dofIdx < numElem; ++ dofIdx) {
615 permeability_[dofIdx] = 0.0;
616 permeability_[dofIdx][0][0] = permxData[dofIdx];
617 permeability_[dofIdx][1][1] = permyData[dofIdx];
618 permeability_[dofIdx][2][2] = permzData[dofIdx];
619 }
620
621 // for now we don't care about non-diagonal entries
622
623 }
624 else
625 throw std::logic_error("Can't read the intrinsic permeability from the ecl state. "
626 "(The PERM{X,Y,Z} keywords are missing)");
627}
628
629template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
631extractPermeability_(const std::function<unsigned int(unsigned int)>& map)
632{
633 unsigned numElem = gridView_.size(/*codim=*/0);
634 permeability_.resize(numElem);
635
636 // read the intrinsic permeabilities from the eclState. Note that all arrays
637 // provided by eclState are one-per-cell of "uncompressed" grid, whereas the
638 // simulation grid might remove a few elements. (e.g. because it is distributed
639 // over several processes.)
640 const auto& fp = eclState_.fieldProps();
641 if (fp.has_double("PERMX")) {
642 const std::vector<double>& permxData = this-> lookUpData_.assignFieldPropsDoubleOnLeaf(fp,"PERMX");
643
644 std::vector<double> permyData;
645 if (fp.has_double("PERMY"))
646 permyData = this-> lookUpData_.assignFieldPropsDoubleOnLeaf(fp,"PERMY");
647 else
648 permyData = permxData;
649
650 std::vector<double> permzData;
651 if (fp.has_double("PERMZ"))
652 permzData = this-> lookUpData_.assignFieldPropsDoubleOnLeaf(fp,"PERMZ");
653 else
654 permzData = permxData;
655
656 for (std::size_t dofIdx = 0; dofIdx < numElem; ++ dofIdx) {
657 permeability_[dofIdx] = 0.0;
658 std::size_t inputDofIdx = map(dofIdx);
659 permeability_[dofIdx][0][0] = permxData[inputDofIdx];
660 permeability_[dofIdx][1][1] = permyData[inputDofIdx];
661 permeability_[dofIdx][2][2] = permzData[inputDofIdx];
662 }
663
664 // for now we don't care about non-diagonal entries
665
666 }
667 else
668 throw std::logic_error("Can't read the intrinsic permeability from the ecl state. "
669 "(The PERM{X,Y,Z} keywords are missing)");
670}
671
672template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
675{
676 // read the intrinsic porosity from the eclState. Note that all arrays
677 // provided by eclState are one-per-cell of "uncompressed" grid, whereas the
678 // simulation grid might remove a few elements. (e.g. because it is distributed
679 // over several processes.)
680 const auto& fp = eclState_.fieldProps();
681 if (fp.has_double("PORO")) {
682 porosity_ = this-> lookUpData_.assignFieldPropsDoubleOnLeaf(fp,"PORO");
683 }
684 else
685 throw std::logic_error("Can't read the porosityfrom the ecl state. "
686 "(The PORO keywords are missing)");
687}
688
689template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
692{
693 if (!enableDispersivity_) {
694 throw std::runtime_error("Dispersion disabled at compile time, but the deck "
695 "contains the DISPERC keyword.");
696 }
697 const auto& fp = eclState_.fieldProps();
698 dispersion_ = this-> lookUpData_.assignFieldPropsDoubleOnLeaf(fp,"DISPERC");
699}
700
701template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
704{
705 const auto& cartDims = cartMapper_.cartesianDimensions();
706 for (auto&& trans: trans_) {
707 //either remove all NNC transmissibilities or those less than the threshold (by default 1e-6 in the deck's unit system)
708 if (removeAll or trans.second < transmissibilityThreshold_) {
709 const auto& id = trans.first;
710 const auto& elements = details::isIdReverse(id);
711 int gc1 = std::min(cartMapper_.cartesianIndex(elements.first), cartMapper_.cartesianIndex(elements.second));
712 int gc2 = std::max(cartMapper_.cartesianIndex(elements.first), cartMapper_.cartesianIndex(elements.second));
713
714 // only adjust the NNCs
715 // When LGRs, all neighbors in the LGR are cartesian neighbours on the level grid representing the LGR.
716 // When elements on the leaf grid view have the same parent cell, gc1 and gc2 coincide.
717 if (gc2 - gc1 == 1 || gc2 - gc1 == cartDims[0] || gc2 - gc1 == cartDims[0]*cartDims[1] || gc2 - gc1 == 0)
718 continue;
719
720 trans.second = 0.0;
721 }
722 }
723}
724
725template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
727applyAllZMultipliers_(Scalar& trans,
728 unsigned insideFaceIdx,
729 unsigned outsideFaceIdx,
730 unsigned insideCartElemIdx,
731 unsigned outsideCartElemIdx,
732 const TransMult& transMult,
733 const std::array<int, dimWorld>& cartDims,
734 bool pinchTop)
735{
736 if(grid_.maxLevel()> 0) {
737 OPM_THROW(std::invalid_argument, "MULTZ not support with LGRS, yet.");
738 }
739 if (insideFaceIdx > 3) { // top or or bottom
740 assert(insideFaceIdx==5); // as insideCartElemIdx < outsideCartElemIdx holds for the Z column
741 // For CpGrid with LGRs, insideCartElemIdx == outsideCartElemIdx when cells on the leaf have the same parent cell on level zero.
742 assert(outsideCartElemIdx >= insideCartElemIdx);
743 unsigned lastCartElemIdx;
744 if (outsideCartElemIdx == insideCartElemIdx) {
745 lastCartElemIdx = outsideCartElemIdx;
746 }
747 else {
748 lastCartElemIdx = outsideCartElemIdx - cartDims[0]*cartDims[1];
749 }
750 // Last multiplier using (Z+)*(Z-)
751 Scalar mult = transMult.getMultiplier(lastCartElemIdx , FaceDir::ZPlus) *
752 transMult.getMultiplier(outsideCartElemIdx , FaceDir::ZMinus);
753
754 if ( !pinchTop )
755 {
756 // pick the smallest multiplier using (Z+)*(Z-) while looking down
757 // the pillar until reaching the other end of the connection
758 for(auto cartElemIdx = insideCartElemIdx; cartElemIdx < lastCartElemIdx;)
759 {
760 auto multiplier = transMult.getMultiplier(cartElemIdx, FaceDir::ZPlus);
761 cartElemIdx += cartDims[0]*cartDims[1];
762 multiplier *= transMult.getMultiplier(cartElemIdx, FaceDir::ZMinus);
763 mult = std::min(mult, multiplier);
764 }
765 }
766
767 trans *= mult;
768 }
769 else
770 {
771 applyMultipliers_(trans, insideFaceIdx, insideCartElemIdx, transMult);
772 applyMultipliers_(trans, outsideFaceIdx, outsideCartElemIdx, transMult);
773 }
774}
775
776template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
778updateFromEclState_(bool global)
779{
780 const FieldPropsManager* fp =
781 (global) ? &(eclState_.fieldProps()) :
782 &(eclState_.globalFieldProps());
783
784 std::array<bool,3> is_tran {fp->tran_active("TRANX"),
785 fp->tran_active("TRANY"),
786 fp->tran_active("TRANZ")};
787
788 if( !(is_tran[0] ||is_tran[1] || is_tran[2]) )
789 {
790 // Skip unneeded expensive traversals
791 return;
792 }
793
794 std::array<std::string, 3> keywords {"TRANX", "TRANY", "TRANZ"};
795 std::array<std::vector<double>,3> trans = createTransmissibilityArrays_(is_tran);
796 auto key = keywords.begin();
797 auto perform = is_tran.begin();
798
799 for (auto it = trans.begin(); it != trans.end(); ++it, ++key, ++perform)
800 {
801 if(*perform) {
802 if(grid_.maxLevel()>0) {
803 OPM_THROW(std::invalid_argument, "Calculations on TRANX/TRANY/TRANZ arrays are not support with LGRS, yet.");
804 }
805 fp->apply_tran(*key, *it);
806 }
807 }
808
809 resetTransmissibilityFromArrays_(is_tran, trans);
810}
811
812template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
813std::array<std::vector<double>,3>
815createTransmissibilityArrays_(const std::array<bool,3>& is_tran)
816{
817 const auto& cartDims = cartMapper_.cartesianDimensions();
818 ElementMapper elemMapper(gridView_, Dune::mcmgElementLayout());
819
820 auto numElem = gridView_.size(/*codim=*/0);
821 std::array<std::vector<double>,3> trans =
822 { std::vector<double>(is_tran[0] ? numElem : 0, 0),
823 std::vector<double>(is_tran[1] ? numElem : 0, 0),
824 std::vector<double>(is_tran[2] ? numElem : 0, 0)};
825
826 // compute the transmissibilities for all intersections
827 for (const auto& elem : elements(gridView_)) {
828 for (const auto& intersection : intersections(gridView_, elem)) {
829 // store intersection, this might be costly
830 if (!intersection.neighbor())
831 continue; // intersection is on the domain boundary
832
833 // In the EclState TRANX[c1] is transmissibility in X+
834 // direction. Ordering of compressed (c1,c2) and cartesian index
835 // (gc1, gc2) is coherent (c1 < c2 <=> gc1 < gc2). This also
836 // holds for the global grid. While distributing changes the
837 // order of the local indices, the transmissibilities are still
838 // stored at the cell with the lower global cartesian index as
839 // the fieldprops are communicated by the grid.
840 unsigned c1 = elemMapper.index(intersection.inside());
841 unsigned c2 = elemMapper.index(intersection.outside());
842 int gc1 = cartMapper_.cartesianIndex(c1);
843 int gc2 = cartMapper_.cartesianIndex(c2);
844
845 if (c1 > c2)
846 continue; // we only need to handle each connection once, thank you.
847
848 auto isID = details::isId(c1, c2);
849
850 // For CpGrid with LGRs, when leaf grid view cells with indices c1 and c2
851 // have the same parent cell on level zero, then gc2 - gc1 == 0. In that case,
852 // 'intersection.indexInSIde()' needed to be checked to determine the direction, i.e.
853 // add in the if/else-if 'gc2 == gc1 && intersection.indexInInside() == ... '
854 if ((gc2 - gc1 == 1 || (gc2 == gc1 && (intersection.indexInInside() == 0 || intersection.indexInInside() == 1)))
855 && cartDims[0] > 1) {
856 if (is_tran[0])
857 // set simulator internal transmissibilities to values from inputTranx
858 trans[0][c1] = trans_[isID];
859 }
860 else if ((gc2 - gc1 == cartDims[0] || (gc2 == gc1 && (intersection.indexInInside() == 2 || intersection.indexInInside() == 3)))
861 && cartDims[1] > 1) {
862 if (is_tran[1])
863 // set simulator internal transmissibilities to values from inputTrany
864 trans[1][c1] = trans_[isID];
865 }
866 else if (gc2 - gc1 == cartDims[0]*cartDims[1] ||
867 (gc2 == gc1 && (intersection.indexInInside() == 4 || intersection.indexInInside() == 5))) {
868 if (is_tran[2])
869 // set simulator internal transmissibilities to values from inputTranz
870 trans[2][c1] = trans_[isID];
871 }
872
873 //else.. We don't support modification of NNC at the moment.
874 }
875 }
876
877 return trans;
878}
879
880template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
882resetTransmissibilityFromArrays_(const std::array<bool,3>& is_tran,
883 const std::array<std::vector<double>,3>& trans)
884{
885 const auto& cartDims = cartMapper_.cartesianDimensions();
886 ElementMapper elemMapper(gridView_, Dune::mcmgElementLayout());
887
888 // compute the transmissibilities for all intersections
889 for (const auto& elem : elements(gridView_)) {
890 for (const auto& intersection : intersections(gridView_, elem)) {
891 if (!intersection.neighbor())
892 continue; // intersection is on the domain boundary
893
894 // In the EclState TRANX[c1] is transmissibility in X+
895 // direction. Ordering of compressed (c1,c2) and cartesian index
896 // (gc1, gc2) is coherent (c1 < c2 <=> gc1 < gc2). This also
897 // holds for the global grid. While distributing changes the
898 // order of the local indices, the transmissibilities are still
899 // stored at the cell with the lower global cartesian index as
900 // the fieldprops are communicated by the grid.
904 unsigned c1 = elemMapper.index(intersection.inside());
905 unsigned c2 = elemMapper.index(intersection.outside());
906 int gc1 = cartMapper_.cartesianIndex(c1);
907 int gc2 = cartMapper_.cartesianIndex(c2);
908 if (c1 > c2)
909 continue; // we only need to handle each connection once, thank you.
910
911 auto isID = details::isId(c1, c2);
912
913 // For CpGrid with LGRs, when leaf grid view cells with indices c1 and c2
914 // have the same parent cell on level zero, then gc2 - gc1 == 0. In that case,
915 // 'intersection.indexInSIde()' needed to be checked to determine the direction, i.e.
916 // add in the if/else-if 'gc2 == gc1 && intersection.indexInInside() == ... '
917 if ((gc2 - gc1 == 1 || (gc2 == gc1 && (intersection.indexInInside() == 0 || intersection.indexInInside() == 1)))
918 && cartDims[0] > 1) {
919 if (is_tran[0])
920 // set simulator internal transmissibilities to values from inputTranx
921 trans_[isID] = trans[0][c1];
922 }
923 else if ((gc2 - gc1 == cartDims[0] || (gc2 == gc1 && (intersection.indexInInside() == 2|| intersection.indexInInside() == 3)))
924 && cartDims[1] > 1) {
925 if (is_tran[1])
926 // set simulator internal transmissibilities to values from inputTrany
927 trans_[isID] = trans[1][c1];
928 }
929 else if (gc2 - gc1 == cartDims[0]*cartDims[1] ||
930 (gc2 == gc1 && (intersection.indexInInside() == 4 || intersection.indexInInside() == 5))) {
931 if (is_tran[2])
932 // set simulator internal transmissibilities to values from inputTranz
933 trans_[isID] = trans[2][c1];
934 }
935
936 //else.. We don't support modification of NNC at the moment.
937 }
938 }
939}
940
941template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
942template<class Intersection>
944computeFaceProperties(const Intersection& intersection,
945 const int,
946 const int,
947 const int,
948 const int,
949 DimVector& faceCenterInside,
950 DimVector& faceCenterOutside,
951 DimVector& faceAreaNormal,
952 /*isCpGrid=*/std::false_type) const
953{
954 // default implementation for DUNE grids
955 const auto& geometry = intersection.geometry();
956 faceCenterInside = geometry.center();
957 faceCenterOutside = faceCenterInside;
958
959 faceAreaNormal = intersection.centerUnitOuterNormal();
960 faceAreaNormal *= geometry.volume();
961}
962
963template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
964template<class Intersection>
966computeFaceProperties(const Intersection& intersection,
967 const int insideElemIdx,
968 const int insideFaceIdx,
969 const int outsideElemIdx,
970 const int outsideFaceIdx,
971 DimVector& faceCenterInside,
972 DimVector& faceCenterOutside,
973 DimVector& faceAreaNormal,
974 /*isCpGrid=*/std::true_type) const
975{
976 int faceIdx = intersection.id();
977
978 if(grid_.maxLevel() == 0) {
979 faceCenterInside = grid_.faceCenterEcl(insideElemIdx, insideFaceIdx, intersection);
980 faceCenterOutside = grid_.faceCenterEcl(outsideElemIdx, outsideFaceIdx, intersection);
981 faceAreaNormal = grid_.faceAreaNormalEcl(faceIdx);
982 }
983 else {
984 if ((intersection.inside().level() != intersection.outside().level())) {
985
986 // For CpGrid with LGRs, intersection laying on the boundary of an LGR, having two neighboring cells:
987 // one coarse neighboring cell and one refined neighboring cell, we get the corresponding parent
988 // intersection (from level 0), and use the center of the parent intersection for the coarse
989 // neighboring cell.
990
991 // Get parent intersection and its geometry
992 const auto& parentIntersection = grid_.getParentIntersectionFromLgrBoundaryFace(intersection);
993 const auto& parentIntersectionGeometry = parentIntersection.geometry();
994
995 // For the coarse neighboring cell, take the center of the parent intersection.
996 // For the refined neighboring cell, take the 'usual' center.
997 faceCenterInside = (intersection.inside().level() == 0) ? parentIntersectionGeometry.center() :
998 grid_.faceCenterEcl(insideElemIdx, insideFaceIdx, intersection);
999 faceCenterOutside = (intersection.outside().level() == 0) ? parentIntersectionGeometry.center() :
1000 grid_.faceCenterEcl(outsideElemIdx, outsideFaceIdx, intersection);
1001
1002 // For some computations, it seems to be benefitial to replace the actual area of the refined face, by
1003 // the area of its parent face.
1004 // faceAreaNormal = parentIntersection.centerUnitOuterNormal();
1005 // faceAreaNormal *= parentIntersectionGeometry.volume();
1006
1008 faceAreaNormal = intersection.centerUnitOuterNormal();
1009 faceAreaNormal *= intersection.geometry().volume();
1010 }
1011 else {
1012 assert(intersection.inside().level() == intersection.outside().level());
1013
1014 faceCenterInside = grid_.faceCenterEcl(insideElemIdx, insideFaceIdx, intersection);
1015 faceCenterOutside = grid_.faceCenterEcl(outsideElemIdx, outsideFaceIdx, intersection);
1016
1017 // When the CpGrid has LGRs, we compute the face area normal differently.
1018 if (intersection.inside().level() > 0) { // remove intersection.inside().level() > 0
1019 faceAreaNormal = intersection.centerUnitOuterNormal();
1020 faceAreaNormal *= intersection.geometry().volume();
1021 }
1022 else {
1023 faceAreaNormal = grid_.faceAreaNormalEcl(faceIdx);
1024 }
1025 }
1026 }
1027}
1028
1029template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
1030void
1032applyNncToGridTrans_(const std::unordered_map<std::size_t,int>& cartesianToCompressed)
1033{
1034 // First scale NNCs with EDITNNC.
1035 const auto& nnc_input = eclState_.getInputNNC().input();
1036
1037 for (const auto& nncEntry : nnc_input) {
1038 auto c1 = nncEntry.cell1;
1039 auto c2 = nncEntry.cell2;
1040 auto lowIt = cartesianToCompressed.find(c1);
1041 auto highIt = cartesianToCompressed.find(c2);
1042 int low = (lowIt == cartesianToCompressed.end())? -1 : lowIt->second;
1043 int high = (highIt == cartesianToCompressed.end())? -1 : highIt->second;
1044
1045 if (low > high)
1046 std::swap(low, high);
1047
1048 if (low == -1 && high == -1)
1049 // Silently discard as it is not between active cells
1050 continue;
1051
1052 if (low == -1 || high == -1) {
1053 // Discard the NNC if it is between active cell and inactive cell
1054 std::ostringstream sstr;
1055 sstr << "NNC between active and inactive cells ("
1056 << low << " -> " << high << ") with globalcell is (" << c1 << "->" << c2 <<")";
1057 OpmLog::warning(sstr.str());
1058 continue;
1059 }
1060
1061 {
1062 auto candidate = trans_.find(details::isId(low, high));
1063 if (candidate != trans_.end()) {
1064 // NNC is represented by the grid and might be a neighboring connection
1065 // In this case the transmissibilty is added to the value already
1066 // set or computed.
1067 candidate->second += nncEntry.trans;
1068 }
1069 }
1070 // if (enableEnergy_) {
1071 // auto candidate = thermalHalfTrans_.find(details::directionalIsId(low, high));
1072 // if (candidate != trans_.end()) {
1073 // // NNC is represented by the grid and might be a neighboring connection
1074 // // In this case the transmissibilty is added to the value already
1075 // // set or computed.
1076 // candidate->second += nncEntry.transEnergy1;
1077 // }
1078 // auto candidate = thermalHalfTrans_.find(details::directionalIsId(high, low));
1079 // if (candidate != trans_.end()) {
1080 // // NNC is represented by the grid and might be a neighboring connection
1081 // // In this case the transmissibilty is added to the value already
1082 // // set or computed.
1083 // candidate->second += nncEntry.transEnergy2;
1084 // }
1085 // }
1086 // if (enableDiffusivity_) {
1087 // auto candidate = diffusivity_.find(details::isId(low, high));
1088 // if (candidate != trans_.end()) {
1089 // // NNC is represented by the grid and might be a neighboring connection
1090 // // In this case the transmissibilty is added to the value already
1091 // // set or computed.
1092 // candidate->second += nncEntry.transDiffusion;
1093 // }
1094 // }
1095 }
1096}
1097
1098template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
1100applyEditNncToGridTrans_(const std::unordered_map<std::size_t,int>& globalToLocal)
1101{
1102 const auto& input = eclState_.getInputNNC();
1103 applyEditNncToGridTransHelper_(globalToLocal, "EDITNNC",
1104 input.edit(),
1105 [&input](const NNCdata& nnc){
1106 return input.edit_location(nnc);},
1107 // Multiply transmissibility with EDITNNC value
1108 [](double& trans, const double& rhs){ trans *= rhs;});
1109}
1110
1111template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
1113applyEditNncrToGridTrans_(const std::unordered_map<std::size_t,int>& globalToLocal)
1114{
1115 const auto& input = eclState_.getInputNNC();
1116 applyEditNncToGridTransHelper_(globalToLocal, "EDITNNCR",
1117 input.editr(),
1118 [&input](const NNCdata& nnc){
1119 return input.editr_location(nnc);},
1120 // Replace Transmissibility with EDITNNCR value
1121 [](double& trans, const double& rhs){ trans = rhs;});
1122}
1123
1124template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
1126applyEditNncToGridTransHelper_(const std::unordered_map<std::size_t,int>& globalToLocal,
1127 const std::string& keyword,
1128 const std::vector<NNCdata>& nncs,
1129 const std::function<KeywordLocation(const NNCdata&)>& getLocation,
1130 const std::function<void(double&, const double&)>& apply)
1131{
1132 if (nncs.empty())
1133 return;
1134 const auto& cartDims = cartMapper_.cartesianDimensions();
1135
1136 auto format_ijk = [&cartDims](std::size_t cell) -> std::string {
1137 auto i = cell % cartDims[0]; cell /= cartDims[0];
1138 auto j = cell % cartDims[1];
1139 auto k = cell / cartDims[1];
1140
1141 return fmt::format("({},{},{})", i + 1,j + 1,k + 1);
1142 };
1143
1144 auto print_warning = [&format_ijk, &getLocation, &keyword] (const NNCdata& nnc) {
1145 const auto& location = getLocation( nnc );
1146 auto warning = fmt::format("Problem with {} keyword\n"
1147 "In {} line {} \n"
1148 "No NNC defined for connection {} -> {}", keyword, location.filename,
1149 location.lineno, format_ijk(nnc.cell1), format_ijk(nnc.cell2));
1150 OpmLog::warning(keyword, warning);
1151 };
1152
1153 // editNnc is supposed to only reference non-neighboring connections and not
1154 // neighboring connections. Use all entries for scaling if there is an NNC.
1155 // variable nnc incremented in loop body.
1156 auto nnc = nncs.begin();
1157 auto end = nncs.end();
1158 std::size_t warning_count = 0;
1159 while (nnc != end) {
1160 auto c1 = nnc->cell1;
1161 auto c2 = nnc->cell2;
1162 auto lowIt = globalToLocal.find(c1);
1163 auto highIt = globalToLocal.find(c2);
1164
1165 if (lowIt == globalToLocal.end() || highIt == globalToLocal.end()) {
1166 print_warning(*nnc);
1167 ++nnc;
1168 warning_count++;
1169 continue;
1170 }
1171
1172 auto low = lowIt->second, high = highIt->second;
1173
1174 if (low > high)
1175 std::swap(low, high);
1176
1177 auto candidate = trans_.find(details::isId(low, high));
1178 if (candidate == trans_.end()) {
1179 print_warning(*nnc);
1180 ++nnc;
1181 warning_count++;
1182 }
1183 else {
1184 // NNC exists
1185 while (nnc!= end && c1==nnc->cell1 && c2==nnc->cell2) {
1186 apply(candidate->second, nnc->trans);
1187 ++nnc;
1188 }
1189 }
1190 }
1191
1192 if (warning_count > 0) {
1193 auto warning = fmt::format("Problems with {} keyword\n"
1194 "A total of {} connections not defined in grid", keyword, warning_count);
1195 OpmLog::warning(warning);
1196 }
1197}
1198
1199template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
1200void
1202applyNncMultreg_(const std::unordered_map<std::size_t,int>& cartesianToCompressed)
1203{
1204 const auto& inputNNC = this->eclState_.getInputNNC();
1205 const auto& transMult = this->eclState_.getTransMult();
1206
1207 auto compressedIdx = [&cartesianToCompressed](const std::size_t globIdx)
1208 {
1209 auto ixPos = cartesianToCompressed.find(globIdx);
1210 return (ixPos == cartesianToCompressed.end()) ? -1 : ixPos->second;
1211 };
1212
1213 // Apply region-based transmissibility multipliers (i.e., the MULTREGT
1214 // keyword) to those transmissibilities that are directly assigned from
1215 // the input.
1216 //
1217 // * NNC::input() covers the NNC keyword and any numerical aquifers
1218 // * NNC::editr() covers the EDITNNCR keyword
1219 //
1220 // Note: We do not apply MULTREGT to the entries in NNC::edit() since
1221 // those act as regular multipliers and have already been fully
1222 // accounted for in the multiplier part of the main loop of update() and
1223 // the applyEditNncToGridTrans_() member function.
1224 for (const auto& nncList : { &NNC::input, &NNC::editr }) {
1225 for (const auto& nncEntry : (inputNNC.*nncList)()) {
1226 const auto c1 = nncEntry.cell1;
1227 const auto c2 = nncEntry.cell2;
1228
1229 auto low = compressedIdx(c1);
1230 auto high = compressedIdx(c2);
1231
1232 if ((low == -1) || (high == -1)) {
1233 continue;
1234 }
1235
1236 if (low > high) {
1237 std::swap(low, high);
1238 }
1239
1240 auto candidate = this->trans_.find(details::isId(low, high));
1241 if (candidate != this->trans_.end()) {
1242 candidate->second *= transMult.getRegionMultiplierNNC(c1, c2);
1243 }
1244 }
1245 }
1246}
1247
1248template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
1250computeHalfTrans_(Scalar& halfTrans,
1251 const DimVector& areaNormal,
1252 int faceIdx, // in the reference element that contains the intersection
1253 const DimVector& distance,
1254 const DimMatrix& perm) const
1255{
1256 assert(faceIdx >= 0);
1257 unsigned dimIdx = faceIdx/2;
1258 assert(dimIdx < dimWorld);
1259 halfTrans = perm[dimIdx][dimIdx];
1260
1261 Scalar val = 0;
1262 for (unsigned i = 0; i < areaNormal.size(); ++i)
1263 val += areaNormal[i]*distance[i];
1264
1265 halfTrans *= std::abs(val);
1266 halfTrans /= distance.two_norm2();
1267}
1268
1269template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
1271computeHalfDiffusivity_(Scalar& halfDiff,
1272 const DimVector& areaNormal,
1273 const DimVector& distance,
1274 const Scalar& poro) const
1275{
1276 halfDiff = poro;
1277 Scalar val = 0;
1278 for (unsigned i = 0; i < areaNormal.size(); ++i)
1279 val += areaNormal[i]*distance[i];
1280
1281 halfDiff *= std::abs(val);
1282 halfDiff /= distance.two_norm2();
1283}
1284
1285template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
1288distanceVector_(const DimVector& center,
1289 int faceIdx, // in the reference element that contains the intersection
1290 unsigned elemIdx,
1291 const std::array<std::vector<DimVector>, dimWorld>& axisCentroids) const
1292{
1293 assert(faceIdx >= 0);
1294 unsigned dimIdx = faceIdx/2;
1295 assert(dimIdx < dimWorld);
1296 DimVector x = center;
1297 x -= axisCentroids[dimIdx][elemIdx];
1298
1299 return x;
1300}
1301
1302template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
1304applyMultipliers_(Scalar& trans,
1305 unsigned faceIdx,
1306 unsigned cartElemIdx,
1307 const TransMult& transMult) const
1308{
1309 // apply multiplyer for the transmissibility of the face. (the
1310 // face index is the index of the reference-element face which
1311 // contains the intersection of interest.)
1312 switch (faceIdx) {
1313 case 0: // left
1314 trans *= transMult.getMultiplier(cartElemIdx, FaceDir::XMinus);
1315 break;
1316 case 1: // right
1317 trans *= transMult.getMultiplier(cartElemIdx, FaceDir::XPlus);
1318 break;
1319
1320 case 2: // front
1321 trans *= transMult.getMultiplier(cartElemIdx, FaceDir::YMinus);
1322 break;
1323 case 3: // back
1324 trans *= transMult.getMultiplier(cartElemIdx, FaceDir::YPlus);
1325 break;
1326
1327 case 4: // bottom
1328 trans *= transMult.getMultiplier(cartElemIdx, FaceDir::ZMinus);
1329 break;
1330 case 5: // top
1331 trans *= transMult.getMultiplier(cartElemIdx, FaceDir::ZPlus);
1332 break;
1333 }
1334}
1335
1336template<class Grid, class GridView, class ElementMapper, class CartesianIndexMapper, class Scalar>
1338applyNtg_(Scalar& trans,
1339 unsigned faceIdx,
1340 unsigned elemIdx,
1341 const std::vector<double>& ntg) const
1342{
1343 // apply multiplyer for the transmissibility of the face. (the
1344 // face index is the index of the reference-element face which
1345 // contains the intersection of interest.)
1346 switch (faceIdx) {
1347 case 0: // left
1348 trans *= ntg[elemIdx];
1349 break;
1350 case 1: // right
1351 trans *= ntg[elemIdx];
1352 break;
1353
1354 case 2: // front
1355 trans *= ntg[elemIdx];
1356 break;
1357 case 3: // back
1358 trans *= ntg[elemIdx];
1359 break;
1360
1361 // NTG does not apply to top and bottom faces
1362 }
1363}
1364
1365} // namespace Opm
1366
1367#endif // OPM_TRANSMISSIBILITY_IMPL_HPP
void applyMultipliers_(Scalar &trans, unsigned faceIdx, unsigned cartElemIdx, const TransMult &transMult) const
Definition: Transmissibility_impl.hpp:1304
Scalar diffusivity(unsigned elemIdx1, unsigned elemIdx2) const
Return the diffusivity for the intersection between two elements.
Definition: Transmissibility_impl.hpp:142
void computeFaceProperties(const Intersection &intersection, const int, const int, const int, const int, DimVector &faceCenterInside, DimVector &faceCenterOutside, DimVector &faceAreaNormal, std::false_type) const
Definition: Transmissibility_impl.hpp:944
void applyEditNncToGridTrans_(const std::unordered_map< std::size_t, int > &globalToLocal)
Multiplies the grid transmissibilities according to EDITNNC.
Definition: Transmissibility_impl.hpp:1100
Dune::FieldMatrix< Scalar, dimWorld, dimWorld > DimMatrix
Definition: Transmissibility.hpp:59
void applyNncMultreg_(const std::unordered_map< std::size_t, int > &globalToLocal)
Definition: Transmissibility_impl.hpp:1202
void applyEditNncrToGridTrans_(const std::unordered_map< std::size_t, int > &globalToLocal)
Resets the grid transmissibilities according to EDITNNCR.
Definition: Transmissibility_impl.hpp:1113
void extractDispersion_()
Definition: Transmissibility_impl.hpp:691
Scalar thermalHalfTrans(unsigned insideElemIdx, unsigned outsideElemIdx) const
Return the thermal "half transmissibility" for the intersection between two elements.
Definition: Transmissibility_impl.hpp:128
void computeHalfDiffusivity_(Scalar &halfDiff, const DimVector &areaNormal, const DimVector &distance, const Scalar &poro) const
Definition: Transmissibility_impl.hpp:1271
void applyAllZMultipliers_(Scalar &trans, unsigned insideFaceIdx, unsigned outsideFaceIdx, unsigned insideCartElemIdx, unsigned outsideCartElemIdx, const TransMult &transMult, const std::array< int, dimWorld > &cartDims, bool pinchTop)
Apply the Multipliers for the case PINCH(4)==TOPBOT.
Definition: Transmissibility_impl.hpp:727
Scalar transmissibilityBoundary(unsigned elemIdx, unsigned boundaryFaceIdx) const
Return the transmissibility for a given boundary segment.
Definition: Transmissibility_impl.hpp:121
void computeHalfTrans_(Scalar &halfTrans, const DimVector &areaNormal, int faceIdx, const DimVector &distance, const DimMatrix &perm) const
Definition: Transmissibility_impl.hpp:1250
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:815
void updateFromEclState_(bool global)
Definition: Transmissibility_impl.hpp:778
void update(bool global, const std::function< unsigned int(unsigned int)> &map={}, bool applyNncMultRegT=false)
Compute all transmissibilities.
Definition: Transmissibility_impl.hpp:162
void extractPermeability_()
Definition: Transmissibility_impl.hpp:589
void resetTransmissibilityFromArrays_(const std::array< bool, 3 > &is_tran, const std::array< std::vector< double >, 3 > &trans)
overwrites calculated transmissibilities
Definition: Transmissibility_impl.hpp:882
Dune::FieldVector< Scalar, dimWorld > DimVector
Definition: Transmissibility.hpp:60
Scalar transmissibilityThreshold_
Definition: Transmissibility.hpp:275
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(double &, const double &)> &apply)
Definition: Transmissibility_impl.hpp:1126
Scalar transmissibility(unsigned elemIdx1, unsigned elemIdx2) const
Return the transmissibility for the intersection between two elements.
Definition: Transmissibility_impl.hpp:114
const EclipseState & eclState_
Definition: Transmissibility.hpp:270
Scalar dispersivity(unsigned elemIdx1, unsigned elemIdx2) const
Return the dispersivity for the intersection between two elements.
Definition: Transmissibility_impl.hpp:152
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:89
void applyNncToGridTrans_(const std::unordered_map< std::size_t, int > &cartesianToCompressed)
Definition: Transmissibility_impl.hpp:1032
void applyNtg_(Scalar &trans, unsigned faceIdx, unsigned elemIdx, const std::vector< double > &ntg) const
Definition: Transmissibility_impl.hpp:1338
DimVector distanceVector_(const DimVector &center, int faceIdx, unsigned elemIdx, const std::array< std::vector< DimVector >, dimWorld > &axisCentroids) const
Definition: Transmissibility_impl.hpp:1288
void extractPorosity_()
Definition: Transmissibility_impl.hpp:674
void removeNonCartesianTransmissibilities_(bool removeAll)
Definition: Transmissibility_impl.hpp:703
Scalar thermalHalfTransBoundary(unsigned insideElemIdx, unsigned boundaryFaceIdx) const
Definition: Transmissibility_impl.hpp:135
constexpr unsigned elemIdxShift
Definition: Transmissibility_impl.hpp:60
std::pair< std::uint32_t, std::uint32_t > isIdReverse(const std::uint64_t &id)
Definition: Transmissibility_impl.hpp:70
std::uint64_t directionalIsId(std::uint32_t elemIdx1, std::uint32_t elemIdx2)
Definition: Transmissibility_impl.hpp:81
std::uint64_t isId(std::uint32_t elemIdx1, std::uint32_t elemIdx2)
Definition: Transmissibility_impl.hpp:62
Definition: BlackoilPhases.hpp:27