]> Creatis software - creaMaracasVisu.git/blob - lib/maracasVisuLib/src/kernel/include/vector.cxx
BUG macOs
[creaMaracasVisu.git] / lib / maracasVisuLib / src / kernel / include / vector.cxx
1 /*=========================================================================
2
3   Program:   wxMaracas
4   Module:    $RCSfile: vector.cxx,v $
5   Language:  C++
6   Date:      $Date: 2009/05/14 13:55:08 $
7   Version:   $Revision: 1.1 $
8
9   Copyright: (c) 2002, 2003
10   License:
11   
12      This software is distributed WITHOUT ANY WARRANTY; without even 
13      the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR 
14      PURPOSE.  See the above copyright notice for more information.
15
16 =========================================================================*/
17
18 #include "gslobj.hxx"
19 #include <string>
20 #include <gsl/gsl_blas.h>
21
22 // ---------------------------------------------------------------------
23 std::ostream& operator<<( std::ostream& os, const gslobj_vector& v )
24 {
25     for( int i = 0; i < v.size( ); i++ ) os << " " << v( i );
26     return( os );
27 }
28
29 // ---------------------------------------------------------------------
30 gslobj_vector::gslobj_vector( size_t s )
31     : _copy( true )
32 {
33     _gsl = gsl_vector_alloc( s );
34 }
35
36 // ---------------------------------------------------------------------
37 gslobj_vector::gslobj_vector( const gslobj_vector& v )
38     : _copy( true )
39 {
40     _gsl = gsl_vector_alloc( v._gsl->size );
41     gsl_vector_memcpy( _gsl, v._gsl );
42 }
43
44 // ---------------------------------------------------------------------
45 gslobj_vector::gslobj_vector( double* v, size_t s )
46     : _copy( false )
47 {
48     _view = gsl_vector_view_array( v, s );
49     _gsl = &( _view.vector );
50 }
51
52 // ---------------------------------------------------------------------
53 gslobj_vector::gslobj_vector( gsl_vector* v )
54     : _copy( false )
55 {
56     _gsl = v;
57 }
58
59 // ---------------------------------------------------------------------
60 /*==PS========================================
61 void gslobj_vector::resize( size_t s )
62 {
63     if( _copy ) {
64
65         if( s != size( ) )
66             gsl_vector_free( _gsl );
67         _gsl = gsl_vector_alloc( s );
68
69     } // fi
70 }
71 ==PS========================================*/
72 // ---------------------------------------------------------------------
73 gslobj_vector& gslobj_vector::operator=( const gslobj_vector& o )
74 {
75     gsl_vector_memcpy( _gsl, o._gsl );
76     return( *this );
77 }
78
79 // ---------------------------------------------------------------------
80 gslobj_vector& gslobj_vector::operator=( double o )
81 {
82     gsl_vector_set_all( _gsl, o );
83     return( *this );
84 }
85
86 // ---------------------------------------------------------------------
87 gslobj_vector& gslobj_vector::operator=( double* o )
88 {
89     memcpy( _gsl->data, o, _gsl->size * sizeof( double ) );
90     return( *this );
91 }
92
93 // ---------------------------------------------------------------------
94 gslobj_vector& gslobj_vector::operator=( gsl_vector* o )
95 {
96     gsl_vector_memcpy( _gsl, o );
97     return( *this );
98 }
99
100 // ---------------------------------------------------------------------
101 bool gslobj_vector::operator==( const gslobj_vector& o ) const
102 {
103     bool ret = true;
104
105     if( _gsl->size != o._gsl->size )
106         return( false );
107
108     for( int i = 0; i < _gsl->size && ret; i++ )
109         ret = ret && ( _gsl->data[ i ] != o._gsl->data[ i ] );
110     return( ret );
111 }
112
113 // ---------------------------------------------------------------------
114 gslobj_vector gslobj_vector::operator+( const gslobj_vector& o )
115 {
116     gslobj_vector r( *this );
117     gsl_vector_add( r._gsl, o._gsl );
118     return( r );
119 }
120
121 // ---------------------------------------------------------------------
122 gslobj_vector gslobj_vector::operator+( double o )
123 {
124     gslobj_vector r( *this );
125     gsl_vector_add_constant( r._gsl, o );
126     return( r );
127 }
128
129 // ---------------------------------------------------------------------
130 gslobj_vector gslobj_vector::operator+( double* o )
131 {
132     gslobj_vector r( *this );
133     gsl_vector_view vw = gsl_vector_view_array( o, _gsl->size );
134     gsl_vector_add( r._gsl, &( vw.vector ) );
135     return( r );
136 }
137
138 // ---------------------------------------------------------------------
139 gslobj_vector gslobj_vector::operator+( gsl_vector* o )
140 {
141     gslobj_vector r( *this );
142     gsl_vector_add( r._gsl, o );
143     return( r );
144 }
145
146 // ---------------------------------------------------------------------
147 gslobj_vector& gslobj_vector::operator+=( const gslobj_vector& o )
148 {
149     gsl_vector_add( _gsl, o._gsl );
150     return( *this );
151 }
152
153 // ---------------------------------------------------------------------
154 gslobj_vector& gslobj_vector::operator+=( double o )
155 {
156     gsl_vector_add_constant( _gsl, o );
157     return( *this );
158 }
159
160 // ---------------------------------------------------------------------
161 gslobj_vector& gslobj_vector::operator+=( double* o )
162 {
163     gsl_vector_view vw = gsl_vector_view_array( o, _gsl->size );
164     gsl_vector_add( _gsl, &( vw.vector ) );
165     return( *this );
166 }
167
168 // ---------------------------------------------------------------------
169 gslobj_vector& gslobj_vector::operator+=( gsl_vector* o )
170 {
171     gsl_vector_add( _gsl, o );
172     return( *this );
173 }
174
175 // ---------------------------------------------------------------------
176 gslobj_vector gslobj_vector::operator-( const gslobj_vector& o )
177 {
178     gslobj_vector r( *this );
179     gsl_vector_sub( r._gsl, o._gsl );
180     return( r );
181 }
182
183 // ---------------------------------------------------------------------
184 gslobj_vector gslobj_vector::operator-( double o )
185 {
186     gslobj_vector r( *this );
187     gsl_vector_add_constant( r._gsl, -o );
188     return( r );
189 }
190
191 // ---------------------------------------------------------------------
192 gslobj_vector gslobj_vector::operator-( double* o )
193 {
194     gslobj_vector r( *this );
195     gsl_vector_view vw = gsl_vector_view_array( o, _gsl->size );
196     gsl_vector_sub( r._gsl, &( vw.vector ) );
197     return( r );
198 }
199
200 // ---------------------------------------------------------------------
201 gslobj_vector gslobj_vector::operator-( gsl_vector* o )
202 {
203     gslobj_vector r( *this );
204     gsl_vector_sub( r._gsl, o );
205     return( r );
206 }
207
208 // ---------------------------------------------------------------------
209 gslobj_vector& gslobj_vector::operator-=( const gslobj_vector& o )
210 {
211     gsl_vector_sub( _gsl, o._gsl );
212     return( *this );
213 }
214
215 // ---------------------------------------------------------------------
216 gslobj_vector& gslobj_vector::operator-=( double o )
217 {
218     gsl_vector_add_constant( _gsl, -o );
219     return( *this );
220 }
221
222 // ---------------------------------------------------------------------
223 gslobj_vector& gslobj_vector::operator-=( double* o )
224 {
225     gsl_vector_view vw = gsl_vector_view_array( o, _gsl->size );
226     gsl_vector_sub( _gsl, &( vw.vector ) );
227     return( *this );
228 }
229
230 // ---------------------------------------------------------------------
231 gslobj_vector& gslobj_vector::operator-=( gsl_vector* o )
232 {
233     gsl_vector_sub( _gsl, o );
234     return( *this );
235 }
236
237 // ---------------------------------------------------------------------
238 gslobj_vector gslobj_vector::operator*( double o )
239 {
240     gslobj_vector r( *this );
241     gsl_vector_scale( r._gsl, o );
242     return( r );
243 }
244
245 // ---------------------------------------------------------------------
246 gslobj_vector& gslobj_vector::operator*=( double o )
247 {
248     gsl_vector_scale( _gsl, o );
249     return( *this );
250 }
251
252 // ---------------------------------------------------------------------
253 double gslobj_vector::dot( const gslobj_vector& o )
254 {
255     double r;
256     gsl_blas_ddot( _gsl, o._gsl, &r );
257     return( r );
258 }
259
260 // ---------------------------------------------------------------------
261 double gslobj_vector::dot( double* o )
262 {
263     double r;
264     gsl_vector_view vw = gsl_vector_view_array( o, _gsl->size );
265     gsl_blas_ddot( _gsl, &( vw.vector ), &r );
266     return( r );
267 }
268
269 // ---------------------------------------------------------------------
270 double gslobj_vector::dot( gsl_vector* o )
271 {
272     double r;
273     gsl_blas_ddot( _gsl, o, &r );
274     return( r );
275 }
276
277 // ---------------------------------------------------------------------
278 gslobj_vector gslobj_vector::cross( const gslobj_vector& o )
279 {
280     return( cross( o._gsl->data ) );
281 }
282
283 // ---------------------------------------------------------------------
284 gslobj_vector gslobj_vector::cross( double* o )
285 {
286     gslobj_vector r( *this );
287     double S, s;
288         
289     // 1st element
290     S = _gsl->data[ 1 ] * o[ 2 ];
291     s = _gsl->data[ 2 ] * o[ 1 ];
292     r._gsl->data[ 0 ] = S - s;
293                 
294     // 2nd element
295     S = _gsl->data[ 2 ] * o[ 0 ];
296     s = _gsl->data[ 0 ] * o[ 2 ];
297     r._gsl->data[ 1 ] = S - s;
298                 
299     // 3rd element
300     S = _gsl->data[ 0 ] * o[ 1 ];
301     s = _gsl->data[ 1 ] * o[ 0 ];
302     r._gsl->data[ 2 ] = S - s;
303     return( r );
304 }
305
306 // ---------------------------------------------------------------------
307 gslobj_vector gslobj_vector::cross( gsl_vector* o )
308 {
309     return( cross( o->data ) );
310 }
311
312 // ---------------------------------------------------------------------
313 int gslobj_vector::scross( const gslobj_vector& o )
314 {
315     *this = cross( o );
316     return( 1 );
317 }
318
319 // ---------------------------------------------------------------------
320 int gslobj_vector::scross( double* o )
321 {
322     *this = cross( o );
323     return( 1 );
324 }
325
326 // ---------------------------------------------------------------------
327 int gslobj_vector::scross( gsl_vector* o )
328 {
329     *this = cross( o );
330     return( 1 );
331 }
332
333 // ---------------------------------------------------------------------
334 gslobj_vector gslobj_vector::operator/( double o )
335 {
336     gslobj_vector r( *this );
337     gsl_vector_scale( r._gsl, 1.0 / o );
338     return( r );
339 }
340
341 // ---------------------------------------------------------------------
342 gslobj_vector& gslobj_vector::operator/=( double o )
343 {
344     gsl_vector_scale( _gsl, 1.0 / o );
345     return( *this );
346 }
347
348 // ---------------------------------------------------------------------
349 /*==PS========================================
350 double gslobj_vector::norm1( )
351 {
352     return( gsl_blas_dasum( _gsl ) );
353 }
354 ==PS========================================*/
355 // ---------------------------------------------------------------------
356 double gslobj_vector::norm2( )
357 {
358     return( gsl_blas_dnrm2( _gsl ) );
359 }
360
361 // ---------------------------------------------------------------------
362 gslobj_vector gslobj_vector::normalize( )
363 {
364     gslobj_vector r( *this );
365     gsl_vector_scale( r._gsl, 1.0 / gsl_blas_dnrm2( _gsl ) );
366     return( r );
367 }
368
369 // ---------------------------------------------------------------------
370 int gslobj_vector::snormalize( )
371 {
372     return( gsl_vector_scale( _gsl, 1.0 / gsl_blas_dnrm2( _gsl ) ) );
373 }
374
375 // eof - gslobj_vector.cxx