aboutsummaryrefslogtreecommitdiff
path: root/Source/Diagnostics/ComputeDiagFunctors
diff options
context:
space:
mode:
Diffstat (limited to '')
-rw-r--r--Source/Diagnostics/ComputeDiagFunctors/CMakeLists.txt1
-rw-r--r--Source/Diagnostics/ComputeDiagFunctors/CellCenterFunctor.cpp35
-rw-r--r--Source/Diagnostics/ComputeDiagFunctors/ComputeDiagFunctor.H35
-rw-r--r--Source/Diagnostics/ComputeDiagFunctors/JFunctor.H50
-rw-r--r--Source/Diagnostics/ComputeDiagFunctors/JFunctor.cpp59
-rw-r--r--Source/Diagnostics/ComputeDiagFunctors/Make.package1
-rw-r--r--Source/Diagnostics/ComputeDiagFunctors/RhoFunctor.cpp30
7 files changed, 152 insertions, 59 deletions
diff --git a/Source/Diagnostics/ComputeDiagFunctors/CMakeLists.txt b/Source/Diagnostics/ComputeDiagFunctors/CMakeLists.txt
index 01f780a57..2ef6af16b 100644
--- a/Source/Diagnostics/ComputeDiagFunctors/CMakeLists.txt
+++ b/Source/Diagnostics/ComputeDiagFunctors/CMakeLists.txt
@@ -5,6 +5,7 @@ foreach(D IN LISTS WarpX_DIMS)
CellCenterFunctor.cpp
DivBFunctor.cpp
DivEFunctor.cpp
+ JFunctor.cpp
RhoFunctor.cpp
PartPerCellFunctor.cpp
PartPerGridFunctor.cpp
diff --git a/Source/Diagnostics/ComputeDiagFunctors/CellCenterFunctor.cpp b/Source/Diagnostics/ComputeDiagFunctors/CellCenterFunctor.cpp
index f43714a9c..fe5037131 100644
--- a/Source/Diagnostics/ComputeDiagFunctors/CellCenterFunctor.cpp
+++ b/Source/Diagnostics/ComputeDiagFunctors/CellCenterFunctor.cpp
@@ -1,11 +1,6 @@
#include "CellCenterFunctor.H"
-#include "Utils/TextMsg.H"
-#ifdef WARPX_DIM_RZ
-# include "WarpX.H"
-#endif
-
-#include <ablastr/coarsen/sample.H>
+#include "WarpX.H"
#include <AMReX.H>
#include <AMReX_IntVect.H>
@@ -21,29 +16,7 @@ CellCenterFunctor::CellCenterFunctor(amrex::MultiFab const * mf_src, int lev,
void
CellCenterFunctor::operator()(amrex::MultiFab& mf_dst, int dcomp, const int /*i_buffer*/) const
{
-#ifdef WARPX_DIM_RZ
- if (m_convertRZmodes2cartesian) {
- // In cylindrical geometry, sum real part of all modes of m_mf_src in
- // temporary multifab mf_dst_stag, and cell-center it to mf_dst.
- WARPX_ALWAYS_ASSERT_WITH_MESSAGE(
- nComp()==1,
- "The RZ averaging over modes must write into 1 single component");
- auto& warpx = WarpX::GetInstance();
- amrex::MultiFab mf_dst_stag(m_mf_src->boxArray(), warpx.DistributionMap(m_lev), 1, m_mf_src->nGrowVect());
- // Mode 0
- amrex::MultiFab::Copy(mf_dst_stag, *m_mf_src, 0, 0, 1, m_mf_src->nGrowVect());
- for (int ic=1 ; ic < m_mf_src->nComp() ; ic += 2) {
- // All modes > 0
- amrex::MultiFab::Add(mf_dst_stag, *m_mf_src, ic, 0, 1, m_mf_src->nGrowVect());
- }
- ablastr::coarsen::sample::Coarsen( mf_dst, mf_dst_stag, dcomp, 0, nComp(), 0, m_crse_ratio);
- } else {
- ablastr::coarsen::sample::Coarsen( mf_dst, *m_mf_src, dcomp, 0, nComp(), 0, m_crse_ratio);
- }
-#else
- // In cartesian geometry, coarsen and interpolate from simulation MultiFab, m_mf_src,
- // to output diagnostic MultiFab, mf_dst.
- ablastr::coarsen::sample::Coarsen(mf_dst, *m_mf_src, dcomp, 0, nComp(), mf_dst.nGrowVect(), m_crse_ratio);
- amrex::ignore_unused(m_lev, m_convertRZmodes2cartesian);
-#endif
+ auto& warpx = WarpX::GetInstance();
+ InterpolateMFForDiag(mf_dst, *m_mf_src, dcomp, warpx.DistributionMap(m_lev),
+ m_convertRZmodes2cartesian);
}
diff --git a/Source/Diagnostics/ComputeDiagFunctors/ComputeDiagFunctor.H b/Source/Diagnostics/ComputeDiagFunctors/ComputeDiagFunctor.H
index 1f71732d4..53a9aff8f 100644
--- a/Source/Diagnostics/ComputeDiagFunctors/ComputeDiagFunctor.H
+++ b/Source/Diagnostics/ComputeDiagFunctors/ComputeDiagFunctor.H
@@ -2,6 +2,9 @@
#define WARPX_COMPUTEDIAGFUNCTOR_H_
#include "ComputeDiagFunctor_fwd.H"
+#include "Utils/TextMsg.H"
+
+#include <ablastr/coarsen/sample.H>
#include <AMReX.H>
#include <AMReX_MultiFab.H>
@@ -57,9 +60,41 @@ public:
k_index_zlab, max_box_size, snapshot_full);
}
virtual void InitData() {}
+
+ void InterpolateMFForDiag (
+ amrex::MultiFab& mf_dst, const amrex::MultiFab& mf_src, int dcomp,
+ const amrex::DistributionMapping& dm, bool convertRZmodes2cartesian ) const
+ {
+#ifdef WARPX_DIM_RZ
+ if (convertRZmodes2cartesian) {
+ // In cylindrical geometry, sum real part of all modes of mf_src in
+ // temporary MultiFab mf_dst_stag, and cell-center it to mf_dst
+ WARPX_ALWAYS_ASSERT_WITH_MESSAGE(
+ nComp()==1,
+ "The RZ averaging over modes must write into one single component");
+ amrex::MultiFab mf_dst_stag( mf_src.boxArray(), dm, 1, mf_src.nGrowVect() );
+ // Mode 0
+ amrex::MultiFab::Copy( mf_dst_stag, mf_src, 0, 0, 1, mf_src.nGrowVect() );
+ for (int ic=1 ; ic < mf_src.nComp() ; ic += 2) {
+ // Real part of all modes > 0
+ amrex::MultiFab::Add( mf_dst_stag, mf_src, ic, 0, 1, mf_src.nGrowVect() );
+ }
+ ablastr::coarsen::sample::Coarsen( mf_dst, mf_dst_stag, dcomp, 0, nComp(), 0, m_crse_ratio );
+ } else {
+ ablastr::coarsen::sample::Coarsen( mf_dst, mf_src, dcomp, 0, nComp(), 0, m_crse_ratio );
+ }
+#else
+ // In Cartesian geometry, coarsen and interpolate from temporary MultiFab mf_src
+ // to output diagnostic MultiFab mf_dst
+ ablastr::coarsen::sample::Coarsen(mf_dst, mf_src, dcomp, 0, nComp(), mf_dst.nGrowVect(), m_crse_ratio );
+ amrex::ignore_unused(convertRZmodes2cartesian, dm);
+#endif
+ }
+
private:
/** Number of components of mf_dst that this functor updates. */
int m_ncomp;
+
protected:
/** Coarsening ratio used to interpolate fields from simulation MultiFabs to output MultiFab. */
amrex::IntVect m_crse_ratio;
diff --git a/Source/Diagnostics/ComputeDiagFunctors/JFunctor.H b/Source/Diagnostics/ComputeDiagFunctors/JFunctor.H
new file mode 100644
index 000000000..c3be88bdb
--- /dev/null
+++ b/Source/Diagnostics/ComputeDiagFunctors/JFunctor.H
@@ -0,0 +1,50 @@
+#ifndef WARPX_JFUNCTOR_H_
+#define WARPX_JFUNCTOR_H_
+
+#include "ComputeDiagFunctor.H"
+
+#include <AMReX_BaseFwd.H>
+
+/**
+ * \brief Functor to cell-center MF for current density and store result in mf_out.
+ */
+class
+JFunctor final : public ComputeDiagFunctor
+{
+public:
+ /** Constructor.
+ *
+ * \param[in] dir direction of vector field to operate on
+ * \param[in] lev level of multifab. Used for averaging in rz.
+ * \param[in] crse_ratio for interpolating field values from the simulation MultiFab, src_mf,
+ to the output diagnostic MultiFab, mf_dst.
+ * \param[in] convertRZmodes2cartesian (in cylindrical) whether to
+ * sum all modes in mf_src before cell-centering into dst multifab.
+ * \param[in] ncomp Number of component of mf_src to cell-center in dst multifab.
+ */
+ JFunctor (const int dir, const int lev,
+ const amrex::IntVect crse_ratio,
+ bool convertRZmodes2cartesian=true,
+ bool deposit_current=false, int ncomp=1);
+ /** \brief Cell-center m_mf_src and write the result in mf_dst.
+ *
+ * In cylindrical geometry, by default this functor average all components
+ * of a MultiFab and writes into one single component.
+ *
+ * \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, const int /*i_buffer=0*/) const override;
+private:
+ /** direction of the current density to save */
+ const int m_dir;
+ /** level on which mf_src is defined */
+ int m_lev;
+ /** (for cylindrical) whether to average all modes into 1 comp */
+ bool m_convertRZmodes2cartesian;
+ /** whether to deposit current density before saving */
+ bool m_deposit_current;
+};
+
+#endif // WARPX_JFUNCTOR_H_
diff --git a/Source/Diagnostics/ComputeDiagFunctors/JFunctor.cpp b/Source/Diagnostics/ComputeDiagFunctors/JFunctor.cpp
new file mode 100644
index 000000000..4663a5e1a
--- /dev/null
+++ b/Source/Diagnostics/ComputeDiagFunctors/JFunctor.cpp
@@ -0,0 +1,59 @@
+/* This file is part of WarpX.
+ *
+ * Authors: Roelof Groenewald
+ * License: BSD-3-Clause-LBNL
+ */
+
+#include "JFunctor.H"
+
+#include "Particles/MultiParticleContainer.H"
+#include "WarpX.H"
+
+#include <AMReX.H>
+#include <AMReX_IntVect.H>
+#include <AMReX_MultiFab.H>
+
+JFunctor::JFunctor (const int dir, int lev,
+ amrex::IntVect crse_ratio,
+ bool convertRZmodes2cartesian,
+ bool deposit_current, int ncomp)
+ : ComputeDiagFunctor(ncomp, crse_ratio), m_dir(dir), m_lev(lev),
+ m_convertRZmodes2cartesian(convertRZmodes2cartesian),
+ m_deposit_current(deposit_current)
+{ }
+
+void
+JFunctor::operator() (amrex::MultiFab& mf_dst, int dcomp, const int /*i_buffer*/) const
+{
+ auto& warpx = WarpX::GetInstance();
+ /** pointer to source multifab (can be multi-component) */
+ amrex::MultiFab* m_mf_src = warpx.get_pointer_current_fp(m_lev, m_dir);
+
+ // Deposit current if no solver or the electrostatic solver is being used
+ if (m_deposit_current)
+ {
+ // allocate temporary multifab to deposit current density into
+ amrex::Vector<std::array< std::unique_ptr<amrex::MultiFab>, 3 > > current_fp_temp;
+ current_fp_temp.resize(1);
+
+ auto& current_fp_x = warpx.getcurrent_fp(m_lev,0);
+ current_fp_temp[0][0] = std::make_unique<amrex::MultiFab>(
+ current_fp_x, amrex::make_alias, 0, current_fp_x.nComp()
+ );
+
+ auto& current_fp_y = warpx.getcurrent_fp(m_lev,1);
+ current_fp_temp[0][1] = std::make_unique<amrex::MultiFab>(
+ current_fp_y, amrex::make_alias, 0, current_fp_y.nComp()
+ );
+ auto& current_fp_z = warpx.getcurrent_fp(m_lev,2);
+ current_fp_temp[0][2] = std::make_unique<amrex::MultiFab>(
+ current_fp_z, amrex::make_alias, 0, current_fp_z.nComp()
+ );
+
+ auto& mypc = warpx.GetPartContainer();
+ mypc.DepositCurrent(current_fp_temp, warpx.getdt(m_lev), 0.0);
+ }
+
+ InterpolateMFForDiag(mf_dst, *m_mf_src, dcomp, warpx.DistributionMap(m_lev),
+ m_convertRZmodes2cartesian);
+}
diff --git a/Source/Diagnostics/ComputeDiagFunctors/Make.package b/Source/Diagnostics/ComputeDiagFunctors/Make.package
index 611ccb60d..9c329263b 100644
--- a/Source/Diagnostics/ComputeDiagFunctors/Make.package
+++ b/Source/Diagnostics/ComputeDiagFunctors/Make.package
@@ -3,6 +3,7 @@ CEXE_sources += PartPerCellFunctor.cpp
CEXE_sources += PartPerGridFunctor.cpp
CEXE_sources += DivBFunctor.cpp
CEXE_sources += DivEFunctor.cpp
+CEXE_sources += JFunctor.cpp
CEXE_sources += RhoFunctor.cpp
CEXE_sources += BackTransformFunctor.cpp
CEXE_sources += BackTransformParticleFunctor.cpp
diff --git a/Source/Diagnostics/ComputeDiagFunctors/RhoFunctor.cpp b/Source/Diagnostics/ComputeDiagFunctors/RhoFunctor.cpp
index 6f4fcda4d..73feef4ed 100644
--- a/Source/Diagnostics/ComputeDiagFunctors/RhoFunctor.cpp
+++ b/Source/Diagnostics/ComputeDiagFunctors/RhoFunctor.cpp
@@ -8,11 +8,8 @@
#endif
#include "Particles/MultiParticleContainer.H"
#include "Particles/WarpXParticleContainer.H"
-#include "Utils/TextMsg.H"
#include "WarpX.H"
-#include <ablastr/coarsen/sample.H>
-
#include <AMReX.H>
#include <AMReX_IntVect.H>
#include <AMReX_MultiFab.H>
@@ -69,29 +66,6 @@ RhoFunctor::operator() ( amrex::MultiFab& mf_dst, const int dcomp, const int /*i
}
#endif
-
-#ifdef WARPX_DIM_RZ
- if (m_convertRZmodes2cartesian) {
- // In cylindrical geometry, sum real part of all modes of rho in
- // temporary MultiFab mf_dst_stag, and cell-center it to mf_dst
- WARPX_ALWAYS_ASSERT_WITH_MESSAGE(
- nComp()==1,
- "The RZ averaging over modes must write into one single component");
- amrex::MultiFab mf_dst_stag( rho->boxArray(), warpx.DistributionMap(m_lev), 1, rho->nGrowVect() );
- // Mode 0
- amrex::MultiFab::Copy( mf_dst_stag, *rho, 0, 0, 1, rho->nGrowVect() );
- for (int ic=1 ; ic < rho->nComp() ; ic += 2) {
- // Real part of all modes > 0
- amrex::MultiFab::Add( mf_dst_stag, *rho, ic, 0, 1, rho->nGrowVect() );
- }
- ablastr::coarsen::sample::Coarsen( mf_dst, mf_dst_stag, dcomp, 0, nComp(), 0, m_crse_ratio );
- } else {
- ablastr::coarsen::sample::Coarsen( mf_dst, *rho, dcomp, 0, nComp(), 0, m_crse_ratio );
- }
-#else
- // In Cartesian geometry, coarsen and interpolate from temporary MultiFab rho
- // to output diagnostic MultiFab mf_dst
- ablastr::coarsen::sample::Coarsen(mf_dst, *rho, dcomp, 0, nComp(), mf_dst.nGrowVect(), m_crse_ratio );
- amrex::ignore_unused(m_convertRZmodes2cartesian);
-#endif
+ InterpolateMFForDiag(mf_dst, *rho, dcomp, warpx.DistributionMap(m_lev),
+ m_convertRZmodes2cartesian);
}