Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
AMReX_ParticleInit.H
Go to the documentation of this file.
1#ifndef AMREX_PARTICLEINIT_H
2#define AMREX_PARTICLEINIT_H
3#include <AMReX_Config.H>
4
5#include <iterator>
6#include <utility>
7
8namespace amrex {
9
10/*
11 \brief Initialize particles from an Ascii file in the following format:
12
13 8
14 0.1 0.1 16.2 1000.0 4.0 1.0 6.0
15 8.1 0.1 16.2 1000.0 -5 0.0 -7.0
16 16.1 0.1 16.2 1000.0 6.0 -8.0 2.0
17 24.1 0.1 16.2 1000.0 9.0 4.0 8.0
18 0.1 8.1 16.2 1000.0 -8.0 -3.0 -10.0
19 8.1 8.1 16.2 1000.0 2.0 1.0 0.0
20 16.1 8.1 16.2 1000.0 0.0 2.0 3.0
21 24.1 8.1 16.2 1000.0 -9.0 7.0 5.0
22
23 The first line is the number of particles. The remaining lines give the particle
24 data to read, one particle on each line. The first AMREX_SPACEDIM components are
25 positions. The next extradata components are the additional real data to read in.
26
27 Integer data is not currently supported by this function.
28
29 Parameters:
30
31 file: the name of the Ascii file with the particles
32 extradata: the number of real components to read in, beyond the AMREX_SPACEDIM
33 positions
34 Nrep: pointer to IntVect that lets you replicate the incoming particles
35 across the domain so that you only need to specify a sub-volume of
36 them. By default particles are not replicated.
37 */
38template <typename ParticleType, int NArrayReal, int NArrayInt,
39 template<class> class Allocator, class CellAssignor>
40void
41ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt, Allocator, CellAssignor>
42::InitFromAsciiFile (const std::string& file, int extradata, const IntVect* Nrep)
43{
44 BL_PROFILE("ParticleContainer<NSR, NSI, NAR, NAI>::InitFromAsciiFile()");
45 AMREX_ASSERT(!file.empty());
46 AMREX_ASSERT(extradata <= NStructReal + NumRealComps());
47
48 const int MyProc = ParallelDescriptor::MyProc();
49 const int NProcs = ParallelDescriptor::NProcs();
50 const auto strttime = amrex::second();
51
52 int NReaders = MaxReaders(); // num ranks to read with
53 int NRedist = 1; // number of times to redistribute while reading
54
55 // try to use some sensible heuristics
56 if (NProcs <= 1024)
57 {
58 if (NReaders > 1)
59 {
60 NRedist = 2;
61 }
62 }
63 else if (NProcs <= 4096)
64 {
65 NReaders = std::max(NReaders,128);
66 NRedist = 4;
67 }
68 else if (NProcs <= 8192)
69 {
70 NReaders = std::max(NReaders,384);
71 NRedist = 32;
72 }
73 else if (NProcs <= 16384)
74 {
75 NReaders = std::max(NReaders,512);
76 NRedist = 48;
77 }
78
79 resizeData();
80
81 IntVect lNrep(AMREX_D_DECL(1,1,1));
82
83 if (Nrep != nullptr) {
84 lNrep = *Nrep;
85 }
86
87 Long how_many = 0;
88 Long how_many_read = 0;
89
92 if (extradata > NStructReal) { nreals.resize(extradata - NStructReal); }
93
94 if (MyProc < NReaders)
95 {
97
98 std::ifstream ifs;
99
100 ifs.rdbuf()->pubsetbuf(io_buffer.dataPtr(), io_buffer.size());
101
102 ifs.open(file.c_str(), std::ios::in);
103
104 if (!ifs.good())
105 {
107 }
108
109 int cnt = 0;
110
111 ifs >> cnt >> std::ws;
112
113 ParticleLocData pld;
114 ParticleType p, p_rep;
115 Vector<Real> r;
116
117 if (extradata > NStructReal) { r.resize(extradata - NStructReal); }
118
119 const int Chunk = cnt / NReaders;
120
121 for (int i = 0; i < MyProc*Chunk; i++)
122 {
123 ifs.ignore(std::numeric_limits<std::streamsize>::max(), '\n');
124 ifs >> std::ws; // Eat newline.
125 }
126
127 if (!ifs.good())
128 {
129 std::string msg("ParticleContainer::InitFromAsciiFile(");
130 msg += file;
131 msg += ") failed @ 1";
132 amrex::Error(msg.c_str());
133 }
134
135 int MyCnt = Chunk;
136
137 if (MyProc == (NReaders - 1))
138 {
139 MyCnt += cnt % NReaders;
140 }
141
142 const Geometry& geom = Geom(0);
143
144 for (int i = 0; i < MyCnt; i++)
145 {
146 AMREX_D_TERM(ifs >> p.pos(0);,
147 ifs >> p.pos(1);,
148 ifs >> p.pos(2););
149
150 for (int n = 0; n < extradata; n++)
151 {
152 if (n < NStructReal)
153 {
154 ifs >> p.rdata(n);
155 }
156 else
157 {
158 ifs >> r[n - NStructReal];
159 }
160 }
161
162 if (!ifs.good())
163 {
164 std::string msg("ParticleContainer::InitFromAsciiFile(");
165 msg += file; msg += ") failed @ 2";
166 amrex::Error(msg.c_str());
167 }
168
169 if (!Where(p, pld))
170 {
171 PeriodicShift(p);
172
173 if (!Where(p, pld))
174 {
175 if (m_verbose) {
176 amrex::AllPrint() << "BAD PARTICLE ID WOULD BE " << ParticleType::NextID() << '\n'
177 << "BAD PARTICLE POS "
178 << AMREX_D_TERM( p.pos(0),
179 << p.pos(1),
180 << p.pos(2))
181 << "\n";
182 }
183 amrex::Abort("ParticleContainer::InitFromAsciiFile(): invalid particle");
184 }
185 }
186
187 // set these rather than reading them in
188 p.id() = ParticleType::NextID();
189 p.cpu() = MyProc;
190
191 nparticles.push_back(p);
192 if(nreals.size() > extradata - NStructReal) {
193 for (int n = NStructReal; n < extradata; n++)
194 {
195 nreals[n-NStructReal].push_back(r[n-NStructReal]);
196 }
197 }
198
199 how_many++;
200 how_many_read++;
201
202 const Real DomSize[AMREX_SPACEDIM] =
203 { AMREX_D_DECL((geom.ProbHi(0)-geom.ProbLo(0))/lNrep[0],
204 (geom.ProbHi(1)-geom.ProbLo(1))/lNrep[1],
205 (geom.ProbHi(2)-geom.ProbLo(2))/lNrep[2]) };
206 int rep[AMREX_SPACEDIM];
207
208#if AMREX_SPACEDIM > 2
209 for (rep[2] = 1; rep[2] <= lNrep[2]; rep[2]++)
210 {
211#endif
212#if AMREX_SPACEDIM > 1
213 for (rep[1] = 1; rep[1] <= lNrep[1]; rep[1]++)
214 {
215#endif
216 for (rep[0] = 1; rep[0] <= lNrep[0]; rep[0]++)
217 {
218 if (!(AMREX_D_TERM( (rep[0] == 1), && (rep[1] == 1), && (rep[2] == 1) ) ) )
219 {
220 // Shift the position.
221 for (int d=0; d<AMREX_SPACEDIM; ++d)
222 {
223 p_rep.pos(d) = static_cast<ParticleReal>(p.pos(d) + Real(rep[d]-1)*DomSize[d]);
224 }
225
226 // Copy the extra data
227 for (int n = 0; n < extradata; n++)
228 {
229 if (n < NStructReal)
230 {
231 p_rep.rdata(n) = p.rdata(n);
232 }
233 }
234
235 if (!Where(p_rep, pld))
236 {
237 PeriodicShift(p_rep);
238 if (!Where(p_rep, pld))
239 {
240 if (m_verbose) {
241 amrex::AllPrint() << "BAD REPLICATED PARTICLE ID WOULD BE " << ParticleType::NextID() << "\n";
242 }
243 amrex::Abort("ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt>::InitFromAsciiFile(): invalid replicated particle");
244 }
245 }
246
247 p_rep.id() = ParticleType::NextID();
248 p_rep.cpu() = MyProc;
249
250 nparticles.push_back(p_rep);
251
252 if (nreals.size() > extradata - NStructReal) {
253 for (int n = NStructReal; n < extradata; n++)
254 {
255 nreals[n-NStructReal].push_back(r[n-NStructReal]);
256 }
257 }
258
259 how_many++;
260 }
261 }
262#if AMREX_SPACEDIM > 1
263 }
264#endif
265#if AMREX_SPACEDIM > 2
266 }
267#endif
268
269 }
270
271 }
272
273 // Now Redistribute() each chunk separately to minimize memory bloat.
274 int NRedist_chunk = NReaders / NRedist;
275
276 ParticleLocData pld;
277
278 for (int nr = 0; nr < NRedist; nr++)
279 {
281 host_particles.reserve(15);
282 host_particles.resize(finestLevel()+1);
283
284 Vector<std::map<std::pair<int, int>, std::array<Gpu::HostVector<ParticleReal>, NArrayReal > > > host_real_attribs;
285 host_real_attribs.reserve(15);
286 host_real_attribs.resize(finestLevel()+1);
287
288 if (m_verbose > 0) {
289 amrex::Print() << "Redistributing from processor "
290 << nr*NRedist_chunk << " to "
291 << (nr+1)*NRedist_chunk-1 << '\n';
292 }
293
294 for (int which = nr*NRedist_chunk; which < (nr+1)*NRedist_chunk; which++)
295 {
296 if (which == MyProc)
297 {
298 while (!nparticles.empty())
299 {
300 ParticleType& p = nparticles.back();
301 Where(p, pld);
302
303 host_particles[pld.m_lev][std::make_pair(pld.m_grid, pld.m_tile)].push_back(p);
304
305 if (nreals.size() > extradata - NStructReal && NArrayReal > 0)
306 {
307 for (int n = NStructReal; n < extradata; n++)
308 {
309 Real rdata = nreals[n-NStructReal].back();
310 host_real_attribs[pld.m_lev][std::make_pair(pld.m_grid, pld.m_tile)][n-NStructReal].push_back(rdata);
311 }
312 }
313
314 nparticles.pop_back();
315
316 if (nreals.size() > extradata - NStructReal)
317 {
318 for (int n = NStructReal; n < extradata; n++)
319 {
320 nreals[n-NStructReal].pop_back();
321 }
322 }
323 }
324 }
325 }
326
327 for (int lev = 0; lev < std::ssize(host_particles); ++lev)
328 {
329 for (auto& kv : host_particles[lev])
330 {
331
332 auto grid = kv.first.first;
333 auto tile = kv.first.second;
334 const auto& src_tile = kv.second;
335
336 auto& dst_tile = GetParticles(lev)[std::make_pair(grid,tile)];
337 auto old_size = dst_tile.GetArrayOfStructs().size();
338 auto new_size = old_size + src_tile.size();
339 dst_tile.resize(new_size);
340
341 Gpu::copyAsync(Gpu::hostToDevice, src_tile.begin(), src_tile.end(),
342 dst_tile.GetArrayOfStructs().begin() + old_size);
343
344 if (std::ssize(host_real_attribs[lev][std::make_pair(grid, tile)]) > NArrayReal) {
345 for (int i = 0; i < NArrayReal; ++i) {
347 host_real_attribs[lev][std::make_pair(grid,tile)][i].begin(),
348 host_real_attribs[lev][std::make_pair(grid,tile)][i].end(),
349 dst_tile.GetStructOfArrays().GetRealData(i).begin() + old_size);
350 }
351 }
352 }
353 }
355
356 Redistribute();
357 }
358
359 if (m_verbose > 0)
360 {
361 if (NRedist*NRedist_chunk < NReaders) {
362 amrex::Print() << "Redistributing from processor "
363 << NRedist*NRedist_chunk << " to "
364 << NReaders << '\n';
365 }
366 }
367 for (int which = NRedist*NRedist_chunk; which < NReaders; which++)
368 {
370 host_particles.reserve(15);
371 host_particles.resize(finestLevel()+1);
372
373 Vector<std::map<std::pair<int, int>, std::array<Gpu::HostVector<ParticleReal>, NArrayReal > > > host_real_attribs;
374 host_real_attribs.reserve(15);
375 host_real_attribs.resize(finestLevel()+1);
376
377 if (which == MyProc)
378 {
379 while (!nparticles.empty())
380 {
381 ParticleType& p = nparticles.back();
382 Where(p, pld);
383
384 host_particles[pld.m_lev][std::make_pair(pld.m_grid, pld.m_tile)].push_back(p);
385 if (std::ssize(host_real_attribs[pld.m_lev][std::make_pair(pld.m_grid, pld.m_tile)]) >
386 extradata - NStructReal) {
387 for (int n = NStructReal; n < extradata; n++)
388 {
389 Real rdata = nreals[n-NStructReal].back();
390 host_real_attribs[pld.m_lev][std::make_pair(pld.m_grid, pld.m_tile)][n-NStructReal].push_back(rdata);
391 }
392 }
393
394 nparticles.pop_back();
396 if (nreals.size() > extradata - NStructReal) {
397 for (int n = NStructReal; n < extradata; n++)
399 nreals[n-NStructReal].pop_back();
401 }
402 }
403 }
404
405 for (int lev = 0; lev < std::ssize(host_particles); ++lev)
406 {
407 for (auto& kv : host_particles[lev])
408 {
409 auto grid = kv.first.first;
410 auto tile = kv.first.second;
411 const auto& src_tile = kv.second;
412
413 auto& dst_tile = GetParticles(lev)[std::make_pair(grid,tile)];
414 auto old_size = dst_tile.GetArrayOfStructs().size();
415 auto new_size = old_size + src_tile.size();
416 dst_tile.resize(new_size);
417
418 Gpu::copyAsync(Gpu::hostToDevice, src_tile.begin(), src_tile.end(),
419 dst_tile.GetArrayOfStructs().begin() + old_size);
420
421 for (int i = 0; i < NArrayReal; ++i) {
423 host_real_attribs[lev][std::make_pair(grid,tile)][i].begin(),
424 host_real_attribs[lev][std::make_pair(grid,tile)][i].end(),
425 dst_tile.GetStructOfArrays().GetRealData(i).begin() + old_size);
426 }
427 }
428 }
430
431 Redistribute();
432
433 }
434
435 if (m_verbose > 0)
436 {
437 const int IOProcNumber = ParallelDescriptor::IOProcessorNumber();
439 Long num_particles = how_many;
440
441 ParallelDescriptor::ReduceLongSum(num_particles, IOProcNumber);
442
443 if (AMREX_D_TERM(lNrep[0] == 1, && lNrep[1] == 1, && lNrep[2] == 1))
444 {
445 amrex::Print() << "Total number of particles: " << num_particles << '\n';
446 }
447 else
448 {
449 Long num_particles_read = how_many_read;
450
451 ParallelDescriptor::ReduceLongSum(num_particles_read, IOProcNumber);
452
453 amrex::Print() << "Replication the domain with vector "
454 << AMREX_D_TERM(lNrep[0] << " ", << lNrep[1] << " ", << lNrep[2]) << "\n"
455 << "Total number of particles read in : " << num_particles_read << '\n'
456 << "Total number of particles after replication: " << num_particles << '\n';
457 }
458 }
459
460 AMREX_ASSERT(OK());
461
462 if (m_verbose > 1)
463 {
464 ByteSpread();
465
466 auto runtime = amrex::second() - strttime;
467
469
470 amrex::Print() << "InitFromAsciiFile() time: " << runtime << '\n';
471 }
472}
473
474//
475// The format of a binary particle init file:
476//
477// NP -- The number of particles in the file. A "Long".
478// DM -- Our dimension. Either 1, 2, or 3. A "int".
479// NX -- The amount of "extra" data. A "int".
480// NP*(DM+NX) native floating-point numbers. A "float" or "double".
481//
482// Note that there is nothing separating all these values.
483// They're packed into the binary file like sardines.
484//
485template <typename ParticleType, int NArrayReal, int NArrayInt,
486 template<class> class Allocator, class CellAssignor>
487void
489InitFromBinaryFile (const std::string& file,
490 int extradata)
491{
492 BL_PROFILE("ParticleContainer<NSR, NSI, NAR, NAI>::InitFromBinaryFile()");
493 AMREX_ASSERT(!file.empty());
494 AMREX_ASSERT(extradata <= NStructReal);
495
496 const int MyProc = ParallelDescriptor::MyProc();
497 const int NProcs = ParallelDescriptor::NProcs();
498 const int IOProc = ParallelDescriptor::IOProcessorNumber();
499 const auto strttime = amrex::second();
500 //
501 // The number of MPI processes that read from the file.
502 // We limit this to a rather small number since there's a limit
503 // to the number of independent I/O channels on most filesystems.
504 //
505 const int NReaders = MaxReaders();
506
507 AMREX_ASSERT(NReaders <= NProcs);
508 //
509 // How many particles each NReaders reads before redistributing.
510 //
511 const Long NPartPerRedist = MaxParticlesPerRead();
512
513 if (m_verbose > 0)
514 {
515 amrex::Print() << "Reading with " << NReaders << " readers\n"
516 << "Redistributing after every " << NPartPerRedist << " particles for each reader\n";
517 }
518 //
519 // tmp_particles should mirror how m_particles is built.
520 // At the end of this routine it'll be the new m_particles.
521 //
522 Vector<ParticleLevel> tmp_particles;
523
524 tmp_particles.reserve(15); // So we don't ever have to do any copying on a resize.
525 tmp_particles.resize(finestLevel()+1);
526
527 resizeData();
528
529 //
530 // All the processors need to participate in Redistribute() though.
531 //
532 int NX = 0;
533 Long NP = 0;
534
536
537 std::ifstream ifs;
538
539 //
540 // The "set" of MPI processor numbers of the readers. We want
541 // to be able to easily check whether or not we're a reader.
542 //
543 std::set<int> readers;
544 //
545 // The same set but in ascending order.
546 //
547 Vector<int> rprocs(NReaders);
548
549 if (NReaders == NProcs)
550 {
551 //
552 // Just set'm.
553 //
554 for (int i = 0; i < NProcs; i++) {
555 rprocs[i] = i;
556 }
557 }
558 else
559 {
560 //
561 // The I/O Proc builds a set of NReader integers in the range: [0,NProcs-1].
562 //
563 // It then broadcast'm to all MPI procs.
564 //
565 // We want these to be as evenly distributed over the full set of
566 // [0,NProcs-1] MPI processors as possible, so that when reading we
567 // minimize the number of readers per Node, and hence can use more
568 // of the available Node memory for reading.
569 //
571 {
572 do
573 {
574 auto n = int(amrex::Random() * Real(NProcs-1));
575
576 AMREX_ASSERT(n >= 0);
577 AMREX_ASSERT(n < NProcs);
578
579 readers.insert(n);
580 }
581 while (std::ssize(readers) < NReaders);
582
583 AMREX_ASSERT(std::ssize(readers) == rprocs.size());
584
585 int i = 0;
586
587 for (auto it = readers.cbegin(), End = readers.cend();
588 it != End;
589 ++it, ++i)
590 {
591 rprocs[i] = *it;
592 }
593 }
594
595 ParallelDescriptor::Bcast(rprocs.dataPtr(), rprocs.size(), IOProc);
596 }
597
598 if (readers.empty())
599 {
600 //
601 // Set readers for non I/O procs.
602 //
603 readers.insert(rprocs.begin(), rprocs.end());
604
605 AMREX_ASSERT(std::ssize(readers) == rprocs.size());
606
607 AMREX_ASSERT(rprocs.size() == NReaders);
608 }
609
610 int RealSizeInFile = 0;
611
612 if (readers.contains(MyProc))
613 {
614 int DM = 0;
615
616 ifs.rdbuf()->pubsetbuf(io_buffer.dataPtr(), io_buffer.size());
617
618 ifs.open(file.c_str(), std::ios::in|std::ios::binary);
619
620 if (!ifs.good()) {
622 }
623
624 ifs.read((char*)&NP, sizeof(NP));
625 ifs.read((char*)&DM, sizeof(DM));
626 ifs.read((char*)&NX, sizeof(NX));
627 //
628 // NP MUST be positive!
629 //
630 if (NP <= 0) {
631 amrex::Abort("ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt>::InitFromBinaryFile(): NP <= 0");
632 }
633 //
634 // DM must equal AMREX_SPACEDIM.
635 //
636 if (DM != AMREX_SPACEDIM) {
637 amrex::Abort("ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt>::InitFromBinaryFile(): DM != AMREX_SPACEDIM");
638 }
639 //
640 // NX MUST be in [0,N].
641 //
642 if (NX < 0 || NX > NStructReal) {
643 amrex::Abort("ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt>::InitFromBinaryFile(): NX < 0 || NX > N");
644 }
645 //
646 // Can't ask for more data than exists in the file!
647 //
648 if (extradata > NX) {
649 amrex::Abort("ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt>::InitFromBinaryFile(): extradata > NX");
650 }
651 //
652 // Figure out whether we're dealing with floats or doubles.
653 //
654 // First get our current position.
655 //
656 const std::streamoff CURPOS = ifs.tellg();
657 //
658 // Seek to end of file.
659 //
660 ifs.seekg(0,std::ios::end);
661 //
662 // ENDPOS - CURPOS should bracket the particle data.
663 //
664 const std::streamoff ENDPOS = ifs.tellg();
665
666 RealSizeInFile = int((ENDPOS - CURPOS) / (NP*(DM+NX)));
667
668 AMREX_ASSERT(RealSizeInFile == sizeof(float) || RealSizeInFile == sizeof(double));
669 //
670 // Now set stream back to earlier position.
671 //
672 ifs.seekg(CURPOS, std::ios::beg);
673 //
674 // Skip to our place in the file.
675 //
676 int id = 0;
677 for ( ; id < NReaders; id++) {
678 if (rprocs[id] == MyProc) {
679 break;
680 }
681 }
682
683 AMREX_ASSERT(id >= 0 && id < NReaders);
684
685 const std::streamoff NSKIP = id * (NP/NReaders) * (DM+NX) * RealSizeInFile;
686
687 if (NSKIP > 0)
688 {
689 ifs.seekg(NSKIP, std::ios::cur);
690 }
691
692 if (!ifs.good())
693 {
694 std::string msg("ParticleContainer::InitFromBinaryFile(");
695 msg += file;
696 msg += ") failed @ 1";
697 amrex::Error(msg.c_str());
698 }
699 }
700 //
701 // Everyone needs to know NP -- the number of particles in the file.
702 //
704 //
705 // How many particles each reader gets to read.
706 //
707 Long MyCnt = NP / NReaders;
708
709 if (MyProc == rprocs[NReaders-1]) {
710 //
711 // Give any remainder to the last reader.
712 //
713 MyCnt += NP % NReaders;
714 }
715
716 Long how_many_redists = NP / (NPartPerRedist*NReaders), how_many_read = 0;
717
718 if (NP % (NPartPerRedist*NReaders)) { how_many_redists++; }
719
720 Vector<float> fxtra, fignore;
721 Vector<double> dxtra, dignore;
722
723 if (extradata > 0)
724 {
725 fxtra.resize(extradata);
726 dxtra.resize(extradata);
727 }
728
729 if ((NX-extradata) > 0)
730 {
731 fignore.resize(NX-extradata);
732 dignore.resize(NX-extradata);
733 }
734
735 ParticleLocData pld;
736
737 for (int j = 0; j < how_many_redists; j++)
738 {
739
741 host_particles.reserve(15);
742 host_particles.resize(finestLevel()+1);
743
744 if (readers.contains(MyProc))
745 {
746 ParticleType p;
747
748 AMREX_ASSERT(MyCnt > how_many_read);
749
750 const Long NRead = std::min((MyCnt-how_many_read), NPartPerRedist);
751
752 for (Long i = 0; i < NRead; i++)
753 {
754 //
755 // We don't read in idata.id or idata.cpu. We'll set those later
756 // in a manner to guarantee the global uniqueness of the pair.
757 //
758 if (RealSizeInFile == sizeof(float))
759 {
760 float fpos[AMREX_SPACEDIM];
761
762 ifs.read((char*)&fpos[0], AMREX_SPACEDIM*sizeof(float));
763
764 AMREX_D_TERM(p.pos(0) = fpos[0];,
765 p.pos(1) = fpos[1];,
766 p.pos(2) = fpos[2];);
767
768 }
769 else if (RealSizeInFile == sizeof(double))
770 {
771 double dpos[AMREX_SPACEDIM];
772
773 ifs.read((char*)&dpos[0], AMREX_SPACEDIM*sizeof(double));
774
775 AMREX_D_TERM(p.pos(0) = static_cast<ParticleReal>(dpos[0]);,
776 p.pos(1) = static_cast<ParticleReal>(dpos[1]);,
777 p.pos(2) = static_cast<ParticleReal>(dpos[2]););
778 }
779
780 //
781 // Read in any "extradata".
782 //
783 if (extradata > 0)
784 {
785 if (RealSizeInFile == sizeof(float))
786 {
787 ifs.read((char*)fxtra.data(), std::streamsize(extradata*sizeof(float)));
788
789 for (int ii = 0; ii < extradata; ii++) {
790 p.rdata(ii) = static_cast<ParticleReal>(fxtra[ii]);
791 }
792 }
793 else if (RealSizeInFile == sizeof(double))
794 {
795 ifs.read((char*)dxtra.data(), std::streamsize(extradata*sizeof(double)));
796
797 for (int ii = 0; ii < extradata; ii++) {
798 p.rdata(ii) = static_cast<ParticleReal>(dxtra[ii]);
799 }
800 }
801 }
802 //
803 // Read any remaining data for this particle.
804 //
805 if ((NX-extradata) > 0)
806 {
807 if (RealSizeInFile == sizeof(float))
808 {
809 ifs.read((char*)fignore.data(), std::streamsize((NX-extradata)*sizeof(float)));
810 }
811 else if (RealSizeInFile == sizeof(double))
812 {
813 ifs.read((char*)dignore.data(), std::streamsize((NX-extradata)*sizeof(double)));
814 }
815 }
816
817 if (!ifs.good())
818 {
819 std::string msg("ParticleContainer::InitFromBinaryFile(");
820 msg += file;
821 msg += ") failed @ 2";
822 amrex::Error(msg.c_str());
823 }
824
825 if (!Where(p, pld))
826 {
827 PeriodicShift(p);
828
829 if (!Where(p, pld))
830 {
831 if (m_verbose) {
832 amrex::AllPrint() << "BAD PARTICLE ID WOULD BE " << ParticleType::NextID() << '\n'
833 << "BAD PARTICLE POS "
834 << AMREX_D_TERM( p.pos(0),
835 << p.pos(1),
836 << p.pos(2))
837 << "\n";
838 }
839 amrex::Abort("ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt>::InitFromBinaryFile(): invalid particle");
840 }
841 }
842
843 p.id() = ParticleType::NextID();
844 p.cpu() = MyProc;
845
846 host_particles[pld.m_lev][std::make_pair(pld.m_grid, pld.m_tile)].push_back(p);
847 }
848
849 how_many_read += NRead;
850 }
851
852 for (int host_lev = 0; host_lev < std::ssize(host_particles); ++host_lev)
853 {
854 for (auto& kv : host_particles[host_lev]) {
855 auto grid = kv.first.first;
856 auto tile = kv.first.second;
857 const auto& src_tile = kv.second;
858
859 auto& dst_tile = GetParticles(host_lev)[std::make_pair(grid,tile)];
860 auto old_size = dst_tile.GetArrayOfStructs().size();
861 auto new_size = old_size + src_tile.size();
862 dst_tile.resize(new_size);
863
864 Gpu::copyAsync(Gpu::hostToDevice, src_tile.begin(), src_tile.end(),
865 dst_tile.GetArrayOfStructs().begin() + old_size);
866 }
867 }
869
870 Redistribute();
871
872 //
873 // Move particles in m_particles into tmp_particles so that
874 // we don't keep trying to redistribute particles that have
875 // already been redistributed correctly.
876 //
877 for (int lev = 0; lev < std::ssize(m_particles); lev++)
878 {
879 auto& pmap = m_particles[lev];
880 auto& tmp_pmap = tmp_particles[lev];
881
882 for (auto& kv : pmap) {
883 auto& aos = kv.second.GetArrayOfStructs()();
884 auto& tmp_aos = tmp_pmap[kv.first].GetArrayOfStructs()();
885
886 tmp_aos.insert(tmp_aos.end(), aos.begin(), aos.end());
887 ParticleVector().swap(aos);
888 }
889
890 ParticleLevel().swap(pmap);
891 }
892 }
893 //
894 // Make tmp_particles the new m_particles.
895 //
896 tmp_particles.swap(m_particles);
897 //
898 // Add up all the particles read in to get the total number of particles.
899 //
900 if (m_verbose > 0)
901 {
902 Long num_particles_read = how_many_read;
903
905
906 amrex::Print() << "\nTotal number of particles: " << num_particles_read << '\n';
907 }
908
909 AMREX_ASSERT(OK());
910
911 if (m_verbose > 1)
912 {
913 ByteSpread();
914
915 auto runtime = amrex::second() - strttime;
916
918
919 amrex::Print() << "InitFromBinaryFile() time: " << runtime << '\n';
920 }
921
923}
924
925//
926// This function expects to read a file containing the pathnames of
927// binary particles files needing to be read in for input. It expects
928// one file name per line.
929//
930
931template <typename ParticleType, int NArrayReal, int NArrayInt,
932 template<class> class Allocator, class CellAssignor>
933void
935InitFromBinaryMetaFile (const std::string& metafile,
936 int extradata)
937{
938 BL_PROFILE("ParticleContainer<NSR, NSI, NAR, NAI>::InitFromBinaryMetaFile()");
939 const auto strttime = amrex::second();
940
941 std::ifstream ifs(metafile.c_str(), std::ios::in);
942
943 std::string file;
944
945 for (;;)
946 {
947 std::getline(ifs,file);
948
949 if (!ifs.good()) { break; }
950
951 if (m_verbose > 1) {
952 amrex::Print() << "InitFromBinaryMetaFile: processing file: " << file << '\n';
953 }
954
955 InitFromBinaryFile(file, extradata);
956 }
957
958 if (m_verbose > 1)
959 {
960 ByteSpread();
961
962 auto runtime = amrex::second() - strttime;
963
965
966 amrex::Print() << "InitFromBinaryMetaFile() time: " << runtime << '\n';
967 }
968}
969
970template <typename ParticleType, int NArrayReal, int NArrayInt,
971 template<class> class Allocator, class CellAssignor>
972void
974InitRandom (Long icount,
975 ULong iseed,
976 const ParticleInitData& pdata,
977 bool serialize,
978 RealBox containing_bx)
979{
980 BL_PROFILE("ParticleContainer<NSR, NSI, NAR, NAI>::InitRandom()");
981 AMREX_ASSERT(iseed > 0);
982 AMREX_ASSERT(icount > 0);
983
984 AMREX_ASSERT(m_gdb != nullptr);
985
986 const int MyProc = ParallelDescriptor::MyProc();
987 const int NProcs = ParallelDescriptor::NProcs();
988 const int IOProc = ParallelDescriptor::IOProcessorNumber();
989 const auto strttime = amrex::second();
990 const Geometry& geom = Geom(0);
991
992 Real r, x, len[AMREX_SPACEDIM] = { AMREX_D_DECL(geom.ProbLength(0),
993 geom.ProbLength(1),
994 geom.ProbLength(2)) };
995
996 // We will enforce that the particles are within the containing_bx.
997 // If containing_bx is not passed in, it defaults to the full domain.
998 if (!containing_bx.ok()) { containing_bx = geom.ProbDomain(); }
999
1000 // containing_bx is assumed to lie within the domain.
1001 if (!geom.ProbDomain().contains(containing_bx))
1002 {
1003 containing_bx.setLo(geom.ProbLo());
1004 containing_bx.setHi(geom.ProbHi());
1005 }
1006
1007 const Real* xlo = containing_bx.lo();
1008 const Real* xhi = containing_bx.hi();
1009
1010 amrex::InitRandom(iseed+MyProc);
1011
1012 if (serialize)
1013 {
1014 if(icount*AMREX_SPACEDIM >= std::numeric_limits<int>::max())
1015 {
1017 "InitRandom has serialize=true, but this would cause too much "
1018 "particle data to be sent from IOProc. Set serialize=false, "
1019 "or use fewer than " +
1020 std::to_string(
1021 amrex::Math::ceil(
1022 amrex::Real(std::numeric_limits<int>::max()) /
1023 amrex::Real(AMREX_SPACEDIM)
1024 )
1025 ) + " particles."
1026 );
1027 }
1028 //
1029 // We'll let IOProc generate the particles so we get the same
1030 // positions no matter how many CPUs we have. This is here
1031 // mainly for debugging purposes. It's not really useful for
1032 // very large numbers of particles.
1033 //
1034 //
1035 Vector<typename ParticleType::RealType> pos(icount*AMREX_SPACEDIM);
1036
1038 {
1039 for (Long j = 0; j < icount; j++)
1040 {
1041 for (int i = 0; i < AMREX_SPACEDIM; i++)
1042 {
1043 do
1044 {
1045 r = amrex::Random();
1046 x = geom.ProbLo(i) + (r * len[i]);
1047 }
1048 while (static_cast<ParticleReal>(x) < static_cast<ParticleReal>(xlo[i]) || static_cast<ParticleReal>(x) >= static_cast<ParticleReal>(xhi[i]));
1049
1050 pos[j*AMREX_SPACEDIM + i] = static_cast<ParticleReal>(x);
1051 }
1052 }
1053 }
1054
1055 ParallelDescriptor::Bcast(pos.dataPtr(), icount*AMREX_SPACEDIM, IOProc);
1056
1057 ParticleLocData pld;
1058
1060 host_particles.reserve(15);
1061 host_particles.resize(finestLevel()+1);
1062
1063 Vector<std::map<std::pair<int, int>, std::array<Gpu::HostVector<ParticleReal>, NArrayReal > > > host_real_attribs;
1064 host_real_attribs.reserve(15);
1065 host_real_attribs.resize(finestLevel()+1);
1066
1067 Vector<std::map<std::pair<int, int>, std::array<Gpu::HostVector<int>, NArrayInt > > > host_int_attribs;
1068 host_int_attribs.reserve(15);
1069 host_int_attribs.resize(finestLevel()+1);
1070
1072 host_idcpu.reserve(15);
1073 host_idcpu.resize(finestLevel()+1);
1074
1075 for (Long j = 0; j < icount; j++)
1076 {
1077 Particle<0, 0> ptest;
1078
1079 for (int i = 0; i < AMREX_SPACEDIM; i++) {
1080 ptest.pos(i) = pos[j*AMREX_SPACEDIM + i];
1081 }
1082
1083 if (!Where(ptest, pld)) {
1084 amrex::Abort("ParticleContainer::InitRandom(): invalid particle");
1085 }
1086
1087 AMREX_ASSERT(pld.m_lev >= 0 && pld.m_lev <= finestLevel());
1088 std::pair<int, int> ind(pld.m_grid, pld.m_tile);
1089
1090 const int who = ParticleDistributionMap(pld.m_lev)[pld.m_grid];
1091
1092 if (who == MyProc) {
1093
1094 if constexpr(!ParticleType::is_soa_particle)
1095 {
1096 ParticleType p;
1097 for (int i = 0; i < AMREX_SPACEDIM; i++) {
1098 p.pos(i) = pos[j*AMREX_SPACEDIM + i];
1099 }
1100
1101 // We own it. Add it at the appropriate level.
1102 p.id() = ParticleType::NextID();
1103 p.cpu() = MyProc;
1104
1105 for (int i = 0; i < NStructInt; i++) {
1106 p.idata(i) = pdata.int_struct_data[i];
1107 }
1108
1109 for (int i = 0; i < NStructReal; i++) {
1110 p.rdata(i) = static_cast<ParticleReal>(pdata.real_struct_data[i]);
1111 }
1112
1113 // add the struct
1114 host_particles[pld.m_lev][ind].push_back(p);
1115
1116 // add the real...
1117 for (int i = 0; i < NArrayReal; i++) {
1118 host_real_attribs[pld.m_lev][ind][i].push_back(static_cast<ParticleReal>(pdata.real_array_data[i]));
1119 }
1120
1121 // ... and int array data
1122 for (int i = 0; i < NArrayInt; i++) {
1123 host_int_attribs[pld.m_lev][ind][i].push_back(pdata.int_array_data[i]);
1124 }
1125 } else {
1126 for (int i = 0; i < AMREX_SPACEDIM; i++) {
1127 host_real_attribs[pld.m_lev][ind][i].push_back(pos[j*AMREX_SPACEDIM+i]);
1128 }
1129
1130 host_idcpu[pld.m_lev][ind].push_back(0);
1131 ParticleIDWrapper(host_idcpu[pld.m_lev][ind].back()) = ParticleType::NextID();
1132 ParticleCPUWrapper(host_idcpu[pld.m_lev][ind].back()) = ParallelDescriptor::MyProc();
1133
1134 host_particles[pld.m_lev][ind];
1135
1136 // add the real...
1137 for (int i = AMREX_SPACEDIM; i < NArrayReal; i++) {
1138 host_real_attribs[pld.m_lev][ind][i].push_back(static_cast<ParticleReal>(pdata.real_array_data[i]));
1139 }
1140
1141 // ... and int array data
1142 for (int i = 2; i < NArrayInt; i++) {
1143 host_int_attribs[pld.m_lev][ind][i].push_back(pdata.int_array_data[i]);
1144 }
1145 }
1146 }
1147 }
1148
1149 for (int host_lev = 0; host_lev < std::ssize(host_particles); ++host_lev)
1150 {
1151 for (auto& kv : host_particles[host_lev]) {
1152 auto grid = kv.first.first;
1153 auto tile = kv.first.second;
1154 const auto& src_tile = kv.second;
1155
1156 auto& dst_tile = GetParticles(host_lev)[std::make_pair(grid,tile)];
1157 auto old_size = dst_tile.size();
1158 auto new_size = old_size;
1159 if constexpr(!ParticleType::is_soa_particle)
1160 {
1161 new_size += src_tile.size();
1162 } else {
1163 new_size += host_real_attribs[host_lev][std::make_pair(grid,tile)][0].size();
1164 }
1165 dst_tile.resize(new_size);
1166
1167 if constexpr(!ParticleType::is_soa_particle)
1168 {
1169 Gpu::copyAsync(Gpu::hostToDevice, src_tile.begin(), src_tile.end(),
1170 dst_tile.GetArrayOfStructs().begin() + old_size);
1171 } else {
1173 host_idcpu[host_lev][std::make_pair(grid,tile)].begin(),
1174 host_idcpu[host_lev][std::make_pair(grid,tile)].end(),
1175 dst_tile.GetStructOfArrays().GetIdCPUData().begin() + old_size);
1176 }
1177
1178 for (int i = 0; i < NArrayReal; ++i) { // NOLINT(readability-misleading-indentation)
1180 host_real_attribs[host_lev][std::make_pair(grid,tile)][i].begin(),
1181 host_real_attribs[host_lev][std::make_pair(grid,tile)][i].end(),
1182 dst_tile.GetStructOfArrays().GetRealData(i).begin() + old_size);
1183 }
1184
1185 for (int i = 0; i < NArrayInt; ++i) {
1187 host_int_attribs[host_lev][std::make_pair(grid,tile)][i].begin(),
1188 host_int_attribs[host_lev][std::make_pair(grid,tile)][i].end(),
1189 dst_tile.GetStructOfArrays().GetIntData(i).begin() + old_size);
1190 }
1191 }
1192 }
1194
1195 AMREX_ASSERT(OK());
1196 }
1197 else {
1198 // We'll generate the particles in parallel.
1199 // Each CPU will key off the given seed to get independent streams of random numbers.
1200 Long M = icount / NProcs;
1201 // Processor 0 will get the slop.
1202 if (MyProc == 0) {
1203 M += (icount % NProcs);
1204 }
1205
1206 ParticleLocData pld;
1207
1209 host_particles.reserve(15);
1210 host_particles.resize(finestLevel()+1);
1211
1212 Vector<std::map<std::pair<int, int>, std::array<Gpu::HostVector<ParticleReal>, NArrayReal > > > host_real_attribs;
1213 host_real_attribs.reserve(15);
1214 host_real_attribs.resize(finestLevel()+1);
1215
1216 Vector<std::map<std::pair<int, int>, std::array<Gpu::HostVector<int>, NArrayInt > > > host_int_attribs;
1217 host_int_attribs.reserve(15);
1218 host_int_attribs.resize(finestLevel()+1);
1219
1221 host_idcpu.reserve(15);
1222 host_idcpu.resize(finestLevel()+1);
1223
1224 for (Long icnt = 0; icnt < M; icnt++) {
1225 Particle<0, 0> ptest;
1226 for (int i = 0; i < AMREX_SPACEDIM; i++) {
1227 do {
1228 r = amrex::Random();
1229 x = geom.ProbLo(i) + (r * len[i]);
1230 }
1231 while (static_cast<ParticleReal>(x) < static_cast<ParticleReal>(xlo[i]) || static_cast<ParticleReal>(x) >= static_cast<ParticleReal>(xhi[i]));
1232
1233 ptest.pos(i) = static_cast<ParticleReal>(x);
1234
1235 AMREX_ASSERT(ptest.pos(i) < geom.ProbHi(i));
1236 }
1237
1238 ptest.id() = ParticleType::NextID();
1239 ptest.cpu() = ParallelDescriptor::MyProc();
1240
1241 // locate the particle
1242 if (!Where(ptest, pld))
1243 {
1244 amrex::Abort("ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt>::InitRandom(): invalid particle");
1245 }
1246 AMREX_ASSERT(pld.m_lev >= 0 && pld.m_lev <= finestLevel());
1247 std::pair<int, int> ind(pld.m_grid, pld.m_tile);
1248
1249 if constexpr(!ParticleType::is_soa_particle)
1250 {
1251 ParticleType p;
1252 for (int i = 0; i < AMREX_SPACEDIM; i++) {
1253 p.pos(i) = ptest.pos(i);;
1254 }
1255
1256 p.id() = ptest.id();
1257 p.cpu() = ptest.cpu();
1258
1259 for (int i = 0; i < NStructReal; i++) {
1260 p.rdata(i) = static_cast<ParticleReal>(pdata.real_struct_data[i]);
1261 }
1262
1263 for (int i = 0; i < NStructInt; i++) {
1264 p.idata(i) = pdata.int_struct_data[i];
1265 }
1266
1267 // add the struct
1268 host_particles[pld.m_lev][ind].push_back(p);
1269
1270 // add the real...
1271 for (int i = 0; i < NArrayReal; i++) {
1272 host_real_attribs[pld.m_lev][ind][i].push_back(static_cast<ParticleReal>(pdata.real_array_data[i]));
1273 }
1274
1275 // ... and int array data
1276 for (int i = 0; i < NArrayInt; i++) {
1277 host_int_attribs[pld.m_lev][ind][i].push_back(pdata.int_array_data[i]);
1278 }
1279 } else {
1280 for (int i = 0; i < AMREX_SPACEDIM; i++) {
1281 host_real_attribs[pld.m_lev][ind][i].push_back(ptest.pos(i));
1282 }
1283
1284 host_idcpu[pld.m_lev][ind].push_back(0);
1285 ParticleIDWrapper(host_idcpu[pld.m_lev][ind].back()) = ParticleType::NextID();
1286 ParticleCPUWrapper(host_idcpu[pld.m_lev][ind].back()) = ParallelDescriptor::MyProc();
1287
1288 host_particles[pld.m_lev][ind];
1289
1290 // add the real...
1291 for (int i = AMREX_SPACEDIM; i < NArrayReal; i++) {
1292 host_real_attribs[pld.m_lev][ind][i].push_back(static_cast<ParticleReal>(pdata.real_array_data[i]));
1293 }
1294
1295 // ... and int array data
1296 for (int i = 2; i < NArrayInt; i++) {
1297 host_int_attribs[pld.m_lev][ind][i].push_back(pdata.int_array_data[i]);
1298 }
1299 }
1300 }
1301
1302 for (int host_lev = 0; host_lev < std::ssize(host_particles); ++host_lev)
1303 {
1304 for (auto& kv : host_particles[host_lev]) {
1305 auto grid = kv.first.first;
1306 auto tile = kv.first.second;
1307 const auto& src_tile = kv.second;
1308
1309 auto& dst_tile = GetParticles(host_lev)[std::make_pair(grid,tile)];
1310 auto old_size = dst_tile.size();
1311 auto new_size = old_size;
1312 if constexpr(!ParticleType::is_soa_particle)
1313 {
1314 new_size += src_tile.size();
1315 } else {
1316 new_size += host_real_attribs[host_lev][std::make_pair(grid,tile)][0].size();
1317 }
1318 dst_tile.resize(new_size);
1319
1320 if constexpr(!ParticleType::is_soa_particle)
1321 {
1322 Gpu::copyAsync(Gpu::hostToDevice, src_tile.begin(), src_tile.end(),
1323 dst_tile.GetArrayOfStructs().begin() + old_size);
1324 } else {
1326 host_idcpu[host_lev][std::make_pair(grid,tile)].begin(),
1327 host_idcpu[host_lev][std::make_pair(grid,tile)].end(),
1328 dst_tile.GetStructOfArrays().GetIdCPUData().begin() + old_size);
1329 }
1330
1331 for (int i = 0; i < NArrayReal; ++i) { // NOLINT(readability-misleading-indentation)
1333 host_real_attribs[host_lev][std::make_pair(grid,tile)][i].begin(),
1334 host_real_attribs[host_lev][std::make_pair(grid,tile)][i].end(),
1335 dst_tile.GetStructOfArrays().GetRealData(i).begin() + old_size);
1336 }
1337
1338 for (int i = 0; i < NArrayInt; ++i) {
1340 host_int_attribs[host_lev][std::make_pair(grid,tile)][i].begin(),
1341 host_int_attribs[host_lev][std::make_pair(grid,tile)][i].end(),
1342 dst_tile.GetStructOfArrays().GetIntData(i).begin() + old_size);
1343 }
1344 }
1345 }
1347
1348 // Let Redistribute() sort out where the particles belong.
1349 Redistribute();
1350 }
1351
1352 if (m_verbose > 1)
1353 {
1354 auto stoptime = amrex::second() - strttime;
1355
1356 ParallelDescriptor::ReduceRealMax(stoptime,IOProc);
1357
1358 amrex::Print() << "ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt>::InitRandom() time: " << stoptime << '\n';
1359 }
1360
1362}
1363
1364template <typename ParticleType, int NArrayReal, int NArrayInt,
1365 template<class> class Allocator, class CellAssignor>
1366void
1368::InitRandomPerBox (Long icount_per_box,
1369 ULong iseed,
1370 const ParticleInitData& pdata)
1371{
1372 BL_PROFILE("ParticleContainer<NSR, NSI, NAR, NAI>::InitRandomPerBox()");
1373 AMREX_ASSERT(iseed > 0);
1374 AMREX_ASSERT(icount_per_box > 0);
1375
1376 AMREX_ASSERT(m_gdb != nullptr);
1377
1378 const int IOProc = ParallelDescriptor::IOProcessorNumber();
1379 const auto strttime = amrex::second();
1380 const Geometry& geom = Geom(0);
1381
1382 ParticleLocData pld;
1383 ParticleType p;
1384
1385 // This assumes level 0 since geom = m_gdb->Geom(0)
1386 const Real* dx = geom.CellSize();
1387
1388 // We use exactly the same seed for every grid
1389 std::mt19937 mt(iseed);
1390 std::uniform_real_distribution<double> dist(0.0, 1.0);
1391
1392 m_particles.resize(m_gdb->finestLevel()+1);
1393
1394 for (int lev = 0; lev < std::ssize(m_particles); lev++)
1395 {
1396 AMREX_ASSERT(m_particles[lev].empty());
1397 }
1398
1399 // We'll generate the particles in parallel -- but no tiling here.
1400 for (MFIter mfi(*m_dummy_mf[0], false); mfi.isValid(); ++mfi)
1401 {
1402 Box grid = m_gdb->ParticleBoxArray(0)[mfi.index()];
1403 RealBox grid_box = RealBox(grid,dx,geom.ProbLo());
1404
1405 for (Long icnt = 0; icnt < icount_per_box; icnt++) {
1406 for (Long jcnt = 0; jcnt < icount_per_box; jcnt++) {
1407 for (Long kcnt = 0; kcnt < icount_per_box; kcnt++)
1408 {
1410 p.pos(0) = static_cast<ParticleReal>(grid_box.lo(0) + (dist(mt) + double(icnt)) / double(icount_per_box) * grid_box.length(0));,
1411 p.pos(1) = static_cast<ParticleReal>(grid_box.lo(1) + (dist(mt) + double(jcnt)) / double(icount_per_box) * grid_box.length(1));,
1412 p.pos(2) = static_cast<ParticleReal>(grid_box.lo(2) + (dist(mt) + double(kcnt)) / double(icount_per_box) * grid_box.length(2));
1413 );
1414
1415 for (int i = 0; i < AMREX_SPACEDIM; i++) {
1416 AMREX_ASSERT(p.pos(i) < grid_box.hi(i));
1417 }
1418
1419 // the real struct data
1420 for (int i = 0; i < NStructReal; i++) {
1421 p.rdata(i) = static_cast<ParticleReal>(pdata.real_struct_data[i]);
1422 }
1423
1424 // the int struct data
1425 p.id() = ParticleType::NextID();
1426 p.cpu() = ParallelDescriptor::MyProc();
1427
1428 for (int i = 0; i < NStructInt; i++) {
1429 p.idata(i) = pdata.int_struct_data[i];
1430 }
1431
1432 // locate the particle
1433 if (!Where(p, pld)) {
1434 amrex::Abort("ParticleContainer::InitRandomPerBox(): invalid particle");
1435 }
1436 AMREX_ASSERT(pld.m_lev >= 0 && pld.m_lev <= finestLevel());
1437 std::pair<int, int> ind(pld.m_grid, pld.m_tile);
1438
1439 // add the struct
1440 m_particles[pld.m_lev][ind].push_back(p);
1441
1442 // add the real...
1443 for (int i = 0; i < NArrayReal; i++) {
1444 m_particles[pld.m_lev][ind].push_back_real(i, static_cast<ParticleReal>(pdata.real_array_data[i]));
1445 }
1446
1447 // ... and int array data
1448 for (int i = 0; i < NArrayInt; i++) {
1449 m_particles[pld.m_lev][ind].push_back_int(i, pdata.int_array_data[i]);
1450 }
1451
1452 } } }
1453 }
1454
1455 if (m_verbose > 1)
1456 {
1457 auto stoptime = amrex::second() - strttime;
1458
1459 ParallelDescriptor::ReduceRealMax(stoptime,IOProc);
1460
1461 amrex::Print() << "ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt>::InitRandomPerBox() time: " << stoptime << '\n';
1462 }
1463}
1464
1465template <typename ParticleType, int NArrayReal, int NArrayInt,
1466 template<class> class Allocator, class CellAssignor>
1467void
1469InitOnePerCell (Real x_off, Real y_off, Real z_off, const ParticleInitData& pdata)
1470{
1471 amrex::ignore_unused(y_off,z_off);
1472
1473 BL_PROFILE("ParticleContainer<NSR, NSI, NAR, NAI>::InitOnePerCell()");
1474
1475 AMREX_ASSERT(m_gdb != nullptr);
1476
1477 // Note that x_off, y_off and z_off are the offsets from the lower left corner
1478 // in each cell as measured in units of dx, so they must be in [0,1]
1479 AMREX_ASSERT(x_off >= 0. && y_off >= 0. && z_off >= 0.);
1480 AMREX_ASSERT(x_off <= 1. && y_off <= 1. && z_off <= 1.);
1481
1482 const int IOProc = ParallelDescriptor::IOProcessorNumber();
1483 const auto strttime = amrex::second();
1484 const Geometry& geom = Geom(0); // generates particles on level 0;
1485
1486 const Real* dx = geom.CellSize();
1487
1488 ParticleType p;
1489
1490 // We'll generate the particles in parallel -- but no tiling of the grid here.
1491 for (MFIter mfi(*m_dummy_mf[0], false); mfi.isValid(); ++mfi) {
1492 Box grid = ParticleBoxArray(0)[mfi.index()];
1493 auto ind = std::make_pair(mfi.index(), mfi.LocalTileIndex());
1494 RealBox grid_box (grid,dx,geom.ProbLo());
1496 for (IntVect beg = grid.smallEnd(), end=grid.bigEnd(), cell = grid.smallEnd(); cell <= end; grid.next(cell))
1497 {
1498 // the real struct data
1499 AMREX_D_TERM(p.pos(0) = static_cast<ParticleReal>(grid_box.lo(0) + (x_off + cell[0]-beg[0])*dx[0]);,
1500 p.pos(1) = static_cast<ParticleReal>(grid_box.lo(1) + (y_off + cell[1]-beg[1])*dx[1]);,
1501 p.pos(2) = static_cast<ParticleReal>(grid_box.lo(2) + (z_off + cell[2]-beg[2])*dx[2]););
1502
1503 for (int d = 0; d < AMREX_SPACEDIM; ++d) {
1504 AMREX_ASSERT(p.pos(d) < grid_box.hi(d));
1505 }
1506
1507 for (int i = 0; i < NStructReal; i++) {
1508 p.rdata(i) = static_cast<ParticleReal>(pdata.real_struct_data[i]);
1509 }
1510
1511 // the int struct data
1512 p.id() = ParticleType::NextID();
1513 p.cpu() = ParallelDescriptor::MyProc();
1514
1515 for (int i = 0; i < NStructInt; i++) {
1516 p.idata(i) = pdata.int_struct_data[i];
1517 }
1518
1519 // add the struct
1520 ptile_tmp.push_back(p);
1521
1522 // add the real...
1523 for (int i = 0; i < NArrayReal; i++) {
1524 ptile_tmp.push_back_real(i, static_cast<ParticleReal>(pdata.real_array_data[i]));
1525 }
1526
1527 // ... and int array data
1528 for (int i = 0; i < NArrayInt; i++) {
1529 ptile_tmp.push_back_int(i, pdata.int_array_data[i]);
1530 }
1531 }
1532
1533 m_particles[0][ind].resize(ptile_tmp.numParticles());
1534 amrex::copyParticles(m_particles[0][ind], ptile_tmp);
1536 }
1537
1538 Redistribute();
1539
1540 if (m_verbose > 1) {
1541 auto stoptime = amrex::second() - strttime;
1542
1543 ParallelDescriptor::ReduceRealMax(stoptime,IOProc);
1544
1545 amrex::Print() << "ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt>::InitOnePerCell() time: " << stoptime << '\n';
1546 }
1547}
1548
1549template <typename ParticleType, int NArrayReal, int NArrayInt,
1550 template<class> class Allocator, class CellAssignor>
1551void
1553InitNRandomPerCell (int n_per_cell, const ParticleInitData& pdata)
1554{
1555 BL_PROFILE("ParticleContainer<NSR, NSI, NAR, NAI>::InitNRandomPerCell()");
1556
1557 AMREX_ASSERT(m_gdb != nullptr);
1558
1559 const int IOProc = ParallelDescriptor::IOProcessorNumber();
1560 const auto strttime = amrex::second();
1561 const Geometry& geom = Geom(0);
1562
1563 // This assumes level 0 since geom = Geom(0)
1564 const Real* dx = geom.CellSize();
1565
1566 resizeData();
1567
1568 for (int lev = 0; lev < std::ssize(m_particles); lev++) {
1569 AMREX_ASSERT(m_particles[lev].empty());
1570 }
1571
1572 ParticleLocData pld;
1573 ParticleType p;
1574 Real r;
1575
1576 // We'll generate the particles in parallel -- but no tiling here.
1577 for (MFIter mfi(*m_dummy_mf[0], false); mfi.isValid(); ++mfi)
1578 {
1579 Box grid = ParticleBoxArray(0)[mfi.index()];
1580 RealBox grid_box (grid,dx,geom.ProbLo());
1581
1583 host_particles.reserve(15);
1584 host_particles.resize(finestLevel()+1);
1585
1586 Vector<std::map<std::pair<int, int>, std::array<Gpu::HostVector<ParticleReal>, NArrayReal > > > host_real_attribs;
1587 host_real_attribs.reserve(15);
1588 host_real_attribs.resize(finestLevel()+1);
1589
1590 Vector<std::map<std::pair<int, int>, std::array<Gpu::HostVector<int>, NArrayInt > > > host_int_attribs;
1591 host_int_attribs.reserve(15);
1592 host_int_attribs.resize(finestLevel()+1);
1593
1594 for (IntVect beg = grid.smallEnd(), end=grid.bigEnd(),
1595 cell = grid.smallEnd(); cell <= end; grid.next(cell)) {
1596
1597 for (int n = 0; n < n_per_cell; n++)
1598 {
1599 // the real struct data
1600 for (int i = 0; i < AMREX_SPACEDIM; i++) {
1601 constexpr int max_iter = 10;
1602 int iter = 0;
1603 while (iter < max_iter) {
1604 r = amrex::Random();
1605 p.pos(i) = static_cast<ParticleReal>(grid_box.lo(i) + (r + Real(cell[i]-beg[i]))*dx[i]);
1606 if (p.pos(i) < grid_box.hi(i)) { break; }
1607 iter++;
1608 }
1609 AMREX_ASSERT(p.pos(i) < grid_box.hi(i));
1610 }
1611
1612 for (int i = 0; i < NStructReal; i++) {
1613 p.rdata(i) = static_cast<ParticleReal>(pdata.real_struct_data[i]);
1614 }
1615
1616 // the int struct data
1617 p.id() = ParticleType::NextID();
1618 p.cpu() = ParallelDescriptor::MyProc();
1619
1620 for (int i = 0; i < NStructInt; i++) {
1621 p.idata(i) = pdata.int_struct_data[i];
1622 }
1623
1624 // locate the particle
1625 if (!Where(p, pld)) {
1626 amrex::Abort("ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt>::InitNRandomPerCell(): invalid particle");
1627 }
1628 AMREX_ASSERT(pld.m_lev >= 0 && pld.m_lev <= finestLevel());
1629 std::pair<int, int> ind(pld.m_grid, pld.m_tile);
1630
1631 // add the struct
1632 host_particles[pld.m_lev][ind].push_back(p);
1633
1634 // add the real...
1635 for (int i = 0; i < NArrayReal; i++) {
1636 host_real_attribs[pld.m_lev][ind][i].push_back(pdata.real_array_data[i]);
1637 }
1638
1639 // ... and int array data
1640 for (int i = 0; i < NArrayInt; i++) {
1641 host_int_attribs[pld.m_lev][ind][i].push_back(pdata.int_array_data[i]);
1642 }
1643 }
1644 }
1645
1646 for (int host_lev = 0; host_lev < std::ssize(host_particles); ++host_lev)
1647 {
1648 for (auto& kv : host_particles[host_lev]) {
1649 auto gid = kv.first.first;
1650 auto tid = kv.first.second;
1651 const auto& src_tid = kv.second;
1652
1653 auto& dst_tile = GetParticles(host_lev)[std::make_pair(gid,tid)];
1654 auto old_size = dst_tile.GetArrayOfStructs().size();
1655 auto new_size = old_size + src_tid.size();
1656 dst_tile.resize(new_size);
1657
1658 Gpu::copyAsync(Gpu::hostToDevice, src_tid.begin(), src_tid.end(),
1659 dst_tile.GetArrayOfStructs().begin() + old_size);
1660
1661 for (int i = 0; i < NArrayReal; ++i)
1662 {
1664 host_real_attribs[host_lev][std::make_pair(gid,tid)][i].begin(),
1665 host_real_attribs[host_lev][std::make_pair(gid,tid)][i].end(),
1666 dst_tile.GetStructOfArrays().GetRealData(i).begin() + old_size);
1667 }
1668
1669 for (int i = 0; i < NArrayInt; ++i)
1670 {
1672 host_int_attribs[host_lev][std::make_pair(gid,tid)][i].begin(),
1673 host_int_attribs[host_lev][std::make_pair(gid,tid)][i].end(),
1674 dst_tile.GetStructOfArrays().GetIntData(i).begin() + old_size);
1675 }
1676 }
1677 }
1679 }
1680
1681 if (m_verbose > 1)
1682 {
1683 auto stoptime = amrex::second() - strttime;
1684
1685 ParallelDescriptor::ReduceRealMax(stoptime,IOProc);
1686
1687 amrex::Print() << "ParticleContainer_impl<ParticleType, NArrayReal, NArrayInt>::InitNRandomPerCell() time: " << stoptime << '\n';
1688 }
1689
1691}
1692
1693}
1694
1695#endif /*AMREX_PARTICLEINIT_H*/
#define BL_PROFILE(a)
Definition AMReX_BLProfiler.H:551
#define AMREX_ASSERT(EX)
Definition AMReX_BLassert.H:38
#define AMREX_D_TERM(a, b, c)
Definition AMReX_SPACE.H:172
#define AMREX_D_DECL(a, b, c)
Definition AMReX_SPACE.H:171
Print on all processors of the default communicator.
Definition AMReX_Print.H:113
__host__ __device__ const IntVectND< dim > & bigEnd() const &noexcept
Return the inclusive upper bound of the box.
Definition AMReX_Box.H:123
__host__ __device__ void next(IntVectND< dim > &) const noexcept
Step through the rectangle. It is a runtime error to give a point not inside rectangle....
Definition AMReX_Box.H:1121
__host__ __device__ const IntVectND< dim > & smallEnd() const &noexcept
Return the inclusive lower bound of the box.
Definition AMReX_Box.H:111
const Real * CellSize() const noexcept
Returns the cellsize for each coordinate direction.
Definition AMReX_CoordSys.H:71
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:75
const Real * ProbHi() const noexcept
Returns the hi end of the problem domain in each dimension.
Definition AMReX_Geometry.H:186
const RealBox & ProbDomain() const noexcept
Returns the problem domain.
Definition AMReX_Geometry.H:176
Real ProbLength(int dir) const noexcept
Returns length of problem domain in specified dimension.
Definition AMReX_Geometry.H:214
const Real * ProbLo() const noexcept
Returns the lo end of the problem domain in each dimension.
Definition AMReX_Geometry.H:184
static void streamSynchronize() noexcept
Definition AMReX_GpuDevice.cpp:856
Iterator for looping ever tiles and boxes of amrex::FabArray based containers.
Definition AMReX_MFIter.H:88
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition AMReX_MFIter.H:172
Dynamically allocated vector for trivially copyable data.
Definition AMReX_PODVector.H:308
void pop_back() noexcept
Definition AMReX_PODVector.H:644
T & back() noexcept
Definition AMReX_PODVector.H:662
bool empty() const noexcept
Definition AMReX_PODVector.H:652
void push_back(const T &a_value)
Definition AMReX_PODVector.H:629
A distributed container for Particles sorted onto the levels, grids, and tiles of a block-structured ...
Definition AMReX_ParticleContainer.H:149
std::map< std::pair< int, int >, ParticleTileType > ParticleLevel
Definition AMReX_ParticleContainer.H:196
void InitNRandomPerCell(int n_per_cell, const ParticleInitData &pdata)
This initializes the particle container with n_per_cell randomly distributed particles per cell,...
Definition AMReX_ParticleInit.H:1553
typename AoS::ParticleVector ParticleVector
Definition AMReX_ParticleContainer.H:204
void InitOnePerCell(Real x_off, Real y_off, Real z_off, const ParticleInitData &pdata)
This initializes the particle container with one particle per cell, where the other particle data and...
Definition AMReX_ParticleInit.H:1469
void InitRandom(Long icount, ULong iseed, const ParticleInitData &pdata, bool serialize=false, RealBox bx=RealBox())
This initializes the particle container with icount randomly distributed particles....
Definition AMReX_ParticleInit.H:974
void InitFromBinaryMetaFile(const std::string &file, int extradata)
Definition AMReX_ParticleInit.H:935
void InitFromBinaryFile(const std::string &file, int extradata)
Definition AMReX_ParticleInit.H:489
T_ParticleType ParticleType
Definition AMReX_ParticleContainer.H:151
This class provides the user with a few print options.
Definition AMReX_Print.H:35
A Box with real dimensions.
Definition AMReX_RealBox.H:28
__host__ __device__ bool contains(const Real *point, Real eps=0.0) const noexcept
Is the specified point contained in the RealBox?
Definition AMReX_RealBox.H:109
void setLo(const Real *a_lo) noexcept
Sets lo side.
Definition AMReX_RealBox.H:71
__host__ __device__ const Real * lo() const &noexcept
Returns lo side.
Definition AMReX_RealBox.H:53
__host__ __device__ const Real * hi() const &noexcept
Returns hide side.
Definition AMReX_RealBox.H:58
__host__ __device__ Real length(int dir) const noexcept
Returns length in specified direction.
Definition AMReX_RealBox.H:69
void setHi(const Real *a_hi) noexcept
Sets hi side.
Definition AMReX_RealBox.H:79
__host__ __device__ bool ok() const noexcept
Is the RealBox OK; i.e. does it have non-negative volume?
Definition AMReX_RealBox.H:88
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
amrex_real Real
Floating Point Type for Fields.
Definition AMReX_REAL.H:79
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
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
void InitRandom(ULong cpu_seed, int nprocs, ULong gpu_seed)
Set the seed of the random number generator.
Definition AMReX_Random.cpp:99
Real Random()
Generate a psuedo-random real from uniform distribution.
Definition AMReX_Random.cpp:133
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
void Bcast(void *, int, MPI_Datatype, int, MPI_Comm)
Definition AMReX_ParallelDescriptor.cpp:1295
void ReduceRealMax(Vector< std::reference_wrapper< Real > > const &)
Definition AMReX_ParallelDescriptor.cpp:1228
Definition AMReX_Amr.cpp:50
__host__ __device__ void ignore_unused(const Ts &...)
This shuts up the compiler about unused variables.
Definition AMReX.H:139
void FileOpenFailed(const std::string &file)
Output a message and abort when couldn't open the file.
Definition AMReX_Utility.cpp:137
void 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
__host__ __device__ Dim3 begin(BoxND< dim > const &box) noexcept
Definition AMReX_Box.H:2006
double second() noexcept
Definition AMReX_Utility.cpp:940
void Error(const std::string &msg)
Print out message to cerr and exit via amrex::Abort().
Definition AMReX.cpp:235
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:2015
Definition AMReX_Particle.H:259
Definition AMReX_Particle.H:154
A struct used to pass initial data into the various Init methods. This struct is used to pass initial...
Definition AMReX_ParticleContainer.H:119
std::array< int, NStructInt > int_struct_data
Definition AMReX_ParticleContainer.H:121
std::array< int, NArrayInt > int_array_data
Definition AMReX_ParticleContainer.H:123
std::array< double, NArrayReal > real_array_data
Definition AMReX_ParticleContainer.H:122
std::array< double, NStructReal > real_struct_data
Definition AMReX_ParticleContainer.H:120
A struct used for storing a particle's position in the AMR hierarchy.
Definition AMReX_ParticleContainer.H:93
int m_grid
Definition AMReX_ParticleContainer.H:95
int m_tile
Definition AMReX_ParticleContainer.H:96
int m_lev
Definition AMReX_ParticleContainer.H:94
Definition AMReX_ParticleTile.H:764
void push_back(const ParticleType &p)
Definition AMReX_ParticleTile.H:1036
int numParticles() const
Returns the number of real particles (excluding neighbors)
Definition AMReX_ParticleTile.H:947
void push_back_int(int comp, int v)
Definition AMReX_ParticleTile.H:1140
void push_back_real(int comp, ParticleReal v)
Definition AMReX_ParticleTile.H:1085
The struct used to store particles.
Definition AMReX_Particle.H:405
__host__ __device__ RealVect pos() const &
Definition AMReX_Particle.H:456
__host__ __device__ ParticleCPUWrapper cpu() &
Definition AMReX_Particle.H:424
__host__ __device__ ParticleIDWrapper id() &
Definition AMReX_Particle.H:427