kraichnan_field.cpp 15.5 KB
Newer Older
Lukas Bentkamp's avatar
Lukas Bentkamp committed
1
2
3
4
/******************************************************************************
*                                                                             *
*  Copyright 2019 Max Planck Institute for Dynamics and Self-Organization     *
*                                                                             *
5
*  This file is part of TurTLE.                                               *
Lukas Bentkamp's avatar
Lukas Bentkamp committed
6
*                                                                             *
7
*  TurTLE is free software: you can redistribute it and/or modify             *
Lukas Bentkamp's avatar
Lukas Bentkamp committed
8
9
10
11
*  it under the terms of the GNU General Public License as published          *
*  by the Free Software Foundation, either version 3 of the License,          *
*  or (at your option) any later version.                                     *
*                                                                             *
12
*  TurTLE is distributed in the hope that it will be useful,                  *
Lukas Bentkamp's avatar
Lukas Bentkamp committed
13
14
15
16
17
*  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          *
18
*  along with TurTLE.  If not, see <http://www.gnu.org/licenses/>             *
Lukas Bentkamp's avatar
Lukas Bentkamp committed
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
*                                                                             *
* Contact: Cristian.Lalescu@ds.mpg.de                                         *
*                                                                             *
******************************************************************************/



#include <string>
#include <cmath>
#include "kraichnan_field.hpp"
#include "scope_timer.hpp"
#include "fftw_tools.hpp"


template <typename rnumber>
int kraichnan_field<rnumber>::initialize(void)
{
    TIMEZONE("kraichnan_file::initialize");
    this->read_iteration();
    this->read_parameters();
    if (this->myrank == 0)
    {
        // set caching parameters
        hid_t fapl = H5Pcreate(H5P_FILE_ACCESS);
        herr_t cache_err = H5Pset_cache(fapl, 0, 521, 134217728, 1.0);
        variable_used_only_in_assert(cache_err);
        DEBUG_MSG("when setting stat_file cache I got %d\n", cache_err);
        this->stat_file = H5Fopen(
                (this->simname + ".h5").c_str(),
                H5F_ACC_RDWR,
                fapl);
    }
    int data_file_problem;
    if (this->myrank == 0)
        data_file_problem = this->grow_file_datasets();
    MPI_Bcast(&data_file_problem, 1, MPI_INT, 0, this->comm);
    if (data_file_problem > 0)
    {
        std::cerr <<
            data_file_problem <<
            " problems growing file datasets.\ntrying to exit now." <<
            std::endl;
        return EXIT_FAILURE;
    }

    this->velocity = new field<rnumber, FFTW, THREE>(
            nx, ny, nz,
            this->comm,
67
	        fftw_planner_string_to_flag[this->fftw_plan_rigor]);
68
    this->velocity->real_space_representation = true;
Lukas Bentkamp's avatar
Lukas Bentkamp committed
69

70
    // the velocity layout should be the same as the vorticity one
Lukas Bentkamp's avatar
Lukas Bentkamp committed
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
    this->kk = new kspace<FFTW, SMOOTH>(
            this->velocity->clayout, this->dkx, this->dky, this->dkz);

    // initialize particles
    this->ps = particles_system_builder(
                this->velocity,              // (field object)
                this->kk,                     // (kspace object, contains dkx, dky, dkz)
                tracers0_integration_steps, // to check coherency between parameters and hdf input file (nb rhs)
                (long long int)nparticles,  // to check coherency between parameters and hdf input file
	            this->get_current_fname(),    // particles input filename
                std::string("/tracers0/state/") + std::to_string(this->iteration), // dataset name for initial input
                std::string("/tracers0/rhs/")  + std::to_string(this->iteration),  // dataset name for initial input
                tracers0_neighbours,        // parameter (interpolation no neighbours)
                tracers0_smoothness,        // parameter
                this->comm,
                this->iteration+1);
    // initialize output objects
    this->particles_output_writer_mpi = new particles_output_hdf5<
        long long int, double, 3>(
                MPI_COMM_WORLD,
                "tracers0",
                nparticles,
                tracers0_integration_steps);
    this->particles_output_writer_mpi->setParticleFileLayout(this->ps->getParticleFileLayout());
    this->particles_sample_writer_mpi = new particles_output_sampling_hdf5<
        long long int, double, 3>(
                MPI_COMM_WORLD,
                this->ps->getGlobalNbParticles(),
                (this->simname + "_particles.h5"),
                "tracers0",
                "position/0");
    this->particles_sample_writer_mpi->setParticleFileLayout(this->ps->getParticleFileLayout());

104
105
    DEBUG_MSG("Coefficient: %g\n",
            this->spectrum_coefficient);
Lukas Bentkamp's avatar
Lukas Bentkamp committed
106

107
    this->generate_random_velocity();
108
109
110
111
112
113
114
115
116
117
118
    // is this the first iteration?
    if ((this->iteration == 0) and (this->output_velocity == 1))
    {
        // if yes, generate random field and save it
        this->velocity->io(
                this->get_current_fname(),
                "velocity",
                this->iteration,
                false);
    }

Lukas Bentkamp's avatar
Lukas Bentkamp committed
119
120
121
122
    return EXIT_SUCCESS;
}

