Chi-Tech
quadrature.cc
Go to the documentation of this file.
1#include "quadrature.h"
2
4
6
7namespace chi_math
8{
9
12
13RegisterWrapperFunction(/*namespace_in_lua=*/chi_math,
14 /*name_in_lua=*/Get1DQuadratureData,
15 /*syntax_function=*/GetSyntax_Get1DQuadratureData,
16 /*actual_function=*/Get1DQuadratureData);
17
19{
21
23 "Lua wrapper function for getting the data from a quadrature.");
24 params.SetDocGroup("LuaQuadrature");
25
26 params.AddRequiredParameter<size_t>("arg0", "Handle to the quadrature");
27
28 return params;
29}
30
32{
34
35 const size_t handle = params.GetParamValue<size_t>("arg0");
36
37 auto& quad = Chi::GetStackItem<chi_math::Quadrature>(
38 Chi::object_stack, handle, __FUNCTION__);
39
40 chi::ParameterBlock qpoints_block("qpoints");
41 chi::ParameterBlock weights_block("weights");
42 {
43 size_t k = 0;
44 for (const auto& qpointXYZ : quad.qpoints_)
45 {
46 qpoints_block.AddParameter(
47 std::to_string(k++),
48 std::vector<double>{qpointXYZ.x, qpointXYZ.y, qpointXYZ.z});
49 }
50 k = 0;
51 for (const double w : quad.weights_)
52 weights_block.AddParameter(std::to_string(k++), w);
53 }
54 qpoints_block.ChangeToArray();
55 weights_block.ChangeToArray();
56
57 output.AddParameter(qpoints_block);
58 output.AddParameter(weights_block);
59
60 chi::ParameterBlock output_as_table;
61 output_as_table.AddParameter(output);
62
63 return output_as_table;
64}
65
67{
69
70 params.SetGeneralDescription("\\defgroup chi_math__Quadrature\n"
71 "\\ingroup LuaQuadrature\n"
72 "Base class for 1D quadratures");
73
74 params.AddRequiredParameter<int>("order", "Quadrature order.");
75 params.AddOptionalParameter("verbose", false, "Enables verbose operations");
76
77 using namespace chi_data_types;
78 params.ConstrainParameterRange("order",
79 AllowableRangeLowHighLimit::New(0, 43));
80
81 return params;
82}
83
85 : ChiObject(params),
86 order_(static_cast<QuadratureOrder>(params.GetParamValue<int>("order"))),
87 verbose_(params.GetParamValue<bool>("verbose"))
88{
89}
90
91void Quadrature::SetRange(const std::pair<double, double>& in_range)
92{
93 const auto& old_range = range_;
94 const auto& new_range = in_range;
95
96 const double h_new = new_range.second - new_range.first;
97 const double h_old = old_range.second - old_range.first;
98
99 ChiInvalidArgumentIf(h_new <= 0.0 or h_old <= 0.0,
100 "Called with negative or zero ranges.");
101
103 "Called with no abscissae initialized.");
104
105 const double scale_factor = h_new / h_old;
106
107 for (unsigned int i = 0; i < qpoints_.size(); ++i)
108 {
109 qpoints_[i](0) =
110 new_range.first + (qpoints_[i][0] - old_range.first) * scale_factor;
111
112 weights_[i] *= scale_factor;
113 }
114
115 range_ = in_range;
116}
117
118} // namespace chi_math
#define ChiInvalidArgumentIf(condition, message)
static std::vector< ChiObjectPtr > object_stack
Definition: chi_runtime.h:96
static chi::InputParameters GetInputParameters()
Definition: ChiObject.cc:4
void SetDocGroup(const std::string &doc_group)
void AddRequiredParameter(const std::string &name, const std::string &doc_string)
void ConstrainParameterRange(const std::string &param_name, AllowableRangePtr allowable_range)
void AddOptionalParameter(const std::string &name, T value, const std::string &doc_string)
void SetGeneralDescription(const std::string &description)
void AddParameter(ParameterBlock block)
T GetParamValue(const std::string &param_name) const
static chi::InputParameters GetInputParameters()
Definition: quadrature.cc:66
void SetRange(const std::pair< double, double > &in_range)
Definition: quadrature.cc:91
std::vector< chi_math::QuadraturePointXYZ > qpoints_
Definition: quadrature.h:37
std::pair< double, double > range_
Definition: quadrature.h:45
Quadrature(const chi::InputParameters &params)
Definition: quadrature.cc:84
std::vector< double > weights_
Definition: quadrature.h:38
chi::ParameterBlock Get1DQuadratureData(const chi::InputParameters &params)
Definition: quadrature.cc:31
QuadratureOrder
Definition: quadrature.h:12
chi::InputParameters GetSyntax_Get1DQuadratureData()
Definition: quadrature.cc:18
RegisterWrapperFunction(chi_math, Get1DQuadratureData, GetSyntax_Get1DQuadratureData, Get1DQuadratureData)