gstatrandomfield.hh 7.48 KB
Newer Older
1
2
3
4
5
6
7
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*****************************************************************************
 *   See the file COPYING for full copying permissions.                      *
 *                                                                           *
 *   This program is free software: you can redistribute it and/or modify    *
 *   it under the terms of the GNU General Public License as published by    *
8
 *   the Free Software Foundation, either version 3 of the License, or       *
9
10
11
12
13
14
15
16
17
18
19
20
 *   (at your option) any later version.                                     *
 *                                                                           *
 *   This program 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 this program.  If not, see <http://www.gnu.org/licenses/>.   *
 *****************************************************************************/
/*!
 * \file
21
 * \ingroup SpatialParameters
22
23
24
25
26
27
28
29
 * \brief Creating random fields using gstat
 */
#ifndef DUMUX_GSTAT_RANDOM_FIELD_HH
#define DUMUX_GSTAT_RANDOM_FIELD_HH

#include <iostream>

#include <dune/grid/common/mcmgmapper.hh>
30
#include <dune/grid/io/file/vtk.hh>
31

32
namespace Dumux {
33
34

/*!
35
 * \ingroup SpatialParameters
36
37
38
 * \brief Creating random fields using gstat
 *
 * gstat is an open source software tool which can (among other things) generate
39
40
 * geostatistical random fields (see <a href="http://www.gstat.org">http://www.gstat.org</a>
 * or \cite Pebesma1998a).
41
 *
42
 * To use this class, execute the installexternal.py from your DuMuX root
43
44
45
 * directory or donwload, unpack and install the tarball from the gstat-website.
 * Then rerun cmake (in the second case set GSTAT_ROOT in your input file to the
 * path where gstat is installed).
46
47
48
49
 */
template<class GridView, class Scalar>
class GstatRandomField
{
50
51
    enum { dim = GridView::dimension };
    enum { dimWorld = GridView::dimensionworld };
52

53
54
    using DataVector = std::vector<Scalar>;
    using Element = typename GridView::Traits::template Codim<0>::Entity;
55
    using ElementMapper = Dune::MultipleCodimMultipleGeomTypeMapper<GridView>;
56
57

public:
58
59
60
    // Add field types if you want to implement e.g. tensor permeabilities.
    enum FieldType { scalar, log10 };

61
    /*!
Bilal's avatar
Bilal committed
62
     * \brief Constructor
63
64
     *
     * \param gridView the used gridView
65
     * \param elementMapper	Maps elements of the given grid view
66
     */
67
68
69
70
71
    GstatRandomField(const GridView& gridView, const ElementMapper& elementMapper)
    : gridView_(gridView)
    , elementMapper_(elementMapper)
    , data_(gridView.size(0))
    {}
72

73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
    /*!
     * \brief Creates a new field with random variables, if desired.
     * Otherwise creates a data field from already available data.
     * For the random field generation three files are necessary.
     *
     * A \a gstatControlFile in which all commands and in/output files for gstat are specified.
     * A \a gstatInputFile contains all coordinates (cell centers) of the grid, so that
     * gstat can perform its random realization. The filename must be same as in the gstatControlFile.
     * A \a gstatOutputFile in which gstat writes the random values to this file.
     * The filename must be the same as in the gstatControlFile.
     * \param fieldType
     * \param gstatControlFile name of control file for gstat
     * \param gstatInputFile name of input file for gstat
     * \param gstatOutputFile name of the gstat output file
     * \param createNew set true to create a new field
     */
89
90
91
92
93
    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)
94
    {
95
        fieldType_ = fieldType;
96
97
        if (createNew)
        {
98
99
100
101
102
#if !HAVE_GSTAT
            DUNE_THROW(Dune::InvalidStateException, "Requested data field generation with gstat"
              << " but gstat was not found on your system. Set GSTAT_ROOT to the path where gstat "
              << " is installed and pass it to CMake, e.g. through an opts file.");
#else
103
            std::ofstream gstatInput(gstatInputFile);
104
            for (const auto& element : elements(gridView_))
105
                // get global coordinates of cell centers
106
                gstatInput << element.geometry().center() << std::endl;
107
108
109
110
111
112
113
114
115
116

            gstatInput.close();

            std::string syscom;
            syscom = GSTAT_EXECUTABLE;
            syscom += " ";
            syscom += gstatControlFile;

            if (!gstatInput.good())
            {
117
118
                DUNE_THROW(Dune::IOError, "Reading the gstat control file: "
                             << gstatControlFile << " failed." << std::endl);
119
120
121
122
123
124
125
            }

            if (system(syscom.c_str()))
            {
                DUNE_THROW(Dune::IOError, "Executing gstat failed.");
            }
#endif
126
        }
127
128
129
130

        std::ifstream gstatOutput(gstatOutputFile);
        if (!gstatOutput.good())
        {
131
132
            DUNE_THROW(Dune::IOError, "Reading from file: "
                         << gstatOutputFile << " failed." << std::endl);
133
134
135
136
137
        }

        std::string line;
        std::getline(gstatOutput, line);

138
139
        Scalar trash, dataValue;
        for (const auto& element : elements(gridView_))
140
141
142
143
        {
            std::getline(gstatOutput, line);
            std::istringstream curLine(line);
            if (dim == 1)
144
                curLine >> trash >> dataValue;
145
            else if (dim == 2)
146
147
148
149
150
151
152
                curLine >> trash >> trash >> dataValue;
            else if (dim == 3)
                curLine >> trash >> trash >> trash >> dataValue;
            else
                DUNE_THROW(Dune::InvalidStateException, "Invalid dimension " << dim);

            data_[elementMapper_.index(element)] = dataValue;
153
154
155
        }
        gstatOutput.close();

156
        // post processing
157
        using std::pow;
158
        if (fieldType_ == FieldType::log10)
159
          std::for_each(data_.begin(), data_.end(), [](Scalar& s){ s = pow(10.0, s); });
160
161
162
    }

    //! \brief Return an entry of the data vector
163
    Scalar data(const Element& e) const
164
165
166
167
    {
        return data_[elementMapper_.index(e)];
    }

168
169
170
171
172
173
    //! \brief Return the data vector for analysis or external vtk output
    const DataVector& data() const
    {
        return data_;
    }

174
    //! \brief Write the data to a vtk file
175
176
    void writeVtk(const std::string& vtkName,
                  const std::string& dataName = "data") const
177
178
179
180
    {
        Dune::VTKWriter<GridView> vtkwriter(gridView_);
        vtkwriter.addCellData(data_, dataName);

181
182
        DataVector logPerm;
        if (fieldType_ == FieldType::log10)
183
        {
184
            logPerm = data_;
185
186
            using std::log10;
            std::for_each(logPerm.begin(), logPerm.end(), [](Scalar& s){ s = log10(s); });
187
            vtkwriter.addCellData(logPerm, "log10 of " + dataName);
188
189
190
191
192
        }
        vtkwriter.write(vtkName, Dune::VTK::OutputType::ascii);
    }

private:
193
194
    const GridView gridView_;
    const ElementMapper& elementMapper_;
195
196
    DataVector data_;
    FieldType fieldType_;
197
198
};

199
} // end namespace Dumux
200
201

#endif