WellOperators.hpp
Go to the documentation of this file.
1/*
2 Copyright 2016 IRIS AS
3 Copyright 2019, 2020 Equinor ASA
4 Copyright 2020 SINTEF
5
6 This file is part of the Open Porous Media project (OPM).
7
8 OPM is free software: you can redistribute it and/or modify
9 it under the terms of the GNU General Public License as published by
10 the Free Software Foundation, either version 3 of the License, or
11 (at your option) any later version.
12
13 OPM is distributed in the hope that it will be useful,
14 but WITHOUT ANY WARRANTY; without even the implied warranty of
15 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16 GNU General Public License for more details.
17
18 You should have received a copy of the GNU General Public License
19 along with OPM. If not, see <http://www.gnu.org/licenses/>.
20*/
21
22#ifndef OPM_WELLOPERATORS_HEADER_INCLUDED
23#define OPM_WELLOPERATORS_HEADER_INCLUDED
24
25#include <dune/common/parallel/communication.hh>
26#include <dune/istl/operators.hh>
27#include <dune/istl/bcrsmatrix.hh>
28
29#include <opm/common/TimingMacros.hpp>
30
32#include <dune/common/shared_ptr.hh>
33#include <dune/istl/paamg/smoother.hh>
34
35#include <cstddef>
36
37namespace Opm {
38
39//=====================================================================
40// Implementation for ISTL-matrix based operators
41// Note: the classes WellModelMatrixAdapter and
42// WellModelGhostLastMatrixAdapter were moved from ISTLSolver.hpp
43// and subsequently modified.
44//=====================================================================
45
55template <class X, class Y>
56class LinearOperatorExtra : public Dune::LinearOperator<X, Y>
57{
58public:
59 using field_type = typename X::field_type;
60 using PressureMatrix = Dune::BCRSMatrix<MatrixBlock<field_type, 1, 1>>;
62 const X& weights,
63 const bool use_well_weights) const = 0;
64 virtual void addWellPressureEquationsStruct(PressureMatrix& jacobian) const = 0;
65 virtual int getNumberOfExtraEquations() const = 0;
66};
67
68template <class WellModel, class X, class Y>
70{
71public:
73 using field_type = typename Base::field_type;
75 explicit WellModelAsLinearOperator(const WellModel& wm)
76 : wellMod_(wm)
77 {
78 }
79
84 void apply(const X& x, Y& y) const override
85 {
86 OPM_TIMEBLOCK(apply);
87 for (const auto& well : this->wellMod_) {
88 this->applySingleWell(x, y, well, well->cells());
89 }
90 }
91
93 void applyscaleadd(field_type alpha, const X& x, Y& y) const override
94 {
95 OPM_TIMEBLOCK(applyscaleadd);
96 if (this->wellMod_.empty()) {
97 return;
98 }
99
100 if (scaleAddRes_.size() != y.size()) {
101 scaleAddRes_.resize(y.size());
102 }
103
104 scaleAddRes_ = 0.0;
105 // scaleAddRes_ = - C D^-1 B x
107 // Ax = Ax + alpha * scaleAddRes_
108 y.axpy(alpha, scaleAddRes_);
109 }
110
116 Dune::SolverCategory::Category category() const override
117 {
118 return Dune::SolverCategory::sequential;
119 }
120
122 const X& weights,
123 const bool use_well_weights) const override
124 {
125 OPM_TIMEBLOCK(addWellPressureEquations);
126 wellMod_.addWellPressureEquations(jacobian, weights, use_well_weights);
127 }
128
129 void addWellPressureEquationsStruct(PressureMatrix& jacobian) const override
130 {
131 OPM_TIMEBLOCK(addWellPressureEquationsStruct);
132 wellMod_.addWellPressureEquationsStruct(jacobian);
133 }
134
135 int getNumberOfExtraEquations() const override
136 {
137 return wellMod_.numLocalWellsEnd();
138 }
139
140protected:
141 const WellModel& wellMod_;
142
143 template<class WellType, class ArrayType>
144 void applySingleWell(const X& x, Y& y,
145 const WellType& well,
146 const ArrayType& cells) const
147 {
148 // Well equations B and C uses only the perforated cells, so need to apply on local vectors
149 x_local_.resize(cells.size());
150 Ax_local_.resize(cells.size());
151
152 for (size_t i = 0; i < cells.size(); ++i) {
153 x_local_[i] = x[cells[i]];
154 Ax_local_[i] = y[cells[i]];
155 }
156
157 well->apply(x_local_, Ax_local_);
158
159 for (size_t i = 0; i < cells.size(); ++i) {
160 // only need to update Ax
161 y[cells[i]] = Ax_local_[i];
162 }
163
164 }
165
166 // These members are used to avoid reallocation.
167 // Their state is not relevant between function calls, so they can
168 // (and must) be mutable, as the functions using them are const.
169 mutable X x_local_{};
170 mutable Y Ax_local_{};
171 mutable Y scaleAddRes_{};
172};
173
174template <class WellModel, class X, class Y>
176{
177public:
179 using WBase::WBase; // inherit all constructors from the base class
182
183 void setDomainIndex(int index) { domainIndex_ = index; }
184
185 void apply(const X& x, Y& y) const override
186 {
187 OPM_TIMEBLOCK(apply);
188 std::size_t well_index = 0;
189 for (const auto& well : this->wellMod_) {
190 if (this->wellMod_.well_domain().at(well->name()) == domainIndex_) {
191 this->applySingleWell(x, y, well,
192 this->wellMod_.well_local_cells()[well_index]);
193 }
194 ++well_index;
195 }
196 }
197
199 const X& weights,
200 const bool use_well_weights) const override
201 {
202 OPM_TIMEBLOCK(addWellPressureEquations);
203 this->wellMod_.addWellPressureEquationsDomain(jacobian,
204 weights,
205 use_well_weights,
206 domainIndex_);
207 }
208
209private:
210 int domainIndex_ = -1;
211};
212
223template<class M, class X, class Y>
224class WellModelMatrixAdapter : public Dune::AssembledLinearOperator<M,X,Y>
225{
226public:
227 using matrix_type = M;
228 using domain_type = X;
229 using range_type = Y;
230 using field_type = typename X::field_type;
231 using PressureMatrix = Dune::BCRSMatrix<MatrixBlock<field_type, 1, 1>>;
232
233 Dune::SolverCategory::Category category() const override
234 {
235 return Dune::SolverCategory::sequential;
236 }
237
240 const LinearOperatorExtra<X, Y>& wellOper)
241 : A_( A ), wellOper_( wellOper )
242 {}
243
244 void apply( const X& x, Y& y ) const override
245 {
246 OPM_TIMEBLOCK(apply);
247 A_.mv(x, y);
248
249 // add well model modification to y
250 wellOper_.apply(x, y);
251 }
252
253 // y += \alpha * A * x
254 void applyscaleadd (field_type alpha, const X& x, Y& y) const override
255 {
256 OPM_TIMEBLOCK(applyscaleadd);
257 A_.usmv(alpha, x, y);
258
259 // add scaled well model modification to y
260 wellOper_.applyscaleadd(alpha, x, y);
261 }
262
263 const matrix_type& getmat() const override { return A_; }
264
266 const X& weights,
267 const bool use_well_weights) const
268 {
269 OPM_TIMEBLOCK(addWellPressureEquations);
270 wellOper_.addWellPressureEquations(jacobian, weights, use_well_weights);
271 }
272
274 {
275 OPM_TIMEBLOCK(addWellPressureEquations);
276 wellOper_.addWellPressureEquationsStruct(jacobian);
277 }
278
280 {
281 return wellOper_.getNumberOfExtraEquations();
282 }
283
284protected:
287};
288
297template<class M, class X, class Y, bool overlapping >
298class WellModelGhostLastMatrixAdapter : public Dune::AssembledLinearOperator<M,X,Y>
299{
300public:
301 using matrix_type = M;
302 using domain_type = X;
303 using range_type = Y;
304 using field_type = typename X::field_type;
305 using PressureMatrix = Dune::BCRSMatrix<MatrixBlock<field_type, 1, 1>>;
306#if HAVE_MPI
307 using communication_type = Dune::OwnerOverlapCopyCommunication<int,int>;
308#else
309 using communication_type = Dune::Communication<int>;
310#endif
311
312 Dune::SolverCategory::Category category() const override
313 {
314 return overlapping ?
315 Dune::SolverCategory::overlapping : Dune::SolverCategory::sequential;
316 }
317
320 const LinearOperatorExtra<X, Y>& wellOper,
321 const std::size_t interiorSize )
322 : A_( A ), wellOper_( wellOper ), interiorSize_(interiorSize)
323 {}
324
325 void apply(const X& x, Y& y) const override
326 {
327 OPM_TIMEBLOCK(apply);
328 for (auto row = A_.begin(); row.index() < interiorSize_; ++row)
329 {
330 y[row.index()]=0;
331 auto endc = (*row).end();
332 for (auto col = (*row).begin(); col != endc; ++col)
333 (*col).umv(x[col.index()], y[row.index()]);
334 }
335
336 // add well model modification to y
337 wellOper_.apply(x, y);
338
340 }
341
342 // y += \alpha * A * x
343 void applyscaleadd (field_type alpha, const X& x, Y& y) const override
344 {
345 OPM_TIMEBLOCK(applyscaleadd);
346 for (auto row = A_.begin(); row.index() < interiorSize_; ++row)
347 {
348 auto endc = (*row).end();
349 for (auto col = (*row).begin(); col != endc; ++col)
350 (*col).usmv(alpha, x[col.index()], y[row.index()]);
351 }
352 // add scaled well model modification to y
353 wellOper_.applyscaleadd(alpha, x, y);
354
356 }
357
358 const matrix_type& getmat() const override { return A_; }
359
361 const X& weights,
362 const bool use_well_weights) const
363 {
364 OPM_TIMEBLOCK(addWellPressureEquations);
365 wellOper_.addWellPressureEquations(jacobian, weights, use_well_weights);
366 }
367
369 {
370 OPM_TIMEBLOCK(addWellPressureEquationsStruct);
371 wellOper_.addWellPressureEquationsStruct(jacobian);
372 }
373
375 {
376 return wellOper_.getNumberOfExtraEquations();
377 }
378
379protected:
380 void ghostLastProject(Y& y) const
381 {
382 std::size_t end = y.size();
383 for (std::size_t i = interiorSize_; i < end; ++i)
384 y[i] = 0;
385 }
386
389 std::size_t interiorSize_;
390};
391
400template<class M, class X, class Y, class C>
401class GhostLastMatrixAdapter : public Dune::AssembledLinearOperator<M,X,Y>
402{
403public:
404 typedef M matrix_type;
405 typedef X domain_type;
406 typedef Y range_type;
407 typedef typename X::field_type field_type;
408
409
411
412 Dune::SolverCategory::Category category() const override
413 {
414 return Dune::SolverCategory::overlapping;
415 }
416
419 const communication_type& comm)
420 : A_( Dune::stackobject_to_shared_ptr(A) ), comm_(comm)
421 {
422 interiorSize_ = setInteriorSize(comm_);
423 }
424
425 GhostLastMatrixAdapter (const std::shared_ptr<M> A,
426 const communication_type& comm)
427 : A_( A ), comm_(comm)
428 {
429 interiorSize_ = setInteriorSize(comm_);
430 }
431
432 virtual void apply( const X& x, Y& y ) const override
433 {
434 for (auto row = A_->begin(); row.index() < interiorSize_; ++row)
435 {
436 y[row.index()]=0;
437 auto endc = (*row).end();
438 for (auto col = (*row).begin(); col != endc; ++col)
439 (*col).umv(x[col.index()], y[row.index()]);
440 }
441
442 ghostLastProject( y );
443 }
444
445 // y += \alpha * A * x
446 virtual void applyscaleadd (field_type alpha, const X& x, Y& y) const override
447 {
448 for (auto row = A_->begin(); row.index() < interiorSize_; ++row)
449 {
450 auto endc = (*row).end();
451 for (auto col = (*row).begin(); col != endc; ++col)
452 (*col).usmv(alpha, x[col.index()], y[row.index()]);
453 }
454
455 ghostLastProject( y );
456 }
457
458 virtual const matrix_type& getmat() const override { return *A_; }
459
460 size_t getInteriorSize() const { return interiorSize_;}
461
462private:
463 void ghostLastProject(Y& y) const
464 {
465 size_t end = y.size();
466 for (size_t i = interiorSize_; i < end; ++i)
467 y[i] = 0; //project to interiorsize, i.e. ignore ghost
468 }
469
470 size_t setInteriorSize(const communication_type& comm) const
471 {
472 auto indexSet = comm.indexSet();
473
474 size_t is = 0;
475 // Loop over index set
476 for (auto idx = indexSet.begin(); idx!=indexSet.end(); ++idx) {
477 //Only take "owner" indices
478 if (idx->local().attribute()==1) {
479 //get local index
480 auto loc = idx->local().local();
481 // if loc is higher than "old interior size", update it
482 if (loc > is) {
483 is = loc;
484 }
485 }
486 }
487 return is + 1; //size is plus 1 since we start at 0
488 }
489 const std::shared_ptr<const matrix_type> A_ ;
490 const communication_type& comm_;
491 size_t interiorSize_;
492};
493
494} // namespace Opm
495
496namespace Dune {
497namespace Amg {
498
499template<class M, class X, class Y, class C>
500class ConstructionTraits<Opm::GhostLastMatrixAdapter<M,X,Y,C> >
501{
502public:
503 typedef ParallelOperatorArgs<M,C> Arguments;
504
505 static inline std::shared_ptr<Opm::GhostLastMatrixAdapter<M,X,Y,C>> construct(const Arguments& args)
506 {
507 return std::make_shared<Opm::GhostLastMatrixAdapter<M,X,Y,C>>
508 (args.matrix_, args.comm_);
509 }
510};
511
512} // end namespace Amg
513} // end namespace Dune
514
515#endif // OPM_WELLOPERATORS_HEADER_INCLUDED
516
static std::shared_ptr< Opm::GhostLastMatrixAdapter< M, X, Y, C > > construct(const Arguments &args)
Definition: WellOperators.hpp:505
ParallelOperatorArgs< M, C > Arguments
Definition: WellOperators.hpp:503
Definition: WellOperators.hpp:176
void addWellPressureEquations(PressureMatrix &jacobian, const X &weights, const bool use_well_weights) const override
Definition: WellOperators.hpp:198
void setDomainIndex(int index)
Definition: WellOperators.hpp:183
void apply(const X &x, Y &y) const override
Definition: WellOperators.hpp:185
Dune linear operator that assumes ghost rows are ordered after interior rows. Avoids some computation...
Definition: WellOperators.hpp:402
C communication_type
Definition: WellOperators.hpp:410
M matrix_type
Definition: WellOperators.hpp:404
Dune::SolverCategory::Category category() const override
Definition: WellOperators.hpp:412
X::field_type field_type
Definition: WellOperators.hpp:407
X domain_type
Definition: WellOperators.hpp:405
virtual const matrix_type & getmat() const override
Definition: WellOperators.hpp:458
GhostLastMatrixAdapter(const std::shared_ptr< M > A, const communication_type &comm)
Definition: WellOperators.hpp:425
GhostLastMatrixAdapter(const M &A, const communication_type &comm)
constructor: just store a reference to a matrix
Definition: WellOperators.hpp:418
virtual void applyscaleadd(field_type alpha, const X &x, Y &y) const override
Definition: WellOperators.hpp:446
Y range_type
Definition: WellOperators.hpp:406
virtual void apply(const X &x, Y &y) const override
Definition: WellOperators.hpp:432
size_t getInteriorSize() const
Definition: WellOperators.hpp:460
Definition: WellOperators.hpp:57
virtual void addWellPressureEquationsStruct(PressureMatrix &jacobian) const =0
virtual void addWellPressureEquations(PressureMatrix &jacobian, const X &weights, const bool use_well_weights) const =0
Dune::BCRSMatrix< MatrixBlock< field_type, 1, 1 > > PressureMatrix
Definition: WellOperators.hpp:60
typename X::field_type field_type
Definition: WellOperators.hpp:59
virtual int getNumberOfExtraEquations() const =0
Definition: WellOperators.hpp:70
typename Base::PressureMatrix PressureMatrix
Definition: WellOperators.hpp:74
const WellModel & wellMod_
Definition: WellOperators.hpp:141
Y Ax_local_
Definition: WellOperators.hpp:170
void apply(const X &x, Y &y) const override
apply operator to x: The input vector is consistent and the output must also be consistent on the in...
Definition: WellOperators.hpp:84
WellModelAsLinearOperator(const WellModel &wm)
Definition: WellOperators.hpp:75
X x_local_
Definition: WellOperators.hpp:169
typename Base::field_type field_type
Definition: WellOperators.hpp:73
void addWellPressureEquations(PressureMatrix &jacobian, const X &weights, const bool use_well_weights) const override
Definition: WellOperators.hpp:121
void applySingleWell(const X &x, Y &y, const WellType &well, const ArrayType &cells) const
Definition: WellOperators.hpp:144
void applyscaleadd(field_type alpha, const X &x, Y &y) const override
apply operator to x, scale and add:
Definition: WellOperators.hpp:93
Dune::SolverCategory::Category category() const override
Definition: WellOperators.hpp:116
int getNumberOfExtraEquations() const override
Definition: WellOperators.hpp:135
void addWellPressureEquationsStruct(PressureMatrix &jacobian) const override
Definition: WellOperators.hpp:129
Y scaleAddRes_
Definition: WellOperators.hpp:171
Adapter to combine a matrix and another linear operator into a combined linear operator.
Definition: WellOperators.hpp:299
void addWellPressureEquationsStruct(PressureMatrix &jacobian) const
Definition: WellOperators.hpp:368
typename X::field_type field_type
Definition: WellOperators.hpp:304
void addWellPressureEquations(PressureMatrix &jacobian, const X &weights, const bool use_well_weights) const
Definition: WellOperators.hpp:360
const matrix_type & A_
Definition: WellOperators.hpp:387
const matrix_type & getmat() const override
Definition: WellOperators.hpp:358
Dune::SolverCategory::Category category() const override
Definition: WellOperators.hpp:312
void apply(const X &x, Y &y) const override
Definition: WellOperators.hpp:325
WellModelGhostLastMatrixAdapter(const M &A, const LinearOperatorExtra< X, Y > &wellOper, const std::size_t interiorSize)
constructor: just store a reference to a matrix
Definition: WellOperators.hpp:319
Dune::OwnerOverlapCopyCommunication< int, int > communication_type
Definition: WellOperators.hpp:307
Dune::BCRSMatrix< MatrixBlock< field_type, 1, 1 > > PressureMatrix
Definition: WellOperators.hpp:305
X domain_type
Definition: WellOperators.hpp:302
void applyscaleadd(field_type alpha, const X &x, Y &y) const override
Definition: WellOperators.hpp:343
int getNumberOfExtraEquations() const
Definition: WellOperators.hpp:374
Y range_type
Definition: WellOperators.hpp:303
std::size_t interiorSize_
Definition: WellOperators.hpp:389
const LinearOperatorExtra< X, Y > & wellOper_
Definition: WellOperators.hpp:388
M matrix_type
Definition: WellOperators.hpp:301
void ghostLastProject(Y &y) const
Definition: WellOperators.hpp:380
Adapter to combine a matrix and another linear operator into a combined linear operator.
Definition: WellOperators.hpp:225
Y range_type
Definition: WellOperators.hpp:229
const matrix_type & A_
Definition: WellOperators.hpp:285
void addWellPressureEquationsStruct(PressureMatrix &jacobian) const
Definition: WellOperators.hpp:273
int getNumberOfExtraEquations() const
Definition: WellOperators.hpp:279
void apply(const X &x, Y &y) const override
Definition: WellOperators.hpp:244
void addWellPressureEquations(PressureMatrix &jacobian, const X &weights, const bool use_well_weights) const
Definition: WellOperators.hpp:265
M matrix_type
Definition: WellOperators.hpp:227
typename X::field_type field_type
Definition: WellOperators.hpp:230
Dune::BCRSMatrix< MatrixBlock< field_type, 1, 1 > > PressureMatrix
Definition: WellOperators.hpp:231
const matrix_type & getmat() const override
Definition: WellOperators.hpp:263
WellModelMatrixAdapter(const M &A, const LinearOperatorExtra< X, Y > &wellOper)
constructor: just store a reference to a matrix
Definition: WellOperators.hpp:239
const LinearOperatorExtra< X, Y > & wellOper_
Definition: WellOperators.hpp:286
void applyscaleadd(field_type alpha, const X &x, Y &y) const override
Definition: WellOperators.hpp:254
X domain_type
Definition: WellOperators.hpp:228
Dune::SolverCategory::Category category() const override
Definition: WellOperators.hpp:233
Definition: fvbaseprimaryvariables.hh:141
Definition: blackoilboundaryratevector.hh:39