/home/runner/work/amr-wind/amr-wind/amr-wind/immersed_boundary/bluff_body/box_ops.H Source File

AMR-Wind API: /home/runner/work/amr-wind/amr-wind/amr-wind/immersed_boundary/bluff_body/box_ops.H Source File
AMR-Wind API v0.1.0
CFD solver for wind plant simulations
Loading...
Searching...
No Matches
box_ops.H
Go to the documentation of this file.
1#ifndef BOX_OPS_H
2#define BOX_OPS_H
3
8#include "AMReX_REAL.H"
9
10using namespace amrex::literals;
11
12namespace amr_wind::ib::ops {
13
14template <>
16{
17 void
18 operator()(Box::DataType& data, const ::amr_wind::utils::MultiParser& pp)
19 {
20 auto& wdata = data.meta();
21 auto& info = data.info();
22
23 bluff_body::read_inputs(wdata, info, pp);
24
25 pp.get("center", wdata.center_loc);
26 pp.get("length", wdata.length);
27 pp.get("width", wdata.width);
28 pp.get("height", wdata.height);
29
30 amrex::Real search_radius =
31 2.0_rt * amrex::max<amrex::Real>(
32 wdata.length,
33 amrex::max<amrex::Real>(wdata.width, wdata.height));
34
35 // clang-format off
36 const auto& origin=wdata.center_loc;
37 info.bound_box = amrex::RealBox(
38 origin[0] - search_radius,
39 origin[1] - search_radius,
40 origin[2] - search_radius,
41 origin[0] + search_radius,
42 origin[1] + search_radius,
43 origin[2] + search_radius);
44 // clang-format on
45 }
46};
47
48template <>
50{
52 {
53 const auto& wdata = data.meta();
54
55 auto& sim = data.sim();
56 // cppcheck-suppress constVariableReference
57 auto& mask_node = sim.repo().get_int_field("mask_node");
58 // cppcheck-suppress constVariableReference
59 auto& levelset = sim.repo().get_field("ib_levelset");
60
61 auto nlevels = sim.repo().num_active_levels();
62 auto geom = sim.mesh().Geom();
63
64 for (int lev = 0; lev < nlevels; ++lev) {
65 const auto& problo = geom[lev].ProbLoArray();
66 const auto& dx = geom[lev].CellSizeArray();
67 for (amrex::MFIter mfi(levelset(lev)); mfi.isValid(); ++mfi) {
68 const auto& bx = mfi.growntilebox();
69 const auto& epsilon_node = mask_node(lev).array(mfi);
70 const auto& phi = levelset(lev).array(mfi);
71
72 const amrex::Real x0 = wdata.center_loc[0];
73 const amrex::Real y0 = wdata.center_loc[1];
74 const amrex::Real z0 = wdata.center_loc[2];
75 const amrex::Real l = wdata.length;
76 const amrex::Real w = wdata.width;
77 const amrex::Real h = wdata.height;
78 amrex::ParallelFor(
79 bx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
80 const amrex::Real x = problo[0] + (i + 0.5_rt) * dx[0];
81 const amrex::Real y = problo[1] + (j + 0.5_rt) * dx[1];
82 const amrex::Real z = problo[2] + (k + 0.5_rt) * dx[2];
83
84 amrex::Real phi_loc = 0;
85 if ((std::abs(x - x0) <= 0.5_rt * l) &&
86 (std::abs(y - y0) <= 0.5_rt * w) &&
87 (std::abs(z - z0) <= 0.5_rt * h)) {
88 phi_loc = amrex::max<amrex::Real>(
89 std::abs(x - x0) - 0.5_rt * l,
90 amrex::max<amrex::Real>(
91 std::abs(y - y0) - 0.5_rt * w,
92 std::abs(z - z0) - 0.5_rt * h));
93 } else {
94 phi_loc = amrex::min<amrex::Real>(
95 std::abs(std::abs(x - x0) - 0.5_rt * l),
96 amrex::min<amrex::Real>(
97 std::abs(std::abs(y - y0) - 0.5_rt * w),
98 std::abs(std::abs(z - z0) - 0.5_rt * h)));
99 }
100 phi(i, j, k) =
101 amrex::min<amrex::Real>(phi_loc, phi(i, j, k));
102 });
103 const auto& nbx = mfi.nodaltilebox();
104 amrex::ParallelFor(
105 nbx, [=] AMREX_GPU_DEVICE(int i, int j, int k) noexcept {
106 const amrex::Real x = problo[0] + i * dx[0];
107 const amrex::Real y = problo[1] + j * dx[1];
108 const amrex::Real z = problo[2] + k * dx[2];
109
110 if ((std::abs(x - x0) <= 0.5_rt * l + 0.5_rt * dx[0]) &&
111 (std::abs(y - y0) <= 0.5_rt * w + 0.5_rt * dx[1]) &&
112 (std::abs(z - z0) <= 0.5_rt * h + 0.5_rt * dx[2])) {
113 epsilon_node(i, j, k) = 0;
114 }
115 });
116 }
117 }
118 }
119};
120
121} // namespace amr_wind::ib::ops
122
123#endif /* BOX_OPS_H */
FieldRepo & repo()
Return the field repository.
Definition CFDSim.H:75
FieldRepo & repo() const
FieldRepo instance that manages this field.
Definition Field.H:159
Field & get_field(const std::string &name, const FieldState fstate=FieldState::New) const
Definition FieldRepo.cpp:152
int num_active_levels() const noexcept
Total number of levels currently active in the AMR mesh.
Definition FieldRepo.H:361
IntField & get_int_field(const std::string &name, const FieldState fstate=FieldState::New) const
Return a reference to an integer field.
Definition FieldRepo.cpp:279
const FieldRepo & repo() const
Reference to the FieldRepo that holds the fabs.
Definition IntField.H:43
IBTrait::MetaType & meta()
Definition IBTypes.H:92
IBTrait::InfoType & info()
Definition IBTypes.H:89
CFDSim & sim()
Definition IBTypes.H:86
void read_inputs(BluffBodyBaseData &wdata, IBInfo &, const ::amr_wind::utils::MultiParser &pp)
Definition bluff_body_ops.cpp:19
Definition bluff_body_ops.H:49
Definition Box.H:17
IBDataHolder< Box > DataType
Definition Box.H:20
void operator()(Box::DataType &data)
Definition box_ops.H:51
Definition IBOps.H:32
void operator()(Box::DataType &data, const ::amr_wind::utils::MultiParser &pp)
Definition box_ops.H:18
Definition IBOps.H:19