opm-simulators
Loading...
Searching...
No Matches
getQuasiImpesWeights.hpp
1/*
2 Copyright 2019 SINTEF Digital, Mathematics and Cybernetics.
3
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 3 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18*/
19
20#ifndef OPM_GET_QUASI_IMPES_WEIGHTS_HEADER_INCLUDED
21#define OPM_GET_QUASI_IMPES_WEIGHTS_HEADER_INCLUDED
22
23#include <dune/common/fvector.hh>
24
25#include <opm/grid/utility/ElementChunks.hpp>
26#include <opm/simulators/utils/DeferredLoggingErrorHelpers.hpp>
27#include <opm/material/common/MathToolbox.hpp>
29#include <algorithm>
30#include <cmath>
31
32#if HAVE_CUDA
33#if USE_HIP
34#include <opm/simulators/linalg/gpuistl_hip/detail/cpr_amg_operations.hpp>
35#else
36#include <opm/simulators/linalg/gpuistl/detail/cpr_amg_operations.hpp>
37#endif
38#endif
39
40
41namespace Opm
42{
43
44namespace Details
45{
46 template <class DenseMatrix>
47 DenseMatrix transposeDenseMatrix(const DenseMatrix& M)
48 {
49 DenseMatrix tmp;
50 for (int i = 0; i < M.rows; ++i)
51 for (int j = 0; j < M.cols; ++j)
52 tmp[j][i] = M[i][j];
53
54 return tmp;
55 }
56} // namespace Details
57
58namespace Amg
59{
60 template <class Matrix, class Vector>
61 void getQuasiImpesWeights(const Matrix& matrix, const int pressureVarIndex, const bool transpose, Vector& weights)
62 {
63 using VectorBlockType = typename Vector::block_type;
64 using MatrixBlockType = typename Matrix::block_type;
65 const Matrix& A = matrix;
66
67 VectorBlockType rhs(0.0);
68 rhs[pressureVarIndex] = 1.0;
69
70 // Declare variables outside the loop to avoid repetitive allocation
71 MatrixBlockType diag_block;
72 VectorBlockType bweights;
73 MatrixBlockType diag_block_transpose;
74
75 // Use OpenMP to parallelize over matrix rows
76 #ifdef _OPENMP
77 #pragma omp parallel for private(diag_block, bweights, diag_block_transpose)
78 #endif
79 for (int row_idx = 0; row_idx < static_cast<int>(A.N()); ++row_idx) {
80 diag_block = MatrixBlockType(0.0);
81 // Find diagonal block for this row
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()) {
86 diag_block = (*j);
87 break;
88 }
89 }
90 if (transpose) {
91 diag_block.solve(bweights, rhs);
92 } else {
93 diag_block_transpose = Details::transposeDenseMatrix(diag_block);
95 }
96
97 double abs_max = *std::max_element(
98 bweights.begin(), bweights.end(), [](double a, double b) { return std::fabs(a) < std::fabs(b); });
99 bweights /= std::fabs(abs_max);
100 weights[row_idx] = bweights;
101 }
102 }
103
104 template <class Matrix, class Vector>
105 Vector getQuasiImpesWeights(const Matrix& matrix, const int pressureVarIndex, const bool transpose)
106 {
107 Vector weights(matrix.N());
108 getQuasiImpesWeights(matrix, pressureVarIndex, transpose, weights);
109 return weights;
110 }
111
112#if HAVE_CUDA
113 template <typename T>
114 std::vector<int> precomputeDiagonalIndices(const gpuistl::GpuSparseMatrixWrapper<T>& matrix) {
115 std::vector<int> diagonalIndices(matrix.N(), -1);
116 const auto rowIndices = matrix.getRowIndices().asStdVector();
117 const auto colIndices = matrix.getColumnIndices().asStdVector();
118
119 for (auto row = 0; row < Opm::gpuistl::detail::to_int(matrix.N()); ++row) {
120 for (auto i = rowIndices[row]; i < rowIndices[row+1]; ++i) {
121 if (colIndices[i] == row) {
122 diagonalIndices[row] = i;
123 break;
124 }
125 }
126 }
127 return diagonalIndices;
128 }
129
130 // GPU version that delegates to the GPU implementation
131 template <typename T, bool transpose>
132 void getQuasiImpesWeights(const gpuistl::GpuSparseMatrixWrapper<T>& matrix,
133 const int pressureVarIndex,
134 gpuistl::GpuVector<T>& weights,
135 const gpuistl::GpuVector<int>& diagonalIndices)
136 {
137 gpuistl::detail::getQuasiImpesWeights<T, transpose>(matrix, pressureVarIndex, weights, diagonalIndices);
138 }
139
140 template <typename T, bool transpose>
141 gpuistl::GpuVector<T> getQuasiImpesWeights(const gpuistl::GpuSparseMatrixWrapper<T>& matrix,
142 const int pressureVarIndex,
143 const gpuistl::GpuVector<int>& diagonalIndices)
144 {
145 gpuistl::GpuVector<T> weights(matrix.N() * matrix.blockSize());
146 getQuasiImpesWeights<T, transpose>(matrix, pressureVarIndex, weights, diagonalIndices);
147 return weights;
148 }
149#endif
150
151 template<class Vector, class ElementContext, class Model, class ElementChunksType>
152 void getTrueImpesWeights(int pressureVarIndex, Vector& weights,
153 const ElementContext& elemCtx,
154 const Model& model,
155 const ElementChunksType& element_chunks)
156 {
157 using VectorBlockType = typename Vector::block_type;
158 using Matrix = typename std::decay_t<decltype(model.linearizer().jacobian())>;
159 using MatrixBlockType = typename Matrix::MatrixBlock;
160 constexpr int numEq = VectorBlockType::size();
161 using Evaluation = typename std::decay_t<decltype(model.localLinearizer(ThreadManager::threadId()).localResidual().residual(0))>
162 ::block_type;
163
164 VectorBlockType rhs(0.0);
165 rhs[pressureVarIndex] = 1.0;
166
167 // Declare variables outside the loop to avoid repetitive allocation
168 MatrixBlockType block;
169 VectorBlockType bweights;
170 MatrixBlockType block_transpose;
171 Dune::FieldVector<Evaluation, numEq> storage;
172
173 OPM_BEGIN_PARALLEL_TRY_CATCH();
174 #ifdef _OPENMP
175 #pragma omp parallel for private(block, bweights, block_transpose, storage)
176 #endif
177 for (const auto& chunk : element_chunks) {
178 const std::size_t thread_id = ThreadManager::threadId();
179 ElementContext localElemCtx(elemCtx.simulator());
180
181 for (const auto& elem : chunk) {
182 localElemCtx.updatePrimaryStencil(elem);
183 localElemCtx.updatePrimaryIntensiveQuantities(/*timeIdx=*/0);
184
185 model.localLinearizer(thread_id).localResidual().computeStorage(storage, localElemCtx, /*spaceIdx=*/0, /*timeIdx=*/0);
186
187 auto extrusionFactor = localElemCtx.intensiveQuantities(0, /*timeIdx=*/0).extrusionFactor();
188 auto scvVolume = localElemCtx.stencil(/*timeIdx=*/0).subControlVolume(0).volume() * extrusionFactor;
189 auto storage_scale = scvVolume / localElemCtx.simulator().timeStepSize();
190 const double pressure_scale = 50e5;
191
192 // Build the transposed matrix directly to avoid separate transpose step
193 for (int ii = 0; ii < numEq; ++ii) {
194 for (int jj = 0; jj < numEq; ++jj) {
195 block_transpose[jj][ii] = storage[ii].derivative(jj)/storage_scale;
196 if (jj == pressureVarIndex) {
198 }
199 }
200 }
201 block_transpose.solve(bweights, rhs);
202
203 double abs_max = *std::max_element(
204 bweights.begin(), bweights.end(), [](double a, double b) { return std::fabs(a) < std::fabs(b); });
205 // probably a scaling which could give approximately total compressibility would be better
206 bweights /= std::fabs(abs_max); // given normal densities this scales weights to about 1.
207
208 const auto index = localElemCtx.globalSpaceIndex(/*spaceIdx=*/0, /*timeIdx=*/0);
209 weights[index] = bweights;
210 }
211 }
212 OPM_END_PARALLEL_TRY_CATCH("getTrueImpesWeights() failed: ", elemCtx.simulator().vanguard().grid().comm());
213 }
214
215 template <class Vector, class ElementContext, class Model, class ElementChunksType>
216 void getTrueImpesWeightsAnalytic(int /*pressureVarIndex*/,
217 Vector& weights,
218 const ElementContext& elemCtx,
219 const Model& model,
220 const ElementChunksType& element_chunks)
221 {
222 // The sequential residual is a linear combination of the
223 // mass balance residuals, with coefficients equal to (for
224 // water, oil, gas):
225 // 1/bw,
226 // (1/bo - rs/bg)/(1-rs*rv)
227 // (1/bg - rv/bo)/(1-rs*rv)
228 // These coefficients must be applied for both the residual and
229 // Jacobian.
230 using FluidSystem = typename Model::FluidSystem;
231 using LhsEval = double;
232
233 using PrimaryVariables = typename Model::PrimaryVariables;
234 using VectorBlockType = typename Vector::block_type;
235 using Evaluation =
236 typename std::decay_t<decltype(model.localLinearizer(ThreadManager::threadId()).localResidual().residual(0))>::block_type;
237 using Toolbox = MathToolbox<Evaluation>;
238
239 const auto& solution = model.solution(/*timeIdx*/ 0);
240 VectorBlockType bweights;
241
242 // Use OpenMP to parallelize over element chunks
243 OPM_BEGIN_PARALLEL_TRY_CATCH();
244 #ifdef _OPENMP
245 #pragma omp parallel for private(bweights)
246 #endif
247 for (const auto& chunk : element_chunks) {
248
249 // Each thread gets a unique copy of elemCtx
250 ElementContext localElemCtx(elemCtx.simulator());
251
252 for (const auto& elem : chunk) {
253 localElemCtx.updatePrimaryStencil(elem);
254 localElemCtx.updatePrimaryIntensiveQuantities(/*timeIdx=*/0);
255
256 const auto index = localElemCtx.globalSpaceIndex(/*spaceIdx=*/0, /*timeIdx=*/0);
257 const auto& intQuants = localElemCtx.intensiveQuantities(/*spaceIdx=*/0, /*timeIdx=*/0);
258 const auto& fs = intQuants.fluidState();
259
260 if (FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx)) {
261 const unsigned activeCompIdx = FluidSystem::canonicalToActiveCompIdx(
262 FluidSystem::solventComponentIndex(FluidSystem::waterPhaseIdx));
264 = Toolbox::template decay<LhsEval>(1 / fs.invB(FluidSystem::waterPhaseIdx));
265 }
266
267 double denominator = 1.0;
268 double rs = Toolbox::template decay<double>(fs.Rs());
269 double rv = Toolbox::template decay<double>(fs.Rv());
270 const auto& priVars = solution[index];
271 if (priVars.primaryVarsMeaningGas() == PrimaryVariables::GasMeaning::Rv) {
272 rs = 0.0;
273 }
274 if (priVars.primaryVarsMeaningGas() == PrimaryVariables::GasMeaning::Rs) {
275 rv = 0.0;
276 }
277 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)
278 && FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
279 denominator = Toolbox::template decay<LhsEval>(1 - rs * rv);
280 }
281
282 if (FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx)) {
283 const unsigned activeCompIdx = FluidSystem::canonicalToActiveCompIdx(
284 FluidSystem::solventComponentIndex(FluidSystem::oilPhaseIdx));
285 bweights[activeCompIdx] = Toolbox::template decay<LhsEval>(
286 (1 / fs.invB(FluidSystem::oilPhaseIdx) - rs / fs.invB(FluidSystem::gasPhaseIdx))
287 / denominator);
288 }
289 if (FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx)) {
290 const unsigned activeCompIdx = FluidSystem::canonicalToActiveCompIdx(
291 FluidSystem::solventComponentIndex(FluidSystem::gasPhaseIdx));
292 bweights[activeCompIdx] = Toolbox::template decay<LhsEval>(
293 (1 / fs.invB(FluidSystem::gasPhaseIdx) - rv / fs.invB(FluidSystem::oilPhaseIdx))
294 / denominator);
295 }
296
297 weights[index] = bweights;
298 }
299 }
300 OPM_END_PARALLEL_TRY_CATCH("getTrueImpesAnalyticWeights() failed: ", elemCtx.simulator().vanguard().grid().comm());
301 }
302} // namespace Amg
303
304} // namespace Opm
305
306#endif // OPM_GET_QUASI_IMPES_WEIGHTS_HEADER_INCLUDED
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.