vtkmultiwriter.hh 16.4 KB
Newer Older
1
2
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
3
/*****************************************************************************
4
 *   Copyright (C) 2008-2011 by Andreas Lauser                               *
Andreas Lauser's avatar
Andreas Lauser committed
5
 *   Institute for Modelling Hydraulic and Environmental Systems             *
6
7
8
 *   University of Stuttgart, Germany                                        *
 *   email: <givenname>.<name>@iws.uni-stuttgart.de                          *
 *                                                                           *
9
 *   This program is free software: you can redistribute it and/or modify    *
10
 *   it under the terms of the GNU General Public License as published by    *
11
12
 *   the Free Software Foundation, either version 2 of the License, or       *
 *   (at your option) any later version.                                     *
13
 *                                                                           *
14
15
16
17
18
19
20
 *   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/>.   *
21
22
23
 *****************************************************************************/
/*!
 * \file
24
 * \brief Simplifies writing multi-file VTK datasets.
25
26
27
28
 */
#ifndef VTK_MULTI_WRITER_HH
#define VTK_MULTI_WRITER_HH

29
30
#include "vtknestedfunction.hh"

31
32
33
#include <dune/common/fvector.hh>
#include <dune/istl/bvector.hh>

Andreas Lauser's avatar
Andreas Lauser committed
34
35
36
37
#include <dune/grid/io/file/vtk/vtkwriter.hh>

#include <dumux/common/valgrind.hh>

38
39
40
41
#if HAVE_MPI
#include <mpi.h>
#endif

42
43
44
45
#include <list>
#include <iostream>
#include <string>

46
47
#include <limits>

48
49
namespace Dumux {
/*!
50
 * \brief Simplifies writing multi-file VTK datasets.
51
52
53
 *
 * This class automatically keeps the meta file up to date and
 * simplifies writing datasets consisting of multiple files. (i.e.
54
 * multiple time steps or grid refinements within a time step.)
55
56
57
58
 */
template<class GridView>
class VtkMultiWriter
{
59
60
61
62
63
64
    enum { dim = GridView::dimension };

    typedef Dune::MultipleCodimMultipleGeomTypeMapper<GridView, Dune::MCMGVertexLayout> VertexMapper;
    typedef Dune::MultipleCodimMultipleGeomTypeMapper<GridView, Dune::MCMGElementLayout> ElementMapper;

public:
65
    typedef Dune::VTKWriter<GridView> VtkWriter;
Andreas Lauser's avatar
Andreas Lauser committed
66
67
    VtkMultiWriter(const GridView &gridView,
                   const std::string &simName = "",
68
69
70
71
                   std::string multiFileName = "")
        : gridView_(gridView)
        , elementMapper_(gridView)
        , vertexMapper_(gridView)
72
73
74
75
76
77
78
79
    {
        simName_ = (simName.empty())?"sim":simName;
        multiFileName_ = multiFileName;
        if (multiFileName_.empty()) {
            multiFileName_ = simName_;
            multiFileName_ += ".pvd";
        }

80
        curWriterNum_ = 0;
81

82
83
        commRank_ = 0;
        commSize_ = 1;
84
85
86
87
#if HAVE_MPI
        MPI_Comm_rank(MPI_COMM_WORLD, &commRank_);
        MPI_Comm_size(MPI_COMM_WORLD, &commSize_);
#endif
88
89
90
91
    }

