-
Notifications
You must be signed in to change notification settings - Fork 638
/
fields_dump.cpp
287 lines (257 loc) · 10.1 KB
/
fields_dump.cpp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
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
67
68
69
70
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
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
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
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
241
242
243
244
245
246
247
248
249
250
251
252
253
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
/* Copyright (C) 2005-2022 Massachusetts Institute of Technology
%
% 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
% the Free Software Foundation; either version 2, or (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, write to the Free Software Foundation,
% Inc., 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA.
*/
// Dump/load raw fields data to/from an HDF5 file. Only
// works if the number of processors/chunks is the same.
#include <math.h>
#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include <cassert>
#include "meep.hpp"
#include "meep_internals.hpp"
namespace meep {
void fields::dump_fields_chunk_field(h5file *h5f, bool single_parallel_file,
const std::string &field_name,
FieldPtrGetter field_ptr_getter) {
/*
* make/save a num_chunks x NUM_FIELD_COMPONENTS x 2 array counting
* the number of entries in the 'field_name' array for each chunk.
*
* When 'single_parallel_file' is true, we are creating a single block of data
* for ALL chunks (that are merged using MPI). Otherwise, we are just
* making a copy of just the chunks that are ours.
*/
int my_num_chunks = 0;
for (int i = 0; i < num_chunks; i++) {
my_num_chunks += (single_parallel_file || chunks[i]->is_mine());
}
size_t num_f_size = my_num_chunks * NUM_FIELD_COMPONENTS * 2;
std::vector<size_t> num_f_(num_f_size);
size_t my_ntot = 0;
for (int i = 0, chunk_i = 0; i < num_chunks; i++) {
if (chunks[i]->is_mine()) {
FOR_FIELD_TYPES(ft) {
if (chunks[i]->pol[ft] != NULL)
meep::abort("non-null polarization_state in fields::dump (unsupported)");
}
size_t ntot = chunks[i]->gv.ntot();
for (int c = 0; c < NUM_FIELD_COMPONENTS; ++c) {
for (int d = 0; d < 2; ++d) {
realnum **f = field_ptr_getter(chunks[i], c, d);
if (*f) {
my_ntot += (num_f_[(chunk_i * NUM_FIELD_COMPONENTS + c) * 2 + d] = ntot);
}
}
}
}
chunk_i += (chunks[i]->is_mine() || single_parallel_file);
}
std::vector<size_t> num_f;
if (single_parallel_file) {
num_f.resize(num_f_size);
sum_to_master(num_f_.data(), num_f.data(), num_f_size);
} else {
num_f = std::move(num_f_);
}
/* determine total dataset size and offset of this process's data */
size_t my_start = 0;
size_t ntotal = my_ntot;
if (single_parallel_file) {
my_start = partial_sum_to_all(my_ntot) - my_ntot;
ntotal = sum_to_all(my_ntot);
}
size_t dims[3] = {(size_t)my_num_chunks, NUM_FIELD_COMPONENTS, 2};
size_t start[3] = {0, 0, 0};
std::string num_f_name = std::string("num_") + field_name;
h5f->create_data(num_f_name.c_str(), 3, dims);
if (am_master() || !single_parallel_file) {
h5f->write_chunk(3, start, dims, num_f.data());
}
/* write the data */
h5f->create_data(field_name.c_str(), 1, &ntotal, false /* append_data */, false /* single_precision */);
for (int i = 0; i < num_chunks; i++) {
if (chunks[i]->is_mine()) {
size_t ntot = chunks[i]->gv.ntot();
for (int c = 0; c < NUM_FIELD_COMPONENTS; ++c) {
for (int d = 0; d < 2; ++d) {
realnum **f = field_ptr_getter(chunks[i], c, d);
if (*f) {
h5f->write_chunk(1, &my_start, &ntot, *f);
my_start += ntot;
}
}
}
}
}
}
void fields::dump(const char *filename, bool single_parallel_file) {
if (verbosity > 0) {
printf("creating fields output file \"%s\" (%d)...\n", filename, single_parallel_file);
}
h5file file(filename, h5file::WRITE, single_parallel_file, !single_parallel_file);
// Write out the current time 't'
size_t dims[1] = {1};
size_t start[1] = {0};
size_t _t[1] = {(size_t)t};
file.create_data("t", 1, dims);
if (am_master() || !single_parallel_file) file.write_chunk(1, start, dims, _t);
dump_fields_chunk_field(
&file, single_parallel_file, "f",
[](fields_chunk *chunk, int c, int d) { return &(chunk->f[c][d]); });
dump_fields_chunk_field(
&file, single_parallel_file, "f_u",
[](fields_chunk *chunk, int c, int d) { return &(chunk->f_u[c][d]); });
dump_fields_chunk_field(
&file, single_parallel_file, "f_w",
[](fields_chunk *chunk, int c, int d) { return &(chunk->f_w[c][d]); });
dump_fields_chunk_field(
&file, single_parallel_file, "f_cond",
[](fields_chunk *chunk, int c, int d) { return &(chunk->f_cond[c][d]); });
dump_fields_chunk_field(
&file, single_parallel_file, "f_w_prev",
[](fields_chunk *chunk, int c, int d) { return &(chunk->f_w_prev[c][d]); });
// Dump DFT chunks.
for (int i = 0; i < num_chunks; i++) {
if (single_parallel_file || chunks[i]->is_mine()) {
char dataname[1024];
snprintf(dataname, 1024, "chunk%02d", i);
save_dft_hdf5(chunks[i]->dft_chunks, dataname, &file, 0, single_parallel_file);
}
}
}
void fields::load_fields_chunk_field(h5file *h5f, bool single_parallel_file,
const std::string &field_name,
FieldPtrGetter field_ptr_getter) {
int my_num_chunks = 0;
for (int i = 0; i < num_chunks; i++) {
my_num_chunks += (single_parallel_file || chunks[i]->is_mine());
}
size_t num_f_size = my_num_chunks * NUM_FIELD_COMPONENTS * 2;
std::vector<size_t> num_f(num_f_size);
int rank;
size_t dims[3], _dims[3] = {(size_t)my_num_chunks, NUM_FIELD_COMPONENTS, 2};
size_t start[3] = {0, 0, 0};
std::string num_f_name = std::string("num_") + field_name;
h5f->read_size(num_f_name.c_str(), &rank, dims, 3);
if (rank != 3 || _dims[0] != dims[0] || _dims[1] != dims[1] || _dims[2] != dims[2])
meep::abort("chunk mismatch in fields::load");
if (am_master() || !single_parallel_file) h5f->read_chunk(3, start, dims, num_f.data());
if (single_parallel_file) {
h5f->prevent_deadlock();
broadcast(0, num_f.data(), dims[0] * dims[1] * dims[2]);
}
/* allocate data as needed and check sizes */
size_t my_ntot = 0;
for (int i = 0, chunk_i = 0; i < num_chunks; i++) {
if (chunks[i]->is_mine()) {
size_t ntot = chunks[i]->gv.ntot();
for (int c = 0; c < NUM_FIELD_COMPONENTS; ++c) {
for (int d = 0; d < 2; ++d) {
size_t n = num_f[(chunk_i * NUM_FIELD_COMPONENTS + c) * 2 + d];
realnum **f = field_ptr_getter(chunks[i], c, d);
if (n == 0) {
delete[] *f;
*f = NULL;
} else {
if (n != ntot)
meep::abort("grid size mismatch %zd vs %zd in fields::load", n, ntot);
// here we need to allocate the fields array for H in the PML region
// because of H = B in fields_chunk::alloc_f whereby H is lazily
// allocated in fields_chunk::update_eh during the first timestep
const direction d_c = component_direction(c);
if (!(*f) || (*f && is_magnetic(component(c)) && chunks[i]->s->sigsize[d_c] > 1))
*f = new realnum[ntot];
my_ntot += ntot;
}
}
}
}
chunk_i += (chunks[i]->is_mine() || single_parallel_file);
}
/* determine total dataset size and offset of this process's data */
size_t my_start = 0;
size_t ntotal = my_ntot;
if (single_parallel_file) {
my_start = partial_sum_to_all(my_ntot) - my_ntot;
ntotal = sum_to_all(my_ntot);
}
/* read the data */
h5f->read_size(field_name.c_str(), &rank, dims, 1);
if (rank != 1 || dims[0] != ntotal) {
meep::abort(
"inconsistent data size for '%s' in fields::load (rank, dims[0]): "
"(%d, %zu) != (1, %zu)",
field_name.c_str(), rank, dims[0], ntotal);
}
for (int i = 0; i < num_chunks; i++) {
if (chunks[i]->is_mine()) {
size_t ntot = chunks[i]->gv.ntot();
for (int c = 0; c < NUM_FIELD_COMPONENTS; ++c) {
for (int d = 0; d < 2; ++d) {
realnum **f = field_ptr_getter(chunks[i], c, d);
if (*f) {
h5f->read_chunk(1, &my_start, &ntot, *f);
my_start += ntot;
}
}
}
}
}
}
void fields::load(const char *filename, bool single_parallel_file) {
if (verbosity > 0)
printf("reading fields from file \"%s\" (%d)...\n", filename, single_parallel_file);
h5file file(filename, h5file::READONLY, single_parallel_file, !single_parallel_file);
// Read in the current time 't'
int rank;
size_t dims[1] = {1};
size_t start[1] = {0};
size_t _t[1];
file.read_size("t", &rank, dims, 1);
if (rank != 1 || dims[0] != 1) meep::abort("time size mismatch in fields::load");
if (am_master() || !single_parallel_file) file.read_chunk(1, start, dims, _t);
if (single_parallel_file) {
file.prevent_deadlock();
broadcast(0, _t, dims[0]);
}
t = static_cast<int>(_t[0]);
calc_sources(time());
load_fields_chunk_field(
&file, single_parallel_file, "f",
[](fields_chunk *chunk, int c, int d) { return &(chunk->f[c][d]); });
load_fields_chunk_field(
&file, single_parallel_file, "f_u",
[](fields_chunk *chunk, int c, int d) { return &(chunk->f_u[c][d]); });
load_fields_chunk_field(
&file, single_parallel_file, "f_w",
[](fields_chunk *chunk, int c, int d) { return &(chunk->f_w[c][d]); });
load_fields_chunk_field(
&file, single_parallel_file, "f_cond",
[](fields_chunk *chunk, int c, int d) { return &(chunk->f_cond[c][d]); });
load_fields_chunk_field(
&file, single_parallel_file, "f_w_prev",
[](fields_chunk *chunk, int c, int d) { return &(chunk->f_w_prev[c][d]); });
// Load DFT chunks.
for (int i = 0; i < num_chunks; i++) {
if (single_parallel_file || chunks[i]->is_mine()) {
char dataname[1024];
snprintf(dataname, 1024, "chunk%02d", i);
load_dft_hdf5(chunks[i]->dft_chunks, dataname, &file, 0, single_parallel_file);
}
}
}
} // namespace meep