mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
Adapt to changes in EclipseGrids
getCellCenter() in EclipseGrid no longer returns tuble<double,double,double> but array<double,3> This simplifies its usage. Variable name is changed from cellCentroid to cellCenter to emphasize a difference between the value in the unstructured grid and in EclispeGrid
This commit is contained in:
parent
2a13275943
commit
67eadc30e5
@ -392,23 +392,11 @@ namespace Opm
|
||||
#endif
|
||||
|
||||
int cartesianCellIdx = AutoDiffGrid::globalCell(grid)[cellIdx];
|
||||
auto cellCentroid = eclGrid->getCellCenter(cartesianCellIdx);
|
||||
auto cellCenter = eclGrid->getCellCenter(cartesianCellIdx);
|
||||
|
||||
for (int indx = 0; indx < dim; ++indx) {
|
||||
|
||||
double Ci = Opm::UgGridHelpers::faceCentroid(grid, faceIdx)[indx];
|
||||
switch (indx) {
|
||||
case 0:
|
||||
Ci -= std::get<0>(cellCentroid);
|
||||
break;
|
||||
case 1:
|
||||
Ci -= std::get<1>(cellCentroid);
|
||||
break;
|
||||
case 2:
|
||||
Ci -= std::get<2>(cellCentroid);
|
||||
break;
|
||||
|
||||
}
|
||||
const double Ci = Opm::UgGridHelpers::faceCentroid(grid, faceIdx)[indx] - cellCenter[indx];
|
||||
dist += Ci*Ci;
|
||||
cn += sgn * Ci * scaledFaceNormal[ indx ];
|
||||
}
|
||||
|
Loading…
Reference in New Issue
Block a user