1#ifndef AMREX_PARTICLEHDF5_H
2#define AMREX_PARTICLEHDF5_H
4#include <AMReX_Config.H>
11#ifdef AMREX_USE_HDF5_ZFP
12#include "H5Zzfp_lib.h"
13#include "H5Zzfp_props.h"
16#ifdef AMREX_USE_HDF5_SZ
22template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
23 template<
class>
class Allocator,
class CellAssignor>
25ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
26::CheckpointHDF5 (
const std::string& dir,
27 const std::string& name,
bool ,
28 const Vector<std::string>& real_comp_names,
29 const Vector<std::string>& int_comp_names,
30 const std::string& compression)
const
32 Vector<int> write_real_comp;
33 Vector<std::string> tmp_real_comp_names;
36 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
37 for (
int i = first_rcomp; i < NStructReal + NumRealComps(); ++i )
39 write_real_comp.push_back(1);
40 if (real_comp_names.empty())
42 tmp_real_comp_names.push_back(getDefaultCompNameReal<ParticleType>(i));
46 tmp_real_comp_names.push_back(real_comp_names[i-first_rcomp]);
50 Vector<int> write_int_comp;
51 Vector<std::string> tmp_int_comp_names;
52 for (
int i = 0; i < NStructInt + NumIntComps(); ++i )
54 write_int_comp.push_back(1);
55 if (int_comp_names.empty())
57 tmp_int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
61 tmp_int_comp_names.push_back(int_comp_names[i]);
66 tmp_real_comp_names, tmp_int_comp_names,
74template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
75 template<
class>
class Allocator,
class CellAssignor>
77ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
78::CheckpointHDF5 (
const std::string& dir,
const std::string& name,
79 const std::string& compression)
const
81 Vector<int> write_real_comp;
82 Vector<std::string> real_comp_names;
84 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
85 for (
int i = first_rcomp; i < NStructReal + NumRealComps(); ++i )
87 write_real_comp.push_back(1);
88 real_comp_names.push_back(getDefaultCompNameReal<ParticleType>(i));
91 Vector<int> write_int_comp;
92 Vector<std::string> int_comp_names;
93 for (
int i = 0; i < NStructInt + NumIntComps(); ++i )
95 write_int_comp.push_back(1);
96 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
100 real_comp_names, int_comp_names,
108template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
109 template<
class>
class Allocator,
class CellAssignor>
111ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
112::WritePlotFileHDF5 (
const std::string& dir,
const std::string& name,
113 const std::string& compression)
const
115 Vector<int> write_real_comp;
116 Vector<std::string> real_comp_names;
118 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
119 for (
int i = first_rcomp; i < NStructReal + NumRealComps(); ++i )
121 write_real_comp.push_back(1);
122 real_comp_names.push_back(getDefaultCompNameReal<ParticleType>(i));
125 Vector<int> write_int_comp;
126 Vector<std::string> int_comp_names;
127 for (
int i = 0; i < NStructInt + NumIntComps(); ++i )
129 write_int_comp.push_back(1);
130 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
134 real_comp_names, int_comp_names,
142template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
143 template<
class>
class Allocator,
class CellAssignor>
145ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
146::WritePlotFileHDF5 (
const std::string& dir,
const std::string& name,
147 const Vector<std::string>& real_comp_names,
148 const Vector<std::string>& int_comp_names,
149 const std::string& compression)
const
152 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
153 AMREX_ALWAYS_ASSERT(
static_cast<int>(real_comp_names.size()) == NStructReal + NumRealComps() - first_rcomp);
154 AMREX_ALWAYS_ASSERT(
static_cast<int>(int_comp_names.size()) == NStructInt + NumIntComps() );
156 Vector<int> write_real_comp;
157 for (
int i = first_rcomp; i < NStructReal + NumRealComps(); ++i) { write_real_comp.push_back(1); }
159 Vector<int> write_int_comp;
160 for (
int i = 0; i < NStructInt + NumIntComps(); ++i) { write_int_comp.push_back(1); }
163 write_real_comp, write_int_comp,
164 real_comp_names, int_comp_names,
172template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
173 template<
class>
class Allocator,
class CellAssignor>
175ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
176::WritePlotFileHDF5 (
const std::string& dir,
const std::string& name,
177 const Vector<std::string>& real_comp_names,
178 const std::string& compression)
const
181 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
182 AMREX_ALWAYS_ASSERT(
static_cast<int>(real_comp_names.size()) == NStructReal + NumRealComps() - first_rcomp);
184 Vector<int> write_real_comp;
185 for (
int i = first_rcomp; i < NStructReal + NumRealComps(); ++i) write_real_comp.push_back(1);
187 Vector<int> write_int_comp;
188 for (
int i = 0; i < NStructInt + NumIntComps(); ++i) write_int_comp.push_back(1);
190 Vector<std::string> int_comp_names;
191 for (
int i = 0; i < NStructInt + NumIntComps(); ++i )
193 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
197 write_real_comp, write_int_comp,
198 real_comp_names, int_comp_names,
206template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
207 template<
class>
class Allocator,
class CellAssignor>
209ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
210::WritePlotFileHDF5 (
const std::string& dir,
211 const std::string& name,
212 const Vector<int>& write_real_comp,
213 const Vector<int>& write_int_comp,
214 const std::string& compression)
const
217 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
218 AMREX_ALWAYS_ASSERT(
static_cast<int>(write_real_comp.size()) == NStructReal + NumRealComps() - first_rcomp);
219 AMREX_ALWAYS_ASSERT(
static_cast<int>(write_int_comp.size()) == NStructInt + NumIntComps() );
221 Vector<std::string> real_comp_names;
222 for (
int i = first_rcomp; i < NStructReal + NumRealComps(); ++i )
224 real_comp_names.push_back(getDefaultCompNameReal<ParticleType>(i));
227 Vector<std::string> int_comp_names;
228 for (
int i = 0; i < NStructInt + NumIntComps(); ++i )
230 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
234 real_comp_names, int_comp_names,
242template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
243 template<
class>
class Allocator,
class CellAssignor>
245ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>::
246WritePlotFileHDF5 (
const std::string& dir,
const std::string& name,
247 const Vector<int>& write_real_comp,
248 const Vector<int>& write_int_comp,
249 const Vector<std::string>& real_comp_names,
250 const Vector<std::string>& int_comp_names,
251 const std::string& compression)
const
253 BL_PROFILE(
"ParticleContainer::WritePlotFile()");
256 write_real_comp, write_int_comp,
257 real_comp_names, int_comp_names,
265template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
266 template<
class>
class Allocator,
class CellAssignor>
267template <
class F, std::enable_if_t<!std::is_same_v<F, Vector<std::
string>>>*>
269ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
270::WritePlotFileHDF5 (
const std::string& dir,
const std::string& name,
271 const std::string& compression,
F&& f)
const
273 Vector<int> write_real_comp;
274 Vector<std::string> real_comp_names;
276 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
277 for (
int i = first_rcomp; i < NStructReal + NumRealComps(); ++i )
279 write_real_comp.push_back(1);
280 real_comp_names.push_back(getDefaultCompNameReal<ParticleType>(i));
283 Vector<int> write_int_comp;
284 Vector<std::string> int_comp_names;
285 for (
int i = 0; i < NStructInt + NumIntComps(); ++i )
287 write_int_comp.push_back(1);
288 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
292 real_comp_names, int_comp_names, compression,
296template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
297 template<
class>
class Allocator,
class CellAssignor>
300ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
301::WritePlotFileHDF5 (
const std::string& dir,
const std::string& name,
302 const Vector<std::string>& real_comp_names,
303 const Vector<std::string>& int_comp_names,
304 const std::string& compression,
F&& f)
const
307 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
308 AMREX_ALWAYS_ASSERT(
static_cast<int>(real_comp_names.size()) == NStructReal + NumRealComps() - first_rcomp);
309 AMREX_ALWAYS_ASSERT(
static_cast<int>(int_comp_names.size()) == NStructInt + NumIntComps() );
311 Vector<int> write_real_comp;
312 for (
int i = first_rcomp; i < NStructReal + NumRealComps(); ++i) { write_real_comp.push_back(1); }
314 Vector<int> write_int_comp;
315 for (
int i = 0; i < NStructInt + NumIntComps(); ++i) { write_int_comp.push_back(1); }
318 write_real_comp, write_int_comp,
319 real_comp_names, int_comp_names,
320 compression, std::forward<F>(f));
323template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
324 template<
class>
class Allocator,
class CellAssignor>
325template <
class F, std::enable_if_t<!std::is_same_v<F, Vector<std::
string>>>*>
327ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
328::WritePlotFileHDF5 (
const std::string& dir,
const std::string& name,
329 const Vector<std::string>& real_comp_names,
330 const std::string& compression,
F&& f)
const
333 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
334 AMREX_ALWAYS_ASSERT(
static_cast<int>(real_comp_names.size()) == NStructReal + NumRealComps() - first_rcomp);
336 Vector<int> write_real_comp;
337 for (
int i = first_rcomp; i < NStructReal + NumRealComps(); ++i) { write_real_comp.push_back(1); }
339 Vector<int> write_int_comp;
340 for (
int i = 0; i < NStructInt + NumIntComps(); ++i) { write_int_comp.push_back(1); }
342 Vector<std::string> int_comp_names;
343 for (
int i = 0; i < NStructInt + NumIntComps(); ++i )
345 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
349 write_real_comp, write_int_comp,
350 real_comp_names, int_comp_names,
351 compression, std::forward<F>(f));
354template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
355 template<
class>
class Allocator,
class CellAssignor>
358ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
359::WritePlotFileHDF5 (
const std::string& dir,
360 const std::string& name,
361 const Vector<int>& write_real_comp,
362 const Vector<int>& write_int_comp,
363 const std::string& compression,
F&& f)
const
366 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
367 AMREX_ALWAYS_ASSERT(
static_cast<int>(write_real_comp.size()) == NStructReal + NumRealComps() - first_rcomp);
368 AMREX_ALWAYS_ASSERT(
static_cast<int>(write_int_comp.size()) == NStructInt + NumIntComps() );
370 Vector<std::string> real_comp_names;
371 for (
int i = first_rcomp; i < NStructReal + NumRealComps(); ++i )
373 real_comp_names.push_back(getDefaultCompNameReal<ParticleType>(i));
376 Vector<std::string> int_comp_names;
377 for (
int i = 0; i < NStructInt + NumIntComps(); ++i )
379 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
383 real_comp_names, int_comp_names,
384 compression, std::forward<F>(f));
387template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
388 template<
class>
class Allocator,
class CellAssignor>
391ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>::
392WritePlotFileHDF5 (
const std::string& dir,
const std::string& name,
393 const Vector<int>& write_real_comp,
394 const Vector<int>& write_int_comp,
395 const Vector<std::string>& real_comp_names,
396 const Vector<std::string>& int_comp_names,
397 const std::string& compression,
F&& f)
const
399 BL_PROFILE(
"ParticleContainer::WritePlotFile()");
402 write_real_comp, write_int_comp,
403 real_comp_names, int_comp_names,
404 compression, std::forward<F>(f));
407template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
408 template<
class>
class Allocator,
class CellAssignor>
411ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
412::WriteHDF5ParticleData (
const std::string& dir,
const std::string& name,
413 const Vector<int>& write_real_comp,
414 const Vector<int>& write_int_comp,
415 const Vector<std::string>& real_comp_names,
416 const Vector<std::string>& int_comp_names,
417 const std::string& compression,
418 F&& f,
bool is_checkpoint)
const
428 write_real_comp, write_int_comp,
429 real_comp_names, int_comp_names,
431 std::forward<F>(f), is_checkpoint);
435template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
436 template<
class>
class Allocator,
class CellAssignor>
438ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
439::CheckpointPreHDF5 ()
445 BL_PROFILE(
"ParticleContainer::CheckpointPre()");
449 Long maxnextid = ParticleType::NextID();
451 for (
int lev = 0; lev < m_particles.size(); lev++) {
452 const auto& pmap = m_particles[lev];
453 for (
const auto& kv : pmap) {
454 const auto& aos = kv.second.GetArrayOfStructs();
455 for (
int k = 0; k < aos.numParticles(); ++k) {
456 const ParticleType& p = aos[k];
457 if (p.id().is_valid()) {
468 ParticleType::NextID(maxnextid);
471 nparticlesPrePost = nparticles;
472 maxnextidPrePost = maxnextid;
474 nParticlesAtLevelPrePost.clear();
475 nParticlesAtLevelPrePost.resize(finestLevel() + 1, 0);
476 for(
int lev(0); lev <= finestLevel(); ++lev) {
477 nParticlesAtLevelPrePost[lev] = NumberOfParticlesAtLevel(lev);
480 whichPrePost.clear();
481 whichPrePost.resize(finestLevel() + 1);
482 countPrePost.clear();
483 countPrePost.resize(finestLevel() + 1);
484 wherePrePost.clear();
485 wherePrePost.resize(finestLevel() + 1);
487 filePrefixPrePost.clear();
488 filePrefixPrePost.resize(finestLevel() + 1);
492template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
493 template<
class>
class Allocator,
class CellAssignor>
495ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
496::CheckpointPostHDF5 ()
502 BL_PROFILE(
"ParticleContainer::CheckpointPostHDF5()");
505 std::ofstream HdrFile;
506 HdrFile.open(HdrFileNamePrePost.c_str(), std::ios::out | std::ios::app);
508 for(
int lev(0); lev <= finestLevel(); ++lev) {
515 for(
int j(0); j < whichPrePost[lev].size(); ++j) {
516 HdrFile << whichPrePost[lev][j] <<
' ' << countPrePost[lev][j] <<
' ' << wherePrePost[lev][j] <<
'\n';
519 const bool gotsome = (nParticlesAtLevelPrePost[lev] > 0);
520 if(gotsome && doUnlink) {
523 Vector<Long> cnt(nOutFilesPrePost,0);
525 for(
int i(0), N = countPrePost[lev].size(); i < N; ++i) {
526 cnt[whichPrePost[lev][i]] += countPrePost[lev][i];
529 for(
int i(0), N = cnt.size(); i < N; ++i) {
542 if( ! HdrFile.good()) {
543 amrex::Abort(
"ParticleContainer::CheckpointPostHDF5(): problem writing HdrFile");
548template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
549 template<
class>
class Allocator,
class CellAssignor>
551ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
552::WritePlotFilePreHDF5 ()
558template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
559 template<
class>
class Allocator,
class CellAssignor>
561ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
562::WritePlotFilePostHDF5 ()
568template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
569 template<
class>
class Allocator,
class CellAssignor>
571ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
572::WriteParticlesHDF5 (
int lev, hid_t grp,
573 Vector<int>& , Vector<int>& count, Vector<Long>& ,
574 const Vector<int>& write_real_comp,
575 const Vector<int>& write_int_comp,
576 const std::string& compression,
577 const Vector<std::map<std::pair<int, int>, IntVector>>& particle_io_flags,
578 bool is_checkpoint)
const
580 BL_PROFILE(
"ParticleContainer::WriteParticlesHDF5()");
583 std::map<int, Vector<int> > tile_map;
586 hid_t dxpl_col, dxpl_ind, dcpl_int, dcpl_real;
587 dxpl_col = H5Pcreate(H5P_DATASET_XFER);
588 dxpl_ind = H5Pcreate(H5P_DATASET_XFER);
590 H5Pset_dxpl_mpio(dxpl_col, H5FD_MPIO_COLLECTIVE);
592 dcpl_int = H5Pcreate(H5P_DATASET_CREATE);
593 dcpl_real = H5Pcreate(H5P_DATASET_CREATE);
595 std::string mode_env, value_env;
596 double comp_value = -1.0;
597 std::string::size_type pos = compression.find(
'@');
598 if (pos != std::string::npos) {
599 mode_env = compression.substr(0, pos);
600 value_env = compression.substr(pos+1);
601 if (!value_env.empty()) {
602 comp_value = atof(value_env.c_str());
606 H5Pset_alloc_time(dcpl_int, H5D_ALLOC_TIME_LATE);
607 H5Pset_alloc_time(dcpl_real, H5D_ALLOC_TIME_LATE);
609 if (!mode_env.empty() && mode_env !=
"None") {
610 const char *chunk_env = NULL;
611 hsize_t chunk_dim = 1024;
612 chunk_env = getenv(
"HDF5_CHUNK_SIZE");
613 if (chunk_env != NULL) {
614 chunk_dim = atoi(chunk_env);
617 H5Pset_chunk(dcpl_int, 1, &chunk_dim);
618 H5Pset_chunk(dcpl_real, 1, &chunk_dim);
620#ifdef AMREX_USE_HDF5_ZFP
621 pos = mode_env.find(
"ZFP");
622 if (pos != std::string::npos) {
623 ret = H5Z_zfp_initialize();
624 if (ret < 0) {
amrex::Abort(
"ZFP initialize failed!"); }
628 if (mode_env ==
"ZLIB") {
629 H5Pset_shuffle(dcpl_int);
630 H5Pset_shuffle(dcpl_real);
631 H5Pset_deflate(dcpl_int, (
int)comp_value);
632 H5Pset_deflate(dcpl_real, (
int)comp_value);
634#ifdef AMREX_USE_HDF5_SZ
635 else if (mode_env ==
"SZ") {
636 ret = H5Z_SZ_Init((
char*)value_env.c_str());
638 std::cout <<
"SZ config file:" << value_env.c_str() << std::endl;
639 amrex::Abort(
"SZ initialize failed, check SZ config file!");
643#ifdef AMREX_USE_HDF5_ZFP
644 else if (mode_env ==
"ZFP_RATE") {
645 H5Pset_zfp_rate(dcpl_int, comp_value);
646 H5Pset_zfp_rate(dcpl_real, comp_value);
648 else if (mode_env ==
"ZFP_PRECISION") {
649 H5Pset_zfp_precision(dcpl_int, (
unsigned int)comp_value);
650 H5Pset_zfp_precision(dcpl_real, (
unsigned int)comp_value);
652 else if (mode_env ==
"ZFP_ACCURACY") {
653 H5Pset_zfp_accuracy(dcpl_int, comp_value);
654 H5Pset_zfp_accuracy(dcpl_real, comp_value);
656 else if (mode_env ==
"ZFP_REVERSIBLE") {
657 H5Pset_zfp_reversible(dcpl_int);
658 H5Pset_zfp_reversible(dcpl_real);
662 std::cout <<
"\nHDF5 particle checkpoint using " << mode_env <<
", " <<
663 value_env <<
", " << chunk_dim << std::endl;
667 for (
const auto& kv : m_particles[lev])
669 const int grid = kv.first.first;
670 const int tile = kv.first.second;
671 tile_map[grid].push_back(tile);
672 const auto& pflags = particle_io_flags[lev].at(kv.first);
675 count[grid] += particle_detail::countFlags(pflags);
679 info.SetAlloc(
false);
680 MultiFab state(ParticleBoxArray(lev), ParticleDistributionMap(lev), 1,0,info);
683 ULong my_mfi_int_total_size = 0, my_mfi_real_total_size = 0, int_size, real_size;
685 Vector<int> my_mfi_real_size;
686 Vector<int> my_mfi_int_size;
687 Vector<int> my_nparticles;
690 hid_t offset_id, offset_space, real_mem_space, real_dset_space, real_dset_id;
691 hid_t int_mem_space, int_dset_id, int_dset_space;
692 hsize_t total_mfi = 0, total_real_size = 0, total_int_size = 0, real_file_offset = 0, int_file_offset = 0;
693 hsize_t my_int_offset, my_int_count, my_real_offset, my_real_count;
697 int real_comp_count = AMREX_SPACEDIM;
699 for (
int i = 0; i < static_cast<int>(write_real_comp.size()); ++i ) {
700 if (write_real_comp[i]) { ++real_comp_count; }
703 int int_comp_count = 2;
705 for (
int i = 0; i < NStructInt + NumIntComps(); ++i ) {
706 if (write_int_comp[i]) { ++int_comp_count; }
710 for (MFIter mfi(state); mfi.isValid(); ++mfi) {
711 const int grid = mfi.index();
712 if (count[grid] == 0)
715 int_size = count[grid] * int_comp_count;
716 my_mfi_int_size.push_back(int_size);
717 my_nparticles.push_back(count[grid]);
718 my_mfi_int_total_size += int_size;
721 real_size = count[grid] * real_comp_count;
722 my_mfi_real_size.push_back(real_size);
723 my_mfi_real_total_size += real_size;
733 total_mfi += all_mfi_cnt[i];
739 all_mfi_cnt[0] = my_mfi_cnt;
740 all_mfi_int_total_size[0] = my_mfi_int_total_size;
745 int_file_offset += all_mfi_int_total_size[i];
746 my_int_offset = int_file_offset;
750 total_int_size += all_mfi_int_total_size[i];
764 if (H5Pget_layout(dcpl_int) == H5D_CHUNKED) {
765 if (H5Pget_chunk(dcpl_int, 1, &chunk_size) > -1) {
766 if (chunk_size > total_int_size) {
767 H5Pset_chunk(dcpl_int, 1, &total_int_size);
772 int_dset_space = H5Screate_simple(1, &total_int_size, NULL);
773#ifdef AMREX_USE_HDF5_ASYNC
774 int_dset_id = H5Dcreate_async(grp,
"data:datatype=0", H5T_NATIVE_INT, int_dset_space, H5P_DEFAULT, dcpl_int, H5P_DEFAULT, es_par_g);
776 int_dset_id = H5Dcreate(grp,
"data:datatype=0", H5T_NATIVE_INT, int_dset_space, H5P_DEFAULT, dcpl_int, H5P_DEFAULT);
779 H5Sclose(int_dset_space);
786 all_mfi_real_total_size[0] = my_mfi_real_total_size;
790 total_real_size += all_mfi_real_total_size[i];
792#ifdef AMREX_USE_HDF5_SZ
793 if (mode_env ==
"SZ") {
795 unsigned int* cd_values = NULL;
796 unsigned filter_config;
797 if (
sizeof(
typename ParticleType::RealType) == 4) {
798 SZ_metaDataToCdArray(&cd_nelmts, &cd_values, SZ_FLOAT, 0, 0, 0, 0, total_real_size);
799 H5Pset_filter(dcpl_real, H5Z_FILTER_SZ, H5Z_FLAG_MANDATORY, cd_nelmts, cd_values);
802 SZ_metaDataToCdArray(&cd_nelmts, &cd_values, SZ_DOUBLE, 0, 0, 0, 0, total_real_size);
803 H5Pset_filter(dcpl_real, H5Z_FILTER_SZ, H5Z_FLAG_MANDATORY, cd_nelmts, cd_values);
808 if (H5Pget_layout(dcpl_real) == H5D_CHUNKED) {
809 if (H5Pget_chunk(dcpl_real, 1, &chunk_size) > -1) {
810 if (chunk_size > total_real_size) {
811 H5Pset_chunk(dcpl_real, 1, &total_real_size);
816 real_dset_space = H5Screate_simple(1, &total_real_size, NULL);
817 if (
sizeof(
typename ParticleType::RealType) == 4) {
818#ifdef AMREX_USE_HDF5_ASYNC
819 real_dset_id = H5Dcreate_async(grp,
"data:datatype=1", H5T_NATIVE_FLOAT, real_dset_space,
820 H5P_DEFAULT, dcpl_real, H5P_DEFAULT, es_par_g);
822 real_dset_id = H5Dcreate(grp,
"data:datatype=1", H5T_NATIVE_FLOAT, real_dset_space,
823 H5P_DEFAULT, dcpl_real, H5P_DEFAULT);
827#ifdef AMREX_USE_HDF5_ASYNC
828 real_dset_id = H5Dcreate_async(grp,
"data:datatype=1", H5T_NATIVE_DOUBLE, real_dset_space,
829 H5P_DEFAULT, dcpl_real, H5P_DEFAULT, es_par_g);
831 real_dset_id = H5Dcreate(grp,
"data:datatype=1", H5T_NATIVE_DOUBLE, real_dset_space,
832 H5P_DEFAULT, dcpl_real, H5P_DEFAULT);
835 H5Sclose(real_dset_space);
837 real_file_offset = 0;
839 real_file_offset += all_mfi_real_total_size[i];
841 my_real_offset = real_file_offset;
844 int max_mfi_count = 0, write_count = 0;
846 if (max_mfi_count < all_mfi_cnt[i]) {
847 max_mfi_count = all_mfi_cnt[i];
852 for (MFIter mfi(state); mfi.isValid(); ++mfi)
854 const int grid = mfi.index();
856 if (count[grid] == 0) {
continue; }
859 Vector<ParticleReal> rstuff;
860 particle_detail::packIOData(istuff, rstuff, *
this, lev, grid,
861 write_real_comp, write_int_comp,
862 particle_io_flags, tile_map[grid], count[grid],
865 my_int_count = istuff.size();
866 int_mem_space = H5Screate_simple(1, &my_int_count, NULL);
869 int_dset_space = H5Screate_simple(1, &total_int_size, NULL);
870 H5Sselect_hyperslab (int_dset_space, H5S_SELECT_SET, &my_int_offset, NULL, &my_int_count, NULL);
872#ifdef AMREX_USE_HDF5_ASYNC
873 ret = H5Dwrite_async(int_dset_id, H5T_NATIVE_INT, int_mem_space, int_dset_space, dxpl_col, istuff.dataPtr(), es_par_g);
875 ret = H5Dwrite(int_dset_id, H5T_NATIVE_INT, int_mem_space, int_dset_space, dxpl_col, istuff.dataPtr());
877 if (ret < 0) {
amrex::Abort(
"H5Dwrite int_dset failed!"); }
879 H5Sclose(int_dset_space);
880 H5Sclose(int_mem_space);
882 my_int_offset += my_int_count;
884 my_real_count = rstuff.size();
885 real_mem_space = H5Screate_simple(1, &my_real_count, NULL);
888 real_dset_space = H5Screate_simple(1, &total_real_size, NULL);
889 H5Sselect_hyperslab (real_dset_space, H5S_SELECT_SET, &my_real_offset, NULL, &my_real_count, NULL);
890 if (
sizeof(
typename ParticleType::RealType) == 4) {
891#ifdef AMREX_USE_HDF5_ASYNC
892 ret = H5Dwrite_async(real_dset_id, H5T_NATIVE_FLOAT, real_mem_space, real_dset_space, dxpl_col, rstuff.dataPtr(), es_par_g);
894 ret = H5Dwrite(real_dset_id, H5T_NATIVE_FLOAT, real_mem_space, real_dset_space, dxpl_col, rstuff.dataPtr());
897#ifdef AMREX_USE_HDF5_ASYNC
898 ret = H5Dwrite_async(real_dset_id, H5T_NATIVE_DOUBLE, real_mem_space, real_dset_space, dxpl_col, rstuff.dataPtr(), es_par_g);
900 ret = H5Dwrite(real_dset_id, H5T_NATIVE_DOUBLE, real_mem_space, real_dset_space, dxpl_col, rstuff.dataPtr());
904 if (ret < 0) {
amrex::Abort(
"H5Dwrite real_dset failed!"); }
906 H5Sclose(real_mem_space);
907 H5Sclose(real_dset_space);
909 my_real_offset += my_real_count;
915 while (write_count < max_mfi_count) {
916 int_dset_space = H5Screate_simple(1, &total_int_size, NULL);
917 real_dset_space = H5Screate_simple(1, &total_real_size, NULL);
918 H5Sselect_none(int_dset_space);
919 H5Sselect_none(real_dset_space);
921#ifdef AMREX_USE_HDF5_ASYNC
922 H5Dwrite_async(int_dset_id, H5T_NATIVE_INT, int_dset_space, int_dset_space, dxpl_col, NULL, es_par_g);
923 if (
sizeof(
typename ParticleType::RealType) == 4) {
924 H5Dwrite_async(real_dset_id, H5T_NATIVE_FLOAT, real_dset_space, real_dset_space, dxpl_col, NULL, es_par_g);
926 H5Dwrite_async(real_dset_id, H5T_NATIVE_DOUBLE, real_dset_space, real_dset_space, dxpl_col, NULL, es_par_g);
929 H5Dwrite(int_dset_id, H5T_NATIVE_INT, int_dset_space, int_dset_space, dxpl_col, NULL);
930 if (
sizeof(
typename ParticleType::RealType) == 4) {
931 H5Dwrite(real_dset_id, H5T_NATIVE_FLOAT, real_dset_space, real_dset_space, dxpl_col, NULL);
933 H5Dwrite(real_dset_id, H5T_NATIVE_DOUBLE, real_dset_space, real_dset_space, dxpl_col, NULL);
937 H5Sclose(int_dset_space);
938 H5Sclose(real_dset_space);
943#ifdef AMREX_USE_HDF5_ASYNC
944 H5Dclose_async(real_dset_id, es_par_g);
945 H5Dclose_async(int_dset_id, es_par_g);
947 H5Dclose(real_dset_id);
948 H5Dclose(int_dset_id);
952 offset_space = H5Screate_simple(1, &total_mfi, NULL);
953#ifdef AMREX_USE_HDF5_ASYNC
954 offset_id = H5Dcreate_async(grp,
"nparticles_grid", H5T_NATIVE_INT, offset_space, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT, es_par_g);
956 offset_id = H5Dcreate(grp,
"nparticles_grid", H5T_NATIVE_INT, offset_space, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
961 my_int_offset += all_mfi_cnt[i];
963 my_int_count = my_mfi_cnt;
964 int_mem_space = H5Screate_simple(1, &my_int_count, NULL);
967 H5Sselect_hyperslab (offset_space, H5S_SELECT_SET, &my_int_offset, NULL, &my_int_count, NULL);
969#ifdef AMREX_USE_HDF5_ASYNC
970 ret = H5Dwrite_async(offset_id, H5T_NATIVE_INT, int_mem_space, offset_space, dxpl_col, my_nparticles.dataPtr(), es_par_g);
972 ret = H5Dwrite(offset_id, H5T_NATIVE_INT, int_mem_space, offset_space, dxpl_col, my_nparticles.dataPtr());
974 if (ret < 0) {
amrex::Abort(
"H5Dwrite offset failed!"); }
981 H5Sclose(int_mem_space);
982 H5Sclose(offset_space);
984#ifdef AMREX_USE_HDF5_ASYNC
985 H5Dclose_async(offset_id, es_par_g);
994template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
995 template<
class>
class Allocator,
class CellAssignor>
997ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
998::RestartHDF5 (
const std::string& dir,
const std::string& file,
bool )
1003template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
1004 template<
class>
class Allocator,
class CellAssignor>
1006ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
1007::RestartHDF5 (
const std::string& dir,
const std::string& file)
1009 BL_PROFILE(
"ParticleContainer::RestartHDF5()");
1015 std::string fullname = dir;
1016 if (!fullname.empty() && fullname[fullname.size()-1] !=
'/') {
1023 hid_t fid, dset, grp, fapl, attr, atype, dspace, int_dset, real_dset;
1026 fapl = H5Pcreate (H5P_FILE_ACCESS);
1033 fid = H5Fopen(fullname.c_str(), H5F_ACC_RDONLY, fapl);
1035 std::string msg(
"ParticleContainer::RestartHDF5(): unable to open file: ");
1040 attr = H5Aopen(fid,
"version_name", H5P_DEFAULT);
1042 std::string msg(
"ParticleContainer::RestartHDF5(): unable to open version attribute ");
1046 atype = H5Aget_type(attr);
1048 std::string msg(
"ParticleContainer::RestartHDF5(): unable to get type of attribute ");
1052 size_t attr_len = H5Tget_size(atype);
1053 if (attr_len == 0) {
1054 std::string msg(
"ParticleContainer::RestartHDF5(): unable to get size of attribute ");
1058 std::string version;
1059 version.resize(attr_len+1);
1061 ret = H5Aread(attr, atype, &version[0]);
1076 bool convert_ids =
false;
1077 if (version.find(
"Version_Two_Dot_One") != std::string::npos) {
1080 if (version.find(
"_single") != std::string::npos) {
1083 else if (version.find(
"_double") != std::string::npos) {
1087 std::string msg(
"ParticleContainer::Restart(): bad version string: ");
1092 std::string gname =
"Chombo_global";
1093 std::string aname =
"SpaceDim";
1095 grp = H5Gopen(fid, gname.c_str(), H5P_DEFAULT);
1097 std::string msg(
"ParticleContainer::RestartHDF5(): unable to open group ");
1101 ret =
ReadHDF5Attr(grp, aname.c_str(), &dm, H5T_NATIVE_INT);
1103 std::string msg(
"ParticleContainer::RestartHDF5(): unable to read attribute ");
1109 aname =
"num_component_real";
1111 ret =
ReadHDF5Attr(fid, aname.c_str(), &nr, H5T_NATIVE_INT);
1113 std::string msg(
"ParticleContainer::RestartHDF5(): unable to read attribute ");
1117 int nrc = ParticleType::is_soa_particle ? NStructReal + NumRealComps() - AMREX_SPACEDIM : NStructReal + NumRealComps();
1119 amrex::Abort(
"ParticleContainer::RestartHDF5(): nr not the expected value");
1122 aname =
"num_component_int";
1124 ret =
ReadHDF5Attr(fid, aname.c_str(), &ni, H5T_NATIVE_INT);
1126 std::string msg(
"ParticleContainer::RestartHDF5(): unable to read attribute ");
1130 if (ni != NStructInt + NumIntComps()) {
1131 amrex::Abort(
"ParticleContainer::Restart(): ni != NStructInt");
1134 aname =
"nparticles";
1136 ret =
ReadHDF5Attr(fid, aname.c_str(), &nparticles, H5T_NATIVE_LLONG);
1138 std::string msg(
"ParticleContainer::RestartHDF5(): unable to read attribute ");
1144 aname =
"maxnextid";
1146 ret =
ReadHDF5Attr(fid, aname.c_str(), &maxnextid, H5T_NATIVE_LLONG);
1148 std::string msg(
"ParticleContainer::RestartHDF5(): unable to read attribute ");
1153 ParticleType::NextID(maxnextid);
1155 aname =
"finest_level";
1156 int finest_level_in_file;
1157 ret =
ReadHDF5Attr(fid, aname.c_str(), &finest_level_in_file, H5T_NATIVE_INT);
1159 std::string msg(
"ParticleContainer::RestartHDF5(): unable to read attribute ");
1165 hid_t comp_dtype = H5Tcreate (H5T_COMPOUND, 2 * AMREX_SPACEDIM *
sizeof(
int));
1166 if (1 == AMREX_SPACEDIM) {
1167 H5Tinsert (comp_dtype,
"lo_i", 0 *
sizeof(
int), H5T_NATIVE_INT);
1168 H5Tinsert (comp_dtype,
"hi_i", 1 *
sizeof(
int), H5T_NATIVE_INT);
1170 else if (2 == AMREX_SPACEDIM) {
1171 H5Tinsert (comp_dtype,
"lo_i", 0 *
sizeof(
int), H5T_NATIVE_INT);
1172 H5Tinsert (comp_dtype,
"lo_j", 1 *
sizeof(
int), H5T_NATIVE_INT);
1173 H5Tinsert (comp_dtype,
"hi_i", 2 *
sizeof(
int), H5T_NATIVE_INT);
1174 H5Tinsert (comp_dtype,
"hi_j", 3 *
sizeof(
int), H5T_NATIVE_INT);
1176 else if (3 == AMREX_SPACEDIM) {
1177 H5Tinsert (comp_dtype,
"lo_i", 0 *
sizeof(
int), H5T_NATIVE_INT);
1178 H5Tinsert (comp_dtype,
"lo_j", 1 *
sizeof(
int), H5T_NATIVE_INT);
1179 H5Tinsert (comp_dtype,
"lo_k", 2 *
sizeof(
int), H5T_NATIVE_INT);
1180 H5Tinsert (comp_dtype,
"hi_i", 3 *
sizeof(
int), H5T_NATIVE_INT);
1181 H5Tinsert (comp_dtype,
"hi_j", 4 *
sizeof(
int), H5T_NATIVE_INT);
1182 H5Tinsert (comp_dtype,
"hi_k", 5 *
sizeof(
int), H5T_NATIVE_INT);
1186 Vector<BoxArray> particle_box_arrays(finest_level_in_file + 1);
1187 bool dual_grid =
false;
1189 for (
int lev = 0; lev <= finest_level_in_file; lev++)
1191 if (lev > finestLevel())
1197 gname =
"level_" + std::to_string(lev);
1198 grp = H5Gopen(fid, gname.c_str(), H5P_DEFAULT);
1200 std::string msg(
"ParticleContainer::RestartHDF5(): unable to open group ");
1205 dset = H5Dopen(grp,
"boxes", H5P_DEFAULT);
1206 dspace = H5Dget_space(dset);
1208 H5Sget_simple_extent_dims(dspace, &ngrid, NULL);
1210 Vector<int> boxes(ngrid*AMREX_SPACEDIM*2);
1211 ret = H5Dread(dset, comp_dtype, H5S_ALL, H5S_ALL, H5P_DEFAULT, &(boxes[0]));
1213 std::string msg(
"ParticleContainer::RestartHDF5(): unable to read nparticles_grid dataset");
1227 for (
int i = 0; i < (
int)ngrid; i++) {
1228 const int base = i * 2 * AMREX_SPACEDIM;
1233 boxes[base+AMREX_SPACEDIM+1],
1234 boxes[base+AMREX_SPACEDIM+2])});
1237 particle_box_arrays[lev] = BoxArray(std::move(bl));
1240 if (! particle_box_arrays[lev].CellEqual(ParticleBoxArray(lev)))
1245 for (
int lev = 0; lev <= finestLevel(); lev++) {
1246 SetParticleBoxArray(lev, particle_box_arrays[lev]);
1247 DistributionMapping pdm(particle_box_arrays[lev]);
1248 SetParticleDistributionMap(lev, pdm);
1252 Vector<int> ngrids(finest_level_in_file+1);
1253 for (
int lev = 0; lev <= finest_level_in_file; lev++) {
1254 gname =
"level_" + std::to_string(lev);
1255 grp = H5Gopen(fid, gname.c_str(), H5P_DEFAULT);
1257 std::string msg(
"ParticleContainer::RestartHDF5(): unable to open group ");
1263 ret =
ReadHDF5Attr(grp, aname.c_str(), &ngrids[lev], H5T_NATIVE_INT);
1265 std::string msg(
"ParticleContainer::RestartHDF5(): unable to read attribute ");
1271 if (lev <= finestLevel()) {
1272 AMREX_ASSERT(ngrids[lev] ==
int(ParticleBoxArray(lev).size()));
1280 if (finest_level_in_file > finestLevel()) {
1281 m_particles.resize(finest_level_in_file+1);
1284 for (
int lev = 0; lev <= finest_level_in_file; lev++) {
1286 gname =
"level_" + std::to_string(lev);
1287 grp = H5Gopen(fid, gname.c_str(), H5P_DEFAULT);
1289 std::string msg(
"ParticleContainer::RestartHDF5(): unable to open group ");
1294 dset = H5Dopen(grp,
"nparticles_grid", H5P_DEFAULT);
1295 Vector<int> count(ngrids[lev]);
1296 ret = H5Dread(dset, H5T_NATIVE_INT, H5S_ALL, H5S_ALL, H5P_DEFAULT, &(count[0]));
1298 std::string msg(
"ParticleContainer::RestartHDF5(): unable to read nparticles_grid dataset");
1303 Vector<hsize_t>
offset(ngrids[lev]);
1305 for (
int i = 1; i < ngrids[lev]; i++) {
1309 int_dset = H5Dopen(grp,
"data:datatype=0", H5P_DEFAULT);
1311 std::string msg(
"ParticleContainer::RestartHDF5(): unable to open int dataset");
1314 real_dset = H5Dopen(grp,
"data:datatype=1", H5P_DEFAULT);
1315 if (real_dset < 0) {
1316 std::string msg(
"ParticleContainer::RestartHDF5(): unable to open int dataset");
1320 Vector<int> grids_to_read;
1321 if (lev <= finestLevel()) {
1322 for (MFIter mfi(*m_dummy_mf[lev]); mfi.isValid(); ++mfi) {
1323 grids_to_read.push_back(mfi.index());
1331 const int NReaders = MaxReaders();
1332 if (rank >= NReaders) {
1334 H5Dclose(real_dset);
1339 const int Navg = ngrids[lev] / NReaders;
1340 const int Nleft = ngrids[lev] - Navg * NReaders;
1344 lo = rank*(Navg + 1);
1348 lo = rank * Navg + Nleft;
1352 for (
int i = lo; i < hi; ++i) {
1353 grids_to_read.push_back(i);
1357 for(
int igrid = 0; igrid < static_cast<int>(grids_to_read.size()); ++igrid) {
1358 const int grid = grids_to_read[igrid];
1360 if (how ==
"single") {
1361 ReadParticlesHDF5<float>(
offset[grid], count[grid], grid, lev, int_dset, real_dset, finest_level_in_file, convert_ids);
1363 else if (how ==
"double") {
1364 ReadParticlesHDF5<double>(
offset[grid], count[grid], grid, lev, int_dset, real_dset, finest_level_in_file, convert_ids);
1367 std::string msg(
"ParticleContainer::Restart(): bad parameter: ");
1374 H5Dclose(real_dset);
1378 H5Tclose(comp_dtype);
1386 if (m_verbose > 1) {
1389 amrex::Print() <<
"ParticleContainer::Restart() time: " << stoptime <<
'\n';
1394template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
1395 template<
class>
class Allocator,
class CellAssignor>
1396template <
class RTYPE>
1398ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
1399::ReadParticlesHDF5 (hsize_t
offset, hsize_t cnt,
int grd,
int lev,
1400 hid_t int_dset, hid_t real_dset,
int finest_level_in_file,
1403 BL_PROFILE(
"ParticleContainer::ReadParticlesHDF5()");
1407 hid_t int_dspace, int_fspace, real_dspace, real_fspace;
1412 const int iChunkSize = 2 + NStructInt + NumIntComps();
1413 Vector<int> istuff(cnt*iChunkSize);
1414 int_fspace = H5Dget_space(int_dset);
1415 hsize_t int_cnt = cnt*iChunkSize;
1416 hsize_t int_offset =
offset*iChunkSize;
1417 int_dspace = H5Screate_simple(1, &int_cnt, NULL);
1418 H5Sselect_hyperslab (int_fspace, H5S_SELECT_SET, &int_offset, NULL, &int_cnt, NULL);
1419 H5Dread(int_dset, H5T_NATIVE_INT, int_dspace, int_fspace, H5P_DEFAULT, istuff.dataPtr());
1421 H5Sclose(int_fspace);
1422 H5Sclose(int_dspace);
1427 const int rChunkSize = ParticleType::is_soa_particle ?
1428 NStructReal + NumRealComps() : AMREX_SPACEDIM + NStructReal + NumRealComps();
1429 Vector<RTYPE> rstuff(cnt*rChunkSize);
1430 real_fspace = H5Dget_space(real_dset);
1431 hsize_t real_cnt = cnt*rChunkSize;
1432 hsize_t real_offset =
offset*rChunkSize;
1433 real_dspace = H5Screate_simple(1, &real_cnt, NULL);
1434 H5Sselect_hyperslab (real_fspace, H5S_SELECT_SET, &real_offset, NULL, &real_cnt, NULL);
1435 if (
sizeof(RTYPE) == 4) {
1436 H5Dread(real_dset, H5T_NATIVE_FLOAT, real_dspace, real_fspace, H5P_DEFAULT, rstuff.dataPtr());
1438 H5Dread(real_dset, H5T_NATIVE_DOUBLE, real_dspace, real_fspace, H5P_DEFAULT, rstuff.dataPtr());
1441 H5Sclose(real_fspace);
1442 H5Sclose(real_dspace);
1445 int* iptr = istuff.dataPtr();
1446 RTYPE* rptr = rstuff.dataPtr();
1450 Particle<NStructReal, NStructInt> ptemp;
1451 ParticleLocData pld;
1453 Vector<std::map<std::pair<int, int>, Gpu::HostVector<ParticleType> > > host_particles;
1454 host_particles.reserve(15);
1455 host_particles.resize(finest_level_in_file+1);
1457 Vector<std::map<std::pair<int, int>,
1458 std::vector<Gpu::HostVector<ParticleReal> > > > host_real_attribs;
1459 host_real_attribs.reserve(15);
1460 host_real_attribs.resize(finest_level_in_file+1);
1462 Vector<std::map<std::pair<int, int>,
1463 std::vector<Gpu::HostVector<int> > > > host_int_attribs;
1464 host_int_attribs.reserve(15);
1465 host_int_attribs.resize(finest_level_in_file+1);
1467 Vector<std::map<std::pair<int, int>, Gpu::HostVector<uint64_t> > > host_idcpu;
1468 host_idcpu.reserve(15);
1469 host_idcpu.resize(finest_level_in_file+1);
1471 for (hsize_t i = 0; i < cnt; i++) {
1475 std::int32_t xi, yi;
1476 std::uint32_t xu, yu;
1479 std::memcpy(&xu, &xi,
sizeof(xi));
1480 std::memcpy(&yu, &yi,
sizeof(yi));
1481 ptemp.m_idcpu = ((std::uint64_t)xu) << 32 | yu;
1483 ptemp.id() = iptr[0];
1484 ptemp.cpu() = iptr[1];
1488 for (
int j = 0; j < NStructInt; j++)
1490 ptemp.idata(j) = *iptr;
1500 rptr += AMREX_SPACEDIM;
1502 for (
int j = 0; j < NStructReal; j++)
1508 locateParticle(ptemp, pld, 0, finestLevel(), 0);
1510 std::pair<int, int> ind(grd, pld.m_tile);
1512 host_real_attribs[lev][ind].resize(NumRealComps());
1513 host_int_attribs[lev][ind].resize(NumIntComps());
1515 if constexpr (!ParticleType::is_soa_particle)
1518 host_particles[lev][ind].push_back(ptemp);
1521 for (
int icomp = 0; icomp < NumRealComps(); icomp++) {
1522 host_real_attribs[lev][ind][icomp].push_back(
ParticleReal(*rptr));
1527 for (
int icomp = 0; icomp < NumIntComps(); icomp++) {
1528 host_int_attribs[lev][ind][icomp].push_back(*iptr);
1533 host_particles[lev][ind];
1535 for (
int j = 0; j < AMREX_SPACEDIM; j++) {
1536 host_real_attribs[lev][ind][j].push_back(ptemp.pos(j));
1539 host_idcpu[lev][ind].push_back(ptemp.m_idcpu);
1542 for (
int icomp = AMREX_SPACEDIM; icomp < NumRealComps(); icomp++) {
1543 host_real_attribs[lev][ind][icomp].push_back(
ParticleReal(*rptr));
1548 for (
int icomp = 0; icomp < NumIntComps(); icomp++) {
1549 host_int_attribs[lev][ind][icomp].push_back(*iptr);
1555 for (
int host_lev = 0; host_lev < static_cast<int>(host_particles.size()); ++host_lev)
1557 for (
auto& kv : host_particles[host_lev]) {
1558 auto grid = kv.first.first;
1559 auto tile = kv.first.second;
1560 const auto& src_tile = kv.second;
1562 auto& dst_tile = DefineAndReturnParticleTile(host_lev, grid, tile);
1563 auto old_size = dst_tile.size();
1564 auto new_size = old_size;
1565 if constexpr (!ParticleType::is_soa_particle)
1567 new_size += src_tile.size();
1570 new_size += host_real_attribs[host_lev][std::make_pair(grid,tile)][0].size();
1572 dst_tile.resize(new_size);
1574 if constexpr (!ParticleType::is_soa_particle)
1577 dst_tile.GetArrayOfStructs().begin() + old_size);
1580 host_idcpu[host_lev][std::make_pair(grid,tile)].
begin(),
1581 host_idcpu[host_lev][std::make_pair(grid,tile)].
end(),
1582 dst_tile.GetStructOfArrays().GetIdCPUData().begin() + old_size);
1585 for (
int i = 0; i < NumRealComps(); ++i) {
1587 host_real_attribs[host_lev][std::make_pair(grid,tile)][i].
begin(),
1588 host_real_attribs[host_lev][std::make_pair(grid,tile)][i].
end(),
1589 dst_tile.GetStructOfArrays().GetRealData(i).begin() + old_size);
1592 for (
int i = 0; i < NumIntComps(); ++i) {
1594 host_int_attribs[host_lev][std::make_pair(grid,tile)][i].
begin(),
1595 host_int_attribs[host_lev][std::make_pair(grid,tile)][i].
end(),
1596 dst_tile.GetStructOfArrays().GetIntData(i).begin() + old_size);
#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_HOST_DEVICE
Definition AMReX_GpuQualifiers.H:20
Array4< int const > offset
Definition AMReX_HypreMLABecLap.cpp:1139
void CheckpointPreHDF5()
Prepare metadata used by deferred HDF5 checkpoint writes.
void RestartHDF5(const std::string &dir, const std::string &file)
Restart from checkpoint.
void WriteHDF5ParticleData(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, const std::string &compression, F &&f, bool is_checkpoint=false) const
Writes particle data to disk in the HDF5 format.
void CheckpointPostHDF5()
Finalize a checkpoint header handed off by CheckpointPreHDF5().
#define AMREX_D_TERM(a, b, c)
Definition AMReX_SPACE.H:172
#define AMREX_D_DECL(a, b, c)
Definition AMReX_SPACE.H:171
const std::string & FileName() const
Definition AMReX_NFiles.H:160
This class provides the user with a few print options.
Definition AMReX_Print.H:35
amrex_ulong ULong
Unsigned integer type guaranteed to be wider than unsigned int.
Definition AMReX_INT.H:32
amrex_particle_real ParticleReal
Floating Point Type for Particles.
Definition AMReX_REAL.H:90
amrex_long Long
Definition AMReX_INT.H:30
int MyProc() noexcept
Definition AMReX_ParallelDescriptor.H:128
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
bool Remove(std::string const &filename)
Definition AMReX_FileSystem.cpp:190
void copyAsync(HostToDevice, InIter begin, InIter end, OutIter result) noexcept
A host-to-device copy routine. Note this is just a wrapper around memcpy, so it assumes contiguous st...
Definition AMReX_GpuContainers.H:228
static constexpr HostToDevice hostToDevice
Definition AMReX_GpuContainers.H:105
void streamSynchronize() noexcept
Definition AMReX_GpuDevice.H:310
MPI_Comm Communicator() noexcept
Definition AMReX_ParallelDescriptor.H:223
void ReduceRealMax(Vector< std::reference_wrapper< Real > > const &)
Definition AMReX_ParallelDescriptor.cpp:1228
__host__ __device__ bool is_valid(const uint64_t idcpu) noexcept
Definition AMReX_Particle.H:146
Definition AMReX_Amr.cpp:49
__host__ __device__ void ignore_unused(const Ts &...)
This shuts up the compiler about unused variables.
Definition AMReX.H:139
BoxND< 3 > Box
Box is an alias for amrex::BoxND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:30
__host__ __device__ Dim3 begin(BoxND< dim > const &box) noexcept
Definition AMReX_Box.H:2006
double second() noexcept
Definition AMReX_Utility.cpp:940
IntVectND< 3 > IntVect
IntVect is an alias for amrex::IntVectND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:33
void Error(const std::string &msg)
Print out message to cerr and exit via amrex::Abort().
Definition AMReX.cpp:234
static void SetHDF5fapl(hid_t fapl, MPI_Comm comm)
Definition AMReX_PlotFileUtilHDF5.cpp:111
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:240
const int[]
Definition AMReX_BLProfiler.cpp:1664
void WriteHDF5ParticleDataSync(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, const std::string &compression, F &&f, bool is_checkpoint)
Synchronous HDF5 particle writer shared by checkpoint and plotfile routines.
Definition AMReX_WriteBinaryParticleDataHDF5.H:182
__host__ __device__ Dim3 end(BoxND< dim > const &box) noexcept
Definition AMReX_Box.H:2015
static int ReadHDF5Attr(hid_t loc, const char *name, void *data, hid_t dtype)
Read an attribute from an HDF5 object into the supplied buffer.
Definition AMReX_WriteBinaryParticleDataHDF5.H:98
static MPI_Datatype type()