Skip to content

Commit

Permalink
add diagnostic for processor number
Browse files Browse the repository at this point in the history
  • Loading branch information
hklion authored and RevathiJambunathan committed Mar 20, 2024
1 parent 29a0a56 commit 3feb26a
Show file tree
Hide file tree
Showing 4 changed files with 86 additions and 0 deletions.
1 change: 1 addition & 0 deletions Source/Diagnostics/ComputeDiagFunctors/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -12,5 +12,6 @@ foreach(D IN LISTS WarpX_DIMS)
BackTransformFunctor.cpp
BackTransformParticleFunctor.cpp
ParticleReductionFunctor.cpp
ProcessorMapFunctor.cpp
)
endforeach()
45 changes: 45 additions & 0 deletions Source/Diagnostics/ComputeDiagFunctors/ProcessorMapFunctor.H
Original file line number Diff line number Diff line change
@@ -0,0 +1,45 @@
#ifndef WARPX_PROCESSORMAPFUNCTOR_H_
#define WARPX_PROCESSORMAPFUNCTOR_H_

#include "ComputeDiagFunctor.H"

#include <AMReX_BaseFwd.H>

/**
* \brief Functor to cell-center MF and store result in mf_out.
*/
class ProcessorMapFunctor final : public ComputeDiagFunctor
{
public:
/** Constructor.
*
* \param[in] mf_src source multifab.
* \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.
*/
ProcessorMapFunctor(const amrex::MultiFab * mf_src, int lev,
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
* 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
*/
void operator()(amrex::MultiFab& mf_dst, int dcomp, int /*i_buffer=0*/) const override;
private:
/** pointer to source multifab (can be multi-component) */
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;
};

#endif // WARPX_CELLCENTERFUNCTOR_H_
37 changes: 37 additions & 0 deletions Source/Diagnostics/ComputeDiagFunctors/ProcessorMapFunctor.cpp
Original file line number Diff line number Diff line change
@@ -0,0 +1,37 @@
#include "ProcessorMapFunctor.H"

#include "WarpX.H"

#include <AMReX.H>
#include <AMReX_IntVect.H>
#include <AMReX_MultiFab.H>

ProcessorMapFunctor::ProcessorMapFunctor(amrex::MultiFab const * mf_src, int lev,
amrex::IntVect crse_ratio,
bool convertRZmodes2cartesian, int ncomp)
: ComputeDiagFunctor(ncomp, crse_ratio), m_mf_src(mf_src), m_lev(lev),
m_convertRZmodes2cartesian(convertRZmodes2cartesian)
{}

void
ProcessorMapFunctor::operator()(amrex::MultiFab& mf_dst, int dcomp, const int /*i_buffer*/) const
{
std::unique_ptr<amrex::MultiFab> tmp;
tmp = std::make_unique<amrex::MultiFab>(m_mf_src->boxArray(), m_mf_src->DistributionMap(), 1, m_mf_src->nGrowVect());

auto& warpx = WarpX::GetInstance();
const amrex::DistributionMapping& dm = warpx.DistributionMap(m_lev);
for (amrex::MFIter mfi(*tmp, false); mfi.isValid(); ++mfi)
{
auto bx = mfi.growntilebox(tmp->nGrowVect());
auto arr = tmp->array(mfi);
int iproc = dm[mfi.index()];
amrex::ParallelFor(bx,
[=] AMREX_GPU_DEVICE(int i, int j, int k)
{
arr(i,j,k) = iproc;
});
}
InterpolateMFForDiag(mf_dst, *tmp, dcomp, warpx.DistributionMap(m_lev),
m_convertRZmodes2cartesian);
}
3 changes: 3 additions & 0 deletions Source/Diagnostics/FullDiagnostics.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -7,6 +7,7 @@
#include "ComputeDiagFunctors/PartPerGridFunctor.H"
#include "ComputeDiagFunctors/ParticleReductionFunctor.H"
#include "ComputeDiagFunctors/RhoFunctor.H"
#include "ComputeDiagFunctors/ProcessorMapFunctor.H"
#include "Diagnostics/Diagnostics.H"
#include "Diagnostics/ParticleDiag/ParticleDiag.H"
#include "FlushFormats/FlushFormat.H"
Expand Down Expand Up @@ -659,6 +660,8 @@ FullDiagnostics::InitializeFieldFunctors (int lev)
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);
} else if ( m_varnames[comp] == "proc_number" ){
m_all_field_functors[lev][comp] = std::make_unique<ProcessorMapFunctor>(warpx.get_pointer_Efield_aux(lev,0), lev, m_crse_ratio);
}
else {

Expand Down

0 comments on commit 3feb26a

Please sign in to comment.