/* Copyright 2015 SINTEF ICT, Applied Mathematics. This file is part of the Open Porous Media project (OPM). OPM is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version. OPM is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with OPM. If not, see . */ #include #include #include #include #include #include namespace Opm { template BlackoilModelParameters::BlackoilModelParameters() { dbhp_max_rel_ = Parameters::Get>(); dwell_fraction_max_ = Parameters::Get>(); max_residual_allowed_ = Parameters::Get>(); relaxed_max_pv_fraction_ = Parameters::Get>(); tolerance_mb_ = Parameters::Get>(); tolerance_mb_relaxed_ = std::max(tolerance_mb_, Parameters::Get>()); tolerance_energy_balance_ = Parameters::Get>(); tolerance_energy_balance_relaxed_ = std::max(tolerance_energy_balance_, Parameters::Get>()); tolerance_cnv_ = Parameters::Get>(); tolerance_cnv_relaxed_ = std::max(tolerance_cnv_, Parameters::Get>()); tolerance_cnv_energy_ = Parameters::Get>(); tolerance_cnv_energy_relaxed_ = std::max(tolerance_cnv_energy_, Parameters::Get>()); tolerance_wells_ = Parameters::Get>(); tolerance_well_control_ = Parameters::Get>(); max_welleq_iter_ = Parameters::Get(); use_multisegment_well_ = Parameters::Get(); tolerance_pressure_ms_wells_ = Parameters::Get>(); relaxed_tolerance_flow_well_ = Parameters::Get>(); relaxed_tolerance_pressure_ms_well_ = Parameters::Get>(); max_pressure_change_ms_wells_ = Parameters::Get>(); max_inner_iter_ms_wells_ = Parameters::Get(); strict_inner_iter_wells_ = Parameters::Get(); strict_outer_iter_wells_ = Parameters::Get(); regularization_factor_wells_ = Parameters::Get>(); max_niter_inner_well_iter_ = Parameters::Get(); shut_unsolvable_wells_ = Parameters::Get(); max_inner_iter_wells_ = Parameters::Get(); maxSinglePrecisionTimeStep_ = Parameters::Get>() * 24 * 60 * 60; min_strict_cnv_iter_ = Parameters::Get(); min_strict_mb_iter_ = Parameters::Get(); solve_welleq_initially_ = Parameters::Get(); update_equations_scaling_ = Parameters::Get(); use_update_stabilization_ = Parameters::Get(); matrix_add_well_contributions_ = Parameters::Get(); check_well_operability_ = Parameters::Get(); check_well_operability_iter_ = Parameters::Get(); max_number_of_well_switches_ = Parameters::Get(); use_average_density_ms_wells_ = Parameters::Get(); local_well_solver_control_switching_ = Parameters::Get(); use_implicit_ipr_ = Parameters::Get(); nonlinear_solver_ = Parameters::Get(); const auto approach = Parameters::Get(); if (approach == "jacobi") { local_solve_approach_ = DomainSolveApproach::Jacobi; } else if (approach == "gauss-seidel") { local_solve_approach_ = DomainSolveApproach::GaussSeidel; } else { throw std::runtime_error("Invalid domain solver approach '" + approach + "' specified."); } max_local_solve_iterations_ = Parameters::Get(); local_tolerance_scaling_mb_ = Parameters::Get>(); local_tolerance_scaling_cnv_ = Parameters::Get>(); nldd_num_initial_newton_iter_ = Parameters::Get(); num_local_domains_ = Parameters::Get(); local_domain_partition_imbalance_ = std::max(Scalar{1.0}, Parameters::Get>()); local_domain_partition_method_ = Parameters::Get(); deck_file_name_ = Parameters::Get(); network_max_strict_iterations_ = Parameters::Get(); network_max_iterations_ = Parameters::Get(); local_domain_ordering_ = domainOrderingMeasureFromString(Parameters::Get()); write_partitions_ = Parameters::Get(); } template void BlackoilModelParameters::registerParameters() { Parameters::Register> ("Maximum relative change of the bottom-hole pressure in a single iteration"); Parameters::Register> ("Maximum absolute change of a well's volume fraction in a single iteration"); Parameters::Register> ("Absolute maximum tolerated for residuals without cutting the time step size"); Parameters::Register> ("The fraction of the pore volume of the reservoir " "where the volumetric error (CNV) may be voilated " "during strict Newton iterations."); Parameters::Register> ("Tolerated mass balance error relative to total mass present"); Parameters::Register> ("Relaxed tolerated mass balance error that applies for iterations " "after the iterations with the strict tolerance"); Parameters::Register> ("Tolerated energy balance error relative to (scaled) total energy present"); Parameters::Register> ("Relaxed tolerated energy balance error that applies for iterations " "after the iterations with the strict tolerance"); Parameters::Register> ("Local convergence tolerance (Maximum of local saturation errors)"); Parameters::Register> ("Relaxed local convergence tolerance that applies for iterations " "after the iterations with the strict tolerance"); Parameters::Register> ("Local energy convergence tolerance (Maximum of local energy errors)"); Parameters::Register> ("Relaxed local energy convergence tolerance that applies for iterations " "after the iterations with the strict tolerance"); Parameters::Register> ("Well convergence tolerance"); Parameters::Register> ("Tolerance for the well control equations"); Parameters::Register ("Maximum number of iterations to determine solution the well equations"); Parameters::Register ("Use the well model for multi-segment wells instead of the " "one for single-segment wells"); Parameters::Register> ("Tolerance for the pressure equations for multi-segment wells"); Parameters::Register> ("Relaxed tolerance for the well flow residual"); Parameters::Register> ("Relaxed tolerance for the MSW pressure solution"); Parameters::Register> ("Maximum relative pressure change for a single iteration " "of the multi-segment well model"); Parameters::Register ("Maximum number of inner iterations for multi-segment wells"); Parameters::Register ("Number of inner well iterations with strict tolerance"); Parameters::Register ("Number of newton iterations for which wells are checked with strict tolerance"); Parameters::Register ("Maximum newton iterations with inner well iterations"); Parameters::Register ("Shut unsolvable wells"); Parameters::Register ("Maximum number of inner iterations for standard wells"); Parameters::Register ("Use alternative well rate initialization procedure"); Parameters::Register> ("Regularization factor for wells"); Parameters::Register> ("Maximum time step size where single precision floating point " "arithmetic can be used solving for the linear systems of equations"); Parameters::Register ("Minimum number of Newton iterations before relaxed tolerances " "can be used for the CNV convergence criterion"); Parameters::Register ("Minimum number of Newton iterations before relaxed tolerances " "can be used for the MB convergence criterion. " "Default -1 means that the relaxed tolerance is used when maximum " "number of Newton iterations are reached."); Parameters::Register ("Fully solve the well equations before each iteration of the reservoir model"); Parameters::Register ("Update scaling factors for mass balance equations during the run"); Parameters::Register ("Try to detect and correct oscillations or stagnation during the Newton method"); Parameters::Register ("Explicitly specify the influences of wells between cells in " "the Jacobian and preconditioner matrices"); Parameters::Register ("Enable the well operability checking"); Parameters::Register ("Enable the well operability checking during iterations"); Parameters::Register ("Maximum number of times a well can switch to the same control"); Parameters::Register ("Approximate segment densitities by averaging over segment and its outlet"); Parameters::Register ("Allow control switching during local well solutions"); Parameters::Register ("Compute implict IPR for stability checks and stable solution search"); Parameters::Register ("Maximum iterations in network solver before relaxing tolerance"); Parameters::Register ("Maximum number of iterations in the network solver before giving up"); Parameters::Register ("Choose nonlinear solver. Valid choices are newton or nldd."); Parameters::Register ("Choose local solve approach. Valid choices are jacobi and gauss-seidel"); Parameters::Register ("Max iterations for local solves with NLDD nonlinear solver."); Parameters::Register> ("Set lower than 1.0 to use stricter convergence tolerance for local solves."); Parameters::Register> ("Set lower than 1.0 to use stricter convergence tolerance for local solves."); Parameters::Register ("Number of initial global Newton iterations when running the NLDD nonlinear solver."); Parameters::Register ("Number of local domains for NLDD nonlinear solver."); Parameters::Register> ("Subdomain partitioning imbalance tolerance. 1.03 is 3 percent imbalance."); Parameters::Register ("Subdomain partitioning method. Allowed values are " "'zoltan', " "'simple', " "and the name of a partition file ending with '.partition'."); Parameters::Register ("Subdomain ordering measure. Allowed values are " "'maxpressure', " "'averagepressure' " "and 'residual'."); Parameters::Register ("Whether or not to emit cell partitions as a debugging aid."); Parameters::Hide(); // if openMP is available, determine the number threads per process automatically. #if _OPENMP Parameters::SetDefault(-1); #endif } template struct BlackoilModelParameters; #if FLOW_INSTANTIATE_FLOAT template struct BlackoilModelParameters; #endif } // namespace Opm