Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
AMReX_WriteBinaryParticleData.H
Go to the documentation of this file.
1#ifndef AMREX_WRITE_BINARY_PARTICLE_DATA_H
2#define AMREX_WRITE_BINARY_PARTICLE_DATA_H
3#include <AMReX_Config.H>
4
5#include <AMReX_TypeTraits.H>
7#include <AMReX_GpuDevice.H>
8#include <AMReX_VisMF.H>
9
10namespace amrex {
11
13{
14 template <typename SrcData>
16 int operator() (const SrcData& src, int i) const noexcept
17 {
18 return (src.id(i).is_valid());
19 }
20};
21
23namespace particle_detail {
24
25template <typename ParticleReal>
26std::size_t PSizeInFile (const Vector<int>& wrc, const Vector<int>& wic)
27{
28 std::size_t rsize = sizeof(ParticleReal)*std::accumulate(wrc.begin(), wrc.end(), 0);
29 std::size_t isize = sizeof(int)*std::accumulate(wic.begin(), wic.end(), 0);
30 return rsize + isize + AMREX_SPACEDIM*sizeof(ParticleReal) + 2*sizeof(int);
31}
32
33template <class Container,
34 class PTile,
35 class F>
36void
37fillFlagsGpu (Container& pflags, const PTile& ptile, F const& f)
38{
39 const auto& ptd = ptile.getConstParticleTileData();
40 const auto np = ptile.numParticles();
41 pflags.resize(np, 0);
42 auto flag_ptr = pflags.data();
44 [=] AMREX_GPU_DEVICE (int k, amrex::RandomEngine const& engine) noexcept
45 {
46 const auto p = ptd.getSuperParticle(k);
47 amrex::ignore_unused(flag_ptr, f, engine);
48 if constexpr (IsCallable<F,decltype(p),RandomEngine>::value) {
49 flag_ptr[k] = f(p,engine);
50 } else if constexpr (IsCallable<F,decltype(p)>::value) {
51 flag_ptr[k] = f(p);
52 } else if constexpr (IsCallable<F,decltype(ptd),int,RandomEngine>::value) {
53 flag_ptr[k] = f(ptd,k,engine);
54 } else {
55 flag_ptr[k] = f(ptd,k);
56 }
57 });
58}
59
60template <class Container,
61 class PTile,
62 class F>
63void
64fillFlagsCpu (Container& pflags, const PTile& ptile, F const& f)
65{
66 const auto& ptd = ptile.getConstParticleTileData();
67 const auto np = ptile.numParticles();
68 pflags.resize(np, 0);
69 auto flag_ptr = pflags.data();
70 for (int k = 0; k < np; ++k) {
71 const auto p = ptd.getSuperParticle(k);
72 if constexpr (IsCallable<F,decltype(p),RandomEngine>::value) {
73 flag_ptr[k] = f(p,getInvalidRandomEngine());
74 } else if constexpr (IsCallable<F,decltype(p)>::value) {
75 flag_ptr[k] = f(p);
76 } else if constexpr (IsCallable<F,decltype(ptd),int,RandomEngine>::value) {
77 flag_ptr[k] = f(ptd,k,getInvalidRandomEngine());
78 } else {
79 flag_ptr[k] = f(ptd,k);
80 }
81 }
82}
83
84template <template <class, class> class Container,
85 class Allocator,
86 class PTile,
87 class F>
88void
89fillFlags (Container<int, Allocator>& pflags, const PTile& ptile, F const& f)
90{
91 if constexpr (IsPolymorphicArenaAllocator<Allocator>::value) {
92 if (pflags.arena()->isManaged() || pflags.arena()->isDevice()) {
93 fillFlagsGpu(pflags, ptile, f);
94 } else {
95 fillFlagsCpu(pflags, ptile, f);
96 }
97 } else {
98 if constexpr (RunOnGpu<Allocator>::value) {
99 fillFlagsGpu(pflags, ptile, f);
100 } else {
101 fillFlagsCpu(pflags, ptile, f);
102 }
103 }
104}
105
106template <class Container, class PC>
108countFlagsGpu (const Vector<std::map<std::pair<int,int>,Container>>& particle_io_flags, const PC& pc)
109{
110 ReduceOps<ReduceOpSum> reduce_op;
111 ReduceData<Long> reduce_data(reduce_op);
112 using ReduceTuple = typename decltype(reduce_data)::Type;
113
114 for (int lev = 0; lev < pc.GetParticles().size(); lev++)
115 {
116 const auto& pmap = pc.GetParticles(lev);
117 for (const auto& kv : pmap)
118 {
119 const auto& pflags = particle_io_flags[lev].at(kv.first);
120 const auto flag_ptr = pflags.data();
121 reduce_op.eval(pflags.size(), reduce_data,
122 [=] AMREX_GPU_DEVICE (const int i) -> ReduceTuple
123 {
124 return flag_ptr[i] ? 1 : 0;
125 });
126 }
127 }
128 ReduceTuple hv = reduce_data.value(reduce_op);
129 return amrex::get<0>(hv);
130}
131
132template <class Container>
134countFlagsGpu (const Container& pflags)
135{
136 ReduceOps<ReduceOpSum> reduce_op;
137 ReduceData<Long> reduce_data(reduce_op);
138 using ReduceTuple = typename decltype(reduce_data)::Type;
139
140 const auto flag_ptr = pflags.data();
141 reduce_op.eval(pflags.size(), reduce_data,
142 [=] AMREX_GPU_DEVICE (const amrex::Long i) -> ReduceTuple
143 {
144 return flag_ptr[i] ? 1 : 0;
145 });
146 ReduceTuple hv = reduce_data.value(reduce_op);
147 return amrex::get<0>(hv);
148}
149
150template <class Container, class PC>
152countFlagsCpu (const Vector<std::map<std::pair<int,int>,Container>>& particle_io_flags, const PC& pc)
153{
154 amrex::Long nparticles = 0;
155 for (int lev = 0; lev < pc.GetParticles().size(); lev++)
156 {
157 const auto& pmap = pc.GetParticles(lev);
158 for (const auto& kv : pmap)
159 {
160 const auto& pflags = particle_io_flags[lev].at(kv.first);
161 for (int k = 0; k < kv.second.numParticles(); ++k)
162 {
163 if (pflags[k]) { nparticles++; }
164 }
165 }
166 }
167 return nparticles;
168}
169
170template <class Container>
172countFlagsCpu (const Container& pflags)
173{
174 amrex::Long nparticles = 0;
175 for (std::size_t k = 0; k < pflags.size(); ++k)
176 {
177 if (pflags[k]) { nparticles++; }
178 }
179 return nparticles;
180}
181
182template <template <class, class> class Container, class Allocator, class PC>
184countFlags (const Vector<std::map<std::pair<int,int>,Container<int,Allocator>>>& particle_io_flags, const PC& pc)
185{
186 if constexpr (IsPolymorphicArenaAllocator<Allocator>::value) {
187 if (pc.arena()->isManaged() || pc.arena()->isDevice()) {
188 return countFlagsGpu(particle_io_flags, pc);
189 } else {
190 return countFlagsCpu(particle_io_flags, pc);
191 }
192 } else {
193 if constexpr (RunOnGpu<Allocator>::value) {
194 return countFlagsGpu(particle_io_flags, pc);
195 } else {
196 return countFlagsCpu(particle_io_flags, pc);
197 }
198 }
199}
200
201template <template <class, class> class Container, class Allocator>
203countFlags (const Container<int,Allocator>& pflags)
204{
205 if constexpr (IsPolymorphicArenaAllocator<Allocator>::value) {
206 if (pflags.arena()->isManaged() || pflags.arena()->isDevice()) {
207 return countFlagsGpu(pflags);
208 } else {
209 return countFlagsCpu(pflags);
210 }
211 } else {
212 if constexpr (RunOnGpu<Allocator>::value) {
213 return countFlagsGpu(pflags);
214 } else {
215 return countFlagsCpu(pflags);
216 }
217 }
218}
219
220template <typename I>
222void packParticleIDs (I* idata, const std::uint64_t idcpu, bool is_checkpoint) noexcept
223{
224 if (is_checkpoint) {
225 std::int32_t xi, yi;
226 std::uint32_t xu, yu;
227 xu = (std::uint32_t)((idcpu & 0xFFFFFFFF00000000LL) >> 32);
228 yu = (std::uint32_t)( idcpu & 0xFFFFFFFFLL);
229 amrex::Gpu::memcpy(&xi, &xu, sizeof(xu));
230 amrex::Gpu::memcpy(&yi, &yu, sizeof(yu));
231 idata[0] = xi;
232 idata[1] = yi;
233 } else {
234 idata[0] = ConstParticleIDWrapper{idcpu};
235 idata[1] = ConstParticleCPUWrapper{idcpu};
236 }
237}
238
239template<class PTD>
241rPackParticleData (const PTD& ptd, int idx, std::remove_const_t<typename PTD::RealType> * rdata_ptr,
242 const int * write_real_comp)
243{
244 std::size_t rout_index = 0;
245
246 if constexpr (!PTD::ParticleType::is_rtsoa_particle) {
247 for (int j = 0; j < AMREX_SPACEDIM; ++j) {
248 rdata_ptr[rout_index] = ptd.pos(j, idx);
249 rout_index++;
250 }
251
252 if constexpr (!PTD::ParticleType::is_soa_particle) {
253 const auto& p = ptd[idx];
254
255 for (int j = 0; j < PTD::ParticleType::NReal; ++j) {
256 if (write_real_comp[j]) {
257 rdata_ptr[rout_index] = p.rdata(j);
258 rout_index++;
259 }
260 }
261 }
262
263 constexpr int real_start_offset = PTD::ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
264
265 if constexpr (PTD::NAR > 0) {
266 for (int j = real_start_offset; j < PTD::NAR; ++j) {
267 if (write_real_comp[PTD::ParticleType::NReal + j - real_start_offset]) {
268 rdata_ptr[rout_index] = ptd.rdata(j)[idx];
269 rout_index++;
270 }
271 }
272 }
273
274 for (int j = 0; j < ptd.m_num_runtime_real; ++j) {
275 if (write_real_comp[PTD::ParticleType::NReal + PTD::NAR + j - real_start_offset]) {
276 rdata_ptr[rout_index] = ptd.m_runtime_rdata[j][idx];
277 rout_index++;
278 }
279 }
280 } else {
281 for (int j = 0; j < ptd.m_n_real; ++j) {
282 if ((j < AMREX_SPACEDIM) || write_real_comp[j - AMREX_SPACEDIM]) {
283 rdata_ptr[rout_index] = ptd.rdata(j)[idx];
284 rout_index++;
285 }
286 }
287 }
288}
289
290template<class PTD>
292iPackParticleData (const PTD& ptd, int idx, std::remove_const_t<typename PTD::IntType> * idata_ptr,
293 const int * write_int_comp, bool is_checkpoint)
294{
295 std::size_t iout_index = 0;
296
297 packParticleIDs(&idata_ptr[iout_index], ptd.idcpu(idx), is_checkpoint);
298 iout_index += 2;
299
300 if constexpr (!PTD::ParticleType::is_rtsoa_particle) {
301 if constexpr (!PTD::ParticleType::is_soa_particle) {
302 const auto& p = ptd[idx];
303
304 for (int j = 0; j < PTD::ParticleType::NInt; ++j) {
305 if (write_int_comp[j]) {
306 idata_ptr[iout_index] = p.idata(j);
307 iout_index++;
308 }
309 }
310 }
311
312 if constexpr (PTD::NAI > 0) {
313 for (int j = 0; j < PTD::NAI; ++j) {
314 if (write_int_comp[PTD::ParticleType::NInt + j]) {
315 idata_ptr[iout_index] = ptd.idata(j)[idx];
316 iout_index++;
317 }
318 }
319 }
320
321 for (int j = 0; j < ptd.m_num_runtime_int; ++j) {
322 if (write_int_comp[PTD::ParticleType::NInt + PTD::NAI + j]) {
323 idata_ptr[iout_index] = ptd.m_runtime_idata[j][idx];
324 iout_index++;
325 }
326 }
327 } else {
328 for (int j = 0; j < ptd.m_n_int; ++j) {
329 if (write_int_comp[j]) {
330 idata_ptr[iout_index] = ptd.idata(j)[idx];
331 iout_index++;
332 }
333 }
334 }
335}
336
337template <class PC>
338void
339packIODataGpu (Vector<int>& idata, Vector<ParticleReal>& rdata, const PC& pc, int lev, int grid,
340 const Vector<int>& write_real_comp, const Vector<int>& write_int_comp,
341 const Vector<std::map<std::pair<int, int>, typename PC::FlagsVector>>& particle_io_flags,
342 const Vector<int>& tiles, int np, bool is_checkpoint)
343{
344 int num_output_int = 0;
345 for (int i = 0; i < pc.NumIntComps() + PC::NStructInt; ++i) {
346 if (write_int_comp[i]) { ++num_output_int; }
347 }
348
349 const Long iChunkSize = 2 + num_output_int;
350 idata.resize(np*iChunkSize);
351
352 int num_output_real = 0;
353 for (int i : write_real_comp) {
354 if (i) { ++num_output_real; }
355 }
356
357 const Long rChunkSize = AMREX_SPACEDIM + num_output_real;
358 rdata.resize(np*rChunkSize);
359
360 Gpu::DeviceVector<int> write_int_comp_d(write_int_comp.size());
361 Gpu::DeviceVector<int> write_real_comp_d(write_real_comp.size());
362 Gpu::copyAsync(Gpu::hostToDevice, write_int_comp.begin(), write_int_comp.end(),
363 write_int_comp_d.begin());
364 Gpu::copyAsync(Gpu::hostToDevice, write_real_comp.begin(), write_real_comp.end(),
365 write_real_comp_d.begin());
366
367 const auto* write_int_comp_d_ptr = write_int_comp_d.data();
368 const auto* write_real_comp_d_ptr = write_real_comp_d.data();
369
370 std::size_t poffset = 0;
371 for (int tile : tiles) {
372 const auto& ptile = pc.ParticlesAt(lev, grid, tile);
373 const auto& pflags = particle_io_flags[lev].at(std::make_pair(grid, tile));
374 int np_tile = ptile.numParticles();
375 Gpu::DeviceVector<int> offsets(np_tile);
376 int num_copies = Scan::ExclusiveSum(np_tile, pflags.begin(), offsets.begin(), Scan::retSum);
377
378 Gpu::DeviceVector<int> idata_d(num_copies*iChunkSize);
379 Gpu::DeviceVector<ParticleReal> rdata_d(num_copies*rChunkSize);
380
381 const auto* flag_ptr = pflags.data();
382 const auto* offset_ptr = offsets.data();
383
384 auto* idata_d_ptr = idata_d.data();
385 auto* rdata_d_ptr = rdata_d.data();
386
387 const auto& ptd = ptile.getConstParticleTileData();
388 amrex::ParallelFor(ptile.numParticles(),
389 [=] AMREX_GPU_DEVICE (int pindex) noexcept
390 {
391 if (flag_ptr[pindex]) {
392 const int out_indx = offset_ptr[pindex];
393 iPackParticleData(ptd, pindex, idata_d_ptr + out_indx * iChunkSize,
394 write_int_comp_d_ptr, is_checkpoint);
395
396 rPackParticleData(ptd, pindex, rdata_d_ptr + out_indx * rChunkSize,
397 write_real_comp_d_ptr);
398 }
399 });
400
401 Gpu::copyAsync(Gpu::deviceToHost, idata_d.begin(), idata_d.end(),
402 idata.begin() + static_cast<Long>(poffset * iChunkSize));
403 Gpu::copyAsync(Gpu::deviceToHost, rdata_d.begin(), rdata_d.end(),
404 rdata.begin() + static_cast<Long>(poffset * rChunkSize));
405 Gpu::Device::streamSynchronize();
406
407 poffset += num_copies;
408 }
409}
410
411template <class PC>
412void
413packIODataCpu (Vector<int>& idata, Vector<ParticleReal>& rdata, const PC& pc, int lev, int grid,
414 const Vector<int>& write_real_comp, const Vector<int>& write_int_comp,
415 const Vector<std::map<std::pair<int, int>, typename PC::FlagsVector>>& particle_io_flags,
416 const Vector<int>& tiles, int np, bool is_checkpoint)
417{
418 int num_output_int = 0;
419 for (int i = 0; i < pc.NumIntComps() + PC::NStructInt; ++i) {
420 if (write_int_comp[i]) { ++num_output_int; }
421 }
422
423 const Long iChunkSize = 2 + num_output_int;
424 idata.resize(np*iChunkSize);
425
426 int num_output_real = 0;
427 for (int i : write_real_comp) {
428 if (i) { ++num_output_real; }
429 }
430
431 const Long rChunkSize = AMREX_SPACEDIM + num_output_real;
432 rdata.resize(np*rChunkSize);
433
434 int* iptr = idata.dataPtr();
435 ParticleReal* rptr = rdata.dataPtr();
436 for (int tile : tiles) {
437 const auto& ptile = pc.ParticlesAt(lev, grid, tile);
438 const auto& pflags = particle_io_flags[lev].at(std::make_pair(grid, tile));
439 const auto& ptd = ptile.getConstParticleTileData();
440
441 for (int pindex = 0; pindex < ptile.numParticles(); ++pindex) {
442 if (pflags[pindex]) {
443 iPackParticleData(ptd, pindex, iptr,
444 write_int_comp.dataPtr(), is_checkpoint);
445 iptr += iChunkSize;
446
447 rPackParticleData(ptd, pindex, rptr,
448 write_real_comp.dataPtr());
449 rptr += rChunkSize;
450 }
451 }
452 }
453}
454
455template <class PC>
456void
457packIOData (Vector<int>& idata, Vector<ParticleReal>& rdata, const PC& pc, int lev, int grid,
458 const Vector<int>& write_real_comp, const Vector<int>& write_int_comp,
459 const Vector<std::map<std::pair<int, int>, typename PC::FlagsVector>>& particle_io_flags,
460 const Vector<int>& tiles, int np, bool is_checkpoint)
461{
462 if constexpr (IsPolymorphicArenaAllocator<typename PC::FlagsVector::allocator_type>::value) {
463 if (pc.arena()->isManaged() || pc.arena()->isDevice()) {
464 packIODataGpu(idata, rdata, pc, lev, grid, write_real_comp, write_int_comp,
465 particle_io_flags, tiles, np, is_checkpoint);
466 } else {
467 packIODataCpu(idata, rdata, pc, lev, grid, write_real_comp, write_int_comp,
468 particle_io_flags, tiles, np, is_checkpoint);
469 }
470 } else {
471 if constexpr (RunOnGpu<typename PC::FlagsVector::allocator_type>::value) {
472 packIODataGpu(idata, rdata, pc, lev, grid, write_real_comp, write_int_comp,
473 particle_io_flags, tiles, np, is_checkpoint);
474 } else {
475 packIODataCpu(idata, rdata, pc, lev, grid, write_real_comp, write_int_comp,
476 particle_io_flags, tiles, np, is_checkpoint);
477 }
478 }
479}
480
481}
483
484template <class PC, class F>
485requires (IsParticleContainer<PC>::value)
487 const std::string& dir, const std::string& name,
488 const Vector<int>& write_real_comp,
489 const Vector<int>& write_int_comp,
490 const Vector<std::string>& real_comp_names,
491 const Vector<std::string>& int_comp_names,
492 F const& f, bool is_checkpoint)
493{
494 BL_PROFILE("WriteBinaryParticleData()");
495 AMREX_ASSERT(pc.OK());
496
497 AMREX_ASSERT(sizeof(typename PC::ParticleType::RealType) == 4 ||
498 sizeof(typename PC::ParticleType::RealType) == 8);
499
500 constexpr int NStructReal = PC::NStructReal;
501 constexpr int NStructInt = PC::NStructInt;
502
503 const int NProcs = ParallelDescriptor::NProcs();
504 const int IOProcNumber = ParallelDescriptor::IOProcessorNumber();
505
506 if constexpr(PC::ParticleType::is_soa_particle) {
507 AMREX_ALWAYS_ASSERT(real_comp_names.size() == pc.NumRealComps() + NStructReal - AMREX_SPACEDIM); // pure SoA: skip positions
508 } else {
509 AMREX_ALWAYS_ASSERT(real_comp_names.size() == pc.NumRealComps() + NStructReal);
510 }
511 AMREX_ALWAYS_ASSERT( int_comp_names.size() == pc.NumIntComps() + NStructInt);
512
513 std::string pdir = dir;
514 if ( ! pdir.empty() && pdir[pdir.size()-1] != '/') { pdir += '/'; }
515 pdir += name;
516
517 if ( ! pc.GetLevelDirectoriesCreated()) {
519 {
520 if ( ! amrex::UtilCreateDirectory(pdir, 0755))
521 {
523 }
524 }
526 }
527
528 std::ofstream HdrFile;
529
530 Long nparticles = 0;
531 Long maxnextid;
532
533 // evaluate f for every particle to determine which ones to output
534 Vector<std::map<std::pair<int, int>, typename PC::FlagsVector > >
535 particle_io_flags(pc.GetParticles().size());
536
537 for (int lev = 0; lev < pc.GetParticles().size(); lev++)
538 {
539 const auto& pmap = pc.GetParticles(lev);
540 for (const auto& kv : pmap)
541 {
542 auto& flags = particle_io_flags[lev][kv.first];
543 if constexpr (PC::has_polymorphic_allocator) {
544 flags.setArena(pc.arena());
545 }
546 particle_detail::fillFlags(flags, kv.second, f);
547 }
548 }
549
551
552 if(pc.GetUsePrePost())
553 {
554 nparticles = pc.GetNParticlesPrePost();
555 maxnextid = pc.GetMaxNextIDPrePost();
556 }
557 else
558 {
559 nparticles = particle_detail::countFlags(particle_io_flags, pc);
560 maxnextid = PC::ParticleType::NextID();
561 ParallelDescriptor::ReduceLongSum(nparticles, IOProcNumber);
562 PC::ParticleType::NextID(maxnextid);
563 ParallelDescriptor::ReduceLongMax(maxnextid, IOProcNumber);
564 }
565
567 {
568 std::string HdrFileName = pdir;
569
570 if ( ! HdrFileName.empty() && HdrFileName[HdrFileName.size()-1] != '/') {
571 HdrFileName += '/';
572 }
573
574 HdrFileName += "Header";
575 pc.HdrFileNamePrePost = HdrFileName;
576
577 HdrFile.open(HdrFileName.c_str(), std::ios::out|std::ios::trunc);
578
579 if ( ! HdrFile.good()) { amrex::FileOpenFailed(HdrFileName); }
580
581 //
582 // First thing written is our version string.
583 // We append "_single" or "_double" to the version string indicating
584 // whether we're using "float" or "double" floating point data.
585 //
586 std::string version_string = is_checkpoint ? PC::CheckpointVersion() : PC::PlotfileVersion();
587 if (sizeof(typename PC::ParticleType::RealType) == 4)
588 {
589 HdrFile << version_string << "_single" << '\n';
590 }
591 else
592 {
593 HdrFile << version_string << "_double" << '\n';
594 }
595
596 int num_output_real = 0;
597 for (int i : write_real_comp) {
598 if (i) { ++num_output_real; }
599 }
600
601 int num_output_int = 0;
602 for (int i = 0; i < pc.NumIntComps() + NStructInt; ++i) {
603 if (write_int_comp[i]) { ++num_output_int; }
604 }
605
606 // AMREX_SPACEDIM and N for sanity checking.
607 HdrFile << AMREX_SPACEDIM << '\n';
608
609 // The number of extra real parameters
610 HdrFile << num_output_real << '\n';
611
612 // Real component names
613 for (int i = 0; i < (int) real_comp_names.size(); ++i ) {
614 if (write_real_comp[i]) { HdrFile << real_comp_names[i] << '\n'; }
615 }
616
617 // The number of extra int parameters
618 HdrFile << num_output_int << '\n';
619
620 // int component names
621 for (int i = 0; i < NStructInt + pc.NumIntComps(); ++i ) {
622 if (write_int_comp[i]) { HdrFile << int_comp_names[i] << '\n'; }
623 }
624
625 bool is_checkpoint_legacy = true; // legacy
626 HdrFile << is_checkpoint_legacy << '\n';
627
628 // The total number of particles.
629 HdrFile << nparticles << '\n';
630
631 // The value of nextid that we need to restore on restart.
632 HdrFile << maxnextid << '\n';
633
634 // Then the finest level of the AMR hierarchy.
635 HdrFile << pc.finestLevel() << '\n';
636
637 // Then the number of grids at each level.
638 for (int lev = 0; lev <= pc.finestLevel(); lev++) {
639 HdrFile << pc.ParticleBoxArray(lev).size() << '\n';
640 }
641 }
642
643 // We want to write the data out in parallel.
644 // We'll allow up to nOutFiles active writers at a time.
645 int nOutFiles(256);
646
647 ParmParse pp("particles");
648 pp.queryAdd("particles_nfiles",nOutFiles);
649 if(nOutFiles == -1) { nOutFiles = NProcs; }
650 nOutFiles = std::max(1, std::min(nOutFiles,NProcs));
651 pc.nOutFilesPrePost = nOutFiles;
652
653 for (int lev = 0; lev <= pc.finestLevel(); lev++)
654 {
655 bool gotsome;
656 if(pc.usePrePost)
657 {
658 gotsome = (pc.nParticlesAtLevelPrePost[lev] > 0);
659 }
660 else
661 {
662 gotsome = (pc.NumberOfParticlesAtLevel(lev) > 0);
663 }
664
665 // We store the particles at each level in their own subdirectory.
666 std::string LevelDir = pdir;
667
668 if (gotsome)
669 {
670 if ( ! LevelDir.empty() && LevelDir[LevelDir.size()-1] != '/') { LevelDir += '/'; }
671
672 LevelDir = amrex::Concatenate(LevelDir.append("Level_"), lev, 1);
673
674 if ( ! pc.GetLevelDirectoriesCreated())
675 {
677 if ( ! amrex::UtilCreateDirectory(LevelDir, 0755)) {
679 }
680 }
682 }
683 }
684
685 // Write out the header for each particle
686 if (gotsome && ParallelDescriptor::IOProcessor()) {
687 std::string HeaderFileName = LevelDir;
688 HeaderFileName += "/Particle_H";
689 std::ofstream ParticleHeader(HeaderFileName);
690
691 pc.ParticleBoxArray(lev).writeOn(ParticleHeader);
692 ParticleHeader << '\n';
693
694 ParticleHeader.flush();
695 ParticleHeader.close();
696 }
697
698 MFInfo info;
699 info.SetAlloc(false);
700 MultiFab state(pc.ParticleBoxArray(lev),
701 pc.ParticleDistributionMap(lev),
702 1,0,info);
703
704 // We eventually want to write out the file name and the offset
705 // into that file into which each grid of particles is written.
706 Vector<int> which(state.size(),0);
707 Vector<int > count(state.size(),0);
708 Vector<Long> where(state.size(),0);
709
710 std::string filePrefix(LevelDir);
711 filePrefix += '/';
712 filePrefix += PC::DataPrefix();
713 if(pc.usePrePost) {
714 pc.filePrefixPrePost[lev] = filePrefix;
715 }
716 bool groupSets(false), setBuf(true);
717
718 if (gotsome)
719 {
720 for(NFilesIter nfi(nOutFiles, filePrefix, groupSets, setBuf); nfi.ReadyToWrite(); ++nfi)
721 {
722 auto& myStream = (std::ofstream&) nfi.Stream();
723 pc.WriteParticles(lev, myStream, nfi.FileNumber(), which, count, where,
724 write_real_comp, write_int_comp, particle_io_flags, is_checkpoint);
725 }
726
727 if(pc.usePrePost) {
728 pc.whichPrePost[lev] = which;
729 pc.countPrePost[lev] = count;
730 pc.wherePrePost[lev] = where;
731 } else {
732 ParallelDescriptor::ReduceIntSum (which.dataPtr(), static_cast<int>(which.size()), IOProcNumber);
733 ParallelDescriptor::ReduceIntSum (count.dataPtr(), static_cast<int>(count.size()), IOProcNumber);
734 ParallelDescriptor::ReduceLongSum(where.dataPtr(), static_cast<int>(where.size()), IOProcNumber);
735 }
736 }
737
739 {
740 if(pc.GetUsePrePost()) {
741 // ---- write to the header and unlink in CheckpointPost
742 } else {
743 for (int j = 0; j < state.size(); j++)
744 {
745 HdrFile << which[j] << ' ' << count[j] << ' ' << where[j] << '\n';
746 }
747
748 if (gotsome && pc.doUnlink)
749 {
750 // Unlink any zero-length data files.
751 Vector<Long> cnt(nOutFiles,0);
752
753 for (int i = 0, N=static_cast<int>(count.size()); i < N; i++) {
754 cnt[which[i]] += count[i];
755 }
756
757 for (int i = 0, N=static_cast<int>(cnt.size()); i < N; i++)
758 {
759 if (cnt[i] == 0)
760 {
761 std::string FullFileName = NFilesIter::FileName(i, filePrefix);
762 FileSystem::Remove(FullFileName);
763 }
764 }
765 }
766 }
767 }
770 }
771 }
772
774 {
775 HdrFile.flush();
776 HdrFile.close();
777 if ( ! HdrFile.good())
778 {
779 amrex::Abort("amrex::WriteBinaryParticleDataSync(): problem writing HdrFile");
780 }
781 }
782}
783
784template <class PC>
785requires (IsParticleContainer<PC>::value)
787 const std::string& dir, const std::string& name,
788 const Vector<int>& write_real_comp,
789 const Vector<int>& write_int_comp,
790 const Vector<std::string>& real_comp_names,
791 const Vector<std::string>& int_comp_names, bool is_checkpoint)
792{
793 BL_PROFILE("WriteBinaryParticleDataAsync");
794 AMREX_ASSERT(pc.OK());
795
796 AMREX_ASSERT(sizeof(typename PC::ParticleType::RealType) == 4 ||
797 sizeof(typename PC::ParticleType::RealType) == 8);
798
799 constexpr int NStructReal = PC::NStructReal;
800 constexpr int NStructInt = PC::NStructInt;
801 constexpr int NArrayReal = PC::NArrayReal;
802 constexpr int NArrayInt = PC::NArrayInt;
803
804 const int MyProc = ParallelDescriptor::MyProc();
805 const int NProcs = ParallelDescriptor::NProcs();
806 const int IOProcNumber = NProcs - 1;
807
808 if constexpr(PC::ParticleType::is_soa_particle) {
809 AMREX_ALWAYS_ASSERT(real_comp_names.size() == pc.NumRealComps() + NStructReal - AMREX_SPACEDIM); // pure SoA: skip positions
810 } else {
811 AMREX_ALWAYS_ASSERT(real_comp_names.size() == pc.NumRealComps() + NStructReal);
812 }
813 AMREX_ALWAYS_ASSERT( int_comp_names.size() == pc.NumIntComps() + NStructInt);
814
815 Vector<LayoutData<Long> > np_per_grid_local(pc.finestLevel()+1);
816 for (int lev = 0; lev <= pc.finestLevel(); lev++)
817 {
818 np_per_grid_local[lev].define(pc.ParticleBoxArray(lev), pc.ParticleDistributionMap(lev));
819 using ParIter = typename PC::ParConstIterType;
820 for (ParIter pti(pc, lev); pti.isValid(); ++pti)
821 {
822 int gid = pti.index();
823 const auto& ptile = pc.ParticlesAt(lev, pti);
824 const auto& ptd = ptile.getConstParticleTileData();
825 const int np = ptile.numParticles();
826
827 ReduceOps<ReduceOpSum> reduce_op;
828 ReduceData<int> reduce_data(reduce_op);
829 using ReduceTuple = typename decltype(reduce_data)::Type;
830
831 reduce_op.eval(np, reduce_data,
832 [=] AMREX_GPU_DEVICE (int i) -> ReduceTuple
833 {
834 return (ptd.id(i).is_valid()) ? 1 : 0;
835 });
836
837 int np_valid = amrex::get<0>(reduce_data.value(reduce_op));
838 np_per_grid_local[lev][gid] += np_valid;
839 }
840 }
841
842 Vector<Vector<Long> > np_per_grid_global(pc.finestLevel()+1);
843 Long total_np = 0;
844 Vector<Long> np_per_level(pc.finestLevel()+1);
845 for (int lev = 0; lev <= pc.finestLevel(); lev++)
846 {
847 np_per_grid_global[lev].resize(np_per_grid_local[lev].size());
849 np_per_grid_global[lev],
850 IOProcNumber);
851 np_per_level[lev] = std::accumulate(np_per_grid_global[lev].begin(),
852 np_per_grid_global[lev].end(), 0L);
853 total_np += np_per_level[lev];
854 }
855
856 std::string pdir = dir;
857 if ( ! pdir.empty() && pdir[pdir.size()-1] != '/') { pdir += '/'; }
858 pdir += name;
859
860 if (MyProc == IOProcNumber)
861 {
862 if ( ! pc.GetLevelDirectoriesCreated())
863 {
864 if ( ! amrex::UtilCreateDirectory(pdir, 0755))
865 {
867 }
868 }
869
870 for (int lev = 0; lev <= pc.finestLevel(); lev++)
871 {
872 std::string LevelDir = pdir;
873 bool gotsome = np_per_level[lev];
874
875 if (gotsome)
876 {
877 if ( ! LevelDir.empty() && LevelDir[LevelDir.size()-1] != '/') { LevelDir += '/'; }
878
879 LevelDir = amrex::Concatenate(LevelDir.append("Level_"), lev, 1);
880
881 if ( ! pc.GetLevelDirectoriesCreated())
882 {
883 if ( ! amrex::UtilCreateDirectory(LevelDir, 0755))
884 {
886 }
887 }
888
889 std::string HeaderFileName = LevelDir;
890 HeaderFileName += "/Particle_H";
891 std::ofstream ParticleHeader(HeaderFileName);
892
893 pc.ParticleBoxArray(lev).writeOn(ParticleHeader);
894 ParticleHeader << '\n';
895
896 ParticleHeader.flush();
897 ParticleHeader.close();
898 }
899 }
900 }
902
903 Long maxnextid = PC::ParticleType::NextID();
904 ParallelDescriptor::ReduceLongMax(maxnextid, IOProcNumber);
905
906 Vector<Long> np_on_rank(NProcs, 0L);
907 std::size_t psize = particle_detail::PSizeInFile<ParticleReal>(write_real_comp, write_int_comp);
908 Vector<int64_t> rank_start_offset(NProcs);
909 if (MyProc == IOProcNumber)
910 {
911 for (int lev = 0; lev <= pc.finestLevel(); lev++)
912 {
913 for (int k = 0; k < pc.ParticleBoxArray(lev).size(); ++k)
914 {
915 int rank = pc.ParticleDistributionMap(lev)[k];
916 np_on_rank[rank] += np_per_grid_global[lev][k];
917 }
918 }
919
920 for (int ip = 0; ip < NProcs; ++ip)
921 {
922 auto info = AsyncOut::GetWriteInfo(ip);
923 rank_start_offset[ip] = (info.ispot == 0) ? 0 : static_cast<int64_t>(rank_start_offset[ip-1] + np_on_rank[ip-1]*psize);
924 }
925 }
926
927 // make tmp particle tiles in pinned memory to write
928 using PinnedPTile = std::conditional_t<
929 PC::is_rtsoa_pc,
930 typename PC::ParticleTileType,
932 >;
933 auto myptiles = std::make_shared<Vector<std::map<std::pair<int, int>,PinnedPTile> > >();
934 myptiles->resize(pc.finestLevel()+1);
935 for (int lev = 0; lev <= pc.finestLevel(); lev++)
936 {
937 for (MFIter mfi = pc.MakeMFIter(lev); mfi.isValid(); ++mfi)
938 {
939 auto& new_ptile = (*myptiles)[lev][std::make_pair(mfi.index(),
940 mfi.LocalTileIndex())];
941
942 if (np_per_grid_local[lev][mfi.index()] > 0)
943 {
944 const auto& ptile = pc.ParticlesAt(lev, mfi);
945
946 const auto np = np_per_grid_local[lev][mfi.index()];
947
948 const auto runtime_real_comps = ptile.NumRuntimeRealComps();
949 const auto runtime_int_comps = ptile.NumRuntimeIntComps();
950
951 new_ptile.define(runtime_real_comps, runtime_int_comps,
952 nullptr, nullptr, The_Pinned_Arena());
953
954 new_ptile.resize(np);
955
956 amrex::filterParticles(new_ptile, ptile, KeepValidFilter());
957 }
958 }
959 }
960
961 int finest_level = pc.finestLevel();
964 for (int lev = 0; lev <= pc.finestLevel(); lev++)
965 {
966 bas.push_back(pc.ParticleBoxArray(lev));
967 dms.push_back(pc.ParticleDistributionMap(lev));
968 }
969
970 int nic = pc.NumIntComps();
971 int rnames_size = (int) real_comp_names.size();
972
973 auto RD = pc.ParticleRealDescriptor;
974
975 AsyncOut::Submit([=] ()
976#if defined(__GNUC__) && (__GNUC__ == 8) && (__GNUC_MINOR__ == 1)
977 mutable // workaround for bug in gcc 8.1
978#endif
979 {
980 if (MyProc == IOProcNumber)
981 {
982 std::string HdrFileName = pdir;
983 std::ofstream HdrFile;
984
985 if ( ! HdrFileName.empty() && HdrFileName[HdrFileName.size()-1] != '/') {
986 HdrFileName += '/';
987 }
988
989 HdrFileName += "Header";
990
991 HdrFile.open(HdrFileName.c_str(), std::ios::out|std::ios::trunc);
992
993 if ( ! HdrFile.good()) { amrex::FileOpenFailed(HdrFileName); }
994
995 std::string version_string = is_checkpoint ? PC::CheckpointVersion() : PC::PlotfileVersion();
996 if (sizeof(typename PC::ParticleType::RealType) == 4)
997 {
998 HdrFile << version_string << "_single" << '\n';
999 }
1000 else
1001 {
1002 HdrFile << version_string << "_double" << '\n';
1003 }
1004
1005 int num_output_real = 0;
1006 for (int i = 0; i < rnames_size; ++i) {
1007 if (write_real_comp[i]) { ++num_output_real; }
1008 }
1009
1010 int num_output_int = 0;
1011 for (int i = 0; i < nic + NStructInt; ++i) {
1012 if (write_int_comp[i]) { ++num_output_int; }
1013 }
1014
1015 // AMREX_SPACEDIM and N for sanity checking.
1016 HdrFile << AMREX_SPACEDIM << '\n';
1017
1018 // The number of extra real parameters
1019 HdrFile << num_output_real << '\n';
1020
1021 // Real component names
1022 for (int i = 0; i < rnames_size; ++i ) {
1023 if (write_real_comp[i]) { HdrFile << real_comp_names[i] << '\n'; }
1024 }
1025
1026 // The number of extra int parameters
1027 HdrFile << num_output_int << '\n';
1028
1029 // int component names
1030 for (int i = 0; i < NStructInt + nic; ++i ) {
1031 if (write_int_comp[i]) { HdrFile << int_comp_names[i] << '\n'; }
1032 }
1033
1034 bool is_checkpoint_legacy = true; // legacy
1035 HdrFile << is_checkpoint_legacy << '\n';
1036
1037 // The total number of particles.
1038 HdrFile << total_np << '\n';
1039
1040 // The value of nextid that we need to restore on restart.
1041 HdrFile << maxnextid << '\n';
1042
1043 // Then the finest level of the AMR hierarchy.
1044 HdrFile << finest_level << '\n';
1045
1046 // Then the number of grids at each level.
1047 for (int lev = 0; lev <= finest_level; lev++) {
1048 HdrFile << dms[lev].size() << '\n';
1049 }
1050
1051 for (int lev = 0; lev <= finest_level; lev++)
1052 {
1053 Vector<int64_t> grid_offset(NProcs, 0);
1054 for (int k = 0; k < bas[lev].size(); ++k)
1055 {
1056 int rank = dms[lev][k];
1057 auto info = AsyncOut::GetWriteInfo(rank);
1058 HdrFile << info.ifile << ' '
1059 << np_per_grid_global[lev][k] << ' '
1060 << grid_offset[rank] + rank_start_offset[rank] << '\n';
1061 grid_offset[rank] += static_cast<int64_t>(np_per_grid_global[lev][k]*psize);
1062 }
1063 }
1064
1065 HdrFile.flush();
1066 HdrFile.close();
1067 if ( ! HdrFile.good())
1068 {
1069 amrex::Abort("amrex::WriteBinaryParticleDataAsync(): problem writing HdrFile");
1070 }
1071 }
1072
1073 AsyncOut::Wait(); // Wait for my turn
1074
1075 for (int lev = 0; lev <= finest_level; lev++)
1076 {
1077 // For a each grid, the tiles it contains
1078 std::map<int, Vector<int> > tile_map;
1079
1080 for (const auto& kv : (*myptiles)[lev])
1081 {
1082 const int grid = kv.first.first;
1083 const int tile = kv.first.second;
1084 tile_map[grid].push_back(tile);
1085 }
1086
1087 std::string LevelDir = pdir;
1088 if ( ! LevelDir.empty() && LevelDir[LevelDir.size()-1] != '/') { LevelDir += '/'; }
1089 LevelDir = amrex::Concatenate(LevelDir.append("Level_"), lev, 1);
1090 std::string filePrefix(LevelDir);
1091 filePrefix += '/';
1092 filePrefix += PC::DataPrefix();
1093 auto info = AsyncOut::GetWriteInfo(MyProc);
1094 std::string file_name = amrex::Concatenate(filePrefix, info.ifile, 5);
1095 std::ofstream ofs;
1096 ofs.open(file_name.c_str(), (info.ispot == 0) ? (std::ios::binary | std::ios::trunc)
1097 : (std::ios::binary | std::ios::app));
1098
1099 for (int k = 0; k < bas[lev].size(); ++k)
1100 {
1101 int rank = dms[lev][k];
1102 if (rank != MyProc) { continue; }
1103 const int grid = k;
1104 if (np_per_grid_local[lev][grid] == 0) { continue; }
1105
1106 // First write out the integer data in binary.
1107 int num_output_int = 0;
1108 for (int i = 0; i < nic + NStructInt; ++i) {
1109 if (write_int_comp[i]) { ++num_output_int; }
1110 }
1111
1112 const Long iChunkSize = 2 + num_output_int;
1113 Vector<int> istuff(np_per_grid_local[lev][grid]*iChunkSize);
1114 int* iptr = istuff.dataPtr();
1115
1116 for (unsigned i = 0; i < tile_map[grid].size(); i++) {
1117 auto ptile_index = std::make_pair(grid, tile_map[grid][i]);
1118 const auto& pbox = (*myptiles)[lev][ptile_index];
1119 const auto& ptd = pbox.getConstParticleTileData();
1120 for (int pindex = 0; pindex < pbox.numParticles(); ++pindex)
1121 {
1122 if (!ptd.id(pindex).is_valid()) { continue; }
1123
1124 particle_detail::iPackParticleData(ptd, pindex, iptr,
1125 write_int_comp.dataPtr(), is_checkpoint);
1126 iptr += iChunkSize;
1127 }
1128 }
1129
1130 writeIntData(istuff.dataPtr(), istuff.size(), ofs);
1131 ofs.flush(); // Some systems require this flush() (probably due to a bug)
1132
1133 // Write the Real data in binary.
1134 int num_output_real = 0;
1135 for (int i = 0; i < rnames_size; ++i) {
1136 if (write_real_comp[i]) { ++num_output_real; }
1137 }
1138
1139 const Long rChunkSize = AMREX_SPACEDIM + num_output_real;
1140 Vector<typename PC::ParticleType::RealType> rstuff(np_per_grid_local[lev][grid]*rChunkSize);
1141 typename PC::ParticleType::RealType* rptr = rstuff.dataPtr();
1142
1143 for (unsigned i = 0; i < tile_map[grid].size(); i++) {
1144 auto ptile_index = std::make_pair(grid, tile_map[grid][i]);
1145 const auto& pbox = (*myptiles)[lev][ptile_index];
1146 const auto& ptd = pbox.getConstParticleTileData();
1147 for (int pindex = 0; pindex < pbox.numParticles(); ++pindex)
1148 {
1149 if (!ptd.id(pindex).is_valid()) { continue; }
1150
1151 particle_detail::rPackParticleData(ptd, pindex, rptr,
1152 write_real_comp.dataPtr());
1153 rptr += rChunkSize;
1154 }
1155 }
1156
1157 if (sizeof(typename PC::ParticleType::RealType) == 4) {
1158 writeFloatData((float*) rstuff.dataPtr(), rstuff.size(), ofs, RD);
1159 }
1160 else if (sizeof(typename PC::ParticleType::RealType) == 8) {
1161 writeDoubleData((double*) rstuff.dataPtr(), rstuff.size(), ofs, RD);
1162 }
1163
1164 ofs.flush(); // Some systems require this flush() (probably due to a bug)
1165 }
1166 }
1167 AsyncOut::Notify(); // Notify others I am done
1168 });
1169}
1170
1171}
1172
1173#ifdef AMREX_USE_HDF5
1175#endif
1176
1177#endif /*AMREX_WRITE_BINARY_PARTICLE_DATA_H*/
#define BL_PROFILE(a)
Definition AMReX_BLProfiler.H:551
#define AMREX_ASSERT(EX)
Definition AMReX_BLassert.H:38
#define AMREX_ALWAYS_ASSERT(EX)
Definition AMReX_BLassert.H:50
#define AMREX_GPU_DEVICE
Definition AMReX_GpuQualifiers.H:18
#define AMREX_GPU_HOST_DEVICE
Definition AMReX_GpuQualifiers.H:20
amrex::ParmParse pp
Input file parser instance for the given namespace.
Definition AMReX_HypreIJIface.cpp:15
int size() const noexcept
Return the number of FABs in the FabArray.
Definition AMReX_FabArrayBase.H:110
static void streamSynchronize() noexcept
Definition AMReX_GpuDevice.cpp:856
Iterator for looping ever tiles and boxes of amrex::FabArray based containers.
Definition AMReX_MFIter.H:88
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition AMReX_MFIter.H:172
int index() const noexcept
The index into the underlying BoxArray of the current FAB.
Definition AMReX_MFIter.H:175
A collection (stored as an array) of FArrayBox objects.
Definition AMReX_MultiFab.H:40
This class encapsulates writing to nfiles.
Definition AMReX_NFiles.H:27
bool ReadyToWrite(bool appendFirst=false)
if appendFirst is true, the first set for this iterator will open the files in append mode
Definition AMReX_NFiles.cpp:205
const std::string & FileName() const
Definition AMReX_NFiles.H:160
Definition AMReX_ParIter.H:118
Parse Parameters From Command Line and Input Files.
Definition AMReX_ParmParse.H:351
int queryAdd(std::string_view name, T &ref)
If name is found, the value in the ParmParse database will be stored in the ref argument....
Definition AMReX_ParmParse.H:1047
Definition AMReX_Reduce.H:438
Type value()
Definition AMReX_Reduce.H:473
Definition AMReX_Reduce.H:597
void eval(MF const &mf, IntVect const &nghost, D &reduce_data, F &&f)
Definition AMReX_Reduce.H:731
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:29
T * dataPtr() noexcept
get access to the underlying data pointer
Definition AMReX_Vector.H:50
Long size() const noexcept
Definition AMReX_Vector.H:54
static bool GetBarrierAfterLevel()
Definition AMReX_VisMF.H:287
amrex_particle_real ParticleReal
Floating Point Type for Particles.
Definition AMReX_REAL.H:90
amrex_long Long
Definition AMReX_INT.H:30
Arena * The_Pinned_Arena()
Definition AMReX_Arena.cpp:860
int MyProc() noexcept
Definition AMReX_ParallelDescriptor.H:128
void Barrier(const std::string &)
Definition AMReX_ParallelDescriptor.cpp:1215
void ReduceIntSum(int &)
Definition AMReX_ParallelDescriptor.cpp:1265
int NProcs() noexcept
Definition AMReX_ParallelDescriptor.H:255
void ReduceLongSum(Long &)
Definition AMReX_ParallelDescriptor.cpp:1236
int IOProcessorNumber() noexcept
The MPI rank number of the I/O Processor (probably rank 0). This rank is usually used to write to std...
Definition AMReX_ParallelDescriptor.H:279
void ReduceLongMax(Long &)
Definition AMReX_ParallelDescriptor.cpp:1237
bool IOProcessor() noexcept
Is this CPU the I/O Processor? To get the rank number, call IOProcessorNumber()
Definition AMReX_ParallelDescriptor.H:289
WriteInfo GetWriteInfo(int rank)
Definition AMReX_AsyncOut.cpp:72
void Wait()
Definition AMReX_AsyncOut.cpp:112
void Notify()
Definition AMReX_AsyncOut.cpp:127
void Submit(std::function< void()> &&a_f)
Definition AMReX_AsyncOut.cpp:95
bool Remove(std::string const &filename)
Remove a file, symbolic link, or empty directory.
Definition AMReX_FileSystem.cpp:41
__host__ __device__ void * memcpy(void *dest, const void *src, std::size_t count)
Definition AMReX_GpuUtility.H:226
void GatherLayoutDataToVector(const LayoutData< T > &sendbuf, Vector< T > &recvbuf, int root)
Gather LayoutData values to a vector on root.
Definition AMReX_ParallelDescriptor.H:1277
Definition AMReX_Amr.cpp:50
void writeIntData(const From *data, std::size_t size, std::ostream &os, const amrex::IntDescriptor &id)
Definition AMReX_IntConv.H:23
__host__ __device__ void ignore_unused(const Ts &...)
This shuts up the compiler about unused variables.
Definition AMReX.H:139
void FileOpenFailed(const std::string &file)
Output a message and abort when couldn't open the file.
Definition AMReX_Utility.cpp:137
void WriteBinaryParticleDataAsync(PC const &pc, const std::string &dir, const std::string &name, const Vector< int > &write_real_comp, const Vector< int > &write_int_comp, const Vector< std::string > &real_comp_names, const Vector< std::string > &int_comp_names, bool is_checkpoint)
Definition AMReX_WriteBinaryParticleData.H:786
__host__ __device__ Dim3 begin(BoxND< dim > const &box) noexcept
Definition AMReX_Box.H:2018
void writeFloatData(const float *data, std::size_t size, std::ostream &os, const RealDescriptor &rd=FPC::Native32RealDescriptor())
Definition AMReX_VectorIO.cpp:130
void ParallelFor(TypeList< CTOs... > ctos, std::array< int, sizeof...(CTOs)> const &runtime_options, T N, F &&f)
Definition AMReX_CTOParallelForImpl.H:202
void WriteBinaryParticleDataSync(PC const &pc, const std::string &dir, const std::string &name, const Vector< int > &write_real_comp, const Vector< int > &write_int_comp, const Vector< std::string > &real_comp_names, const Vector< std::string > &int_comp_names, F const &f, bool is_checkpoint)
Definition AMReX_WriteBinaryParticleData.H:486
std::string Concatenate(const std::string &root, int num, int mindigits)
Returns rootNNNN where NNNN == num.
Definition AMReX_String.cpp:35
void CreateDirectoryFailed(const std::string &dir)
Output a message and abort when couldn't create the directory.
Definition AMReX_Utility.cpp:129
RandomEngine getInvalidRandomEngine()
Definition AMReX_RandomEngine.H:86
bool UtilCreateDirectory(const std::string &path, mode_t mode, bool verbose=false)
Creates the specified directories. path may be either a full pathname or a relative pathname....
Definition AMReX_Utility.cpp:116
void writeDoubleData(const double *data, std::size_t size, std::ostream &os, const RealDescriptor &rd=FPC::Native64RealDescriptor())
Definition AMReX_VectorIO.cpp:142
Index filterParticles(DstTile &dst, const SrcTile &src, const Index *mask) noexcept
Conditionally copy particles from src to dst based on the value of mask.
Definition AMReX_ParticleTransformation.H:389
AMREX_ATTRIBUTE_FLATTEN_FOR void ParallelForRNG(T n, L const &f) noexcept
Definition AMReX_GpuLaunchFunctsC.H:1151
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:241
const int[]
Definition AMReX_BLProfiler.cpp:1664
__host__ __device__ Dim3 end(BoxND< dim > const &box) noexcept
Definition AMReX_Box.H:2028
Definition AMReX_WriteBinaryParticleData.H:13
__host__ __device__ int operator()(const SrcData &src, int i) const noexcept
Definition AMReX_WriteBinaryParticleData.H:16
FabArray memory allocation information.
Definition AMReX_FabArray.H:68
MFInfo & SetAlloc(bool a) noexcept
Definition AMReX_FabArray.H:75
Metadata prefix of an AMReX particle "Header" file.
Definition AMReX_ParticleHeader.H:33
Definition AMReX_ParticleTile.H:764
void resize(std::size_t count, GrowthStrategy strategy=GrowthStrategy::Poisson)
Definition AMReX_ParticleTile.H:1014
Definition AMReX_RandomEngine.H:72