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 wellMod_.apply(x, y);
88 }
89
91 void applyscaleadd(field_type alpha, const X& x, Y& y) const override
92 {
93 OPM_TIMEBLOCK(applyscaleadd);
94 wellMod_.applyScaleAdd(alpha, x, y);
95 }
96
102 Dune::SolverCategory::Category category() const override
103 {
104 return Dune::SolverCategory::sequential;
105 }
106
108 const X& weights,
109 const bool use_well_weights) const override
110 {
111 OPM_TIMEBLOCK(addWellPressureEquations);
112 wellMod_.addWellPressureEquations(jacobian, weights, use_well_weights);
113 }
114
115 void addWellPressureEquationsStruct(PressureMatrix& jacobian) const override
116 {
117 OPM_TIMEBLOCK(addWellPressureEquationsStruct);
118 wellMod_.addWellPressureEquationsStruct(jacobian);
119 }
120
121 int getNumberOfExtraEquations() const override
122 {
123 return wellMod_.numLocalWellsEnd();
124 }
125
126private:
127 const WellModel& wellMod_;
128};
129
141template<class M, class X, class Y, bool overlapping >
142class WellModelMatrixAdapter : public Dune::AssembledLinearOperator<M,X,Y>
143{
144public:
145 using matrix_type = M;
146 using domain_type = X;
147 using range_type = Y;
148 using field_type = typename X::field_type;
149 using PressureMatrix = Dune::BCRSMatrix<MatrixBlock<field_type, 1, 1>>;
150#if HAVE_MPI
151 using communication_type = Dune::OwnerOverlapCopyCommunication<int,int>;
152#else
153 using communication_type = Dune::Communication<int>;
154#endif
155
156 Dune::SolverCategory::Category category() const override
157 {
158 return overlapping ?
159 Dune::SolverCategory::overlapping : Dune::SolverCategory::sequential;
160 }
161
164 const LinearOperatorExtra<X, Y>& wellOper,
165 const std::shared_ptr<communication_type>& comm = {})
166 : A_( A ), wellOper_( wellOper ), comm_(comm)
167 {}
168
169 void apply( const X& x, Y& y ) const override
170 {
171 OPM_TIMEBLOCK(apply);
172 A_.mv(x, y);
173
174 // add well model modification to y
175 wellOper_.apply(x, y);
176
177 #if HAVE_MPI
178 if (comm_) {
179 comm_->project(y);
180 }
181 #endif
182 }
183
184 // y += \alpha * A * x
185 void applyscaleadd (field_type alpha, const X& x, Y& y) const override
186 {
187 OPM_TIMEBLOCK(applyscaleadd);
188 A_.usmv(alpha, x, y);
189
190 // add scaled well model modification to y
191 wellOper_.applyscaleadd(alpha, x, y);
192
193 #if HAVE_MPI
194 if (comm_) {
195 comm_->project( y );
196 }
197 #endif
198 }
199
200 const matrix_type& getmat() const override { return A_; }
201
203 const X& weights,
204 const bool use_well_weights) const
205 {
206 OPM_TIMEBLOCK(addWellPressureEquations);
207 wellOper_.addWellPressureEquations(jacobian, weights, use_well_weights);
208 }
209
211 {
212 OPM_TIMEBLOCK(addWellPressureEquations);
213 wellOper_.addWellPressureEquationsStruct(jacobian);
214 }
215
217 {
218 return wellOper_.getNumberOfExtraEquations();
219 }
220
221protected:
224 std::shared_ptr<communication_type> comm_;
225};
226
235template<class M, class X, class Y, bool overlapping >
236class WellModelGhostLastMatrixAdapter : public Dune::AssembledLinearOperator<M,X,Y>
237{
238public:
239 using matrix_type = M;
240 using domain_type = X;
241 using range_type = Y;
242 using field_type = typename X::field_type;
243 using PressureMatrix = Dune::BCRSMatrix<MatrixBlock<field_type, 1, 1>>;
244#if HAVE_MPI
245 using communication_type = Dune::OwnerOverlapCopyCommunication<int,int>;
246#else
247 using communication_type = Dune::Communication<int>;
248#endif
249
250 Dune::SolverCategory::Category category() const override
251 {
252 return overlapping ?
253 Dune::SolverCategory::overlapping : Dune::SolverCategory::sequential;
254 }
255
258 const LinearOperatorExtra<X, Y>& wellOper,
259 const std::size_t interiorSize )
260 : A_( A ), wellOper_( wellOper ), interiorSize_(interiorSize)
261 {}
262
263 void apply(const X& x, Y& y) const override
264 {
265 OPM_TIMEBLOCK(apply);
266 for (auto row = A_.begin(); row.index() < interiorSize_; ++row)
267 {
268 y[row.index()]=0;
269 auto endc = (*row).end();
270 for (auto col = (*row).begin(); col != endc; ++col)
271 (*col).umv(x[col.index()], y[row.index()]);
272 }
273
274 // add well model modification to y
275 wellOper_.apply(x, y);
276
278 }
279
280 // y += \alpha * A * x
281 void applyscaleadd (field_type alpha, const X& x, Y& y) const override
282 {
283 OPM_TIMEBLOCK(applyscaleadd);
284 for (auto row = A_.begin(); row.index() < interiorSize_; ++row)
285 {
286 auto endc = (*row).end();
287 for (auto col = (*row).begin(); col != endc; ++col)
288 (*col).usmv(alpha, x[col.index()], y[row.index()]);
289 }
290 // add scaled well model modification to y
291 wellOper_.applyscaleadd(alpha, x, y);
292
294 }
295
296 const matrix_type& getmat() const override { return A_; }
297
299 const X& weights,
300 const bool use_well_weights) const
301 {
302 OPM_TIMEBLOCK(addWellPressureEquations);
303 wellOper_.addWellPressureEquations(jacobian, weights, use_well_weights);
304 }
305
307 {
308 OPM_TIMEBLOCK(addWellPressureEquationsStruct);
309 wellOper_.addWellPressureEquationsStruct(jacobian);
310 }
311
313 {
314 return wellOper_.getNumberOfExtraEquations();
315 }
316
317protected:
318 void ghostLastProject(Y& y) const
319 {
320 std::size_t end = y.size();
321 for (std::size_t i = interiorSize_; i < end; ++i)
322 y[i] = 0;
323 }
324
327 std::size_t interiorSize_;
328};
329
338template<class M, class X, class Y, class C>
339class GhostLastMatrixAdapter : public Dune::AssembledLinearOperator<M,X,Y>
340{
341public:
342 typedef M matrix_type;
343 typedef X domain_type;
344 typedef Y range_type;
345 typedef typename X::field_type field_type;
346
347
349
350 Dune::SolverCategory::Category category() const override
351 {
352 return Dune::SolverCategory::overlapping;
353 }
354
357 const communication_type& comm)
358 : A_( Dune::stackobject_to_shared_ptr(A) ), comm_(comm)
359 {
360 interiorSize_ = setInteriorSize(comm_);
361 }
362
363 GhostLastMatrixAdapter (const std::shared_ptr<M> A,
364 const communication_type& comm)
365 : A_( A ), comm_(comm)
366 {
367 interiorSize_ = setInteriorSize(comm_);
368 }
369
370 virtual void apply( const X& x, Y& y ) const override
371 {
372 for (auto row = A_->begin(); row.index() < interiorSize_; ++row)
373 {
374 y[row.index()]=0;
375 auto endc = (*row).end();
376 for (auto col = (*row).begin(); col != endc; ++col)
377 (*col).umv(x[col.index()], y[row.index()]);
378 }
379
380 ghostLastProject( y );
381 }
382
383 // y += \alpha * A * x
384 virtual void applyscaleadd (field_type alpha, const X& x, Y& y) const override
385 {
386 for (auto row = A_->begin(); row.index() < interiorSize_; ++row)
387 {
388 auto endc = (*row).end();
389 for (auto col = (*row).begin(); col != endc; ++col)
390 (*col).usmv(alpha, x[col.index()], y[row.index()]);
391 }
392
393 ghostLastProject( y );
394 }
395
396 virtual const matrix_type& getmat() const override { return *A_; }
397
398 size_t getInteriorSize() const { return interiorSize_;}
399
400private:
401 void ghostLastProject(Y& y) const
402 {
403 size_t end = y.size();
404 for (size_t i = interiorSize_; i < end; ++i)
405 y[i] = 0; //project to interiorsize, i.e. ignore ghost
406 }
407
408 size_t setInteriorSize(const communication_type& comm) const
409 {
410 auto indexSet = comm.indexSet();
411
412 size_t is = 0;
413 // Loop over index set
414 for (auto idx = indexSet.begin(); idx!=indexSet.end(); ++idx) {
415 //Only take "owner" indices
416 if (idx->local().attribute()==1) {
417 //get local index
418 auto loc = idx->local().local();
419 // if loc is higher than "old interior size", update it
420 if (loc > is) {
421 is = loc;
422 }
423 }
424 }
425 return is + 1; //size is plus 1 since we start at 0
426 }
427 const std::shared_ptr<const matrix_type> A_ ;
428 const communication_type& comm_;
429 size_t interiorSize_;
430};
431
432} // namespace Opm
433
434namespace Dune {
435namespace Amg {
436
437template<class M, class X, class Y, class C>
438class ConstructionTraits<Opm::GhostLastMatrixAdapter<M,X,Y,C> >
439{
440public:
441 typedef ParallelOperatorArgs<M,C> Arguments;
442
443 static inline std::shared_ptr<Opm::GhostLastMatrixAdapter<M,X,Y,C>> construct(const Arguments& args)
444 {
445 return std::make_shared<Opm::GhostLastMatrixAdapter<M,X,Y,C>>
446 (args.matrix_, args.comm_);
447 }
448};
449
450} // end namespace Amg
451} // end namespace Dune
452
453#endif // OPM_WELLOPERATORS_HEADER_INCLUDED
454
static std::shared_ptr< Opm::GhostLastMatrixAdapter< M, X, Y, C > > construct(const Arguments &args)
Definition: WellOperators.hpp:443
ParallelOperatorArgs< M, C > Arguments
Definition: WellOperators.hpp:441
Dune linear operator that assumes ghost rows are ordered after interior rows. Avoids some computation...
Definition: WellOperators.hpp:340
C communication_type
Definition: WellOperators.hpp:348
M matrix_type
Definition: WellOperators.hpp:342
Dune::SolverCategory::Category category() const override
Definition: WellOperators.hpp:350
X::field_type field_type
Definition: WellOperators.hpp:345
X domain_type
Definition: WellOperators.hpp:343
virtual const matrix_type & getmat() const override
Definition: WellOperators.hpp:396
GhostLastMatrixAdapter(const std::shared_ptr< M > A, const communication_type &comm)
Definition: WellOperators.hpp:363
GhostLastMatrixAdapter(const M &A, const communication_type &comm)
constructor: just store a reference to a matrix
Definition: WellOperators.hpp:356
virtual void applyscaleadd(field_type alpha, const X &x, Y &y) const override
Definition: WellOperators.hpp:384
Y range_type
Definition: WellOperators.hpp:344
virtual void apply(const X &x, Y &y) const override
Definition: WellOperators.hpp:370
size_t getInteriorSize() const
Definition: WellOperators.hpp:398
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
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
void addWellPressureEquations(PressureMatrix &jacobian, const X &weights, const bool use_well_weights) const override
Definition: WellOperators.hpp:107
void applyscaleadd(field_type alpha, const X &x, Y &y) const override
apply operator to x, scale and add:
Definition: WellOperators.hpp:91
Dune::SolverCategory::Category category() const override
Definition: WellOperators.hpp:102
int getNumberOfExtraEquations() const override
Definition: WellOperators.hpp:121
void addWellPressureEquationsStruct(PressureMatrix &jacobian) const override
Definition: WellOperators.hpp:115
Adapter to combine a matrix and another linear operator into a combined linear operator.
Definition: WellOperators.hpp:237
void addWellPressureEquationsStruct(PressureMatrix &jacobian) const
Definition: WellOperators.hpp:306
typename X::field_type field_type
Definition: WellOperators.hpp:242
void addWellPressureEquations(PressureMatrix &jacobian, const X &weights, const bool use_well_weights) const
Definition: WellOperators.hpp:298
const matrix_type & A_
Definition: WellOperators.hpp:325
const matrix_type & getmat() const override
Definition: WellOperators.hpp:296
Dune::SolverCategory::Category category() const override
Definition: WellOperators.hpp:250
void apply(const X &x, Y &y) const override
Definition: WellOperators.hpp:263
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:257
Dune::OwnerOverlapCopyCommunication< int, int > communication_type
Definition: WellOperators.hpp:245
Dune::BCRSMatrix< MatrixBlock< field_type, 1, 1 > > PressureMatrix
Definition: WellOperators.hpp:243
X domain_type
Definition: WellOperators.hpp:240
void applyscaleadd(field_type alpha, const X &x, Y &y) const override
Definition: WellOperators.hpp:281
int getNumberOfExtraEquations() const
Definition: WellOperators.hpp:312
Y range_type
Definition: WellOperators.hpp:241
std::size_t interiorSize_
Definition: WellOperators.hpp:327
const LinearOperatorExtra< X, Y > & wellOper_
Definition: WellOperators.hpp:326
M matrix_type
Definition: WellOperators.hpp:239
void ghostLastProject(Y &y) const
Definition: WellOperators.hpp:318
Adapter to combine a matrix and another linear operator into a combined linear operator.
Definition: WellOperators.hpp:143
void apply(const X &x, Y &y) const override
Definition: WellOperators.hpp:169
const matrix_type & A_
Definition: WellOperators.hpp:222
typename X::field_type field_type
Definition: WellOperators.hpp:148
int getNumberOfExtraEquations() const
Definition: WellOperators.hpp:216
void addWellPressureEquationsStruct(PressureMatrix &jacobian) const
Definition: WellOperators.hpp:210
X domain_type
Definition: WellOperators.hpp:146
M matrix_type
Definition: WellOperators.hpp:145
void applyscaleadd(field_type alpha, const X &x, Y &y) const override
Definition: WellOperators.hpp:185
WellModelMatrixAdapter(const M &A, const LinearOperatorExtra< X, Y > &wellOper, const std::shared_ptr< communication_type > &comm={})
constructor: just store a reference to a matrix
Definition: WellOperators.hpp:163
Dune::SolverCategory::Category category() const override
Definition: WellOperators.hpp:156
const matrix_type & getmat() const override
Definition: WellOperators.hpp:200
Dune::OwnerOverlapCopyCommunication< int, int > communication_type
Definition: WellOperators.hpp:151
Dune::BCRSMatrix< MatrixBlock< field_type, 1, 1 > > PressureMatrix
Definition: WellOperators.hpp:149
std::shared_ptr< communication_type > comm_
Definition: WellOperators.hpp:224
const LinearOperatorExtra< X, Y > & wellOper_
Definition: WellOperators.hpp:223
void addWellPressureEquations(PressureMatrix &jacobian, const X &weights, const bool use_well_weights) const
Definition: WellOperators.hpp:202
Y range_type
Definition: WellOperators.hpp:147
Definition: fvbaseprimaryvariables.hh:141
Definition: blackoilboundaryratevector.hh:37