BlackoilModelParameters.hpp
Go to the documentation of this file.
1/*
2 Copyright 2015 SINTEF ICT, Applied Mathematics.
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_BLACKOILMODELPARAMETERS_HEADER_INCLUDED
21#define OPM_BLACKOILMODELPARAMETERS_HEADER_INCLUDED
22
23#include <opm/models/discretization/common/fvbaseproperties.hh>
24
25#include <opm/models/utils/basicproperties.hh>
26#include <opm/models/utils/parametersystem.hh>
27#include <opm/models/utils/propertysystem.hh>
28
30
31#include <stdexcept>
32#include <string>
33
34namespace Opm::Properties {
35
36namespace TTag {
38}
39
40template<class TypeTag, class MyTypeTag>
42 using type = UndefinedProperty;
43};
44template<class TypeTag, class MyTypeTag>
45struct DbhpMaxRel {
46 using type = UndefinedProperty;
47};
48template<class TypeTag, class MyTypeTag>
50 using type = UndefinedProperty;
51};
52template<class TypeTag, class MyTypeTag>
54 using type = UndefinedProperty;
55};
56template<class TypeTag, class MyTypeTag>
58 using type = UndefinedProperty;
59};
60template<class TypeTag, class MyTypeTag>
62 using type = UndefinedProperty;
63};
64template<class TypeTag, class MyTypeTag>
66 using type = UndefinedProperty;
67};
68template<class TypeTag, class MyTypeTag>
70 using type = UndefinedProperty;
71};
72template<class TypeTag, class MyTypeTag>
74 using type = UndefinedProperty;
75};
76template<class TypeTag, class MyTypeTag>
78 using type = UndefinedProperty;
79};
80template<class TypeTag, class MyTypeTag>
82 using type = UndefinedProperty;
83};
84template<class TypeTag, class MyTypeTag>
86 using type = UndefinedProperty;
87};
88template<class TypeTag, class MyTypeTag>
90 using type = UndefinedProperty;
91};
92template<class TypeTag, class MyTypeTag>
94 using type = UndefinedProperty;
95};
96template<class TypeTag, class MyTypeTag>
98 using type = UndefinedProperty;
99};
100template<class TypeTag, class MyTypeTag>
102 using type = UndefinedProperty;
103};
104template<class TypeTag, class MyTypeTag>
106 using type = UndefinedProperty;
107};
108template<class TypeTag, class MyTypeTag>
110 using type = UndefinedProperty;
111};
112template<class TypeTag, class MyTypeTag>
114 using type = UndefinedProperty;
115};
116template<class TypeTag, class MyTypeTag>
118 using type = UndefinedProperty;
119};
120template<class TypeTag, class MyTypeTag>
122 using type = UndefinedProperty;
123};
124template<class TypeTag, class MyTypeTag>
126 using type = UndefinedProperty;
127};
128template<class TypeTag, class MyTypeTag>
130 using type = UndefinedProperty;
131};
132// parameters for multisegment wells
133template<class TypeTag, class MyTypeTag>
135 using type = UndefinedProperty;
136};
137template<class TypeTag, class MyTypeTag>
139 using type = UndefinedProperty;
140};
141template<class TypeTag, class MyTypeTag>
143 using type = UndefinedProperty;
144};
145template<class TypeTag, class MyTypeTag>
147 using type = UndefinedProperty;
148};
149template<class TypeTag, class MyTypeTag>
151 using type = UndefinedProperty;
152};
153template<class TypeTag, class MyTypeTag>
155 using type = UndefinedProperty;
156};
157template<class TypeTag, class MyTypeTag>
159 using type = UndefinedProperty;
160};
161template<class TypeTag, class MyTypeTag>
163 using type = UndefinedProperty;
164};
165template<class TypeTag, class MyTypeTag>
167 using type = UndefinedProperty;
168};
169template<class TypeTag, class MyTypeTag>
171 using type = UndefinedProperty;
172};
173template<class TypeTag, class MyTypeTag>
175 using type = UndefinedProperty;
176};
177template<class TypeTag, class MyTypeTag>
179 using type = UndefinedProperty;
180};
181template<class TypeTag, class MyTypeTag>
183 using type = UndefinedProperty;
184};
185template<class TypeTag, class MyTypeTag>
187 using type = UndefinedProperty;
188};
189template<class TypeTag, class MyTypeTag>
191 using type = UndefinedProperty;
192};
193template<class TypeTag, class MyTypeTag>
195 using type = UndefinedProperty;
196};
197// Network solver parameters
198template<class TypeTag, class MyTypeTag>
200 using type = UndefinedProperty;
201};
202template<class TypeTag, class MyTypeTag>
204 using type = UndefinedProperty;
205};
206template<class TypeTag, class MyTypeTag>
208 using type = UndefinedProperty;
209};
210template<class TypeTag, class MyTypeTag>
212 using type = UndefinedProperty;
213};
214template<class TypeTag, class MyTypeTag>
216 using type = UndefinedProperty;
217};
218template<class TypeTag, class MyTypeTag>
220 using type = UndefinedProperty;
221};
222template<class TypeTag, class MyTypeTag>
224 using type = UndefinedProperty;
225};
226template<class TypeTag, class MyTypeTag>
228 using type = UndefinedProperty;
229};
230template<class TypeTag, class MyTypeTag>
232 using type = UndefinedProperty;
233};
234template<class TypeTag, class MyTypeTag>
236 using type = UndefinedProperty;
237};
238template<class TypeTag, class MyTypeTag>
240 using type = UndefinedProperty;
241};
242template<class TypeTag, class MyTypeTag>
244 using type = UndefinedProperty;
245};
246template<class TypeTag>
247struct DbhpMaxRel<TypeTag, TTag::FlowModelParameters> {
248 using type = GetPropType<TypeTag, Scalar>;
249 static constexpr type value = 1.0;
250};
251template<class TypeTag>
252struct DwellFractionMax<TypeTag, TTag::FlowModelParameters> {
253 using type = GetPropType<TypeTag, Scalar>;
254 static constexpr type value = 0.2;
255};
256template<class TypeTag>
257struct MaxResidualAllowed<TypeTag, TTag::FlowModelParameters> {
258 using type = GetPropType<TypeTag, Scalar>;
259 static constexpr type value = 1e7;
260};
261template<class TypeTag>
262struct RelaxedMaxPvFraction<TypeTag, TTag::FlowModelParameters> {
263 using type = GetPropType<TypeTag, Scalar>;
264 static constexpr type value = 0.03;
265};
266template<class TypeTag>
267struct ToleranceMb<TypeTag, TTag::FlowModelParameters> {
268 using type = GetPropType<TypeTag, Scalar>;
269 static constexpr type value = 1e-6;
270};
271template<class TypeTag>
272struct ToleranceMbRelaxed<TypeTag, TTag::FlowModelParameters> {
273 using type = GetPropType<TypeTag, Scalar>;
274 static constexpr type value = 1e-6;
275};
276template<class TypeTag>
277struct ToleranceCnv<TypeTag, TTag::FlowModelParameters> {
278 using type = GetPropType<TypeTag, Scalar>;
279 static constexpr type value = 1e-2;
280};
281template<class TypeTag>
282struct ToleranceCnvRelaxed<TypeTag, TTag::FlowModelParameters> {
283 using type = GetPropType<TypeTag, Scalar>;
284 static constexpr type value = 1;
285};
286template<class TypeTag>
287struct ToleranceWells<TypeTag, TTag::FlowModelParameters> {
288 using type = GetPropType<TypeTag, Scalar>;
289 static constexpr type value = 1e-4;
290};
291template<class TypeTag>
292struct ToleranceWellControl<TypeTag, TTag::FlowModelParameters> {
293 using type = GetPropType<TypeTag, Scalar>;
294 static constexpr type value = 1e-7;
295};
296template<class TypeTag>
297struct MaxWelleqIter<TypeTag, TTag::FlowModelParameters> {
298 static constexpr int value = 30;
299};
300template<class TypeTag>
301struct UseMultisegmentWell<TypeTag, TTag::FlowModelParameters> {
302 static constexpr bool value = true;
303};
304template<class TypeTag>
305struct MaxSinglePrecisionDays<TypeTag, TTag::FlowModelParameters> {
306 using type = GetPropType<TypeTag, Scalar>;
307 static constexpr type value = 20.0;
308};
309template<class TypeTag>
310struct MinStrictCnvIter<TypeTag, TTag::FlowModelParameters> {
311 static constexpr int value = 0;
312};
313template<class TypeTag>
314struct MinStrictMbIter<TypeTag, TTag::FlowModelParameters> {
315 static constexpr int value = -1;
316};
317template<class TypeTag>
318struct SolveWelleqInitially<TypeTag, TTag::FlowModelParameters> {
319 static constexpr bool value = true;
320};
321template<class TypeTag>
322struct UpdateEquationsScaling<TypeTag, TTag::FlowModelParameters> {
323 static constexpr bool value = false;
324};
325template<class TypeTag>
326struct UseUpdateStabilization<TypeTag, TTag::FlowModelParameters> {
327 static constexpr bool value = true;
328};
329template<class TypeTag>
330struct MatrixAddWellContributions<TypeTag, TTag::FlowModelParameters> {
331 static constexpr bool value = false;
332};
333template<class TypeTag>
334struct TolerancePressureMsWells<TypeTag, TTag::FlowModelParameters> {
335 using type = GetPropType<TypeTag, Scalar>;
336 static constexpr type value = 0.01*1e5;
337};
338template<class TypeTag>
339struct MaxPressureChangeMsWells<TypeTag, TTag::FlowModelParameters> {
340 using type = GetPropType<TypeTag, Scalar>;
341 static constexpr type value = 10*1e5;
342};
343template<class TypeTag>
344struct MaxNewtonIterationsWithInnerWellIterations<TypeTag, TTag::FlowModelParameters> {
345 static constexpr int value = 8;
346};
347template<class TypeTag>
348struct MaxInnerIterMsWells<TypeTag, TTag::FlowModelParameters> {
349 static constexpr int value = 100;
350};
351template<class TypeTag>
352struct MaxInnerIterWells<TypeTag, TTag::FlowModelParameters> {
353 static constexpr int value = 50;
354};
355template<class TypeTag>
356struct ShutUnsolvableWells<TypeTag, TTag::FlowModelParameters> {
357 static constexpr bool value = true;
358};
359template<class TypeTag>
360struct AlternativeWellRateInit<TypeTag, TTag::FlowModelParameters> {
361 static constexpr bool value = true;
362};
363template<class TypeTag>
364struct StrictOuterIterWells<TypeTag, TTag::FlowModelParameters> {
365 static constexpr int value = 6;
366};
367template<class TypeTag>
368struct StrictInnerIterWells<TypeTag, TTag::FlowModelParameters> {
369 static constexpr int value = 40;
370};
371template<class TypeTag>
372struct RegularizationFactorWells<TypeTag, TTag::FlowModelParameters> {
373 using type = GetPropType<TypeTag, Scalar>;
374 static constexpr type value = 100;
375};
376template<class TypeTag>
377struct EnableWellOperabilityCheck<TypeTag, TTag::FlowModelParameters> {
378 static constexpr bool value = true;
379};
380template<class TypeTag>
381struct EnableWellOperabilityCheckIter<TypeTag, TTag::FlowModelParameters> {
382 static constexpr bool value = false;
383};
384template<class TypeTag>
385struct DebugEmitCellPartition<TypeTag, TTag::FlowModelParameters> {
386 static constexpr bool value = false;
387};
388template<class TypeTag>
389struct RelaxedWellFlowTol<TypeTag, TTag::FlowModelParameters> {
390 using type = GetPropType<TypeTag, Scalar>;
391 static constexpr type value = 1e-3;
392};
393template<class TypeTag>
394struct RelaxedPressureTolMsw<TypeTag, TTag::FlowModelParameters> {
395 using type = GetPropType<TypeTag, Scalar>;
396 static constexpr type value = 1.0e4;
397};
398template<class TypeTag>
399struct MaximumNumberOfWellSwitches<TypeTag, TTag::FlowModelParameters> {
400 static constexpr int value = 3;
401};
402template<class TypeTag>
403struct UseAverageDensityMsWells<TypeTag, TTag::FlowModelParameters> {
404 static constexpr bool value = false;
405};
406template<class TypeTag>
407struct LocalWellSolveControlSwitching<TypeTag, TTag::FlowModelParameters> {
408 static constexpr bool value = false;
409};
410template<class TypeTag>
411struct UseImplicitIpr<TypeTag, TTag::FlowModelParameters> {
412 static constexpr bool value = false;
413};
414
415// Network solver parameters
416template<class TypeTag>
417struct NetworkMaxStrictIterations<TypeTag, TTag::FlowModelParameters> {
418 static constexpr int value = 100;
419};
420template<class TypeTag>
421struct NetworkMaxIterations<TypeTag, TTag::FlowModelParameters> {
422 static constexpr int value = 200;
423};
424template<class TypeTag>
425struct NonlinearSolver<TypeTag, TTag::FlowModelParameters> {
426 static constexpr auto value = "newton";
427};
428template<class TypeTag>
429struct LocalSolveApproach<TypeTag, TTag::FlowModelParameters> {
430 static constexpr auto value = "gauss-seidel";
431};
432template<class TypeTag>
433struct MaxLocalSolveIterations<TypeTag, TTag::FlowModelParameters> {
434 static constexpr int value = 20;
435};
436template<class TypeTag>
437struct LocalToleranceScalingMb<TypeTag, TTag::FlowModelParameters> {
438 using type = GetPropType<TypeTag, Scalar>;
439 static constexpr type value = 1.0;
440};
441template<class TypeTag>
442struct LocalToleranceScalingCnv<TypeTag, TTag::FlowModelParameters> {
443 using type = GetPropType<TypeTag, Scalar>;
444 static constexpr type value = 0.1;
445};
446template<class TypeTag>
447struct NlddNumInitialNewtonIter<TypeTag, TTag::FlowModelParameters> {
448 using type = int;
449 static constexpr auto value = type{1};
450};
451template<class TypeTag>
452struct NumLocalDomains<TypeTag, TTag::FlowModelParameters> {
453 using type = int;
454 static constexpr auto value = 0;
455};
456template<class TypeTag>
457struct LocalDomainsPartitioningImbalance<TypeTag, TTag::FlowModelParameters> {
458 using type = GetPropType<TypeTag, Scalar>;
459 static constexpr auto value = type{1.03};
460};
461template<class TypeTag>
462struct LocalDomainsPartitioningMethod<TypeTag, TTag::FlowModelParameters> {
463 static constexpr auto value = "zoltan";
464};
465template<class TypeTag>
466struct LocalDomainsOrderingMeasure<TypeTag, TTag::FlowModelParameters> {
467 static constexpr auto value = "maxpressure";
468};
469// if openMP is available, determine the number threads per process automatically.
470#if _OPENMP
471template<class TypeTag>
472struct ThreadsPerProcess<TypeTag, TTag::FlowModelParameters> {
473 static constexpr int value = -1;
474};
475#endif
476
477} // namespace Opm::Properties
478
479namespace Opm
480{
482 template <class TypeTag>
484 {
485 private:
486 using Scalar = GetPropType<TypeTag, Properties::Scalar>;
487
488 public:
509 // TODO: it might need to distinguish between rate control and pressure control later
515
518
521
524
527
530
533
536
539
542
545
549
552
555
558
561
564
571
573 std::string deck_file_name_;
574
577
582
585
588
591
594
597
600
602 std::string nonlinear_solver_;
605
607
610
616
617 bool write_partitions_{false};
618
621 {
622 dbhp_max_rel_= Parameters::get<TypeTag, Properties::DbhpMaxRel>();
623 dwell_fraction_max_ = Parameters::get<TypeTag, Properties::DwellFractionMax>();
624 max_residual_allowed_ = Parameters::get<TypeTag, Properties::MaxResidualAllowed>();
625 relaxed_max_pv_fraction_ = Parameters::get<TypeTag, Properties::RelaxedMaxPvFraction>();
626 tolerance_mb_ = Parameters::get<TypeTag, Properties::ToleranceMb>();
627 tolerance_mb_relaxed_ = Parameters::get<TypeTag, Properties::ToleranceMbRelaxed>();
628 tolerance_cnv_ = Parameters::get<TypeTag, Properties::ToleranceCnv>();
629 tolerance_cnv_relaxed_ = Parameters::get<TypeTag, Properties::ToleranceCnvRelaxed>();
630 tolerance_wells_ = Parameters::get<TypeTag, Properties::ToleranceWells>();
631 tolerance_well_control_ = Parameters::get<TypeTag, Properties::ToleranceWellControl>();
632 max_welleq_iter_ = Parameters::get<TypeTag, Properties::MaxWelleqIter>();
633 use_multisegment_well_ = Parameters::get<TypeTag, Properties::UseMultisegmentWell>();
634 tolerance_pressure_ms_wells_ = Parameters::get<TypeTag, Properties::TolerancePressureMsWells>();
635 relaxed_tolerance_flow_well_ = Parameters::get<TypeTag, Properties::RelaxedWellFlowTol>();
636 relaxed_tolerance_pressure_ms_well_ = Parameters::get<TypeTag, Properties::RelaxedPressureTolMsw>();
637 max_pressure_change_ms_wells_ = Parameters::get<TypeTag, Properties::MaxPressureChangeMsWells>();
638 max_inner_iter_ms_wells_ = Parameters::get<TypeTag, Properties::MaxInnerIterMsWells>();
639 strict_inner_iter_wells_ = Parameters::get<TypeTag, Properties::StrictInnerIterWells>();
640 strict_outer_iter_wells_ = Parameters::get<TypeTag, Properties::StrictOuterIterWells>();
641 regularization_factor_wells_ = Parameters::get<TypeTag, Properties::RegularizationFactorWells>();
642 max_niter_inner_well_iter_ = Parameters::get<TypeTag, Properties::MaxNewtonIterationsWithInnerWellIterations>();
643 shut_unsolvable_wells_ = Parameters::get<TypeTag, Properties::ShutUnsolvableWells>();
644 max_inner_iter_wells_ = Parameters::get<TypeTag, Properties::MaxInnerIterWells>();
645 maxSinglePrecisionTimeStep_ = Parameters::get<TypeTag, Properties::MaxSinglePrecisionDays>() * 24 * 60 * 60;
646 min_strict_cnv_iter_ = Parameters::get<TypeTag, Properties::MinStrictCnvIter>();
647 min_strict_mb_iter_ = Parameters::get<TypeTag, Properties::MinStrictMbIter>();
648 solve_welleq_initially_ = Parameters::get<TypeTag, Properties::SolveWelleqInitially>();
649 update_equations_scaling_ = Parameters::get<TypeTag, Properties::UpdateEquationsScaling>();
650 use_update_stabilization_ = Parameters::get<TypeTag, Properties::UseUpdateStabilization>();
651 matrix_add_well_contributions_ = Parameters::get<TypeTag, Properties::MatrixAddWellContributions>();
652 check_well_operability_ = Parameters::get<TypeTag, Properties::EnableWellOperabilityCheck>();
653 check_well_operability_iter_ = Parameters::get<TypeTag, Properties::EnableWellOperabilityCheckIter>();
654 max_number_of_well_switches_ = Parameters::get<TypeTag, Properties::MaximumNumberOfWellSwitches>();
655 use_average_density_ms_wells_ = Parameters::get<TypeTag, Properties::UseAverageDensityMsWells>();
656 local_well_solver_control_switching_ = Parameters::get<TypeTag, Properties::LocalWellSolveControlSwitching>();
657 use_implicit_ipr_ = Parameters::get<TypeTag, Properties::UseImplicitIpr>();
658 nonlinear_solver_ = Parameters::get<TypeTag, Properties::NonlinearSolver>();
659 const auto approach = Parameters::get<TypeTag, Properties::LocalSolveApproach>();
660 if (approach == "jacobi") {
662 } else if (approach == "gauss-seidel") {
664 } else {
665 throw std::runtime_error("Invalid domain solver approach '" + approach + "' specified.");
666 }
667
668 max_local_solve_iterations_ = Parameters::get<TypeTag, Properties::MaxLocalSolveIterations>();
669 local_tolerance_scaling_mb_ = Parameters::get<TypeTag, Properties::LocalToleranceScalingMb>();
670 local_tolerance_scaling_cnv_ = Parameters::get<TypeTag, Properties::LocalToleranceScalingCnv>();
671 nldd_num_initial_newton_iter_ = Parameters::get<TypeTag, Properties::NlddNumInitialNewtonIter>();
672 num_local_domains_ = Parameters::get<TypeTag, Properties::NumLocalDomains>();
673 local_domain_partition_imbalance_ = std::max(1.0, Parameters::get<TypeTag, Properties::LocalDomainsPartitioningImbalance>());
674 local_domain_partition_method_ = Parameters::get<TypeTag, Properties::LocalDomainsPartitioningMethod>();
675 deck_file_name_ = Parameters::get<TypeTag, Properties::EclDeckFileName>();
676 network_max_strict_iterations_ = Parameters::get<TypeTag, Properties::NetworkMaxStrictIterations>();
677 network_max_iterations_ = Parameters::get<TypeTag, Properties::NetworkMaxIterations>();
678 local_domain_ordering_ = domainOrderingMeasureFromString(Parameters::get<TypeTag, Properties::LocalDomainsOrderingMeasure>());
679 write_partitions_ = Parameters::get<TypeTag, Properties::DebugEmitCellPartition>();
680 }
681
682 static void registerParameters()
683 {
684 Parameters::registerParam<TypeTag, Properties::DbhpMaxRel>
685 ("Maximum relative change of the bottom-hole pressure in a single iteration");
686 Parameters::registerParam<TypeTag, Properties::DwellFractionMax>
687 ("Maximum absolute change of a well's volume fraction in a single iteration");
688 Parameters::registerParam<TypeTag, Properties::MaxResidualAllowed>
689 ("Absolute maximum tolerated for residuals without cutting the time step size");
690 Parameters::registerParam<TypeTag, Properties::RelaxedMaxPvFraction>
691 ("The fraction of the pore volume of the reservoir "
692 "where the volumetric error (CNV) may be voilated "
693 "during strict Newton iterations.");
694 Parameters::registerParam<TypeTag, Properties::ToleranceMb>
695 ("Tolerated mass balance error relative to total mass present");
696 Parameters::registerParam<TypeTag, Properties::ToleranceMbRelaxed>
697 ("Relaxed tolerated mass balance error that applies for iterations "
698 "after the iterations with the strict tolerance");
699 Parameters::registerParam<TypeTag, Properties::ToleranceCnv>
700 ("Local convergence tolerance (Maximum of local saturation errors)");
701 Parameters::registerParam<TypeTag, Properties::ToleranceCnvRelaxed>
702 ("Relaxed local convergence tolerance that applies for iterations "
703 "after the iterations with the strict tolerance");
704 Parameters::registerParam<TypeTag, Properties::ToleranceWells>
705 ("Well convergence tolerance");
706 Parameters::registerParam<TypeTag, Properties::ToleranceWellControl>
707 ("Tolerance for the well control equations");
708 Parameters::registerParam<TypeTag, Properties::MaxWelleqIter>
709 ("Maximum number of iterations to determine solution the well equations");
710 Parameters::registerParam<TypeTag, Properties::UseMultisegmentWell>
711 ("Use the well model for multi-segment wells instead of the "
712 "one for single-segment wells");
713 Parameters::registerParam<TypeTag, Properties::TolerancePressureMsWells>
714 ("Tolerance for the pressure equations for multi-segment wells");
715 Parameters::registerParam<TypeTag, Properties::RelaxedWellFlowTol>
716 ("Relaxed tolerance for the well flow residual");
717 Parameters::registerParam<TypeTag, Properties::RelaxedPressureTolMsw>
718 ("Relaxed tolerance for the MSW pressure solution");
719 Parameters::registerParam<TypeTag, Properties::MaxPressureChangeMsWells>
720 ("Maximum relative pressure change for a single iteration "
721 "of the multi-segment well model");
722 Parameters::registerParam<TypeTag, Properties::MaxInnerIterMsWells>
723 ("Maximum number of inner iterations for multi-segment wells");
724 Parameters::registerParam<TypeTag, Properties::StrictInnerIterWells>
725 ("Number of inner well iterations with strict tolerance");
726 Parameters::registerParam<TypeTag, Properties::StrictOuterIterWells>
727 ("Number of newton iterations for which wells are checked with strict tolerance");
728 Parameters::registerParam<TypeTag, Properties::MaxNewtonIterationsWithInnerWellIterations>
729 ("Maximum newton iterations with inner well iterations");
730 Parameters::registerParam<TypeTag, Properties::ShutUnsolvableWells>
731 ("Shut unsolvable wells");
732 Parameters::registerParam<TypeTag, Properties::MaxInnerIterWells>
733 ("Maximum number of inner iterations for standard wells");
734 Parameters::registerParam<TypeTag, Properties::AlternativeWellRateInit>
735 ("Use alternative well rate initialization procedure");
736 Parameters::registerParam<TypeTag, Properties::RegularizationFactorWells>
737 ("Regularization factor for wells");
738 Parameters::registerParam<TypeTag, Properties::MaxSinglePrecisionDays>
739 ("Maximum time step size where single precision floating point "
740 "arithmetic can be used solving for the linear systems of equations");
741 Parameters::registerParam<TypeTag, Properties::MinStrictCnvIter>
742 ("Minimum number of Newton iterations before relaxed tolerances "
743 "can be used for the CNV convergence criterion");
744 Parameters::registerParam<TypeTag, Properties::MinStrictMbIter>
745 ("Minimum number of Newton iterations before relaxed tolerances "
746 "can be used for the MB convergence criterion. "
747 "Default -1 means that the relaxed tolerance is used when maximum "
748 "number of Newton iterations are reached.");
749 Parameters::registerParam<TypeTag, Properties::SolveWelleqInitially>
750 ("Fully solve the well equations before each iteration of the reservoir model");
751 Parameters::registerParam<TypeTag, Properties::UpdateEquationsScaling>
752 ("Update scaling factors for mass balance equations during the run");
753 Parameters::registerParam<TypeTag, Properties::UseUpdateStabilization>
754 ("Try to detect and correct oscillations or stagnation during the Newton method");
755 Parameters::registerParam<TypeTag, Properties::MatrixAddWellContributions>
756 ("Explicitly specify the influences of wells between cells in "
757 "the Jacobian and preconditioner matrices");
758 Parameters::registerParam<TypeTag, Properties::EnableWellOperabilityCheck>
759 ("Enable the well operability checking");
760 Parameters::registerParam<TypeTag, Properties::EnableWellOperabilityCheckIter>
761 ("Enable the well operability checking during iterations");
762 Parameters::registerParam<TypeTag, Properties::MaximumNumberOfWellSwitches>
763 ("Maximum number of times a well can switch to the same control");
764 Parameters::registerParam<TypeTag, Properties::UseAverageDensityMsWells>
765 ("Approximate segment densitities by averaging over segment and its outlet");
766 Parameters::registerParam<TypeTag, Properties::LocalWellSolveControlSwitching>
767 ("Allow control switching during local well solutions");
768 Parameters::registerParam<TypeTag, Properties::UseImplicitIpr>
769 ("Compute implict IPR for stability checks and stable solution search");
770 Parameters::registerParam<TypeTag, Properties::NetworkMaxStrictIterations>
771 ("Maximum iterations in network solver before relaxing tolerance");
772 Parameters::registerParam<TypeTag, Properties::NetworkMaxIterations>
773 ("Maximum number of iterations in the network solver before giving up");
774 Parameters::registerParam<TypeTag, Properties::NonlinearSolver>
775 ("Choose nonlinear solver. Valid choices are newton or nldd.");
776 Parameters::registerParam<TypeTag, Properties::LocalSolveApproach>
777 ("Choose local solve approach. Valid choices are jacobi and gauss-seidel");
778 Parameters::registerParam<TypeTag, Properties::MaxLocalSolveIterations>
779 ("Max iterations for local solves with NLDD nonlinear solver.");
780 Parameters::registerParam<TypeTag, Properties::LocalToleranceScalingMb>
781 ("Set lower than 1.0 to use stricter convergence tolerance for local solves.");
782 Parameters::registerParam<TypeTag, Properties::LocalToleranceScalingCnv>
783 ("Set lower than 1.0 to use stricter convergence tolerance for local solves.");
784 Parameters::registerParam<TypeTag, Properties::NlddNumInitialNewtonIter>
785 ("Number of initial global Newton iterations when running the NLDD nonlinear solver.");
786 Parameters::registerParam<TypeTag, Properties::NumLocalDomains>
787 ("Number of local domains for NLDD nonlinear solver.");
788 Parameters::registerParam<TypeTag, Properties::LocalDomainsPartitioningImbalance>
789 ("Subdomain partitioning imbalance tolerance. 1.03 is 3 percent imbalance.");
790 Parameters::registerParam<TypeTag, Properties::LocalDomainsPartitioningMethod>
791 ("Subdomain partitioning method. Allowed values are "
792 "'zoltan', "
793 "'simple', "
794 "and the name of a partition file ending with '.partition'.");
795 Parameters::registerParam<TypeTag, Properties::LocalDomainsOrderingMeasure>
796 ("Subdomain ordering measure. Allowed values are "
797 "'maxpressure', "
798 "'averagepressure' "
799 "and 'residual'.");
800 Parameters::registerParam<TypeTag, Properties::DebugEmitCellPartition>
801 ("Whether or not to emit cell partitions as a debugging aid.");
802
803 Parameters::hideParam<TypeTag, Properties::DebugEmitCellPartition>();
804 }
805 };
806} // namespace Opm
807
808#endif // OPM_BLACKOILMODELPARAMETERS_HEADER_INCLUDED
Definition: AluGridVanguard.hpp:57
Definition: BlackoilPhases.hpp:27
DomainOrderingMeasure
Measure to use for domain ordering.
Definition: SubDomain.hpp:39
DomainOrderingMeasure domainOrderingMeasureFromString(const std::string_view measure)
Definition: SubDomain.hpp:45
DomainSolveApproach
Solver approach for NLDD.
Definition: SubDomain.hpp:33
Solver parameters for the BlackoilModel.
Definition: BlackoilModelParameters.hpp:484
bool shut_unsolvable_wells_
Whether to shut unsolvable well.
Definition: BlackoilModelParameters.hpp:538
int num_local_domains_
Definition: BlackoilModelParameters.hpp:612
double local_tolerance_scaling_mb_
Definition: BlackoilModelParameters.hpp:608
bool check_well_operability_
Whether to check well operability.
Definition: BlackoilModelParameters.hpp:579
double relaxed_max_pv_fraction_
Definition: BlackoilModelParameters.hpp:497
bool update_equations_scaling_
Update scaling factors for mass balance equations.
Definition: BlackoilModelParameters.hpp:560
double maxSinglePrecisionTimeStep_
Definition: BlackoilModelParameters.hpp:548
double max_pressure_change_ms_wells_
Maximum pressure change over an iteratio for ms wells.
Definition: BlackoilModelParameters.hpp:520
double tolerance_cnv_relaxed_
Relaxed local convergence tolerance (can be used when iter >= min_strict_cnv_iter_ && cnvViolatedPV <...
Definition: BlackoilModelParameters.hpp:505
std::string local_domain_partition_method_
Definition: BlackoilModelParameters.hpp:614
bool solve_welleq_initially_
Solve well equation initially.
Definition: BlackoilModelParameters.hpp:557
bool write_partitions_
Definition: BlackoilModelParameters.hpp:617
double max_residual_allowed_
Absolute max limit for residuals.
Definition: BlackoilModelParameters.hpp:494
bool use_implicit_ipr_
Whether to use implicit IPR for thp stability checks and solution search.
Definition: BlackoilModelParameters.hpp:593
bool use_average_density_ms_wells_
Whether to approximate segment densities by averaging over segment and its outlet.
Definition: BlackoilModelParameters.hpp:587
int min_strict_mb_iter_
Minimum number of Newton iterations before we can use relaxed MB convergence criterion.
Definition: BlackoilModelParameters.hpp:554
double local_tolerance_scaling_cnv_
Definition: BlackoilModelParameters.hpp:609
double tolerance_wells_
Well convergence tolerance.
Definition: BlackoilModelParameters.hpp:507
double tolerance_cnv_
Local convergence tolerance (max of local saturation errors).
Definition: BlackoilModelParameters.hpp:503
double tolerance_pressure_ms_wells_
Tolerance for the pressure equations for multisegment wells.
Definition: BlackoilModelParameters.hpp:512
int network_max_iterations_
Maximum number of iterations in the network solver before giving up.
Definition: BlackoilModelParameters.hpp:599
int max_welleq_iter_
Maximum iteration number of the well equation solution.
Definition: BlackoilModelParameters.hpp:544
bool matrix_add_well_contributions_
Whether to add influences of wells between cells to the matrix and preconditioner matrix.
Definition: BlackoilModelParameters.hpp:576
DomainOrderingMeasure local_domain_ordering_
Definition: BlackoilModelParameters.hpp:615
double dbhp_max_rel_
Max relative change in bhp in single iteration.
Definition: BlackoilModelParameters.hpp:490
bool local_well_solver_control_switching_
Whether to allow control switching during local well solutions.
Definition: BlackoilModelParameters.hpp:590
int strict_outer_iter_wells_
Newton iteration where wells are stricly convergent.
Definition: BlackoilModelParameters.hpp:529
int nldd_num_initial_newton_iter_
Definition: BlackoilModelParameters.hpp:611
std::string deck_file_name_
The file name of the deck.
Definition: BlackoilModelParameters.hpp:573
int max_number_of_well_switches_
Maximum number of times a well can switch to the same controt.
Definition: BlackoilModelParameters.hpp:584
DomainSolveApproach local_solve_approach_
'jacobi' and 'gauss-seidel' supported.
Definition: BlackoilModelParameters.hpp:604
int max_local_solve_iterations_
Definition: BlackoilModelParameters.hpp:606
bool use_update_stabilization_
Try to detect oscillation or stagnation.
Definition: BlackoilModelParameters.hpp:563
int max_niter_inner_well_iter_
Maximum newton iterations with inner well iterations.
Definition: BlackoilModelParameters.hpp:535
double dwell_fraction_max_
Max absolute change in well volume fraction in single iteration.
Definition: BlackoilModelParameters.hpp:492
double relaxed_tolerance_pressure_ms_well_
Relaxed tolerance for the MSW pressure solution.
Definition: BlackoilModelParameters.hpp:517
std::string nonlinear_solver_
Nonlinear solver type: newton or nldd.
Definition: BlackoilModelParameters.hpp:602
double tolerance_well_control_
Tolerance for the well control equations.
Definition: BlackoilModelParameters.hpp:510
double regularization_factor_wells_
Regularization factor for wells.
Definition: BlackoilModelParameters.hpp:532
double tolerance_mb_relaxed_
Relaxed mass balance tolerance (can be used when iter >= min_strict_mb_iter_).
Definition: BlackoilModelParameters.hpp:501
bool use_multisegment_well_
Definition: BlackoilModelParameters.hpp:570
int network_max_strict_iterations_
Maximum number of iterations in the network solver before relaxing tolerance.
Definition: BlackoilModelParameters.hpp:596
int max_inner_iter_ms_wells_
Maximum inner iteration number for ms wells.
Definition: BlackoilModelParameters.hpp:523
double local_domain_partition_imbalance_
Definition: BlackoilModelParameters.hpp:613
int max_inner_iter_wells_
Maximum inner iteration number for standard wells.
Definition: BlackoilModelParameters.hpp:541
double relaxed_tolerance_flow_well_
Relaxed tolerance for for the well flow residual.
Definition: BlackoilModelParameters.hpp:514
static void registerParameters()
Definition: BlackoilModelParameters.hpp:682
int min_strict_cnv_iter_
Minimum number of Newton iterations before we can use relaxed CNV convergence criterion.
Definition: BlackoilModelParameters.hpp:551
BlackoilModelParameters()
Construct from user parameters or defaults.
Definition: BlackoilModelParameters.hpp:620
bool check_well_operability_iter_
Whether to check well operability during iterations.
Definition: BlackoilModelParameters.hpp:581
int strict_inner_iter_wells_
Strict inner iteration number for wells.
Definition: BlackoilModelParameters.hpp:526
double tolerance_mb_
Relative mass balance tolerance (total mass balance error).
Definition: BlackoilModelParameters.hpp:499
Definition: BlackoilModelParameters.hpp:178
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:179
GetPropType< TypeTag, Scalar > type
Definition: BlackoilModelParameters.hpp:248
Definition: BlackoilModelParameters.hpp:45
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:46
Definition: BlackoilModelParameters.hpp:129
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:130
GetPropType< TypeTag, Scalar > type
Definition: BlackoilModelParameters.hpp:253
Definition: BlackoilModelParameters.hpp:49
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:50
Definition: BlackoilModelParameters.hpp:41
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:42
Definition: BlackoilModelParameters.hpp:125
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:126
Definition: BlackoilModelParameters.hpp:121
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:122
Definition: BlackoilModelParameters.hpp:243
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:244
GetPropType< TypeTag, Scalar > type
Definition: BlackoilModelParameters.hpp:458
Definition: BlackoilModelParameters.hpp:235
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:236
Definition: BlackoilModelParameters.hpp:239
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:240
Definition: BlackoilModelParameters.hpp:211
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:212
GetPropType< TypeTag, Scalar > type
Definition: BlackoilModelParameters.hpp:443
Definition: BlackoilModelParameters.hpp:223
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:224
GetPropType< TypeTag, Scalar > type
Definition: BlackoilModelParameters.hpp:438
Definition: BlackoilModelParameters.hpp:219
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:220
Definition: BlackoilModelParameters.hpp:190
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:191
Definition: BlackoilModelParameters.hpp:117
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:118
Definition: BlackoilModelParameters.hpp:142
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:143
Definition: BlackoilModelParameters.hpp:174
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:175
Definition: BlackoilModelParameters.hpp:215
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:216
Definition: BlackoilModelParameters.hpp:166
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:167
GetPropType< TypeTag, Scalar > type
Definition: BlackoilModelParameters.hpp:340
Definition: BlackoilModelParameters.hpp:138
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:139
GetPropType< TypeTag, Scalar > type
Definition: BlackoilModelParameters.hpp:258
Definition: BlackoilModelParameters.hpp:53
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:54
GetPropType< TypeTag, Scalar > type
Definition: BlackoilModelParameters.hpp:306
Definition: BlackoilModelParameters.hpp:93
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:94
Definition: BlackoilModelParameters.hpp:85
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:86
Definition: BlackoilModelParameters.hpp:182
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:183
Definition: BlackoilModelParameters.hpp:97
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:98
Definition: BlackoilModelParameters.hpp:101
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:102
Definition: BlackoilModelParameters.hpp:203
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:204
Definition: BlackoilModelParameters.hpp:199
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:200
Definition: BlackoilModelParameters.hpp:227
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:228
Definition: BlackoilModelParameters.hpp:207
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:208
int type
Definition: BlackoilModelParameters.hpp:453
Definition: BlackoilModelParameters.hpp:231
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:232
GetPropType< TypeTag, Scalar > type
Definition: BlackoilModelParameters.hpp:373
Definition: BlackoilModelParameters.hpp:162
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:163
GetPropType< TypeTag, Scalar > type
Definition: BlackoilModelParameters.hpp:263
Definition: BlackoilModelParameters.hpp:57
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:58
GetPropType< TypeTag, Scalar > type
Definition: BlackoilModelParameters.hpp:395
Definition: BlackoilModelParameters.hpp:158
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:159
GetPropType< TypeTag, Scalar > type
Definition: BlackoilModelParameters.hpp:390
Definition: BlackoilModelParameters.hpp:150
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:151
Definition: BlackoilModelParameters.hpp:170
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:171
Definition: BlackoilModelParameters.hpp:105
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:106
Definition: BlackoilModelParameters.hpp:146
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:147
Definition: BlackoilModelParameters.hpp:154
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:155
Definition: BlackoilModelParameters.hpp:37
GetPropType< TypeTag, Scalar > type
Definition: BlackoilModelParameters.hpp:278
Definition: BlackoilModelParameters.hpp:69
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:70
GetPropType< TypeTag, Scalar > type
Definition: BlackoilModelParameters.hpp:283
Definition: BlackoilModelParameters.hpp:73
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:74
GetPropType< TypeTag, Scalar > type
Definition: BlackoilModelParameters.hpp:268
Definition: BlackoilModelParameters.hpp:61
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:62
GetPropType< TypeTag, Scalar > type
Definition: BlackoilModelParameters.hpp:273
Definition: BlackoilModelParameters.hpp:65
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:66
GetPropType< TypeTag, Scalar > type
Definition: BlackoilModelParameters.hpp:335
Definition: BlackoilModelParameters.hpp:134
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:135
GetPropType< TypeTag, Scalar > type
Definition: BlackoilModelParameters.hpp:293
Definition: BlackoilModelParameters.hpp:81
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:82
GetPropType< TypeTag, Scalar > type
Definition: BlackoilModelParameters.hpp:288
Definition: BlackoilModelParameters.hpp:77
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:78
Definition: BlackoilModelParameters.hpp:109
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:110
Definition: BlackoilModelParameters.hpp:186
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:187
Definition: BlackoilModelParameters.hpp:194
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:195
Definition: BlackoilModelParameters.hpp:89
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:90
Definition: BlackoilModelParameters.hpp:113
UndefinedProperty type
Definition: BlackoilModelParameters.hpp:114