Skip to content
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

[WIP] Add collisional diagnostic buffers #5543

Draft
wants to merge 10 commits into
base: development
Choose a base branch
from

Conversation

budjensen
Copy link
Contributor

Stemming out of conversation #5381 with @roelof-groenewald, included is initial work on adding multifab buffers to track collisional diagnostics in a MCC-type simulation.

The work on two buffers is in progress--one for ionization and one for total collisional energy loss. Both multifabs have been included in Fields.H and initialized in WarpX.cpp. Relevant flags to turn on MCC tracking have been created in WarpX.H and initialized in WarpX.cpp. The energy loss multifab has been added into the collision protocol, and future commits will add the ionization buffer.

In the end, I aim to set up Python access to these multifabs, with similar functionality to the particle_containers.ParticleBoundaryBufferWrapper() wrappers.

@dpgrote
Copy link
Member

dpgrote commented Jan 8, 2025

This looks good. A comment - these new MultiFabs will automatically be accessible from Python since they are stored in the MultiFab map. You can also add convenience wrappers for them by adding the appropriate wrapper routines in Python/pywarpx/fields.py.

@budjensen
Copy link
Contributor Author

@dpgrote, thanks for the tip!

const auto [ii, jj, kk] = amrex::getParticleCell(p , plo, dxi).dim3();

// store energy lost in the buffer
coll_ionization_arr(ii, jj, kk, 0) += w[ip + np_elec];
Copy link
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I've added everything that I think I'll need for this diagnostic, but am running into the following compilation error:

/home/bj8080/src/warpx/Source/Particles/Collision/BackgroundMCC/BackgroundMCCCollision.cpp:603:67: error:
assignment of read-only location 'coll_ionization_arr.amrex::Array4<const double>::operator()<>(((int)ii), ((int)jj), ((int)kk), 0)'
  603 |                                coll_ionization_arr(ii, jj, kk, 0) += w[ip + np_elec];
      |                                ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~^~~~~~~~~~~~~~~~~~

since coll_ionization_tracking is a constant. When I correct for this by removing const from the initialization of coll_ionization_arr and the earlier line initializing coll_ionization, I get an error:

/home/bj8080/src/warpx/Source/Particles/Collision/BackgroundMCC/BackgroundMCCCollision.cpp:586:67:
error: cannot bind non-const lvalue reference of type 'amrex::Array4<double>&' to an rvalue of type 'amrex::Array4<double>'
  586 |             auto& coll_ionization_arr = coll_ionization[pti].array();
      |                                         ~~~~~~~~~~~~~~~~~~~~~~~~~~^~

I think to get around this I should create a multifab to sum up the ionization events, and then add that multifab to coll_ionization_arr. How would I do this? Or, is there a better approach?

Copy link
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Pinging this with @roelof-groenewald, since we've talked about this feature in discussion 5381. Would you be able to take a look?

Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I think you should do

const auto& coll_ionization_arr = coll_ionization.array(pti);

Does that help?

Copy link
Member

@roelof-groenewald roelof-groenewald left a comment

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Sorry for the slow response (again). Things have been quite hectic since the start of the year 😃

Comment on lines +85 to +86
collision_energy_change, /**< Stores collisional energy transfer data for Monte Carlo Collisions */
collision_ionization /**< Stores ionization data for Monte Carlo Collisions */
Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

For clarity it would be good to add mcc to these multifab names (to make it clear they are tied to the MCC routine). Maybe something like mcc_collision_energy_change, etc.

Comment on lines +438 to +450
if (WarpX::do_MCC_energy_tracking) {
// Calculate the energy lost for tracking
v_lost2 -= ux[ip]*ux[ip] + uy[ip]*uy[ip] + uz[ip]*uz[ip];
ParticleUtils::getEnergy(v_lost2, m, E_lost);

// store the energy (add E_lost * w[ip] to the buffer of energy sent to the background)
// get the particle's cell index (method used here is similar to TemperatureFunctor.cpp)
const auto p = WarpXParticleContainer::ParticleType(ptd, ip);
const auto [ii, jj, kk] = amrex::getParticleCell(p , plo, dxi).dim3();

// store energy lost in the buffer
coll_E_change_arr(ii, jj, kk, 0) += E_lost * w[ip];
}
Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

It would be nice if this code does not need to be duplicated as it is currently. Perhaps the break can be removed, and a else used for the non-charge exchange processes?

Comment on lines +870 to +873
// Control whether Monte Carlo Collision statistics are collected
static bool do_MCC_energy_tracking;
static bool do_MCC_ionization_tracking;

Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

These should instead be made member variables of the BackgroundMCCCollision class.

@@ -2224,6 +2274,25 @@ WarpX::AllocLevelMFs (int lev, const BoxArray& ba, const DistributionMapping& dm
m_fields.alloc_init(FieldType::current_fp, Direction{1}, lev, amrex::convert(ba, jy_nodal_flag), dm, ncomps, ngJ, 0.0_rt);
m_fields.alloc_init(FieldType::current_fp, Direction{2}, lev, amrex::convert(ba, jz_nodal_flag), dm, ncomps, ngJ, 0.0_rt);

// Allocate extra multifabs needed for MCC collision tracking
if (do_MCC_ionization_tracking || do_MCC_energy_tracking) {
Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

These multifabs could also be initialized in the constructor for BackgroundMCCCollision.

const auto& coll_ionization_arr = coll_ionization[pti].array();

// get Struct-Of-Array particle data
// as long as w is set after the copy, will it include newly created particles (be of length np_elec + num_added)?
Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Yes, at this point the length will include the new particles (since this happens after setNewParticleIDs.

const auto p = WarpXParticleContainer::ParticleType(ptd, ip + np_elec);
const auto [ii, jj, kk] = amrex::getParticleCell(p , plo, dxi).dim3();

// store energy lost in the buffer
Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

This buffer is tracking number of electrons created, not energy, right?

const auto [ii, jj, kk] = amrex::getParticleCell(p , plo, dxi).dim3();

// store energy lost in the buffer
coll_ionization_arr(ii, jj, kk, 0) += w[ip + np_elec];
Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I think you should do

const auto& coll_ionization_arr = coll_ionization.array(pti);

Does that help?

Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
Labels
None yet
Projects
None yet
Development

Successfully merging this pull request may close these issues.

3 participants