-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathDepthStretchTransform_8H_source.html
344 lines (342 loc) · 39.7 KB
/
DepthStretchTransform_8H_source.html
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
<!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>REMORA: /home/runner/work/REMORA/REMORA/Source/Utils/DepthStretchTransform.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">REMORA
</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('DepthStretchTransform_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">DepthStretchTransform.H</div> </div>
</div><!--header-->
<div class="contents">
<a href="DepthStretchTransform_8H.html">Go to the documentation of this file.</a><div class="fragment"><div class="line"><a name="l00001"></a><span class="lineno"> 1</span> <span class="preprocessor">#ifndef STRETCH_H_</span></div>
<div class="line"><a name="l00002"></a><span class="lineno"> 2</span> <span class="preprocessor">#define STRETCH_H_</span></div>
<div class="line"><a name="l00003"></a><span class="lineno"> 3</span>  </div>
<div class="line"><a name="l00004"></a><span class="lineno"> 4</span> <span class="preprocessor">#include <cmath></span></div>
<div class="line"><a name="l00005"></a><span class="lineno"> 5</span> <span class="preprocessor">#include <<a class="code" href="DataStruct_8H.html">DataStruct.H</a>></span></div>
<div class="line"><a name="l00006"></a><span class="lineno"> 6</span> <span class="preprocessor">#include <<a class="code" href="REMORA_8H.html">REMORA.H</a>></span></div>
<div class="line"><a name="l00007"></a><span class="lineno"> 7</span> <span class="preprocessor">#include <<a class="code" href="prob__common_8H.html">prob_common.H</a>></span></div>
<div class="line"><a name="l00008"></a><span class="lineno"> 8</span>  </div>
<div class="line"><a name="l00009"></a><span class="lineno"> 9</span> <span class="keyword">using namespace </span><a class="code" href="DataStruct_8H.html#a4036964f680afc43cfce078311111a9eabbe17e876cafcbc06d3ea0d2a23854d3">amrex</a>;</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="keywordtype">void</span></div>
<div class="line"><a name="l00012"></a><span class="lineno"><a class="line" href="classREMORA.html#a20f12090c7b81ca46d953ccf1319c0f5"> 12</a></span> <a class="code" href="classREMORA.html#a20f12090c7b81ca46d953ccf1319c0f5">REMORA::stretch_transform</a> (<span class="keywordtype">int</span> lev)</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>  std::unique_ptr<MultiFab>& mf_z_w = vec_z_w[lev];</div>
<div class="line"><a name="l00015"></a><span class="lineno"> 15</span>  std::unique_ptr<MultiFab>& mf_z_r = vec_z_r[lev];</div>
<div class="line"><a name="l00016"></a><span class="lineno"> 16</span>  std::unique_ptr<MultiFab>& mf_s_r = vec_s_r[lev];</div>
<div class="line"><a name="l00017"></a><span class="lineno"> 17</span>  std::unique_ptr<MultiFab>& mf_Hz = vec_Hz[lev];</div>
<div class="line"><a name="l00018"></a><span class="lineno"> 18</span>  std::unique_ptr<MultiFab>& mf_h = vec_hOfTheConfusingName[lev];</div>
<div class="line"><a name="l00019"></a><span class="lineno"> 19</span>  std::unique_ptr<MultiFab>& mf_Zt_avg1 = vec_Zt_avg1[lev];</div>
<div class="line"><a name="l00020"></a><span class="lineno"> 20</span>  std::unique_ptr<MultiFab>& mf_z_phys_nd = vec_z_phys_nd[lev];</div>
<div class="line"><a name="l00021"></a><span class="lineno"> 21</span>  <span class="keyword">auto</span> N_loop = Geom(lev).Domain().size()[2]-1; <span class="comment">// Number of vertical "levs" aka, NZ</span></div>
<div class="line"><a name="l00022"></a><span class="lineno"> 22</span>  </div>
<div class="line"><a name="l00023"></a><span class="lineno"> 23</span>  <span class="keywordflow">for</span> ( MFIter mfi(*cons_new[lev], TilingIfNotGPU()); mfi.isValid(); ++mfi )</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>  Array4<Real> <span class="keyword">const</span>& z_w = (mf_z_w)->array(mfi);</div>
<div class="line"><a name="l00026"></a><span class="lineno"> 26</span>  Array4<Real> <span class="keyword">const</span>& z_r = (mf_z_r)->array(mfi);</div>
<div class="line"><a name="l00027"></a><span class="lineno"> 27</span>  Array4<Real> <span class="keyword">const</span>& s_r = (mf_s_r)->array(mfi);</div>
<div class="line"><a name="l00028"></a><span class="lineno"> 28</span>  Array4<Real> <span class="keyword">const</span>& Hz = (mf_Hz)->array(mfi);</div>
<div class="line"><a name="l00029"></a><span class="lineno"> 29</span>  Array4<Real> <span class="keyword">const</span>& h = (mf_h)->array(mfi);</div>
<div class="line"><a name="l00030"></a><span class="lineno"> 30</span>  Array4<Real> <span class="keyword">const</span>& Zt_avg1 = (mf_Zt_avg1)->array(mfi);</div>
<div class="line"><a name="l00031"></a><span class="lineno"> 31</span>  Box bx = mfi.tilebox();</div>
<div class="line"><a name="l00032"></a><span class="lineno"> 32</span>  Box gbx2 = bx;</div>
<div class="line"><a name="l00033"></a><span class="lineno"> 33</span>  gbx2.grow(IntVect(<a class="code" href="IndexDefines_8H.html#a124da9e47fb1b1f3debaf49368f42e38">NGROW</a>,<a class="code" href="IndexDefines_8H.html#a124da9e47fb1b1f3debaf49368f42e38">NGROW</a>,0));</div>
<div class="line"><a name="l00034"></a><span class="lineno"> 34</span>  Box gbx3 = bx;</div>
<div class="line"><a name="l00035"></a><span class="lineno"> 35</span>  gbx3.grow(IntVect(<a class="code" href="IndexDefines_8H.html#a124da9e47fb1b1f3debaf49368f42e38">NGROW</a>+1,<a class="code" href="IndexDefines_8H.html#a124da9e47fb1b1f3debaf49368f42e38">NGROW</a>+1,0));</div>
<div class="line"><a name="l00036"></a><span class="lineno"> 36</span>  Box gbx2D = gbx2;</div>
<div class="line"><a name="l00037"></a><span class="lineno"> 37</span>  gbx2D.makeSlab(2,0);</div>
<div class="line"><a name="l00038"></a><span class="lineno"> 38</span>  Box gbx3D = gbx3;</div>
<div class="line"><a name="l00039"></a><span class="lineno"> 39</span>  gbx3D.makeSlab(2,0);</div>
<div class="line"><a name="l00040"></a><span class="lineno"> 40</span>  Box wgbx3 = gbx3.surroundingNodes(2);</div>
<div class="line"><a name="l00041"></a><span class="lineno"> 41</span>  </div>
<div class="line"><a name="l00042"></a><span class="lineno"> 42</span>  <span class="keyword">const</span> <span class="keyword">auto</span> & geomdata = Geom(lev).data();</div>
<div class="line"><a name="l00043"></a><span class="lineno"> 43</span>  </div>
<div class="line"><a name="l00044"></a><span class="lineno"> 44</span>  <span class="keywordtype">int</span> nz = geom[lev].Domain().length(2);</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> N = nz; <span class="comment">// Number of vertical "levels" aka, NZ</span></div>
<div class="line"><a name="l00047"></a><span class="lineno"> 47</span>  <span class="comment">//forcing tcline to be the same as probhi for now, one in DataStruct.H other in inputs</span></div>
<div class="line"><a name="l00048"></a><span class="lineno"> 48</span>  <a class="code" href="DataStruct_8H.html#aa307030d99f8c225e1a1adb20224ed36a7f80fcc452c2f1ed2bb51b39d0864df1">Real</a> hc=-min(geomdata.ProbHi(2),-solverChoice.tcline); <span class="comment">// Do we need to enforce min here?</span></div>
<div class="line"><a name="l00049"></a><span class="lineno"> 49</span>  <span class="keyword">const</span> <span class="keyword">auto</span> local_theta_s = solverChoice.theta_s;</div>
<div class="line"><a name="l00050"></a><span class="lineno"> 50</span>  <span class="keyword">const</span> <span class="keyword">auto</span> local_theta_b = solverChoice.theta_b;</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>  amrex::ParallelFor(wgbx3, [=] AMREX_GPU_DEVICE (<span class="keywordtype">int</span> i, <span class="keywordtype">int</span> j, <span class="keywordtype">int</span> k)</div>
<div class="line"><a name="l00053"></a><span class="lineno"> 53</span>  {</div>
<div class="line"><a name="l00054"></a><span class="lineno"> 54</span>  z_w(i,j,k) = h(i,j,0);</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>  </div>
<div class="line"><a name="l00057"></a><span class="lineno"> 57</span>  <span class="comment">// ROMS Transform 2</span></div>
<div class="line"><a name="l00058"></a><span class="lineno"> 58</span>  Gpu::streamSynchronize();</div>
<div class="line"><a name="l00059"></a><span class="lineno"> 59</span>  </div>
<div class="line"><a name="l00060"></a><span class="lineno"> 60</span>  amrex::ParallelFor(gbx3D, [=] AMREX_GPU_DEVICE (<span class="keywordtype">int</span> i, <span class="keywordtype">int</span> j, <span class="keywordtype">int</span> )</div>
<div class="line"><a name="l00061"></a><span class="lineno"> 61</span>  {</div>
<div class="line"><a name="l00062"></a><span class="lineno"> 62</span>  <span class="keywordflow">for</span> (<span class="keywordtype">int</span> k=-1; k<=N_loop; k++) {</div>
<div class="line"><a name="l00063"></a><span class="lineno"> 63</span>  <span class="comment">// const Real z = prob_lo[2] + (k + 0.5) * dx[2];</span></div>
<div class="line"><a name="l00064"></a><span class="lineno"> 64</span>  <span class="comment">// const auto prob_lo = geomdata.ProbLo();</span></div>
<div class="line"><a name="l00065"></a><span class="lineno"> 65</span>  <span class="comment">// const auto dx = geomdata.CellSize();</span></div>
<div class="line"><a name="l00066"></a><span class="lineno"> 66</span>  <span class="comment">// This is the z for the bottom of the cell this k corresponds to</span></div>
<div class="line"><a name="l00067"></a><span class="lineno"> 67</span>  <span class="comment">// if we weren't stretching and transforming</span></div>
<div class="line"><a name="l00068"></a><span class="lineno"> 68</span>  <span class="comment">// const Real z = prob_lo[2] + (k) * dx[2];</span></div>
<div class="line"><a name="l00069"></a><span class="lineno"> 69</span>  <span class="comment">// h(i,j,0) = -prob_lo[2]; // conceptually</span></div>
<div class="line"><a name="l00070"></a><span class="lineno"> 70</span> <span class="comment"></span> </div>
<div class="line"><a name="l00071"></a><span class="lineno"> 71</span> <span class="comment"> ////////////////////////////////////////////////////////////////////</span></div>
<div class="line"><a name="l00072"></a><span class="lineno"> 72</span> <span class="comment"></span> <span class="comment">//ROMS Stretching 4</span></div>
<div class="line"><a name="l00073"></a><span class="lineno"> 73</span>  <span class="comment">// Move this block to it's own function for maintainability if needed</span></div>
<div class="line"><a name="l00074"></a><span class="lineno"> 74</span>  <span class="comment">// Information about the problem dimension would need to be added</span></div>
<div class="line"><a name="l00075"></a><span class="lineno"> 75</span>  <span class="comment">// This file would need a k dependent function to return the</span></div>
<div class="line"><a name="l00076"></a><span class="lineno"> 76</span>  <span class="comment">// stretching scalars, or access to 4 vectors of length prob_length(2)</span><span class="comment"></span></div>
<div class="line"><a name="l00077"></a><span class="lineno"> 77</span> <span class="comment"> /////////////////////////////////////////////////////////////////////</span></div>
<div class="line"><a name="l00078"></a><span class="lineno"> 78</span> <span class="comment"></span> <a class="code" href="DataStruct_8H.html#aa307030d99f8c225e1a1adb20224ed36a7f80fcc452c2f1ed2bb51b39d0864df1">Real</a> ds = 1.0_rt / <a class="code" href="DataStruct_8H.html#aa307030d99f8c225e1a1adb20224ed36a7f80fcc452c2f1ed2bb51b39d0864df1">Real</a>(N);</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>  <a class="code" href="DataStruct_8H.html#aa307030d99f8c225e1a1adb20224ed36a7f80fcc452c2f1ed2bb51b39d0864df1">Real</a> cff_r, cff_w, cff1_r, cff1_w, cff2_r, cff2_w, Csur, Cbot;</div>
<div class="line"><a name="l00081"></a><span class="lineno"> 81</span>  <a class="code" href="DataStruct_8H.html#aa307030d99f8c225e1a1adb20224ed36a7f80fcc452c2f1ed2bb51b39d0864df1">Real</a> sc_r,sc_w,Cs_r,Cs_w;</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="keywordflow">if</span> (k==N) <span class="comment">// end of array // pretend we're storing 0?</span></div>
<div class="line"><a name="l00084"></a><span class="lineno"> 84</span>  {</div>
<div class="line"><a name="l00085"></a><span class="lineno"> 85</span>  sc_w=0.0; <span class="comment">//sc_w / hc</span></div>
<div class="line"><a name="l00086"></a><span class="lineno"> 86</span>  Cs_w=0.0; <span class="comment">//Cs_w</span></div>
<div class="line"><a name="l00087"></a><span class="lineno"> 87</span>  }</div>
<div class="line"><a name="l00088"></a><span class="lineno"> 88</span>  <span class="keywordflow">else</span> <span class="keywordflow">if</span> (k==0) <span class="comment">// beginning of array</span></div>
<div class="line"><a name="l00089"></a><span class="lineno"> 89</span>  {</div>
<div class="line"><a name="l00090"></a><span class="lineno"> 90</span>  sc_w=-1.0; <span class="comment">//sc_w / hc</span></div>
<div class="line"><a name="l00091"></a><span class="lineno"> 91</span>  Cs_w=-1.0; <span class="comment">//Cs_w</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>  <span class="keywordflow">else</span></div>
<div class="line"><a name="l00094"></a><span class="lineno"> 94</span>  {</div>
<div class="line"><a name="l00095"></a><span class="lineno"> 95</span>  sc_w=ds*(k-N);</div>
<div class="line"><a name="l00096"></a><span class="lineno"> 96</span>  </div>
<div class="line"><a name="l00097"></a><span class="lineno"> 97</span>  <span class="keywordflow">if</span> (local_theta_s > 0.0_rt) {</div>
<div class="line"><a name="l00098"></a><span class="lineno"> 98</span>  Csur=(1.0_rt-std::cosh(local_theta_s*sc_w))/</div>
<div class="line"><a name="l00099"></a><span class="lineno"> 99</span>  (std::cosh(local_theta_s)-1.0_rt);</div>
<div class="line"><a name="l00100"></a><span class="lineno"> 100</span>  } <span class="keywordflow">else</span> {</div>
<div class="line"><a name="l00101"></a><span class="lineno"> 101</span>  Csur=-sc_w*sc_w;</div>
<div class="line"><a name="l00102"></a><span class="lineno"> 102</span>  }</div>
<div class="line"><a name="l00103"></a><span class="lineno"> 103</span>  </div>
<div class="line"><a name="l00104"></a><span class="lineno"> 104</span>  <span class="keywordflow">if</span> (local_theta_b > 0.0_rt) {</div>
<div class="line"><a name="l00105"></a><span class="lineno"> 105</span>  Cbot=(std::exp(local_theta_b*Csur)-1.0_rt)/</div>
<div class="line"><a name="l00106"></a><span class="lineno"> 106</span>  (1.0_rt-std::exp(-local_theta_b));</div>
<div class="line"><a name="l00107"></a><span class="lineno"> 107</span>  Cs_w=Cbot;</div>
<div class="line"><a name="l00108"></a><span class="lineno"> 108</span>  } <span class="keywordflow">else</span> {</div>
<div class="line"><a name="l00109"></a><span class="lineno"> 109</span>  Cs_w=Csur;</div>
<div class="line"><a name="l00110"></a><span class="lineno"> 110</span>  }</div>
<div class="line"><a name="l00111"></a><span class="lineno"> 111</span>  } <span class="comment">// k test</span></div>
<div class="line"><a name="l00112"></a><span class="lineno"> 112</span>  </div>
<div class="line"><a name="l00113"></a><span class="lineno"> 113</span>  cff_w=hc*sc_w;</div>
<div class="line"><a name="l00114"></a><span class="lineno"> 114</span>  cff1_w=Cs_w;</div>
<div class="line"><a name="l00115"></a><span class="lineno"> 115</span>  </div>
<div class="line"><a name="l00116"></a><span class="lineno"> 116</span>  <span class="comment">//cff_r => sc_r *hc</span></div>
<div class="line"><a name="l00117"></a><span class="lineno"> 117</span>  <span class="comment">//cff1_r => Cs_r</span></div>
<div class="line"><a name="l00118"></a><span class="lineno"> 118</span>  <span class="comment">//Don't do anything special for first/last index</span></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>  sc_r=ds*(k-N+0.5_rt);</div>
<div class="line"><a name="l00121"></a><span class="lineno"> 121</span>  </div>
<div class="line"><a name="l00122"></a><span class="lineno"> 122</span>  <span class="keywordflow">if</span> (local_theta_s > 0.0_rt) {</div>
<div class="line"><a name="l00123"></a><span class="lineno"> 123</span>  Csur=(1.0_rt-std::cosh(local_theta_s*sc_r))/</div>
<div class="line"><a name="l00124"></a><span class="lineno"> 124</span>  (std::cosh(local_theta_s)-1.0_rt);</div>
<div class="line"><a name="l00125"></a><span class="lineno"> 125</span>  } <span class="keywordflow">else</span> {</div>
<div class="line"><a name="l00126"></a><span class="lineno"> 126</span>  Csur=-sc_r*sc_r;</div>
<div class="line"><a name="l00127"></a><span class="lineno"> 127</span>  }</div>
<div class="line"><a name="l00128"></a><span class="lineno"> 128</span>  </div>
<div class="line"><a name="l00129"></a><span class="lineno"> 129</span>  <span class="keywordflow">if</span> (local_theta_b > 0.0_rt) {</div>
<div class="line"><a name="l00130"></a><span class="lineno"> 130</span>  Cbot=(std::exp(local_theta_b*Csur)-1.0_rt)/</div>
<div class="line"><a name="l00131"></a><span class="lineno"> 131</span>  (1.0_rt-std::exp(-local_theta_b));</div>
<div class="line"><a name="l00132"></a><span class="lineno"> 132</span>  Cs_r=Cbot;</div>
<div class="line"><a name="l00133"></a><span class="lineno"> 133</span>  } <span class="keywordflow">else</span> {</div>
<div class="line"><a name="l00134"></a><span class="lineno"> 134</span>  Cs_r=Csur;</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>  }</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="keywordflow">if</span> (i==0&&j==0&&k<N&&k>=0) {</div>
<div class="line"><a name="l00139"></a><span class="lineno"> 139</span>  s_r(0,0,k) = sc_r;</div>
<div class="line"><a name="l00140"></a><span class="lineno"> 140</span>  }</div>
<div class="line"><a name="l00141"></a><span class="lineno"> 141</span>  </div>
<div class="line"><a name="l00142"></a><span class="lineno"> 142</span>  cff_r=hc*sc_r;</div>
<div class="line"><a name="l00143"></a><span class="lineno"> 143</span>  cff1_r=Cs_r;</div>
<div class="line"><a name="l00144"></a><span class="lineno"> 144</span> <span class="comment"></span> </div>
<div class="line"><a name="l00145"></a><span class="lineno"> 145</span> <span class="comment"> ////////////////////////////////////////////////////////////////////</span></div>
<div class="line"><a name="l00146"></a><span class="lineno"> 146</span> <span class="comment"></span> <a class="code" href="DataStruct_8H.html#aa307030d99f8c225e1a1adb20224ed36a7f80fcc452c2f1ed2bb51b39d0864df1">Real</a> hwater=h(i,j,0);</div>
<div class="line"><a name="l00147"></a><span class="lineno"> 147</span>  <span class="comment">//</span></div>
<div class="line"><a name="l00148"></a><span class="lineno"> 148</span>  <span class="comment">// if (k==0) //extra guess added (maybe not actually defined in ROMS)</span></div>
<div class="line"><a name="l00149"></a><span class="lineno"> 149</span>  <span class="comment">// {</span></div>
<div class="line"><a name="l00150"></a><span class="lineno"> 150</span>  <span class="comment">// Real hinv=1.0_rt/(hc+hwater);</span></div>
<div class="line"><a name="l00151"></a><span class="lineno"> 151</span>  <span class="comment">// cff2_r=(cff_r+cff1_r*hwater)*hinv;</span></div>
<div class="line"><a name="l00152"></a><span class="lineno"> 152</span>  <span class="comment">// // z_w(i,j,k-2) = hwater;</span></div>
<div class="line"><a name="l00153"></a><span class="lineno"> 153</span>  <span class="comment">// // z_w(i,j,k-1)= -hwater;</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>  <span class="comment">// z_r(i,j,k) = Zt_avg1(i,j,0)+(Zt_avg1(i,j,0)+hwater)*cff2_r;</span></div>
<div class="line"><a name="l00156"></a><span class="lineno"> 156</span>  <span class="comment">// Hz(i,j,k)=z_w(i,j,k)+hwater;//-z_w(i,j,k-1);</span></div>
<div class="line"><a name="l00157"></a><span class="lineno"> 157</span>  <span class="comment">// } else</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="comment">//Note, we are not supporting ICESHELF flag</span></div>
<div class="line"><a name="l00160"></a><span class="lineno"> 160</span>  </div>
<div class="line"><a name="l00161"></a><span class="lineno"> 161</span>  <a class="code" href="DataStruct_8H.html#aa307030d99f8c225e1a1adb20224ed36a7f80fcc452c2f1ed2bb51b39d0864df1">Real</a> hinv=1.0_rt/(hc+hwater);</div>
<div class="line"><a name="l00162"></a><span class="lineno"> 162</span>  cff2_r=(cff_r+cff1_r*hwater)*hinv;</div>
<div class="line"><a name="l00163"></a><span class="lineno"> 163</span>  cff2_w=(cff_w+cff1_w*hwater)*hinv;</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>(k==0) {</div>
<div class="line"><a name="l00166"></a><span class="lineno"> 166</span>  <span class="comment">// HACK: should actually be the normal expression with coeffs evaluated at k=N-1</span></div>
<div class="line"><a name="l00167"></a><span class="lineno"> 167</span>  z_w(i,j,N-1)=Zt_avg1(i,j,0);</div>
<div class="line"><a name="l00168"></a><span class="lineno"> 168</span>  </div>
<div class="line"><a name="l00169"></a><span class="lineno"> 169</span>  } <span class="keywordflow">else</span> <span class="keywordflow">if</span> (k==-1) {</div>
<div class="line"><a name="l00170"></a><span class="lineno"> 170</span>  h(i,j,0,1) = Zt_avg1(i,j,0)+(Zt_avg1(i,j,0)+hwater)*cff2_w;</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="keywordflow">else</span> {</div>
<div class="line"><a name="l00173"></a><span class="lineno"> 173</span>  z_w(i,j,k-1)=Zt_avg1(i,j,0)+(Zt_avg1(i,j,0)+hwater)*cff2_w;</div>
<div class="line"><a name="l00174"></a><span class="lineno"> 174</span>  </div>
<div class="line"><a name="l00175"></a><span class="lineno"> 175</span>  }</div>
<div class="line"><a name="l00176"></a><span class="lineno"> 176</span>  </div>
<div class="line"><a name="l00177"></a><span class="lineno"> 177</span>  <span class="keywordflow">if</span>(k!=-1) {</div>
<div class="line"><a name="l00178"></a><span class="lineno"> 178</span>  z_r(i,j,k)=Zt_avg1(i,j,0)+(Zt_avg1(i,j,0)+hwater)*cff2_r;</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="comment">// k</span></div>
<div class="line"><a name="l00181"></a><span class="lineno"> 181</span>  });</div>
<div class="line"><a name="l00182"></a><span class="lineno"> 182</span>  </div>
<div class="line"><a name="l00183"></a><span class="lineno"> 183</span>  Gpu::streamSynchronize();</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>  amrex::ParallelFor(gbx3, [=] AMREX_GPU_DEVICE (<span class="keywordtype">int</span> i, <span class="keywordtype">int</span> j, <span class="keywordtype">int</span> k)</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="keywordflow">if</span> (k==0) {</div>
<div class="line"><a name="l00188"></a><span class="lineno"> 188</span>  Hz(i,j,k)=z_w(i,j,k)+h(i,j,0);</div>
<div class="line"><a name="l00189"></a><span class="lineno"> 189</span>  } <span class="keywordflow">else</span> {</div>
<div class="line"><a name="l00190"></a><span class="lineno"> 190</span>  Hz(i,j,k)=z_w(i,j,k)-z_w(i,j,k-1);</div>
<div class="line"><a name="l00191"></a><span class="lineno"> 191</span>  }</div>
<div class="line"><a name="l00192"></a><span class="lineno"> 192</span>  });</div>
<div class="line"><a name="l00193"></a><span class="lineno"> 193</span>  } <span class="comment">// mfi</span></div>
<div class="line"><a name="l00194"></a><span class="lineno"> 194</span>  </div>
<div class="line"><a name="l00195"></a><span class="lineno"> 195</span>  <a class="code" href="DataStruct_8H.html#aa307030d99f8c225e1a1adb20224ed36a7f80fcc452c2f1ed2bb51b39d0864df1">Real</a> time = t_new[lev];</div>
<div class="line"><a name="l00196"></a><span class="lineno"> 196</span>  FillPatch(lev, time, *vec_z_w[lev], GetVecOfPtrs(vec_z_w));</div>
<div class="line"><a name="l00197"></a><span class="lineno"> 197</span>  FillPatch(lev, time, *vec_z_r[lev], GetVecOfPtrs(vec_z_r));</div>
<div class="line"><a name="l00198"></a><span class="lineno"> 198</span>  FillPatch(lev, time, *vec_s_r[lev], GetVecOfPtrs(vec_s_r));</div>
<div class="line"><a name="l00199"></a><span class="lineno"> 199</span>  FillPatch(lev, time, *vec_Hz[lev] , GetVecOfPtrs(vec_Hz));</div>
<div class="line"><a name="l00200"></a><span class="lineno"> 200</span>  </div>
<div class="line"><a name="l00201"></a><span class="lineno"> 201</span>  <span class="comment">// Define nodal z as average of z on w-faces</span></div>
<div class="line"><a name="l00202"></a><span class="lineno"> 202</span>  <span class="keywordflow">for</span> ( MFIter mfi(*cons_new[lev], TilingIfNotGPU()); mfi.isValid(); ++mfi )</div>
<div class="line"><a name="l00203"></a><span class="lineno"> 203</span>  {</div>
<div class="line"><a name="l00204"></a><span class="lineno"> 204</span>  Array4<Real> <span class="keyword">const</span>& z_w = (mf_z_w)->array(mfi);</div>
<div class="line"><a name="l00205"></a><span class="lineno"> 205</span>  Array4<Real> <span class="keyword">const</span>& z_phys_nd = (mf_z_phys_nd)->array(mfi);</div>
<div class="line"><a name="l00206"></a><span class="lineno"> 206</span>  </div>
<div class="line"><a name="l00207"></a><span class="lineno"> 207</span>  Box z_w_box = Box(z_w);</div>
<div class="line"><a name="l00208"></a><span class="lineno"> 208</span>  <span class="keyword">auto</span> <span class="keyword">const</span> lo = amrex::lbound(z_w_box);</div>
<div class="line"><a name="l00209"></a><span class="lineno"> 209</span>  <span class="keyword">auto</span> <span class="keyword">const</span> hi = amrex::ubound(z_w_box);</div>
<div class="line"><a name="l00210"></a><span class="lineno"> 210</span>  </div>
<div class="line"><a name="l00211"></a><span class="lineno"> 211</span>  <span class="comment">//</span></div>
<div class="line"><a name="l00212"></a><span class="lineno"> 212</span>  <span class="comment">// WARNING: z_w(i,j,k) refers to the face on the HIGH side of cell (i,j,k)</span></div>
<div class="line"><a name="l00213"></a><span class="lineno"> 213</span>  <span class="comment">// WARNING: z_phys_nd(i,j,k) refers to the node on the LOW side of cell (i,j,k)</span></div>
<div class="line"><a name="l00214"></a><span class="lineno"> 214</span>  <span class="comment">//</span></div>
<div class="line"><a name="l00215"></a><span class="lineno"> 215</span>  ParallelFor(Box(z_phys_nd), [=] AMREX_GPU_DEVICE (<span class="keywordtype">int</span> i, <span class="keywordtype">int</span> j, <span class="keywordtype">int</span> k)</div>
<div class="line"><a name="l00216"></a><span class="lineno"> 216</span>  {</div>
<div class="line"><a name="l00217"></a><span class="lineno"> 217</span>  <span class="comment">// For now assume all boundaries are constant height --</span></div>
<div class="line"><a name="l00218"></a><span class="lineno"> 218</span>  <span class="comment">// we will enforce periodicity below</span></div>
<div class="line"><a name="l00219"></a><span class="lineno"> 219</span>  <span class="keywordtype">int</span> kk = (k == 0) ? hi.z : k-1;</div>
<div class="line"><a name="l00220"></a><span class="lineno"> 220</span>  if ( i >= lo.x && i <= hi.x-1 && j >= lo.y && j <= hi.y-1 )</div>
<div class="line"><a name="l00221"></a><span class="lineno"> 221</span>  {</div>
<div class="line"><a name="l00222"></a><span class="lineno"> 222</span>  z_phys_nd(i,j,k)=0.25*( z_w(i,j ,kk) + z_w(i+1,j ,kk) +</div>
<div class="line"><a name="l00223"></a><span class="lineno"> 223</span>  z_w(i,j+1,kk) + z_w(i+1,j+1,kk) );</div>
<div class="line"><a name="l00224"></a><span class="lineno"> 224</span>  } <span class="keywordflow">else</span> {</div>
<div class="line"><a name="l00225"></a><span class="lineno"> 225</span>  <span class="keywordtype">int</span> ii = std::min(std::max(i, lo.x), hi.x);</div>
<div class="line"><a name="l00226"></a><span class="lineno"> 226</span>  <span class="keywordtype">int</span> jj = std::min(std::max(j, lo.y), hi.y);</div>
<div class="line"><a name="l00227"></a><span class="lineno"> 227</span>  z_phys_nd(i,j,k) = z_w(ii,jj,kk);</div>
<div class="line"><a name="l00228"></a><span class="lineno"> 228</span>  }</div>
<div class="line"><a name="l00229"></a><span class="lineno"> 229</span>  <span class="keywordflow">if</span> (k == 0) z_phys_nd(i,j,k) *= -1.;</div>
<div class="line"><a name="l00230"></a><span class="lineno"> 230</span>  });</div>
<div class="line"><a name="l00231"></a><span class="lineno"> 231</span>  } <span class="comment">// mf</span></div>
<div class="line"><a name="l00232"></a><span class="lineno"> 232</span>  </div>
<div class="line"><a name="l00233"></a><span class="lineno"> 233</span>  <span class="comment">// Note that we do *not* want to do a multilevel fill here -- we have</span></div>
<div class="line"><a name="l00234"></a><span class="lineno"> 234</span>  <span class="comment">// already filled z_phys_nd on the grown boxes, but we enforce periodicity just in case</span></div>
<div class="line"><a name="l00235"></a><span class="lineno"> 235</span>  vec_z_phys_nd[lev]->FillBoundary(geom[lev].periodicity());</div>
<div class="line"><a name="l00236"></a><span class="lineno"> 236</span> }</div>
<div class="line"><a name="l00237"></a><span class="lineno"> 237</span>  </div>
<div class="line"><a name="l00238"></a><span class="lineno"> 238</span> <span class="preprocessor">#endif</span></div>
<div class="ttc" id="aDataStruct_8H_html"><div class="ttname"><a href="DataStruct_8H.html">DataStruct.H</a></div></div>
<div class="ttc" id="aDataStruct_8H_html_a4036964f680afc43cfce078311111a9eabbe17e876cafcbc06d3ea0d2a23854d3"><div class="ttname"><a href="DataStruct_8H.html#a4036964f680afc43cfce078311111a9eabbe17e876cafcbc06d3ea0d2a23854d3">PlotfileType::amrex</a></div><div class="ttdeci">@ amrex</div></div>
<div class="ttc" id="aDataStruct_8H_html_aa307030d99f8c225e1a1adb20224ed36a7f80fcc452c2f1ed2bb51b39d0864df1"><div class="ttname"><a href="DataStruct_8H.html#aa307030d99f8c225e1a1adb20224ed36a7f80fcc452c2f1ed2bb51b39d0864df1">IC_BC_Type::Real</a></div><div class="ttdeci">@ Real</div></div>
<div class="ttc" id="aIndexDefines_8H_html_a124da9e47fb1b1f3debaf49368f42e38"><div class="ttname"><a href="IndexDefines_8H.html#a124da9e47fb1b1f3debaf49368f42e38">NGROW</a></div><div class="ttdeci">#define NGROW</div><div class="ttdef"><b>Definition:</b> IndexDefines.H:13</div></div>
<div class="ttc" id="aREMORA_8H_html"><div class="ttname"><a href="REMORA_8H.html">REMORA.H</a></div></div>
<div class="ttc" id="aclassREMORA_html_a20f12090c7b81ca46d953ccf1319c0f5"><div class="ttname"><a href="classREMORA.html#a20f12090c7b81ca46d953ccf1319c0f5">REMORA::stretch_transform</a></div><div class="ttdeci">void stretch_transform(int lev)</div><div class="ttdef"><b>Definition:</b> DepthStretchTransform.H:12</div></div>
<div class="ttc" id="aprob__common_8H_html"><div class="ttname"><a href="prob__common_8H.html">prob_common.H</a></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_5c09e96eccedf512ae411d636afd2712.html">Utils</a></li><li class="navelem"><a class="el" href="DepthStretchTransform_8H.html">DepthStretchTransform.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>