ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
Utils.H File Reference
#include <AMReX.H>
#include <AMReX_MultiFab.H>
#include <AMReX_BCRec.H>
#include <DataStruct.H>
#include <IndexDefines.H>
#include <ABLMost.H>
#include <ERF_FillPatcher.H>
Include dependency graph for Utils.H:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Functions

void make_J (const amrex::Geometry &geom, amrex::MultiFab &z_phys_nd, amrex::MultiFab &detJ_cc)
 
void make_areas (const amrex::Geometry &geom, amrex::MultiFab &z_phys_nd, amrex::MultiFab &ax, amrex::MultiFab &ay, amrex::MultiFab &az)
 
void make_zcc (const amrex::Geometry &geom, amrex::MultiFab &z_phys_nd, amrex::MultiFab &z_phys_cc)
 
void MomentumToVelocity (amrex::MultiFab &xvel_out, amrex::MultiFab &yvel_out, amrex::MultiFab &zvel_out, const amrex::MultiFab &cons_in, const amrex::MultiFab &xmom_in, const amrex::MultiFab &ymom_in, const amrex::MultiFab &zmom_in, const amrex::Box &domain, const amrex::Vector< amrex::BCRec > &domain_bcs_type_h)
 
void VelocityToMomentum (const amrex::MultiFab &xvel_in, const amrex::IntVect &xvel_ngrow, const amrex::MultiFab &yvel_in, const amrex::IntVect &yvel_ngrow, const amrex::MultiFab &zvel_in, const amrex::IntVect &zvel_ngrow, const amrex::MultiFab &cons_in, amrex::MultiFab &xmom_out, amrex::MultiFab &ymom_out, amrex::MultiFab &zmom_out, const amrex::Box &domain, const amrex::Vector< amrex::BCRec > &domain_bcs_type_h)
 
void compute_interior_ghost_bxs_xy (const amrex::Box &bx, const amrex::Box &domain, const int &width, const int &set_width, amrex::Box &bx_xlo, amrex::Box &bx_xhi, amrex::Box &bx_ylo, amrex::Box &bx_yhi, const amrex::IntVect &ng_vect=amrex::IntVect(0, 0, 0), const bool get_int_ng=false)
 
void realbdy_compute_interior_ghost_rhs (const std::string &init_type, const amrex::Real &bdy_time_interval, const amrex::Real &start_bdy_time, const amrex::Real &time, const amrex::Real &delta_t, int width, int set_width, const amrex::Geometry &geom, amrex::Vector< amrex::MultiFab > &S_rhs, amrex::Vector< amrex::MultiFab > &S_old_data, amrex::Vector< amrex::MultiFab > &S_cur_data, amrex::Vector< amrex::Vector< amrex::FArrayBox >> &bdy_data_xlo, amrex::Vector< amrex::Vector< amrex::FArrayBox >> &bdy_data_xhi, amrex::Vector< amrex::Vector< amrex::FArrayBox >> &bdy_data_ylo, amrex::Vector< amrex::Vector< amrex::FArrayBox >> &bdy_data_yhi)
 
void fine_compute_interior_ghost_rhs (const amrex::Real &time, const amrex::Real &delta_t, const int &width, const int &set_width, const amrex::Geometry &geom, ERFFillPatcher *FPr_c, ERFFillPatcher *FPr_u, ERFFillPatcher *FPr_v, ERFFillPatcher *FPr_w, amrex::Vector< amrex::BCRec > &domain_bcs_type, amrex::Vector< amrex::MultiFab > &S_rhs_f, amrex::Vector< amrex::MultiFab > &S_data_f)
 
void Time_Avg_Vel_atCC (const amrex::Real &dt, amrex::Real &t_avg_cnt, amrex::MultiFab *vel_t_avg, amrex::MultiFab &xvel, amrex::MultiFab &yvel, amrex::MultiFab &zvel)
 
