REMORA
Regional Modeling of Oceans Refined Adaptively
Loading...
Searching...
No Matches
REMORA_init_from_netcdf.cpp
Go to the documentation of this file.
1/**
2 * \file REMORA_init_from_netcdf.cpp
3 */
4
5#include <REMORA.H>
6#include <REMORA_Constants.H>
8#include <REMORA_DataStruct.H>
9
10using namespace amrex;
11
12#ifdef REMORA_USE_NETCDF
13
14/** \brief helper function for reading in initial state data from netcdf */
15void
16read_data_from_netcdf (int /*lev*/, const Box& domain, const std::string& fname,
17 FArrayBox& NC_temp_fab, FArrayBox& NC_salt_fab,
18 FArrayBox& NC_xvel_fab, FArrayBox& NC_yvel_fab);
19
20/** \brief helper function for reading in full domain high-resolution initial state data from netcdf */
21void
22read_data_full_domain_from_netcdf (int /*lev*/, const Box& domain, const std::string& fname,
23 FArrayBox& NC_temp_fab, FArrayBox& NC_salt_fab,
24 FArrayBox& NC_xvel_fab, FArrayBox& NC_yvel_fab,
25 IntVect ngrow);
26
27/** \brief helper function for reading in land-sea masks from netcdf */
28void
29read_masks_from_netcdf (int /*lev*/, const Box& domain, const std::string& fname,
30 FArrayBox& NC_mskr_fab, FArrayBox& NC_msku_fab,
31 FArrayBox& NC_mskv_fab);
32
33/** \brief helper function for reading boundary data from netcdf */
34Real
35read_bdry_from_netcdf (const Box& domain, const std::string& fname,
36 Vector<Vector<FArrayBox>>& bdy_data_xlo,
37 Vector<Vector<FArrayBox>>& bdy_data_xhi,
38 Vector<Vector<FArrayBox>>& bdy_data_ylo,
39 Vector<Vector<FArrayBox>>& bdy_data_yhi,
40 int& width, amrex::Real& start_bdy_time,
41 std::string bdry_time_varname,
42 amrex::GpuArray<amrex::GpuArray<bool, AMREX_SPACEDIM*2>,BdyVars::NumTypes+1>&);
43
44/** \brief helper function to initialize state from netcdf */
45void
47 FArrayBox& temp_fab, FArrayBox& salt_fab,
48 FArrayBox& x_vel_fab, FArrayBox& y_vel_fab,
49 const Vector<FArrayBox>& NC_temp_fab,
50 const Vector<FArrayBox>& NC_salt_fab,
51 const Vector<FArrayBox>& NC_xvel_fab,
52 const Vector<FArrayBox>& NC_yvel_fab);
53
54/** \brief helper function to read bathymetry from netcdf */
55void
56read_bathymetry_from_netcdf (int lev, const Box& domain, const std::string& fname,
57 FArrayBox& NC_h_fab);
58
59/** \brief helper function to read grid variables from netcdf */
60void
61read_grid_vars_from_netcdf (int lev, const Box& domain, const std::string& fname,
62 FArrayBox& NC_pm_fab, FArrayBox& NC_pn_fab,
63 FArrayBox& NC_xr_fab, FArrayBox& NC_yr_fab,
64 FArrayBox& NC_xu_fab, FArrayBox& NC_yu_fab,
65 FArrayBox& NC_xv_fab, FArrayBox& NC_yv_fab,
66 FArrayBox& NC_xp_fab, FArrayBox& NC_yp_fab);
67
68/** \brief helper function to read full-domain high resolution bathymetry from netcdf */
69void
70read_bathymetry_full_domain_from_netcdf (const Box& domain, const std::string& fname,
71 FArrayBox& NC_h_fab, IntVect ngrow);
72
73/** \brief helper function to read full-domain high resolution grid variables from netcdf */
74void
75read_grid_vars_full_domain_from_netcdf (const Box& domain, const std::string& fname,
76 FArrayBox& NC_pm_fab, FArrayBox& NC_pn_fab,
77 IntVect ngrow);
78
79/** \brief helper function to read coriolis factor from netcdf */
80void
81read_coriolis_from_netcdf (int lev, const Box& domain, const std::string& fname, FArrayBox& NC_fcor_fab);
82
83/** \brief helper function to read sea surface height from netcdf */
84void
85read_zeta_from_netcdf (int lev, const Box& domain, const std::string& fname,
86 FArrayBox& NC_zeta_fab);
87
88/** \brief helper function to read high-resolution full-domain sea surface height from netcdf */
89void
90read_zeta_full_domain_from_netcdf (int lev, const Box& domain, const std::string& fname,
91 FArrayBox& NC_zeta_fab, IntVect ngrow);
92
93/** \brief helper function to read climatology nudging from netcdf */
94void
95read_clim_nudg_coeff_from_netcdf (int lev, const Box& domain, const std::string& fname,
96 bool do_m2_clim_nudg,
97 bool do_m3_clim_nudg,
98 bool do_temp_clim_nudg,
99 bool do_salt_clim_nudg,
100 FArrayBox& NC_M2NC_fab,
101 FArrayBox& NC_M3NC_fab,
102 FArrayBox& NC_TempNC_fab,
103 FArrayBox& NC_SaltNC_fab);
104
105/** \brief helper function to read in vector of data from netcdf */
106void read_vec_from_netcdf (int lev, const amrex::Vector<std::string>& fnames, const std::string& field_name, amrex::Vector<int>& vec_dat);
107
108/**
109 * @param lev Integer specifying the current level
110 */
111void
113{
114 // *** FArrayBox's at this level for holding the INITIAL data
115 Vector<FArrayBox> NC_temp_fab ; NC_temp_fab.resize(num_boxes_at_level[lev]);
116 Vector<FArrayBox> NC_salt_fab ; NC_salt_fab.resize(num_boxes_at_level[lev]);
117 Vector<FArrayBox> NC_xvel_fab ; NC_xvel_fab.resize(num_boxes_at_level[lev]);
118 Vector<FArrayBox> NC_yvel_fab ; NC_yvel_fab.resize(num_boxes_at_level[lev]);
119
120 for (int idx = 0; idx < num_boxes_at_level[lev]; idx++)
121 {
122 read_data_from_netcdf(lev, boxes_at_level[lev][idx], nc_init_file[lev][idx],
123 NC_temp_fab[idx], NC_salt_fab[idx],
124 NC_xvel_fab[idx], NC_yvel_fab[idx]);
125 }
126
127
128 MultiFab mf_temp(*cons_new[lev], make_alias, Temp_comp, 1);
129 MultiFab mf_salt(*cons_new[lev], make_alias, Salt_comp, 1);
130
131#ifdef _OPENMP
132#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
133#endif
134 {
135 // Don't tile this since we are operating on full FABs in this routine
136 for ( MFIter mfi(*cons_new[lev], false); mfi.isValid(); ++mfi )
137 {
138 // Define fabs for holding the initial data
139 FArrayBox &temp_fab = mf_temp[mfi];
140 FArrayBox &salt_fab = mf_salt[mfi];
141 FArrayBox &xvel_fab = (*xvel_new[lev])[mfi];
142 FArrayBox &yvel_fab = (*yvel_new[lev])[mfi];
143
144 init_state_from_netcdf(lev, temp_fab, salt_fab,
145 xvel_fab, yvel_fab,
146 NC_temp_fab, NC_salt_fab,
147 NC_xvel_fab, NC_yvel_fab);
148 } // mf
149 } // omp
150
151 if (nscalar > 1) {
152 cons_new[lev]->setVal(0.0_rt, Tracer_comp + 1, nscalar - 1, cons_new[lev]->nGrowVect());
153 }
154}
155
156void
158{
159 // *** FArrayBox's at this level for holding the INITIAL data
160 Vector<FArrayBox> NC_temp_fab ; NC_temp_fab.resize(1);
161 Vector<FArrayBox> NC_salt_fab ; NC_salt_fab.resize(1);
162 Vector<FArrayBox> NC_xvel_fab ; NC_xvel_fab.resize(1);
163 Vector<FArrayBox> NC_yvel_fab ; NC_yvel_fab.resize(1);
164
166 NC_temp_fab[0], NC_salt_fab[0],
167 NC_xvel_fab[0], NC_yvel_fab[0],
169
170 MultiFab mf_temp(*vec_cons_full_domain[hires_init_level], make_alias, Temp_comp, 1);
171 MultiFab mf_salt(*vec_cons_full_domain[hires_init_level], make_alias, Salt_comp, 1);
172
173#ifdef _OPENMP
174#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
175#endif
176 {
177 // Don't tile this since we are operating on full FABs in this routine
178 for ( MFIter mfi(mf_temp, false); mfi.isValid(); ++mfi )
179 {
180 // Define fabs for holding the initial data
181 FArrayBox &temp_fab = mf_temp[mfi];
182 FArrayBox &salt_fab = mf_salt[mfi];
183 FArrayBox &xvel_fab = (*vec_xvel_full_domain[hires_init_level])[mfi];
184 FArrayBox &yvel_fab = (*vec_yvel_full_domain[hires_init_level])[mfi];
185
186 init_state_from_netcdf(hires_init_level, temp_fab, salt_fab,
187 xvel_fab, yvel_fab,
188 NC_temp_fab, NC_salt_fab,
189 NC_xvel_fab, NC_yvel_fab);
190 } // mf
191 } // omp
192
194 if (nscalar > 1) {
196 }
197
198 // Average down to fill levels below hires_grid_level. Use a special average_down so
199 // grow cells get populated by averaged down fine data
200 for (int lev=hires_init_level-1; lev >= 0; lev--) {
204 }
205}
206
207/**
208 * @param lev Integer specifying the current level
209 */
210void
212{
213 // *** FArrayBox's at this level for holding the INITIAL data
214 Vector<FArrayBox> NC_zeta_fab ; NC_zeta_fab.resize(num_boxes_at_level[lev]);
215
216 for (int idx = 0; idx < num_boxes_at_level[lev]; idx++)
217 {
218 read_zeta_from_netcdf(lev,boxes_at_level[lev][idx], nc_init_file[lev][idx],
219 NC_zeta_fab[idx]);
220
221#ifdef _OPENMP
222#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
223#endif
224 {
225 // Don't tile this since we are operating on full FABs in this routine
226 for ( MFIter mfi(*cons_new[lev], false); mfi.isValid(); ++mfi )
227 {
228 FArrayBox &zeta_fab = (*vec_zeta[lev])[mfi];
229
230 //
231 // FArrayBox to FArrayBox copy does "copy on intersection"
232 // This only works here because we have broadcast the FArrayBox of data from the netcdf file to all ranks
233 //
234
235 zeta_fab.template copy<RunOn::Device>(NC_zeta_fab[idx],0,0,1);
236 } // mf
237 } // omp
238 } // idx
239
240 vec_zeta[lev]->FillBoundary(geom[lev].periodicity());
241 (*physbcs[lev])(*vec_zeta[lev],*vec_mskr[lev].get(),0,1,vec_zeta[lev]->nGrowVect(),t_new[lev],zeta_bc(),0,*vec_zeta[lev],*vec_msku[lev],*vec_mskv[lev]);
242// (*physbcs[lev])(*vec_zeta[lev],*vec_mskr[lev].get(),1,1,vec_zeta[lev]->nGrowVect(),t_new[lev],BCVars::zeta_bc);
243// (*physbcs[lev])(*vec_zeta[lev],*vec_mskr[lev].get(),2,1,vec_zeta[lev]->nGrowVect(),t_new[lev],BCVars::zeta_bc);
244
246 Real told = t_new[lev];
247 fill_from_bdyfiles(lev, *vec_zeta[lev], *vec_mskr[lev], told, zeta_bc(),BdyVars::zeta,0,0,
248 *vec_zeta[lev]);
249 }
250 if (lev>0) {
251 FillPatch(lev, t_old[lev], *vec_zeta[lev], GetVecOfPtrs(vec_zeta), zeta_bc(), BdyVars::zeta,
252 0, false,false,0,0,0.0,*vec_zeta[lev]);
253 }
254// fill_from_bdyfiles(lev, *vec_zeta[lev], *vec_mskr[lev], told, BCVars::zeta_bc,BdyVars::zeta,1,1);
255// fill_from_bdyfiles(lev, *vec_zeta[lev], *vec_mskr[lev], told, BCVars::zeta_bc,BdyVars::zeta,2,2);
256}
257
258void
260{
261 // *** FArrayBox's at this level for holding the INITIAL data
262 Vector<FArrayBox> NC_zeta_fab ; NC_zeta_fab.resize(1);
263
265 NC_zeta_fab[0], cum_ref_ratios[hires_init_level]);
266
267#ifdef _OPENMP
268#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
269#endif
270 {
271 // Don't tile this since we are operating on full FABs in this routine
272 for ( MFIter mfi(*vec_zeta_full_domain[hires_init_level], false); mfi.isValid(); ++mfi )
273 {
274 FArrayBox &zeta_fab = (*vec_zeta_full_domain[hires_init_level])[mfi];
275
276 //
277 // FArrayBox to FArrayBox copy does "copy on intersection"
278 // This only works here because we have broadcast the FArrayBox of data from the netcdf file to all ranks
279 //
280
281 zeta_fab.template copy<RunOn::Device>(NC_zeta_fab[0],0,0,1);
282 } // mf
283 } // omp
284
285 // Average down to fill levels below hires_grid_level. Use a special average_down so
286 // grow cells get populated by averaged down fine data
287 for (int lev=hires_init_level-1; lev >= 0; lev--) {
289 }
290}
291
292
293/**
294 * @param lev Integer specifying the current level
295 */
296void
298{
299 // *** FArrayBox's at this level for holding the INITIAL data
300 Vector<FArrayBox> NC_pm_fab ; NC_pm_fab.resize(num_boxes_at_level[lev]);
301 Vector<FArrayBox> NC_pn_fab ; NC_pn_fab.resize(num_boxes_at_level[lev]);
302
303 Vector<FArrayBox> NC_xr_fab ; NC_xr_fab.resize(num_boxes_at_level[lev]);
304 Vector<FArrayBox> NC_yr_fab ; NC_yr_fab.resize(num_boxes_at_level[lev]);
305 Vector<FArrayBox> NC_xu_fab ; NC_xu_fab.resize(num_boxes_at_level[lev]);
306 Vector<FArrayBox> NC_yu_fab ; NC_yu_fab.resize(num_boxes_at_level[lev]);
307 Vector<FArrayBox> NC_xv_fab ; NC_xv_fab.resize(num_boxes_at_level[lev]);
308 Vector<FArrayBox> NC_yv_fab ; NC_yv_fab.resize(num_boxes_at_level[lev]);
309 Vector<FArrayBox> NC_xp_fab ; NC_xp_fab.resize(num_boxes_at_level[lev]);
310 Vector<FArrayBox> NC_yp_fab ; NC_yp_fab.resize(num_boxes_at_level[lev]);
311
312 for (int idx = 0; idx < num_boxes_at_level[lev]; idx++)
313 {
315 NC_pm_fab[idx], NC_pn_fab[idx],
316 NC_xr_fab[idx], NC_yr_fab[idx],
317 NC_xu_fab[idx], NC_yu_fab[idx],
318 NC_xv_fab[idx], NC_yv_fab[idx],
319 NC_xp_fab[idx], NC_yp_fab[idx]);
320
321#ifdef _OPENMP
322#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
323#endif
324 {
325 // Don't tile this since we are operating on full FABs in this routine
326 for ( MFIter mfi(*cons_new[lev], false); mfi.isValid(); ++mfi )
327 {
328 FArrayBox &pm_fab = (*vec_pm[lev])[mfi];
329 FArrayBox &pn_fab = (*vec_pn[lev])[mfi];
330 FArrayBox &xr_fab = (*vec_xr[lev])[mfi];
331 FArrayBox &yr_fab = (*vec_yr[lev])[mfi];
332 FArrayBox &xu_fab = (*vec_xu[lev])[mfi];
333 FArrayBox &yu_fab = (*vec_yu[lev])[mfi];
334 FArrayBox &xv_fab = (*vec_xv[lev])[mfi];
335 FArrayBox &yv_fab = (*vec_yv[lev])[mfi];
336 FArrayBox &xp_fab = (*vec_xp[lev])[mfi];
337 FArrayBox &yp_fab = (*vec_yp[lev])[mfi];
338
339 //
340 // FArrayBox to FArrayBox copy does "copy on intersection"
341 // This only works here because we have broadcast the FArrayBox of data from the netcdf file to all ranks
342 //
343
344 pm_fab.template copy<RunOn::Device>(NC_pm_fab[idx]);
345 pn_fab.template copy<RunOn::Device>(NC_pn_fab[idx]);
346
347 xr_fab.template copy<RunOn::Device>(NC_xr_fab[idx]);
348 yr_fab.template copy<RunOn::Device>(NC_yr_fab[idx]);
349 xu_fab.template copy<RunOn::Device>(NC_xu_fab[idx]);
350 yu_fab.template copy<RunOn::Device>(NC_yu_fab[idx]);
351 xv_fab.template copy<RunOn::Device>(NC_xv_fab[idx]);
352 yv_fab.template copy<RunOn::Device>(NC_yv_fab[idx]);
353 xp_fab.template copy<RunOn::Device>(NC_xp_fab[idx]);
354 yp_fab.template copy<RunOn::Device>(NC_yp_fab[idx]);
355 } // mf
356 } // omp
357 } // idx
358
359 Real dummy_time = 0.0_rt;
360 if (lev > 0) {
361 FillPatch(lev,dummy_time,*vec_pm[lev],GetVecOfPtrs(vec_pm),
363 BdyVars::null,0,false,true);
364 FillPatch(lev,dummy_time,*vec_pn[lev],GetVecOfPtrs(vec_pn),
366 BdyVars::null,0,false,true);
367 }
368
369
370 int ng = vec_pm[lev]->nGrow();
371
372 const auto& dom_lo = amrex::lbound(geom[lev].Domain());
373 const auto& dom_hi = amrex::ubound(geom[lev].Domain());
374
375 //
376 // We need values of pm and pn outside the domain so we fill
377 // them here with foextrap
378 //
379 // We first fill interior ghost cells because we will need to extrapolate
380 // from ghost cells inside the domain to ghost cells outside the domain
381 //
382 vec_pm[lev]->FillBoundary(geom[lev].periodicity());
383 vec_pn[lev]->FillBoundary(geom[lev].periodicity());
384
385 vec_xr[lev]->FillBoundary(geom[lev].periodicity());
386 vec_yr[lev]->FillBoundary(geom[lev].periodicity());
387 vec_xu[lev]->FillBoundary(geom[lev].periodicity());
388 vec_yu[lev]->FillBoundary(geom[lev].periodicity());
389 vec_xv[lev]->FillBoundary(geom[lev].periodicity());
390 vec_yv[lev]->FillBoundary(geom[lev].periodicity());
391 vec_xp[lev]->FillBoundary(geom[lev].periodicity());
392 vec_yp[lev]->FillBoundary(geom[lev].periodicity());
393
394 for ( MFIter mfi(*vec_pm[lev]); mfi.isValid(); ++mfi )
395 {
396 Box bx = mfi.tilebox();
397
398 auto pm_fab = vec_pm[lev]->array(mfi);
399 auto pn_fab = vec_pn[lev]->array(mfi);
400
401 Box gbx_lox = adjCellLo(bx,0,ng); gbx_lox.grow(1,ng); gbx_lox.setBig (0,dom_lo.x-2);
402 Box gbx_hix = adjCellHi(bx,0,ng); gbx_hix.grow(1,ng); gbx_hix.setSmall(0,dom_hi.x+2);
403 Box gbx_loy = adjCellLo(bx,1,ng); gbx_loy.grow(0,ng); gbx_loy.setBig (1,dom_lo.y-2);
404 Box gbx_hiy = adjCellHi(bx,1,ng); gbx_hiy.grow(0,ng); gbx_hiy.setSmall(1,dom_hi.y+2);
405
406 // if (gbx_lox.ok()) amrex::AllPrint() << "GBX_XLO " << gbx_lox << std::endl;
407 // if (gbx_hix.ok()) amrex::AllPrint() << "GBX_XHI " << gbx_hix << std::endl;
408 // if (gbx_loy.ok()) amrex::AllPrint() << "GBX_YLO " << gbx_loy << std::endl;
409 // if (gbx_hiy.ok()) amrex::AllPrint() << "GBX_YHI " << gbx_hiy << std::endl;
410
411 if (gbx_lox.ok()) {
412 ParallelFor(gbx_lox, [=] AMREX_GPU_DEVICE (int i, int j, int k)
413 {
414 pm_fab(i,j,k,0) = pm_fab(dom_lo.x-1,j,k,0);
415 pn_fab(i,j,k,0) = pn_fab(dom_lo.x-1,j,k,0);
416 });
417 }
418 if (gbx_hix.ok()) {
419 ParallelFor(gbx_hix, [=] AMREX_GPU_DEVICE (int i, int j, int k)
420 {
421 pm_fab(i,j,k,0) = pm_fab(dom_hi.x+1,j,k,0);
422 pn_fab(i,j,k,0) = pn_fab(dom_hi.x+1,j,k,0);
423 });
424 }
425 if (gbx_loy.ok()) {
426 ParallelFor(gbx_loy, [=] AMREX_GPU_DEVICE (int i, int j, int k)
427 {
428 pm_fab(i,j,k,0) = pm_fab(i,dom_lo.y-1,k,0);
429 pn_fab(i,j,k,0) = pn_fab(i,dom_lo.y-1,k,0);
430 });
431 }
432 if (gbx_hiy.ok()) {
433 ParallelFor(gbx_hiy, [=] AMREX_GPU_DEVICE (int i, int j, int k)
434 {
435 pm_fab(i,j,k,0) = pm_fab(i,dom_hi.y+1,k,0);
436 pn_fab(i,j,k,0) = pn_fab(i,dom_hi.y+1,k,0);
437 });
438 }
439 } // mfi
440}
441/**
442 * @param lev Integer specifying the current level
443 */
444void
446{
447 // *** FArrayBox's at this level for holding the INITIAL data
448 Vector<FArrayBox> NC_h_fab ; NC_h_fab.resize(num_boxes_at_level[lev]);
449 for (int idx = 0; idx < num_boxes_at_level[lev]; idx++)
450 {
452 NC_h_fab[idx]);
453
454#ifdef _OPENMP
455#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
456#endif
457 {
458 // Don't tile this since we are operating on full FABs in this routine
459 for ( MFIter mfi(*cons_new[lev], false); mfi.isValid(); ++mfi )
460 {
461 FArrayBox &h_fab = (*vec_h[lev])[mfi];
462
463 //
464 // FArrayBox to FArrayBox copy does "copy on intersection"
465 // This only works here because we have broadcast the FArrayBox of data from the netcdf file to all ranks
466 //
467
468 // Copy into both components of h
469 h_fab.template copy<RunOn::Device>(NC_h_fab[idx],0,0,1);
470 h_fab.template copy<RunOn::Device>(NC_h_fab[idx],0,1,1);
471 } // mf
472 } // omp
473 } // idx
474
475 const double dummy_time = 0.0_rt;
476 // Unconditional foextrap will overwrite periodicity, but EnforcePeriodicity will
477 // be called on h afterwards
478 FillPatch(lev,dummy_time,*vec_h[lev],GetVecOfPtrs(vec_h),
480 BdyVars::null,0,false,true,1);
481 FillPatch(lev,dummy_time,*vec_h[lev],GetVecOfPtrs(vec_h),
483 BdyVars::null,1,false,true,1);
484
485 vec_h[lev]->FillBoundary(geom[lev].periodicity());
486}
487
488/**
489 * @param lev Integer specifying the current level
490 */
491void
493{
494 // *** FArrayBox's at this level for holding the INITIAL data
495 Vector<FArrayBox> NC_fcor_fab ; NC_fcor_fab.resize(num_boxes_at_level[lev]);
496
497 for (int idx = 0; idx < num_boxes_at_level[lev]; idx++)
498 {
500 NC_fcor_fab[idx]);
501
502#ifdef _OPENMP
503#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
504#endif
505 {
506 // Don't tile this since we are operating on full FABs in this routine
507 for ( MFIter mfi(*cons_new[lev], false); mfi.isValid(); ++mfi )
508 {
509 FArrayBox &fcor_fab = (*vec_fcor[lev])[mfi];
510
511 //
512 // FArrayBox to FArrayBox copy does "copy on intersection"
513 // This only works here because we have broadcast the FArrayBox of data from the netcdf file to all ranks
514 //
515
516 fcor_fab.template copy<RunOn::Device>(NC_fcor_fab[idx]);
517 } // mf
518 } // omp
519 } // idx
520 vec_fcor[lev]->FillBoundary(geom[lev].periodicity());
521}
522
523/**
524 * @param lev Integer specifying the current level
525 */
526void
528{
529 // *** FArrayBox's at this level for holding the INITIAL data
530 Vector<FArrayBox> NC_mskr_fab ; NC_mskr_fab.resize(num_boxes_at_level[lev]);
531 Vector<FArrayBox> NC_msku_fab ; NC_msku_fab.resize(num_boxes_at_level[lev]);
532 Vector<FArrayBox> NC_mskv_fab ; NC_mskv_fab.resize(num_boxes_at_level[lev]);
533
534 for (int idx = 0; idx < num_boxes_at_level[lev]; idx++)
535 {
536 read_masks_from_netcdf(lev,boxes_at_level[lev][idx], nc_grid_file[lev][idx],
537 NC_mskr_fab[idx],NC_msku_fab[idx],
538 NC_mskv_fab[idx]);
539
540#ifdef _OPENMP
541#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
542#endif
543 {
544 // Don't tile this since we are operating on full FABs in this routine
545 for ( MFIter mfi(*cons_new[lev], false); mfi.isValid(); ++mfi )
546 {
547 FArrayBox &mskr_fab = (*vec_mskr[lev])[mfi];
548 FArrayBox &msku_fab = (*vec_msku[lev])[mfi];
549 FArrayBox &mskv_fab = (*vec_mskv[lev])[mfi];
550
551 //
552 // FArrayBox to FArrayBox copy does "copy on intersection"
553 // This only works here because we have broadcast the FArrayBox of data from the netcdf file to all ranks
554 //
555
556 mskr_fab.template copy<RunOn::Device>(NC_mskr_fab[idx]);
557 msku_fab.template copy<RunOn::Device>(NC_msku_fab[idx]);
558 mskv_fab.template copy<RunOn::Device>(NC_mskv_fab[idx]);
559 } // mf
560 } // omp
561 } // idx
562
563 update_mskp(lev);
564 vec_mskr[lev]->FillBoundary(geom[lev].periodicity());
565 vec_msku[lev]->FillBoundary(geom[lev].periodicity());
566 vec_mskv[lev]->FillBoundary(geom[lev].periodicity());
567 vec_mskp[lev]->FillBoundary(geom[lev].periodicity());
568}
569
570/**
571 * @param lev Integer specifying the current level
572 */
573void
575{
576 if (nc_bdry_file.empty() || nc_bdry_file[0].empty()) {
577 amrex::Error("NetCDF boundary file name must be provided via input");
578 }
579
580 amrex::Vector<std::string> field_name = {"u", "v", "temp", "salt", "ubar", "vbar", "zeta"};
581 amrex::Vector<IntVect > index_types = {IntVect(1,0,0), IntVect(0,1,0),
582 IntVect(0,0,0), IntVect(0,0,0),
583 IntVect(1,0,0), IntVect(0,1,0),
584 IntVect(0,0,0)};
585 std::vector<bool > is_2d = {false, false, false, false, true, true, true};
586
587 amrex::Print() << "DOING INIT AT LEVEL " << lev << std::endl;
588 int rx = 1; int ry = 1;
589 if (lev > 0) {
590 for (int k = lev-1; k >= 0; k--) {
591 rx *= ref_ratio[k][0];
592 ry *= ref_ratio[k][1];
593 }
594 }
595 for (int ivar = 0; ivar < BdyVars::NumTypes; ivar++) {
596 boundary_series[lev].push_back(std::unique_ptr<NCTimeSeriesBoundary>(new NCTimeSeriesBoundary(lev, geom, nc_bdry_file, field_name[ivar],
598 index_types[ivar],
599 &phys_bc_need_data[ivar], is_2d[ivar], rx, ry)));
600 boundary_series[lev][ivar]->Initialize();
601 }
602}
603
604/**
605 * \brief Helper function to initialize state and velocity data in a Fab.
606 *
607 * @param lev Integer specifying current level
608 * @param state_fab FArrayBox object holding the state data we initialize
609 * @param temp_fab FArrayBox object holding the temperature data we initialize
610 * @param salt_fab FArrayBox object holding the salt data we initialize
611 * @param x_vel_fab FArrayBox object holding the x-velocity data we initialize
612 * @param y_vel_fab FArrayBox object holding the y-velocity data we initialize
613 * @param NC_temp_fab Vector of FArrayBox objects with the REMORA dataset specifying temperature
614 * @param NC_salt_fab Vector of FArrayBox objects with the REMORA dataset specifying salinity
615 * @param NC_xvel_fab Vector of FArrayBox objects with the REMORA dataset specifying x-velocity
616 * @param NC_yvel_fab Vector of FArrayBox objects with the REMORA dataset specifying y-velocity
617 */
618void
620 FArrayBox& temp_fab, FArrayBox& salt_fab,
621 FArrayBox& x_vel_fab, FArrayBox& y_vel_fab,
622 const Vector<FArrayBox>& NC_temp_fab,
623 const Vector<FArrayBox>& NC_salt_fab,
624 const Vector<FArrayBox>& NC_xvel_fab,
625 const Vector<FArrayBox>& NC_yvel_fab)
626{
627 int nboxes = NC_xvel_fab.size();
628 for (int idx = 0; idx < nboxes; idx++)
629 {
630 //
631 // FArrayBox to FArrayBox copy does "copy on intersection"
632 // This only works here because we have broadcast the FArrayBox of data from the netcdf file to all ranks
633 //
634 temp_fab.template copy<RunOn::Device>(NC_temp_fab[idx]);
635 salt_fab.template copy<RunOn::Device>(NC_salt_fab[idx]);
636 x_vel_fab.template copy<RunOn::Device>(NC_xvel_fab[idx]);
637 y_vel_fab.template copy<RunOn::Device>(NC_yvel_fab[idx]);
638 } // idx
639}
640
641/**
642 * @param lev Integer specifying the current level
643 */
644void
646{
647 // *** FArrayBox's at this level for holding the INITIAL data
648 Vector<FArrayBox> NC_M2NC_fab ; NC_M2NC_fab.resize(num_boxes_at_level[lev]);
649 Vector<FArrayBox> NC_M3NC_fab ; NC_M3NC_fab.resize(num_boxes_at_level[lev]);
650 Vector<FArrayBox> NC_TempNC_fab ; NC_TempNC_fab.resize(num_boxes_at_level[lev]);
651 Vector<FArrayBox> NC_SaltNC_fab ; NC_SaltNC_fab.resize(num_boxes_at_level[lev]);
652
653 for (int idx = 0; idx < num_boxes_at_level[lev]; idx++)
654 {
660 NC_M2NC_fab[idx],NC_M3NC_fab[idx],
661 NC_TempNC_fab[idx],NC_SaltNC_fab[idx]);
662
663#ifdef _OPENMP
664#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
665#endif
666 {
667 // Don't tile this since we are operating on full FABs in this routine
668 for ( MFIter mfi(*cons_new[lev], false); mfi.isValid(); ++mfi )
669 {
671 FArrayBox &ubarNC_fab = (*vec_nudg_coeff[BdyVars::ubar][lev])[mfi];
672 ubarNC_fab.template copy<RunOn::Device>(NC_M2NC_fab[idx]);
673 FArrayBox &vbarNC_fab = (*vec_nudg_coeff[BdyVars::vbar][lev])[mfi];
674 vbarNC_fab.template copy<RunOn::Device>(NC_M2NC_fab[idx]);
675 }
677 FArrayBox &uNC_fab = (*vec_nudg_coeff[BdyVars::u][lev])[mfi];
678 uNC_fab.template copy<RunOn::Device>(NC_M3NC_fab[idx]);
679 FArrayBox &vNC_fab = (*vec_nudg_coeff[BdyVars::v][lev])[mfi];
680 vNC_fab.template copy<RunOn::Device>(NC_M3NC_fab[idx]);
681 }
683 FArrayBox &TempNC_fab = (*vec_nudg_coeff[BdyVars::t][lev])[mfi];
684 TempNC_fab.template copy<RunOn::Device>(NC_TempNC_fab[idx]);
685 }
687 FArrayBox &SaltNC_fab = (*vec_nudg_coeff[BdyVars::s][lev])[mfi];
688 SaltNC_fab.template copy<RunOn::Device>(NC_SaltNC_fab[idx]);
689 }
690
691 } // mf
692 } // omp
693 } // idx
694
696 vec_nudg_coeff[BdyVars::ubar][lev]->FillBoundary(geom[lev].periodicity());
697 vec_nudg_coeff[BdyVars::vbar][lev]->FillBoundary(geom[lev].periodicity());
700 }
702 vec_nudg_coeff[BdyVars::u][lev]->FillBoundary(geom[lev].periodicity());
703 vec_nudg_coeff[BdyVars::v][lev]->FillBoundary(geom[lev].periodicity());
706 }
708 vec_nudg_coeff[BdyVars::t][lev]->FillBoundary(geom[lev].periodicity());
710 }
712 vec_nudg_coeff[BdyVars::s][lev]->FillBoundary(geom[lev].periodicity());
714 }
715}
716
717/**
718 * @param[in ] lev level to read in river data
719 */
720void
722{
723 amrex::Vector<int> river_pos_x;
724 amrex::Vector<int> river_pos_y;
725 amrex::Vector<int> river_direction_tmp;
726
727 std::string river_x_name = "river_Xposition";
728 std::string river_y_name = "river_Eposition";
729 std::string river_dir_name = "river_direction";
730
731 read_vec_from_netcdf(lev, nc_riv_file, river_x_name, river_pos_x);
732 read_vec_from_netcdf(lev, nc_riv_file, river_y_name, river_pos_y);
733 read_vec_from_netcdf(lev, nc_riv_file, river_dir_name, river_direction_tmp);
734
735 int nriv = river_pos_x.size();
736 amrex::Gpu::DeviceVector<int> xpos_d(nriv);
737 amrex::Gpu::DeviceVector<int> ypos_d(nriv);
738 river_direction.resize(nriv);
739#ifdef AMREX_USE_GPU
740 Gpu::htod_memcpy(xpos_d.data(), river_pos_x.data(), sizeof(int)*nriv);
741 Gpu::htod_memcpy(ypos_d.data(), river_pos_y.data(), sizeof(int)*nriv);
742 Gpu::htod_memcpy(river_direction.data(), river_direction_tmp.data(), sizeof(int)*nriv);
743#else
744 std::memcpy(xpos_d.data(), river_pos_x.data(), sizeof(int)*nriv);
745 std::memcpy(ypos_d.data(), river_pos_y.data(), sizeof(int)*nriv);
746 std::memcpy(river_direction.data(), river_direction_tmp.data(), sizeof(int)*nriv);
747#endif
748 const int* xpos_ptr = xpos_d.data();
749 const int* ypos_ptr = ypos_d.data();
750
751 for (amrex::MFIter mfi(*(vec_river_position[lev]).get(),true); mfi.isValid(); ++mfi) {
752 amrex::Box bx = mfi.growntilebox(amrex::IntVect(NGROW,NGROW,0));
753 auto river_pos = vec_river_position[lev]->array(mfi);
754 ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int ) {
755 for (int iriv=0; iriv < nriv; iriv++) {
756 int xriv = xpos_ptr[iriv]-1;
757 int yriv = ypos_ptr[iriv]-1;
758 if (i==xriv && j==yriv) {
759 river_pos(i,j,0) = iriv;
760 }
761 }
762 });
763 }
764}
765
766/**
767 * @param[inout] mf multifab of data to convert
768 */
769void
771 Real inv_days_to_inv_s = 1.0_rt / (3600._rt * 24._rt);
772
773 for ( MFIter mfi(*mf, TilingIfNotGPU()); mfi.isValid(); ++mfi )
774 {
775 Array4<Real> const& arr = mf->array(mfi);
776 Box bx = mfi.growntilebox(IntVect(NGROW,NGROW,0));
777 ParallelFor(bx, [=] AMREX_GPU_DEVICE (int i, int j, int k) {
778 arr(i,j,k) *= inv_days_to_inv_s;
779 });
780 }
781
782}
783
784void
786{
787 if (nc_grid_file_hires.empty()) {
788 Abort("Must specify high-resolution grid file when initializing from NetCDF and hires_grid_level > 0");
789 }
790 Vector<FArrayBox> NC_h_fab ; NC_h_fab.resize(1);
792
793 // Don't tile this since we are operating on full FABs in this routine
794 for ( MFIter mfi(*vec_h_full_domain[hires_grid_level], false); mfi.isValid(); ++mfi )
795 {
796 FArrayBox &h_fab = (*vec_h_full_domain[hires_grid_level])[mfi];
797 h_fab.template copy<RunOn::Device>(NC_h_fab[0]);
798 }
799
800 // Average down to fill levels below hires_grid_level. Use a special average_down so
801 // grow cells get populated by averaged down fine data
802 for (int lev=hires_grid_level-1; lev >= 0; lev--) {
804 }
805}
806
807void
809{
810 if (nc_grid_file_hires.empty()) {
811 Abort("Must specify high-resolution grid file when initializing from NetCDF and hires_grid_level > 0");
812 }
813 Vector<FArrayBox> NC_pm_fab ; NC_pm_fab.resize(1);
814 Vector<FArrayBox> NC_pn_fab ; NC_pn_fab.resize(1);
815
817 NC_pm_fab[0], NC_pn_fab[0],
819
820 // Don't tile this since we are operating on full FABs in this routine
821 for ( MFIter mfi(*vec_h_full_domain[hires_grid_level], false); mfi.isValid(); ++mfi )
822 {
823 FArrayBox &pm_fab = (*vec_pm_full_domain[hires_grid_level])[mfi];
824 FArrayBox &pn_fab = (*vec_pn_full_domain[hires_grid_level])[mfi];
825 pm_fab.template copy<RunOn::Device>(NC_pm_fab[0]);
826 pn_fab.template copy<RunOn::Device>(NC_pn_fab[0]);
827 }
828
829 // Average down to fill levels below hires_grid_level. Use a special average_down so
830 // grow cells get populated by averaged down fine data
831 for (int lev=hires_grid_level-1; lev >= 0; lev--) {
834
835 int rrx = ref_ratio[lev][0];
836 int rry = ref_ratio[lev][1];
837 // pm and pn need to be rescaled by the refinement ratio
838 for ( MFIter mfi(*vec_h_full_domain[lev], TilingIfNotGPU()); mfi.isValid(); ++mfi )
839 {
840 Array4<Real> const& pm = vec_pm_full_domain[lev]->array(mfi);
841 Array4<Real> const& pn = vec_pn_full_domain[lev]->array(mfi);
842 Box ubx = mfi.growntilebox(cum_ref_ratios[lev] - IntVect(1,0,0));;
843 Box vbx = mfi.growntilebox(cum_ref_ratios[lev] - IntVect(0,1,0));;
844 ParallelFor(makeSlab(ubx,2,0), [=] AMREX_GPU_DEVICE (int i, int j, int ) {
845 pm(i,j,0) = pm(i,j,0) / Real(rrx);
846 });
847 ParallelFor(makeSlab(vbx,2,0), [=] AMREX_GPU_DEVICE (int i, int j, int ) {
848 pn(i,j,0) = pn(i,j,0) / Real(rry);
849 });
850 }
851 }
852
853 for (int lev=0; lev<=hires_grid_level; lev++) {
854 int ng = vec_pm_full_domain[lev]->nGrow();
855
856 const auto& dom_lo = amrex::lbound(geom[lev].Domain());
857 const auto& dom_hi = amrex::ubound(geom[lev].Domain());
858
859 for ( MFIter mfi(*vec_pm_full_domain[lev]); mfi.isValid(); ++mfi )
860 {
861 Box bx = mfi.tilebox();
862
863 auto pm_fab = vec_pm_full_domain[lev]->array(mfi);
864 auto pn_fab = vec_pn_full_domain[lev]->array(mfi);
865
866 Box gbx_lox = adjCellLo(bx,0,ng); gbx_lox.grow(1,ng); gbx_lox.setBig (0,dom_lo.x-2);
867 Box gbx_hix = adjCellHi(bx,0,ng); gbx_hix.grow(1,ng); gbx_hix.setSmall(0,dom_hi.x+2);
868 Box gbx_loy = adjCellLo(bx,1,ng); gbx_loy.grow(0,ng); gbx_loy.setBig (1,dom_lo.y-2);
869 Box gbx_hiy = adjCellHi(bx,1,ng); gbx_hiy.grow(0,ng); gbx_hiy.setSmall(1,dom_hi.y+2);
870
871 // if (gbx_lox.ok()) amrex::AllPrint() << "GBX_XLO " << gbx_lox << std::endl;
872 // if (gbx_hix.ok()) amrex::AllPrint() << "GBX_XHI " << gbx_hix << std::endl;
873 // if (gbx_loy.ok()) amrex::AllPrint() << "GBX_YLO " << gbx_loy << std::endl;
874 // if (gbx_hiy.ok()) amrex::AllPrint() << "GBX_YHI " << gbx_hiy << std::endl;
875
876 if (gbx_lox.ok()) {
877 ParallelFor(gbx_lox, [=] AMREX_GPU_DEVICE (int i, int j, int k)
878 {
879 pm_fab(i,j,k,0) = pm_fab(dom_lo.x-1,j,k,0);
880 pn_fab(i,j,k,0) = pn_fab(dom_lo.x-1,j,k,0);
881 });
882 }
883 if (gbx_hix.ok()) {
884 ParallelFor(gbx_hix, [=] AMREX_GPU_DEVICE (int i, int j, int k)
885 {
886 pm_fab(i,j,k,0) = pm_fab(dom_hi.x+1,j,k,0);
887 pn_fab(i,j,k,0) = pn_fab(dom_hi.x+1,j,k,0);
888 });
889 }
890 if (gbx_loy.ok()) {
891 ParallelFor(gbx_loy, [=] AMREX_GPU_DEVICE (int i, int j, int k)
892 {
893 pm_fab(i,j,k,0) = pm_fab(i,dom_lo.y-1,k,0);
894 pn_fab(i,j,k,0) = pn_fab(i,dom_lo.y-1,k,0);
895 });
896 }
897 if (gbx_hiy.ok()) {
898 ParallelFor(gbx_hiy, [=] AMREX_GPU_DEVICE (int i, int j, int k)
899 {
900 pm_fab(i,j,k,0) = pm_fab(i,dom_hi.y+1,k,0);
901 pn_fab(i,j,k,0) = pn_fab(i,dom_hi.y+1,k,0);
902 });
903 }
904 } // mfi
905 }
906}
907
908#endif // REMORA_USE_NETCDF
#define NGROW
#define Temp_comp
#define Tracer_comp
#define Salt_comp
void read_vec_from_netcdf(int lev, const amrex::Vector< std::string > &fnames, const std::string &field_name, amrex::Vector< int > &vec_dat)
helper function to read in vector of data from netcdf
void read_coriolis_from_netcdf(int lev, const Box &domain, const std::string &fname, FArrayBox &NC_fcor_fab)
helper function to read coriolis factor from netcdf
void read_grid_vars_from_netcdf(int lev, const Box &domain, const std::string &fname, FArrayBox &NC_pm_fab, FArrayBox &NC_pn_fab, FArrayBox &NC_xr_fab, FArrayBox &NC_yr_fab, FArrayBox &NC_xu_fab, FArrayBox &NC_yu_fab, FArrayBox &NC_xv_fab, FArrayBox &NC_yv_fab, FArrayBox &NC_xp_fab, FArrayBox &NC_yp_fab)
helper function to read grid variables from netcdf
void read_zeta_from_netcdf(int lev, const Box &domain, const std::string &fname, FArrayBox &NC_zeta_fab)
helper function to read sea surface height from netcdf
Real read_bdry_from_netcdf(const Box &domain, const std::string &fname, Vector< Vector< FArrayBox > > &bdy_data_xlo, Vector< Vector< FArrayBox > > &bdy_data_xhi, Vector< Vector< FArrayBox > > &bdy_data_ylo, Vector< Vector< FArrayBox > > &bdy_data_yhi, int &width, amrex::Real &start_bdy_time, std::string bdry_time_varname, amrex::GpuArray< amrex::GpuArray< bool, AMREX_SPACEDIM *2 >, BdyVars::NumTypes+1 > &)
helper function for reading boundary data from netcdf
void read_bathymetry_from_netcdf(int lev, const Box &domain, const std::string &fname, FArrayBox &NC_h_fab)
helper function to read bathymetry from netcdf
void read_data_full_domain_from_netcdf(int, const Box &domain, const std::string &fname, FArrayBox &NC_temp_fab, FArrayBox &NC_salt_fab, FArrayBox &NC_xvel_fab, FArrayBox &NC_yvel_fab, IntVect ngrow)
helper function for reading in full domain high-resolution initial state data from netcdf
void init_state_from_netcdf(int lev, FArrayBox &temp_fab, FArrayBox &salt_fab, FArrayBox &x_vel_fab, FArrayBox &y_vel_fab, const Vector< FArrayBox > &NC_temp_fab, const Vector< FArrayBox > &NC_salt_fab, const Vector< FArrayBox > &NC_xvel_fab, const Vector< FArrayBox > &NC_yvel_fab)
helper function to initialize state from netcdf
void read_bathymetry_full_domain_from_netcdf(const Box &domain, const std::string &fname, FArrayBox &NC_h_fab, IntVect ngrow)
helper function to read full-domain high resolution bathymetry from netcdf
void read_masks_from_netcdf(int, const Box &domain, const std::string &fname, FArrayBox &NC_mskr_fab, FArrayBox &NC_msku_fab, FArrayBox &NC_mskv_fab)
helper function for reading in land-sea masks from netcdf
void read_data_from_netcdf(int, const Box &domain, const std::string &fname, FArrayBox &NC_temp_fab, FArrayBox &NC_salt_fab, FArrayBox &NC_xvel_fab, FArrayBox &NC_yvel_fab)
helper function for reading in initial state data from netcdf
void read_zeta_full_domain_from_netcdf(int lev, const Box &domain, const std::string &fname, FArrayBox &NC_zeta_fab, IntVect ngrow)
helper function to read high-resolution full-domain sea surface height from netcdf
void read_grid_vars_full_domain_from_netcdf(const Box &domain, const std::string &fname, FArrayBox &NC_pm_fab, FArrayBox &NC_pn_fab, IntVect ngrow)
helper function to read full-domain high resolution grid variables from netcdf
void read_clim_nudg_coeff_from_netcdf(int lev, const Box &domain, const std::string &fname, bool do_m2_clim_nudg, bool do_m3_clim_nudg, bool do_temp_clim_nudg, bool do_salt_clim_nudg, FArrayBox &NC_M2NC_fab, FArrayBox &NC_M3NC_fab, FArrayBox &NC_TempNC_fab, FArrayBox &NC_SaltNC_fab)
helper function to read climatology nudging from netcdf
A class to hold and interpolate time series data read from a NetCDF file.
std::string nc_grid_file_hires
Grid file for high resolution bathymetry.
Definition REMORA.H:1556
amrex::Vector< std::string > nc_riv_file
NetCDF river file(s)
Definition REMORA.H:1573
int foextrap_periodic_bc() const noexcept
Definition REMORA.H:1145
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_zeta_full_domain
high resolution initial free surface height (2D)
Definition REMORA.H:456
std::string nc_init_file_hires
Init file for high resolution.
Definition REMORA.H:1563
void init_bathymetry_from_netcdf(int lev)
Bathymetry data initialization from NetCDF file.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_fcor
coriolis factor (2D)
Definition REMORA.H:480
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_xvel_full_domain
multilevel data container for high res initial x velocities (u in ROMS)
Definition REMORA.H:305
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_h
multilevel data container for current step's z velocities (largely unused; W stored separately)
Definition REMORA.H:314
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_pm
horizontal scaling factor: 1 / dx (2D)
Definition REMORA.H:470
void init_zeta_full_domain_from_netcdf()
Full-domain high res sea-surface height data initialization from NetCDF file.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_yv
y_grid on v-points (2D)
Definition REMORA.H:495
amrex::Vector< amrex::MultiFab * > cons_new
multilevel data container for current step's scalar data: temperature, salinity, passive tracer
Definition REMORA.H:294
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_mskr
land/sea mask at cell centers (2D)
Definition REMORA.H:459
void init_grid_vars_from_netcdf(int lev)
Grid variable initialization from NetCDF file.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_yp
y_grid on psi-points (2D)
Definition REMORA.H:500
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_xr
x_grid on rho points (2D)
Definition REMORA.H:483
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_xv
x_grid on v-points (2D)
Definition REMORA.H:493
amrex::Vector< amrex::Vector< amrex::Box > > boxes_at_level
the boxes specified at each level by tagging criteria
Definition REMORA.H:1357
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_msku
land/sea mask at x-faces (2D)
Definition REMORA.H:461
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_pm_full_domain
horizontal scaling factor: 1 / dx (2D) on whole domain
Definition REMORA.H:474
void init_data_from_netcdf(int lev)
Problem initialization from NetCDF file.
void init_masks_from_netcdf(int lev)
Mask data initialization from NetCDF file.
amrex::Vector< amrex::MultiFab * > yvel_new
multilevel data container for current step's y velocities (v in ROMS)
Definition REMORA.H:298
amrex::Box nc_hires_init_box
Box for the full domain at nc_hires_init_level.
Definition REMORA.H:1565
int zeta_bc() const noexcept
Definition REMORA.H:1143
amrex::Vector< amrex::IntVect > cum_ref_ratios
Cumulative refinement ratio between level 0 and level i.
Definition REMORA.H:1568
amrex::Vector< int > num_boxes_at_level
how many boxes specified at each level by tagging criteria
Definition REMORA.H:1353
amrex::Vector< amrex::MultiFab * > xvel_new
multilevel data container for current step's x velocities (u in ROMS)
Definition REMORA.H:296
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_mskp
land/sea mask at cell corners (2D)
Definition REMORA.H:465
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_mskv
land/sea mask at y-faces (2D)
Definition REMORA.H:463
amrex::Vector< std::unique_ptr< REMORAPhysBCFunct > > physbcs
Vector (over level) of functors to apply physical boundary conditions.
Definition REMORA.H:1374
void FillPatch(int lev, amrex::Real time, amrex::MultiFab &mf_to_be_filled, amrex::Vector< amrex::MultiFab * > const &mfs, const int bccomp, const int bdy_var_type=BdyVars::null, const int icomp=0, const bool fill_all=true, const bool fill_set=true, const int n_not_fill=0, const int icomp_calc=0, const amrex::Real dt=amrex::Real(0.0), const amrex::MultiFab &mf_calc=amrex::MultiFab())
Fill a new MultiFab by copying in phi from valid region and filling ghost cells.
std::string nc_clim_coeff_file
NetCDF climatology coefficient file.
Definition REMORA.H:1578
amrex::Vector< std::string > bdry_time_name_byvar
Name of time fields for boundary data.
Definition REMORA.H:1583
void init_riv_pos_from_netcdf(int lev)
static amrex::Vector< std::string > nc_bdry_file
NetCDF boundary data.
Definition REMORA.H:51
void update_mskp(int lev)
Set psi-point mask to be consistent with rho-point mask.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_yvel_full_domain
multilevel data container for high res initial y velocities (v in ROMS)
Definition REMORA.H:307
void init_zeta_from_netcdf(int lev)
Sea-surface height data initialization from NetCDF file.
void init_coriolis_from_netcdf(int lev)
Coriolis parameter data initialization from NetCDF file.
void fill_from_bdyfiles(int lev, amrex::MultiFab &mf_to_fill, const amrex::MultiFab &mf_mask, const amrex::Real time, const int bccomp, const int bdy_var_type, const int icomp_to_fill, const int icomp_calc=0, const amrex::MultiFab &mf_calc=amrex::MultiFab(), const amrex::Real=amrex::Real(0.0))
Fill boundary data from netcdf file.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_h_full_domain
Bathymetry data on the whole domain at each potential level.
Definition REMORA.H:317
void init_data_full_domain_from_netcdf()
High resolution roblem initialization from NetCDF file.
int hires_init_level
Which level the high resolution initialization data is at.
Definition REMORA.H:1561
int nscalar
Number of passive scalars carried in the state.
Definition REMORA.H:1435
void average_down_with_grow_cells(int lev, amrex::Vector< std::unique_ptr< amrex::MultiFab > > &mf)
Average down from level lev+1 to lev in mf, including grow cells.
Definition REMORA.cpp:1831
amrex::Vector< amrex::Real > t_new
new time at each level
Definition REMORA.H:1364
void init_bdry_from_netcdf(int lev)
Boundary data initialization from NetCDF file.
static SolverChoice solverChoice
Container for algorithmic choices.
Definition REMORA.H:1494
amrex::Vector< std::unique_ptr< amrex::iMultiFab > > vec_river_position
iMultiFab for river positions; contents are indices of rivers
Definition REMORA.H:1302
amrex::Vector< amrex::Vector< std::unique_ptr< NCTimeSeriesBoundary > > > boundary_series
Vector over BdyVars of boundary series data containers.
Definition REMORA.H:1298
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_xp
x_grid on psi-points (2D)
Definition REMORA.H:498
static amrex::Vector< amrex::Vector< std::string > > nc_grid_file
NetCDF grid file.
Definition REMORA.H:53
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_zeta
free surface height (2D)
Definition REMORA.H:453
amrex::Gpu::DeviceVector< int > river_direction
Vector over rivers of river direction: 0: u-face; 1: v-face; 2: w-face.
Definition REMORA.H:1304
amrex::Box nc_hires_grid_box
Box for the full domain at nc_hires_grid_level.
Definition REMORA.H:1558
void init_clim_nudg_coeff_from_netcdf(int lev)
Climatology nudging coefficient initialization from NetCDF file.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_cons_full_domain
multilevel data container for high res initial data: temperature, salinity, passive tracer
Definition REMORA.H:303
void init_bathymetry_full_domain_from_netcdf()
Full domain high-res bathymetry data initialization from NetCDF file.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_yu
y_grid on u-points (2D)
Definition REMORA.H:490
void init_grid_vars_full_domain_from_netcdf()
Full domain high-res grid variable initialization from NetCDF file.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_xu
x_grid on u-points (2D)
Definition REMORA.H:488
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_pn_full_domain
horizontal scaling factor: 1 / dy (2D) on whole domain
Definition REMORA.H:476
amrex::GpuArray< amrex::GpuArray< bool, AMREX_SPACEDIM *2 >, BdyVars::NumTypes+1 > phys_bc_need_data
These are flags that indicate whether we need to read in boundary data from file.
Definition REMORA.H:1400
amrex::Vector< amrex::Vector< std::unique_ptr< amrex::MultiFab > > > vec_nudg_coeff
Climatology nudging coefficients.
Definition REMORA.H:534
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_pn
horizontal scaling factor: 1 / dy (2D)
Definition REMORA.H:472
static amrex::Vector< amrex::Vector< std::string > > nc_init_file
NetCDF initialization file.
Definition REMORA.H:52
amrex::Vector< amrex::Real > t_old
old time at each level
Definition REMORA.H:1366
int hires_grid_level
Which level the high resolution bathymetry is at.
Definition REMORA.H:1554
void convert_inv_days_to_inv_s(amrex::MultiFab *)
Convert data in a multifab from inverse days to inverse seconds.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_yr
y_grid on rho points (2D)
Definition REMORA.H:485