File: ortho_poly.h

package info (click to toggle)
lorene 0.0.0~cvs20161116%2Bdfsg-1
  • links: PTS, VCS
  • area: main
  • in suites: bookworm, bullseye, buster, stretch
  • size: 26,444 kB
  • ctags: 13,953
  • sloc: cpp: 212,946; fortran: 21,645; makefile: 1,750; sh: 4
file content (435 lines) | stat: -rw-r--r-- 15,369 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
/*
 * Declaration of class Ortho_poly and derived classes
 */

/*
 *   Copyright (c) 2005 Eric Gourgoulhon
 *
 *   This file is part of LORENE.
 *
 *   LORENE is free software; you can redistribute it and/or modify
 *   it under the terms of the GNU General Public License as published by
 *   the Free Software Foundation; either version 2 of the License, or
 *   (at your option) any later version.
 *
 *   LORENE is distributed in the hope that it will be useful,
 *   but WITHOUT ANY WARRANTY; without even the implied warranty of
 *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
 *   GNU General Public License for more details.
 *
 *   You should have received a copy of the GNU General Public License
 *   along with LORENE; if not, write to the Free Software
 *   Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
 *
 */

/*
 * $Id: ortho_poly.h,v 1.1 2005/11/14 01:57:00 e_gourgoulhon Exp $
 * $Log: ortho_poly.h,v $
 * Revision 1.1  2005/11/14 01:57:00  e_gourgoulhon
 * First version
 *
 *
 * $Header: /cvsroot/Lorene/School05/Monday/ortho_poly.h,v 1.1 2005/11/14 01:57:00 e_gourgoulhon Exp $
 *
 */

#ifndef __ORTHO_POLY_H_ 
#define __ORTHO_POLY_H_ 

class Grid ; 

            //----------------------------------//
            //          class Ortho_poly        //
            //----------------------------------//

/** Base class for orthogonal polynomials on [-1,1]. \ingroup (poly)
 * Storage of polynomials of an orthogonal family \f$(p_i)\f$ up to a 
 * given degree \e N.
 * The class \c Ortho_poly is abstract: it cannot be instanciated. 
 *
 */
class Ortho_poly {

    // Data : 
    // -----
    protected:
        const int nn ; ///< \e N = maximum degree of the polynomials
        
        /** Pointer on the Gauss nodes
         */
        mutable const Grid* p_gauss_nodes ; 

        /** Pointer on the Gauss weights
         */
        mutable double* p_gauss_weights ; 

        /** Pointer on the gamma factors (squares of the polynomials with
         *  respect to the discrete scalar product associated with the
         *   Gauss nodes)
         */
        mutable double* p_gauss_gamma ; 

        /** Pointer on the Gauss-Lobatto nodes
         */
        mutable const Grid* p_gauss_lobatto_nodes ; 
        
        /** Pointer on the Gauss-Lobatto weights
         */
        mutable double* p_gauss_lobatto_weights ; 

        /** Pointer on the gamma factors (square of the polynomials with
         *  respect to the discrete scalar product associated with the
         *   Gauss-Lobatto nodes)
         */
        mutable double* p_gauss_lobatto_gamma ; 

        
    // Constructors - Destructor
    // -------------------------
    protected:
	    Ortho_poly(const Ortho_poly& ) ; ///< Copy constructor
        
        /** Constructor to be used only by derived classes (hence protected)
         *  @param ni maximum degree \e N of the polynomials 
         */
        Ortho_poly(int ni) ;
         
    public:
	    virtual ~Ortho_poly() ;			///< Destructor
 
    // Mutators / assignment
    // ---------------------
    public:
	    /// Assignment to another Ortho_poly
	    void operator=(const Ortho_poly&)  ;	
	
    // Accessors
    // ---------
    public:
        int n() const ; ///< returns \e N, i.e. the maximum degree of the polynomials
        
    // Computation
    // -----------
    public: 
        /** Weight function \f$w\f$.
         *
         * @param x point in [-1,1] where the weight function is to be evaluated
         * @return value of  \f$w(x)\f$
         */
        virtual double weight(double x) const = 0 ; 

        /** Value of the polynomial \f$p_i\f$
         *
         * @param i no. of the polynomial
         * @param x point in [-1,1] where the polynomial is to be evaluated
         * @return value of  \f$p_i(x)\f$
         *
         */
        virtual double operator()(int i, double x) const = 0 ; 
        
        /// Gauss nodes
        virtual const Grid& gauss_nodes() const = 0 ;
        
        /** Gauss weights
         * 
         * @param i index of the weight
         * @return weight \f$w_i\f$
         */
        virtual double gauss_weight(int i) const = 0 ;
        
        /** Gamma factor \f$\gamma_i\f$ (square of the polynomial 
         * number \e i with
         *  respect to the discrete scalar product associated with the
         *   Gauss nodes)
         */
        virtual double gauss_gamma(int i) const = 0 ;

