-
Notifications
You must be signed in to change notification settings - Fork 332
/
fextract.cpp
403 lines (369 loc) · 16.2 KB
/
fextract.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
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
#include <AMReX.H>
#include <AMReX_Print.H>
#include <AMReX_PlotFileUtil.H>
#include <AMReX_MultiFabUtil.H>
#include <AMReX_ParallelDescriptor.H>
#include <limits>
#include <iterator>
#include <fstream>
using namespace amrex;
void main_main()
{
const int narg = amrex::command_argument_count();
std::string slicefile;
std::string pltfile;
int idir = 0;
std::string varnames_arg;
bool center = true;
int coarse_level = 0;
int fine_level = -1; // This will be fixed later
Real xcoord = std::numeric_limits<Real>::lowest();
Real ycoord = std::numeric_limits<Real>::lowest();
Real zcoord = std::numeric_limits<Real>::lowest();
bool scientific = false;
bool csv = false;
int precision = 17;
Real tolerance = std::numeric_limits<Real>::lowest();
bool print_info = false;
int farg = 1;
while (farg <= narg) {
const std::string& name = amrex::get_command_argument(farg);
if (name == "-s" || name == "--slicefile") {
slicefile = amrex::get_command_argument(++farg);
} else if (name == "-d" || name == "--direction") {
idir = std::stoi(amrex::get_command_argument(++farg));
} else if (name == "-v" || name == "--variable") {
varnames_arg = amrex::get_command_argument(++farg);
} else if (name == "-x") {
xcoord = Real(std::stod(amrex::get_command_argument(++farg)));
} else if (name == "-y") {
ycoord = Real(std::stod(amrex::get_command_argument(++farg)));
} else if (name == "-z") {
zcoord = Real(std::stod(amrex::get_command_argument(++farg)));
} else if (name == "-l" || name == "--lower_left") {
center = false;
} else if (name == "-c" || name == "--coarse_level") {
coarse_level = std::stoi(amrex::get_command_argument(++farg));
} else if (name == "-f" || name == "--fine_level") {
fine_level = std::stoi(amrex::get_command_argument(++farg));
} else if (name == "-e" || name == "--scientific") {
scientific = true;
} else if (name == "-csv" || name == "--csv") {
csv = true;
} else if (name == "-p" || name == "--precision") {
precision = std::stoi(amrex::get_command_argument(++farg));
} else if (name == "-t" || name == "--tolerance") {
tolerance = Real(std::stod(amrex::get_command_argument(++farg)));
} else if (name == "-i" || name == "--info") {
print_info = true;
} else {
break;
}
++farg;
}
if (pltfile.empty() && farg <= narg) {
pltfile = amrex::get_command_argument(farg);
}
if (pltfile.empty()) {
amrex::Print()
<< "\n"
<< " Extract at 1D slice through a plotfile in any coordinate direction.\n"
<< " Works with 1-, 2-, or 3-d datasets.\n"
<< "\n"
<< " Usage:\n"
<< " fextract [-s outfile] [-d dir] [-v variable] [-y ycood] [-c coarse_level] [-f fine_level] plotfile\n"
<< "\n"
<< " args [-s|--slicefile] slice file : slice file (optional)\n"
<< " [-d|--direction] idir : slice direction {0 (default), 1, or 2}\n"
<< " [-v|--variable] varname(s) : only output the values of variable\n"
<< " varname (space separated string for\n"
<< " multiple variables)\n"
<< " [-l|--lower_left] : slice through lower left corner\n"
<< " instead of center\n"
<< " [-x][-y][-z] : (x,y,z)-coordinate to pass through\n"
<< " (overrides center/lower-left)\n"
<< " [-c|--coarse_level] coarse level : coarsest level to extract from\n"
<< " [-f|--fine_level] fine level : finest level to extract from\n"
<< " [-e|--scientific] : output data in scientific notation\n"
<< " [-p|--precision] precision : decimal precision {17 (default)}\n"
<< " [-t|--tolerance] tolerance : set to 0 any value lower than tolerance\n"
<< " [-i|--info] : output job info\n"
<< "\n"
<< " If a job_info file is present in the plotfile, that information is made\n"
<< " available at the end of the slice file (commented out), for reference.\n"
<< '\n';
return;
}
// slicefile not defined, default to plotfile.slice
if (slicefile.empty()) {
slicefile = pltfile;
if (slicefile.back() == '/') {
slicefile.pop_back();
}
slicefile += ".slice";
}
PlotFileData pf(pltfile);
const Vector<std::string>& var_names_pf = pf.varNames();
Vector<std::string> var_names;
if (varnames_arg.empty()) {
var_names = var_names_pf;
} else {
std::istringstream is(varnames_arg);
var_names.assign(std::istream_iterator<std::string>{is},
std::istream_iterator<std::string>{ });
var_names.erase(std::remove_if(var_names.begin(), var_names.end(),
[&](std::string const& x) {
return var_names_pf.end() ==
std::find(var_names_pf.begin(), var_names_pf.end(), x);
}),
var_names.end());
if (var_names.empty()) {
amrex::Abort("ERROR: no valid variable names");
}
}
Array<Real,AMREX_SPACEDIM> problo = pf.probLo();
Array<Real,AMREX_SPACEDIM> dx0 = pf.cellSize(0);
Box probdom0 = pf.probDomain(0);
const auto lo0 = amrex::lbound(probdom0);
const auto hi0 = amrex::ubound(probdom0);
// compute the index of the center or lower left of the domain on the
// coarse grid. These are used to set the position of the slice in
// the transverse direction.
int iloc = 0, jloc = 0, kloc = 0;
if (center) {
iloc = (hi0.x-lo0.x+1)/2 + lo0.x;
jloc = (hi0.y-lo0.y+1)/2 + lo0.y;
kloc = (hi0.z-lo0.z+1)/2 + lo0.z;
}
if (xcoord > -1.e36 && AMREX_SPACEDIM >= 1) {
// we specified the x value to pass through
iloc = hi0.x;
for (int i = lo0.x; i <= hi0.x; ++i) {
amrex::Real xc = problo[0] + (Real(i)+Real(0.5))*dx0[0];
if (xc > xcoord) {
iloc = i;
break;
}
}
}
if (ycoord > -1.e36 && AMREX_SPACEDIM >= 2) {
// we specified the y value to pass through
jloc = hi0.y;
for (int j = lo0.y; j <= hi0.y; ++j) {
amrex::Real yc = problo[1] + (Real(j)+Real(0.5))*dx0[1];
if (yc > ycoord) {
jloc = j;
break;
}
}
}
if (zcoord > -1.e36 && AMREX_SPACEDIM == 3) {
// we specified the z value to pass through
kloc = hi0.z;
for (int k = lo0.z; k <= hi0.z; ++k) {
amrex::Real zc = problo[2] + (Real(k)+Real(0.5))*dx0[2];
if (zc > zcoord) {
kloc = k;
break;
}
}
}
const int dim = pf.spaceDim();
if (idir < 0 || idir >= dim) {
amrex::Print() << " invalid direction\n";
return;
} else if (idir == 0) {
amrex::Print() << " slicing along x-direction at coarse grid (j,k)=(" << jloc << "," << kloc << ") and output to " << slicefile << "\n";
} else if (idir == 1) {
amrex::Print() << " slicing along y-direction at coarse grid (i,k)=(" << iloc << "," << kloc << ") and output to " << slicefile << "\n";
} else if (idir == 2) {
amrex::Print() << " slicing along z-direction at coarse grid (i,j)=(" << iloc << "," << jloc << ") and output to " << slicefile << "\n";
}
const IntVect ivloc{AMREX_D_DECL(iloc,jloc,kloc)};
if (fine_level < 0) { fine_level = pf.finestLevel(); }
// sanity check on valid selected levels
if (fine_level > pf.finestLevel() || coarse_level < 0 || coarse_level > fine_level) {
amrex::Abort("Invalid level selection");
}
Vector<Real> pos;
Vector<Vector<Real> > data(var_names.size());
IntVect rr{1};
for (int ilev = coarse_level; ilev <= fine_level; ++ilev) {
Box slice_box(ivloc*rr,ivloc*rr);
slice_box.setSmall(idir, std::numeric_limits<int>::lowest());
slice_box.setBig(idir, std::numeric_limits<int>::max());
Array<Real,AMREX_SPACEDIM> dx = pf.cellSize(ilev);
if (ilev < fine_level) {
IntVect ratio{pf.refRatio(ilev)};
for (int idim = dim; idim < AMREX_SPACEDIM; ++idim) {
ratio[idim] = 1;
}
const iMultiFab mask = makeFineMask(pf.boxArray(ilev), pf.DistributionMap(ilev),
pf.boxArray(ilev+1), ratio);
for (int ivar = 0; ivar < var_names.size(); ++ivar) {
const MultiFab& mf = pf.get(ilev, var_names[ivar]);
for (MFIter mfi(mf); mfi.isValid(); ++mfi) {
const Box& bx = mfi.validbox() & slice_box;
if (bx.ok()) {
const auto& m = mask.array(mfi);
const auto& fab = mf.array(mfi);
const auto lo = amrex::lbound(bx);
const auto hi = amrex::ubound(bx);
for (int k = lo.z; k <= hi.z; ++k) {
for (int j = lo.y; j <= hi.y; ++j) {
for (int i = lo.x; i <= hi.x; ++i) {
if (m(i,j,k) == 0) { // not covered by fine
if (pos.size() == data[ivar].size()) {
Array<Real,AMREX_SPACEDIM> p
= {AMREX_D_DECL(problo[0]+static_cast<Real>(i+0.5)*dx[0],
problo[1]+static_cast<Real>(j+0.5)*dx[1],
problo[2]+static_cast<Real>(k+0.5)*dx[2])};
pos.push_back(p[idir]);
}
data[ivar].push_back(fab(i,j,k));
}
}
}
}
}
}
}
rr *= ratio;
} else {
for (int ivar = 0; ivar < var_names.size(); ++ivar) {
const MultiFab& mf = pf.get(ilev, var_names[ivar]);
for (MFIter mfi(mf); mfi.isValid(); ++mfi) {
const Box& bx = mfi.validbox() & slice_box;
if (bx.ok()) {
const auto& fab = mf.array(mfi);
const auto lo = amrex::lbound(bx);
const auto hi = amrex::ubound(bx);
for (int k = lo.z; k <= hi.z; ++k) {
for (int j = lo.y; j <= hi.y; ++j) {
for (int i = lo.x; i <= hi.x; ++i) {
if (pos.size() == data[ivar].size()) {
Array<Real,AMREX_SPACEDIM> p
= {AMREX_D_DECL(problo[0]+static_cast<Real>(i+0.5)*dx[0],
problo[1]+static_cast<Real>(j+0.5)*dx[1],
problo[2]+static_cast<Real>(k+0.5)*dx[2])};
pos.push_back(p[idir]);
}
data[ivar].push_back(fab(i,j,k));
}
}
}
}
}
}
}
}
#ifdef BL_USE_MPI
{
const auto numpts = int(pos.size());
auto numpts_vec = ParallelDescriptor::Gather(numpts,
ParallelDescriptor::IOProcessorNumber());
Vector<int> recvcnt, disp;
Vector<Real> allpos;
Vector<Vector<Real> > alldata(data.size());
if (ParallelDescriptor::IOProcessor()) {
recvcnt.resize(numpts_vec.size());
disp.resize(numpts_vec.size());
int ntot = 0;
disp[0] = 0;
for (int i = 0, N = int(numpts_vec.size()); i < N; ++i) {
ntot += numpts_vec[i];
recvcnt[i] = numpts_vec[i];
if (i+1 < N) {
disp[i+1] = disp[i] + numpts_vec[i];
}
}
allpos.resize(ntot);
alldata.resize(data.size());
for (auto& v : alldata) {
v.resize(ntot);
}
} else {
recvcnt.resize(1);
disp.resize(1);
allpos.resize(1);
for (auto& v: alldata) {
v.resize(1);
}
}
ParallelDescriptor::Gatherv(pos.data(), numpts, allpos.data(), recvcnt, disp,
ParallelDescriptor::IOProcessorNumber());
for (int i = 0; i < data.size(); ++i) {
ParallelDescriptor::Gatherv(data[i].data(), numpts, alldata[i].data(), recvcnt, disp,
ParallelDescriptor::IOProcessorNumber());
}
if (ParallelDescriptor::IOProcessor()) {
pos = std::move(allpos);
data = std::move(alldata);
}
}
#endif
if (ParallelDescriptor::IOProcessor()) {
Vector<std::pair<Real,int> > posidx;
posidx.reserve(pos.size());
for (int i = 0; i < pos.size(); ++i) {
posidx.emplace_back(pos[i],i);
}
std::sort(posidx.begin(), posidx.end());
const std::string dirstr = (idir == 0) ? "x" : ((idir == 1) ? "y" : "z");
std::ofstream ofs(slicefile, std::ios::trunc);
if (scientific) {
ofs << std::scientific;
}
if (csv)
{
ofs << std::setw(24) << std::left << dirstr;
for (auto const& vname : var_names) {
ofs << ", " << std::setw(24) << std::left << vname;
}
ofs << "\n";
for (auto const& pos_idx : posidx) {
ofs << std::setw(25) << std::left << std::setprecision(17) << pos_idx.first;
for (auto const& vr : data) {
ofs << ", " << std::setw(25) << std::left << std::setprecision(17) << vr[pos_idx.second];
}
ofs << "\n";
}
}
else
{
ofs << "# 1-d slice in " << dirstr << "-direction, file: " << pltfile << "\n";
ofs << "# time = " << std::setw(20) << std::setprecision(precision) << pf.time() << "\n";
ofs << "#" << std::setw(24) << dirstr;
for (auto const& vname : var_names) {
ofs << " " << std::setw(24) << std::right << vname;
}
ofs << "\n";
for (auto const& pos_idx : posidx) {
ofs << std::setw(25) << std::right << std::setprecision(precision) << pos_idx.first;
for (auto& vr : data) {
if (std::abs(vr[pos_idx.second])< tolerance ) { vr[pos_idx.second] = 0.; }
ofs << std::setw(25) << std::right << std::setprecision(precision) << vr[pos_idx.second];
}
ofs << "\n";
}
if (print_info) {
// job_info? if so write it out to the slice file end
std::ifstream jobinfo(pltfile+"/job_info");
if (jobinfo.good()) {
ofs << "\n";
std::string s;
while (std::getline(jobinfo, s)) {
ofs << "#" << s << "\n";
}
}
}
}
}
}
int main (int argc, char* argv[])
{
amrex::SetVerbose(0);
amrex::Initialize(argc, argv, false);
main_main();
amrex::Finalize();
}