24#ifndef DUMUX_GSTAT_RANDOM_FIELD_HH
25#define DUMUX_GSTAT_RANDOM_FIELD_HH
29#include <dune/grid/common/mcmgmapper.hh>
30#include <dune/grid/io/file/vtk.hh>
47template<
class Gr
idView,
class Scalar>
50 enum { dim = GridView::dimension };
51 enum { dimWorld = GridView::dimensionworld };
53 using DataVector = std::vector<Scalar>;
54 using Element =
typename GridView::Traits::template Codim<0>::Entity;
55 using ElementMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
69 , elementMapper_(elementMapper)
70 , data_(gridView.size(0))
89 void create(
const std::string& gstatControlFile,
90 const std::string& gstatInputFile =
"gstatInput.txt",
91 const std::string& gstatOutputFile =
"permeab.dat",
93 bool createNew =
true)
95 fieldType_ = fieldType;
99 DUNE_THROW(Dune::InvalidStateException,
"Requested data field generation with gstat"
100 <<
" but gstat was not found on your system. Set GSTAT_ROOT to the path where gstat "
101 <<
" is installed and pass it to CMake, e.g. through an opts file.");
103 std::ofstream gstatInput(gstatInputFile);
104 for (
const auto& element : elements(gridView_))
106 gstatInput << element.geometry().center() << std::endl;
111 syscom = GSTAT_EXECUTABLE;
113 syscom += gstatControlFile;
115 if (!gstatInput.good())
117 DUNE_THROW(Dune::IOError,
"Reading the gstat control file: "
118 << gstatControlFile <<
" failed." << std::endl);
121 if (system(syscom.c_str()))
123 DUNE_THROW(Dune::IOError,
"Executing gstat failed.");
128 std::ifstream gstatOutput(gstatOutputFile);
129 if (!gstatOutput.good())
131 DUNE_THROW(Dune::IOError,
"Reading from file: "
132 << gstatOutputFile <<
" failed." << std::endl);
136 std::getline(gstatOutput, line);
138 Scalar trash, dataValue;
139 for (
const auto& element : elements(gridView_))
141 std::getline(gstatOutput, line);
142 std::istringstream curLine(line);
144 curLine >> trash >> dataValue;
146 curLine >> trash >> trash >> dataValue;
148 curLine >> trash >> trash >> trash >> dataValue;
150 DUNE_THROW(Dune::InvalidStateException,
"Invalid dimension " << dim);
152 data_[elementMapper_.index(element)] = dataValue;
158 if (fieldType_ == FieldType::log10)
159 std::for_each(data_.begin(), data_.end(), [](Scalar& s){ s = pow(10.0, s); });
163 Scalar
data(
const Element& e)
const
165 return data_[elementMapper_.index(e)];
176 const std::string& dataName =
"data")
const
178 Dune::VTKWriter<GridView> vtkwriter(gridView_);
179 vtkwriter.addCellData(data_, dataName);
182 if (fieldType_ == FieldType::log10)
186 std::for_each(logPerm.begin(), logPerm.end(), [](Scalar& s){ s = log10(s); });
187 vtkwriter.addCellData(logPerm,
"log10 of " + dataName);
189 vtkwriter.write(vtkName, Dune::VTK::OutputType::ascii);
193 const GridView gridView_;
194 const ElementMapper& elementMapper_;
Creating random fields using gstat.
Definition: gstatrandomfield.hh:49
void writeVtk(const std::string &vtkName, const std::string &dataName="data") const
Write the data to a vtk file.
Definition: gstatrandomfield.hh:175
Scalar data(const Element &e) const
Return an entry of the data vector.
Definition: gstatrandomfield.hh:163
void create(const std::string &gstatControlFile, const std::string &gstatInputFile="gstatInput.txt", const std::string &gstatOutputFile="permeab.dat", FieldType fieldType=FieldType::scalar, bool createNew=true)
Creates a new field with random variables, if desired. Otherwise creates a data field from already av...
Definition: gstatrandomfield.hh:89
GstatRandomField(const GridView &gridView, const ElementMapper &elementMapper)
Constructor.
Definition: gstatrandomfield.hh:67
const DataVector & data() const
Return the data vector for analysis or external vtk output.
Definition: gstatrandomfield.hh:169
FieldType
Definition: gstatrandomfield.hh:59
@ scalar
Definition: gstatrandomfield.hh:59
@ log10
Definition: gstatrandomfield.hh:59