WarpX
Loading...
Searching...
No Matches
MultiFabRegister.H
Go to the documentation of this file.
1/* Copyright 2024 The ABLASTR Community
2 *
3 * This file is part of ABLASTR.
4 *
5 * License: BSD-3-Clause-LBNL
6 * Authors: Axel Huebl
7 */
8#ifndef ABLASTR_FIELDS_MF_REGISTER_H
9#define ABLASTR_FIELDS_MF_REGISTER_H
10
11#include <AMReX_BaseFwd.H>
12#include <AMReX_Enum.H>
13#include <AMReX_Extension.H>
14#include <AMReX_IntVect.H>
15#include <AMReX_MultiFab.H>
16#include <AMReX_REAL.H>
17#include <AMReX_Vector.H>
18
19#include <array>
20#include <map>
21#include <memory>
22#include <optional>
23#include <stdexcept>
24#include <string>
25#include <type_traits>
26#include <utility>
27#include <vector>
28
29
30namespace
31{
32 // type trait helpers in lieu of an amrex::is_amrex_enum
33 template <typename, typename = std::void_t<>>
34 struct is_castable_to_string : std::false_type {};
35
36 template <typename T>
37 struct is_castable_to_string<T, std::void_t<decltype(static_cast<std::string>(std::declval<T>()))>> : std::true_type {};
38
40 template<typename T>
41 std::string getExtractedName (T name)
42 {
43 if constexpr(is_castable_to_string<T>())
44 {
45 // already a unique string key
46 return std::string(name);
47 } else
48 {
49 // user-defined AMREX_ENUM or compile error
50 return amrex::getEnumNameString(name);
51 }
52 }
53}
54
55namespace ablastr::fields
56{
71 {
72 int dir = 0;
73
74 public:
75 constexpr explicit Direction (int d) : dir(d) {}
76
77#if defined(WARPX_DIM_RZ) || defined(WARPX_DIM_RCYLINDER) || defined(WARPX_DIM_RSPHERE)
78 static const Direction r;
79 static const Direction theta;
80#endif
81#if defined(WARPX_DIM_RSPHERE)
82 static const Direction phi;
83#endif
84#if defined(WARPX_DIM_3D) || defined(WARPX_DIM_XZ) || defined(WARPX_DIM_1D_Z)
85 static const Direction x;
86 static const Direction y;
87#endif
88#if !defined(WARPX_DIM_RSPHERE)
89 static const Direction z;
90#endif
91
92 bool operator<(const Direction& other) const
93 {
94 return other.dir < this->dir;
95 }
96
97 operator std::string() const
98 {
99#if defined(WARPX_DIM_RZ) || defined(WARPX_DIM_RCYLINDER) || defined(WARPX_DIM_RSPHERE)
100 if (dir == r) { return "r"; }
101 if (dir == theta) { return "theta"; }
102#endif
103#if defined(WARPX_DIM_RSPHERE)
104 if (dir == phi) { return "phi"; }
105#endif
106#if defined(WARPX_DIM_3D) || defined(WARPX_DIM_XZ) || defined(WARPX_DIM_1D_Z)
107 if (dir == x) { return "x"; }
108 if (dir == y) { return "y"; }
109#endif
110#if !defined(WARPX_DIM_RSPHERE)
111 if (dir == z) { return "z"; }
112#endif
113 throw std::runtime_error("invalid direction: " + std::to_string(dir));
114 return std::to_string(dir);
115 }
116
117 Direction (std::string const & s)
118 {
119 dir = -1;
120#if defined(WARPX_DIM_RZ) || defined(WARPX_DIM_RCYLINDER) || defined(WARPX_DIM_RSPHERE)
121 if (s.compare(r) == 0) { dir = r; }
122 if (s.compare(theta) == 0) { dir = theta; }
123#endif
124#if defined(WARPX_DIM_RSPHERE)
125 if (s.compare(phi) == 0) { dir = phi; }
126#endif
127#if defined(WARPX_DIM_3D) || defined(WARPX_DIM_XZ) || defined(WARPX_DIM_1D_Z)
128 if (s.compare(x) == 0) { dir = x; }
129 if (s.compare(y) == 0) { dir = y; }
130#endif
131#if !defined(WARPX_DIM_RSPHERE)
132 if (s.compare(z) == 0) { dir = z; }
133#endif
134
135 if (dir == -1) {
136 throw std::runtime_error("invalid direction: " + s);
137 }
138 }
139
140 // e.g., "x"
141 Direction (char const * c) : dir{Direction{c}} {}
142
143 // e.g., 'x'
144 Direction (char const c) : dir{Direction{c}} {}
145
146 // rule of five: define default constructors and the destructor
147 // because we defined a special one above for char
148 ~Direction () = default;
149 Direction (const Direction&) = default;
150 Direction& operator= (const Direction&) = default;
151 Direction (Direction&&) = default;
153
154 /* TODO: just temporary int compatibility */
155 operator int() const { return dir; }
156
157 };
158
159#if defined(WARPX_DIM_RZ) || defined(WARPX_DIM_RCYLINDER) || defined(WARPX_DIM_RSPHERE)
160 inline constexpr Direction Direction::r = Direction{0};
161 inline constexpr Direction Direction::theta = Direction{1};
162#endif
163#if defined(WARPX_DIM_RSPHERE)
164 inline constexpr Direction Direction::phi = Direction{2};
165#endif
166
167#if defined(WARPX_DIM_3D) || defined(WARPX_DIM_XZ) || defined(WARPX_DIM_1D_Z)
168 inline constexpr Direction Direction::x = Direction{0};
169 inline constexpr Direction Direction::y = Direction{1};
170#endif
171
172#if defined(WARPX_DIM_3D) || defined(WARPX_DIM_XZ) || defined(WARPX_DIM_1D_Z) || defined(WARPX_DIM_RZ) || defined(WARPX_DIM_RCYLINDER)
173 inline constexpr Direction Direction::z = Direction{2};
174#endif
175
181
187
190 //using VectorField = ablastr::utils::ConstMap<Direction, amrex::MultiFab *>;
191 using VectorField = std::array<amrex::MultiFab *, 3>;
192
195 //using VectorField = ablastr::utils::ConstMap<Direction, amrex::MultiFab const *>;
196 using ConstVectorField = std::array<amrex::MultiFab const *, 3>;
197
201
205
209
213
221 {
222 // TODO: also add iMultiFab via std::variant
223
226
228 std::optional<Direction> m_dir = std::nullopt;
229
231 int m_level = 0;
232
234 bool m_remake = true;
235
238
240 std::string m_owner;
241
244 bool
245 is_vector () const { return m_dir.has_value(); }
246
252 bool
253 is_alias () const { return !m_owner.empty(); }
254 };
255
262 {
263 // Avoid accidental copies when passing to member functions
264 MultiFabRegister() = default;
269 ~MultiFabRegister() = default;
270
287 template<typename T>
290 T name,
291 int level,
292 amrex::BoxArray const & ba,
294 int ncomp,
295 amrex::IntVect const & ngrow,
296 std::optional<amrex::Real const> initial_value = std::nullopt,
297 bool remake = true,
298 bool redistribute_on_remake = true
299 )
300 {
301 return internal_alloc_init(
302 getExtractedName(name),
303 level,
304 ba,
305 dm,
306 ncomp,
307 ngrow,
308 initial_value,
309 remake,
310 redistribute_on_remake
311 );
312 }
313
331 template<typename T>
334 T name,
335 Direction dir,
336 int level,
337 amrex::BoxArray const & ba,
339 int ncomp,
340 amrex::IntVect const & ngrow,
341 std::optional<amrex::Real const> initial_value = std::nullopt,
342 bool remake = true,
343 bool redistribute_on_remake = true
344 )
345 {
346 return internal_alloc_init(
347 getExtractedName(name),
348 dir,
349 level,
350 ba,
351 dm,
352 ncomp,
353 ngrow,
354 initial_value,
355 remake,
356 redistribute_on_remake
357 );
358 }
359
371 template<typename N, typename A>
374 N new_name,
375 A alias_name,
376 int level,
377 std::optional<amrex::Real const> initial_value = std::nullopt
378 )
379 {
380 return internal_alias_init(
381 getExtractedName(new_name),
382 getExtractedName(alias_name),
383 level,
384 initial_value
385 );
386 }
387
400 template<typename N, typename A>
403 N new_name,
404 A alias_name,
405 Direction dir,
406 int level,
407 std::optional<amrex::Real const> initial_value = std::nullopt
408 )
409 {
410 return internal_alias_init(
411 getExtractedName(new_name),
412 getExtractedName(alias_name),
413 dir,
414 level,
415 initial_value
416 );
417 }
418
425 template<typename T>
426 [[nodiscard]] bool
428 T name,
429 int level
430 ) const
431 {
432 return internal_has(
433 getExtractedName(name),
434 level
435 );
436 }
437
445 template<typename T>
446 [[nodiscard]] bool
448 T name,
449 Direction dir,
450 int level
451 ) const
452 {
453 return internal_has(
454 getExtractedName(name),
455 dir,
456 level
457 );
458 }
459
466 template<typename T>
467 [[nodiscard]] bool
469 T name,
470 int level
471 ) const
472 {
473 return internal_has_vector(
474 getExtractedName(name),
475 level
476 );
477 }
478
487 template<typename T>
488 [[nodiscard]] amrex::MultiFab*
490 T name,
491 int level
492 )
493 {
494 return internal_get(
495 getExtractedName(name),
496 level
497 );
498 }
499
509 template<typename T>
510 [[nodiscard]] amrex::MultiFab*
512 T name,
513 Direction dir,
514 int level
515 )
516 {
517 return internal_get(
518 getExtractedName(name),
519 dir,
520 level
521 );
522 }
523
532 template<typename T>
533 [[nodiscard]] amrex::MultiFab const *
535 T name,
536 int level
537 ) const
538 {
539 return internal_get(
540 getExtractedName(name),
541 level
542 );
543 }
544
554 template<typename T>
555 [[nodiscard]] amrex::MultiFab const *
557 T name,
558 Direction dir,
559 int level
560 ) const
561 {
562 return internal_get(
563 getExtractedName(name),
564 dir,
565 level
566 );
567 }
568
579 template<typename T>
580 [[nodiscard]] MultiLevelScalarField
582 T name,
583 int finest_level,
584 bool skip_level_0=false
585 )
586 {
588 getExtractedName(name),
589 finest_level,
590 skip_level_0
591 );
592 }
593 template<typename T>
594 [[nodiscard]] ConstMultiLevelScalarField
596 T name,
597 int finest_level,
598 bool skip_level_0=false
599 ) const
600 {
602 getExtractedName(name),
603 finest_level,
604 skip_level_0
605 );
606 }
607
608
618 template<typename T>
619 [[nodiscard]] VectorField
621 T name,
622 int level
623 )
624 {
626 getExtractedName(name),
627 level
628 );
629 }
630 template<typename T>
631 [[nodiscard]] ConstVectorField
633 T name,
634 int level
635 ) const
636 {
638 getExtractedName(name),
639 level
640 );
641 }
642
643
655 template<typename T>
656 [[nodiscard]] MultiLevelVectorField
658 T name,
659 int finest_level,
660 bool skip_level_0=false
661 )
662 {
664 getExtractedName(name),
665 finest_level,
666 skip_level_0
667 );
668 }
669 template<typename T>
670 [[nodiscard]] ConstMultiLevelVectorField
672 T name,
673 int finest_level,
674 bool skip_level_0=false
675 ) const
676 {
678 getExtractedName(name),
679 finest_level,
680 skip_level_0
681 );
682 }
683
684
689 [[nodiscard]] std::vector<std::string>
690 list () const;
691
697 template<typename T>
698 void
700 T name,
701 int level
702 )
703 {
704 internal_erase(getExtractedName(name), level);
705 }
706
713 template<typename T>
714 void
716 T name,
717 Direction dir,
718 int level
719 )
720 {
721 internal_erase(getExtractedName(name), dir, level);
722 }
723
730 void
732 int level
733 );
734
742 void
744 int other_level,
745 amrex::DistributionMapping const & new_dm
746 );
747
754 [[nodiscard]] std::string
755 mf_name (
756 std::string name,
757 int level
758 ) const;
759
767 [[nodiscard]] std::string
768 mf_name (
769 std::string name,
770 Direction dir,
771 int level
772 ) const;
773
775 [[nodiscard]] bool
777 std::string const & internal_name
778 );
779 [[nodiscard]] amrex::MultiFab *
781 std::string const & internal_name
782 );
783
784 private:
785
786 [[nodiscard]] amrex::MultiFab const *
788 std::string const & internal_name
789 ) const;
790
793 std::string const & name,
794 int level,
795 amrex::BoxArray const & ba,
797 int ncomp,
798 amrex::IntVect const & ngrow,
799 std::optional<amrex::Real const> initial_value = std::nullopt,
800 bool remake = true,
801 bool redistribute_on_remake = true
802 );
805 std::string const & name,
806 Direction dir,
807 int level,
808 amrex::BoxArray const & ba,
810 int ncomp,
811 amrex::IntVect const & ngrow,
812 std::optional<amrex::Real const> initial_value = std::nullopt,
813 bool remake = true,
814 bool redistribute_on_remake = true
815 );
816
819 std::string const & new_name,
820 std::string const & alias_name,
821 int level,
822 std::optional<amrex::Real const> initial_value = std::nullopt
823 );
826 std::string const & new_name,
827 std::string const & alias_name,
828 Direction dir,
829 int level,
830 std::optional<amrex::Real const> initial_value = std::nullopt
831 );
832
833 [[nodiscard]] bool
835 std::string const & name,
836 int level
837 ) const;
838 [[nodiscard]] bool
840 std::string const & name,
841 Direction dir,
842 int level
843 ) const;
844 [[nodiscard]] bool
846 std::string const & name,
847 int level
848 ) const;
849
850 [[nodiscard]] amrex::MultiFab *
852 std::string const & name,
853 int level
854 );
855 [[nodiscard]] amrex::MultiFab const *
857 std::string const & name,
858 int level
859 ) const;
860 [[nodiscard]] amrex::MultiFab *
862 std::string const & name,
863 Direction dir,
864 int level
865 );
866 [[nodiscard]] amrex::MultiFab const *
868 std::string const & name,
869 Direction dir,
870 int level
871 ) const;
872 [[nodiscard]] MultiLevelScalarField
874 std::string const & name,
875 int finest_level,
876 bool skip_level_0
877 );
878 [[nodiscard]] ConstMultiLevelScalarField
880 std::string const & name,
881 int finest_level,
882 bool skip_level_0
883 ) const;
884 [[nodiscard]] VectorField
886 std::string const & name,
887 int level
888 );
889 [[nodiscard]] ConstVectorField
891 std::string const & name,
892 int level
893 ) const;
894 [[nodiscard]] MultiLevelVectorField
896 std::string const & name,
897 int finest_level,
898 bool skip_level_0
899 );
900 [[nodiscard]] ConstMultiLevelVectorField
902 std::string const & name,
903 int finest_level,
904 bool skip_level_0
905 ) const;
906
907 void
909 std::string const & name,
910 int level
911 );
912 void
914 std::string const & name,
915 Direction dir,
916 int level
917 );
918
920 std::map<
921 std::string,
924
925 public:
934 static inline std::vector<Direction> m_all_dirs =
935#if defined(WARPX_DIM_3D) || defined(WARPX_DIM_XZ) || defined(WARPX_DIM_1D_Z)
937#elif defined(WARPX_DIM_RZ) || defined(WARPX_DIM_RCYLINDER)
939#elif defined(WARPX_DIM_RSPHERE)
940 {Direction::r, Direction::theta, Direction::phi};
941#endif
942 };
943
949 a2m (
950 std::array< std::unique_ptr<amrex::MultiFab>, 3 > const & old_vectorfield
951 );
952
953} // namespace ablastr::fields
954
955#endif // ABLASTR_FIELDS_MF_REGISTER_H
#define AMREX_INLINE
Definition MultiFabRegister.H:71
static const Direction r
Definition MultiFabRegister.H:160
static const Direction z
Definition MultiFabRegister.H:173
static const Direction y
Definition MultiFabRegister.H:169
constexpr Direction(int d)
Definition MultiFabRegister.H:75
static const Direction x
Definition MultiFabRegister.H:168
static const Direction theta
Definition MultiFabRegister.H:161
Definition MultiFabRegister.H:71
static const Direction r
Definition MultiFabRegister.H:78
int dir
Definition MultiFabRegister.H:72
Direction(Direction &&)=default
static const Direction z
Definition MultiFabRegister.H:89
Direction & operator=(const Direction &)=default
Direction(std::string const &s)
Definition MultiFabRegister.H:117
static const Direction y
Definition MultiFabRegister.H:86
constexpr Direction(int d)
Definition MultiFabRegister.H:75
static const Direction x
Definition MultiFabRegister.H:85
Direction(char const c)
Definition MultiFabRegister.H:144
bool operator<(const Direction &other) const
Definition MultiFabRegister.H:92
Direction(char const *c)
Definition MultiFabRegister.H:141
Direction(const Direction &)=default
static const Direction theta
Definition MultiFabRegister.H:79
Definition EffectivePotentialPoissonSolver.H:63
VectorField a2m(std::array< std::unique_ptr< amrex::MultiFab >, 3 > const &old_vectorfield)
Definition MultiFabRegister.cpp:645
amrex::Vector< ConstVectorField > ConstMultiLevelVectorField
Definition MultiFabRegister.H:212
std::array< amrex::MultiFab const *, 3 > ConstVectorField
Definition MultiFabRegister.H:196
std::array< amrex::MultiFab *, 3 > VectorField
Definition MultiFabRegister.H:191
amrex::Vector< ScalarField > MultiLevelScalarField
Definition MultiFabRegister.H:200
amrex::MultiFab * ScalarField
Definition MultiFabRegister.H:180
amrex::Vector< ConstScalarField > ConstMultiLevelScalarField
Definition MultiFabRegister.H:204
amrex::MultiFab const * ConstScalarField
Definition MultiFabRegister.H:186
amrex::Vector< VectorField > MultiLevelVectorField
Definition MultiFabRegister.H:208
std::string getEnumNameString(T const &v)
IntVectND< 3 > IntVect
Definition MultiFabRegister.H:221
amrex::MultiFab m_mf
Definition MultiFabRegister.H:225
AMREX_INLINE bool is_alias() const
Definition MultiFabRegister.H:253
int m_level
Definition MultiFabRegister.H:231
std::optional< Direction > m_dir
Definition MultiFabRegister.H:228
AMREX_INLINE bool is_vector() const
Definition MultiFabRegister.H:245
bool m_redistribute_on_remake
Definition MultiFabRegister.H:237
std::string m_owner
Definition MultiFabRegister.H:240
bool m_remake
Definition MultiFabRegister.H:234
amrex::MultiFab * alias_init(N new_name, A alias_name, int level, std::optional< amrex::Real const > initial_value=std::nullopt)
Definition MultiFabRegister.H:373
bool internal_has(std::string const &internal_name)
Definition MultiFabRegister.cpp:342
ConstMultiLevelScalarField get_mr_levels(T name, int finest_level, bool skip_level_0=false) const
Definition MultiFabRegister.H:595
MultiFabRegister & operator=(MultiFabRegister &&)=delete
MultiLevelScalarField get_mr_levels(T name, int finest_level, bool skip_level_0=false)
Definition MultiFabRegister.H:581
amrex::MultiFab const * get(T name, Direction dir, int level) const
Definition MultiFabRegister.H:556
amrex::MultiFab * alloc_init(T name, int level, amrex::BoxArray const &ba, amrex::DistributionMapping const &dm, int ncomp, amrex::IntVect const &ngrow, std::optional< amrex::Real const > initial_value=std::nullopt, bool remake=true, bool redistribute_on_remake=true)
Definition MultiFabRegister.H:289
std::vector< std::string > list() const
Definition MultiFabRegister.cpp:560
std::string mf_name(std::string name, int level) const
Definition MultiFabRegister.cpp:615
MultiFabRegister(MultiFabRegister &&)=delete
VectorField get_alldirs(T name, int level)
Definition MultiFabRegister.H:620
ConstVectorField get_alldirs(T name, int level) const
Definition MultiFabRegister.H:632
bool internal_has_vector(std::string const &name, int level) const
Definition MultiFabRegister.cpp:326
void remake_level(int other_level, amrex::DistributionMapping const &new_dm)
Definition MultiFabRegister.cpp:247
void clear_level(int level)
Definition MultiFabRegister.cpp:599
void internal_erase(std::string const &name, int level)
Definition MultiFabRegister.cpp:570
static std::vector< Direction > m_all_dirs
Definition MultiFabRegister.H:934
MultiLevelScalarField internal_get_mr_levels(std::string const &name, int finest_level, bool skip_level_0)
Definition MultiFabRegister.cpp:418
amrex::MultiFab * alloc_init(T name, Direction dir, int level, amrex::BoxArray const &ba, amrex::DistributionMapping const &dm, int ncomp, amrex::IntVect const &ngrow, std::optional< amrex::Real const > initial_value=std::nullopt, bool remake=true, bool redistribute_on_remake=true)
Definition MultiFabRegister.H:333
std::map< std::string, MultiFabOwner > m_mf_register
Definition MultiFabRegister.H:923
amrex::MultiFab * get(T name, Direction dir, int level)
Definition MultiFabRegister.H:511
amrex::MultiFab * internal_get(std::string const &internal_name)
Definition MultiFabRegister.cpp:350
amrex::MultiFab * get(T name, int level)
Definition MultiFabRegister.H:489
MultiLevelVectorField internal_get_mr_levels_alldirs(std::string const &name, int finest_level, bool skip_level_0)
Definition MultiFabRegister.cpp:498
ConstMultiLevelVectorField get_mr_levels_alldirs(T name, int finest_level, bool skip_level_0=false) const
Definition MultiFabRegister.H:671
bool has_vector(T name, int level) const
Definition MultiFabRegister.H:468
MultiLevelVectorField get_mr_levels_alldirs(T name, int finest_level, bool skip_level_0=false)
Definition MultiFabRegister.H:657
amrex::MultiFab const * get(T name, int level) const
Definition MultiFabRegister.H:534
amrex::MultiFab * alias_init(N new_name, A alias_name, Direction dir, int level, std::optional< amrex::Real const > initial_value=std::nullopt)
Definition MultiFabRegister.H:402
amrex::MultiFab * internal_alloc_init(std::string const &name, int level, amrex::BoxArray const &ba, amrex::DistributionMapping const &dm, int ncomp, amrex::IntVect const &ngrow, std::optional< amrex::Real const > initial_value=std::nullopt, bool remake=true, bool redistribute_on_remake=true)
Definition MultiFabRegister.cpp:26
MultiFabRegister & operator=(MultiFabRegister const &)=delete
bool has(T name, int level) const
Definition MultiFabRegister.H:427
amrex::MultiFab * internal_alias_init(std::string const &new_name, std::string const &alias_name, int level, std::optional< amrex::Real const > initial_value=std::nullopt)
Definition MultiFabRegister.cpp:129
void erase(T name, int level)
Definition MultiFabRegister.H:699
VectorField internal_get_alldirs(std::string const &name, int level)
Definition MultiFabRegister.cpp:464
bool has(T name, Direction dir, int level) const
Definition MultiFabRegister.H:447
void erase(T name, Direction dir, int level)
Definition MultiFabRegister.H:715
MultiFabRegister(MultiFabRegister const &)=delete