ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_InitForEnsemble.cpp File Reference
#include <ERF.H>
#include <ERF_TileNoZ.H>
#include <AMReX_PlotFileUtil.H>
#include <filesystem>
Include dependency graph for ERF_InitForEnsemble.cpp:

Functions

void NormalizeMultiFabRMS_PerComponent (MultiFab &mf_cc_pert)
 
void ApplyNeumannBCs (const Geometry &geom, MultiFab &mf_cc)
 
void ReadCustomDataFile (const std::string &filename_custom, int &nx, int &ny, int &nz, int &ng, int &ncomp, std::array< Real, 3 > &problo_ext, std::array< Real, 3 > &probhi_ext, Vector< Real > &data_rho, Vector< Real > &data_theta, Vector< Real > &data_xvel, Vector< Real > &data_yvel, Vector< Real > &data_zvel)
 
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE int idx (int i, int j, int k, int nx, int ny)
 
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Real interp_trilinear (const Real *f, int i, int j, int k, Real tx, Real ty, Real tz, int nx, int ny, int nz)
 
void InterpolateToFineMF (const Vector< Real > &data_rho, const Vector< Real > &data_theta, const Vector< Real > &data_xvel, const Vector< Real > &data_yvel, const Vector< Real > &data_zvel, int nx, int ny, int nz, const std::array< Real, 3 > &problo, const std::array< Real, 3 > &probhi, MultiFab &mf_fine, const Geometry &geom_fine)
 
void MakeFinalMultiFabs (const MultiFab &mf_cc_fine, MultiFab &cons_pert, MultiFab &xvel_pert, MultiFab &yvel_pert, MultiFab &zvel_pert)
 
void AddPertToBckgnd (MultiFab &mf_cc_fine, const MultiFab &mf_cc_pert)
 
std::string get_last_plotfile (const std::string &plotfile_dir)
 
void read_plot_file (PlotFileData &pf, const std::vector< std::string > varnames, MultiFab &mf)
 

Function Documentation

◆ AddPertToBckgnd()

void AddPertToBckgnd ( MultiFab &  mf_cc_fine,
const MultiFab &  mf_cc_pert 
)
492 {
493  const int ncomp = mf_cc_fine.nComp();
494 
495  // Optional safety check (recommended)
496  AMREX_ALWAYS_ASSERT(mf_cc_pert.nComp() == ncomp);
497 
498  for (MFIter mfi(mf_cc_fine, TilingIfNotGPU()); mfi.isValid(); ++mfi)
499  {
500  const Box& bx = mfi.tilebox();
501 
502  auto const& bg = mf_cc_fine.array(mfi);
503  auto const& pert = mf_cc_pert.const_array(mfi);
504 
505  amrex::ParallelFor(bx, ncomp,
506  [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
507  {
508  Real ens_amp = 0.02*std::abs(bg(i,j,k,n));
509  bg(i,j,k,n) += ens_amp*pert(i,j,k,n);
510  });
511  }
512 }
AMREX_ALWAYS_ASSERT(bx.length()[2]==khi+1)
ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept { const Real *dx=geomdata.CellSize();const Real x=(i+0.5) *dx[0];const Real y=(j+0.5) *dx[1];const Real Omg=erf_vortex_Gaussian(x, y, xc, yc, R, beta, sigma);const Real deltaT=-(gamma - 1.0)/(2.0 *sigma *sigma) *Omg *Omg;const Real rho_norm=std::pow(1.0+deltaT, inv_gm1);const Real T=(1.0+deltaT) *T_inf;const Real p=std::pow(rho_norm, Gamma)/Gamma *rho_0 *a_inf *a_inf;const Real rho_theta=rho_0 *rho_norm *(T *std::pow(p_0/p, rdOcp));state_pert(i, j, k, RhoTheta_comp)=rho_theta - getRhoThetagivenP(p_hse(i, j, k));const Real r2d_xy=std::sqrt((x-xc) *(x-xc)+(y-yc) *(y-yc));state_pert(i, j, k, RhoScalar_comp)=0.25 *(1.0+std::cos(PI *std::min(r2d_xy, R)/R));})
amrex::Real Real
Definition: ERF_ShocInterface.H:19
real(c_double), private bg
Definition: ERF_module_mp_morr_two_moment.F90:182

