1#ifndef AMREX_PARTICLEHDF5_H
2#define AMREX_PARTICLEHDF5_H
4#include <AMReX_Config.H>
13#ifdef AMREX_USE_HDF5_ZFP
14#include "H5Zzfp_lib.h"
15#include "H5Zzfp_props.h"
18#ifdef AMREX_USE_HDF5_SZ
24template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
25 template<
class>
class Allocator,
class CellAssignor>
27ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
28::CheckpointHDF5 (
const std::string& dir,
29 const std::string& name,
bool ,
30 const Vector<std::string>& real_comp_names,
31 const Vector<std::string>& int_comp_names,
32 const std::string& compression)
const
34 Vector<int> write_real_comp;
35 Vector<std::string> tmp_real_comp_names;
38 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
39 for (
int i = first_rcomp; i < NStructReal + NumRealComps(); ++i )
41 write_real_comp.push_back(1);
42 if (real_comp_names.empty())
44 tmp_real_comp_names.push_back(getDefaultCompNameReal<ParticleType>(i));
48 tmp_real_comp_names.push_back(real_comp_names[i-first_rcomp]);
52 Vector<int> write_int_comp;
53 Vector<std::string> tmp_int_comp_names;
54 for (
int i = 0; i < NStructInt + NumIntComps(); ++i )
56 write_int_comp.push_back(1);
57 if (int_comp_names.empty())
59 tmp_int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
63 tmp_int_comp_names.push_back(int_comp_names[i]);
68 tmp_real_comp_names, tmp_int_comp_names,
76template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
77 template<
class>
class Allocator,
class CellAssignor>
79ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
80::CheckpointHDF5 (
const std::string& dir,
const std::string& name,
81 const std::string& compression)
const
83 Vector<int> write_real_comp;
84 Vector<std::string> real_comp_names;
86 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
87 for (
int i = first_rcomp; i < NStructReal + NumRealComps(); ++i )
89 write_real_comp.push_back(1);
90 real_comp_names.push_back(getDefaultCompNameReal<ParticleType>(i));
93 Vector<int> write_int_comp;
94 Vector<std::string> int_comp_names;
95 for (
int i = 0; i < NStructInt + NumIntComps(); ++i )
97 write_int_comp.push_back(1);
98 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
102 real_comp_names, int_comp_names,
110template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
111 template<
class>
class Allocator,
class CellAssignor>
113ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
114::WritePlotFileHDF5 (
const std::string& dir,
const std::string& name,
115 const std::string& compression)
const
117 Vector<int> write_real_comp;
118 Vector<std::string> real_comp_names;
120 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
121 for (
int i = first_rcomp; i < NStructReal + NumRealComps(); ++i )
123 write_real_comp.push_back(1);
124 real_comp_names.push_back(getDefaultCompNameReal<ParticleType>(i));
127 Vector<int> write_int_comp;
128 Vector<std::string> int_comp_names;
129 for (
int i = 0; i < NStructInt + NumIntComps(); ++i )
131 write_int_comp.push_back(1);
132 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
136 real_comp_names, int_comp_names,
144template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
145 template<
class>
class Allocator,
class CellAssignor>
147ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
148::WritePlotFileHDF5 (
const std::string& dir,
const std::string& name,
149 const Vector<std::string>& real_comp_names,
150 const Vector<std::string>& int_comp_names,
151 const std::string& compression)
const
154 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
155 AMREX_ALWAYS_ASSERT(std::ssize(real_comp_names) == NStructReal + NumRealComps() - first_rcomp);
158 Vector<int> write_real_comp;
159 for (
int i = first_rcomp; i < NStructReal + NumRealComps(); ++i) { write_real_comp.push_back(1); }
161 Vector<int> write_int_comp;
162 for (
int i = 0; i < NStructInt + NumIntComps(); ++i) { write_int_comp.push_back(1); }
165 write_real_comp, write_int_comp,
166 real_comp_names, int_comp_names,
174template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
175 template<
class>
class Allocator,
class CellAssignor>
177ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
178::WritePlotFileHDF5 (
const std::string& dir,
const std::string& name,
179 const Vector<std::string>& real_comp_names,
180 const std::string& compression)
const
183 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
184 AMREX_ALWAYS_ASSERT(std::ssize(real_comp_names) == NStructReal + NumRealComps() - first_rcomp);
186 Vector<int> write_real_comp;
187 for (
int i = first_rcomp; i < NStructReal + NumRealComps(); ++i) write_real_comp.push_back(1);
189 Vector<int> write_int_comp;
190 for (
int i = 0; i < NStructInt + NumIntComps(); ++i) write_int_comp.push_back(1);
192 Vector<std::string> int_comp_names;
193 for (
int i = 0; i < NStructInt + NumIntComps(); ++i )
195 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
199 write_real_comp, write_int_comp,
200 real_comp_names, int_comp_names,
208template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
209 template<
class>
class Allocator,
class CellAssignor>
211ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
212::WritePlotFileHDF5 (
const std::string& dir,
213 const std::string& name,
214 const Vector<int>& write_real_comp,
215 const Vector<int>& write_int_comp,
216 const std::string& compression)
const
219 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
220 AMREX_ALWAYS_ASSERT(std::ssize(write_real_comp) == NStructReal + NumRealComps() - first_rcomp);
223 Vector<std::string> real_comp_names;
224 for (
int i = first_rcomp; i < NStructReal + NumRealComps(); ++i )
226 real_comp_names.push_back(getDefaultCompNameReal<ParticleType>(i));
229 Vector<std::string> int_comp_names;
230 for (
int i = 0; i < NStructInt + NumIntComps(); ++i )
232 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
236 real_comp_names, int_comp_names,
244template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
245 template<
class>
class Allocator,
class CellAssignor>
247ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>::
248WritePlotFileHDF5 (
const std::string& dir,
const std::string& name,
249 const Vector<int>& write_real_comp,
250 const Vector<int>& write_int_comp,
251 const Vector<std::string>& real_comp_names,
252 const Vector<std::string>& int_comp_names,
253 const std::string& compression)
const
255 BL_PROFILE(
"ParticleContainer::WritePlotFile()");
258 write_real_comp, write_int_comp,
259 real_comp_names, int_comp_names,
267template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
268 template<
class>
class Allocator,
class CellAssignor>
270requires (!std::same_as<F, Vector<std::string>>)
272ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
274 const std::string& compression,
F&& f)
const
276 Vector<int> write_real_comp;
277 Vector<std::string> real_comp_names;
279 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
280 for (
int i = first_rcomp; i < NStructReal + NumRealComps(); ++i )
282 write_real_comp.push_back(1);
283 real_comp_names.push_back(getDefaultCompNameReal<ParticleType>(i));
286 Vector<int> write_int_comp;
287 Vector<std::string> int_comp_names;
288 for (
int i = 0; i < NStructInt + NumIntComps(); ++i )
290 write_int_comp.push_back(1);
291 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
295 real_comp_names, int_comp_names, compression,
299template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
300 template<
class>
class Allocator,
class CellAssignor>
303ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
304::WritePlotFileHDF5 (
const std::string& dir,
const std::string& name,
305 const Vector<std::string>& real_comp_names,
306 const Vector<std::string>& int_comp_names,
307 const std::string& compression,
F&& f)
const
310 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
311 AMREX_ALWAYS_ASSERT(std::ssize(real_comp_names) == NStructReal + NumRealComps() - first_rcomp);
314 Vector<int> write_real_comp;
315 for (
int i = first_rcomp; i < NStructReal + NumRealComps(); ++i) { write_real_comp.push_back(1); }
317 Vector<int> write_int_comp;
318 for (
int i = 0; i < NStructInt + NumIntComps(); ++i) { write_int_comp.push_back(1); }
321 write_real_comp, write_int_comp,
322 real_comp_names, int_comp_names,
323 compression, std::forward<F>(f));
326template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
327 template<
class>
class Allocator,
class CellAssignor>
329requires (!std::same_as<F, Vector<std::string>>)
331ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
333 const Vector<std::string>& real_comp_names,
334 const std::string& compression,
F&& f)
const
337 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
338 AMREX_ALWAYS_ASSERT(std::ssize(real_comp_names) == NStructReal + NumRealComps() - first_rcomp);
340 Vector<int> write_real_comp;
341 for (
int i = first_rcomp; i < NStructReal + NumRealComps(); ++i) { write_real_comp.push_back(1); }
343 Vector<int> write_int_comp;
344 for (
int i = 0; i < NStructInt + NumIntComps(); ++i) { write_int_comp.push_back(1); }
346 Vector<std::string> int_comp_names;
347 for (
int i = 0; i < NStructInt + NumIntComps(); ++i )
349 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
353 write_real_comp, write_int_comp,
354 real_comp_names, int_comp_names,
355 compression, std::forward<F>(f));
358template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
359 template<
class>
class Allocator,
class CellAssignor>
362ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
363::WritePlotFileHDF5 (
const std::string& dir,
364 const std::string& name,
365 const Vector<int>& write_real_comp,
366 const Vector<int>& write_int_comp,
367 const std::string& compression,
F&& f)
const
370 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
371 AMREX_ALWAYS_ASSERT(std::ssize(write_real_comp) == NStructReal + NumRealComps() - first_rcomp);
374 Vector<std::string> real_comp_names;
375 for (
int i = first_rcomp; i < NStructReal + NumRealComps(); ++i )
377 real_comp_names.push_back(getDefaultCompNameReal<ParticleType>(i));
380 Vector<std::string> int_comp_names;
381 for (
int i = 0; i < NStructInt + NumIntComps(); ++i )
383 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
387 real_comp_names, int_comp_names,
388 compression, std::forward<F>(f));
391template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
392 template<
class>
class Allocator,
class CellAssignor>
395ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>::
396WritePlotFileHDF5 (
const std::string& dir,
const std::string& name,
397 const Vector<int>& write_real_comp,
398 const Vector<int>& write_int_comp,
399 const Vector<std::string>& real_comp_names,
400 const Vector<std::string>& int_comp_names,
401 const std::string& compression,
F&& f)
const
403 BL_PROFILE(
"ParticleContainer::WritePlotFile()");
406 write_real_comp, write_int_comp,
407 real_comp_names, int_comp_names,
408 compression, std::forward<F>(f));
411template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
412 template<
class>
class Allocator,
class CellAssignor>
415ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
416::WriteHDF5ParticleData (
const std::string& dir,
const std::string& name,
417 const Vector<int>& write_real_comp,
418 const Vector<int>& write_int_comp,
419 const Vector<std::string>& real_comp_names,
420 const Vector<std::string>& int_comp_names,
421 const std::string& compression,
422 F&& f,
bool is_checkpoint)
const
432 write_real_comp, write_int_comp,
433 real_comp_names, int_comp_names,
435 std::forward<F>(f), is_checkpoint);
439template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
440 template<
class>
class Allocator,
class CellAssignor>
442ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
443::CheckpointPreHDF5 ()
449 BL_PROFILE(
"ParticleContainer::CheckpointPre()");
453 Long maxnextid = ParticleType::NextID();
455 for (
int lev = 0; lev < std::ssize(m_particles); lev++) {
456 const auto& pmap = m_particles[lev];
457 for (
const auto& kv : pmap) {
458 const auto& aos = kv.second.GetArrayOfStructs();
459 for (
int k = 0; k < aos.numParticles(); ++k) {
460 const ParticleType& p = aos[k];
461 if (p.id().is_valid()) {
472 ParticleType::NextID(maxnextid);
475 nparticlesPrePost = nparticles;
476 maxnextidPrePost = maxnextid;
478 nParticlesAtLevelPrePost.clear();
479 nParticlesAtLevelPrePost.resize(finestLevel() + 1, 0);
480 for(
int lev(0); lev <= finestLevel(); ++lev) {
481 nParticlesAtLevelPrePost[lev] = NumberOfParticlesAtLevel(lev);
484 whichPrePost.clear();
485 whichPrePost.resize(finestLevel() + 1);
486 countPrePost.clear();
487 countPrePost.resize(finestLevel() + 1);
488 wherePrePost.clear();
489 wherePrePost.resize(finestLevel() + 1);
491 filePrefixPrePost.clear();
492 filePrefixPrePost.resize(finestLevel() + 1);
496template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
497 template<
class>
class Allocator,
class CellAssignor>
499ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
500::CheckpointPostHDF5 ()
506 BL_PROFILE(
"ParticleContainer::CheckpointPostHDF5()");
509 std::ofstream HdrFile;
510 HdrFile.open(HdrFileNamePrePost.c_str(), std::ios::out | std::ios::app);
512 for(
int lev(0); lev <= finestLevel(); ++lev) {
519 for(
int j(0); j < whichPrePost[lev].size(); ++j) {
520 HdrFile << whichPrePost[lev][j] <<
' ' << countPrePost[lev][j] <<
' ' << wherePrePost[lev][j] <<
'\n';
523 const bool gotsome = (nParticlesAtLevelPrePost[lev] > 0);
524 if(gotsome && doUnlink) {
527 Vector<Long> cnt(nOutFilesPrePost,0);
529 for(
int i(0), N = countPrePost[lev].size(); i < N; ++i) {
530 cnt[whichPrePost[lev][i]] += countPrePost[lev][i];
533 for(
int i(0), N = cnt.size(); i < N; ++i) {
546 if( ! HdrFile.good()) {
547 amrex::Abort(
"ParticleContainer::CheckpointPostHDF5(): problem writing HdrFile");
552template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
553 template<
class>
class Allocator,
class CellAssignor>
555ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
556::WritePlotFilePreHDF5 ()
562template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
563 template<
class>
class Allocator,
class CellAssignor>
565ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
566::WritePlotFilePostHDF5 ()
572template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
573 template<
class>
class Allocator,
class CellAssignor>
575ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
576::WriteParticlesHDF5 (
int lev, hid_t grp,
577 Vector<int>& , Vector<int>& count, Vector<Long>& ,
578 const Vector<int>& write_real_comp,
579 const Vector<int>& write_int_comp,
580 const std::string& compression,
581 const Vector<std::map<std::pair<int, int>, IntVector>>& particle_io_flags,
582 bool is_checkpoint)
const
584 BL_PROFILE(
"ParticleContainer::WriteParticlesHDF5()");
587 std::map<int, Vector<int> > tile_map;
590 hid_t dxpl_col, dxpl_ind, dcpl_int, dcpl_real;
591 dxpl_col = H5Pcreate(H5P_DATASET_XFER);
592 dxpl_ind = H5Pcreate(H5P_DATASET_XFER);
594 H5Pset_dxpl_mpio(dxpl_col, H5FD_MPIO_COLLECTIVE);
596 dcpl_int = H5Pcreate(H5P_DATASET_CREATE);
597 dcpl_real = H5Pcreate(H5P_DATASET_CREATE);
599 std::string mode_env, value_env;
600 double comp_value = -1.0;
601 std::string::size_type pos = compression.find(
'@');
602 if (pos != std::string::npos) {
603 mode_env = compression.substr(0, pos);
604 value_env = compression.substr(pos+1);
605 if (!value_env.empty()) {
606 comp_value = atof(value_env.c_str());
610 H5Pset_alloc_time(dcpl_int, H5D_ALLOC_TIME_LATE);
611 H5Pset_alloc_time(dcpl_real, H5D_ALLOC_TIME_LATE);
613 if (!mode_env.empty() && mode_env !=
"None") {
614 const char *chunk_env = NULL;
615 hsize_t chunk_dim = 1024;
616 chunk_env = getenv(
"HDF5_CHUNK_SIZE");
617 if (chunk_env != NULL) {
618 chunk_dim = atoi(chunk_env);
621 H5Pset_chunk(dcpl_int, 1, &chunk_dim);
622 H5Pset_chunk(dcpl_real, 1, &chunk_dim);
624#ifdef AMREX_USE_HDF5_ZFP
625 pos = mode_env.find(
"ZFP");
626 if (pos != std::string::npos) {
627 ret = H5Z_zfp_initialize();
628 if (ret < 0) {
amrex::Abort(
"ZFP initialize failed!"); }
632 if (mode_env ==
"ZLIB") {
633 H5Pset_shuffle(dcpl_int);
634 H5Pset_shuffle(dcpl_real);
635 H5Pset_deflate(dcpl_int, (
int)comp_value);
636 H5Pset_deflate(dcpl_real, (
int)comp_value);
638#ifdef AMREX_USE_HDF5_SZ
639 else if (mode_env ==
"SZ") {
640 ret = H5Z_SZ_Init((
char*)value_env.c_str());
642 std::cout <<
"SZ config file:" << value_env.c_str() << std::endl;
643 amrex::Abort(
"SZ initialize failed, check SZ config file!");
647#ifdef AMREX_USE_HDF5_ZFP
648 else if (mode_env ==
"ZFP_RATE") {
649 H5Pset_zfp_rate(dcpl_int, comp_value);
650 H5Pset_zfp_rate(dcpl_real, comp_value);
652 else if (mode_env ==
"ZFP_PRECISION") {
653 H5Pset_zfp_precision(dcpl_int, (
unsigned int)comp_value);
654 H5Pset_zfp_precision(dcpl_real, (
unsigned int)comp_value);
656 else if (mode_env ==
"ZFP_ACCURACY") {
657 H5Pset_zfp_accuracy(dcpl_int, comp_value);
658 H5Pset_zfp_accuracy(dcpl_real, comp_value);
660 else if (mode_env ==
"ZFP_REVERSIBLE") {
661 H5Pset_zfp_reversible(dcpl_int);
662 H5Pset_zfp_reversible(dcpl_real);
666 std::cout <<
"\nHDF5 particle checkpoint using " << mode_env <<
", " <<
667 value_env <<
", " << chunk_dim << std::endl;
671 for (
const auto& kv : m_particles[lev])
673 const int grid = kv.first.first;
674 const int tile = kv.first.second;
675 tile_map[grid].push_back(tile);
676 const auto& pflags = particle_io_flags[lev].at(kv.first);
679 count[grid] += particle_detail::countFlags(pflags);
683 info.SetAlloc(
false);
684 MultiFab state(ParticleBoxArray(lev), ParticleDistributionMap(lev), 1,0,info);
687 ULong my_mfi_int_total_size = 0, my_mfi_real_total_size = 0, int_size, real_size;
689 Vector<int> my_mfi_real_size;
690 Vector<int> my_mfi_int_size;
691 Vector<int> my_nparticles;
694 hid_t offset_id, offset_space, real_mem_space, real_dset_space, real_dset_id;
695 hid_t int_mem_space, int_dset_id, int_dset_space;
696 hsize_t total_mfi = 0, total_real_size = 0, total_int_size = 0, real_file_offset = 0, int_file_offset = 0;
697 hsize_t my_int_offset, my_int_count, my_real_offset, my_real_count;
701 int real_comp_count = AMREX_SPACEDIM;
703 for (
int i = 0; i < std::ssize(write_real_comp); ++i ) {
704 if (write_real_comp[i]) { ++real_comp_count; }
707 int int_comp_count = 2;
709 for (
int i = 0; i < NStructInt + NumIntComps(); ++i ) {
710 if (write_int_comp[i]) { ++int_comp_count; }
714 for (MFIter mfi(state); mfi.isValid(); ++mfi) {
715 const int grid = mfi.index();
716 if (count[grid] == 0)
719 int_size = count[grid] * int_comp_count;
720 my_mfi_int_size.push_back(int_size);
721 my_nparticles.push_back(count[grid]);
722 my_mfi_int_total_size += int_size;
725 real_size = count[grid] * real_comp_count;
726 my_mfi_real_size.push_back(real_size);
727 my_mfi_real_total_size += real_size;
737 total_mfi += all_mfi_cnt[i];
743 all_mfi_cnt[0] = my_mfi_cnt;
744 all_mfi_int_total_size[0] = my_mfi_int_total_size;
749 int_file_offset += all_mfi_int_total_size[i];
750 my_int_offset = int_file_offset;
754 total_int_size += all_mfi_int_total_size[i];
768 if (H5Pget_layout(dcpl_int) == H5D_CHUNKED) {
769 if (H5Pget_chunk(dcpl_int, 1, &chunk_size) > -1) {
770 if (chunk_size > total_int_size) {
771 H5Pset_chunk(dcpl_int, 1, &total_int_size);
776 int_dset_space = H5Screate_simple(1, &total_int_size, NULL);
777#ifdef AMREX_USE_HDF5_ASYNC
778 int_dset_id = H5Dcreate_async(grp,
"data:datatype=0", H5T_NATIVE_INT, int_dset_space, H5P_DEFAULT, dcpl_int, H5P_DEFAULT, es_par_g);
780 int_dset_id = H5Dcreate(grp,
"data:datatype=0", H5T_NATIVE_INT, int_dset_space, H5P_DEFAULT, dcpl_int, H5P_DEFAULT);
783 H5Sclose(int_dset_space);
790 all_mfi_real_total_size[0] = my_mfi_real_total_size;
794 total_real_size += all_mfi_real_total_size[i];
796#ifdef AMREX_USE_HDF5_SZ
797 if (mode_env ==
"SZ") {
799 unsigned int* cd_values = NULL;
800 unsigned filter_config;
801 if (
sizeof(
typename ParticleType::RealType) == 4) {
802 SZ_metaDataToCdArray(&cd_nelmts, &cd_values, SZ_FLOAT, 0, 0, 0, 0, total_real_size);
803 H5Pset_filter(dcpl_real, H5Z_FILTER_SZ, H5Z_FLAG_MANDATORY, cd_nelmts, cd_values);
806 SZ_metaDataToCdArray(&cd_nelmts, &cd_values, SZ_DOUBLE, 0, 0, 0, 0, total_real_size);
807 H5Pset_filter(dcpl_real, H5Z_FILTER_SZ, H5Z_FLAG_MANDATORY, cd_nelmts, cd_values);
812 if (H5Pget_layout(dcpl_real) == H5D_CHUNKED) {
813 if (H5Pget_chunk(dcpl_real, 1, &chunk_size) > -1) {
814 if (chunk_size > total_real_size) {
815 H5Pset_chunk(dcpl_real, 1, &total_real_size);
820 real_dset_space = H5Screate_simple(1, &total_real_size, NULL);
821 if (
sizeof(
typename ParticleType::RealType) == 4) {
822#ifdef AMREX_USE_HDF5_ASYNC
823 real_dset_id = H5Dcreate_async(grp,
"data:datatype=1", H5T_NATIVE_FLOAT, real_dset_space,
824 H5P_DEFAULT, dcpl_real, H5P_DEFAULT, es_par_g);
826 real_dset_id = H5Dcreate(grp,
"data:datatype=1", H5T_NATIVE_FLOAT, real_dset_space,
827 H5P_DEFAULT, dcpl_real, H5P_DEFAULT);
831#ifdef AMREX_USE_HDF5_ASYNC
832 real_dset_id = H5Dcreate_async(grp,
"data:datatype=1", H5T_NATIVE_DOUBLE, real_dset_space,
833 H5P_DEFAULT, dcpl_real, H5P_DEFAULT, es_par_g);
835 real_dset_id = H5Dcreate(grp,
"data:datatype=1", H5T_NATIVE_DOUBLE, real_dset_space,
836 H5P_DEFAULT, dcpl_real, H5P_DEFAULT);
839 H5Sclose(real_dset_space);
841 real_file_offset = 0;
843 real_file_offset += all_mfi_real_total_size[i];
845 my_real_offset = real_file_offset;
848 int max_mfi_count = 0, write_count = 0;
850 if (max_mfi_count < all_mfi_cnt[i]) {
851 max_mfi_count = all_mfi_cnt[i];
856 for (MFIter mfi(state); mfi.isValid(); ++mfi)
858 const int grid = mfi.index();
860 if (count[grid] == 0) {
continue; }
863 Vector<ParticleReal> rstuff;
864 particle_detail::packIOData(istuff, rstuff, *
this, lev, grid,
865 write_real_comp, write_int_comp,
866 particle_io_flags, tile_map[grid], count[grid],
869 my_int_count = istuff.size();
870 int_mem_space = H5Screate_simple(1, &my_int_count, NULL);
873 int_dset_space = H5Screate_simple(1, &total_int_size, NULL);
874 H5Sselect_hyperslab (int_dset_space, H5S_SELECT_SET, &my_int_offset, NULL, &my_int_count, NULL);
876#ifdef AMREX_USE_HDF5_ASYNC
877 ret = H5Dwrite_async(int_dset_id, H5T_NATIVE_INT, int_mem_space, int_dset_space, dxpl_col, istuff.dataPtr(), es_par_g);
879 ret = H5Dwrite(int_dset_id, H5T_NATIVE_INT, int_mem_space, int_dset_space, dxpl_col, istuff.dataPtr());
881 if (ret < 0) {
amrex::Abort(
"H5Dwrite int_dset failed!"); }
883 H5Sclose(int_dset_space);
884 H5Sclose(int_mem_space);
886 my_int_offset += my_int_count;
888 my_real_count = rstuff.size();
889 real_mem_space = H5Screate_simple(1, &my_real_count, NULL);
892 real_dset_space = H5Screate_simple(1, &total_real_size, NULL);
893 H5Sselect_hyperslab (real_dset_space, H5S_SELECT_SET, &my_real_offset, NULL, &my_real_count, NULL);
894 if (
sizeof(
typename ParticleType::RealType) == 4) {
895#ifdef AMREX_USE_HDF5_ASYNC
896 ret = H5Dwrite_async(real_dset_id, H5T_NATIVE_FLOAT, real_mem_space, real_dset_space, dxpl_col, rstuff.dataPtr(), es_par_g);
898 ret = H5Dwrite(real_dset_id, H5T_NATIVE_FLOAT, real_mem_space, real_dset_space, dxpl_col, rstuff.dataPtr());
901#ifdef AMREX_USE_HDF5_ASYNC
902 ret = H5Dwrite_async(real_dset_id, H5T_NATIVE_DOUBLE, real_mem_space, real_dset_space, dxpl_col, rstuff.dataPtr(), es_par_g);
904 ret = H5Dwrite(real_dset_id, H5T_NATIVE_DOUBLE, real_mem_space, real_dset_space, dxpl_col, rstuff.dataPtr());
908 if (ret < 0) {
amrex::Abort(
"H5Dwrite real_dset failed!"); }
910 H5Sclose(real_mem_space);
911 H5Sclose(real_dset_space);
913 my_real_offset += my_real_count;
919 while (write_count < max_mfi_count) {
920 int_dset_space = H5Screate_simple(1, &total_int_size, NULL);
921 real_dset_space = H5Screate_simple(1, &total_real_size, NULL);
922 H5Sselect_none(int_dset_space);
923 H5Sselect_none(real_dset_space);
925#ifdef AMREX_USE_HDF5_ASYNC
926 H5Dwrite_async(int_dset_id, H5T_NATIVE_INT, int_dset_space, int_dset_space, dxpl_col, NULL, es_par_g);
927 if (
sizeof(
typename ParticleType::RealType) == 4) {
928 H5Dwrite_async(real_dset_id, H5T_NATIVE_FLOAT, real_dset_space, real_dset_space, dxpl_col, NULL, es_par_g);
930 H5Dwrite_async(real_dset_id, H5T_NATIVE_DOUBLE, real_dset_space, real_dset_space, dxpl_col, NULL, es_par_g);
933 H5Dwrite(int_dset_id, H5T_NATIVE_INT, int_dset_space, int_dset_space, dxpl_col, NULL);
934 if (
sizeof(
typename ParticleType::RealType) == 4) {
935 H5Dwrite(real_dset_id, H5T_NATIVE_FLOAT, real_dset_space, real_dset_space, dxpl_col, NULL);
937 H5Dwrite(real_dset_id, H5T_NATIVE_DOUBLE, real_dset_space, real_dset_space, dxpl_col, NULL);
941 H5Sclose(int_dset_space);
942 H5Sclose(real_dset_space);
947#ifdef AMREX_USE_HDF5_ASYNC
948 H5Dclose_async(real_dset_id, es_par_g);
949 H5Dclose_async(int_dset_id, es_par_g);
951 H5Dclose(real_dset_id);
952 H5Dclose(int_dset_id);
956 offset_space = H5Screate_simple(1, &total_mfi, NULL);
957#ifdef AMREX_USE_HDF5_ASYNC
958 offset_id = H5Dcreate_async(grp,
"nparticles_grid", H5T_NATIVE_INT, offset_space, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT, es_par_g);
960 offset_id = H5Dcreate(grp,
"nparticles_grid", H5T_NATIVE_INT, offset_space, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
965 my_int_offset += all_mfi_cnt[i];
967 my_int_count = my_mfi_cnt;
968 int_mem_space = H5Screate_simple(1, &my_int_count, NULL);
971 H5Sselect_hyperslab (offset_space, H5S_SELECT_SET, &my_int_offset, NULL, &my_int_count, NULL);
973#ifdef AMREX_USE_HDF5_ASYNC
974 ret = H5Dwrite_async(offset_id, H5T_NATIVE_INT, int_mem_space, offset_space, dxpl_col, my_nparticles.dataPtr(), es_par_g);
976 ret = H5Dwrite(offset_id, H5T_NATIVE_INT, int_mem_space, offset_space, dxpl_col, my_nparticles.dataPtr());
978 if (ret < 0) {
amrex::Abort(
"H5Dwrite offset failed!"); }
985 H5Sclose(int_mem_space);
986 H5Sclose(offset_space);
988#ifdef AMREX_USE_HDF5_ASYNC
989 H5Dclose_async(offset_id, es_par_g);
998template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
999 template<
class>
class Allocator,
class CellAssignor>
1001ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
1002::RestartHDF5 (
const std::string& dir,
const std::string& file,
bool )
1007template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
1008 template<
class>
class Allocator,
class CellAssignor>
1010ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
1011::RestartHDF5 (
const std::string& dir,
const std::string& file)
1013 BL_PROFILE(
"ParticleContainer::RestartHDF5()");
1019 std::string fullname = dir;
1020 if (!fullname.empty() && fullname[fullname.size()-1] !=
'/') {
1027 hid_t fid, dset, grp, fapl, attr, atype, dspace, int_dset, real_dset;
1030 fapl = H5Pcreate (H5P_FILE_ACCESS);
1037 fid = H5Fopen(fullname.c_str(), H5F_ACC_RDONLY, fapl);
1039 std::string msg(
"ParticleContainer::RestartHDF5(): unable to open file: ");
1044 attr = H5Aopen(fid,
"version_name", H5P_DEFAULT);
1046 std::string msg(
"ParticleContainer::RestartHDF5(): unable to open version attribute ");
1050 atype = H5Aget_type(attr);
1052 std::string msg(
"ParticleContainer::RestartHDF5(): unable to get type of attribute ");
1056 size_t attr_len = H5Tget_size(atype);
1057 if (attr_len == 0) {
1058 std::string msg(
"ParticleContainer::RestartHDF5(): unable to get size of attribute ");
1062 std::string version;
1063 version.resize(attr_len+1);
1065 ret = H5Aread(attr, atype, &version[0]);
1080 bool convert_ids =
false;
1081 if (version.find(
"Version_Two_Dot_One") != std::string::npos) {
1084 if (version.find(
"_single") != std::string::npos) {
1087 else if (version.find(
"_double") != std::string::npos) {
1091 std::string msg(
"ParticleContainer::Restart(): bad version string: ");
1096 std::string gname =
"Chombo_global";
1097 std::string aname =
"SpaceDim";
1099 grp = H5Gopen(fid, gname.c_str(), H5P_DEFAULT);
1101 std::string msg(
"ParticleContainer::RestartHDF5(): unable to open group ");
1105 ret =
ReadHDF5Attr(grp, aname.c_str(), &dm, H5T_NATIVE_INT);
1107 std::string msg(
"ParticleContainer::RestartHDF5(): unable to read attribute ");
1113 aname =
"num_component_real";
1115 ret =
ReadHDF5Attr(fid, aname.c_str(), &nr, H5T_NATIVE_INT);
1117 std::string msg(
"ParticleContainer::RestartHDF5(): unable to read attribute ");
1121 int nrc = ParticleType::is_soa_particle ? NStructReal + NumRealComps() - AMREX_SPACEDIM : NStructReal + NumRealComps();
1123 amrex::Abort(
"ParticleContainer::RestartHDF5(): nr not the expected value");
1126 aname =
"num_component_int";
1128 ret =
ReadHDF5Attr(fid, aname.c_str(), &ni, H5T_NATIVE_INT);
1130 std::string msg(
"ParticleContainer::RestartHDF5(): unable to read attribute ");
1134 if (ni != NStructInt + NumIntComps()) {
1135 amrex::Abort(
"ParticleContainer::Restart(): ni != NStructInt");
1138 aname =
"nparticles";
1140 ret =
ReadHDF5Attr(fid, aname.c_str(), &nparticles, H5T_NATIVE_LLONG);
1142 std::string msg(
"ParticleContainer::RestartHDF5(): unable to read attribute ");
1148 aname =
"maxnextid";
1150 ret =
ReadHDF5Attr(fid, aname.c_str(), &maxnextid, H5T_NATIVE_LLONG);
1152 std::string msg(
"ParticleContainer::RestartHDF5(): unable to read attribute ");
1157 ParticleType::NextID(maxnextid);
1159 aname =
"finest_level";
1160 int finest_level_in_file;
1161 ret =
ReadHDF5Attr(fid, aname.c_str(), &finest_level_in_file, H5T_NATIVE_INT);
1163 std::string msg(
"ParticleContainer::RestartHDF5(): unable to read attribute ");
1169 hid_t comp_dtype = H5Tcreate (H5T_COMPOUND, 2 * AMREX_SPACEDIM *
sizeof(
int));
1170 if (1 == AMREX_SPACEDIM) {
1171 H5Tinsert (comp_dtype,
"lo_i", 0 *
sizeof(
int), H5T_NATIVE_INT);
1172 H5Tinsert (comp_dtype,
"hi_i", 1 *
sizeof(
int), H5T_NATIVE_INT);
1174 else if (2 == AMREX_SPACEDIM) {
1175 H5Tinsert (comp_dtype,
"lo_i", 0 *
sizeof(
int), H5T_NATIVE_INT);
1176 H5Tinsert (comp_dtype,
"lo_j", 1 *
sizeof(
int), H5T_NATIVE_INT);
1177 H5Tinsert (comp_dtype,
"hi_i", 2 *
sizeof(
int), H5T_NATIVE_INT);
1178 H5Tinsert (comp_dtype,
"hi_j", 3 *
sizeof(
int), H5T_NATIVE_INT);
1180 else if (3 == AMREX_SPACEDIM) {
1181 H5Tinsert (comp_dtype,
"lo_i", 0 *
sizeof(
int), H5T_NATIVE_INT);
1182 H5Tinsert (comp_dtype,
"lo_j", 1 *
sizeof(
int), H5T_NATIVE_INT);
1183 H5Tinsert (comp_dtype,
"lo_k", 2 *
sizeof(
int), H5T_NATIVE_INT);
1184 H5Tinsert (comp_dtype,
"hi_i", 3 *
sizeof(
int), H5T_NATIVE_INT);
1185 H5Tinsert (comp_dtype,
"hi_j", 4 *
sizeof(
int), H5T_NATIVE_INT);
1186 H5Tinsert (comp_dtype,
"hi_k", 5 *
sizeof(
int), H5T_NATIVE_INT);
1190 Vector<BoxArray> particle_box_arrays(finest_level_in_file + 1);
1191 bool dual_grid =
false;
1193 for (
int lev = 0; lev <= finest_level_in_file; lev++)
1195 if (lev > finestLevel())
1201 gname =
"level_" + std::to_string(lev);
1202 grp = H5Gopen(fid, gname.c_str(), H5P_DEFAULT);
1204 std::string msg(
"ParticleContainer::RestartHDF5(): unable to open group ");
1209 dset = H5Dopen(grp,
"boxes", H5P_DEFAULT);
1210 dspace = H5Dget_space(dset);
1212 H5Sget_simple_extent_dims(dspace, &ngrid, NULL);
1214 Vector<int> boxes(ngrid*AMREX_SPACEDIM*2);
1215 ret = H5Dread(dset, comp_dtype, H5S_ALL, H5S_ALL, H5P_DEFAULT, &(boxes[0]));
1217 std::string msg(
"ParticleContainer::RestartHDF5(): unable to read nparticles_grid dataset");
1231 for (
int i = 0; i < (
int)ngrid; i++) {
1232 const int base = i * 2 * AMREX_SPACEDIM;
1237 boxes[base+AMREX_SPACEDIM+1],
1238 boxes[base+AMREX_SPACEDIM+2])});
1241 particle_box_arrays[lev] = BoxArray(std::move(bl));
1244 if (! particle_box_arrays[lev].CellEqual(ParticleBoxArray(lev)))
1249 for (
int lev = 0; lev <= finestLevel(); lev++) {
1250 SetParticleBoxArray(lev, particle_box_arrays[lev]);
1251 DistributionMapping pdm(particle_box_arrays[lev]);
1252 SetParticleDistributionMap(lev, pdm);
1256 Vector<int> ngrids(finest_level_in_file+1);
1257 for (
int lev = 0; lev <= finest_level_in_file; lev++) {
1258 gname =
"level_" + std::to_string(lev);
1259 grp = H5Gopen(fid, gname.c_str(), H5P_DEFAULT);
1261 std::string msg(
"ParticleContainer::RestartHDF5(): unable to open group ");
1267 ret =
ReadHDF5Attr(grp, aname.c_str(), &ngrids[lev], H5T_NATIVE_INT);
1269 std::string msg(
"ParticleContainer::RestartHDF5(): unable to read attribute ");
1275 if (lev <= finestLevel()) {
1276 AMREX_ASSERT(ngrids[lev] ==
int(ParticleBoxArray(lev).size()));
1284 if (finest_level_in_file > finestLevel()) {
1285 m_particles.resize(finest_level_in_file+1);
1288 for (
int lev = 0; lev <= finest_level_in_file; lev++) {
1290 gname =
"level_" + std::to_string(lev);
1291 grp = H5Gopen(fid, gname.c_str(), H5P_DEFAULT);
1293 std::string msg(
"ParticleContainer::RestartHDF5(): unable to open group ");
1298 dset = H5Dopen(grp,
"nparticles_grid", H5P_DEFAULT);
1299 Vector<int> count(ngrids[lev]);
1300 ret = H5Dread(dset, H5T_NATIVE_INT, H5S_ALL, H5S_ALL, H5P_DEFAULT, &(count[0]));
1302 std::string msg(
"ParticleContainer::RestartHDF5(): unable to read nparticles_grid dataset");
1307 Vector<hsize_t>
offset(ngrids[lev]);
1309 for (
int i = 1; i < ngrids[lev]; i++) {
1313 int_dset = H5Dopen(grp,
"data:datatype=0", H5P_DEFAULT);
1315 std::string msg(
"ParticleContainer::RestartHDF5(): unable to open int dataset");
1318 real_dset = H5Dopen(grp,
"data:datatype=1", H5P_DEFAULT);
1319 if (real_dset < 0) {
1320 std::string msg(
"ParticleContainer::RestartHDF5(): unable to open int dataset");
1324 Vector<int> grids_to_read;
1325 if (lev <= finestLevel()) {
1326 for (MFIter mfi(*m_dummy_mf[lev]); mfi.isValid(); ++mfi) {
1327 grids_to_read.push_back(mfi.index());
1335 const int NReaders = MaxReaders();
1336 if (rank >= NReaders) {
1338 H5Dclose(real_dset);
1343 const int Navg = ngrids[lev] / NReaders;
1344 const int Nleft = ngrids[lev] - Navg * NReaders;
1348 lo = rank*(Navg + 1);
1352 lo = rank * Navg + Nleft;
1356 for (
int i = lo; i < hi; ++i) {
1357 grids_to_read.push_back(i);
1361 for(
int igrid = 0; igrid < std::ssize(grids_to_read); ++igrid) {
1362 const int grid = grids_to_read[igrid];
1364 if (how ==
"single") {
1365 ReadParticlesHDF5<float>(
offset[grid], count[grid], grid, lev, int_dset, real_dset, finest_level_in_file, convert_ids);
1367 else if (how ==
"double") {
1368 ReadParticlesHDF5<double>(
offset[grid], count[grid], grid, lev, int_dset, real_dset, finest_level_in_file, convert_ids);
1371 std::string msg(
"ParticleContainer::Restart(): bad parameter: ");
1378 H5Dclose(real_dset);
1382 H5Tclose(comp_dtype);
1390 if (m_verbose > 1) {
1393 amrex::Print() <<
"ParticleContainer::Restart() time: " << stoptime <<
'\n';
1398template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
1399 template<
class>
class Allocator,
class CellAssignor>
1400template <
class RTYPE>
1402ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
1403::ReadParticlesHDF5 (hsize_t
offset, hsize_t cnt,
int grd,
int lev,
1404 hid_t int_dset, hid_t real_dset,
int finest_level_in_file,
1407 BL_PROFILE(
"ParticleContainer::ReadParticlesHDF5()");
1411 hid_t int_dspace, int_fspace, real_dspace, real_fspace;
1416 const int iChunkSize = 2 + NStructInt + NumIntComps();
1417 Vector<int> istuff(cnt*iChunkSize);
1418 int_fspace = H5Dget_space(int_dset);
1419 hsize_t int_cnt = cnt*iChunkSize;
1420 hsize_t int_offset =
offset*iChunkSize;
1421 int_dspace = H5Screate_simple(1, &int_cnt, NULL);
1422 H5Sselect_hyperslab (int_fspace, H5S_SELECT_SET, &int_offset, NULL, &int_cnt, NULL);
1423 H5Dread(int_dset, H5T_NATIVE_INT, int_dspace, int_fspace, H5P_DEFAULT, istuff.dataPtr());
1425 H5Sclose(int_fspace);
1426 H5Sclose(int_dspace);
1431 const int rChunkSize = ParticleType::is_soa_particle ?
1432 NStructReal + NumRealComps() : AMREX_SPACEDIM + NStructReal + NumRealComps();
1433 Vector<RTYPE> rstuff(cnt*rChunkSize);
1434 real_fspace = H5Dget_space(real_dset);
1435 hsize_t real_cnt = cnt*rChunkSize;
1436 hsize_t real_offset =
offset*rChunkSize;
1437 real_dspace = H5Screate_simple(1, &real_cnt, NULL);
1438 H5Sselect_hyperslab (real_fspace, H5S_SELECT_SET, &real_offset, NULL, &real_cnt, NULL);
1439 if (
sizeof(RTYPE) == 4) {
1440 H5Dread(real_dset, H5T_NATIVE_FLOAT, real_dspace, real_fspace, H5P_DEFAULT, rstuff.dataPtr());
1442 H5Dread(real_dset, H5T_NATIVE_DOUBLE, real_dspace, real_fspace, H5P_DEFAULT, rstuff.dataPtr());
1445 H5Sclose(real_fspace);
1446 H5Sclose(real_dspace);
1449 int* iptr = istuff.dataPtr();
1450 RTYPE* rptr = rstuff.dataPtr();
1454 Particle<NStructReal, NStructInt> ptemp;
1455 ParticleLocData pld;
1457 Vector<std::map<std::pair<int, int>, Gpu::HostVector<ParticleType> > > host_particles;
1458 host_particles.reserve(15);
1459 host_particles.resize(finest_level_in_file+1);
1461 Vector<std::map<std::pair<int, int>,
1462 std::vector<Gpu::HostVector<ParticleReal> > > > host_real_attribs;
1463 host_real_attribs.reserve(15);
1464 host_real_attribs.resize(finest_level_in_file+1);
1466 Vector<std::map<std::pair<int, int>,
1467 std::vector<Gpu::HostVector<int> > > > host_int_attribs;
1468 host_int_attribs.reserve(15);
1469 host_int_attribs.resize(finest_level_in_file+1);
1471 Vector<std::map<std::pair<int, int>, Gpu::HostVector<uint64_t> > > host_idcpu;
1472 host_idcpu.reserve(15);
1473 host_idcpu.resize(finest_level_in_file+1);
1475 for (hsize_t i = 0; i < cnt; i++) {
1479 std::int32_t xi, yi;
1480 std::uint32_t xu, yu;
1483 std::memcpy(&xu, &xi,
sizeof(xi));
1484 std::memcpy(&yu, &yi,
sizeof(yi));
1485 ptemp.m_idcpu = ((std::uint64_t)xu) << 32 | yu;
1487 ptemp.id() = iptr[0];
1488 ptemp.cpu() = iptr[1];
1492 for (
int j = 0; j < NStructInt; j++)
1494 ptemp.idata(j) = *iptr;
1504 rptr += AMREX_SPACEDIM;
1506 for (
int j = 0; j < NStructReal; j++)
1512 locateParticle(ptemp, pld, 0, finestLevel(), 0);
1514 std::pair<int, int> ind(grd, pld.m_tile);
1516 host_real_attribs[lev][ind].resize(NumRealComps());
1517 host_int_attribs[lev][ind].resize(NumIntComps());
1519 if constexpr (!ParticleType::is_soa_particle)
1522 host_particles[lev][ind].push_back(ptemp);
1525 for (
int icomp = 0; icomp < NumRealComps(); icomp++) {
1526 host_real_attribs[lev][ind][icomp].push_back(
ParticleReal(*rptr));
1531 for (
int icomp = 0; icomp < NumIntComps(); icomp++) {
1532 host_int_attribs[lev][ind][icomp].push_back(*iptr);
1537 host_particles[lev][ind];
1539 for (
int j = 0; j < AMREX_SPACEDIM; j++) {
1540 host_real_attribs[lev][ind][j].push_back(ptemp.pos(j));
1543 host_idcpu[lev][ind].push_back(ptemp.m_idcpu);
1546 for (
int icomp = AMREX_SPACEDIM; icomp < NumRealComps(); icomp++) {
1547 host_real_attribs[lev][ind][icomp].push_back(
ParticleReal(*rptr));
1552 for (
int icomp = 0; icomp < NumIntComps(); icomp++) {
1553 host_int_attribs[lev][ind][icomp].push_back(*iptr);
1559 for (
int host_lev = 0; host_lev < std::ssize(host_particles); ++host_lev)
1561 for (
auto& kv : host_particles[host_lev]) {
1562 auto grid = kv.first.first;
1563 auto tile = kv.first.second;
1564 const auto& src_tile = kv.second;
1566 auto& dst_tile = DefineAndReturnParticleTile(host_lev, grid, tile);
1567 auto old_size = dst_tile.size();
1568 auto new_size = old_size;
1569 if constexpr (!ParticleType::is_soa_particle)
1571 new_size += src_tile.size();
1574 new_size += host_real_attribs[host_lev][std::make_pair(grid,tile)][0].size();
1576 dst_tile.resize(new_size);
1578 if constexpr (!ParticleType::is_soa_particle)
1581 dst_tile.GetArrayOfStructs().begin() + old_size);
1584 host_idcpu[host_lev][std::make_pair(grid,tile)].
begin(),
1585 host_idcpu[host_lev][std::make_pair(grid,tile)].
end(),
1586 dst_tile.GetStructOfArrays().GetIdCPUData().begin() + old_size);
1589 for (
int i = 0; i < NumRealComps(); ++i) {
1591 host_real_attribs[host_lev][std::make_pair(grid,tile)][i].
begin(),
1592 host_real_attribs[host_lev][std::make_pair(grid,tile)][i].
end(),
1593 dst_tile.GetStructOfArrays().GetRealData(i).begin() + old_size);
1596 for (
int i = 0; i < NumIntComps(); ++i) {
1598 host_int_attribs[host_lev][std::make_pair(grid,tile)][i].
begin(),
1599 host_int_attribs[host_lev][std::make_pair(grid,tile)][i].
end(),
1600 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 WritePlotFileHDF5(const std::string &dir, const std::string &name, const std::string &compression) const
This version of WritePlotFile writes all components and assigns component names.
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)
Remove a file, symbolic link, or empty directory.
Definition AMReX_FileSystem.cpp:41
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:50
__host__ __device__ void ignore_unused(const Ts &...)
This shuts up the compiler about unused variables.
Definition AMReX.H:139
__host__ __device__ Dim3 begin(BoxND< dim > const &box) noexcept
Definition AMReX_Box.H:2018
BoxND< 3 > Box
Box is an alias for amrex::BoxND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:30
double second() noexcept
Definition AMReX_Utility.cpp:940
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:185
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:235
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:241
const int[]
Definition AMReX_BLProfiler.cpp:1664
__host__ __device__ Dim3 end(BoxND< dim > const &box) noexcept
Definition AMReX_Box.H:2028
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:100
static MPI_Datatype type()