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

Functions

void ComputeStressConsVisc_EB (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< const Real > &vfrac, Array4< Real > &tau13i, Array4< Real > &tau23i, Array4< Real > &tau33i)
 
void ComputeStressVarVisc_EB (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< const Real > &vfrac, Array4< Real > &tau13i, Array4< Real > &tau23i, Array4< Real > &tau33i)
 

Function Documentation

◆ ComputeStressConsVisc_EB()

void ComputeStressConsVisc_EB ( 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< const Real > &  vfrac,
Array4< Real > &  tau13i,
Array4< Real > &  tau23i,
Array4< Real > &  tau33i 
)

Function for computing the stress with constant viscosity for EB.

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
35 {
36  Real OneThird = (1./3.);
37 
38  // NOTE: mu_eff includes factor of 2
39 
40  if (cell_data)
41  // constant alpha (stored in mu_eff)
42  {
43  // Cell centered strains
44  ParallelFor(bxcc, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
45  {
46  Real rhoAlpha = cell_data(i, j, k, Rho_comp) * mu_eff;
47  if (tau33i) tau33i(i,j,k) = -rhoAlpha * tau33(i,j,k);
48  tau11(i,j,k) = -rhoAlpha * ( tau11(i,j,k) - OneThird*er_arr(i,j,k) );
49  tau22(i,j,k) = -rhoAlpha * ( tau22(i,j,k) - OneThird*er_arr(i,j,k) );
50  tau33(i,j,k) = -rhoAlpha * ( tau33(i,j,k) - OneThird*er_arr(i,j,k) );
51  });
52 
53  // Off-diagonal strains
54  ParallelFor(tbxxy,tbxxz,tbxyz,
55  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
56  Real vol_sum = vfrac(i,j,k) + vfrac(i-1,j,k) + vfrac(i,j-1,k) + vfrac(i-1,j-1,k);
57  Real rho_bar = 0.0;
58  if (vol_sum > 1.e-16) {
59  rho_bar = ( vfrac(i-1,j,k) * cell_data(i-1, j , k, Rho_comp)
60  + vfrac(i,j,k) * cell_data(i, j , k, Rho_comp)
61  + vfrac(i-1,j-1,k) * cell_data(i-1, j-1, k, Rho_comp)
62  + vfrac(i,j-1,k) * cell_data(i, j-1, k, Rho_comp) ) / vol_sum;
63  } else {
64  rho_bar = 0.25*( cell_data(i-1, j , k, Rho_comp) + cell_data(i, j , k, Rho_comp)
65  + cell_data(i-1, j-1, k, Rho_comp) + cell_data(i, j-1, k, Rho_comp) );
66  }
67  tau12(i,j,k) *= -rho_bar * mu_eff;
68  },
69  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
70  Real vol_sum = vfrac(i,j,k) + vfrac(i-1,j,k) + vfrac(i,j,k-1) + vfrac(i-1,j,k-1);
71  Real rho_bar = 0.0;
72  if (vol_sum > 1.e-16) {
73  rho_bar = ( vfrac(i-1,j,k) * cell_data(i-1, j, k , Rho_comp)
74  + vfrac(i,j,k) * cell_data(i, j, k , Rho_comp)
75  + vfrac(i-1,j,k-1) * cell_data(i-1, j, k-1, Rho_comp)
76  + vfrac(i,j,k-1) * cell_data(i, j, k-1, Rho_comp) )/ vol_sum;
77  } else {
78  rho_bar = 0.25*( cell_data(i-1, j, k , Rho_comp) + cell_data(i, j, k , Rho_comp)
79  + cell_data(i-1, j, k-1, Rho_comp) + cell_data(i, j, k-1, Rho_comp) );
80  }
81  tau13(i,j,k) *= -rho_bar * mu_eff;
82 
83  if (tau13i) tau13i(i,j,k) *= -rho_bar * mu_eff;
84  },
85  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
86  Real vol_sum = vfrac(i,j,k) + vfrac(i,j-1,k) + vfrac(i,j,k-1) + vfrac(i,j-1,k-1);
87  Real rho_bar = 0.0;
88  if (vol_sum > 1.e-16) {
89  rho_bar = ( vfrac(i,j-1,k) * cell_data(i, j-1, k , Rho_comp)
90  + vfrac(i,j,k) * cell_data(i, j, k , Rho_comp)
91  + vfrac(i,j-1,k-1) * cell_data(i, j-1, k-1, Rho_comp)
92  + vfrac(i,j,k-1) * cell_data(i, j, k-1, Rho_comp) ) / vol_sum;
93  } else {
94  rho_bar = 0.25*( cell_data(i, j-1, k , Rho_comp) + cell_data(i, j, k , Rho_comp)
95  + cell_data(i, j-1, k-1, Rho_comp) + cell_data(i, j, k-1, Rho_comp) );
96  }
97  tau23(i,j,k) *= -rho_bar * mu_eff;
98 
99  if (tau23i) tau23i(i,j,k) *= -rho_bar * mu_eff;
100  });
101  }
102  else
103  // constant mu_eff
104  {
105  // Cell centered strains
106  ParallelFor(bxcc, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept
107  {
108  if (tau33i) tau33i(i,j,k) = -mu_eff * tau33(i,j,k);
109  tau11(i,j,k) = -mu_eff * ( tau11(i,j,k) - OneThird*er_arr(i,j,k) );
110  tau22(i,j,k) = -mu_eff * ( tau22(i,j,k) - OneThird*er_arr(i,j,k) );
111  tau33(i,j,k) = -mu_eff * ( tau33(i,j,k) - OneThird*er_arr(i,j,k) );
112  });
113 
114  // Off-diagonal strains
115  ParallelFor(tbxxy,tbxxz,tbxyz,
116  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
117  tau12(i,j,k) *= -mu_eff;
118  },
119  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
120  tau13(i,j,k) *= -mu_eff;
121 
122  if (tau13i) tau13i(i,j,k) *= -mu_eff;
123  },
124  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
125  tau23(i,j,k) *= -mu_eff;
126 
127  if (tau23i) tau23i(i,j,k) *= -mu_eff;
128  });
129  }
130 }
@ tau12
Definition: ERF_DataStruct.H:32
@ tau23
Definition: ERF_DataStruct.H:32
@ tau33
Definition: ERF_DataStruct.H:32
@ tau22
Definition: ERF_DataStruct.H:32
@ tau11
Definition: ERF_DataStruct.H:32
@ tau13
Definition: ERF_DataStruct.H:32
#define Rho_comp
Definition: ERF_IndexDefines.H:36
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

