1 #ifndef SUPERDROPLET_PC_H_
2 #define SUPERDROPLET_PC_H_
4 #ifdef ERF_USE_PARTICLES
11 #include <AMReX_StructOfArrays.H>
13 class SuperDropletPC :
public ERFPC
15 using MFPtr = std::unique_ptr<amrex::MultiFab>;
16 using BCTypeArr = amrex::GpuArray<ERF_BC, AMREX_SPACEDIM*2>;
21 SuperDropletPC ( amrex::ParGDBBase* a_gdb,
22 const std::vector<Species::Name>& a_species_mat,
23 const std::vector<Species::Name>& a_aerosol_mat,
25 const std::string& a_name =
"super_droplets" )
26 : ERFPC (a_gdb, a_name)
28 define( a_species_mat,
30 a_gdb->ParticleBoxArray(m_lev),
31 a_gdb->ParticleDistributionMap(m_lev),
36 SuperDropletPC (
const amrex::Geometry& a_geom,
37 const amrex::DistributionMapping& a_dmap,
38 const amrex::BoxArray& a_ba,
39 const std::vector<Species::Name>& a_species_mat,
40 const std::vector<Species::Name>& a_aerosol_mat,
42 const std::string& a_name =
"super_droplets" )
43 : ERFPC (a_geom, a_dmap, a_ba, a_name)
45 define(a_species_mat, a_aerosol_mat, a_ba, a_dmap, a_dt);
51 if (m_mass_change_logging) {
52 fclose(m_mass_change_log);
58 virtual void setSpeciesMaterial (
const Species::Name& a_name )
60 if (a_name == Species::Name::H2O) { m_idx_w = m_species_mat.size(); }
61 m_species_mat.push_back(std::make_unique<MaterialProperties>(a_name));
62 m_device_props_initialized =
false;
67 virtual void setSpeciesMaterial (
const std::vector<Species::Name>& a_names )
69 for (
auto& name : a_names) { setSpeciesMaterial(name); }
74 virtual const MaterialProperties& getSpeciesMaterial(
const Species::Name& a_name)
const
76 for (
auto& species : m_species_mat) {
77 if (species->m_name == a_name) {
return *species; }
79 amrex::Abort(
"SuperDropletPC::getSpeciesMaterial() - species not found");
80 return *m_species_mat[0];
85 virtual void setAerosolMaterial (
const Species::Name& a_name )
87 m_aerosol_mat.push_back(std::make_unique<MaterialProperties>(a_name));
88 m_device_props_initialized =
false;
93 virtual void setAerosolMaterial (
const std::vector<Species::Name>& a_names )
95 for (
auto& name : a_names) { setAerosolMaterial(name); }
99 void updateDeviceProperties();
102 template<
typename SOAType>
103 void setupMassPointers(SOAType& soa, SDPCDefn::SDSpeciesMassArr& sp_mass_ptrs,
104 SDPCDefn::SDAerosolMassArr& ae_mass_ptrs)
const
107 sp_mass_ptrs[i] = soa.GetRealData(idx_s(i, m_num_aerosols,
m_num_species)).data();
109 for (
int i = 0; i < m_num_aerosols; i++) {
110 ae_mass_ptrs[i] = soa.GetRealData(idx_a(i, m_num_aerosols,
m_num_species)).data();
115 AMREX_GPU_DEVICE AMREX_FORCE_INLINE
116 static void updateParticleAttributes(
118 amrex::ParticleReal* radius_ptr,
119 amrex::ParticleReal* mass_ptr,
121 amrex::ParticleReal rho_w,
122 int num_sp,
int num_ae,
123 const int* sp_sol_arr,
124 const int* ae_sol_arr,
125 const SDPCDefn::SDSpeciesMassArr sp_mass_ptrs,
126 const SDPCDefn::SDAerosolMassArr ae_mass_ptrs,
127 const amrex::ParticleReal* sp_rho_arr,
128 const amrex::ParticleReal* ae_rho_arr)
131 radius_ptr[particle_idx] = SD_effective_radius(
132 particle_idx, idx_w, rho_w,
134 sp_sol_arr, ae_sol_arr,
135 sp_mass_ptrs, ae_mass_ptrs,
136 sp_rho_arr, ae_rho_arr);
139 mass_ptr[particle_idx] = SD_total_mass(
140 particle_idx, num_sp, num_ae,
141 sp_mass_ptrs, ae_mass_ptrs);
145 [[nodiscard]]
virtual amrex::Vector<std::string> varNames ()
const override;
148 [[nodiscard]]
virtual amrex::Vector<std::string> meshPlotVarNames ()
const override;
151 virtual void computeMeshVar(
const std::string&,
153 const amrex::MultiFab&,
154 const int )
const override;
157 virtual void InitializeParticles (
const amrex::Real,
const MFPtr& a_ptr)
override;
166 void setNumSDBoxDistribution( amrex::iMultiFab&,
169 const amrex::RealBox&,
173 void setNumSDBubbleDistribution( amrex::iMultiFab&,
176 const amrex::RealBox&,
182 virtual void SetAttributes ( amrex::MultiFab& a_mf );
187 virtual void DensityScaling (
const amrex::MultiFab& a_mf );
190 virtual void EvolveParticles (
int,
192 amrex::Vector<amrex::Vector<amrex::MultiFab>>&,
193 const amrex::Vector<MFPtr>& )
override
195 amrex::Abort(
"SuperDropletPC::EvolveParticles() is intentionally disabled.");
199 virtual void AdvectParticles (
int a_lev,
202 const amrex::MultiFab*
const a_flow_vel,
203 const amrex::MultiFab& a_density,
204 const amrex::MultiFab& a_pressure,
205 const amrex::MultiFab& a_temperature,
206 const amrex::Vector<MFPtr>& a_z_phys_nd,
207 const BCTypeArr& a_bctypes,
208 const bool a_recycle);
211 virtual void MassChange (
int a_lev,
213 const Species::Name& a_vap_name,
214 const amrex::MultiFab& a_temperature,
215 const amrex::MultiFab& a_pressure,
216 const amrex::MultiFab& a_sat_pressure,
217 const amrex::MultiFab& a_sat_ratio,
218 const amrex::Vector<MFPtr>& a_z_phys_nd,
219 const bool a_is_water);
222 virtual void Coalescence (
int a_lev,
224 const amrex::MultiFab& a_pressure,
225 const amrex::MultiFab& a_temperature);
228 virtual void Recycle (
const int a_lev,
229 const amrex::Vector<MFPtr>& a_z_phys_nd,
232 const bool a_recycle);
236 [[nodiscard]]
inline virtual amrex::Long NumSuperDroplets ()
238 return ERFPC::TotalNumberOfParticles();
242 [[nodiscard]]
virtual amrex::Real TotalNumberOfParticles ();
245 [[nodiscard]]
virtual amrex::Long NumSDDeactivated ();
250 virtual void ComputeDistributions (
int,
252 amrex::ParticleReal );
253 #ifdef ERF_USE_ML_UPHYS_DIAGNOSTICS
255 virtual void ComputeBinnedDistributions (
int );
257 virtual void ComputeBinnedDistributionsCell (
int,
amrex::Real );
261 virtual void SDNumberDensity ( amrex::MultiFab&,
const amrex::MultiFab& a_z_phys_nd,
const int a_comp=0 )
const;
263 virtual void numberDensity ( amrex::MultiFab&,
const amrex::MultiFab& a_z_phys_nd,
const int a_comp=0 )
const;
265 virtual void massDensity ( amrex::MultiFab&,
const amrex::MultiFab& a_z_phys_nd,
const int a_comp=0 )
const;
267 virtual void massFlux ( amrex::MultiFab&,
const amrex::MultiFab& a_z_phys_nd,
const int,
const int a_comp=0 )
const;
270 virtual void speciesMassDensity ( amrex::MultiFab&,
271 const amrex::MultiFab& a_z_phys_nd,
273 const int a_comp = 0)
const;
276 virtual void cloudRainDensity ( amrex::MultiFab&,
277 const amrex::MultiFab& a_z_phys_nd,
280 const int a_comp = 0)
const;
283 virtual void speciesMassFlux ( amrex::MultiFab&,
const amrex::MultiFab& a_z_phys_nd,
const int,
const int,
const int a_comp=0 )
const;
286 virtual void aerosolMassDensity ( amrex::MultiFab&,
const amrex::MultiFab& a_z_phys_nd,
const int,
const int a_comp=0 )
const;
288 virtual void aerosolMassFlux ( amrex::MultiFab&,
const amrex::MultiFab& a_z_phys_nd,
const int,
const int,
const int a_comp=0 )
const;
291 virtual void effectiveRadius ( amrex::MultiFab&,
const amrex::MultiFab& a_z_phys_nd,
const int a_comp=0 )
const;
294 virtual void applyBoundaryTreatment (
int,
const amrex::Vector<MFPtr>&,
const BCTypeArr&,
const bool );
300 SDCoalescenceKernelType m_coalescence_kernel;
301 bool m_include_brownian_coalescence;
303 SDTerminalVelocityType m_term_vel_type_w;
305 int m_num_sd_per_cell;
306 bool m_density_scaling;
307 bool m_nucleate_particles;
308 bool m_prescribed_advection;
311 std::vector<std::unique_ptr<MaterialProperties>> m_species_mat;
313 std::vector<std::unique_ptr<MaterialProperties>> m_aerosol_mat;
322 bool m_mass_change_logging;
323 FILE* m_mass_change_log;
324 std::string m_mass_change_log_fname;
327 SDMassChangeTIMethod m_mass_change_ti;
328 long m_num_unconverged_particles;
330 amrex::IntVect m_coalescence_bin_size;
332 int m_distribution_grid_size;
334 amrex::MultiFab m_mf_buf;
336 #ifdef ERF_USE_ML_UPHYS_DIAGNOSTICS
341 amrex::MultiFab m_mass_ln_R_mf;
343 amrex::MultiFab m_num_ln_R_mf;
350 int m_num_initializations = 1;
352 std::vector< std::unique_ptr<SDInitialization> > m_initializations;
355 int m_num_injections = 0;
357 std::vector< std::unique_ptr<SDInjection> > m_injections;
363 bool m_place_randomly_in_cells;
366 std::mt19937 m_rndeng;
375 bool m_save_inactive;
378 amrex::Gpu::DeviceVector<amrex::ParticleReal> m_sp_density;
379 amrex::Gpu::DeviceVector<int> m_sp_solubility;
380 amrex::Gpu::DeviceVector<amrex::ParticleReal> m_sp_ionization;
381 amrex::Gpu::DeviceVector<amrex::ParticleReal> m_sp_mol_weight;
383 amrex::Gpu::DeviceVector<amrex::ParticleReal> m_ae_density;
384 amrex::Gpu::DeviceVector<int> m_ae_solubility;
385 amrex::Gpu::DeviceVector<amrex::ParticleReal> m_ae_ionization;
386 amrex::Gpu::DeviceVector<amrex::ParticleReal> m_ae_mol_weight;
389 bool m_device_props_initialized =
false;
400 void initializeDeviceProperties();
403 SDProcess::ProcessContext buildProcessContext(
int a_lev)
const
405 SDProcess::ProcessContext ctx;
406 const amrex::Geometry& geom = m_gdb->Geom(a_lev);
407 ctx.plo = geom.ProbLoArray();
408 ctx.phi = geom.ProbHiArray();
409 ctx.dxi = geom.InvCellSizeArray();
410 ctx.dx = geom.CellSizeArray();
411 ctx.domain = geom.Domain();
412 for (
int d = 0; d < AMREX_SPACEDIM; d++) {
413 ctx.is_periodic[d] = geom.isPeriodic(d) ? 1 : 0;
415 const auto cell_size = geom.CellSize();
416 ctx.cell_volume = AMREX_D_TERM(cell_size[0], *cell_size[1], *cell_size[2]);
418 ctx.num_aerosols = m_num_aerosols;
419 ctx.idx_water = m_idx_w;
420 ctx.rho_water = m_species_mat[m_idx_w]->m_density;
425 template<
typename SOAType,
typename AOSType>
426 void setupParticlePointers(
429 SDProcess::ParticlePointers& ptrs)
const
431 using namespace SDPCDefn;
433 constexpr
int rtoff_i = SuperDropletsIntIdxSoA::ncomps;
434 constexpr
int rtoff_r = SuperDropletsRealIdxSoA::ncomps;
436 ptrs.num_particles = aos.numParticles();
437 ptrs.mass_ptr = soa.GetRealData(SuperDropletsRealIdxSoA::mass).data();
438 ptrs.radius_ptr = soa.GetRealData(rtoff_r + SuperDropletsRealIdxSoA_RT::radius).data();
439 ptrs.active_ptr = soa.GetIntData(rtoff_i + SuperDropletsIntIdxSoA_RT::active).data();
440 ptrs.v_ptr[0] = soa.GetRealData(SuperDropletsRealIdxSoA::vx).data();
441 ptrs.v_ptr[1] = soa.GetRealData(SuperDropletsRealIdxSoA::vy).data();
442 ptrs.v_ptr[2] = soa.GetRealData(SuperDropletsRealIdxSoA::vz).data();
443 ptrs.vterm_ptr = soa.GetRealData(rtoff_r + SuperDropletsRealIdxSoA_RT::term_vel).data();
444 ptrs.mult_ptr = soa.GetRealData(rtoff_r + SuperDropletsRealIdxSoA_RT::multiplicity).data();
445 setupMassPointers(soa, ptrs.sp_mass_ptrs, ptrs.ae_mass_ptrs);
446 if (!m_device_props_initialized) {
447 const_cast<SuperDropletPC*
>(
this)->initializeDeviceProperties();
449 ptrs.sp_rho_arr = m_sp_density.data();
450 ptrs.sp_sol_arr = m_sp_solubility.data();
451 ptrs.sp_ion_arr = m_sp_ionization.data();
452 ptrs.sp_mw_arr = m_sp_mol_weight.data();
453 ptrs.ae_rho_arr = m_ae_density.data();
454 ptrs.ae_sol_arr = m_ae_solubility.data();
455 ptrs.ae_ion_arr = m_ae_ionization.data();
456 ptrs.ae_mw_arr = m_ae_mol_weight.data();
460 template<
typename TileFunc>
461 void forEachParticleTileBody(
462 ParIterType& pti,
int a_lev,
463 const SDProcess::ProcessContext& ctx,
466 int grid = pti.index();
467 auto& ptile = ParticlesAt(a_lev, pti);
468 auto& aos = ptile.GetArrayOfStructs();
469 auto& soa = ptile.GetStructOfArrays();
470 if (aos.numParticles() == 0) {
return; }
471 auto* p_pbox = aos().data();
472 SDProcess::ParticlePointers ptrs;
473 setupParticlePointers(soa, aos, ptrs);
474 func(pti, grid, p_pbox, ptrs, ctx);
478 template<
typename TileFunc>
480 void forEachParticleTile(
482 const SDProcess::ProcessContext& ctx,
486 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
488 for (ParIterType pti(*
this, a_lev); pti.isValid(); ++pti) {
489 forEachParticleTileBody(pti, a_lev, ctx, std::forward<TileFunc>(func));
496 template<
typename TileFunc>
498 void forEachParticleTileSerial(
int a_lev,
const SDProcess::ProcessContext& ctx, TileFunc&& func)
500 for (ParIterType pti(*
this, a_lev); pti.isValid(); ++pti) {
501 forEachParticleTileBody(pti, a_lev, ctx, std::forward<TileFunc>(func));
506 template<
typename TileFunc>
507 void forEachParticleTile(
int a_lev, TileFunc&& func)
510 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
512 for (ParIterType pti(*
this, a_lev); pti.isValid(); ++pti) {
513 int grid = pti.index();
514 auto& ptile = ParticlesAt(a_lev, pti);
515 auto& aos = ptile.GetArrayOfStructs();
516 auto& soa = ptile.GetStructOfArrays();
517 const int num_particles = aos.numParticles();
518 if (num_particles == 0) {
continue; }
519 auto* p_pbox = aos().data();
520 SDProcess::ParticlePointers ptrs;
521 setupParticlePointers(soa, aos, ptrs);
522 func(pti, grid, p_pbox, ptrs, num_particles);
527 virtual void readInputs ()
override
529 amrex::Abort(
"SuperDropletPC::readInputs(): Do not use this interface.");
536 void initializeParticlesNull (
const MFPtr&) { }
541 void define (
const std::vector<Species::Name>&,
542 const std::vector<Species::Name>&,
543 const amrex::BoxArray&,
544 const amrex::DistributionMapping&,
548 void add_superdroplet_attributes();
551 using ERFPC::massDensity;
int m_num_species
Definition: ERF_InitCustomPert_MultiSpeciesBubble.H:28
amrex::Real Real
Definition: ERF_ShocInterface.H:19
Common data structures for SuperDroplet physical processes.
Super-droplets initial properties.
Definition: ERF_SDInitialization.H:150
Definition: ERF_MaterialProperties.H:160