LORENE
donne_lm.C
1 /*
2  * Copyright (c) 1999-2001 Philippe Grandclement
3  * Copyright (c) 2000-2001 Eric Gourgoulhon
4  * Copyright (c) 2000-2001 Jerome Novak
5  *
6  * This file is part of LORENE.
7  *
8  * LORENE is free software; you can redistribute it and/or modify
9  * it under the terms of the GNU General Public License as published by
10  * the Free Software Foundation; either version 2 of the License, or
11  * (at your option) any later version.
12  *
13  * LORENE is distributed in the hope that it will be useful,
14  * but WITHOUT ANY WARRANTY; without even the implied warranty of
15  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16  * GNU General Public License for more details.
17  *
18  * You should have received a copy of the GNU General Public License
19  * along with LORENE; if not, write to the Free Software
20  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
21  *
22  */
23 
24 
25 
26 
27 
28 /*
29  * $Id: donne_lm.C,v 1.11 2016/12/05 16:18:02 j_novak Exp $
30  * $Log: donne_lm.C,v $
31  * Revision 1.11 2016/12/05 16:18:02 j_novak
32  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
33  *
34  * Revision 1.10 2014/10/13 08:53:12 j_novak
35  * Lorene classes and functions now belong to the namespace Lorene.
36  *
37  * Revision 1.9 2014/10/06 15:16:02 j_novak
38  * Modified #include directives to use c++ syntax.
39  *
40  * Revision 1.8 2009/10/26 10:48:37 j_novak
41  * Completed the T_LEG_MI case.
42  *
43  * Revision 1.7 2009/10/23 12:54:47 j_novak
44  * New base T_LEG_MI
45  *
46  * Revision 1.6 2009/10/13 19:45:01 j_novak
47  * New base T_LEG_MP.
48  *
49  * Revision 1.5 2005/02/18 13:14:13 j_novak
50  * Changing of malloc/free to new/delete + suppression of some unused variables
51  * (trying to avoid compilation warnings).
52  *
53  * Revision 1.4 2004/11/23 15:13:50 m_forot
54  * Added the bases for the cases without any equatorial symmetry
55  * (T_COSSIN_C, T_COSSIN_S, T_LEG, R_CHEBPI_P, R_CHEBPI_I).
56  *
57  * Revision 1.3 2003/09/16 12:11:59 j_novak
58  * Added the base T_LEG_II.
59  *
60  * Revision 1.2 2002/10/16 14:36:54 j_novak
61  * Reorganization of #include instructions of standard C++, in order to
62  * use experimental version 3 of gcc.
63  *
64  * Revision 1.1.1.1 2001/11/20 15:19:29 e_gourgoulhon
65  * LORENE
66  *
67  * Revision 3.0 2000/10/09 09:15:01 novak
68  * correction pour les cas de bases en r alternees
69  *
70  * Revision 2.8 2000/10/04 14:55:33 eric
71  * Ajout des bases T_LEG_IP et T_LEG_PI.
72  *
73  * Revision 2.7 1999/12/16 16:41:33 phil
74  * *** empty log message ***
75  *
76  * Revision 2.6 1999/12/16 16:23:39 phil
77  * vire un assert
78  *
79  * Revision 2.5 1999/12/16 16:21:38 phil
80  * correction cas nt = 1
81  *
82  * Revision 2.4 1999/09/16 12:06:11 phil
83  * correction des cas antisymetriques en z=0
84  *
85  * Revision 2.3 1999/09/14 17:52:59 phil
86  * *** empty log message ***
87  *
88  * Revision 2.2 1999/09/14 17:44:16 phil
89  * *** empty log message ***
90  *
91  * Revision 2.1 1999/04/13 13:50:01 phil
92  * *** empty log message ***
93  *
94  * Revision 2.0 1999/04/13 13:31:01 phil
95  * *** empty log message ***
96  *
97  * Revision 1.1 1999/04/13 13:30:28 phil
98  * Initial revision
99  *
100  *
101  * $Header: /cvsroot/Lorene/C++/Source/Non_class_members/Coef/donne_lm.C,v 1.11 2016/12/05 16:18:02 j_novak Exp $
102  *
103  */
104 
105 // Entetes C
106 #include <cstdlib>
107 
108 // Entete Lorene
109 #include "headcpp.h"
110 #include "type_parite.h"
111 #include "base_val.h"
112 
113 /*
114  * Fonction affection les nombres l_quant, m_quant et la base en r
115  *
116  * ENTREES : nz : le nombre de zones
117  * zone : la zone de travail
118  * j et k : indices en theta et phi, respectivement
119  * base : la base de developpement
120  *
121  * SORTIES : les variables m_quant, l_quant et base_r.
122  *
123  */
124 
125 //-----------------------------------------------------------------
126 // Developpement en P_COSSIN pour phi et T_LEG en theta
127 //-------------------------------------------------------------------
128 
129 namespace Lorene {
130 void donne_lm_nonsymTP (int j, int k, int &m_quant, int &l_quant) {
131 
132  m_quant = (k%2 == 0) ? k/2 : (k-1)/2;
133  l_quant = j ;
134 
135 }
136 
137 
138  //-----------------------------------------------------------------
139  // Developpement en P_COSSIN pour phi et T_LEG_P en theta
140  //-------------------------------------------------------------------
141 
142 void donne_lm_nonsym (int j, int k, int &m_quant, int &l_quant) {
143 
144  m_quant = (k%2 == 0) ? k/2 : (k-1)/2;
145  l_quant = (m_quant%2 == 0) ? 2*j : 2*j+1 ;
146 
147 }
148 
149  //-----------------------------------------------------------------
150  // Developpement en P_COSSIN pour phi et T_LEG_I en theta
151  //-------------------------------------------------------------------
152 
153 void donne_lm_nonsym_anti (int j, int k, int &m_quant, int &l_quant) {
154 
155  m_quant = (k%2 == 0) ? k/2 : (k-1)/2;
156  l_quant = (m_quant%2 == 1) ? 2*j : 2*j+1 ;
157 
158 }
159 
160  //------------------------------------------------------
161  // Developpement en P_COSSIN_P pour phi et T_LEG_PP en theta
162  //-------------------------------------------------------
163 
164 void donne_lm_sym (int j, int k, int &m_quant, int &l_quant) {
165 
166  m_quant = (k%2 == 0) ? k : k-1;
167  l_quant = 2*j ;
168 
169 }
170 
171 
172  //-------------------------------------------------------
173  // Developpement en P_COSSIN_P pour phi et T_LEG_IP en theta
174  //---------------------------------------------------------
175 
176 void donne_lm_t_leg_ip (int j, int k, int &m_quant, int &l_quant) {
177 
178  m_quant = (k%2 == 0) ? k : k-1 ;
179  l_quant = 2*j+1 ;
180 
181 }
182 
183 
184  //----------------------------------------------------------
185  // Developpement en P_COSSIN_P pour phi et T_LEG_MP en theta
186  //------------------------------------------------------------
187 
188 void donne_lm_t_leg_mp (int j, int k, int &m_quant, int &l_quant) {
189 
190  m_quant = (k%2 == 0) ? k : k-1;
191  l_quant = j ;
192 
193 }
194 
195  //----------------------------------------------------------
196  // Developpement en P_COSSIN_I pour phi et T_LEG_MI en theta
197  //------------------------------------------------------------
198 
199 void donne_lm_t_leg_mi (int j, int k, int &m_quant, int &l_quant) {
200 
201  m_quant = 2*((k-1)/2 ) + 1 ;
202  l_quant = j ;
203 
204 }
205 
206  //-------------------------------------------------------
207  // Developpement en P_COSSIN_I pour phi et T_LEG_PI en theta
208  //---------------------------------------------------------
209 
210 void donne_lm_t_leg_pi (int j, int k, int &m_quant, int &l_quant) {
211 
212  if (k<=2) {
213  m_quant = 1 ;
214  }
215  else{
216  m_quant = (k%2 == 0) ? k-1 : k ;
217  }
218 
219  l_quant = 2*j+1 ;
220 
221 }
222 
223  //-------------------------------------------------------
224  // Developpement en P_COSSIN_I pour phi et T_LEG_II en theta
225  //---------------------------------------------------------
226 
227 void donne_lm_t_leg_ii (int j, int k, int &m_quant, int &l_quant) {
228 
229  if (k<=2) {
230  m_quant = 1 ;
231  }
232  else{
233  m_quant = (k%2 == 0) ? k-1 : k ;
234  }
235 
236  l_quant = 2*j ;
237 
238 }
239 
240 
241 
242  //-----------------------------
243  // La fonction
244  //-------------------------------
245 
246 void donne_lm (int nz, int zone, int j, int k, Base_val base,
247  int &m_quant, int &l_quant, int& base_r) {
248 
249  //verifications :
250  assert (zone >= 0) ;
251  assert (zone < nz) ;
252 
253  int base_t = (base.b[zone] & MSQ_T) ;
254  int base_p = (base.b[zone] & MSQ_P) ;
255  base_r = (base.b[zone] & MSQ_R) ;
256 
257  switch (base_p) {
258  case P_COSSIN :
259  // cas sym ou antisym en z=0 ...
260  switch (base_t) {
261 
262  case T_LEG :
263  donne_lm_nonsymTP (j, k, m_quant, l_quant) ;
264  break ;
265 
266  case T_LEG_P :
267  donne_lm_nonsym (j, k, m_quant, l_quant) ;
268  break ;
269 
270  case T_LEG_I :
271  donne_lm_nonsym_anti (j, k, m_quant, l_quant) ;
272  break ;
273 
274  default :
275  cout << "donne_lm : cas inconnu ..." << endl ;
276  abort() ;
277  break ;
278  }
279  break ;
280 
281  case P_COSSIN_P :
282  switch (base_t) {
283 
284  case T_LEG_PP :
285  donne_lm_sym (j, k, m_quant, l_quant) ;
286  break ;
287 
288  case T_LEG_MP :
289  donne_lm_t_leg_mp (j, k, m_quant, l_quant) ;
290  break ;
291 
292  case T_LEG_IP :
293  donne_lm_t_leg_ip (j, k, m_quant, l_quant);
294  break ;
295 
296  default :
297  cout << "donne_lm : cas inconnu ..." << endl ;
298  abort() ;
299  break ;
300  }
301  break ;
302 
303  case P_COSSIN_I :
304  switch (base_t) {
305 
306  case T_LEG_PI :
307  donne_lm_t_leg_pi (j, k, m_quant, l_quant) ;
308  break ;
309 
310  case T_LEG_II :
311  donne_lm_t_leg_ii (j, k, m_quant, l_quant) ;
312  break ;
313 
314  case T_LEG_MI :
315  donne_lm_t_leg_mp (j, k, m_quant, l_quant) ;
316  break ;
317 
318  default :
319  cout << "donne_lm : cas inconnu ..." << endl ;
320  abort() ;
321  break ;
322 
323  }
324  break ;
325 
326  default :
327  cout << "donne_lm : cas inconnu ..." << endl ;
328  cout << nz << endl ; // to avoid compilation warnings ...
329  abort() ;
330  break ;
331  }
332  switch (base_r) {
333 
334  case R_CHEBPI_P :
335  base_r = (l_quant%2 == 0) ? R_CHEBP : R_CHEBI ;
336  break ;
337 
338  case R_CHEBPI_I :
339  base_r = (l_quant%2 == 1) ? R_CHEBP : R_CHEBI ;
340  break ;
341 
342  case R_CHEBPIM_P :
343  base_r = (m_quant%2 == 0) ? R_CHEBP : R_CHEBI ;
344  break ;
345 
346  case R_CHEBPIM_I :
347  base_r = (m_quant%2 == 1) ? R_CHEBP : R_CHEBI ;
348  break ;
349 
350  }
351 }
352 }
#define T_LEG
fct. de Legendre associees
Definition: type_parite.h:236
#define T_LEG_MP
fct. de Legendre associees avec m pair
Definition: type_parite.h:238
#define P_COSSIN
dev. standart
Definition: type_parite.h:245
#define T_LEG_PI
fct. de Legendre associees paires avec m impair
Definition: type_parite.h:224
#define R_CHEBPI_I
Cheb. pair-impair suivant l impair pour l=0.
Definition: type_parite.h:174
#define T_LEG_MI
fct. de Legendre associees avec m impair
Definition: type_parite.h:240
Lorene prototypes.
Definition: app_hor.h:67
#define MSQ_P
Extraction de l&#39;info sur Phi.
Definition: type_parite.h:156
#define T_LEG_I
fct. de Legendre associees impaires
Definition: type_parite.h:220
#define R_CHEBI
base de Cheb. impaire (rare) seulement
Definition: type_parite.h:170
#define R_CHEBP
base de Cheb. paire (rare) seulement
Definition: type_parite.h:168
#define MSQ_T
Extraction de l&#39;info sur Theta.
Definition: type_parite.h:154
#define MSQ_R
Extraction de l&#39;info sur R.
Definition: type_parite.h:152
#define T_LEG_IP
fct. de Legendre associees impaires avec m pair
Definition: type_parite.h:222
#define R_CHEBPIM_I
Cheb. pair-impair suivant m, impair pour m=0.
Definition: type_parite.h:178
#define R_CHEBPIM_P
Cheb. pair-impair suivant m, pair pour m=0.
Definition: type_parite.h:176
#define T_LEG_P
fct. de Legendre associees paires
Definition: type_parite.h:216
#define R_CHEBPI_P
Cheb. pair-impair suivant l pair pour l=0.
Definition: type_parite.h:172
#define P_COSSIN_I
dev. sur Phi = 2*phi, freq. impaires
Definition: type_parite.h:249
#define P_COSSIN_P
dev. sur Phi = 2*phi, freq. paires
Definition: type_parite.h:247
#define T_LEG_II
fct. de Legendre associees impaires avec m impair
Definition: type_parite.h:226
#define T_LEG_PP
fct. de Legendre associees paires avec m pair
Definition: type_parite.h:218