ReservoirPropertyCommon_impl.hpp
Go to the documentation of this file.
1 //===========================================================================
2 //
3 // File: ReservoirPropertyCommon_impl.hpp
4 //
5 // Created: Mon Oct 26 08:29:09 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_RESERVOIRPROPERTYCOMMON_IMPL_HEADER
37 #define OPENRS_RESERVOIRPROPERTYCOMMON_IMPL_HEADER
38 
39 #include <opm/core/io/eclipse/EclipseGridInspector.hpp>
40 #include <opm/parser/eclipse/Deck/Deck.hpp>
41 
42 #include <fstream>
43 #include <array>
44 
45 namespace Opm
46 {
47 
48  namespace {
49 
69  PermeabilityKind classifyPermeability(Opm::DeckConstPtr deck)
70  {
71  const bool xx = deck->hasKeyword("PERMX" );
72  const bool xy = deck->hasKeyword("PERMXY");
73  const bool xz = deck->hasKeyword("PERMXZ");
74 
75  const bool yx = deck->hasKeyword("PERMYX");
76  const bool yy = deck->hasKeyword("PERMY" );
77  const bool yz = deck->hasKeyword("PERMYZ");
78 
79  const bool zx = deck->hasKeyword("PERMZX");
80  const bool zy = deck->hasKeyword("PERMZY");
81  const bool zz = deck->hasKeyword("PERMZ" );
82 
83  int num_cross_comp = xy + xz + yx + yz + zx + zy;
84  int num_comp = xx + yy + zz + num_cross_comp;
85  PermeabilityKind retval = None;
86  if (num_cross_comp > 0) {
87  retval = TensorPerm;
88  } else {
89  if (num_comp == 1) {
90  retval = ScalarPerm;
91  } else if (num_comp >= 2) {
92  retval = DiagonalPerm;
93  }
94  }
95 
96  bool ok = true;
97  if (num_comp > 0) {
98  // At least one tensor component specified on input.
99  // Verify that any remaining components are OK from a
100  // structural point of view. In particular, there
101  // must not be any cross-components (e.g., k_{xy})
102  // unless the corresponding diagonal component (e.g.,
103  // k_{xx}) is present as well...
104  //
105  ok = xx || !(xy || xz || yx || zx) ;
106  ok = ok && (yy || !(yx || yz || xy || zy));
107  ok = ok && (zz || !(zx || zy || xz || yz));
108  }
109  if (!ok) {
110  retval = Invalid;
111  }
112 
113  return retval;
114  }
115 
116 
137  void setScalarPermIfNeeded(std::array<int,9>& kmap,
138  int i, int j, int k)
139  {
140  if (kmap[j] == 0) { kmap[j] = kmap[i]; }
141  if (kmap[k] == 0) { kmap[k] = kmap[i]; }
142  }
143 
144 
177  PermeabilityKind fillTensor(Opm::DeckConstPtr deck,
178  std::vector<const std::vector<double>*>& tensor,
179  std::array<int,9>& kmap)
180  {
181  PermeabilityKind kind = classifyPermeability(deck);
182  if (kind == Invalid) {
183  OPM_THROW(std::runtime_error, "Invalid set of permeability fields given.");
184  }
185  assert (tensor.size() == 1);
186  for (int i = 0; i < 9; ++i) { kmap[i] = 0; }
187 
188  enum { xx, xy, xz, // 0, 1, 2
189  yx, yy, yz, // 3, 4, 5
190  zx, zy, zz }; // 6, 7, 8
191 
192  // -----------------------------------------------------------
193  // 1st row: [kxx, kxy, kxz]
194  if (deck->hasKeyword("PERMX" )) {
195  kmap[xx] = tensor.size();
196  tensor.push_back(&deck->getKeyword("PERMX")->getSIDoubleData());
197 
198  setScalarPermIfNeeded(kmap, xx, yy, zz);
199  }
200  if (deck->hasKeyword("PERMXY")) {
201  kmap[xy] = kmap[yx] = tensor.size(); // Enforce symmetry.
202  tensor.push_back(&deck->getKeyword("PERMXY")->getSIDoubleData());
203  }
204  if (deck->hasKeyword("PERMXZ")) {
205  kmap[xz] = kmap[zx] = tensor.size(); // Enforce symmetry.
206  tensor.push_back(&deck->getKeyword("PERMXZ")->getSIDoubleData());
207  }
208 
209  // -----------------------------------------------------------
210  // 2nd row: [kyx, kyy, kyz]
211  if (deck->hasKeyword("PERMYX")) {
212  kmap[yx] = kmap[xy] = tensor.size(); // Enforce symmetry.
213  tensor.push_back(&deck->getKeyword("PERMYX")->getSIDoubleData());
214  }
215  if (deck->hasKeyword("PERMY" )) {
216  kmap[yy] = tensor.size();
217  tensor.push_back(&deck->getKeyword("PERMY")->getSIDoubleData());
218 
219  setScalarPermIfNeeded(kmap, yy, zz, xx);
220  }
221  if (deck->hasKeyword("PERMYZ")) {
222  kmap[yz] = kmap[zy] = tensor.size(); // Enforce symmetry.
223  tensor.push_back(&deck->getKeyword("PERMYZ")->getSIDoubleData());
224  }
225 
226  // -----------------------------------------------------------
227  // 3rd row: [kzx, kzy, kzz]
228  if (deck->hasKeyword("PERMZX")) {
229  kmap[zx] = kmap[xz] = tensor.size(); // Enforce symmetry.
230  tensor.push_back(&deck->getKeyword("PERMZX")->getSIDoubleData());
231  }
232  if (deck->hasKeyword("PERMZY")) {
233  kmap[zy] = kmap[yz] = tensor.size(); // Enforce symmetry.
234  tensor.push_back(&deck->getKeyword("PERMZY")->getSIDoubleData());
235  }
236  if (deck->hasKeyword("PERMZ" )) {
237  kmap[zz] = tensor.size();
238  tensor.push_back(&deck->getKeyword("PERMZ")->getSIDoubleData());
239 
240  setScalarPermIfNeeded(kmap, zz, xx, yy);
241  }
242  return kind;
243  }
244 
245  } // anonymous namespace
246 
247 
248 
249 
250  // ----- Methods of ReservoirPropertyCommon -----
251 
252 
253 
254 
255  template <int dim, class RPImpl, class RockType>
257 #if 1
258  : density1_ (1013.9*Opm::unit::kilogram/Opm::unit::cubic(Opm::unit::meter)),
259  density2_ ( 834.7*Opm::unit::kilogram/Opm::unit::cubic(Opm::unit::meter)),
260  viscosity1_( 1.0*Opm::prefix::centi*Opm::unit::Poise),
261  viscosity2_( 3.0*Opm::prefix::centi*Opm::unit::Poise),
262 #else
263  : density1_ (1000.0*Opm::unit::kilogram/Opm::unit::cubic(Opm::unit::meter)),
264  density2_ (1000.0*Opm::unit::kilogram/Opm::unit::cubic(Opm::unit::meter)),
265  viscosity1_( 1000.0*Opm::prefix::centi*Opm::unit::Poise),
266  viscosity2_( 1000.0*Opm::prefix::centi*Opm::unit::Poise),
267 #endif
268  permeability_kind_(Invalid)
269  {
270  }
271 
272 
273  template <int dim, class RPImpl, class RockType>
275  const std::vector<int>& global_cell,
276  const double perm_threshold,
277  const std::string* rock_list_filename,
278  const bool use_jfunction_scaling,
279  const double sigma,
280  const double theta)
281  {
282  // This code is mostly copied from ReservoirPropertyCommon::init(...).
283  static_assert(dim == 3, "");
284 
285  permfield_valid_.assign(global_cell.size(),
286  std::vector<unsigned char>::value_type(0));
287 
288  assignPorosity (deck, global_cell);
289  assignNTG (deck, global_cell);
290  assignSWCR (deck, global_cell);
291  assignSOWCR (deck, global_cell);
292  assignPermeability(deck, global_cell, perm_threshold);
293  assignRockTable (deck, global_cell);
294 
295  if (rock_list_filename) {
296  readRocks(*rock_list_filename);
297  }
298 
299  // Added section. This is a hack, because not all rock classes
300  // may care about J-scaling. They still have to implement
301  // setUseJfunctionScaling() and setSigmaAndTheta(), though the
302  // latter may throw if called for a rock where it does not make sense.
303  int num_rocks = rock_.size();
304  for (int i = 0; i < num_rocks; ++i) {
305  rock_[i].setUseJfunctionScaling(use_jfunction_scaling);
306  if (use_jfunction_scaling) {
307  rock_[i].setSigmaAndTheta(sigma, theta);
308  }
309  }
310  // End of added section.
311 
312  asImpl().computeCflFactors();
313  }
314 
315 
316 
317  template <int dim, class RPImpl, class RockType>
319  const double uniform_poro,
320  const double uniform_perm)
321  {
322  permfield_valid_.assign(num_cells, std::vector<unsigned char>::value_type(1));
323  porosity_.assign(num_cells, uniform_poro);
324  permeability_.assign(dim*dim*num_cells, 0.0);
325  for (int i = 0; i < num_cells; ++i) {
326  SharedPermTensor K = permeabilityModifiable(i);
327  for (int dd = 0; dd < dim; ++dd) {
328  K(dd, dd) = uniform_perm;
329  }
330  }
331  cell_to_rock_.assign(num_cells, 0);
332  asImpl().computeCflFactors();
333  }
334 
335  template <int dim, class RPImpl, class RockType>
337  {
338  viscosity1_ = v1;
339  viscosity2_ = v2;
340  }
341 
342  template <int dim, class RPImpl, class RockType>
344  {
345  density1_ = d1;
346  density2_ = d2;
347  }
348 
349  template <int dim, class RPImpl, class RockType>
351  {
352  return viscosity1_;
353  }
354 
355 
356  template <int dim, class RPImpl, class RockType>
358  {
359  return viscosity2_;
360  }
361 
362 
363  template <int dim, class RPImpl, class RockType>
365  {
366  return density1_;
367  }
368 
369 
370  template <int dim, class RPImpl, class RockType>
372  {
373  return density2_;
374  }
375 
376 
377  template <int dim, class RPImpl, class RockType>
379  {
380  return porosity_[cell_index];
381  }
382 
383  template <int dim, class RPImpl, class RockType>
385  {
386  return ntg_[cell_index];
387  }
388 
389  template <int dim, class RPImpl, class RockType>
391  {
392  return swcr_[cell_index];
393  }
394 
395  template <int dim, class RPImpl, class RockType>
397  {
398  return sowcr_[cell_index];
399  }
400 
401 
402  template <int dim, class RPImpl, class RockType>
405  {
406  assert (permfield_valid_[cell_index]);
407 
408  const PermTensor K(dim, dim, &permeability_[dim*dim*cell_index]);
409  return K;
410  }
411 
412 
413  template <int dim, class RPImpl, class RockType>
416  {
417  // Typically only used for assigning synthetic perm values.
418  SharedPermTensor K(dim, dim, &permeability_[dim*dim*cell_index]);
419 
420  // Trust caller!
421  permfield_valid_[cell_index] = std::vector<unsigned char>::value_type(1);
422 
423  return K;
424  }
425 
426 
427  template <int dim, class RPImpl, class RockType>
428  template<class Vector>
429  void ReservoirPropertyCommon<dim, RPImpl, RockType>::phaseDensities(int /*cell_index*/, Vector& density) const
430  {
431  assert (density.size() >= NumberOfPhases);
432  density[0] = densityFirstPhase();
433  density[1] = densitySecondPhase();
434  }
435 
436 
437  template <int dim, class RPImpl, class RockType>
439  {
440  return density1_ - density2_;
441  }
442 
443 
444  template <int dim, class RPImpl, class RockType>
446  {
447  return cfl_factor_;
448  }
449 
450 
451  template <int dim, class RPImpl, class RockType>
453  {
454  return cfl_factor_gravity_;
455  }
456 
457 
458  template <int dim, class RPImpl, class RockType>
460  {
461  return cfl_factor_capillary_;
462  }
463 
464  template <int dim, class RPImpl, class RockType>
465  double ReservoirPropertyCommon<dim, RPImpl, RockType>::capillaryPressure(int cell_index, double saturation) const
466  {
467  if (rock_.size() > 0) {
468  int r = cell_to_rock_[cell_index];
469  return rock_[r].capPress(permeability(cell_index), porosity(cell_index), saturation);
470  } else {
471  // HACK ALERT!
472  // Use zero capillary pressure if no known rock table exists.
473  return 1e5*(1-saturation);
474  }
475  }
476 
477  template <int dim, class RPImpl, class RockType>
478  double ReservoirPropertyCommon<dim, RPImpl, RockType>::capillaryPressureDeriv(int cell_index, double saturation) const
479  {
480  if (rock_.size() > 0) {
481  int r = cell_to_rock_[cell_index];
482  double dpc = rock_[r].capPressDeriv(permeability(cell_index), porosity(cell_index), saturation);
483  return dpc;
484  } else {
485  // HACK ALERT!
486  // Use zero capillary pressure if no known rock table exists.
487  return -1.0e5;
488  }
489  }
490  template <int dim, class RPImpl, class RockType>
492  {
493  if (rock_.size() > 0) {
494  int r = cell_to_rock_[cell_index];
495  return rock_[r].s_min();
496  } else {
497  // HACK ALERT!
498  // Use zero as minimum saturation if no known rock table exists.
499  return 0;
500  }
501  }
502 
503  template <int dim, class RPImpl, class RockType>
505  {
506  if (rock_.size() > 0) {
507  int r = cell_to_rock_[cell_index];
508  return rock_[r].s_max();
509  } else {
510  // HACK ALERT!
511  // Use 1 as maximum saturation if no known rock table exists.
512  return 1;
513  }
514  }
515 
516  template <int dim, class RPImpl, class RockType>
518  {
519  if (rock_.size() > 0) {
520  int r = cell_to_rock_[cell_index];
521  return rock_[r].satFromCapPress(permeability(cell_index), porosity(cell_index), cap_press);
522  } else {
523  // HACK ALERT!
524  // Use a zero saturation if no known rock table exists.
525  return 0.0;
526  }
527  }
528 
529 
530  template <int dim, class RPImpl, class RockType>
532  {
533  int num_cells = porosity_.size();
534  // Write porosity.
535  {
536  std::string filename = grid_prefix + "-poro.dat";
537  std::ofstream file(filename.c_str());
538  if (!file) {
539  OPM_THROW(std::runtime_error, "Could not open file " << filename);
540  }
541  file << num_cells << '\n';
542  std::copy(porosity_.begin(), porosity_.end(), std::ostream_iterator<double>(file, "\n"));
543  }
544  // Write permeability.
545  {
546  std::string filename = grid_prefix + "-perm.dat";
547  std::ofstream file(filename.c_str());
548  if (!file) {
549  OPM_THROW(std::runtime_error, "Could not open file " << filename);
550  }
551  file << num_cells << '\n';
552  switch (permeability_kind_) {
553  case TensorPerm:
554  std::copy(permeability_.begin(), permeability_.end(), std::ostream_iterator<double>(file, "\n"));
555  break;
556  case DiagonalPerm:
557  for (int c = 0; c < num_cells; ++c) {
558  int index = c*dim*dim;
559  for (int dd = 0; dd < dim; ++dd) {
560  file << permeability_[index + (dim + 1)*dd] << ' ';
561  }
562  file << '\n';
563  }
564  break;
565  case ScalarPerm:
566  case None: // Treated like a scalar permeability.
567  for (int c = 0; c < num_cells; ++c) {
568  file << permeability_[c*dim*dim] << '\n';
569  }
570  break;
571  default:
572  OPM_THROW(std::runtime_error, "Cannot write invalid permeability.");
573  }
574  }
575  }
576 
577 
578  // ------ Private methods ------
579 
580 
581  template <int dim, class RPImpl, class RockType>
583  {
584  return static_cast<RPImpl&>(*this);
585  }
586 
587 
588 
589 
590  template <int dim, class RPImpl, class RockType>
592  const std::vector<int>& global_cell)
593  {
594  porosity_.assign(global_cell.size(), 1.0);
595 
596  if (deck->hasKeyword("PORO")) {
597  Opm::EclipseGridInspector insp(deck);
598  std::array<int, 3> dims = insp.gridSize();
599  int num_global_cells = dims[0]*dims[1]*dims[2];
600  const std::vector<double>& poro = deck->getKeyword("PORO")->getSIDoubleData();
601  if (int(poro.size()) != num_global_cells) {
602  OPM_THROW(std::runtime_error, "PORO field must have the same size as the "
603  "logical cartesian size of the grid: "
604  << poro.size() << " != " << num_global_cells);
605  }
606  for (int c = 0; c < int(porosity_.size()); ++c) {
607  porosity_[c] = poro[global_cell[c]];
608  }
609  }
610  }
611 
612  template <int dim, class RPImpl, class RockType>
614  const std::vector<int>& global_cell)
615  {
616  ntg_.assign(global_cell.size(), 1.0);
617 
618  if (deck->hasKeyword("NTG")) {
619  Opm::EclipseGridInspector insp(deck);
620  std::array<int, 3> dims = insp.gridSize();
621  int num_global_cells = dims[0]*dims[1]*dims[2];
622  const std::vector<double>& ntg = deck->getKeyword("NTG")->getSIDoubleData();
623  if (int(ntg.size()) != num_global_cells) {
624  OPM_THROW(std::runtime_error, "NTG field must have the same size as the "
625  "logical cartesian size of the grid: "
626  << ntg.size() << " != " << num_global_cells);
627  }
628  for (int c = 0; c < int(ntg_.size()); ++c) {
629  ntg_[c] = ntg[global_cell[c]];
630  }
631  }
632  }
633 
634  template <int dim, class RPImpl, class RockType>
636  const std::vector<int>& global_cell)
637  {
638  swcr_.assign(global_cell.size(), 0.0);
639 
640  if (deck->hasKeyword("SWCR")) {
641  Opm::EclipseGridInspector insp(deck);
642  std::array<int, 3> dims = insp.gridSize();
643  int num_global_cells = dims[0]*dims[1]*dims[2];
644  const std::vector<double>& swcr =
645  deck->getKeyword("SWCR")->getSIDoubleData();
646  if (int(swcr.size()) != num_global_cells) {
647  OPM_THROW(std::runtime_error, "SWCR field must have the same size as the "
648  "logical cartesian size of the grid: "
649  << swcr.size() << " != " << num_global_cells);
650  }
651  for (int c = 0; c < int(swcr_.size()); ++c) {
652  swcr_[c] = swcr[global_cell[c]];
653  }
654  }
655  }
656 
657  template <int dim, class RPImpl, class RockType>
659  const std::vector<int>& global_cell)
660  {
661  sowcr_.assign(global_cell.size(), 0.0);
662 
663  if (deck->hasKeyword("SOWCR")) {
664  Opm::EclipseGridInspector insp(deck);
665  std::array<int, 3> dims = insp.gridSize();
666  int num_global_cells = dims[0]*dims[1]*dims[2];
667  const std::vector<double>& sowcr =
668  deck->getKeyword("SOWCR")->getSIDoubleData();
669  if (int(sowcr.size()) != num_global_cells) {
670  OPM_THROW(std::runtime_error, "SOWCR field must have the same size as the "
671  "logical cartesian size of the grid: "
672  << sowcr.size() << " != " << num_global_cells);
673  }
674  for (int c = 0; c < int(sowcr_.size()); ++c) {
675  sowcr_[c] = sowcr[global_cell[c]];
676  }
677  }
678  }
679 
680  template <int dim, class RPImpl, class RockType>
682  const std::vector<int>& global_cell,
683  double perm_threshold)
684  {
685  Opm::EclipseGridInspector insp(deck);
686  std::array<int, 3> dims = insp.gridSize();
687  int num_global_cells = dims[0]*dims[1]*dims[2];
688  assert (num_global_cells > 0);
689 
690  permeability_.assign(dim * dim * global_cell.size(), 0.0);
691 
692  std::vector<const std::vector<double>*> tensor;
693  tensor.reserve(10);
694 
695  const std::vector<double> zero(num_global_cells, 0.0);
696  tensor.push_back(&zero);
697 
698  static_assert(dim == 3, "");
699  std::array<int,9> kmap;
700  permeability_kind_ = fillTensor(deck, tensor, kmap);
701  for (int i = 1; i < int(tensor.size()); ++i) {
702  if (int(tensor[i]->size()) != num_global_cells) {
703  OPM_THROW(std::runtime_error, "All permeability fields must have the same size as the "
704  "logical cartesian size of the grid: "
705  << (tensor[i]->size()) << " != " << num_global_cells);
706  }
707  }
708 
709  // Assign permeability values only if such values are
710  // given in the input deck represented by 'deck'. In
711  // other words: Don't set any (arbitrary) default values.
712  // It is infinitely better to experience a reproducible
713  // crash than subtle errors resulting from a (poorly
714  // chosen) default value...
715  //
716  if (tensor.size() > 1) {
717  const int nc = global_cell.size();
718  int off = 0;
719 
720  for (int c = 0; c < nc; ++c, off += dim*dim) {
721  SharedPermTensor K(dim, dim, &permeability_[off]);
722  int kix = 0;
723  const int glob = global_cell[c];
724 
725  for (int i = 0; i < dim; ++i) {
726  for (int j = 0; j < dim; ++j, ++kix) {
727  K(i,j) = (*tensor[kmap[kix]])[glob];
728  }
729  K(i,i) = std::max(K(i,i), perm_threshold);
730  }
731 
732  permfield_valid_[c] = std::vector<unsigned char>::value_type(1);
733  }
734  }
735  }
736 
737 
738 
739 
740  template <int dim, class RPImpl, class RockType>
742  const std::vector<int>& global_cell)
743  {
744  const int nc = global_cell.size();
745 
746  cell_to_rock_.assign(nc, 0);
747 
748  if (deck->hasKeyword("SATNUM")) {
749  Opm::EclipseGridInspector insp(deck);
750  std::array<int, 3> dims = insp.gridSize();
751  int num_global_cells = dims[0]*dims[1]*dims[2];
752  const std::vector<int>& satnum = deck->getKeyword("SATNUM")->getIntData();
753  if (int(satnum.size()) != num_global_cells) {
754  OPM_THROW(std::runtime_error, "SATNUM field must have the same size as the "
755  "logical cartesian size of the grid: "
756  << satnum.size() << " != " << num_global_cells);
757  }
758  for (int c = 0; c < nc; ++c) {
759  // Note: SATNUM is FORTRANish, ranging from 1 to n, therefore we subtract one.
760  cell_to_rock_[c] = satnum[global_cell[c]] - 1;
761  }
762  }
763  else if (deck->hasKeyword("ROCKTYPE")) {
764  Opm::EclipseGridInspector insp(deck);
765  std::array<int, 3> dims = insp.gridSize();
766  int num_global_cells = dims[0]*dims[1]*dims[2];
767  const std::vector<int>& satnum = deck->getKeyword("ROCKTYPE")->getIntData();
768  if (int(satnum.size()) != num_global_cells) {
769  OPM_THROW(std::runtime_error, "ROCKTYPE field must have the same size as the "
770  "logical cartesian size of the grid: "
771  << satnum.size() << " != " << num_global_cells);
772  }
773  for (int c = 0; c < nc; ++c) {
774  // Note: ROCKTYPE is FORTRANish, ranging from 1 to n, therefore we subtract one.
775  cell_to_rock_[c] = satnum[global_cell[c]] - 1;
776  }
777  }
778  }
779 
780 
781 
782 
783  template <int dim, class RPImpl, class RockType>
784  void ReservoirPropertyCommon<dim, RPImpl, RockType>::readRocks(const std::string& rock_list_file)
785  {
786  std::ifstream rl(rock_list_file.c_str());
787  if (!rl) {
788  OPM_THROW(std::runtime_error, "Could not open file " << rock_list_file);
789  }
790  int num_rocks = -1;
791  rl >> num_rocks;
792  assert(num_rocks >= 1);
793  rock_.resize(num_rocks);
794  std::string dir(rock_list_file.begin(), rock_list_file.begin() + rock_list_file.find_last_of('/') + 1);
795  for (int i = 0; i < num_rocks; ++i) {
796  std::string spec;
797  while (spec.empty()) {
798  std::getline(rl, spec);
799  }
800  rock_[i].read(dir, spec);
801  }
802  }
803 
804 
805 
806 
807 } // namespace Opm
808 
809 
810 #endif // OPENRS_RESERVOIRPROPERTYCOMMON_IMPL_HEADER
double ntg(int cell_index) const
Read-access to ntg.
Definition: ReservoirPropertyCommon_impl.hpp:384
void setDensities(double d1, double d2)
Set densitities of both faces.
Definition: ReservoirPropertyCommon_impl.hpp:343
double cflFactorGravity() const
A factor useful in gravity cfl computations.
Definition: ReservoirPropertyCommon_impl.hpp:452
double viscositySecondPhase() const
Viscosity of second (oil) phase.
Definition: ReservoirPropertyCommon_impl.hpp:357
PermTensor permeability(int cell_index) const
Read-access to permeability.
Definition: ReservoirPropertyCommon_impl.hpp:404
void assignRockTable(Opm::DeckConstPtr deck, const std::vector< int > &global_cell)
Definition: ReservoirPropertyCommon_impl.hpp:741
double densitySecondPhase() const
Density of second (oil) phase.
Definition: ReservoirPropertyCommon_impl.hpp:371
Definition: BlackoilFluid.hpp:31
double porosity(int cell_index) const
Read-access to porosity.
Definition: ReservoirPropertyCommon_impl.hpp:378
double swcr(int cell_index) const
Read-access to swcr.
Definition: ReservoirPropertyCommon_impl.hpp:390
double cflFactorCapillary() const
A factor useful in gravity cfl computations.
Definition: ReservoirPropertyCommon_impl.hpp:459
void readRocks(const std::string &rock_list_file)
Definition: ReservoirPropertyCommon_impl.hpp:784
double viscosityFirstPhase() const
Viscosity of first (water) phase.
Definition: ReservoirPropertyCommon_impl.hpp:350
SharedCMatrix SharedPermTensor
Tensor type for read and write access to permeability.
Definition: ReservoirPropertyCommon.hpp:66
double s_min(int c) const
Definition: ReservoirPropertyCommon_impl.hpp:491
PermeabilityKind
Enum for the kind of permeability field originally retrieved.
Definition: ReservoirPropertyCommon.hpp:50
double capillaryPressure(int cell_index, double saturation) const
Capillary pressure.
Definition: ReservoirPropertyCommon_impl.hpp:465
void init(Opm::DeckConstPtr deck, const std::vector< int > &global_cell, const double perm_threshold=0.0, const std::string *rock_list_filename=0, const bool use_jfunction_scaling=true, const double sigma=1.0, const double theta=0.0)
Initialize from a grdecl file.
Definition: ReservoirPropertyCommon_impl.hpp:274
double densityDifference() const
Difference of densities.
Definition: ReservoirPropertyCommon_impl.hpp:438
double densityFirstPhase() const
Density of first (water) phase.
Definition: ReservoirPropertyCommon_impl.hpp:364
void assignPermeability(Opm::DeckConstPtr deck, const std::vector< int > &global_cell, const double perm_threshold)
Definition: ReservoirPropertyCommon_impl.hpp:681
Definition: ReservoirPropertyCommon.hpp:50
double s_max(int c) const
Definition: ReservoirPropertyCommon_impl.hpp:504
double saturationFromCapillaryPressure(int cell_index, double cap_press) const
Inverse of the capillary pressure function.
Definition: ReservoirPropertyCommon_impl.hpp:517
Definition: ReservoirPropertyCommon.hpp:50
double cflFactor() const
A factor useful in cfl computations.
Definition: ReservoirPropertyCommon_impl.hpp:445
RPImpl & asImpl()
Definition: ReservoirPropertyCommon_impl.hpp:582
Definition: ReservoirPropertyCommon.hpp:50
void assignPorosity(Opm::DeckConstPtr deck, const std::vector< int > &global_cell)
Definition: ReservoirPropertyCommon_impl.hpp:591
void phaseDensities(int, Vector &density) const
Densities for both phases.
Definition: ReservoirPropertyCommon_impl.hpp:429
void assignSOWCR(Opm::DeckConstPtr deck, const std::vector< int > &global_cell)
Definition: ReservoirPropertyCommon_impl.hpp:658
ImmutableCMatrix PermTensor
Tensor type for read-only access to permeability.
Definition: ReservoirPropertyCommon.hpp:62
ReservoirPropertyCommon()
Default constructor.
Definition: ReservoirPropertyCommon_impl.hpp:256
Definition: ReservoirPropertyCommon.hpp:50
void assignNTG(Opm::DeckConstPtr deck, const std::vector< int > &global_cell)
Definition: ReservoirPropertyCommon_impl.hpp:613
void assignSWCR(Opm::DeckConstPtr deck, const std::vector< int > &global_cell)
Definition: ReservoirPropertyCommon_impl.hpp:635
void setViscosities(double v1, double v2)
Set viscosities of both faces.
Definition: ReservoirPropertyCommon_impl.hpp:336
Definition: ReservoirPropertyCommon.hpp:50
double sowcr(int cell_index) const
Read-access to sowcr.
Definition: ReservoirPropertyCommon_impl.hpp:396
void zero(Matrix &A)
Zero-fill a.
Definition: Matrix.hpp:603
SharedPermTensor permeabilityModifiable(int cell_index)
Read- and write-access to permeability. Use with caution.
Definition: ReservoirPropertyCommon_impl.hpp:415
void writeSintefLegacyFormat(const std::string &grid_prefix) const
Write permeability and porosity in the Sintef legacy format.
Definition: ReservoirPropertyCommon_impl.hpp:531
double capillaryPressureDeriv(int c, double s) const
Derivative of Capillary pressure.
Definition: ReservoirPropertyCommon_impl.hpp:478