ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_InteriorGhostCells.cpp File Reference
#include <ERF_Utils.H>
Include dependency graph for ERF_InteriorGhostCells.cpp:

Functions

void realbdy_interior_bxs_xy (const Box &bx, const Box &domain, const int &width, Box &bx_xlo, Box &bx_xhi, Box &bx_ylo, Box &bx_yhi, const IntVect &ng_vect, const bool get_int_ng)
 
void realbdy_compute_interior_ghost_rhs (const Real &time, const Real &delta_t, const Real &start_bdy_time, const Real &final_bdy_time, const Real &bdy_time_interval, const Real &nudge_factor, int width, bool do_upwind, const Geometry &geom, Vector< MultiFab > &S_rhs, Vector< MultiFab > &S_cur_data, Vector< Vector< FArrayBox >> &bdy_data_xlo, Vector< Vector< FArrayBox >> &bdy_data_xhi, Vector< Vector< FArrayBox >> &bdy_data_ylo, Vector< Vector< FArrayBox >> &bdy_data_yhi)
 
void fine_compute_interior_ghost_rhs (const Real &time, const Real &delta_t, const int &width, const int &set_width, const Geometry &geom, ERFFillPatcher *FPr_c, ERFFillPatcher *FPr_u, ERFFillPatcher *FPr_v, ERFFillPatcher *FPr_w, Vector< BCRec > &domain_bcs_type, Vector< MultiFab > &S_rhs_f, Vector< MultiFab > &S_data_f)
 

Variables

PhysBCFunctNoOp void_bc
 

Function Documentation

◆ fine_compute_interior_ghost_rhs()

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

Compute the RHS in the fine relaxation zone