    ~VtkMultiWriter()
    {
92
        finishMultiFile_();
93
94
95
96
97

        if (commRank_ == 0)
            multiFile_.close();
    }

98
99
100
101
102
103
    /*!
     * \brief Returns the number of the current VTK file.
     */
    int curWriterNum() const
    { return curWriterNum_; }

104
105
106
107
108
109
110
111
112
113
114
115
116
    /*!
     * \brief Updates the internal data structures after mesh
     *        refinement.
     *
     * If the grid changes between two calls of beginWrite(), this
     * method _must_ be called before the second beginWrite()!
     */
    void gridChanged()
    {
        elementMapper_.update();
        vertexMapper_.update();
    }

117
    /*!
118
     * \brief Called when ever a new time step or a new grid
119
120
     *        must be written.
     */
121
    void beginWrite(double t)
122
123
    {
        if (!multiFile_.is_open()) {
124
            startMultiFile_(multiFileName_);
125
126
127
        }


128
129
        curWriter_ = new VtkWriter(gridView_,
                                   Dune::VTKOptions::conforming);
130
        ++curWriterNum_;
131
132
133
134

        curTime_ = t;
        curOutFileName_ = fileName_();
    }
Andreas Lauser's avatar
Andreas Lauser committed
135

136
    /*!
137
138
139
140
     * \brief Allocate a managed buffer for a scalar field
     *
     * The buffer will be deleted automatically after the data has
     * been written by to disk.
141
     */
142
    template <class Scalar, int nComp>
143
    Dune::BlockVector<Dune::FieldVector<Scalar, nComp> > *allocateManagedBuffer(int nEntities)
144
145
    {
        typedef Dune::BlockVector<Dune::FieldVector<Scalar, nComp> > VectorField;
Andreas Lauser's avatar
Andreas Lauser committed
146

147
148
149
        ManagedVectorField_<VectorField> *vfs =
            new ManagedVectorField_<VectorField>(nEntities);
        managedObjects_.push_back(vfs);
150
151
152
        return &(vfs->vf);
    }

153
154
155
156
157
158
159
160
161
    // todo: remove these two functions as soon as we depend on a
    //       contemporary compilers which support default template
    //       arguments for function templates
    template <class Scalar>
    Dune::BlockVector<Dune::FieldVector<Scalar, 1> > *allocateManagedBuffer(int nEntities)
    { return allocateManagedBuffer<Scalar, 1>(nEntities); }
    Dune::BlockVector<Dune::FieldVector<double, 1> > *allocateManagedBuffer(int nEntities)
    { return allocateManagedBuffer<double, 1>(nEntities); }

162
163
    /*!
     * \brief Add a finished vertex centered vector field to the
Andreas Lauser's avatar
Andreas Lauser committed
164
165
     *        output.
     * \brief Add a vertex-centered quantity to the output.
166
167
168
169
     *
     * If the buffer is managed by the VtkMultiWriter, it must have
     * been created using createField() and may not be used by
     * anywhere after calling this method. After the data is written
Andreas Lauser's avatar
Andreas Lauser committed
170
     * to disk, it will be deleted automatically.
171
172
173
     *
     * If the buffer is not managed by the MultiWriter, the buffer
     * must exist at least until the call to endWrite()
Andreas Lauser's avatar
Andreas Lauser committed
174
     * finishes.
175
176
177
     *
     * In both cases, modifying the buffer between the call to this
     * method and endWrite() results in _undefined behaviour_.
178
     */
179
    template <class DataBuffer>
180
    void attachVertexData(DataBuffer &buf, std::string name, int nComps=1)
181
    {
182
183
        sanitizeBuffer_(buf, nComps);

184
        typedef typename VtkWriter::VTKFunctionPtr FunctionPtr;
185
        typedef Dumux::VtkNestedFunction<GridView, VertexMapper, DataBuffer> VtkFn;
186
        FunctionPtr fnPtr(new VtkFn(name,
187
                                    gridView_,
188
189
190
191
192
                                    vertexMapper_,
                                    buf,
                                    /*codim=*/dim,
                                    nComps));
        curWriter_->addVertexData(fnPtr);
193
    }
194

195
    /*!
196
197
198
199
200
     * \brief Add a cell centered quantity to the output.
     *
     * If the buffer is managed by the VtkMultiWriter, it must have
     * been created using createField() and may not be used by
     * anywhere after calling this method. After the data is written
Andreas Lauser's avatar
Andreas Lauser committed
201
     * to disk, it will be deleted automatically.
202
203
204
     *
     * If the buffer is not managed by the MultiWriter, the buffer
     * must exist at least until the call to endWrite()
Andreas Lauser's avatar
Andreas Lauser committed
205
     * finishes.
206
207
208
     *
     * In both cases, modifying the buffer between the call to this
     * method and endWrite() results in _undefined behaviour_.
209
     */
210
    template <class DataBuffer>
211
    void attachCellData(DataBuffer &buf, std::string name, int nComps = 1)
Andreas Lauser's avatar
Andreas Lauser committed
212
    {
213
214
        sanitizeBuffer_(buf, nComps);

215
        typedef typename VtkWriter::VTKFunctionPtr FunctionPtr;
216
        typedef Dumux::VtkNestedFunction<GridView, ElementMapper, DataBuffer> VtkFn;
217
        FunctionPtr fnPtr(new VtkFn(name,
218
                                    gridView_,
219
220
221
222
223
                                    elementMapper_,
                                    buf,
                                    /*codim=*/0,
                                    nComps));
        curWriter_->addCellData(fnPtr);
224
225
226
227
228
    }

    /*!
     * \brief Finalizes the current writer.
     *
229
230
231
     * This means that everything will be written to disk, except if
     * the onlyDiscard argument is true. In this case only all managed
     * buffers are deleted, but no output is written.
232
     */
233
    void endWrite(bool onlyDiscard = false)
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
    {
        if (!onlyDiscard) {
            curWriter_->write(curOutFileName_.c_str(),
                              Dune::VTKOptions::ascii);

            // determine name to write into the multi-file for the
            // current time step
            std::string fileName;
            std::string suffix = fileSuffix_();
            if (commSize_ == 1) {
                fileName = curOutFileName_;
                multiFile_.precision(16);
                multiFile_ << "   <DataSet timestep=\""
                           << curTime_
                           << "\" file=\""
                           << fileName << "." << suffix << "\"/>\n";
            }
251
            if (commSize_ > 1 && commRank_ == 0) {
252
253
254
255
256
257
258
259
260
261
262
263
264
                // only the first process updates the multi-file
                for (int part=0; part < commSize_; ++part) {
                    fileName = fileName_(part);
                    multiFile_.precision(16);
                    multiFile_ << "   <DataSet "
                               << " part=\"" << part << "\""
                               << " timestep=\"" << curTime_ << "\""
                               << " file=\"" << fileName << "." << suffix << "\"/>\n";
                }
            }

        }
        else
265
            -- curWriterNum_;
266

267
        // discard managed objects and the current VTK writer
268
        delete curWriter_;
269
270
271
        while (managedObjects_.begin() != managedObjects_.end()) {
            delete managedObjects_.front();
            managedObjects_.pop_front();
272
273
        }

274
275
276
277
        // temporarily write the closing XML mumbo-jumbo to the mashup
        // file so that the data set can be loaded even if the
        // simulation is aborted (or not yet finished)
        finishMultiFile_();
278
279
280
281
282
283
284
285
    };

    /*!
     * \brief Write the multi-writer's state to a restart file.
     */
    template <class Restarter>
    void serialize(Restarter &res)
    {
286
        res.serializeSectionBegin("VTKMultiWriter");
287
        res.serializeStream() << curWriterNum_ << "\n";
288

289
        if (commRank_ == 0) {
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
            size_t fileLen = 0;
            size_t filePos = 0;
            if (multiFile_.is_open()) {
                // write the meta file into the restart file
                filePos = multiFile_.tellp();
                multiFile_.seekp(0, std::ios::end);
                fileLen = multiFile_.tellp();
                multiFile_.seekp(filePos);
            }

            res.serializeStream() << fileLen << "  " << filePos << "\n";

            if (fileLen > 0) {
                std::ifstream multiFileIn(multiFileName_.c_str());
                char *tmp = new char[fileLen];
                multiFileIn.read(tmp, fileLen);
                res.serializeStream().write(tmp, fileLen);
                delete[] tmp;
            }
309
        }
310
311

        res.serializeSectionEnd();
312
313
314
315
316
317
318
319
    }

    /*!
     * \brief Read the multi-writer's state from a restart file.
     */
    template <class Restarter>
    void deserialize(Restarter &res)
    {
320
        res.deserializeSectionBegin("VTKMultiWriter");
321
        res.deserializeStream() >> curWriterNum_;
Andreas Lauser's avatar
Andreas Lauser committed
322

323
        if (commRank_ == 0) {
324
325
326
            std::string dummy;
            std::getline(res.deserializeStream(), dummy);

327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
            // recreate the meta file from the restart file
            size_t filePos, fileLen;
            res.deserializeStream() >> fileLen >> filePos;
            std::getline(res.deserializeStream(), dummy);
            if (multiFile_.is_open())
                multiFile_.close();

            if (fileLen > 0) {
                multiFile_.open(multiFileName_.c_str());

                char *tmp = new char[fileLen];
                res.deserializeStream().read(tmp, fileLen);
                multiFile_.write(tmp, fileLen);
                delete[] tmp;
            }

            multiFile_.seekp(filePos);
344
        }
345
346
347
348
        else {
            std::string tmp;
            std::getline(res.deserializeStream(), tmp);
        };
349
        res.deserializeSectionEnd();
350
351
352
353
354
    }

private:
    std::string fileName_()
    {
355
356
357
        std::ostringstream oss;
        oss << simName_ << "-" << std::setw(5) << std::setfill('0') << curWriterNum_;
        return oss.str();
358
359
360
361
362
    }

