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>
50 for (
int i = 0; i < M.rows; ++i)
51 for (
int j = 0; j < M.cols; ++j)
60 template <
class Matrix,
class Vector>
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;
77#pragma omp parallel for private(diag_block, bweights, diag_block_transpose) if(enable_thread_parallel)
83 const auto endj = (*row_it).end();
84 for (
auto j = (*row_it).begin(); j !=
endj; ++j) {
85 if (
row_it.index() == j.index()) {
97 double abs_max = *std::max_element(
98 bweights.begin(),
bweights.end(), [](
double a,
double b) { return std::fabs(a) < std::fabs(b); });
104 template <
class Matrix,
class Vector>
107 Vector weights(matrix.N());
113 template <
typename T>
116 const auto rowIndices = matrix.getRowIndices().asStdVector();
117 const auto colIndices = matrix.getColumnIndices().asStdVector();
120 for (
auto i = rowIndices[
row]; i < rowIndices[
row+1]; ++i) {
121 if (colIndices[i] ==
row) {
131 template <
typename T,
bool transpose>
132 void getQuasiImpesWeights(
const gpuistl::GpuSparseMatrixWrapper<T>& matrix,
133 const int pressureVarIndex,
134 gpuistl::GpuVector<T>& weights,
137 gpuistl::detail::getQuasiImpesWeights<T, transpose>(matrix, pressureVarIndex, weights,
diagonalIndices);
140 template <
typename T,
bool transpose>
141 gpuistl::GpuVector<T> getQuasiImpesWeights(
const gpuistl::GpuSparseMatrixWrapper<T>& matrix,
142 const int pressureVarIndex,
145 gpuistl::GpuVector<T> weights(matrix.N() * matrix.blockSize());
151 template<
class Vector,
class ElementContext,
class Model,
class ElementChunksType>
152 void getTrueImpesWeights(
int pressureVarIndex, Vector& weights,
153 const ElementContext& elemCtx,
158 using VectorBlockType =
typename Vector::block_type;
159 using Matrix =
typename std::decay_t<
decltype(model.linearizer().jacobian())>;
160 using MatrixBlockType =
typename Matrix::MatrixBlock;
161 constexpr int numEq = VectorBlockType::size();
162 using Evaluation =
typename std::decay_t<
decltype(model.localLinearizer(
ThreadManager::threadId()).localResidual().residual(0))>
165 VectorBlockType rhs(0.0);
166 rhs[pressureVarIndex] = 1.0;
169 MatrixBlockType block;
172 Dune::FieldVector<Evaluation, numEq>
storage;
174 OPM_BEGIN_PARALLEL_TRY_CATCH();
176#pragma omp parallel for private(block, bweights, block_transpose, storage) if(enable_thread_parallel)
188 auto extrusionFactor =
localElemCtx.intensiveQuantities(0, 0).extrusionFactor();
194 for (
int ii = 0;
ii < numEq; ++
ii) {
195 for (
int jj = 0;
jj < numEq; ++
jj) {
197 if (
jj == pressureVarIndex) {
204 double abs_max = *std::max_element(
205 bweights.begin(),
bweights.end(), [](
double a,
double b) { return std::fabs(a) < std::fabs(b); });
213 OPM_END_PARALLEL_TRY_CATCH(
"getTrueImpesWeights() failed: ", elemCtx.simulator().vanguard().grid().comm());
216 template <
class Vector,
class ElementContext,
class Model,
class ElementChunksType>
217 void getTrueImpesWeightsAnalytic(
int ,
219 const ElementContext& elemCtx,
232 using FluidSystem =
typename Model::FluidSystem;
235 using PrimaryVariables =
typename Model::PrimaryVariables;
236 using VectorBlockType =
typename Vector::block_type;
238 typename std::decay_t<
decltype(model.localLinearizer(
ThreadManager::threadId()).localResidual().residual(0))>::block_type;
241 const auto& solution = model.solution( 0);
245 OPM_BEGIN_PARALLEL_TRY_CATCH();
247#pragma omp parallel for private(bweights) if(enable_thread_parallel)
259 const auto& intQuants =
localElemCtx.intensiveQuantities(0, 0);
260 const auto& fs = intQuants.fluidState();
262 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
263 const unsigned activeCompIdx = FluidSystem::canonicalToActiveCompIdx(
264 FluidSystem::solventComponentIndex(FluidSystem::waterPhaseIdx));
266 = Toolbox::template
decay<LhsEval>(1 / fs.invB(FluidSystem::waterPhaseIdx));
272 const auto& priVars = solution[index];
273 if (priVars.primaryVarsMeaningGas() == PrimaryVariables::GasMeaning::Rv) {
276 if (priVars.primaryVarsMeaningGas() == PrimaryVariables::GasMeaning::Rs) {
279 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)
280 && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
284 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
285 const unsigned activeCompIdx = FluidSystem::canonicalToActiveCompIdx(
286 FluidSystem::solventComponentIndex(FluidSystem::oilPhaseIdx));
288 (1 / fs.invB(FluidSystem::oilPhaseIdx) - rs / fs.invB(FluidSystem::gasPhaseIdx))
291 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
292 const unsigned activeCompIdx = FluidSystem::canonicalToActiveCompIdx(
293 FluidSystem::solventComponentIndex(FluidSystem::gasPhaseIdx));
295 (1 / fs.invB(FluidSystem::gasPhaseIdx) - rv / fs.invB(FluidSystem::oilPhaseIdx))
302 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
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:43
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240
Simplifies multi-threaded capabilities.