-
Notifications
You must be signed in to change notification settings - Fork 201
New issue
Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.
By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.
Already on GitHub? Sign in to your account
Improve rho diagnostic #1099
Improve rho diagnostic #1099
Changes from 10 commits
f5d8e8c
c3889df
a23633d
6fbd0cf
0b1b54a
a169cdd
702ca97
789c9e9
8836a80
beaaa1f
09295a4
1c3ba57
6e337dc
d1666b2
40666d2
File filter
Filter by extension
Conversations
Jump to
Diff view
Diff view
There are no files selected for viewing
Original file line number | Diff line number | Diff line change |
---|---|---|
|
@@ -29,4 +29,4 @@ | |
"particle_position_y": 2.6214399663009544, | ||
"particle_position_z": 2.6214399657276033 | ||
} | ||
} | ||
} |
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,47 @@ | ||
#ifndef WARPX_RHOFUNCTOR_H_ | ||
#define WARPX_RHOFUNCTOR_H_ | ||
|
||
#include "ComputeDiagFunctor.H" | ||
|
||
/** | ||
* \brief Functor to compute charge density rho into mf_out | ||
*/ | ||
class | ||
RhoFunctor final : public ComputeDiagFunctor | ||
{ | ||
|
||
public: | ||
|
||
/** | ||
* \brief Constructor | ||
* | ||
* \param[in] mf_src source MultiFab (must be nullptr as no source MF is | ||
* needed to compute the charge density) | ||
* \param[in] lev level of MultiFab | ||
* \param[in] crse_ratio coarsening ratio for interpolation of field values | ||
* from simulation MultiFabs to the output MultiFab mf_dst | ||
* \param[in] convertRZmodes2cartesian | ||
* \param[in] ncomp optional number of component of source MultiFab mf_src | ||
* to be cell-centered in output MultiFab mf_dst | ||
*/ | ||
RhoFunctor ( const int lev, const amrex::IntVect crse_ratio, | ||
bool convertRZmodes2cartesian=true, const int ncomp=1 ); | ||
|
||
/** | ||
* \brief Compute rho 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 are stored | ||
*/ | ||
virtual void operator() ( amrex::MultiFab& mf_dst, const int dcomp ) const override; | ||
|
||
private: | ||
|
||
// Level on which source MultiFab mf_src is defined in RZ geometry | ||
int const m_lev; | ||
|
||
// Whether to average all modes into 1 component in RZ geometry | ||
bool m_convertRZmodes2cartesian; | ||
}; | ||
|
||
#endif // WARPX_RHOFUNCTOR_H_ |
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,43 @@ | ||
#include "WarpX.H" | ||
#include "RhoFunctor.H" | ||
#include "Utils/CoarsenIO.H" | ||
|
||
RhoFunctor::RhoFunctor ( const int lev, const amrex::IntVect crse_ratio, | ||
bool convertRZmodes2cartesian, const int ncomp ) | ||
: ComputeDiagFunctor(ncomp, crse_ratio), m_lev(lev), | ||
m_convertRZmodes2cartesian(convertRZmodes2cartesian) | ||
{} | ||
|
||
void | ||
RhoFunctor::operator() ( amrex::MultiFab& mf_dst, const int dcomp ) const | ||
{ | ||
auto& warpx = WarpX::GetInstance(); | ||
auto& mypc = warpx.GetPartContainer(); | ||
|
||
// Deposit charge density | ||
std::unique_ptr<amrex::MultiFab> rho = mypc.GetChargeDensity( m_lev ); | ||
|
||
#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 | ||
AMREX_ALWAYS_ASSERT_WITH_MESSAGE( | ||
nComp()==1, | ||
"The RZ averaging over modes must write into 1 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() ); | ||
} | ||
CoarsenIO::Coarsen( mf_dst, mf_dst_stag, dcomp, 0, nComp(), 0, m_crse_ratio ); | ||
} else { | ||
CoarsenIO::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 | ||
CoarsenIO::Coarsen( mf_dst, *rho, dcomp, 0, nComp(), 0, m_crse_ratio ); | ||
#endif | ||
} |
Original file line number | Diff line number | Diff line change |
---|---|---|
|
@@ -6,6 +6,7 @@ | |
#include "ComputeDiagFunctors/PartPerGridFunctor.H" | ||
#include "ComputeDiagFunctors/DivBFunctor.H" | ||
#include "ComputeDiagFunctors/DivEFunctor.H" | ||
#include "ComputeDiagFunctors/RhoFunctor.H" | ||
#include "FlushFormats/FlushFormat.H" | ||
#include "FlushFormats/FlushFormatPlotfile.H" | ||
#include "FlushFormats/FlushFormatCheckpoint.H" | ||
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. Some code needs to be added to the routine FullDiagnostics::AddRZModesToDiags to add the writing of the rho components. Duplicate the code that is there for divE, but writing rho. Note that the check for rho_requested can be simplified to by There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. 👍 looks addressed below now. |
||
|
@@ -48,7 +49,6 @@ FullDiagnostics::ReadParameters () | |
{ | ||
// Read list of full diagnostics fields requested by the user. | ||
bool checkpoint_compatibility = BaseReadParameters(); | ||
auto & warpx = WarpX::GetInstance(); | ||
amrex::ParmParse pp(m_diag_name); | ||
AMREX_ALWAYS_ASSERT_WITH_MESSAGE( | ||
m_format == "plotfile" || m_format == "openpmd" || | ||
|
@@ -146,6 +146,10 @@ FullDiagnostics::AddRZModesToDiags (int lev) | |
} | ||
} | ||
|
||
// Check if rho is requested | ||
// If so, all components will be written out | ||
bool rho_requested = WarpXUtilStr::is_in( m_varnames, "rho" ); | ||
|
||
// First index of m_all_field_functors[lev] where RZ modes are stored | ||
int icomp = m_all_field_functors[0].size(); | ||
const std::array<std::string, 3> coord {"r", "theta", "z"}; | ||
|
@@ -156,6 +160,9 @@ FullDiagnostics::AddRZModesToDiags (int lev) | |
if (divE_requested) { | ||
n_new_fields += 1; | ||
} | ||
if (rho_requested) { | ||
n_new_fields += 1; | ||
} | ||
m_all_field_functors[lev].resize( m_all_field_functors[0].size() + n_new_fields ); | ||
// E | ||
for (int dim=0; dim<3; dim++){ | ||
|
@@ -194,6 +201,12 @@ FullDiagnostics::AddRZModesToDiags (int lev) | |
icomp += 1; | ||
AddRZModesToOutputNames(std::string("divE"), ncomp_multimodefab); | ||
} | ||
// rho | ||
if (rho_requested) { | ||
m_all_field_functors[lev][icomp] = std::make_unique<RhoFunctor>(lev, m_crse_ratio, false, ncomp_multimodefab); | ||
icomp += 1; | ||
AddRZModesToOutputNames(std::string("rho"), ncomp_multimodefab); | ||
} | ||
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. Thanks for adding this! |
||
// Sum the number of components in input vector m_all_field_functors | ||
// and check that it corresponds to the number of components in m_varnames | ||
// and m_mf_output | ||
|
@@ -342,13 +355,7 @@ FullDiagnostics::InitializeFieldFunctors (int lev) | |
} 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 | ||
amrex::MultiFab* rho_new = new amrex::MultiFab(*warpx.get_pointer_rho_fp(lev), amrex::make_alias, 1, 1); | ||
m_all_field_functors[lev][comp] = std::make_unique<CellCenterFunctor>(rho_new, 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 | ||
m_all_field_functors[lev][comp] = std::make_unique<RhoFunctor>(lev, m_crse_ratio); | ||
} 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" ){ | ||
|
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
this params lacks a description :)