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 <string>
15#include <type_traits>
16#include <vector>
17
18
19namespace amrex {
20
21// Forward Declaration
22template <class RType, class IType>
23struct RTSoAParticle;
24
26
27template <class T>
28struct ArrayView {
29 using size_type = Long;
30
31 T* m_data = nullptr;
33
35 T& operator[] (const size_type i) const {
37 return m_data[i];
38 }
39
41 T* data () const {
42 return m_data;
43 }
44
46 T* dataPtr () const {
47 return m_data;
48 }
49
51 T* begin () const {
52 return m_data;
53 }
54
56 T* end () const {
57 return m_data + m_capacity;
58 }
59
61 explicit operator T* () const {
62 return m_data;
63 }
64};
65
66template <class T>
67ArrayView(T* data, std::size_t capacity) -> ArrayView<T>;
68
69template <class RType, class IType>
71{
73 using size_type = Long;
74
76 using RealType = RType;
77 using IntType = IType;
78
79 static constexpr bool is_particle_tile_data = true;
80 static constexpr bool is_const = std::is_const_v<RType>;
81
83 std::conditional_t<is_const,
84 const uint64_t* AMREX_RESTRICT,
85 uint64_t* AMREX_RESTRICT> m_idcpu = nullptr;
86 RType* AMREX_RESTRICT m_rdata = nullptr;
87 IType* AMREX_RESTRICT m_idata = nullptr;
88 int m_n_real = 0;
89 int m_n_int = 0;
90
92 constexpr ParticleTileDataRT () = default;
93
95 constexpr ParticleTileDataRT (size_type a_capacity,
96 decltype(m_idcpu) a_idcpu, decltype(m_rdata) a_rdata,
97 decltype(m_idata) a_idata, int a_n_real, int a_n_int) noexcept :
98 m_capacity{a_capacity},
99 m_idcpu{a_idcpu},
100 m_rdata{a_rdata},
101 m_idata{a_idata},
102 m_n_real{a_n_real},
103 m_n_int{a_n_int} {}
104
105 template <class aRType, class aIType>
108 m_capacity{rhs.m_capacity},
109 m_idcpu{rhs.m_idcpu},
110 m_rdata{rhs.m_rdata},
111 m_idata{rhs.m_idata},
112 m_n_real{rhs.m_n_real},
113 m_n_int{rhs.m_n_int} {}
114
116 RType& pos (const int dir, const size_type index) const &
117 {
118 AMREX_ASSERT(index < m_capacity && 0 <= dir && dir < m_n_real);
119 return this->m_rdata[dir * m_capacity + index];
120 }
121
123 decltype(auto) id (const size_type index) const &
124 {
125 AMREX_ASSERT(index < m_capacity);
126 if constexpr (is_const) {
127 return ConstParticleIDWrapper(this->m_idcpu[index]);
128 } else {
129 return ParticleIDWrapper(this->m_idcpu[index]);
130 }
131 }
132
134 decltype(auto) cpu (const size_type index) const &
135 {
136 AMREX_ASSERT(index < m_capacity);
137 if constexpr (is_const) {
138 return ConstParticleCPUWrapper(this->m_idcpu[index]);
139 } else {
140 return ParticleCPUWrapper(this->m_idcpu[index]);
141 }
142 }
143
145 decltype(auto) idcpu (const size_type index) const &
146 {
147 AMREX_ASSERT(index < m_capacity);
148 return this->m_idcpu[index];
149 }
150
152 RType * rdata (const int comp_index) const &
153 {
154 AMREX_ASSERT(0 <= comp_index && comp_index < m_n_real);
155 return this->m_rdata + comp_index * m_capacity;
156 }
157
159 IType * idata (const int comp_index) const &
160 {
161 AMREX_ASSERT(0 <= comp_index && comp_index < m_n_int);
162 return this->m_idata + comp_index * m_capacity;
163 }
164
166 RType& rdata (const int comp_index, const size_type partice_index) const &
167 {
168 AMREX_ASSERT(partice_index < m_capacity && 0 <= comp_index && comp_index < m_n_real);
169 return this->m_rdata[comp_index * m_capacity + partice_index];
170 }
171
173 IType& idata (const int comp_index, const size_type partice_index) const &
174 {
175 AMREX_ASSERT(partice_index < m_capacity && 0 <= comp_index && comp_index < m_n_int);
176 return this->m_idata[comp_index * m_capacity + partice_index];
177 }
178
180 decltype(auto) operator[] (const size_type index) const
181 {
182 AMREX_ASSERT(index < m_capacity);
183 return RTSoAParticle<RType, IType>(*this, index);
184 }
185
187 void packParticleData (char* buffer, size_type src_index, std::size_t dst_offset,
188 const int* comm_real, const int * comm_int) const noexcept
189 {
190 AMREX_ASSERT(src_index < m_capacity);
191 auto* dst = buffer + dst_offset;
192
193 memcpy(dst, m_idcpu + src_index, sizeof(uint64_t));
194 dst += sizeof(uint64_t);
195
196 for (int i = 0; i < m_n_real; ++i) {
197 if (comm_real[i]) {
198 memcpy(dst, m_rdata + i * m_capacity + src_index, sizeof(RType));
199 dst += sizeof(RType);
200 }
201 }
202
203 for (int i = 0; i < m_n_int; ++i) {
204 if (comm_int[i]) {
205 memcpy(dst, m_idata + i * m_capacity + src_index, sizeof(IType));
206 dst += sizeof(IType);
207 }
208 }
209 }
210
212 void unpackParticleData (const char* buffer, Long src_offset, size_type dst_index,
213 const int* comm_real, const int* comm_int) const noexcept
214 {
215 AMREX_ASSERT(dst_index < m_capacity);
216 const auto* src = buffer + src_offset;
217
218 memcpy(m_idcpu + dst_index, src, sizeof(uint64_t));
219 src += sizeof(uint64_t);
220
221 for (int i = 0; i < m_n_real; ++i) {
222 if (comm_real[i]) {
223 memcpy(m_rdata + i * m_capacity + dst_index, src, sizeof(RType));
224 src += sizeof(RType);
225 }
226 }
227
228 for (int i = 0; i < m_n_int; ++i) {
229 if (comm_int[i]) {
230 memcpy(m_idata + i * m_capacity + dst_index, src, sizeof(IType));
231 src += sizeof(IType);
232 }
233 }
234 }
235
238 {
239 return {{
240 { AMREX_D_DECL( pos(0, index), pos(1, index), pos(2, index) ) },
241 idcpu(index)
242 }};
243 }
244};
245
246// SOA Particle Structure
247template <class RType, class IType>
249{
251 using ConstType = RTSoAParticle<std::add_const_t<RType>, std::add_const_t<IType>>;
252 static constexpr bool is_constsoa_particle = PTD::is_const;
253 static constexpr int NReal = 0;
254 static constexpr int NInt = 0;
255 static constexpr bool is_soa_particle = true;
256 static constexpr bool is_rtsoa_particle = true;
257 using size_type = typename PTD::size_type;
258 using RealType = RType;
259 using IntType = IType;
260
262 RTSoAParticle (PTD const& ptd, size_type i) noexcept :
263 m_particle_tile_data(ptd),
264 m_index(i) {}
265
266 template <class aRType, class aIType>
267 friend struct RTSoAParticle;
268
269 template <class aRType, class aIType>
272 m_particle_tile_data(rhs.m_particle_tile_data),
273 m_index(rhs.m_index) {}
274
275 // functions to get id and cpu in the SOA data
276
278 decltype(auto) cpu () const & { return this->m_particle_tile_data.cpu(m_index); }
279
281 decltype(auto) id () const & { return this->m_particle_tile_data.id(m_index); }
282
284 decltype(auto) idcpu () const & { return this->m_particle_tile_data.idcpu(m_index); }
285
287 RType& rdata (const int comp_index) const & {
288 return this->m_particle_tile_data.rdata(comp_index, m_index);
289 }
290
292 IType& idata (const int comp_index) const & {
293 return this->m_particle_tile_data.idata(comp_index, m_index);
294 }
295
296 // functions to get positions of the particle in the SOA data
297
299 RealVect pos () const & {
300 return RealVect(AMREX_D_DECL(
301 this->m_particle_tile_data.pos(0, m_index),
302 this->m_particle_tile_data.pos(1, m_index),
303 this->m_particle_tile_data.pos(2, m_index)
304 ));
305 }
306
308 RType& pos (int position_index) const & {
309 return this->m_particle_tile_data.pos(position_index, m_index);
310 }
311
312 static Long NextID ();
313
317 static Long UnprotectedNextID ();
318
324 static void NextID (Long nextid);
325
326private:
327
328 static_assert(std::is_trivially_copyable<PTD>(), "ParticleTileData is not trivially copyable");
329
330 PTD m_particle_tile_data;
331 size_type m_index;
332};
333
335 inline static Long the_next_id = 1;
336};
337
338template <class RType, class IType>
339Long
341{
342 Long next;
343// we should be able to test on _OPENMP < 201107 for capture (version 3.1)
344// but we must work around a bug in gcc < 4.9
345#if defined(AMREX_USE_OMP) && defined(_OPENMP) && _OPENMP < 201307
346#pragma omp critical (amrex_particle_nextid)
347#elif defined(AMREX_USE_OMP)
348#pragma omp atomic capture
349#endif
351
353 amrex::Abort("RTSoAParticle<RType, IType>::NextID() -- too many particles");
354 }
355
356 return next;
357}
358
359template <class RType, class IType>
360Long
362{
365 amrex::Abort("RTSoAParticle<RType, IType>::NextID() -- too many particles");
366 }
367 return next;
368}
369
370template <class RType, class IType>
371void
376
377template<class RType, class IType>
379
380template <class RType=ParticleReal, class IType=int>
382{
384 using RealType = RType;
385 using IntType = IType;
386
389
391
394
397
398 ParticleTileRT () = default;
399
400#ifndef _WIN32 // workaround windows compiler bug
401 ~ParticleTileRT () = default;
402
404 ParticleTileRT (ParticleTileRT &&) noexcept = default;
405
406 ParticleTileRT& operator= (ParticleTileRT const&) = delete;
407 ParticleTileRT& operator= (ParticleTileRT &&) noexcept = default;
408#endif
409
410 void define (
411 int a_num_real_comps,
412 int a_num_int_comps,
413 std::vector<std::string>* a_rdata_names=nullptr,
414 std::vector<std::string>* a_idata_names=nullptr,
415 Arena* a_arena=nullptr
416 )
417 {
418 if (m_defined) {
419 AMREX_ALWAYS_ASSERT_WITH_MESSAGE(a_arena == m_arena,
420 "ParticleTileRT redefined with different memory arena");
421 if (a_num_real_comps != m_n_real || a_num_int_comps != m_n_int) {
422 realloc_and_move(m_size, m_capacity, a_num_real_comps, a_num_int_comps);
423 }
424 } else {
425 AMREX_ALWAYS_ASSERT_WITH_MESSAGE(a_arena != nullptr,
426 "ParticleTileRT defined with no memory arena! "
427 "Make sure to call setArena() on the ParticleContainer before initialization or "
428 "to pass an Arena to ParticleTile::define()");
429 m_defined = true;
430 m_arena = a_arena;
431 m_n_real = a_num_real_comps;
432 m_n_int = a_num_int_comps;
433 }
434
435 if (a_rdata_names != nullptr) {
436 m_real_names = a_rdata_names;
437 }
438 if (a_idata_names != nullptr) {
439 m_int_names = a_idata_names;
440 }
441 }
442
443 [[nodiscard]] bool empty () const
444 {
445 AMREX_ALWAYS_ASSERT(m_defined);
446 return size() == 0;
447 }
448
452 [[nodiscard]] size_type size () const
453 {
454 AMREX_ALWAYS_ASSERT(m_defined);
455 return m_size;
456 }
457
461 [[nodiscard]] size_type numParticles () const
462 {
463 AMREX_ALWAYS_ASSERT(m_defined);
464 return m_size;
465 }
466
470 [[nodiscard]] size_type numRealParticles () const
471 {
472 AMREX_ALWAYS_ASSERT(m_defined);
473 return m_size;
474 }
475
480 [[nodiscard]] size_type numNeighborParticles () const
481 {
482 AMREX_ALWAYS_ASSERT(m_defined);
483 return 0;
484 }
485
489 [[nodiscard]] size_type numTotalParticles () const
490 {
491 AMREX_ALWAYS_ASSERT(m_defined);
492 return m_size;
493 }
494
499 [[nodiscard]] size_type getNumNeighbors () const
500 {
501 AMREX_ALWAYS_ASSERT(m_defined);
502 return 0;
503 }
504
505 [[nodiscard]] int NumRealComps () const noexcept
506 {
507 AMREX_ALWAYS_ASSERT(m_defined);
508 return m_n_real;
509 }
510
511 [[nodiscard]] int NumIntComps () const noexcept
512 {
513 AMREX_ALWAYS_ASSERT(m_defined);
514 return m_n_int;
515 }
516
517 [[nodiscard]] int NumRuntimeRealComps () const noexcept
518 {
519 AMREX_ALWAYS_ASSERT(m_defined);
520 return m_n_real;
521 }
522
523 [[nodiscard]] int NumRuntimeIntComps () const noexcept
524 {
525 AMREX_ALWAYS_ASSERT(m_defined);
526 return m_n_int;
527 }
528
530 [[nodiscard]] std::vector<std::string> GetRealNames () const
531 {
532 AMREX_ALWAYS_ASSERT(m_defined);
533 if (m_real_names != nullptr) {
534 return *m_real_names;
535 } else {
536 return std::vector<std::string>();
537 }
538 }
539
541 [[nodiscard]] std::vector<std::string> GetIntNames () const
542 {
543 AMREX_ALWAYS_ASSERT(m_defined);
544 if (m_int_names != nullptr) {
545 return *m_int_names;
546 } else {
547 return std::vector<std::string>();
548 }
549 }
550
553 AMREX_ALWAYS_ASSERT(m_defined);
554 return {m_idcpu_data.data(), m_capacity};
555 }
556
559 AMREX_ALWAYS_ASSERT(m_defined);
560 return {m_idcpu_data.data(), m_capacity};
561 }
562
567 [[nodiscard]] ArrayView<RType> GetRealData (int comp_index) {
568 AMREX_ALWAYS_ASSERT(m_defined);
569 AMREX_ALWAYS_ASSERT(0 <= comp_index && comp_index < m_n_real);
570 return {m_real_data.data() + comp_index * m_capacity, m_capacity};
571 }
572
577 [[nodiscard]] ArrayView<const RType> GetRealData (int comp_index) const {
578 AMREX_ALWAYS_ASSERT(m_defined);
579 AMREX_ALWAYS_ASSERT(0 <= comp_index && comp_index < m_n_real);
580 return {m_real_data.data() + comp_index * m_capacity, m_capacity};
581 }
582
587 [[nodiscard]] ArrayView<IType> GetIntData (int comp_index) {
588 AMREX_ALWAYS_ASSERT(m_defined);
589 AMREX_ALWAYS_ASSERT(0 <= comp_index && comp_index < m_n_int);
590 return {m_int_data.data() + comp_index * m_capacity, m_capacity};
591 }
592
597 [[nodiscard]] ArrayView<const IType> GetIntData (int comp_index) const {
598 AMREX_ALWAYS_ASSERT(m_defined);
599 AMREX_ALWAYS_ASSERT(0 <= comp_index && comp_index < m_n_int);
600 return {m_int_data.data() + comp_index * m_capacity, m_capacity};
601 }
602
607 [[nodiscard]] ArrayView<RType> GetRealData (std::string const & name) {
608 return GetRealData(get_idx_from_str(name, m_real_names));
609 }
610
615 [[nodiscard]] ArrayView<const RType> GetRealData (std::string const & name) const {
616 return GetRealData(get_idx_from_str(name, m_real_names));
617 }
618
623 [[nodiscard]] ArrayView<IType> GetIntData (std::string const & name) {
624 return GetIntData(get_idx_from_str(name, m_int_names));
625 }
626
631 [[nodiscard]] ArrayView<const IType> GetIntData (std::string const & name) const {
632 return GetIntData(get_idx_from_str(name, m_int_names));
633 }
634
636 {
637 AMREX_ALWAYS_ASSERT(m_defined);
638
639 if (m_capacity < new_size) {
640 size_type new_grown_capacity =
641 grow_podvector_capacity(strategy, new_size, m_capacity, sizeof(RType));
642 realloc_and_move(new_size, new_grown_capacity, m_n_real, m_n_int);
643 } else {
644 m_size = new_size;
645 }
646 }
647
649 {
650 AMREX_ALWAYS_ASSERT(m_defined);
651
652 if (m_capacity < new_capacity) {
653 size_type new_grown_capacity =
654 grow_podvector_capacity(strategy, new_capacity, m_capacity, sizeof(RType));
655 realloc_and_move(m_size, new_grown_capacity, m_n_real, m_n_int);
656 }
657 }
658
659 static void reserve (std::map<ParticleTileRT<RType, IType>*, int> const& addsizes)
660 {
661 // we don't do anything special here for now
662 for (auto [p,s] : addsizes) {
663 p->reserve(p->size()+s);
664 }
665 }
666
667 void realloc_and_move (size_type new_size, size_type new_capacity,
668 int new_n_real, int new_n_int)
669 {
670 AMREX_ALWAYS_ASSERT(m_defined);
671
672 if (new_capacity <= new_size) {
673 new_capacity = new_size;
674 }
675
676 decltype(m_idcpu_data) new_idcpu_data;
677 decltype(m_real_data) new_real_data;
678 decltype(m_int_data) new_int_data;
679
680 new_idcpu_data.setArena(m_arena);
681 new_real_data.setArena(m_arena);
682 new_int_data.setArena(m_arena);
683
684 new_idcpu_data.resize(new_capacity, GrowthStrategy::Exact);
685 new_real_data.resize(new_capacity * new_n_real, GrowthStrategy::Exact);
686 new_int_data.resize(new_capacity * new_n_int, GrowthStrategy::Exact);
687
688 auto old_ptd = getParticleTileData();
689 auto copy_size = std::min(m_size, new_size);
690 auto copy_n_real = std::min(m_n_real, new_n_real);
691 auto copy_n_int = std::min(m_n_int, new_n_int);
692
693 m_size = new_size;
694 m_capacity = new_capacity;
695 m_n_real = new_n_real;
696 m_n_int = new_n_int;
697 m_idcpu_data.swap(new_idcpu_data);
698 m_real_data.swap(new_real_data);
699 m_int_data.swap(new_int_data);
700
701 auto new_ptd = getParticleTileData();
702
703 if (new_size == 0) {
704 return;
705 }
706
707 if (m_arena->isDeviceAccessible()) {
708 ParallelFor(copy_size,
709 [=] AMREX_GPU_DEVICE (size_type i) noexcept {
710 new_ptd.idcpu(i) = old_ptd.idcpu(i);
711
712 for (int n = 0; n < copy_n_real; ++n) {
713 new_ptd.rdata(n, i) = old_ptd.rdata(n, i);
714 }
715
716 for (int n = 0; n < copy_n_int; ++n) {
717 new_ptd.idata(n, i) = old_ptd.idata(n, i);
718 }
719 });
721 } else {
722 for (size_type i = 0; i < copy_size; ++i) {
723 new_ptd.idcpu(i) = old_ptd.idcpu(i);
724
725 for (int n = 0; n < copy_n_real; ++n) {
726 new_ptd.rdata(n, i) = old_ptd.rdata(n, i);
727 }
728
729 for (int n = 0; n < copy_n_int; ++n) {
730 new_ptd.idata(n, i) = old_ptd.idata(n, i);
731 }
732 }
733 }
734 }
735
737 {
738 AMREX_ALWAYS_ASSERT(m_defined);
739
740 if (m_size != m_capacity) {
741 realloc_and_move(m_size, m_size, m_n_real, m_n_int);
742 }
743 }
744
745 [[nodiscard]] size_type capacity () const
746 {
747 AMREX_ALWAYS_ASSERT(m_defined);
748 return m_capacity * (sizeof(uint64_t) + m_n_real * sizeof(RType) + m_n_int * sizeof(IType));
749 }
750
751 void swap (ParticleTileRT<RType, IType>& other) noexcept
752 {
753 std::swap(m_defined, other.m_defined);
754 std::swap(m_arena, other.m_arena);
755 std::swap(m_capacity, other.m_capacity);
756 std::swap(m_size, other.m_size);
757 std::swap(m_n_real, other.m_n_real);
758 std::swap(m_n_int, other.m_n_int);
759 m_idcpu_data.swap(other.m_idcpu_data);
760 m_real_data.swap(other.m_real_data);
761 m_int_data.swap(other.m_int_data);
762 std::swap(m_real_names, other.m_real_names);
763 std::swap(m_int_names, other.m_int_names);
764 }
765
767 {
768 AMREX_ALWAYS_ASSERT(m_defined);
770 m_capacity,
771 m_idcpu_data.data(),
772 m_real_data.data(),
773 m_int_data.data(),
774 m_n_real,
775 m_n_int
776 };
777 }
778
780 {
781 AMREX_ALWAYS_ASSERT(m_defined);
783 m_capacity,
784 m_idcpu_data.data(),
785 m_real_data.data(),
786 m_int_data.data(),
787 m_n_real,
788 m_n_int
789 };
790 }
791
792 [[nodiscard]] auto& GetStructOfArrays () { return *this; }
793 [[nodiscard]] const auto& GetStructOfArrays () const { return *this; }
794
795 [[nodiscard]] Arena* arena() const {
796 return m_arena;
797 }
798
799private:
800
801 static int get_idx_from_str (std::string const & name, std::vector<std::string>* name_list) {
802 AMREX_ALWAYS_ASSERT(name_list != nullptr);
803 auto const pos = std::find(name_list->begin(), name_list->end(), name);
804 AMREX_ALWAYS_ASSERT(pos != name_list->end());
805 return static_cast<int>(std::distance(name_list->begin(), pos));
806 }
807
808 bool m_defined = false;
809
810 Arena* m_arena = nullptr;
811
812 size_type m_capacity = 0;
813 size_type m_size = 0;
814
815 int m_n_real = 0;
816 int m_n_int = 0;
817
821
822 std::vector<std::string>* m_real_names = nullptr;
823 std::vector<std::string>* m_int_names = nullptr;
824};
825
826} // namespace amrex
827
828#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
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:49
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:35
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:240
Definition AMReX_ParticleTileRT.H:28
__host__ __device__ T * end() const
Definition AMReX_ParticleTileRT.H:56
__host__ __device__ T * data() const
Definition AMReX_ParticleTileRT.H:41
__host__ __device__ T * dataPtr() const
Definition AMReX_ParticleTileRT.H:46
__host__ __device__ T & operator[](const size_type i) const
Definition AMReX_ParticleTileRT.H:35
size_type m_capacity
Definition AMReX_ParticleTileRT.H:32
Long size_type
Definition AMReX_ParticleTileRT.H:29
T * m_data
Definition AMReX_ParticleTileRT.H:31
__host__ __device__ T * begin() const
Definition AMReX_ParticleTileRT.H:51
Definition AMReX_Particle.H:327
Definition AMReX_Particle.H:301
Definition AMReX_ParticleTileRT.H:334
static Long the_next_id
Definition AMReX_ParticleTileRT.H:335
Definition AMReX_Particle.H:259
Definition AMReX_Particle.H:154
Definition AMReX_ParticleTileRT.H:71
IType IntType
Definition AMReX_ParticleTileRT.H:77
RType *__restrict__ m_rdata
Definition AMReX_ParticleTileRT.H:86
__host__ __device__ RType * rdata(const int comp_index) const &
Definition AMReX_ParticleTileRT.H:152
__host__ __device__ constexpr ParticleTileDataRT()=default
static constexpr bool is_const
Definition AMReX_ParticleTileRT.H:80
RType RealType
Definition AMReX_ParticleTileRT.H:76
__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:95
Long size_type
Definition AMReX_ParticleTileRT.H:73
__host__ __device__ IType * idata(const int comp_index) const &
Definition AMReX_ParticleTileRT.H:159
__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:187
__host__ __device__ decltype(auto) idcpu(const size_type index) const &
Definition AMReX_ParticleTileRT.H:145
__host__ __device__ IType & idata(const int comp_index, const size_type partice_index) const &
Definition AMReX_ParticleTileRT.H:173
static constexpr bool is_particle_tile_data
Definition AMReX_ParticleTileRT.H:79
__host__ __device__ RType & rdata(const int comp_index, const size_type partice_index) const &
Definition AMReX_ParticleTileRT.H:166
__host__ __device__ decltype(auto) id(const size_type index) const &
Definition AMReX_ParticleTileRT.H:123
IType *__restrict__ m_idata
Definition AMReX_ParticleTileRT.H:87
size_type m_capacity
Definition AMReX_ParticleTileRT.H:82
__host__ __device__ RType & pos(const int dir, const size_type index) const &
Definition AMReX_ParticleTileRT.H:116
__host__ __device__ constexpr ParticleTileDataRT(ParticleTileDataRT< aRType, aIType > const &rhs) noexcept
Definition AMReX_ParticleTileRT.H:107
int m_n_real
Definition AMReX_ParticleTileRT.H:88
int m_n_int
Definition AMReX_ParticleTileRT.H:89
__host__ __device__ Particle< 0, 0 > getSuperParticle(size_type index) const noexcept
Definition AMReX_ParticleTileRT.H:237
std::conditional_t< is_const, const uint64_t *__restrict__, uint64_t *__restrict__ > m_idcpu
Definition AMReX_ParticleTileRT.H:85
__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:212
__host__ __device__ decltype(auto) cpu(const size_type index) const &
Definition AMReX_ParticleTileRT.H:134
Definition AMReX_ParticleTileRT.H:382
const auto & GetStructOfArrays() const
Definition AMReX_ParticleTileRT.H:793
int NumRuntimeIntComps() const noexcept
Definition AMReX_ParticleTileRT.H:523
size_type numParticles() const
Returns the number of particles.
Definition AMReX_ParticleTileRT.H:461
ArrayView< RType > GetRealData(std::string const &name)
Definition AMReX_ParticleTileRT.H:607
int NumRealComps() const noexcept
Definition AMReX_ParticleTileRT.H:505
ArrayView< RType > GetRealData(int comp_index)
Definition AMReX_ParticleTileRT.H:567
ArrayView< IType > GetIntData(std::string const &name)
Definition AMReX_ParticleTileRT.H:623
void swap(ParticleTileRT< RType, IType > &other) noexcept
Definition AMReX_ParticleTileRT.H:751
typename ParticleTileDataType::size_type size_type
Definition AMReX_ParticleTileRT.H:390
RType RealType
Definition AMReX_ParticleTileRT.H:384
void realloc_and_move(size_type new_size, size_type new_capacity, int new_n_real, int new_n_int)
Definition AMReX_ParticleTileRT.H:667
bool empty() const
Definition AMReX_ParticleTileRT.H:443
void shrink_to_fit()
Definition AMReX_ParticleTileRT.H:736
size_type getNumNeighbors() const
Returns the number of neighbor particles. For ParticleTileRT this is always zero.
Definition AMReX_ParticleTileRT.H:499
std::vector< std::string > GetIntNames() const
Definition AMReX_ParticleTileRT.H:541
ArrayView< const RType > GetRealData(std::string const &name) const
Definition AMReX_ParticleTileRT.H:615
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:410
ConstParticleTileDataType getConstParticleTileData() const
Definition AMReX_ParticleTileRT.H:779
int NumRuntimeRealComps() const noexcept
Definition AMReX_ParticleTileRT.H:517
ParticleTileDataType getParticleTileData()
Definition AMReX_ParticleTileRT.H:766
size_type size() const
Returns the number of particles.
Definition AMReX_ParticleTileRT.H:452
size_type numRealParticles() const
Returns the number of particles.
Definition AMReX_ParticleTileRT.H:470
ParticleTileRT(ParticleTileRT &&) noexcept=default
Arena * arena() const
Definition AMReX_ParticleTileRT.H:795
int NumIntComps() const noexcept
Definition AMReX_ParticleTileRT.H:511
ArrayView< uint64_t > GetIdCPUData()
Definition AMReX_ParticleTileRT.H:552
ArrayView< const IType > GetIntData(std::string const &name) const
Definition AMReX_ParticleTileRT.H:631
std::vector< std::string > GetRealNames() const
Definition AMReX_ParticleTileRT.H:530
ArrayView< IType > GetIntData(int comp_index)
Definition AMReX_ParticleTileRT.H:587
ArrayView< const RType > GetRealData(int comp_index) const
Definition AMReX_ParticleTileRT.H:577
ArrayView< const uint64_t > GetIdCPUData() const
Definition AMReX_ParticleTileRT.H:558
void resize(size_type new_size, GrowthStrategy strategy=GrowthStrategy::Poisson)
Definition AMReX_ParticleTileRT.H:635
size_type numTotalParticles() const
Returns the number of particles.
Definition AMReX_ParticleTileRT.H:489
ParticleTileRT(ParticleTileRT const &)=delete
ArrayView< const IType > GetIntData(int comp_index) const
Definition AMReX_ParticleTileRT.H:597
size_type capacity() const
Definition AMReX_ParticleTileRT.H:745
auto & GetStructOfArrays()
Definition AMReX_ParticleTileRT.H:792
size_type numNeighborParticles() const
Returns the number of neighbor particles. For ParticleTileRT this is always zero.
Definition AMReX_ParticleTileRT.H:480
static void reserve(std::map< ParticleTileRT< RType, IType > *, int > const &addsizes)
Definition AMReX_ParticleTileRT.H:659
IType IntType
Definition AMReX_ParticleTileRT.H:385
void reserve(size_type new_capacity, GrowthStrategy strategy=GrowthStrategy::Poisson)
Definition AMReX_ParticleTileRT.H:648
The struct used to store particles.
Definition AMReX_Particle.H:405
Definition AMReX_ParticleTileRT.H:249
typename PTD::size_type size_type
Definition AMReX_ParticleTileRT.H:257
static Long NextID()
Definition AMReX_ParticleTileRT.H:340
__host__ __device__ IType & idata(const int comp_index) const &
Definition AMReX_ParticleTileRT.H:292
static constexpr int NReal
Definition AMReX_ParticleTileRT.H:253
__host__ __device__ RTSoAParticle(PTD const &ptd, size_type i) noexcept
Definition AMReX_ParticleTileRT.H:262
__host__ __device__ RType & pos(int position_index) const &
Definition AMReX_ParticleTileRT.H:308
__host__ __device__ RType & rdata(const int comp_index) const &
Definition AMReX_ParticleTileRT.H:287
__host__ __device__ decltype(auto) idcpu() const &
Definition AMReX_ParticleTileRT.H:284
__host__ __device__ RealVect pos() const &
Definition AMReX_ParticleTileRT.H:299
__host__ __device__ RTSoAParticle(RTSoAParticle< aRType, aIType > const &rhs) noexcept
Definition AMReX_ParticleTileRT.H:271
ParticleTileDataRT< RType, IType > PTD
Definition AMReX_ParticleTileRT.H:250
IType IntType
Definition AMReX_ParticleTileRT.H:259
__host__ __device__ decltype(auto) cpu() const &
Definition AMReX_ParticleTileRT.H:278
static constexpr bool is_rtsoa_particle
Definition AMReX_ParticleTileRT.H:256
RType RealType
Definition AMReX_ParticleTileRT.H:258
static Long UnprotectedNextID()
This version can only be used inside omp critical.
Definition AMReX_ParticleTileRT.H:361
static constexpr int NInt
Definition AMReX_ParticleTileRT.H:254
static constexpr bool is_constsoa_particle
Definition AMReX_ParticleTileRT.H:252
Definition AMReX_ParticleTile.H:716
Definition AMReX_MakeParticle.H:13