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();
46template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
47 template<
class>
class Allocator,
class CellAssignor>
49ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
50::Checkpoint (
const std::string& dir,
51 const std::string& name,
bool ,
58 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
59 for (
int i = first_rcomp; i < NStructReal + NumRealComps(); ++i )
61 write_real_comp.push_back(1);
62 if (real_comp_names.empty())
64 tmp_real_comp_names.push_back(getDefaultCompNameReal<ParticleType>(i));
68 tmp_real_comp_names.push_back(real_comp_names[i-first_rcomp]);
74 for (
int i = 0; i < NStructInt + NumIntComps(); ++i )
76 write_int_comp.push_back(1);
77 if (int_comp_names.empty())
79 tmp_int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
83 tmp_int_comp_names.push_back(int_comp_names[i]);
87 WriteBinaryParticleData(dir, name, write_real_comp, write_int_comp,
88 tmp_real_comp_names, tmp_int_comp_names,
92template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
93 template<
class>
class Allocator,
class CellAssignor>
102 WriteBinaryParticleData(dir, name, write_real_comp, write_int_comp,
103 real_comp_names, int_comp_names,
107template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
108 template<
class>
class Allocator,
class CellAssignor>
116 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
117 for (
int i = first_rcomp; i < NStructReal + NumRealComps(); ++i )
119 write_real_comp.push_back(1);
120 real_comp_names.push_back(getDefaultCompNameReal<ParticleType>(i));
125 for (
int i = 0; i < NStructInt + NumIntComps(); ++i )
127 write_int_comp.push_back(1);
128 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
131 WriteBinaryParticleData(dir, name, write_real_comp, write_int_comp,
132 real_comp_names, int_comp_names,
136template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
137 template<
class>
class Allocator,
class CellAssignor>
144 if constexpr(ParticleType::is_soa_particle) {
152 int nrc = ParticleType::is_soa_particle ? NStructReal + NumRealComps() - AMREX_SPACEDIM : NStructReal + NumRealComps();
153 for (
int i = 0; i < nrc; ++i) {
154 write_real_comp.push_back(1);
158 for (
int i = 0; i < NStructInt + NumIntComps(); ++i) {
159 write_int_comp.push_back(1);
162 WriteBinaryParticleData(dir, name,
163 write_real_comp, write_int_comp,
164 real_comp_names, int_comp_names,
168template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
169 template<
class>
class Allocator,
class CellAssignor>
175 if constexpr(ParticleType::is_soa_particle) {
182 int nrc = ParticleType::is_soa_particle ? NStructReal + NumRealComps() - AMREX_SPACEDIM : NStructReal + NumRealComps();
183 for (
int i = 0; i < nrc; ++i) {
184 write_real_comp.push_back(1);
188 for (
int i = 0; i < NStructInt + NumIntComps(); ++i) {
189 write_int_comp.push_back(1);
193 for (
int i = 0; i < NStructInt + NumIntComps(); ++i )
195 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
198 WriteBinaryParticleData(dir, name,
199 write_real_comp, write_int_comp,
200 real_comp_names, int_comp_names,
204template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
205 template<
class>
class Allocator,
class CellAssignor>
209 const std::string& name,
214 if constexpr(ParticleType::is_soa_particle) {
222 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
223 for (
int i = first_rcomp; i < NStructReal + NumRealComps(); ++i )
225 real_comp_names.push_back(getDefaultCompNameReal<ParticleType>(i));
229 for (
int i = 0; i < NStructInt + NumIntComps(); ++i )
231 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
234 WriteBinaryParticleData(dir, name, write_real_comp, write_int_comp,
235 real_comp_names, int_comp_names,
239template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
240 template<
class>
class Allocator,
class CellAssignor>
243WritePlotFile (
const std::string& dir,
const std::string& name,
249 BL_PROFILE(
"ParticleContainer::WritePlotFile()");
251 WriteBinaryParticleData(dir, name,
252 write_real_comp, write_int_comp,
253 real_comp_names, int_comp_names,
257template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
258 template<
class>
class Allocator,
class CellAssignor>
259template <
class F, std::enable_if_t<!std::is_same_v<std::decay_t<F>, Vector<std::
string>>>*>
267 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
268 for (
int i = first_rcomp; i < NStructReal + NumRealComps(); ++i )
270 write_real_comp.push_back(1);
271 real_comp_names.push_back(getDefaultCompNameReal<ParticleType>(i));
276 for (
int i = 0; i < NStructInt + NumIntComps(); ++i )
278 write_int_comp.push_back(1);
279 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
282 WriteBinaryParticleData(dir, name, write_real_comp, write_int_comp,
283 real_comp_names, int_comp_names,
287template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
288 template<
class>
class Allocator,
class CellAssignor>
296 if constexpr(ParticleType::is_soa_particle) {
304 int nrc = ParticleType::is_soa_particle ? NStructReal + NumRealComps() - AMREX_SPACEDIM : NStructReal + NumRealComps();
305 for (
int i = 0; i < nrc; ++i) {
306 write_real_comp.push_back(1);
310 for (
int i = 0; i < NStructInt + NumIntComps(); ++i) {
311 write_int_comp.push_back(1);
314 WriteBinaryParticleData(dir, name,
315 write_real_comp, write_int_comp,
316 real_comp_names, int_comp_names,
320template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
321 template<
class>
class Allocator,
class CellAssignor>
322template <
class F, std::enable_if_t<!std::is_same_v<std::decay_t<F>, Vector<std::
string>>>*>
328 if constexpr(ParticleType::is_soa_particle) {
335 int nrc = ParticleType::is_soa_particle ? NStructReal + NumRealComps() - AMREX_SPACEDIM : NStructReal + NumRealComps();
336 for (
int i = 0; i < nrc; ++i) {
337 write_real_comp.push_back(1);
341 for (
int i = 0; i < NStructInt + NumIntComps(); ++i) {
342 write_int_comp.push_back(1);
346 for (
int i = 0; i < NStructInt + NumIntComps(); ++i )
348 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
351 WriteBinaryParticleData(dir, name,
352 write_real_comp, write_int_comp,
353 real_comp_names, int_comp_names,
357template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
358 template<
class>
class Allocator,
class CellAssignor>
363 const std::string& name,
367 if constexpr(ParticleType::is_soa_particle) {
375 int first_rcomp = ParticleType::is_soa_particle ? AMREX_SPACEDIM : 0;
376 for (
int i = first_rcomp; i < NStructReal + NumRealComps(); ++i )
378 real_comp_names.push_back(getDefaultCompNameReal<ParticleType>(i));
382 for (
int i = 0; i < NStructInt + NumIntComps(); ++i )
384 int_comp_names.push_back(getDefaultCompNameInt<ParticleType>(i));
387 WriteBinaryParticleData(dir, name, write_real_comp, write_int_comp,
388 real_comp_names, int_comp_names,
392template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
393 template<
class>
class Allocator,
class CellAssignor>
397WritePlotFile (
const std::string& dir,
const std::string& name,
404 BL_PROFILE(
"ParticleContainer::WritePlotFile()");
406 WriteBinaryParticleData(dir, name,
407 write_real_comp, write_int_comp,
408 real_comp_names, int_comp_names,
412template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
413 template<
class>
class Allocator,
class CellAssignor>
422 F&& f,
bool is_checkpoint)
const
426 write_real_comp, write_int_comp,
427 real_comp_names, int_comp_names, is_checkpoint);
431 write_real_comp, write_int_comp,
432 real_comp_names, int_comp_names,
433 std::forward<F>(f), is_checkpoint);
437template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
438 template<
class>
class Allocator,
class CellAssignor>
447 BL_PROFILE(
"ParticleContainer::CheckpointPre()");
451 Long maxnextid = ParticleType::NextID();
453 for (
int lev = 0; lev < m_particles.size(); lev++) {
454 const auto& pmap = m_particles[lev];
455 for (
const auto& kv : pmap) {
456 const auto& aos = kv.second.GetArrayOfStructs();
457 for (
int k = 0; k < aos.numParticles(); ++k) {
459 if (p.id().is_valid()) {
470 ParticleType::NextID(maxnextid);
473 nparticlesPrePost = nparticles;
474 maxnextidPrePost = int(maxnextid);
476 nParticlesAtLevelPrePost.clear();
477 nParticlesAtLevelPrePost.resize(finestLevel() + 1, 0);
478 for(
int lev(0); lev <= finestLevel(); ++lev) {
479 nParticlesAtLevelPrePost[lev] = NumberOfParticlesAtLevel(lev);
482 whichPrePost.clear();
483 whichPrePost.resize(finestLevel() + 1);
484 countPrePost.clear();
485 countPrePost.resize(finestLevel() + 1);
486 wherePrePost.clear();
487 wherePrePost.resize(finestLevel() + 1);
489 filePrefixPrePost.clear();
490 filePrefixPrePost.resize(finestLevel() + 1);
494template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
495 template<
class>
class Allocator,
class CellAssignor>
504 BL_PROFILE(
"ParticleContainer::CheckpointPost()");
507 std::ofstream HdrFile;
508 HdrFile.open(HdrFileNamePrePost.c_str(), std::ios::out | std::ios::app);
510 for(
int lev(0); lev <= finestLevel(); ++lev) {
517 for(
int j(0); j < whichPrePost[lev].size(); ++j) {
518 HdrFile << whichPrePost[lev][j] <<
' ' << countPrePost[lev][j] <<
' ' << wherePrePost[lev][j] <<
'\n';
521 const bool gotsome = (nParticlesAtLevelPrePost[lev] > 0);
522 if(gotsome && doUnlink) {
527 for(
int i(0), N = countPrePost[lev].size(); i < N; ++i) {
528 cnt[whichPrePost[lev][i]] += countPrePost[lev][i];
531 for(
int i(0), N =
int(cnt.
size()); i < N; ++i) {
544 if( ! HdrFile.good()) {
545 amrex::Abort(
"ParticleContainer::CheckpointPost(): problem writing HdrFile");
550template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
551 template<
class>
class Allocator,
class CellAssignor>
560template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
561 template<
class>
class Allocator,
class CellAssignor>
570template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
571 template<
class>
class Allocator,
class CellAssignor>
578 const Vector<std::map<std::pair<int, int>,
IntVector>>& particle_io_flags,
579 bool is_checkpoint)
const
581 BL_PROFILE(
"ParticleContainer::WriteParticles()");
584 std::map<int, Vector<int> > tile_map;
586 for (
const auto& kv : m_particles[lev])
588 const int grid = kv.first.first;
589 const int tile = kv.first.second;
590 tile_map[grid].push_back(tile);
591 const auto& pflags = particle_io_flags[lev].at(kv.first);
599 MultiFab state(ParticleBoxArray(lev), ParticleDistributionMap(lev), 1,0,info);
603 const int grid = mfi.index();
608 if (count[grid] == 0) {
continue; }
613 write_real_comp, write_int_comp,
614 particle_io_flags, tile_map[grid], count[grid], is_checkpoint);
621 WriteParticleRealData(rstuff.
dataPtr(), rstuff.
size(), ofs);
629template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
630 template<
class>
class Allocator,
class CellAssignor>
633::Restart (
const std::string& dir,
const std::string& file,
bool )
638template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
639 template<
class>
class Allocator,
class CellAssignor>
642::Restart (
const std::string& dir,
const std::string& file)
650 int DATA_Digits_Read(5);
652 pp.
query(
"datadigits_read",DATA_Digits_Read);
654 std::string fullname = dir;
655 if (!fullname.empty() && fullname[fullname.size()-1] !=
'/') {
659 std::string HdrFileName = fullname;
660 if (!HdrFileName.empty() && HdrFileName[HdrFileName.size()-1] !=
'/') {
663 HdrFileName +=
"Header";
667 std::string fileCharPtrString(fileCharPtr.
dataPtr());
668 std::istringstream HdrFile(fileCharPtrString, std::istringstream::in);
682 bool convert_ids =
false;
683 if (version.find(
"Version_Two_Dot_One") != std::string::npos) {
686 if (version.find(
"Version_One_Dot_Zero") != std::string::npos) {
689 else if (version.find(
"Version_One_Dot_One") != std::string::npos ||
690 version.find(
"Version_Two_Dot_Zero") != std::string::npos ||
691 version.find(
"Version_Two_Dot_One") != std::string::npos) {
692 if (version.find(
"_single") != std::string::npos) {
695 else if (version.find(
"_double") != std::string::npos) {
699 std::string msg(
"ParticleContainer::Restart(): bad version string: ");
705 std::string msg(
"ParticleContainer::Restart(): unknown version string: ");
712 if (dm != AMREX_SPACEDIM) {
713 amrex::Abort(
"ParticleContainer::Restart(): dm != AMREX_SPACEDIM");
718 int nrc = ParticleType::is_soa_particle ? NStructReal + NumRealComps() - AMREX_SPACEDIM : NStructReal + NumRealComps();
720 amrex::Abort(
"ParticleContainer::Restart(): nr not the expected value");
723 std::string comp_name;
724 for (
int i = 0; i < nr; ++i) {
725 HdrFile >> comp_name;
730 if (ni != NStructInt + NumIntComps()) {
731 amrex::Abort(
"ParticleContainer::Restart(): ni != NStructInt");
734 for (
int i = 0; i < ni; ++i) {
735 HdrFile >> comp_name;
739 HdrFile >> checkpoint;
742 HdrFile >> nparticles;
746 HdrFile >> maxnextid;
748 ParticleType::NextID(maxnextid);
750 int finest_level_in_file;
751 HdrFile >> finest_level_in_file;
758 bool dual_grid =
false;
760 bool have_pheaders =
false;
761 for (
int lev = 0; lev <= finest_level_in_file; lev++)
763 std::string phdr_name = fullname;
764 phdr_name +=
"/Level_";
766 phdr_name +=
"/Particle_H";
769 have_pheaders =
true;
776 for (
int lev = 0; lev <= finestLevel(); lev++)
778 old_dms[lev] = ParticleDistributionMap(lev);
779 old_bas[lev] = ParticleBoxArray(lev);
780 std::string phdr_name = fullname;
781 phdr_name +=
"/Level_";
783 phdr_name +=
"/Particle_H";
789 std::string phdr_string(phdr_chars.
dataPtr());
790 std::istringstream phdr_file(phdr_string, std::istringstream::in);
792 particle_box_arrays[lev].readFrom(phdr_file);
793 if (! particle_box_arrays[lev].CellEqual(ParticleBoxArray(lev))) { dual_grid =
true; }
801 for (
int lev = 0; lev <= finestLevel(); lev++) {
803 if (particle_box_arrays[lev].empty()) {
804 particle_box_arrays[lev] = BoxArray(Geom(lev).Domain());
806 SetParticleBoxArray(lev, particle_box_arrays[lev]);
808 SetParticleDistributionMap(lev, pdm);
812 Vector<int> ngrids(finest_level_in_file+1);
813 for (
int lev = 0; lev <= finest_level_in_file; lev++) {
814 HdrFile >> ngrids[lev];
820 if (finest_level_in_file > finestLevel()) {
821 m_particles.resize(finest_level_in_file+1);
824 for (
int lev = 0; lev <= finest_level_in_file; lev++) {
828 for (
int i = 0; i < ngrids[lev]; i++) {
829 HdrFile >> which[i] >> count[i] >> where[i];
833 if (lev <= finestLevel()) {
835 grids_to_read.push_back(mfi.index());
843 const int NReaders = MaxReaders();
844 if (rank >= NReaders) {
return; }
846 const int Navg = ngrids[lev] / NReaders;
847 const int Nleft = ngrids[lev] - Navg * NReaders;
851 lo = rank*(Navg + 1);
855 lo = rank * Navg + Nleft;
859 for (
int i = lo; i < hi; ++i) {
860 grids_to_read.push_back(i);
864 for(
int grid : grids_to_read) {
865 if (count[grid] <= 0) {
continue; }
868 std::string name = fullname;
870 if (!name.empty() && name[name.size()-1] !=
'/') {
877 name += DataPrefix();
880 std::ifstream ParticleFile;
882 ParticleFile.open(name.c_str(), std::ios::in | std::ios::binary);
884 if (!ParticleFile.good()) {
888 ParticleFile.seekg(where[grid], std::ios::beg);
893 if (how ==
"single") {
894 if constexpr (std::is_same_v<ParticleReal, float>) {
895 ReadParticles<float>(count[grid], grid, lev, ParticleFile, finest_level_in_file, convert_ids);
897 amrex::Error(
"File contains single-precision data, while AMReX is compiled with ParticleReal==double");
900 else if (how ==
"double") {
901 if constexpr (std::is_same_v<ParticleReal, double>) {
902 ReadParticles<double>(count[grid], grid, lev, ParticleFile, finest_level_in_file, convert_ids);
904 amrex::Error(
"File contains double-precision data, while AMReX is compiled with ParticleReal==float");
908 std::string msg(
"ParticleContainer::Restart(): bad parameter: ");
913 ParticleFile.close();
915 if (!ParticleFile.good()) {
916 amrex::Abort(
"ParticleContainer::Restart(): problem reading particles");
922 for (
int lev = 0; lev <= finestLevel(); lev++) {
923 SetParticleBoxArray(lev, old_bas[lev]);
924 SetParticleDistributionMap(lev, old_dms[lev]);
935 amrex::Print() <<
"ParticleContainer::Restart() time: " << stoptime <<
'\n';
940template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
941 template<
class>
class Allocator,
class CellAssignor>
942template <
class RTYPE>
944ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
945::ReadParticles (
int cnt,
int grd,
int lev, std::ifstream& ifs,
946 int finest_level_in_file,
bool convert_ids)
948 BL_PROFILE(
"ParticleContainer::ReadParticles()");
955 const int iChunkSize = 2 + NStructInt + NumIntComps();
960 const int rChunkSize = ParticleType::is_soa_particle ? NStructReal + NumRealComps() : AMREX_SPACEDIM + NStructReal + NumRealComps();
962 ReadParticleRealData(rstuff.
dataPtr(), rstuff.
size(), ifs);
966 RTYPE* rptr = rstuff.
dataPtr();
972 host_particles.reserve(15);
973 host_particles.resize(finest_level_in_file+1);
976 std::vector<Gpu::HostVector<RTYPE> > > > host_real_attribs;
977 host_real_attribs.reserve(15);
978 host_real_attribs.resize(finest_level_in_file+1);
981 std::vector<Gpu::HostVector<int> > > > host_int_attribs;
982 host_int_attribs.reserve(15);
983 host_int_attribs.resize(finest_level_in_file+1);
986 host_idcpu.reserve(15);
987 host_idcpu.resize(finestLevel()+1);
989 for (
int i = 0; i < cnt; i++) {
994 std::uint32_t xu, yu;
997 std::memcpy(&xu, &xi,
sizeof(xi));
998 std::memcpy(&yu, &yi,
sizeof(yi));
999 ptemp.
m_idcpu = ((std::uint64_t)xu) << 32 | yu;
1002 ptemp.
cpu() = iptr[1];
1006 for (
int j = 0; j < NStructInt; j++)
1008 ptemp.
idata(j) = *iptr;
1015 ptemp.
pos(1) = ParticleReal(rptr[1]);,
1016 ptemp.
pos(2) = ParticleReal(rptr[2]););
1018 rptr += AMREX_SPACEDIM;
1020 for (
int j = 0; j < NStructReal; j++)
1022 ptemp.
rdata(j) = ParticleReal(*rptr);
1026 locateParticle(ptemp, pld, 0, finestLevel(), 0);
1028 std::pair<int, int> ind(grd, pld.
m_tile);
1030 host_real_attribs[lev][ind].resize(NumRealComps());
1031 host_int_attribs[lev][ind].resize(NumIntComps());
1034 if constexpr(!ParticleType::is_soa_particle)
1036 host_particles[lev][ind].push_back(ptemp);
1039 for (
int icomp = 0; icomp < NumRealComps(); icomp++) {
1040 host_real_attribs[lev][ind][icomp].push_back(*rptr);
1045 for (
int icomp = 0; icomp < NumIntComps(); icomp++) {
1046 host_int_attribs[lev][ind][icomp].push_back(*iptr);
1050 host_particles[lev][ind];
1052 for (
int j = 0; j < AMREX_SPACEDIM; j++) {
1053 host_real_attribs[lev][ind][j].push_back(ptemp.
pos(j));
1056 host_idcpu[lev][ind].push_back(ptemp.
m_idcpu);
1060 for (
int icomp = AMREX_SPACEDIM; icomp < NumRealComps(); icomp++) {
1061 host_real_attribs[lev][ind][icomp].push_back(*rptr);
1066 for (
int icomp = 0; icomp < NumIntComps(); icomp++) {
1067 host_int_attribs[lev][ind][icomp].push_back(*iptr);
1073 for (
int host_lev = 0; host_lev < static_cast<int>(host_particles.
size()); ++host_lev)
1075 for (
auto& kv : host_particles[host_lev]) {
1076 auto grid = kv.first.first;
1077 auto tile = kv.first.second;
1078 const auto& src_tile = kv.second;
1080 auto& dst_tile = DefineAndReturnParticleTile(host_lev, grid, tile);
1081 auto old_size = dst_tile.size();
1082 auto new_size = old_size;
1083 if constexpr(!ParticleType::is_soa_particle)
1085 new_size += src_tile.size();
1088 new_size += host_real_attribs[host_lev][std::make_pair(grid,tile)][0].
size();
1090 dst_tile.resize(new_size);
1092 if constexpr(!ParticleType::is_soa_particle)
1095 dst_tile.GetArrayOfStructs().begin() + old_size);
1098 host_idcpu[host_lev][std::make_pair(grid,tile)].
begin(),
1099 host_idcpu[host_lev][std::make_pair(grid,tile)].
end(),
1100 dst_tile.GetStructOfArrays().GetIdCPUData().begin() + old_size);
1103 for (
int i = 0; i < NumRealComps(); ++i) {
1105 host_real_attribs[host_lev][std::make_pair(grid,tile)][i].
begin(),
1106 host_real_attribs[host_lev][std::make_pair(grid,tile)][i].
end(),
1107 dst_tile.GetStructOfArrays().GetRealData(i).begin() + old_size);
1110 for (
int i = 0; i < NumIntComps(); ++i) {
1112 host_int_attribs[host_lev][std::make_pair(grid,tile)][i].
begin(),
1113 host_int_attribs[host_lev][std::make_pair(grid,tile)][i].
end(),
1114 dst_tile.GetStructOfArrays().GetIntData(i).begin() + old_size);
1122template <
typename ParticleType,
int NArrayReal,
int NArrayInt,
1123 template<
class>
class Allocator,
class CellAssignor>
1125ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
1126::WriteAsciiFile (
const std::string& filename)
1128 BL_PROFILE(
"ParticleContainer::WriteAsciiFile()");
1135 Long nparticles = 0;
1137 for (
int lev = 0; lev < m_particles.size(); lev++) {
1138 auto& pmap = m_particles[lev];
1139 for (
const auto& kv : pmap) {
1140 const auto& aos = kv.second.GetArrayOfStructs();
1141 auto np = aos.numParticles();
1144 for (
int k = 0; k < np; ++k) {
1146 if (p.id().is_valid()) {
1168 File.open(filename.c_str(), std::ios::out|std::ios::trunc);
1174 File << nparticles <<
'\n';
1175 File << NStructReal <<
'\n';
1176 File << NStructInt <<
'\n';
1177 File << NumRealComps() <<
'\n';
1178 File << NumIntComps() <<
'\n';
1185 amrex::Abort(
"ParticleContainer::WriteAsciiFile(): problem writing file");
1204 File.rdbuf()->pubsetbuf(io_buffer.
dataPtr(), io_buffer.
size());
1206 File.open(filename.c_str(), std::ios::out|std::ios::app);
1214 for (
int lev = 0; lev < m_particles.size(); lev++) {
1215 auto& pmap = m_particles[lev];
1216 for (
const auto& kv : pmap) {
1219 pinned_ptile.define(NumRuntimeRealComps(), NumRuntimeIntComps(),
1221 pinned_ptile.resize(kv.second.numParticles());
1223 const auto& host_aos = pinned_ptile.GetArrayOfStructs();
1224 const auto& host_soa = pinned_ptile.GetStructOfArrays();
1226 auto np = host_aos.numParticles();
1227 for (
int index = 0; index < np; ++index) {
1229 if (it->id().is_valid()) {
1233 << it->pos(1) <<
' ',
1234 << it->pos(2) <<
' ');
1236 for (
int i = 0; i < NStructReal; i++) {
1237 File << it->rdata(i) <<
' ';
1240 File << it->id() <<
' ';
1241 File << it->cpu() <<
' ';
1243 for (
int i = 0; i < NStructInt; i++) {
1244 File << it->idata(i) <<
' ';
1248 for (
int i = 0; i < NumRealComps(); i++) {
1249 File << host_soa.GetRealData(i)[index] <<
' ';
1252 for (
int i = 0; i < NumIntComps(); i++) {
1253 File << host_soa.GetIntData(i)[index] <<
' ';
1267 amrex::Abort(
"ParticleContainer::WriteAsciiFile(): problem writing file");
1281 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
Calculates the distribution of FABs to MPI processes.
Definition AMReX_DistributionMapping.H:41
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
Definition AMReX_MFIter.H:57
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition AMReX_MFIter.H:141
A collection (stored as an array) of FArrayBox objects.
Definition AMReX_MultiFab.H:38
const std::string & FileName() const
Definition AMReX_NFiles.H:160
Definition AMReX_PODVector.H:297
iterator begin() noexcept
Definition AMReX_PODVector.H:663
Parse Parameters From Command Line and Input Files.
Definition AMReX_ParmParse.H:343
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:1393
A distributed container for Particles sorted onto the levels, grids, and tiles of a block-structured ...
Definition AMReX_ParticleContainer.H:147
typename SoA::IntVector IntVector
Definition AMReX_ParticleContainer.H:195
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:111
T_ParticleType ParticleType
Definition AMReX_ParticleContainer.H:149
Definition AMReX_GpuAllocators.H:144
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
@ IO_Buffer_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:282
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:121
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:221
static constexpr DeviceToHost deviceToHost
Definition AMReX_GpuContainers.H:99
static constexpr HostToDevice hostToDevice
Definition AMReX_GpuContainers.H:98
void streamSynchronize() noexcept
Definition AMReX_GpuDevice.H:260
void ReduceIntSum(int &)
Integer sum reduction.
Definition AMReX_ParallelDescriptor.cpp:1261
int MyProc() noexcept
return the rank number local to the current Parallel Context
Definition AMReX_ParallelDescriptor.H:125
void ReadAndBcastFile(const std::string &filename, Vector< char > &charBuf, bool bExitOnError, const MPI_Comm &comm)
Definition AMReX_ParallelDescriptor.cpp:1491
void ReduceLongMax(Long &)
Long max reduction.
Definition AMReX_ParallelDescriptor.cpp:1233
void ReduceLongSum(Long &)
Long sum reduction.
Definition AMReX_ParallelDescriptor.cpp:1232
int NProcs() noexcept
return the number of MPI ranks local to the current Parallel Context
Definition AMReX_ParallelDescriptor.H:243
int IOProcessorNumber() noexcept
Definition AMReX_ParallelDescriptor.H:266
void Barrier(const std::string &)
Definition AMReX_ParallelDescriptor.cpp:1211
bool IOProcessor() noexcept
Is this CPU the I/O Processor? To get the rank number, call IOProcessorNumber()
Definition AMReX_ParallelDescriptor.H:275
void ReduceRealMax(Vector< std::reference_wrapper< Real > > const &)
Definition AMReX_ParallelDescriptor.cpp:1224
void packIOData(Vector< int > &idata, Vector< ParticleReal > &rdata, const PC &pc, int lev, int grid, const Vector< int > &write_real_comp, const Vector< int > &write_int_comp, const Vector< std::map< std::pair< int, int >, typename PC::IntVector > > &particle_io_flags, const Vector< int > &tiles, int np, bool is_checkpoint)
Definition AMReX_WriteBinaryParticleData.H:438
amrex::Long countFlags(const Vector< std::map< std::pair< int, int >, Container< int, Allocator > > > &particle_io_flags, const PC &pc)
Definition AMReX_WriteBinaryParticleData.H:183
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:138
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:465
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:161
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:1899
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:224
Arena * The_Pinned_Arena()
Definition AMReX_Arena.cpp:745
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:230
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:763
__host__ __device__ Dim3 end(BoxND< dim > const &box) noexcept
Definition AMReX_Box.H:1908
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:92
int m_tile
Definition AMReX_ParticleContainer.H:95
Definition AMReX_ParticleTile.H:720
The struct used to store particles.
Definition AMReX_Particle.H:402
__host__ __device__ RealVect pos() const &
Definition AMReX_Particle.H:452
__host__ __device__ int & idata(int index) &
Definition AMReX_Particle.H:541
__host__ __device__ ParticleCPUWrapper cpu() &
Definition AMReX_Particle.H:420
__host__ __device__ ParticleIDWrapper id() &
Definition AMReX_Particle.H:423
__host__ __device__ RealType & rdata(int index) &
Definition AMReX_Particle.H:470