template <typename rnumber>
123
int kraichnan_field<rnumber>::generate_random_velocity(void)
Lukas Bentkamp's avatar
Lukas Bentkamp committed
124
{
125
126
    TIMEZONE("kraichnan_file::make_random_field");
    // generate Gaussian random field
Lukas Bentkamp's avatar
Lukas Bentkamp committed
127
128
129
    make_gaussian_random_field(
        this->kk,
        this->velocity,
130
        this->iteration, // not an ideal choice because resulting field sequence depends on MPI/OpenMP configuration. See note below
Lukas Bentkamp's avatar
Lukas Bentkamp committed
131
132
        this->spectrum_slope,
        this->spectrum_k_cutoff,
133
        this->spectrum_coefficient * 3./2.); // incompressibility projection factor
134
135
    // this->velocity is now in Fourier space
    // project divfree, requires field in Fourier space
136
137
138
139
140
141
142
    // Note on the choice of random seed:
    // If in the future the simulation will continue with a smaller number of total threads (number of processes times number of threads per process),
    // then during that run some of the threads will be seeded with a seed that has already been used for a previous iteration.
    // So some sequences of Fourier modes will be identical to sequences of Fourier modes that occured in the past.
    // Also see implementation of "make_gaussian_random_field".
    // One work-around would be to multiply "this->iteration" with 10 or so ---
    // it is unlikely the simulation will be continued with less than 0.1 of the initial total number of threads.
143
    DEBUG_MSG("L2Norm before: %g\n",
Lukas Bentkamp's avatar
Lukas Bentkamp committed
144
            this->velocity->L2norm(this->kk));
Lukas Bentkamp's avatar
Lukas Bentkamp committed
145
    this->kk->template project_divfree<rnumber>(this->velocity->get_cdata());
146
    DEBUG_MSG("L2Norm after: %g\n",
Lukas Bentkamp's avatar
Lukas Bentkamp committed
147
            this->velocity->L2norm(this->kk));
148
149
150
151
    // transform velocity to real space representation
    this->velocity->ift();
    return EXIT_SUCCESS;
}
152

153
154
155
156
template <typename rnumber>
int kraichnan_field<rnumber>::step(void)
{
    TIMEZONE("kraichnan_file::step");
Lukas Bentkamp's avatar
Lukas Bentkamp committed
157
    this->ps->completeLoop(sqrt(this->dt));
158
    this->generate_random_velocity();
Lukas Bentkamp's avatar
Lukas Bentkamp committed
159
160
161
162
163
164
165
166
    this->iteration++;
    return EXIT_SUCCESS;
}

template <typename rnumber>
int kraichnan_field<rnumber>::write_checkpoint(void)
{
    TIMEZONE("kraichnan_file::write_checkpoint");
167
168
169
170
171
172
173
174
175
176
177
178
179
    this->update_checkpoint();
    // output field
    if (this->output_velocity == 1)
    {
        assert(this->velocity->real_space_representation);
        std::string fname = this->get_current_fname();
        this->velocity->io(
                fname,
                "velocity",
                this->iteration,
                false);
    }
    // output particles
Lukas Bentkamp's avatar
Lukas Bentkamp committed
180
181
182
183
184
185
186
187
    this->particles_output_writer_mpi->open_file(this->get_current_fname());
    this->particles_output_writer_mpi->template save<3>(
            this->ps->getParticlesState(),
            this->ps->getParticlesRhs(),
            this->ps->getParticlesIndexes(),
            this->ps->getLocalNbParticles(),
            this->iteration);
    this->particles_output_writer_mpi->close_file();
188
    this->write_iteration();
Lukas Bentkamp's avatar
Lukas Bentkamp committed
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
    return EXIT_SUCCESS;
}

