blob: f6a297741b9490e1f5f2c92c799a1803898f82f4 (
plain) (
blame)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
|
#ifndef WARPX_DIVBFUNCTOR_H_
#define WARPX_DIVBFUNCTOR_H_
#include "ComputeDiagFunctor.H"
#include <AMReX_BaseFwd.H>
#include <array>
/**
* \brief Functor to compute divB into mf_out.
*/
class
DivBFunctor final : public ComputeDiagFunctor
{
public:
/** Constructor.
* \param[in] arr_mf_src source multifabs (3 elements for x y z).
* \param[in] lev level of multifab.
* \param[in] crse_ratio for interpolating field values from simulation MultiFabs
to the output diagnostic MultiFab, mf_dst.
* \param[in] convertRZmodes2cartesian whether to generate the result in Cartesian coordinates
* (summing over modes)
* \param[in] ncomp Number of component of mf_src to cell-center in dst multifab.
*/
DivBFunctor(std::array<const amrex::MultiFab* const, 3> arr_mf_src, int lev, amrex::IntVect crse_ratio,
bool convertRZmodes2cartesian=true, int ncomp=1);
/** \brief Compute DivB directly into mf_dst.
*
* \param[out] mf_dst output MultiFab where the result is written
* \param[in] dcomp first component of mf_dst in which cell-centered
* data is stored
*/
virtual void operator()(amrex::MultiFab& mf_dst, int dcomp, int /*i_buffer*/) const override;
private:
/** Vector of pointer to source multifab Bx, By, Bz */
std::array<const amrex::MultiFab * const, 3> m_arr_mf_src;
int const m_lev; /**< level on which mf_src is defined (used in cylindrical) */
/**< (for cylindrical) whether to average all modes into 1 comp */
bool m_convertRZmodes2cartesian;
};
#endif // WARPX_DIVBFUNCTOR_H_
|