REMORA
Regional Modeling of Oceans Refined Adaptively
Loading...
Searching...
No Matches
REMORA_Tagging.cpp
Go to the documentation of this file.
1#include <REMORA.H>
2#include <REMORA_Derive.H>
3
4using namespace amrex;
5
6/**
7 * Function to tag cells for refinement -- this overrides the pure virtual function in AmrCore
8 *
9 * @param[in] levc level of refinement (0 is coarsest level)
10 * @param[out] tags array of tagged cells
11 * @param[in] time current time
12 * @param[in] ngrow number of grow cells
13*/
14void
15REMORA::ErrorEst (int levc, TagBoxArray& tags, Real time, int /*ngrow*/)
16{
17 const int clearval = TagBox::CLEAR;
18 const int tagval = TagBox::SET;
19 for (int j=0; j < ref_tags.size(); ++j)
20 {
21 std::unique_ptr<MultiFab> mf;
22
23 // This allows dynamic refinement based on the value of the scalar
24 if (ref_tags[j].Field() == "scalar")
25 {
26 mf = std::make_unique<MultiFab>(grids[levc], dmap[levc], 1, 0);
27 MultiFab::Copy(*mf,*cons_new[levc],Scalar_comp,0,1,0);
28 } else if (ref_tags[j].Field() == "temp") {
29 mf = std::make_unique<MultiFab>(grids[levc], dmap[levc], 1, 0);
30 MultiFab::Copy(*mf,*cons_new[levc],Temp_comp,0,1,0);
31 } else if (ref_tags[j].Field() == "salt") {
32 mf = std::make_unique<MultiFab>(grids[levc], dmap[levc], 1, 0);
33 MultiFab::Copy(*mf,*cons_new[levc],Salt_comp,0,1,0);
34 } else if (ref_tags[j].Field() == "x_velocity") {
35 mf = std::make_unique<MultiFab>(grids[levc], dmap[levc], 1, 0);
36 MultiFab::Copy(*mf,*xvel_new[levc],0,0,1,0);
37 } else if (ref_tags[j].Field() == "y_velocity") {
38 mf = std::make_unique<MultiFab>(grids[levc], dmap[levc], 1, 0);
39 MultiFab::Copy(*mf,*yvel_new[levc],0,0,1,0);
40 } else if (ref_tags[j].Field() == "z_velocity") {
41 mf = std::make_unique<MultiFab>(grids[levc], dmap[levc], 1, 0);
42 MultiFab::Copy(*mf,*zvel_new[levc],0,0,1,0);
43 } else if (ref_tags[j].Field() == "vorticity") {
44 mf = std::make_unique<MultiFab>(grids[levc], dmap[levc], 1, 0);
45 MultiFab mf_cc_vel(grids[levc],dmap[levc],3,1);
46 average_face_to_cellcenter(mf_cc_vel,0,
47 Array<const MultiFab*,3>{xvel_new[levc],
48 yvel_new[levc],
49 zvel_new[levc]});
50 // Impose bc's at domain boundaries at all levels
51 FillBdyCCVels(levc, mf_cc_vel);
52
53#ifdef _OPENMP
54#pragma omp parallel if (Gpu::notInLaunchRegion())
55#endif
56 for (MFIter mfi(*mf, TilingIfNotGPU()); mfi.isValid(); ++mfi)
57 {
58 const Box& bx = mfi.tilebox();
59 auto& dfab = (*mf)[mfi];
60 auto& sfab = mf_cc_vel[mfi];
61 auto pm = vec_pm[levc]->const_array(mfi);
62 auto pn = vec_pn[levc]->const_array(mfi);
63 derived::remora_dervort(bx, dfab, 0, 1, sfab, pm, pn, Geom(levc), time, nullptr, levc);
64 } // mfi
65
66#ifdef REMORA_USE_PARTICLES
67 } else {
68 //
69 // This allows dynamic refinement based on the number of particles per cell
70 //
71 // Note that we must count all the particles in levels both at and above the current,
72 // since otherwise, e.g., if the particles are all at level 1, counting particles at
73 // level 0 will not trigger refinement when regridding so level 1 will disappear,
74 // then come back at the next regridding
75 //
76 const auto& particles_namelist( particleData.getNames() );
77 mf->setVal(0.0);
78 for (ParticlesNamesVector::size_type i = 0; i < particles_namelist.size(); i++)
79 {
80 std::string tmp_string(particles_namelist[i]+"_count");
81 IntVect rr = IntVect::TheUnitVector();
82 if (ref_tags[j].Field() == tmp_string) {
83 for (int lev = levc; lev <= finest_level; lev++)
84 {
85 MultiFab temp_dat(grids[lev], dmap[lev], 1, 0); temp_dat.setVal(0);
86 particleData[particles_namelist[i]]->IncrementWithTotal(temp_dat, lev);
87
88 MultiFab temp_dat_crse(grids[levc], dmap[levc], 1, 0); temp_dat_crse.setVal(0);
89
90 if (lev == levc) {
91 MultiFab::Copy(*mf, temp_dat, 0, 0, 1, 0);
92 } else {
93 for (int d = 0; d < AMREX_SPACEDIM; d++) {
94 rr[d] *= ref_ratio[levc][d];
95 }
96 average_down(temp_dat, temp_dat_crse, 0, 1, rr);
97 MultiFab::Add(*mf, temp_dat_crse, 0, 0, 1, 0);
98 }
99 }
100 }
101 }
102
103
104#endif
105 }
106
107 ref_tags[j](tags,mf.get(),clearval,tagval,time,levc,geom[levc]);
108 }
109}
110
111/**
112 * Function to define the refinement criteria based on user input
113*/
114void
116{
117 if (max_level > 0)
118 {
119 ParmParse pp(pp_prefix);
120 Vector<std::string> refinement_indicators;
121 pp.queryarr("refinement_indicators",refinement_indicators,0,pp.countval("refinement_indicators"));
122
123 for (int i=0; i<refinement_indicators.size(); ++i)
124 {
125 std::string ref_prefix = pp_prefix + "." + refinement_indicators[i];
126
127 ParmParse ppr(ref_prefix);
128 RealBox realbox;
129 int lev_for_box;
130 if (ppr.countval("in_box_lo")) {
131 std::vector<Real> box_lo(3), box_hi(3);
132 ppr.get("max_level",lev_for_box);
133 if (lev_for_box <= max_level)
134 {
135 ppr.getarr("in_box_lo",box_lo,0,2);
136 ppr.getarr("in_box_hi",box_hi,0,2);
137 box_lo[2] = geom[0].ProbLo(2);
138 box_hi[2] = geom[0].ProbHi(2);
139 realbox = RealBox(&(box_lo[0]),&(box_hi[0]));
140
141 amrex::Print() << "Reading " << realbox << " at level " << lev_for_box << std::endl;
142 num_boxes_at_level[lev_for_box] += 1;
143
144 const auto* dx = geom[lev_for_box].CellSize();
145 const Real* plo = geom[lev_for_box].ProbLo();
146 int ilo = static_cast<int>((box_lo[0] - plo[0])/dx[0]);
147 int jlo = static_cast<int>((box_lo[1] - plo[1])/dx[1]);
148 int klo = static_cast<int>((box_lo[2] - plo[2])/dx[2]);
149 int ihi = static_cast<int>((box_hi[0] - plo[0])/dx[0]-1);
150 int jhi = static_cast<int>((box_hi[1] - plo[1])/dx[1]-1);
151 int khi = static_cast<int>((box_hi[2] - plo[2])/dx[2]-1);
152 Box bx(IntVect(ilo,jlo,klo),IntVect(ihi,jhi,khi));
153 if ( (ilo%ref_ratio[lev_for_box-1][0] != 0) || ((ihi+1)%ref_ratio[lev_for_box-1][0] != 0) ||
154 (jlo%ref_ratio[lev_for_box-1][1] != 0) || ((jhi+1)%ref_ratio[lev_for_box-1][1] != 0) )
155 amrex::Error("Fine box is not legit with this ref_ratio");
156 boxes_at_level[lev_for_box].push_back(bx);
157 amrex::Print() << "Saving in 'boxes at level' as " << bx << std::endl;
158 } // lev
159 }
160
161 AMRErrorTagInfo info;
162
163 if (realbox.ok()) {
164 info.SetRealBox(realbox);
165 }
166 if (ppr.countval("start_time") > 0) {
167 Real ref_min_time; ppr.get("start_time",ref_min_time);
168 info.SetMinTime(ref_min_time);
169 }
170 if (ppr.countval("end_time") > 0) {
171 Real ref_max_time; ppr.get("end_time",ref_max_time);
172 info.SetMaxTime(ref_max_time);
173 }
174 if (ppr.countval("max_level") > 0) {
175 int ref_max_level; ppr.get("max_level",ref_max_level);
176 info.SetMaxLevel(ref_max_level);
177 }
178
179 if (ppr.countval("value_greater")) {
180 int num_val = ppr.countval("value_greater");
181 Vector<Real> value(num_val);
182 ppr.getarr("value_greater",value,0,num_val);
183 std::string field; ppr.get("field_name",field);
184 ref_tags.push_back(AMRErrorTag(value,AMRErrorTag::GREATER,field,info));
185 }
186 else if (ppr.countval("value_less")) {
187 int num_val = ppr.countval("value_less");
188 Vector<Real> value(num_val);
189 ppr.getarr("value_less",value,0,num_val);
190 std::string field; ppr.get("field_name",field);
191 ref_tags.push_back(AMRErrorTag(value,AMRErrorTag::LESS,field,info));
192 }
193 else if (ppr.countval("adjacent_difference_greater")) {
194 int num_val = ppr.countval("adjacent_difference_greater");
195 Vector<Real> value(num_val);
196 ppr.getarr("adjacent_difference_greater",value,0,num_val);
197 std::string field; ppr.get("field_name",field);
198 ref_tags.push_back(AMRErrorTag(value,AMRErrorTag::GRAD,field,info));
199 }
200 else if (realbox.ok())
201 {
202 ref_tags.push_back(AMRErrorTag(info));
203 } else {
204 Abort(std::string("Unrecognized refinement indicator for " + refinement_indicators[i]).c_str());
205 }
206 } // loop over criteria
207 } // if max_level > 0
208}
209
#define Temp_comp
#define Scalar_comp
#define Salt_comp
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_pm
horizontal scaling factor: 1 / dx (2D)
Definition REMORA.H:355
amrex::Vector< amrex::MultiFab * > cons_new
multilevel data container for current step's scalar data: temperature, salinity, passive scalar
Definition REMORA.H:217
amrex::Vector< amrex::MultiFab * > zvel_new
multilevel data container for current step's z velocities (largely unused; W stored separately)
Definition REMORA.H:223
amrex::Vector< amrex::Vector< amrex::Box > > boxes_at_level
the boxes specified at each level by tagging criteria
Definition REMORA.H:1091
static amrex::Vector< amrex::AMRErrorTag > ref_tags
Holds info for dynamically generated tagging criteria.
Definition REMORA.H:1301
amrex::Vector< amrex::MultiFab * > yvel_new
multilevel data container for current step's y velocities (v in ROMS)
Definition REMORA.H:221
amrex::Vector< int > num_boxes_at_level
how many boxes specified at each level by tagging criteria
Definition REMORA.H:1087
amrex::Vector< amrex::MultiFab * > xvel_new
multilevel data container for current step's x velocities (u in ROMS)
Definition REMORA.H:219
void refinement_criteria_setup()
Set refinement criteria.
virtual void ErrorEst(int lev, amrex::TagBoxArray &tags, amrex::Real time, int ngrow) override
Tag cells for refinement.
void FillBdyCCVels(int lev, amrex::MultiFab &mf_cc_vel)
Fill the physical boundary conditions for cell-centered velocity (diagnostic only)
std::string pp_prefix
default prefix for input file parameters
Definition REMORA.H:205
amrex::Vector< std::unique_ptr< amrex::MultiFab > > vec_pn
horizontal scaling factor: 1 / dy (2D)
Definition REMORA.H:357
void remora_dervort(const amrex::Box &bx, amrex::FArrayBox &derfab, int dcomp, int ncomp, const amrex::FArrayBox &datfab, const amrex::Array4< const amrex::Real > &pm, const amrex::Array4< const amrex::Real > &pn, const amrex::Geometry &, amrex::Real, const int *, const int)