ReFRACtor
fp_gsl_matrix.h
Go to the documentation of this file.
1 #ifndef FP_GSL_MATRIX_H
2 #define FP_GSL_MATRIX_H
3 #include <gsl/gsl_matrix.h>
4 #include <blitz/array.h>
5 #include <boost/utility.hpp>
6 
7 namespace FullPhysics {
8 
9 /****************************************************************/
20 class GslMatrix : public boost::noncopyable {
21 public:
22 //-----------------------------------------------------------------------
24 //-----------------------------------------------------------------------
25 
27 
28 //-----------------------------------------------------------------------
33 //-----------------------------------------------------------------------
34 
35  GslMatrix(gsl_matrix* M, bool Owned = true) : gsl_matrix_(0)
36  {reset(M, Owned);}
37 
38 //-----------------------------------------------------------------------
47 //-----------------------------------------------------------------------
48 
49  GslMatrix(blitz::Array<double, 2>& M) : gsl_matrix_(0) {reset(M);}
50  ~GslMatrix();
51  void reset(gsl_matrix* M, bool Owned = true);
52  void reset(blitz::Array<double, 2>& M);
53 
54 //-----------------------------------------------------------------------
56 //-----------------------------------------------------------------------
57 
58  const gsl_matrix* gsl() const { return gsl_matrix_; }
59 
60 //-----------------------------------------------------------------------
62 //-----------------------------------------------------------------------
63 
64  gsl_matrix* gsl() { return gsl_matrix_; }
65 
66 //-----------------------------------------------------------------------
68 //-----------------------------------------------------------------------
69 
70  const blitz::Array<double, 2>& blitz_array() const {return blitz_array_;}
71 
72 //-----------------------------------------------------------------------
74 //-----------------------------------------------------------------------
75 
76  blitz::Array<double, 2>& blitz_array() {return blitz_array_;}
77 public:
78  blitz::Array<double, 2> blitz_array_;
79  gsl_matrix* gsl_matrix_;
80  bool owned_; // This only applies to gsl_matrix_,
81  // blitz_array manages itself.
83 };
84 
85 /****************************************************************/
96 class GslVector : public boost::noncopyable {
97 public:
98 //-----------------------------------------------------------------------
100 //-----------------------------------------------------------------------
101 
102  GslVector() : gsl_vector_(0) {}
103 
104 //-----------------------------------------------------------------------
109 //-----------------------------------------------------------------------
110 
111  GslVector(gsl_vector* M, bool Owned = true) :gsl_vector_(0)
112  { reset(M, Owned);}
113 
114 //-----------------------------------------------------------------------
123 //-----------------------------------------------------------------------
124 
125  GslVector(blitz::Array<double, 1>& M)
126  :gsl_vector_(0) { reset(M);}
127  ~GslVector();
128  void reset(gsl_vector* M, bool Owned = true);
129  void reset(blitz::Array<double, 1>& M);
130 
131 //-----------------------------------------------------------------------
133 //-----------------------------------------------------------------------
134 
135  const gsl_vector* gsl() const { return gsl_vector_; }
136 
137 //-----------------------------------------------------------------------
139 //-----------------------------------------------------------------------
140 
141  gsl_vector* gsl() { return gsl_vector_; }
142 
143 //-----------------------------------------------------------------------
145 //-----------------------------------------------------------------------
146 
147  const blitz::Array<double, 1>& blitz_array() const {return blitz_array_;}
148 
149 //-----------------------------------------------------------------------
151 //-----------------------------------------------------------------------
152 
153  blitz::Array<double, 1>& blitz_array() {return blitz_array_;}
154 public:
155  blitz::Array<double, 1> blitz_array_;
156  gsl_vector* gsl_vector_;
157  bool owned_; // This only applies to gsl_vector_,
158  // blitz_array manages itself.
160 };
161 }
162 #endif
blitz::Array< double, 2 > & blitz_array()
Return blitz::Array look at data.
Definition: fp_gsl_matrix.h:76
GslMatrix()
Default constructor.
Definition: fp_gsl_matrix.h:26
GslVector(gsl_vector *M, bool Owned=true)
Use data owned by gsl_vector.
gsl_matrix * gsl_matrix_
Definition: fp_gsl_matrix.h:79
const blitz::Array< double, 2 > & blitz_array() const
Return blitz::Array look at data.
Definition: fp_gsl_matrix.h:70
GslMatrix(blitz::Array< double, 2 > &M)
Use data owned by blitz::Array.
Definition: fp_gsl_matrix.h:49
blitz::Array< double, 1 > blitz_array_
This provides thin wrapper around the GNU Scientific Library gsl_matrix.
Definition: fp_gsl_matrix.h:20
GslMatrix(gsl_matrix *M, bool Owned=true)
Use data owned by gsl_matrix.
Definition: fp_gsl_matrix.h:35
const blitz::Array< double, 1 > & blitz_array() const
Return blitz::Array look at data.
~GslMatrix()
Destructor.
blitz::Array< double, 2 > blitz_array_
Definition: fp_gsl_matrix.h:78
gsl_vector * gsl()
Return gsl_vector look at data.
Contains classes to abstract away details in various Spurr Radiative Transfer software.
Definition: doxygen_python.h:1
const gsl_matrix * gsl() const
Return gsl_matrix look at data.
Definition: fp_gsl_matrix.h:58
GslVector()
Default constructor.
void reset(gsl_matrix *M, bool Owned=true)
Reset class to point to new matrix.
blitz::Array< double, 1 > & blitz_array()
Return blitz::Array look at data.
gsl_vector * gsl_vector_
gsl_matrix * gsl()
Return gsl_matrix look at data.
Definition: fp_gsl_matrix.h:64
const gsl_vector * gsl() const
Return gsl_vector look at data.
This provides thin wrapper around the GNU Scientific Library gsl_vector.
Definition: fp_gsl_matrix.h:96
GslVector(blitz::Array< double, 1 > &M)
Use data owned by blitz::Array.

Copyright © 2017, California Institute of Technology.
ALL RIGHTS RESERVED.
U.S. Government Sponsorship acknowledged.
Generated Fri Aug 24 2018 15:44:11