DiFfRG
Loading...
Searching...
No Matches
flow_equations.hh
Go to the documentation of this file.
1#pragma once
2
3// DiFfRG
7
8// external libraries
9#include <autodiff/forward/dual.hpp>
10#include <autodiff/forward/real.hpp>
11#include <deal.II/base/quadrature_lib.h>
12
13namespace DiFfRG
14{
15 using namespace dealii;
16
18 {
19 protected:
20 FlowEquations(const JSONValue &json, std::function<double(double)> optimize_x);
21
22 public:
23 void set_k(const double k);
24
29 void print_parameters(const std::string logname) const;
30
31 template <typename NT, int d, typename FUN> NT integrate(const FUN &fun) const
32 {
33 if constexpr (std::is_same_v<NT, autodiff::real> || std::is_same_v<NT, autodiff::dual>)
35 else
37 }
38 template <typename NT, int d, typename FUN> NT angle_integrate(const FUN &fun) const
39 {
40 if constexpr (std::is_same_v<NT, autodiff::real> || std::is_same_v<NT, autodiff::dual>)
42 else
44 }
45 template <typename NT, int d, typename FUN> NT two_angle_integrate(const FUN &fun) const
46 {
47 if constexpr (std::is_same_v<NT, autodiff::real> || std::is_same_v<NT, autodiff::dual>)
49 else
51 }
52 template <typename NT, int d, typename FUN> NT three_angle_integrate(const FUN &fun) const
53 {
54 if constexpr (std::is_same_v<NT, autodiff::real> || std::is_same_v<NT, autodiff::dual>)
56 else
58 }
59
60 template <typename NT, int d, typename FUN> NT integrate(const FUN &fun, const double k) const
61 {
62 if constexpr (std::is_same_v<NT, autodiff::real> || std::is_same_v<NT, autodiff::dual>)
64 else
66 }
67 template <typename NT, int d, typename FUN> NT angle_integrate(const FUN &fun, const double k) const
68 {
69 if constexpr (std::is_same_v<NT, autodiff::real> || std::is_same_v<NT, autodiff::dual>)
71 else
73 }
74 template <typename NT, int d, typename FUN> NT two_angle_integrate(const FUN &fun, const double k) const
75 {
76 if constexpr (std::is_same_v<NT, autodiff::real> || std::is_same_v<NT, autodiff::dual>)
78 else
80 }
81 template <typename NT, int d, typename FUN> NT three_angle_integrate(const FUN &fun, const double k) const
82 {
83 if constexpr (std::is_same_v<NT, autodiff::real> || std::is_same_v<NT, autodiff::dual>)
85 else
87 }
88
90
91 protected:
93
95 const QGauss<1> x_quadrature;
96 double x_extent;
97
99 const QGauss<1> angle_quadrature;
100
102
107
108 const double x_extent_tolerance;
109
110 const std::function<double(double)> optimize_x;
111
112 double k;
114 const int verbosity;
115 };
116
118 {
119 protected:
120 FlowEquationsFiniteT(const JSONValue &json, const double T, std::function<double(double)> optimize_x,
121 std::function<double(double)> optimize_x0, std::function<double(double)> optimize_q0);
122
123 public:
124 void set_k(const double k);
125
130 void print_parameters(const std::string logname) const;
131
132 template <typename NT, int d, typename FUN> NT integrate(const FUN &fun) const
133 {
134 if constexpr (std::is_same_v<NT, autodiff::real> || std::is_same_v<NT, autodiff::dual>)
136 else
138 }
139 template <typename NT, int d, typename FUN> NT angle_integrate(const FUN &fun) const
140 {
141 if constexpr (std::is_same_v<NT, autodiff::real> || std::is_same_v<NT, autodiff::dual>)
143 else
145 }
146 template <typename NT, int d, typename FUN> NT two_angle_integrate(const FUN &fun) const
147 {
148 if constexpr (std::is_same_v<NT, autodiff::real> || std::is_same_v<NT, autodiff::dual>)
150 else
152 }
153 template <typename NT, int d, typename FUN> NT three_angle_integrate(const FUN &fun) const
154 {
155 if constexpr (std::is_same_v<NT, autodiff::real> || std::is_same_v<NT, autodiff::dual>)
157 else
159 }
160 template <typename NT, int d, typename FUN> NT spatial_integrate_and_integrate(const FUN &fun) const
161 {
162 if constexpr (std::is_same_v<NT, autodiff::real> || std::is_same_v<NT, autodiff::dual>)
165 else
167 q0_extent);
168 }
169 template <typename NT, int d, typename FUN> NT spatial_angle_integrate_and_integrate(const FUN &fun) const
170 {
171 if constexpr (std::is_same_v<NT, autodiff::real> || std::is_same_v<NT, autodiff::dual>)
174 else
177 }
178 template <typename NT, int d, typename FUN> NT spatial_integrate_and_sum(const FUN &fun, const double T) const
179 {
180 if constexpr (std::is_same_v<NT, autodiff::real> || std::is_same_v<NT, autodiff::dual>)
183 else
186 }
187 template <typename NT, int d, typename FUN> NT spatial_angle_integrate_and_sum(const FUN &fun, const double T) const
188 {
189 if constexpr (std::is_same_v<NT, autodiff::real> || std::is_same_v<NT, autodiff::dual>)
192 else
195 }
196 template <typename NT, typename FUN> NT sum(const FUN &fun, const double T) const
197 {
198 if constexpr (std::is_same_v<NT, autodiff::real> || std::is_same_v<NT, autodiff::dual>)
200 else
202 }
203
204 template <typename NT, int d, typename FUN> NT integrate(const FUN &fun, const double k) const
205 {
206 if constexpr (std::is_same_v<NT, autodiff::real> || std::is_same_v<NT, autodiff::dual>)
208 else
210 }
211 template <typename NT, int d, typename FUN> NT angle_integrate(const FUN &fun, const double k) const
212 {
213 if constexpr (std::is_same_v<NT, autodiff::real> || std::is_same_v<NT, autodiff::dual>)
215 else
217 }
218 template <typename NT, int d, typename FUN> NT two_angle_integrate(const FUN &fun, const double k) const
219 {
220 if constexpr (std::is_same_v<NT, autodiff::real> || std::is_same_v<NT, autodiff::dual>)
222 else
224 }
225 template <typename NT, int d, typename FUN> NT three_angle_integrate(const FUN &fun, const double k) const
226 {
227 if constexpr (std::is_same_v<NT, autodiff::real> || std::is_same_v<NT, autodiff::dual>)
229 else
231 }
232 template <typename NT, int d, typename FUN> NT spatial_integrate_and_integrate(const FUN &fun, const double k) const
233 {
234 if constexpr (std::is_same_v<NT, autodiff::real> || std::is_same_v<NT, autodiff::dual>)
237 else
239 q0_extent);
240 }
241 template <typename NT, int d, typename FUN>
242 NT spatial_angle_integrate_and_integrate(const FUN &fun, const double k) const
243 {
244 if constexpr (std::is_same_v<NT, autodiff::real> || std::is_same_v<NT, autodiff::dual>)
247 else
250 }
251 template <typename NT, int d, typename FUN>
252 NT spatial_integrate_and_sum(const FUN &fun, const double T, const double k) const
253 {
254 if constexpr (std::is_same_v<NT, autodiff::real> || std::is_same_v<NT, autodiff::dual>)
257 else
260 }
261 template <typename NT, int d, typename FUN>
262 NT spatial_angle_integrate_and_sum(const FUN &fun, const double T, const double k) const
263 {
264 if constexpr (std::is_same_v<NT, autodiff::real> || std::is_same_v<NT, autodiff::dual>)
267 else
270 }
271 template <typename NT, typename FUN> NT sum(const FUN &fun, const double T, const double k) const
272 {
273 if constexpr (std::is_same_v<NT, autodiff::real> || std::is_same_v<NT, autodiff::dual>)
275 else
277 }
278
280
281 protected:
285
287 const QGauss<1> x_quadrature;
288 double x_extent;
289
292 const QGauss<1> x0_quadrature;
293 double x0_extent;
294
297 const QGauss<1> q0_quadrature;
298 double q0_extent;
299
301 const QGauss<1> angle_quadrature;
302
304
313
314 const double x_extent_tolerance;
317
318 const double T;
319
320 const std::function<double(double)> optimize_x;
321 const std::function<double(double)> optimize_x0;
322 const std::function<double(double)> optimize_q0;
323
324 double k;
326 const int verbosity;
327 };
328} // namespace DiFfRG
Definition flow_equations.hh:118
NT spatial_integrate_and_integrate(const FUN &fun) const
Definition flow_equations.hh:160
NT angle_integrate(const FUN &fun, const double k) const
Definition flow_equations.hh:211
void print_parameters(const std::string logname) const
Print all deduced extents to the given spdlog.
NT spatial_angle_integrate_and_sum(const FUN &fun, const double T, const double k) const
Definition flow_equations.hh:262
NT spatial_angle_integrate_and_sum(const FUN &fun, const double T) const
Definition flow_equations.hh:187
NT two_angle_integrate(const FUN &fun, const double k) const
Definition flow_equations.hh:218
uint jac_q0_quadrature_order
Definition flow_equations.hh:309
NT spatial_integrate_and_integrate(const FUN &fun, const double k) const
Definition flow_equations.hh:232
NT spatial_integrate_and_sum(const FUN &fun, const double T) const
Definition flow_equations.hh:178
NT spatial_angle_integrate_and_integrate(const FUN &fun) const
Definition flow_equations.hh:169
const uint q0_quadrature_order
Definition flow_equations.hh:296
const double x0_extent_tolerance
Definition flow_equations.hh:315
const std::function< double(double)> optimize_x0
Definition flow_equations.hh:321
double k
Definition flow_equations.hh:324
const QGauss< 1 > x0_quadrature
Definition flow_equations.hh:292
double x0_extent
Definition flow_equations.hh:293
const QGauss< 1 > angle_quadrature
Definition flow_equations.hh:301
NT three_angle_integrate(const FUN &fun, const double k) const
Definition flow_equations.hh:225
const uint x0_quadrature_order
Definition flow_equations.hh:291
NT sum(const FUN &fun, const double T) const
Definition flow_equations.hh:196
double x_extent
Definition flow_equations.hh:288
void set_k(const double k)
double q0_extent
Definition flow_equations.hh:298
NT integrate(const FUN &fun) const
Definition flow_equations.hh:132
NT spatial_angle_integrate_and_integrate(const FUN &fun, const double k) const
Definition flow_equations.hh:242
FlowEquationsFiniteT(const JSONValue &json, const double T, std::function< double(double)> optimize_x, std::function< double(double)> optimize_x0, std::function< double(double)> optimize_q0)
const uint x_quadrature_order
Definition flow_equations.hh:286
const int verbosity
Definition flow_equations.hh:326
QGauss< 1 > jac_x_quadrature
Definition flow_equations.hh:306
NT three_angle_integrate(const FUN &fun) const
Definition flow_equations.hh:153
QGauss< 1 > jac_x0_quadrature
Definition flow_equations.hh:308
uint jac_angle_quadrature_order
Definition flow_equations.hh:311
NT spatial_integrate_and_sum(const FUN &fun, const double T, const double k) const
Definition flow_equations.hh:252
NT sum(const FUN &fun, const double T, const double k) const
Definition flow_equations.hh:271
uint jac_x0_quadrature_order
Definition flow_equations.hh:307
const double q0_extent_tolerance
Definition flow_equations.hh:316
const uint x0_summands
Definition flow_equations.hh:290
QGauss< 1 > jac_q0_quadrature
Definition flow_equations.hh:310
uint jac_x_quadrature_order
Definition flow_equations.hh:305
const double T
Definition flow_equations.hh:318
NT two_angle_integrate(const FUN &fun) const
Definition flow_equations.hh:146
QGauss< 1 > jac_angle_quadrature
Definition flow_equations.hh:312
const QGauss< 1 > q0_quadrature
Definition flow_equations.hh:297
NT integrate(const FUN &fun, const double k) const
Definition flow_equations.hh:204
void set_jacobian_quadrature_factor(const double jacobian_quadrature_factor)
const uint q0_summands
Definition flow_equations.hh:295
bool unoptimized
Definition flow_equations.hh:325
double jacobian_quadrature_factor
Definition flow_equations.hh:303
const QGauss< 1 > x_quadrature
Definition flow_equations.hh:287
const std::function< double(double)> optimize_x
Definition flow_equations.hh:320
NT angle_integrate(const FUN &fun) const
Definition flow_equations.hh:139
const std::function< double(double)> optimize_q0
Definition flow_equations.hh:322
const double x_extent_tolerance
Definition flow_equations.hh:314
const uint angle_quadrature_order
Definition flow_equations.hh:300
Definition flow_equations.hh:18
NT integrate(const FUN &fun) const
Definition flow_equations.hh:31
uint jac_angle_quadrature_order
Definition flow_equations.hh:105
NT angle_integrate(const FUN &fun) const
Definition flow_equations.hh:38
NT three_angle_integrate(const FUN &fun, const double k) const
Definition flow_equations.hh:81
FlowEquations(const JSONValue &json, std::function< double(double)> optimize_x)
double k
Definition flow_equations.hh:112
const uint angle_quadrature_order
Definition flow_equations.hh:98
NT two_angle_integrate(const FUN &fun) const
Definition flow_equations.hh:45
NT three_angle_integrate(const FUN &fun) const
Definition flow_equations.hh:52
NT integrate(const FUN &fun, const double k) const
Definition flow_equations.hh:60
NT two_angle_integrate(const FUN &fun, const double k) const
Definition flow_equations.hh:74
void set_k(const double k)
QGauss< 1 > jac_x_quadrature
Definition flow_equations.hh:104
uint jac_x_quadrature_order
Definition flow_equations.hh:103
void set_jacobian_quadrature_factor(const double jacobian_quadrature_factor)
bool unoptimized
Definition flow_equations.hh:113
const QGauss< 1 > x_quadrature
Definition flow_equations.hh:95
const int verbosity
Definition flow_equations.hh:114
const std::function< double(double)> optimize_x
Definition flow_equations.hh:110
const double x_extent_tolerance
Definition flow_equations.hh:108
NT angle_integrate(const FUN &fun, const double k) const
Definition flow_equations.hh:67
const uint x_quadrature_order
Definition flow_equations.hh:94
double jacobian_quadrature_factor
Definition flow_equations.hh:101
double x_extent
Definition flow_equations.hh:96
void print_parameters(const std::string logname) const
Print all deduced extents to the given spdlog.
QGauss< 1 > jac_angle_quadrature
Definition flow_equations.hh:106
const QGauss< 1 > angle_quadrature
Definition flow_equations.hh:99
A wrapper around the boost json value class.
Definition json.hh:19
NT spatial_angle_integrate_and_sum(const FUN &fun, const QGauss< 1 > &x_quadrature, const double x_extent, const double k, const QGauss< 1 > &cos_quadrature, const int m_order, const QGauss< 1 > &m_quadrature, const double m_extent, const double T)
Performs the integral.
Definition loop_integrals.hh:529
NT two_angle_integrate(const FUN &fun, const QGauss< 1 > &x_quadrature, const double x_extent, const double k, const QGauss< 1 > &cos_quadrature)
Performs the integral.
Definition loop_integrals.hh:145
NT spatial_integrate_and_integrate(const FUN &fun, const QGauss< 1 > &x_quadrature, const double x_extent, const double k, const QGauss< 1 > &m_quadrature, const double m_extent)
Performs the integral.
Definition loop_integrals.hh:296
NT spatial_angle_integrate_and_integrate(const FUN &fun, const QGauss< 1 > &x_quadrature, const double x_extent, const double k, const QGauss< 1 > &cos_quadrature, const QGauss< 1 > &m_quadrature, const double m_extent)
Performs the integral.
Definition loop_integrals.hh:360
NT angle_integrate(const FUN &fun, const QGauss< 1 > &x_quadrature, const double x_extent, const double k, const QGauss< 1 > &cos_quadrature)
Performs the integral.
Definition loop_integrals.hh:84
NT three_angle_integrate(const FUN &fun, const QGauss< 1 > &x_quadrature, const double x_extent, const double k, const QGauss< 1 > &cos_quadrature)
Performs the integral.
Definition loop_integrals.hh:215
NT integrate(const FUN &fun, const QGauss< 1 > &x_quadrature, const double x_extent, const double k)
Performs the integral.
Definition loop_integrals.hh:35
NT sum(const FUN &fun, const int m_order, const QGauss< 1 > &m_quadrature, const double m_extent, const double T)
Definition loop_integrals.hh:619
NT spatial_integrate_and_sum(const FUN &fun, const QGauss< 1 > &x_quadrature, const double x_extent, const double k, const int m_order, const QGauss< 1 > &m_quadrature, const double m_extent, const double T)
Performs the integral.
Definition loop_integrals.hh:436
Definition complex_math.hh:14
unsigned int uint
Definition utils.hh:22