Parameters
[in]timecurrent (elapsed) time
[in]delta_ttimestep
[in]widthnumber of cells in (relaxation+specified) zone
[in]set_widthnumber of cells in (specified) zone
[in]FPr_ccons fine patch container
[in]FPr_uuvel fine patch container
[in]FPr_vvvel fine patch container
[in]FPr_wwvel fine patch container
[in]boxes_at_levelboxes at current level
[in]domain_bcs_typeboundary condition types
[out]S_rhsRHS to be computed here
[in]S_datacurrent value of the solution
469 {
470  BL_PROFILE_REGION("fine_compute_interior_ghost_RHS()");
471 
472  // Relaxation constants
473  Real F1 = 1./(10.*delta_t);
474  Real F2 = 1./(50.*delta_t);
475 
476  // Vector of MFs to hold data (dm differs w/ fine patch)
477  Vector<MultiFab> fmf_p_v;
478 
479  // Loop over the variables
480  for (int ivar_idx = 0; ivar_idx < IntVars::NumTypes; ++ivar_idx)
481  {
482  // Fine mfs
483  MultiFab& fmf = S_data_f[ivar_idx];
484  MultiFab& rhs = S_rhs_f [ivar_idx];
485 
486  // NOTE: These temporary MFs and copy operations are horrible
487  // for memory usage and efficiency. However, we need to
488  // have access to ghost cells in the cons array to convert
489  // from primitive u/v/w to momentum. Furthermore, the BA
490  // for the fine patches in ERFFillPatcher don't match the
491  // BA for the data/RHS. For this reason, the data is copied
492  // to a vector of MFs (with ghost cells) so the BAs match
493  // the BA of data/RHS and we have access to rho to convert
494  // prim to conserved.
495 
496  // Temp MF on box (distribution map differs w/ fine patch)
497  int num_var = fmf.nComp();
498  fmf_p_v.emplace_back(fmf.boxArray(), fmf.DistributionMap(), num_var, fmf.nGrowVect());
499  MultiFab& fmf_p = fmf_p_v[ivar_idx];
500  MultiFab::Copy(fmf_p,fmf, 0, 0, num_var, fmf.nGrowVect());
501 
502  // Integer mask MF
503  int set_mask_val;
504  int relax_mask_val;
505  iMultiFab* mask;
506 
507  // Fill fine patch on interior halo region
508  //==========================================================
509  if (ivar_idx == IntVars::cons)
510  {
511  FPr_c->FillRelax(fmf_p, time, void_bc, domain_bcs_type);
512  mask = FPr_c->GetMask();
513  set_mask_val = FPr_c->GetSetMaskVal();
514  relax_mask_val = FPr_c->GetRelaxMaskVal();
515  }
516  else if (ivar_idx == IntVars::xmom)
517  {
518  FPr_u->FillRelax(fmf_p, time, void_bc, domain_bcs_type);
519  mask = FPr_u->GetMask();
520  set_mask_val = FPr_u->GetSetMaskVal();
521  relax_mask_val = FPr_u->GetRelaxMaskVal();
522 
523 #ifdef _OPENMP
524 #pragma omp parallel if (Gpu::notInLaunchRegion())
525 #endif
526  for ( MFIter mfi(fmf_p,TilingIfNotGPU()); mfi.isValid(); ++mfi)
527  {
528  Box tbx = mfi.tilebox();
529 
530  const Array4<Real>& prim_arr = fmf_p.array(mfi);
531  const Array4<const Real>& rho_arr = fmf_p_v[0].const_array(mfi);
532  const Array4<const int>& mask_arr = mask->const_array(mfi);
533 
534  ParallelFor(tbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
535  {
536  if (mask_arr(i,j,k) == relax_mask_val) {
537  Real rho_interp = 0.5 * ( rho_arr(i-1,j,k) + rho_arr(i,j,k) );
538  prim_arr(i,j,k) *= rho_interp;
539  }
540  });
541  } // mfi
542  }
543  else if (ivar_idx == IntVars::ymom)
544  {
545  FPr_v->FillRelax(fmf_p, time, void_bc, domain_bcs_type);
546  mask = FPr_v->GetMask();
547  set_mask_val = FPr_v->GetSetMaskVal();
548  relax_mask_val = FPr_v->GetRelaxMaskVal();
549 
550 #ifdef _OPENMP
551 #pragma omp parallel if (Gpu::notInLaunchRegion())
552 #endif
553  for ( MFIter mfi(fmf_p,TilingIfNotGPU()); mfi.isValid(); ++mfi)
554  {
555  Box tbx = mfi.tilebox();
556 
557  const Array4<Real>& prim_arr = fmf_p.array(mfi);
558  const Array4<const Real>& rho_arr = fmf_p_v[0].const_array(mfi);
559  const Array4<const int>& mask_arr = mask->const_array(mfi);
560 
561  ParallelFor(tbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
562  {
563  if (mask_arr(i,j,k) == relax_mask_val) {
564  Real rho_interp = 0.5 * ( rho_arr(i,j-1,k) + rho_arr(i,j,k) );
565  prim_arr(i,j,k) *= rho_interp;
566  }
567  });
568  } // mfi
569  }
570  else if (ivar_idx == IntVars::zmom)
571  {
572  FPr_w->FillRelax(fmf_p, time, void_bc, domain_bcs_type);
573  mask = FPr_w->GetMask();
574  set_mask_val = FPr_w->GetSetMaskVal();
575  relax_mask_val = FPr_w->GetRelaxMaskVal();
576 
577 #ifdef _OPENMP
578 #pragma omp parallel if (Gpu::notInLaunchRegion())
579 #endif
580  for ( MFIter mfi(fmf_p,TilingIfNotGPU()); mfi.isValid(); ++mfi)
581  {
582  Box tbx = mfi.tilebox();
583 
584  const Array4<Real>& prim_arr = fmf_p.array(mfi);
585  const Array4<const Real>& rho_arr = fmf_p_v[0].const_array(mfi);
586  const Array4<const int>& mask_arr = mask->const_array(mfi);
587 
588  ParallelFor(tbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
589  {
590  if (mask_arr(i,j,k) == relax_mask_val) {
591  Real rho_interp = 0.5 * ( rho_arr(i,j,k-1) + rho_arr(i,j,k) );
592  prim_arr(i,j,k) *= rho_interp;
593  }
594  });
595  } // mfi
596  } else {
597  Abort("Dont recognize this variable type in fine_compute_interior_ghost_RHS");
598  }
599 
600 
601  // Zero RHS in set region
602  //==========================================================
603 #ifdef _OPENMP
604 #pragma omp parallel if (Gpu::notInLaunchRegion())
605 #endif
606  for ( MFIter mfi(rhs,TilingIfNotGPU()); mfi.isValid(); ++mfi)
607  {
608  Box tbx = mfi.tilebox();
609  const Array4<Real>& rhs_arr = rhs.array(mfi);
610  const Array4<const int>& mask_arr = mask->const_array(mfi);
611 
612  ParallelFor(tbx, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
613  {
614  if (mask_arr(i,j,k) == set_mask_val) {
615  rhs_arr(i,j,k) = 0.0;
616  }
617  });
618  } // mfi
619 
620  // For Laplacian stencil
621  rhs.FillBoundary(geom.periodicity());
622 
623 
624  // Compute RHS in relaxation region
625  //==========================================================
626 #ifdef _OPENMP
627 #pragma omp parallel if (Gpu::notInLaunchRegion())
628 #endif
629  for ( MFIter mfi(fmf_p,TilingIfNotGPU()); mfi.isValid(); ++mfi)
630  {
631  Box tbx = mfi.tilebox();
632  const Array4<Real>& rhs_arr = rhs.array(mfi);
633  const Array4<const Real>& fine_arr = fmf_p.const_array(mfi);
634  const Array4<const Real>& data_arr = fmf.const_array(mfi);
635  const Array4<const int>& mask_arr = mask->const_array(mfi);
636 
637  Box vbx = mfi.validbox();
638  const auto& vbx_lo = lbound(vbx);
639  const auto& vbx_hi = ubound(vbx);
640 
641  int icomp = 0;
642 
643  int Spec_z = set_width;
644  int Relax_z = width - Spec_z;
645  Real num = Real(Spec_z + Relax_z);
646  Real denom = Real(Relax_z - 1);
647  ParallelFor(tbx, num_var, [=] AMREX_GPU_DEVICE (int i, int j, int k, int n) noexcept
648  {
649  if (mask_arr(i,j,k) == relax_mask_val) {
650 
651  // Indices
652  Real n_ind(-1); // Set to -1 to quiet compiler warning
653  int ii(width-1); int jj(width-1);
654  bool near_x_lo_wall(false); bool near_x_hi_wall(false);
655  bool near_y_lo_wall(false); bool near_y_hi_wall(false);
656  bool mask_x_found(false); bool mask_y_found(false);
657 
658  // Near x-wall
659  if ((i-vbx_lo.x) < width) {
660  near_x_lo_wall = true;
661  ii = i-vbx_lo.x;
662  if (mask_arr(vbx_lo.x,j,k) == 2) mask_x_found = true;
663  } else if ((vbx_hi.x-i) < width) {
664  near_x_hi_wall = true;
665  ii = vbx_hi.x-i;
666  if (mask_arr(vbx_hi.x,j,k) == 2) mask_x_found = true;
667  }
668 
669  // Near y-wall
670  if ((j-vbx_lo.y) < width) {
671  near_y_lo_wall = true;
672  jj = j-vbx_lo.y;
673  if (mask_arr(i,vbx_lo.y,k) == 2) mask_y_found = true;
674  } else if ((vbx_hi.y-j) < width) {
675  near_y_hi_wall = true;
676  jj = vbx_hi.y-j;
677  if (mask_arr(i,vbx_hi.y,k) == 2) mask_y_found = true;
678  }
679 
680  // Found a nearby masked cell (valid n_ind)
681  if (mask_x_found && mask_y_found) {
682  n_ind = std::min(ii,jj) + 1.0;
683  } else if (mask_x_found) {
684  n_ind = ii + 1.0;
685  } else if (mask_y_found) {
686  n_ind = jj + 1.0;
687  // Pesky corner cell
688  } else {
689  if (near_x_lo_wall || near_x_hi_wall) {
690  Real dj_min{width-1.0};
691  int j_lb = std::max(vbx_lo.y,j-width);
692  int j_ub = std::min(vbx_hi.y,j+width);
693  int li = (near_x_lo_wall) ? vbx_lo.x : vbx_hi.x;
694  for (int lj(j_lb); lj<=j_ub; ++lj) {
695  if (mask_arr(li,lj,k) == 2) {
696  mask_y_found = true;
697  dj_min = std::min(dj_min,(Real) std::abs(lj-j));
698  }
699  }
700  if (mask_y_found) {
701  Real mag = sqrt( Real(dj_min*dj_min + ii*ii) );
702  n_ind = std::min(mag,width-1.0) + 1.0;
703  } else {
704  Abort("Mask not found near x wall!");
705  }
706  } else if (near_y_lo_wall || near_y_hi_wall) {
707  Real di_min{width-1.0};
708  int i_lb = std::max(vbx_lo.x,i-width);
709  int i_ub = std::min(vbx_hi.x,i+width);
710  int lj = (near_y_lo_wall) ? vbx_lo.y : vbx_hi.y;
711  for (int li(i_lb); li<=i_ub; ++li) {
712  if (mask_arr(li,lj,k) == 2) {
713  mask_x_found = true;
714  di_min = std::min(di_min,(Real) std::abs(li-i));
715  }
716  }
717  if (mask_x_found) {
718  Real mag = sqrt( Real(di_min*di_min + jj*jj) );
719  n_ind = std::min(mag,width-1.0) + 1.0;
720  } else {
721  Abort("Mask not found near y wall!");
722  }
723  } else {
724  Abort("Relaxation cell must be near a wall!");
725  }
726  }
727 
728  Real Factor = (num - n_ind)/denom;
729  Real d = data_arr(i ,j ,k ,n+icomp) + delta_t*rhs_arr(i , j , k ,n+icomp);
730  Real d_ip1 = data_arr(i+1,j ,k ,n+icomp) + delta_t*rhs_arr(i+1, j , k ,n+icomp);
731  Real d_im1 = data_arr(i-1,j ,k ,n+icomp) + delta_t*rhs_arr(i-1, j , k ,n+icomp);
732  Real d_jp1 = data_arr(i ,j+1,k ,n+icomp) + delta_t*rhs_arr(i , j+1, k ,n+icomp);
733  Real d_jm1 = data_arr(i ,j-1,k ,n+icomp) + delta_t*rhs_arr(i , j-1, k ,n+icomp);
734  Real delta = fine_arr(i ,j ,k,n) - d;
735  Real delta_xp = fine_arr(i+1,j ,k,n) - d_ip1;
736  Real delta_xm = fine_arr(i-1,j ,k,n) - d_im1;
737  Real delta_yp = fine_arr(i ,j+1,k,n) - d_jp1;
738  Real delta_ym = fine_arr(i ,j-1,k,n) - d_jm1;
739  Real Laplacian = delta_xp + delta_xm + delta_yp + delta_ym - 4.0*delta;
740  rhs_arr(i,j,k,n) += (F1*delta - F2*Laplacian) * Factor;
741  }
742  });
743  } // mfi
744  } // ivar_idx
745 }
@ num
Definition: ERF_DataStruct.H:24
ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept { const auto prob_lo=geomdata.ProbLo();const auto dx=geomdata.CellSize();const Real x=(prob_lo[0]+(i+0.5) *dx[0])/mf_m(i, j, 0);const Real z=z_cc(i, j, k);Real L=std::sqrt(std::pow((x - x_c)/x_r, 2)+std::pow((z - z_c)/z_r, 2));if(L<=1.0) { Real dT=T_pert *(std::cos(PI *L)+1.0)/2.0;Real Tbar_hse=p_hse(i, j, k)/(R_d *r_hse(i, j, k));Real theta_perturbed=(Tbar_hse+dT) *std::pow(p_0/p_hse(i, j, k), rdOcp);Real theta_0=(Tbar_hse) *std::pow(p_0/p_hse(i, j, k), rdOcp);if(const_rho) { state_pert(i, j, k, RhoTheta_comp)=r_hse(i, j, k) *(theta_perturbed - theta_0);} else { state_pert(i, j, k, Rho_comp)=getRhoThetagivenP(p_hse(i, j, k))/theta_perturbed - r_hse(i, j, k);} } })
PhysBCFunctNoOp void_bc
Definition: ERF_InteriorGhostCells.cpp:5
amrex::Real Real
Definition: ERF_ShocInterface.H:19
amrex::iMultiFab * GetMask()
Definition: ERF_FillPatcher.H:43
void FillRelax(amrex::MultiFab &mf, amrex::Real time, BC &cbc, amrex::Vector< amrex::BCRec > const &bcs)
Definition: ERF_FillPatcher.H:105
int GetSetMaskVal()
Definition: ERF_FillPatcher.H:39
int GetRelaxMaskVal()
Definition: ERF_FillPatcher.H:41
@ NumTypes
Definition: ERF_IndexDefines.H:162
@ ymom
Definition: ERF_IndexDefines.H:160
@ cons
Definition: ERF_IndexDefines.H:158
@ zmom
Definition: ERF_IndexDefines.H:161
@ xmom
Definition: ERF_IndexDefines.H:159
Here is the call graph for this function:

