12 #include <AMReX_AmrCore.H>
13 #include <AMReX_BCRec.H>
14 #include <AMReX_InterpFaceRegister.H>
16 #include <AMReX_ParallelDescriptor.H>
17 #include <AMReX_ParmParse.H>
18 #include <AMReX_MultiFabUtil.H>
19 #include <AMReX_FillPatchUtil.H>
20 #include <AMReX_VisMF.H>
21 #include <AMReX_PhysBCFunct.H>
22 #include <AMReX_YAFluxRegister.H>
23 #include <AMReX_ErrorList.H>
25 #ifdef AMREX_MEM_PROFILING
26 #include <AMReX_MemProfiler.H>
36 #ifdef REMORA_USE_PARTICLES
40 #ifdef REMORA_USE_NETCDF
47 #include <AMReX_Lazy.H>
62 template<
typename V,
typename T>
64 return std::find(iterable.begin(), iterable.end(), query) != iterable.end();
68 :
public amrex::AmrCore
86 virtual void ErrorEst (
int lev, amrex::TagBoxArray& tags, amrex::Real time,
int ngrow)
override;
92 void init_only (
int lev, amrex::Real time);
98 void post_timestep (
int nstep, amrex::Real time, amrex::Real dt_lev);
108 const amrex::MultiFab& mf,
int comp,
bool local,
bool finemask);
112 int action_interval, amrex::Real action_per);
119 const amrex::DistributionMapping& dm)
override;
125 virtual void RemakeLevel (
int lev, amrex::Real time,
const amrex::BoxArray& ba,
126 const amrex::DistributionMapping& dm)
override;
138 const amrex::DistributionMapping& dm)
override;
148 amrex::MultiFab& source,
149 const amrex::Geometry fine_geom,
150 const amrex::Real
dt,
const amrex::Real time);
156 void WritePlotFile (
int which, amrex::Vector<std::string> plot_var_names);
160 const amrex::Vector<const amrex::MultiFab*> &mf,
161 const amrex::Vector<const amrex::MultiFab*> &mf_nd,
162 const amrex::Vector<std::string> &varnames,
164 const amrex::Vector<int> &level_steps,
165 const std::string &versionName =
"HyperCLaw-V1.1",
166 const std::string &levelPrefix =
"Level_",
167 const std::string &mfPrefix =
"Cell",
168 const amrex::Vector<std::string>& extra_dirs = amrex::Vector<std::string>())
const;
173 const amrex::Vector<amrex::BoxArray> &bArray,
174 const amrex::Vector<std::string> &varnames,
176 const amrex::Vector<int> &level_steps,
177 const std::string &versionName,
178 const std::string &levelPrefix,
179 const std::string &mfPrefix)
const;
200 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_Hz;
202 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_Huon;
204 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_Hvom;
206 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_ru;
208 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_rv;
210 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_rufrc;
212 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_rvfrc;
214 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_Akv;
216 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_Akt;
222 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_diff2;
225 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_x_r;
228 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_y_r;
231 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_z_r;
234 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_z_w;
237 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_s_r;
246 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_sustr;
249 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_svstr;
252 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_rdrag;
255 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_bustr;
258 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_bvstr;
269 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_rubar;
271 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_rvbar;
273 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_rzeta;
275 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_ubar;
277 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_vbar;
279 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_zeta;
282 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_pm;
285 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_pn;
288 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_fcor;
295 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_rhoS;
296 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
vec_rhoA;
302 void Advance (
int lev, amrex::Real time, amrex::Real dt_lev,
int iteration,
int ncycle);
305 void setup_step(
int lev, amrex::Real time, amrex::Real dt_lev);
311 void advance_2d_onestep (
int lev, amrex::Real dt_lev, amrex::Real dtfast_lev,
int my_iif,
int nfast_counter);
315 amrex::MultiFab
const* mf_rhoS,
316 amrex::MultiFab
const* mf_rhoA,
317 amrex::MultiFab * mf_ru,
318 amrex::MultiFab * mf_rv,
319 amrex::MultiFab * mf_rufrc,
320 amrex::MultiFab * mf_rvfrc,
321 amrex::MultiFab * mf_Zt_avg1,
322 std::unique_ptr<amrex::MultiFab>& mf_DU_avg1,
323 std::unique_ptr<amrex::MultiFab>& mf_DU_avg2,
324 std::unique_ptr<amrex::MultiFab>& mf_DV_avg1,
325 std::unique_ptr<amrex::MultiFab>& mf_DV_avg2,
326 std::unique_ptr<amrex::MultiFab>& mf_rubar,
327 std::unique_ptr<amrex::MultiFab>& mf_rvbar,
328 std::unique_ptr<amrex::MultiFab>& mf_rzeta,
329 std::unique_ptr<amrex::MultiFab>& mf_ubar,
330 std::unique_ptr<amrex::MultiFab>& mf_vbar,
331 amrex::MultiFab * mf_zeta,
332 amrex::MultiFab
const* mf_h,
333 amrex::MultiFab
const* mf_pm,
334 amrex::MultiFab
const* mf_pn,
335 amrex::MultiFab
const* mf_fcor,
336 amrex::MultiFab
const* mf_visc2_p,
337 amrex::MultiFab
const* mf_visc2_r,
338 const amrex::Real dtfast_lev,
339 const bool predictor_2d_step,
340 const bool first_2d_step,
int my_iif,
345 amrex::MultiFab& mf_cons,
346 amrex::MultiFab& mf_u , amrex::MultiFab& mf_v,
347 amrex::MultiFab* mf_sstore,
348 amrex::MultiFab* mf_ru , amrex::MultiFab* mf_rv,
349 std::unique_ptr<amrex::MultiFab>& mf_DU_avg1,
350 std::unique_ptr<amrex::MultiFab>& mf_DU_avg2,
351 std::unique_ptr<amrex::MultiFab>& mf_DV_avg1,
352 std::unique_ptr<amrex::MultiFab>& mf_DV_avg2,
353 std::unique_ptr<amrex::MultiFab>& mf_ubar,
354 std::unique_ptr<amrex::MultiFab>& mf_vbar,
355 std::unique_ptr<amrex::MultiFab>& mf_Akv,
356 std::unique_ptr<amrex::MultiFab>& mf_Akt,
357 std::unique_ptr<amrex::MultiFab>& mf_Hz,
358 std::unique_ptr<amrex::MultiFab>& mf_Huon,
359 std::unique_ptr<amrex::MultiFab>& mf_Hvom,
360 std::unique_ptr<amrex::MultiFab>& mf_z_w,
361 amrex::MultiFab
const* mf_h,
362 amrex::MultiFab
const* mf_pm,
363 amrex::MultiFab
const* mf_pn,
365 const amrex::Real dt_lev);
369 amrex::MultiFab& mf_uold, amrex::MultiFab& mf_vold,
370 amrex::MultiFab& mf_u, amrex::MultiFab& mf_v,
371 amrex::MultiFab* mf_ru,
372 amrex::MultiFab* mf_rv,
373 amrex::MultiFab& S_old,
374 amrex::MultiFab& S_new,
375 amrex::MultiFab& mf_W, amrex::MultiFab& mf_DC,
377 const amrex::MultiFab* mf_z_r,
378 const amrex::MultiFab* mf_z_w,
379 const amrex::MultiFab* mf_h,
380 const amrex::MultiFab* mf_pm,
381 const amrex::MultiFab* mf_pn,
382 const amrex::MultiFab* mf_sustr,
383 const amrex::MultiFab* mf_svstr,
384 const amrex::MultiFab* mf_bustr,
385 const amrex::MultiFab* mf_bvstr,
386 const int iic,
const int nfirst,
387 const int nnew,
int nstp,
int nrhs,
388 int N,
const amrex::Real dt_lev);
392 const amrex::Box& gbx,
393 const amrex::Array4<amrex::Real >& tempold,
394 const amrex::Array4<amrex::Real >& tempcache,
395 const amrex::Array4<amrex::Real >& Hz,
396 const amrex::Array4<amrex::Real >& Huon,
397 const amrex::Array4<amrex::Real >& Hvom,
398 const amrex::Array4<amrex::Real >& W,
399 const amrex::Array4<amrex::Real >& DC,
400 const amrex::Array4<amrex::Real >& FC,
401 const amrex::Array4<amrex::Real >& sstore,
402 const amrex::Array4<amrex::Real const>& z_w,
403 const amrex::Array4<amrex::Real const>& h,
404 const amrex::Array4<amrex::Real const>& pm,
405 const amrex::Array4<amrex::Real const>& pn,
406 int iic,
int ntfirst,
int nrhs,
int N,
407 const amrex::Real dt_lev);
409 void rhs_t_3d (
const amrex::Box& bx,
410 const amrex::Box& gbx,
411 const amrex::Array4<amrex::Real >&
t,
412 const amrex::Array4<amrex::Real const>& tempstore,
413 const amrex::Array4<amrex::Real const>& Huon,
414 const amrex::Array4<amrex::Real const>& Hvom,
415 const amrex::Array4<amrex::Real const>& Hz,
416 const amrex::Array4<amrex::Real const>& pn,
417 const amrex::Array4<amrex::Real const>& pm,
418 const amrex::Array4<amrex::Real const>& W,
419 const amrex::Array4<amrex::Real >& FC,
420 int nrhs,
int nnew,
int N,
const amrex::Real dt_lev);
424 const amrex::Box& ybx,
425 const amrex::Array4<amrex::Real const>& uold,
426 const amrex::Array4<amrex::Real const>& vold,
427 const amrex::Array4<amrex::Real >& ru,
428 const amrex::Array4<amrex::Real >& rv,
429 const amrex::Array4<amrex::Real >& rufrc,
430 const amrex::Array4<amrex::Real >& rvfrc,
431 const amrex::Array4<amrex::Real const>& sustr,
432 const amrex::Array4<amrex::Real const>& svstr,
433 const amrex::Array4<amrex::Real const>& bustr,
434 const amrex::Array4<amrex::Real const>& bvstr,
435 const amrex::Array4<amrex::Real const>& Huon,
436 const amrex::Array4<amrex::Real const>& Hvom,
437 const amrex::Array4<amrex::Real const>& pm,
438 const amrex::Array4<amrex::Real const>& pn,
439 const amrex::Array4<amrex::Real const>& W,
440 const amrex::Array4<amrex::Real >& FC,
444 const amrex::Box& ybx,
445 const amrex::Array4<amrex::Real const>& uold,
446 const amrex::Array4<amrex::Real const>& vold,
447 const amrex::Array4<amrex::Real >& ru,
448 const amrex::Array4<amrex::Real >& rv,
449 const amrex::Array4<amrex::Real const>& Duon,
450 const amrex::Array4<amrex::Real const>& Dvom,
453 void rho_eos (
const amrex::Box& bx,
454 const amrex::Array4<amrex::Real const>& state,
455 const amrex::Array4<amrex::Real >& rho,
456 const amrex::Array4<amrex::Real >& rhoA,
457 const amrex::Array4<amrex::Real >& rhoS,
458 const amrex::Array4<amrex::Real const>& Hz,
459 const amrex::Array4<amrex::Real const>& z_w,
460 const amrex::Array4<amrex::Real const>& h,
463 void prsgrd (
const amrex::Box& bx,
464 const amrex::Box& gbx,
465 const amrex::Box& utbx,
466 const amrex::Box& vtbx,
467 const amrex::Array4<amrex::Real >& ru,
468 const amrex::Array4<amrex::Real >& rv,
469 const amrex::Array4<amrex::Real const>& pn,
470 const amrex::Array4<amrex::Real const>& pm,
471 const amrex::Array4<amrex::Real const>& rho,
472 const amrex::Array4<amrex::Real >& FC,
473 const amrex::Array4<amrex::Real const>& Hz,
474 const amrex::Array4<amrex::Real const>& z_r,
475 const amrex::Array4<amrex::Real const>& z_w,
476 const int nrhs,
const int N);
481 const amrex::Box& gbx,
482 const int ioff,
const int joff,
483 const amrex::Array4<amrex::Real >& vel,
484 const amrex::Array4<amrex::Real const>& vel_old,
485 const amrex::Array4<amrex::Real >& rvel,
486 const amrex::Array4<amrex::Real const>& Hz,
487 const amrex::Array4<amrex::Real const>& Akv,
488 const amrex::Array4<amrex::Real >& DC,
489 const amrex::Array4<amrex::Real >& FC,
490 const amrex::Array4<amrex::Real const>& sstr,
491 const amrex::Array4<amrex::Real const>& bstr,
492 const amrex::Array4<amrex::Real const>& z_r,
493 const amrex::Array4<amrex::Real const>& pm,
494 const amrex::Array4<amrex::Real const>& pn,
495 const int iic,
const int ntfirst,
const int nnew,
int nstp,
int nrhs,
int N,
496 const amrex::Real lambda,
const amrex::Real dt_lev);
499 const int ioff,
const int joff,
500 const amrex::Array4<amrex::Real >& phi,
501 const amrex::Array4<amrex::Real const>& Hz,
502 const amrex::Array4<amrex::Real >& Hzk,
503 const amrex::Array4<amrex::Real >& AK,
504 const amrex::Array4<amrex::Real >& Akv,
505 const amrex::Array4<amrex::Real >& BC,
506 const amrex::Array4<amrex::Real >& DC,
507 const amrex::Array4<amrex::Real >& FC,
508 const amrex::Array4<amrex::Real >& CF,
509 const int nnew,
const int N,
510 const amrex::Real dt_lev);
513 const int ioff,
const int joff,
514 const amrex::Array4<amrex::Real >& phi,
515 const amrex::Array4<amrex::Real >& phibar,
516 const amrex::Array4<amrex::Real >& Hphi,
517 const amrex::Array4<amrex::Real const>& Hz,
518 const amrex::Array4<amrex::Real const>& pm_or_pn,
519 const amrex::Array4<amrex::Real const>& Dphi1,
520 const amrex::Array4<amrex::Real const>& Dphi2,
521 const amrex::Array4<amrex::Real >& DC,
522 const amrex::Array4<amrex::Real >& FC,
526 const int ioff,
const int joff,
527 const amrex::Array4<amrex::Real >& phi,
528 const amrex::Array4<amrex::Real const>& Hz,
529 const amrex::Array4<amrex::Real const>& Dphi_avg1,
530 const amrex::Array4<amrex::Real >& DC,
531 const amrex::Array4<amrex::Real >& CF,
532 const amrex::Array4<amrex::Real const>& dxlen,
533 const int nnew,
const int N);
536 void uv3dmix (
const amrex::Box& xbx,
537 const amrex::Box& ybx,
538 const amrex::Array4<amrex::Real >&
u,
539 const amrex::Array4<amrex::Real >&
v,
540 const amrex::Array4<amrex::Real const>& uold,
541 const amrex::Array4<amrex::Real const>& vold,
542 const amrex::Array4<amrex::Real >& rufrc,
543 const amrex::Array4<amrex::Real >& rvfrc,
544 const amrex::Array4<amrex::Real const>& visc2_p,
545 const amrex::Array4<amrex::Real const>& visc2_r,
546 const amrex::Array4<amrex::Real const>& Hz,
547 const amrex::Array4<amrex::Real const>& pm,
548 const amrex::Array4<amrex::Real const>& pn,
550 const amrex::Real dt_lev);
553 void t3dmix (
const amrex::Box& bx,
554 const amrex::Array4<amrex::Real >& state,
555 const amrex::Array4<amrex::Real const>& diff2,
556 const amrex::Array4<amrex::Real const>& Hz,
557 const amrex::Array4<amrex::Real const>& pm,
558 const amrex::Array4<amrex::Real const>& pn,
559 const amrex::Real dt_lev,
562 void coriolis (
const amrex::Box& xbx,
563 const amrex::Box& ybx,
564 const amrex::Array4<amrex::Real const>& uold,
565 const amrex::Array4<amrex::Real const>& vold,
566 const amrex::Array4<amrex::Real >& ru,
567 const amrex::Array4<amrex::Real >& rv,
568 const amrex::Array4<amrex::Real const>& Hz,
569 const amrex::Array4<amrex::Real const>& fomn,
591 void FillPatch (
int lev, amrex::Real time,
592 amrex::MultiFab& mf_to_be_filled,
593 amrex::Vector<amrex::MultiFab*>
const& mfs,
597 const amrex::Real time,
598 const int bdy_var_type);
602 #ifdef REMORA_USE_NETCDF
603 static int total_nc_plot_file_step;
624 void init_stuff (
int lev,
const amrex::BoxArray& ba,
const amrex::DistributionMapping& dm);
634 void FillCoarsePatch (
int lev, amrex::Real time, amrex::MultiFab* mf_fine, amrex::MultiFab* mf_crse);
640 void timeStep (
int lev, amrex::Real time,
int iteration);
643 void timeStepML (amrex::Real time,
int iteration);
658 amrex::Vector<std::string>
PlotFileVarNames (amrex::Vector<std::string> plot_var_names)
const;
661 void setPlotVariables (
const std::string& pp_plot_var_names, amrex::Vector<std::string>& plot_var_names);
664 #ifdef REMORA_USE_NETCDF
666 void WriteNCPlotFile(
int istep)
const;
667 void WriteNCPlotFile_which(
int istep,
int lev,
int which_subdomain,
671 void WriteNCMultiFab (
const amrex::FabArray<amrex::FArrayBox> &fab,
672 const std::string& name,
673 bool set_ghost =
false)
const;
676 void ReadNCMultiFab (amrex::FabArray<amrex::FArrayBox> &fab,
677 const std::string &name,
678 int coordinatorProc = amrex::ParallelDescriptor::IOProcessorNumber(),
679 int allow_empty_mf = 0);
682 amrex::Vector<amrex::Vector<amrex::FArrayBox>> bdy_data_xlo;
683 amrex::Vector<amrex::Vector<amrex::FArrayBox>> bdy_data_xhi;
684 amrex::Vector<amrex::Vector<amrex::FArrayBox>> bdy_data_ylo;
685 amrex::Vector<amrex::Vector<amrex::FArrayBox>> bdy_data_yhi;
687 amrex::Real start_bdy_time;
688 amrex::Real bdy_time_interval;
690 static bool write_history_file;
693 int bdy_set_width{0};
695 void init_data_from_netcdf (
int lev);
696 void init_bdry_from_netcdf ();
697 void init_bathymetry_from_netcdf (
int lev);
698 void init_coriolis_from_netcdf (
int lev);
719 void post_update (amrex::MultiFab& state_mf,
const amrex::Real time,
const amrex::Geometry& geom);
720 void fill_rhs (amrex::MultiFab& rhs_mf,
const amrex::MultiFab& state_mf,
const amrex::Real time,
const amrex::Geometry& geom);
735 amrex::Vector<amrex::Real>
dt;
737 amrex::Vector<std::unique_ptr<REMORAPhysBCFunct>>
physbcs;
741 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
mapfac_m;
742 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
mapfac_u;
743 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
mapfac_v;
745 amrex::Vector<std::unique_ptr<amrex::MultiFab>>
sst;
776 amrex::Real
stop_time = std::numeric_limits<amrex::Real>::max();
811 const amrex::Vector<std::string>
cons_names {
"temp",
"salt",
"scalar"};
814 #ifdef REMORA_USE_PARTICLES
815 const amrex::Vector<std::string>
derived_names {
"dpdx",
"dpdy" ,
"tracer_particle_count"};
824 #ifdef REMORA_USE_PARTICLES
893 switch (spatial_order) {
910 amrex::Error(
"Must specify spatial order to be 2,3,4,5 or 6");
942 int numCores = amrex::ParallelDescriptor::NProcs();
944 numCores = numCores * omp_get_max_threads();
948 numCores * (amrex::ParallelDescriptor::second() -
startCPUTime) +
956 if (amrex::ParallelDescriptor::IOProcessor())
958 datalog[i] = std::make_unique<std::fstream>();
959 datalog[i]->open(filename.c_str(),std::ios::out|std::ios::app);
961 amrex::FileOpenFailed(filename);
964 amrex::ParallelDescriptor::Barrier(
"REMORA::setRecordDataInfo");
967 amrex::Vector<std::unique_ptr<std::fstream> >
datalog;
PlotfileType
Definition: DataStruct.H:37
#define NCONS
Definition: IndexDefines.H:11
bool containerHasElement(const V &iterable, const T &query)
Definition: REMORA.H:63
static PlotfileType plotfile_type
Definition: REMORA.H:839
const amrex::Vector< std::string > cons_names
Definition: REMORA.H:811
static void GotoNextLine(std::istream &is)
Definition: Checkpoint.cpp:8
static int bndry_output_planes_interval
Definition: REMORA.H:861
int nfast
Definition: REMORA.H:790
static int column_interval
Definition: REMORA.H:853
amrex::Vector< std::string > plot_var_names_1
Definition: REMORA.H:809
void rhs_uv_2d(const amrex::Box &xbx, const amrex::Box &ybx, const amrex::Array4< amrex::Real const > &uold, const amrex::Array4< amrex::Real const > &vold, const amrex::Array4< amrex::Real > &ru, const amrex::Array4< amrex::Real > &rv, const amrex::Array4< amrex::Real const > &Duon, const amrex::Array4< amrex::Real const > &Dvom, const int nrhs)
Definition: REMORA_rhs_uv_2d.cpp:20
amrex::Array< amrex::Array< amrex::Real, AMREX_SPACEDIM *2 >, AMREX_SPACEDIM+NCONS > m_bc_extdir_vals
Definition: REMORA.H:760
static amrex::Real fixed_dt
Definition: REMORA.H:787
void prestep_diffusion(const amrex::Box &bx, const amrex::Box &gbx, const int ioff, const int joff, const amrex::Array4< amrex::Real > &vel, const amrex::Array4< amrex::Real const > &vel_old, const amrex::Array4< amrex::Real > &rvel, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &Akv, const amrex::Array4< amrex::Real > &DC, const amrex::Array4< amrex::Real > &FC, const amrex::Array4< amrex::Real const > &sstr, const amrex::Array4< amrex::Real const > &bstr, const amrex::Array4< amrex::Real const > &z_r, const amrex::Array4< amrex::Real const > &pm, const amrex::Array4< amrex::Real const > &pn, const int iic, const int ntfirst, const int nnew, int nstp, int nrhs, int N, const amrex::Real lambda, const amrex::Real dt_lev)
Definition: REMORA_prestep_diffusion.cpp:19
AMREX_FORCE_INLINE amrex::YAFluxRegister * getAdvFluxReg(int lev)
Definition: REMORA.H:917
void init_bcs()
Definition: REMORA_init_bcs.cpp:9
amrex::Vector< amrex::BCRec > domain_bcs_type
Definition: REMORA.H:753
void rhs_t_3d(const amrex::Box &bx, const amrex::Box &gbx, const amrex::Array4< amrex::Real > &t, const amrex::Array4< amrex::Real const > &tempstore, const amrex::Array4< amrex::Real const > &Huon, const amrex::Array4< amrex::Real const > &Hvom, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &pn, const amrex::Array4< amrex::Real const > &pm, const amrex::Array4< amrex::Real const > &W, const amrex::Array4< amrex::Real > &FC, int nrhs, int nnew, int N, const amrex::Real dt_lev)
Definition: REMORA_rhs_t_3d.cpp:25
static amrex::Real previousCPUTimeUsed
Definition: REMORA.H:937
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_fcor
Definition: REMORA.H:288
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_hOfTheConfusingName
Definition: REMORA.H:197
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rubar
Definition: REMORA.H:269
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_pm
Definition: REMORA.H:282
amrex::MultiFab fine_mask
Definition: REMORA.H:886
static amrex::Real init_shrink
Definition: REMORA.H:783
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_DU_avg2
Definition: REMORA.H:263
int plot_int_2
Definition: REMORA.H:803
amrex::Vector< amrex::MultiFab * > cons_new
Definition: REMORA.H:188
virtual void MakeNewLevelFromCoarse(int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) override
Definition: REMORA_make_new_level.cpp:18
void stretch_transform(int lev)
Definition: DepthStretchTransform.H:12
void advance_2d(int lev, amrex::MultiFab const *mf_rhoS, amrex::MultiFab const *mf_rhoA, amrex::MultiFab *mf_ru, amrex::MultiFab *mf_rv, amrex::MultiFab *mf_rufrc, amrex::MultiFab *mf_rvfrc, amrex::MultiFab *mf_Zt_avg1, std::unique_ptr< amrex::MultiFab > &mf_DU_avg1, std::unique_ptr< amrex::MultiFab > &mf_DU_avg2, std::unique_ptr< amrex::MultiFab > &mf_DV_avg1, std::unique_ptr< amrex::MultiFab > &mf_DV_avg2, std::unique_ptr< amrex::MultiFab > &mf_rubar, std::unique_ptr< amrex::MultiFab > &mf_rvbar, std::unique_ptr< amrex::MultiFab > &mf_rzeta, std::unique_ptr< amrex::MultiFab > &mf_ubar, std::unique_ptr< amrex::MultiFab > &mf_vbar, amrex::MultiFab *mf_zeta, amrex::MultiFab const *mf_h, amrex::MultiFab const *mf_pm, amrex::MultiFab const *mf_pn, amrex::MultiFab const *mf_fcor, amrex::MultiFab const *mf_visc2_p, amrex::MultiFab const *mf_visc2_r, const amrex::Real dtfast_lev, const bool predictor_2d_step, const bool first_2d_step, int my_iif, int &next_indx1)
Definition: REMORA_advance_2d.cpp:37
void prestep_t_advection(const amrex::Box &tbx, const amrex::Box &gbx, const amrex::Array4< amrex::Real > &tempold, const amrex::Array4< amrex::Real > &tempcache, const amrex::Array4< amrex::Real > &Hz, const amrex::Array4< amrex::Real > &Huon, const amrex::Array4< amrex::Real > &Hvom, const amrex::Array4< amrex::Real > &W, const amrex::Array4< amrex::Real > &DC, const amrex::Array4< amrex::Real > &FC, const amrex::Array4< amrex::Real > &sstore, const amrex::Array4< amrex::Real const > &z_w, const amrex::Array4< amrex::Real const > &h, const amrex::Array4< amrex::Real const > &pm, const amrex::Array4< amrex::Real const > &pn, int iic, int ntfirst, int nrhs, int N, const amrex::Real dt_lev)
Definition: REMORA_prestep_t_advection.cpp:10
static int sum_interval
Definition: REMORA.H:835
amrex::Real stop_time
Definition: REMORA.H:776
static amrex::Real bndry_output_planes_start_time
Definition: REMORA.H:863
int do_substep
Definition: REMORA.H:793
amrex::Vector< std::string > plot_var_names_2
Definition: REMORA.H:810
void Evolve()
Definition: REMORA.cpp:142
void ReadCheckpointFile()
Definition: Checkpoint.cpp:214
void writeJobInfo(const std::string &dir) const
Definition: writeJobInfo.cpp:7
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_sustr
Definition: REMORA.H:246
virtual void ClearLevel(int lev) override
Definition: REMORA_make_new_level.cpp:327
amrex::Vector< std::string > datalogname
Definition: REMORA.H:968
amrex::Vector< amrex::MultiFab * > zvel_new
Definition: REMORA.H:191
int plot_int_1
Definition: REMORA.H:802
AMREX_FORCE_INLINE int ComputeGhostCells(const int &spatial_order)
Definition: REMORA.H:890
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_sstore
Definition: REMORA.H:291
static amrex::Real column_loc_x
Definition: REMORA.H:855
void init_only(int lev, amrex::Real time)
Definition: REMORA.cpp:447
void set_coriolis(int lev)
Definition: REMORA.cpp:396
amrex::Vector< amrex::Vector< amrex::Box > > boxes_at_level
Definition: REMORA.H:727
static amrex::Vector< amrex::AMRErrorTag > ref_tags
Definition: REMORA.H:880
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Hz
Definition: REMORA.H:200
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Akt
Definition: REMORA.H:216
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rvfrc
Definition: REMORA.H:212
static amrex::Real column_loc_y
Definition: REMORA.H:856
std::string check_file
Definition: REMORA.H:806
static amrex::Real startCPUTime
Definition: REMORA.H:936
void set_smflux(int lev, amrex::Real time)
Definition: REMORA.cpp:438
amrex::Vector< amrex::MultiFab * > xvel_old
Definition: REMORA.H:184
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_x_r
Definition: REMORA.H:225
amrex::Vector< amrex::MultiFab * > yvel_new
Definition: REMORA.H:190
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rufrc
Definition: REMORA.H:210
static amrex::Real fixed_fast_dt
Definition: REMORA.H:788
int regrid_int
Definition: REMORA.H:797
amrex::Vector< std::unique_ptr< amrex::MultiFab > > mapfac_u
Definition: REMORA.H:742
void uv3dmix(const amrex::Box &xbx, const amrex::Box &ybx, const amrex::Array4< amrex::Real > &u, const amrex::Array4< amrex::Real > &v, const amrex::Array4< amrex::Real const > &uold, const amrex::Array4< amrex::Real const > &vold, const amrex::Array4< amrex::Real > &rufrc, const amrex::Array4< amrex::Real > &rvfrc, const amrex::Array4< amrex::Real const > &visc2_p, const amrex::Array4< amrex::Real const > &visc2_r, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &pm, const amrex::Array4< amrex::Real const > &pn, int nrhs, int nnew, const amrex::Real dt_lev)
Definition: REMORA_uv3dmix.cpp:6
static int output_bndry_planes
Definition: REMORA.H:860
amrex::Vector< amrex::Real > vec_weight2
Definition: REMORA.H:299
int last_plot_file_step_2
Definition: REMORA.H:766
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_visc2_p
Definition: REMORA.H:218
amrex::MultiFab & build_fine_mask(int lev)
Definition: REMORA_SumIQ.cpp:115
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rvbar
Definition: REMORA.H:271
amrex::Vector< int > num_files_at_level
Definition: REMORA.H:726
void advance_2d_onestep(int lev, amrex::Real dt_lev, amrex::Real dtfast_lev, int my_iif, int nfast_counter)
Definition: REMORA_advance_2d_onestep.cpp:6
void AverageDownTo(int crse_lev)
Definition: REMORA.cpp:641
void post_timestep(int nstep, amrex::Real time, amrex::Real dt_lev)
Definition: REMORA.cpp:231
int max_step
Definition: REMORA.H:775
amrex::Vector< amrex::MultiFab * > zvel_old
Definition: REMORA.H:186
std::string plot_file_2
Definition: REMORA.H:801
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_z_r
Definition: REMORA.H:231
amrex::Vector< std::unique_ptr< amrex::MultiFab > > mapfac_m
Definition: REMORA.H:741
void WritePlotFile(int which, amrex::Vector< std::string > plot_var_names)
Definition: Plotfile.cpp:99
amrex::Vector< amrex::Real > h_havg_temperature
Definition: REMORA.H:869
amrex::Vector< int > num_boxes_at_level
Definition: REMORA.H:725
amrex::Vector< amrex::MultiFab * > xvel_new
Definition: REMORA.H:189
void refinement_criteria_setup()
Definition: REMORA_Tagging.cpp:39
void remora_advance(int level, amrex::MultiFab &cons_old, amrex::MultiFab &cons_new, amrex::MultiFab &xvel_old, amrex::MultiFab &yvel_old, amrex::MultiFab &zvel_old, amrex::MultiFab &xvel_new, amrex::MultiFab &yvel_new, amrex::MultiFab &zvel_new, amrex::MultiFab &source, const amrex::Geometry fine_geom, const amrex::Real dt, const amrex::Real time)
int last_check_file_step
Definition: REMORA.H:768
void rhs_uv_3d(const amrex::Box &xbx, const amrex::Box &ybx, const amrex::Array4< amrex::Real const > &uold, const amrex::Array4< amrex::Real const > &vold, const amrex::Array4< amrex::Real > &ru, const amrex::Array4< amrex::Real > &rv, const amrex::Array4< amrex::Real > &rufrc, const amrex::Array4< amrex::Real > &rvfrc, const amrex::Array4< amrex::Real const > &sustr, const amrex::Array4< amrex::Real const > &svstr, const amrex::Array4< amrex::Real const > &bustr, const amrex::Array4< amrex::Real const > &bvstr, const amrex::Array4< amrex::Real const > &Huon, const amrex::Array4< amrex::Real const > &Hvom, const amrex::Array4< amrex::Real const > &pm, const amrex::Array4< amrex::Real const > &pn, const amrex::Array4< amrex::Real const > &W, const amrex::Array4< amrex::Real > &FC, int nrhs, int N)
Definition: REMORA_rhs_uv_3d.cpp:31
void init_beta_plane_coriolis(int lev)
Definition: REMORA_init.cpp:54
void set_vmix(int lev)
Definition: REMORA.cpp:418
virtual void RemakeLevel(int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) override
Definition: REMORA_make_new_level.cpp:49
void vert_mean_3d(const amrex::Box &bx, const int ioff, const int joff, const amrex::Array4< amrex::Real > &phi, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &Dphi_avg1, const amrex::Array4< amrex::Real > &DC, const amrex::Array4< amrex::Real > &CF, const amrex::Array4< amrex::Real const > &dxlen, const int nnew, const int N)
Definition: REMORA_vert_mean_3d.cpp:10
amrex::Vector< int > nsubsteps
Definition: REMORA.H:730
amrex::Vector< std::unique_ptr< REMORAPhysBCFunct > > physbcs
Definition: REMORA.H:737
void ComputeDt()
Definition: REMORA_ComputeTimestep.cpp:8
void WriteGenericPlotfileHeaderWithBathymetry(std::ostream &HeaderFile, int nlevels, const amrex::Vector< amrex::BoxArray > &bArray, const amrex::Vector< std::string > &varnames, amrex::Real time, const amrex::Vector< int > &level_steps, const std::string &versionName, const std::string &levelPrefix, const std::string &mfPrefix) const
Definition: Plotfile.cpp:390
static void writeBuildInfo(std::ostream &os)
Definition: writeJobInfo.cpp:135
void init_custom(int lev)
Definition: REMORA_init.cpp:13
virtual void ErrorEst(int lev, amrex::TagBoxArray &tags, amrex::Real time, int ngrow) override
Definition: REMORA_Tagging.cpp:15
std::string plot_file_1
Definition: REMORA.H:800
amrex::Real getCPUTime() const
Definition: REMORA.H:940
void InitData()
Definition: REMORA.cpp:260
amrex::Vector< int > istep
Definition: REMORA.H:729
void t3dmix(const amrex::Box &bx, const amrex::Array4< amrex::Real > &state, const amrex::Array4< amrex::Real const > &diff2, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &pm, const amrex::Array4< amrex::Real const > &pn, const amrex::Real dt_lev, const int ncomp)
Definition: REMORA_t3dmix.cpp:6
void setRecordDataInfo(int i, const std::string &filename)
Definition: REMORA.H:954
void advance_3d_ml(int lev, amrex::Real dt_lev)
Definition: REMORA_advance_3d_ml.cpp:6
const std::string DataLogName(int i) const noexcept
The filename of the ith datalog file.
Definition: REMORA.H:971
amrex::Gpu::DeviceVector< amrex::Real > d_havg_temperature
Definition: REMORA.H:872
amrex::Vector< amrex::Real > h_havg_pressure
Definition: REMORA.H:870
AMREX_FORCE_INLINE std::ostream & DataLog(int i)
Definition: REMORA.H:924
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rhoS
Definition: REMORA.H:295
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_visc2_r
Definition: REMORA.H:220
amrex::Vector< std::unique_ptr< amrex::MultiFab > > sst
Definition: REMORA.H:745
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_svstr
Definition: REMORA.H:249
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Huon
Definition: REMORA.H:202
REMORA()
Definition: REMORA.cpp:61
static amrex::Real change_max
Definition: REMORA.H:784
void WriteMultiLevelPlotfileWithBathymetry(const std::string &plotfilename, int nlevels, const amrex::Vector< const amrex::MultiFab * > &mf, const amrex::Vector< const amrex::MultiFab * > &mf_nd, const amrex::Vector< std::string > &varnames, amrex::Real time, const amrex::Vector< int > &level_steps, const std::string &versionName="HyperCLaw-V1.1", const std::string &levelPrefix="Level_", const std::string &mfPrefix="Cell", const amrex::Vector< std::string > &extra_dirs=amrex::Vector< std::string >()) const
Definition: Plotfile.cpp:305
void FillCoarsePatch(int lev, amrex::Real time, amrex::MultiFab *mf_fine, amrex::MultiFab *mf_crse)
Definition: REMORA_FillPatch.cpp:176
std::string pp_prefix
Definition: REMORA.H:181
void set_bathymetry(int lev)
Definition: REMORA.cpp:372
static bool use_tracer_particles
Definition: REMORA.H:832
void initHSE()
Initialize HSE.
amrex::Vector< amrex::MultiFab * > yvel_old
Definition: REMORA.H:185
void init_stuff(int lev, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm)
Definition: REMORA_make_new_level.cpp:191
TimeInterpolatedData GetDataAtTime(int lev, amrex::Real time)
Definition: REMORA_FillPatch.cpp:115
void Advance(int lev, amrex::Real time, amrex::Real dt_lev, int iteration, int ncycle)
Definition: REMORA_Advance.cpp:7
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rhoA
Definition: REMORA.H:296
amrex::Gpu::DeviceVector< amrex::Real > d_havg_pressure
Definition: REMORA.H:873
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_DV_avg1
Definition: REMORA.H:265
amrex::Vector< amrex::Real > h_havg_density
Definition: REMORA.H:868
int last_plot_file_step_1
Definition: REMORA.H:765
amrex::Vector< std::unique_ptr< amrex::MultiFab > > mapfac_v
Definition: REMORA.H:743
amrex::Vector< amrex::YAFluxRegister * > advflux_reg
Definition: REMORA.H:748
amrex::Vector< amrex::Real > t_new
Definition: REMORA.H:733
amrex::Gpu::DeviceVector< amrex::Real > d_havg_density
Definition: REMORA.H:871
static SolverChoice solverChoice
Definition: REMORA.H:822
void resize_stuff(int lev)
Definition: REMORA_make_new_level.cpp:137
void ReadParameters()
Definition: REMORA.cpp:484
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_visc3d_r
Definition: REMORA.H:293
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ru
Definition: REMORA.H:206
void sum_integrated_quantities(amrex::Real time)
Definition: REMORA_SumIQ.cpp:9
amrex::Vector< amrex::Real > vec_weight1
Definition: REMORA.H:298
void fill_from_bdyfiles(amrex::MultiFab &mf_to_fill, const amrex::Real time, const int bdy_var_type)
void prestep(int lev, amrex::MultiFab &mf_uold, amrex::MultiFab &mf_vold, amrex::MultiFab &mf_u, amrex::MultiFab &mf_v, amrex::MultiFab *mf_ru, amrex::MultiFab *mf_rv, amrex::MultiFab &S_old, amrex::MultiFab &S_new, amrex::MultiFab &mf_W, amrex::MultiFab &mf_DC, const amrex::MultiFab *mf_z_r, const amrex::MultiFab *mf_z_w, const amrex::MultiFab *mf_h, const amrex::MultiFab *mf_pm, const amrex::MultiFab *mf_pn, const amrex::MultiFab *mf_sustr, const amrex::MultiFab *mf_svstr, const amrex::MultiFab *mf_bustr, const amrex::MultiFab *mf_bvstr, const int iic, const int nfirst, const int nnew, int nstp, int nrhs, int N, const amrex::Real dt_lev)
Definition: REMORA_prestep.cpp:34
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_zeta
Definition: REMORA.H:279
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_vbar
Definition: REMORA.H:277
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_y_r
Definition: REMORA.H:228
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_DU_avg1
Definition: REMORA.H:261
int plot_file_on_restart
Definition: REMORA.H:769
void set_2darrays(int lev)
Definition: REMORA_init.cpp:82
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_ubar
Definition: REMORA.H:275
void fill_rhs(amrex::MultiFab &rhs_mf, const amrex::MultiFab &state_mf, const amrex::Real time, const amrex::Geometry &geom)
void InitializeLevelFromData(int lev, const amrex::MultiFab &initial_data)
amrex::Vector< amrex::MultiFab * > cons_old
Definition: REMORA.H:183
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_bustr
Definition: REMORA.H:255
AMREX_FORCE_INLINE int NumDataLogs() noexcept
Definition: REMORA.H:931
amrex::Real volWgtSumMF(int lev, const amrex::MultiFab &mf, int comp, bool local, bool finemask)
Definition: REMORA_SumIQ.cpp:89
void vert_visc_3d(const amrex::Box &bx, const int ioff, const int joff, const amrex::Array4< amrex::Real > &phi, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real > &Hzk, const amrex::Array4< amrex::Real > &AK, const amrex::Array4< amrex::Real > &Akv, const amrex::Array4< amrex::Real > &BC, const amrex::Array4< amrex::Real > &DC, const amrex::Array4< amrex::Real > &FC, const amrex::Array4< amrex::Real > &CF, const int nnew, const int N, const amrex::Real dt_lev)
Definition: REMORA_vert_visc_3d.cpp:10
bool is_it_time_for_action(int nstep, amrex::Real time, amrex::Real dt, int action_interval, amrex::Real action_per)
Definition: REMORA_SumIQ.cpp:143
std::string PlotFileName(int lev) const
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_DV_avg2
Definition: REMORA.H:267
void set_hmixcoef(int lev)
Definition: REMORA.cpp:427
amrex::Array< std::string, 2 *AMREX_SPACEDIM > domain_bc_type
Definition: REMORA.H:757
static std::string column_file_name
Definition: REMORA.H:857
amrex::Real estTimeStep(int lev) const
Definition: REMORA_ComputeTimestep.cpp:40
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_bvstr
Definition: REMORA.H:258
void timeStep(int lev, amrex::Real time, int iteration)
Definition: REMORA_TimeStep.cpp:9
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rzeta
Definition: REMORA.H:273
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_z_phys_nd
Definition: REMORA.H:240
void AverageDown()
Definition: REMORA.cpp:631
static int fixed_ndtfast_ratio
Definition: REMORA.H:789
void update_massflux_3d(const amrex::Box &bx, const int ioff, const int joff, const amrex::Array4< amrex::Real > &phi, const amrex::Array4< amrex::Real > &phibar, const amrex::Array4< amrex::Real > &Hphi, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &pm_or_pn, const amrex::Array4< amrex::Real const > &Dphi1, const amrex::Array4< amrex::Real const > &Dphi2, const amrex::Array4< amrex::Real > &DC, const amrex::Array4< amrex::Real > &FC, const int nnew)
Definition: REMORA_update_massflux_3d.cpp:23
void timeStepML(amrex::Real time, int iteration)
Definition: REMORA_TimeStepML.cpp:9
static int input_bndry_planes
Definition: REMORA.H:866
void set_weights(int lev)
Definition: REMORA_set_weights.cpp:10
static int output_1d_column
Definition: REMORA.H:852
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_pn
Definition: REMORA.H:285
amrex::Vector< std::unique_ptr< std::fstream > > datalog
Definition: REMORA.H:967
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Akv
Definition: REMORA.H:214
static amrex::Real cfl
Definition: REMORA.H:782
void post_update(amrex::MultiFab &state_mf, const amrex::Real time, const amrex::Geometry &geom)
virtual void MakeNewLevelFromScratch(int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) override
Definition: REMORA_make_new_level.cpp:102
static int verbose
Definition: REMORA.H:829
void setup_step(int lev, amrex::Real time, amrex::Real dt_lev)
Definition: REMORA_setup_step.cpp:7
amrex::GpuArray< REMORA_BC, AMREX_SPACEDIM *2 > phys_bc_type
Definition: REMORA.H:763
void initRayleigh()
Initialize Rayleigh damping profiles.
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rdrag
Definition: REMORA.H:252
static std::string nc_bdry_file
Definition: REMORA.H:846
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Zt_avg1
Definition: REMORA.H:243
amrex::Vector< std::string > PlotFileVarNames(amrex::Vector< std::string > plot_var_names) const
Definition: Plotfile.cpp:87
int check_int
Definition: REMORA.H:807
std::string restart_chkfile
Definition: REMORA.H:779
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_rv
Definition: REMORA.H:208
const amrex::Vector< std::string > derived_names
Definition: REMORA.H:817
amrex::Vector< amrex::Real > t_old
Definition: REMORA.H:734
void FillPatch(int lev, amrex::Real time, amrex::MultiFab &mf_to_be_filled, amrex::Vector< amrex::MultiFab * > const &mfs, const int bdy_var_type=BdyVars::null)
Definition: REMORA_FillPatch.cpp:16
void WriteCheckpointFile() const
Definition: Checkpoint.cpp:15
void rho_eos(const amrex::Box &bx, const amrex::Array4< amrex::Real const > &state, const amrex::Array4< amrex::Real > &rho, const amrex::Array4< amrex::Real > &rhoA, const amrex::Array4< amrex::Real > &rhoS, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &z_w, const amrex::Array4< amrex::Real const > &h, const int N)
Definition: REMORA_rho_eos.cpp:20
amrex::Vector< amrex::Real > dt
Definition: REMORA.H:735
static amrex::Real sum_per
Definition: REMORA.H:836
static std::string input_sounding_file
Definition: REMORA.H:849
static amrex::Vector< amrex::Vector< std::string > > nc_grid_file
Definition: REMORA.H:843
void InitializeFromFile()
void coriolis(const amrex::Box &xbx, const amrex::Box &ybx, const amrex::Array4< amrex::Real const > &uold, const amrex::Array4< amrex::Real const > &vold, const amrex::Array4< amrex::Real > &ru, const amrex::Array4< amrex::Real > &rv, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &fomn, int nrhs, int nr)
Definition: REMORA_coriolis.cpp:10
static amrex::Real bndry_output_planes_per
Definition: REMORA.H:862
amrex::Gpu::DeviceVector< amrex::BCRec > domain_bcs_type_d
Definition: REMORA.H:754
virtual ~REMORA()
Definition: REMORA.cpp:136
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_z_w
Definition: REMORA.H:234
static amrex::Vector< amrex::Vector< std::string > > nc_init_file
Definition: REMORA.H:842
void restart()
Definition: REMORA.cpp:363
void prsgrd(const amrex::Box &bx, const amrex::Box &gbx, const amrex::Box &utbx, const amrex::Box &vtbx, const amrex::Array4< amrex::Real > &ru, const amrex::Array4< amrex::Real > &rv, const amrex::Array4< amrex::Real const > &pn, const amrex::Array4< amrex::Real const > &pm, const amrex::Array4< amrex::Real const > &rho, const amrex::Array4< amrex::Real > &FC, const amrex::Array4< amrex::Real const > &Hz, const amrex::Array4< amrex::Real const > &z_r, const amrex::Array4< amrex::Real const > &z_w, const int nrhs, const int N)
Definition: REMORA_prsgrd.cpp:6
static amrex::Real column_per
Definition: REMORA.H:854
void advance_3d(int lev, amrex::MultiFab &mf_cons, amrex::MultiFab &mf_u, amrex::MultiFab &mf_v, amrex::MultiFab *mf_sstore, amrex::MultiFab *mf_ru, amrex::MultiFab *mf_rv, std::unique_ptr< amrex::MultiFab > &mf_DU_avg1, std::unique_ptr< amrex::MultiFab > &mf_DU_avg2, std::unique_ptr< amrex::MultiFab > &mf_DV_avg1, std::unique_ptr< amrex::MultiFab > &mf_DV_avg2, std::unique_ptr< amrex::MultiFab > &mf_ubar, std::unique_ptr< amrex::MultiFab > &mf_vbar, std::unique_ptr< amrex::MultiFab > &mf_Akv, std::unique_ptr< amrex::MultiFab > &mf_Akt, std::unique_ptr< amrex::MultiFab > &mf_Hz, std::unique_ptr< amrex::MultiFab > &mf_Huon, std::unique_ptr< amrex::MultiFab > &mf_Hvom, std::unique_ptr< amrex::MultiFab > &mf_z_w, amrex::MultiFab const *mf_h, amrex::MultiFab const *mf_pm, amrex::MultiFab const *mf_pn, const int N, const amrex::Real dt_lev)
Definition: REMORA_advance_3d.cpp:10
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_s_r
Definition: REMORA.H:237
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_Hvom
Definition: REMORA.H:204
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_diff2
Definition: REMORA.H:222
void setPlotVariables(const std::string &pp_plot_var_names, amrex::Vector< std::string > &plot_var_names)
Definition: Plotfile.cpp:9
Definition: NCInterface.H:285
@ u
Definition: IndexDefines.H:34
@ v
Definition: IndexDefines.H:35
@ null
Definition: IndexDefines.H:33
@ t
Definition: IndexDefines.H:36
@ CellConservativeQuartic
Definition: REMORA.H:58
@ CellBilinear
Definition: REMORA.H:55
@ PCInterp
Definition: REMORA.H:52
@ NodeBilinear
Definition: REMORA.H:53
@ CellConservativeLinear
Definition: REMORA.H:54
@ CellQuadratic
Definition: REMORA.H:56
@ CellConservativeProtected
Definition: REMORA.H:57
Definition: ParticleData.H:19
Definition: DataStruct.H:41
Definition: TimeInterpolatedData.H:8