added: MultisegmentWellEquations::init

this initializes the equation system.
use the new method in the well implementation.
This commit is contained in:
Arne Morten Kvarving 2022-11-11 21:41:24 +01:00
parent ade6d99289
commit 8fe6b3968e
3 changed files with 98 additions and 64 deletions

View File

@ -22,8 +22,88 @@
#include <config.h>
#include <opm/simulators/wells/MultisegmentWellEquations.hpp>
#include <opm/simulators/wells/MultisegmentWellGeneric.hpp>
namespace Opm {
template<class Scalar, int numWellEq, int numEq>
MultisegmentWellEquations<Scalar,numWellEq,numEq>::
MultisegmentWellEquations(const MultisegmentWellGeneric<Scalar>& well)
: well_(well)
{
}
template<class Scalar, int numWellEq, int numEq>
void MultisegmentWellEquations<Scalar,numWellEq,numEq>::
init(const int num_cells,
const int numPerfs,
const std::vector<int>& cells)
{
duneB_.setBuildMode(OffDiagMatWell::row_wise);
duneC_.setBuildMode(OffDiagMatWell::row_wise);
duneD_.setBuildMode(DiagMatWell::row_wise);
// set the size and patterns for all the matrices and vectors
// [A C^T [x = [ res
// B D] x_well] res_well]
// calculating the NNZ for duneD_
// NNZ = number_of_segments + 2 * (number_of_inlets / number_of_outlets)
{
int nnz_d = well_.numberOfSegments();
for (const std::vector<int>& inlets : well_.segmentInlets()) {
nnz_d += 2 * inlets.size();
}
duneD_.setSize(well_.numberOfSegments(), well_.numberOfSegments(), nnz_d);
}
duneB_.setSize(well_.numberOfSegments(), num_cells, numPerfs);
duneC_.setSize(well_.numberOfSegments(), num_cells, numPerfs);
// we need to add the off diagonal ones
for (auto row = duneD_.createbegin(),
end = duneD_.createend(); row != end; ++row) {
// the number of the row corrspnds to the segment now
const int seg = row.index();
// adding the item related to outlet relation
const Segment& segment = well_.segmentSet()[seg];
const int outlet_segment_number = segment.outletSegment();
if (outlet_segment_number > 0) { // if there is a outlet_segment
const int outlet_segment_index = well_.segmentNumberToIndex(outlet_segment_number);
row.insert(outlet_segment_index);
}
// Add nonzeros for diagonal
row.insert(seg);
// insert the item related to its inlets
for (const int& inlet : well_.segmentInlets()[seg]) {
row.insert(inlet);
}
}
// make the C matrix
for (auto row = duneC_.createbegin(),
end = duneC_.createend(); row != end; ++row) {
// the number of the row corresponds to the segment number now.
for (const int& perf : well_.segmentPerforations()[row.index()]) {
const int cell_idx = cells[perf];
row.insert(cell_idx);
}
}
// make the B^T matrix
for (auto row = duneB_.createbegin(),
end = duneB_.createend(); row != end; ++row) {
// the number of the row corresponds to the segment number now.
for (const int& perf : well_.segmentPerforations()[row.index()]) {
const int cell_idx = cells[perf];
row.insert(cell_idx);
}
}
resWell_.resize(well_.numberOfSegments());
}
template<class Scalar, int numWellEq, int numEq>
void MultisegmentWellEquations<Scalar,numWellEq,numEq>::clear()
{

View File

@ -36,6 +36,8 @@ template<class M> class UMFPack;
namespace Opm
{
template<class Scalar> class MultisegmentWellGeneric;
template<class Scalar, int numWellEq, int numEq>
class MultisegmentWellEquations
{
@ -59,6 +61,17 @@ public:
using OffDiagMatrixBlockWellType = Dune::FieldMatrix<Scalar,numWellEq,numEq>;
using OffDiagMatWell = Dune::BCRSMatrix<OffDiagMatrixBlockWellType>;
MultisegmentWellEquations(const MultisegmentWellGeneric<Scalar>& well);
//! \brief Setup sparsity pattern for the matrices.
//! \param num_cells Total number of cells
//! \param numPerfs Number of perforations
//! \param cells Cell indices for perforations
void init(const int num_cells,
const int numPerfs,
const std::vector<int>& cells);
//! \brief Set all coefficients to 0.
void clear();
// TODO, the following should go to a class for computing purpose
@ -75,6 +88,9 @@ public:
// residuals of the well equations
BVectorWell resWell_;
private:
const MultisegmentWellGeneric<Scalar>& well_; //!< Reference to well
};
}

View File

@ -57,6 +57,7 @@ MultisegmentWellEval<FluidSystem,Indices,Scalar>::
MultisegmentWellEval(WellInterfaceIndices<FluidSystem,Indices,Scalar>& baseif)
: MultisegmentWellGeneric<Scalar>(baseif)
, baseif_(baseif)
, linSys_(*this)
, upwinding_segments_(this->numberOfSegments(), 0)
, segment_densities_(this->numberOfSegments(), 0.0)
, segment_mass_rates_(this->numberOfSegments(), 0.0)
@ -74,70 +75,7 @@ void
MultisegmentWellEval<FluidSystem,Indices,Scalar>::
initMatrixAndVectors(const int num_cells)
{
linSys_.duneB_.setBuildMode(Equations::OffDiagMatWell::row_wise);
linSys_.duneC_.setBuildMode(Equations::OffDiagMatWell::row_wise);
linSys_.duneD_.setBuildMode(Equations::DiagMatWell::row_wise);
// set the size and patterns for all the matrices and vectors
// [A C^T [x = [ res
// B D] x_well] res_well]
// calculatiing the NNZ for duneD_
// NNZ = number_of_segments + 2 * (number_of_inlets / number_of_outlets)
{
int nnz_d = this->numberOfSegments();
for (const std::vector<int>& inlets : this->segment_inlets_) {
nnz_d += 2 * inlets.size();
}
linSys_.duneD_.setSize(this->numberOfSegments(), this->numberOfSegments(), nnz_d);
}
linSys_.duneB_.setSize(this->numberOfSegments(), num_cells, baseif_.numPerfs());
linSys_.duneC_.setSize(this->numberOfSegments(), num_cells, baseif_.numPerfs());
// we need to add the off diagonal ones
for (auto row = linSys_.duneD_.createbegin(),
end = linSys_.duneD_.createend(); row != end; ++row) {
// the number of the row corrspnds to the segment now
const int seg = row.index();
// adding the item related to outlet relation
const Segment& segment = this->segmentSet()[seg];
const int outlet_segment_number = segment.outletSegment();
if (outlet_segment_number > 0) { // if there is a outlet_segment
const int outlet_segment_index = this->segmentNumberToIndex(outlet_segment_number);
row.insert(outlet_segment_index);
}
// Add nonzeros for diagonal
row.insert(seg);
// insert the item related to its inlets
for (const int& inlet : this->segment_inlets_[seg]) {
row.insert(inlet);
}
}
// make the C matrix
for (auto row = linSys_.duneC_.createbegin(),
end = linSys_.duneC_.createend(); row != end; ++row) {
// the number of the row corresponds to the segment number now.
for (const int& perf : this->segment_perforations_[row.index()]) {
const int cell_idx = baseif_.cells()[perf];
row.insert(cell_idx);
}
}
// make the B^T matrix
for (auto row = linSys_.duneB_.createbegin(),
end = linSys_.duneB_.createend(); row != end; ++row) {
// the number of the row corresponds to the segment number now.
for (const int& perf : this->segment_perforations_[row.index()]) {
const int cell_idx = baseif_.cells()[perf];
row.insert(cell_idx);
}
}
linSys_.resWell_.resize(this->numberOfSegments());
linSys_.init(num_cells, baseif_.numPerfs(), baseif_.cells());
primary_variables_.resize(this->numberOfSegments());
primary_variables_evaluation_.resize(this->numberOfSegments());
}