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

Functions

void ComputeStressConsVisc_N (Box bxcc, Box tbxxy, Box tbxxz, Box tbxyz, Real mu_eff, const Array4< const Real > &cell_data, Array4< Real > &tau11, Array4< Real > &tau22, Array4< Real > &tau33, Array4< Real > &tau12, Array4< Real > &tau13, Array4< Real > &tau23, const Array4< const Real > &er_arr)
 
void ComputeStressVarVisc_N (Box bxcc, Box tbxxy, Box tbxxz, Box tbxyz, Real mu_eff, const Array4< const Real > &mu_turb, const Array4< const Real > &cell_data, Array4< Real > &tau11, Array4< Real > &tau22, Array4< Real > &tau33, Array4< Real > &tau12, Array4< Real > &tau13, Array4< Real > &tau23, const Array4< const Real > &er_arr)
 

Function Documentation

◆ ComputeStressConsVisc_N()

void ComputeStressConsVisc_N ( Box  bxcc,
Box  tbxxy,
Box  tbxxz,
Box  tbxyz,
Real  mu_eff,
const Array4< const Real > &  cell_data,
Array4< Real > &  tau11,
Array4< Real > &  tau22,
Array4< Real > &  tau33,
Array4< Real > &  tau12,
Array4< Real > &  tau13,
Array4< Real > &  tau23,
const Array4< const Real > &  er_arr 
)

Function for computing the stress with constant viscosity and without terrain.

Parameters
[in]bxcccell center box for tau_ii
[in]tbxxynodal xy box for tau_12
[in]tbxxznodal xz box for tau_13
[in]tbxyznodal yz box for tau_23
[in]mu_effconstant molecular viscosity
[in]cell_datato access rho if ConstantAlpha
[in,out]tau1111 strain -> stress
[in,out]tau2222 strain -> stress
[in,out]tau3333 strain -> stress
[in,out]tau1212 strain -> stress
[in,out]tau1313 strain -> stress
[in,out]tau2323 strain -> stress
[in]er_arrexpansion rate
28 {
29  Real OneThird = (1./3.);
30 
31  if (cell_data)
32  // constant alpha (stored in mu_eff)
33  {
34  // Cell centered strains
35  ParallelFor(bxcc, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
36  {
37  Real rhoAlpha = cell_data(i, j, k, Rho_comp) * mu_eff;
38  tau11(i,j,k) = -rhoAlpha * ( tau11(i,j,k) - OneThird*er_arr(i,j,k) );
39  tau22(i,j,k) = -rhoAlpha * ( tau22(i,j,k) - OneThird*er_arr(i,j,k) );
40  tau33(i,j,k) = -rhoAlpha * ( tau33(i,j,k) - OneThird*er_arr(i,j,k) );
41  });
42 
43  // Off-diagonal strains
44  ParallelFor(tbxxy,tbxxz,tbxyz,
45  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
46  Real rho_bar = 0.25*( cell_data(i-1, j , k, Rho_comp) + cell_data(i, j , k, Rho_comp)
47  + cell_data(i-1, j-1, k, Rho_comp) + cell_data(i, j-1, k, Rho_comp) );
48  tau12(i,j,k) *= -rho_bar * mu_eff;
49  },
50  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
51  Real rho_bar = 0.25*( cell_data(i-1, j, k , Rho_comp) + cell_data(i, j, k , Rho_comp)
52  + cell_data(i-1, j, k-1, Rho_comp) + cell_data(i, j, k-1, Rho_comp) );
53  tau13(i,j,k) *= -rho_bar * mu_eff;
54  },
55  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
56  Real rho_bar = 0.25*( cell_data(i, j-1, k , Rho_comp) + cell_data(i, j, k , Rho_comp)
57  + cell_data(i, j-1, k-1, Rho_comp) + cell_data(i, j, k-1, Rho_comp) );
58  tau23(i,j,k) *= -rho_bar * mu_eff;
59  });
60  }
61  else
62  // constant mu_eff
63  {
64  // Cell centered strains
65  ParallelFor(bxcc, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
66  {
67  tau11(i,j,k) = -mu_eff * ( tau11(i,j,k) - OneThird*er_arr(i,j,k) );
68  tau22(i,j,k) = -mu_eff * ( tau22(i,j,k) - OneThird*er_arr(i,j,k) );
69  tau33(i,j,k) = -mu_eff * ( tau33(i,j,k) - OneThird*er_arr(i,j,k) );
70  });
71 
72  // Off-diagonal strains
73  ParallelFor(tbxxy,tbxxz,tbxyz,
74  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
75  tau12(i,j,k) *= -mu_eff;
76  },
77  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
78  tau13(i,j,k) *= -mu_eff;
79  },
80  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
81  tau23(i,j,k) *= -mu_eff;
82  });
83  }
84 }
#define Rho_comp
Definition: ERF_IndexDefines.H:36

Referenced by erf_make_tau_terms().

Here is the caller graph for this function:

◆ ComputeStressVarVisc_N()

void ComputeStressVarVisc_N ( Box  bxcc,
Box  tbxxy,
Box  tbxxz,
Box  tbxyz,
Real  mu_eff,
const Array4< const Real > &  mu_turb,
const Array4< const Real > &  cell_data,
Array4< Real > &  tau11,
Array4< Real > &  tau22,
Array4< Real > &  tau33,
Array4< Real > &  tau12,
Array4< Real > &  tau13,
Array4< Real > &  tau23,
const Array4< const Real > &  er_arr 
)

Function for computing the stress with variable viscosity and without terrain.

Parameters
[in]bxcccell center box for tau_ii
[in]tbxxynodal xy box for tau_12
[in]tbxxznodal xz box for tau_13
[in]tbxyznodal yz box for tau_23
[in]mu_effconstant molecular viscosity
[in]mu_turbvariable turbulent viscosity
[in]cell_datato access rho if ConstantAlpha
[in,out]tau1111 strain -> stress
[in,out]tau2222 strain -> stress
[in,out]tau3333 strain -> stress
[in,out]tau1212 strain -> stress
[in,out]tau1313 strain -> stress
[in,out]tau2323 strain -> stress
[in]er_arrexpansion rate
111 {
112  Real OneThird = (1./3.);
113 
114  if (cell_data)
115  // constant alpha (stored in mu_eff)
116  {
117  // Cell centered strains
118  ParallelFor(bxcc, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
119  Real rhoAlpha = cell_data(i, j, k, Rho_comp) * mu_eff;
120  Real mu_11 = rhoAlpha + 2.0 * mu_turb(i, j, k, EddyDiff::Mom_h);
121  Real mu_22 = mu_11;
122  Real mu_33 = rhoAlpha + 2.0 * mu_turb(i, j, k, EddyDiff::Mom_v);
123  tau11(i,j,k) = -mu_11 * ( tau11(i,j,k) - OneThird*er_arr(i,j,k) );
124  tau22(i,j,k) = -mu_22 * ( tau22(i,j,k) - OneThird*er_arr(i,j,k) );
125  tau33(i,j,k) = -mu_33 * ( tau33(i,j,k) - OneThird*er_arr(i,j,k) );
126  });
127 
128  // Off-diagonal strains
129  ParallelFor(tbxxy,tbxxz,tbxyz,
130  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
131  Real rho_bar = 0.25*( cell_data(i-1, j , k, Rho_comp) + cell_data(i, j , k, Rho_comp)
132  + cell_data(i-1, j-1, k, Rho_comp) + cell_data(i, j-1, k, Rho_comp) );
133  Real mu_bar = 0.25*( mu_turb(i-1, j , k, EddyDiff::Mom_h) + mu_turb(i, j , k, EddyDiff::Mom_h)
134  + mu_turb(i-1, j-1, k, EddyDiff::Mom_h) + mu_turb(i, j-1, k, EddyDiff::Mom_h) );
135  Real mu_12 = rho_bar*mu_eff + 2.0*mu_bar;
136  tau12(i,j,k) *= -mu_12;
137  },
138  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
139  Real rho_bar = 0.25*( cell_data(i-1, j, k , Rho_comp) + cell_data(i, j, k , Rho_comp)
140  + cell_data(i-1, j, k-1, Rho_comp) + cell_data(i, j, k-1, Rho_comp) );
141  Real mu_bar = 0.25*( mu_turb(i-1, j, k , EddyDiff::Mom_v) + mu_turb(i, j, k , EddyDiff::Mom_v)
142  + mu_turb(i-1, j, k-1, EddyDiff::Mom_v) + mu_turb(i, j, k-1, EddyDiff::Mom_v) );
143  Real mu_13 = rho_bar*mu_eff + 2.0*mu_bar;
144  tau13(i,j,k) *= -mu_13;
145  },
146  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
147  Real rho_bar = 0.25*( cell_data(i, j-1, k , Rho_comp) + cell_data(i, j, k , Rho_comp)
148  + cell_data(i, j-1, k-1, Rho_comp) + cell_data(i, j, k-1, Rho_comp) );
149  Real mu_bar = 0.25*( mu_turb(i, j-1, k , EddyDiff::Mom_v) + mu_turb(i, j, k , EddyDiff::Mom_v)
150  + mu_turb(i, j-1, k-1, EddyDiff::Mom_v) + mu_turb(i, j, k-1, EddyDiff::Mom_v) );
151  Real mu_23 = rho_bar*mu_eff + 2.0*mu_bar;
152  tau23(i,j,k) *= -mu_23;
153  });
154  }
155  else
156  // constant mu_eff
157  {
158  // Cell centered strains
159  ParallelFor(bxcc, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
160  Real mu_11 = mu_eff + 2.0 * mu_turb(i, j, k, EddyDiff::Mom_h);
161  Real mu_22 = mu_11;
162  Real mu_33 = mu_eff + 2.0 * mu_turb(i, j, k, EddyDiff::Mom_v);
163  tau11(i,j,k) = -mu_11 * ( tau11(i,j,k) - OneThird*er_arr(i,j,k) );
164  tau22(i,j,k) = -mu_22 * ( tau22(i,j,k) - OneThird*er_arr(i,j,k) );
165  tau33(i,j,k) = -mu_33 * ( tau33(i,j,k) - OneThird*er_arr(i,j,k) );
166  });
167 
168  // Off-diagonal strains
169  ParallelFor(tbxxy,tbxxz,tbxyz,
170  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
171  Real mu_bar = 0.25*( mu_turb(i-1, j , k, EddyDiff::Mom_h) + mu_turb(i, j , k, EddyDiff::Mom_h)
172  + mu_turb(i-1, j-1, k, EddyDiff::Mom_h) + mu_turb(i, j-1, k, EddyDiff::Mom_h) );
173  Real mu_12 = mu_eff + 2.0*mu_bar;
174  tau12(i,j,k) *= -mu_12;
175  },
176  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
177  Real mu_bar = 0.25*( mu_turb(i-1, j, k , EddyDiff::Mom_v) + mu_turb(i, j, k , EddyDiff::Mom_v)
178  + mu_turb(i-1, j, k-1, EddyDiff::Mom_v) + mu_turb(i, j, k-1, EddyDiff::Mom_v) );
179  Real mu_13 = mu_eff + 2.0*mu_bar;
180  tau13(i,j,k) *= -mu_13;
181  },
182  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
183  Real mu_bar = 0.25*( mu_turb(i, j-1, k , EddyDiff::Mom_v) + mu_turb(i, j, k , EddyDiff::Mom_v)
184  + mu_turb(i, j-1, k-1, EddyDiff::Mom_v) + mu_turb(i, j, k-1, EddyDiff::Mom_v) );
185  Real mu_23 = mu_eff + 2.0*mu_bar;
186  tau23(i,j,k) *= -mu_23;
187  });
188  }
189 }
@ Mom_h
Definition: ERF_IndexDefines.H:151
@ Mom_v
Definition: ERF_IndexDefines.H:156

Referenced by erf_make_tau_terms().

Here is the caller graph for this function: