/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 "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#ifdef AMREX_USE_OMP
55#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
56#endif
57 for (amrex::MFIter mfi(mfab, amrex::TilingIfNotGPU()); mfi.isValid();
58 ++mfi) {
59 const auto& bx = mfi.tilebox();
60 const auto& tag = tags.array(mfi);
61 const auto& farr = mfab.const_array(mfi);
62
63 if (tag_field) {
64 const auto fld_err = m_field_error[level];
65 amrex::ParallelFor(
66 bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
67 const amrex::RealVect coord = {AMREX_D_DECL(
68 prob_lo[0] + (i + 0.5) * dx[0],
69 prob_lo[1] + (j + 0.5) * dx[1],
70 prob_lo[2] + (k + 0.5) * dx[2])};
71
72 if ((farr(i, j, k) > fld_err) &&
73 (tagging_box.contains(coord))) {
74 tag(i, j, k) = amrex::TagBox::SET;
75 }
76 });
77 }
78
79 if (tag_grad) {
80 const auto gerr = m_grad_error[level];
81 amrex::ParallelFor(
82 bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
83 const amrex::RealVect coord = {AMREX_D_DECL(
84 prob_lo[0] + (i + 0.5) * dx[0],
85 prob_lo[1] + (j + 0.5) * dx[1],
86 prob_lo[2] + (k + 0.5) * dx[2])};
87 const auto axp =
88 std::abs(farr(i + 1, j, k) - farr(i, j, k));
89 const auto ayp =
90 std::abs(farr(i, j + 1, k) - farr(i, j, k));
91 const auto azp =
92 std::abs(farr(i, j, k + 1) - farr(i, j, k));
93 const auto axm =
94 std::abs(farr(i - 1, j, k) - farr(i, j, k));
95 const auto aym =
96 std::abs(farr(i, j - 1, k) - farr(i, j, k));
97 const auto azm =
98 std::abs(farr(i, j, k - 1) - farr(i, j, k));
99 const auto ax = amrex::max(axp, axm);
100 const auto ay = amrex::max(ayp, aym);
101 const auto az = amrex::max(azp, azm);
102 if ((amrex::max(ax, ay, az) >= gerr) &&
103 (tagging_box.contains(coord))) {
104 tag(i, j, k) = amrex::TagBox::SET;
105 }
106 });
107 }
108 }
109 }
110
111private:
112 const CFDSim& m_sim;
113
114 Field* m_field{nullptr};
116
117 amrex::Vector<amrex::Real> m_field_error;
118 amrex::Vector<amrex::Real> m_grad_error;
119
122 amrex::RealBox m_tagging_box;
123};
124
125} // namespace amr_wind
126
127#endif /* FIELDREFINEMENT_H */
Definition CFDSim.H:47
FieldRepo & repo()
Return the field repository.
Definition CFDSim.H:62
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:114
int m_max_lev_grad
Definition FieldRefinement.H:121
amrex::Vector< amrex::Real > m_field_error
Definition FieldRefinement.H:117
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:112
amrex::Vector< amrex::Real > m_grad_error
Definition FieldRefinement.H:118
amrex::RealBox m_tagging_box
Definition FieldRefinement.H:122
IntField * m_int_field
Definition FieldRefinement.H:115
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:120
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