Hex  1.0
Hydrogen-electron collision solver
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
bspline.h
Go to the documentation of this file.
1 /* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *\
2  * *
3  * / / / / __ \ \ / / *
4  * / /__ / / / _ \ \ \/ / *
5  * / ___ / | |/_/ / /\ \ *
6  * / / / / \_\ / / \ \ *
7  * *
8  * Jakub Benda (c) 2014 *
9  * Charles University in Prague *
10  * *
11 \* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
12 
13 #ifndef HEX_BSPLINE
14 #define HEX_BSPLINE
15 
16 #include "arrays.h"
17 #include "complex.h"
18 #include "misc.h"
19 
35 class Bspline
36 {
37  public:
38 
48  Bspline (int order, rArrayView const & rknots, double th, rArrayView const & cknots);
49 
70  Complex bspline(int i, int iknot, int k, Complex r) const;
71 
84  Complex dspline(int i, int iknot, int k, Complex r) const;
85 
92  void B(int i, int iknot, int n, const Complex* x, Complex* y) const;
93 
100  void dB(int i, int iknot, int n, const Complex* x, Complex* y) const;
101 
107  inline Complex rotate (double r) const
108  {
109  return (r <= R0_) ? r : R0_ + (r - R0_) * rotation_;
110  };
111 
117  inline double unrotate (Complex z) const
118  {
119  return (z.imag() == 0.) ? z.real() : R0_ + ((z - R0_) / rotation_).real();
120  };
121 
131  cArray zip (const cArrayView coeff, const rArrayView grid) const;
132 
144  cArray zip (const cArrayView coeff, const rArrayView xgrid, const rArrayView ygrid) const;
145 
149  void writeVTK (std::ofstream & out, const cArrayView coeff, const rArrayView xgrid, const rArrayView ygrid) const;
150 
157  int knot (Complex x) const;
158 
169  Complex eval(const cArrayView coeff, double x) const;
170 
181  Complex eval(const cArrayView coeff, double x, double y) const;
182 
183  // getters
184 
186  inline Complex const & t (int i) const { return *(t_ + i); }
187 
189  inline int Nspline() const { return Nspline_; }
190 
192  inline int Nknot() const { return Nknot_; }
193 
195  inline int Nreknot() const { return Nreknot_; }
196 
198  inline int order() const { return order_; }
199 
201  inline double R0() const { return R0_; };
202 
204  inline double Rmax() const { return Rmax_; };
205 
207  inline double ECStheta() const { return theta_; }
208 
210  inline rArray const & rknots() const { return rknots_; }
211 
213  inline rArray const & cknots() const { return cknots_; }
214 
215  private:
216 
218  rArray rknots_;
219 
221  rArray cknots_;
222 
224  Complex * restrict t_;
225 
227  double theta_;
228 
230  Complex rotation_;
231 
233  double R0_;
234 
236  double Rmax_;
237 
239  int Nknot_;
240 
242  int Nreknot_;
243 
245  int Nspline_;
246 
248  int Nintval_;
249 
251  int order_;
252 };
253 
254 #endif
Bspline(int order, rArrayView const &rknots, double th, rArrayView const &cknots)
Constructor.
Definition: bspline.cpp:318
Complex dspline(int i, int iknot, int k, Complex r) const
Evaluate derivative of a B-spline.
Definition: bspline.cpp:39
int Nknot() const
Number of knots.
Definition: bspline.h:192
double ECStheta() const
ECS rotation angle.
Definition: bspline.h:207
A comfortable number array class.
Definition: arrays.h:171
int order() const
B-spline order.
Definition: bspline.h:198
Complex const & t(int i) const
B-spline knot sequence.
Definition: bspline.h:186
rArray const & cknots() const
complex knots
Definition: bspline.h:213
B-spline environment.
Definition: bspline.h:35
#define restrict
Definition: misc.h:88
void B(int i, int iknot, int n, const Complex *x, Complex *y) const
B-spline.
Definition: bspline.cpp:59
Complex rotate(double r) const
Apply the ECS transformation.
Definition: bspline.h:107
int Nspline() const
Number of B-splines.
Definition: bspline.h:189
void dB(int i, int iknot, int n, const Complex *x, Complex *y) const
Derivative of a B-spline.
Definition: bspline.cpp:65
void writeVTK(std::ofstream &out, const cArrayView coeff, const rArrayView xgrid, const rArrayView ygrid) const
Zip 2D expansion to VTK.
Definition: bspline.cpp:268
double Rmax() const
End of complex grid (real, unrotated).
Definition: bspline.h:204
int Nreknot() const
Number of real knots.
Definition: bspline.h:195
Complex bspline(int i, int iknot, int k, Complex r) const
Evaluate B-spline.
Definition: bspline.cpp:25
cArray zip(const cArrayView coeff, const rArrayView grid) const
Zip 1D expansion.
Definition: bspline.cpp:77
double R0() const
End of real grid.
Definition: bspline.h:201
double unrotate(Complex z) const
Apply the inverse ECS transformation.
Definition: bspline.h:117
std::complex< double > Complex
Definition: complex.h:20
int knot(Complex x) const
Get knot index for coordinate.
Definition: bspline.cpp:342
Complex eval(const cArrayView coeff, double x) const
Evaluate 1D B-spline expansion.
Definition: bspline.cpp:362
rArray const & rknots() const
real knots
Definition: bspline.h:210