Referenced by ERF::create_background_state_for_ensemble().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ ApplyNeumannBCs()

void ApplyNeumannBCs ( const Geometry &  geom,
MultiFab &  mf_cc 
)
164 {
165 
166  // -------------------------------------------------
167  // 2. Fill interior + periodic ghost cells
168  // -------------------------------------------------
169  mf_cc.FillBoundary(geom.periodicity());
170  // -------------------------------------------------
171  // 3. Apply FOExtrap (Neumann) at domain boundaries
172  // -------------------------------------------------
173  const Box& domain = geom.Domain();
174 
175  for (MFIter mfi(mf_cc, TilingIfNotGPU()); mfi.isValid(); ++mfi)
176  {
177  const Box& gbx = mfi.growntilebox(); // includes ghost cells
178  const Box& vbx = mfi.validbox();
179 
180  auto const& arr = mf_cc.array(mfi);
181  int ncomp = mf_cc.nComp();
182 
183  ParallelFor(gbx, ncomp,
184  [=] AMREX_GPU_DEVICE (int i, int j, int k, int n)
185  {
186  if (vbx.contains(i,j,k)) return;
187 
188  int ii = i;
189  int jj = j;
190  int kk = k;
191 
192  // Clamp to domain interior (FOExtrap)
193  ii = amrex::max(domain.smallEnd(0),
194  amrex::min(i, domain.bigEnd(0)));
195 
196  jj = amrex::max(domain.smallEnd(1),
197  amrex::min(j, domain.bigEnd(1)));
198 
199  kk = amrex::max(domain.smallEnd(2),
200  amrex::min(k, domain.bigEnd(2)));
201 
202  arr(i,j,k,n) = arr(ii,jj,kk,n);
203  });
204  }
205 }

Referenced by ERF::create_background_state_for_ensemble().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_last_plotfile()

std::string get_last_plotfile ( const std::string &  plotfile_dir)
520 {
521  std::vector<std::string> pltfiles;
522 
523  for (const auto& entry : fs::directory_iterator(plotfile_dir)) {
524  if (entry.is_directory()) {
525  std::string name = entry.path().filename().string();
526  if (name.find("plt") == 0) {
527  pltfiles.push_back(entry.path().string());
528  }
529  }
530  }
531 
532  std::sort(pltfiles.begin(), pltfiles.end());
533  return pltfiles.back(); // last one
534 }

◆ idx()

◆ interp_trilinear()

AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Real interp_trilinear ( const Real f,
int  i,
int  j,
int  k,
Real  tx,
Real  ty,
Real  tz,
int  nx,
int  ny,
int  nz 
)
299 {
300  int i1 = amrex::min(i+1, nx-1);
301  int j1 = amrex::min(j+1, ny-1);
302  int k1 = amrex::min(k+1, nz-1);
303 
304  Real c000 = f[idx(i ,j ,k ,nx,ny)];
305  Real c100 = f[idx(i1,j ,k ,nx,ny)];
306  Real c010 = f[idx(i ,j1,k ,nx,ny)];
307  Real c110 = f[idx(i1,j1,k ,nx,ny)];
308  Real c001 = f[idx(i ,j ,k1,nx,ny)];
309  Real c101 = f[idx(i1,j ,k1,nx,ny)];
310  Real c011 = f[idx(i ,j1,k1,nx,ny)];
311  Real c111 = f[idx(i1,j1,k1,nx,ny)];
312 
313  Real c00 = c000*(1-tx) + c100*tx;
314  Real c10 = c010*(1-tx) + c110*tx;
315  Real c01 = c001*(1-tx) + c101*tx;
316  Real c11 = c011*(1-tx) + c111*tx;
317 
318  Real c0 = c00*(1-ty) + c10*ty;
319  Real c1 = c01*(1-ty) + c11*ty;
320 
321  return c0*(1-tz) + c1*tz;
322 }
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE int idx(int i, int j, int k, int nx, int ny)
Definition: ERF_InitForEnsemble.cpp:287
real(c_double), private k1
Definition: ERF_module_mp_morr_two_moment.F90:213
real(c_double), private c1
Definition: ERF_module_mp_morr_two_moment.F90:212

