Skip to content

Commit

Permalink
remove calls to WarpX:GetInstance() from MacroscopicProperties.cpp (#…
Browse files Browse the repository at this point in the history
  • Loading branch information
lucafedeli88 authored Apr 24, 2024
1 parent 28434e9 commit 7e65117
Show file tree
Hide file tree
Showing 3 changed files with 48 additions and 23 deletions.
Original file line number Diff line number Diff line change
Expand Up @@ -13,8 +13,11 @@
#include "Utils/WarpXConst.H"

#include <AMReX_Array.H>
#include <AMReX_BoxArray.H>
#include <AMReX_DistributionMapping.H>
#include <AMReX_Extension.H>
#include <AMReX_GpuQualifiers.H>
#include <AMReX_IntVect.H>
#include <AMReX_MultiFab.H>
#include <AMReX_Parser.H>
#include <AMReX_REAL.H>
Expand All @@ -33,8 +36,26 @@ public:
MacroscopicProperties (); // constructor
/** Read user-defined macroscopic properties. Called in constructor. */
void ReadParameters ();
/** Initialize multifabs storing macroscopic multifabs */
void InitData ();

/**
* \brief Initialize multifabs storing macroscopic multifabs
*
* @param[in] ba the box array associated to the multifabs E and B
* @param[in] dmap the distribution mapping
* @param[in] ng_EB_alloc guard cells allocated for multifabs E and B
* @param[in] geom the geometry
* @param[in] Ex_stag staggering of the Ex field
* @param[in] Ey_stag staggering of the Ey field
* @param[in] Ez_stag staggering of the Ez field
*/
void InitData (
const amrex::BoxArray& ba,
const amrex::DistributionMapping& dmap,
const amrex::IntVect& ng_EB_alloc,
const amrex::Geometry& geom,
const amrex::IntVect& Ex_stag,
const amrex::IntVect& Ey_stag,
const amrex::IntVect& Ez_stag);

/** return MultiFab, sigma (conductivity) of the medium. */
amrex::MultiFab& getsigma_mf () {return (*m_sigma_mf);}
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -3,18 +3,14 @@
#include "FieldSolver/Fields.H"
#include "Utils/Parser/ParserUtils.H"
#include "Utils/TextMsg.H"
#include "WarpX.H"

#include <ablastr/warn_manager/WarnManager.H>

#include <AMReX_Array4.H>
#include <AMReX_BoxArray.H>
#include <AMReX_Config.H>
#include <AMReX_DistributionMapping.H>
#include <AMReX_Geometry.H>
#include <AMReX_GpuLaunch.H>
#include <AMReX_IndexType.H>
#include <AMReX_IntVect.H>
#include <AMReX_MFIter.H>
#include <AMReX_ParmParse.H>
#include <AMReX_Print.H>
Expand Down Expand Up @@ -123,16 +119,17 @@ MacroscopicProperties::ReadParameters ()
}

void
MacroscopicProperties::InitData ()
MacroscopicProperties::InitData (
const amrex::BoxArray& ba,
const amrex::DistributionMapping& dmap,
const amrex::IntVect& ng_EB_alloc,
const amrex::Geometry& geom,
const amrex::IntVect& Ex_stag,
const amrex::IntVect& Ey_stag,
const amrex::IntVect& Ez_stag)
{
amrex::Print() << Utils::TextMsg::Info("we are in init data of macro");
auto & warpx = WarpX::GetInstance();

// Get BoxArray and DistributionMap of warpx instance.
const int lev = 0;
const amrex::BoxArray ba = warpx.boxArray(lev);
const amrex::DistributionMapping dmap = warpx.DistributionMap(lev);
const amrex::IntVect ng_EB_alloc = warpx.getngEB();
// Define material property multifabs using ba and dmap from WarpX instance
// sigma is cell-centered MultiFab
m_sigma_mf = std::make_unique<amrex::MultiFab>(ba, dmap, 1, ng_EB_alloc);
Expand All @@ -148,7 +145,7 @@ MacroscopicProperties::InitData ()
} else if (m_sigma_s == "parse_sigma_function") {

InitializeMacroMultiFabUsingParser(m_sigma_mf.get(), m_sigma_parser->compile<3>(),
warpx.Geom(lev).CellSizeArray(), warpx.Geom(lev).ProbDomain());
geom.CellSizeArray(), geom.ProbDomain());
}
// Initialize epsilon
if (m_epsilon_s == "constant") {
Expand All @@ -158,7 +155,7 @@ MacroscopicProperties::InitData ()
} else if (m_epsilon_s == "parse_epsilon_function") {

InitializeMacroMultiFabUsingParser(m_eps_mf.get(), m_epsilon_parser->compile<3>(),
warpx.Geom(lev).CellSizeArray(), warpx.Geom(lev).ProbDomain());
geom.CellSizeArray(), geom.ProbDomain());

}
// In the Maxwell solver, `epsilon` is used in the denominator.
Expand All @@ -175,16 +172,14 @@ MacroscopicProperties::InitData ()
} else if (m_mu_s == "parse_mu_function") {

InitializeMacroMultiFabUsingParser(m_mu_mf.get(), m_mu_parser->compile<3>(),
warpx.Geom(lev).CellSizeArray(), warpx.Geom(lev).ProbDomain());
geom.CellSizeArray(), geom.ProbDomain());

}

amrex::IntVect sigma_stag = m_sigma_mf->ixType().toIntVect();
amrex::IntVect epsilon_stag = m_eps_mf->ixType().toIntVect();
amrex::IntVect mu_stag = m_mu_mf->ixType().toIntVect();
amrex::IntVect Ex_stag = warpx.getField(FieldType::Efield_fp, 0,0).ixType().toIntVect();
amrex::IntVect Ey_stag = warpx.getField(FieldType::Efield_fp, 0,1).ixType().toIntVect();
amrex::IntVect Ez_stag = warpx.getField(FieldType::Efield_fp, 0,2).ixType().toIntVect();
const amrex::IntVect sigma_stag = m_sigma_mf->ixType().toIntVect();
const amrex::IntVect epsilon_stag = m_eps_mf->ixType().toIntVect();
const amrex::IntVect mu_stag = m_mu_mf->ixType().toIntVect();


for ( int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
sigma_IndexType[idim] = sigma_stag[idim];
Expand Down
11 changes: 10 additions & 1 deletion Source/Initialization/WarpXInitData.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -474,7 +474,16 @@ WarpX::InitData ()
BuildBufferMasks();

if (WarpX::em_solver_medium==1) {
m_macroscopic_properties->InitData();
const int lev_zero = 0;
m_macroscopic_properties->InitData(
boxArray(lev_zero),
DistributionMap(lev_zero),
getngEB(),
Geom(lev_zero),
getField(FieldType::Efield_fp, lev_zero,0).ixType().toIntVect(),
getField(FieldType::Efield_fp, lev_zero,1).ixType().toIntVect(),
getField(FieldType::Efield_fp, lev_zero,2).ixType().toIntVect()
);
}

if (WarpX::electromagnetic_solver_id == ElectromagneticSolverAlgo::HybridPIC) {
Expand Down

0 comments on commit 7e65117

Please sign in to comment.