aboutsummaryrefslogtreecommitdiff
path: root/Source/Python
diff options
context:
space:
mode:
Diffstat (limited to 'Source/Python')
-rw-r--r--Source/Python/WarpXWrappers.H10
-rw-r--r--Source/Python/WarpXWrappers.cpp31
2 files changed, 41 insertions, 0 deletions
diff --git a/Source/Python/WarpXWrappers.H b/Source/Python/WarpXWrappers.H
index 8d2ee4364..bd47a3ebb 100644
--- a/Source/Python/WarpXWrappers.H
+++ b/Source/Python/WarpXWrappers.H
@@ -110,6 +110,16 @@ extern "C" {
void warpx_clearParticleBoundaryBuffer ();
+ /**
+ * \brief This function is used to deposit a given species' charge density
+ * in the rho_fp multifab which can then be accessed from python via
+ * pywarpx.fields.RhoFPWrapper()
+ *
+ * @param[in] species_name specifying the name of the species to deposit
+ * @param[in] lev mesh refinement level
+ */
+ void warpx_depositChargeDensity (const char* species_name, int lev);
+
void warpx_ComputeDt ();
void warpx_MoveWindow (int step, bool move_j);
diff --git a/Source/Python/WarpXWrappers.cpp b/Source/Python/WarpXWrappers.cpp
index 7350bbbb6..061e80c6a 100644
--- a/Source/Python/WarpXWrappers.cpp
+++ b/Source/Python/WarpXWrappers.cpp
@@ -606,6 +606,37 @@ namespace
particle_buffers.clearParticles();
}
+ void warpx_depositChargeDensity (const char* char_species_name, int lev) {
+ // this function is used to deposit a given species' charge density
+ // in the rho_fp multifab which can then be accessed from python via
+ // pywarpx.fields.RhoFPWrapper()
+ WarpX& warpx = WarpX::GetInstance();
+ const auto & mypc = warpx.GetPartContainer();
+ const std::string species_name(char_species_name);
+ auto & myspc = mypc.GetParticleContainerFromName(species_name);
+ auto * rho_fp = warpx.get_pointer_rho_fp(lev);
+
+ if (rho_fp == nullptr) {
+ warpx.RecordWarning(
+ "WarpXWrappers", "rho_fp is not allocated", WarnPriority::low
+ );
+ return;
+ }
+
+ // reset rho before depositing
+ rho_fp->setVal(0.);
+
+ for (WarpXParIter pti(myspc, lev); pti.isValid(); ++pti)
+ {
+ const long np = pti.numParticles();
+ auto& wp = pti.GetAttribs(PIdx::w);
+ myspc.DepositCharge(pti, wp, nullptr, rho_fp, 0, 0, np, 0, lev, lev);
+ }
+#ifdef WARPX_DIM_RZ
+ warpx.ApplyInverseVolumeScalingToChargeDensity(rho_fp, lev);
+#endif
+ }
+
void warpx_ComputeDt () {
WarpX& warpx = WarpX::GetInstance();
warpx.ComputeDt ();