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
621 {
622  const auto dx = cgeom.CellSizeArray();
623  const auto prob_lo = cgeom.ProbLoArray();
624 
625  for (MFIter mfi(*tags); mfi.isValid(); ++mfi) {
626  TagBox& tag = (*tags)[mfi];
627  auto tag_arr = tag.array(); // Get device-accessible array
628 
629  const Box& tile_box = mfi.tilebox(); // The box for this tile
630 
631  ParallelFor(tile_box, [=] AMREX_GPU_DEVICE(int i, int j, int k) {
632  // Compute cell center coordinates
633  Real x = prob_lo[0] + (i + 0.5) * dx[0];
634  Real y = prob_lo[1] + (j + 0.5) * dx[1];
635 
636  Real dist = std::sqrt((x - eye_x)*(x - eye_x) + (y - eye_y)*(y - eye_y));
637 
638  if (dist < rad_tag) {
639  tag_arr(i,j,k) = TagBox::SET;
640  } else {
641  tag_arr(i,j,k) = TagBox::CLEAR;
642  }
643  });
644  }
645 }
amrex::Real Real
Definition: ERF_ShocInterface.H:19

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

Here is the caller graph for this function: