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

Functions

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

Function Documentation

◆ Time_Avg_Vel_atCC()

void Time_Avg_Vel_atCC ( const Real dt,
Real t_avg_cnt,
MultiFab *  vel_t_avg,
MultiFab &  xvel,
MultiFab &  yvel,
MultiFab &  zvel 
)
15 {
16  // Augment the counter
17  t_avg_cnt += dt;
18 
19 #ifdef _OPENMP
20 #pragma omp parallel if (Gpu::notInLaunchRegion())
21 #endif
22  for ( MFIter mfi(*(vel_t_avg),TilingIfNotGPU()); mfi.isValid(); ++mfi)
23  {
24  // CC tilebox
25  Box tbx = mfi.tilebox();
26 
27  // Velocity on faces
28  const Array4<Real>& velx = xvel.array(mfi);
29  const Array4<Real>& vely = yvel.array(mfi);
30  const Array4<Real>& velz = zvel.array(mfi);
31 
32  // Time average at CC
33  Array4<Real> vel_t_avg_arr = vel_t_avg->array(mfi);
34 
35  ParallelFor(tbx, [=] AMREX_GPU_DEVICE (int i, int j, int k)
36  {
37  Real u_cc = 0.5 * ( velx(i,j,k) + velx(i+1,j ,k ) );
38  Real v_cc = 0.5 * ( vely(i,j,k) + vely(i ,j+1,k ) );
39  Real w_cc = 0.5 * ( velz(i,j,k) + velz(i ,j ,k+1) );
40  Real umag_cc = std::sqrt(u_cc*u_cc + v_cc*v_cc + w_cc*w_cc);
41  vel_t_avg_arr(i,j,k,0) += u_cc * dt;
42  vel_t_avg_arr(i,j,k,1) += v_cc * dt;
43  vel_t_avg_arr(i,j,k,2) += w_cc * dt;
44  vel_t_avg_arr(i,j,k,3) += umag_cc * dt;
45  });
46  }
47 }
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);} } })
amrex::Real Real
Definition: ERF_ShocInterface.H:19
@ xvel
Definition: ERF_IndexDefines.H:141
@ zvel
Definition: ERF_IndexDefines.H:143
@ yvel
Definition: ERF_IndexDefines.H:142

Referenced by ERF::Advance(), and ERF::InitData_post().

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