template <typename rnumber>
int kraichnan_field<rnumber>::finalize(void)
{
    TIMEZONE("kraichnan_field::finalize");
    if (this->myrank == 0)
        H5Fclose(this->stat_file);
    // good practice rule number n+1: always delete in inverse order of allocation
    delete this->ps.release();
    delete this->particles_output_writer_mpi;
    delete this->particles_sample_writer_mpi;
    delete this->kk;
    delete this->velocity;
    return EXIT_SUCCESS;
}

/** \brief Compute statistics.
 *
 */

template <typename rnumber>
int kraichnan_field<rnumber>::do_stats()
{
    TIMEZONE("kraichnan_field::do_stats");
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
    /// compute field statistics
    if (this->iteration % this->niter_stat == 0)
    {
        hid_t stat_group;
        if (this->myrank == 0)
            stat_group = H5Gopen(
                    this->stat_file,
                    "statistics",
                    H5P_DEFAULT);
        else
            stat_group = 0;
        field<rnumber, FFTW, THREE> *tvelocity = new field<rnumber, FFTW, THREE>(
                this->nx, this->ny, this->nz,
                this->comm,
	            fftw_planner_string_to_flag[this->fftw_plan_rigor]);
        *tvelocity = *this->velocity;
        tvelocity->compute_stats(
            this->kk,
            stat_group,
            "velocity",
            this->iteration / this->niter_stat,
            this->max_velocity_estimate/sqrt(3));
        if (this->myrank == 0)
            H5Gclose(stat_group);
        delete tvelocity;
    }
Lukas Bentkamp's avatar
Lukas Bentkamp committed
241
242
243
244
245
246
247
248
249
250
251
252
    /// either one of two conditions suffices to compute statistics:
    /// 1) current iteration is a multiple of niter_part
    /// 2) we are within niter_part_fine_duration/2 of a multiple of niter_part_fine_period
    if (!(this->iteration % this->niter_part == 0 ||
          ((this->iteration + this->niter_part_fine_duration/2) % this->niter_part_fine_period <=
           this->niter_part_fine_duration)))
        return EXIT_SUCCESS;

    // allocate temporary data array
    std::unique_ptr<double[]> pdata(new double[3*this->ps->getLocalNbParticles()]);

    /// copy position data
253

Lukas Bentkamp's avatar
Lukas Bentkamp committed
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
    /// sample position
    std::copy(this->ps->getParticlesState(),
              this->ps->getParticlesState()+3*this->ps->getLocalNbParticles(),
              pdata.get());
    this->particles_sample_writer_mpi->template save_dataset<3>(
            "tracers0",
            "position",
            this->ps->getParticlesState(),
            &pdata,
            this->ps->getParticlesIndexes(),
            this->ps->getLocalNbParticles(),
            this->ps->get_step_idx()-1);

    /// sample velocity
    std::fill_n(pdata.get(), 3*this->ps->getLocalNbParticles(), 0);
    this->ps->sample_compute_field(*this->velocity, pdata.get());
    this->particles_sample_writer_mpi->template save_dataset<3>(
            "tracers0",
            "velocity",
            this->ps->getParticlesState(),
            &pdata,
            this->ps->getParticlesIndexes(),
            this->ps->getLocalNbParticles(),
            this->ps->get_step_idx()-1);

    // deallocate temporary data array
    delete[] pdata.release();

    return EXIT_SUCCESS;
}

template <typename rnumber>
int kraichnan_field<rnumber>::read_parameters(void)
{
    TIMEZONE("kraichnan_field::read_parameters");
    this->direct_numerical_simulation::read_parameters();
    hid_t parameter_file = H5Fopen((this->simname + ".h5").c_str(), H5F_ACC_RDONLY, H5P_DEFAULT);
    this->fftw_plan_rigor = hdf5_tools::read_string(parameter_file, "parameters/fftw_plan_rigor");
    this->dt = hdf5_tools::read_value<double>(parameter_file, "parameters/dt");
    this->niter_part = hdf5_tools::read_value<int>(parameter_file, "parameters/niter_part");
    this->niter_part_fine_period = hdf5_tools::read_value<int>(parameter_file, "parameters/niter_part_fine_period");
    this->niter_part_fine_duration = hdf5_tools::read_value<int>(parameter_file, "parameters/niter_part_fine_duration");
    this->nparticles = hdf5_tools::read_value<int>(parameter_file, "parameters/nparticles");
    this->tracers0_integration_steps = hdf5_tools::read_value<int>(parameter_file, "parameters/tracers0_integration_steps");
    this->tracers0_neighbours = hdf5_tools::read_value<int>(parameter_file, "parameters/tracers0_neighbours");
    this->tracers0_smoothness = hdf5_tools::read_value<int>(parameter_file, "parameters/tracers0_smoothness");
300
301
302
303
    this->spectrum_slope = hdf5_tools::read_value<double>(parameter_file, "parameters/spectrum_slope");
    this->spectrum_k_cutoff = hdf5_tools::read_value<double>(parameter_file, "parameters/spectrum_k_cutoff");
    this->spectrum_coefficient = hdf5_tools::read_value<double>(parameter_file, "parameters/spectrum_coefficient");
    this->field_random_seed = hdf5_tools::read_value<int>(parameter_file, "parameters/field_random_seed");
304
305
    this->output_velocity = hdf5_tools::read_value<int>(parameter_file, "parameters/output_velocity");
    this->max_velocity_estimate = hdf5_tools::read_value<double>(parameter_file, "parameters/max_velocity_estimate");
Lukas Bentkamp's avatar
Lukas Bentkamp committed
306
307
308
309
    H5Fclose(parameter_file);
    return EXIT_SUCCESS;
}

310
311
312
template <class rnumber>
void kraichnan_field<rnumber>::update_checkpoint()
{
313
314
    TIMEZONE("kraichnan_field::update_checkpoint");
    DEBUG_MSG("entered kraichan_field::update_checkpoint\n");
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
    std::string fname = this->get_current_fname();
    if (this->kk->layout->myrank == 0)
    {
        bool file_exists = false;
        {
            struct stat file_buffer;
            file_exists = (stat(fname.c_str(), &file_buffer) == 0);
        }
        if (file_exists)
        {
            // check how many fields there are in the checkpoint file
            // increment checkpoint if needed
            hsize_t fields_stored;
            hid_t fid, group_id;
            fid = H5Fopen(fname.c_str(), H5F_ACC_RDONLY, H5P_DEFAULT);
330
331
332
333
334
            group_id = H5Gopen(fid, "tracers0/state", H5P_DEFAULT);
            bool dset_exists;
            if (group_id > 0)
            {
                H5Gget_num_objs(
335
336
                    group_id,
                    &fields_stored);
337
                dset_exists = H5Lexists(
338
339
340
                    group_id,
                    std::to_string(this->iteration).c_str(),
                    H5P_DEFAULT);
341
342
343
344
345
346
            }
            else
            {
                fields_stored = 0;
                dset_exists = false;
            }
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
            H5Gclose(group_id);
            H5Fclose(fid);
            if ((int(fields_stored) >= this->checkpoints_per_file) &&
                !dset_exists)
                this->checkpoint++;
        }
        else
        {
            // create file, create fields_stored dset
            hid_t fid = H5Fcreate(
                    fname.c_str(),
                    H5F_ACC_EXCL,
                    H5P_DEFAULT,
                    H5P_DEFAULT);
            hid_t gg = H5Gcreate(
                    fid,
363
                    "tracers0",
364
365
366
367
368
                    H5P_DEFAULT,
                    H5P_DEFAULT,
                    H5P_DEFAULT);
            hid_t ggg = H5Gcreate(
                    gg,
369
                    "state",
370
371
372
373
374
375
376
377
378
                    H5P_DEFAULT,
                    H5P_DEFAULT,
                    H5P_DEFAULT);
            H5Gclose(ggg);
            H5Gclose(gg);
            H5Fclose(fid);
        }
    }
    MPI_Bcast(&this->checkpoint, 1, MPI_INT, 0, this->kk->layout->comm);
379
    DEBUG_MSG("exiting kraichan_field::update_checkpoint\n");
380
381
}

Lukas Bentkamp's avatar
Lukas Bentkamp committed
382
383
384
template class kraichnan_field<float>;
template class kraichnan_field<double>;