// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*- // vi: set et ts=4 sw=4 sts=4: /* 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 2 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 . Consult the COPYING file in the top-level source directory of this module for the precise wording of the license and the list of copyright holders. */ /*! * \file * * \brief The common settings for all flowexp variants. */ #ifndef FLOW_EXP_HPP #define FLOW_EXP_HPP #include #include #include #include #include #include #include #include #include #include namespace Opm { template class FlowExpProblem; } namespace Opm::Properties { namespace TTag { struct FlowExpTypeTag { using InheritsFrom = std::tuple; }; } // Set the problem class template struct Problem { using type = FlowExpProblem; }; // Enable experimental features for flowexp: flowexp is the research simulator of the OPM // project. If you're looking for a more stable "production quality" simulator, consider // using `flow` template struct EnableExperiments { static constexpr bool value = true; }; // use flow's well model for now template struct WellModel { using type = BlackoilWellModel; }; template struct NewtonMethod { using type = FlowExpNewtonMethod; }; // flow's well model only works with surface volumes template struct BlackoilConserveSurfaceVolume { static constexpr bool value = true; }; // the values for the residual are for the whole cell instead of for a cubic meter of the cell template struct UseVolumetricResidual { static constexpr bool value = false; }; // by default use flow's aquifer model for now template struct AquiferModel { using type = BlackoilAquiferModel; }; // use flow's linear solver backend for now template struct LinearSolverSplice { using type = TTag::FlowIstlSolver; }; template<> struct LinearSolverBackend { using type = ISTLSolver; }; template struct LinearSolverBackend { using type = ISTLSolver; }; } // namespace Opm::Properties namespace Opm::Parameters { // if openMP is available, set the default the number of threads per process for the main // simulation to 2 (instead of grabbing everything that is available). #if _OPENMP template struct ThreadsPerProcess { static constexpr int value = 2; }; #endif // By default, flowexp accepts the result of the time integration unconditionally if the // smallest time step size is reached. template struct ContinueOnConvergenceError { static constexpr bool value = true; }; template struct EnableTerminalOutput { static constexpr bool value = false; }; // the default for the allowed volumetric error for oil per second template struct NewtonTolerance { using type = GetPropType; static constexpr type value = 1e-1; }; // set the maximum number of Newton iterations to 8 so that we fail quickly (albeit // relatively often) template struct NewtonMaxIterations { static constexpr int value = 8; }; // the maximum volumetric error of a cell in the relaxed region template struct EclNewtonRelaxedTolerance { using type = GetPropType; static constexpr auto baseValue = Parameters::NewtonTolerance::value; static constexpr type value = 1e6 * baseValue; }; // currently, flowexp uses the non-multisegment well model by default to avoid // regressions. the --use-multisegment-well=true|false command line parameter is still // available in flowexp, but hidden from view. template struct UseMultisegmentWell { static constexpr bool value = false; }; // set some properties that are only required by the well model template struct MatrixAddWellContributions { static constexpr bool value = true; }; } // namespace Opm::Parameters namespace Opm { template class FlowExpProblem : public FlowProblem //, public FvBaseProblem { typedef FlowProblem ParentType; using BaseType = ParentType; // GetPropType; public: void writeOutput(bool verbose = true) { OPM_TIMEBLOCK(problemWriteOutput); // use the generic code to prepare the output fields and to // write the desired VTK files. if (Parameters::get() || this->simulator().episodeWillBeOver()) { // \Note: the SimulatorTimer does not carry any useful information, so PRT file (if it gets output) will contain wrong // timing information. BaseType::writeOutput(SimulatorTimer{}, verbose); } } static void registerParameters() { ParentType::registerParameters(); BlackoilModelParameters::registerParameters(); Parameters::registerParam("Do *NOT* use!"); Parameters::hideParam(); Parameters::hideParam(); Parameters::hideParam(); Parameters::hideParam(); Parameters::hideParam(); Parameters::hideParam(); Parameters::hideParam(); Parameters::hideParam(); Parameters::hideParam(); Parameters::hideParam(); Parameters::hideParam(); Parameters::hideParam(); Parameters::hideParam(); Parameters::hideParam(); Parameters::hideParam(); Parameters::hideParam(); Parameters::hideParam(); Parameters::hideParam(); Parameters::hideParam(); Parameters::hideParam(); Parameters::hideParam(); Parameters::hideParam(); Parameters::hideParam(); Parameters::hideParam(); } // inherit the constructors using ParentType::FlowProblem; }; } #endif