2016-01-05 12:32:36 -06:00
|
|
|
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
|
|
|
|
// vi: set et ts=4 sw=4 sts=4:
|
|
|
|
/*
|
|
|
|
This file is part of the Open Porous Media project (OPM).
|
|
|
|
|
|
|
|
OPM is free software: you can redistribute it and/or modify
|
|
|
|
it under the terms of the GNU General Public License as published by
|
|
|
|
the Free Software Foundation, either version 2 of the License, or
|
|
|
|
(at your option) any later version.
|
|
|
|
|
|
|
|
OPM is distributed in the hope that it will be useful,
|
|
|
|
but WITHOUT ANY WARRANTY; without even the implied warranty of
|
|
|
|
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
|
|
|
GNU General Public License for more details.
|
|
|
|
|
|
|
|
You should have received a copy of the GNU General Public License
|
|
|
|
along with OPM. If not, see <http://www.gnu.org/licenses/>.
|
2016-03-14 07:21:47 -05:00
|
|
|
|
|
|
|
Consult the COPYING file in the top-level source directory of this
|
|
|
|
module for the precise wording of the license and the list of
|
|
|
|
copyright holders.
|
2016-01-05 12:32:36 -06:00
|
|
|
*/
|
|
|
|
/*!
|
|
|
|
* \file
|
2019-09-05 10:04:39 -05:00
|
|
|
* \copydoc Opm::EclAluGridVanguard
|
2016-01-05 12:32:36 -06:00
|
|
|
*/
|
2018-02-01 09:26:58 -06:00
|
|
|
#ifndef EWOMS_ECL_ALU_GRID_VANGUARD_HH
|
|
|
|
#define EWOMS_ECL_ALU_GRID_VANGUARD_HH
|
2016-01-05 12:32:36 -06:00
|
|
|
|
2018-02-01 09:26:58 -06:00
|
|
|
#include "eclbasevanguard.hh"
|
2016-01-05 12:32:36 -06:00
|
|
|
#include "alucartesianindexmapper.hh"
|
|
|
|
|
|
|
|
#include <dune/alugrid/grid.hh>
|
|
|
|
#include <dune/alugrid/common/fromtogridfactory.hh>
|
2018-02-08 05:20:17 -06:00
|
|
|
#include <opm/grid/CpGrid.hpp>
|
2016-01-05 12:32:36 -06:00
|
|
|
|
2019-09-05 10:04:39 -05:00
|
|
|
namespace Opm {
|
2016-01-05 12:32:36 -06:00
|
|
|
template <class TypeTag>
|
2018-02-01 09:26:58 -06:00
|
|
|
class EclAluGridVanguard;
|
2016-01-05 12:32:36 -06:00
|
|
|
|
2019-09-05 10:04:39 -05:00
|
|
|
} // namespace Opm
|
2018-06-14 09:06:05 -05:00
|
|
|
|
|
|
|
BEGIN_PROPERTIES
|
|
|
|
|
2018-02-01 09:26:58 -06:00
|
|
|
NEW_TYPE_TAG(EclAluGridVanguard, INHERITS_FROM(EclBaseVanguard));
|
2016-01-05 12:32:36 -06:00
|
|
|
|
|
|
|
// declare the properties
|
2019-09-05 10:04:39 -05:00
|
|
|
SET_TYPE_PROP(EclAluGridVanguard, Vanguard, Opm::EclAluGridVanguard<TypeTag>);
|
2018-02-01 09:26:58 -06:00
|
|
|
SET_TYPE_PROP(EclAluGridVanguard, Grid, Dune::ALUGrid<3, 3, Dune::cube, Dune::nonconforming>);
|
|
|
|
SET_TYPE_PROP(EclAluGridVanguard, EquilGrid, Dune::CpGrid);
|
2018-06-14 09:06:05 -05:00
|
|
|
|
|
|
|
END_PROPERTIES
|
|
|
|
|
2019-09-05 10:04:39 -05:00
|
|
|
namespace Opm {
|
2016-01-05 12:32:36 -06:00
|
|
|
|
|
|
|
/*!
|
|
|
|
* \ingroup EclBlackOilSimulator
|
|
|
|
*
|
|
|
|
* \brief Helper class for grid instantiation of ECL file-format using problems.
|
|
|
|
*
|
|
|
|
* This class uses Dune::ALUGrid as the simulation grid.
|
|
|
|
*/
|
|
|
|
template <class TypeTag>
|
2018-02-01 09:26:58 -06:00
|
|
|
class EclAluGridVanguard : public EclBaseVanguard<TypeTag>
|
2016-01-05 12:32:36 -06:00
|
|
|
{
|
2018-02-01 09:26:58 -06:00
|
|
|
friend class EclBaseVanguard<TypeTag>;
|
|
|
|
typedef EclBaseVanguard<TypeTag> ParentType;
|
2016-01-05 12:32:36 -06:00
|
|
|
|
|
|
|
typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
|
|
|
|
typedef typename GET_PROP_TYPE(TypeTag, Simulator) Simulator;
|
|
|
|
|
|
|
|
public:
|
|
|
|
typedef typename GET_PROP_TYPE(TypeTag, Grid) Grid;
|
|
|
|
typedef typename GET_PROP_TYPE(TypeTag, EquilGrid) EquilGrid;
|
|
|
|
typedef typename GET_PROP_TYPE(TypeTag, GridView) GridView;
|
|
|
|
|
2016-02-26 12:11:05 -06:00
|
|
|
private:
|
2019-09-05 10:04:39 -05:00
|
|
|
typedef Opm::AluCartesianIndexMapper<Grid> CartesianIndexMapper;
|
2016-01-17 14:15:18 -06:00
|
|
|
typedef Dune::CartesianIndexMapper<EquilGrid> EquilCartesianIndexMapper;
|
2016-01-05 12:32:36 -06:00
|
|
|
|
|
|
|
static const int dimension = Grid::dimension;
|
|
|
|
|
|
|
|
public:
|
2019-10-01 14:18:17 -05:00
|
|
|
EclAluGridVanguard(Simulator& simulator)
|
|
|
|
: EclBaseVanguard<TypeTag>(simulator)
|
|
|
|
{
|
|
|
|
this->callImplementationInit();
|
|
|
|
}
|
2016-01-05 12:32:36 -06:00
|
|
|
|
2018-02-01 09:26:58 -06:00
|
|
|
~EclAluGridVanguard()
|
2016-01-05 12:32:36 -06:00
|
|
|
{
|
|
|
|
delete cartesianIndexMapper_;
|
2016-01-17 14:15:18 -06:00
|
|
|
delete equilCartesianIndexMapper_;
|
2016-01-05 12:32:36 -06:00
|
|
|
delete grid_;
|
|
|
|
delete equilGrid_;
|
|
|
|
}
|
|
|
|
|
|
|
|
/*!
|
|
|
|
* \brief Return a reference to the simulation grid.
|
|
|
|
*/
|
|
|
|
Grid& grid()
|
|
|
|
{ return *grid_; }
|
|
|
|
|
|
|
|
/*!
|
|
|
|
* \brief Return a reference to the simulation grid.
|
|
|
|
*/
|
|
|
|
const Grid& grid() const
|
|
|
|
{ return *grid_; }
|
|
|
|
|
|
|
|
/*!
|
|
|
|
* \brief Returns a refefence to the grid which should be used by the EQUIL
|
|
|
|
* initialization code.
|
|
|
|
*
|
|
|
|
* The EQUIL keyword is used to specify the initial condition of the reservoir in
|
|
|
|
* hydrostatic equilibrium. Since the code which does this is not accepting arbitrary
|
|
|
|
* DUNE grids (the code is part of the opm-core module), this is not necessarily the
|
|
|
|
* same as the grid which is used for the actual simulation.
|
|
|
|
*/
|
|
|
|
const EquilGrid& equilGrid() const
|
|
|
|
{ return *equilGrid_; }
|
|
|
|
|
|
|
|
/*!
|
|
|
|
* \brief Indicates that the initial condition has been computed and the memory used
|
|
|
|
* by the EQUIL grid can be released.
|
|
|
|
*
|
|
|
|
* Depending on the implementation, subsequent accesses to the EQUIL grid lead to
|
|
|
|
* crashes.
|
|
|
|
*/
|
|
|
|
void releaseEquilGrid()
|
|
|
|
{
|
2016-01-17 14:15:18 -06:00
|
|
|
delete equilCartesianIndexMapper_;
|
|
|
|
equilCartesianIndexMapper_ = 0;
|
|
|
|
|
2016-01-05 12:32:36 -06:00
|
|
|
delete equilGrid_;
|
|
|
|
equilGrid_ = 0;
|
|
|
|
}
|
|
|
|
|
|
|
|
/*!
|
|
|
|
* \brief Distribute the simulation grid over multiple processes
|
|
|
|
*
|
|
|
|
* (For parallel simulation runs.)
|
|
|
|
*/
|
|
|
|
void loadBalance()
|
|
|
|
{
|
|
|
|
auto gridView = grid().leafGridView();
|
|
|
|
auto dataHandle = cartesianIndexMapper_->dataHandle(gridView);
|
|
|
|
grid().loadBalance(*dataHandle);
|
|
|
|
|
|
|
|
// communicate non-interior cells values
|
|
|
|
grid().communicate(*dataHandle,
|
|
|
|
Dune::InteriorBorder_All_Interface,
|
|
|
|
Dune::ForwardCommunication );
|
|
|
|
}
|
|
|
|
|
|
|
|
/*!
|
|
|
|
* \brief Returns the object which maps a global element index of the simulation grid
|
|
|
|
* to the corresponding element index of the logically Cartesian index.
|
|
|
|
*/
|
|
|
|
const CartesianIndexMapper& cartesianIndexMapper() const
|
|
|
|
{ return *cartesianIndexMapper_; }
|
|
|
|
|
2016-01-17 14:15:18 -06:00
|
|
|
/*!
|
|
|
|
* \brief Returns mapper from compressed to cartesian indices for the EQUIL grid
|
|
|
|
*/
|
|
|
|
const EquilCartesianIndexMapper& equilCartesianIndexMapper() const
|
|
|
|
{ return *equilCartesianIndexMapper_; }
|
|
|
|
|
2016-01-05 12:32:36 -06:00
|
|
|
protected:
|
|
|
|
void createGrids_()
|
|
|
|
{
|
2017-06-22 08:55:01 -05:00
|
|
|
const auto& gridProps = this->eclState().get3DProperties();
|
2016-11-07 08:14:07 -06:00
|
|
|
const std::vector<double>& porv = gridProps.getDoubleGridProperty("PORV").getData();
|
2016-01-05 12:32:36 -06:00
|
|
|
|
|
|
|
// we use separate grid objects: one for the calculation of the initial condition
|
|
|
|
// via EQUIL and one for the actual simulation. The reason is that the EQUIL code
|
|
|
|
// cannot cope with arbitrary Dune grids and is also allergic to distributed
|
|
|
|
// grids.
|
|
|
|
|
|
|
|
/////
|
|
|
|
// create the EQUIL grid
|
|
|
|
/////
|
|
|
|
equilGrid_ = new EquilGrid();
|
2017-06-22 08:55:01 -05:00
|
|
|
equilGrid_->processEclipseFormat(this->eclState().getInputGrid(),
|
2016-01-05 12:32:36 -06:00
|
|
|
/*isPeriodic=*/false,
|
|
|
|
/*flipNormals=*/false,
|
|
|
|
/*clipZ=*/false,
|
|
|
|
porv);
|
|
|
|
|
|
|
|
cartesianCellId_ = equilGrid_->globalCell();
|
|
|
|
|
|
|
|
for (unsigned i = 0; i < dimension; ++i)
|
|
|
|
cartesianDimension_[i] = equilGrid_->logicalCartesianSize()[i];
|
|
|
|
|
2016-01-17 14:15:18 -06:00
|
|
|
equilCartesianIndexMapper_ = new EquilCartesianIndexMapper(*equilGrid_);
|
|
|
|
|
|
|
|
/////
|
|
|
|
// create the simulation grid
|
|
|
|
/////
|
2016-01-05 12:32:36 -06:00
|
|
|
Dune::FromToGridFactory<Grid> factory;
|
|
|
|
grid_ = factory.convert(*equilGrid_, cartesianCellId_);
|
|
|
|
|
|
|
|
cartesianIndexMapper_ =
|
|
|
|
new CartesianIndexMapper(*grid_, cartesianDimension_, cartesianCellId_);
|
|
|
|
}
|
|
|
|
|
2018-06-28 08:49:45 -05:00
|
|
|
void filterConnections_()
|
2018-02-07 04:21:57 -06:00
|
|
|
{
|
|
|
|
// not handling the removal of completions for this type of grid yet.
|
|
|
|
}
|
|
|
|
|
2016-01-17 14:15:18 -06:00
|
|
|
Grid* grid_;
|
|
|
|
EquilGrid* equilGrid_;
|
2016-01-05 12:32:36 -06:00
|
|
|
std::vector<int> cartesianCellId_;
|
|
|
|
std::array<int,dimension> cartesianDimension_;
|
2016-01-17 14:15:18 -06:00
|
|
|
CartesianIndexMapper* cartesianIndexMapper_;
|
|
|
|
EquilCartesianIndexMapper* equilCartesianIndexMapper_;
|
2016-01-05 12:32:36 -06:00
|
|
|
};
|
|
|
|
|
2019-09-05 10:04:39 -05:00
|
|
|
} // namespace Opm
|
2016-01-05 12:32:36 -06:00
|
|
|
|
|
|
|
#endif
|