ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_SuperDropletPC.H
Go to the documentation of this file.
1 #ifndef SUPERDROPLET_PC_H_
2 #define SUPERDROPLET_PC_H_
3 
4 #ifdef ERF_USE_PARTICLES
5 
6 #include <random>
7 #include "ERF_Constants.H"
11 #include <AMReX_StructOfArrays.H>
12 
13 class SuperDropletPC : public ERFPC
14 {
15  using MFPtr = std::unique_ptr<amrex::MultiFab>;
16  using BCTypeArr = amrex::GpuArray<ERF_BC, AMREX_SPACEDIM*2>;
17 
18  public:
19 
20  /*! \brief Constructor */
21  SuperDropletPC ( amrex::ParGDBBase* a_gdb,
22  const std::vector<Species::Name>& a_species_mat,
23  const std::vector<Species::Name>& a_aerosol_mat,
24  const amrex::Real a_dt,
25  const std::string& a_name = "super_droplets" )
26  : ERFPC (a_gdb, a_name)
27  {
28  define( a_species_mat,
29  a_aerosol_mat,
30  a_gdb->ParticleBoxArray(m_lev),
31  a_gdb->ParticleDistributionMap(m_lev),
32  a_dt );
33  }
34 
35  /*! \brief Constructor */
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,
41  const amrex::Real a_dt,
42  const std::string& a_name = "super_droplets" )
43  : ERFPC (a_geom, a_dmap, a_ba, a_name)
44  {
45  define(a_species_mat, a_aerosol_mat, a_ba, a_dmap, a_dt);
46  }
47 
48  /*! \brief Destructor */
49  ~SuperDropletPC()
50  {
51  if (m_mass_change_logging) {
52  fclose(m_mass_change_log);
53  }
54  }
55 
56  /*! \brief Set vapour species material */
57  AMREX_FORCE_INLINE
58  virtual void setSpeciesMaterial ( const Species::Name& a_name )
59  {
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;
63  }
64 
65  /*! \brief Set species material */
66  AMREX_FORCE_INLINE
67  virtual void setSpeciesMaterial ( const std::vector<Species::Name>& a_names )
68  {
69  for (auto& name : a_names) { setSpeciesMaterial(name); }
70  }
71 
72  /*! \brief Get vapour material */
73  AMREX_FORCE_INLINE
74  virtual const MaterialProperties& getSpeciesMaterial(const Species::Name& a_name) const
75  {
76  for (auto& species : m_species_mat) {
77  if (species->m_name == a_name) { return *species; }
78  }
79  amrex::Abort("SuperDropletPC::getSpeciesMaterial() - species not found");
80  return *m_species_mat[0];
81  }
82 
83  /*! \brief Set aerosol material */
84  AMREX_FORCE_INLINE
85  virtual void setAerosolMaterial ( const Species::Name& a_name )
86  {
87  m_aerosol_mat.push_back(std::make_unique<MaterialProperties>(a_name));
88  m_device_props_initialized = false;
89  }
90 
91  /*! \brief Set aerosol material */
92  AMREX_FORCE_INLINE
93  virtual void setAerosolMaterial ( const std::vector<Species::Name>& a_names )
94  {
95  for (auto& name : a_names) { setAerosolMaterial(name); }
96  }
97 
98  /*! \brief Update device properties if material properties change */
99  void updateDeviceProperties();
100 
101  /*! \brief Setup species and aerosol mass pointers from SOA */
102  template<typename SOAType>
103  void setupMassPointers(SOAType& soa, SDPCDefn::SDSpeciesMassArr& sp_mass_ptrs,
104  SDPCDefn::SDAerosolMassArr& ae_mass_ptrs) const
105  {
106  for (int i = 0; i < m_num_species; i++) {
107  sp_mass_ptrs[i] = soa.GetRealData(idx_s(i, m_num_aerosols, m_num_species)).data();
108  }
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();
111  }
112  }
113 
114  /*! \brief Update particle radius and total mass from species/aerosol masses */
115  AMREX_GPU_DEVICE AMREX_FORCE_INLINE
116  static void updateParticleAttributes(
117  int particle_idx,
118  amrex::ParticleReal* radius_ptr,
119  amrex::ParticleReal* mass_ptr,
120  int idx_w,
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)
129  {
130  // Update effective radius
131  radius_ptr[particle_idx] = SD_effective_radius(
132  particle_idx, idx_w, rho_w,
133  num_sp, num_ae,
134  sp_sol_arr, ae_sol_arr,
135  sp_mass_ptrs, ae_mass_ptrs,
136  sp_rho_arr, ae_rho_arr);
137 
138  // Update total mass
139  mass_ptr[particle_idx] = SD_total_mass(
140  particle_idx, num_sp, num_ae,
141  sp_mass_ptrs, ae_mass_ptrs);
142  }
143 
144  /*! \brief Get real-type particle attribute names */
145  [[nodiscard]] virtual amrex::Vector<std::string> varNames () const override;
146 
147  /*! \brief Get Eulerian plot variable names */
148  [[nodiscard]] virtual amrex::Vector<std::string> meshPlotVarNames () const override;
149 
150  /*! \brief Compute mesh variable from particles */
151  virtual void computeMeshVar( const std::string&,
152  amrex::MultiFab&,
153  const amrex::MultiFab&,
154  const int ) const override;
155 
156  /*! \brief Initialize super-droplets in domain */
157  virtual void InitializeParticles (const amrex::Real, const MFPtr& a_ptr) override;
158 
159  /*! \brief Inject super-droplets in domain */
160  virtual void InjectParticles (const amrex::Real, const MFPtr& a_ptr, const amrex::Real);
161 
162  /*! \brief add particles for a given set of initialization parameters */
163  void addParticles ( const MFPtr&, const SDInitProperties& );
164 
165  /*! \brief Set initial number of superdroplets per cell as a box with uniform density */
166  void setNumSDBoxDistribution( amrex::iMultiFab&,
167  const int,
168  const MFPtr&,
169  const amrex::RealBox&,
170  const bool );
171 
172  /*! \brief Set initial number of superdroplets per cell as a bubble with uniform density */
173  void setNumSDBubbleDistribution( amrex::iMultiFab&,
174  const int,
175  const MFPtr&,
176  const amrex::RealBox&,
177  const bool );
178 
179  /*! \brief Set super-droplets attributes from a given condensate mass density
180  * \param[in] a_mf Multifab containing the initial condensate mass density
181  */
182  virtual void SetAttributes ( amrex::MultiFab& a_mf );
183 
184  /*! \brief Scale the SDM number density with air density
185  * \param[in] a_mf Multifab containing the air density field
186  */
187  virtual void DensityScaling ( const amrex::MultiFab& a_mf );
188 
189  /*! \brief Evolve particles for one time step */
190  virtual void EvolveParticles ( int,
191  amrex::Real,
192  amrex::Vector<amrex::Vector<amrex::MultiFab>>&,
193  const amrex::Vector<MFPtr>& ) override
194  {
195  amrex::Abort("SuperDropletPC::EvolveParticles() is intentionally disabled.");
196  }
197 
198  /*! \brief Advect particles for one time step */
199  virtual void AdvectParticles ( int a_lev,
200  amrex::Real a_time,
201  amrex::Real a_dt,
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);
209 
210  /*! \brief Grow/shrink particles for one time step */
211  virtual void MassChange ( int a_lev,
212  amrex::Real a_dt,
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);
220 
221  /*! \brief Coalescence of super-droplets */
222  virtual void Coalescence ( int a_lev,
223  amrex::Real a_dt,
224  const amrex::MultiFab& a_pressure,
225  const amrex::MultiFab& a_temperature);
226 
227  /*! \brief Recycle super-droplets */
228  virtual void Recycle ( const int a_lev,
229  const amrex::Vector<MFPtr>& a_z_phys_nd,
230  const int a_iter,
231  const amrex::Real a_dt,
232  const bool a_recycle);
233 
234  /*! \brief Return the number of super-droplets
235  This returns the number of AMReX particles being used in the simulation. */
236  [[nodiscard]] inline virtual amrex::Long NumSuperDroplets ()
237  {
238  return ERFPC::TotalNumberOfParticles();
239  }
240 
241  /*! \brief Return the total number of physical particles */
242  [[nodiscard]] virtual amrex::Real TotalNumberOfParticles ();
243 
244  /*! \brief Return the number of deactivated super-droplets */
245  [[nodiscard]] virtual amrex::Long NumSDDeactivated ();
246 
247  /*! \brief Compute and print diagnostics */
248  virtual void Diagnostics (int, amrex::Real, bool);
249  /*! \brief Compute mass density distribution */
250  virtual void ComputeDistributions ( int,
251  amrex::ParticleReal,
252  amrex::ParticleReal );
253 #ifdef ERF_USE_ML_UPHYS_DIAGNOSTICS
254  /*! \brief Compute binned size distributions */
255  virtual void ComputeBinnedDistributions ( int );
256  /*! \brief Compute cell-wise binned size distributions */
257  virtual void ComputeBinnedDistributionsCell ( int, amrex::Real );
258 #endif
259 
260  /*! \brief Compute the SD number density on a mesh */
261  virtual void SDNumberDensity ( amrex::MultiFab&, const amrex::MultiFab& a_z_phys_nd, const int a_comp=0 ) const;
262  /*! \brief Compute the particle number density on a mesh */
263  virtual void numberDensity ( amrex::MultiFab&, const amrex::MultiFab& a_z_phys_nd, const int a_comp=0 ) const;
264  /*! \brief Compute the particle mass density on a mesh */
265  virtual void massDensity ( amrex::MultiFab&, const amrex::MultiFab& a_z_phys_nd, const int a_comp=0 ) const;
266  /*! \brief Compute the mass flux component on a mesh */
267  virtual void massFlux ( amrex::MultiFab&, const amrex::MultiFab& a_z_phys_nd, const int, const int a_comp=0 ) const;
268 
269  /*! \brief Compute the species mass density on a mesh */
270  virtual void speciesMassDensity ( amrex::MultiFab&,
271  const amrex::MultiFab& a_z_phys_nd,
272  int,
273  const int a_comp = 0) const;
274 
275  /*! \brief Compute the cloud/rain mass density on a mesh */
276  virtual void cloudRainDensity ( amrex::MultiFab&,
277  const amrex::MultiFab& a_z_phys_nd,
278  amrex::Real,
279  amrex::Real,
280  const int a_comp = 0) const;
281 
282  /*! \brief Compute the species mass flux component on a mesh */
283  virtual void speciesMassFlux ( amrex::MultiFab&, const amrex::MultiFab& a_z_phys_nd, const int, const int, const int a_comp=0 ) const;
284 
285  /*! \brief Compute the aerosol mass density on a mesh */
286  virtual void aerosolMassDensity ( amrex::MultiFab&, const amrex::MultiFab& a_z_phys_nd, const int, const int a_comp=0 ) const;
287  /*! \brief Compute the aerosol mass flux component on a mesh */
288  virtual void aerosolMassFlux ( amrex::MultiFab&, const amrex::MultiFab& a_z_phys_nd, const int, const int, const int a_comp=0 ) const;
289 
290  /*! \brief Compute the particle effective radius on a mesh */
291  virtual void effectiveRadius ( amrex::MultiFab&, const amrex::MultiFab& a_z_phys_nd, const int a_comp=0 ) const;
292 
293  /*! \brief Applies boundary treatment to the particle container */
294  virtual void applyBoundaryTreatment ( int, const amrex::Vector<MFPtr>&, const BCTypeArr&, const bool );
295 
296  protected:
297 
298  const int m_lev = 0; /*!< AMR Level: always 0 for super-droplets PC */
299 
300  SDCoalescenceKernelType m_coalescence_kernel; /*!< Coalescence kernel */
301  bool m_include_brownian_coalescence; /*!< Include Brownian coalescence? */
302 
303  SDTerminalVelocityType m_term_vel_type_w; /*!< Terminal velocity model for water droplets */
304 
305  int m_num_sd_per_cell; /*!< Number of super-droplets per cell */
306  bool m_density_scaling; /*!< Scale initial number density with air density */
307  bool m_nucleate_particles; /*!< nucleate new superdroplets from vapour */
308  bool m_prescribed_advection; /*!< move particles with prescribed vertical velocity */
309 
310  int m_num_species; /*!< Number of vapour/condensate species */
311  std::vector<std::unique_ptr<MaterialProperties>> m_species_mat; /*!< Vapour/condensate material */
312  int m_num_aerosols; /*!< Number of aerosols */
313  std::vector<std::unique_ptr<MaterialProperties>> m_aerosol_mat; /*!< Aerosol materials */
314 
315  /* Mass change equation Newton solver parameters */
316  amrex::Real m_newton_rtol; /*!< Newton solver - relative tolerance */
317  amrex::Real m_newton_atol; /*!< Newton solver - absolute tolerance */
318  amrex::Real m_newton_stol; /*!< Newton solver - step tolerance */
319  int m_newton_maxits; /*!< Newton solver - maximum iterations */
320 
321  /* Mass change equation (un)convergence logs */
322  bool m_mass_change_logging; /*!< Whether to log unconverged particles */
323  FILE* m_mass_change_log; /*!< File handle for log with unconverged info */
324  std::string m_mass_change_log_fname; /*!< Name of unconverged info log file */
325 
326  amrex::Real m_mass_change_cfl; /*!< CFL for phase change equation */
327  SDMassChangeTIMethod m_mass_change_ti; /*!< time integrator for mass change ODE */
328  long m_num_unconverged_particles; /*!< total number of unconverged particles */
329 
330  amrex::IntVect m_coalescence_bin_size; /*!< bin size for coalescence */
331 
332  int m_distribution_grid_size; /*!< Size of 1D grid to compute distributions on */
333 
334  amrex::MultiFab m_mf_buf;
335 
336 #ifdef ERF_USE_ML_UPHYS_DIAGNOSTICS
337  amrex::Real m_bindist_rmin; /*!< min radius for binned distribution */
338  amrex::Real m_bindist_rmax; /*!< max radius for binned distribution */
339 
340  /*! Cell-wise mass distribution */
341  amrex::MultiFab m_mass_ln_R_mf;
342  /*! Cell-wise number distribution */
343  amrex::MultiFab m_num_ln_R_mf;
344 #endif
345 
346  /*! maximum coalescence time scale */
347  amrex::Real m_t_coalescence;
348 
349  /*! Initialization parameters objects */
350  int m_num_initializations = 1;
351  /*! vector of initializations */
352  std::vector< std::unique_ptr<SDInitialization> > m_initializations;
353 
354  /*! Injection parameters objects */
355  int m_num_injections = 0;
356  /*! vector of injections */
357  std::vector< std::unique_ptr<SDInjection> > m_injections;
358 
359  /*! sigma_0 parameter for mass distribution calculation */
360  amrex::Real m_sigma0;
361 
362  /*! place the initial SDs randomly within cells */
363  bool m_place_randomly_in_cells;
364 
365  /*! random engine */
366  std::mt19937 m_rndeng;
367 
368  /*! species index of water */
369  int m_idx_w = -1;
370 
371  /*! inactive particles threshold */
372  amrex::Real m_deac_threshold;
373 
374  /*! save inactive particles to file */
375  bool m_save_inactive;
376 
377  /*! Persistent device vectors for material properties */
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;
382 
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;
387 
388  /*! Flag to track if device properties are initialized */
389  bool m_device_props_initialized = false;
390 
391  /* recycled particle position bounds */
392  amrex::Real m_recyc_xmin;
393  amrex::Real m_recyc_xmax;
394  amrex::Real m_recyc_ymin;
395  amrex::Real m_recyc_ymax;
396  amrex::Real m_recyc_zmin;
397  amrex::Real m_recyc_zmax;
398 
399  /*! Method to initialize device properties */
400  void initializeDeviceProperties();
401 
402  /*! \brief Build ProcessContext with geometry and species info */
403  SDProcess::ProcessContext buildProcessContext(int a_lev) const
404  {
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;
414  }
415  const auto cell_size = geom.CellSize();
416  ctx.cell_volume = AMREX_D_TERM(cell_size[0], *cell_size[1], *cell_size[2]);
417  ctx.num_species = m_num_species;
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;
421  return ctx;
422  }
423 
424  /*! \brief Setup particle attribute pointers for a tile */
425  template<typename SOAType, typename AOSType>
426  void setupParticlePointers(
427  SOAType& soa,
428  AOSType& aos,
429  SDProcess::ParticlePointers& ptrs) const
430  {
431  using namespace SDPCDefn;
432 
433  constexpr int rtoff_i = SuperDropletsIntIdxSoA::ncomps;
434  constexpr int rtoff_r = SuperDropletsRealIdxSoA::ncomps;
435 
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();
448  }
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();
457  }
458 
459  /*! \brief Common body for forEachParticleTile iterations */
460  template<typename TileFunc>
461  void forEachParticleTileBody(
462  ParIterType& pti, int a_lev,
463  const SDProcess::ProcessContext& ctx,
464  TileFunc&& func)
465  {
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);
475  }
476 
477  /*! \brief Iterate over particle tiles with common setup (parallel) */
478  template<typename TileFunc>
479  AMREX_FORCE_INLINE
480  void forEachParticleTile(
481  int a_lev,
482  const SDProcess::ProcessContext& ctx,
483  TileFunc&& func)
484  {
485 #ifdef AMREX_USE_OMP
486 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
487 #endif
488  for (ParIterType pti(*this, a_lev); pti.isValid(); ++pti) {
489  forEachParticleTileBody(pti, a_lev, ctx, std::forward<TileFunc>(func));
490  }
491  }
492 
493  /*! \brief Iterate over particle tiles WITHOUT OpenMP (serial)
494  * Use when per-tile operations are not thread-safe (e.g., DenseBins).
495  */
496  template<typename TileFunc>
497  AMREX_FORCE_INLINE
498  void forEachParticleTileSerial(int a_lev, const SDProcess::ProcessContext& ctx, TileFunc&& func)
499  {
500  for (ParIterType pti(*this, a_lev); pti.isValid(); ++pti) {
501  forEachParticleTileBody(pti, a_lev, ctx, std::forward<TileFunc>(func));
502  }
503  }
504 
505  /*! \brief Lightweight iteration (no ProcessContext) */
506  template<typename TileFunc>
507  void forEachParticleTile(int a_lev, TileFunc&& func)
508  {
509 #ifdef AMREX_USE_OMP
510 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
511 #endif
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);
523  }
524  }
525 
526  /*! \brief read inputs from file */
527  virtual void readInputs () override
528  {
529  amrex::Abort("SuperDropletPC::readInputs(): Do not use this interface.");
530  }
531 
532  /*! \brief read inputs from file */
533  virtual void readInputs (const amrex::Real);
534 
535  /*! \brief Particle initialization - null (no particles are initialized) */
536  void initializeParticlesNull (const MFPtr&) { }
537 
538  private:
539 
540  /*! \brief define super-droplets */
541  void define ( const std::vector<Species::Name>&,
542  const std::vector<Species::Name>&,
543  const amrex::BoxArray&,
544  const amrex::DistributionMapping&,
545  const amrex::Real );
546 
547  /*! \brief add super-droplet method-specific particle attributes */
548  void add_superdroplet_attributes();
549 
550  /* do not use */
551  using ERFPC::massDensity;
552 
553 };
554 
555 #endif
556 #endif
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