/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
FieldRefinement.H
Go to the documentation of this file.
1#ifndef FIELDREFINEMENT_H
2#define FIELDREFINEMENT_H
3
4#include "amr-wind/CFDSim.H"
6
7namespace amr_wind {
8class CFDSim;
9class Field;
10class IntField;
11
26class FieldRefinement : public RefinementCriteria::Register<FieldRefinement>
27{
28public:
29 static std::string identifier() { return "FieldRefinement"; }
30
31 explicit FieldRefinement(const CFDSim& sim);
32
33 ~FieldRefinement() override = default;
34
36 void initialize(const std::string& key) override;
37
38 void operator()(
39 const int level,
40 amrex::TagBoxArray& tags,
41 const amrex::Real time,
42 const int ngrow) override;
43
44 template <typename MF>
45 void tag(const int level, amrex::TagBoxArray& tags, const MF& mfab)
46 {
47 const bool tag_field = level <= m_max_lev_field;
48 const bool tag_grad = level <= m_max_lev_grad;
49 const auto& geom = m_sim.repo().mesh().Geom(level);
50 const auto& prob_lo = geom.ProbLoArray();
51 const auto& dx = geom.CellSizeArray();
52 const auto tagging_box = m_tagging_box;
53
54 const auto& tag_arrs = tags.arrays();
55 const auto& farrs = mfab.const_arrays();
56 if (tag_field) {
57 const auto fld_err = m_field_error[level];
58 amrex::ParallelFor(
59 mfab,
60 [=] AMREX_GPU_DEVICE(int nbx, int i, int j, int k) noexcept {
61 const amrex::RealVect coord = {AMREX_D_DECL(
62 prob_lo[0] + (i + 0.5) * dx[0],
63 prob_lo[1] + (j + 0.5) * dx[1],
64 prob_lo[2] + (k + 0.5) * 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,
77 [=] AMREX_GPU_DEVICE(int nbx, int i, int j, int k) noexcept {
78 const amrex::RealVect coord = {AMREX_D_DECL(
79 prob_lo[0] + (i + 0.5) * dx[0],
80 prob_lo[1] + (j + 0.5) * dx[1],
81 prob_lo[2] + (k + 0.5) * dx[2])};
82 const auto axp =
83 std::abs(farrs[nbx](i + 1, j, k) - farrs[nbx](i, j, k));
84 const auto ayp =
85 std::abs(farrs[nbx](i, j + 1, k) - farrs[nbx](i, j, k));
86 const auto azp =
87 std::abs(farrs[nbx](i, j, k + 1) - farrs[nbx](i, j, k));
88 const auto axm =
89 std::abs(farrs[nbx](i - 1, j, k) - farrs[nbx](i, j, k));
90 const auto aym =
91 std::abs(farrs[nbx](i, j - 1, k) - farrs[nbx](i, j, k));
92 const auto azm =
93 std::abs(farrs[nbx](i, j, k - 1) - farrs[nbx](i, j, k));
94 const auto ax = amrex::max(axp, axm);
95 const auto ay = amrex::max(ayp, aym);
96 const auto az = amrex::max(azp, azm);
97 if ((amrex::max(ax, ay, az) >= gerr) &&
98 (tagging_box.contains(coord))) {
99 tag_arrs[nbx](i, j, k) = amrex::TagBox::SET;
100 }
101 });
102 }
103 }
104
105private:
106 const CFDSim& m_sim;
107
108 Field* m_field{nullptr};
110
111 amrex::Vector<amrex::Real> m_field_error;
112 amrex::Vector<amrex::Real> m_grad_error;
113
116 amrex::RealBox m_tagging_box;
117};
118
119} // namespace amr_wind
120
121#endif /* FIELDREFINEMENT_H */
Definition CFDSim.H:54
FieldRepo & repo()
Return the field repository.
Definition CFDSim.H:75
Definition Field.H:116
Definition FieldRefinement.H:27
void initialize(const std::string &key) override
Read input file and initialize boxarray used to refine each level.
Definition FieldRefinement.cpp:18
FieldRefinement(const CFDSim &sim)
Definition FieldRefinement.cpp:9
~FieldRefinement() override=default
Field * m_field
Definition FieldRefinement.H:108
int m_max_lev_grad
Definition FieldRefinement.H:115
amrex::Vector< amrex::Real > m_field_error
Definition FieldRefinement.H:111
static std::string identifier()
Definition FieldRefinement.H:29
void tag(const int level, amrex::TagBoxArray &tags, const MF &mfab)
Definition FieldRefinement.H:45
const CFDSim & m_sim
Definition FieldRefinement.H:106
amrex::Vector< amrex::Real > m_grad_error
Definition FieldRefinement.H:112
amrex::RealBox m_tagging_box
Definition FieldRefinement.H:116
IntField * m_int_field
Definition FieldRefinement.H:109
void operator()(const int level, amrex::TagBoxArray &tags, const amrex::Real time, const int ngrow) override
Definition FieldRefinement.cpp:77
int m_max_lev_field
Definition FieldRefinement.H:114
const amrex::AmrCore & mesh() const
Return a reference to the underlying AMR mesh instance.
Definition FieldRepo.H:358
Definition IntField.H:20
Definition BCInterface.cpp:7