base_val_quantum.C

00001 /*
00002  *   Copyright (c) 2004 Philippe Grandclement
00003  *
00004  *   This file is part of LORENE.
00005  *
00006  *   LORENE is free software; you can redistribute it and/or modify
00007  *   it under the terms of the GNU General Public License as published by
00008  *   the Free Software Foundation; either version 2 of the License, or
00009  *   (at your option) any later version.
00010  *
00011  *   LORENE is distributed in the hope that it will be useful,
00012  *   but WITHOUT ANY WARRANTY; without even the implied warranty of
00013  *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00014  *   GNU General Public License for more details.
00015  *
00016  *   You should have received a copy of the GNU General Public License
00017  *   along with LORENE; if not, write to the Free Software
00018  *   Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00019  *
00020  */
00021 
00022 
00023 char base_val_quantum_C[] = "$Header: /cvsroot/Lorene/C++/Source/Base_val/base_val_quantum.C,v 1.7 2009/10/23 12:55:16 j_novak Exp $" ;
00024 
00025 /*
00026  * $Id: base_val_quantum.C,v 1.7 2009/10/23 12:55:16 j_novak Exp $
00027  * $Log: base_val_quantum.C,v $
00028  * Revision 1.7  2009/10/23 12:55:16  j_novak
00029  * New base T_LEG_MI
00030  *
00031  * Revision 1.6  2009/10/08 16:20:13  j_novak
00032  * Addition of new bases T_COS and T_SIN.
00033  *
00034  * Revision 1.5  2007/12/11 15:28:09  jl_cornou
00035  * Jacobi(0,2) polynomials partially implemented
00036  *
00037  * Revision 1.4  2005/09/07 13:09:50  j_novak
00038  * New method for determining the highest multipole that can be described on a 3D
00039  * grid.
00040  *
00041  * Revision 1.3  2004/11/23 15:08:01  m_forot
00042  * Added the bases for the cases without any equatorial symmetry
00043  * (T_COSSIN_C, T_COSSIN_S, T_LEG, R_CHEBPI_P, R_CHEBPI_I).
00044  *
00045  * Revision 1.2  2004/08/30 16:27:59  r_prix
00046  * added #include <stdlib.h> (got ERROR 'abort' is undefined without this...)
00047  *
00048  * Revision 1.1  2004/08/24 09:14:41  p_grandclement
00049  * Addition of some new operators, like Poisson in 2d... It now requieres the
00050  * GSL library to work.
00051  *
00052  * Also, the way a variable change is stored by a Param_elliptic is changed and
00053  * no longer uses Change_var but rather 2 Scalars. The codes using that feature
00054  * will requiere some modification. (It should concern only the ones about monopoles)
00055  *
00056  *
00057  * $Header: /cvsroot/Lorene/C++/Source/Base_val/base_val_quantum.C,v 1.7 2009/10/23 12:55:16 j_novak Exp $
00058  *
00059  */
00060 
00061 // Headers C
00062 #include <stdio.h>
00063 #include <assert.h>
00064 #include <stdlib.h>
00065 
00066 // Headers Lorene
00067 #include "grilles.h"
00068 #include "base_val.h"
00069 #include "utilitaires.h"
00070 
00071 void Base_val::give_quant_numbers (int l, int k, int j, 
00072              int& m_quant, int& l_quant, int& base_r_1d) const {
00073 
00074   int base_p = get_base_p(l) ;
00075   int base_t = get_base_t(l) ;
00076   int base_r = get_base_r(l) ;
00077   
00078   switch (base_p) {
00079   case P_COSSIN :
00080       m_quant = k/2 ;
00081     break;
00082     
00083   case P_COSSIN_P :
00084     if (k%2 == 0)
00085       m_quant = k ;
00086     else
00087       m_quant = (k-1) ;
00088     break;
00089 
00090   case P_COSSIN_I :
00091       m_quant = 2*( (k-1) / 2) + 1 ;
00092     break;
00093   default:
00094     cout << "Unknown basis in phi in give_quant_numbers ..." << endl ;
00095     abort() ;
00096     break ;
00097   }
00098 
00099   switch (base_t) {
00100   case T_COS_P :
00101     l_quant = 2*j ;
00102     break;
00103     
00104   case T_SIN_P :
00105     l_quant = 2*j ;
00106     break;
00107 
00108   case T_COS_I :
00109     l_quant = 2*j+1 ;
00110     break;
00111 
00112   case T_SIN_I :
00113     l_quant = 2*j+1 ;
00114     break;
00115  
00116   case T_COSSIN_CP :
00117     if (m_quant%2 == 0)
00118       l_quant = 2*j ;
00119     else
00120       l_quant = 2*j+1 ;
00121     break ;
00122     
00123   case T_COSSIN_SP :
00124     if (m_quant%2 == 0)
00125       l_quant = 2*j ;
00126     else
00127       l_quant = 2*j+1 ;
00128     break ;
00129     
00130   case T_COSSIN_CI :
00131     if (m_quant%2 == 0)
00132       l_quant = 2*j+1 ;
00133     else
00134       l_quant = 2*j ;
00135     break ;
00136     
00137   case T_COSSIN_SI :
00138     if (m_quant%2 == 0)
00139       l_quant = 2*j+1 ;
00140     else
00141       l_quant = 2*j ;
00142     break ;
00143 
00144   case T_COSSIN_C :
00145        l_quant = j ;
00146     break ;
00147     
00148   case T_COSSIN_S :
00149        l_quant = j ;
00150     break ;   
00151 
00152   case T_COS :
00153        l_quant = j ;
00154     break ;
00155     
00156   case T_LEG_P :
00157     if (m_quant%2 == 0)
00158       l_quant = 2*j ;
00159     else
00160       l_quant = 2*j+1 ;
00161     break ;
00162 
00163  case T_LEG_PP :
00164    l_quant = 2*j ;
00165    break ;
00166 
00167   case T_LEG_I :
00168     if (m_quant%2 == 0)
00169       l_quant = 2*j+1 ;
00170     else
00171       l_quant = 2*j ;
00172     break ;
00173    
00174   case T_LEG_IP :
00175     l_quant = 2*j+1 ;
00176     break ;
00177 
00178   case T_LEG_PI :
00179     l_quant = 2*j+1 ;
00180     break ; 
00181 
00182   case T_LEG_II :
00183     l_quant = 2*j ;
00184     break ; 
00185 
00186   case T_LEG :
00187    l_quant = j ;
00188    break ;
00189 
00190   case T_LEG_MP :
00191    l_quant = j ;
00192    break ;
00193 
00194   case T_LEG_MI :
00195    l_quant = j ;
00196    break ;
00197 
00198   case T_CL_COS_P:
00199     l_quant = 2*j ;
00200     break ;
00201   
00202   case T_CL_SIN_P:
00203     l_quant = 2*j ;
00204     break ;
00205    
00206   case T_CL_COS_I:
00207     l_quant = 2*j+1 ;
00208     break ;
00209     
00210   case T_CL_SIN_I:
00211     l_quant = 2*j+1 ;
00212     break ;
00213     
00214   default:
00215     cout << "Unknown basis in theta in give_quant_numbers ..." << endl ;
00216     abort() ;
00217     break ;
00218   }
00219 
00220   switch (base_r) {
00221   case R_CHEB :
00222     base_r_1d = R_CHEB ;
00223     break ;
00224 
00225   case R_JACO02 :
00226     base_r_1d = R_JACO02 ;
00227     break ;
00228     
00229   case R_CHEBP :
00230     base_r_1d = R_CHEBP ;
00231     break ;
00232     
00233    case R_CHEBI :
00234     base_r_1d = R_CHEBI ;
00235     break ;
00236     
00237   case R_CHEBPIM_P :    
00238     if (m_quant%2 == 0) 
00239       base_r_1d = R_CHEBP ;
00240     else
00241       base_r_1d = R_CHEBI ;
00242     break ;
00243   case R_CHEBPIM_I :
00244     if (m_quant%2 == 0)
00245       base_r_1d = R_CHEBI ;
00246     else
00247       base_r_1d = R_CHEBP ;
00248     break ;
00249   case R_CHEBPI_P :    
00250     if (l_quant%2 == 0) 
00251       base_r_1d = R_CHEBP ;
00252     else
00253       base_r_1d = R_CHEBI ;
00254     break ;
00255   case R_CHEBPI_I :
00256     if (l_quant%2 == 0)
00257       base_r_1d = R_CHEBI ;
00258     else
00259       base_r_1d = R_CHEBP ;
00260     break ;   
00261   case R_CHEBU :
00262     base_r_1d = R_CHEBU ;
00263     break ;
00264     
00265   default:
00266     cout << "Unknown basis in r in give_quant_numbers ..." << endl ;
00267     abort() ;
00268     break ;
00269   }
00270 }
00271 
00272 int Base_val::give_lmax(const Mg3d& mgrid, int lz) const {
00273 
00274 #ifndef NDEBUG
00275     int nz = mgrid.get_nzone() ;
00276     assert (lz < nz) ;
00277 #endif
00278 
00279     int ntm1 = mgrid.get_nt(lz) - 1;
00280     int base_t = get_base_t(lz) ;
00281     bool m_odd = (mgrid.get_np(lz) > 2) ;
00282 
00283     int l_max = 0 ;
00284 
00285     switch (base_t) {
00286     case T_COS_P :
00287         l_max = 2*ntm1 ;
00288         break;
00289     
00290     case T_SIN_P :
00291         l_max = 2*ntm1 ;
00292         break;
00293 
00294     case T_COS_I :
00295         l_max = 2*ntm1+1 ;
00296         break;
00297 
00298     case T_SIN_I :
00299         l_max = 2*ntm1+1 ;
00300         break;
00301  
00302     case T_COSSIN_CP :
00303         if (!m_odd)
00304         l_max = 2*ntm1 ;
00305         else
00306         l_max = 2*ntm1+1 ;
00307         break ;
00308     
00309     case T_COSSIN_SP :
00310         if (!m_odd)
00311         l_max = 2*ntm1 ;
00312         else
00313         l_max = 2*ntm1+1 ;
00314         break ;
00315         
00316     case T_COSSIN_CI :
00317         if (!m_odd)
00318         l_max = 2*ntm1+1 ;
00319         else
00320         l_max = 2*ntm1 ;
00321         break ;
00322         
00323     case T_COSSIN_SI :
00324         if (!m_odd)
00325         l_max = 2*ntm1+1 ;
00326         else
00327         l_max = 2*ntm1 ;
00328         break ;
00329         
00330     case T_COSSIN_C :
00331         l_max = ntm1 ;
00332         break ;
00333         
00334     case T_COSSIN_S :
00335         l_max = ntm1 ;
00336         break ;   
00337         
00338     case T_COS :
00339         l_max = ntm1 ;
00340         break ;   
00341         
00342     case T_LEG_P :
00343         if (!m_odd)
00344         l_max = 2*ntm1 ;
00345         else
00346         l_max = 2*ntm1+1 ;
00347         break ;
00348         
00349     case T_LEG_PP :
00350         l_max = 2*ntm1 ;
00351         break ;
00352         
00353     case T_LEG_I :
00354         if (!m_odd)
00355         l_max = 2*ntm1+1 ;
00356         else
00357         l_max = 2*ntm1 ;
00358         break ;
00359         
00360     case T_LEG_IP :
00361         l_max = 2*ntm1+1 ;
00362         break ;
00363         
00364     case T_LEG_PI :
00365         l_max = 2*ntm1+1 ;
00366         break ; 
00367         
00368     case T_LEG_II :
00369         l_max = 2*ntm1 ;
00370         break ; 
00371         
00372     case T_LEG :
00373         l_max = ntm1 ;
00374         break ;
00375         
00376     case T_LEG_MP :
00377         l_max = ntm1 ;
00378         break ;
00379         
00380     case T_LEG_MI :
00381         l_max = ntm1 ;
00382         break ;
00383         
00384     case T_CL_COS_P:
00385         l_max = 2*ntm1 ;
00386         break ;
00387         
00388     case T_CL_SIN_P:
00389         l_max = 2*ntm1 ;
00390         break ;
00391         
00392     case T_CL_COS_I:
00393         l_max = 2*ntm1+1 ;
00394         break ;
00395         
00396     case T_CL_SIN_I:
00397         l_max = 2*ntm1+1 ;
00398         break ;
00399         
00400     default:
00401         cout << "Unknown basis in theta in Base_val::get_lmax ..." 
00402          << endl ;
00403         abort() ;
00404         break ;
00405     }
00406     return l_max ;
00407 }

Generated on Tue Feb 7 01:35:14 2012 for LORENE by  doxygen 1.4.6