Referenced by InterpolateToFineMF().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ InterpolateToFineMF()

void InterpolateToFineMF ( const Vector< Real > &  data_rho,
const Vector< Real > &  data_theta,
const Vector< Real > &  data_xvel,
const Vector< Real > &  data_yvel,
const Vector< Real > &  data_zvel,
int  nx,
int  ny,
int  nz,
const std::array< Real, 3 > &  problo,
const std::array< Real, 3 > &  probhi,
MultiFab &  mf_fine,
const Geometry &  geom_fine 
)
336 {
337  // coarse spacing
338  Real dx_c[3];
339  dx_c[0] = (probhi[0] - problo[0]) / nx;
340  dx_c[1] = (probhi[1] - problo[1]) / ny;
341  dx_c[2] = (probhi[2] - problo[2]) / nz;
342 
343  const auto problo_f = geom_fine.ProbLoArray();
344  const auto dx_f = geom_fine.CellSizeArray();
345 
346  // Step 1: declare device vectors with correct size
347  amrex::Gpu::DeviceVector<Real> d_rho(data_rho.size());
348  amrex::Gpu::DeviceVector<Real> d_theta(data_theta.size());
349  amrex::Gpu::DeviceVector<Real> d_xvel(data_xvel.size());
350  amrex::Gpu::DeviceVector<Real> d_yvel(data_yvel.size());
351  amrex::Gpu::DeviceVector<Real> d_zvel(data_zvel.size());
352 
353  // Step 2: copy data from host to device
354  amrex::Gpu::copyAsync(amrex::Gpu::hostToDevice,
355  data_rho.begin(), data_rho.end(),
356  d_rho.begin());
357 
358  amrex::Gpu::copyAsync(amrex::Gpu::hostToDevice,
359  data_theta.begin(), data_theta.end(),
360  d_theta.begin());
361 
362  amrex::Gpu::copyAsync(amrex::Gpu::hostToDevice,
363  data_xvel.begin(), data_xvel.end(),
364  d_xvel.begin());
365 
366  amrex::Gpu::copyAsync(amrex::Gpu::hostToDevice,
367  data_yvel.begin(), data_yvel.end(),
368  d_yvel.begin());
369 
370  amrex::Gpu::copyAsync(amrex::Gpu::hostToDevice,
371  data_zvel.begin(), data_zvel.end(),
372  d_zvel.begin());
373 
374  const Real* rho_ptr = d_rho.data();
375  const Real* theta_ptr = d_theta.data();
376  const Real* xvel_ptr = d_xvel.data();
377  const Real* yvel_ptr = d_yvel.data();
378  const Real* zvel_ptr = d_zvel.data();
379  // -------------------------------
380  // GPU kernel over MultiFab
381  // -------------------------------
382  for (MFIter mfi(mf_fine); mfi.isValid(); ++mfi)
383  {
384  const Box& bx = mfi.validbox();
385  auto arr = mf_fine.array(mfi);
386 
388  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
389  {
390  // physical location (fine cell center)
391  Real x = problo_f[0] + (i + 0.5) * dx_f[0];
392  Real y = problo_f[1] + (j + 0.5) * dx_f[1];
393  Real z = problo_f[2] + (k + 0.5) * dx_f[2];
394 
395  // map to coarse index space
396  Real rx = (x - problo[0]) / dx_c[0] - 0.5;
397  Real ry = (y - problo[1]) / dx_c[1] - 0.5;
398  Real rz = (z - problo[2]) / dx_c[2] - 0.5;
399 
400  int ic = static_cast<int>(floor(rx));
401  int jc = static_cast<int>(floor(ry));
402  int kc = static_cast<int>(floor(rz));
403 
404  Real tx = rx - ic;
405  Real ty = ry - jc;
406  Real tz = rz - kc;
407 
408  // clamp
409  ic = amrex::max(0, amrex::min(ic, nx-1));
410  jc = amrex::max(0, amrex::min(jc, ny-1));
411  kc = amrex::max(0, amrex::min(kc, nz-1));
412 
413  //printf("The values are x, y, z, rx, ry, rz = %0.15g %0.15g %0.15g %0.15g %0.15g %0.15g\n", x, y, z, rx, ry, rz);
414 
415  // interpolate each component using device trilinear
416  arr(i,j,k,0) = interp_trilinear(rho_ptr, ic,jc,kc, tx,ty,tz, nx,ny,nz);
417  arr(i,j,k,1) = interp_trilinear(theta_ptr, ic,jc,kc, tx,ty,tz, nx,ny,nz);
418  arr(i,j,k,2) = interp_trilinear(xvel_ptr, ic,jc,kc, tx,ty,tz, nx,ny,nz);
419  arr(i,j,k,3) = interp_trilinear(yvel_ptr, ic,jc,kc, tx,ty,tz, nx,ny,nz);
420  arr(i,j,k,4) = interp_trilinear(zvel_ptr, ic,jc,kc, tx,ty,tz, nx,ny,nz);
421  });
422  }
423 }
auto probhi
Definition: ERF_InitCustomPertVels_ABL.H:21
auto problo
Definition: ERF_InitCustomPertVels_ABL.H:20
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE Real interp_trilinear(const Real *f, int i, int j, int k, Real tx, Real ty, Real tz, int nx, int ny, int nz)
Definition: ERF_InitForEnsemble.cpp:294

