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, Array4< Real > &tau13i, Array4< Real > &tau23i, Array4< Real > &tau33i)
 
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, Array4< Real > &tau13i, Array4< Real > &tau23i, Array4< Real > &tau33i)
 

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,
Array4< Real > &  tau13i,
Array4< Real > &  tau23i,
Array4< Real > &  tau33i 
)

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
[in,out]tau13icontribution to stress from du/dz
[in,out]tau23icontribution to stress from dv/dz
[in,out]tau33icontribution to stress from dw/dz
34 {
35  Real OneThird = (1./3.);
36 
37  // NOTE: mu_eff includes factor of 2
38 
39  if (cell_data)
40  // constant alpha (stored in mu_eff)
41  {
42  // Cell centered strains
43  ParallelFor(bxcc, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
44  {
45  Real rhoAlpha = cell_data(i, j, k, Rho_comp) * mu_eff;
46  if (tau33i) tau33i(i,j,k) = -rhoAlpha * tau33(i,j,k);
47  tau11(i,j,k) = -rhoAlpha * ( tau11(i,j,k) - OneThird*er_arr(i,j,k) );
48  tau22(i,j,k) = -rhoAlpha * ( tau22(i,j,k) - OneThird*er_arr(i,j,k) );
49  tau33(i,j,k) = -rhoAlpha * ( tau33(i,j,k) - OneThird*er_arr(i,j,k) );
50  });
51 
52  // Off-diagonal strains
53  ParallelFor(tbxxy,tbxxz,tbxyz,
54  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
55  Real rho_bar = 0.25*( cell_data(i-1, j , k, Rho_comp) + cell_data(i, j , k, Rho_comp)
56  + cell_data(i-1, j-1, k, Rho_comp) + cell_data(i, j-1, k, Rho_comp) );
57  tau12(i,j,k) *= -rho_bar * mu_eff;
58  },
59  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
60  Real rho_bar = 0.25*( cell_data(i-1, j, k , Rho_comp) + cell_data(i, j, k , Rho_comp)
61  + cell_data(i-1, j, k-1, Rho_comp) + cell_data(i, j, k-1, Rho_comp) );
62  tau13(i,j,k) *= -rho_bar * mu_eff;
63 
64  if (tau13i) tau13i(i,j,k) *= -rho_bar * mu_eff;
65  },
66  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
67  Real rho_bar = 0.25*( cell_data(i, j-1, k , Rho_comp) + cell_data(i, j, k , Rho_comp)
68  + cell_data(i, j-1, k-1, Rho_comp) + cell_data(i, j, k-1, Rho_comp) );
69  tau23(i,j,k) *= -rho_bar * mu_eff;
70 
71  if (tau23i) tau23i(i,j,k) *= -rho_bar * mu_eff;
72  });
73  }
74  else
75  // constant mu_eff
76  {
77  // Cell centered strains
78  ParallelFor(bxcc, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
79  {
80  if (tau33i) tau33i(i,j,k) = -mu_eff * tau33(i,j,k);
81  tau11(i,j,k) = -mu_eff * ( tau11(i,j,k) - OneThird*er_arr(i,j,k) );
82  tau22(i,j,k) = -mu_eff * ( tau22(i,j,k) - OneThird*er_arr(i,j,k) );
83  tau33(i,j,k) = -mu_eff * ( tau33(i,j,k) - OneThird*er_arr(i,j,k) );
84  });
85 
86  // Off-diagonal strains
87  ParallelFor(tbxxy,tbxxz,tbxyz,
88  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
89  tau12(i,j,k) *= -mu_eff;
90  },
91  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
92  tau13(i,j,k) *= -mu_eff;
93 
94  if (tau13i) tau13i(i,j,k) *= -mu_eff;
95  },
96  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
97  tau23(i,j,k) *= -mu_eff;
98 
99  if (tau23i) tau23i(i,j,k) *= -mu_eff;
100  });
101  }
102 }
@ tau12
Definition: ERF_DataStruct.H:31
@ tau23
Definition: ERF_DataStruct.H:31
@ tau33
Definition: ERF_DataStruct.H:31
@ tau22
Definition: ERF_DataStruct.H:31
@ tau11
Definition: ERF_DataStruct.H:31
@ tau13
Definition: ERF_DataStruct.H:31
#define Rho_comp
Definition: ERF_IndexDefines.H:36
amrex::Real Real
Definition: ERF_ShocInterface.H:19

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,
Array4< Real > &  tau13i,
Array4< Real > &  tau23i,
Array4< Real > &  tau33i 
)

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
[in,out]tau13icontribution to stress from du/dz
[in,out]tau23icontribution to stress from dv/dz
[in,out]tau33icontribution to stress from dw/dz
135 {
136  Real OneThird = (1./3.);
137 
138  // NOTE: mu_eff includes factor of 2
139 
140  if (cell_data)
141  // constant alpha (stored in mu_eff)
142  {
143  // Cell centered strains
144  ParallelFor(bxcc, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
145  Real rhoAlpha = cell_data(i, j, k, Rho_comp) * mu_eff;
146  Real mu_11 = rhoAlpha + 2.0 * mu_turb(i, j, k, EddyDiff::Mom_h);
147  Real mu_22 = mu_11;
148  Real mu_33 = rhoAlpha + 2.0 * mu_turb(i, j, k, EddyDiff::Mom_v);
149  if (tau33i) tau33i(i,j,k) = -mu_33 * tau33(i,j,k);
150  tau11(i,j,k) = -mu_11 * ( tau11(i,j,k) - OneThird*er_arr(i,j,k) );
151  tau22(i,j,k) = -mu_22 * ( tau22(i,j,k) - OneThird*er_arr(i,j,k) );
152  tau33(i,j,k) = -mu_33 * ( tau33(i,j,k) - OneThird*er_arr(i,j,k) );
153  });
154 
155  // Off-diagonal strains
156  ParallelFor(tbxxy,tbxxz,tbxyz,
157  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
158  Real rho_bar = 0.25*( cell_data(i-1, j , k, Rho_comp) + cell_data(i, j , k, Rho_comp)
159  + cell_data(i-1, j-1, k, Rho_comp) + cell_data(i, j-1, k, Rho_comp) );
160  Real mu_bar = 0.25*( mu_turb(i-1, j , k, EddyDiff::Mom_h) + mu_turb(i, j , k, EddyDiff::Mom_h)
161  + mu_turb(i-1, j-1, k, EddyDiff::Mom_h) + mu_turb(i, j-1, k, EddyDiff::Mom_h) );
162  Real mu_12 = rho_bar*mu_eff + 2.0*mu_bar;
163  tau12(i,j,k) *= -mu_12;
164  },
165  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
166  Real rho_bar = 0.25*( cell_data(i-1, j, k , Rho_comp) + cell_data(i, j, k , Rho_comp)
167  + cell_data(i-1, j, k-1, Rho_comp) + cell_data(i, j, k-1, Rho_comp) );
168  Real mu_bar = 0.25*( mu_turb(i-1, j, k , EddyDiff::Mom_v) + mu_turb(i, j, k , EddyDiff::Mom_v)
169  + mu_turb(i-1, j, k-1, EddyDiff::Mom_v) + mu_turb(i, j, k-1, EddyDiff::Mom_v) );
170  Real mu_13 = rho_bar*mu_eff + 2.0*mu_bar;
171  tau13(i,j,k) *= -mu_13;
172 
173  if (tau13i) tau13i(i,j,k) *= -mu_13;
174  },
175  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
176  Real rho_bar = 0.25*( cell_data(i, j-1, k , Rho_comp) + cell_data(i, j, k , Rho_comp)
177  + cell_data(i, j-1, k-1, Rho_comp) + cell_data(i, j, k-1, Rho_comp) );
178  Real mu_bar = 0.25*( mu_turb(i, j-1, k , EddyDiff::Mom_v) + mu_turb(i, j, k , EddyDiff::Mom_v)
179  + mu_turb(i, j-1, k-1, EddyDiff::Mom_v) + mu_turb(i, j, k-1, EddyDiff::Mom_v) );
180  Real mu_23 = rho_bar*mu_eff + 2.0*mu_bar;
181  tau23(i,j,k) *= -mu_23;
182 
183  if (tau23i) tau23i(i,j,k) *= -mu_23;
184  });
185  }
186  else
187  // constant mu_eff
188  {
189  // Cell centered strains
190  ParallelFor(bxcc, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
191  Real mu_11 = mu_eff + 2.0 * mu_turb(i, j, k, EddyDiff::Mom_h);
192  Real mu_22 = mu_11;
193  Real mu_33 = mu_eff + 2.0 * mu_turb(i, j, k, EddyDiff::Mom_v);
194  if (tau33i) tau33i(i,j,k) = -mu_33 * tau33(i,j,k);
195  tau11(i,j,k) = -mu_11 * ( tau11(i,j,k) - OneThird*er_arr(i,j,k) );
196  tau22(i,j,k) = -mu_22 * ( tau22(i,j,k) - OneThird*er_arr(i,j,k) );
197  tau33(i,j,k) = -mu_33 * ( tau33(i,j,k) - OneThird*er_arr(i,j,k) );
198  });
199 
200  // Off-diagonal strains
201  ParallelFor(tbxxy,tbxxz,tbxyz,
202  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
203  Real mu_bar = 0.25*( mu_turb(i-1, j , k, EddyDiff::Mom_h) + mu_turb(i, j , k, EddyDiff::Mom_h)
204  + mu_turb(i-1, j-1, k, EddyDiff::Mom_h) + mu_turb(i, j-1, k, EddyDiff::Mom_h) );
205  Real mu_12 = mu_eff + 2.0*mu_bar;
206  tau12(i,j,k) *= -mu_12;
207  },
208  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
209  Real mu_bar = 0.25*( mu_turb(i-1, j, k , EddyDiff::Mom_v) + mu_turb(i, j, k , EddyDiff::Mom_v)
210  + mu_turb(i-1, j, k-1, EddyDiff::Mom_v) + mu_turb(i, j, k-1, EddyDiff::Mom_v) );
211  Real mu_13 = mu_eff + 2.0*mu_bar;
212  tau13(i,j,k) *= -mu_13;
213 
214  if (tau13i) tau13i(i,j,k) *= -mu_13;
215  },
216  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
217  Real mu_bar = 0.25*( mu_turb(i, j-1, k , EddyDiff::Mom_v) + mu_turb(i, j, k , EddyDiff::Mom_v)
218  + mu_turb(i, j-1, k-1, EddyDiff::Mom_v) + mu_turb(i, j, k-1, EddyDiff::Mom_v) );
219  Real mu_23 = mu_eff + 2.0*mu_bar;
220  tau23(i,j,k) *= -mu_23;
221 
222  if (tau23i) tau23i(i,j,k) *= -mu_23;
223  });
224  }
225 }
@ Mom_h
Definition: ERF_IndexDefines.H:170
@ Mom_v
Definition: ERF_IndexDefines.H:175

Referenced by erf_make_tau_terms().

Here is the caller graph for this function: