MultisegmentWellContribution: template Scalar type

This commit is contained in:
Arne Morten Kvarving 2024-04-15 16:57:54 +02:00
parent ba1c6db855
commit 216f0bea0d
4 changed files with 59 additions and 42 deletions

View File

@ -29,21 +29,27 @@
namespace Opm
{
MultisegmentWellContribution::MultisegmentWellContribution(unsigned int dim_, unsigned int dim_wells_,
unsigned int Mb_,
std::vector<double> &Bvalues, std::vector<unsigned int> &BcolIndices, std::vector<unsigned int> &BrowPointers,
unsigned int DnumBlocks_, double *Dvalues, UMFPackIndex *DcolPointers, UMFPackIndex *DrowIndices,
std::vector<double> &Cvalues)
:
dim(dim_), // size of blockvectors in vectors x and y, equal to MultisegmentWell::numEq
dim_wells(dim_wells_), // size of blocks in C, B and D, equal to MultisegmentWell::numWellEq
M(Mb_ * dim_wells), // number of rows, M == dim_wells*Mb
Mb(Mb_), // number of blockrows in C, D and B
DnumBlocks(DnumBlocks_), // number of blocks in D
template<class Scalar>
MultisegmentWellContribution<Scalar>::
MultisegmentWellContribution(unsigned int dim_, unsigned int dim_wells_,
unsigned int Mb_,
std::vector<Scalar>& Bvalues,
std::vector<unsigned int>& BcolIndices,
std::vector<unsigned int>& BrowPointers,
unsigned int DnumBlocks_,
Scalar* Dvalues,
UMFPackIndex* DcolPointers,
UMFPackIndex* DrowIndices,
std::vector<Scalar>& Cvalues)
: dim(dim_) // size of blockvectors in vectors x and y, equal to MultisegmentWell::numEq
, dim_wells(dim_wells_) // size of blocks in C, B and D, equal to MultisegmentWell::numWellEq
, M(Mb_ * dim_wells) // number of rows, M == dim_wells*Mb
, Mb(Mb_) // number of blockrows in C, D and B
, DnumBlocks(DnumBlocks_) // number of blocks in D
// copy data for matrix D into vectors to prevent it going out of scope
Dvals(Dvalues, Dvalues + DnumBlocks * dim_wells * dim_wells),
Dcols(DcolPointers, DcolPointers + M + 1),
Drows(DrowIndices, DrowIndices + DnumBlocks * dim_wells * dim_wells)
, Dvals(Dvalues, Dvalues + DnumBlocks * dim_wells * dim_wells)
, Dcols(DcolPointers, DcolPointers + M + 1)
, Drows(DrowIndices, DrowIndices + DnumBlocks * dim_wells * dim_wells)
{
Cvals = std::move(Cvalues);
Bvals = std::move(Bvalues);
@ -57,17 +63,18 @@ MultisegmentWellContribution::MultisegmentWellContribution(unsigned int dim_, un
umfpack_di_numeric(Dcols.data(), Drows.data(), Dvals.data(), UMFPACK_Symbolic, &UMFPACK_Numeric, nullptr, nullptr);
}
MultisegmentWellContribution::~MultisegmentWellContribution()
template<class Scalar>
MultisegmentWellContribution<Scalar>::~MultisegmentWellContribution()
{
umfpack_di_free_symbolic(&UMFPACK_Symbolic);
umfpack_di_free_numeric(&UMFPACK_Numeric);
}
// Apply the MultisegmentWellContribution, similar to MultisegmentWell::apply()
// h_x and h_y reside on host
// y -= (C^T * (D^-1 * (B * x)))
void MultisegmentWellContribution::apply(double *h_x, double *h_y)
template<class Scalar>
void MultisegmentWellContribution<Scalar>::apply(Scalar* h_x, Scalar* h_y)
{
OPM_TIMEBLOCK(apply);
// reset z1 and z2
@ -80,7 +87,7 @@ void MultisegmentWellContribution::apply(double *h_x, double *h_y)
for (unsigned int blockID = Brows[row]; blockID < Brows[row + 1]; ++blockID) {
unsigned int colIdx = Bcols[blockID];
for (unsigned int j = 0; j < dim_wells; ++j) {
double temp = 0.0;
Scalar temp = 0.0;
for (unsigned int k = 0; k < dim; ++k) {
temp += Bvals[blockID * dim * dim_wells + j * dim + k] * h_x[colIdx * dim + k];
}
@ -100,7 +107,7 @@ void MultisegmentWellContribution::apply(double *h_x, double *h_y)
for (unsigned int blockID = Brows[row]; blockID < Brows[row + 1]; ++blockID) {
unsigned int colIdx = Bcols[blockID];
for (unsigned int j = 0; j < dim; ++j) {
double temp = 0.0;
Scalar temp = 0.0;
for (unsigned int k = 0; k < dim_wells; ++k) {
temp += Cvals[blockID * dim * dim_wells + j + k * dim] * z2[row * dim_wells + k];
}
@ -111,11 +118,14 @@ void MultisegmentWellContribution::apply(double *h_x, double *h_y)
}
#if HAVE_CUDA
void MultisegmentWellContribution::setCudaStream(cudaStream_t stream_)
template<class Scalar>
void MultisegmentWellContribution<Scalar>::setCudaStream(cudaStream_t stream_)
{
stream = stream_;
}
#endif
template class MultisegmentWellContribution<double>;
} //namespace Opm

View File

@ -41,6 +41,7 @@ namespace Opm
/// B*x and D*B*x are a vector with M*numWellEq doubles
/// C*D*B*x is a vector with N*numEq doubles.
template<class Scalar>
class MultisegmentWellContribution
{
@ -57,15 +58,15 @@ private:
// C and B are stored in BCRS format, D is stored in CSC format (Dune::UMFPack)
// Sparsity pattern for C is not stored, since it is the same as B
unsigned int DnumBlocks; // number of blocks in D
std::vector<double> Cvals;
std::vector<double> Dvals;
std::vector<double> Bvals;
std::vector<Scalar> Cvals;
std::vector<Scalar> Dvals;
std::vector<Scalar> Bvals;
std::vector<int> Dcols; // Columnpointers, contains M+1 entries
std::vector<unsigned int> Bcols;
std::vector<int> Drows; // Rowindicies, contains DnumBlocks*dim*dim_wells entries
std::vector<unsigned int> Brows;
std::vector<double> z1; // z1 = B * x
std::vector<double> z2; // z2 = D^-1 * B * x
std::vector<Scalar> z1; // z1 = B * x
std::vector<Scalar> z2; // z2 = D^-1 * B * x
void *UMFPACK_Symbolic, *UMFPACK_Numeric;
/// Translate the columnIndex if needed
@ -97,9 +98,14 @@ public:
/// \param[in] Cvalues nonzero values of matrix C
MultisegmentWellContribution(unsigned int dim, unsigned int dim_wells,
unsigned int Mb,
std::vector<double> &Bvalues, std::vector<unsigned int> &BcolIndices, std::vector<unsigned int> &BrowPointers,
unsigned int DnumBlocks, double *Dvalues, UMFPackIndex *DcolPointers,
UMFPackIndex *DrowIndices, std::vector<double> &Cvalues);
std::vector<Scalar>& Bvalues,
std::vector<unsigned int>& BcolIndices,
std::vector<unsigned int>& BrowPointers,
unsigned int DnumBlocks,
Scalar* Dvalues,
UMFPackIndex* DcolPointers,
UMFPackIndex* DrowIndices,
std::vector<Scalar>& Cvalues);
/// Destroy a MultisegmentWellContribution, and free memory
~MultisegmentWellContribution();
@ -108,7 +114,7 @@ public:
/// performs y -= (C^T * (D^-1 * (B*x))) for MultisegmentWell
/// \param[in] h_x vector x, must be on CPU
/// \param[inout] h_y vector y, must be on CPU
void apply(double *h_x, double *h_y);
void apply(Scalar* h_x, Scalar* h_y);
};
} //namespace Opm

View File

@ -157,17 +157,18 @@ void WellContributions::addMultisegmentWellContribution(unsigned int dim_,
std::vector<double>& Cvalues)
{
assert(dim==dim_);
multisegments.push_back(std::make_unique<MultisegmentWellContribution>(dim_,
dim_wells_,
Mb,
Bvalues,
BcolIndices,
BrowPointers,
DnumBlocks,
Dvalues,
DcolPointers,
DrowIndices,
Cvalues));
using MSW = MultisegmentWellContribution<double>;
multisegments.push_back(std::make_unique<MSW>(dim_,
dim_wells_,
Mb,
Bvalues,
BcolIndices,
BrowPointers,
DnumBlocks,
Dvalues,
DcolPointers,
DrowIndices,
Cvalues));
++num_ms_wells;
}

View File

@ -30,7 +30,7 @@
namespace Opm {
class MultisegmentWellContribution;
template<class Scalar> class MultisegmentWellContribution;
/// This class serves to eliminate the need to include the WellContributions into the matrix (with --matrix-add-well-contributions=true) for the cusparseSolver or openclSolver.
/// If the --matrix-add-well-contributions commandline parameter is true, this class should still be used, but be empty.
@ -74,7 +74,7 @@ protected:
unsigned int num_std_wells_so_far = 0; // keep track of where next data is written
std::vector<unsigned int> val_pointers; // val_pointers[wellID] == index of first block for this well in Ccols and Bcols
std::vector<std::unique_ptr<MultisegmentWellContribution>> multisegments;
std::vector<std::unique_ptr<MultisegmentWellContribution<double>>> multisegments;
public:
unsigned int getNumWells(){