LORENE
base_val.h
1 /*
2  * Definition of Lorene class Base_val
3  *
4  */
5 
6 /*
7  * Copyright (c) 1999-2000 Jean-Alain Marck
8  * Copyright (c) 1999-2001 Eric Gourgoulhon
9  * Copyright (c) 1999-2001 Philippe Grandclement
10  * Copyright (c) 2001 Jerome Novak
11  *
12  * This file is part of LORENE.
13  *
14  * LORENE is free software; you can redistribute it and/or modify
15  * it under the terms of the GNU General Public License as published by
16  * the Free Software Foundation; either version 2 of the License, or
17  * (at your option) any later version.
18  *
19  * LORENE is distributed in the hope that it will be useful,
20  * but WITHOUT ANY WARRANTY; without even the implied warranty of
21  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
22  * GNU General Public License for more details.
23  *
24  * You should have received a copy of the GNU General Public License
25  * along with LORENE; if not, write to the Free Software
26  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
27  *
28  */
29 
30 
31 #ifndef __BASE_VAL_H_
32 #define __BASE_VAL_H_
33 
34 /*
35  * $Id: base_val.h,v 1.24 2023/06/28 10:04:32 j_novak Exp $
36  * $Log: base_val.h,v $
37  * Revision 1.24 2023/06/28 10:04:32 j_novak
38  * Use of C++ strings and flows instead of C types.
39  *
40  * Revision 1.23 2014/10/13 08:52:31 j_novak
41  * Lorene classes and functions now belong to the namespace Lorene.
42  *
43  * Revision 1.22 2014/10/06 15:09:39 j_novak
44  * Modified #include directives to use c++ syntax.
45  *
46  * Revision 1.21 2013/06/05 15:10:41 j_novak
47  * Suppression of FINJAC sampling in r. This Jacobi(0,2) base is now
48  * available by setting colloc_r to BASE_JAC02 in the Mg3d constructor.
49  *
50  * Revision 1.20 2013/01/11 08:20:10 j_novak
51  * New radial spectral bases with Legendre polynomials (R_LEG, R_LEGP, R_LEGI).
52  *
53  * Revision 1.19 2012/01/17 10:18:17 j_penner
54  * changed nzone from a private member to a protected member
55  *
56  * Revision 1.18 2009/10/23 12:55:46 j_novak
57  * New base T_LEG_MI
58  *
59  * Revision 1.17 2009/10/08 16:19:32 j_novak
60  * Addition of new bases T_COS and T_SIN.
61  *
62  * Revision 1.16 2008/12/03 15:21:20 j_novak
63  * New method mult_cost.
64  *
65  * Revision 1.15 2007/12/11 15:28:05 jl_cornou
66  * Jacobi(0,2) polynomials partially implemented
67  *
68  * Revision 1.14 2005/09/07 13:09:48 j_novak
69  * New method for determining the highest multipole that can be described on a 3D
70  * grid.
71  *
72  * Revision 1.13 2004/11/23 15:03:14 m_forot
73  * Corrected error in comments.
74  *
75  * Revision 1.12 2004/11/04 15:39:45 e_gourgoulhon
76  * Modif documentation: added description of bases without any symmetry
77  * in theta.
78  *
79  * Revision 1.11 2004/08/24 09:14:40 p_grandclement
80  * Addition of some new operators, like Poisson in 2d... It now requieres the
81  * GSL library to work.
82  *
83  * Also, the way a variable change is stored by a Param_elliptic is changed and
84  * no longer uses Change_var but rather 2 Scalars. The codes using that feature
85  * will requiere some modification. (It should concern only the ones about monopoles)
86  *
87  * Revision 1.10 2004/03/22 13:12:40 j_novak
88  * Modification of comments to use doxygen instead of doc++
89  *
90  * Revision 1.9 2004/01/27 14:31:25 j_novak
91  * New method Base_val::mult_sint()
92  *
93  * Revision 1.8 2004/01/27 14:13:58 j_novak
94  * Added method Base_val::mult_x()
95  *
96  * Revision 1.7 2003/10/20 06:41:43 e_gourgoulhon
97  * Corrected documentation.
98  *
99  * Revision 1.6 2003/10/19 19:42:50 e_gourgoulhon
100  * -- member nzone now private
101  * -- introduced new methods get_nzone() and get_b()
102  * -- introduced new methods name_r, name_theta and name_phi.
103  *
104  * Revision 1.5 2003/09/16 08:53:05 j_novak
105  * Addition of the T_LEG_II base (odd in theta, only for odd m) and the
106  * transformation functions from and to T_SIN_P.
107  *
108  * Revision 1.4 2002/10/16 14:36:28 j_novak
109  * Reorganization of #include instructions of standard C++, in order to
110  * use experimental version 3 of gcc.
111  *
112  * Revision 1.3 2002/09/13 09:17:31 j_novak
113  * Modif. commentaires
114  *
115  * Revision 1.2 2002/06/17 14:05:15 j_novak
116  * friend functions are now also declared outside the class definition
117  *
118  * Revision 1.1.1.1 2001/11/20 15:19:27 e_gourgoulhon
119  * LORENE
120  *
121  * Revision 2.20 2001/10/12 14:56:09 novak
122  * Ajout des methodes dsdx(), dsdt(), etc...
123  *
124  * Revision 2.19 2001/05/29 16:08:45 eric
125  * Modif commentaires (mise en conformite Doc++ 3.4.7).
126  *
127  * Revision 2.18 2000/09/28 10:19:50 eric
128  * Modif commentaires (nouvelles bases T_LEG_IP et T_LEG_PI).
129  *
130  * Revision 2.17 2000/09/08 11:42:55 eric
131  * *** empty log message ***
132  *
133  * Revision 2.16 2000/09/08 10:14:11 eric
134  * *** empty log message ***
135  *
136  * Revision 2.15 2000/09/08 10:11:49 eric
137  * *** empty log message ***
138  *
139  * Revision 2.14 2000/09/08 10:10:24 eric
140  * *** empty log message ***
141  *
142  * Revision 2.13 2000/09/08 10:06:19 eric
143  * Ajout des methodes set_base_r, etc... et get_base_r, etc...
144  *
145  * Revision 2.12 1999/12/29 10:49:12 eric
146  * theta_functions et phi_functions declarees const.
147  *
148  * Revision 2.11 1999/12/29 10:36:58 eric
149  * Modif commentaires.
150  *
151  * Revision 2.10 1999/12/28 12:57:44 eric
152  * Introduction des methodes theta_functions et phi_functions.
153  *
154  * Revision 2.9 1999/11/19 14:53:11 phil
155  * *** empty log message ***
156  *
157  * Revision 2.8 1999/11/18 13:48:35 phil
158  * *** empty log message ***
159  *
160  * Revision 2.7 1999/11/18 12:51:12 phil
161  * commentaire de operator*
162  *
163  * Revision 2.6 1999/10/26 13:23:12 phil
164  * *** empty log message ***
165  *
166  * Revision 2.5 1999/10/26 13:18:06 phil
167  * ajout de l'operator*
168  *
169  * Revision 2.4 1999/10/13 15:49:12 eric
170  * *** empty log message ***
171  *
172  * Revision 2.3 1999/10/12 10:02:17 eric
173  * Amelioration des commentaires: description de toutes les bases.
174  *
175  * Revision 2.2 1999/10/01 15:55:58 eric
176  * Depoussierage.
177  * Documentation
178  *
179  * Revision 2.1 1999/09/13 14:38:08 phil
180  * ajout de l'operateur ==
181  *
182  * Revision 2.0 1999/02/22 15:17:47 hyc
183  * *** empty log message ***
184  *
185  *
186  * $Header: /cvsroot/Lorene/C++/Include/base_val.h,v 1.24 2023/06/28 10:04:32 j_novak Exp $
187  *
188  */
189 
190 #include <cstdio>
191 #include <cassert>
192 #include "headcpp.h"
193 
194 
195 #include "type_parite.h"
196 namespace Lorene {
197 class Tbl ;
198 class Mg3d ;
199 
325 class Base_val {
326 
327  // Data
328  // ----
329  protected:
330  int nzone ;
331 
332  public:
334  int* b ;
335 
336  // Constructors - Destructor
337  // -------------------------
338  public:
339  explicit Base_val(int nb_of_domains) ;
340  Base_val(const Base_val& ) ;
341 
343  explicit Base_val(FILE* ) ;
344 
345  ~Base_val() ;
346 
347 
348  // Mutators / assignment
349  // ---------------------
350  public:
351  void set_base_nondef() ;
352 
362  void set_base_r(int l, int base_r) ;
363 
372  void set_base_t(int base_t) ;
373 
382  void set_base_p(int base_p) ;
383 
384  void operator=(const Base_val& ) ;
385 
386  // Accessors
387  // ---------
388  public:
389  bool operator==(const Base_val& ) const ;
390 
392  int get_b(int l) const {
393  assert( (l>=0) && (l<nzone) ) ;
394  return b[l] ;
395  }
396 
403  int get_base_r(int l) const {
404  assert( (l>=0) && (l<nzone) ) ;
405  return b[l] & MSQ_R ;
406  } ;
407 
414  int get_base_t(int l) const {
415  assert( (l>=0) && (l<nzone) ) ;
416  return b[l] & MSQ_T ;
417  } ;
418 
425  int get_base_p(int l) const {
426  assert( (l>=0) && (l<nzone) ) ;
427  return b[l] & MSQ_P ;
428  } ;
429 
440  void name_r(int l, int k, int j, int i, string& basename) const ;
441 
450  void name_theta(int l, int k, int j, string& basename) const ;
451 
458  void name_phi(int l, int k, string& basename) const ;
459 
460 
483  const Tbl& theta_functions(int l, int nt) const ;
484 
499  const Tbl& phi_functions(int l, int np) const ;
500 
502  int get_nzone() const {return nzone; } ;
503 
504  // Manipulation of basis
505  // ----------------------
506  public:
511  void dsdx() ;
512 
517  void sx() ;
518 
523  void mult_x() ;
524 
529  void dsdt() ;
530 
535  void ssint() ;
536 
541  void mult_sint() ;
542 
547  void mult_cost() ;
548 
553  void ylm() ;
554 
558  void give_quant_numbers (int, int, int,
559  int&, int&, int&) const ;
560 
568  int give_lmax(const Mg3d& mgrid, int lz) const ;
569 
570  // Outputs
571  // -------
572  public:
573  void sauve(FILE *) const ;
574 
575  friend ostream& operator<<(ostream& , const Base_val& ) ;
576  friend Base_val operator*(const Base_val&, const Base_val&) ;
577 };
578 ostream& operator<<(ostream& , const Base_val& ) ;
596 Base_val operator*(const Base_val&, const Base_val&) ;
597 
600 }
601 #endif
void operator=(const Base_val &)
Assignment.
Definition: base_val.C:218
void mult_cost()
The basis is transformed as with a multiplication.
int give_lmax(const Mg3d &mgrid, int lz) const
Returns the highest multipole for a given grid.
void sx()
The basis is transformed as with a multiplication.
void set_base_r(int l, int base_r)
Sets the expansion basis for r ( ) functions in a given domain.
Definition: base_val.C:188
void give_quant_numbers(int, int, int, int &, int &, int &) const
Computes the various quantum numbers and 1d radial base.
Lorene prototypes.
Definition: app_hor.h:67
int get_base_t(int l) const
Returns the expansion basis for functions in the domain of index l (e.g.
Definition: base_val.h:414
Base_val operator*(const Base_val &, const Base_val &)
This operator is used when calling multiplication or division of Valeur .
#define MSQ_P
Extraction de l&#39;info sur Phi.
Definition: type_parite.h:156
friend ostream & operator<<(ostream &, const Base_val &)
Display.
Definition: base_val.C:241
void dsdx()
The basis is transformed as with a operation.
friend Base_val operator*(const Base_val &, const Base_val &)
This operator is used when calling multiplication or division of Valeur .
void ylm()
The basis is transformed as with a transformation to basis.
void name_theta(int l, int k, int j, string &basename) const
Name of the basis function in .
~Base_val()
Destructor.
Definition: base_val.C:180
void set_base_t(int base_t)
Sets the expansion basis for functions in all domains.
Definition: base_val.C:198
int nzone
Number of domains (zones)
Definition: base_val.h:330
int get_base_r(int l) const
Returns the expansion basis for r ( ) functions in the domain of index l (e.g.
Definition: base_val.h:403
#define MSQ_T
Extraction de l&#39;info sur Theta.
Definition: type_parite.h:154
void sauve(FILE *) const
Save in a file.
Definition: base_val.C:231
void set_base_p(int base_p)
Sets the expansion basis for functions in all domains.
Definition: base_val.C:207
int get_nzone() const
Returns the number of domains.
Definition: base_val.h:502
int * b
Array (size: nzone ) of the spectral basis in each domain.
Definition: base_val.h:334
#define MSQ_R
Extraction de l&#39;info sur R.
Definition: type_parite.h:152
int get_base_p(int l) const
Returns the expansion basis for functions in the domain of index l (e.g.
Definition: base_val.h:425
bool operator==(const Base_val &) const
Comparison operator.
Definition: base_val.C:338
void name_r(int l, int k, int j, int i, string &basename) const
Name of the basis function in r ( )
void name_phi(int l, int k, string &basename) const
Name of the basis function in .
void ssint()
The basis is transformed as with a multiplication.
const Tbl & theta_functions(int l, int nt) const
Values of the theta basis functions at the theta collocation points.
Base_val(int nb_of_domains)
Standard constructor.
Definition: base_val.C:153
void dsdt()
The basis is transformed as with a operation.
Multi-domain grid.
Definition: grilles.h:279
Bases of the spectral expansions.
Definition: base_val.h:325
int get_b(int l) const
Returns the code for the expansion basis in domain no. l.
Definition: base_val.h:392
void mult_x()
The basis is transformed as with a multiplication by .
void set_base_nondef()
Sets the spectral bases to NONDEF.
Definition: base_val.C:329
const Tbl & phi_functions(int l, int np) const
Values of the phi basis functions at the phi collocation points.
Basic array class.
Definition: tbl.h:164
void mult_sint()
The basis is transformed as with a multiplication.