Referenced by ERF::create_background_state_for_ensemble().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ MakeFinalMultiFabs()

void MakeFinalMultiFabs ( const MultiFab &  mf_cc_fine,
MultiFab &  cons_pert,
MultiFab &  xvel_pert,
MultiFab &  yvel_pert,
MultiFab &  zvel_pert 
)
431 {
432 
433  for (MFIter mfi(cons_pert, TilingIfNotGPU()); mfi.isValid(); ++mfi)
434  {
435  const Box& bx = mfi.tilebox();
436 
437  auto const& mf_cc_fine_arr = mf_cc_fine.const_array(mfi);
438  auto const& cons_pert_arr = cons_pert.array(mfi);
439 
440  amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
441  {
442  Real tmp_rho = mf_cc_fine_arr(i,j,k,0);
443  Real tmp_theta = mf_cc_fine_arr(i,j,k,1);
444  cons_pert_arr(i,j,k,Rho_comp) = tmp_rho;
445  cons_pert_arr(i,j,k,RhoTheta_comp) = tmp_rho*tmp_theta;
446  });
447  }
448 
449  // --- X-faces (component 2) ---
450  for (MFIter mfi(xvel_pert, TilingIfNotGPU()); mfi.isValid(); ++mfi)
451  {
452  const Box& bx = mfi.tilebox();
453  auto const& uface = xvel_pert.array(mfi);
454  auto const& cc = mf_cc_fine.const_array(mfi);
455 
456  amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
457  {
458  uface(i,j,k) = 0.5 * (cc(i-1,j,k,2) + cc(i,j,k,2));
459  });
460  }
461 
462  // --- Y-faces (component 3) ---
463  for (MFIter mfi(yvel_pert, TilingIfNotGPU()); mfi.isValid(); ++mfi)
464  {
465  const Box& bx = mfi.tilebox();
466  auto const& vface = yvel_pert.array(mfi);
467  auto const& cc = mf_cc_fine.const_array(mfi);
468 
469  amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
470  {
471  vface(i,j,k) = 0.5 * (cc(i,j-1,k,3) + cc(i,j,k,3));
472  });
473  }
474 
475  // --- Z-faces (component 4) ---
476  for (MFIter mfi(zvel_pert, TilingIfNotGPU()); mfi.isValid(); ++mfi)
477  {
478  const Box& bx = mfi.tilebox();
479  auto const& wface = zvel_pert.array(mfi);
480  //auto const& cc = mf_cc_fine.const_array(mfi);
481 
482  amrex::ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
483  {
484  wface(i,j,k) = 0.0;//0.5 * (cc(i,j,k-1,4) + cc(i,j,k,4));
485  });
486  }
487 }
#define Rho_comp
Definition: ERF_IndexDefines.H:36
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37

