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

Functions

Real compute_A (Real z)
 
Real compute_Aijk (Real z_k, Real z_kp1)
 
void fitch_advance (int lev, const Geometry &geom, const Real &dt_advance, MultiFab &cons_in, MultiFab &U_old, MultiFab &V_old, MultiFab &W_old, MultiFab &mf_vars_fitch, const amrex::MultiFab &mf_Nturb)
 
void fitch_update (const Real &dt_advance, MultiFab &cons_in, MultiFab &U_old, MultiFab &V_old, const MultiFab &mf_vars_fitch)
 
void fitch_source_terms_cellcentered (const Geometry &geom, const MultiFab &cons_in, const MultiFab &U_old, const MultiFab &V_old, const MultiFab &W_old, MultiFab &mf_vars_fitch, const amrex::MultiFab &mf_Nturb)
 

Variables

Real R = 30.0
 
Real z_c = 100.0
 
Real C_T = 0.8
 
Real C_TKE = 0.0
 

Function Documentation

◆ compute_A()

Real compute_A ( Real  z)
10 {
11 
12  Real d = std::min(std::fabs(z - z_c), R);
13  Real theta = std::acos(d/R);
14  Real A_s = R*R*theta - d*std::pow(R*R - d*d, 0.5);
15  Real A = M_PI*R*R/2.0 - A_s;
16 
17  return A;
18 }
Real z_c
Definition: AdvanceFitch.cpp:6
Real R
Definition: AdvanceFitch.cpp:5
@ theta
Definition: MM5.H:20

Referenced by compute_Aijk().

Here is the caller graph for this function:

◆ compute_Aijk()

Real compute_Aijk ( Real  z_k,
Real  z_kp1 
)
21 {
22 
23  Real A_k = compute_A(z_k);
24  Real A_kp1 = compute_A(z_kp1);
25 
26  Real check = (z_k - z_c)*(z_kp1 - z_c);
27  Real A_ijk;
28  if(check > 0){
29  A_ijk = std::fabs(A_k -A_kp1);
30  }
31  else{
32  A_ijk = A_k + A_kp1;
33  }
34 
35  return A_ijk;
36 }
Real compute_A(Real z)
Definition: AdvanceFitch.cpp:9

Referenced by fitch_source_terms_cellcentered().

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

◆ fitch_advance()

void fitch_advance ( int  lev,
const Geometry &  geom,
const Real &  dt_advance,
MultiFab &  cons_in,
MultiFab &  U_old,
MultiFab &  V_old,
MultiFab &  W_old,
MultiFab &  mf_vars_fitch,
const amrex::MultiFab &  mf_Nturb 
)
45 {
46  fitch_source_terms_cellcentered(geom, cons_in, U_old, V_old, W_old, mf_vars_fitch, mf_Nturb);
47  fitch_update(dt_advance, cons_in, U_old, V_old, mf_vars_fitch);
48 }
void fitch_source_terms_cellcentered(const Geometry &geom, const MultiFab &cons_in, const MultiFab &U_old, const MultiFab &V_old, const MultiFab &W_old, MultiFab &mf_vars_fitch, const amrex::MultiFab &mf_Nturb)
Definition: AdvanceFitch.cpp:84
void fitch_update(const Real &dt_advance, MultiFab &cons_in, MultiFab &U_old, MultiFab &V_old, const MultiFab &mf_vars_fitch)
Definition: AdvanceFitch.cpp:51

Referenced by ERF::Advance().

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

◆ fitch_source_terms_cellcentered()

