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