/home/runner/work/amr-wind/amr-wind/amr-wind/incflo.H Source File

AMR-Wind API: /home/runner/work/amr-wind/amr-wind/amr-wind/incflo.H Source File
AMR-Wind API v0.1.0
CFD solver for wind plant simulations
Loading...
Searching...
No Matches
incflo.H
Go to the documentation of this file.
1#ifndef INCFLO_H
2#define INCFLO_H
3
4#include <AMReX_AmrCore.H>
5#include <AMReX_MultiFabUtil.H>
6#include <AMReX_ParmParse.H>
7#include <AMReX_iMultiFab.H>
8#include <hydro_NodalProjector.H>
9
11#include "amr-wind/CFDSim.H"
15
17class MultiBlockContainer;
18namespace amr_wind {
19namespace pde {
20class PDEBase;
21}
24} // namespace amr_wind
25
46
51class incflo : public amrex::AmrCore
52{
53public:
54 incflo();
55 ~incflo() override;
56
57 // Initialize multilevel AMR data
58 void InitData();
59
60 // Evolve solution to final time through repeated calls to Advance()
61 void Evolve();
62
63 // Tag cells for refinement
64 void
65 ErrorEst(int lev, amrex::TagBoxArray& tags, amrex::Real time, int ngrow)
66 override;
67
68 // Make a new level from scratch using provided BoxArray and
69 // DistributionMapping Only used during initialization
71 int lev,
72 amrex::Real time,
73 const amrex::BoxArray& new_grids,
74 const amrex::DistributionMapping& new_dmap) override;
75
76 // Make a new level using provided BoxArray and DistributionMapping,
77 // and fill with interpolated coarse level data
79 int lev,
80 amrex::Real time,
81 const amrex::BoxArray& ba,
82 const amrex::DistributionMapping& dm) override;
83
84 // Remake an existing level using provided BoxArray and DistributionMapping,
85 // and fill with existing fine and coarse data
86 void RemakeLevel(
87 int lev,
88 amrex::Real time,
89 const amrex::BoxArray& ba,
90 const amrex::DistributionMapping& dm) override;
91
92 // Delete level data
93 void ClearLevel(int lev) override;
94
95 void init_mesh();
98 bool regrid_and_update();
99 void pre_advance_stage1();
100 void pre_advance_stage2();
101 void prepare_time_step();
102 void do_advance(int fixed_point_iteration);
103 void advance(int fixed_point_iteration);
104 void prescribe_advance();
105 void post_advance_work();
106
108 const amr_wind::SimTime& time() const { return m_time; }
110 const amr_wind::FieldRepo& repo() const { return m_repo; }
111
112 amr_wind::pde::PDEBase& icns() { return m_sim.pde_manager().icns(); }
114 {
115 return m_sim.pde_manager().icns();
116 }
118 {
119 return m_sim.pde_manager().scalar_eqns();
120 }
122 {
123 return m_sim.pde_manager().scalar_eqns();
124 }
125
126 amr_wind::Field& velocity() const { return m_repo.get_field("velocity"); }
127 amr_wind::Field& density() const { return m_repo.get_field("density"); }
129 {
130 return m_repo.get_field("temperature");
131 }
132 amr_wind::Field& pressure() const { return m_repo.get_field("p"); }
133 amr_wind::Field& grad_p() const { return m_repo.get_field("gp"); }
134
135 void compute_dt();
137 void advance_time() { m_time.advance_time(); }
138
139 void ApplyPredictor(
140 bool incremental_projection = false, int fixed_point_iteration = 0);
141 void ApplyCorrector();
142 void ApplyPrescribeStep();
143
144 void ApplyProjection(
145 amrex::Vector<amrex::MultiFab const*> density,
146 amrex::Real time,
147 amrex::Real scaling_factor,
148 bool incremental);
149
152
153 void ReadCheckpointFile();
154
155 void SetMultiBlockPointer(MultiBlockContainer* mbc) { m_sim.set_mbc(mbc); }
157 {
158 m_sim.set_read_erf(std::move(fcn));
159 }
160
161private:
162 //
163 // member variables
164 //
165
170
171 std::unique_ptr<amr_wind::RefineCriteriaManager> m_mesh_refiner;
172
173 // Be verbose?
174 int m_verbose = 0;
175
176 // Initial projection / iterations
177 bool m_do_initial_proj = true;
179
180 bool m_use_godunov = true;
181
182 // Prescribe advection velocity
183 bool m_prescribe_vel = false;
184
185 // Perform a dry run (0 steps, output plotfile)
186 bool m_dry_run = false;
187
188 // Fixed point iterations every timestep
190
192 amrex::Long m_cell_count{-1};
193
195
198
199 //
200 // end of member variables
201 //
202
203 amrex::FabFactory<amrex::FArrayBox> const& Factory(int lev) const
204 {
205 return m_repo.factory(lev);
206 }
207
208 bool need_divtau() const
209 {
211 }
212
213 //
214 // setup
215 //
217
218 void CheckAndSetUpDryRun();
219 void ReadParameters();
220 void InitialProjection();
221 void InitialIterations();
222
224 //
225 // utilities
226 //
228 void PrintMaxValues(const std::string& header);
229 void PrintMaxVelLocations(const std::string& header);
230 void PrintMaxVel(int lev) const;
231 void PrintMaxGp(int lev) const;
232 void CheckForNans(int lev) const;
233};
234
235#endif
std::function< void( const amrex::Real time, amrex::Vector< amrex::Real > &, amr_wind::InletData &, const amrex::Vector< amr_wind::Field * > &, MultiBlockContainer *)> ReadERFFunction
Definition ABLReadERFFunction.H:10
Definition CFDSim.H:54
amrex::Vector< TypePtr > TypeVector
Definition CollMgr.H:25
Definition Field.H:112
Definition FieldRepo.H:86
Definition OversetOps.H:13
Definition RefinementCriteria.H:57
Definition RefinementCriteria.H:33
Definition SimTime.H:33
Definition PDEBase.H:51
void Evolve()
Definition incflo.cpp:283
const amr_wind::pde::PDEMgr::TypeVector & scalar_eqns() const
Definition incflo.H:121
const amr_wind::SimTime & time() const
Definition incflo.H:108
amr_wind::Field & density() const
Definition incflo.H:127
void ApplyPrescribeStep()
Definition incflo_advance.cpp:655
void pre_advance_stage1()
Definition incflo_advance.cpp:15
void ClearLevel(int lev) override
Definition incflo_regrid.cpp:41
void advance(int fixed_point_iteration)
Definition incflo_advance.cpp:57
incflo()
Definition incflo.cpp:15
void ApplyProjection(amrex::Vector< amrex::MultiFab const * > density, amrex::Real time, amrex::Real scaling_factor, bool incremental)
Definition incflo_apply_nodal_projection.cpp:157
amr_wind::Field & temperature() const
Definition incflo.H:128
std::unique_ptr< amr_wind::RefineCriteriaManager > m_mesh_refiner
Definition incflo.H:171
bool m_dry_run
Definition incflo.H:186
bool m_prescribe_vel
Definition incflo.H:183
amr_wind::Field & velocity() const
Definition incflo.H:126
void prepare_for_time_integration()
Definition incflo.cpp:135
void post_advance_work()
Definition incflo.cpp:251
void InitialIterations()
Definition init.cpp:104
void init_amr_wind_modules()
Definition incflo.cpp:98
void ReadParameters()
Definition init.cpp:48
void RemakeLevel(int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) override
Definition incflo_regrid.cpp:24
void do_advance(int fixed_point_iteration)
Definition incflo.cpp:353
amr_wind::CFDSim m_sim
Definition incflo.H:166
int m_initial_iterations
Definition incflo.H:178
void CheckAndSetUpDryRun()
Definition init.cpp:14
int m_fixed_point_iterations
Definition incflo.H:189
void ReadCheckpointFile()
Definition io.cpp:12
bool need_divtau() const
Definition incflo.H:208
void PrintMaxGp(int lev) const
Definition diagnostics.cpp:719
void PrintMaxVelLocations(const std::string &header)
Definition diagnostics.cpp:700
amr_wind::OversetOps m_ovst_ops
Definition incflo.H:169
amr_wind::Field & grad_p() const
Definition incflo.H:133
DiffusionType m_diff_type
Definition incflo.H:194
bool m_use_godunov
Definition incflo.H:180
amr_wind::pde::PDEBase & icns()
Definition incflo.H:112
amr_wind::SimTime & m_time
Definition incflo.H:167
void advance_time()
Definition incflo.H:137
void MakeNewLevelFromScratch(int lev, amrex::Real time, const amrex::BoxArray &new_grids, const amrex::DistributionMapping &new_dmap) override
Definition incflo.cpp:376
amrex::Long m_cell_count
number of cells on all levels including covered cells
Definition incflo.H:192
void SetMultiBlockPointer(MultiBlockContainer *mbc)
Definition incflo.H:155
void compute_prescribe_dt()
Definition incflo_compute_dt.cpp:223
amr_wind::CFDSim & sim()
Definition incflo.H:107
void prepare_time_step()
Definition incflo_advance.cpp:31
void prescribe_advance()
Definition incflo_advance.cpp:646
amr_wind::FieldRepo & m_repo
Definition incflo.H:168
~incflo() override
void ErrorEst(int lev, amrex::TagBoxArray &tags, amrex::Real time, int ngrow) override
Definition incflo.cpp:459
bool m_do_initial_proj
Definition incflo.H:177
void MakeNewLevelFromCoarse(int lev, amrex::Real time, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm) override
Definition incflo_regrid.cpp:6
void init_physics_and_pde()
Initialize Physics instances as well as PDEs (include turbulence models)
Definition incflo.cpp:406
void InitData()
Definition incflo.cpp:171
amr_wind::FieldRepo & repo()
Definition incflo.H:109
void set_read_erf(ReadERFFunction fcn)
Definition incflo.H:156
void InitialProjection()
Definition init.cpp:201
amrex::FabFactory< amrex::FArrayBox > const & Factory(int lev) const
Definition incflo.H:203
bool m_reconstruct_true_pressure
reconstruct true pressure
Definition incflo.H:197
bool regrid_and_update()
Definition incflo.cpp:184
void PrintMaxVel(int lev) const
Definition diagnostics.cpp:708
const amr_wind::pde::PDEBase & icns() const
Definition incflo.H:113
int m_verbose
Definition incflo.H:174
amr_wind::Field & pressure() const
Definition incflo.H:132
void pre_advance_stage2()
Definition incflo_advance.cpp:21
void CheckForNans(int lev) const
Definition diagnostics.cpp:728
void init_mesh()
Definition incflo.cpp:42
void PrintMaxValues(const std::string &header)
Definition diagnostics.cpp:657
amr_wind::pde::PDEMgr::TypeVector & scalar_eqns()
Definition incflo.H:117
void ApplyCorrector()
Definition incflo_advance.cpp:520
const amr_wind::FieldRepo & repo() const
Definition incflo.H:110
void compute_dt()
Definition incflo_compute_dt.cpp:33
void ApplyPredictor(bool incremental_projection=false, int fixed_point_iteration=0)
Definition incflo_advance.cpp:180
DiffusionType
Definition incflo_enums.H:6
@ Implicit
Definition incflo_enums.H:10
@ Crank_Nicolson
Definition incflo_enums.H:9
Definition AdvOp_Godunov.H:21
This test case is intended as an evaluation of the momentum advection scheme.
Definition BCInterface.cpp:10