◆ realbdy_compute_interior_ghost_rhs()

void realbdy_compute_interior_ghost_rhs ( const Real time,
const Real delta_t,
const Real start_bdy_time,
const Real final_bdy_time,
const Real bdy_time_interval,
const Real nudge_factor,
int  width,
bool  do_upwind,
const Geometry &  geom,
Vector< MultiFab > &  S_rhs,
Vector< MultiFab > &  S_cur_data,
Vector< Vector< FArrayBox >> &  bdy_data_xlo,
Vector< Vector< FArrayBox >> &  bdy_data_xhi,
Vector< Vector< FArrayBox >> &  bdy_data_ylo,
Vector< Vector< FArrayBox >> &  bdy_data_yhi 
)

Compute the RHS in the relaxation zone

Parameters
[in]timecurrent (total) time
[in]delta_ttimestep
[in]start_bdy_timefull time of the first time slice of boundary data
[in]final_bdy_timefull time of the last time slice of boundary data
[in]bdy_time_intervaltime interval between boundary condition time stamps
[in]widthnumber of cells in (relaxation+specified) zone
[in]set_widthnumber of cells in (specified) zone
[in]geomcontainer for geometric information
[out]S_rhsRHS to be computed here
[in]S_datacurrent value of the solution
[in]bdy_data_xloboundary data on interior of low x-face
[in]bdy_data_xhiboundary data on interior of high x-face
[in]bdy_data_yloboundary data on interior of low y-face
[in]bdy_data_yhiboundary data on interior of high y-face
111 {
112  BL_PROFILE_REGION("realbdy_compute_interior_ghost_RHS()");
113 
114  //
115  // Note that time (= start_time+old_stage_time) is measured as total time
116  // start_bdy_time and final_bdy_time are also measured as total time
117  //
118 
119  // Relaxation constants
120  Real F1 = 1./(nudge_factor*delta_t);
121 
122  // Time interpolation
123  Real dT = bdy_time_interval;
124 
125  int n_time = static_cast<int>( (time-start_bdy_time) / dT);
126  int n_time_p1 = n_time + 1;
127  Real alpha = ((time-start_bdy_time) - n_time * dT) / dT;
128 
129  // Do not over run the last bdy file
130  if (time >= final_bdy_time) {
131  n_time = static_cast<int>( (final_bdy_time - start_bdy_time)/ dT);
132  n_time_p1 = n_time;
133  alpha = 0.0;
134  }
135 
136  AMREX_ALWAYS_ASSERT( alpha >= 0. && alpha <= 1.0);
137  Real oma = 1.0 - alpha;
138 
139  /*
140  // UNIT TEST DEBUG
141  oma = 1.0; alpha = 0.0;
142  */
143 
144  // Temporary FABs for storage (owned/filled on all ranks)
145  FArrayBox U_xlo, U_xhi, U_ylo, U_yhi;
146  FArrayBox V_xlo, V_xhi, V_ylo, V_yhi;
147  FArrayBox T_xlo, T_xhi, T_ylo, T_yhi;
148 
149  // Variable index map (WRFBdyVars -> Vars)
150  Vector<int> var_map = {Vars::xvel, Vars::yvel, Vars::cons, Vars::cons};
151  Vector<int> ivar_map = {IntVars::xmom, IntVars::ymom, IntVars::cons, IntVars::cons};
152 
153  // Variable icomp map
154  Vector<int> comp_map = {0, 0, RhoTheta_comp};
155 
156  // Indices
157  int ivarU = RealBdyVars::U;
158  int ivarV = RealBdyVars::V;
159  int ivarT = RealBdyVars::T;
160  int BdyEnd = RealBdyVars::NumTypes-1;
161 
162 
163  // NOTE: The sizing of the temporary BDY FABS is
164  // GLOBAL and occurs over the entire BDY region.
165 
166  // Size the FABs
167  //==========================================================
168  for (int ivar(ivarU); ivar < BdyEnd; ivar++) {
169  int ivar_idx = var_map[ivar];
170  Box domain = geom.Domain();
171  auto ixtype = S_cur_data[ivar_idx].boxArray().ixType();
172  domain.convert(ixtype);
173 
174  // NOTE: Ghost cells needed for idx type mismatch between mask and data (do_upwind)
175  IntVect ng_vect(1,1,0);
176  Box gdom(domain); gdom.grow(ng_vect);
177  Box bx_xlo, bx_xhi, bx_ylo, bx_yhi;
178  realbdy_interior_bxs_xy(gdom, domain, width,
179  bx_xlo, bx_xhi,
180  bx_ylo, bx_yhi,
181  ng_vect, true);
182 
183  // Size the FABs
184  if (ivar == ivarU) {
185  U_xlo.resize(bx_xlo,1,The_Async_Arena()); U_xhi.resize(bx_xhi,1,The_Async_Arena());
186  U_ylo.resize(bx_ylo,1,The_Async_Arena()); U_yhi.resize(bx_yhi,1,The_Async_Arena());
187  } else if (ivar == ivarV) {
188  V_xlo.resize(bx_xlo,1,The_Async_Arena()); V_xhi.resize(bx_xhi,1,The_Async_Arena());
189  V_ylo.resize(bx_ylo,1,The_Async_Arena()); V_yhi.resize(bx_yhi,1,The_Async_Arena());
190  } else if (ivar == ivarT){
191  T_xlo.resize(bx_xlo,1,The_Async_Arena()); T_xhi.resize(bx_xhi,1,The_Async_Arena());
192  T_ylo.resize(bx_ylo,1,The_Async_Arena()); T_yhi.resize(bx_yhi,1,The_Async_Arena());
193  } else {
194  continue;
195  }
196  } // ivar
197 
198 
199  // NOTE: These operations use the BDY FABS and RHO. The
200  // use of RHO to go from PRIM -> CONS requires that
201  // these operations be LOCAL. So we have allocated
202  // enough space to do global operations (1 rank) but
203  // will fill a subset of that data that the rank owns.
204 
205  // Populate FABs from bdy interpolation (primitive vars)
206  //==========================================================
207  for (int ivar(ivarU); ivar < BdyEnd; ivar++) {
208  int ivar_idx = var_map[ivar];
209  Box domain = geom.Domain();
210  auto ixtype = S_cur_data[ivar_idx].boxArray().ixType();
211  domain.convert(ixtype);
212  const auto& dom_lo = lbound(domain);
213  const auto& dom_hi = ubound(domain);
214 
215 #ifdef _OPENMP
216 #pragma omp parallel if (Gpu::notInLaunchRegion())
217 #endif
218  for (MFIter mfi(S_cur_data[ivar_idx],TilingIfNotGPU()); mfi.isValid(); ++mfi) {
219  // NOTE: Ghost cells needed for idx type mismatch between mask and data (do_upwind)
220  IntVect ng_vect(1,1,0);
221  Box gtbx = grow(mfi.tilebox(ixtype.toIntVect()),ng_vect);
222  Box tbx_xlo, tbx_xhi, tbx_ylo, tbx_yhi;
223  realbdy_interior_bxs_xy(gtbx, domain, width,
224  tbx_xlo, tbx_xhi,
225  tbx_ylo, tbx_yhi,
226  ng_vect, true);
227 
228  Array4<Real> arr_xlo; Array4<Real> arr_xhi;
229  Array4<Real> arr_ylo; Array4<Real> arr_yhi;
230  if (ivar == ivarU) {
231  arr_xlo = U_xlo.array(); arr_xhi = U_xhi.array();
232  arr_ylo = U_ylo.array(); arr_yhi = U_yhi.array();
233  } else if (ivar == ivarV) {
234  arr_xlo = V_xlo.array(); arr_xhi = V_xhi.array();
235  arr_ylo = V_ylo.array(); arr_yhi = V_yhi.array();
236  } else if (ivar == ivarT){
237  arr_xlo = T_xlo.array(); arr_xhi = T_xhi.array();
238  arr_ylo = T_ylo.array(); arr_yhi = T_yhi.array();
239  } else {
240  continue;
241  }
242 
243  // Boundary data at fixed time intervals
244  const auto& bdatxlo_n = bdy_data_xlo[n_time ][ivar].const_array();
245  const auto& bdatxlo_np1 = bdy_data_xlo[n_time_p1][ivar].const_array();
246  const auto& bdatxhi_n = bdy_data_xhi[n_time ][ivar].const_array();
247  const auto& bdatxhi_np1 = bdy_data_xhi[n_time_p1][ivar].const_array();
248  const auto& bdatylo_n = bdy_data_ylo[n_time ][ivar].const_array();
249  const auto& bdatylo_np1 = bdy_data_ylo[n_time_p1][ivar].const_array();
250  const auto& bdatyhi_n = bdy_data_yhi[n_time ][ivar].const_array();
251  const auto& bdatyhi_np1 = bdy_data_yhi[n_time_p1][ivar].const_array();
252 
253  // Current density to convert to conserved vars
254  Array4<Real> r_arr = S_cur_data[IntVars::cons].array(mfi);
255 
256  // Limiting offset
257  int offset = width - 1;
258 
259  // Populate with interpolation (protect from ghost cells)
260  ParallelFor(tbx_xlo, tbx_xhi,
261  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
262  {
263  int ii = std::max(i , dom_lo.x);
264  ii = std::min(ii, dom_lo.x+offset);
265  int jj = std::max(j , dom_lo.y);
266  jj = std::min(jj, dom_hi.y);
267 
268  Real rho_interp;
269  if (ivar==ivarU) {
270  rho_interp = 0.5 * ( r_arr(i-1,j ,k) + r_arr(i,j,k) );
271  } else if (ivar==ivarV) {
272  rho_interp = 0.5 * ( r_arr(i ,j-1,k) + r_arr(i,j,k) );
273  } else {
274  rho_interp = r_arr(i,j,k);
275  }
276  arr_xlo(i,j,k) = rho_interp * ( oma * bdatxlo_n (ii,jj,k,0)
277  + alpha * bdatxlo_np1(ii,jj,k,0) );
278  },
279  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
280  {
281  int ii = std::max(i , dom_hi.x-offset);
282  ii = std::min(ii, dom_hi.x);
283  int jj = std::max(j , dom_lo.y);
284  jj = std::min(jj, dom_hi.y);
285 
286  Real rho_interp;
287  if (ivar==ivarU) {
288  rho_interp = 0.5 * ( r_arr(i-1,j ,k) + r_arr(i,j,k) );
289  } else if (ivar==ivarV) {
290  rho_interp = 0.5 * ( r_arr(i ,j-1,k) + r_arr(i,j,k) );
291  } else {
292  rho_interp = r_arr(i,j,k);
293  }
294  arr_xhi(i,j,k) = rho_interp * ( oma * bdatxhi_n (ii,jj,k,0)
295  + alpha * bdatxhi_np1(ii,jj,k,0) );
296  });
297 
298  ParallelFor(tbx_ylo, tbx_yhi,
299  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
300  {
301  int ii = std::max(i , dom_lo.x);
302  ii = std::min(ii, dom_hi.x);
303  int jj = std::max(j , dom_lo.y);
304  jj = std::min(jj, dom_lo.y+offset);
305 
306  Real rho_interp;
307  if (ivar==ivarU) {
308  rho_interp = 0.5 * ( r_arr(i-1,j ,k) + r_arr(i,j,k) );
309  } else if (ivar==ivarV) {
310  rho_interp = 0.5 * ( r_arr(i ,j-1,k) + r_arr(i,j,k) );
311  } else {
312  rho_interp = r_arr(i,j,k);
313  }
314  arr_ylo(i,j,k) = rho_interp * ( oma * bdatylo_n (ii,jj,k,0)
315  + alpha * bdatylo_np1(ii,jj,k,0) );
316  },
317  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
318  {
319  int ii = std::max(i , dom_lo.x);
320  ii = std::min(ii, dom_hi.x);
321  int jj = std::max(j , dom_hi.y-offset);
322  jj = std::min(jj, dom_hi.y);
323 
324  Real rho_interp;
325  if (ivar==ivarU) {
326  rho_interp = 0.5 * ( r_arr(i-1,j ,k) + r_arr(i,j,k) );
327  } else if (ivar==ivarV) {
328  rho_interp = 0.5 * ( r_arr(i ,j-1,k) + r_arr(i,j,k) );
329  } else {
330  rho_interp = r_arr(i,j,k);
331  }
332  arr_yhi(i,j,k) = rho_interp * ( oma * bdatyhi_n (ii,jj,k,0)
333  + alpha * bdatyhi_np1(ii,jj,k,0) );
334  });
335  } // mfi
336  } // ivar
337 
338 
339  // Compute RHS in relaxation region
340  //==========================================================
341  auto dx = geom.CellSizeArray();
342  auto ProbLo = geom.ProbLoArray();
343  auto ProbHi = geom.ProbHiArray();
344  for (int ivar(ivarU); ivar < BdyEnd; ivar++) {
345  int ivar_idx = ivar_map[ivar];
346  int icomp = comp_map[ivar];
347 
348  Box domain = geom.Domain();
349  domain.convert(S_cur_data[ivar_idx].boxArray().ixType());
350  IntVect ng_vect(0);
351 
352 #ifdef _OPENMP
353 #pragma omp parallel if (Gpu::notInLaunchRegion())
354 #endif
355  for (MFIter mfi(S_cur_data[ivar_idx],TilingIfNotGPU()); mfi.isValid(); ++mfi) {
356  Box tbx = mfi.tilebox();
357  Box tbx_xlo, tbx_xhi, tbx_ylo, tbx_yhi;
358  realbdy_interior_bxs_xy(tbx, domain, width,
359  tbx_xlo, tbx_xhi,
360  tbx_ylo, tbx_yhi,
361  ng_vect);
362 
363  Array4<Real> rhs_arr; Array4<Real> data_arr;
364  Array4<Real> arr_xlo; Array4<Real> arr_xhi;
365  Array4<Real> arr_ylo; Array4<Real> arr_yhi;
366  if (ivar == ivarU) {
367  arr_xlo = U_xlo.array(); arr_xhi = U_xhi.array();
368  arr_ylo = U_ylo.array(); arr_yhi = U_yhi.array();
369  rhs_arr = S_rhs[IntVars::xmom].array(mfi);
370  data_arr = S_cur_data[IntVars::xmom].array(mfi);
371  } else if (ivar == ivarV) {
372  arr_xlo = V_xlo.array(); arr_xhi = V_xhi.array();
373  arr_ylo = V_ylo.array(); arr_yhi = V_yhi.array();
374  rhs_arr = S_rhs[IntVars::ymom].array(mfi);
375  data_arr = S_cur_data[IntVars::ymom].array(mfi);
376  } else if (ivar == ivarT){
377  arr_xlo = T_xlo.array(); arr_xhi = T_xhi.array();
378  arr_ylo = T_ylo.array(); arr_yhi = T_yhi.array();
379  rhs_arr = S_rhs[IntVars::cons].array(mfi);
380  data_arr = S_cur_data[IntVars::cons].array(mfi);
381  } else {
382  continue;
383  }
384 
385  Array4<Real> u_xlo = U_xlo.array(); Array4<Real> u_xhi = U_xhi.array();
386  Array4<Real> v_xlo = V_xlo.array(); Array4<Real> v_xhi = V_xhi.array();
387  Array4<Real> v_ylo = V_ylo.array(); Array4<Real> v_yhi = V_yhi.array();
388 
390  width, dx, ProbLo, ProbHi, F1, geom.Domain(),
391  tbx_xlo , tbx_xhi , tbx_ylo , tbx_yhi ,
392  arr_xlo , arr_xhi , arr_ylo , arr_yhi ,
393  u_xlo, u_xhi, v_xlo, v_xhi, v_ylo, v_yhi,
394  data_arr, rhs_arr, do_upwind);
395 
396  /*
397  // UNIT TEST DEBUG
398  realbdy_interior_bxs_xy(tbx, domain, width,
399  tbx_xlo, tbx_xhi,
400  tbx_ylo, tbx_yhi);
401  ParallelFor(tbx_xlo, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
402  {
403  if (std::fabs(arr_xlo(i,j,k) - data_arr(i,j,k,icomp)) > 0.01) {
404  Print() << "ERROR XLO: " << ivar << ' ' << icomp << ' ' << IntVect(i,j,k) << "\n";
405  Print() << "DATA: " << data_arr(i,j,k,icomp) << ' ' << arr_xlo(i,j,k) << "\n";
406  exit(0);
407  }
408  });
409  ParallelFor(tbx_xhi, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
410  {
411  if (std::fabs(arr_xhi(i,j,k) - data_arr(i,j,k,icomp)) > 0.01) {
412  Print() << "ERROR XHI: " << ivar << ' ' << icomp << ' ' << IntVect(i,j,k) << "\n";
413  Print() << "DATA: " << data_arr(i,j,k,icomp) << ' ' << arr_xhi(i,j,k) << "\n";
414  exit(0);
415  }
416  });
417  ParallelFor(tbx_ylo, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
418  {
419  if (std::fabs(arr_ylo(i,j,k) - data_arr(i,j,k,icomp)) > 0.01) {
420  Print() << "ERROR YLO: " << ivar << ' ' << icomp << ' ' << IntVect(i,j,k) << "\n";
421  Print() << "DATA: " << data_arr(i,j,k,icomp) << ' ' << arr_ylo(i,j,k) << "\n";
422  exit(0);
423  }
424  });
425  ParallelFor(tbx_yhi, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
426  {
427  if (std::fabs(arr_yhi(i,j,k)-data_arr(i,j,k,icomp)) > 0.01) {
428  Print() << "ERROR YHI: " << ivar << ' ' << icomp << ' ' << IntVect(i,j,k) << "\n";
429  Print() << "DATA: " << data_arr(i,j,k,icomp) << ' ' << arr_yhi(i,j,k) << "\n";
430  exit(0);
431  }
432  });
433  */
434  } // mfi
435  } // ivar
436  //ParallelDescriptor::Barrier();
437  //exit(0);
438 }
#define RhoTheta_comp
Definition: ERF_IndexDefines.H:37
const auto dx
Definition: ERF_InitCustomPertVels_ParticleTests.H:15
amrex::Real alpha
Definition: ERF_InitCustomPert_IsentropicVortex.H:8
AMREX_ALWAYS_ASSERT(bx.length()[2]==khi+1)
void realbdy_interior_bxs_xy(const Box &bx, const Box &domain, const int &width, Box &bx_xlo, Box &bx_xhi, Box &bx_ylo, Box &bx_yhi, const IntVect &ng_vect, const bool get_int_ng)
Definition: ERF_InteriorGhostCells.cpp:23
AMREX_FORCE_INLINE IntVect offset(const int face_dir, const int normal)
Definition: ERF_ReadBndryPlanes.cpp:28
const auto & dom_hi
Definition: ERF_SetupVertDiff.H:2
const auto & dom_lo
Definition: ERF_SetupVertDiff.H:1
AMREX_GPU_HOST AMREX_FORCE_INLINE void realbdy_compute_relaxation(const int &icomp, const int &num_var, const int &width, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &dx, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &ProbLo, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &ProbHi, const amrex::Real &F1, const amrex::Box &domain_cc, 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 > &u_xlo, const amrex::Array4< const amrex::Real > &u_xhi, const amrex::Array4< const amrex::Real > &v_xlo, const amrex::Array4< const amrex::Real > &v_xhi, const amrex::Array4< const amrex::Real > &v_ylo, const amrex::Array4< const amrex::Real > &v_yhi, const amrex::Array4< const amrex::Real > &data_arr, const amrex::Array4< amrex::Real > &rhs_arr, const bool &do_upwind)
Definition: ERF_Utils.H:180
@ U
Definition: ERF_IndexDefines.H:108
@ NumTypes
Definition: ERF_IndexDefines.H:112
@ T
Definition: ERF_IndexDefines.H:110
@ V
Definition: ERF_IndexDefines.H:109
@ xvel
Definition: ERF_IndexDefines.H:141
@ cons
Definition: ERF_IndexDefines.H:140
@ yvel
Definition: ERF_IndexDefines.H:142
Here is the call graph for this function:

◆ realbdy_interior_bxs_xy()

void realbdy_interior_bxs_xy ( const Box &  bx,
const Box &  domain,
const int &  width,
Box &  bx_xlo,
Box &  bx_xhi,
Box &  bx_ylo,
Box &  bx_yhi,
const IntVect &  ng_vect,
const bool  get_int_ng 
)

Get the boxes for looping over interior/exterior ghost cells for use by fillpatch, erf_slow_rhs_pre, and erf_slow_rhs_post.

Parameters
[in]bxbox to intersect with 4 halo regions
[in]domainbox of the whole domain
[in]widthnumber of cells in (relaxation+specified) zone
[in]set_widthnumber of cells in (specified) zone
[out]bx_xlohalo box at x_lo boundary
[out]bx_xhihalo box at x_hi boundary
[out]bx_ylohalo box at y_lo boundary
[out]bx_yhihalo box at y_hi boundary
[in]ng_vectnumber of ghost cells in each direction
[in]get_int_ngflag to get ghost cells inside the domain
32 {
33  AMREX_ALWAYS_ASSERT(bx.ixType() == domain.ixType());
34 
35  //==================================================================
36  // NOTE: X-face boxes take ownership of the overlapping region.
37  // With exterior ghost cells (ng_vect != 0), the x-face
38  // boxes will have exterior ghost cells in both x & y.
39  //==================================================================
40 
41  // Domain bounds without ghost cells
42  const auto& dom_lo = lbound(domain);
43  const auto& dom_hi = ubound(domain);
44 
45  // Four boxes matching the domain
46  Box gdom_xlo(domain); Box gdom_xhi(domain);
47  Box gdom_ylo(domain); Box gdom_yhi(domain);
48 
49  // Trim the boxes to only include internal ghost cells
50  gdom_xlo.setBig(0,dom_lo.x+width-1); gdom_xhi.setSmall(0,dom_hi.x-width+1);
51  gdom_ylo.setBig(1,dom_lo.y+width-1); gdom_yhi.setSmall(1,dom_hi.y-width+1);
52 
53  // Remove overlapping corners from y-face boxes
54  gdom_ylo.setSmall(0,gdom_xlo.bigEnd(0)+1); gdom_ylo.setBig(0,gdom_xhi.smallEnd(0)-1);
55  gdom_yhi.setSmall(0,gdom_xlo.bigEnd(0)+1); gdom_yhi.setBig(0,gdom_xhi.smallEnd(0)-1);
56 
57  // Grow boxes to get external ghost cells only
58  gdom_xlo.growLo(0,ng_vect[0]); gdom_xhi.growHi(0,ng_vect[0]);
59  gdom_xlo.grow (1,ng_vect[1]); gdom_xhi.grow (1,ng_vect[1]);
60  gdom_ylo.growLo(1,ng_vect[1]); gdom_yhi.growHi(1,ng_vect[1]);
61 
62  // Grow boxes to get internal ghost cells
63  if (get_int_ng) {
64  gdom_xlo.growHi(0,ng_vect[0]); gdom_xhi.growLo(0,ng_vect[0]);
65  gdom_ylo.grow (0,ng_vect[0]); gdom_yhi.grow (0,ng_vect[0]);
66  gdom_ylo.growHi(1,ng_vect[1]); gdom_yhi.growLo(1,ng_vect[1]);
67  }
68 
69  // Populate everything
70  bx_xlo = (bx & gdom_xlo);
71  bx_xhi = (bx & gdom_xhi);
72  bx_ylo = (bx & gdom_ylo);
73  bx_yhi = (bx & gdom_yhi);
74 }

Referenced by realbdy_compute_interior_ghost_rhs().

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

Variable Documentation

◆ void_bc

PhysBCFunctNoOp void_bc