LORENE
base_val.C
1/*
2 * Methods of class Base_val
3 *
4 * (see file base_val.h for documentation)
5 *
6 */
7
8/*
9 * Copyright (c) 1999-2000 Jean-Alain Marck
10 * Copyright (c) 1999-2001 Eric Gourgoulhon
11 * Copyright (c) 1999-2001 Philippe Grandclement
12 *
13 * This file is part of LORENE.
14 *
15 * LORENE is free software; you can redistribute it and/or modify
16 * it under the terms of the GNU General Public License as published by
17 * the Free Software Foundation; either version 2 of the License, or
18 * (at your option) any later version.
19 *
20 * LORENE is distributed in the hope that it will be useful,
21 * but WITHOUT ANY WARRANTY; without even the implied warranty of
22 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
23 * GNU General Public License for more details.
24 *
25 * You should have received a copy of the GNU General Public License
26 * along with LORENE; if not, write to the Free Software
27 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
28 *
29 */
30
31
32
33
34/*
35 * $Id: base_val.C,v 1.18 2016/12/05 16:17:44 j_novak Exp $
36 * $Log: base_val.C,v $
37 * Revision 1.18 2016/12/05 16:17:44 j_novak
38 * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
39 *
40 * Revision 1.17 2014/10/13 08:52:38 j_novak
41 * Lorene classes and functions now belong to the namespace Lorene.
42 *
43 * Revision 1.16 2014/10/06 15:12:56 j_novak
44 * Modified #include directives to use c++ syntax.
45 *
46 * Revision 1.15 2013/01/11 08:20:11 j_novak
47 * New radial spectral bases with Legendre polynomials (R_LEG, R_LEGP, R_LEGI).
48 *
49 * Revision 1.14 2012/01/17 14:44:19 j_penner
50 * Modified phi variables to only use 16 integers in arrays
51 *
52 * Revision 1.13 2009/10/23 12:55:16 j_novak
53 * New base T_LEG_MI
54 *
55 * Revision 1.12 2009/10/08 16:20:13 j_novak
56 * Addition of new bases T_COS and T_SIN.
57 *
58 * Revision 1.11 2008/08/19 06:41:59 j_novak
59 * Minor modifications to avoid warnings with gcc 4.3. Most of them concern
60 * cast-type operations, and constant strings that must be defined as const char*
61 *
62 * Revision 1.10 2008/02/18 13:53:38 j_novak
63 * Removal of special indentation instructions.
64 *
65 * Revision 1.9 2007/12/11 15:28:09 jl_cornou
66 * Jacobi(0,2) polynomials partially implemented
67 *
68 * Revision 1.8 2004/12/17 13:35:01 m_forot
69 * Add the case T_LEG
70 *
71 * Revision 1.7 2004/11/04 15:19:02 e_gourgoulhon
72 * operator<< : added names R_CHEBPI_P, R_CHEBPI_I, T_COSSIN_C, T_COSSIN_S.
73 *
74 * Revision 1.6 2004/08/24 09:14:41 p_grandclement
75 * Addition of some new operators, like Poisson in 2d... It now requieres the
76 * GSL library to work.
77 *
78 * Also, the way a variable change is stored by a Param_elliptic is changed and
79 * no longer uses Change_var but rather 2 Scalars. The codes using that feature
80 * will requiere some modification. (It should concern only the ones about monopoles)
81 *
82 * Revision 1.5 2003/09/16 08:54:09 j_novak
83 * Addition of the T_LEG_II base (odd in theta, only for odd m) and the
84 * transformation functions to and from the T_SIN_P base.
85 *
86 * Revision 1.4 2002/11/13 15:05:59 j_novak
87 * Affichage de la base T_COS
88 *
89 * Revision 1.3 2002/10/16 14:36:30 j_novak
90 * Reorganization of #include instructions of standard C++, in order to
91 * use experimental version 3 of gcc.
92 *
93 * Revision 1.2 2001/12/04 21:27:52 e_gourgoulhon
94 *
95 * All writing/reading to a binary file are now performed according to
96 * the big endian convention, whatever the system is big endian or
97 * small endian, thanks to the functions fwrite_be and fread_be
98 *
99 * Revision 1.1.1.1 2001/11/20 15:19:28 e_gourgoulhon
100 * LORENE
101 *
102 * Revision 2.8 2000/09/28 10:20:19 eric
103 * Affichage: nouvelles bases T_LEG_IP et T_LEG_PI.
104 *
105 * Revision 2.7 2000/09/08 10:06:45 eric
106 * Ajout des methodes set_base_r, etc...
107 *
108 * Revision 2.6 1999/12/28 12:57:26 eric
109 * Reorganisation des headers.
110 *
111 * Revision 2.5 1999/10/12 10:02:51 eric
112 * Implementation de sauve().
113 * Modif de << : affichage du nom des bases et non plus de leur numero.
114 *
115 * Revision 2.4 1999/10/01 15:56:21 eric
116 * Depoussierage.
117 * Documentation.
118 *
119 * Revision 2.3 1999/09/13 14:57:06 phil
120 * ajout de l'operateur ==
121 *
122 * Revision 2.2 1999/03/02 15:22:23 eric
123 * Affichage des bases.
124 *
125 * Revision 2.1 1999/03/01 14:54:01 eric
126 * *** empty log message ***
127 *
128 * Revision 2.0 1999/02/22 15:19:20 hyc
129 * *** empty log message ***
130 *
131 *
132 * $Header: /cvsroot/Lorene/C++/Source/Base_val/base_val.C,v 1.18 2016/12/05 16:17:44 j_novak Exp $
133 *
134 */
135
136// Headers C
137#include <cstdio>
138#include <cassert>
139
140// Headers Lorene
141#include "headcpp.h"
142#include "type_parite.h"
143#include "base_val.h"
144#include "utilitaires.h"
145
146
147 //---------------//
148 // Constructeurs //
149 //---------------//
150
151// Constructeur
152namespace Lorene {
154 b = new int[nzone] ;
155 for (int i=0 ; i<nzone ; i++) { // Boucle sur les zones
156 b[i] = NONDEF ;
157 }
158}
159
160// Copie
162 b = new int[nzone] ;
163 for (int i=0 ; i<nzone ; i++) { // Boucle sur les zones
164 b[i] = bi.b[i] ;
165 }
166}
167
168// From file
170 fread_be(&nzone, sizeof(int), 1, fd) ; // nzone
171 b = new int[nzone] ;
172 fread_be(b, sizeof(int), nzone, fd) ; // b[]
173}
174
175 //--------------//
176 // Destructeurs //
177 //--------------//
178
179// Destructeur
181 delete [] b ;
182}
183
184 //-------------//
185 // Affectation //
186 //-------------//
187
188void Base_val::set_base_r(int l, int base_r) {
189
190 assert( (l>=0) && (l<nzone) ) ;
191
192 int base_t = b[l] & MSQ_T ;
193 int base_p = b[l] & MSQ_P ;
194 b[l] = base_p | base_t | base_r ;
195
196}
197
198void Base_val::set_base_t(int base_t) {
199
200 for (int l=0; l<nzone; l++) {
201 int base_r = b[l] & MSQ_R ;
202 int base_p = b[l] & MSQ_P ;
203 b[l] = base_p | base_t | base_r ;
204 }
205}
206
207void Base_val::set_base_p(int base_p) {
208
209 for (int l=0; l<nzone; l++) {
210 int base_r = b[l] & MSQ_R ;
211 int base_t = b[l] & MSQ_T ;
212 b[l] = base_p | base_t | base_r ;
213 }
214}
215
216
217// From Base_val
219 // Protection
220 assert(nzone == bi.nzone) ;
221 for (int i=0 ; i<nzone ; i++) { // Boucle sur les zones
222 b[i] = bi.b[i] ;
223 }
224}
225
226 //------------//
227 // Sauvegarde //
228 //------------//
229
230// Save in a file
231void Base_val::sauve(FILE* fd) const {
232 fwrite_be(&nzone, sizeof(int), 1, fd) ; // nzone
233 fwrite_be(b, sizeof(int), nzone, fd) ; // b[]
234}
235
236 //------------//
237 // Impression //
238 //------------//
239
240// Operateurs <<
241ostream& operator<<(ostream& o, const Base_val & bi) {
242
243 static bool premier_appel = true ;
244 static const char* nom_r[MAX_BASE] ;
245 static const char* nom_t[MAX_BASE] ;
246 static const char* nom_p[MAX_BASE_2] ;
247
248 if (premier_appel) { // First call initializations
249
250 premier_appel = false ;
251
252 for (int i=0; i<MAX_BASE; i++) {
253 nom_r[i] = "UNDEFINED" ;
254 nom_t[i] = "UNDEFINED" ;
255 if(i%2==0){
256 nom_p[i/2] = "UNDEFINED" ; // saves a loop
257 }
258 }
259
260 nom_r[R_CHEB >> TRA_R] = "R_CHEB " ;
261 nom_r[R_CHEBU >> TRA_R] = "R_CHEBU " ;
262 nom_r[R_CHEBP >> TRA_R] = "R_CHEBP " ;
263 nom_r[R_CHEBI >> TRA_R] = "R_CHEBI " ;
264 nom_r[R_CHEBPIM_P >> TRA_R] = "R_CHEBPIM_P" ;
265 nom_r[R_CHEBPIM_I >> TRA_R] = "R_CHEBPIM_I" ;
266 nom_r[R_CHEBPI_P >> TRA_R] = "R_CHEBPI_P " ;
267 nom_r[R_CHEBPI_I >> TRA_R] = "R_CHEBPI_I " ;
268 nom_r[R_LEG >> TRA_R] = "R_LEG " ;
269 nom_r[R_LEGP >> TRA_R] = "R_LEGP " ;
270 nom_r[R_LEGI >> TRA_R] = "R_LEGI " ;
271 nom_r[R_JACO02 >> TRA_R] = "R_JACO02 " ;
272
273 nom_t[T_COS >> TRA_T] = "T_COS " ;
274 nom_t[T_SIN >> TRA_T] = "T_SIN " ;
275 nom_t[T_COS_P >> TRA_T] = "T_COS_P " ;
276 nom_t[T_COS_I >> TRA_T] = "T_COS_I " ;
277 nom_t[T_SIN_P >> TRA_T] = "T_SIN_P " ;
278 nom_t[T_SIN_I >> TRA_T] = "T_SIN_I " ;
279 nom_t[T_COSSIN_CP >> TRA_T] = "T_COSSIN_CP" ;
280 nom_t[T_COSSIN_SI >> TRA_T] = "T_COSSIN_SI" ;
281 nom_t[T_COSSIN_SP >> TRA_T] = "T_COSSIN_SP" ;
282 nom_t[T_COSSIN_CI >> TRA_T] = "T_COSSIN_CI" ;
283 nom_t[T_COSSIN_C >> TRA_T] = "T_COSSIN_C " ;
284 nom_t[T_COSSIN_S >> TRA_T] = "T_COSSIN_S " ;
285 nom_t[T_LEG >> TRA_T] = "T_LEG " ;
286 nom_t[T_LEG_MP >> TRA_T] = "T_LEG_MP " ;
287 nom_t[T_LEG_MI >> TRA_T] = "T_LEG_MI " ;
288 nom_t[T_LEG_P >> TRA_T] = "T_LEG_P " ;
289 nom_t[T_LEG_PP >> TRA_T] = "T_LEG_PP " ;
290 nom_t[T_LEG_I >> TRA_T] = "T_LEG_I " ;
291 nom_t[T_LEG_IP >> TRA_T] = "T_LEG_IP " ;
292 nom_t[T_LEG_PI >> TRA_T] = "T_LEG_PI " ;
293 nom_t[T_LEG_II >> TRA_T] = "T_LEG_II " ;
294 nom_t[T_CL_COS_P >> TRA_T] = "T_CL_COS_P " ;
295 nom_t[T_CL_SIN_P >> TRA_T] = "T_CL_SIN_P " ;
296 nom_t[T_CL_COS_I >> TRA_T] = "T_CL_COS_I " ;
297 nom_t[T_CL_SIN_I >> TRA_T] = "T_CL_SIN_I " ;
298
299 nom_p[P_COSSIN >> TRA_P] = "P_COSSIN " ;
300 nom_p[P_COSSIN_P >> TRA_P] = "P_COSSIN_P " ;
301 nom_p[P_COSSIN_I >> TRA_P] = "P_COSSIN_I " ;
302
303
304 } // End of first call operations
305
306
307 // Intro - Nombre de zones
308 int nzone = bi.nzone ;
309 o << "Bases of spectral expansions: " ;
310 for (int l=0 ; l<nzone ; l++) {
311 int base_r = (bi.b[l] & MSQ_R) >> TRA_R ;
312 int base_t = (bi.b[l] & MSQ_T) >> TRA_T ;
313 int base_p = (bi.b[l] & MSQ_P) >> TRA_P ;
314 o << endl ;
315 o << "Domain #" << l << " : r: " << nom_r[base_r]
316 << ", theta: " << nom_t[base_t]
317 << ", phi: " << nom_p[base_p] ;
318 }
319 o << endl ;
320
321 //Termine
322 return o ;
323}
324
325 //----------------------//
326 // Manipulation de base //
327 //----------------------//
328
330 for (int i=0 ; i<nzone ; i++) {
331 b[i] = NONDEF ;
332 }
333}
334 //----------------------//
335 // operateur logique //
336 //----------------------//
337
338bool Base_val::operator== (const Base_val& c2) const {
339
340 return (*b == *c2.b) ;
341}
342}
Base_val(int nb_of_domains)
Standard constructor.
Definition base_val.C:153
friend ostream & operator<<(ostream &, const Base_val &)
Display.
Definition base_val.C:241
void sauve(FILE *) const
Save in a file.
Definition base_val.C:231
void set_base_t(int base_t)
Sets the expansion basis for functions in all domains.
Definition base_val.C:198
~Base_val()
Destructor.
Definition base_val.C:180
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
bool operator==(const Base_val &) const
Comparison operator.
Definition base_val.C:338
int * b
Array (size: nzone ) of the spectral basis in each domain.
Definition base_val.h:334
int nzone
Number of domains (zones).
Definition base_val.h:330
void set_base_nondef()
Sets the spectral bases to NONDEF.
Definition base_val.C:329
void operator=(const Base_val &)
Assignment.
Definition base_val.C:218
void set_base_p(int base_p)
Sets the expansion basis for functions in all domains.
Definition base_val.C:207
#define MAX_BASE_2
Smaller maximum bases used for phi (and higher dimensions for now).
#define T_CL_COS_I
CL of odd cosines.
#define T_CL_SIN_I
CL of odd sines.
#define TRA_T
Translation en Theta, used for a bitwise shift (in hex).
#define T_LEG_MP
fct. de Legendre associees avec m pair
#define P_COSSIN_P
dev. sur Phi = 2*phi, freq. paires
#define R_LEGP
base de Legendre paire (rare) seulement
#define MAX_BASE
Nombre max. de bases differentes.
#define T_LEG_PI
fct. de Legendre associees paires avec m impair
#define R_CHEBU
base de Chebychev ordinaire (fin), dev. en 1/r
#define R_JACO02
base de Jacobi(0,2) ordinaire (finjac)
#define T_COSSIN_SP
sin pair-cos impair alternes, sin pour m=0
#define R_LEGI
base de Legendre impaire (rare) seulement
#define R_CHEBI
base de Cheb. impaire (rare) seulement
#define T_LEG
fct. de Legendre associees
#define T_SIN_P
dev. sin seulement, harmoniques paires
#define MSQ_R
Extraction de l'info sur R.
#define T_COSSIN_S
dev. cos-sin alternes, sin pour m=0
#define R_CHEBPIM_I
Cheb. pair-impair suivant m, impair pour m=0.
#define T_COSSIN_SI
sin impair-cos pair alternes, sin pour m=0
#define R_CHEBPI_I
Cheb. pair-impair suivant l impair pour l=0.
#define R_LEG
base de Legendre ordinaire (fin)
#define T_CL_SIN_P
CL of even sines.
#define T_LEG_P
fct. de Legendre associees paires
#define T_LEG_IP
fct. de Legendre associees impaires avec m pair
#define NONDEF
base inconnue
#define T_COS_P
dev. cos seulement, harmoniques paires
#define T_COSSIN_CI
cos impair-sin pair alternes, cos pour m=0
#define P_COSSIN
dev. standart
#define P_COSSIN_I
dev. sur Phi = 2*phi, freq. impaires
#define R_CHEBPIM_P
Cheb. pair-impair suivant m, pair pour m=0.
#define TRA_R
Translation en R, used for a bitwise shift (in hex).
#define T_CL_COS_P
CL of even cosines.
#define T_LEG_MI
fct. de Legendre associees avec m impair
#define MSQ_T
Extraction de l'info sur Theta.
#define T_COSSIN_CP
cos pair-sin impair alternes, cos pour m=0
#define T_LEG_II
fct. de Legendre associees impaires avec m impair
#define R_CHEB
base de Chebychev ordinaire (fin)
#define T_LEG_I
fct. de Legendre associees impaires
#define T_SIN_I
dev. sin seulement, harmoniques impaires
#define T_COS
dev. cos seulement
#define R_CHEBP
base de Cheb. paire (rare) seulement
#define MSQ_P
Extraction de l'info sur Phi.
#define T_SIN
dev. sin seulement
#define T_LEG_PP
fct. de Legendre associees paires avec m pair
#define T_COS_I
dev. cos seulement, harmoniques impaires
#define T_COSSIN_C
dev. cos-sin alternes, cos pour m=0
#define TRA_P
Translation en Phi, used for a bitwise shift (in hex).
#define R_CHEBPI_P
Cheb. pair-impair suivant l pair pour l=0.
int fread_be(int *aa, int size, int nb, FILE *fich)
Reads integer(s) from a binary file according to the big endian convention.
Definition fread_be.C:72
int fwrite_be(const int *aa, int size, int nb, FILE *fich)
Writes integer(s) into a binary file according to the big endian convention.
Definition fwrite_be.C:73
Lorene prototypes.
Definition app_hor.h:67