GeoProps: fix the use_local_perm implementation for flow_cp. The problem that the

faceNormals are scaled differently remains.
This commit is contained in:
Robert Kloefkorn 2015-07-08 13:15:00 +02:00
parent 8cd014f937
commit bdb30ec023
2 changed files with 28 additions and 12 deletions

View File

@ -323,15 +323,7 @@ try
bool use_gravity = (gravity[0] != 0.0 || gravity[1] != 0.0 || gravity[2] != 0.0);
const double *grav = use_gravity ? &gravity[0] : 0;
#if USE_DUNE_CORNERPOINTGRID
if(output_cout)
{
std::cout << std::endl << "Warning: use of local perm is not yet implemented for CpGrid!" << std::endl << std::endl;
}
const bool use_local_perm = false;
#else
const bool use_local_perm = param.getDefault("use_local_perm", true);
#endif
DerivedGeology geoprops(grid, new_props, eclipseState, use_local_perm, grav);
boost::any parallel_information;

View File

@ -91,7 +91,7 @@ namespace Opm
}
// get grid from parser.
// Get original grid cell volume.
EclipseGridConstPtr eclgrid = eclState->getEclipseGrid();
// Pore volume.
@ -107,7 +107,7 @@ namespace Opm
pvol_[cellIdx] *= AutoDiffGrid::cellVolume(grid, cellIdx);
} else {
pvol_[cellIdx] *= eclgrid->getCellVolume(cartesianCellIdx);
}
}
}
// Use volume weighted arithmetic average of the NTG values for
// the cells effected by the current OPM cpgrid process algorithm
@ -366,17 +366,41 @@ namespace Opm
// d = 0: XPERM d = 4: YPERM d = 8: ZPERM ignores off-diagonal permeability values.
const int d = std::floor(faceTag/2) * 4;
//std::cout << d << " d value | faceTag " << faceTag << std::endl;
// compute the half transmissibility
double dist = 0.0;
double cn = 0.0;
double sgn = 2.0 * (faceCells(faceIdx, 0) == cellIdx) - 1;
const int dim = Opm::UgGridHelpers::dimensions(grid);
/*
std::cout << "faceNormal ";
Dune::DynamicVector< double > normal( dim );
for (int indx = 0; indx < dim; ++indx) {
const double Ci = Opm::UgGridHelpers::faceCentroid(grid, faceIdx)[indx] -
normal[ indx ] = Opm::UgGridHelpers::faceNormal(grid,faceIdx)[indx];
}
if( ! ( std::abs( normal.two_norm() - Opm::UgGridHelpers::faceArea(grid, faceIdx) ) < 1e-12 ) )
normal *=
std::abort();
std::cout << std::endl;
*/
const double* faceNormal = Opm::UgGridHelpers::faceNormal(grid, faceIdx);
assert( dim <= 3 );
Dune::FieldVector< double, 3 > scaledFaceNormal;
for (int indx = 0; indx < dim; ++indx) {
scaledFaceNormal[ indx ] = faceNormal[ indx ];
}
// compute unit normal incase the normal is already scaled
scaledFaceNormal /= scaledFaceNormal.two_norm();
// compute proper normal scaled with face area
scaledFaceNormal *= Opm::UgGridHelpers::faceArea(grid, faceIdx);
for (int indx = 0; indx < dim; ++indx) {
const double Ci = Opm::UgGridHelpers::faceCentroid(grid, faceIdx)[indx] -
Opm::UgGridHelpers::cellCentroidCoordinate(grid, cellIdx, indx);
dist += Ci*Ci;
cn += sgn * Ci * Opm::UgGridHelpers::faceNormal(grid, faceIdx)[indx];
cn += sgn * Ci * scaledFaceNormal[ indx ]; //Opm::UgGridHelpers::faceNormal(grid, faceIdx)[indx];
}
if (cn < 0){