WarpX
Loading...
Searching...
No Matches
WarpXSolverVec.H
Go to the documentation of this file.
1/* Copyright 2024 Justin Angus
2 *
3 * This file is part of WarpX.
4 *
5 * License: BSD-3-Clause-LBNL
6 */
7#ifndef WarpXSolverVec_H_
8#define WarpXSolverVec_H_
9
10#include "Utils/TextMsg.H"
12#include "Fields.H"
13
17
18#include <AMReX.H>
19#include <AMReX_Array.H>
20#include <AMReX_BLassert.H>
21#include <AMReX_IntVect.H>
22#include <AMReX_LayoutData.H>
23#include <AMReX_MultiFab.H>
24#include <AMReX_iMultiFab.H>
25#include <AMReX_ParmParse.H>
26#include <AMReX_Print.H>
27#include <AMReX_REAL.H>
28#include <AMReX_Utility.H>
29#include <AMReX_Vector.H>
30
31#include <algorithm>
32#include <array>
33#include <memory>
34#include <ostream>
35#include <vector>
36
38
39// forward declaration
40class WarpX;
41
58{
59public:
60
61 WarpXSolverVec() = default;
62
64
66
67 using value_type = amrex::Real;
68 using RT = value_type;
69
70 [[nodiscard]] inline bool IsDefined () const { return m_is_defined; }
71
72 void Define ( WarpX* a_WarpX,
73 const std::string& a_vector_type_name,
74 const std::string& a_scalar_type_name = "none" );
75
76 inline
77 void Define ( const WarpXSolverVec& a_solver_vec )
78 {
79 assertIsDefined( a_solver_vec );
82 a_solver_vec.getVectorType(),
83 a_solver_vec.getScalarType() );
84 }
85
86 [[nodiscard]] RT dotProduct( const WarpXSolverVec& a_X ) const;
87
88 void Copy ( FieldType a_array_type,
89 FieldType a_scalar_type = FieldType::None,
90 bool allow_type_mismatch = false);
91
92 inline
93 void Copy ( const WarpXSolverVec& a_solver_vec )
94 {
95 assertIsDefined( a_solver_vec );
96 if (IsDefined()) { assertSameType( a_solver_vec ); }
97 else { Define(a_solver_vec); }
98
99 for (int lev = 0; lev < m_num_amr_levels; ++lev) {
101 for (int n = 0; n < 3; ++n) {
102 const amrex::MultiFab* this_field = a_solver_vec.getArrayVec()[lev][n];
103 amrex::MultiFab::Copy( *m_array_vec[lev][n], *this_field, 0, 0, m_ncomp,
105 }
106 }
108 const amrex::MultiFab* this_scalar = a_solver_vec.getScalarVec()[lev];
109 amrex::MultiFab::Copy( *m_scalar_vec[lev], *this_scalar, 0, 0, m_ncomp,
111 }
112 }
113 }
114
115 // Prohibit Copy assignment operator
116 WarpXSolverVec& operator= ( const WarpXSolverVec& a_solver_vec ) = delete;
117
118 // Move assignment operator
119 WarpXSolverVec(WarpXSolverVec&&) noexcept = default;
120 WarpXSolverVec& operator= ( WarpXSolverVec&& a_solver_vec ) noexcept
121 {
122 if (this != &a_solver_vec) {
123 m_array_vec = std::move(a_solver_vec.m_array_vec);
124 m_scalar_vec = std::move(a_solver_vec.m_scalar_vec);
125 m_array_type = a_solver_vec.m_array_type;
126 m_scalar_type = a_solver_vec.m_scalar_type;
127 m_is_defined = true;
128 }
129 return *this;
130 }
131
132 inline
133 void operator+= ( const WarpXSolverVec& a_solver_vec )
134 {
135 assertIsDefined( a_solver_vec );
136 assertSameType( a_solver_vec );
137 for (int lev = 0; lev < m_num_amr_levels; ++lev) {
139 m_array_vec[lev][0]->plus(*(a_solver_vec.getArrayVec()[lev][0]), 0, 1, 0);
140 m_array_vec[lev][1]->plus(*(a_solver_vec.getArrayVec()[lev][1]), 0, 1, 0);
141 m_array_vec[lev][2]->plus(*(a_solver_vec.getArrayVec()[lev][2]), 0, 1, 0);
142 }
144 m_scalar_vec[lev]->plus(*(a_solver_vec.getScalarVec()[lev]), 0, 1, 0);
145 }
146 }
147 }
148
149 inline
150 void operator-= (const WarpXSolverVec& a_solver_vec)
151 {
152 assertIsDefined( a_solver_vec );
153 assertSameType( a_solver_vec );
154 for (int lev = 0; lev < m_num_amr_levels; ++lev) {
156 m_array_vec[lev][0]->minus(*(a_solver_vec.getArrayVec()[lev][0]), 0, 1, 0);
157 m_array_vec[lev][1]->minus(*(a_solver_vec.getArrayVec()[lev][1]), 0, 1, 0);
158 m_array_vec[lev][2]->minus(*(a_solver_vec.getArrayVec()[lev][2]), 0, 1, 0);
159 }
161 m_scalar_vec[lev]->minus(*(a_solver_vec.getScalarVec()[lev]), 0, 1, 0);
162 }
163 }
164 }
165
169 inline
170 void linComb (const RT a, const WarpXSolverVec& X, const RT b, const WarpXSolverVec& Y)
171 {
172 assertIsDefined( X );
173 assertIsDefined( Y );
174 assertSameType( X );
175 assertSameType( Y );
176 for (int lev = 0; lev < m_num_amr_levels; ++lev) {
178 for (int n = 0; n < 3; n++) {
179 amrex::MultiFab::LinComb(*m_array_vec[lev][n], a, *X.getArrayVec()[lev][n], 0,
180 b, *Y.getArrayVec()[lev][n], 0,
181 0, 1, 0);
182 }
183 }
186 b, *Y.getScalarVec()[lev], 0,
187 0, 1, 0);
188 }
189 }
190 }
191
195 void increment (const WarpXSolverVec& X, const RT a)
196 {
197 assertIsDefined( X );
198 assertSameType( X );
199 for (int lev = 0; lev < m_num_amr_levels; ++lev) {
201 for (int n = 0; n < 3; n++) {
202 amrex::MultiFab::Saxpy( *m_array_vec[lev][n], a, *X.getArrayVec()[lev][n],
204 }
205 }
209 }
210 }
211 }
212
216 inline
217 void scale (RT a_a)
218 {
220 IsDefined(),
221 "WarpXSolverVec::scale() called on undefined WarpXSolverVec");
222 for (int lev = 0; lev < m_num_amr_levels; ++lev) {
224 m_array_vec[lev][0]->mult(a_a, 0, 1);
225 m_array_vec[lev][1]->mult(a_a, 0, 1);
226 m_array_vec[lev][2]->mult(a_a, 0, 1);
227 }
229 m_scalar_vec[lev]->mult(a_a, 0, 1);
230 }
231 }
232 }
233
234 inline
235 void zero () { setVal(0.0); }
236
237 inline
238 void setVal ( const RT a_val )
239 {
241 IsDefined(),
242 "WarpXSolverVec::setVal() called on undefined WarpXSolverVec");
243 for (int lev = 0; lev < m_num_amr_levels; ++lev) {
245 m_array_vec[lev][0]->setVal(a_val);
246 m_array_vec[lev][1]->setVal(a_val);
247 m_array_vec[lev][2]->setVal(a_val);
248 }
250 m_scalar_vec[lev]->setVal(a_val);
251 }
252 }
253 }
254
255 inline
256 void assertIsDefined( const WarpXSolverVec& a_solver_vec ) const
257 {
259 a_solver_vec.IsDefined(),
260 "WarpXSolverVec::function(X) called with undefined WarpXSolverVec X");
261 }
262
263 inline
264 void assertSameType( const WarpXSolverVec& a_solver_vec ) const
265 {
267 a_solver_vec.getArrayVecType()==m_array_type &&
268 a_solver_vec.getScalarVecType()==m_scalar_type,
269 "WarpXSolverVec::function(X) called with WarpXSolverVec X of different type");
270 }
271
272 [[nodiscard]] inline RT norm2 () const
273 {
274 auto const norm = dotProduct(*this);
275 return std::sqrt(norm);
276 }
277
280
283
284 // solver vector types are type warpx::fields::FieldType
285 [[nodiscard]] warpx::fields::FieldType getArrayVecType () const { return m_array_type; }
286 [[nodiscard]] warpx::fields::FieldType getScalarVecType () const { return m_scalar_type; }
287
288 // solver vector type names
289 [[nodiscard]] std::string getVectorType () const { return m_vector_type_name; }
290 [[nodiscard]] std::string getScalarType () const { return m_scalar_type_name; }
291
292 // vector sizes
293 [[nodiscard]] amrex::Long nDOF_local () const
294 {
296 (m_dofs != nullptr),
297 "WarpXSolverVec::nDOF_local() DOF object is a nullptr");
298 return m_dofs->m_nDoFs_l;
299 }
300 [[nodiscard]] amrex::Long nDOF_global () const
301 {
303 (m_dofs != nullptr),
304 "WarpXSolverVec::nDOF_global() DOF object is a nullptr");
305 return m_dofs->m_nDoFs_g;
306 }
307
308 // copy to/from Real-type arrays
309 void copyFrom (const amrex::Real* const);
310 void copyTo (amrex::Real* const) const;
311
312 // return WarpX pointer
313 [[nodiscard]] auto getWarpX () const { return m_WarpX; }
314
315 // return the number of AMR levels
316 [[nodiscard]] auto numAMRLevels () const { return m_num_amr_levels; }
317
318 // return DOFs object pointer
319 inline const auto& getDOFsObject () const { return m_dofs; }
320
321private:
322
323 bool m_is_defined = false;
324
327
330
331 std::string m_vector_type_name = "none";
332 std::string m_scalar_type_name = "none";
333
334 static constexpr int m_ncomp = 1;
336
337 inline static bool m_warpx_ptr_defined = false;
338 inline static WarpX* m_WarpX = nullptr;
339
340 static std::unique_ptr<WarpXSolverDOF> m_dofs;
341};
342
343#endif
#define WARPX_ALWAYS_ASSERT_WITH_MESSAGE(EX, MSG)
Definition TextMsg.H:13
Definition WarpX.H:85
This is a wrapper class around a Vector of pointers to MultiFabs that contains basic math operators a...
Definition WarpXSolverVec.H:58
void assertSameType(const WarpXSolverVec &a_solver_vec) const
Definition WarpXSolverVec.H:264
void operator-=(const WarpXSolverVec &a_solver_vec)
Definition WarpXSolverVec.H:150
amrex::Long nDOF_global() const
Definition WarpXSolverVec.H:300
std::string m_scalar_type_name
Definition WarpXSolverVec.H:332
WarpXSolverVec(WarpXSolverVec &&) noexcept=default
bool m_is_defined
Definition WarpXSolverVec.H:323
WarpXSolverVec()=default
void Copy(const WarpXSolverVec &a_solver_vec)
Definition WarpXSolverVec.H:93
auto numAMRLevels() const
Definition WarpXSolverVec.H:316
std::string getVectorType() const
Definition WarpXSolverVec.H:289
std::string m_vector_type_name
Definition WarpXSolverVec.H:331
warpx::fields::FieldType m_array_type
Definition WarpXSolverVec.H:328
const auto & getDOFsObject() const
Definition WarpXSolverVec.H:319
static constexpr int m_ncomp
Definition WarpXSolverVec.H:334
void assertIsDefined(const WarpXSolverVec &a_solver_vec) const
Definition WarpXSolverVec.H:256
static bool m_warpx_ptr_defined
Definition WarpXSolverVec.H:337
value_type RT
Definition WarpXSolverVec.H:68
static std::unique_ptr< WarpXSolverDOF > m_dofs
Definition WarpXSolverVec.H:340
void linComb(const RT a, const WarpXSolverVec &X, const RT b, const WarpXSolverVec &Y)
Y = a*X + b*Y.
Definition WarpXSolverVec.H:170
void operator+=(const WarpXSolverVec &a_solver_vec)
Definition WarpXSolverVec.H:133
void Define(WarpX *a_WarpX, const std::string &a_vector_type_name, const std::string &a_scalar_type_name="none")
Definition WarpXSolverVec.cpp:24
static WarpX * m_WarpX
Definition WarpXSolverVec.H:338
ablastr::fields::MultiLevelScalarField m_scalar_vec
Definition WarpXSolverVec.H:326
ablastr::fields::MultiLevelVectorField m_array_vec
Definition WarpXSolverVec.H:325
void Define(const WarpXSolverVec &a_solver_vec)
Definition WarpXSolverVec.H:77
amrex::Long nDOF_local() const
Definition WarpXSolverVec.H:293
warpx::fields::FieldType getArrayVecType() const
Definition WarpXSolverVec.H:285
void scale(RT a_a)
Scale Y by a (Y *= a)
Definition WarpXSolverVec.H:217
int m_num_amr_levels
Definition WarpXSolverVec.H:335
bool IsDefined() const
Definition WarpXSolverVec.H:70
void zero()
Definition WarpXSolverVec.H:235
~WarpXSolverVec()
Definition WarpXSolverVec.cpp:13
warpx::fields::FieldType m_scalar_type
Definition WarpXSolverVec.H:329
ablastr::fields::MultiLevelVectorField & getArrayVec()
Definition WarpXSolverVec.H:279
RT dotProduct(const WarpXSolverVec &a_X) const
Definition WarpXSolverVec.cpp:232
auto getWarpX() const
Definition WarpXSolverVec.H:313
amrex::Real value_type
Definition WarpXSolverVec.H:67
warpx::fields::FieldType getScalarVecType() const
Definition WarpXSolverVec.H:286
const ablastr::fields::MultiLevelScalarField & getScalarVec() const
Definition WarpXSolverVec.H:281
void setVal(const RT a_val)
Definition WarpXSolverVec.H:238
ablastr::fields::MultiLevelScalarField & getScalarVec()
Definition WarpXSolverVec.H:282
void copyTo(amrex::Real *const) const
Definition WarpXSolverVec.cpp:187
std::string getScalarType() const
Definition WarpXSolverVec.H:290
const ablastr::fields::MultiLevelVectorField & getArrayVec() const
Definition WarpXSolverVec.H:278
void increment(const WarpXSolverVec &X, const RT a)
Increment Y by a*X (Y += a*X)
Definition WarpXSolverVec.H:195
void copyFrom(const amrex::Real *const)
Definition WarpXSolverVec.cpp:142
RT norm2() const
Definition WarpXSolverVec.H:272
void Copy(FieldType a_array_type, FieldType a_scalar_type=FieldType::None, bool allow_type_mismatch=false)
Definition WarpXSolverVec.cpp:114
WarpXSolverVec & operator=(const WarpXSolverVec &a_solver_vec)=delete
WarpXSolverVec(const WarpXSolverVec &)=delete
__host__ static __device__ constexpr IntVectND< dim > TheZeroVector() noexcept
static void LinComb(MultiFab &dst, Real a, const MultiFab &x, int xcomp, Real b, const MultiFab &y, int ycomp, int dstcomp, int numcomp, int nghost)
static void Saxpy(MultiFab &dst, Real a, const MultiFab &src, int srccomp, int dstcomp, int numcomp, int nghost)
static void Copy(MultiFab &dst, const MultiFab &src, int srccomp, int dstcomp, int numcomp, int nghost)
amrex::Vector< ScalarField > MultiLevelScalarField
Definition MultiFabRegister.H:200
amrex::Vector< VectorField > MultiLevelVectorField
Definition MultiFabRegister.H:208
FieldType
Definition Fields.H:94
@ None
Definition Fields.H:94