/home/runner/work/amr-wind/amr-wind/amr-wind/utilities/tagging/FieldRefinement.H Source File

AMR-Wind API: /home/runner/work/amr-wind/amr-wind/amr-wind/utilities/tagging/FieldRefinement.H Source File
AMR-Wind API v0.1.0
CFD solver for wind plant simulations
Loading...
Searching...
No Matches
FieldRefinement.H
Go to the documentation of this file.
1#ifndef FIELDREFINEMENT_H
2#define FIELDREFINEMENT_H
3
4#include "AMReX_REAL.H"
5#include "amr-wind/CFDSim.H"
7
8using namespace amrex::literals;
9
10namespace amr_wind {
11class CFDSim;
12class Field;
13class IntField;
14
29class FieldRefinement : public RefinementCriteria::Register<FieldRefinement>
30{
31public:
32 static std::string identifier() { return "FieldRefinement"; }
33
34 explicit FieldRefinement(const CFDSim& sim);
35
36 ~FieldRefinement() override = default;
37
39 void initialize(const std::string& key) override;
40
41 void
42 operator()(int level, amrex::TagBoxArray& tags, amrex::Real time, int ngrow)
43 override;
44
45 template <typename MF>
46 void tag(const int level, amrex::TagBoxArray& tags, const MF& mfab)
47 {
48 const bool tag_field = level <= m_max_lev_field;
49 const bool tag_grad = level <= m_max_lev_grad;
50 const auto& geom = m_sim.repo().mesh().Geom(level);
51 const auto& prob_lo = geom.ProbLoArray();
52 const auto& dx = geom.CellSizeArray();
53 const auto tagging_box = m_tagging_box;
54
55 const auto& tag_arrs = tags.arrays();
56 const auto& farrs = mfab.const_arrays();
57 if (tag_field) {
58 const auto fld_err = m_field_error[level];
59 amrex::ParallelFor(
60 mfab, [=] AMREX_GPU_DEVICE(int nbx, int i, int j, int k) {
61 const amrex::RealVect coord = {AMREX_D_DECL(
62 prob_lo[0] + ((i + 0.5_rt) * dx[0]),
63 prob_lo[1] + ((j + 0.5_rt) * dx[1]),
64 prob_lo[2] + ((k + 0.5_rt) * dx[2]))};
65
66 if ((farrs[nbx](i, j, k) > fld_err) &&
67 (tagging_box.contains(coord))) {
68 tag_arrs[nbx](i, j, k) = amrex::TagBox::SET;
69 }
70 });
71 }
72
73 if (tag_grad) {
74 const auto gerr = m_grad_error[level];
75 amrex::ParallelFor(
76 mfab, [=] AMREX_GPU_DEVICE(int nbx, int i, int j, int k) {
77 const amrex::RealVect coord = {AMREX_D_DECL(
78 prob_lo[0] + ((i + 0.5_rt) * dx[0]),
79 prob_lo[1] + ((j + 0.5_rt) * dx[1]),
80 prob_lo[2] + ((k + 0.5_rt) * dx[2]))};
81 const auto axp =
82 std::abs(farrs[nbx](i + 1, j, k) - farrs[nbx](i, j, k));
83 const auto ayp =
84 std::abs(farrs[nbx](i, j + 1, k) - farrs[nbx](i, j, k));
85 const auto azp =
86 std::abs(farrs[nbx](i, j, k + 1) - farrs[nbx](i, j, k));
87 const auto axm =
88 std::abs(farrs[nbx](i - 1, j, k) - farrs[nbx](i, j, k));
89 const auto aym =
90 std::abs(farrs[nbx](i, j - 1, k) - farrs[nbx](i, j, k));
91 const auto azm =
92 std::abs(farrs[nbx](i, j, k - 1) - farrs[nbx](i, j, k));
93 const auto ax = amrex::max<amrex::Real>(axp, axm);
94 const auto ay = amrex::max<amrex::Real>(ayp, aym);
95 const auto az = amrex::max<amrex::Real>(azp, azm);
96 if ((amrex::max<amrex::Real>(ax, ay, az) >= gerr) &&
97 (tagging_box.contains(coord))) {
98 tag_arrs[nbx](i, j, k) = amrex::TagBox::SET;
99 }
100 });
101 }
102 }
103
104private:
105 const CFDSim& m_sim;
106
107 Field* m_field{nullptr};
109
110 amrex::Vector<amrex::Real> m_field_error;
111 amrex::Vector<amrex::Real> m_grad_error;
112
115 amrex::RealBox m_tagging_box;
116};
117
118} // namespace amr_wind
119
120#endif /* FIELDREFINEMENT_H */
Definition CFDSim.H:54
Definition Field.H:112
void initialize(const std::string &key) override
Read input file and initialize boxarray used to refine each level.
Definition FieldRefinement.cpp:21
FieldRefinement(const CFDSim &sim)
Definition FieldRefinement.cpp:12
~FieldRefinement() override=default
Field * m_field
Definition FieldRefinement.H:107
int m_max_lev_grad
Definition FieldRefinement.H:114
amrex::Vector< amrex::Real > m_field_error
Definition FieldRefinement.H:110
static std::string identifier()
Definition FieldRefinement.H:32
void tag(const int level, amrex::TagBoxArray &tags, const MF &mfab)
Definition FieldRefinement.H:46
void operator()(int level, amrex::TagBoxArray &tags, amrex::Real time, int ngrow) override
Definition FieldRefinement.cpp:80
const CFDSim & m_sim
Definition FieldRefinement.H:105
amrex::Vector< amrex::Real > m_grad_error
Definition FieldRefinement.H:111
amrex::RealBox m_tagging_box
Definition FieldRefinement.H:115
IntField * m_int_field
Definition FieldRefinement.H:108
int m_max_lev_field
Definition FieldRefinement.H:113
Definition IntField.H:20
This test case is intended as an evaluation of the momentum advection scheme.
Definition BCInterface.cpp:10