forked from seahorce-scidac/REMORA
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathREMORA.H
976 lines (788 loc) · 39 KB
/
REMORA.H
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
#ifndef REMORA_H_
#define REMORA_H_
#include <string>
#include <limits>
#include <memory>
#ifdef _OPENMP
#include <omp.h>
#endif
#include <AMReX_AmrCore.H>
#include <AMReX_BCRec.H>
#include <AMReX_InterpFaceRegister.H>
#include <AMReX_ParallelDescriptor.H>
#include <AMReX_ParmParse.H>
#include <AMReX_MultiFabUtil.H>
#include <AMReX_FillPatchUtil.H>
#include <AMReX_VisMF.H>
#include <AMReX_PhysBCFunct.H>
#include <AMReX_YAFluxRegister.H>
#include <AMReX_ErrorList.H>
#ifdef AMREX_MEM_PROFILING
#include <AMReX_MemProfiler.H>
#endif
#include <REMORA_Math.H>
#include <REMORA_PhysBCFunct.H>
#include <REMORA_FillPatcher.H>
#include <IndexDefines.H>
#include <TimeInterpolatedData.H>
#include <DataStruct.H>
#include <Derive.H>
#ifdef REMORA_USE_PARTICLES
#include "ParticleData.H"
#endif
#ifdef REMORA_USE_NETCDF
#include "NCInterface.H"
#endif
#include <iostream>
#ifdef AMREX_LAZY
#include <AMReX_Lazy.H>
#endif
namespace InterpType {
enum {
PCInterp = 0,
NodeBilinear,
CellConservativeLinear,
CellBilinear,
CellQuadratic,
CellConservativeProtected,
CellConservativeQuartic
};
}
template<typename V, typename T>
bool containerHasElement(const V& iterable, const T& query) {
return std::find(iterable.begin(), iterable.end(), query) != iterable.end();
}
class REMORA
: public amrex::AmrCore
{
public:
////////////////
// public member functions
/**
* constructor - reads in parameters from inputs file
* - sizes multilevel arrays and data structures
*/
REMORA ();
virtual ~REMORA();
/** Advance solution to final time */
void Evolve ();
/** Tag cells for refinement */
virtual void ErrorEst (int lev, amrex::TagBoxArray& tags, amrex::Real time, int ngrow) override;
/** Initialize multilevel data */
void InitData ();
/** Init (NOT restart or regrid) */
void init_only (int lev, amrex::Real time);
/** Restart */
void restart ();
/** Called after every level 0 timestep */
void post_timestep (int nstep, amrex::Real time, amrex::Real dt_lev);
// Diagnostics
/** Integrate conserved quantities for diagnostics */
void sum_integrated_quantities(amrex::Real time);
/* Perform the volume-weighted sum */
amrex::Real
volWgtSumMF(int lev,
const amrex::MultiFab& mf, int comp, bool local, bool finemask);
/* Decide if it is time to take an action */
bool is_it_time_for_action(int nstep, amrex::Real time, amrex::Real dt,
int action_interval, amrex::Real action_per);
/** Make a new level using provided BoxArray and DistributionMapping and
* fill with interpolated coarse level data.
* Overrides the pure virtual function in AmrCore
*/
virtual void MakeNewLevelFromCoarse (int lev, amrex::Real time, const amrex::BoxArray& ba,
const amrex::DistributionMapping& dm) override;
/** Remake an existing level using provided BoxArray and DistributionMapping and
* fill with existing fine and coarse data.
* Overrides the pure virtual function in AmrCore
*/
virtual void RemakeLevel (int lev, amrex::Real time, const amrex::BoxArray& ba,
const amrex::DistributionMapping& dm) override;
/** Delete level data
* Overrides the pure virtual function in AmrCore
*/
virtual void ClearLevel (int lev) override;
/** Make a new level from scratch using provided BoxArray and DistributionMapping.
* Only used during initialization.
* Overrides the pure virtual function in AmrCore
*/
virtual void MakeNewLevelFromScratch (int lev, amrex::Real time, const amrex::BoxArray& ba,
const amrex::DistributionMapping& dm) override;
/** compute dt from CFL considerations */
amrex::Real estTimeStep (int lev) const;
/** Interface for advancing the data at one level by one "slow" timestep */
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);
/** Make mask to zero out covered cells*/
amrex::MultiFab& build_fine_mask(int lev);
/** write plotfile */
void WritePlotFile ();
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;
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;
std::string pp_prefix {"remora"};
amrex::Vector<amrex::MultiFab*> cons_old;
amrex::Vector<amrex::MultiFab*> xvel_old;
amrex::Vector<amrex::MultiFab*> yvel_old;
amrex::Vector<amrex::MultiFab*> zvel_old;
amrex::Vector<amrex::MultiFab*> cons_new;
amrex::Vector<amrex::MultiFab*> xvel_new;
amrex::Vector<amrex::MultiFab*> yvel_new;
amrex::Vector<amrex::MultiFab*> zvel_new;
// Program state data is represented by vectors of pointers to AMReX Multifabs.
// There is one pointer per level
/** Bathymetry data (2D, positive valued, h in ROMS) **/
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_hOfTheConfusingName;
/** Width of cells in the vertical (z-) direction (3D, Hz in ROMS) */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_Hz;
/** u-volume flux (3D) */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_Huon;
/** v-volume flux (3D) */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_Hvom;
/** u velocity RHS (3D, includes horizontal and vertical advection) */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_ru;
/** v velocity RHS (3D, includes horizontal and vertical advection) */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_rv;
/** u velocity RHS, integrated, including advection and bottom/surface stresses (2D) */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_rufrc;
/** v velocity RHS, integrated, including advection and bottom/surface stresses (2D) */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_rvfrc;
/** Vertical viscosity coefficient (3D), set in .in file */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_Akv;
/** Vertical diffusion coefficient (3D), set in .in file */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_Akt;
/** Harmonic viscosity devined on the psi points (corners of horizontal grid cells) */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_visc2_p;
/** Harmonic viscosity devined on the rho points (centers) */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_visc2_r;
/** Harmonic diffusivity for temperature / salinity */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_diff2;
/** x coordinates at rho points (cell centers) */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_x_r;
/** y coordinates at rho points (cell centers) */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_y_r;
/** z coordinates at rho points (cell centers) */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_z_r;
/** z coordinates at w points (faces between z-cells) */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_z_w;
/** Scaled vertical coordinate (range [0,1]) that transforms to z, defined at rho points (cell centers) */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_s_r;
/** z coordinates at psi points (cell nodes) **/
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_z_phys_nd;
/** Average of the free surface, zeta (2D) */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_Zt_avg1;
/** Surface stress in the u direction */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_sustr;
/** Surface stress in the v direction */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_svstr;
/** Linear drag coefficient [m/s], defined at rho points */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_rdrag;
/** Bottom stress in the u direction */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_bustr;
/** Bottom stress in the v direction */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_bvstr;
/** time average of barotropic x velocity flux (2D) */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_DU_avg1;
/** correct time average of barotropic x velocity flux for coupling (2D) */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_DU_avg2;
/** time average of barotropic y velocity flux */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_DV_avg1;
/** correct time average of barotropic y velocity flux for coupling (2D) */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_DV_avg2;
/** barotropic x velocity for the RHS (2D) */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_rubar;
/** barotropic y velocity for the RHS (2D) */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_rvbar;
/** free surface height for the RHS (2D) */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_rzeta;
/** barotropic x velocity (2D) */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_ubar;
/** barotropic y velocity (2D) */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_vbar;
/** free surface height (2D) */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_zeta;
/** land/sea mask at cell centers (2D) */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_mskr;
/** land/sea mask at x-faces (2D) */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_msku;
/** land/sea mask at y-faces (2D) */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_mskv;
/** map factor (2D) */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_pm;
/** map factor (2D) */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_pn;
/** coriolis factor (2D) */
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_fcor;
/** saltstore, tempstore, etc*/
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_sstore;
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_visc3d_r;
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_rhoS;
amrex::Vector<std::unique_ptr<amrex::MultiFab>> vec_rhoA;
amrex::Vector<amrex::Real> vec_weight1;
amrex::Vector<amrex::Real> vec_weight2;
/** advance a single level for a single time step */
void Advance (int lev, amrex::Real time, amrex::Real dt_lev, int iteration, int ncycle);
/** Set everything up for a step on a level */
void setup_step(int lev, amrex::Real time, amrex::Real dt_lev);
/** 3D advance on a single level */
void advance_3d_ml (int lev, amrex::Real dt_lev);
/** 2D advance, one predictor/corrector step */
void advance_2d_onestep (int lev, amrex::Real dt_lev, amrex::Real dtfast_lev, int my_iif, int nfast_counter);
/** Perform a 2D predictor (predictor_2d_step=True) or corrector (predictor_2d_step=False) step */
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);
/** Advance the 3D variables */
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);
/** Wrapper function for prestep */
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,
/* MF mf_FC? */
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);
/** Prestep calculations for the tracers */
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);
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);
/** Calculation of the RHS */
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);
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);
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);
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);
/** Update velocities or tracers with diffusion/viscosity
* as the last part of the prestep */
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);
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);
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);
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);
/** Harmonic viscosity */
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);
/** Harmonic diffusivity for tracers */
void t3dmix (const amrex::Box& bx,
const amrex::Array4<amrex::Real >& state,
const amrex::Array4<amrex::Real >& state_rhs,
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);
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);
void set_2darrays (int lev);
void set_bathymetry (int lev);
void set_coriolis (int lev);
void stretch_transform (int lev);
void set_vmix (int lev);
void set_smflux (int lev, amrex::Real time);
void set_hmixcoef (int lev);
void set_drag (int lev);
void set_weights (int lev);
// Fill a new MultiFab by copying in phi from valid region and filling ghost cells
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,
const int icomp=0,
const bool fill_all=true);
void FillPatchNoBC (int lev, amrex::Real time,
amrex::MultiFab& mf_to_be_filled,
amrex::Vector<amrex::MultiFab*> const& mfs,
const int bdy_var_type = BdyVars::null,
const int icomp=0,
const bool fill_all=true);
void fill_from_bdyfiles (amrex::MultiFab& mf_to_fill,
const amrex::Real time,
const int bdy_var_type);
void init_beta_plane_coriolis (int lev);
void init_data_from_netcdf (int lev);
void init_bdry_from_netcdf ();
void init_bathymetry_from_netcdf (int lev);
void init_coriolis_from_netcdf (int lev);
#ifdef REMORA_USE_NETCDF
static int total_nc_plot_file_step;
#endif
private:
///////////////////////////
// private member functions
///////////////////////////
/** read in some parameters from inputs file */
void ReadParameters();
/** set covered coarse cells to be the average of overlying fine cells */
void AverageDown ();
void init_bcs();
void init_custom(int lev);
void init_stuff (int lev, const amrex::BoxArray& ba, const amrex::DistributionMapping& dm);
void resize_stuff (int lev);
/** more flexible version of AverageDown() that lets you average down across multiple levels */
void AverageDownTo (int crse_lev);
void Construct_REMORAFillPatchers (int lev);
void Define_REMORAFillPatchers (int lev);
/** fill an entire multifab by interpolating from the coarser level
* this comes into play when a new level of refinement appears
*/
void FillCoarsePatch (int lev, amrex::Real time, amrex::MultiFab* mf_fine, amrex::MultiFab* mf_crse);
/** utility to copy in data from old and/or new state into another multifab */
TimeInterpolatedData GetDataAtTime (int lev, amrex::Real time);
/** advance a level by dt, includes a recursive call for finer levels */
void timeStep (int lev, amrex::Real time, int iteration);
/** advance all levels by dt, loops over finer levels */
void timeStepML (amrex::Real time, int iteration);
/** a wrapper for estTimeStep() */
void ComputeDt ();
/** get plotfile name */
std::string PlotFileName (int lev) const;
/** set which variables and derived quantities go into plotfiles */
void setPlotVariables (const std::string& pp_plot_var_names);
#ifdef REMORA_USE_NETCDF
//! Write plotfile using NETCDF
void WriteNCPlotFile (int istep) const;
void WriteNCPlotFile_which (int lev, int which_subdomain,
bool write_header, ncutils::NCFile& ncf) const;
//! Write MultiFab in NetCDF format
void WriteNCMultiFab (const amrex::FabArray<amrex::FArrayBox> &fab,
const std::string& name,
bool set_ghost = false) const;
//! Read MultiFab in NetCDF format
void ReadNCMultiFab (amrex::FabArray<amrex::FArrayBox> &fab,
const std::string &name,
int coordinatorProc = amrex::ParallelDescriptor::IOProcessorNumber(),
int allow_empty_mf = 0);
// Vectors (over time) of Vector (over variables) of FArrayBoxs for holding the data read from the wrfbdy NetCDF file
amrex::Vector<amrex::Vector<amrex::FArrayBox>> bdy_data_xlo;
amrex::Vector<amrex::Vector<amrex::FArrayBox>> bdy_data_xhi;
amrex::Vector<amrex::Vector<amrex::FArrayBox>> bdy_data_ylo;
amrex::Vector<amrex::Vector<amrex::FArrayBox>> bdy_data_yhi;
amrex::Real start_bdy_time;
amrex::Real bdy_time_interval;
static bool write_history_file;
int bdy_width{0};
int bdy_set_width{0};
#endif // REMORA_USE_NETCDF
// Fillpatcher classes for coarse-fine boundaries
int cf_width{0};
int cf_set_width{0};
amrex::Vector<REMORAFillPatcher> FPr_c;
amrex::Vector<REMORAFillPatcher> FPr_u;
amrex::Vector<REMORAFillPatcher> FPr_v;
amrex::Vector<REMORAFillPatcher> FPr_w;
// write checkpoint file to disk
void WriteCheckpointFile () const;
// read checkpoint file from disk
void ReadCheckpointFile ();
// Read the file passed to amr.restart and use it as an initial condition for
// the current simulation. Supports a different number of components and
// ghost cells.
void InitializeFromFile ();
// Initialize the new-time data at a level from the initial_data MultiFab
void InitializeLevelFromData (int lev, const amrex::MultiFab& initial_data);
// utility to skip to next line in Header
static void GotoNextLine (std::istream& is);
// Single level functions called by advance()
void post_update (amrex::MultiFab& state_mf, const amrex::Real time, const amrex::Geometry& geom);
void fill_rhs (amrex::MultiFab& rhs_mf, const amrex::MultiFab& state_mf, const amrex::Real time, const amrex::Geometry& geom);
////////////////
// private data members
amrex::Vector<int> num_boxes_at_level; // how many boxes specified at each level by tagging criteria
amrex::Vector<int> num_files_at_level; // how many wrfinput files specified at each level
amrex::Vector<amrex::Vector<amrex::Box>> boxes_at_level; // the boxes specified at each level by tagging criteria
amrex::Vector<int> istep; // which step?
amrex::Vector<int> nsubsteps; // how many substeps on each level?
// keep track of old time, new time, and time step at each level
amrex::Vector<amrex::Real> t_new;
amrex::Vector<amrex::Real> t_old;
amrex::Vector<amrex::Real> dt;
amrex::Vector<std::unique_ptr<REMORAPhysBCFunct>> physbcs;
// array of multifabs to store the solution at each level of refinement
// after advancing a level we use "swap".
amrex::Vector<std::unique_ptr<amrex::MultiFab>> mapfac_m;
amrex::Vector<std::unique_ptr<amrex::MultiFab>> mapfac_u;
amrex::Vector<std::unique_ptr<amrex::MultiFab>> mapfac_v;
amrex::Vector<std::unique_ptr<amrex::MultiFab>> sst;
// array of flux registers
amrex::Vector<amrex::YAFluxRegister*> advflux_reg;
// A BCRec is essentially a 2*DIM integer array storing the boundary
// condition type at each lo/hi walls in each direction. We have one BCRec
// for each component of the cell-centered variables and each velocity component.
amrex::Vector <amrex::BCRec> domain_bcs_type;
amrex::Gpu::DeviceVector<amrex::BCRec> domain_bcs_type_d;
// We store these so that we can print them out in the job_info file
amrex::Array<std::string,2*AMREX_SPACEDIM> domain_bc_type;
// These hold the Dirichlet values at walls which need them ...
amrex::Array<amrex::Array<amrex::Real, AMREX_SPACEDIM*2>,AMREX_SPACEDIM+NCONS> m_bc_extdir_vals;
// These are the "physical" boundary condition types (e.g. "inflow")
amrex::GpuArray<REMORA_BC, AMREX_SPACEDIM*2> phys_bc_type;
int last_plot_file_step;
int last_check_file_step;
int plot_file_on_restart = 1;
////////////////
// runtime parameters
// maximum number of steps and stop time
int max_step = std::numeric_limits<int>::max();
amrex::Real stop_time = std::numeric_limits<amrex::Real>::max();
// if >= 0 we restart from a checkpoint
std::string restart_chkfile = "";
// Time step controls
static amrex::Real cfl;
static amrex::Real init_shrink;
static amrex::Real change_max;
// Fixed dt for level 0 timesteps (only used if positive)
static amrex::Real fixed_dt;
static amrex::Real fixed_fast_dt;
static int fixed_ndtfast_ratio;
int nfast;
// Whether to substep fine levels in time
int do_substep;
// how often each level regrids the higher levels of refinement
// (after a level advances that many time steps)
int regrid_int = 2;
// plotfile prefix and frequency
std::string plot_file_name {"plt_"};
int plot_int = -1;
// Checkpoint type, prefix and frequency
std::string check_file {"chk"};
int check_int = -1;
amrex::Vector<std::string> plot_var_names;
const amrex::Vector<std::string> cons_names {"temp", "salt", "scalar"};
// Note that the order of variable names here must match the order in Derive.cpp
#ifdef REMORA_USE_PARTICLES
const amrex::Vector<std::string> derived_names {"dpdx", "dpdy" ,"tracer_particle_count"};
#else
const amrex::Vector<std::string> derived_names {"dpdx", "dpdy"};
#endif
// algorithm choices
static SolverChoice solverChoice;
#ifdef REMORA_USE_PARTICLES
// Particle info
static ParticleData particleData;
#endif
static int verbose;
// flag to turn tracer particle generation on/off
static bool use_tracer_particles;
// Diagnostic output interval
static int sum_interval;
static amrex::Real sum_per;
// Native or NetCDF
static PlotfileType plotfile_type;
// NetCDF initialization files for solution and grid
static amrex::Vector<amrex::Vector<std::string>> nc_init_file;
static amrex::Vector<amrex::Vector<std::string>> nc_grid_file;
// NetCDF file of time series of boundary data
static std::string nc_bdry_file;
void refinement_criteria_setup();
//
// Holds info for dynamically generated tagging criteria
//
static amrex::Vector<amrex::AMRErrorTag> ref_tags;
//
// Build a mask that zeroes out values on a coarse level underlying
// grids on the next finest level
//
amrex::MultiFab fine_mask;
AMREX_FORCE_INLINE
int
ComputeGhostCells(const int& spatial_order) {
int nGhostCells = 0;
switch (spatial_order) {
case 2:
nGhostCells = 2; // We need this many to compute the eddy viscosity in the ghost cells
break;
case 3:
nGhostCells = 2;
break;
case 4:
nGhostCells = 2;
break;
case 5:
nGhostCells = 3;
break;
case 6:
nGhostCells = 3;
break;
default:
amrex::Error("Must specify spatial order to be 2,3,4,5 or 6");
}
return nGhostCells;
}
AMREX_FORCE_INLINE
amrex::YAFluxRegister* getAdvFluxReg (int lev)
{
return advflux_reg[lev];
}
AMREX_FORCE_INLINE
std::ostream&
DataLog (int i)
{
return *datalog[i];
}
AMREX_FORCE_INLINE
int
NumDataLogs () noexcept
{
return datalog.size();
}
static amrex::Real startCPUTime;
static amrex::Real previousCPUTimeUsed;
amrex::Real
getCPUTime() const
{
int numCores = amrex::ParallelDescriptor::NProcs();
#ifdef _OPENMP
numCores = numCores * omp_get_max_threads();
#endif
amrex::Real T =
numCores * (amrex::Real(amrex::ParallelDescriptor::second()) - startCPUTime) +
previousCPUTimeUsed;
return T;
}
void setRecordDataInfo (int i, const std::string& filename)
{
if (amrex::ParallelDescriptor::IOProcessor())
{
datalog[i] = std::make_unique<std::fstream>();
datalog[i]->open(filename.c_str(),std::ios::out|std::ios::app);
if (!datalog[i]->good()) {
amrex::FileOpenFailed(filename);
}
}
amrex::ParallelDescriptor::Barrier("REMORA::setRecordDataInfo");
}
amrex::Vector<std::unique_ptr<std::fstream> > datalog;
amrex::Vector<std::string> datalogname;
//! The filename of the ith datalog file.
const std::string DataLogName (int i) const noexcept { return datalogname[i]; }
public:
void writeJobInfo (const std::string& dir) const;
static void writeBuildInfo (std::ostream& os);
static void print_banner(MPI_Comm /*comm*/, std::ostream& /*out*/);
static void print_usage(MPI_Comm /*comm*/, std::ostream& /*out*/);
static void print_error(MPI_Comm /*comm*/, const std::string& msg);
static void print_summary(std::ostream&);
static void print_tpls(std::ostream& /*out*/);
};
#endif