-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathERF__AdvectionSrcForState_8cpp.html
816 lines (811 loc) · 82.9 KB
/
ERF__AdvectionSrcForState_8cpp.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
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
<!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/Advection/ERF_AdvectionSrcForState.cpp File Reference</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__AdvectionSrcForState_8cpp.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="summary">
<a href="#func-members">Functions</a> </div>
<div class="headertitle">
<div class="title">ERF_AdvectionSrcForState.cpp File Reference</div> </div>
</div><!--header-->
<div class="contents">
<div class="textblock"><code>#include <<a class="el" href="ERF__IndexDefines_8H_source.html">ERF_IndexDefines.H</a>></code><br />
<code>#include <<a class="el" href="ERF__TerrainMetrics_8H_source.html">ERF_TerrainMetrics.H</a>></code><br />
<code>#include <<a class="el" href="ERF__Advection_8H_source.html">ERF_Advection.H</a>></code><br />
<code>#include <<a class="el" href="ERF__AdvectionSrcForScalars_8H_source.html">ERF_AdvectionSrcForScalars.H</a>></code><br />
</div><div class="textblock"><div class="dynheader">
Include dependency graph for ERF_AdvectionSrcForState.cpp:</div>
<div class="dyncontent">
<div class="center"><img src="ERF__AdvectionSrcForState_8cpp__incl.png" border="0" usemap="#aSource_2Advection_2ERF__AdvectionSrcForState_8cpp" alt=""/></div>
<map name="aSource_2Advection_2ERF__AdvectionSrcForState_8cpp" id="aSource_2Advection_2ERF__AdvectionSrcForState_8cpp">
<area shape="rect" title=" " alt="" coords="1559,5,1753,47"/>
<area shape="rect" href="ERF__IndexDefines_8H.html" title=" " alt="" coords="1915,468,2069,495"/>
<area shape="rect" href="ERF__TerrainMetrics_8H.html" title=" " alt="" coords="1463,319,1623,345"/>
<area shape="rect" href="ERF__Advection_8H.html" title=" " alt="" coords="1721,95,1855,121"/>
<area shape="rect" href="ERF__AdvectionSrcForScalars_8H.html" title=" " alt="" coords="966,95,1186,121"/>
<area shape="rect" title=" " alt="" coords="2297,543,2420,569"/>
<area shape="rect" title=" " alt="" coords="1928,543,2056,569"/>
<area shape="rect" title=" " alt="" coords="1641,543,1724,569"/>
<area shape="rect" title=" " alt="" coords="933,468,1086,495"/>
<area shape="rect" title=" " alt="" coords="1913,393,2058,420"/>
<area shape="rect" title=" " alt="" coords="1925,169,2057,196"/>
<area shape="rect" title=" " alt="" coords="1513,169,1700,196"/>
<area shape="rect" href="ERF__DataStruct_8H.html" title=" " alt="" coords="706,319,843,345"/>
<area shape="rect" href="ERF__ABLMost_8H.html" title=" " alt="" coords="1775,169,1900,196"/>
<area shape="rect" title=" " alt="" coords="423,468,481,495"/>
<area shape="rect" title=" " alt="" coords="505,468,585,495"/>
<area shape="rect" title=" " alt="" coords="750,468,909,495"/>
<area shape="rect" title=" " alt="" coords="277,468,397,495"/>
<area shape="rect" title=" " alt="" coords="610,468,726,495"/>
<area shape="rect" href="ERF__Constants_8H.html" title=" " alt="" coords="2407,468,2540,495"/>
<area shape="rect" href="ERF__AdvStruct_8H.html" title=" " alt="" coords="705,393,836,420"/>
<area shape="rect" href="ERF__DiffStruct_8H.html" title=" " alt="" coords="271,393,401,420"/>
<area shape="rect" href="ERF__SpongeStruct_8H.html" title=" " alt="" coords="476,393,631,420"/>
<area shape="rect" href="ERF__TurbStruct_8H.html" title=" " alt="" coords="12,393,145,420"/>
<area shape="rect" href="ERF__TurbPertStruct_8H.html" title=" " alt="" coords="1215,393,1374,420"/>
<area shape="rect" href="ERF__MYNNStruct_8H.html" title=" " alt="" coords="5,468,152,495"/>
<area shape="rect" title=" " alt="" coords="1199,468,1367,495"/>
<area shape="rect" href="ERF__TileNoZ_8H.html" title=" " alt="" coords="1391,468,1510,495"/>
<area shape="rect" title=" " alt="" coords="1111,468,1175,495"/>
<area shape="rect" title=" " alt="" coords="2075,319,2231,345"/>
<area shape="rect" title=" " alt="" coords="1748,319,1897,345"/>
<area shape="rect" title=" " alt="" coords="1967,244,2132,271"/>
<area shape="rect" href="ERF__MOSTAverage_8H.html" title=" " alt="" coords="1581,244,1739,271"/>
<area shape="rect" href="ERF__MOSTStress_8H.html" title=" " alt="" coords="2409,319,2556,345"/>
<area shape="rect" href="ERF__PBLHeight_8H.html" title=" " alt="" coords="1398,393,1533,420"/>
<area shape="rect" href="ERF__MOSTRoughness_8H.html" title=" " alt="" coords="2641,393,2818,420"/>
<area shape="rect" href="ERF__Wstar_8H.html" title=" " alt="" coords="2511,393,2617,420"/>
<area shape="rect" href="ERF__Thetav_8H.html" title=" " alt="" coords="1533,468,1645,495"/>
<area shape="rect" href="ERF__Interpolation_8H.html" title=" " alt="" coords="703,169,854,196"/>
<area shape="rect" href="ERF__Interpolation__UPW_8H.html" title=" " alt="" coords="237,244,425,271"/>
<area shape="rect" href="ERF__Interpolation__WENO_8H.html" title=" " alt="" coords="449,244,647,271"/>
<area shape="rect" href="ERF__Interpolation__WENO__Z_8H.html" title=" " alt="" coords="671,244,886,271"/>
</map>
</div>
</div><table class="memberdecls">
<tr class="heading"><td colspan="2"><h2 class="groupheader"><a name="func-members"></a>
Functions</h2></td></tr>
<tr class="memitem:a69b1270b4619e0ecfd16ab625cba13e1"><td class="memItemLeft" align="right" valign="top">void </td><td class="memItemRight" valign="bottom"><a class="el" href="ERF__AdvectionSrcForState_8cpp.html#a69b1270b4619e0ecfd16ab625cba13e1">AdvectionSrcForRho</a> (const Box &bx, const Array4< Real > &advectionSrc, const Array4< const Real > &rho_u, const Array4< const Real > &rho_v, const Array4< const Real > &Omega, const Array4< Real > &avg_xmom, const Array4< Real > &avg_ymom, const Array4< Real > &avg_zmom, const Array4< const Real > &ax_arr, const Array4< const Real > &ay_arr, const Array4< const Real > &az_arr, const Array4< const Real > &detJ, const GpuArray< Real, AMREX_SPACEDIM > &cellSizeInv, const Array4< const Real > &mf_m, const Array4< const Real > &mf_u, const Array4< const Real > &mf_v, const GpuArray< const Array4< Real >, AMREX_SPACEDIM > &flx_arr, const bool fixed_rho)</td></tr>
<tr class="separator:a69b1270b4619e0ecfd16ab625cba13e1"><td class="memSeparator" colspan="2"> </td></tr>
<tr class="memitem:a244d2d95f771960db0f0ec982d3595f8"><td class="memItemLeft" align="right" valign="top">void </td><td class="memItemRight" valign="bottom"><a class="el" href="ERF__AdvectionSrcForState_8cpp.html#a244d2d95f771960db0f0ec982d3595f8">AdvectionSrcForScalars</a> (const Real &dt, const Box &bx, const int icomp, const int ncomp, const Array4< const Real > &avg_xmom, const Array4< const Real > &avg_ymom, const Array4< const Real > &avg_zmom, const Array4< const Real > &cur_cons, const Array4< const Real > &cell_prim, const Array4< Real > &advectionSrc, const bool &use_mono_adv, Real *max_s_ptr, Real *min_s_ptr, const Array4< const Real > &detJ, const GpuArray< Real, AMREX_SPACEDIM > &cellSizeInv, const Array4< const Real > &mf_m, const <a class="el" href="ERF__IndexDefines_8H.html#ada2d84e0dbbfb8d748defbf018748a70">AdvType</a> horiz_adv_type, const <a class="el" href="ERF__IndexDefines_8H.html#ada2d84e0dbbfb8d748defbf018748a70">AdvType</a> vert_adv_type, const Real horiz_upw_frac, const Real vert_upw_frac, const GpuArray< const Array4< Real >, AMREX_SPACEDIM > &flx_arr, const GpuArray< Array4< Real >, AMREX_SPACEDIM > &flx_tmp_arr, const Box &domain, const BCRec *bc_ptr_h)</td></tr>
<tr class="separator:a244d2d95f771960db0f0ec982d3595f8"><td class="memSeparator" colspan="2"> </td></tr>
</table>
<h2 class="groupheader">Function Documentation</h2>
<a id="a69b1270b4619e0ecfd16ab625cba13e1"></a>
<h2 class="memtitle"><span class="permalink"><a href="#a69b1270b4619e0ecfd16ab625cba13e1">◆ </a></span>AdvectionSrcForRho()</h2>
<div class="memitem">
<div class="memproto">
<table class="memname">
<tr>
<td class="memname">void AdvectionSrcForRho </td>
<td>(</td>
<td class="paramtype">const Box & </td>
<td class="paramname"><em>bx</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const Array4< Real > & </td>
<td class="paramname"><em>advectionSrc</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const Array4< const Real > & </td>
<td class="paramname"><em>rho_u</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const Array4< const Real > & </td>
<td class="paramname"><em>rho_v</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const Array4< const Real > & </td>
<td class="paramname"><em>Omega</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const Array4< Real > & </td>
<td class="paramname"><em>avg_xmom</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const Array4< Real > & </td>
<td class="paramname"><em>avg_ymom</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const Array4< Real > & </td>
<td class="paramname"><em>avg_zmom</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const Array4< const Real > & </td>
<td class="paramname"><em>ax_arr</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const Array4< const Real > & </td>
<td class="paramname"><em>ay_arr</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const Array4< const Real > & </td>
<td class="paramname"><em>az_arr</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const Array4< const Real > & </td>
<td class="paramname"><em>detJ</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const GpuArray< Real, AMREX_SPACEDIM > & </td>
<td class="paramname"><em>cellSizeInv</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const Array4< const Real > & </td>
<td class="paramname"><em>mf_m</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const Array4< const Real > & </td>
<td class="paramname"><em>mf_u</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const Array4< const Real > & </td>
<td class="paramname"><em>mf_v</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const GpuArray< const Array4< Real >, AMREX_SPACEDIM > & </td>
<td class="paramname"><em>flx_arr</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const bool </td>
<td class="paramname"><em>fixed_rho</em> </td>
</tr>
<tr>
<td></td>
<td>)</td>
<td></td><td></td>
</tr>
</table>
</div><div class="memdoc">
<p>Function for computing the advective tendency for the update equations for rho and (rho theta) This routine has explicit expressions for all cases (terrain or not) when the horizontal and vertical spatial orders are <= 2, and calls more specialized functions when either (or both) spatial order(s) is greater than 2.</p>
<dl class="params"><dt>Parameters</dt><dd>
<table class="params">
<tr><td class="paramdir">[in]</td><td class="paramname">bx</td><td>box over which the scalars are updated </td></tr>
<tr><td class="paramdir">[out]</td><td class="paramname">advectionSrc</td><td>tendency for the scalar update equation </td></tr>
<tr><td class="paramdir">[in]</td><td class="paramname">rho_u</td><td>x-component of momentum </td></tr>
<tr><td class="paramdir">[in]</td><td class="paramname">rho_v</td><td>y-component of momentum </td></tr>
<tr><td class="paramdir">[in]</td><td class="paramname">Omega</td><td>component of momentum normal to the z-coordinate surface </td></tr>
<tr><td class="paramdir">[out]</td><td class="paramname">avg_xmom</td><td>x-component of time-averaged momentum defined in this routine </td></tr>
<tr><td class="paramdir">[out]</td><td class="paramname">avg_ymom</td><td>y-component of time-averaged momentum defined in this routine </td></tr>
<tr><td class="paramdir">[out]</td><td class="paramname">avg_zmom</td><td>z-component of time-averaged momentum defined in this routine </td></tr>
<tr><td class="paramdir">[in]</td><td class="paramname">detJ</td><td>Jacobian of the metric transformation (= 1 if use_terrain is false) </td></tr>
<tr><td class="paramdir">[in]</td><td class="paramname">cellSizeInv</td><td>inverse of the mesh spacing </td></tr>
<tr><td class="paramdir">[in]</td><td class="paramname">mf_m</td><td>map factor at cell centers </td></tr>
<tr><td class="paramdir">[in]</td><td class="paramname">mf_u</td><td>map factor at x-faces </td></tr>
<tr><td class="paramdir">[in]</td><td class="paramname">mf_v</td><td>map factor at y-faces </td></tr>
</table>
</dd>
</dl>
<div class="fragment"><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>  BL_PROFILE_VAR(<span class="stringliteral">"AdvectionSrcForRho"</span>, <a class="code" href="ERF__AdvectionSrcForState_8cpp.html#a69b1270b4619e0ecfd16ab625cba13e1">AdvectionSrcForRho</a>);</div>
<div class="line"><a name="l00050"></a><span class="lineno"> 50</span>  <span class="keyword">auto</span> dxInv = cellSizeInv[0], dyInv = cellSizeInv[1], dzInv = cellSizeInv[2];</div>
<div class="line"><a name="l00051"></a><span class="lineno"> 51</span>  </div>
<div class="line"><a name="l00052"></a><span class="lineno"> 52</span>  <span class="keyword">const</span> Box xbx = surroundingNodes(bx,0);</div>
<div class="line"><a name="l00053"></a><span class="lineno"> 53</span>  <span class="keyword">const</span> Box ybx = surroundingNodes(bx,1);</div>
<div class="line"><a name="l00054"></a><span class="lineno"> 54</span>  <span class="keyword">const</span> Box zbx = surroundingNodes(bx,2);</div>
<div class="line"><a name="l00055"></a><span class="lineno"> 55</span>  </div>
<div class="line"><a name="l00056"></a><span class="lineno"> 56</span>  ParallelFor(xbx, [=] AMREX_GPU_DEVICE (<span class="keywordtype">int</span> i, <span class="keywordtype">int</span> j, <span class="keywordtype">int</span> k) noexcept</div>
<div class="line"><a name="l00057"></a><span class="lineno"> 57</span>  {</div>
<div class="line"><a name="l00058"></a><span class="lineno"> 58</span>  (flx_arr[0])(i,j,k,0) = ax_arr(i,j,k) * rho_u(i,j,k) / mf_u(i,j,0);</div>
<div class="line"><a name="l00059"></a><span class="lineno"> 59</span>  avg_xmom(i,j,k) = (flx_arr[0])(i,j,k,0);</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>  ParallelFor(ybx, [=] AMREX_GPU_DEVICE (<span class="keywordtype">int</span> i, <span class="keywordtype">int</span> j, <span class="keywordtype">int</span> k) noexcept</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>  (flx_arr[1])(i,j,k,0) = ay_arr(i,j,k) * rho_v(i,j,k) / mf_v(i,j,0);</div>
<div class="line"><a name="l00064"></a><span class="lineno"> 64</span>  avg_ymom(i,j,k) = (flx_arr[1])(i,j,k,0);</div>
<div class="line"><a name="l00065"></a><span class="lineno"> 65</span>  });</div>
<div class="line"><a name="l00066"></a><span class="lineno"> 66</span>  ParallelFor(zbx, [=] AMREX_GPU_DEVICE (<span class="keywordtype">int</span> i, <span class="keywordtype">int</span> j, <span class="keywordtype">int</span> k) noexcept</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>  Real mfsq = mf_m(i,j,0) * mf_m(i,j,0);</div>
<div class="line"><a name="l00069"></a><span class="lineno"> 69</span>  (flx_arr[2])(i,j,k,0) = az_arr(i,j,k) * Omega(i,j,k) / mfsq;</div>
<div class="line"><a name="l00070"></a><span class="lineno"> 70</span>  avg_zmom(i,j,k) = (flx_arr[2])(i,j,k,0);</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>  </div>
<div class="line"><a name="l00073"></a><span class="lineno"> 73</span>  <span class="keywordflow">if</span> (fixed_rho) {</div>
<div class="line"><a name="l00074"></a><span class="lineno"> 74</span>  ParallelFor(bx, [=] AMREX_GPU_DEVICE (<span class="keywordtype">int</span> i, <span class="keywordtype">int</span> j, <span class="keywordtype">int</span> k) noexcept</div>
<div class="line"><a name="l00075"></a><span class="lineno"> 75</span>  {</div>
<div class="line"><a name="l00076"></a><span class="lineno"> 76</span>  advectionSrc(i,j,k,0) = 0.0;</div>
<div class="line"><a name="l00077"></a><span class="lineno"> 77</span>  });</div>
<div class="line"><a name="l00078"></a><span class="lineno"> 78</span>  } <span class="keywordflow">else</span></div>
<div class="line"><a name="l00079"></a><span class="lineno"> 79</span>  {</div>
<div class="line"><a name="l00080"></a><span class="lineno"> 80</span>  ParallelFor(bx, [=] AMREX_GPU_DEVICE (<span class="keywordtype">int</span> i, <span class="keywordtype">int</span> j, <span class="keywordtype">int</span> k) noexcept</div>
<div class="line"><a name="l00081"></a><span class="lineno"> 81</span>  {</div>
<div class="line"><a name="l00082"></a><span class="lineno"> 82</span>  <span class="keywordflow">if</span> (detJ(i,j,k) > 0.) {</div>
<div class="line"><a name="l00083"></a><span class="lineno"> 83</span>  Real mfsq = mf_m(i,j,0) * mf_m(i,j,0);</div>
<div class="line"><a name="l00084"></a><span class="lineno"> 84</span>  advectionSrc(i,j,k,0) = - mfsq / detJ(i,j,k) * (</div>
<div class="line"><a name="l00085"></a><span class="lineno"> 85</span>  ( (flx_arr[0])(i+1,j,k,0) - (flx_arr[0])(i ,j,k,0) ) * dxInv +</div>
<div class="line"><a name="l00086"></a><span class="lineno"> 86</span>  ( (flx_arr[1])(i,j+1,k,0) - (flx_arr[1])(i,j ,k,0) ) * dyInv +</div>
<div class="line"><a name="l00087"></a><span class="lineno"> 87</span>  ( (flx_arr[2])(i,j,k+1,0) - (flx_arr[2])(i,j,k ,0) ) * dzInv );</div>
<div class="line"><a name="l00088"></a><span class="lineno"> 88</span>  } <span class="keywordflow">else</span> {</div>
<div class="line"><a name="l00089"></a><span class="lineno"> 89</span>  advectionSrc(i,j,k,0) = 0.;</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>  }</div>
<div class="line"><a name="l00093"></a><span class="lineno"> 93</span> }</div>
<div class="ttc" id="aERF__AdvectionSrcForState_8cpp_html_a69b1270b4619e0ecfd16ab625cba13e1"><div class="ttname"><a href="ERF__AdvectionSrcForState_8cpp.html#a69b1270b4619e0ecfd16ab625cba13e1">AdvectionSrcForRho</a></div><div class="ttdeci">void AdvectionSrcForRho(const Box &bx, const Array4< Real > &advectionSrc, const Array4< const Real > &rho_u, const Array4< const Real > &rho_v, const Array4< const Real > &Omega, const Array4< Real > &avg_xmom, const Array4< Real > &avg_ymom, const Array4< Real > &avg_zmom, const Array4< const Real > &ax_arr, const Array4< const Real > &ay_arr, const Array4< const Real > &az_arr, const Array4< const Real > &detJ, const GpuArray< Real, AMREX_SPACEDIM > &cellSizeInv, const Array4< const Real > &mf_m, const Array4< const Real > &mf_u, const Array4< const Real > &mf_v, const GpuArray< const Array4< Real >, AMREX_SPACEDIM > &flx_arr, const bool fixed_rho)</div><div class="ttdef"><b>Definition:</b> ERF_AdvectionSrcForState.cpp:30</div></div>
</div><!-- fragment -->
</div>
</div>
<a id="a244d2d95f771960db0f0ec982d3595f8"></a>
<h2 class="memtitle"><span class="permalink"><a href="#a244d2d95f771960db0f0ec982d3595f8">◆ </a></span>AdvectionSrcForScalars()</h2>
<div class="memitem">
<div class="memproto">
<table class="memname">
<tr>
<td class="memname">void AdvectionSrcForScalars </td>
<td>(</td>
<td class="paramtype">const Real & </td>
<td class="paramname"><em>dt</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const Box & </td>
<td class="paramname"><em>bx</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const int </td>
<td class="paramname"><em>icomp</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const int </td>
<td class="paramname"><em>ncomp</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const Array4< const Real > & </td>
<td class="paramname"><em>avg_xmom</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const Array4< const Real > & </td>
<td class="paramname"><em>avg_ymom</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const Array4< const Real > & </td>
<td class="paramname"><em>avg_zmom</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const Array4< const Real > & </td>
<td class="paramname"><em>cur_cons</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const Array4< const Real > & </td>
<td class="paramname"><em>cell_prim</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const Array4< Real > & </td>
<td class="paramname"><em>advectionSrc</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const bool & </td>
<td class="paramname"><em>use_mono_adv</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">Real * </td>
<td class="paramname"><em>max_s_ptr</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">Real * </td>
<td class="paramname"><em>min_s_ptr</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const Array4< const Real > & </td>
<td class="paramname"><em>detJ</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const GpuArray< Real, AMREX_SPACEDIM > & </td>
<td class="paramname"><em>cellSizeInv</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const Array4< const Real > & </td>
<td class="paramname"><em>mf_m</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const <a class="el" href="ERF__IndexDefines_8H.html#ada2d84e0dbbfb8d748defbf018748a70">AdvType</a> </td>
<td class="paramname"><em>horiz_adv_type</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const <a class="el" href="ERF__IndexDefines_8H.html#ada2d84e0dbbfb8d748defbf018748a70">AdvType</a> </td>
<td class="paramname"><em>vert_adv_type</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const Real </td>
<td class="paramname"><em>horiz_upw_frac</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const Real </td>
<td class="paramname"><em>vert_upw_frac</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const GpuArray< const Array4< Real >, AMREX_SPACEDIM > & </td>
<td class="paramname"><em>flx_arr</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const GpuArray< Array4< Real >, AMREX_SPACEDIM > & </td>
<td class="paramname"><em>flx_tmp_arr</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const Box & </td>
<td class="paramname"><em>domain</em>, </td>
</tr>
<tr>
<td class="paramkey"></td>
<td></td>
<td class="paramtype">const BCRec * </td>
<td class="paramname"><em>bc_ptr_h</em> </td>
</tr>
<tr>
<td></td>
<td>)</td>
<td></td><td></td>
</tr>
</table>
</div><div class="memdoc">
<p>Function for computing the advective tendency for the update equations for all scalars other than rho and (rho theta) This routine has explicit expressions for all cases (terrain or not) when the horizontal and vertical spatial orders are <= 2, and calls more specialized functions when either (or both) spatial order(s) is greater than 2.</p>
<dl class="params"><dt>Parameters</dt><dd>
<table class="params">
<tr><td class="paramdir">[in]</td><td class="paramname">bx</td><td>box over which the scalars are updated if no external boundary conditions </td></tr>
<tr><td class="paramdir">[in]</td><td class="paramname">icomp</td><td>component of first scalar to be updated </td></tr>
<tr><td class="paramdir">[in]</td><td class="paramname">ncomp</td><td>number of components to be updated </td></tr>
<tr><td class="paramdir">[in]</td><td class="paramname">avg_xmom</td><td>x-component of time-averaged momentum defined in this routine </td></tr>
<tr><td class="paramdir">[in]</td><td class="paramname">avg_ymom</td><td>y-component of time-averaged momentum defined in this routine </td></tr>
<tr><td class="paramdir">[in]</td><td class="paramname">avg_zmom</td><td>z-component of time-averaged momentum defined in this routine </td></tr>
<tr><td class="paramdir">[in]</td><td class="paramname">cell_prim</td><td>primitive form of scalar variables, here only potential temperature theta </td></tr>
<tr><td class="paramdir">[out]</td><td class="paramname">advectionSrc</td><td>tendency for the scalar update equation </td></tr>
<tr><td class="paramdir">[in]</td><td class="paramname">detJ</td><td>Jacobian of the metric transformation (= 1 if use_terrain is false) </td></tr>
<tr><td class="paramdir">[in]</td><td class="paramname">cellSizeInv</td><td>inverse of the mesh spacing </td></tr>
<tr><td class="paramdir">[in]</td><td class="paramname">mf_m</td><td>map factor at cell centers </td></tr>
<tr><td class="paramdir">[in]</td><td class="paramname">horiz_adv_type</td><td>advection scheme to be used in horiz. directions for dry scalars </td></tr>
<tr><td class="paramdir">[in]</td><td class="paramname">vert_adv_type</td><td>advection scheme to be used in vert. directions for dry scalars </td></tr>
<tr><td class="paramdir">[in]</td><td class="paramname">horiz_upw_frac</td><td>upwinding fraction to be used in horiz. directions for dry scalars (for Blended schemes only) </td></tr>
<tr><td class="paramdir">[in]</td><td class="paramname">vert_upw_frac</td><td>upwinding fraction to be used in vert. directions for dry scalars (for Blended schemes only) </td></tr>
</table>
</dd>
</dl>
<div class="fragment"><div class="line"><a name="l00143"></a><span class="lineno"> 143</span> {</div>
<div class="line"><a name="l00144"></a><span class="lineno"> 144</span>  BL_PROFILE_VAR(<span class="stringliteral">"AdvectionSrcForScalars"</span>, <a class="code" href="ERF__AdvectionSrcForState_8cpp.html#a244d2d95f771960db0f0ec982d3595f8">AdvectionSrcForScalars</a>);</div>
<div class="line"><a name="l00145"></a><span class="lineno"> 145</span>  <span class="keyword">auto</span> dxInv = cellSizeInv[0], dyInv = cellSizeInv[1], dzInv = cellSizeInv[2];</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="keyword">const</span> Box xbx = surroundingNodes(bx,0);</div>
<div class="line"><a name="l00148"></a><span class="lineno"> 148</span>  <span class="keyword">const</span> Box ybx = surroundingNodes(bx,1);</div>
<div class="line"><a name="l00149"></a><span class="lineno"> 149</span>  <span class="keyword">const</span> Box zbx = surroundingNodes(bx,2);</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="comment">// Open bc will be imposed upon all vars (we only access cons here for simplicity)</span></div>
<div class="line"><a name="l00152"></a><span class="lineno"> 152</span>  <span class="keyword">const</span> <span class="keywordtype">bool</span> xlo_open = (bc_ptr_h[<a class="code" href="namespaceBCVars.html#a00367571247ba920989360b94c156f4aaaf601346b01655b32acf66d07b856aec">BCVars::cons_bc</a>].lo(0) == <a class="code" href="namespaceERFBCType.html#a40fbdd5a379a92debb4bcaaccac6e01ca2b4e1a73027449d9af83a25ed67ae914">ERFBCType::open</a>);</div>
<div class="line"><a name="l00153"></a><span class="lineno"> 153</span>  <span class="keyword">const</span> <span class="keywordtype">bool</span> xhi_open = (bc_ptr_h[<a class="code" href="namespaceBCVars.html#a00367571247ba920989360b94c156f4aaaf601346b01655b32acf66d07b856aec">BCVars::cons_bc</a>].hi(0) == <a class="code" href="namespaceERFBCType.html#a40fbdd5a379a92debb4bcaaccac6e01ca2b4e1a73027449d9af83a25ed67ae914">ERFBCType::open</a>);</div>
<div class="line"><a name="l00154"></a><span class="lineno"> 154</span>  <span class="keyword">const</span> <span class="keywordtype">bool</span> ylo_open = (bc_ptr_h[<a class="code" href="namespaceBCVars.html#a00367571247ba920989360b94c156f4aaaf601346b01655b32acf66d07b856aec">BCVars::cons_bc</a>].lo(1) == <a class="code" href="namespaceERFBCType.html#a40fbdd5a379a92debb4bcaaccac6e01ca2b4e1a73027449d9af83a25ed67ae914">ERFBCType::open</a>);</div>
<div class="line"><a name="l00155"></a><span class="lineno"> 155</span>  <span class="keyword">const</span> <span class="keywordtype">bool</span> yhi_open = (bc_ptr_h[<a class="code" href="namespaceBCVars.html#a00367571247ba920989360b94c156f4aaaf601346b01655b32acf66d07b856aec">BCVars::cons_bc</a>].hi(1) == <a class="code" href="namespaceERFBCType.html#a40fbdd5a379a92debb4bcaaccac6e01ca2b4e1a73027449d9af83a25ed67ae914">ERFBCType::open</a>);</div>
<div class="line"><a name="l00156"></a><span class="lineno"> 156</span>  </div>
<div class="line"><a name="l00157"></a><span class="lineno"> 157</span>  <span class="comment">// Only advection operations in bndry normal direction with OPEN BC</span></div>
<div class="line"><a name="l00158"></a><span class="lineno"> 158</span>  Box bx_xlo, bx_xhi, bx_ylo, bx_yhi;</div>
<div class="line"><a name="l00159"></a><span class="lineno"> 159</span>  <span class="keywordflow">if</span> (xlo_open) {</div>
<div class="line"><a name="l00160"></a><span class="lineno"> 160</span>  <span class="keywordflow">if</span> ( bx.smallEnd(0) == domain.smallEnd(0)) { bx_xlo = makeSlab( bx,0,domain.smallEnd(0));}</div>
<div class="line"><a name="l00161"></a><span class="lineno"> 161</span>  }</div>
<div class="line"><a name="l00162"></a><span class="lineno"> 162</span>  <span class="keywordflow">if</span> (xhi_open) {</div>
<div class="line"><a name="l00163"></a><span class="lineno"> 163</span>  <span class="keywordflow">if</span> ( bx.bigEnd(0) == domain.bigEnd(0)) { bx_xhi = makeSlab( bx,0,domain.bigEnd(0) );}</div>
<div class="line"><a name="l00164"></a><span class="lineno"> 164</span>  }</div>
<div class="line"><a name="l00165"></a><span class="lineno"> 165</span>  <span class="keywordflow">if</span> (ylo_open) {</div>
<div class="line"><a name="l00166"></a><span class="lineno"> 166</span>  <span class="keywordflow">if</span> ( bx.smallEnd(1) == domain.smallEnd(1)) { bx_ylo = makeSlab( bx,1,domain.smallEnd(1));}</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>  <span class="keywordflow">if</span> (yhi_open) {</div>
<div class="line"><a name="l00169"></a><span class="lineno"> 169</span>  <span class="keywordflow">if</span> ( bx.bigEnd(1) == domain.bigEnd(1)) { bx_yhi = makeSlab( bx,1,domain.bigEnd(1) );}</div>
<div class="line"><a name="l00170"></a><span class="lineno"> 170</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="comment">// Inline with 2nd order for efficiency</span></div>
<div class="line"><a name="l00173"></a><span class="lineno"> 173</span>  <span class="comment">// NOTE: we don't need to weight avg_xmom, avg_ymom, avg_zmom with terrain metrics</span></div>
<div class="line"><a name="l00174"></a><span class="lineno"> 174</span>  <span class="comment">// (or with EB area fractions)</span></div>
<div class="line"><a name="l00175"></a><span class="lineno"> 175</span>  <span class="comment">// because that was done when they were constructed in AdvectionSrcForRhoAndTheta</span></div>
<div class="line"><a name="l00176"></a><span class="lineno"> 176</span>  <span class="keywordflow">if</span> (horiz_adv_type == <a class="code" href="ERF__IndexDefines_8H.html#ada2d84e0dbbfb8d748defbf018748a70ada570868919bafaab1776170229a7049">AdvType::Centered_2nd</a> && vert_adv_type == <a class="code" href="ERF__IndexDefines_8H.html#ada2d84e0dbbfb8d748defbf018748a70ada570868919bafaab1776170229a7049">AdvType::Centered_2nd</a>)</div>
<div class="line"><a name="l00177"></a><span class="lineno"> 177</span>  {</div>
<div class="line"><a name="l00178"></a><span class="lineno"> 178</span>  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="l00179"></a><span class="lineno"> 179</span>  {</div>
<div class="line"><a name="l00180"></a><span class="lineno"> 180</span>  <span class="keyword">const</span> <span class="keywordtype">int</span> cons_index = icomp + n;</div>
<div class="line"><a name="l00181"></a><span class="lineno"> 181</span>  <span class="keyword">const</span> <span class="keywordtype">int</span> prim_index = cons_index - 1;</div>
<div class="line"><a name="l00182"></a><span class="lineno"> 182</span>  <span class="keyword">const</span> Real prim_on_face = 0.5 * (cell_prim(i,j,k,prim_index) + cell_prim(i-1,j,k,prim_index));</div>
<div class="line"><a name="l00183"></a><span class="lineno"> 183</span>  (flx_arr[0])(i,j,k,cons_index) = avg_xmom(i,j,k) * prim_on_face;</div>
<div class="line"><a name="l00184"></a><span class="lineno"> 184</span>  });</div>
<div class="line"><a name="l00185"></a><span class="lineno"> 185</span>  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="l00186"></a><span class="lineno"> 186</span>  {</div>
<div class="line"><a name="l00187"></a><span class="lineno"> 187</span>  <span class="keyword">const</span> <span class="keywordtype">int</span> cons_index = icomp + n;</div>
<div class="line"><a name="l00188"></a><span class="lineno"> 188</span>  <span class="keyword">const</span> <span class="keywordtype">int</span> prim_index = cons_index - 1;</div>
<div class="line"><a name="l00189"></a><span class="lineno"> 189</span>  <span class="keyword">const</span> Real prim_on_face = 0.5 * (cell_prim(i,j,k,prim_index) + cell_prim(i,j-1,k,prim_index));</div>
<div class="line"><a name="l00190"></a><span class="lineno"> 190</span>  (flx_arr[1])(i,j,k,cons_index) = avg_ymom(i,j,k) * prim_on_face;</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>  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="l00193"></a><span class="lineno"> 193</span>  {</div>
<div class="line"><a name="l00194"></a><span class="lineno"> 194</span>  <span class="keyword">const</span> <span class="keywordtype">int</span> cons_index = icomp + n;</div>
<div class="line"><a name="l00195"></a><span class="lineno"> 195</span>  <span class="keyword">const</span> <span class="keywordtype">int</span> prim_index = cons_index - 1;</div>
<div class="line"><a name="l00196"></a><span class="lineno"> 196</span>  <span class="keyword">const</span> Real prim_on_face = 0.5 * (cell_prim(i,j,k,prim_index) + cell_prim(i,j,k-1,prim_index));</div>
<div class="line"><a name="l00197"></a><span class="lineno"> 197</span>  (flx_arr[2])(i,j,k,cons_index) = avg_zmom(i,j,k) * prim_on_face;</div>
<div class="line"><a name="l00198"></a><span class="lineno"> 198</span>  });</div>
<div class="line"><a name="l00199"></a><span class="lineno"> 199</span>  </div>
<div class="line"><a name="l00200"></a><span class="lineno"> 200</span>  <span class="comment">// Template higher order methods (horizontal first)</span></div>
<div class="line"><a name="l00201"></a><span class="lineno"> 201</span>  } <span class="keywordflow">else</span> {</div>
<div class="line"><a name="l00202"></a><span class="lineno"> 202</span>  <span class="keywordflow">switch</span>(horiz_adv_type) {</div>
<div class="line"><a name="l00203"></a><span class="lineno"> 203</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="l00204"></a><span class="lineno"> 204</span>  AdvectionSrcForScalarsVert<CENTERED2>(bx, ncomp, icomp, flx_arr, cell_prim,</div>
<div class="line"><a name="l00205"></a><span class="lineno"> 205</span>  avg_xmom, avg_ymom, avg_zmom,</div>
<div class="line"><a name="l00206"></a><span class="lineno"> 206</span>  horiz_upw_frac, vert_upw_frac, vert_adv_type);</div>
<div class="line"><a name="l00207"></a><span class="lineno"> 207</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00208"></a><span class="lineno"> 208</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="l00209"></a><span class="lineno"> 209</span>  AdvectionSrcForScalarsVert<UPWIND3>(bx, ncomp, icomp, flx_arr, cell_prim,</div>
<div class="line"><a name="l00210"></a><span class="lineno"> 210</span>  avg_xmom, avg_ymom, avg_zmom,</div>
<div class="line"><a name="l00211"></a><span class="lineno"> 211</span>  horiz_upw_frac, vert_upw_frac, vert_adv_type);</div>
<div class="line"><a name="l00212"></a><span class="lineno"> 212</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00213"></a><span class="lineno"> 213</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="l00214"></a><span class="lineno"> 214</span>  AdvectionSrcForScalarsVert<CENTERED4>(bx, ncomp, icomp, flx_arr, cell_prim,</div>
<div class="line"><a name="l00215"></a><span class="lineno"> 215</span>  avg_xmom, avg_ymom, avg_zmom,</div>
<div class="line"><a name="l00216"></a><span class="lineno"> 216</span>  horiz_upw_frac, vert_upw_frac, vert_adv_type);</div>
<div class="line"><a name="l00217"></a><span class="lineno"> 217</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00218"></a><span class="lineno"> 218</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="l00219"></a><span class="lineno"> 219</span>  AdvectionSrcForScalarsVert<UPWIND5>(bx, ncomp, icomp, flx_arr, cell_prim,</div>
<div class="line"><a name="l00220"></a><span class="lineno"> 220</span>  avg_xmom, avg_ymom, avg_zmom,</div>
<div class="line"><a name="l00221"></a><span class="lineno"> 221</span>  horiz_upw_frac, vert_upw_frac, vert_adv_type);</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#ada2d84e0dbbfb8d748defbf018748a70a8c9475eb6bfde28e5bcb5230ac2e22a9">AdvType::Centered_6th</a>:</div>
<div class="line"><a name="l00224"></a><span class="lineno"> 224</span>  AdvectionSrcForScalarsVert<CENTERED6>(bx, ncomp, icomp, flx_arr, cell_prim,</div>
<div class="line"><a name="l00225"></a><span class="lineno"> 225</span>  avg_xmom, avg_ymom, avg_zmom,</div>
<div class="line"><a name="l00226"></a><span class="lineno"> 226</span>  horiz_upw_frac, vert_upw_frac, vert_adv_type);</div>
<div class="line"><a name="l00227"></a><span class="lineno"> 227</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00228"></a><span class="lineno"> 228</span>  <span class="keywordflow">case</span> <a class="code" href="ERF__IndexDefines_8H.html#ada2d84e0dbbfb8d748defbf018748a70ada35597e5f00fa867fe6e39d651a2195">AdvType::Weno_3</a>:</div>
<div class="line"><a name="l00229"></a><span class="lineno"> 229</span>  AdvectionSrcForScalarsWrapper<WENO3,WENO3>(bx, ncomp, icomp, flx_arr, cell_prim,</div>
<div class="line"><a name="l00230"></a><span class="lineno"> 230</span>  avg_xmom, avg_ymom, avg_zmom,</div>
<div class="line"><a name="l00231"></a><span class="lineno"> 231</span>  horiz_upw_frac, vert_upw_frac);</div>
<div class="line"><a name="l00232"></a><span class="lineno"> 232</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00233"></a><span class="lineno"> 233</span>  <span class="keywordflow">case</span> <a class="code" href="ERF__IndexDefines_8H.html#ada2d84e0dbbfb8d748defbf018748a70af2429e8de544f3ebd7c9668e3dbebaa8">AdvType::Weno_5</a>:</div>
<div class="line"><a name="l00234"></a><span class="lineno"> 234</span>  AdvectionSrcForScalarsWrapper<WENO5,WENO5>(bx, ncomp, icomp, flx_arr, cell_prim,</div>
<div class="line"><a name="l00235"></a><span class="lineno"> 235</span>  avg_xmom, avg_ymom, avg_zmom,</div>
<div class="line"><a name="l00236"></a><span class="lineno"> 236</span>  horiz_upw_frac, vert_upw_frac);</div>
<div class="line"><a name="l00237"></a><span class="lineno"> 237</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00238"></a><span class="lineno"> 238</span>  <span class="keywordflow">case</span> <a class="code" href="ERF__IndexDefines_8H.html#ada2d84e0dbbfb8d748defbf018748a70a3095d588151e351b421893f62411d15e">AdvType::Weno_7</a>:</div>
<div class="line"><a name="l00239"></a><span class="lineno"> 239</span>  AdvectionSrcForScalarsWrapper<WENO7,WENO7>(bx, ncomp, icomp, 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>  horiz_upw_frac, vert_upw_frac);</div>
<div class="line"><a name="l00242"></a><span class="lineno"> 242</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00243"></a><span class="lineno"> 243</span>  <span class="keywordflow">case</span> <a class="code" href="ERF__IndexDefines_8H.html#ada2d84e0dbbfb8d748defbf018748a70a9c29b64a7627155706a1cd93682b4351">AdvType::Weno_3Z</a>:</div>
<div class="line"><a name="l00244"></a><span class="lineno"> 244</span>  AdvectionSrcForScalarsWrapper<WENO_Z3,WENO_Z3>(bx, ncomp, icomp, flx_arr, cell_prim,</div>
<div class="line"><a name="l00245"></a><span class="lineno"> 245</span>  avg_xmom, avg_ymom, avg_zmom,</div>
<div class="line"><a name="l00246"></a><span class="lineno"> 246</span>  horiz_upw_frac, vert_upw_frac);</div>
<div class="line"><a name="l00247"></a><span class="lineno"> 247</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00248"></a><span class="lineno"> 248</span>  <span class="keywordflow">case</span> <a class="code" href="ERF__IndexDefines_8H.html#ada2d84e0dbbfb8d748defbf018748a70af8b942d20ed08b25df71ea1c87244210">AdvType::Weno_3MZQ</a>:</div>
<div class="line"><a name="l00249"></a><span class="lineno"> 249</span>  AdvectionSrcForScalarsWrapper<WENO_MZQ3,WENO_MZQ3>(bx, ncomp, icomp, flx_arr, cell_prim,</div>
<div class="line"><a name="l00250"></a><span class="lineno"> 250</span>  avg_xmom, avg_ymom, avg_zmom,</div>
<div class="line"><a name="l00251"></a><span class="lineno"> 251</span>  horiz_upw_frac, vert_upw_frac);</div>
<div class="line"><a name="l00252"></a><span class="lineno"> 252</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00253"></a><span class="lineno"> 253</span>  <span class="keywordflow">case</span> <a class="code" href="ERF__IndexDefines_8H.html#ada2d84e0dbbfb8d748defbf018748a70aeb75f5943c7837da80e2048f2b344efc">AdvType::Weno_5Z</a>:</div>
<div class="line"><a name="l00254"></a><span class="lineno"> 254</span>  AdvectionSrcForScalarsWrapper<WENO_Z5,WENO_Z5>(bx, ncomp, icomp, flx_arr, cell_prim,</div>
<div class="line"><a name="l00255"></a><span class="lineno"> 255</span>  avg_xmom, avg_ymom, avg_zmom,</div>
<div class="line"><a name="l00256"></a><span class="lineno"> 256</span>  horiz_upw_frac, vert_upw_frac);</div>
<div class="line"><a name="l00257"></a><span class="lineno"> 257</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00258"></a><span class="lineno"> 258</span>  <span class="keywordflow">case</span> <a class="code" href="ERF__IndexDefines_8H.html#ada2d84e0dbbfb8d748defbf018748a70a6edd20e2cefe2061bec0e1ebbf5c8573">AdvType::Weno_7Z</a>:</div>
<div class="line"><a name="l00259"></a><span class="lineno"> 259</span>  AdvectionSrcForScalarsWrapper<WENO_Z7,WENO_Z7>(bx, ncomp, icomp, flx_arr, cell_prim,</div>
<div class="line"><a name="l00260"></a><span class="lineno"> 260</span>  avg_xmom, avg_ymom, avg_zmom,</div>
<div class="line"><a name="l00261"></a><span class="lineno"> 261</span>  horiz_upw_frac, vert_upw_frac);</div>
<div class="line"><a name="l00262"></a><span class="lineno"> 262</span>  <span class="keywordflow">break</span>;</div>
<div class="line"><a name="l00263"></a><span class="lineno"> 263</span>  <span class="keywordflow">default</span>:</div>
<div class="line"><a name="l00264"></a><span class="lineno"> 264</span>  AMREX_ASSERT_WITH_MESSAGE(<span class="keyword">false</span>, <span class="stringliteral">"Unknown advection scheme!"</span>);</div>
<div class="line"><a name="l00265"></a><span class="lineno"> 265</span>  }</div>
<div class="line"><a name="l00266"></a><span class="lineno"> 266</span>  }</div>
<div class="line"><a name="l00267"></a><span class="lineno"> 267</span>  </div>
<div class="line"><a name="l00268"></a><span class="lineno"> 268</span>  <span class="comment">/* =======================================================================</span></div>
<div class="line"><a name="l00269"></a><span class="lineno"> 269</span> <span class="comment"> Monotonicity preserving order reduction for scalars (0-th upwind).</span></div>
<div class="line"><a name="l00270"></a><span class="lineno"> 270</span> <span class="comment"></span> </div>
<div class="line"><a name="l00271"></a><span class="lineno"> 271</span> <span class="comment"> NOTE:</span></div>
<div class="line"><a name="l00272"></a><span class="lineno"> 272</span> <span class="comment"> The following order reduction operations for scalar advection have been</span></div>
<div class="line"><a name="l00273"></a><span class="lineno"> 273</span> <span class="comment"> adapted from the PINACLES code developed at PPNL by K. Pressel et al.;</span></div>
<div class="line"><a name="l00274"></a><span class="lineno"> 274</span> <span class="comment"> see https://github.com/pnnl/pinacles.git and This version utilizes global</span></div>
<div class="line"><a name="l00275"></a><span class="lineno"> 275</span> <span class="comment"> min/max values for simplicity rather than the local average around each</span></div>
<div class="line"><a name="l00276"></a><span class="lineno"> 276</span> <span class="comment"> cell. Motivation for this modification is the compressible dycore in ERF</span></div>
<div class="line"><a name="l00277"></a><span class="lineno"> 277</span> <span class="comment"> as opposed to incompressible/analestic in PINACLES.</span></div>
<div class="line"><a name="l00278"></a><span class="lineno"> 278</span> <span class="comment"> ======================================================================= */</span></div>
<div class="line"><a name="l00279"></a><span class="lineno"> 279</span>  <span class="keywordflow">if</span> (use_mono_adv) {</div>
<div class="line"><a name="l00280"></a><span class="lineno"> 280</span>  <span class="comment">// Copy flux data to flx_arr to avoid race condition on GPU</span></div>
<div class="line"><a name="l00281"></a><span class="lineno"> 281</span>  ParallelFor(bx, 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="l00282"></a><span class="lineno"> 282</span>  {</div>
<div class="line"><a name="l00283"></a><span class="lineno"> 283</span>  <span class="keyword">const</span> <span class="keywordtype">int</span> cons_index = icomp + n;</div>
<div class="line"><a name="l00284"></a><span class="lineno"> 284</span>  (flx_tmp_arr[0])(i,j,k,cons_index) = (flx_arr[0])(i,j,k,cons_index);</div>
<div class="line"><a name="l00285"></a><span class="lineno"> 285</span>  (flx_tmp_arr[1])(i,j,k,cons_index) = (flx_arr[1])(i,j,k,cons_index);</div>
<div class="line"><a name="l00286"></a><span class="lineno"> 286</span>  (flx_tmp_arr[2])(i,j,k,cons_index) = (flx_arr[2])(i,j,k,cons_index);</div>
<div class="line"><a name="l00287"></a><span class="lineno"> 287</span>  });</div>
<div class="line"><a name="l00288"></a><span class="lineno"> 288</span>  </div>
<div class="line"><a name="l00289"></a><span class="lineno"> 289</span>  <span class="comment">// Mono limiting</span></div>
<div class="line"><a name="l00290"></a><span class="lineno"> 290</span>  ParallelFor(bx, 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="l00291"></a><span class="lineno"> 291</span>  {</div>
<div class="line"><a name="l00292"></a><span class="lineno"> 292</span>  <span class="keyword">const</span> <span class="keywordtype">int</span> cons_index = icomp + n;</div>
<div class="line"><a name="l00293"></a><span class="lineno"> 293</span>  <span class="keyword">const</span> <span class="keywordtype">int</span> prim_index = cons_index - 1;</div>
<div class="line"><a name="l00294"></a><span class="lineno"> 294</span>  </div>
<div class="line"><a name="l00295"></a><span class="lineno"> 295</span>  Real max_val = max_s_ptr[cons_index];</div>
<div class="line"><a name="l00296"></a><span class="lineno"> 296</span>  Real min_val = min_s_ptr[cons_index];</div>
<div class="line"><a name="l00297"></a><span class="lineno"> 297</span>  </div>
<div class="line"><a name="l00298"></a><span class="lineno"> 298</span>  Real invdetJ = (detJ(i,j,k) > 0.) ? 1. / detJ(i,j,k) : 1.;</div>
<div class="line"><a name="l00299"></a><span class="lineno"> 299</span>  Real mfsq = mf_m(i,j,0) * mf_m(i,j,0);</div>
<div class="line"><a name="l00300"></a><span class="lineno"> 300</span>  </div>
<div class="line"><a name="l00301"></a><span class="lineno"> 301</span>  Real RHS = - invdetJ * mfsq * (</div>
<div class="line"><a name="l00302"></a><span class="lineno"> 302</span>  ( (flx_arr[0])(i+1,j ,k ,cons_index) - (flx_arr[0])(i,j,k,cons_index) ) * dxInv +</div>
<div class="line"><a name="l00303"></a><span class="lineno"> 303</span>  ( (flx_arr[1])(i ,j+1,k ,cons_index) - (flx_arr[1])(i,j,k,cons_index) ) * dyInv +</div>
<div class="line"><a name="l00304"></a><span class="lineno"> 304</span>  ( (flx_arr[2])(i ,j ,k+1,cons_index) - (flx_arr[2])(i,j,k,cons_index) ) * dzInv );</div>
<div class="line"><a name="l00305"></a><span class="lineno"> 305</span>  </div>
<div class="line"><a name="l00306"></a><span class="lineno"> 306</span>  <span class="comment">// NOTE: This forward prediction uses the `cur_cons` as opposed to `old_cons` since</span></div>
<div class="line"><a name="l00307"></a><span class="lineno"> 307</span>  <span class="comment">// source terms may cause increase/decrease in a variable each RK stage. We</span></div>
<div class="line"><a name="l00308"></a><span class="lineno"> 308</span>  <span class="comment">// want to ensure the advection operator does not induce over/under-shoot</span></div>
<div class="line"><a name="l00309"></a><span class="lineno"> 309</span>  <span class="comment">// from the current state. If the `old_cons` is used and significant forcing is</span></div>
<div class="line"><a name="l00310"></a><span class="lineno"> 310</span>  <span class="comment">// present, we could trip an order reduction just due to the source terms.</span></div>
<div class="line"><a name="l00311"></a><span class="lineno"> 311</span>  Real tmp_upd = cur_cons(i,j,k,cons_index) + RHS*dt;</div>
<div class="line"><a name="l00312"></a><span class="lineno"> 312</span>  <span class="keywordflow">if</span> (tmp_upd<min_val || tmp_upd>max_val) {</div>
<div class="line"><a name="l00313"></a><span class="lineno"> 313</span>  <span class="comment">// HI</span></div>
<div class="line"><a name="l00314"></a><span class="lineno"> 314</span>  <span class="keywordflow">if</span> (avg_xmom(i+1,j,k)>0.0) {</div>
<div class="line"><a name="l00315"></a><span class="lineno"> 315</span>  (flx_tmp_arr[0])(i+1,j,k,cons_index) = avg_xmom(i+1,j,k) * cell_prim(i ,j,k,prim_index);</div>
<div class="line"><a name="l00316"></a><span class="lineno"> 316</span>  } <span class="keywordflow">else</span> {</div>
<div class="line"><a name="l00317"></a><span class="lineno"> 317</span>  (flx_tmp_arr[0])(i+1,j,k,cons_index) = avg_xmom(i+1,j,k) * cell_prim(i+1,j,k,prim_index);</div>
<div class="line"><a name="l00318"></a><span class="lineno"> 318</span>  }</div>
<div class="line"><a name="l00319"></a><span class="lineno"> 319</span>  <span class="keywordflow">if</span> (avg_ymom(i,j+1,k)>0.0) {</div>
<div class="line"><a name="l00320"></a><span class="lineno"> 320</span>  (flx_tmp_arr[1])(i,j+1,k,cons_index) = avg_ymom(i,j+1,k) * cell_prim(i,j ,k,prim_index);</div>
<div class="line"><a name="l00321"></a><span class="lineno"> 321</span>  } <span class="keywordflow">else</span> {</div>
<div class="line"><a name="l00322"></a><span class="lineno"> 322</span>  (flx_tmp_arr[1])(i,j+1,k,cons_index) = avg_ymom(i,j+1,k) * cell_prim(i,j+1,k,prim_index);</div>
<div class="line"><a name="l00323"></a><span class="lineno"> 323</span>  }</div>
<div class="line"><a name="l00324"></a><span class="lineno"> 324</span>  <span class="keywordflow">if</span> (avg_zmom(i,j,k+1)>0.0) {</div>
<div class="line"><a name="l00325"></a><span class="lineno"> 325</span>  (flx_tmp_arr[2])(i,j,k+1,cons_index) = avg_zmom(i,j,k+1) * cell_prim(i,j,k ,prim_index);</div>
<div class="line"><a name="l00326"></a><span class="lineno"> 326</span>  } <span class="keywordflow">else</span> {</div>
<div class="line"><a name="l00327"></a><span class="lineno"> 327</span>  (flx_tmp_arr[2])(i,j,k+1,cons_index) = avg_zmom(i,j,k+1) * cell_prim(i,j,k+1,prim_index);</div>
<div class="line"><a name="l00328"></a><span class="lineno"> 328</span>  }</div>
<div class="line"><a name="l00329"></a><span class="lineno"> 329</span>  </div>
<div class="line"><a name="l00330"></a><span class="lineno"> 330</span>  <span class="comment">// LO</span></div>
<div class="line"><a name="l00331"></a><span class="lineno"> 331</span>  <span class="keywordflow">if</span> (avg_xmom(i,j,k)>0.0) {</div>
<div class="line"><a name="l00332"></a><span class="lineno"> 332</span>  (flx_tmp_arr[0])(i,j,k,cons_index) = avg_xmom(i,j,k) * cell_prim(i-1,j,k,prim_index);</div>
<div class="line"><a name="l00333"></a><span class="lineno"> 333</span>  } <span class="keywordflow">else</span> {</div>
<div class="line"><a name="l00334"></a><span class="lineno"> 334</span>  (flx_tmp_arr[0])(i,j,k,cons_index) = avg_xmom(i,j,k) * cell_prim(i ,j,k,prim_index);</div>
<div class="line"><a name="l00335"></a><span class="lineno"> 335</span>  }</div>
<div class="line"><a name="l00336"></a><span class="lineno"> 336</span>  <span class="keywordflow">if</span> (avg_ymom(i,j,k)>0.0) {</div>
<div class="line"><a name="l00337"></a><span class="lineno"> 337</span>  (flx_tmp_arr[1])(i,j,k,cons_index) = avg_ymom(i,j,k) * cell_prim(i,j-1,k,prim_index);</div>
<div class="line"><a name="l00338"></a><span class="lineno"> 338</span>  } <span class="keywordflow">else</span> {</div>
<div class="line"><a name="l00339"></a><span class="lineno"> 339</span>  (flx_tmp_arr[1])(i,j,k,cons_index) = avg_ymom(i,j,k) * cell_prim(i,j ,k,prim_index);</div>
<div class="line"><a name="l00340"></a><span class="lineno"> 340</span>  }</div>
<div class="line"><a name="l00341"></a><span class="lineno"> 341</span>  <span class="keywordflow">if</span> (avg_zmom(i,j,k)>0.0) {</div>
<div class="line"><a name="l00342"></a><span class="lineno"> 342</span>  (flx_tmp_arr[2])(i,j,k,cons_index) = avg_zmom(i,j,k) * cell_prim(i,j,k-1,prim_index);</div>
<div class="line"><a name="l00343"></a><span class="lineno"> 343</span>  } <span class="keywordflow">else</span> {</div>
<div class="line"><a name="l00344"></a><span class="lineno"> 344</span>  (flx_tmp_arr[2])(i,j,k,cons_index) = avg_zmom(i,j,k) * cell_prim(i,j,k ,prim_index);</div>
<div class="line"><a name="l00345"></a><span class="lineno"> 345</span>  }</div>
<div class="line"><a name="l00346"></a><span class="lineno"> 346</span>  }</div>
<div class="line"><a name="l00347"></a><span class="lineno"> 347</span>  });</div>
<div class="line"><a name="l00348"></a><span class="lineno"> 348</span>  </div>
<div class="line"><a name="l00349"></a><span class="lineno"> 349</span>  <span class="comment">// Copy back to flx_arr to avoid race condition on GPU</span></div>
<div class="line"><a name="l00350"></a><span class="lineno"> 350</span>  ParallelFor(bx, 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="l00351"></a><span class="lineno"> 351</span>  {</div>
<div class="line"><a name="l00352"></a><span class="lineno"> 352</span>  <span class="keyword">const</span> <span class="keywordtype">int</span> cons_index = icomp + n;</div>
<div class="line"><a name="l00353"></a><span class="lineno"> 353</span>  (flx_arr[0])(i,j,k,cons_index) = (flx_tmp_arr[0])(i,j,k,cons_index);</div>
<div class="line"><a name="l00354"></a><span class="lineno"> 354</span>  (flx_arr[1])(i,j,k,cons_index) = (flx_tmp_arr[1])(i,j,k,cons_index);</div>
<div class="line"><a name="l00355"></a><span class="lineno"> 355</span>  (flx_arr[2])(i,j,k,cons_index) = (flx_tmp_arr[2])(i,j,k,cons_index);</div>
<div class="line"><a name="l00356"></a><span class="lineno"> 356</span>  });</div>
<div class="line"><a name="l00357"></a><span class="lineno"> 357</span>  }</div>
<div class="line"><a name="l00358"></a><span class="lineno"> 358</span>  </div>
<div class="line"><a name="l00359"></a><span class="lineno"> 359</span>  ParallelFor(bx, 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="l00360"></a><span class="lineno"> 360</span>  {</div>
<div class="line"><a name="l00361"></a><span class="lineno"> 361</span>  <span class="keyword">const</span> <span class="keywordtype">int</span> cons_index = icomp + n;</div>
<div class="line"><a name="l00362"></a><span class="lineno"> 362</span>  <span class="keywordflow">if</span> (detJ(i,j,k) > 0.)</div>
<div class="line"><a name="l00363"></a><span class="lineno"> 363</span>  {</div>
<div class="line"><a name="l00364"></a><span class="lineno"> 364</span>  Real invdetJ = 1.0 / detJ(i,j,k);</div>
<div class="line"><a name="l00365"></a><span class="lineno"> 365</span>  Real mfsq = mf_m(i,j,0) * mf_m(i,j,0);</div>
<div class="line"><a name="l00366"></a><span class="lineno"> 366</span>  </div>
<div class="line"><a name="l00367"></a><span class="lineno"> 367</span>  advectionSrc(i,j,k,cons_index) = - invdetJ * mfsq * (</div>
<div class="line"><a name="l00368"></a><span class="lineno"> 368</span>  ( (flx_arr[0])(i+1,j,k,cons_index) - (flx_arr[0])(i ,j,k,cons_index) ) * dxInv +</div>
<div class="line"><a name="l00369"></a><span class="lineno"> 369</span>  ( (flx_arr[1])(i,j+1,k,cons_index) - (flx_arr[1])(i,j ,k,cons_index) ) * dyInv +</div>
<div class="line"><a name="l00370"></a><span class="lineno"> 370</span>  ( (flx_arr[2])(i,j,k+1,cons_index) - (flx_arr[2])(i,j,k ,cons_index) ) * dzInv );</div>
<div class="line"><a name="l00371"></a><span class="lineno"> 371</span>  } <span class="keywordflow">else</span> {</div>
<div class="line"><a name="l00372"></a><span class="lineno"> 372</span>  advectionSrc(i,j,k,cons_index) = 0.;</div>
<div class="line"><a name="l00373"></a><span class="lineno"> 373</span>  }</div>
<div class="line"><a name="l00374"></a><span class="lineno"> 374</span>  });</div>
<div class="line"><a name="l00375"></a><span class="lineno"> 375</span>  </div>
<div class="line"><a name="l00376"></a><span class="lineno"> 376</span>  <span class="comment">// Special advection operator for open BC (bndry tangent operations)</span></div>
<div class="line"><a name="l00377"></a><span class="lineno"> 377</span>  <span class="keywordflow">if</span> (xlo_open) {</div>
<div class="line"><a name="l00378"></a><span class="lineno"> 378</span>  <span class="keywordtype">bool</span> do_lo = <span class="keyword">true</span>;</div>
<div class="line"><a name="l00379"></a><span class="lineno"> 379</span>  <a class="code" href="ERF__Advection_8H.html#abd84d0909db1f87c22b5a92aae03e919">AdvectionSrcForOpenBC_Tangent_Cons</a>(bx_xlo, 0, icomp, ncomp, advectionSrc, cell_prim,</div>
<div class="line"><a name="l00380"></a><span class="lineno"> 380</span>  avg_xmom, avg_ymom, avg_zmom,</div>
<div class="line"><a name="l00381"></a><span class="lineno"> 381</span>  detJ, cellSizeInv, do_lo);</div>
<div class="line"><a name="l00382"></a><span class="lineno"> 382</span>  }</div>
<div class="line"><a name="l00383"></a><span class="lineno"> 383</span>  <span class="keywordflow">if</span> (xhi_open) {</div>
<div class="line"><a name="l00384"></a><span class="lineno"> 384</span>  <a class="code" href="ERF__Advection_8H.html#abd84d0909db1f87c22b5a92aae03e919">AdvectionSrcForOpenBC_Tangent_Cons</a>(bx_xhi, 0, icomp, ncomp, advectionSrc, cell_prim,</div>
<div class="line"><a name="l00385"></a><span class="lineno"> 385</span>  avg_xmom, avg_ymom, avg_zmom,</div>
<div class="line"><a name="l00386"></a><span class="lineno"> 386</span>  detJ, cellSizeInv);</div>
<div class="line"><a name="l00387"></a><span class="lineno"> 387</span>  }</div>
<div class="line"><a name="l00388"></a><span class="lineno"> 388</span>  <span class="keywordflow">if</span> (ylo_open) {</div>
<div class="line"><a name="l00389"></a><span class="lineno"> 389</span>  <span class="keywordtype">bool</span> do_lo = <span class="keyword">true</span>;</div>
<div class="line"><a name="l00390"></a><span class="lineno"> 390</span>  <a class="code" href="ERF__Advection_8H.html#abd84d0909db1f87c22b5a92aae03e919">AdvectionSrcForOpenBC_Tangent_Cons</a>(bx_ylo, 1, icomp, ncomp, advectionSrc, cell_prim,</div>
<div class="line"><a name="l00391"></a><span class="lineno"> 391</span>  avg_xmom, avg_ymom, avg_zmom,</div>
<div class="line"><a name="l00392"></a><span class="lineno"> 392</span>  detJ, cellSizeInv, do_lo);</div>
<div class="line"><a name="l00393"></a><span class="lineno"> 393</span>  }</div>
<div class="line"><a name="l00394"></a><span class="lineno"> 394</span>  <span class="keywordflow">if</span> (yhi_open) {</div>
<div class="line"><a name="l00395"></a><span class="lineno"> 395</span>  <a class="code" href="ERF__Advection_8H.html#abd84d0909db1f87c22b5a92aae03e919">AdvectionSrcForOpenBC_Tangent_Cons</a>(bx_yhi, 1, icomp, ncomp, advectionSrc, cell_prim,</div>
<div class="line"><a name="l00396"></a><span class="lineno"> 396</span>  avg_xmom, avg_ymom, avg_zmom,</div>
<div class="line"><a name="l00397"></a><span class="lineno"> 397</span>  detJ, cellSizeInv);</div>
<div class="line"><a name="l00398"></a><span class="lineno"> 398</span>  }</div>
<div class="line"><a name="l00399"></a><span class="lineno"> 399</span> }</div>
<div class="ttc" id="aERF__AdvectionSrcForState_8cpp_html_a244d2d95f771960db0f0ec982d3595f8"><div class="ttname"><a href="ERF__AdvectionSrcForState_8cpp.html#a244d2d95f771960db0f0ec982d3595f8">AdvectionSrcForScalars</a></div><div class="ttdeci">void AdvectionSrcForScalars(const Real &dt, const Box &bx, const int icomp, const int ncomp, const Array4< const Real > &avg_xmom, const Array4< const Real > &avg_ymom, const Array4< const Real > &avg_zmom, const Array4< const Real > &cur_cons, const Array4< const Real > &cell_prim, const Array4< Real > &advectionSrc, const bool &use_mono_adv, Real *max_s_ptr, Real *min_s_ptr, const Array4< const Real > &detJ, const GpuArray< Real, AMREX_SPACEDIM > &cellSizeInv, const Array4< const Real > &mf_m, const AdvType horiz_adv_type, const AdvType vert_adv_type, const Real horiz_upw_frac, const Real vert_upw_frac, const GpuArray< const Array4< Real >, AMREX_SPACEDIM > &flx_arr, const GpuArray< Array4< Real >, AMREX_SPACEDIM > &flx_tmp_arr, const Box &domain, const BCRec *bc_ptr_h)</div><div class="ttdef"><b>Definition:</b> ERF_AdvectionSrcForState.cpp:119</div></div>
<div class="ttc" id="aERF__Advection_8H_html_abd84d0909db1f87c22b5a92aae03e919"><div class="ttname"><a href="ERF__Advection_8H.html#abd84d0909db1f87c22b5a92aae03e919">AdvectionSrcForOpenBC_Tangent_Cons</a></div><div class="ttdeci">void AdvectionSrcForOpenBC_Tangent_Cons(const amrex::Box &bx, const int &dir, const int &icomp, const int &ncomp, const amrex::Array4< amrex::Real > &cell_rhs, 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::Real > &detJ, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const bool do_lo=false)</div></div>
<div class="ttc" id="aERF__IndexDefines_8H_html_ada2d84e0dbbfb8d748defbf018748a70a3095d588151e351b421893f62411d15e"><div class="ttname"><a href="ERF__IndexDefines_8H.html#ada2d84e0dbbfb8d748defbf018748a70a3095d588151e351b421893f62411d15e">AdvType::Weno_7</a></div><div class="ttdeci">@ Weno_7</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_ada2d84e0dbbfb8d748defbf018748a70a6edd20e2cefe2061bec0e1ebbf5c8573"><div class="ttname"><a href="ERF__IndexDefines_8H.html#ada2d84e0dbbfb8d748defbf018748a70a6edd20e2cefe2061bec0e1ebbf5c8573">AdvType::Weno_7Z</a></div><div class="ttdeci">@ Weno_7Z</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_ada2d84e0dbbfb8d748defbf018748a70a9c29b64a7627155706a1cd93682b4351"><div class="ttname"><a href="ERF__IndexDefines_8H.html#ada2d84e0dbbfb8d748defbf018748a70a9c29b64a7627155706a1cd93682b4351">AdvType::Weno_3Z</a></div><div class="ttdeci">@ Weno_3Z</div></div>
<div class="ttc" id="aERF__IndexDefines_8H_html_ada2d84e0dbbfb8d748defbf018748a70ada35597e5f00fa867fe6e39d651a2195"><div class="ttname"><a href="ERF__IndexDefines_8H.html#ada2d84e0dbbfb8d748defbf018748a70ada35597e5f00fa867fe6e39d651a2195">AdvType::Weno_3</a></div><div class="ttdeci">@ Weno_3</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__IndexDefines_8H_html_ada2d84e0dbbfb8d748defbf018748a70aeb75f5943c7837da80e2048f2b344efc"><div class="ttname"><a href="ERF__IndexDefines_8H.html#ada2d84e0dbbfb8d748defbf018748a70aeb75f5943c7837da80e2048f2b344efc">AdvType::Weno_5Z</a></div><div class="ttdeci">@ Weno_5Z</div></div>
<div class="ttc" id="aERF__IndexDefines_8H_html_ada2d84e0dbbfb8d748defbf018748a70af2429e8de544f3ebd7c9668e3dbebaa8"><div class="ttname"><a href="ERF__IndexDefines_8H.html#ada2d84e0dbbfb8d748defbf018748a70af2429e8de544f3ebd7c9668e3dbebaa8">AdvType::Weno_5</a></div><div class="ttdeci">@ Weno_5</div></div>
<div class="ttc" id="aERF__IndexDefines_8H_html_ada2d84e0dbbfb8d748defbf018748a70af8b942d20ed08b25df71ea1c87244210"><div class="ttname"><a href="ERF__IndexDefines_8H.html#ada2d84e0dbbfb8d748defbf018748a70af8b942d20ed08b25df71ea1c87244210">AdvType::Weno_3MZQ</a></div><div class="ttdeci">@ Weno_3MZQ</div></div>
<div class="ttc" id="anamespaceBCVars_html_a00367571247ba920989360b94c156f4aaaf601346b01655b32acf66d07b856aec"><div class="ttname"><a href="namespaceBCVars.html#a00367571247ba920989360b94c156f4aaaf601346b01655b32acf66d07b856aec">BCVars::cons_bc</a></div><div class="ttdeci">@ cons_bc</div><div class="ttdef"><b>Definition:</b> ERF_IndexDefines.H:76</div></div>
<div class="ttc" id="anamespaceERFBCType_html_a40fbdd5a379a92debb4bcaaccac6e01ca2b4e1a73027449d9af83a25ed67ae914"><div class="ttname"><a href="namespaceERFBCType.html#a40fbdd5a379a92debb4bcaaccac6e01ca2b4e1a73027449d9af83a25ed67ae914">ERFBCType::open</a></div><div class="ttdeci">@ open</div><div class="ttdef"><b>Definition:</b> ERF_IndexDefines.H:197</div></div>
</div><!-- fragment --><div class="dynheader">
Here is the call graph for this function:</div>
<div class="dyncontent">
<div class="center"><img src="ERF__AdvectionSrcForState_8cpp_a244d2d95f771960db0f0ec982d3595f8_cgraph.png" border="0" usemap="#aERF__AdvectionSrcForState_8cpp_a244d2d95f771960db0f0ec982d3595f8_cgraph" alt=""/></div>
<map name="aERF__AdvectionSrcForState_8cpp_a244d2d95f771960db0f0ec982d3595f8_cgraph" id="aERF__AdvectionSrcForState_8cpp_a244d2d95f771960db0f0ec982d3595f8_cgraph">
<area shape="rect" title=" " alt="" coords="5,13,180,39"/>
<area shape="rect" href="ERF__Advection_8H.html#abd84d0909db1f87c22b5a92aae03e919" title=" " alt="" coords="228,5,407,47"/>
</map>
</div>
</div>
</div>
</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_87c27e56fd01e6f1a6e2085b6fe8a1a5.html">Advection</a></li><li class="navelem"><a class="el" href="ERF__AdvectionSrcForState_8cpp.html">ERF_AdvectionSrcForState.cpp</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>