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

Functions

void tag_on_distance_from_eye (const Geometry &cgeom, TagBoxArray *tags, const Real eye_x, const Real eye_y, const Real rad_tag)
 

Function Documentation

◆ tag_on_distance_from_eye()

void tag_on_distance_from_eye ( const Geometry &  cgeom,
TagBoxArray *  tags,
const Real  eye_x,
const Real  eye_y,
const Real  rad_tag 
)

Function to tag cells for refinement – this overrides the pure virtual function in AmrCore

Parameters
[in]levclevel of refinement at which we tag cells (0 is coarsest level)
[out]tagsarray of tagged cells
[in]timecurrent time
745 {
746  const auto dx = cgeom.CellSizeArray();
747  const auto prob_lo = cgeom.ProbLoArray();
748 
749  for (MFIter mfi(*tags); mfi.isValid(); ++mfi) {
750  TagBox& tag = (*tags)[mfi];
751  auto tag_arr = tag.array(); // Get device-accessible array
752 
753  const Box& tile_box = mfi.tilebox(); // The box for this tile
754 
755  ParallelFor(tile_box, [=] AMREX_GPU_DEVICE(int i, int j, int k) {
756  // Compute cell center coordinates
757  Real x = prob_lo[0] + (i + 0.5) * dx[0];
758  Real y = prob_lo[1] + (j + 0.5) * dx[1];
759 
760  Real dist = std::sqrt((x - eye_x)*(x - eye_x) + (y - eye_y)*(y - eye_y));
761 
762  if (dist < rad_tag) {
763  tag_arr(i,j,k) = TagBox::SET;
764  } else {
765  tag_arr(i,j,k) = TagBox::CLEAR;
766  }
767  });
768  }
769 }
const auto dx
Definition: ERF_InitCustomPertVels_ParticleTests.H:15
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);} } })
const amrex::Real * prob_lo
Definition: ERF_InitCustomPert_IsentropicVortex.H:16
amrex::Real Real
Definition: ERF_ShocInterface.H:19

Referenced by ERF::ErrorEst(), and ERF::HurricaneTracker().

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