/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 operator()(
42 const int level,
43 amrex::TagBoxArray& tags,
44 const amrex::Real time,
45 const int ngrow) override;
46
47 template <typename MF>
48 void tag(const int level, amrex::TagBoxArray& tags, const MF& mfab)
49 {
50 const bool tag_field = level <= m_max_lev_field;
51 const bool tag_grad = level <= m_max_lev_grad;
52 const auto& geom = m_sim.repo().mesh().Geom(level);
53 const auto& prob_lo = geom.ProbLoArray();
54 const auto& dx = geom.CellSizeArray();
55 const auto tagging_box = m_tagging_box;
56
57 const auto& tag_arrs = tags.arrays();
58 const auto& farrs = mfab.const_arrays();
59 if (tag_field) {
60 const auto fld_err = m_field_error[level];
61 amrex::ParallelFor(
62 mfab,
63 [=] AMREX_GPU_DEVICE(int nbx, int i, int j, int k) noexcept {
64 const amrex::RealVect coord = {AMREX_D_DECL(
65 prob_lo[0] + (i + 0.5_rt) * dx[0],
66 prob_lo[1] + (j + 0.5_rt) * dx[1],
67 prob_lo[2] + (k + 0.5_rt) * dx[2])};
68
69 if ((farrs[nbx](i, j, k) > fld_err) &&
70 (tagging_box.contains(coord))) {
71 tag_arrs[nbx](i, j, k) = amrex::TagBox::SET;
72 }
73 });
74 }
75
76 if (tag_grad) {
77 const auto gerr = m_grad_error[level];
78 amrex::ParallelFor(
79 mfab,
80 [=] AMREX_GPU_DEVICE(int nbx, int i, int j, int k) noexcept {
81 const amrex::RealVect coord = {AMREX_D_DECL(
82 prob_lo[0] + (i + 0.5_rt) * dx[0],
83 prob_lo[1] + (j + 0.5_rt) * dx[1],
84 prob_lo[2] + (k + 0.5_rt) * dx[2])};
85 const auto axp =
86 std::abs(farrs[nbx](i + 1, j, k) - farrs[nbx](i, j, k));
87 const auto ayp =
88 std::abs(farrs[nbx](i, j + 1, k) - farrs[nbx](i, j, k));
89 const auto azp =
90 std::abs(farrs[nbx](i, j, k + 1) - farrs[nbx](i, j, k));
91 const auto axm =
92 std::abs(farrs[nbx](i - 1, j, k) - farrs[nbx](i, j, k));
93 const auto aym =
94 std::abs(farrs[nbx](i, j - 1, k) - farrs[nbx](i, j, k));
95 const auto azm =
96 std::abs(farrs[nbx](i, j, k - 1) - farrs[nbx](i, j, k));
97 const auto ax = amrex::max<amrex::Real>(axp, axm);
98 const auto ay = amrex::max<amrex::Real>(ayp, aym);
99 const auto az = amrex::max<amrex::Real>(azp, azm);
100 if ((amrex::max<amrex::Real>(ax, ay, az) >= gerr) &&
101 (tagging_box.contains(coord))) {
102 tag_arrs[nbx](i, j, k) = amrex::TagBox::SET;
103 }
104 });
105 }
106 }
107
108private:
109 const CFDSim& m_sim;
110
111 Field* m_field{nullptr};
113
114 amrex::Vector<amrex::Real> m_field_error;
115 amrex::Vector<amrex::Real> m_grad_error;
116
119 amrex::RealBox m_tagging_box;
120};
121
122} // namespace amr_wind
123
124#endif /* FIELDREFINEMENT_H */
Definition CFDSim.H:54
Definition Field.H:116
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:111
int m_max_lev_grad
Definition FieldRefinement.H:118
amrex::Vector< amrex::Real > m_field_error
Definition FieldRefinement.H:114
static std::string identifier()
Definition FieldRefinement.H:32
void tag(const int level, amrex::TagBoxArray &tags, const MF &mfab)
Definition FieldRefinement.H:48
const CFDSim & m_sim
Definition FieldRefinement.H:109
amrex::Vector< amrex::Real > m_grad_error
Definition FieldRefinement.H:115
amrex::RealBox m_tagging_box
Definition FieldRefinement.H:119
IntField * m_int_field
Definition FieldRefinement.H:112
void operator()(const int level, amrex::TagBoxArray &tags, const amrex::Real time, const int ngrow) override
Definition FieldRefinement.cpp:80
int m_max_lev_field
Definition FieldRefinement.H:117
Definition IntField.H:20
This test case is intended as an evaluation of the momentum advection scheme.
Definition BCInterface.cpp:10