ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
Time_Avg_Vel.cpp File Reference
#include <Utils.H>
Include dependency graph for Time_Avg_Vel.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;
42  vel_t_avg_arr(i,j,k,1) += v_cc;
43  vel_t_avg_arr(i,j,k,2) += w_cc;
44  vel_t_avg_arr(i,j,k,3) += umag_cc;
45  });
46  }
47 }
@ xvel
Definition: IndexDefines.H:100
@ zvel
Definition: IndexDefines.H:102
@ yvel
Definition: IndexDefines.H:101

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

Here is the caller graph for this function: