Files
opm-simulators/opm/polymer/fullyimplicit/SimulatorFullyImplicitCompressiblePolymer.hpp
Andreas Lauser ef2a560fb3 flow_ebos: print statistics about failed time steps
the performance summary at the end of a Norne run which are printed by
`flow_ebos` now looks like this on my machine:

```
Total time (seconds):         773.757
Solver time (seconds):        753.349
 Assembly time (seconds):     377.218 (Failed: 23.537; 6.23965%)
 Linear solve time (seconds): 352.022 (Failed: 23.2757; 6.61201%)
 Update time (seconds):       16.3658 (Failed: 1.13149; 6.91375%)
 Output write time (seconds): 22.5991
Overall Well Iterations:      870 (Failed: 35; 4.02299%)
Overall Linearizations:       2098 (Failed: 136; 6.48236%)
Overall Newton Iterations:    1756 (Failed: 136; 7.74487%)
Overall Linear Iterations:    26572 (Failed: 1786; 6.72136%)
```

for the flow_legacy family, nothing changes.
2017-04-11 11:12:11 +02:00

150 lines
6.0 KiB
C++

/*
Copyright 2014 SINTEF ICT, Applied Mathematics.
Copyright 2014 STATOIL ASA.
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 <http://www.gnu.org/licenses/>.
*/
#ifndef OPM_SIMULATORFULLYIMPLICITCOMPRESSIBLEPOLYMER_HEADER_INCLUDED
#define OPM_SIMULATORFULLYIMPLICITCOMPRESSIBLEPOLYMER_HEADER_INCLUDED
#include <opm/core/simulator/SimulatorReport.hpp>
#include <opm/core/utility/parameters/ParameterGroup.hpp>
#include <opm/common/ErrorMacros.hpp>
#include <opm/autodiff/GeoProps.hpp>
#include <opm/autodiff/BlackoilPropsAdFromDeck.hpp>
#include <opm/autodiff/WellStateFullyImplicitBlackoil.hpp>
#include <opm/autodiff/SimulatorBase.hpp>
#include <opm/polymer/fullyimplicit/FullyImplicitCompressiblePolymerSolver.hpp>
#include <opm/polymer/fullyimplicit/BlackoilPolymerModel.hpp>
#include <opm/core/grid.h>
#include <opm/core/wells.h>
#include <opm/core/pressure/flow_bc.h>
#include <opm/core/simulator/SimulatorReport.hpp>
#include <opm/simulators/timestepping/SimulatorTimer.hpp>
#include <opm/core/utility/StopWatch.hpp>
#include <opm/output/eclipse/EclipseIO.hpp>
#include <opm/core/utility/miscUtilities.hpp>
#include <opm/core/utility/miscUtilitiesBlackoil.hpp>
#include <opm/core/wells/WellsManager.hpp>
#include <opm/core/props/rock/RockCompressibility.hpp>
#include <opm/core/grid/ColumnExtract.hpp>
#include <opm/polymer/PolymerBlackoilState.hpp>
#include <opm/polymer/PolymerInflow.hpp>
#include <opm/core/simulator/WellState.hpp>
#include <opm/core/transport/reorder/TransportSolverCompressibleTwophaseReorder.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/Schedule.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/ScheduleEnums.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/Well.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/WellProductionProperties.hpp>
#include <opm/parser/eclipse/Deck/Deck.hpp>
#include <boost/filesystem.hpp>
#include <boost/scoped_ptr.hpp>
#include <boost/lexical_cast.hpp>
#include <numeric>
#include <fstream>
#include <iostream>
namespace Opm
{
template <class GridT>
class SimulatorFullyImplicitCompressiblePolymer;
class StandardWells;
template <class GridT>
struct SimulatorTraits<SimulatorFullyImplicitCompressiblePolymer<GridT> >
{
typedef PolymerBlackoilState ReservoirState;
typedef WellStateFullyImplicitBlackoilPolymer WellState;
typedef BlackoilOutputWriter OutputWriter;
typedef GridT Grid;
typedef FullyImplicitCompressiblePolymerSolver Solver;
typedef StandardWells WellModel;
/// Dummy class, this Solver does not use a Model.
struct Model
{
typedef parameter::ParameterGroup ModelParameters;
};
};
/// Class collecting all necessary components for a two-phase simulation.
template <class GridT>
class SimulatorFullyImplicitCompressiblePolymer
: public SimulatorBase<SimulatorFullyImplicitCompressiblePolymer<GridT> >
{
typedef SimulatorFullyImplicitCompressiblePolymer ThisType;
typedef SimulatorBase<ThisType> BaseType;
typedef typename BaseType::Solver Solver;
typedef typename BaseType::WellModel WellModel;
public:
/// Initialise from parameters and objects to observe.
SimulatorFullyImplicitCompressiblePolymer(const parameter::ParameterGroup& param,
const GridT& grid,
DerivedGeology& geo,
BlackoilPropsAdFromDeck& props,
const PolymerPropsAd& polymer_props,
const RockCompressibility* rock_comp_props,
std::shared_ptr<EclipseState> eclipse_state,
BlackoilOutputWriter& output_writer,
const Deck& deck,
NewtonIterationBlackoilInterface& linsolver,
const double* gravity);
std::unique_ptr<Solver> createSolver(const WellModel& well_model);
void handleAdditionalWellInflow(SimulatorTimer& timer,
WellsManager& wells_manager,
typename BaseType::WellState& well_state,
const Wells* wells);
void updateListEconLimited(const std::unique_ptr<Solver>& solver,
const Schedule& schedule,
const int current_step,
const Wells* wells,
const WellState& well_state,
DynamicListEconLimited& list_econ_limited) const;
/// return the statistics if the nonlinearIteration() method failed.
///
/// NOTE: for the flow_legacy simulator family this method is a stub, i.e. the
/// failure report object will *not* contain any meaningful data.
const SimulatorReport& failureReport() const
{ return failureReport_; }
private:
SimulatorReport failureReport_;
const Deck& deck_;
const PolymerPropsAd& polymer_props_;
};
} // namespace Opm
#include "SimulatorFullyImplicitCompressiblePolymer_impl.hpp"
#endif // OPM_SIMULATORFULLYIMPLICITCOMPRESSIBLEPOLYMER_HEADER_INCLUDED