1#ifndef AMREX_PARTICLEIO_H
2#define AMREX_PARTICLEIO_H
3#include <AMReX_Config.H>
10template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
11 template<
class>
class Allocator,
class CellAssignor>
13ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
14::WriteParticleRealData (
void* data,
size_t size, std::ostream& os)
const
16 if (
sizeof(
typename ParticleType::RealType) == 4) {
19 else if (
sizeof(
typename ParticleType::RealType) == 8) {
24template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
25 template<
class>
class Allocator,
class CellAssignor>
30 if (
sizeof(
typename ParticleType::RealType) == 4) {
31 readFloatData((
float*) data, size, is, ParticleRealDescriptor);
33 else if (
sizeof(
typename ParticleType::RealType) == 8) {
42 return p.id().is_valid();
47template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
48 template<
class>
class Allocator,
class CellAssignor>
50ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
51::Checkpoint (
const std::string& dir,
52 const std::string& name,
bool ,
59 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
60 for (
int i = first_rcomp; i < NStructReal + NumRealComps(); ++i )
62 write_real_comp.push_back(1);
63 if (real_comp_names.empty())
65 tmp_real_comp_names.push_back(getDefaultCompNameReal<ParticleType>(i));
69 tmp_real_comp_names.push_back(real_comp_names[i-first_rcomp]);
75 for (
int i = 0; i < NStructInt + NumIntComps(); ++i )
77 write_int_comp.push_back(1);
78 if (int_comp_names.empty())
80 tmp_int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
84 tmp_int_comp_names.push_back(int_comp_names[i]);
88 WriteBinaryParticleData(dir, name, write_real_comp, write_int_comp,
89 tmp_real_comp_names, tmp_int_comp_names,
93template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
94 template<
class>
class Allocator,
class CellAssignor>
103 WriteBinaryParticleData(dir, name, write_real_comp, write_int_comp,
104 real_comp_names, int_comp_names,
109template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
110 template<
class>
class Allocator,
class CellAssignor>
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));
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));
133 WriteBinaryParticleData(dir, name, write_real_comp, write_int_comp,
134 real_comp_names, int_comp_names,
139template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
140 template<
class>
class Allocator,
class CellAssignor>
147 if constexpr(ParticleType::is_soa_particle) {
155 int nrc = ParticleType::is_soa_particle ? NStructReal + NumRealComps() - AMREX_SPACEDIM : NStructReal + NumRealComps();
156 for (
int i = 0; i < nrc; ++i) {
157 write_real_comp.push_back(1);
161 for (
int i = 0; i < NStructInt + NumIntComps(); ++i) {
162 write_int_comp.push_back(1);
165 WriteBinaryParticleData(dir, name,
166 write_real_comp, write_int_comp,
167 real_comp_names, int_comp_names,
172template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
173 template<
class>
class Allocator,
class CellAssignor>
179 if constexpr(ParticleType::is_soa_particle) {
186 int nrc = ParticleType::is_soa_particle ? NStructReal + NumRealComps() - AMREX_SPACEDIM : NStructReal + NumRealComps();
187 for (
int i = 0; i < nrc; ++i) {
188 write_real_comp.push_back(1);
192 for (
int i = 0; i < NStructInt + NumIntComps(); ++i) {
193 write_int_comp.push_back(1);
197 for (
int i = 0; i < NStructInt + NumIntComps(); ++i )
199 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
202 WriteBinaryParticleData(dir, name,
203 write_real_comp, write_int_comp,
204 real_comp_names, int_comp_names,
209template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
210 template<
class>
class Allocator,
class CellAssignor>
214 const std::string& name,
219 if constexpr(ParticleType::is_soa_particle) {
227 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
228 for (
int i = first_rcomp; i < NStructReal + NumRealComps(); ++i )
230 real_comp_names.push_back(getDefaultCompNameReal<ParticleType>(i));
234 for (
int i = 0; i < NStructInt + NumIntComps(); ++i )
236 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
239 WriteBinaryParticleData(dir, name, write_real_comp, write_int_comp,
240 real_comp_names, int_comp_names,
245template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
246 template<
class>
class Allocator,
class CellAssignor>
249WritePlotFile (
const std::string& dir,
const std::string& name,
255 BL_PROFILE(
"ParticleContainer::WritePlotFile()");
257 WriteBinaryParticleData(dir, name,
258 write_real_comp, write_int_comp,
259 real_comp_names, int_comp_names,
264template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
265 template<
class>
class Allocator,
class CellAssignor>
274 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
275 for (
int i = first_rcomp; i < NStructReal + NumRealComps(); ++i )
277 write_real_comp.push_back(1);
278 real_comp_names.push_back(getDefaultCompNameReal<ParticleType>(i));
283 for (
int i = 0; i < NStructInt + NumIntComps(); ++i )
285 write_int_comp.push_back(1);
286 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
289 WriteBinaryParticleData(dir, name, write_real_comp, write_int_comp,
290 real_comp_names, int_comp_names, f);
294template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
295 template<
class>
class Allocator,
class CellAssignor>
303 if constexpr(ParticleType::is_soa_particle) {
311 int nrc = ParticleType::is_soa_particle ? NStructReal + NumRealComps() - AMREX_SPACEDIM : NStructReal + NumRealComps();
312 for (
int i = 0; i < nrc; ++i) {
313 write_real_comp.push_back(1);
317 for (
int i = 0; i < NStructInt + NumIntComps(); ++i) {
318 write_int_comp.push_back(1);
321 WriteBinaryParticleData(dir, name,
322 write_real_comp, write_int_comp,
323 real_comp_names, int_comp_names, f);
327template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
328 template<
class>
class Allocator,
class CellAssignor>
335 if constexpr(ParticleType::is_soa_particle) {
342 int nrc = ParticleType::is_soa_particle ? NStructReal + NumRealComps() - AMREX_SPACEDIM : NStructReal + NumRealComps();
343 for (
int i = 0; i < nrc; ++i) {
344 write_real_comp.push_back(1);
348 for (
int i = 0; i < NStructInt + NumIntComps(); ++i) {
349 write_int_comp.push_back(1);
353 for (
int i = 0; i < NStructInt + NumIntComps(); ++i )
355 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
358 WriteBinaryParticleData(dir, name,
359 write_real_comp, write_int_comp,
360 real_comp_names, int_comp_names, f);
364template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
365 template<
class>
class Allocator,
class CellAssignor>
370 const std::string& name,
374 if constexpr(ParticleType::is_soa_particle) {
382 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
383 for (
int i = first_rcomp; i < NStructReal + NumRealComps(); ++i )
385 real_comp_names.push_back(getDefaultCompNameReal<ParticleType>(i));
389 for (
int i = 0; i < NStructInt + NumIntComps(); ++i )
391 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
394 WriteBinaryParticleData(dir, name, write_real_comp, write_int_comp,
395 real_comp_names, int_comp_names, f);
399template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
400 template<
class>
class Allocator,
class CellAssignor>
404WritePlotFile (
const std::string& dir,
const std::string& name,
411 BL_PROFILE(
"ParticleContainer::WritePlotFile()");
413 WriteBinaryParticleData(dir, name,
414 write_real_comp, write_int_comp,
415 real_comp_names, int_comp_names, f);
418template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
419 template<
class>
class Allocator,
class CellAssignor>
428 F const& f,
bool is_checkpoint)
const
432 write_real_comp, write_int_comp,
433 real_comp_names, int_comp_names, is_checkpoint);
437 write_real_comp, write_int_comp,
438 real_comp_names, int_comp_names,
443template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
444 template<
class>
class Allocator,
class CellAssignor>
453 BL_PROFILE(
"ParticleContainer::CheckpointPre()");
457 Long maxnextid = ParticleType::NextID();
459 for (
int lev = 0; lev < m_particles.size(); lev++) {
460 const auto& pmap = m_particles[lev];
461 for (
const auto& kv : pmap) {
462 const auto& aos = kv.second.GetArrayOfStructs();
463 for (
int k = 0; k < aos.numParticles(); ++k) {
465 if (p.id().is_valid()) {
476 ParticleType::NextID(maxnextid);
479 nparticlesPrePost = nparticles;
480 maxnextidPrePost = maxnextid;
482 nParticlesAtLevelPrePost.clear();
483 nParticlesAtLevelPrePost.resize(finestLevel() + 1, 0);
484 for(
int lev(0); lev <= finestLevel(); ++lev) {
485 nParticlesAtLevelPrePost[lev] = NumberOfParticlesAtLevel(lev);
488 whichPrePost.clear();
489 whichPrePost.resize(finestLevel() + 1);
490 countPrePost.clear();
491 countPrePost.resize(finestLevel() + 1);
492 wherePrePost.clear();
493 wherePrePost.resize(finestLevel() + 1);
495 filePrefixPrePost.clear();
496 filePrefixPrePost.resize(finestLevel() + 1);
500template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
501 template<
class>
class Allocator,
class CellAssignor>
510 BL_PROFILE(
"ParticleContainer::CheckpointPost()");
513 std::ofstream HdrFile;
514 HdrFile.open(HdrFileNamePrePost.c_str(), std::ios::out | std::ios::app);
516 for(
int lev(0); lev <= finestLevel(); ++lev) {
523 for(
int j(0); j < whichPrePost[lev].size(); ++j) {
524 HdrFile << whichPrePost[lev][j] <<
' ' << countPrePost[lev][j] <<
' ' << wherePrePost[lev][j] <<
'\n';
527 const bool gotsome = (nParticlesAtLevelPrePost[lev] > 0);
528 if(gotsome && doUnlink) {
533 for(
int i(0), N = countPrePost[lev].size(); i < N; ++i) {
534 cnt[whichPrePost[lev][i]] += countPrePost[lev][i];
537 for(
int i(0), N =
int(cnt.
size()); i < N; ++i) {
550 if( ! HdrFile.good()) {
551 amrex::Abort(
"ParticleContainer::CheckpointPost(): problem writing HdrFile");
556template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
557 template<
class>
class Allocator,
class CellAssignor>
566template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
567 template<
class>
class Allocator,
class CellAssignor>
576template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
577 template<
class>
class Allocator,
class CellAssignor>
585 bool is_checkpoint)
const
587 BL_PROFILE(
"ParticleContainer::WriteParticles()");
590 std::map<int, Vector<int> > tile_map;
592 for (
const auto& kv : m_particles[lev])
594 const int grid = kv.first.first;
595 const int tile = kv.first.second;
596 tile_map[grid].push_back(tile);
597 const auto& pflags = particle_io_flags[lev].at(kv.first);
600 count[grid] += particle_detail::countFlags(pflags);
605 MultiFab state(ParticleBoxArray(lev), ParticleDistributionMap(lev), 1,0,info);
609 const int grid = mfi.index();
614 if (count[grid] <= 0) {
continue; }
618 particle_detail::packIOData(istuff, rstuff, *
this, lev, grid,
619 write_real_comp, write_int_comp,
620 particle_io_flags, tile_map[grid], count[grid], is_checkpoint);
627 WriteParticleRealData(rstuff.
dataPtr(), rstuff.
size(), ofs);
635template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
636 template<
class>
class Allocator,
class CellAssignor>
639::Restart (
const std::string& dir,
const std::string& file,
bool )
644template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
645 template<
class>
class Allocator,
class CellAssignor>
648::Restart (
const std::string& dir,
const std::string& file)
656 int DATA_Digits_Read(5);
658 pp.
query(
"datadigits_read",DATA_Digits_Read);
660 std::string fullname = dir;
661 if (!fullname.empty() && fullname[fullname.size()-1] !=
'/') {
665 std::string HdrFileName = fullname;
666 if (!HdrFileName.empty() && HdrFileName[HdrFileName.size()-1] !=
'/') {
669 HdrFileName +=
"Header";
673 std::string fileCharPtrString(fileCharPtr.
dataPtr());
674 std::istringstream HdrFile(fileCharPtrString, std::istringstream::in);
688 bool convert_ids =
false;
689 if (version.find(
"Version_Two_Dot_One") != std::string::npos) {
692 if (version.find(
"Version_One_Dot_Zero") != std::string::npos) {
695 else if (version.find(
"Version_One_Dot_One") != std::string::npos ||
696 version.find(
"Version_Two_Dot_Zero") != std::string::npos ||
697 version.find(
"Version_Two_Dot_One") != std::string::npos) {
698 if (version.find(
"_single") != std::string::npos) {
701 else if (version.find(
"_double") != std::string::npos) {
705 std::string msg(
"ParticleContainer::Restart(): bad version string: ");
711 std::string msg(
"ParticleContainer::Restart(): unknown version string: ");
718 if (dm != AMREX_SPACEDIM) {
719 amrex::Abort(
"ParticleContainer::Restart(): dm != AMREX_SPACEDIM");
724 int nrc = ParticleType::is_soa_particle ? NStructReal + NumRealComps() - AMREX_SPACEDIM : NStructReal + NumRealComps();
726 amrex::Abort(
"ParticleContainer::Restart(): nr not the expected value");
729 std::string comp_name;
730 for (
int i = 0; i < nr; ++i) {
731 HdrFile >> comp_name;
736 if (ni != NStructInt + NumIntComps()) {
737 amrex::Abort(
"ParticleContainer::Restart(): ni != NStructInt");
740 for (
int i = 0; i < ni; ++i) {
741 HdrFile >> comp_name;
745 HdrFile >> checkpoint;
748 HdrFile >> nparticles;
752 HdrFile >> maxnextid;
754 ParticleType::NextID(maxnextid);
756 int finest_level_in_file;
757 HdrFile >> finest_level_in_file;
764 int const max_level = std::max(finest_level_in_file, finestLevel());
768 bool dual_grid =
false;
770 bool have_pheaders =
false;
771 for (
int lev = 0; lev <= finest_level_in_file; lev++)
773 std::string phdr_name = fullname;
774 phdr_name +=
"/Level_";
776 phdr_name +=
"/Particle_H";
779 have_pheaders =
true;
786 for (
int lev = 0; lev <= finestLevel(); lev++)
788 old_dms[lev] = ParticleDistributionMap(lev);
789 old_bas[lev] = ParticleBoxArray(lev);
790 std::string phdr_name = fullname;
791 phdr_name +=
"/Level_";
793 phdr_name +=
"/Particle_H";
799 std::string phdr_string(phdr_chars.
dataPtr());
800 std::istringstream phdr_file(phdr_string, std::istringstream::in);
802 particle_box_arrays[lev].readFrom(phdr_file);
803 if (! particle_box_arrays[lev].CellEqual(ParticleBoxArray(lev))) { dual_grid =
true; }
811 for (
int lev = 0; lev <= finestLevel(); lev++) {
813 if (particle_box_arrays[lev].empty()) {
814 particle_box_arrays[lev] = BoxArray(Geom(lev).Domain());
816 SetParticleBoxArray(lev, particle_box_arrays[lev]);
817 DistributionMapping pdm(particle_box_arrays[lev]);
818 SetParticleDistributionMap(lev, pdm);
822 Vector<int> ngrids(finest_level_in_file+1);
823 for (
int lev = 0; lev <= finest_level_in_file; lev++) {
824 HdrFile >> ngrids[lev];
830 if (finest_level_in_file > finestLevel()) {
831 m_particles.resize(finest_level_in_file+1);
834 for (
int lev = 0; lev <= finest_level_in_file; lev++) {
838 for (
int i = 0; i < ngrids[lev]; i++) {
839 HdrFile >> which[i] >> count[i] >> where[i];
843 if (lev <= finestLevel()) {
845 grids_to_read.push_back(mfi.index());
853 const int NReaders = MaxReaders();
854 if (rank >= NReaders) {
continue; }
856 const int Navg = ngrids[lev] / NReaders;
857 const int Nleft = ngrids[lev] - Navg * NReaders;
861 lo = rank*(Navg + 1);
865 lo = rank * Navg + Nleft;
869 for (
int i = lo; i < hi; ++i) {
870 grids_to_read.push_back(i);
874 for(
int grid : grids_to_read) {
875 if (count[grid] <= 0) {
continue; }
878 std::string name = fullname;
880 if (!name.empty() && name[name.size()-1] !=
'/') {
887 name += DataPrefix();
890 std::ifstream ParticleFile;
892 ParticleFile.open(name.c_str(), std::ios::in | std::ios::binary);
894 if (!ParticleFile.good()) {
898 ParticleFile.seekg(where[grid], std::ios::beg);
903 if (how ==
"single") {
904 if constexpr (std::is_same_v<ParticleReal, float>) {
905 ReadParticles<float>(count[grid], grid, lev, ParticleFile, finest_level_in_file, convert_ids);
907 amrex::Error(
"File contains single-precision data, while AMReX is compiled with ParticleReal==double");
910 else if (how ==
"double") {
911 if constexpr (std::is_same_v<ParticleReal, double>) {
912 ReadParticles<double>(count[grid], grid, lev, ParticleFile, finest_level_in_file, convert_ids);
914 amrex::Error(
"File contains double-precision data, while AMReX is compiled with ParticleReal==float");
918 std::string msg(
"ParticleContainer::Restart(): bad parameter: ");
923 ParticleFile.close();
925 if (!ParticleFile.good()) {
926 amrex::Abort(
"ParticleContainer::Restart(): problem reading particles");
932 for (
int lev = 0; lev <= finestLevel(); lev++) {
933 SetParticleBoxArray(lev, old_bas[lev]);
934 SetParticleDistributionMap(lev, old_dms[lev]);
945 amrex::Print() <<
"ParticleContainer::Restart() time: " << stoptime <<
'\n';
950template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
951 template<
class>
class Allocator,
class CellAssignor>
952template <
class RTYPE>
954ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
955::ReadParticles (
int cnt,
int grd,
int lev, std::ifstream& ifs,
956 int finest_level_in_file,
bool convert_ids)
958 BL_PROFILE(
"ParticleContainer::ReadParticles()");
965 const int iChunkSize = 2 + NStructInt + NumIntComps();
970 const int rChunkSize = ParticleType::is_soa_particle ? NStructReal + NumRealComps() : AMREX_SPACEDIM + NStructReal + NumRealComps();
972 ReadParticleRealData(rstuff.
dataPtr(), rstuff.
size(), ifs);
976 RTYPE* rptr = rstuff.
dataPtr();
982 host_particles.reserve(15);
983 host_particles.resize(finest_level_in_file+1);
986 std::vector<Gpu::HostVector<RTYPE> > > > host_real_attribs;
987 host_real_attribs.reserve(15);
988 host_real_attribs.resize(finest_level_in_file+1);
991 std::vector<Gpu::HostVector<int> > > > host_int_attribs;
992 host_int_attribs.reserve(15);
993 host_int_attribs.resize(finest_level_in_file+1);
996 host_idcpu.reserve(15);
997 host_idcpu.resize(finest_level_in_file+1);
999 for (
int i = 0; i < cnt; i++) {
1003 std::int32_t xi, yi;
1004 std::uint32_t xu, yu;
1007 std::memcpy(&xu, &xi,
sizeof(xi));
1008 std::memcpy(&yu, &yi,
sizeof(yi));
1009 ptemp.
m_idcpu = ((std::uint64_t)xu) << 32 | yu;
1011 ptemp.
id() = iptr[0];
1016 for (
int j = 0; j < NStructInt; j++)
1018 ptemp.
idata(j) = *iptr;
1028 rptr += AMREX_SPACEDIM;
1030 for (
int j = 0; j < NStructReal; j++)
1036 locateParticle(ptemp, pld, 0, finestLevel(), 0);
1038 std::pair<int, int> ind(grd, pld.
m_tile);
1040 host_real_attribs[lev][ind].resize(NumRealComps());
1041 host_int_attribs[lev][ind].resize(NumIntComps());
1044 if constexpr(!ParticleType::is_soa_particle)
1046 host_particles[lev][ind].push_back(ptemp);
1049 for (
int icomp = 0; icomp < NumRealComps(); icomp++) {
1050 host_real_attribs[lev][ind][icomp].push_back(*rptr);
1055 for (
int icomp = 0; icomp < NumIntComps(); icomp++) {
1056 host_int_attribs[lev][ind][icomp].push_back(*iptr);
1060 host_particles[lev][ind];
1062 for (
int j = 0; j < AMREX_SPACEDIM; j++) {
1063 host_real_attribs[lev][ind][j].push_back(ptemp.
pos(j));
1066 host_idcpu[lev][ind].push_back(ptemp.
m_idcpu);
1070 for (
int icomp = AMREX_SPACEDIM; icomp < NumRealComps(); icomp++) {
1071 host_real_attribs[lev][ind][icomp].push_back(*rptr);
1076 for (
int icomp = 0; icomp < NumIntComps(); icomp++) {
1077 host_int_attribs[lev][ind][icomp].push_back(*iptr);
1083 for (
int host_lev = 0; host_lev < static_cast<int>(host_particles.
size()); ++host_lev)
1085 for (
auto& kv : host_particles[host_lev]) {
1086 auto grid = kv.first.first;
1087 auto tile = kv.first.second;
1088 const auto& src_tile = kv.second;
1090 auto& dst_tile = DefineAndReturnParticleTile(host_lev, grid, tile);
1091 auto old_size = dst_tile.size();
1092 auto new_size = old_size;
1093 if constexpr(!ParticleType::is_soa_particle)
1095 new_size += src_tile.size();
1098 new_size += host_real_attribs[host_lev][std::make_pair(grid,tile)][0].
size();
1100 dst_tile.resize(new_size);
1102 if constexpr(!ParticleType::is_soa_particle)
1105 dst_tile.GetArrayOfStructs().begin() + old_size);
1108 host_idcpu[host_lev][std::make_pair(grid,tile)].
begin(),
1109 host_idcpu[host_lev][std::make_pair(grid,tile)].
end(),
1110 dst_tile.GetStructOfArrays().GetIdCPUData().begin() + old_size);
1113 for (
int i = 0; i < NumRealComps(); ++i) {
1115 host_real_attribs[host_lev][std::make_pair(grid,tile)][i].
begin(),
1116 host_real_attribs[host_lev][std::make_pair(grid,tile)][i].
end(),
1117 dst_tile.GetStructOfArrays().GetRealData(i).begin() + old_size);
1120 for (
int i = 0; i < NumIntComps(); ++i) {
1122 host_int_attribs[host_lev][std::make_pair(grid,tile)][i].
begin(),
1123 host_int_attribs[host_lev][std::make_pair(grid,tile)][i].
end(),
1124 dst_tile.GetStructOfArrays().GetIntData(i).begin() + old_size);
1132template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
1133 template<
class>
class Allocator,
class CellAssignor>
1135ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
1136::WriteAsciiFile (
const std::string& filename)
1138 BL_PROFILE(
"ParticleContainer::WriteAsciiFile()");
1145 Long nparticles = 0;
1147 for (
int lev = 0; lev < m_particles.size(); lev++) {
1148 auto& pmap = m_particles[lev];
1149 for (
const auto& kv : pmap) {
1150 const auto& aos = kv.second.GetArrayOfStructs();
1151 auto np = aos.numParticles();
1154 for (
int k = 0; k < np; ++k) {
1156 if (p.id().is_valid()) {
1178 File.open(filename.c_str(), std::ios::out|std::ios::trunc);
1184 File << nparticles <<
'\n';
1185 File << NStructReal <<
'\n';
1186 File << NStructInt <<
'\n';
1187 File << NumRealComps() <<
'\n';
1188 File << NumIntComps() <<
'\n';
1195 amrex::Abort(
"ParticleContainer::WriteAsciiFile(): problem writing file");
1214 File.rdbuf()->pubsetbuf(io_buffer.
dataPtr(), io_buffer.
size());
1216 File.open(filename.c_str(), std::ios::out|std::ios::app);
1224 for (
int lev = 0; lev < m_particles.size(); lev++) {
1225 auto& pmap = m_particles[lev];
1226 for (
const auto& kv : pmap) {
1229 pinned_ptile.define(NumRuntimeRealComps(), NumRuntimeIntComps(),
1231 pinned_ptile.resize(kv.second.numParticles());
1233 const auto& host_aos = pinned_ptile.GetArrayOfStructs();
1234 const auto& host_soa = pinned_ptile.GetStructOfArrays();
1236 auto np = host_aos.numParticles();
1237 for (
int index = 0; index < np; ++index) {
1239 if (it->id().is_valid()) {
1243 << it->pos(1) <<
' ',
1244 << it->pos(2) <<
' ');
1246 for (
int i = 0; i < NStructReal; i++) {
1247 File << it->rdata(i) <<
' ';
1250 File << it->id() <<
' ';
1251 File << it->cpu() <<
' ';
1253 for (
int i = 0; i < NStructInt; i++) {
1254 File << it->idata(i) <<
' ';
1258 for (
int i = 0; i < NumRealComps(); i++) {
1259 File << host_soa.GetRealData(i)[index] <<
' ';
1262 for (
int i = 0; i < NumIntComps(); i++) {
1263 File << host_soa.GetIntData(i)[index] <<
' ';
1277 amrex::Abort(
"ParticleContainer::WriteAsciiFile(): problem writing file");
1291 amrex::Print() <<
"ParticleContainer::WriteAsciiFile() time: " << stoptime <<
'\n';
#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
amrex::ParmParse pp
Input file parser instance for the given namespace.
Definition AMReX_HypreIJIface.cpp:15
#define AMREX_D_TERM(a, b, c)
Definition AMReX_SPACE.H:172
static const IntDescriptor & NativeIntDescriptor()
Returns a constant reference to an IntDescriptor describing the native "int" under which AMReX was co...
Definition AMReX_FPC.cpp:76
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
A collection (stored as an array) of FArrayBox objects.
Definition AMReX_MultiFab.H:40
const std::string & FileName() const
Definition AMReX_NFiles.H:160
Dynamically allocated vector for trivially copyable data.
Definition AMReX_PODVector.H:308
iterator begin() noexcept
Definition AMReX_PODVector.H:674
Parse Parameters From Command Line and Input Files.
Definition AMReX_ParmParse.H:349
int query(std::string_view name, bool &ref, int ival=FIRST) const
Same as querykth() but searches for the last occurrence of name.
Definition AMReX_ParmParse.cpp:1946
A distributed container for Particles sorted onto the levels, grids, and tiles of a block-structured ...
Definition AMReX_ParticleContainer.H:149
T_ParticleType ParticleType
Definition AMReX_ParticleContainer.H:151
Definition AMReX_GpuAllocators.H:156
This class provides the user with a few print options.
Definition AMReX_Print.H:35
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:28
T * dataPtr() noexcept
get access to the underlying data pointer
Definition AMReX_Vector.H:49
Long size() const noexcept
Definition AMReX_Vector.H:53
static constexpr int IO_Buffer_Size
We try to do I/O with buffers of this size.
Definition AMReX_VisMFBuffer.H:16
static Long FileOffset(std::ostream &os)
The file offset of the passed ostream.
Definition AMReX_VisMF.cpp:579
static bool GetNoFlushAfterWrite()
Definition AMReX_VisMF.H:283
amrex_particle_real ParticleReal
Floating Point Type for Particles.
Definition AMReX_REAL.H:90
amrex_long Long
Definition AMReX_INT.H:30
void WritePlotFile(const std::string &dir, const std::string &name) const
This version of WritePlotFile writes all components and assigns component names.
Definition AMReX_ParticleIO.H:113
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
void ReadAndBcastFile(const std::string &filename, Vector< char > &charBuf, bool bExitOnError, const MPI_Comm &comm)
Definition AMReX_ParallelDescriptor.cpp:1495
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 UseAsyncOut()
Definition AMReX_AsyncOut.cpp:70
bool Remove(std::string const &filename)
Definition AMReX_FileSystem.cpp:190
void copy(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:128
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 DeviceToHost deviceToHost
Definition AMReX_GpuContainers.H:106
static constexpr HostToDevice hostToDevice
Definition AMReX_GpuContainers.H:105
void streamSynchronize() noexcept
Definition AMReX_GpuDevice.H:310
void ReduceRealMax(Vector< std::reference_wrapper< Real > > const &)
Definition AMReX_ParallelDescriptor.cpp:1228
Definition AMReX_Amr.cpp:49
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 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:485
void FileOpenFailed(const std::string &file)
Output a message and abort when couldn't open the file.
Definition AMReX_Utility.cpp:137
void copyParticles(DstTile &dst, const SrcTile &src) noexcept
Copy particles from src to dst. This version copies all the particles, writing them to the beginning ...
Definition AMReX_ParticleTransformation.H:222
void readFloatData(float *data, std::size_t size, std::istream &is, const RealDescriptor &rd)
Definition AMReX_VectorIO.cpp:120
void writeFloatData(const float *data, std::size_t size, std::ostream &os, const RealDescriptor &rd=FPC::Native32RealDescriptor())
Definition AMReX_VectorIO.cpp:114
bool FileExists(const std::string &filename)
Check if a file already exists. Return true if the filename is an existing file, directory,...
Definition AMReX_Utility.cpp:145
__host__ __device__ Dim3 begin(BoxND< dim > const &box) noexcept
Definition AMReX_Box.H:2006
std::string Concatenate(const std::string &root, int num, int mindigits)
Returns rootNNNN where NNNN == num.
Definition AMReX_String.cpp:34
double second() noexcept
Definition AMReX_Utility.cpp:940
void writeDoubleData(const double *data, std::size_t size, std::ostream &os, const RealDescriptor &rd=FPC::Native64RealDescriptor())
Definition AMReX_VectorIO.cpp:126
void Error(const std::string &msg)
Print out message to cerr and exit via amrex::Abort().
Definition AMReX.cpp:234
void readDoubleData(double *data, std::size_t size, std::istream &is, const RealDescriptor &rd)
Definition AMReX_VectorIO.cpp:132
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:240
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:784
__host__ __device__ Dim3 end(BoxND< dim > const &box) noexcept
Definition AMReX_Box.H:2015
void readIntData(To *data, std::size_t size, std::istream &is, const amrex::IntDescriptor &id)
Definition AMReX_IntConv.H:36
Definition AMReX_ParticleIO.H:38
__host__ __device__ int operator()(const P &p) const
Definition AMReX_ParticleIO.H:41
FabArray memory allocation information.
Definition AMReX_FabArray.H:66
MFInfo & SetAlloc(bool a) noexcept
Definition AMReX_FabArray.H:73
uint64_t m_idcpu
Definition AMReX_Particle.H:359
__host__ __device__ bool is_valid() const noexcept
Definition AMReX_Particle.H:252
A struct used for storing a particle's position in the AMR hierarchy.
Definition AMReX_ParticleContainer.H:93
int m_tile
Definition AMReX_ParticleContainer.H:96
Definition AMReX_ParticleTile.H:723
The struct used to store particles.
Definition AMReX_Particle.H:405
__host__ __device__ RealVect pos() const &
Definition AMReX_Particle.H:456
__host__ __device__ int & idata(int index) &
Definition AMReX_Particle.H:545
__host__ __device__ ParticleCPUWrapper cpu() &
Definition AMReX_Particle.H:424
__host__ __device__ ParticleIDWrapper id() &
Definition AMReX_Particle.H:427
__host__ __device__ RealType & rdata(int index) &
Definition AMReX_Particle.H:474