AMREX_GPU_HOST AMREX_FORCE_INLINE void wrfbdy_set_rhs_in_spec_region (const amrex::Real &dt, const int &icomp, const int &num_var, const int &width, const int &set_width, const amrex::Dim3 &dom_lo, const amrex::Dim3 &dom_hi, const amrex::Box &bx_xlo, const amrex::Box &bx_xhi, const amrex::Box &bx_ylo, const amrex::Box &bx_yhi, const amrex::Array4< const amrex::Real > &arr_xlo, const amrex::Array4< const amrex::Real > &arr_xhi, const amrex::Array4< const amrex::Real > &arr_ylo, const amrex::Array4< const amrex::Real > &arr_yhi, const amrex::Array4< const amrex::Real > &data_arr, const amrex::Array4< amrex::Real > &rhs_arr)
 
AMREX_GPU_HOST AMREX_FORCE_INLINE void wrfbdy_compute_laplacian_relaxation (const int &icomp, const int &num_var, const int &width, const int &set_width, const amrex::Dim3 &dom_lo, const amrex::Dim3 &dom_hi, const amrex::Real &F1, const amrex::Real &F2, const amrex::Box &bx_xlo, const amrex::Box &bx_xhi, const amrex::Box &bx_ylo, const amrex::Box &bx_yhi, const amrex::Array4< const amrex::Real > &arr_xlo, const amrex::Array4< const amrex::Real > &arr_xhi, const amrex::Array4< const amrex::Real > &arr_ylo, const amrex::Array4< const amrex::Real > &arr_yhi, const amrex::Array4< const amrex::Real > &data_arr, const amrex::Array4< amrex::Real > &rhs_arr)
 
AMREX_GPU_HOST AMREX_FORCE_INLINE void ApplyMask (amrex::MultiFab &dst, const amrex::iMultiFab &imask, const int nghost=0)
 
AMREX_GPU_HOST AMREX_FORCE_INLINE void ApplyInvertedMask (amrex::MultiFab &dst, const amrex::iMultiFab &imask, const int nghost=0)
 

Function Documentation

◆ ApplyInvertedMask()

