-
Notifications
You must be signed in to change notification settings - Fork 1
Expand file tree
/
Copy pathERFPCParticleToMesh_8H_source.html
More file actions
269 lines (267 loc) · 38.3 KB
/
ERFPCParticleToMesh_8H_source.html
File metadata and controls
269 lines (267 loc) · 38.3 KB
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
<!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/Particles/ERFPCParticleToMesh.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('ERFPCParticleToMesh_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">ERFPCParticleToMesh.H</div> </div>
</div><!--header-->
<div class="contents">
<a href="ERFPCParticleToMesh_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="comment">/*! \file ERFPCParticleToMesh.H</span></div>
<div class="line"><a name="l00002"></a><span class="lineno"> 2</span> <span class="comment"> * \brief Template implementation of ERFPC::ERFPCParticleToMesh</span></div>
<div class="line"><a name="l00003"></a><span class="lineno"> 3</span> <span class="comment"> *</span></div>
<div class="line"><a name="l00004"></a><span class="lineno"> 4</span> <span class="comment"> * Separated from ERFPC.H to keep the header lightweight.</span></div>
<div class="line"><a name="l00005"></a><span class="lineno"> 5</span> <span class="comment"> * Include this in .cpp files that call ERFPCParticleToMesh().</span></div>
<div class="line"><a name="l00006"></a><span class="lineno"> 6</span> <span class="comment"> */</span></div>
<div class="line"><a name="l00007"></a><span class="lineno"> 7</span>  </div>
<div class="line"><a name="l00008"></a><span class="lineno"> 8</span> <span class="preprocessor">#ifndef ERFPC_PARTICLE_TO_MESH_H_</span></div>
<div class="line"><a name="l00009"></a><span class="lineno"> 9</span> <span class="preprocessor">#define ERFPC_PARTICLE_TO_MESH_H_</span></div>
<div class="line"><a name="l00010"></a><span class="lineno"> 10</span>  </div>
<div class="line"><a name="l00011"></a><span class="lineno"> 11</span> <span class="preprocessor">#include <<a class="code" href="ERFPC_8H.html">ERFPC.H</a>></span></div>
<div class="line"><a name="l00012"></a><span class="lineno"> 12</span> <span class="preprocessor">#include <<a class="code" href="ERF__TerrainConversion_8H.html">ERF_TerrainConversion.H</a>></span></div>
<div class="line"><a name="l00013"></a><span class="lineno"> 13</span>  </div>
<div class="line"><a name="l00014"></a><span class="lineno"> 14</span> <span class="keyword">template</span><<span class="keyword">typename</span> ValueFunc></div>
<div class="line"><a name="l00015"></a><span class="lineno"> 15</span> <span class="keywordtype">void</span> ERFPC::ERFPCParticleToMesh(amrex::MultiFab& a_mf,</div>
<div class="line"><a name="l00016"></a><span class="lineno"> 16</span>  <span class="keyword">const</span> amrex::MultiFab& a_z_phys_nd,</div>
<div class="line"><a name="l00017"></a><span class="lineno"> 17</span>  <span class="keywordtype">int</span> a_lev, <span class="keywordtype">int</span> a_comp,</div>
<div class="line"><a name="l00018"></a><span class="lineno"> 18</span>  ValueFunc&& value_func)<span class="keyword"> const</span></div>
<div class="line"><a name="l00019"></a><span class="lineno"> 19</span> <span class="keyword"></span>{</div>
<div class="line"><a name="l00020"></a><span class="lineno"> 20</span>  <span class="keyword">using namespace </span><a class="code" href="namespaceamrex.html">amrex</a>;</div>
<div class="line"><a name="l00021"></a><span class="lineno"> 21</span>  </div>
<div class="line"><a name="l00022"></a><span class="lineno"> 22</span>  AMREX_ASSERT(OK());</div>
<div class="line"><a name="l00023"></a><span class="lineno"> 23</span>  AMREX_ASSERT(numParticlesOutOfRange(*<span class="keyword">this</span>, 0) == 0);</div>
<div class="line"><a name="l00024"></a><span class="lineno"> 24</span>  </div>
<div class="line"><a name="l00025"></a><span class="lineno"> 25</span>  <span class="keyword">const</span> <span class="keyword">auto</span>& geom = Geom(a_lev);</div>
<div class="line"><a name="l00026"></a><span class="lineno"> 26</span>  <span class="keyword">const</span> <span class="keyword">auto</span> plo = geom.ProbLoArray();</div>
<div class="line"><a name="l00027"></a><span class="lineno"> 27</span>  <span class="keyword">const</span> <span class="keyword">auto</span> dxi = geom.InvCellSizeArray();</div>
<div class="line"><a name="l00028"></a><span class="lineno"> 28</span>  <span class="keyword">const</span> <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a> <a class="code" href="ERF__InitCustomPert__ABL_8H.html#af9913ed40bcb0ccf6257c6727b2960b0">dx</a> = geom.CellSize(0);</div>
<div class="line"><a name="l00029"></a><span class="lineno"> 29</span>  <span class="keyword">const</span> <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a> <a class="code" href="ERF__InitCustomPert__ABL_8H.html#a3d480c71a8d13543a4a79124b1b694e5">dy</a> = geom.CellSize(1);</div>
<div class="line"><a name="l00030"></a><span class="lineno"> 30</span>  <span class="keyword">const</span> Box domain = geom.Domain();</div>
<div class="line"><a name="l00031"></a><span class="lineno"> 31</span>  <span class="keyword">const</span> <span class="keywordtype">int</span> k_max = domain.bigEnd(AMREX_SPACEDIM-1);</div>
<div class="line"><a name="l00032"></a><span class="lineno"> 32</span>  </div>
<div class="line"><a name="l00033"></a><span class="lineno"> 33</span>  a_mf.setVal(0.0);</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="preprocessor">#ifdef AMREX_USE_OMP</span></div>
<div class="line"><a name="l00036"></a><span class="lineno"> 36</span> <span class="preprocessor">#pragma omp parallel if (Gpu::notInLaunchRegion())</span></div>
<div class="line"><a name="l00037"></a><span class="lineno"> 37</span> <span class="preprocessor">#endif</span></div>
<div class="line"><a name="l00038"></a><span class="lineno"> 38</span>  <span class="keywordflow">for</span> (ParConstIterType pti(*<span class="keyword">this</span>, a_lev); pti.isValid(); ++pti)</div>
<div class="line"><a name="l00039"></a><span class="lineno"> 39</span>  {</div>
<div class="line"><a name="l00040"></a><span class="lineno"> 40</span>  <span class="keyword">const</span> <span class="keywordtype">int</span> grid = pti.index();</div>
<div class="line"><a name="l00041"></a><span class="lineno"> 41</span>  <span class="keyword">const</span> <span class="keyword">auto</span>& ptile = ParticlesAt(a_lev, pti);</div>
<div class="line"><a name="l00042"></a><span class="lineno"> 42</span>  <span class="keyword">const</span> <span class="keyword">auto</span>& aos = ptile.GetArrayOfStructs();</div>
<div class="line"><a name="l00043"></a><span class="lineno"> 43</span>  <span class="keyword">const</span> <span class="keywordtype">int</span> np = aos.numParticles();</div>
<div class="line"><a name="l00044"></a><span class="lineno"> 44</span>  <span class="keywordflow">if</span> (np == 0) { <span class="keywordflow">continue</span>; }</div>
<div class="line"><a name="l00045"></a><span class="lineno"> 45</span>  </div>
<div class="line"><a name="l00046"></a><span class="lineno"> 46</span>  <span class="keyword">auto</span> <a class="code" href="ERF__InitCustomPert__Bubble_8H.html#ab8ec92cc3ea8422c9349409bae98d2a0">rho</a> = a_mf[grid].array();</div>
<div class="line"><a name="l00047"></a><span class="lineno"> 47</span>  <span class="keyword">auto</span> zheight = a_z_phys_nd[grid].array();</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="keyword">const</span> <span class="keyword">auto</span> ptd = ptile.getConstParticleTileData();</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>  <a class="code" href="ERF__InitCustomPert__DataAssimilation__ISV_8H.html#a84e204b057acb38e8a4c184a5c036b77">ParallelFor</a>(np, [=] AMREX_GPU_DEVICE (<span class="keywordtype">int</span> i)</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="keyword">const</span> <span class="keyword">auto</span>& <a class="code" href="ERF__InitCustomPert__SquallLine_8H.html#af83c1918e681461dc4cc7ed4c0828b2c">p</a> = ptd.m_aos[i];</div>
<div class="line"><a name="l00054"></a><span class="lineno"> 54</span>  <span class="keywordflow">if</span> (<a class="code" href="ERF__InitCustomPert__SquallLine_8H.html#af83c1918e681461dc4cc7ed4c0828b2c">p</a>.id() <= 0) { return; }</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>  <span class="keyword">const</span> <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a> pval = value_func(ptd, i);</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>  <span class="keyword">const</span> <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a> lx = (<a class="code" href="ERF__InitCustomPert__SquallLine_8H.html#af83c1918e681461dc4cc7ed4c0828b2c">p</a>.pos(0) - plo[0]) * dxi[0] + <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a>(0.5);</div>
<div class="line"><a name="l00059"></a><span class="lineno"> 59</span>  <span class="keyword">const</span> <span class="keywordtype">int</span> ix = <span class="keyword">static_cast<</span><span class="keywordtype">int</span><span class="keyword">></span>(Math::floor(lx)) - 1;</div>
<div class="line"><a name="l00060"></a><span class="lineno"> 60</span>  <span class="keyword">const</span> <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a> wx1 = lx - <span class="keyword">static_cast<</span><a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a><span class="keyword">></span>(ix + 1);</div>
<div class="line"><a name="l00061"></a><span class="lineno"> 61</span>  <span class="keyword">const</span> <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a> wx0 = <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a>(1.0) - wx1;</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="keyword">const</span> <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a> ly = (<a class="code" href="ERF__InitCustomPert__SquallLine_8H.html#af83c1918e681461dc4cc7ed4c0828b2c">p</a>.pos(1) - plo[1]) * dxi[1] + <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a>(0.5);</div>
<div class="line"><a name="l00064"></a><span class="lineno"> 64</span>  <span class="keyword">const</span> <span class="keywordtype">int</span> iy = <span class="keyword">static_cast<</span><span class="keywordtype">int</span><span class="keyword">></span>(Math::floor(ly)) - 1;</div>
<div class="line"><a name="l00065"></a><span class="lineno"> 65</span>  <span class="keyword">const</span> <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a> wy1 = ly - <span class="keyword">static_cast<</span><a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a><span class="keyword">></span>(iy + 1);</div>
<div class="line"><a name="l00066"></a><span class="lineno"> 66</span>  <span class="keyword">const</span> <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a> wy0 = <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a>(1.0) - wy1;</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>  <span class="keyword">const</span> <span class="keywordtype">int</span> kz = <span class="keyword">static_cast<</span><span class="keywordtype">int</span><span class="keyword">></span>(amrex::Math::floor(</div>
<div class="line"><a name="l00069"></a><span class="lineno"> 69</span>  (<a class="code" href="ERF__InitCustomPert__SquallLine_8H.html#af83c1918e681461dc4cc7ed4c0828b2c">p</a>.pos(AMREX_SPACEDIM-1) - plo[AMREX_SPACEDIM-1])</div>
<div class="line"><a name="l00070"></a><span class="lineno"> 70</span>  * dxi[AMREX_SPACEDIM-1]));</div>
<div class="line"><a name="l00071"></a><span class="lineno"> 71</span>  <span class="keyword">const</span> <span class="keywordtype">int</span> ix_p = <span class="keyword">static_cast<</span><span class="keywordtype">int</span><span class="keyword">></span>(Math::floor((<a class="code" href="ERF__InitCustomPert__SquallLine_8H.html#af83c1918e681461dc4cc7ed4c0828b2c">p</a>.pos(0) - plo[0]) * dxi[0]));</div>
<div class="line"><a name="l00072"></a><span class="lineno"> 72</span>  <span class="keyword">const</span> <span class="keywordtype">int</span> iy_p = <span class="keyword">static_cast<</span><span class="keywordtype">int</span><span class="keyword">></span>(Math::floor((<a class="code" href="ERF__InitCustomPert__SquallLine_8H.html#af83c1918e681461dc4cc7ed4c0828b2c">p</a>.pos(1) - plo[1]) * dxi[1]));</div>
<div class="line"><a name="l00073"></a><span class="lineno"> 73</span>  <span class="keyword">const</span> <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a> fx = (<a class="code" href="ERF__InitCustomPert__SquallLine_8H.html#af83c1918e681461dc4cc7ed4c0828b2c">p</a>.pos(0) - plo[0]) * dxi[0] - <span class="keyword">static_cast<</span><a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a><span class="keyword">></span>(ix_p);</div>
<div class="line"><a name="l00074"></a><span class="lineno"> 74</span>  <span class="keyword">const</span> <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a> fy = (<a class="code" href="ERF__InitCustomPert__SquallLine_8H.html#af83c1918e681461dc4cc7ed4c0828b2c">p</a>.pos(1) - plo[1]) * dxi[1] - <span class="keyword">static_cast<</span><a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a><span class="keyword">></span>(iy_p);</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>  <span class="keyword">auto</span> zface = [=] AMREX_GPU_DEVICE (<span class="keywordtype">int</span> k_face) -> <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a> {</div>
<div class="line"><a name="l00077"></a><span class="lineno"> 77</span>  <span class="keywordflow">return</span> zheight(ix_p , iy_p , k_face) * (<a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a>(1.0)-fx) * (<a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a>(1.0)-fy)</div>
<div class="line"><a name="l00078"></a><span class="lineno"> 78</span>  + zheight(ix_p+1, iy_p , k_face) * fx * (<a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a>(1.0)-fy)</div>
<div class="line"><a name="l00079"></a><span class="lineno"> 79</span>  + zheight(ix_p , iy_p+1, k_face) * (<a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a>(1.0)-fx) * fy</div>
<div class="line"><a name="l00080"></a><span class="lineno"> 80</span>  + zheight(ix_p+1, iy_p+1, k_face) * fx * fy;</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>  </div>
<div class="line"><a name="l00083"></a><span class="lineno"> 83</span>  <span class="keyword">const</span> <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a> z_lo_k = zface(kz);</div>
<div class="line"><a name="l00084"></a><span class="lineno"> 84</span>  <span class="keyword">const</span> <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a> z_hi_k = zface(kz + 1);</div>
<div class="line"><a name="l00085"></a><span class="lineno"> 85</span>  <span class="keyword">const</span> <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a> z_ctr_k = <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a>(0.5) * (z_lo_k + z_hi_k);</div>
<div class="line"><a name="l00086"></a><span class="lineno"> 86</span>  </div>
<div class="line"><a name="l00087"></a><span class="lineno"> 87</span>  <span class="comment">// Particle pos(2) is computational zeta; convert to physical z so</span></div>
<div class="line"><a name="l00088"></a><span class="lineno"> 88</span>  <span class="comment">// the cell-center comparison and weight span are in the same units</span></div>
<div class="line"><a name="l00089"></a><span class="lineno"> 89</span>  <span class="comment">// as zface() (which is built from the height array).</span></div>
<div class="line"><a name="l00090"></a><span class="lineno"> 90</span>  <span class="keyword">const</span> <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a> p_z = <a class="code" href="namespaceERF_1_1ParticlePos.html#a00ae39f3fe6835cd448231222f4ad632">ERF::ParticlePos::z_from_zeta</a>(</div>
<div class="line"><a name="l00091"></a><span class="lineno"> 91</span>  <a class="code" href="ERF__InitCustomPert__SquallLine_8H.html#af83c1918e681461dc4cc7ed4c0828b2c">p</a>.pos(0), <a class="code" href="ERF__InitCustomPert__SquallLine_8H.html#af83c1918e681461dc4cc7ed4c0828b2c">p</a>.pos(1), <a class="code" href="ERF__InitCustomPert__SquallLine_8H.html#af83c1918e681461dc4cc7ed4c0828b2c">p</a>.pos(AMREX_SPACEDIM-1),</div>
<div class="line"><a name="l00092"></a><span class="lineno"> 92</span>  plo, dxi, zheight);</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="keywordtype">int</span> kz_lo, kz_hi;</div>
<div class="line"><a name="l00095"></a><span class="lineno"> 95</span>  <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a> wz_lo, wz_hi;</div>
<div class="line"><a name="l00096"></a><span class="lineno"> 96</span>  <span class="keywordflow">if</span> (p_z >= z_ctr_k) {</div>
<div class="line"><a name="l00097"></a><span class="lineno"> 97</span>  kz_lo = kz;</div>
<div class="line"><a name="l00098"></a><span class="lineno"> 98</span>  kz_hi = kz + 1;</div>
<div class="line"><a name="l00099"></a><span class="lineno"> 99</span>  <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a> span;</div>
<div class="line"><a name="l00100"></a><span class="lineno"> 100</span>  <span class="keywordflow">if</span> (kz < k_max) {</div>
<div class="line"><a name="l00101"></a><span class="lineno"> 101</span>  <span class="keyword">const</span> <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a> z_hi2 = zface(kz + 2);</div>
<div class="line"><a name="l00102"></a><span class="lineno"> 102</span>  <span class="keyword">const</span> <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a> z_ctr_hi = <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a>(0.5) * (z_hi_k + z_hi2);</div>
<div class="line"><a name="l00103"></a><span class="lineno"> 103</span>  span = z_ctr_hi - z_ctr_k;</div>
<div class="line"><a name="l00104"></a><span class="lineno"> 104</span>  } <span class="keywordflow">else</span> {</div>
<div class="line"><a name="l00105"></a><span class="lineno"> 105</span>  span = z_hi_k - z_lo_k;</div>
<div class="line"><a name="l00106"></a><span class="lineno"> 106</span>  }</div>
<div class="line"><a name="l00107"></a><span class="lineno"> 107</span>  wz_hi = (p_z - z_ctr_k) / span;</div>
<div class="line"><a name="l00108"></a><span class="lineno"> 108</span>  wz_lo = <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a>(1.0) - wz_hi;</div>
<div class="line"><a name="l00109"></a><span class="lineno"> 109</span>  } <span class="keywordflow">else</span> {</div>
<div class="line"><a name="l00110"></a><span class="lineno"> 110</span>  kz_hi = kz;</div>
<div class="line"><a name="l00111"></a><span class="lineno"> 111</span>  kz_lo = kz - 1;</div>
<div class="line"><a name="l00112"></a><span class="lineno"> 112</span>  <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a> span;</div>
<div class="line"><a name="l00113"></a><span class="lineno"> 113</span>  <span class="keywordflow">if</span> (kz > 0) {</div>
<div class="line"><a name="l00114"></a><span class="lineno"> 114</span>  <span class="keyword">const</span> <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a> z_lo2 = zface(kz - 1);</div>
<div class="line"><a name="l00115"></a><span class="lineno"> 115</span>  <span class="keyword">const</span> <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a> z_ctr_lo = <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a>(0.5) * (z_lo2 + z_lo_k);</div>
<div class="line"><a name="l00116"></a><span class="lineno"> 116</span>  span = z_ctr_k - z_ctr_lo;</div>
<div class="line"><a name="l00117"></a><span class="lineno"> 117</span>  } <span class="keywordflow">else</span> {</div>
<div class="line"><a name="l00118"></a><span class="lineno"> 118</span>  span = z_hi_k - z_lo_k;</div>
<div class="line"><a name="l00119"></a><span class="lineno"> 119</span>  }</div>
<div class="line"><a name="l00120"></a><span class="lineno"> 120</span>  wz_lo = (z_ctr_k - p_z) / span;</div>
<div class="line"><a name="l00121"></a><span class="lineno"> 121</span>  wz_hi = <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a>(1.0) - wz_lo;</div>
<div class="line"><a name="l00122"></a><span class="lineno"> 122</span>  }</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="comment">// Per-corner node-averaged cell thickness for volume normalization.</span></div>
<div class="line"><a name="l00125"></a><span class="lineno"> 125</span>  <span class="keyword">auto</span> dz_cell = [=] AMREX_GPU_DEVICE (<span class="keywordtype">int</span> ic, <span class="keywordtype">int</span> jc, <span class="keywordtype">int</span> kc) -> <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a> {</div>
<div class="line"><a name="l00126"></a><span class="lineno"> 126</span>  <span class="keyword">const</span> <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a> <a class="code" href="ERF__InitCustomPert__Bubble_8H.html#abc695cd190273d99450d0b8f595c5c6b">dz</a> = <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a>(0.25) * (</div>
<div class="line"><a name="l00127"></a><span class="lineno"> 127</span>  zheight(ic , jc , kc+1) - zheight(ic , jc , kc)</div>
<div class="line"><a name="l00128"></a><span class="lineno"> 128</span>  + zheight(ic+1, jc , kc+1) - zheight(ic+1, jc , kc)</div>
<div class="line"><a name="l00129"></a><span class="lineno"> 129</span>  + zheight(ic , jc+1, kc+1) - zheight(ic , jc+1, kc)</div>
<div class="line"><a name="l00130"></a><span class="lineno"> 130</span>  + zheight(ic+1, jc+1, kc+1) - zheight(ic+1, jc+1, kc) );</div>
<div class="line"><a name="l00131"></a><span class="lineno"> 131</span>  <span class="keywordflow">return</span> <a class="code" href="ERF__InitCustomPert__Bubble_8H.html#abc695cd190273d99450d0b8f595c5c6b">dz</a>;</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>  </div>
<div class="line"><a name="l00134"></a><span class="lineno"> 134</span>  <span class="keyword">const</span> <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a> inv_dxdy = <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a>(1.0) / (<a class="code" href="ERF__InitCustomPert__ABL_8H.html#af9913ed40bcb0ccf6257c6727b2960b0">dx</a> * <a class="code" href="ERF__InitCustomPert__ABL_8H.html#a3d480c71a8d13543a4a79124b1b694e5">dy</a>);</div>
<div class="line"><a name="l00135"></a><span class="lineno"> 135</span>  </div>
<div class="line"><a name="l00136"></a><span class="lineno"> 136</span>  <span class="keyword">auto</span> deposit = [=] AMREX_GPU_DEVICE (<span class="keywordtype">int</span> ic, <span class="keywordtype">int</span> jc, <span class="keywordtype">int</span> kc, <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a> w_horiz, <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a> w_z)</div>
<div class="line"><a name="l00137"></a><span class="lineno"> 137</span>  {</div>
<div class="line"><a name="l00138"></a><span class="lineno"> 138</span>  <span class="comment">// Clamp kc only for dz lookup; still write to ghost rho cell.</span></div>
<div class="line"><a name="l00139"></a><span class="lineno"> 139</span>  <span class="keyword">const</span> <span class="keywordtype">int</span> kc_dz = amrex::max(0, amrex::min(k_max, kc));</div>
<div class="line"><a name="l00140"></a><span class="lineno"> 140</span>  <span class="keyword">const</span> <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a> <a class="code" href="ERF__InitCustomPert__Bubble_8H.html#abc695cd190273d99450d0b8f595c5c6b">dz</a> = dz_cell(ic, jc, kc_dz);</div>
<div class="line"><a name="l00141"></a><span class="lineno"> 141</span>  <span class="keyword">const</span> <a class="code" href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a> contrib = pval * w_horiz * w_z * inv_dxdy / <a class="code" href="ERF__InitCustomPert__Bubble_8H.html#abc695cd190273d99450d0b8f595c5c6b">dz</a>;</div>
<div class="line"><a name="l00142"></a><span class="lineno"> 142</span>  Gpu::Atomic::AddNoRet(&<a class="code" href="ERF__InitCustomPert__Bubble_8H.html#ab8ec92cc3ea8422c9349409bae98d2a0">rho</a>(ic, jc, kc, a_comp), contrib);</div>
<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>  </div>
<div class="line"><a name="l00145"></a><span class="lineno"> 145</span>  deposit(ix , iy , kz_lo, wx0*wy0, wz_lo);</div>
<div class="line"><a name="l00146"></a><span class="lineno"> 146</span>  deposit(ix+1, iy , kz_lo, wx1*wy0, wz_lo);</div>
<div class="line"><a name="l00147"></a><span class="lineno"> 147</span>  deposit(ix , iy+1, kz_lo, wx0*wy1, wz_lo);</div>
<div class="line"><a name="l00148"></a><span class="lineno"> 148</span>  deposit(ix+1, iy+1, kz_lo, wx1*wy1, wz_lo);</div>
<div class="line"><a name="l00149"></a><span class="lineno"> 149</span>  deposit(ix , iy , kz_hi, wx0*wy0, wz_hi);</div>
<div class="line"><a name="l00150"></a><span class="lineno"> 150</span>  deposit(ix+1, iy , kz_hi, wx1*wy0, wz_hi);</div>
<div class="line"><a name="l00151"></a><span class="lineno"> 151</span>  deposit(ix , iy+1, kz_hi, wx0*wy1, wz_hi);</div>
<div class="line"><a name="l00152"></a><span class="lineno"> 152</span>  deposit(ix+1, iy+1, kz_hi, wx1*wy1, wz_hi);</div>
<div class="line"><a name="l00153"></a><span class="lineno"> 153</span>  });</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>  </div>
<div class="line"><a name="l00156"></a><span class="lineno"> 156</span>  a_mf.SumBoundary(geom.periodicity());</div>
<div class="line"><a name="l00157"></a><span class="lineno"> 157</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> <span class="preprocessor">#endif</span></div>
<div class="ttc" id="aERFPC_8H_html"><div class="ttname"><a href="ERFPC_8H.html">ERFPC.H</a></div></div>
<div class="ttc" id="aERF__InitCustomPert__ABL_8H_html_a3d480c71a8d13543a4a79124b1b694e5"><div class="ttname"><a href="ERF__InitCustomPert__ABL_8H.html#a3d480c71a8d13543a4a79124b1b694e5">dy</a></div><div class="ttdeci">const Real dy</div><div class="ttdef"><b>Definition:</b> ERF_InitCustomPert_ABL.H:24</div></div>
<div class="ttc" id="aERF__InitCustomPert__ABL_8H_html_af9913ed40bcb0ccf6257c6727b2960b0"><div class="ttname"><a href="ERF__InitCustomPert__ABL_8H.html#af9913ed40bcb0ccf6257c6727b2960b0">dx</a></div><div class="ttdeci">const Real dx</div><div class="ttdef"><b>Definition:</b> ERF_InitCustomPert_ABL.H:23</div></div>
<div class="ttc" id="aERF__InitCustomPert__Bubble_8H_html_ab8ec92cc3ea8422c9349409bae98d2a0"><div class="ttname"><a href="ERF__InitCustomPert__Bubble_8H.html#ab8ec92cc3ea8422c9349409bae98d2a0">rho</a></div><div class="ttdeci">rho</div><div class="ttdef"><b>Definition:</b> ERF_InitCustomPert_Bubble.H:106</div></div>
<div class="ttc" id="aERF__InitCustomPert__Bubble_8H_html_abc695cd190273d99450d0b8f595c5c6b"><div class="ttname"><a href="ERF__InitCustomPert__Bubble_8H.html#abc695cd190273d99450d0b8f595c5c6b">dz</a></div><div class="ttdeci">const Real dz</div><div class="ttdef"><b>Definition:</b> ERF_InitCustomPert_Bubble.H:25</div></div>
<div class="ttc" id="aERF__InitCustomPert__DataAssimilation__ISV_8H_html_a84e204b057acb38e8a4c184a5c036b77"><div class="ttname"><a href="ERF__InitCustomPert__DataAssimilation__ISV_8H.html#a84e204b057acb38e8a4c184a5c036b77">ParallelFor</a></div><div class="ttdeci">ParallelFor(bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept { const Real *dx=geomdata.CellSize();const Real x=(i+myhalf) *dx[0];const Real y=(j+myhalf) *dx[1];const Real Omg=erf_vortex_Gaussian(x, y, xc, yc, R, beta, sigma);const Real deltaT=-(gamma - one)/(two *sigma *sigma) *Omg *Omg;const Real rho_norm=std::pow(one+deltaT, inv_gm1);const Real T=(one+deltaT) *T_inf;const Real p=std::pow(rho_norm, Gamma)/Gamma *rho_0 *a_inf *a_inf;const Real rho_theta=rho_0 *rho_norm *(T *std::pow(p_0/p, rdOcp));state_pert(i, j, k, RhoTheta_comp)=rho_theta - getRhoThetagivenP(p_hse(i, j, k));const Real r2d_xy=std::sqrt((x-xc) *(x-xc)+(y-yc) *(y-yc));state_pert(i, j, k, RhoScalar_comp)=fourth *(one+std::cos(PI *std::min(r2d_xy, R)/R));})</div></div>
<div class="ttc" id="aERF__InitCustomPert__SquallLine_8H_html_af83c1918e681461dc4cc7ed4c0828b2c"><div class="ttname"><a href="ERF__InitCustomPert__SquallLine_8H.html#af83c1918e681461dc4cc7ed4c0828b2c">p</a></div><div class="ttdeci">Real * p</div><div class="ttdef"><b>Definition:</b> ERF_InitCustomPert_SquallLine.H:61</div></div>
<div class="ttc" id="aERF__ShocInterface_8H_html_ab8a9d2a7cbf2084043f890c3d0a68f57"><div class="ttname"><a href="ERF__ShocInterface_8H.html#ab8a9d2a7cbf2084043f890c3d0a68f57">Real</a></div><div class="ttdeci">amrex::Real Real</div><div class="ttdef"><b>Definition:</b> ERF_ShocInterface.H:19</div></div>
<div class="ttc" id="aERF__TerrainConversion_8H_html"><div class="ttname"><a href="ERF__TerrainConversion_8H.html">ERF_TerrainConversion.H</a></div></div>
<div class="ttc" id="anamespaceERF_1_1ParticlePos_html_a00ae39f3fe6835cd448231222f4ad632"><div class="ttname"><a href="namespaceERF_1_1ParticlePos.html#a00ae39f3fe6835cd448231222f4ad632">ERF::ParticlePos::z_from_zeta</a></div><div class="ttdeci">AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real z_from_zeta(amrex::Real x, amrex::Real y, amrex::Real zeta, amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > const &plo, amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > const &dxi, amrex::Array4< amrex::Real const > const &height_arr) noexcept</div><div class="ttdef"><b>Definition:</b> ERF_TerrainConversion.H:51</div></div>
<div class="ttc" id="anamespaceamrex_html"><div class="ttname"><a href="namespaceamrex.html">amrex</a></div><div class="ttdef"><b>Definition:</b> ERF_ConsoleIO.cpp:12</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_fbd11baa4baa1a8b78c4a3d08373cbc6.html">Particles</a></li><li class="navelem"><a class="el" href="ERFPCParticleToMesh_8H.html">ERFPCParticleToMesh.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>