-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathERF__EBAdvectionSrcForScalars_8H_source.html
373 lines (371 loc) · 50.5 KB
/
ERF__EBAdvectionSrcForScalars_8H_source.html
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
<!DOCTYPE html PUBLIC "-//W3C//DTD XHTML 1.0 Transitional//EN" "https://www.w3.org/TR/xhtml1/DTD/xhtml1-transitional.dtd">
<html xmlns="http://www.w3.org/1999/xhtml">
<head>
<meta http-equiv="Content-Type" content="text/xhtml;charset=UTF-8"/>
<meta http-equiv="X-UA-Compatible" content="IE=9"/>
<meta name="generator" content="Doxygen 1.9.1"/>
<meta name="viewport" content="width=device-width, initial-scale=1"/>
<title>ERF: Source/EB/ERF_EBAdvectionSrcForScalars.H Source File</title>
<link href="tabs.css" rel="stylesheet" type="text/css"/>
<script type="text/javascript" src="jquery.js"></script>
<script type="text/javascript" src="dynsections.js"></script>
<link href="navtree.css" rel="stylesheet" type="text/css"/>
<script type="text/javascript" src="resize.js"></script>
<script type="text/javascript" src="navtreedata.js"></script>
<script type="text/javascript" src="navtree.js"></script>
<link href="search/search.css" rel="stylesheet" type="text/css"/>
<script type="text/javascript" src="search/searchdata.js"></script>
<script type="text/javascript" src="search/search.js"></script>
<link href="doxygen.css" rel="stylesheet" type="text/css" />
</head>
<body>
<div id="top"><!-- do not remove this div, it is closed by doxygen! -->
<div id="titlearea">
<table cellspacing="0" cellpadding="0">
<tbody>
<tr style="height: 56px;">
<td id="projectalign" style="padding-left: 0.5em;">
<div id="projectname">ERF
</div>
<div id="projectbrief">Energy Research and Forecasting: An Atmospheric Modeling Code</div>
</td>
</tr>
</tbody>
</table>
</div>
<!-- end header part -->
<!-- Generated by Doxygen 1.9.1 -->
<script type="text/javascript">
/* @license magnet:?xt=urn:btih:cf05388f2679ee054f2beb29a391d25f4e673ac3&dn=gpl-2.0.txt GPL-v2 */
var searchBox = new SearchBox("searchBox", "search",false,'Search','.html');
/* @license-end */
</script>
<script type="text/javascript" src="menudata.js"></script>
<script type="text/javascript" src="menu.js"></script>
<script type="text/javascript">
/* @license magnet:?xt=urn:btih:cf05388f2679ee054f2beb29a391d25f4e673ac3&dn=gpl-2.0.txt GPL-v2 */
$(function() {
initMenu('',true,false,'search.php','Search');
$(document).ready(function() { init_search(); });
});
/* @license-end */</script>
<div id="main-nav"></div>
</div><!-- top -->
<div id="side-nav" class="ui-resizable side-nav-resizable">
<div id="nav-tree">
<div id="nav-tree-contents">
<div id="nav-sync" class="sync"></div>
</div>
</div>
<div id="splitbar" style="-moz-user-select:none;"
class="ui-resizable-handle">
</div>
</div>
<script type="text/javascript">
/* @license magnet:?xt=urn:btih:cf05388f2679ee054f2beb29a391d25f4e673ac3&dn=gpl-2.0.txt GPL-v2 */
$(document).ready(function(){initNavTree('ERF__EBAdvectionSrcForScalars_8H_source.html',''); initResizable(); });
/* @license-end */
</script>
<div id="doc-content">
<!-- window showing the filter options -->
<div id="MSearchSelectWindow"
onmouseover="return searchBox.OnSearchSelectShow()"
onmouseout="return searchBox.OnSearchSelectHide()"
onkeydown="return searchBox.OnSearchSelectKey(event)">
</div>
<!-- iframe showing the search results (closed by default) -->
<div id="MSearchResultsWindow">
<iframe src="javascript:void(0)" frameborder="0"
name="MSearchResults" id="MSearchResults">
</iframe>
</div>
<div class="header">
<div class="headertitle">
<div class="title">ERF_EBAdvectionSrcForScalars.H</div> </div>
</div><!--header-->
<div class="contents">
<a href="ERF__EBAdvectionSrcForScalars_8H.html">Go to the documentation of this file.</a><div class="fragment"><div class="line"><a name="l00001"></a><span class="lineno"> 1</span> <span class="preprocessor">#include <<a class="code" href="ERF__IndexDefines_8H.html">ERF_IndexDefines.H</a>></span></div>
<div class="line"><a name="l00002"></a><span class="lineno"> 2</span> <span class="preprocessor">#include <<a class="code" href="ERF__Interpolation_8H.html">ERF_Interpolation.H</a>></span></div>
<div class="line"><a name="l00003"></a><span class="lineno"> 3</span>  </div>
<div class="line"><a name="l00004"></a><span class="lineno"> 4</span> <span class="preprocessor">#include <iostream></span></div>
<div class="line"><a name="l00005"></a><span class="lineno"> 5</span> <span class="preprocessor">#include <fstream></span></div>
<div class="line"><a name="l00006"></a><span class="lineno"> 6</span> <span class="preprocessor">#include <AMReX_Vector.H></span></div>
<div class="line"><a name="l00007"></a><span class="lineno"> 7</span> <span class="comment"></span> </div>
<div class="line"><a name="l00008"></a><span class="lineno"> 8</span> <span class="comment">/**</span></div>
<div class="line"><a name="l00009"></a><span class="lineno"> 9</span> <span class="comment"> * When terrain_type == TerreinType::EB,</span></div>
<div class="line"><a name="l00010"></a><span class="lineno"> 10</span> <span class="comment"> * Wrapper function for computing the advective tendency w/ spatial order > 2</span></div>
<div class="line"><a name="l00011"></a><span class="lineno"> 11</span> <span class="comment"> */</span></div>
<div class="line"><a name="l00012"></a><span class="lineno"> 12</span> <span class="keyword">template</span><<span class="keyword">typename</span> InterpType_H, <span class="keyword">typename</span> InterpType_V></div>
<div class="line"><a name="l00013"></a><span class="lineno"> 13</span> <span class="keywordtype">void</span></div>
<div class="line"><a name="l00014"></a><span class="lineno"><a class="line" href="ERF__EBAdvectionSrcForScalars_8H.html#add5efd63d823caa0afcb44c57c9abb9f"> 14</a></span> <a class="code" href="ERF__EBAdvectionSrcForScalars_8H.html#add5efd63d823caa0afcb44c57c9abb9f">EBAdvectionSrcForScalarsWrapper</a> (<span class="keyword">const</span> amrex::Box& bx,</div>
<div class="line"><a name="l00015"></a><span class="lineno"> 15</span>  <span class="keyword">const</span> <span class="keywordtype">int</span>& ncomp, <span class="keyword">const</span> <span class="keywordtype">int</span>& icomp,</div>
<div class="line"><a name="l00016"></a><span class="lineno"> 16</span>  <span class="keyword">const</span> amrex::GpuArray<<span class="keyword">const</span> amrex::Array4<amrex::Real>, AMREX_SPACEDIM> flx_arr,</div>
<div class="line"><a name="l00017"></a><span class="lineno"> 17</span>  <span class="keyword">const</span> amrex::Array4<const amrex::Real>& cell_prim,</div>
<div class="line"><a name="l00018"></a><span class="lineno"> 18</span>  <span class="keyword">const</span> amrex::Array4<const amrex::Real>& avg_xmom,</div>
<div class="line"><a name="l00019"></a><span class="lineno"> 19</span>  <span class="keyword">const</span> amrex::Array4<const amrex::Real>& avg_ymom,</div>
<div class="line"><a name="l00020"></a><span class="lineno"> 20</span>  <span class="keyword">const</span> amrex::Array4<const amrex::Real>& avg_zmom,</div>
<div class="line"><a name="l00021"></a><span class="lineno"> 21</span>  <span class="keyword">const</span> amrex::Array4<const amrex::EBCellFlag>& cellflag,</div>
<div class="line"><a name="l00022"></a><span class="lineno"> 22</span>  <span class="keyword">const</span> amrex::Array4<const amrex::Real>& ax_arr,</div>
<div class="line"><a name="l00023"></a><span class="lineno"> 23</span>  <span class="keyword">const</span> amrex::Array4<const amrex::Real>& ay_arr,</div>
<div class="line"><a name="l00024"></a><span class="lineno"> 24</span>  <span class="keyword">const</span> amrex::Array4<const amrex::Real>& az_arr,</div>
<div class="line"><a name="l00025"></a><span class="lineno"> 25</span>  <span class="keyword">const</span> amrex::Real horiz_upw_frac,</div>
<div class="line"><a name="l00026"></a><span class="lineno"> 26</span>  <span class="keyword">const</span> amrex::Real vert_upw_frac)</div>
<div class="line"><a name="l00027"></a><span class="lineno"> 27</span> {</div>
<div class="line"><a name="l00028"></a><span class="lineno"> 28</span>  <span class="comment">// Instantiate structs for vert/horiz interp</span></div>
<div class="line"><a name="l00029"></a><span class="lineno"> 29</span>  InterpType_H interp_prim_h(cell_prim, horiz_upw_frac);</div>
<div class="line"><a name="l00030"></a><span class="lineno"> 30</span>  InterpType_V interp_prim_v(cell_prim, vert_upw_frac);</div>
<div class="line"><a name="l00031"></a><span class="lineno"> 31</span>  </div>
<div class="line"><a name="l00032"></a><span class="lineno"> 32</span>  <span class="keyword">const</span> <span class="keywordtype">int</span> ncells_h = interp_prim_h.GetUpwindCellNumber();</div>
<div class="line"><a name="l00033"></a><span class="lineno"> 33</span>  <span class="keyword">const</span> <span class="keywordtype">int</span> ncells_v = interp_prim_v.GetUpwindCellNumber();</div>
<div class="line"><a name="l00034"></a><span class="lineno"> 34</span>  </div>
<div class="line"><a name="l00035"></a><span class="lineno"> 35</span>  <span class="comment">// Instantiate structs for lower-order vert/horiz interp</span></div>
<div class="line"><a name="l00036"></a><span class="lineno"> 36</span>  <a class="code" href="structCENTERED2.html">CENTERED2</a> interp_prim_CEN2(cell_prim, 0);</div>
<div class="line"><a name="l00037"></a><span class="lineno"> 37</span>  <a class="code" href="structUPWIND3.html">UPWIND3</a> interp_prim_UPW3(cell_prim, horiz_upw_frac);</div>
<div class="line"><a name="l00038"></a><span class="lineno"> 38</span>  </div>
<div class="line"><a name="l00039"></a><span class="lineno"> 39</span>  <span class="keyword">const</span> amrex::Box xbx = amrex::surroundingNodes(bx,0);</div>
<div class="line"><a name="l00040"></a><span class="lineno"> 40</span>  <span class="keyword">const</span> amrex::Box ybx = amrex::surroundingNodes(bx,1);</div>
<div class="line"><a name="l00041"></a><span class="lineno"> 41</span>  <span class="keyword">const</span> amrex::Box zbx = amrex::surroundingNodes(bx,2);</div>
<div class="line"><a name="l00042"></a><span class="lineno"> 42</span>  </div>
<div class="line"><a name="l00043"></a><span class="lineno"> 43</span>  amrex::ParallelFor(xbx, ncomp,[=] AMREX_GPU_DEVICE (<span class="keywordtype">int</span> i, <span class="keywordtype">int</span> j, <span class="keywordtype">int</span> k, <span class="keywordtype">int</span> n) noexcept</div>
<div class="line"><a name="l00044"></a><span class="lineno"> 44</span>  {</div>
<div class="line"><a name="l00045"></a><span class="lineno"> 45</span>  <span class="keyword">const</span> <span class="keywordtype">int</span> cons_index = icomp + n;</div>
<div class="line"><a name="l00046"></a><span class="lineno"> 46</span>  </div>
<div class="line"><a name="l00047"></a><span class="lineno"> 47</span>  <span class="keywordflow">if</span> ( ax_arr(i,j,k) > 0. ) <span class="comment">// Do we need a threshold value for the area fraction to avoid infinitesimal fractions?</span></div>
<div class="line"><a name="l00048"></a><span class="lineno"> 48</span>  {</div>
<div class="line"><a name="l00049"></a><span class="lineno"> 49</span>  <span class="comment">// Find the highest upwind order based on the cell flag</span></div>
<div class="line"><a name="l00050"></a><span class="lineno"> 50</span>  </div>
<div class="line"><a name="l00051"></a><span class="lineno"> 51</span>  <span class="keywordtype">int</span> icell = 0;</div>
<div class="line"><a name="l00052"></a><span class="lineno"> 52</span>  </div>
<div class="line"><a name="l00053"></a><span class="lineno"> 53</span>  <span class="keywordflow">if</span> (avg_xmom(i,j,k)>0.)</div>
<div class="line"><a name="l00054"></a><span class="lineno"> 54</span>  {</div>
<div class="line"><a name="l00055"></a><span class="lineno"> 55</span>  <span class="keywordflow">for</span> (icell=0; icell<ncells_h; ++icell) {</div>
<div class="line"><a name="l00056"></a><span class="lineno"> 56</span>  <span class="keywordflow">if</span> (cellflag(i-icell-1,j,k).isCovered()) {</div>
<div class="line"><a name="l00057"></a><span class="lineno"> 57</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00058"></a><span class="lineno"> 58</span>  }</div>
<div class="line"><a name="l00059"></a><span class="lineno"> 59</span>  }</div>
<div class="line"><a name="l00060"></a><span class="lineno"> 60</span>  }</div>
<div class="line"><a name="l00061"></a><span class="lineno"> 61</span>  <span class="keywordflow">else</span> <span class="keywordflow">if</span> (avg_xmom(i,j,k)<0.)</div>
<div class="line"><a name="l00062"></a><span class="lineno"> 62</span>  {</div>
<div class="line"><a name="l00063"></a><span class="lineno"> 63</span>  <span class="keywordflow">for</span> (icell=0; icell<ncells_h; ++icell) {</div>
<div class="line"><a name="l00064"></a><span class="lineno"> 64</span>  <span class="keywordflow">if</span> (cellflag(i+icell,j,k).isCovered()) {</div>
<div class="line"><a name="l00065"></a><span class="lineno"> 65</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00066"></a><span class="lineno"> 66</span>  }</div>
<div class="line"><a name="l00067"></a><span class="lineno"> 67</span>  }</div>
<div class="line"><a name="l00068"></a><span class="lineno"> 68</span>  }</div>
<div class="line"><a name="l00069"></a><span class="lineno"> 69</span>  </div>
<div class="line"><a name="l00070"></a><span class="lineno"> 70</span>  <span class="comment">// Interpolate the scalar variable (cell_prim) using the highest-order scheme</span></div>
<div class="line"><a name="l00071"></a><span class="lineno"> 71</span>  </div>
<div class="line"><a name="l00072"></a><span class="lineno"> 72</span>  <span class="keyword">const</span> <span class="keywordtype">int</span> prim_index = cons_index - 1;</div>
<div class="line"><a name="l00073"></a><span class="lineno"> 73</span>  amrex::Real interpx(0.);</div>
<div class="line"><a name="l00074"></a><span class="lineno"> 74</span>  </div>
<div class="line"><a name="l00075"></a><span class="lineno"> 75</span>  <span class="keywordflow">if</span> (icell==ncells_h-1) {</div>
<div class="line"><a name="l00076"></a><span class="lineno"> 76</span>  interp_prim_h.InterpolateInX(i,j,k,prim_index,interpx,avg_xmom(i,j,k));</div>
<div class="line"><a name="l00077"></a><span class="lineno"> 77</span>  } <span class="keywordflow">else</span> {</div>
<div class="line"><a name="l00078"></a><span class="lineno"> 78</span>  <span class="keywordflow">if</span> (icell==1) {</div>
<div class="line"><a name="l00079"></a><span class="lineno"> 79</span>  interp_prim_CEN2.<a class="code" href="structCENTERED2.html#adb613a538275649aec33db89afd0e822">InterpolateInX</a>(i,j,k,prim_index,interpx,avg_xmom(i,j,k));</div>
<div class="line"><a name="l00080"></a><span class="lineno"> 80</span>  } <span class="keywordflow">else</span> <span class="keywordflow">if</span> (icell==2) {</div>
<div class="line"><a name="l00081"></a><span class="lineno"> 81</span>  interp_prim_UPW3.<a class="code" href="structUPWIND3.html#afd23b4297df2e3daf4a1b1e40cf40fe2">InterpolateInX</a>(i,j,k,prim_index,interpx,avg_xmom(i,j,k));</div>
<div class="line"><a name="l00082"></a><span class="lineno"> 82</span>  }</div>
<div class="line"><a name="l00083"></a><span class="lineno"> 83</span>  }</div>
<div class="line"><a name="l00084"></a><span class="lineno"> 84</span>  (flx_arr[0])(i,j,k,cons_index) = avg_xmom(i,j,k) * interpx;</div>
<div class="line"><a name="l00085"></a><span class="lineno"> 85</span>  }</div>
<div class="line"><a name="l00086"></a><span class="lineno"> 86</span>  <span class="keywordflow">else</span></div>
<div class="line"><a name="l00087"></a><span class="lineno"> 87</span>  {</div>
<div class="line"><a name="l00088"></a><span class="lineno"> 88</span>  (flx_arr[0])(i,j,k,cons_index) = 0.;</div>
<div class="line"><a name="l00089"></a><span class="lineno"> 89</span>  }</div>
<div class="line"><a name="l00090"></a><span class="lineno"> 90</span>  });</div>
<div class="line"><a name="l00091"></a><span class="lineno"> 91</span>  </div>
<div class="line"><a name="l00092"></a><span class="lineno"> 92</span>  amrex::ParallelFor(ybx, ncomp,[=] AMREX_GPU_DEVICE (<span class="keywordtype">int</span> i, <span class="keywordtype">int</span> j, <span class="keywordtype">int</span> k, <span class="keywordtype">int</span> n) noexcept</div>
<div class="line"><a name="l00093"></a><span class="lineno"> 93</span>  {</div>
<div class="line"><a name="l00094"></a><span class="lineno"> 94</span>  <span class="keyword">const</span> <span class="keywordtype">int</span> cons_index = icomp + n;</div>
<div class="line"><a name="l00095"></a><span class="lineno"> 95</span>  </div>
<div class="line"><a name="l00096"></a><span class="lineno"> 96</span>  <span class="keywordflow">if</span> ( ay_arr(i,j,k) > 0. )</div>
<div class="line"><a name="l00097"></a><span class="lineno"> 97</span>  {</div>
<div class="line"><a name="l00098"></a><span class="lineno"> 98</span>  <span class="comment">// Find the highest upwind order based on the cell flag</span></div>
<div class="line"><a name="l00099"></a><span class="lineno"> 99</span>  </div>
<div class="line"><a name="l00100"></a><span class="lineno"> 100</span>  <span class="keywordtype">int</span> jcell = 0;</div>
<div class="line"><a name="l00101"></a><span class="lineno"> 101</span>  </div>
<div class="line"><a name="l00102"></a><span class="lineno"> 102</span>  <span class="keywordflow">if</span> (avg_ymom(i,j,k)>0.)</div>
<div class="line"><a name="l00103"></a><span class="lineno"> 103</span>  {</div>
<div class="line"><a name="l00104"></a><span class="lineno"> 104</span>  <span class="keywordflow">for</span> (jcell=0; jcell<ncells_h; ++jcell) {</div>
<div class="line"><a name="l00105"></a><span class="lineno"> 105</span>  <span class="keywordflow">if</span> (cellflag(i,j-jcell-1,k).isCovered()) {</div>
<div class="line"><a name="l00106"></a><span class="lineno"> 106</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00107"></a><span class="lineno"> 107</span>  }</div>
<div class="line"><a name="l00108"></a><span class="lineno"> 108</span>  }</div>
<div class="line"><a name="l00109"></a><span class="lineno"> 109</span>  }</div>
<div class="line"><a name="l00110"></a><span class="lineno"> 110</span>  <span class="keywordflow">else</span> <span class="keywordflow">if</span> (avg_ymom(i,j,k)<0.)</div>
<div class="line"><a name="l00111"></a><span class="lineno"> 111</span>  {</div>
<div class="line"><a name="l00112"></a><span class="lineno"> 112</span>  <span class="keywordflow">for</span> (jcell=0; jcell<ncells_h; ++jcell) {</div>
<div class="line"><a name="l00113"></a><span class="lineno"> 113</span>  <span class="keywordflow">if</span> (cellflag(i,j+jcell,k).isCovered()) {</div>
<div class="line"><a name="l00114"></a><span class="lineno"> 114</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00115"></a><span class="lineno"> 115</span>  }</div>
<div class="line"><a name="l00116"></a><span class="lineno"> 116</span>  }</div>
<div class="line"><a name="l00117"></a><span class="lineno"> 117</span>  }</div>
<div class="line"><a name="l00118"></a><span class="lineno"> 118</span>  </div>
<div class="line"><a name="l00119"></a><span class="lineno"> 119</span>  <span class="comment">// Interpolate the scalar variable (cell_prim) using the highest-order scheme</span></div>
<div class="line"><a name="l00120"></a><span class="lineno"> 120</span>  </div>
<div class="line"><a name="l00121"></a><span class="lineno"> 121</span>  <span class="keyword">const</span> <span class="keywordtype">int</span> prim_index = cons_index - 1;</div>
<div class="line"><a name="l00122"></a><span class="lineno"> 122</span>  amrex::Real interpy(0.);</div>
<div class="line"><a name="l00123"></a><span class="lineno"> 123</span>  </div>
<div class="line"><a name="l00124"></a><span class="lineno"> 124</span>  <span class="keywordflow">if</span> (jcell==ncells_h-1) {</div>
<div class="line"><a name="l00125"></a><span class="lineno"> 125</span>  interp_prim_h.InterpolateInY(i,j,k,prim_index,interpy,avg_ymom(i,j,k));</div>
<div class="line"><a name="l00126"></a><span class="lineno"> 126</span>  } <span class="keywordflow">else</span> {</div>
<div class="line"><a name="l00127"></a><span class="lineno"> 127</span>  <span class="keywordflow">if</span> (jcell==1) {</div>
<div class="line"><a name="l00128"></a><span class="lineno"> 128</span>  interp_prim_CEN2.<a class="code" href="structCENTERED2.html#a2b0e42ec56f22f97bb65ed044f6f77a0">InterpolateInY</a>(i,j,k,prim_index,interpy,avg_ymom(i,j,k));</div>
<div class="line"><a name="l00129"></a><span class="lineno"> 129</span>  } <span class="keywordflow">else</span> <span class="keywordflow">if</span> (jcell==2) {</div>
<div class="line"><a name="l00130"></a><span class="lineno"> 130</span>  interp_prim_UPW3.<a class="code" href="structUPWIND3.html#a6719e6c075ae163173c16f37b9e9c4b2">InterpolateInY</a>(i,j,k,prim_index,interpy,avg_ymom(i,j,k));</div>
<div class="line"><a name="l00131"></a><span class="lineno"> 131</span>  }</div>
<div class="line"><a name="l00132"></a><span class="lineno"> 132</span>  }</div>
<div class="line"><a name="l00133"></a><span class="lineno"> 133</span>  (flx_arr[1])(i,j,k,cons_index) = avg_ymom(i,j,k) * interpy;</div>
<div class="line"><a name="l00134"></a><span class="lineno"> 134</span>  }</div>
<div class="line"><a name="l00135"></a><span class="lineno"> 135</span>  <span class="keywordflow">else</span></div>
<div class="line"><a name="l00136"></a><span class="lineno"> 136</span>  {</div>
<div class="line"><a name="l00137"></a><span class="lineno"> 137</span>  (flx_arr[1])(i,j,k,cons_index) = 0.;</div>
<div class="line"><a name="l00138"></a><span class="lineno"> 138</span>  }</div>
<div class="line"><a name="l00139"></a><span class="lineno"> 139</span>  });</div>
<div class="line"><a name="l00140"></a><span class="lineno"> 140</span>  </div>
<div class="line"><a name="l00141"></a><span class="lineno"> 141</span>  interp_prim_UPW3.<a class="code" href="structUPWIND3.html#a424f3d6e525d4e195377c034deb8a345">SetUpwinding</a>(vert_upw_frac);</div>
<div class="line"><a name="l00142"></a><span class="lineno"> 142</span>  </div>
<div class="line"><a name="l00143"></a><span class="lineno"> 143</span>  amrex::ParallelFor(zbx, ncomp,[=] AMREX_GPU_DEVICE (<span class="keywordtype">int</span> i, <span class="keywordtype">int</span> j, <span class="keywordtype">int</span> k, <span class="keywordtype">int</span> n) noexcept</div>
<div class="line"><a name="l00144"></a><span class="lineno"> 144</span>  {</div>
<div class="line"><a name="l00145"></a><span class="lineno"> 145</span>  <span class="keyword">const</span> <span class="keywordtype">int</span> cons_index = icomp + n;</div>
<div class="line"><a name="l00146"></a><span class="lineno"> 146</span>  </div>
<div class="line"><a name="l00147"></a><span class="lineno"> 147</span>  <span class="keywordflow">if</span> ( az_arr(i,j,k) > 0. )</div>
<div class="line"><a name="l00148"></a><span class="lineno"> 148</span>  {</div>
<div class="line"><a name="l00149"></a><span class="lineno"> 149</span>  <span class="comment">// Find the highest upwind order based on the cell flag</span></div>
<div class="line"><a name="l00150"></a><span class="lineno"> 150</span>  </div>
<div class="line"><a name="l00151"></a><span class="lineno"> 151</span>  <span class="keywordtype">int</span> kcell = 0;</div>
<div class="line"><a name="l00152"></a><span class="lineno"> 152</span>  </div>
<div class="line"><a name="l00153"></a><span class="lineno"> 153</span>  <span class="keywordflow">if</span> (avg_zmom(i,j,k)>0.)</div>
<div class="line"><a name="l00154"></a><span class="lineno"> 154</span>  {</div>
<div class="line"><a name="l00155"></a><span class="lineno"> 155</span>  <span class="keywordflow">for</span> (kcell=0; kcell<ncells_v; ++kcell) {</div>
<div class="line"><a name="l00156"></a><span class="lineno"> 156</span>  <span class="keywordflow">if</span> (cellflag(i,j,k-kcell-1).isCovered()) {</div>
<div class="line"><a name="l00157"></a><span class="lineno"> 157</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00158"></a><span class="lineno"> 158</span>  }</div>
<div class="line"><a name="l00159"></a><span class="lineno"> 159</span>  }</div>
<div class="line"><a name="l00160"></a><span class="lineno"> 160</span>  }</div>
<div class="line"><a name="l00161"></a><span class="lineno"> 161</span>  <span class="keywordflow">else</span> <span class="keywordflow">if</span> (avg_zmom(i,j,k)<0.)</div>
<div class="line"><a name="l00162"></a><span class="lineno"> 162</span>  {</div>
<div class="line"><a name="l00163"></a><span class="lineno"> 163</span>  <span class="keywordflow">for</span> (kcell=0; kcell<ncells_v; ++kcell) {</div>
<div class="line"><a name="l00164"></a><span class="lineno"> 164</span>  <span class="keywordflow">if</span> (cellflag(i,j,k+kcell).isCovered()) {</div>
<div class="line"><a name="l00165"></a><span class="lineno"> 165</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00166"></a><span class="lineno"> 166</span>  }</div>
<div class="line"><a name="l00167"></a><span class="lineno"> 167</span>  }</div>
<div class="line"><a name="l00168"></a><span class="lineno"> 168</span>  }</div>
<div class="line"><a name="l00169"></a><span class="lineno"> 169</span>  </div>
<div class="line"><a name="l00170"></a><span class="lineno"> 170</span>  <span class="comment">// Interpolate the scalar variable (cell_prim) using the highest-order scheme</span></div>
<div class="line"><a name="l00171"></a><span class="lineno"> 171</span>  </div>
<div class="line"><a name="l00172"></a><span class="lineno"> 172</span>  <span class="keyword">const</span> <span class="keywordtype">int</span> prim_index = cons_index - 1;</div>
<div class="line"><a name="l00173"></a><span class="lineno"> 173</span>  amrex::Real interpz(0.);</div>
<div class="line"><a name="l00174"></a><span class="lineno"> 174</span>  </div>
<div class="line"><a name="l00175"></a><span class="lineno"> 175</span>  <span class="keywordflow">if</span> (kcell==ncells_v-1) {</div>
<div class="line"><a name="l00176"></a><span class="lineno"> 176</span>  interp_prim_v.InterpolateInZ(i,j,k,prim_index,interpz,avg_zmom(i,j,k));</div>
<div class="line"><a name="l00177"></a><span class="lineno"> 177</span>  } <span class="keywordflow">else</span> {</div>
<div class="line"><a name="l00178"></a><span class="lineno"> 178</span>  <span class="keywordflow">if</span> (kcell==1) {</div>
<div class="line"><a name="l00179"></a><span class="lineno"> 179</span>  interp_prim_CEN2.<a class="code" href="structCENTERED2.html#a66de83f17cc1fa457cb89af60afb0041">InterpolateInZ</a>(i,j,k,prim_index,interpz,avg_zmom(i,j,k));</div>
<div class="line"><a name="l00180"></a><span class="lineno"> 180</span>  } <span class="keywordflow">else</span> <span class="keywordflow">if</span> (kcell==2) {</div>
<div class="line"><a name="l00181"></a><span class="lineno"> 181</span>  interp_prim_UPW3.<a class="code" href="structUPWIND3.html#a420b9e841d549ce4f116205928f1743c">InterpolateInZ</a>(i,j,k,prim_index,interpz,avg_zmom(i,j,k));</div>
<div class="line"><a name="l00182"></a><span class="lineno"> 182</span>  }</div>
<div class="line"><a name="l00183"></a><span class="lineno"> 183</span>  }</div>
<div class="line"><a name="l00184"></a><span class="lineno"> 184</span>  (flx_arr[2])(i,j,k,cons_index) = avg_zmom(i,j,k) * interpz;</div>
<div class="line"><a name="l00185"></a><span class="lineno"> 185</span>  }</div>
<div class="line"><a name="l00186"></a><span class="lineno"> 186</span>  <span class="keywordflow">else</span></div>
<div class="line"><a name="l00187"></a><span class="lineno"> 187</span>  {</div>
<div class="line"><a name="l00188"></a><span class="lineno"> 188</span>  (flx_arr[2])(i,j,k,cons_index) = 0.;</div>
<div class="line"><a name="l00189"></a><span class="lineno"> 189</span>  }</div>
<div class="line"><a name="l00190"></a><span class="lineno"> 190</span>  });</div>
<div class="line"><a name="l00191"></a><span class="lineno"> 191</span> }</div>
<div class="line"><a name="l00192"></a><span class="lineno"> 192</span>  </div>
<div class="line"><a name="l00193"></a><span class="lineno"> 193</span> <span class="comment"></span> </div>
<div class="line"><a name="l00194"></a><span class="lineno"> 194</span> <span class="comment">/**</span></div>
<div class="line"><a name="l00195"></a><span class="lineno"> 195</span> <span class="comment"> * When terrain_type == TerreinType::EB,</span></div>
<div class="line"><a name="l00196"></a><span class="lineno"> 196</span> <span class="comment"> * Wrapper function for templating the vertical advective tendency w/ spatial order > 2.</span></div>
<div class="line"><a name="l00197"></a><span class="lineno"> 197</span> <span class="comment"> */</span></div>
<div class="line"><a name="l00198"></a><span class="lineno"> 198</span> <span class="keyword">template</span><<span class="keyword">typename</span> InterpType_H></div>
<div class="line"><a name="l00199"></a><span class="lineno"> 199</span> <span class="keywordtype">void</span></div>
<div class="line"><a name="l00200"></a><span class="lineno"><a class="line" href="ERF__EBAdvectionSrcForScalars_8H.html#a5c491c40d5d9bb1e7a6e507d4591a387"> 200</a></span> <a class="code" href="ERF__EBAdvectionSrcForScalars_8H.html#a5c491c40d5d9bb1e7a6e507d4591a387">EBAdvectionSrcForScalarsVert</a> (<span class="keyword">const</span> amrex::Box& bx,</div>
<div class="line"><a name="l00201"></a><span class="lineno"> 201</span>  <span class="keyword">const</span> <span class="keywordtype">int</span>& ncomp, <span class="keyword">const</span> <span class="keywordtype">int</span>& icomp,</div>
<div class="line"><a name="l00202"></a><span class="lineno"> 202</span>  <span class="keyword">const</span> amrex::GpuArray<<span class="keyword">const</span> amrex::Array4<amrex::Real>, AMREX_SPACEDIM> flx_arr,</div>
<div class="line"><a name="l00203"></a><span class="lineno"> 203</span>  <span class="keyword">const</span> amrex::Array4<const amrex::Real>& cell_prim,</div>
<div class="line"><a name="l00204"></a><span class="lineno"> 204</span>  <span class="keyword">const</span> amrex::Array4<const amrex::Real>& avg_xmom,</div>
<div class="line"><a name="l00205"></a><span class="lineno"> 205</span>  <span class="keyword">const</span> amrex::Array4<const amrex::Real>& avg_ymom,</div>
<div class="line"><a name="l00206"></a><span class="lineno"> 206</span>  <span class="keyword">const</span> amrex::Array4<const amrex::Real>& avg_zmom,</div>
<div class="line"><a name="l00207"></a><span class="lineno"> 207</span>  <span class="keyword">const</span> amrex::Array4<const amrex::EBCellFlag>& cellflag,</div>
<div class="line"><a name="l00208"></a><span class="lineno"> 208</span>  <span class="keyword">const</span> amrex::Array4<const amrex::Real>& ax_arr,</div>
<div class="line"><a name="l00209"></a><span class="lineno"> 209</span>  <span class="keyword">const</span> amrex::Array4<const amrex::Real>& ay_arr,</div>
<div class="line"><a name="l00210"></a><span class="lineno"> 210</span>  <span class="keyword">const</span> amrex::Array4<const amrex::Real>& az_arr,</div>
<div class="line"><a name="l00211"></a><span class="lineno"> 211</span>  <span class="keyword">const</span> amrex::Real horiz_upw_frac,</div>
<div class="line"><a name="l00212"></a><span class="lineno"> 212</span>  <span class="keyword">const</span> amrex::Real vert_upw_frac,</div>
<div class="line"><a name="l00213"></a><span class="lineno"> 213</span>  <span class="keyword">const</span> <a class="code" href="ERF__IndexDefines_8H.html#ada2d84e0dbbfb8d748defbf018748a70">AdvType</a> vert_adv_type)</div>
<div class="line"><a name="l00214"></a><span class="lineno"> 214</span> {</div>
<div class="line"><a name="l00215"></a><span class="lineno"> 215</span>  <span class="keywordflow">switch</span>(vert_adv_type) {</div>
<div class="line"><a name="l00216"></a><span class="lineno"> 216</span>  <span class="keywordflow">case</span> <a class="code" href="ERF__IndexDefines_8H.html#ada2d84e0dbbfb8d748defbf018748a70ada570868919bafaab1776170229a7049">AdvType::Centered_2nd</a>:</div>
<div class="line"><a name="l00217"></a><span class="lineno"> 217</span>  EBAdvectionSrcForScalarsWrapper<InterpType_H,CENTERED2>(bx, ncomp, icomp,</div>
<div class="line"><a name="l00218"></a><span class="lineno"> 218</span>  flx_arr, cell_prim,</div>
<div class="line"><a name="l00219"></a><span class="lineno"> 219</span>  avg_xmom, avg_ymom, avg_zmom,</div>
<div class="line"><a name="l00220"></a><span class="lineno"> 220</span>  cellflag, ax_arr, ay_arr, az_arr,</div>
<div class="line"><a name="l00221"></a><span class="lineno"> 221</span>  horiz_upw_frac, vert_upw_frac);</div>
<div class="line"><a name="l00222"></a><span class="lineno"> 222</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00223"></a><span class="lineno"> 223</span>  <span class="keywordflow">case</span> <a class="code" href="ERF__IndexDefines_8H.html#ada2d84e0dbbfb8d748defbf018748a70a40adf6565c1a31b2c72d8072412d0b9f">AdvType::Upwind_3rd</a>:</div>
<div class="line"><a name="l00224"></a><span class="lineno"> 224</span>  EBAdvectionSrcForScalarsWrapper<InterpType_H,UPWIND3>(bx, ncomp, icomp,</div>
<div class="line"><a name="l00225"></a><span class="lineno"> 225</span>  flx_arr, cell_prim,</div>
<div class="line"><a name="l00226"></a><span class="lineno"> 226</span>  avg_xmom, avg_ymom, avg_zmom,</div>
<div class="line"><a name="l00227"></a><span class="lineno"> 227</span>  cellflag, ax_arr, ay_arr, az_arr,</div>
<div class="line"><a name="l00228"></a><span class="lineno"> 228</span>  horiz_upw_frac, vert_upw_frac);</div>
<div class="line"><a name="l00229"></a><span class="lineno"> 229</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00230"></a><span class="lineno"> 230</span>  <span class="keywordflow">case</span> <a class="code" href="ERF__IndexDefines_8H.html#ada2d84e0dbbfb8d748defbf018748a70a778e0ed47937713eb1c4bedef6b00161">AdvType::Centered_4th</a>:</div>
<div class="line"><a name="l00231"></a><span class="lineno"> 231</span>  EBAdvectionSrcForScalarsWrapper<InterpType_H,CENTERED4>(bx, ncomp, icomp,</div>
<div class="line"><a name="l00232"></a><span class="lineno"> 232</span>  flx_arr, cell_prim,</div>
<div class="line"><a name="l00233"></a><span class="lineno"> 233</span>  avg_xmom, avg_ymom, avg_zmom,</div>
<div class="line"><a name="l00234"></a><span class="lineno"> 234</span>  cellflag, ax_arr, ay_arr, az_arr,</div>
<div class="line"><a name="l00235"></a><span class="lineno"> 235</span>  horiz_upw_frac, vert_upw_frac);</div>
<div class="line"><a name="l00236"></a><span class="lineno"> 236</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00237"></a><span class="lineno"> 237</span>  <span class="keywordflow">case</span> <a class="code" href="ERF__IndexDefines_8H.html#ada2d84e0dbbfb8d748defbf018748a70a83aec89012c35752597047c6ebdd983a">AdvType::Upwind_5th</a>:</div>
<div class="line"><a name="l00238"></a><span class="lineno"> 238</span>  EBAdvectionSrcForScalarsWrapper<InterpType_H,UPWIND5>(bx, ncomp, icomp,</div>
<div class="line"><a name="l00239"></a><span class="lineno"> 239</span>  flx_arr, cell_prim,</div>
<div class="line"><a name="l00240"></a><span class="lineno"> 240</span>  avg_xmom, avg_ymom, avg_zmom,</div>
<div class="line"><a name="l00241"></a><span class="lineno"> 241</span>  cellflag, ax_arr, ay_arr, az_arr,</div>
<div class="line"><a name="l00242"></a><span class="lineno"> 242</span>  horiz_upw_frac, vert_upw_frac);</div>
<div class="line"><a name="l00243"></a><span class="lineno"> 243</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00244"></a><span class="lineno"> 244</span>  <span class="keywordflow">case</span> <a class="code" href="ERF__IndexDefines_8H.html#ada2d84e0dbbfb8d748defbf018748a70a8c9475eb6bfde28e5bcb5230ac2e22a9">AdvType::Centered_6th</a>:</div>
<div class="line"><a name="l00245"></a><span class="lineno"> 245</span>  EBAdvectionSrcForScalarsWrapper<InterpType_H,CENTERED6>(bx, ncomp, icomp,</div>
<div class="line"><a name="l00246"></a><span class="lineno"> 246</span>  flx_arr, cell_prim,</div>
<div class="line"><a name="l00247"></a><span class="lineno"> 247</span>  avg_xmom, avg_ymom, avg_zmom,</div>
<div class="line"><a name="l00248"></a><span class="lineno"> 248</span>  cellflag, ax_arr, ay_arr, az_arr,</div>
<div class="line"><a name="l00249"></a><span class="lineno"> 249</span>  horiz_upw_frac, vert_upw_frac);</div>
<div class="line"><a name="l00250"></a><span class="lineno"> 250</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00251"></a><span class="lineno"> 251</span>  <span class="keywordflow">default</span>:</div>
<div class="line"><a name="l00252"></a><span class="lineno"> 252</span>  AMREX_ASSERT_WITH_MESSAGE(<span class="keyword">false</span>, <span class="stringliteral">"Unknown vertical advection scheme!"</span>);</div>
<div class="line"><a name="l00253"></a><span class="lineno"> 253</span>  }</div>
<div class="line"><a name="l00254"></a><span class="lineno"> 254</span> }</div>
<div class="line"><a name="l00255"></a><span class="lineno"> 255</span>  </div>
<div class="ttc" id="aERF__EBAdvectionSrcForScalars_8H_html_a5c491c40d5d9bb1e7a6e507d4591a387"><div class="ttname"><a href="ERF__EBAdvectionSrcForScalars_8H.html#a5c491c40d5d9bb1e7a6e507d4591a387">EBAdvectionSrcForScalarsVert</a></div><div class="ttdeci">void EBAdvectionSrcForScalarsVert(const amrex::Box &bx, const int &ncomp, const int &icomp, const amrex::GpuArray< const amrex::Array4< amrex::Real >, AMREX_SPACEDIM > flx_arr, const amrex::Array4< const amrex::Real > &cell_prim, const amrex::Array4< const amrex::Real > &avg_xmom, const amrex::Array4< const amrex::Real > &avg_ymom, const amrex::Array4< const amrex::Real > &avg_zmom, const amrex::Array4< const amrex::EBCellFlag > &cellflag, const amrex::Array4< const amrex::Real > &ax_arr, const amrex::Array4< const amrex::Real > &ay_arr, const amrex::Array4< const amrex::Real > &az_arr, const amrex::Real horiz_upw_frac, const amrex::Real vert_upw_frac, const AdvType vert_adv_type)</div><div class="ttdef"><b>Definition:</b> ERF_EBAdvectionSrcForScalars.H:200</div></div>
<div class="ttc" id="aERF__EBAdvectionSrcForScalars_8H_html_add5efd63d823caa0afcb44c57c9abb9f"><div class="ttname"><a href="ERF__EBAdvectionSrcForScalars_8H.html#add5efd63d823caa0afcb44c57c9abb9f">EBAdvectionSrcForScalarsWrapper</a></div><div class="ttdeci">void EBAdvectionSrcForScalarsWrapper(const amrex::Box &bx, const int &ncomp, const int &icomp, const amrex::GpuArray< const amrex::Array4< amrex::Real >, AMREX_SPACEDIM > flx_arr, const amrex::Array4< const amrex::Real > &cell_prim, const amrex::Array4< const amrex::Real > &avg_xmom, const amrex::Array4< const amrex::Real > &avg_ymom, const amrex::Array4< const amrex::Real > &avg_zmom, const amrex::Array4< const amrex::EBCellFlag > &cellflag, const amrex::Array4< const amrex::Real > &ax_arr, const amrex::Array4< const amrex::Real > &ay_arr, const amrex::Array4< const amrex::Real > &az_arr, const amrex::Real horiz_upw_frac, const amrex::Real vert_upw_frac)</div><div class="ttdef"><b>Definition:</b> ERF_EBAdvectionSrcForScalars.H:14</div></div>
<div class="ttc" id="aERF__IndexDefines_8H_html"><div class="ttname"><a href="ERF__IndexDefines_8H.html">ERF_IndexDefines.H</a></div></div>
<div class="ttc" id="aERF__IndexDefines_8H_html_ada2d84e0dbbfb8d748defbf018748a70"><div class="ttname"><a href="ERF__IndexDefines_8H.html#ada2d84e0dbbfb8d748defbf018748a70">AdvType</a></div><div class="ttdeci">AdvType</div><div class="ttdef"><b>Definition:</b> ERF_IndexDefines.H:203</div></div>
<div class="ttc" id="aERF__IndexDefines_8H_html_ada2d84e0dbbfb8d748defbf018748a70a40adf6565c1a31b2c72d8072412d0b9f"><div class="ttname"><a href="ERF__IndexDefines_8H.html#ada2d84e0dbbfb8d748defbf018748a70a40adf6565c1a31b2c72d8072412d0b9f">AdvType::Upwind_3rd</a></div><div class="ttdeci">@ Upwind_3rd</div></div>
<div class="ttc" id="aERF__IndexDefines_8H_html_ada2d84e0dbbfb8d748defbf018748a70a778e0ed47937713eb1c4bedef6b00161"><div class="ttname"><a href="ERF__IndexDefines_8H.html#ada2d84e0dbbfb8d748defbf018748a70a778e0ed47937713eb1c4bedef6b00161">AdvType::Centered_4th</a></div><div class="ttdeci">@ Centered_4th</div></div>
<div class="ttc" id="aERF__IndexDefines_8H_html_ada2d84e0dbbfb8d748defbf018748a70a83aec89012c35752597047c6ebdd983a"><div class="ttname"><a href="ERF__IndexDefines_8H.html#ada2d84e0dbbfb8d748defbf018748a70a83aec89012c35752597047c6ebdd983a">AdvType::Upwind_5th</a></div><div class="ttdeci">@ Upwind_5th</div></div>
<div class="ttc" id="aERF__IndexDefines_8H_html_ada2d84e0dbbfb8d748defbf018748a70a8c9475eb6bfde28e5bcb5230ac2e22a9"><div class="ttname"><a href="ERF__IndexDefines_8H.html#ada2d84e0dbbfb8d748defbf018748a70a8c9475eb6bfde28e5bcb5230ac2e22a9">AdvType::Centered_6th</a></div><div class="ttdeci">@ Centered_6th</div></div>
<div class="ttc" id="aERF__IndexDefines_8H_html_ada2d84e0dbbfb8d748defbf018748a70ada570868919bafaab1776170229a7049"><div class="ttname"><a href="ERF__IndexDefines_8H.html#ada2d84e0dbbfb8d748defbf018748a70ada570868919bafaab1776170229a7049">AdvType::Centered_2nd</a></div><div class="ttdeci">@ Centered_2nd</div></div>
<div class="ttc" id="aERF__Interpolation_8H_html"><div class="ttname"><a href="ERF__Interpolation_8H.html">ERF_Interpolation.H</a></div></div>
<div class="ttc" id="astructCENTERED2_html"><div class="ttname"><a href="structCENTERED2.html">CENTERED2</a></div><div class="ttdef"><b>Definition:</b> ERF_Interpolation_UPW.H:10</div></div>
<div class="ttc" id="astructCENTERED2_html_a2b0e42ec56f22f97bb65ed044f6f77a0"><div class="ttname"><a href="structCENTERED2.html#a2b0e42ec56f22f97bb65ed044f6f77a0">CENTERED2::InterpolateInY</a></div><div class="ttdeci">AMREX_GPU_DEVICE AMREX_FORCE_INLINE void InterpolateInY(const int &i, const int &j, const int &k, const int &qty_index, amrex::Real &val_lo, amrex::Real) const</div><div class="ttdef"><b>Definition:</b> ERF_Interpolation_UPW.H:36</div></div>
<div class="ttc" id="astructCENTERED2_html_a66de83f17cc1fa457cb89af60afb0041"><div class="ttname"><a href="structCENTERED2.html#a66de83f17cc1fa457cb89af60afb0041">CENTERED2::InterpolateInZ</a></div><div class="ttdeci">AMREX_GPU_DEVICE AMREX_FORCE_INLINE void InterpolateInZ(const int &i, const int &j, const int &k, const int &qty_index, amrex::Real &val_lo, amrex::Real) const</div><div class="ttdef"><b>Definition:</b> ERF_Interpolation_UPW.H:54</div></div>
<div class="ttc" id="astructCENTERED2_html_adb613a538275649aec33db89afd0e822"><div class="ttname"><a href="structCENTERED2.html#adb613a538275649aec33db89afd0e822">CENTERED2::InterpolateInX</a></div><div class="ttdeci">AMREX_GPU_DEVICE AMREX_FORCE_INLINE void InterpolateInX(const int &i, const int &j, const int &k, const int &qty_index, amrex::Real &val_lo, amrex::Real) const</div><div class="ttdef"><b>Definition:</b> ERF_Interpolation_UPW.H:18</div></div>
<div class="ttc" id="astructUPWIND3_html"><div class="ttname"><a href="structUPWIND3.html">UPWIND3</a></div><div class="ttdef"><b>Definition:</b> ERF_Interpolation_UPW.H:93</div></div>
<div class="ttc" id="astructUPWIND3_html_a420b9e841d549ce4f116205928f1743c"><div class="ttname"><a href="structUPWIND3.html#a420b9e841d549ce4f116205928f1743c">UPWIND3::InterpolateInZ</a></div><div class="ttdeci">AMREX_GPU_DEVICE AMREX_FORCE_INLINE void InterpolateInZ(const int &i, const int &j, const int &k, const int &qty_index, amrex::Real &val_lo, amrex::Real upw_lo) const</div><div class="ttdef"><b>Definition:</b> ERF_Interpolation_UPW.H:155</div></div>
<div class="ttc" id="astructUPWIND3_html_a424f3d6e525d4e195377c034deb8a345"><div class="ttname"><a href="structUPWIND3.html#a424f3d6e525d4e195377c034deb8a345">UPWIND3::SetUpwinding</a></div><div class="ttdeci">void SetUpwinding(amrex::Real upw_frac)</div><div class="ttdef"><b>Definition:</b> ERF_Interpolation_UPW.H:199</div></div>
<div class="ttc" id="astructUPWIND3_html_a6719e6c075ae163173c16f37b9e9c4b2"><div class="ttname"><a href="structUPWIND3.html#a6719e6c075ae163173c16f37b9e9c4b2">UPWIND3::InterpolateInY</a></div><div class="ttdeci">AMREX_GPU_DEVICE AMREX_FORCE_INLINE void InterpolateInY(const int &i, const int &j, const int &k, const int &qty_index, amrex::Real &val_lo, amrex::Real upw_lo) const</div><div class="ttdef"><b>Definition:</b> ERF_Interpolation_UPW.H:129</div></div>
<div class="ttc" id="astructUPWIND3_html_afd23b4297df2e3daf4a1b1e40cf40fe2"><div class="ttname"><a href="structUPWIND3.html#afd23b4297df2e3daf4a1b1e40cf40fe2">UPWIND3::InterpolateInX</a></div><div class="ttdeci">AMREX_GPU_DEVICE AMREX_FORCE_INLINE void InterpolateInX(const int &i, const int &j, const int &k, const int &qty_index, amrex::Real &val_lo, amrex::Real upw_lo) const</div><div class="ttdef"><b>Definition:</b> ERF_Interpolation_UPW.H:103</div></div>
</div><!-- fragment --></div><!-- contents -->
</div><!-- doc-content -->
<!-- start footer part -->
<div id="nav-path" class="navpath"><!-- id is needed for treeview function! -->
<ul>
<li class="navelem"><a class="el" href="dir_74389ed8173ad57b461b9d623a1f3867.html">Source</a></li><li class="navelem"><a class="el" href="dir_9c08f4b6b7cd879dc853cee646ab55d3.html">EB</a></li><li class="navelem"><a class="el" href="ERF__EBAdvectionSrcForScalars_8H.html">ERF_EBAdvectionSrcForScalars.H</a></li>
<li class="footer">Generated by <a href="https://www.doxygen.org/index.html"><img class="footer" src="doxygen.svg" width="104" height="31" alt="doxygen"/></a> 1.9.1 </li>
</ul>
</div>
</body>
</html>