Referenced by erf_make_tau_terms().

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

◆ ComputeStressVarVisc_EB()

void ComputeStressVarVisc_EB ( 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< const Real > &  vfrac,
Array4< Real > &  tau13i,
Array4< Real > &  tau23i,
Array4< Real > &  tau33i 
)

Function for computing the stress with variable viscosity for EB. rho_bar at off-diagonal locations uses vfrac-weighted averaging as in ComputeStressConsVisc_EB; mu_turb is averaged with the same vfrac weights.

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_turbturbulent viscosity (cell-centered)
[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]vfracvolume fractions
[in,out]tau13icontribution to stress from du/dz
[in,out]tau23icontribution to stress from dv/dz
[in,out]tau33icontribution to stress from dw/dz
167 {
168  Real OneThird = (1./3.);
169 
170  // NOTE: mu_eff includes factor of 2
171 
172  if (cell_data)
173  // constant alpha (stored in mu_eff)
174  {
175  // Cell centered strains (no EB correction needed — cell-centered quantities)
176  ParallelFor(bxcc, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
177  if (vfrac(i,j,k) > 0.0) {
178  Real rhoAlpha = cell_data(i, j, k, Rho_comp) * mu_eff;
179  Real mu_11 = rhoAlpha + 2.0 * mu_turb(i, j, k, EddyDiff::Mom_h);
180  Real mu_22 = mu_11;
181  Real mu_33 = rhoAlpha + 2.0 * mu_turb(i, j, k, EddyDiff::Mom_v);
182  if (tau33i) tau33i(i,j,k) = -mu_33 * tau33(i,j,k);
183  tau11(i,j,k) = -mu_11 * ( tau11(i,j,k) - OneThird*er_arr(i,j,k) );
184  tau22(i,j,k) = -mu_22 * ( tau22(i,j,k) - OneThird*er_arr(i,j,k) );
185  tau33(i,j,k) = -mu_33 * ( tau33(i,j,k) - OneThird*er_arr(i,j,k) );
186  } else {
187  if (tau33i) tau33i(i,j,k) = 0.0;
188  tau11(i,j,k) = 0.0;
189  tau22(i,j,k) = 0.0;
190  tau33(i,j,k) = 0.0;
191  }
192  });
193 
194  // Off-diagonal strains: vfrac-weighted rho_bar and mu_bar
195  ParallelFor(tbxxy,tbxxz,tbxyz,
196  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
197  Real vol_sum = vfrac(i,j,k) + vfrac(i-1,j,k) + vfrac(i,j-1,k) + vfrac(i-1,j-1,k);
198  Real rho_bar, mu_bar;
199  if (vol_sum > 1.e-16) {
200  rho_bar = ( vfrac(i-1,j ,k) * cell_data(i-1, j , k, Rho_comp)
201  + vfrac(i ,j ,k) * cell_data(i , j , k, Rho_comp)
202  + vfrac(i-1,j-1,k) * cell_data(i-1, j-1, k, Rho_comp)
203  + vfrac(i ,j-1,k) * cell_data(i , j-1, k, Rho_comp) ) / vol_sum;
204  mu_bar = ( vfrac(i-1,j ,k) * mu_turb(i-1, j , k, EddyDiff::Mom_h)
205  + vfrac(i ,j ,k) * mu_turb(i , j , k, EddyDiff::Mom_h)
206  + vfrac(i-1,j-1,k) * mu_turb(i-1, j-1, k, EddyDiff::Mom_h)
207  + vfrac(i ,j-1,k) * mu_turb(i , j-1, k, EddyDiff::Mom_h) ) / vol_sum;
208  } else {
209  rho_bar = 0.25*( cell_data(i-1, j , k, Rho_comp) + cell_data(i , j , k, Rho_comp)
210  + cell_data(i-1, j-1, k, Rho_comp) + cell_data(i , j-1, k, Rho_comp) );
211  mu_bar = 0.25*( mu_turb(i-1, j , k, EddyDiff::Mom_h) + mu_turb(i , j , k, EddyDiff::Mom_h)
212  + mu_turb(i-1, j-1, k, EddyDiff::Mom_h) + mu_turb(i , j-1, k, EddyDiff::Mom_h) );
213  }
214  Real mu_12 = rho_bar*mu_eff + 2.0*mu_bar;
215  tau12(i,j,k) *= -mu_12;
216  },
217  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
218  Real vol_sum = vfrac(i,j,k) + vfrac(i-1,j,k) + vfrac(i,j,k-1) + vfrac(i-1,j,k-1);
219  Real rho_bar, mu_bar;
220  if (vol_sum > 1.e-16) {
221  rho_bar = ( vfrac(i-1,j,k ) * cell_data(i-1, j, k , Rho_comp)
222  + vfrac(i ,j,k ) * cell_data(i , j, k , Rho_comp)
223  + vfrac(i-1,j,k-1) * cell_data(i-1, j, k-1, Rho_comp)
224  + vfrac(i ,j,k-1) * cell_data(i , j, k-1, Rho_comp) ) / vol_sum;
225  mu_bar = ( vfrac(i-1,j,k ) * mu_turb(i-1, j, k , EddyDiff::Mom_v)
226  + vfrac(i ,j,k ) * mu_turb(i , j, k , EddyDiff::Mom_v)
227  + vfrac(i-1,j,k-1) * mu_turb(i-1, j, k-1, EddyDiff::Mom_v)
228  + vfrac(i ,j,k-1) * mu_turb(i , j, k-1, EddyDiff::Mom_v) ) / vol_sum;
229  } else {
230  rho_bar = 0.25*( cell_data(i-1, j, k , Rho_comp) + cell_data(i , j, k , Rho_comp)
231  + cell_data(i-1, j, k-1, Rho_comp) + cell_data(i , j, k-1, Rho_comp) );
232  mu_bar = 0.25*( mu_turb(i-1, j, k , EddyDiff::Mom_v) + mu_turb(i , j, k , EddyDiff::Mom_v)
233  + mu_turb(i-1, j, k-1, EddyDiff::Mom_v) + mu_turb(i , j, k-1, EddyDiff::Mom_v) );
234  }
235  Real mu_13 = rho_bar*mu_eff + 2.0*mu_bar;
236  tau13(i,j,k) *= -mu_13;
237  if (tau13i) tau13i(i,j,k) *= -mu_13;
238  },
239  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
240  Real vol_sum = vfrac(i,j,k) + vfrac(i,j-1,k) + vfrac(i,j,k-1) + vfrac(i,j-1,k-1);
241  Real rho_bar, mu_bar;
242  if (vol_sum > 1.e-16) {
243  rho_bar = ( vfrac(i,j-1,k ) * cell_data(i, j-1, k , Rho_comp)
244  + vfrac(i,j ,k ) * cell_data(i, j , k , Rho_comp)
245  + vfrac(i,j-1,k-1) * cell_data(i, j-1, k-1, Rho_comp)
246  + vfrac(i,j ,k-1) * cell_data(i, j , k-1, Rho_comp) ) / vol_sum;
247  mu_bar = ( vfrac(i,j-1,k ) * mu_turb(i, j-1, k , EddyDiff::Mom_v)
248  + vfrac(i,j ,k ) * mu_turb(i, j , k , EddyDiff::Mom_v)
249  + vfrac(i,j-1,k-1) * mu_turb(i, j-1, k-1, EddyDiff::Mom_v)
250  + vfrac(i,j ,k-1) * mu_turb(i, j , k-1, EddyDiff::Mom_v) ) / vol_sum;
251  } else {
252  rho_bar = 0.25*( cell_data(i, j-1, k , Rho_comp) + cell_data(i, j , k , Rho_comp)
253  + cell_data(i, j-1, k-1, Rho_comp) + cell_data(i, j , k-1, Rho_comp) );
254  mu_bar = 0.25*( mu_turb(i, j-1, k , EddyDiff::Mom_v) + mu_turb(i, j , k , EddyDiff::Mom_v)
255  + mu_turb(i, j-1, k-1, EddyDiff::Mom_v) + mu_turb(i, j , k-1, EddyDiff::Mom_v) );
256  }
257  Real mu_23 = rho_bar*mu_eff + 2.0*mu_bar;
258  tau23(i,j,k) *= -mu_23;
259  if (tau23i) tau23i(i,j,k) *= -mu_23;
260  });
261  }
262  else
263  // constant mu_eff
264  {
265  // Cell centered strains
266  ParallelFor(bxcc, [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
267  if (vfrac(i,j,k) > 0.0) {
268  Real mu_11 = mu_eff + 2.0 * mu_turb(i, j, k, EddyDiff::Mom_h);
269  Real mu_22 = mu_11;
270  Real mu_33 = mu_eff + 2.0 * mu_turb(i, j, k, EddyDiff::Mom_v);
271  if (tau33i) tau33i(i,j,k) = -mu_33 * tau33(i,j,k);
272  tau11(i,j,k) = -mu_11 * ( tau11(i,j,k) - OneThird*er_arr(i,j,k) );
273  tau22(i,j,k) = -mu_22 * ( tau22(i,j,k) - OneThird*er_arr(i,j,k) );
274  tau33(i,j,k) = -mu_33 * ( tau33(i,j,k) - OneThird*er_arr(i,j,k) );
275  } else {
276  if (tau33i) tau33i(i,j,k) = 0.0;
277  tau11(i,j,k) = 0.0;
278  tau22(i,j,k) = 0.0;
279  tau33(i,j,k) = 0.0;
280  }
281  });
282 
283  // Off-diagonal strains: vfrac-weighted mu_bar
284  ParallelFor(tbxxy,tbxxz,tbxyz,
285  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
286  Real vol_sum = vfrac(i,j,k) + vfrac(i-1,j,k) + vfrac(i,j-1,k) + vfrac(i-1,j-1,k);
287  Real mu_bar;
288  if (vol_sum > 1.e-16) {
289  mu_bar = ( vfrac(i-1,j ,k) * mu_turb(i-1, j , k, EddyDiff::Mom_h)
290  + vfrac(i ,j ,k) * mu_turb(i , j , k, EddyDiff::Mom_h)
291  + vfrac(i-1,j-1,k) * mu_turb(i-1, j-1, k, EddyDiff::Mom_h)
292  + vfrac(i ,j-1,k) * mu_turb(i , j-1, k, EddyDiff::Mom_h) ) / vol_sum;
293  } else {
294  mu_bar = 0.25*( mu_turb(i-1, j , k, EddyDiff::Mom_h) + mu_turb(i , j , k, EddyDiff::Mom_h)
295  + mu_turb(i-1, j-1, k, EddyDiff::Mom_h) + mu_turb(i , j-1, k, EddyDiff::Mom_h) );
296  }
297  Real mu_12 = mu_eff + 2.0*mu_bar;
298  tau12(i,j,k) *= -mu_12;
299  },
300  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
301  Real vol_sum = vfrac(i,j,k) + vfrac(i-1,j,k) + vfrac(i,j,k-1) + vfrac(i-1,j,k-1);
302  Real mu_bar;
303  if (vol_sum > 1.e-16) {
304  mu_bar = ( vfrac(i-1,j,k ) * mu_turb(i-1, j, k , EddyDiff::Mom_v)
305  + vfrac(i ,j,k ) * mu_turb(i , j, k , EddyDiff::Mom_v)
306  + vfrac(i-1,j,k-1) * mu_turb(i-1, j, k-1, EddyDiff::Mom_v)
307  + vfrac(i ,j,k-1) * mu_turb(i , j, k-1, EddyDiff::Mom_v) ) / vol_sum;
308  } else {
309  mu_bar = 0.25*( mu_turb(i-1, j, k , EddyDiff::Mom_v) + mu_turb(i , j, k , EddyDiff::Mom_v)
310  + mu_turb(i-1, j, k-1, EddyDiff::Mom_v) + mu_turb(i , j, k-1, EddyDiff::Mom_v) );
311  }
312  Real mu_13 = mu_eff + 2.0*mu_bar;
313  tau13(i,j,k) *= -mu_13;
314  if (tau13i) tau13i(i,j,k) *= -mu_13;
315  },
316  [=] AMREX_GPU_DEVICE (int i, int j, int k) noexcept {
317  Real vol_sum = vfrac(i,j,k) + vfrac(i,j-1,k) + vfrac(i,j,k-1) + vfrac(i,j-1,k-1);
318  Real mu_bar;
319  if (vol_sum > 1.e-16) {
320  mu_bar = ( vfrac(i,j-1,k ) * mu_turb(i, j-1, k , EddyDiff::Mom_v)
321  + vfrac(i,j ,k ) * mu_turb(i, j , k , EddyDiff::Mom_v)
322  + vfrac(i,j-1,k-1) * mu_turb(i, j-1, k-1, EddyDiff::Mom_v)
323  + vfrac(i,j ,k-1) * mu_turb(i, j , k-1, EddyDiff::Mom_v) ) / vol_sum;
324  } else {
325  mu_bar = 0.25*( mu_turb(i, j-1, k , EddyDiff::Mom_v) + mu_turb(i, j , k , EddyDiff::Mom_v)
326  + mu_turb(i, j-1, k-1, EddyDiff::Mom_v) + mu_turb(i, j , k-1, EddyDiff::Mom_v) );
327  }
328  Real mu_23 = mu_eff + 2.0*mu_bar;
329  tau23(i,j,k) *= -mu_23;
330  if (tau23i) tau23i(i,j,k) *= -mu_23;
331  });
332  }
333 }
@ Mom_h
Definition: ERF_IndexDefines.H:186
@ Mom_v
Definition: ERF_IndexDefines.H:191

Referenced by erf_make_tau_terms().

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