SimulatorTesterFlexibleBC.hpp
Go to the documentation of this file.
1//===========================================================================
2//
3// File: SimulatorTesterFlexibleBC.hpp
4//
5// Created: Wed Mar 24 11:14:12 2010
6//
7// Author(s): Atgeirr F Rasmussen <atgeirr@sintef.no>
8// Jostein R Natvig <jostein.r.natvig@sintef.no>
9//
10// $Date$
11//
12// $Revision$
13//
14//===========================================================================
15
16/*
17 Copyright 2010 SINTEF ICT, Applied Mathematics.
18 Copyright 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_SIMULATORTESTERFLEXIBLEBC_HEADER
37#define OPENRS_SIMULATORTESTERFLEXIBLEBC_HEADER
38
39#include "SimulatorTester.hpp"
40
41#include <iostream>
42
43namespace Opm
44{
45
46 template <class SimTraits>
48 {
49 protected:
51 typedef typename Super::GridInterface GI;
52 typedef typename Super::Vector Vector;
53
54 virtual void initSources(const Opm::parameter::ParameterGroup& param)
55 {
56 // Zero-initializing first.
57 int nc = this->ginterf_.numberOfCells();
58 this->injection_rates_ = Opm::SparseVector<double>(nc);
59 this->injection_rates_psolver_.resize(nc, 0.0);
60
61// this->injection_rates_.addElement(1.0, 0);
62// this->injection_rates_psolver_[0] = 1.0;
63// this->injection_rates_.addElement(-1.0, nc-1);
64// this->injection_rates_psolver_[nc-1] = -1.0;
65
66 // Initializing blocks.
67 double total_inj = 0.0;
68 bool has_inj_block = param.getDefault("has_inj_block", false);
69 if (has_inj_block) {
70 Vector low;
71 low[0] = param.getDefault("inj_block_low_x", 0.0);
72 low[1] = param.getDefault("inj_block_low_y", 0.0);
73 low[2] = param.getDefault("inj_block_low_z", 0.0);
74 Vector high;
75 high[0] = param.getDefault("inj_block_high_x", 1.0);
76 high[1] = param.getDefault("inj_block_high_y", 1.0);
77 high[2] = param.getDefault("inj_block_high_z", 1.0);
78 double inj_block_density = param.get<double>("inj_block_density");
79 total_inj = setSourceBlock(low, high, inj_block_density, true);
80 }
81 double total_prod = 0.0;
82 bool has_prod_block = param.getDefault("has_prod_block", false);
83 if (has_prod_block) {
84 Vector low;
85 low[0] = param.getDefault("prod_block_low_x", 0.0);
86 low[1] = param.getDefault("prod_block_low_y", 0.0);
87 low[2] = param.getDefault("prod_block_low_z", 0.0);
88 Vector high;
89 high[0] = param.getDefault("prod_block_high_x", 1.0);
90 high[1] = param.getDefault("prod_block_high_y", 1.0);
91 high[2] = param.getDefault("prod_block_high_z", 1.0);
92 double prod_block_density = param.get<double>("prod_block_density");
93 total_prod = setSourceBlock(low, high, prod_block_density, false);
94 }
95 if (has_inj_block || has_prod_block) {
96 std::cout << "Initialized injectors with total rate: " << total_inj
97 << "\nInitialized producers with total rate: " << total_prod
98 << std::endl;
99 }
100 }
101
102 virtual void initBoundaryConditions(const Opm::parameter::ParameterGroup& param)
103 {
104 setupBoundaryConditions(param, this->ginterf_, this->bcond_);
105 }
106
107 private:
108 bool isInside(const Vector& low, const Vector& high, const Vector& pt)
109 {
110 return low[0] < pt[0]
111 && low[1] < pt[1]
112 && low[2] < pt[2]
113 && high[0] > pt[0]
114 && high[1] > pt[1]
115 && high[2] > pt[2];
116 }
117
118 double setSourceBlock(const Vector& low, const Vector& high, double density, bool is_injection)
119 {
120 typedef typename GI::CellIterator CI;
121 // Iterate over all cells, if the centroid of a cell is in the
122 // domain given, set a source term. Accumulate total source terms.
123 double total_rate = 0.0;
124 for (CI c = this->ginterf_.cellbegin(); c != this->ginterf_.cellend(); ++c) {
125 if (isInside(low, high, c->centroid())) {
126 int index = c->index();
127 double rate = density*c->volume();
128 if (!is_injection) {
129 rate = -rate;
130 }
131 total_rate += rate;
132 this->injection_rates_.addElement(rate, index);
133 this->injection_rates_psolver_[index] = rate;
134 }
135 }
136 return total_rate;
137 }
138
139 };
140
141} // namespace Opm
142
143#endif // OPENRS_SIMULATORTESTERFLEXIBLEBC_HEADER
GIE::CellIterator< InterfaceType > CellIterator
Definition: GridInterfaceEuler.hpp:417
CellIterator cellbegin() const
Definition: GridInterfaceEuler.hpp:445
int numberOfCells() const
Definition: GridInterfaceEuler.hpp:453
Opm::SparseVector< double > injection_rates_
Definition: SimulatorBase.hpp:137
std::vector< double > injection_rates_psolver_
Definition: SimulatorBase.hpp:138
GridInterface ginterf_
Definition: SimulatorBase.hpp:134
BCs bcond_
Definition: SimulatorBase.hpp:136
Dune::FieldVector< double, Dimension > Vector
Definition: SimulatorBase.hpp:116
Definition: SimulatorTesterFlexibleBC.hpp:48
virtual void initBoundaryConditions(const Opm::parameter::ParameterGroup &param)
Definition: SimulatorTesterFlexibleBC.hpp:102
Super::GridInterface GI
Definition: SimulatorTesterFlexibleBC.hpp:51
SimulatorTester< SimTraits > Super
Definition: SimulatorTesterFlexibleBC.hpp:50
virtual void initSources(const Opm::parameter::ParameterGroup &param)
Definition: SimulatorTesterFlexibleBC.hpp:54
Super::Vector Vector
Definition: SimulatorTesterFlexibleBC.hpp:52
Definition: SimulatorTester.hpp:54
Definition: BlackoilFluid.hpp:32
void setupBoundaryConditions(const Opm::parameter::ParameterGroup &param, const GridInterface &g, BCs &bcs)
Setup boundary conditions for a simulation. It is assumed that the boundary ids are 1-6,...
Definition: setupBoundaryConditions.hpp:51