Referenced by ERF::create_background_state_for_ensemble().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ NormalizeMultiFabRMS_PerComponent()

void NormalizeMultiFabRMS_PerComponent ( MultiFab &  mf_cc_pert)
41 {
42  const int ncomp = mf_cc_pert.nComp();
43 
44  for (int n = 0; n < ncomp; ++n)
45  {
46  // 1. Set up AMReX reduction (sum of squares + count)
47  ReduceOps<ReduceOpSum, ReduceOpSum> reduce_op;
48  ReduceData<Real, Long> reduce_data(reduce_op);
49  using ReduceTuple = typename decltype(reduce_data)::Type;
50 
51  // 2. Loop over tiles and accumulate
52  for (MFIter mfi(mf_cc_pert, TilingIfNotGPU()); mfi.isValid(); ++mfi)
53  {
54  const Box& bx = mfi.tilebox();
55  auto const& arr = mf_cc_pert.const_array(mfi);
56 
57  reduce_op.eval(bx, reduce_data,
58  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept -> ReduceTuple
59  {
60  Real v = arr(i, j, k, n);
61  return { v * v, 1L };
62  });
63  }
64 
65  // 3. Retrieve results (includes implicit GPU sync + host copy)
66  auto rv = reduce_data.value(reduce_op);
67  Real h_sumsq = amrex::get<0>(rv);
68  Long h_count = amrex::get<1>(rv);
69 
70  // 4. Sum across MPI ranks
71  ParallelDescriptor::ReduceRealSum(h_sumsq);
72  ParallelDescriptor::ReduceLongSum(h_count);
73 
74  // 5. Compute RMS and normalize
75  if (h_count > 0)
76  {
77  Real rms = std::sqrt(h_sumsq / static_cast<Real>(h_count));
78  if (rms > 0.0) {
79  mf_cc_pert.mult(1.0 / rms, n, 1);
80  }
81  }
82  }
83 }

Referenced by ERF::apply_gaussian_smoothing_to_perturbations().

Here is the caller graph for this function:

◆ read_plot_file()

void read_plot_file ( PlotFileData &  pf,
const std::vector< std::string >  varnames,
MultiFab &  mf 
)
543 {
544  // ------------------------------------------------------------
545  // Open plotfile
546  // ------------------------------------------------------------
547  const std::vector<std::string>& var_names_pf = pf.varNames();
548 
549  // ------------------------------------------------------------
550  // Validate requested variables
551  // ------------------------------------------------------------
552  for (auto const& v : varnames) {
553  bool found = false;
554  for (auto const& vpf : var_names_pf) {
555  if (v == vpf) {
556  found = true;
557  break;
558  }
559  }
560  if (!found) {
561  Abort("read_plot_file: invalid variable name: " + v);
562  }
563  }
564 
565  // ------------------------------------------------------------
566  // Define destination MultiFab (single level only)
567  // ------------------------------------------------------------
568  const int level = 0;
569 
570  BoxArray ba = pf.boxArray(level);
571  DistributionMapping dm(ba);
572 
573  int ncomp = varnames.size();
574 
575  mf.define(ba, dm, ncomp, 0);
576 
577  // ------------------------------------------------------------
578  // Copy plotfile data → mf
579  // ------------------------------------------------------------
580  for (int comp = 0; comp < ncomp; ++comp)
581  {
582  const MultiFab& src = pf.get(level, varnames[comp]);
583  MultiFab::Copy(mf, src, 0, comp, 1, 0);
584  }
585 }

