mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
added: MultisegmentWellEquations::clear
this zeros the equation system. use the new method in the well implementation.
This commit is contained in:
parent
57f09050fc
commit
ac245a2e17
@ -22,8 +22,29 @@
|
||||
#include <config.h>
|
||||
#include <opm/simulators/wells/MultisegmentWellEquations.hpp>
|
||||
|
||||
namespace Opm {
|
||||
|
||||
namespace Opm
|
||||
template<class Scalar, int numWellEq, int numEq>
|
||||
void MultisegmentWellEquations<Scalar,numWellEq,numEq>::clear()
|
||||
{
|
||||
duneB_ = 0.0;
|
||||
duneC_ = 0.0;
|
||||
duneD_ = 0.0;
|
||||
resWell_ = 0.0;
|
||||
duneDSolver_.reset();
|
||||
}
|
||||
|
||||
#define INSTANCE(numWellEq, numEq) \
|
||||
template class MultisegmentWellEquations<double,numWellEq,numEq>;
|
||||
|
||||
INSTANCE(2,1)
|
||||
INSTANCE(2,2)
|
||||
INSTANCE(2,6)
|
||||
INSTANCE(3,2)
|
||||
INSTANCE(3,3)
|
||||
INSTANCE(3,4)
|
||||
INSTANCE(4,3)
|
||||
INSTANCE(4,4)
|
||||
INSTANCE(4,5)
|
||||
|
||||
}
|
||||
|
@ -59,6 +59,8 @@ public:
|
||||
using OffDiagMatrixBlockWellType = Dune::FieldMatrix<Scalar,numWellEq,numEq>;
|
||||
using OffDiagMatWell = Dune::BCRSMatrix<OffDiagMatrixBlockWellType>;
|
||||
|
||||
void clear();
|
||||
|
||||
// TODO, the following should go to a class for computing purpose
|
||||
// two off-diagonal matrices
|
||||
OffDiagMatWell duneB_;
|
||||
|
@ -1634,13 +1634,7 @@ namespace Opm
|
||||
computeSegmentFluidProperties(ebosSimulator, deferred_logger);
|
||||
|
||||
// clear all entries
|
||||
this->linSys_.duneB_ = 0.0;
|
||||
this->linSys_.duneC_ = 0.0;
|
||||
|
||||
this->linSys_.duneD_ = 0.0;
|
||||
this->linSys_.resWell_ = 0.0;
|
||||
|
||||
this->linSys_.duneDSolver_.reset();
|
||||
this->linSys_.clear();
|
||||
|
||||
auto& ws = well_state.well(this->index_of_well_);
|
||||
ws.dissolved_gas_rate = 0;
|
||||
|
Loading…
Reference in New Issue
Block a user