        /// Gauss-Lobatto nodes
        virtual const Grid& gauss_lobatto_nodes() const = 0 ;
        
        /** Gauss-Lobatto weights
         * 
         * @param i index of the weight
         * @return weight \f$w_i\f$
         */
        virtual double gauss_lobatto_weight(int i) const = 0 ;
        
        /** Gamma factor \f$\gamma_i\f$ (square of the polynomial number
         *  \e i with
         *  respect to the discrete scalar product associated with the
         *   Gauss-Lobatto nodes)
         */
        virtual double gauss_lobatto_gamma(int i) const = 0 ;

        /** Coefficients of the interpolant polynomial through the 
         *  Gauss nodes.
         *
         *  @param (*f)(double) function to be interpolated         
         *  @param cf [output] array containing the \e N+1 coefficients
         *          of the expansion of the interpolant on the basis polynomials.
         *          This array must be of size at least \e N+1 and must have
         *          been allocated by the user prior to the call of this method.       
         */
        void coef_interpolant_Gauss(double (*f)(double), double* cf) const ;  

        /** Coefficients of the interpolant polynomial through the 
         *  Gauss-Lobatto nodes.
         *
         *  @param (*f)(double) function to be interpolated         
         *  @param cf [output] array containing the \e N+1 coefficients
         *          of the expansion of the interpolant on the basis polynomials.
         *          This array must be of size at least \e N+1 and must have
         *          been allocated by the user prior to the call of this method.       
         */
        void coef_interpolant_GL(double (*f)(double), double* cf) const ;  

        /** Coefficients of the expansion over the basis polynomials of the 
         * orthogonal projection on the space of polynomials of maximum 
         * degree \e N. 
         *
         *  @param (*f)(double) function to be projected      
         *  @param cf [output] array containing the \e N+1 coefficients of the 
         *   expansion over the basis polynomials of the 
         *      orthogonal projection on the space of polynomials of maximum 
         *      degree \e N.          
         *          This array must be of size at least \e N+1 and must have
         *          been allocated by the user prior to the call of this method.       
         */
        virtual void coef_projection(double (*f)(double), double* cf) const = 0 ;  

        /** Value of some expansion over the polynomials:
         * \f$\sum_{i=0}^N a_i\, p_i(x)\f$
         *
         *  @param a array of size at least \e N+1 and storing the \e N+1 
         *          coefficients \f$a_i\f$ of the expansion on the basis 
         *          polynomials \f$p_i(x)\f$.
         *          For an interpolant polynomial, \e a is the array \c cf
         *          computed by the method
         *          \c coef_interpolant_Gauss (Gauss nodes) or
         *          \c coef_interpolant_GL (Gauss-Lobatto nodes) 
         *  @param x point in [-1,1] where the series is to be evaluated
         *  @return value of \f$\sum_{i=0}^N a_i\, p_i(x)\f$
         */
        double series(const double* a, double x) const ;  

};
    
            //----------------------------------//
            //      class Chebyshev_poly        //
            //----------------------------------//

/** Chebyshev polynomials. \ingroup (poly)
 * Each instance of this class represents a family of Chebyshev polynomials 
 * \f$(T_i)\f$ up to a given degree \e N.
 *
 */
class Chebyshev_poly : public Ortho_poly {
        
    // Constructors - Destructor
    // -------------------------
    public:
        /** Standard constructor
         *  @param ni maximum degree of the polynomials 
         */
        Chebyshev_poly(int ni ) ; 
        
	    Chebyshev_poly(const Chebyshev_poly& ) ; ///< Copy constructor
        
    public:
	    virtual ~Chebyshev_poly() ;			///< Destructor
 
    // Mutators / assignment
    // ---------------------
    public:
	    /// Assignment to another Chebyshev_poly
	    void operator=(const Chebyshev_poly&) ;	
	
            
    // Computation
    // -----------
    public: 
        /** Weight function \f$w\f$.
         *
         * @param x point in [-1,1] where the weight function is to be evaluated
         * @return \f$w(x)=1/\sqrt{1-x^2}\f$
         */
        virtual double weight(double x) const ; 

        /** Value of the polynomial \f$p_i\f$
         *
         * @param i no. of the polynomial
         * @param x point in [-1,1] where the polynomial is to be evaluated
         * @return value of  \f$p_i(x)\f$
         *
         */
        virtual double operator()(int i, double x) const ; 

        /// Gauss nodes
        virtual const Grid& gauss_nodes() const ;
        
        /** Gauss weights
         * 
         * @param i index of the weight
         * @return weight \f$w_i\f$
         */
        virtual double gauss_weight(int i) const ;
        
        /** Gamma factor \f$\gamma_i\f$ (square of the polynomial no. \e i with
         *  respect to the discrete scalar product associated with the
         *   Gauss nodes)
         */
        virtual double gauss_gamma(int i) const ;

        /// Gauss-Lobatto nodes
        virtual const Grid& gauss_lobatto_nodes() const ;
        