Referenced by ERF::ComputeAndWriteEnsemblePerturbations().

Here is the caller graph for this function:

◆ ReadCustomDataFile()

void ReadCustomDataFile ( const std::string &  filename_custom,
int &  nx,
int &  ny,
int &  nz,
int &  ng,
int &  ncomp,
std::array< Real, 3 > &  problo_ext,
std::array< Real, 3 > &  probhi_ext,
Vector< Real > &  data_rho,
Vector< Real > &  data_theta,
Vector< Real > &  data_xvel,
Vector< Real > &  data_yvel,
Vector< Real > &  data_zvel 
)
218 {
219  std::ifstream ifs(filename_custom, std::ios::binary);
220  if (!ifs.is_open()) {
221  Abort("Failed to open file " + filename_custom + " for reading");
222  }
223 
224  // ----------------------------
225  // Read header
226  // ----------------------------
227  ifs.read(reinterpret_cast<char*>(&nx), sizeof(int));
228  ifs.read(reinterpret_cast<char*>(&ny), sizeof(int));
229  ifs.read(reinterpret_cast<char*>(&nz), sizeof(int));
230 
231  ifs.read(reinterpret_cast<char*>(&ng), sizeof(int));
232  ifs.read(reinterpret_cast<char*>(&ncomp), sizeof(int));
233 
234  ifs.read(reinterpret_cast<char*>(&problo_ext[0]), sizeof(Real));
235  ifs.read(reinterpret_cast<char*>(&problo_ext[1]), sizeof(Real));
236  ifs.read(reinterpret_cast<char*>(&problo_ext[2]), sizeof(Real));
237 
238  ifs.read(reinterpret_cast<char*>(&probhi_ext[0]), sizeof(Real));
239  ifs.read(reinterpret_cast<char*>(&probhi_ext[1]), sizeof(Real));
240  ifs.read(reinterpret_cast<char*>(&probhi_ext[2]), sizeof(Real));
241 
242  const std::size_t ncell = static_cast<std::size_t>(nx) * ny * nz;
243 
244  // ----------------------------
245  // Allocate storage
246  // ----------------------------
247  data_rho.resize(ncell);
248  data_theta.resize(ncell);
249  data_xvel.resize(ncell);
250  data_yvel.resize(ncell);
251  data_zvel.resize(ncell);
252 
253  // ----------------------------
254  // Read data
255  // ----------------------------
256  std::size_t idx = 0;
257 
258  for (int k = 0; k < nz; ++k)
259  {
260  for (int j = 0; j < ny; ++j)
261  {
262  for (int i = 0; i < nx; ++i)
263  {
264  // Skip coordinates
265  Real x, y, z;
266  ifs.read(reinterpret_cast<char*>(&x), sizeof(Real));
267  ifs.read(reinterpret_cast<char*>(&y), sizeof(Real));
268  ifs.read(reinterpret_cast<char*>(&z), sizeof(Real));
269 
270  // Read components (fixed order)
271  ifs.read(reinterpret_cast<char*>(&data_rho[idx]), sizeof(Real));
272  ifs.read(reinterpret_cast<char*>(&data_theta[idx]), sizeof(Real));
273  ifs.read(reinterpret_cast<char*>(&data_xvel[idx]), sizeof(Real));
274  ifs.read(reinterpret_cast<char*>(&data_yvel[idx]), sizeof(Real));
275  ifs.read(reinterpret_cast<char*>(&data_zvel[idx]), sizeof(Real));
276 
277  ++idx;
278  }
279  }
280  }
281 
282  ifs.close();
283 }
@ ng
Definition: ERF_Morrison.H:48

Referenced by ERF::create_background_state_for_ensemble().

Here is the call graph for this function:
Here is the caller graph for this function: