mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
Added computeSurfacevol() method.
The method is called by the reordering transport solver after computing new saturations in order to update the surface volumes.
This commit is contained in:
@@ -24,6 +24,7 @@
|
||||
#include <opm/core/transport/reorder/reordersequence.h>
|
||||
#include <opm/core/utility/RootFinders.hpp>
|
||||
#include <opm/core/utility/miscUtilities.hpp>
|
||||
#include <opm/core/utility/miscUtilitiesBlackoil.hpp>
|
||||
#include <opm/core/pressure/tpfa/trans_tpfa.h>
|
||||
|
||||
#include <fstream>
|
||||
@@ -112,15 +113,9 @@ namespace Opm
|
||||
&ia_downw_[0], &ja_downw_[0]);
|
||||
reorderAndTransport(grid_, darcyflux);
|
||||
toBothSat(saturation_, saturation);
|
||||
|
||||
|
||||
// Compute surface volume as a postprocessing step from saturation and A_
|
||||
surfacevol = saturation;
|
||||
const int np = props_.numPhases();
|
||||
for (int cell = 0; cell < grid_.number_of_cells; ++cell) {
|
||||
for (int phase = 0; phase < np; ++phase) {
|
||||
surfacevol[np*cell + phase] *= A_[np*np*cell + np*phase + phase];
|
||||
}
|
||||
}
|
||||
computeSurfacevol(grid_.number_of_cells, props_.numPhases(), &A_[0], &saturation[0], &surfacevol[0]);
|
||||
}
|
||||
|
||||
// Residual function r(s) for a single-cell implicit Euler transport
|
||||
|
||||
@@ -220,5 +220,35 @@ namespace Opm
|
||||
}
|
||||
}
|
||||
|
||||
/// Computes the surface volume densities from saturations by the formula
|
||||
/// z = A s
|
||||
/// for a number of data points, where z is the surface volume density,
|
||||
/// s is the saturation (both as column vectors) and A is the
|
||||
/// phase-to-component relation matrix.
|
||||
/// @param[in] n number of data points
|
||||
/// @param[in] np number of phases, must be 2 or 3
|
||||
/// @param[in] A array containing n square matrices of size num_phases^2,
|
||||
/// in Fortran ordering, typically the output of a call
|
||||
/// to the matrix() method of a BlackoilProperties* class.
|
||||
/// @param[in] saturation concatenated saturation values (for all P phases)
|
||||
/// @param[out] surfacevol concatenated surface-volume values (for all P phases)
|
||||
void computeSurfacevol(const int n,
|
||||
const int np,
|
||||
const double* A,
|
||||
const double* saturation,
|
||||
double* surfacevol)
|
||||
{
|
||||
// Note: since this is a simple matrix-vector product, it can
|
||||
// be done by a BLAS call, but then we have to reorder the A
|
||||
// matrix data.
|
||||
std::fill(surfacevol, surfacevol + n*np, 0.0);
|
||||
for (int i = 0; i < n; ++i) {
|
||||
for (int row = 0; row < np; ++row) {
|
||||
for (int col = 0; col < np; ++col) {
|
||||
surfacevol[i*np + row] += A[i*np*np + row*np + col] * saturation[i*np + col];
|
||||
}
|
||||
}
|
||||
}
|
||||
}
|
||||
|
||||
} // namespace Opm
|
||||
|
||||
@@ -112,6 +112,25 @@ namespace Opm
|
||||
const std::vector<double>& s,
|
||||
std::vector<double>& fractional_flows);
|
||||
|
||||
|
||||
/// Computes the surface volume densities from saturations by the formula
|
||||
/// z = A s
|
||||
/// for a number of data points, where z is the surface volume density,
|
||||
/// s is the saturation (both as column vectors) and A is the
|
||||
/// phase-to-component relation matrix.
|
||||
/// @param[in] n number of data points
|
||||
/// @param[in] np number of phases, must be 2 or 3
|
||||
/// @param[in] A array containing n square matrices of size num_phases,
|
||||
/// in Fortran ordering, typically the output of a call
|
||||
/// to the matrix() method of a BlackoilProperties* class.
|
||||
/// @param[in] saturation concatenated saturation values (for all P phases)
|
||||
/// @param[out] surfacevol concatenated surface-volume values (for all P phases)
|
||||
void computeSurfacevol(const int n,
|
||||
const int np,
|
||||
const double* A,
|
||||
const double* saturation,
|
||||
double* surfacevol);
|
||||
|
||||
} // namespace Opm
|
||||
|
||||
#endif // OPM_MISCUTILITIESBLACKOIL_HEADER_INCLUDED
|
||||
|
||||
Reference in New Issue
Block a user