void fitch_source_terms_cellcentered ( const Geometry &  geom,
const MultiFab &  cons_in,
const MultiFab &  U_old,
const MultiFab &  V_old,
const MultiFab &  W_old,
MultiFab &  mf_vars_fitch,
const amrex::MultiFab &  mf_Nturb 
)
88 {
89 
90  auto dx = geom.CellSizeArray();
91  auto ProbHiArr = geom.ProbHiArray();
92  auto ProbLoArr = geom.ProbLoArray();
93 
94  // Domain valid box
95  const amrex::Box& domain = geom.Domain();
96  int domlo_x = domain.smallEnd(0);
97  int domhi_x = domain.bigEnd(0) + 1;
98  int domlo_y = domain.smallEnd(1);
99  int domhi_y = domain.bigEnd(1) + 1;
100  int domlo_z = domain.smallEnd(2);
101  int domhi_z = domain.bigEnd(2) + 1;
102 
103 
104  Real sum = 0.0;
105  Real *sum_area = &sum;
106 
107  // The order of variables are - Vabs dVabsdt, dudt, dvdt, dTKEdt
108  mf_vars_fitch.setVal(0.0);
109 
110  for ( MFIter mfi(cons_in,TilingIfNotGPU()); mfi.isValid(); ++mfi) {
111 
112  const Box& bx = mfi.tilebox();
113  const Box& gbx = mfi.growntilebox(1);
114  auto cons_array = cons_in.array(mfi);
115  auto fitch_array = mf_vars_fitch.array(mfi);
116  auto Nturb_array = mf_Nturb.array(mfi);
117  auto u_vel = U_old.array(mfi);
118  auto v_vel = V_old.array(mfi);
119  auto w_vel = W_old.array(mfi);
120 
121  amrex::IntVect lo = bx.smallEnd();
122 
123  ParallelFor(gbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
124  int ii = amrex::min(amrex::max(i, domlo_x), domhi_x);
125  int jj = amrex::min(amrex::max(j, domlo_y), domhi_y);
126  int kk = amrex::min(amrex::max(k, domlo_z), domhi_z);
127 
128 
129  Real x = (ii+0.5) * dx[0];
130  Real y = (jj+0.5) * dx[1];
131  Real z = (kk+0.5) * dx[2];
132 
133  Real z_k = kk*dx[2];
134  Real z_kp1 = (kk+1)*dx[2];
135 
136  Real A_ijk = compute_Aijk(z_k, z_kp1);
137 
138  // Compute Fitch source terms
139 
140  Real Vabs = std::pow(u_vel(i,j,k)*u_vel(i,j,k) +
141  v_vel(i,j,k)*v_vel(i,j,k) +
142  w_vel(i,j,kk)*w_vel(i,j,kk), 0.5);
143 
144  fitch_array(i,j,k,0) = Vabs;
145  fitch_array(i,j,k,1) = -0.5*Nturb_array(i,j,k)/(dx[0]*dx[1])*C_T*Vabs*Vabs*A_ijk/(z_kp1 - z_k);
146  fitch_array(i,j,k,2) = u_vel(i,j,k)/Vabs*fitch_array(i,j,k,1);
147  fitch_array(i,j,k,3) = v_vel(i,j,k)/Vabs*fitch_array(i,j,k,1);
148  fitch_array(i,j,k,4) = 0.5*Nturb_array(i,j,k)/(dx[0]*dx[1])*C_TKE*std::pow(Vabs,3)*A_ijk/(z_kp1 - z_k);
149 
150  //amrex::Gpu::Atomic::Add(sum_area, A_ijk);
151  });
152  }
153  //std::cout << "Checking sum here...." <<"\n";
154  //printf("%0.15g, %0.15g\n", *sum_area , M_PI*R*R);
155  //exit(0);
156 }
Real C_T
Definition: AdvanceFitch.cpp:7
Real C_TKE
Definition: AdvanceFitch.cpp:7
Real compute_Aijk(Real z_k, Real z_kp1)
Definition: AdvanceFitch.cpp:20

Referenced by fitch_advance().

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

◆ fitch_update()

void fitch_update ( const Real &  dt_advance,
MultiFab &  cons_in,
MultiFab &  U_old,
MultiFab &  V_old,
const MultiFab &  mf_vars_fitch 
)
55 {
56 
57  for ( MFIter mfi(cons_in,TilingIfNotGPU()); mfi.isValid(); ++mfi) {
58 
59  Box bx = mfi.tilebox();
60  Box tbx = mfi.nodaltilebox(0);
61  Box tby = mfi.nodaltilebox(1);
62 
63  auto cons_array = cons_in.array(mfi);
64  auto fitch_array = mf_vars_fitch.array(mfi);
65  auto u_vel = U_old.array(mfi);
66  auto v_vel = V_old.array(mfi);
67 
68  ParallelFor(tbx, tby, bx,
69  [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
70  {
71  u_vel(i,j,k) = u_vel(i,j,k) + (fitch_array(i-1,j,k,2) + fitch_array(i,j,k,2))/2.0*dt_advance;
72  },
73  [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
74  {
75  v_vel(i,j,k) = v_vel(i,j,k) + (fitch_array(i,j-1,k,3) + fitch_array(i,j,k,3))/2.0*dt_advance;
76  },
77  [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept
78  {
79  cons_array(i,j,k,RhoQKE_comp) = cons_array(i,j,k,RhoQKE_comp) + fitch_array(i,j,k,4)*dt_advance;
80  });
81  }
82 }
#define RhoQKE_comp
Definition: IndexDefines.H:15

Referenced by fitch_advance().

Here is the caller graph for this function:

Variable Documentation

◆ C_T

Real C_T = 0.8

◆ C_TKE

Real C_TKE = 0.0

◆ R

Real R = 30.0

Referenced by compute_A().

◆ z_c

Real z_c = 100.0

Referenced by compute_A(), and compute_Aijk().