aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
-rw-r--r--Source/Diagnostics/ComputeDiagFunctors/CellCenterFunctor.H3
-rw-r--r--Source/Diagnostics/ComputeDiagFunctors/ComputeDiagFunctor.H2
-rw-r--r--Source/Diagnostics/Diagnostics.H7
-rw-r--r--Source/Diagnostics/Diagnostics.cpp117
-rw-r--r--Source/Diagnostics/MultiDiagnostics.H5
-rw-r--r--Source/Diagnostics/MultiDiagnostics.cpp9
-rw-r--r--Source/Parallelization/WarpXRegrid.cpp2
7 files changed, 89 insertions, 56 deletions
diff --git a/Source/Diagnostics/ComputeDiagFunctors/CellCenterFunctor.H b/Source/Diagnostics/ComputeDiagFunctors/CellCenterFunctor.H
index bb153517b..cb4cd7a41 100644
--- a/Source/Diagnostics/ComputeDiagFunctors/CellCenterFunctor.H
+++ b/Source/Diagnostics/ComputeDiagFunctors/CellCenterFunctor.H
@@ -22,7 +22,6 @@ public:
CellCenterFunctor(const amrex::MultiFab * const mf_src, const int lev,
const amrex::IntVect crse_ratio,
bool convertRZmodes2cartesian=true, 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
@@ -35,7 +34,7 @@ public:
virtual void operator()(amrex::MultiFab& mf_dst, int dcomp) const override;
private:
/** pointer to source multifab (can be multi-component) */
- amrex::MultiFab const * const m_mf_src;
+ amrex::MultiFab const * const m_mf_src = nullptr;
int 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;
diff --git a/Source/Diagnostics/ComputeDiagFunctors/ComputeDiagFunctor.H b/Source/Diagnostics/ComputeDiagFunctors/ComputeDiagFunctor.H
index 2bef79a7f..2cacf1e54 100644
--- a/Source/Diagnostics/ComputeDiagFunctors/ComputeDiagFunctor.H
+++ b/Source/Diagnostics/ComputeDiagFunctors/ComputeDiagFunctor.H
@@ -13,6 +13,8 @@ ComputeDiagFunctor
public:
ComputeDiagFunctor( int ncomp, amrex::IntVect crse_ratio) :
m_ncomp(ncomp), m_crse_ratio(crse_ratio) {};
+ //** Virtual Destructor to handle clean destruction of derived classes */
+ virtual ~ComputeDiagFunctor() = default;
/** Compute a field and store the result in mf_dst
* \param[out] mf_dst output MultiFab where the result is written
* \param[in] dcomp first component of mf_dst in which the result is written
diff --git a/Source/Diagnostics/Diagnostics.H b/Source/Diagnostics/Diagnostics.H
index e8d8f6e8c..d028aca1a 100644
--- a/Source/Diagnostics/Diagnostics.H
+++ b/Source/Diagnostics/Diagnostics.H
@@ -37,6 +37,10 @@ public:
* \return bool, whether to flush
*/
bool DoDump (int step, bool force_flush=false);
+ /** Initialize functors that store pointers to the fields requested by the user.
+ * \param[in] lev level on which the vector of unique_ptrs to field functors is initialized.
+ */
+ void InitializeFieldFunctors (int lev);
private:
void ReadParameters ();
/** Append varnames with names for all modes of a field
@@ -66,10 +70,11 @@ private:
* cell-center for standard EB fields) as well as more involved operations
* (back-transformed diagnostics, filtering, reconstructing cartesian
* fields in cylindrical). */
- amrex::Vector< amrex::Vector <ComputeDiagFunctor const *> > m_all_field_functors;
+ amrex::Vector< amrex::Vector <std::unique_ptr<ComputeDiagFunctor const> > > m_all_field_functors;
/** output multifab, where all fields are cell-centered and stacked */
amrex::Vector< amrex::MultiFab > m_mf_output;
int nlev; /**< number of levels to output */
+ int nmax_lev; /**< max_level to allocate output multifab and vector of field functors. */
/** This class is responsible for flushing the data to file */
FlushFormat* m_flush_format;
/** Whether to plot raw (i.e., NOT cell-centered) fields */
diff --git a/Source/Diagnostics/Diagnostics.cpp b/Source/Diagnostics/Diagnostics.cpp
index a3a1c2889..26d737895 100644
--- a/Source/Diagnostics/Diagnostics.cpp
+++ b/Source/Diagnostics/Diagnostics.cpp
@@ -90,56 +90,15 @@ Diagnostics::InitData ()
{
Print()<<"Diagnostics::InitData\n";
auto & warpx = WarpX::GetInstance();
+ // Number of levels
nlev = warpx.finestLevel() + 1;
- // Initialize vector of pointers to the fields requested by the user.
- m_all_field_functors.resize( nlev );
- m_mf_output.resize( nlev );
+ // Maximum number of levels that will be allocated in the simulation
+ nmax_lev = warpx.maxLevel() + 1;
+ m_mf_output.resize( nmax_lev );
+ m_all_field_functors.resize( nmax_lev );
for ( int lev=0; lev<nlev; lev++ ){
- m_all_field_functors[lev].resize( m_varnames.size() );
- // Fill vector of functors for all components except individual cylindrical modes.
- for (int comp=0, n=m_all_field_functors[lev].size(); comp<n; comp++){
- if ( m_varnames[comp] == "Ex" ){
- m_all_field_functors[lev][comp] = new CellCenterFunctor(warpx.get_pointer_Efield_aux(lev, 0), lev, m_crse_ratio);
- } else if ( m_varnames[comp] == "Ey" ){
- m_all_field_functors[lev][comp] = new CellCenterFunctor(warpx.get_pointer_Efield_aux(lev, 1), lev, m_crse_ratio);
- } else if ( m_varnames[comp] == "Ez" ){
- m_all_field_functors[lev][comp] = new CellCenterFunctor(warpx.get_pointer_Efield_aux(lev, 2), lev, m_crse_ratio);
- } else if ( m_varnames[comp] == "Bx" ){
- m_all_field_functors[lev][comp] = new CellCenterFunctor(warpx.get_pointer_Bfield_aux(lev, 0), lev, m_crse_ratio);
- } else if ( m_varnames[comp] == "By" ){
- m_all_field_functors[lev][comp] = new CellCenterFunctor(warpx.get_pointer_Bfield_aux(lev, 1), lev, m_crse_ratio);
- } else if ( m_varnames[comp] == "Bz" ){
- m_all_field_functors[lev][comp] = new CellCenterFunctor(warpx.get_pointer_Bfield_aux(lev, 2), lev, m_crse_ratio);
- } else if ( m_varnames[comp] == "jx" ){
- m_all_field_functors[lev][comp] = new CellCenterFunctor(warpx.get_pointer_current_fp(lev, 0), lev, m_crse_ratio);
- } else if ( m_varnames[comp] == "jy" ){
- m_all_field_functors[lev][comp] = new CellCenterFunctor(warpx.get_pointer_current_fp(lev, 1), lev, m_crse_ratio);
- } else if ( m_varnames[comp] == "jz" ){
- m_all_field_functors[lev][comp] = new CellCenterFunctor(warpx.get_pointer_current_fp(lev, 2), lev, m_crse_ratio);
- } else if ( m_varnames[comp] == "rho" ){
- // rho_new is stored in component 1 of rho_fp when using PSATD
-#ifdef WARPX_USE_PSATD
- MultiFab* rho_new = new MultiFab(*warpx.get_pointer_rho_fp(lev), amrex::make_alias, 1, 1);
- m_all_field_functors[lev][comp] = new CellCenterFunctor(rho_new, lev, m_crse_ratio);
-#else
- m_all_field_functors[lev][comp] = new CellCenterFunctor(warpx.get_pointer_rho_fp(lev), lev, m_crse_ratio);
-#endif
- } else if ( m_varnames[comp] == "F" ){
- m_all_field_functors[lev][comp] = new CellCenterFunctor(warpx.get_pointer_F_fp(lev), lev, m_crse_ratio);
- } else if ( m_varnames[comp] == "part_per_cell" ){
- m_all_field_functors[lev][comp] = new PartPerCellFunctor(nullptr, lev, m_crse_ratio);
- } else if ( m_varnames[comp] == "part_per_grid" ){
- m_all_field_functors[lev][comp] = new PartPerGridFunctor(nullptr, lev, m_crse_ratio);
- } else if ( m_varnames[comp] == "divB" ){
- m_all_field_functors[lev][comp] = new DivBFunctor(warpx.get_array_Bfield_aux(lev), lev, m_crse_ratio);
- } else if ( m_varnames[comp] == "divE" ){
- m_all_field_functors[lev][comp] = new DivEFunctor(warpx.get_array_Efield_aux(lev), lev, m_crse_ratio);
- }
- }
-
- AddRZModesToDiags( lev );
-
+ InitializeFieldFunctors( lev );
// At this point, m_varnames.size() >= m_all_field_functors[0].size()
// Initialize member variable m_mf_output depending on m_crse_ratio, m_lo and m_hi
@@ -248,8 +207,8 @@ Diagnostics::AddRZModesToDiags (int lev)
// E
for (int dim=0; dim<3; dim++){
// 3 components, r theta z
- m_all_field_functors[lev][icomp] = new
- CellCenterFunctor(warpx.get_pointer_Efield_aux(lev, dim), lev,
+ m_all_field_functors[lev][icomp] =
+ std::make_unique<CellCenterFunctor>(warpx.get_pointer_Efield_aux(lev, dim), lev,
m_crse_ratio, false, ncomp_multimodefab);
AddRZModesToOutputNames(std::string("E") + coord[dim],
warpx.get_pointer_Efield_aux(0, 0)->nComp());
@@ -258,8 +217,8 @@ Diagnostics::AddRZModesToDiags (int lev)
// B
for (int dim=0; dim<3; dim++){
// 3 components, r theta z
- m_all_field_functors[lev][icomp] = new
- CellCenterFunctor(warpx.get_pointer_Bfield_aux(lev, dim), lev,
+ m_all_field_functors[lev][icomp] =
+ std::make_unique<CellCenterFunctor>(warpx.get_pointer_Bfield_aux(lev, dim), lev,
m_crse_ratio, false, ncomp_multimodefab);
AddRZModesToOutputNames(std::string("B") + coord[dim],
warpx.get_pointer_Bfield_aux(0, 0)->nComp());
@@ -268,8 +227,8 @@ Diagnostics::AddRZModesToDiags (int lev)
// j
for (int dim=0; dim<3; dim++){
// 3 components, r theta z
- m_all_field_functors[lev][icomp] = new
- CellCenterFunctor(warpx.get_pointer_current_fp(lev, dim), lev,
+ m_all_field_functors[lev][icomp] =
+ std::make_unique<CellCenterFunctor>(warpx.get_pointer_current_fp(lev, dim), lev,
m_crse_ratio, false, ncomp_multimodefab);
icomp += 1;
AddRZModesToOutputNames(std::string("J") + coord[dim],
@@ -388,3 +347,55 @@ Diagnostics::DefineDiagMultiFab ( int lev ) {
// Allocate output MultiFab for diagnostics. The data will be stored at cell-centers.
m_mf_output[lev] = MultiFab(ba, dmap, m_varnames.size(), 0);
}
+
+
+void
+Diagnostics::InitializeFieldFunctors (int lev)
+{
+ auto & warpx = WarpX::GetInstance();
+ // Clear any pre-existing vector to release stored data.
+ m_all_field_functors[lev].clear();
+
+ m_all_field_functors[lev].resize( m_varnames.size() );
+ // Fill vector of functors for all components except individual cylindrical modes.
+ for (int comp=0, n=m_all_field_functors[lev].size(); comp<n; comp++){
+ if ( m_varnames[comp] == "Ex" ){
+ m_all_field_functors[lev][comp] = std::make_unique<CellCenterFunctor>(warpx.get_pointer_Efield_aux(lev, 0), lev, m_crse_ratio);
+ } else if ( m_varnames[comp] == "Ey" ){
+ m_all_field_functors[lev][comp] = std::make_unique<CellCenterFunctor>(warpx.get_pointer_Efield_aux(lev, 1), lev, m_crse_ratio);
+ } else if ( m_varnames[comp] == "Ez" ){
+ m_all_field_functors[lev][comp] = std::make_unique<CellCenterFunctor>(warpx.get_pointer_Efield_aux(lev, 2), lev, m_crse_ratio);
+ } else if ( m_varnames[comp] == "Bx" ){
+ m_all_field_functors[lev][comp] = std::make_unique<CellCenterFunctor>(warpx.get_pointer_Bfield_aux(lev, 0), lev, m_crse_ratio);
+ } else if ( m_varnames[comp] == "By" ){
+ m_all_field_functors[lev][comp] = std::make_unique<CellCenterFunctor>(warpx.get_pointer_Bfield_aux(lev, 1), lev, m_crse_ratio);
+ } else if ( m_varnames[comp] == "Bz" ){
+ m_all_field_functors[lev][comp] = std::make_unique<CellCenterFunctor>(warpx.get_pointer_Bfield_aux(lev, 2), lev, m_crse_ratio);
+ } else if ( m_varnames[comp] == "jx" ){
+ m_all_field_functors[lev][comp] = std::make_unique<CellCenterFunctor>(warpx.get_pointer_current_fp(lev, 0), lev, m_crse_ratio);
+ } else if ( m_varnames[comp] == "jy" ){
+ m_all_field_functors[lev][comp] = std::make_unique<CellCenterFunctor>(warpx.get_pointer_current_fp(lev, 1), lev, m_crse_ratio);
+ } else if ( m_varnames[comp] == "jz" ){
+ m_all_field_functors[lev][comp] = std::make_unique<CellCenterFunctor>(warpx.get_pointer_current_fp(lev, 2), lev, m_crse_ratio);
+ } else if ( m_varnames[comp] == "rho" ){
+ // rho_new is stored in component 1 of rho_fp when using PSATD
+#ifdef WARPX_USE_PSATD
+ std::unique_ptr<MultiFab> rho_new = std::make_unique<MultiFab>(*warpx.get_pointer_rho_fp(lev), amrex::make_alias, 1, 1);
+ m_all_field_functors[lev][comp] = std::make_unique<CellCenterFunctor>(rho_new.get(), lev, m_crse_ratio);
+#else
+ m_all_field_functors[lev][comp] = std::make_unique<CellCenterFunctor>(warpx.get_pointer_rho_fp(lev), lev, m_crse_ratio);
+#endif
+ } else if ( m_varnames[comp] == "F" ){
+ m_all_field_functors[lev][comp] = std::make_unique<CellCenterFunctor>(warpx.get_pointer_F_fp(lev), lev, m_crse_ratio);
+ } else if ( m_varnames[comp] == "part_per_cell" ){
+ m_all_field_functors[lev][comp] = std::make_unique<PartPerCellFunctor>(nullptr, lev, m_crse_ratio);
+ } else if ( m_varnames[comp] == "part_per_grid" ){
+ m_all_field_functors[lev][comp] = std::make_unique<PartPerGridFunctor>(nullptr, lev, m_crse_ratio);
+ } else if ( m_varnames[comp] == "divB" ){
+ m_all_field_functors[lev][comp] = std::make_unique<DivBFunctor>(warpx.get_array_Bfield_aux(lev), lev, m_crse_ratio);
+ } else if ( m_varnames[comp] == "divE" ){
+ m_all_field_functors[lev][comp] = std::make_unique<DivEFunctor>(warpx.get_array_Efield_aux(lev), lev, m_crse_ratio);
+ }
+ }
+ AddRZModesToDiags( lev );
+}
diff --git a/Source/Diagnostics/MultiDiagnostics.H b/Source/Diagnostics/MultiDiagnostics.H
index 8313167da..5ac2172a8 100644
--- a/Source/Diagnostics/MultiDiagnostics.H
+++ b/Source/Diagnostics/MultiDiagnostics.H
@@ -19,6 +19,11 @@ public:
void InitData ();
/** \brief Called at each iteration. Compute diags and flush. */
void FilterComputePackFlush (int step, bool force_flush=false);
+ /** \brief Loop over diags in all diags and call their InitializeFieldFunctors.
+ Called when a new partitioning is generated at level, lev.
+ * \param[in] lev level at this the field functors are initialized.
+ */
+ void InitializeFieldFunctors (int lev);
private:
/** Vector of pointers to all diagnostics */
amrex::Vector<std::unique_ptr<Diagnostics> > alldiags;
diff --git a/Source/Diagnostics/MultiDiagnostics.cpp b/Source/Diagnostics/MultiDiagnostics.cpp
index d1100d871..2a66c9b25 100644
--- a/Source/Diagnostics/MultiDiagnostics.cpp
+++ b/Source/Diagnostics/MultiDiagnostics.cpp
@@ -28,6 +28,15 @@ MultiDiagnostics::InitData ()
}
void
+MultiDiagnostics::InitializeFieldFunctors ( int lev )
+{
+ for( auto& diag : alldiags ){
+ // Initialize functors to store pointers to fields.
+ diag->InitializeFieldFunctors( lev );
+ }
+}
+
+void
MultiDiagnostics::ReadParameters ()
{
ParmParse pp("diagnostics");
diff --git a/Source/Parallelization/WarpXRegrid.cpp b/Source/Parallelization/WarpXRegrid.cpp
index bb5e5d9d8..ac8c81934 100644
--- a/Source/Parallelization/WarpXRegrid.cpp
+++ b/Source/Parallelization/WarpXRegrid.cpp
@@ -268,6 +268,8 @@ WarpX::RemakeLevel (int lev, Real /*time*/, const BoxArray& ba, const Distributi
{
amrex::Abort("RemakeLevel: to be implemented");
}
+ // Re-initialize diagnostic functors that stores pointers to the user-requested fields at level, lev.
+ multi_diags->InitializeFieldFunctors( lev );
}
void