        /** Gauss-Lobatto weights
         * 
         * @param i index of the weight
         * @return weight \f$w_i\f$
         */
        virtual double gauss_lobatto_weight(int i) const ;
        
        /** Gamma factor \f$\gamma_i\f$ (square of the polynomial no. \e i with
         *  respect to the discrete scalar product associated with the
         *   Gauss-Lobatto nodes)
         */
        virtual double gauss_lobatto_gamma(int i) const ;

        /** Coefficients of the interpolant polynomial through the 
         *  Gauss-Lobatto nodes via a FFT.
         *
         *  @param (*f)(double) function to be interpolated         
         *  @param cf [output] array of size containing the \e N+1 coefficients
         *          of the expansion of the interpolant on the Chebyshev 
         *          polynomials.
         *          This array must be of size at least \e N+1 and must have
         *          allocated by the user prior to the call of this method.       
         */
        void coef_interpolant_GL_FFT(double (*f)(double), double* cf) const ;  

        /** Coefficients of the expansion over the Chebyshev polynomials of the 
         * orthogonal projection on the space of polynomials of maximum 
         * degree \e N. 
         *
         *  @param (*f)(double) function to be projected      
         *  @param cf [output] array containing the \e N+1 coefficients of the 
         *      expansion over the Chebyshev polynomials of the 
         *      orthogonal projection on the space of polynomials of maximum 
         *      degree \e N.          
         *          This array must be of size at least \e N+1 and must have
         *          been allocated by the user prior to the call of this method.       
         */
        virtual void coef_projection(double (*f)(double), double* cf) const  ;  

};
    
/// Display of an object of class \c Chebyshev_poly
ostream& operator<<(ostream& , const Chebyshev_poly& ) ;	


            //----------------------------------//
            //      class Legendre_poly        //
            //----------------------------------//

/** Legendre polynomials. \ingroup (poly)
 * Each instance of this class represents a family of Legendre polynomials 
 * \f$(T_i)\f$ up to a given degree \e N.
 *
 */
class Legendre_poly : public Ortho_poly {
        
    // Constructors - Destructor
    // -------------------------
    public:
        /** Standard constructor
         *  @param ni maximum degree of the polynomials 
         */
        Legendre_poly(int ni ) ; 
        
	    Legendre_poly(const Legendre_poly& ) ; ///< Copy constructor
        
    public:
	    virtual ~Legendre_poly() ;			///< Destructor
 
    // Mutators / assignment
    // ---------------------
    public:
	    /// Assignment to another Legendre_poly
	    void operator=(const Legendre_poly&) ;	
	
            
    // Computation
    // -----------
    public: 
        /** Weight function \f$w\f$.
         *
         * @param x point in [-1,1] where the weight function is to be evaluated
         * @return \f$w(x)=1\f$ (= constant value 1 for Legendre polynomials)
         */
        virtual double weight(double x) const ; 

        /** Value of the polynomial \f$p_i\f$
         *
         * @param i no. of the polynomial
         * @param x point in [-1,1] where the polynomial is to be evaluated
         * @return value of  \f$p_i(x)\f$
         *
         */
        virtual double operator()(int i, double x) const ; 

        /// Gauss nodes
        virtual const Grid& gauss_nodes() const ;
        
        /** Gauss weights
         * 
         * @param i index of the weight
         * @return weight \f$w_i\f$
         */
        virtual double gauss_weight(int i) const ;
        
        /** Gamma factor \f$\gamma_i\f$ (square of the polynomial no. \e i with
         *  respect to the discrete scalar product associated with the
         *   Gauss nodes)
         */
        virtual double gauss_gamma(int i) const ;

        /// Gauss-Lobatto nodes
        virtual const Grid& gauss_lobatto_nodes() const ;
        
        /** Gauss-Lobatto weights
         * 
         * @param i index of the weight
         * @return weight \f$w_i\f$
         */
        virtual double gauss_lobatto_weight(int i) const ;
        
        /** Gamma factor \f$\gamma_i\f$ (square of the polynomial no. \e i with
         *  respect to the discrete scalar product associated with the
         *   Gauss-Lobatto nodes)
         */
        virtual double gauss_lobatto_gamma(int i) const ;

        /** Coefficients of the expansion over the Legendre polynomials of the 
         * orthogonal projection on the space of polynomials of maximum 
         * degree \e N. 
         *
         *  @param (*f)(double) function to be projected      
         *  @param cf [output] array containing the \e N+1 coefficients of the 
         *      expansion over the Legendre polynomials of the 
         *      orthogonal projection on the space of polynomials of maximum 
         *      degree \e N.          
         *          This array must be of size at least \e N+1 and must have
         *          been allocated by the user prior to the call of this method.       
         */
        virtual void coef_projection(double (*f)(double), double* cf) const  ;  

};
    
/// Display of an object of class \c Legendre_poly
ostream& operator<<(ostream& , const Legendre_poly& ) ;	


#endif