/home/runner/work/kynema-sgf/kynema-sgf/src/utilities/FieldPlaneAveraging.H Source File

Kynema-SGF API: /home/runner/work/kynema-sgf/kynema-sgf/src/utilities/FieldPlaneAveraging.H Source File
Kynema-SGF API v0.1.0
CFD solver for wind plant simulations
Loading...
Searching...
No Matches
FieldPlaneAveraging.H
Go to the documentation of this file.
1#ifndef FIELDPLANEAVERAGING_H
2#define FIELDPLANEAVERAGING_H
3
5#include "src/CFDSim.H"
6#include "src/core/Field.H"
7#include "src/core/SimTime.H"
8
18
19namespace kynema_sgf {
20
28template <typename FType>
30{
31public:
42 const FType& field_in,
43 const kynema_sgf::SimTime& time,
44 int axis_in,
45 int max_level,
46 bool compute_deriv = false);
47
48 virtual ~FPlaneAveraging() = default;
49
50 virtual void operator()();
51
52 void convert_x_to_ind(amrex::Real x, int& ind, amrex::Real& c) const;
53
55 [[nodiscard]] amrex::Real
56 line_average_interpolated(amrex::Real x, int comp) const;
58 [[nodiscard]] amrex::Real line_average_cell(int ind, int comp) const;
59
62 [[nodiscard]] amrex::Real
63 line_derivative_interpolated(amrex::Real x, int comp) const;
66 [[nodiscard]] amrex::Real
67 line_derivative_of_average_cell(int ind, int comp) const;
68
70 const std::string& filename, int step, amrex::Real time);
71 void output_line_average_ascii(int step, amrex::Real time);
72
74 void set_precision(int p) { m_precision = p; };
75
76 [[nodiscard]] amrex::Real dx() const { return m_dx; };
77 [[nodiscard]] amrex::Real xlo() const { return m_xlo; };
78
79 [[nodiscard]] int axis() const { return m_axis; };
80 [[nodiscard]] int level() const { return m_max_level; };
81 [[nodiscard]] int ncomp() const { return m_ncomp; };
82 [[nodiscard]] int ncell_plane() const { return m_ncell_plane; };
83 [[nodiscard]] int ncell_line() const { return m_ncell_line; };
84 [[nodiscard]] int last_updated_index() const
85 {
87 };
88
89 [[nodiscard]] const amrex::Vector<amrex::Real>& line_average() const
90 {
91 return m_line_average;
92 };
93 [[nodiscard]] const amrex::Vector<amrex::Real>& line_deriv() const
94 {
95 return m_line_deriv;
96 };
97 void line_average(int comp, amrex::Vector<amrex::Real>& l_vec);
98 [[nodiscard]] const amrex::Vector<amrex::Real>& line_centroids() const
99 {
100 return m_line_xcentroid;
101 };
102
103 [[nodiscard]] const FType& field() const { return m_field; };
104
105protected:
107
108 amrex::Vector<amrex::Real>
112 amrex::Vector<amrex::Real> m_line_deriv;
113
114 amrex::Vector<amrex::Real> m_line_xcentroid;
116
117 amrex::Real m_dx;
118 amrex::Real m_xlo;
119 amrex::Real m_xhi;
120
123
124 int m_precision = 4;
127
128 const FType& m_field;
130 const int m_axis;
132 const bool m_comp_deriv;
133
134public: // public for GPU
136 template <typename IndexSelector>
137 void compute_averages(const IndexSelector& idxOp);
138
141};
142
145
150{
151public:
152 VelPlaneAveraging(CFDSim& sim, int axis_in, int max_level);
153
154 ~VelPlaneAveraging() override = default;
155
156 void operator()() override;
157
158private:
159 amrex::Vector<amrex::Real>
162
163 amrex::Vector<amrex::Real>
166
167 amrex::Vector<amrex::Real>
170
171public: // public for GPU
173 template <typename IndexSelector>
174 void compute_hvelmag_averages(const IndexSelector& idxOp);
175
177 const amrex::Vector<amrex::Real>& line_hvelmag_average()
178 {
180 };
181
184 [[nodiscard]] amrex::Real
185 line_hvelmag_average_interpolated(amrex::Real x) const;
188 [[nodiscard]] amrex::Real line_su_average_interpolated(amrex::Real x) const;
191 [[nodiscard]] amrex::Real line_sv_average_interpolated(amrex::Real x) const;
192};
193
194} // namespace kynema_sgf
195
196#endif /* FieldPlaneAveraging_H */
Definition CFDSim.H:55
Definition FieldPlaneAveraging.H:30
void compute_line_derivatives()
Definition FieldPlaneAveraging.cpp:426
amrex::Real m_xlo
Definition FieldPlaneAveraging.H:118
amrex::Real line_average_interpolated(amrex::Real x, int comp) const
Definition FieldPlaneAveraging.cpp:151
const amrex::Vector< amrex::Real > & line_average() const
Definition FieldPlaneAveraging.H:89
int m_precision
Definition FieldPlaneAveraging.H:124
int m_ncomp
Definition FieldPlaneAveraging.H:106
amrex::Real line_derivative_interpolated(amrex::Real x, int comp) const
Definition FieldPlaneAveraging.cpp:467
int ncomp() const
Definition FieldPlaneAveraging.H:81
amrex::Vector< amrex::Real > m_line_average
Definition FieldPlaneAveraging.H:109
void line_average(int comp, amrex::Vector< amrex::Real > &l_vec)
Definition FieldPlaneAveraging.cpp:167
void output_line_average_ascii(const std::string &filename, int step, amrex::Real time)
Definition FieldPlaneAveraging.cpp:99
FPlaneAveraging(const FType &field_in, const kynema_sgf::SimTime &time, int axis_in, int max_level, bool compute_deriv=false)
Definition FieldPlaneAveraging.cpp:13
int ncell_line() const
Definition FieldPlaneAveraging.H:83
amrex::Real dx() const
Definition FieldPlaneAveraging.H:76
const int m_axis
Definition FieldPlaneAveraging.H:130
amrex::Real line_average_cell(int ind, int comp) const
Definition FieldPlaneAveraging.cpp:180
amrex::Real m_dx
Definition FieldPlaneAveraging.H:117
amrex::Real xlo() const
Definition FieldPlaneAveraging.H:77
amrex::Real line_derivative_of_average_cell(int ind, int comp) const
Definition FieldPlaneAveraging.cpp:439
amrex::Vector< amrex::Real > m_line_xcentroid
Definition FieldPlaneAveraging.H:114
amrex::Real m_xhi
Definition FieldPlaneAveraging.H:119
void convert_x_to_ind(amrex::Real x, int &ind, amrex::Real &c) const
Definition FieldPlaneAveraging.cpp:78
const SimTime & m_time
Definition FieldPlaneAveraging.H:129
int last_updated_index() const
Definition FieldPlaneAveraging.H:84
virtual ~FPlaneAveraging()=default
void set_precision(int p)
Definition FieldPlaneAveraging.H:74
const FType & field() const
Definition FieldPlaneAveraging.H:103
const Field & m_field
Definition FieldPlaneAveraging.H:128
int level() const
Definition FieldPlaneAveraging.H:80
void compute_averages(const IndexSelector &idxOp)
Definition FieldPlaneAveraging.cpp:221
int m_ncell_plane
Definition FieldPlaneAveraging.H:121
int m_max_level
Definition FieldPlaneAveraging.H:131
const amrex::Vector< amrex::Real > & line_centroids() const
Definition FieldPlaneAveraging.H:98
int ncell_plane() const
Definition FieldPlaneAveraging.H:82
const bool m_comp_deriv
Definition FieldPlaneAveraging.H:132
void output_line_average_ascii(int step, amrex::Real time)
Definition FieldPlaneAveraging.cpp:142
virtual void operator()()
Definition FieldPlaneAveraging.cpp:191
int m_ncell_line
Definition FieldPlaneAveraging.H:122
const amrex::Vector< amrex::Real > & line_deriv() const
Definition FieldPlaneAveraging.H:93
int axis() const
Definition FieldPlaneAveraging.H:79
int m_last_updated_index
Definition FieldPlaneAveraging.H:125
amrex::Vector< amrex::Real > m_line_deriv
Definition FieldPlaneAveraging.H:112
Definition SimTime.H:33
amrex::Vector< amrex::Real > m_line_Su_average
Definition FieldPlaneAveraging.H:164
amrex::Vector< amrex::Real > m_line_hvelmag_average
Definition FieldPlaneAveraging.H:160
amrex::Real line_hvelmag_average_interpolated(amrex::Real x) const
Definition FieldPlaneAveraging.cpp:702
amrex::Real line_su_average_interpolated(amrex::Real x) const
Definition FieldPlaneAveraging.cpp:714
amrex::Vector< amrex::Real > m_line_Sv_average
Definition FieldPlaneAveraging.H:168
void operator()() override
Definition FieldPlaneAveraging.cpp:499
const amrex::Vector< amrex::Real > & line_hvelmag_average()
Definition FieldPlaneAveraging.H:177
~VelPlaneAveraging() override=default
void compute_hvelmag_averages(const IndexSelector &idxOp)
Definition FieldPlaneAveraging.cpp:529
VelPlaneAveraging(CFDSim &sim, int axis_in, int max_level)
Definition FieldPlaneAveraging.cpp:485
amrex::Real line_sv_average_interpolated(amrex::Real x) const
Definition FieldPlaneAveraging.cpp:725
This test case is intended as an evaluation of the momentum advection scheme.
Definition BCInterface.cpp:10
FPlaneAveraging< Field > FieldPlaneAveraging
Definition FieldPlaneAveraging.H:143
FPlaneAveraging< ScratchField > ScratchFieldPlaneAveraging
Definition FieldPlaneAveraging.H:144