Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
AMReX_ParticleTileRT.H
Go to the documentation of this file.
1#ifndef AMREX_PARTICLETILERT_H_
2#define AMREX_PARTICLETILERT_H_
3#include <AMReX_Config.H>
4
5#include <AMReX_Extension.H>
6#include <AMReX_Particle.H>
9#include <AMReX_Vector.H>
10#include <AMReX_REAL.H>
11#include <AMReX_RealVect.H>
12
13#include <array>
14#include <numeric>
15#include <string>
16#include <type_traits>
17#include <vector>
18
19
20namespace amrex {
21
22// Forward Declaration
23template <class RType, class IType>
24struct RTSoAParticle;
25
27
28template <class T>
29struct ArrayView {
30 using size_type = Long;
31
32 T* m_data = nullptr;
34
36 T& operator[] (const size_type i) const {
38 return m_data[i];
39 }
40
42 T* data () const {
43 return m_data;
44 }
45
47 T* dataPtr () const {
48 return m_data;
49 }
50
52 T* begin () const {
53 return m_data;
54 }
55
57 T* end () const {
58 return m_data + m_capacity;
59 }
60
62 explicit operator T* () const {
63 return m_data;
64 }
65};
66
67template <class T>
68ArrayView(T* data, std::size_t capacity) -> ArrayView<T>;
69
70template <class RType, class IType>
72{
74 using size_type = Long;
75
77 using RealType = RType;
78 using IntType = IType;
79
80 static constexpr bool is_particle_tile_data = true;
81 static constexpr bool is_const = std::is_const_v<RType>;
82
84 std::conditional_t<is_const,
85 const uint64_t* AMREX_RESTRICT,
86 uint64_t* AMREX_RESTRICT> m_idcpu = nullptr;
87 RType* AMREX_RESTRICT m_rdata = nullptr;
88 IType* AMREX_RESTRICT m_idata = nullptr;
89 int m_n_real = 0;
90 int m_n_int = 0;
91
93 constexpr ParticleTileDataRT () = default;
94
96 constexpr ParticleTileDataRT (size_type a_capacity,
97 decltype(m_idcpu) a_idcpu, decltype(m_rdata) a_rdata,
98 decltype(m_idata) a_idata, int a_n_real, int a_n_int) noexcept :
99 m_capacity{a_capacity},
100 m_idcpu{a_idcpu},
101 m_rdata{a_rdata},
102 m_idata{a_idata},
103 m_n_real{a_n_real},
104 m_n_int{a_n_int} {}
105
106 template <class aRType, class aIType>
109 m_capacity{rhs.m_capacity},
110 m_idcpu{rhs.m_idcpu},
111 m_rdata{rhs.m_rdata},
112 m_idata{rhs.m_idata},
113 m_n_real{rhs.m_n_real},
114 m_n_int{rhs.m_n_int} {}
115
117 RType& pos (const int dir, const size_type index) const &
118 {
119 AMREX_ASSERT(index < m_capacity && 0 <= dir && dir < m_n_real);
120 return this->m_rdata[dir * m_capacity + index];
121 }
122
124 decltype(auto) id (const size_type index) const &
125 {
126 AMREX_ASSERT(index < m_capacity);
127 if constexpr (is_const) {
128 return ConstParticleIDWrapper(this->m_idcpu[index]);
129 } else {
130 return ParticleIDWrapper(this->m_idcpu[index]);
131 }
132 }
133
135 decltype(auto) cpu (const size_type index) const &
136 {
137 AMREX_ASSERT(index < m_capacity);
138 if constexpr (is_const) {
139 return ConstParticleCPUWrapper(this->m_idcpu[index]);
140 } else {
141 return ParticleCPUWrapper(this->m_idcpu[index]);
142 }
143 }
144
146 decltype(auto) idcpu (const size_type index) const &
147 {
148 AMREX_ASSERT(index < m_capacity);
149 return this->m_idcpu[index];
150 }
151
153 RType * rdata (const int comp_index) const &
154 {
155 AMREX_ASSERT(0 <= comp_index && comp_index < m_n_real);
156 return this->m_rdata + comp_index * m_capacity;
157 }
158
160 IType * idata (const int comp_index) const &
161 {
162 AMREX_ASSERT(0 <= comp_index && comp_index < m_n_int);
163 return this->m_idata + comp_index * m_capacity;
164 }
165
167 RType& rdata (const int comp_index, const size_type partice_index) const &
168 {
169 AMREX_ASSERT(partice_index < m_capacity && 0 <= comp_index && comp_index < m_n_real);
170 return this->m_rdata[comp_index * m_capacity + partice_index];
171 }
172
174 IType& idata (const int comp_index, const size_type partice_index) const &
175 {
176 AMREX_ASSERT(partice_index < m_capacity && 0 <= comp_index && comp_index < m_n_int);
177 return this->m_idata[comp_index * m_capacity + partice_index];
178 }
179
181 decltype(auto) operator[] (const size_type index) const
182 {
183 AMREX_ASSERT(index < m_capacity);
184 return RTSoAParticle<RType, IType>(*this, index);
185 }
186
188 void packParticleData (char* buffer, size_type src_index, std::size_t dst_offset,
189 const int* comm_real, const int * comm_int) const noexcept
190 {
191 AMREX_ASSERT(src_index < m_capacity);
192 auto* dst = buffer + dst_offset;
193
194 memcpy(dst, m_idcpu + src_index, sizeof(uint64_t));
195 dst += sizeof(uint64_t);
196
197 for (int i = 0; i < m_n_real; ++i) {
198 if (comm_real[i]) {
199 memcpy(dst, m_rdata + i * m_capacity + src_index, sizeof(RType));
200 dst += sizeof(RType);
201 }
202 }
203
204 for (int i = 0; i < m_n_int; ++i) {
205 if (comm_int[i]) {
206 memcpy(dst, m_idata + i * m_capacity + src_index, sizeof(IType));
207 dst += sizeof(IType);
208 }
209 }
210 }
211
213 void unpackParticleData (const char* buffer, Long src_offset, size_type dst_index,
214 const int* comm_real, const int* comm_int) const noexcept
215 {
216 AMREX_ASSERT(dst_index < m_capacity);
217 const auto* src = buffer + src_offset;
218
219 memcpy(m_idcpu + dst_index, src, sizeof(uint64_t));
220 src += sizeof(uint64_t);
221
222 for (int i = 0; i < m_n_real; ++i) {
223 if (comm_real[i]) {
224 memcpy(m_rdata + i * m_capacity + dst_index, src, sizeof(RType));
225 src += sizeof(RType);
226 }
227 }
228
229 for (int i = 0; i < m_n_int; ++i) {
230 if (comm_int[i]) {
231 memcpy(m_idata + i * m_capacity + dst_index, src, sizeof(IType));
232 src += sizeof(IType);
233 }
234 }
235 }
236
239 {
240 return {{
241 { AMREX_D_DECL( pos(0, index), pos(1, index), pos(2, index) ) },
242 idcpu(index)
243 }};
244 }
245};
246
247// SOA Particle Structure
248template <class RType, class IType>
250{
252 using ConstType = RTSoAParticle<std::add_const_t<RType>, std::add_const_t<IType>>;
253 static constexpr bool is_constsoa_particle = PTD::is_const;
254 static constexpr int NReal = 0;
255 static constexpr int NInt = 0;
256 static constexpr bool is_soa_particle = true;
257 static constexpr bool is_rtsoa_particle = true;
258 using size_type = typename PTD::size_type;
259 using RealType = RType;
260 using IntType = IType;
261
263 RTSoAParticle (PTD const& ptd, size_type i) noexcept :
264 m_particle_tile_data(ptd),
265 m_index(i) {}
266
267 template <class aRType, class aIType>
268 friend struct RTSoAParticle;
269
270 template <class aRType, class aIType>
273 m_particle_tile_data(rhs.m_particle_tile_data),
274 m_index(rhs.m_index) {}
275
276 // functions to get id and cpu in the SOA data
277
279 decltype(auto) cpu () const & { return this->m_particle_tile_data.cpu(m_index); }
280
282 decltype(auto) id () const & { return this->m_particle_tile_data.id(m_index); }
283
285 decltype(auto) idcpu () const & { return this->m_particle_tile_data.idcpu(m_index); }
286
288 RType& rdata (const int comp_index) const & {
289 return this->m_particle_tile_data.rdata(comp_index, m_index);
290 }
291
293 IType& idata (const int comp_index) const & {
294 return this->m_particle_tile_data.idata(comp_index, m_index);
295 }
296
297 // functions to get positions of the particle in the SOA data
298
300 RealVect pos () const & {
301 return RealVect(AMREX_D_DECL(
302 this->m_particle_tile_data.pos(0, m_index),
303 this->m_particle_tile_data.pos(1, m_index),
304 this->m_particle_tile_data.pos(2, m_index)
305 ));
306 }
307
309 RType& pos (int position_index) const & {
310 return this->m_particle_tile_data.pos(position_index, m_index);
311 }
312
313 static Long NextID ();
314
318 static Long UnprotectedNextID ();
319
325 static void NextID (Long nextid);
326
327private:
328
329 static_assert(std::is_trivially_copyable<PTD>(), "ParticleTileData is not trivially copyable");
330
331 PTD m_particle_tile_data;
332 size_type m_index;
333};
334
336 inline static Long the_next_id = 1;
337};
338
339template <class RType, class IType>
340Long
342{
343 Long next;
344// we should be able to test on _OPENMP < 201107 for capture (version 3.1)
345// but we must work around a bug in gcc < 4.9
346#if defined(AMREX_USE_OMP) && defined(_OPENMP) && _OPENMP < 201307
347#pragma omp critical (amrex_particle_nextid)
348#elif defined(AMREX_USE_OMP)
349#pragma omp atomic capture
350#endif
352
354 amrex::Abort("RTSoAParticle<RType, IType>::NextID() -- too many particles");
355 }
356
357 return next;
358}
359
360template <class RType, class IType>
361Long
363{
366 amrex::Abort("RTSoAParticle<RType, IType>::NextID() -- too many particles");
367 }
368 return next;
369}
370
371template <class RType, class IType>
372void
377
378template<class RType, class IType>
380
381template <class RType=ParticleReal, class IType=int>
383{
385 using RealType = RType;
386 using IntType = IType;
387
390
392
395
398
399 ParticleTileRT () = default;
400
401#ifndef _WIN32 // workaround windows compiler bug
402 ~ParticleTileRT () = default;
403
405 ParticleTileRT (ParticleTileRT &&) noexcept = default;
406
407 ParticleTileRT& operator= (ParticleTileRT const&) = delete;
408 ParticleTileRT& operator= (ParticleTileRT &&) noexcept = default;
409#endif
410
411 void define (
412 int a_num_real_comps,
413 int a_num_int_comps,
414 std::vector<std::string>* a_rdata_names=nullptr,
415 std::vector<std::string>* a_idata_names=nullptr,
416 Arena* a_arena=nullptr
417 )
418 {
419 if (m_defined) {
420 AMREX_ALWAYS_ASSERT_WITH_MESSAGE(a_arena == m_arena,
421 "ParticleTileRT redefined with different memory arena");
422 if (a_num_real_comps != m_n_real || a_num_int_comps != m_n_int) {
423 realloc_and_move(m_size, m_capacity, a_num_real_comps, a_num_int_comps);
424 }
425 } else {
426 AMREX_ALWAYS_ASSERT_WITH_MESSAGE(a_arena != nullptr,
427 "ParticleTileRT defined with no memory arena! "
428 "Make sure to call setArena() on the ParticleContainer before initialization or "
429 "to pass an Arena to ParticleTile::define()");
430 m_defined = true;
431 m_arena = a_arena;
432 m_n_real = a_num_real_comps;
433 m_n_int = a_num_int_comps;
434 }
435
436 if (a_rdata_names != nullptr) {
437 m_real_names = a_rdata_names;
438 }
439 if (a_idata_names != nullptr) {
440 m_int_names = a_idata_names;
441 }
442 }
443
444 [[nodiscard]] bool empty () const
445 {
446 AMREX_ALWAYS_ASSERT(m_defined);
447 return size() == 0;
448 }
449
453 [[nodiscard]] size_type size () const
454 {
455 AMREX_ALWAYS_ASSERT(m_defined);
456 return m_size;
457 }
458
462 [[nodiscard]] size_type numParticles () const
463 {
464 AMREX_ALWAYS_ASSERT(m_defined);
465 return m_size;
466 }
467
471 [[nodiscard]] size_type numRealParticles () const
472 {
473 AMREX_ALWAYS_ASSERT(m_defined);
474 return m_size;
475 }
476
481 [[nodiscard]] size_type numNeighborParticles () const
482 {
483 AMREX_ALWAYS_ASSERT(m_defined);
484 return 0;
485 }
486
490 [[nodiscard]] size_type numTotalParticles () const
491 {
492 AMREX_ALWAYS_ASSERT(m_defined);
493 return m_size;
494 }
495
500 [[nodiscard]] size_type getNumNeighbors () const
501 {
502 AMREX_ALWAYS_ASSERT(m_defined);
503 return 0;
504 }
505
506 [[nodiscard]] int NumRealComps () const noexcept
507 {
508 AMREX_ALWAYS_ASSERT(m_defined);
509 return m_n_real;
510 }
511
512 [[nodiscard]] int NumIntComps () const noexcept
513 {
514 AMREX_ALWAYS_ASSERT(m_defined);
515 return m_n_int;
516 }
517
518 [[nodiscard]] int NumRuntimeRealComps () const noexcept
519 {
520 AMREX_ALWAYS_ASSERT(m_defined);
521 return m_n_real;
522 }
523
524 [[nodiscard]] int NumRuntimeIntComps () const noexcept
525 {
526 AMREX_ALWAYS_ASSERT(m_defined);
527 return m_n_int;
528 }
529
531 [[nodiscard]] std::vector<std::string> GetRealNames () const
532 {
533 AMREX_ALWAYS_ASSERT(m_defined);
534 if (m_real_names != nullptr) {
535 return *m_real_names;
536 } else {
537 return std::vector<std::string>();
538 }
539 }
540
542 [[nodiscard]] std::vector<std::string> GetIntNames () const
543 {
544 AMREX_ALWAYS_ASSERT(m_defined);
545 if (m_int_names != nullptr) {
546 return *m_int_names;
547 } else {
548 return std::vector<std::string>();
549 }
550 }
551
554 AMREX_ALWAYS_ASSERT(m_defined);
555 return {m_idcpu_data.data(), m_capacity};
556 }
557
560 AMREX_ALWAYS_ASSERT(m_defined);
561 return {m_idcpu_data.data(), m_capacity};
562 }
563
568 [[nodiscard]] ArrayView<RType> GetRealData (int comp_index) {
569 AMREX_ALWAYS_ASSERT(m_defined);
570 AMREX_ALWAYS_ASSERT(0 <= comp_index && comp_index < m_n_real);
571 return {m_real_data.data() + comp_index * m_capacity, m_capacity};
572 }
573
578 [[nodiscard]] ArrayView<const RType> GetRealData (int comp_index) const {
579 AMREX_ALWAYS_ASSERT(m_defined);
580 AMREX_ALWAYS_ASSERT(0 <= comp_index && comp_index < m_n_real);
581 return {m_real_data.data() + comp_index * m_capacity, m_capacity};
582 }
583
588 [[nodiscard]] ArrayView<IType> GetIntData (int comp_index) {
589 AMREX_ALWAYS_ASSERT(m_defined);
590 AMREX_ALWAYS_ASSERT(0 <= comp_index && comp_index < m_n_int);
591 return {m_int_data.data() + comp_index * m_capacity, m_capacity};
592 }
593
598 [[nodiscard]] ArrayView<const IType> GetIntData (int comp_index) const {
599 AMREX_ALWAYS_ASSERT(m_defined);
600 AMREX_ALWAYS_ASSERT(0 <= comp_index && comp_index < m_n_int);
601 return {m_int_data.data() + comp_index * m_capacity, m_capacity};
602 }
603
608 [[nodiscard]] ArrayView<RType> GetRealData (std::string const & name) {
609 return GetRealData(get_idx_from_str(name, m_real_names));
610 }
611
616 [[nodiscard]] ArrayView<const RType> GetRealData (std::string const & name) const {
617 return GetRealData(get_idx_from_str(name, m_real_names));
618 }
619
624 [[nodiscard]] ArrayView<IType> GetIntData (std::string const & name) {
625 return GetIntData(get_idx_from_str(name, m_int_names));
626 }
627
632 [[nodiscard]] ArrayView<const IType> GetIntData (std::string const & name) const {
633 return GetIntData(get_idx_from_str(name, m_int_names));
634 }
635
637 {
638 AMREX_ALWAYS_ASSERT(m_defined);
639
640 if (m_capacity < new_size) {
641 size_type new_grown_capacity =
642 grow_podvector_capacity(strategy, new_size, m_capacity, sizeof(RType));
643 new_grown_capacity = align_capacity(new_grown_capacity);
644 realloc_and_move(new_size, new_grown_capacity, m_n_real, m_n_int);
645 } else {
646 m_size = new_size;
647 }
648 }
649
651 {
652 AMREX_ALWAYS_ASSERT(m_defined);
653
654 if (m_capacity < new_capacity) {
655 size_type new_grown_capacity =
656 grow_podvector_capacity(strategy, new_capacity, m_capacity, sizeof(RType));
657 new_grown_capacity = align_capacity(new_grown_capacity);
658 realloc_and_move(m_size, new_grown_capacity, m_n_real, m_n_int);
659 }
660 }
661
662 static void reserve (std::map<ParticleTileRT<RType, IType>*, int> const& addsizes)
663 {
664 // we don't do anything special here for now
665 for (auto [p,s] : addsizes) {
666 p->reserve(p->size()+s);
667 }
668 }
669
670 void realloc_and_move (size_type new_size, size_type new_capacity,
671 int new_n_real, int new_n_int)
672 {
673 AMREX_ALWAYS_ASSERT(m_defined);
674
675 if (new_capacity <= new_size) {
676 new_capacity = new_size;
677 }
678
679 decltype(m_idcpu_data) new_idcpu_data;
680 decltype(m_real_data) new_real_data;
681 decltype(m_int_data) new_int_data;
682
683 new_idcpu_data.setArena(m_arena);
684 new_real_data.setArena(m_arena);
685 new_int_data.setArena(m_arena);
686
687 new_idcpu_data.resize(new_capacity, GrowthStrategy::Exact);
688 new_real_data.resize(new_capacity * new_n_real, GrowthStrategy::Exact);
689 new_int_data.resize(new_capacity * new_n_int, GrowthStrategy::Exact);
690
691 auto old_ptd = getParticleTileData();
692 auto copy_size = std::min(m_size, new_size);
693 auto copy_n_real = std::min(m_n_real, new_n_real);
694 auto copy_n_int = std::min(m_n_int, new_n_int);
695
696 m_size = new_size;
697 m_capacity = new_capacity;
698 m_n_real = new_n_real;
699 m_n_int = new_n_int;
700 m_idcpu_data.swap(new_idcpu_data);
701 m_real_data.swap(new_real_data);
702 m_int_data.swap(new_int_data);
703
704 auto new_ptd = getParticleTileData();
705
706 if (new_size == 0) {
707 return;
708 }
709
710 if (m_arena->isDeviceAccessible()) {
711 ParallelFor(copy_size,
712 [=] AMREX_GPU_DEVICE (size_type i) noexcept {
713 new_ptd.idcpu(i) = old_ptd.idcpu(i);
714
715 for (int n = 0; n < copy_n_real; ++n) {
716 new_ptd.rdata(n, i) = old_ptd.rdata(n, i);
717 }
718
719 for (int n = 0; n < copy_n_int; ++n) {
720 new_ptd.idata(n, i) = old_ptd.idata(n, i);
721 }
722 });
724 } else {
725 for (size_type i = 0; i < copy_size; ++i) {
726 new_ptd.idcpu(i) = old_ptd.idcpu(i);
727
728 for (int n = 0; n < copy_n_real; ++n) {
729 new_ptd.rdata(n, i) = old_ptd.rdata(n, i);
730 }
731
732 for (int n = 0; n < copy_n_int; ++n) {
733 new_ptd.idata(n, i) = old_ptd.idata(n, i);
734 }
735 }
736 }
737 }
738
740 {
741 AMREX_ALWAYS_ASSERT(m_defined);
742
743 size_type aligned_size = align_capacity(m_size);
744 if (aligned_size < m_capacity) {
745 realloc_and_move(m_size, aligned_size, m_n_real, m_n_int);
746 }
747 }
748
749 [[nodiscard]] size_type capacity () const
750 {
751 AMREX_ALWAYS_ASSERT(m_defined);
752 return m_capacity * (sizeof(uint64_t) + m_n_real * sizeof(RType) + m_n_int * sizeof(IType));
753 }
754
755 void swap (ParticleTileRT<RType, IType>& other) noexcept
756 {
757 std::swap(m_defined, other.m_defined);
758 std::swap(m_arena, other.m_arena);
759 std::swap(m_capacity, other.m_capacity);
760 std::swap(m_size, other.m_size);
761 std::swap(m_n_real, other.m_n_real);
762 std::swap(m_n_int, other.m_n_int);
763 m_idcpu_data.swap(other.m_idcpu_data);
764 m_real_data.swap(other.m_real_data);
765 m_int_data.swap(other.m_int_data);
766 std::swap(m_real_names, other.m_real_names);
767 std::swap(m_int_names, other.m_int_names);
768 }
769
771 {
772 AMREX_ALWAYS_ASSERT(m_defined);
774 m_capacity,
775 m_idcpu_data.data(),
776 m_real_data.data(),
777 m_int_data.data(),
778 m_n_real,
779 m_n_int
780 };
781 }
782
784 {
785 AMREX_ALWAYS_ASSERT(m_defined);
787 m_capacity,
788 m_idcpu_data.data(),
789 m_real_data.data(),
790 m_int_data.data(),
791 m_n_real,
792 m_n_int
793 };
794 }
795
796 [[nodiscard]] auto& GetStructOfArrays () { return *this; }
797 [[nodiscard]] const auto& GetStructOfArrays () const { return *this; }
798
799 [[nodiscard]] Arena* arena() const {
800 return m_arena;
801 }
802
803private:
804
805 static int get_idx_from_str (std::string const & name, std::vector<std::string>* name_list) {
806 AMREX_ALWAYS_ASSERT(name_list != nullptr);
807 auto const pos = std::find(name_list->begin(), name_list->end(), name);
808 AMREX_ALWAYS_ASSERT(pos != name_list->end());
809 return static_cast<int>(std::distance(name_list->begin(), pos));
810 }
811
812 static size_type align_capacity (size_type capacity) {
813 constexpr auto aligment_in_bytes = Arena::align_size;
814 constexpr auto aligment_in_elements = aligment_in_bytes / std::gcd(
815 std::gcd(aligment_in_bytes, sizeof(uint64_t)),
816 std::gcd(sizeof(RType), sizeof(IType))
817 );
818 return amrex::aligned_size(aligment_in_elements, capacity);
819 }
820
821 bool m_defined = false;
822
823 Arena* m_arena = nullptr;
824
825 size_type m_capacity = 0;
826 size_type m_size = 0;
827
828 int m_n_real = 0;
829 int m_n_int = 0;
830
834
835 std::vector<std::string>* m_real_names = nullptr;
836 std::vector<std::string>* m_int_names = nullptr;
837};
838
839} // namespace amrex
840
841#endif // AMREX_PARTICLETILERT_H_
#define AMREX_ALWAYS_ASSERT_WITH_MESSAGE(EX, MSG)
Definition AMReX_BLassert.H:49
#define AMREX_ASSERT(EX)
Definition AMReX_BLassert.H:38
#define AMREX_ALWAYS_ASSERT(EX)
Definition AMReX_BLassert.H:50
#define AMREX_FORCE_INLINE
Definition AMReX_Extension.H:119
#define AMREX_RESTRICT
Definition AMReX_Extension.H:32
#define AMREX_GPU_DEVICE
Definition AMReX_GpuQualifiers.H:18
#define AMREX_GPU_HOST_DEVICE
Definition AMReX_GpuQualifiers.H:20
#define AMREX_D_DECL(a, b, c)
Definition AMReX_SPACE.H:171
A virtual base class for objects that manage their own dynamic memory allocation.
Definition AMReX_Arena.H:132
static const std::size_t align_size
Definition AMReX_Arena.H:317
virtual bool isDeviceAccessible() const
Definition AMReX_Arena.cpp:66
Dynamically allocated vector for trivially copyable data.
Definition AMReX_PODVector.H:308
amrex_long Long
Definition AMReX_INT.H:30
void streamSynchronize() noexcept
Definition AMReX_GpuDevice.H:310
constexpr Long LastParticleID
Definition AMReX_Particle.H:21
Definition AMReX_Amr.cpp:50
std::enable_if_t< std::is_integral_v< T > > ParallelFor(TypeList< CTOs... > ctos, std::array< int, sizeof...(CTOs)> const &runtime_options, T N, F &&f)
Definition AMReX_CTOParallelForImpl.H:193
GrowthStrategy
Definition AMReX_PODVector.H:250
RealVectND< 3 > RealVect
Definition AMReX_ParmParse.H:36
std::size_t grow_podvector_capacity(GrowthStrategy strategy, std::size_t new_size, std::size_t old_capacity, std::size_t sizeof_T)
Definition AMReX_PODVector.H:268
void Abort(const std::string &msg)
Print out message to cerr and exit via abort().
Definition AMReX.cpp:241
std::size_t aligned_size(std::size_t align_requirement, std::size_t size) noexcept
Given a minimum required size in bytes, this returns the smallest size greater or equal to size that ...
Definition AMReX_Arena.H:33
Definition AMReX_ParticleTileRT.H:29
__host__ __device__ T * end() const
Definition AMReX_ParticleTileRT.H:57
__host__ __device__ T * data() const
Definition AMReX_ParticleTileRT.H:42
__host__ __device__ T * dataPtr() const
Definition AMReX_ParticleTileRT.H:47
__host__ __device__ T & operator[](const size_type i) const
Definition AMReX_ParticleTileRT.H:36
size_type m_capacity
Definition AMReX_ParticleTileRT.H:33
Long size_type
Definition AMReX_ParticleTileRT.H:30
T * m_data
Definition AMReX_ParticleTileRT.H:32
__host__ __device__ T * begin() const
Definition AMReX_ParticleTileRT.H:52
Definition AMReX_Particle.H:327
Definition AMReX_Particle.H:301
Definition AMReX_ParticleTileRT.H:335
static Long the_next_id
Definition AMReX_ParticleTileRT.H:336
Definition AMReX_Particle.H:259
Definition AMReX_Particle.H:154
Definition AMReX_ParticleTileRT.H:72
IType IntType
Definition AMReX_ParticleTileRT.H:78
RType *__restrict__ m_rdata
Definition AMReX_ParticleTileRT.H:87
__host__ __device__ RType * rdata(const int comp_index) const &
Definition AMReX_ParticleTileRT.H:153
__host__ __device__ constexpr ParticleTileDataRT()=default
static constexpr bool is_const
Definition AMReX_ParticleTileRT.H:81
RType RealType
Definition AMReX_ParticleTileRT.H:77
__host__ __device__ constexpr ParticleTileDataRT(size_type a_capacity, decltype(m_idcpu) a_idcpu, decltype(m_rdata) a_rdata, decltype(m_idata) a_idata, int a_n_real, int a_n_int) noexcept
Definition AMReX_ParticleTileRT.H:96
Long size_type
Definition AMReX_ParticleTileRT.H:74
__host__ __device__ IType * idata(const int comp_index) const &
Definition AMReX_ParticleTileRT.H:160
__host__ __device__ void packParticleData(char *buffer, size_type src_index, std::size_t dst_offset, const int *comm_real, const int *comm_int) const noexcept
Definition AMReX_ParticleTileRT.H:188
__host__ __device__ decltype(auto) idcpu(const size_type index) const &
Definition AMReX_ParticleTileRT.H:146
__host__ __device__ IType & idata(const int comp_index, const size_type partice_index) const &
Definition AMReX_ParticleTileRT.H:174
static constexpr bool is_particle_tile_data
Definition AMReX_ParticleTileRT.H:80
__host__ __device__ RType & rdata(const int comp_index, const size_type partice_index) const &
Definition AMReX_ParticleTileRT.H:167
__host__ __device__ decltype(auto) id(const size_type index) const &
Definition AMReX_ParticleTileRT.H:124
IType *__restrict__ m_idata
Definition AMReX_ParticleTileRT.H:88
size_type m_capacity
Definition AMReX_ParticleTileRT.H:83
__host__ __device__ RType & pos(const int dir, const size_type index) const &
Definition AMReX_ParticleTileRT.H:117
__host__ __device__ constexpr ParticleTileDataRT(ParticleTileDataRT< aRType, aIType > const &rhs) noexcept
Definition AMReX_ParticleTileRT.H:108
int m_n_real
Definition AMReX_ParticleTileRT.H:89
int m_n_int
Definition AMReX_ParticleTileRT.H:90
__host__ __device__ Particle< 0, 0 > getSuperParticle(size_type index) const noexcept
Definition AMReX_ParticleTileRT.H:238
std::conditional_t< is_const, const uint64_t *__restrict__, uint64_t *__restrict__ > m_idcpu
Definition AMReX_ParticleTileRT.H:86
__host__ __device__ void unpackParticleData(const char *buffer, Long src_offset, size_type dst_index, const int *comm_real, const int *comm_int) const noexcept
Definition AMReX_ParticleTileRT.H:213
__host__ __device__ decltype(auto) cpu(const size_type index) const &
Definition AMReX_ParticleTileRT.H:135
Definition AMReX_ParticleTileRT.H:383
const auto & GetStructOfArrays() const
Definition AMReX_ParticleTileRT.H:797
int NumRuntimeIntComps() const noexcept
Definition AMReX_ParticleTileRT.H:524
size_type numParticles() const
Returns the number of particles.
Definition AMReX_ParticleTileRT.H:462
ArrayView< RType > GetRealData(std::string const &name)
Definition AMReX_ParticleTileRT.H:608
int NumRealComps() const noexcept
Definition AMReX_ParticleTileRT.H:506
ArrayView< RType > GetRealData(int comp_index)
Definition AMReX_ParticleTileRT.H:568
ArrayView< IType > GetIntData(std::string const &name)
Definition AMReX_ParticleTileRT.H:624
void swap(ParticleTileRT< RType, IType > &other) noexcept
Definition AMReX_ParticleTileRT.H:755
typename ParticleTileDataType::size_type size_type
Definition AMReX_ParticleTileRT.H:391
RType RealType
Definition AMReX_ParticleTileRT.H:385
void realloc_and_move(size_type new_size, size_type new_capacity, int new_n_real, int new_n_int)
Definition AMReX_ParticleTileRT.H:670
bool empty() const
Definition AMReX_ParticleTileRT.H:444
void shrink_to_fit()
Definition AMReX_ParticleTileRT.H:739
size_type getNumNeighbors() const
Returns the number of neighbor particles. For ParticleTileRT this is always zero.
Definition AMReX_ParticleTileRT.H:500
std::vector< std::string > GetIntNames() const
Definition AMReX_ParticleTileRT.H:542
ArrayView< const RType > GetRealData(std::string const &name) const
Definition AMReX_ParticleTileRT.H:616
void define(int a_num_real_comps, int a_num_int_comps, std::vector< std::string > *a_rdata_names=nullptr, std::vector< std::string > *a_idata_names=nullptr, Arena *a_arena=nullptr)
Definition AMReX_ParticleTileRT.H:411
ConstParticleTileDataType getConstParticleTileData() const
Definition AMReX_ParticleTileRT.H:783
int NumRuntimeRealComps() const noexcept
Definition AMReX_ParticleTileRT.H:518
ParticleTileDataType getParticleTileData()
Definition AMReX_ParticleTileRT.H:770
size_type size() const
Returns the number of particles.
Definition AMReX_ParticleTileRT.H:453
size_type numRealParticles() const
Returns the number of particles.
Definition AMReX_ParticleTileRT.H:471
ParticleTileRT(ParticleTileRT &&) noexcept=default
Arena * arena() const
Definition AMReX_ParticleTileRT.H:799
int NumIntComps() const noexcept
Definition AMReX_ParticleTileRT.H:512
ArrayView< uint64_t > GetIdCPUData()
Definition AMReX_ParticleTileRT.H:553
ArrayView< const IType > GetIntData(std::string const &name) const
Definition AMReX_ParticleTileRT.H:632
std::vector< std::string > GetRealNames() const
Definition AMReX_ParticleTileRT.H:531
ArrayView< IType > GetIntData(int comp_index)
Definition AMReX_ParticleTileRT.H:588
ArrayView< const RType > GetRealData(int comp_index) const
Definition AMReX_ParticleTileRT.H:578
ArrayView< const uint64_t > GetIdCPUData() const
Definition AMReX_ParticleTileRT.H:559
void resize(size_type new_size, GrowthStrategy strategy=GrowthStrategy::Poisson)
Definition AMReX_ParticleTileRT.H:636
size_type numTotalParticles() const
Returns the number of particles.
Definition AMReX_ParticleTileRT.H:490
ParticleTileRT(ParticleTileRT const &)=delete
ArrayView< const IType > GetIntData(int comp_index) const
Definition AMReX_ParticleTileRT.H:598
size_type capacity() const
Definition AMReX_ParticleTileRT.H:749
auto & GetStructOfArrays()
Definition AMReX_ParticleTileRT.H:796
size_type numNeighborParticles() const
Returns the number of neighbor particles. For ParticleTileRT this is always zero.
Definition AMReX_ParticleTileRT.H:481
static void reserve(std::map< ParticleTileRT< RType, IType > *, int > const &addsizes)
Definition AMReX_ParticleTileRT.H:662
IType IntType
Definition AMReX_ParticleTileRT.H:386
void reserve(size_type new_capacity, GrowthStrategy strategy=GrowthStrategy::Poisson)
Definition AMReX_ParticleTileRT.H:650
The struct used to store particles.
Definition AMReX_Particle.H:405
Definition AMReX_ParticleTileRT.H:250
typename PTD::size_type size_type
Definition AMReX_ParticleTileRT.H:258
static Long NextID()
Definition AMReX_ParticleTileRT.H:341
__host__ __device__ IType & idata(const int comp_index) const &
Definition AMReX_ParticleTileRT.H:293
static constexpr int NReal
Definition AMReX_ParticleTileRT.H:254
__host__ __device__ RTSoAParticle(PTD const &ptd, size_type i) noexcept
Definition AMReX_ParticleTileRT.H:263
__host__ __device__ RType & pos(int position_index) const &
Definition AMReX_ParticleTileRT.H:309
__host__ __device__ RType & rdata(const int comp_index) const &
Definition AMReX_ParticleTileRT.H:288
__host__ __device__ decltype(auto) idcpu() const &
Definition AMReX_ParticleTileRT.H:285
__host__ __device__ RealVect pos() const &
Definition AMReX_ParticleTileRT.H:300
__host__ __device__ RTSoAParticle(RTSoAParticle< aRType, aIType > const &rhs) noexcept
Definition AMReX_ParticleTileRT.H:272
ParticleTileDataRT< RType, IType > PTD
Definition AMReX_ParticleTileRT.H:251
IType IntType
Definition AMReX_ParticleTileRT.H:260
__host__ __device__ decltype(auto) cpu() const &
Definition AMReX_ParticleTileRT.H:279
static constexpr bool is_rtsoa_particle
Definition AMReX_ParticleTileRT.H:257
RType RealType
Definition AMReX_ParticleTileRT.H:259
static Long UnprotectedNextID()
This version can only be used inside omp critical.
Definition AMReX_ParticleTileRT.H:362
static constexpr int NInt
Definition AMReX_ParticleTileRT.H:255
static constexpr bool is_constsoa_particle
Definition AMReX_ParticleTileRT.H:253
Definition AMReX_ParticleTile.H:717
Definition AMReX_MakeParticle.H:13