20#ifndef OPM_GET_QUASI_IMPES_WEIGHTS_HEADER_INCLUDED
21#define OPM_GET_QUASI_IMPES_WEIGHTS_HEADER_INCLUDED
23#include <dune/common/fvector.hh>
25#include <opm/grid/utility/ElementChunks.hpp>
26#include <opm/simulators/utils/DeferredLoggingErrorHelpers.hpp>
27#include <opm/material/common/MathToolbox.hpp>
34#include <opm/simulators/linalg/gpuistl_hip/detail/cpr_amg_operations.hpp>
36#include <opm/simulators/linalg/gpuistl/detail/cpr_amg_operations.hpp>
46 template <
class DenseMatrix>
47 DenseMatrix transposeDenseMatrix(
const DenseMatrix& M)
50 for (
int i = 0; i < M.rows; ++i)
51 for (
int j = 0; j < M.cols; ++j)
60 template <
class Matrix,
class Vector>
61 void getQuasiImpesWeights(
const Matrix& matrix,
const int pressureVarIndex,
const bool transpose, Vector& weights,
bool enable_thread_parallel)
63 using VectorBlockType =
typename Vector::block_type;
64 using MatrixBlockType =
typename Matrix::block_type;
65 const Matrix& A = matrix;
67 VectorBlockType rhs(0.0);
68 rhs[pressureVarIndex] = 1.0;
71 MatrixBlockType diag_block;
72 VectorBlockType bweights;
73 MatrixBlockType diag_block_transpose;
77#pragma omp parallel for private(diag_block, bweights, diag_block_transpose) if(enable_thread_parallel)
79 for (
int row_idx = 0; row_idx < static_cast<int>(A.N()); ++row_idx) {
80 diag_block = MatrixBlockType(0.0);
82 const auto row_it = A.begin() + row_idx;
83 const auto endj = (*row_it).end();
84 for (
auto j = (*row_it).begin(); j != endj; ++j) {
85 if (row_it.index() == j.index()) {
91 diag_block.solve(bweights, rhs);
93 diag_block_transpose = Details::transposeDenseMatrix(diag_block);
94 diag_block_transpose.solve(bweights, rhs);
97 double abs_max = std::fabs(*std::max_element(
98 bweights.begin(), bweights.end(),
99 [](
double a,
double b) { return std::fabs(a) < std::fabs(b); }));
103 weights[row_idx] = bweights;
107 template <
class Matrix,
class Vector>
108 Vector getQuasiImpesWeights(
const Matrix& matrix,
const int pressureVarIndex,
const bool transpose,
bool enable_thread_parallel)
110 Vector weights(matrix.N());
111 getQuasiImpesWeights(matrix, pressureVarIndex, transpose, weights, enable_thread_parallel);
116 template <
typename T>
117 std::vector<int> precomputeDiagonalIndices(
const gpuistl::GpuSparseMatrixWrapper<T>& matrix) {
118 std::vector<int> diagonalIndices(matrix.N(), -1);
119 const auto rowIndices = matrix.getRowIndices().asStdVector();
120 const auto colIndices = matrix.getColumnIndices().asStdVector();
123 for (
auto i = rowIndices[row]; i < rowIndices[row+1]; ++i) {
124 if (colIndices[i] == row) {
125 diagonalIndices[row] = i;
130 return diagonalIndices;
134 template <
typename T,
bool transpose>
135 void getQuasiImpesWeights(
const gpuistl::GpuSparseMatrixWrapper<T>& matrix,
136 const int pressureVarIndex,
137 gpuistl::GpuVector<T>& weights,
138 const gpuistl::GpuVector<int>& diagonalIndices)
143 template <
typename T,
bool transpose>
144 gpuistl::GpuVector<T> getQuasiImpesWeights(
const gpuistl::GpuSparseMatrixWrapper<T>& matrix,
145 const int pressureVarIndex,
146 const gpuistl::GpuVector<int>& diagonalIndices)
148 gpuistl::GpuVector<T> weights(matrix.N() * matrix.blockSize());
149 getQuasiImpesWeights<T, transpose>(matrix, pressureVarIndex, weights, diagonalIndices);
154 template<
class Vector,
class ElementContext,
class Model,
class ElementChunksType>
155 void getTrueImpesWeights(
int pressureVarIndex, Vector& weights,
156 const ElementContext& elemCtx,
158 const ElementChunksType& element_chunks,
159 bool enable_thread_parallel)
161 using VectorBlockType =
typename Vector::block_type;
162 using Matrix =
typename std::decay_t<
decltype(model.linearizer().jacobian())>;
163 using MatrixBlockType =
typename Matrix::MatrixBlock;
164 constexpr int numEq = VectorBlockType::size();
165 using Evaluation =
typename std::decay_t<
decltype(model.localLinearizer(
ThreadManager::threadId()).localResidual().residual(0))>
168 VectorBlockType rhs(0.0);
169 rhs[pressureVarIndex] = 1.0;
172 MatrixBlockType block;
173 VectorBlockType bweights;
174 MatrixBlockType block_transpose;
175 Dune::FieldVector<Evaluation, numEq> storage;
177 OPM_BEGIN_PARALLEL_TRY_CATCH();
179#pragma omp parallel for private(block, bweights, block_transpose, storage) if(enable_thread_parallel)
181 for (
const auto& chunk : element_chunks) {
183 ElementContext localElemCtx(elemCtx.simulator());
185 for (
const auto& elem : chunk) {
186 localElemCtx.updatePrimaryStencil(elem);
187 localElemCtx.updatePrimaryIntensiveQuantities(0);
189 model.localLinearizer(thread_id).localResidual().computeStorage(storage, localElemCtx, 0, 0);
191 auto extrusionFactor = localElemCtx.intensiveQuantities(0, 0).extrusionFactor();
192 auto scvVolume = localElemCtx.stencil(0).subControlVolume(0).volume() * extrusionFactor;
193 auto storage_scale = scvVolume / localElemCtx.simulator().timeStepSize();
194 const double pressure_scale = 50e5;
197 for (
int ii = 0; ii < numEq; ++ii) {
198 for (
int jj = 0; jj < numEq; ++jj) {
199 block_transpose[jj][ii] = storage[ii].derivative(jj)/storage_scale;
200 if (jj == pressureVarIndex) {
201 block_transpose[jj][ii] *= pressure_scale;
205 block_transpose.solve(bweights, rhs);
207 double abs_max = std::fabs(*std::max_element(
208 bweights.begin(), bweights.end(),
209 [](
double a,
double b) { return std::fabs(a) < std::fabs(b); }));
216 const auto index = localElemCtx.globalSpaceIndex(0, 0);
217 weights[index] = bweights;
220 OPM_END_PARALLEL_TRY_CATCH(
"getTrueImpesWeights() failed: ", elemCtx.simulator().vanguard().grid().comm());
223 template <
class Vector,
class ElementContext,
class Model,
class ElementChunksType>
224 void getTrueImpesWeightsAnalytic(
int ,
226 const ElementContext& elemCtx,
228 const ElementChunksType& element_chunks,
229 bool enable_thread_parallel)
239 using FluidSystem =
typename Model::FluidSystem;
240 using LhsEval = double;
242 using PrimaryVariables =
typename Model::PrimaryVariables;
243 using VectorBlockType =
typename Vector::block_type;
245 typename std::decay_t<
decltype(model.localLinearizer(
ThreadManager::threadId()).localResidual().residual(0))>::block_type;
246 using Toolbox = MathToolbox<Evaluation>;
248 const auto& solution = model.solution( 0);
249 VectorBlockType bweights;
252 OPM_BEGIN_PARALLEL_TRY_CATCH();
254#pragma omp parallel for private(bweights) if(enable_thread_parallel)
256 for (
const auto& chunk : element_chunks) {
259 ElementContext localElemCtx(elemCtx.simulator());
261 for (
const auto& elem : chunk) {
262 localElemCtx.updatePrimaryStencil(elem);
263 localElemCtx.updatePrimaryIntensiveQuantities(0);
265 const auto index = localElemCtx.globalSpaceIndex(0, 0);
266 const auto& intQuants = localElemCtx.intensiveQuantities(0, 0);
267 const auto& fs = intQuants.fluidState();
269 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
270 const unsigned activeCompIdx = FluidSystem::canonicalToActiveCompIdx(
271 FluidSystem::solventComponentIndex(FluidSystem::waterPhaseIdx));
272 bweights[activeCompIdx]
273 = Toolbox::template decay<LhsEval>(1 / fs.invB(FluidSystem::waterPhaseIdx));
276 double denominator = 1.0;
277 double rs = Toolbox::template decay<double>(fs.Rs());
278 double rv = Toolbox::template decay<double>(fs.Rv());
279 const auto& priVars = solution[index];
280 if (priVars.primaryVarsMeaningGas() == PrimaryVariables::GasMeaning::Rv) {
283 if (priVars.primaryVarsMeaningGas() == PrimaryVariables::GasMeaning::Rs) {
286 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)
287 && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
288 denominator = Toolbox::template decay<LhsEval>(1 - rs * rv);
291 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
292 const unsigned activeCompIdx = FluidSystem::canonicalToActiveCompIdx(
293 FluidSystem::solventComponentIndex(FluidSystem::oilPhaseIdx));
294 bweights[activeCompIdx] = Toolbox::template decay<LhsEval>(
295 (1 / fs.invB(FluidSystem::oilPhaseIdx) - rs / fs.invB(FluidSystem::gasPhaseIdx))
298 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
299 const unsigned activeCompIdx = FluidSystem::canonicalToActiveCompIdx(
300 FluidSystem::solventComponentIndex(FluidSystem::gasPhaseIdx));
301 bweights[activeCompIdx] = Toolbox::template decay<LhsEval>(
302 (1 / fs.invB(FluidSystem::gasPhaseIdx) - rv / fs.invB(FluidSystem::oilPhaseIdx))
306 weights[index] = bweights;
309 OPM_END_PARALLEL_TRY_CATCH(
"getTrueImpesAnalyticWeights() failed: ", elemCtx.simulator().vanguard().grid().comm());
static unsigned threadId()
Return the index of the current OpenMP thread.
Definition threadmanager.cpp:84
int to_int(std::size_t s)
to_int converts a (on most relevant platforms) 64 bits unsigned size_t to a signed 32 bits signed int
Definition safe_conversion.hpp:52
void getQuasiImpesWeights(const GpuSparseMatrixWrapper< T > &matrix, std::size_t pressureVarIndex, GpuVector< T > &weights, const GpuVector< int > &diagonalIndices)
Calculates quasi-IMPES weights for CPR preconditioner on GPU.
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:43
Simplifies multi-threaded capabilities.