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
|
/*
* Copyright (c) 2004 Philippe Grandclement
*
* 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
*
*/
char base_val_quantum_C[] = "$Header: /cvsroot/Lorene/C++/Source/Base_val/base_val_quantum.C,v 1.10 2014/10/13 08:52:39 j_novak Exp $" ;
/*
* $Id: base_val_quantum.C,v 1.10 2014/10/13 08:52:39 j_novak Exp $
* $Log: base_val_quantum.C,v $
* Revision 1.10 2014/10/13 08:52:39 j_novak
* Lorene classes and functions now belong to the namespace Lorene.
*
* Revision 1.9 2014/10/06 15:12:57 j_novak
* Modified #include directives to use c++ syntax.
*
* Revision 1.8 2013/01/11 08:20:11 j_novak
* New radial spectral bases with Legendre polynomials (R_LEG, R_LEGP, R_LEGI).
*
* Revision 1.7 2009/10/23 12:55:16 j_novak
* New base T_LEG_MI
*
* Revision 1.6 2009/10/08 16:20:13 j_novak
* Addition of new bases T_COS and T_SIN.
*
* Revision 1.5 2007/12/11 15:28:09 jl_cornou
* Jacobi(0,2) polynomials partially implemented
*
* Revision 1.4 2005/09/07 13:09:50 j_novak
* New method for determining the highest multipole that can be described on a 3D
* grid.
*
* Revision 1.3 2004/11/23 15:08:01 m_forot
* Added the bases for the cases without any equatorial symmetry
* (T_COSSIN_C, T_COSSIN_S, T_LEG, R_CHEBPI_P, R_CHEBPI_I).
*
* Revision 1.2 2004/08/30 16:27:59 r_prix
* added #include <stdlib.h> (got ERROR 'abort' is undefined without this...)
*
* Revision 1.1 2004/08/24 09:14:41 p_grandclement
* Addition of some new operators, like Poisson in 2d... It now requieres the
* GSL library to work.
*
* Also, the way a variable change is stored by a Param_elliptic is changed and
* no longer uses Change_var but rather 2 Scalars. The codes using that feature
* will requiere some modification. (It should concern only the ones about monopoles)
*
*
* $Header: /cvsroot/Lorene/C++/Source/Base_val/base_val_quantum.C,v 1.10 2014/10/13 08:52:39 j_novak Exp $
*
*/
// Headers C
#include <cstdio>
#include <cassert>
#include <cstdlib>
// Headers Lorene
#include "grilles.h"
#include "base_val.h"
#include "utilitaires.h"
namespace Lorene {
void Base_val::give_quant_numbers (int l, int k, int j,
int& m_quant, int& l_quant, int& base_r_1d) const {
int base_p = get_base_p(l) ;
int base_t = get_base_t(l) ;
int base_r = get_base_r(l) ;
switch (base_p) {
case P_COSSIN :
m_quant = k/2 ;
break;
case P_COSSIN_P :
if (k%2 == 0)
m_quant = k ;
else
m_quant = (k-1) ;
break;
case P_COSSIN_I :
m_quant = 2*( (k-1) / 2) + 1 ;
break;
default:
cout << "Unknown basis in phi in give_quant_numbers ..." << endl ;
abort() ;
break ;
}
switch (base_t) {
case T_COS_P :
l_quant = 2*j ;
break;
case T_SIN_P :
l_quant = 2*j ;
break;
case T_COS_I :
l_quant = 2*j+1 ;
break;
case T_SIN_I :
l_quant = 2*j+1 ;
break;
case T_COSSIN_CP :
if (m_quant%2 == 0)
l_quant = 2*j ;
else
l_quant = 2*j+1 ;
break ;
case T_COSSIN_SP :
if (m_quant%2 == 0)
l_quant = 2*j ;
else
l_quant = 2*j+1 ;
break ;
case T_COSSIN_CI :
if (m_quant%2 == 0)
l_quant = 2*j+1 ;
else
l_quant = 2*j ;
break ;
case T_COSSIN_SI :
if (m_quant%2 == 0)
l_quant = 2*j+1 ;
else
l_quant = 2*j ;
break ;
case T_COSSIN_C :
l_quant = j ;
break ;
case T_COSSIN_S :
l_quant = j ;
break ;
case T_COS :
l_quant = j ;
break ;
case T_LEG_P :
if (m_quant%2 == 0)
l_quant = 2*j ;
else
l_quant = 2*j+1 ;
break ;
case T_LEG_PP :
l_quant = 2*j ;
break ;
case T_LEG_I :
if (m_quant%2 == 0)
l_quant = 2*j+1 ;
else
l_quant = 2*j ;
break ;
case T_LEG_IP :
l_quant = 2*j+1 ;
break ;
case T_LEG_PI :
l_quant = 2*j+1 ;
break ;
case T_LEG_II :
l_quant = 2*j ;
break ;
case T_LEG :
l_quant = j ;
break ;
case T_LEG_MP :
l_quant = j ;
break ;
case T_LEG_MI :
l_quant = j ;
break ;
case T_CL_COS_P:
l_quant = 2*j ;
break ;
case T_CL_SIN_P:
l_quant = 2*j ;
break ;
case T_CL_COS_I:
l_quant = 2*j+1 ;
break ;
case T_CL_SIN_I:
l_quant = 2*j+1 ;
break ;
default:
cout << "Unknown basis in theta in give_quant_numbers ..." << endl ;
abort() ;
break ;
}
switch (base_r) {
case R_CHEB :
base_r_1d = R_CHEB ;
break ;
case R_CHEBP :
base_r_1d = R_CHEBP ;
break ;
case R_CHEBI :
base_r_1d = R_CHEBI ;
break ;
case R_LEG :
base_r_1d = R_LEG ;
break ;
case R_LEGP :
base_r_1d = R_LEGP ;
break ;
case R_LEGI :
base_r_1d = R_LEGI ;
break ;
case R_JACO02 :
base_r_1d = R_JACO02 ;
break ;
case R_CHEBPIM_P :
if (m_quant%2 == 0)
base_r_1d = R_CHEBP ;
else
base_r_1d = R_CHEBI ;
break ;
case R_CHEBPIM_I :
if (m_quant%2 == 0)
base_r_1d = R_CHEBI ;
else
base_r_1d = R_CHEBP ;
break ;
case R_CHEBPI_P :
if (l_quant%2 == 0)
base_r_1d = R_CHEBP ;
else
base_r_1d = R_CHEBI ;
break ;
case R_CHEBPI_I :
if (l_quant%2 == 0)
base_r_1d = R_CHEBI ;
else
base_r_1d = R_CHEBP ;
break ;
case R_CHEBU :
base_r_1d = R_CHEBU ;
break ;
default:
cout << "Unknown basis in r in give_quant_numbers ..." << endl ;
abort() ;
break ;
}
}
int Base_val::give_lmax(const Mg3d& mgrid, int lz) const {
#ifndef NDEBUG
int nz = mgrid.get_nzone() ;
assert (lz < nz) ;
#endif
int ntm1 = mgrid.get_nt(lz) - 1;
int base_t = get_base_t(lz) ;
bool m_odd = (mgrid.get_np(lz) > 2) ;
int l_max = 0 ;
switch (base_t) {
case T_COS_P :
l_max = 2*ntm1 ;
break;
case T_SIN_P :
l_max = 2*ntm1 ;
break;
case T_COS_I :
l_max = 2*ntm1+1 ;
break;
case T_SIN_I :
l_max = 2*ntm1+1 ;
break;
case T_COSSIN_CP :
if (!m_odd)
l_max = 2*ntm1 ;
else
l_max = 2*ntm1+1 ;
break ;
case T_COSSIN_SP :
if (!m_odd)
l_max = 2*ntm1 ;
else
l_max = 2*ntm1+1 ;
break ;
case T_COSSIN_CI :
if (!m_odd)
l_max = 2*ntm1+1 ;
else
l_max = 2*ntm1 ;
break ;
case T_COSSIN_SI :
if (!m_odd)
l_max = 2*ntm1+1 ;
else
l_max = 2*ntm1 ;
break ;
case T_COSSIN_C :
l_max = ntm1 ;
break ;
case T_COSSIN_S :
l_max = ntm1 ;
break ;
case T_COS :
l_max = ntm1 ;
break ;
case T_LEG_P :
if (!m_odd)
l_max = 2*ntm1 ;
else
l_max = 2*ntm1+1 ;
break ;
case T_LEG_PP :
l_max = 2*ntm1 ;
break ;
case T_LEG_I :
if (!m_odd)
l_max = 2*ntm1+1 ;
else
l_max = 2*ntm1 ;
break ;
case T_LEG_IP :
l_max = 2*ntm1+1 ;
break ;
case T_LEG_PI :
l_max = 2*ntm1+1 ;
break ;
case T_LEG_II :
l_max = 2*ntm1 ;
break ;
case T_LEG :
l_max = ntm1 ;
break ;
case T_LEG_MP :
l_max = ntm1 ;
break ;
case T_LEG_MI :
l_max = ntm1 ;
break ;
case T_CL_COS_P:
l_max = 2*ntm1 ;
break ;
case T_CL_SIN_P:
l_max = 2*ntm1 ;
break ;
case T_CL_COS_I:
l_max = 2*ntm1+1 ;
break ;
case T_CL_SIN_I:
l_max = 2*ntm1+1 ;
break ;
default:
cout << "Unknown basis in theta in Base_val::get_lmax ..."
<< endl ;
abort() ;
break ;
}
return l_max ;
}
}
|