blackoilmicpmodules.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
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 2 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 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
28#ifndef EWOMS_BLACK_OIL_MICP_MODULE_HH
29#define EWOMS_BLACK_OIL_MICP_MODULE_HH
30
31#include <dune/common/fvector.hh>
32
35
37
38#include <cstddef>
39#include <stdexcept>
40
41namespace Opm {
47template <class TypeTag, bool enableMICPV = getPropValue<TypeTag, Properties::EnableMICP>()>
49{
62
63 using Toolbox = MathToolbox<Evaluation>;
64
65 static constexpr unsigned microbialConcentrationIdx = Indices::microbialConcentrationIdx;
66 static constexpr unsigned oxygenConcentrationIdx = Indices::oxygenConcentrationIdx;
67 static constexpr unsigned ureaConcentrationIdx = Indices::ureaConcentrationIdx;
68 static constexpr unsigned biofilmConcentrationIdx = Indices::biofilmConcentrationIdx;
69 static constexpr unsigned calciteConcentrationIdx = Indices::calciteConcentrationIdx;
70 static constexpr unsigned contiMicrobialEqIdx = Indices::contiMicrobialEqIdx;
71 static constexpr unsigned contiOxygenEqIdx = Indices::contiOxygenEqIdx;
72 static constexpr unsigned contiUreaEqIdx = Indices::contiUreaEqIdx;
73 static constexpr unsigned contiBiofilmEqIdx = Indices::contiBiofilmEqIdx;
74 static constexpr unsigned contiCalciteEqIdx = Indices::contiCalciteEqIdx;
75 static constexpr unsigned waterPhaseIdx = FluidSystem::waterPhaseIdx;
76
77 static constexpr unsigned enableMICP = enableMICPV;
78
79 static constexpr unsigned numEq = getPropValue<TypeTag, Properties::NumEq>();
80
81public:
84 {
85 params_ = params;
86 }
87
95 static void checkCloggingMICP(const Model& model, const Scalar phi, unsigned dofIdx)
96 {
97 const PrimaryVariables& priVars = model.solution(/*timeIdx=*/1)[dofIdx];
98 if (phi - priVars[biofilmConcentrationIdx] - priVars[calciteConcentrationIdx] < toleranceBeforeClogging())
99 throw std::logic_error("Clogging has been (almost) reached in at least one cell\n");
100 }
101
105 static void registerParameters()
106 {
107 if (!enableMICP)
108 // MICP has been disabled at compile time
109 return;
110
112 }
113
117 static void registerOutputModules(Model& model,
118 Simulator& simulator)
119 {
120 if (!enableMICP)
121 // MICP has been disabled at compile time
122 return;
123
124 model.addOutputModule(new VtkBlackOilMICPModule<TypeTag>(simulator));
125 }
126
127 static bool eqApplies(unsigned eqIdx)
128 {
129 if (!enableMICP)
130 return false;
131
132 // All MICP components are true here
133 return eqIdx == contiMicrobialEqIdx || eqIdx == contiOxygenEqIdx || eqIdx == contiUreaEqIdx || eqIdx == contiBiofilmEqIdx || eqIdx == contiCalciteEqIdx;
134 }
135
136 static Scalar eqWeight([[maybe_unused]] unsigned eqIdx)
137 {
138 assert(eqApplies(eqIdx));
139
140 // TODO: it may be beneficial to chose this differently.
141 return static_cast<Scalar>(1.0);
142 }
143
144 // must be called after water storage is computed
145 template <class LhsEval>
146 static void addStorage(Dune::FieldVector<LhsEval, numEq>& storage,
147 const IntensiveQuantities& intQuants)
148 {
149 if (!enableMICP)
150 return;
151
152 LhsEval surfaceVolumeWater = Toolbox::template decay<LhsEval>(intQuants.porosity());
153 // avoid singular matrix if no water is present.
154 surfaceVolumeWater = max(surfaceVolumeWater, 1e-10);
155 // Suspended microbes in water phase
156 const LhsEval massMicrobes = surfaceVolumeWater * Toolbox::template decay<LhsEval>(intQuants.microbialConcentration());
157 LhsEval accumulationMicrobes = massMicrobes;
158 storage[contiMicrobialEqIdx] += accumulationMicrobes;
159 // Oxygen in water phase
160 const LhsEval massOxygen = surfaceVolumeWater * Toolbox::template decay<LhsEval>(intQuants.oxygenConcentration());
161 LhsEval accumulationOxygen = massOxygen;
162 storage[contiOxygenEqIdx] += accumulationOxygen;
163 // Urea in water phase
164 const LhsEval massUrea = surfaceVolumeWater * Toolbox::template decay<LhsEval>(intQuants.ureaConcentration());
165 LhsEval accumulationUrea = massUrea;
166 storage[contiUreaEqIdx] += accumulationUrea;
167 // Biofilm
168 const LhsEval massBiofilm = Toolbox::template decay<LhsEval>(intQuants.biofilmConcentration());
169 LhsEval accumulationBiofilm = massBiofilm;
170 storage[contiBiofilmEqIdx] += accumulationBiofilm;
171 // Calcite
172 const LhsEval massCalcite = Toolbox::template decay<LhsEval>(intQuants.calciteConcentration());
173 LhsEval accumulationCalcite = massCalcite;
174 storage[contiCalciteEqIdx] += accumulationCalcite;
175 }
176
177 static void computeFlux(RateVector& flux,
178 const ElementContext& elemCtx,
179 unsigned scvfIdx,
180 unsigned timeIdx)
181
182 {
183 if (!enableMICP)
184 return;
185
186 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
187
188 const unsigned upIdx = extQuants.upstreamIndex(waterPhaseIdx);
189 const unsigned inIdx = extQuants.interiorIndex();
190 const auto& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
191
192
193 if (upIdx == inIdx) {
194 flux[contiMicrobialEqIdx] = extQuants.volumeFlux(waterPhaseIdx) * up.microbialConcentration();
195 flux[contiOxygenEqIdx] = extQuants.volumeFlux(waterPhaseIdx) * up.oxygenConcentration();
196 flux[contiUreaEqIdx] = extQuants.volumeFlux(waterPhaseIdx) * up.ureaConcentration();
197 }
198 else {
199 flux[contiMicrobialEqIdx] = extQuants.volumeFlux(waterPhaseIdx) * decay<Scalar>(up.microbialConcentration());
200 flux[contiOxygenEqIdx] = extQuants.volumeFlux(waterPhaseIdx) * decay<Scalar>(up.oxygenConcentration());
201 flux[contiUreaEqIdx] = extQuants.volumeFlux(waterPhaseIdx) * decay<Scalar>(up.ureaConcentration());
202 }
203 }
204
205 // See https://doi.org/10.1016/j.ijggc.2021.103256 for the micp processes in the model.
206 static void addSource(RateVector& source,
207 const ElementContext& elemCtx,
208 unsigned dofIdx,
209 unsigned timeIdx)
210
211 {
212 if (!enableMICP)
213 return;
214
215 // compute dpW (max norm of the pressure gradient in the cell center)
216 const IntensiveQuantities& intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
217 const auto& K = elemCtx.problem().intrinsicPermeability(elemCtx, dofIdx, 0);
218 size_t numInteriorFaces = elemCtx.numInteriorFaces(timeIdx);
219 Evaluation dpW = 0;
220 for (unsigned scvfIdx = 0; scvfIdx < numInteriorFaces; scvfIdx++) {
221 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
222 unsigned upIdx = extQuants.upstreamIndex(waterPhaseIdx);
223 const auto& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
224 const Evaluation& mobWater = up.mobility(waterPhaseIdx);
225
226 // compute water velocity from flux
227 Evaluation waterVolumeVelocity = extQuants.volumeFlux(waterPhaseIdx) / (K[0][0] * mobWater);
228 dpW = std::max(dpW, abs(waterVolumeVelocity));
229 }
230
231 // get the model parameters
232 Scalar k_a = microbialAttachmentRate();
233 Scalar k_d = microbialDeathRate();
234 Scalar rho_b = densityBiofilm();
235 Scalar rho_c = densityCalcite();
236 Scalar k_str = detachmentRate();
237 Scalar k_o = halfVelocityOxygen();
238 Scalar k_u = halfVelocityUrea() / 10.0;//Dividing by scaling factor 10 (see WellInterface_impl.hpp)
239 Scalar mu = maximumGrowthRate();
240 Scalar mu_u = maximumUreaUtilization() / 10.0;//Dividing by scaling factor 10 (see WellInterface_impl.hpp)
241 Scalar Y_sb = yieldGrowthCoefficient();
242 Scalar F = oxygenConsumptionFactor();
243 Scalar Y_uc = 1.67 * 10; //Multiplying by scaling factor 10 (see WellInterface_impl.hpp)
244
245 // compute the processes
246 source[Indices::contiMicrobialEqIdx] += intQuants.microbialConcentration() * intQuants.porosity() *
247 (Y_sb * mu * intQuants.oxygenConcentration() / (k_o + intQuants.oxygenConcentration()) - k_d - k_a)
248 + rho_b * intQuants.biofilmConcentration() * k_str * pow(intQuants.porosity() * dpW, 0.58);
249
250 source[Indices::contiOxygenEqIdx] -= (intQuants.microbialConcentration() * intQuants.porosity() + rho_b * intQuants.biofilmConcentration()) *
251 F * mu * intQuants.oxygenConcentration() / (k_o + intQuants.oxygenConcentration());
252
253 source[Indices::contiUreaEqIdx] -= rho_b * intQuants.biofilmConcentration() * mu_u * intQuants.ureaConcentration() / (k_u + intQuants.ureaConcentration());
254
255 source[Indices::contiBiofilmEqIdx] += intQuants.biofilmConcentration() * (Y_sb * mu * intQuants.oxygenConcentration() / (k_o + intQuants.oxygenConcentration()) - k_d
256 - k_str * pow(intQuants.porosity() * dpW, 0.58) - Y_uc * (rho_b / rho_c) * intQuants.biofilmConcentration() * mu_u *
257 (intQuants.ureaConcentration() / (k_u + intQuants.ureaConcentration())) / (intQuants.porosity() + intQuants.biofilmConcentration()))
258 + k_a * intQuants.microbialConcentration() * intQuants.porosity() / rho_b;
259
260 source[Indices::contiCalciteEqIdx] += (rho_b / rho_c) * intQuants.biofilmConcentration() * Y_uc * mu_u * intQuants.ureaConcentration() / (k_u + intQuants.ureaConcentration());
261 }
262
263 static const Scalar densityBiofilm()
264 {
265 return params_.densityBiofilm_;
266 }
267
268 static const Scalar densityCalcite()
269 {
270 return params_.densityCalcite_;
271 }
272
273 static const Scalar detachmentRate()
274 {
275 return params_.detachmentRate_;
276 }
277
278 static const Scalar criticalPorosity()
279 {
280 return params_.criticalPorosity_;
281 }
282
283 static const Scalar fittingFactor()
284 {
285 return params_.fittingFactor_;
286 }
287
288 static const Scalar halfVelocityOxygen()
289 {
290 return params_.halfVelocityOxygen_;
291 }
292
293 static const Scalar halfVelocityUrea()
294 {
295 return params_.halfVelocityUrea_;
296 }
297
298 static const Scalar maximumGrowthRate()
299 {
300 return params_.maximumGrowthRate_;
301 }
302
303 static const Scalar maximumOxygenConcentration()
304 {
305 return params_.maximumOxygenConcentration_;
306 }
307
308 static const Scalar maximumUreaConcentration()
309 {
310 return params_.maximumUreaConcentration_ / 10.0;//Dividing by scaling factor 10 (see WellInterface_impl.hpp);
311 }
312
313 static const Scalar maximumUreaUtilization()
314 {
315 return params_.maximumUreaUtilization_;
316 }
317
318 static const Scalar microbialAttachmentRate()
319 {
320 return params_.microbialAttachmentRate_;
321 }
322
323 static const Scalar microbialDeathRate()
324 {
325 return params_.microbialDeathRate_;
326 }
327
328 static const Scalar minimumPermeability()
329 {
330 return params_.minimumPermeability_;
331 }
332
333 static const Scalar oxygenConsumptionFactor()
334 {
335 return params_.oxygenConsumptionFactor_;
336 }
337
338 static const Scalar toleranceBeforeClogging()
339 {
340 return params_.toleranceBeforeClogging_;
341 }
342
343 static const Scalar yieldGrowthCoefficient()
344 {
345 return params_.yieldGrowthCoefficient_;
346 }
347
348 static const std::vector<Scalar> phi()
349 {
350 return params_.phi_;
351 }
352
353private:
354 static BlackOilMICPParams<Scalar> params_;
355};
356
357
358template <class TypeTag, bool enableMICPV>
359BlackOilMICPParams<typename BlackOilMICPModule<TypeTag, enableMICPV>::Scalar>
360BlackOilMICPModule<TypeTag, enableMICPV>::params_;
361
369template <class TypeTag, bool enableMICPV = getPropValue<TypeTag, Properties::EnableMICP>()>
371{
373
380
382
383 static constexpr int microbialConcentrationIdx = Indices::microbialConcentrationIdx;
384 static constexpr int oxygenConcentrationIdx = Indices::oxygenConcentrationIdx;
385 static constexpr int ureaConcentrationIdx = Indices::ureaConcentrationIdx;
386 static constexpr int biofilmConcentrationIdx = Indices::biofilmConcentrationIdx;
387 static constexpr int calciteConcentrationIdx = Indices::calciteConcentrationIdx;
388 static constexpr int waterPhaseIdx = FluidSystem::waterPhaseIdx;
389
390
391public:
392
398 void MICPPropertiesUpdate_(const ElementContext& elemCtx,
399 unsigned dofIdx,
400 unsigned timeIdx)
401 {
402 const auto linearizationType = elemCtx.linearizationType();
403 const PrimaryVariables& priVars = elemCtx.primaryVars(dofIdx, timeIdx);
404 const auto& intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
405 const auto& K = elemCtx.problem().intrinsicPermeability(elemCtx, dofIdx, timeIdx);
406 Scalar referencePorosity_ = elemCtx.problem().porosity(elemCtx, dofIdx, timeIdx);
407 Scalar eta = MICPModule::fittingFactor();
408 Scalar k_min = MICPModule::minimumPermeability();
409 Scalar phi_crit = MICPModule::criticalPorosity();
410
411 microbialConcentration_ = priVars.makeEvaluation(microbialConcentrationIdx, timeIdx, linearizationType);
412 oxygenConcentration_ = priVars.makeEvaluation(oxygenConcentrationIdx, timeIdx, linearizationType);
413 ureaConcentration_ = priVars.makeEvaluation(ureaConcentrationIdx, timeIdx, linearizationType);
414 biofilmConcentration_ = priVars.makeEvaluation(biofilmConcentrationIdx, timeIdx, linearizationType);
415 calciteConcentration_ = priVars.makeEvaluation(calciteConcentrationIdx, timeIdx, linearizationType);
416
417 // Permeability reduction due to MICP, by adjusting the water mobility
418 asImp_().mobility_[waterPhaseIdx] *= max((pow((intQuants.porosity() - phi_crit) / (referencePorosity_ - phi_crit), eta) + k_min / K[0][0])/(1. + k_min / K[0][0]), k_min / K[0][0]);
419
420 }
421
422 const Evaluation& microbialConcentration() const
423 { return microbialConcentration_; }
424
425 const Evaluation& oxygenConcentration() const
426 { return oxygenConcentration_; }
427
428 const Evaluation& ureaConcentration() const
429 { return ureaConcentration_; }
430
431 const Evaluation& biofilmConcentration() const
432 { return biofilmConcentration_; }
433
434 const Evaluation& calciteConcentration() const
435 { return calciteConcentration_; }
436
437
438protected:
439 Implementation& asImp_()
440 { return *static_cast<Implementation*>(this); }
441
447
448};
449
450template <class TypeTag>
452{
456
457public:
458 void MICPPropertiesUpdate_(const ElementContext&,
459 unsigned,
460 unsigned)
461 { }
462
463 const Evaluation& microbialConcentration() const
464 { throw std::logic_error("microbialConcentration() called but MICP is disabled"); }
465
466 const Evaluation& oxygenConcentration() const
467 { throw std::logic_error("oxygenConcentration() called but MICP is disabled"); }
468
469 const Evaluation& ureaConcentration() const
470 { throw std::logic_error("ureaConcentration() called but MICP is disabled"); }
471
472 const Evaluation& biofilmConcentration() const
473 { throw std::logic_error("biofilmConcentration() called but MICP is disabled"); }
474
475 const Evaluation& calciteConcentration() const
476 { throw std::logic_error("calciteConcentration() called but MICP is disabled"); }
477};
478
486template <class TypeTag, bool enableMICPV = getPropValue<TypeTag, Properties::EnableMICP>()>
488{
490
491private:
492 Implementation& asImp_()
493 { return *static_cast<Implementation*>(this); }
494
495};
496
497template <class TypeTag>
498class BlackOilMICPExtensiveQuantities<TypeTag, false>{};
499
500} // namespace Opm
501
502#endif
Contains the parameters required to extend the black-oil model by MICP.
Declares the properties required by the black oil model.
Provides the MICP specific extensive quantities to the generic black-oil module's extensive quantitie...
Definition: blackoilmicpmodules.hh:488
void MICPPropertiesUpdate_(const ElementContext &, unsigned, unsigned)
Definition: blackoilmicpmodules.hh:458
const Evaluation & ureaConcentration() const
Definition: blackoilmicpmodules.hh:469
const Evaluation & oxygenConcentration() const
Definition: blackoilmicpmodules.hh:466
const Evaluation & microbialConcentration() const
Definition: blackoilmicpmodules.hh:463
const Evaluation & calciteConcentration() const
Definition: blackoilmicpmodules.hh:475
const Evaluation & biofilmConcentration() const
Definition: blackoilmicpmodules.hh:472
Provides the volumetric quantities required for the equations needed by the MICP extension of the bla...
Definition: blackoilmicpmodules.hh:371
void MICPPropertiesUpdate_(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Update the intensive properties needed to handle MICP from the primary variables.
Definition: blackoilmicpmodules.hh:398
Evaluation biofilmConcentration_
Definition: blackoilmicpmodules.hh:445
const Evaluation & biofilmConcentration() const
Definition: blackoilmicpmodules.hh:431
Implementation & asImp_()
Definition: blackoilmicpmodules.hh:439
Evaluation oxygenConcentration_
Definition: blackoilmicpmodules.hh:443
const Evaluation & calciteConcentration() const
Definition: blackoilmicpmodules.hh:434
Evaluation ureaConcentration_
Definition: blackoilmicpmodules.hh:444
const Evaluation & microbialConcentration() const
Definition: blackoilmicpmodules.hh:422
const Evaluation & ureaConcentration() const
Definition: blackoilmicpmodules.hh:428
Evaluation microbialConcentration_
Definition: blackoilmicpmodules.hh:442
Evaluation calciteConcentration_
Definition: blackoilmicpmodules.hh:446
const Evaluation & oxygenConcentration() const
Definition: blackoilmicpmodules.hh:425
Contains the high level supplements required to extend the black oil model by MICP.
Definition: blackoilmicpmodules.hh:49
static void checkCloggingMICP(const Model &model, const Scalar phi, unsigned dofIdx)
The simulator stops if "clogging" has been (almost) reached in any of the cells.
Definition: blackoilmicpmodules.hh:95
static void registerParameters()
Register all run-time parameters for the black-oil MICP module.
Definition: blackoilmicpmodules.hh:105
static const Scalar maximumGrowthRate()
Definition: blackoilmicpmodules.hh:298
static const Scalar yieldGrowthCoefficient()
Definition: blackoilmicpmodules.hh:343
static const Scalar detachmentRate()
Definition: blackoilmicpmodules.hh:273
static Scalar eqWeight(unsigned eqIdx)
Definition: blackoilmicpmodules.hh:136
static void computeFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Definition: blackoilmicpmodules.hh:177
static void addSource(RateVector &source, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx)
Definition: blackoilmicpmodules.hh:206
static const Scalar halfVelocityOxygen()
Definition: blackoilmicpmodules.hh:288
static const Scalar densityBiofilm()
Definition: blackoilmicpmodules.hh:263
static bool eqApplies(unsigned eqIdx)
Definition: blackoilmicpmodules.hh:127
static const Scalar toleranceBeforeClogging()
Definition: blackoilmicpmodules.hh:338
static const Scalar fittingFactor()
Definition: blackoilmicpmodules.hh:283
static const Scalar microbialDeathRate()
Definition: blackoilmicpmodules.hh:323
static const Scalar maximumUreaConcentration()
Definition: blackoilmicpmodules.hh:308
static void registerOutputModules(Model &model, Simulator &simulator)
Register all MICP specific VTK and ECL output modules.
Definition: blackoilmicpmodules.hh:117
static const Scalar maximumUreaUtilization()
Definition: blackoilmicpmodules.hh:313
static const std::vector< Scalar > phi()
Definition: blackoilmicpmodules.hh:348
static const Scalar densityCalcite()
Definition: blackoilmicpmodules.hh:268
static void addStorage(Dune::FieldVector< LhsEval, numEq > &storage, const IntensiveQuantities &intQuants)
Definition: blackoilmicpmodules.hh:146
static void setParams(BlackOilMICPParams< Scalar > &&params)
Set parameters.
Definition: blackoilmicpmodules.hh:83
static const Scalar minimumPermeability()
Definition: blackoilmicpmodules.hh:328
static const Scalar microbialAttachmentRate()
Definition: blackoilmicpmodules.hh:318
static const Scalar criticalPorosity()
Definition: blackoilmicpmodules.hh:278
static const Scalar maximumOxygenConcentration()
Definition: blackoilmicpmodules.hh:303
static const Scalar oxygenConsumptionFactor()
Definition: blackoilmicpmodules.hh:333
static const Scalar halfVelocityUrea()
Definition: blackoilmicpmodules.hh:293
VTK output module for the MICP model's related quantities.
Definition: vtkblackoilmicpmodule.hpp:53
static void registerParameters()
Register all run-time parameters for the multi-phase VTK output module.
Definition: vtkblackoilmicpmodule.hpp:82
Definition: blackoilboundaryratevector.hh:37
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:235
Struct holding the parameters for the BlackOilMICPModule class.
Definition: blackoilmicpparams.hpp:42