AMREX_GPU_HOST AMREX_FORCE_INLINE void ApplyInvertedMask ( amrex::MultiFab &  dst,
const amrex::iMultiFab &  imask,
const int  nghost = 0 
)
375 {
376  for (amrex::MFIter mfi(dst,amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
377  {
378  const amrex::Box& bx = mfi.growntilebox(nghost);
379  if (bx.ok())
380  {
381  auto dstFab = dst.array(mfi);
382  const auto maskFab = imask.const_array(mfi);
383  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
384  {
385  dstFab(i,j,k) *= (1-maskFab(i,j,k));
386  });
387  }
388  }
389 }

Referenced by add_thin_body_sources().

Here is the caller graph for this function:

◆ ApplyMask()

AMREX_GPU_HOST AMREX_FORCE_INLINE void ApplyMask ( amrex::MultiFab &  dst,
const amrex::iMultiFab &  imask,
const int  nghost = 0 
)
354 {
355  for (amrex::MFIter mfi(dst,amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
356  {
357  const amrex::Box& bx = mfi.growntilebox(nghost);
358  if (bx.ok())
359  {
360  auto dstFab = dst.array(mfi);
361  const auto maskFab = imask.const_array(mfi);
362  ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
363  {
364  dstFab(i,j,k) *= maskFab(i,j,k);
365  });
366  }
367  }
368 }

Referenced by ERF::FillIntermediatePatch().

Here is the caller graph for this function:

◆ compute_interior_ghost_bxs_xy()

void compute_interior_ghost_bxs_xy ( const amrex::Box &  bx,
const amrex::Box &  domain,
const int &  width,
const int &  set_width,
amrex::Box &  bx_xlo,
amrex::Box &  bx_xhi,
amrex::Box &  bx_ylo,
amrex::Box &  bx_yhi,
const amrex::IntVect &  ng_vect = amrex::IntVect(0, 0, 0),
const bool  get_int_ng = false 
)

◆ fine_compute_interior_ghost_rhs()

void fine_compute_interior_ghost_rhs ( const amrex::Real &  time,
const amrex::Real &  delta_t,
const int &  width,
const int &  set_width,
const amrex::Geometry &  geom,
ERFFillPatcher FPr_c,
ERFFillPatcher FPr_u,
ERFFillPatcher FPr_v,
ERFFillPatcher FPr_w,
amrex::Vector< amrex::BCRec > &  domain_bcs_type,
amrex::Vector< amrex::MultiFab > &  S_rhs_f,
amrex::Vector< amrex::MultiFab > &  S_data_f 
)

◆ make_areas()

void make_areas ( const amrex::Geometry &  geom,
amrex::MultiFab &  z_phys_nd,
amrex::MultiFab &  ax,
amrex::MultiFab &  ay,
amrex::MultiFab &  az 
)

◆ make_J()

void make_J ( const amrex::Geometry &  geom,
amrex::MultiFab &  z_phys_nd,
amrex::MultiFab &  detJ_cc 
)

◆ make_zcc()

void make_zcc ( const amrex::Geometry &  geom,
amrex::MultiFab &  z_phys_nd,
amrex::MultiFab &  z_phys_cc 
)

◆ MomentumToVelocity()

void MomentumToVelocity ( amrex::MultiFab &  xvel_out,
amrex::MultiFab &  yvel_out,
amrex::MultiFab &  zvel_out,
const amrex::MultiFab &  cons_in,
const amrex::MultiFab &  xmom_in,
const amrex::MultiFab &  ymom_in,
const amrex::MultiFab &  zmom_in,
const amrex::Box &  domain,
const amrex::Vector< amrex::BCRec > &  domain_bcs_type_h 
)

◆ realbdy_compute_interior_ghost_rhs()

void realbdy_compute_interior_ghost_rhs ( const std::string &  init_type,
const amrex::Real &  bdy_time_interval,
const amrex::Real &  start_bdy_time,
const amrex::Real &  time,
const amrex::Real &  delta_t,
int  width,
int  set_width,
const amrex::Geometry &  geom,
amrex::Vector< amrex::MultiFab > &  S_rhs,
amrex::Vector< amrex::MultiFab > &  S_old_data,
amrex::Vector< amrex::MultiFab > &  S_cur_data,
amrex::Vector< amrex::Vector< amrex::FArrayBox >> &  bdy_data_xlo,
amrex::Vector< amrex::Vector< amrex::FArrayBox >> &  bdy_data_xhi,
amrex::Vector< amrex::Vector< amrex::FArrayBox >> &  bdy_data_ylo,
amrex::Vector< amrex::Vector< amrex::FArrayBox >> &  bdy_data_yhi 
)

◆ Time_Avg_Vel_atCC()

void Time_Avg_Vel_atCC ( const amrex::Real &  dt,
amrex::Real &  t_avg_cnt,
amrex::MultiFab *  vel_t_avg,
amrex::MultiFab &  xvel,
amrex::MultiFab &  yvel,
amrex::MultiFab &  zvel 
)

◆ VelocityToMomentum()

void VelocityToMomentum ( const amrex::MultiFab &  xvel_in,
const amrex::IntVect &  xvel_ngrow,
const amrex::MultiFab &  yvel_in,
const amrex::IntVect &  yvel_ngrow,
const amrex::MultiFab &  zvel_in,
const amrex::IntVect &  zvel_ngrow,
const amrex::MultiFab &  cons_in,
amrex::MultiFab &  xmom_out,
amrex::MultiFab &  ymom_out,
amrex::MultiFab &  zmom_out,
const amrex::Box &  domain,
const amrex::Vector< amrex::BCRec > &  domain_bcs_type_h 
)

Referenced by ERF::advance_dycore(), ERF::AverageDownTo(), ERF::FillCoarsePatch(), ERF::FillIntermediatePatch(), and ERF::FillPatch().

Here is the caller graph for this function:

◆ wrfbdy_compute_laplacian_relaxation()

AMREX_GPU_HOST AMREX_FORCE_INLINE void wrfbdy_compute_laplacian_relaxation ( const int &  icomp,
const int &  num_var,
const int &  width,
const int &  set_width,
const amrex::Dim3 &  dom_lo,
const amrex::Dim3 &  dom_hi,
const amrex::Real &  F1,
const amrex::Real &  F2,
const amrex::Box &  bx_xlo,
const amrex::Box &  bx_xhi,
const amrex::Box &  bx_ylo,
const amrex::Box &  bx_yhi,
const amrex::Array4< const amrex::Real > &  arr_xlo,
const amrex::Array4< const amrex::Real > &  arr_xhi,
const amrex::Array4< const amrex::Real > &  arr_ylo,
const amrex::Array4< const amrex::Real > &  arr_yhi,
const amrex::Array4< const amrex::Real > &  data_arr,
const amrex::Array4< amrex::Real > &  rhs_arr 
)

Compute the Laplacian RHS in the relaxation zone

Parameters
[in]delta_ttime step
[in]icompcomponent offset
[in]num_varnumber of variables to loop
[in]widthwidth of wrf bdy file
[in]set_widthwidth the set region
[in]dom_lolow bound of domain
[in]dom_hihigh bound of domain
[in]F1drift relaxation parameter
[in]F2Laplacian relaxation parameter
[in]bx_xlobox for low x relaxation
[in]bx_xhibox for high x relaxation
[in]bx_ylobox for low y relaxation
[in]bx_yhibox for high y relaxation
[in]arr_xloarray for low x relaxation
[in]arr_xhiarray for high x relaxation
[in]arr_yloarray for low y relaxation
[in]arr_yhiarray for high y relaxation
[in]data_arrdata array
[out]rhs_arrRHS array
245 {
246  // RHS computation
247  int Spec_z = set_width;
248  int Relax_z = width - Spec_z + 1;
249  amrex::Real num = amrex::Real(Spec_z + Relax_z);
250  amrex::Real denom = amrex::Real(Relax_z - 1);
251  amrex::Real SpecExp = -std::log(0.01) / amrex::Real(width - Spec_z);
252  amrex::ParallelFor(bx_xlo, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
253  {
254  // Corners with x boxes
255  int j_lo = std::min(j-dom_lo.y,width-1);
256  int j_hi = std::min(dom_hi.y-j,width-1);
257  int jj = std::min(j_lo,j_hi);
258  int n_ind = std::min(i-dom_lo.x,jj) + 1;
259  AMREX_ASSERT(n_ind > Spec_z);
260  amrex::Real Factor = (num - amrex::Real(n_ind))/denom
261  * std::exp(-SpecExp * amrex::Real(n_ind - Spec_z));
262  amrex::Real d = data_arr(i ,j ,k ,n+icomp);
263  amrex::Real d_ip1 = data_arr(i+1,j ,k ,n+icomp);
264  amrex::Real d_im1 = data_arr(i-1,j ,k ,n+icomp);
265  amrex::Real d_jp1 = data_arr(i ,j+1,k ,n+icomp);
266  amrex::Real d_jm1 = data_arr(i ,j-1,k ,n+icomp);
267  amrex::Real delta = arr_xlo(i ,j ,k,n) - d;
268  amrex::Real delta_xp = arr_xlo(i+1,j ,k,n) - d_ip1;
269  amrex::Real delta_xm = arr_xlo(i-1,j ,k,n) - d_im1;
270  amrex::Real delta_yp = arr_xlo(i ,j+1,k,n) - d_jp1;
271  amrex::Real delta_ym = arr_xlo(i ,j-1,k,n) - d_jm1;
272  amrex::Real Laplacian = delta_xp + delta_xm + delta_yp + delta_ym - 4.0*delta;
273  amrex::Real Temp = (F1*delta - F2*Laplacian) * Factor;
274  rhs_arr(i,j,k,n+icomp) += Temp;
275  },
276  bx_xhi, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
277  {
278  // Corners with x boxes
279  int j_lo = std::min(j-dom_lo.y,width-1);
280  int j_hi = std::min(dom_hi.y-j,width-1);
281  int jj = std::min(j_lo,j_hi);
282  int n_ind = std::min(dom_hi.x-i,jj) + 1;
283  AMREX_ASSERT(n_ind > Spec_z);
284  amrex::Real Factor = (num - amrex::Real(n_ind))/denom
285  * std::exp(-SpecExp * amrex::Real(n_ind - Spec_z));
286  amrex::Real d = data_arr(i ,j ,k ,n+icomp);
287  amrex::Real d_ip1 = data_arr(i+1,j ,k ,n+icomp);
288  amrex::Real d_im1 = data_arr(i-1,j ,k ,n+icomp);
289  amrex::Real d_jp1 = data_arr(i ,j+1,k ,n+icomp);
290  amrex::Real d_jm1 = data_arr(i ,j-1,k ,n+icomp);
291  amrex::Real delta = arr_xhi(i ,j ,k,n) - d;
292  amrex::Real delta_xp = arr_xhi(i+1,j ,k,n) - d_ip1;
293  amrex::Real delta_xm = arr_xhi(i-1,j ,k,n) - d_im1;
294  amrex::Real delta_yp = arr_xhi(i ,j+1,k,n) - d_jp1;
295  amrex::Real delta_ym = arr_xhi(i ,j-1,k,n) - d_jm1;
296  amrex::Real Laplacian = delta_xp + delta_xm + delta_yp + delta_ym - 4.0*delta;
297  amrex::Real Temp = (F1*delta - F2*Laplacian) * Factor;
298  rhs_arr(i,j,k,n+icomp) += Temp;
299  });
300 
301  amrex::ParallelFor(bx_ylo, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
302  {
303  // No corners for y boxes
304  int n_ind = j - dom_lo.y + 1;
305  AMREX_ASSERT(n_ind > Spec_z);
306  amrex::Real Factor = (num - amrex::Real(n_ind))/denom
307  * std::exp(-SpecExp * amrex::Real(n_ind - Spec_z));
308  amrex::Real d = data_arr(i ,j ,k ,n+icomp);
309  amrex::Real d_ip1 = data_arr(i+1,j ,k ,n+icomp);
310  amrex::Real d_im1 = data_arr(i-1,j ,k ,n+icomp);
311  amrex::Real d_jp1 = data_arr(i ,j+1,k ,n+icomp);
312  amrex::Real d_jm1 = data_arr(i ,j-1,k ,n+icomp);
313  amrex::Real delta = arr_ylo(i ,j ,k,n) - d;
314  amrex::Real delta_xp = arr_ylo(i+1,j ,k,n) - d_ip1;
315  amrex::Real delta_xm = arr_ylo(i-1,j ,k,n) - d_im1;
316  amrex::Real delta_yp = arr_ylo(i ,j+1,k,n) - d_jp1;
317  amrex::Real delta_ym = arr_ylo(i ,j-1,k,n) - d_jm1;
318  amrex::Real Laplacian = delta_xp + delta_xm + delta_yp + delta_ym - 4.0*delta;
319  amrex::Real Temp = (F1*delta - F2*Laplacian) * Factor;
320  rhs_arr(i,j,k,n+icomp) += Temp;
321  },
322  bx_yhi, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
323  {
324  // No corners for y boxes
325  int n_ind = dom_hi.y - j + 1;
326  AMREX_ASSERT(n_ind > Spec_z);
327  amrex::Real Factor = (num - amrex::Real(n_ind))/denom
328  * std::exp(-SpecExp * amrex::Real(n_ind - Spec_z));
329  amrex::Real d = data_arr(i ,j ,k ,n+icomp);
330  amrex::Real d_ip1 = data_arr(i+1,j ,k ,n+icomp);
331  amrex::Real d_im1 = data_arr(i-1,j ,k ,n+icomp);
332  amrex::Real d_jp1 = data_arr(i ,j+1,k ,n+icomp);
333  amrex::Real d_jm1 = data_arr(i ,j-1,k ,n+icomp);
334  amrex::Real delta = arr_yhi(i ,j ,k,n) - d;
335  amrex::Real delta_xp = arr_yhi(i+1,j ,k,n) - d_ip1;
336  amrex::Real delta_xm = arr_yhi(i-1,j ,k,n) - d_im1;
337  amrex::Real delta_yp = arr_yhi(i ,j+1,k,n) - d_jp1;
338  amrex::Real delta_ym = arr_yhi(i ,j-1,k,n) - d_jm1;
339  amrex::Real Laplacian = delta_xp + delta_xm + delta_yp + delta_ym - 4.0*delta;
340  amrex::Real Temp = (F1*delta - F2*Laplacian) * Factor;
341  rhs_arr(i,j,k,n+icomp) += Temp;
342  });
343 }

Referenced by realbdy_compute_interior_ghost_rhs().

Here is the caller graph for this function:

◆ wrfbdy_set_rhs_in_spec_region()

AMREX_GPU_HOST AMREX_FORCE_INLINE void wrfbdy_set_rhs_in_spec_region ( const amrex::Real &  dt,
const int &  icomp,
const int &  num_var,
const int &  width,
const int &  set_width,
const amrex::Dim3 &  dom_lo,
const amrex::Dim3 &  dom_hi,
const amrex::Box &  bx_xlo,
const amrex::Box &  bx_xhi,
const amrex::Box &  bx_ylo,
const amrex::Box &  bx_yhi,
const amrex::Array4< const amrex::Real > &  arr_xlo,
const amrex::Array4< const amrex::Real > &  arr_xhi,
const amrex::Array4< const amrex::Real > &  arr_ylo,
const amrex::Array4< const amrex::Real > &  arr_yhi,
const amrex::Array4< const amrex::Real > &  data_arr,
const amrex::Array4< amrex::Real > &  rhs_arr 
)

Zero RHS in the set region

Parameters
[in]icompcomponent offset
[in]num_varnumber of variables to loop
[in]bx_xlobox for low x relaxation
[in]bx_xhibox for high x relaxation
[in]bx_ylobox for low y relaxation
[in]bx_yhibox for high y relaxation
[out]rhs_arrRHS array
154 {
155  int Spec_z = set_width;
156  amrex::ParallelFor(bx_xlo, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
157  {
158  // Corners with x boxes
159  int j_lo = std::min(j-dom_lo.y,width-1);
160  int j_hi = std::min(dom_hi.y-j,width-1);
161  int jj = std::min(j_lo,j_hi);
162  int n_ind = std::min(i-dom_lo.x,jj) + 1;
163  if (n_ind <= Spec_z) {
164  amrex::Real tend = ( arr_xlo(i,j,k) - data_arr(i,j,k,n+icomp) ) / dt;
165  rhs_arr(i,j,k,n+icomp) = tend;
166  }
167  },
168  bx_xhi, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
169  {
170  // Corners with x boxes
171  int j_lo = std::min(j-dom_lo.y,width-1);
172  int j_hi = std::min(dom_hi.y-j,width-1);
173  int jj = std::min(j_lo,j_hi);
174  int n_ind = std::min(dom_hi.x-i,jj) + 1;
175  if (n_ind <= Spec_z) {
176  amrex::Real tend = ( arr_xhi(i,j,k) - data_arr(i,j,k,n+icomp) ) / dt;
177  rhs_arr(i,j,k,n+icomp) = tend;
178  }
179  });
180 
181  amrex::ParallelFor(bx_ylo, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
182  {
183  // No corners for y boxes
184  int n_ind = j - dom_lo.y + 1;
185  if (n_ind <= Spec_z) {
186  amrex::Real tend = ( arr_ylo(i,j,k) - data_arr(i,j,k,n+icomp) ) / dt;
187  rhs_arr(i,j,k,n+icomp) = tend;
188  }
189  },
190  bx_yhi, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
191  {
192  // No corners for y boxes
193  int n_ind = dom_hi.y - j + 1;
194  if (n_ind <= Spec_z) {
195  amrex::Real tend = ( arr_yhi(i,j,k) - data_arr(i,j,k,n+icomp) ) / dt;
196  rhs_arr(i,j,k,n+icomp) = tend;
197  }
198  });
199 }

Referenced by realbdy_compute_interior_ghost_rhs().

Here is the caller graph for this function: