MimeticIPAnisoRelpermEvaluator.hpp
Go to the documentation of this file.
1 //===========================================================================
2 //
3 // File: MimeticIPAnisoRelpermEvaluator.hpp
4 //
5 // Created: Mon Oct 19 10:22:22 2009
6 //
7 // Author(s): Atgeirr F Rasmussen <atgeirr@sintef.no>
8 // B�rd Skaflestad <bard.skaflestad@sintef.no>
9 //
10 // $Date$
11 //
12 // $Revision$
13 //
14 //===========================================================================
15 
16 /*
17  Copyright 2009, 2010 SINTEF ICT, Applied Mathematics.
18  Copyright 2009, 2010 Statoil ASA.
19 
20  This file is part of The Open Reservoir Simulator Project (OpenRS).
21 
22  OpenRS is free software: you can redistribute it and/or modify
23  it under the terms of the GNU General Public License as published by
24  the Free Software Foundation, either version 3 of the License, or
25  (at your option) any later version.
26 
27  OpenRS is distributed in the hope that it will be useful,
28  but WITHOUT ANY WARRANTY; without even the implied warranty of
29  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
30  GNU General Public License for more details.
31 
32  You should have received a copy of the GNU General Public License
33  along with OpenRS. If not, see <http://www.gnu.org/licenses/>.
34 */
35 
36 #ifndef OPENRS_MIMETICIPANISORELPERMEVALUATOR_HEADER
37 #define OPENRS_MIMETICIPANISORELPERMEVALUATOR_HEADER
38 
39 
40 #include <algorithm>
41 #include <vector>
42 
43 #include <boost/bind.hpp>
44 
45 #include <opm/common/ErrorMacros.hpp>
46 #include <opm/core/utility/SparseTable.hpp>
47 
51 
52 namespace Opm {
84  template<class GridInterface, class RockInterface>
86  {
87  public:
90  enum { dim = GridInterface::Dimension };
93  typedef typename GridInterface::CellIterator CellIter;
99  typedef typename CellIter::Scalar Scalar;
100 
101 
104  : max_nf_(-1),
105  prock_(0)
106  {}
107 
108 
117  : max_nf_ (max_nf ),
118  fa_ (max_nf * max_nf),
119  t1_ (max_nf * dim ),
120  t2_ (max_nf * dim ),
121  second_term_ ( ),
122  n_ ( ),
123  Kg_ ( ),
124  prock_ ( 0 )
125  {}
126 
127 
135  void init(const int max_nf)
136  {
137  max_nf_ = max_nf;
138  std::vector<double>(max_nf * max_nf).swap(fa_);
139  std::vector<double>(max_nf * dim ).swap(t1_);
140  std::vector<double>(max_nf * dim ).swap(t2_);
141  }
142 
143 
157  template<class Vector>
158  void reserveMatrices(const Vector& sz)
159  {
160  typedef typename Vector::value_type vt;
161 
162  Vector sz2(sz.size());
163 
164  std::transform(sz.begin(), sz.end(), sz2.begin(),
165  boost::bind(std::multiplies<vt>(), _1, _1));
166 
167  second_term_.allocate(sz2.begin(), sz2.end());
168 
169  std::transform(sz.begin(), sz.end(), sz2.begin(),
170  boost::bind(std::multiplies<vt>(), _1, int(dim)));
171 
172  n_.allocate(sz2.begin(), sz2.end());
173 
174  std::fill(sz2.begin(), sz2.end(), vt(dim));
175  Kg_.allocate(sz2.begin(), sz2.end());
176  }
177 
178 
207  void buildStaticContrib(const CellIter& c,
208  const RockInterface& r,
209  const typename CellIter::Vector& grav,
210  const int nf)
211  {
212  // Binv = (N*lambda*K*N' + t*diag(A)*(I - Q*Q')*diag(A))/vol
213  // ^ ^^^^^^^^^^^^^^^^^^^^^^^^^^
214  // precompute: n_ precompute: second_term_
215  // t = 6/dim * trace(lambda*K)
216 
217  typedef typename CellIter::FaceIterator FI;
218  typedef typename CellIter::Vector CV;
219  typedef typename FI ::Vector FV;
220 
221  // Now we need to remember the rocks, since we will need
222  // the permeability for dynamic assembly.
223  prock_ = &r;
224 
225  const int ci = c->index();
226 
227  static_assert (FV::dimension == int(dim), "");
228  assert (int(t1_.size()) >= nf * dim);
229  assert (int(t2_.size()) >= nf * dim);
230  assert (int(fa_.size()) >= nf * nf);
231 
232  SharedFortranMatrix T2 (nf, dim, &t2_ [0]);
233  SharedFortranMatrix fa (nf, nf , &fa_ [0]);
234  SharedFortranMatrix second_term(nf, nf, &second_term_[ci][0]);
235  SharedFortranMatrix n(nf, dim, &n_[ci][0]);
236 
237  // Clear matrices of any residual data.
238  zero(second_term); zero(n); zero(T2); zero(fa);
239 
240  // Setup: second_term <- I, n <- N, T2 <- C
241  const CV cc = c->centroid();
242  int i = 0;
243  for (FI f = c->facebegin(); f != c->faceend(); ++f, ++i) {
244  second_term(i,i) = Scalar(1.0);
245  fa(i,i) = f->area();
246 
247  FV fc = f->centroid(); fc -= cc; fc *= fa(i,i);
248  FV fn = f->normal (); fn *= fa(i,i);
249 
250  for (int j = 0; j < dim; ++j) {
251  n (i,j) = fn[j];
252  T2(i,j) = fc[j];
253  }
254  }
255  assert (i == nf);
256 
257  // T2 <- orth(T2)
258  if (orthogonalizeColumns(T2) != 0) {
259  assert (false);
260  }
261 
262  // second_term <- second_term - T2*T2' == I - Q*Q'
263  symmetricUpdate(Scalar(-1.0), T2, Scalar(1.0), second_term);
264 
265  // second_term <- diag(A) * second_term * diag(A)
266  symmetricUpdate(fa, second_term);
267 
268  // Gravity term: Kg_ = K * grav
269  vecMulAdd_N(Scalar(1.0), r.permeability(ci), &grav[0],
270  Scalar(0.0), &Kg_[ci][0]);
271  }
272 
273 
298  template<class FluidInterface, class Sat>
299  void computeDynamicParams(const CellIter& c,
300  const FluidInterface& fl,
301  const std::vector<Sat>& s)
302  {
303  const int ci = c->index();
304 
305  std::array<Scalar, dim * dim> lambda_t;
306  std::array<Scalar, dim * dim> pmob_data;
307 
308  SharedFortranMatrix pmob(dim, dim, &pmob_data[0]);
309  SharedFortranMatrix Kg (dim, 1 , &Kg_[ci][0]);
310 
311  std::array<Scalar, FluidInterface::NumberOfPhases> rho;
312  fl.phaseDensities(ci, rho);
313 
314  std::fill(dyn_Kg_.begin(), dyn_Kg_.end(), Scalar(0.0));
315  std::fill(lambda_t.begin(), lambda_t.end(), 0.0);
316 
317  for (int phase = 0; phase < FluidInterface::NumberOfPhases; ++phase) {
318  fl.phaseMobility(phase, ci, s[ci], pmob);
319 
320  // dyn_Kg_ += (\rho_phase \lambda_phase) Kg
321  vecMulAdd_N(rho[phase], pmob, Kg.data(), Scalar(1.0), dyn_Kg_.data());
322 
323  // \lambda_t += \lambda_phase
324  std::transform(lambda_t.begin(), lambda_t.end(), pmob_data.begin(),
325  lambda_t.begin(),
326  std::plus<Scalar>());
327  }
328 
329  // lambdaK_ = (\sum_i \lambda_i) K
330  SharedFortranMatrix lambdaT(dim, dim, lambda_t.data());
331  SharedFortranMatrix lambdaK(dim, dim, lambdaK_.data());
332  prod(lambdaT, prock_->permeability(ci), lambdaK);
333  }
334 
335 
362  template<template<typename> class SP>
363  void getInverseMatrix(const CellIter& c,
364  FullMatrix<Scalar,SP,FortranOrdering>& Binv) const
365  {
366  // Binv = (N*lambda*K*N' + t*diag(A)*(I - Q*Q')*diag(A))/vol
367  // ^ ^^^^^^^^^^^^^^^^^^^^^^^^^^
368  // precomputed: n_ precomputed: second_term_
369  // t = 6/dim * trace(lambda*K)
370  int ci = c->index();
371  int nf = Binv.numRows();
372  ImmutableFortranMatrix n(nf, dim, &n_[ci][0]);
373  ImmutableFortranMatrix t2(nf, nf, &second_term_[ci][0]);
374  Binv = t2;
375  ImmutableFortranMatrix lambdaK(dim, dim, lambdaK_.data());
376  SharedFortranMatrix T2(nf, dim, &t2_[0]);
377 
378  // T2 <- N*lambda*K
379  matMulAdd_NN(Scalar(1.0), n, lambdaK, Scalar(0.0), T2);
380 
381  // Binv <- (T2*N' + t*Binv) / vol(c)
382  // == (N*lambda*K*N' + t*(diag(A) * (I - Q*Q') * diag(A))) / vol(c)
383  //
384  // where t = 6/d * TRACE(lambda*K) (== 2*TRACE(lambda*K) for 3D).
385  //
386  Scalar t = Scalar(6.0) * trace(lambdaK) / dim;
387  matMulAdd_NT(Scalar(1.0) / c->volume(), T2, n,
388  t / c->volume(), Binv );
389  }
390 
391 
405  template<class Vector>
406  void gravityFlux(const CellIter& c,
407  Vector& gflux) const
408  {
409  const int ci = c->index();
410  const int nf = n_.rowSize(ci) / dim;
411 
412  ImmutableFortranMatrix N(nf, dim, &n_[ci][0]);
413 
414  // gflux = N (\sum_i \rho_i \lambda_i) Kg
415  vecMulAdd_N(Scalar(1.0), N, &dyn_Kg_[0],
416  Scalar(0.0), &gflux[0]);
417  }
418 
419  private:
420  int max_nf_ ;
421  mutable std::vector<Scalar> fa_, t1_, t2_;
422  Opm::SparseTable<Scalar> second_term_ ;
423  Opm::SparseTable<Scalar> n_ ;
424  Opm::SparseTable<Scalar> Kg_ ;
425  std::array<Scalar, dim> dyn_Kg_ ;
426  std::array<double, dim*dim> lambdaK_ ;
427  const RockInterface* prock_ ;
428  };
429 } // namespace Opm
430 
431 
432 #endif // OPENRS_MIMETICIPANISORELPERMEVALUATOR_HEADER
void symmetricUpdate(const T &a1, const FullMatrix< T, StoragePolicy, FortranOrdering > &A, const T &a2, FullMatrix< T, StoragePolicy, FortranOrdering > &C)
Symmetric, rank update of symmetric matrix. Specifically, .
Definition: Matrix.hpp:830
void matMulAdd_NT(const T &a1, const FullMatrix< T, SP1, FortranOrdering > &A, const FullMatrix< T, SP2, FortranOrdering > &B, const T &a2, FullMatrix< T, SP3, FortranOrdering > &C)
GEneral Matrix-Matrix product update of other matrix. Specificlly, .
Definition: Matrix.hpp:1195
void getInverseMatrix(const CellIter &c, FullMatrix< Scalar, SP, FortranOrdering > &Binv) const
Retrieve the dynamic (mobility updated) inverse mimetic inner product matrix for specific cell...
Definition: MimeticIPAnisoRelpermEvaluator.hpp:363
void reserveMatrices(const Vector &sz)
Reserve internal space for storing values of (static) IP contributions for given set of cells...
Definition: MimeticIPAnisoRelpermEvaluator.hpp:158
void computeDynamicParams(const CellIter &c, const FluidInterface &fl, const std::vector< Sat > &s)
Evaluate dynamic (saturation dependent) properties in single cell.
Definition: MimeticIPAnisoRelpermEvaluator.hpp:299
Definition: BlackoilFluid.hpp:31
const FullMatrix< double, ImmutableSharedData, FortranOrdering > ImmutableFortranMatrix
Definition: Matrix.hpp:591
MimeticIPAnisoRelpermEvaluator()
Default constructor.
Definition: MimeticIPAnisoRelpermEvaluator.hpp:103
GridInterface::CellIterator CellIter
The iterator type for iterating over grid cells.
Definition: MimeticIPAnisoRelpermEvaluator.hpp:93
void gravityFlux(const CellIter &c, Vector &gflux) const
Compute gravity flux for all faces of single cell.
Definition: MimeticIPAnisoRelpermEvaluator.hpp:406
Dune::FieldVector< typename Matrix::value_type, rows > prod(const Matrix &A, const Dune::FieldVector< typename Matrix::value_type, rows > &x)
Matrix applied to a vector.
Definition: Matrix.hpp:669
FullMatrix< double, SharedData, FortranOrdering > SharedFortranMatrix
Definition: Matrix.hpp:590
Definition: MimeticIPAnisoRelpermEvaluator.hpp:85
void init(const int max_nf)
Initialization routine.
Definition: MimeticIPAnisoRelpermEvaluator.hpp:135
void matMulAdd_NN(const T &a1, const FullMatrix< T, SP1, FortranOrdering > &A, const FullMatrix< T, SP2, FortranOrdering > &B, const T &a2, FullMatrix< T, SP3, FortranOrdering > &C)
GEneral Matrix-Matrix product update of other matrix. Specificlly, .
Definition: Matrix.hpp:1137
int orthogonalizeColumns(FullMatrix< T, StoragePolicy, FortranOrdering > &A)
Construct orthonormal basis for matrix range (i.e., column space). Based on a QR factorization of the...
Definition: Matrix.hpp:730
Definition: MimeticIPAnisoRelpermEvaluator.hpp:90
CellIter::Scalar Scalar
The element type of the matrix representation of the mimetic inner product. Assumed to be a floating ...
Definition: MimeticIPAnisoRelpermEvaluator.hpp:99
Matrix::value_type trace(const Matrix &A)
Compute matrix trace (i.e., sum(diag(A))).
Definition: Matrix.hpp:639
MimeticIPAnisoRelpermEvaluator(const int max_nf)
Constructor.
Definition: MimeticIPAnisoRelpermEvaluator.hpp:116
void buildStaticContrib(const CellIter &c, const RockInterface &r, const typename CellIter::Vector &grav, const int nf)
Main evaluation routine. Computes the inverse of the matrix representation of the mimetic inner produ...
Definition: MimeticIPAnisoRelpermEvaluator.hpp:207
void vecMulAdd_N(const T &a1, const FullMatrix< T, SP, FortranOrdering > &A, const std::vector< T > &x, const T &a2, std::vector< T > &y)
GEneral Matrix-Vector product (GAXPY operation). Specifically, .
Definition: Matrix.hpp:914
void zero(Matrix &A)
Zero-fill a.
Definition: Matrix.hpp:603