    std::string fileName_(int rank)
    {
        if (commSize_ > 1) {
363
364
365
366
367
368
369
370
            std::ostringstream oss;
            oss << "s" << std::setw(4) << std::setfill('0') << commSize_
                << ":p" << rank
                << ":" << simName_ << "-"
                << std::setw(5) << curWriterNum_;
            return oss.str();

            return oss.str();
371
372
        }
        else {
373
374
375
            std::ostringstream oss;
            oss << simName_ << "-" << std::setw(5) << std::setfill('0') << curWriterNum_;
            return oss.str();
376
377
378
379
380
381
382
383
384
        }
    }

    std::string fileSuffix_()
    {
        return (GridView::dimension == 1)?"vtp":"vtu";
    }


385
    void startMultiFile_(const std::string &multiFileName)
386
387
388
    {
        // only the first process writes to the multi-file
        if (commRank_ == 0) {
389
            // generate one meta vtk-file holding the individual time steps
390
391
392
393
394
395
396
397
398
399
400
            multiFile_.open(multiFileName.c_str());
            multiFile_ <<
                "<?xml version=\"1.0\"?>\n"
                "<VTKFile type=\"Collection\"\n"
                "         version=\"0.1\"\n"
                "         byte_order=\"LittleEndian\"\n"
                "         compressor=\"vtkZLibDataCompressor\">\n"
                " <Collection>\n";
        }
    }

401
    void finishMultiFile_()
402
403
404
405
406
407
408
409
410
411
412
413
414
    {
        // only the first process writes to the multi-file
        if (commRank_ == 0) {
            // make sure that we always have a working meta file
            std::ofstream::pos_type pos = multiFile_.tellp();
            multiFile_ <<
                " </Collection>\n"
                "</VTKFile>\n";
            multiFile_.seekp(pos);
            multiFile_.flush();
        }
    }

415
416
    // make sure the field is well defined if running under valgrind
    // and make sure that all values can be displayed by paraview
417
418
    template <class DataBuffer>
    void sanitizeBuffer_(DataBuffer &b, int nComps)
419
    {
420
        for (unsigned int i = 0; i < b.size(); ++i) {
421
422
423
424
425
            for (int j = 0; j < nComps; ++j) {
                Valgrind::CheckDefined(b[i][j]);

                // set values which are too small to 0 to avoid
                // problems with paraview
Andreas Lauser's avatar
Andreas Lauser committed
426
                if (std::abs(b[i][j])
427
428
429
430
431
432
433
434
                    < std::numeric_limits<float>::min())
                {
                    b[i][j] = 0.0;
                }
            }
        }
    }

435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
    //////////////////////////////
    // HACK: when ever we attach some data we need to copy the
    //       vector field (that's because Dune::VTKWriter is not
    //       able to write fields one at a time and using
    //       VTKWriter::add*Data doesn't copy the data's
    //       representation so that once we arrive at the point
    //       where we want to write the data to disk, it might
    //       exist anymore). The problem we encounter there is
    //       that add*Data accepts arbitrary types as vector
    //       fields, but we need a single type for the linked list
    //       which keeps track of the data added. The trick we use
    //       here is to define a non-template base class with a
    //       virtual destructor for the type given to the linked
    //       list and a derived template class which actually
    //       knows the type of the vector field it must delete.

    /** \todo Please doc me! */

453
    class ManagedObject_
454
455
    {
    public:
456
        virtual ~ManagedObject_()
457
458
459
460
461
462
        {}
    };

    /** \todo Please doc me! */

    template <class VF>
463
    class ManagedVectorField_ : public ManagedObject_
464
465
    {
    public:
466
        ManagedVectorField_(int size)
467
468
469
470
471
472
473
            : vf(size)
        { }
        VF vf;
    };
    // end hack
    ////////////////////////////////////

474
475
476
    const GridView gridView_;
    ElementMapper elementMapper_;
    VertexMapper vertexMapper_;
477

478
479
480
    std::string simName_;
    std::ofstream multiFile_;
    std::string multiFileName_;
481
482
483
484

    int commSize_; // number of processes in the communicator
    int commRank_; // rank of the current process in the communicator

485
    VtkWriter *curWriter_;
486
487
    double curTime_;
    std::string curOutFileName_;
488
    int curWriterNum_;
489

490
    std::list<ManagedObject_*> managedObjects_;
491
492
493
494
};
}

#endif