scalar_manip.C

00001 /*
00002  *   Copyright (c) 2003 Eric Gourgoulhon & Jerome Novak
00003  *   
00004  *   Copyright (c) 2000-2001 Philippe Grandclement (Cmp version)
00005  *
00006  *   This file is part of LORENE.
00007  *
00008  *   LORENE is free software; you can redistribute it and/or modify
00009  *   it under the terms of the GNU General Public License as published by
00010  *   the Free Software Foundation; either version 2 of the License, or
00011  *   (at your option) any later version.
00012  *
00013  *   LORENE is distributed in the hope that it will be useful,
00014  *   but WITHOUT ANY WARRANTY; without even the implied warranty of
00015  *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00016  *   GNU General Public License for more details.
00017  *
00018  *   You should have received a copy of the GNU General Public License
00019  *   along with LORENE; if not, write to the Free Software
00020  *   Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00021  *
00022  */
00023 
00024 
00025 char scalar_manip_C[] = "$Header: /cvsroot/Lorene/C++/Source/Tensor/Scalar/scalar_manip.C,v 1.17 2008/10/03 09:03:52 j_novak Exp $" ;
00026 
00027 /*
00028  * $Id: scalar_manip.C,v 1.17 2008/10/03 09:03:52 j_novak Exp $
00029  * $Log: scalar_manip.C,v $
00030  * Revision 1.17  2008/10/03 09:03:52  j_novak
00031  * Correction of yet another mistake (the array values in physical space was not
00032  * destroyed).
00033  *
00034  * Revision 1.16  2008/09/30 08:35:18  j_novak
00035  * Correction of forgotten call to coef()
00036  *
00037  * Revision 1.15  2008/09/29 13:23:51  j_novak
00038  * Implementation of the angular mapping associated with an affine
00039  * mapping. Things must be improved to take into account the domain index.
00040  *
00041  * Revision 1.14  2008/09/22 19:08:01  j_novak
00042  * New methods to deal with boundary conditions
00043  *
00044  * Revision 1.13  2007/06/21 20:00:00  k_taniguchi
00045  * Addition of the method filtre_r (int n, int nz)
00046  *
00047  * Revision 1.12  2006/06/28 07:46:41  j_novak
00048  * Better treatment in the case of a domain set to zero.
00049  *
00050  * Revision 1.11  2005/09/07 13:39:10  j_novak
00051  * *** empty log message ***
00052  *
00053  * Revision 1.10  2005/09/07 13:10:48  j_novak
00054  * Added a filter setting to zero all mulitpoles in a given range.
00055  *
00056  * Revision 1.9  2004/11/23 12:47:44  f_limousin
00057  * Add function filtre_tp(int nn, int nz1, int nz2).
00058  *
00059  * Revision 1.8  2004/05/07 11:24:54  f_limousin
00060  * Implement new method filtre_r (int* nn)
00061  *
00062  * Revision 1.7  2004/02/27 09:47:26  f_limousin
00063  * New methods filtre_phi(int) and filtre_theta(int).
00064  *
00065  * Revision 1.6  2004/01/23 13:26:28  e_gourgoulhon
00066  *  Added methods set_inner_boundary and set_outer_boundary.
00067  *  Methods set_val_inf and set_val_hor, which are particular cases of
00068  *  the above, have been suppressed.
00069  *
00070  * Revision 1.5  2003/11/13 13:43:55  p_grandclement
00071  * Addition of things needed for Bhole::update_metric (const Etoile_bin&, double, double)
00072  *
00073  * Revision 1.4  2003/10/11 14:47:17  e_gourgoulhon
00074  * Lines 112 and 145 : replaced "0" by "double(0)" in the comparison
00075  * statement.
00076  *
00077  * Revision 1.3  2003/10/10 15:57:29  j_novak
00078  * Added the state one (ETATUN) to the class Scalar
00079  *
00080  * Revision 1.2  2003/10/08 14:24:09  j_novak
00081  * replaced mult_r_zec with mult_r_ced
00082  *
00083  * Revision 1.1  2003/09/25 09:33:36  j_novak
00084  * Added methods for integral calculation and various manipulations
00085  *
00086  *
00087  * $Header: /cvsroot/Lorene/C++/Source/Tensor/Scalar/scalar_manip.C,v 1.17 2008/10/03 09:03:52 j_novak Exp $
00088  *
00089  */
00090 
00091 //standard
00092 #include <stdlib.h>
00093 #include <math.h>
00094 
00095 // Lorene
00096 #include "tensor.h"
00097 #include "proto.h"
00098 #include "utilitaires.h"
00099 
00100 /*
00101  * Annule tous les l entre l_min et l_max (compris)
00102  */
00103 
00104 void Scalar::annule_l (int l_min, int l_max, bool ylm_output) {
00105 
00106     assert (etat != ETATNONDEF) ;
00107     assert (l_min <= l_max) ;
00108     assert (l_min >= 0) ;
00109     if (etat == ETATZERO )
00110     return ;
00111 
00112     if (etat == ETATUN) {
00113     if (l_min == 0) set_etat_zero() ;
00114     else return ;
00115     }
00116 
00117     va.ylm() ;
00118     Mtbl_cf& m_coef = *va.c_cf ;
00119     const Base_val& base = va.base ; 
00120     int l_q, m_q, base_r ;
00121     for (int lz=0; lz<mp->get_mg()->get_nzone(); lz++)
00122     if (m_coef(lz).get_etat() != ETATZERO)
00123         for (int k=0; k<mp->get_mg()->get_np(lz)+1; k++)
00124         for (int j=0; j<mp->get_mg()->get_nt(lz); j++)
00125             for (int i=0; i<mp->get_mg()->get_nr(lz); i++) {
00126             base.give_quant_numbers(lz, k, j, m_q, l_q, base_r) ;
00127             if ((l_min <= l_q) && (l_q<= l_max)) 
00128                 m_coef.set(lz, k, j, i) = 0 ;
00129             }
00130     if (va.c != 0x0) {
00131     delete va.c ;
00132     va.c = 0x0 ;
00133     }
00134     if (!ylm_output)
00135     va.ylm_i() ;
00136     
00137     return ;
00138 }
00139 
00140 /*
00141  * Annule les n derniers coefficients en r dans la derniere zone
00142  */
00143  
00144 void Scalar::filtre (int n) {
00145     
00146     assert (etat != ETATNONDEF) ;
00147     if ( (etat == ETATZERO) || (etat == ETATUN) )
00148     return ;
00149     
00150     int nz = mp->get_mg()->get_nzone() ;
00151     int np = mp->get_mg()->get_np(nz-1) ;
00152     int nt = mp->get_mg()->get_nt(nz-1) ;
00153     int nr = mp->get_mg()->get_nr(nz-1) ;
00154     
00155     del_deriv() ;
00156     
00157     va.coef() ;
00158     va.set_etat_cf_qcq() ;
00159 
00160     
00161     for (int k=0 ; k<np+1 ; k++)
00162     if (k!=1)
00163         for (int j=0 ; j<nt ; j++)
00164         for (int i=nr-1 ; i>nr-1-n ; i--)
00165             va.c_cf->set(nz-1, k, j, i) = 0 ;
00166 }
00167 
00168 
00169 /*
00170  * Annule les n derniers coefficients en r dans toutes les zones
00171  */
00172  
00173 void Scalar::filtre_r (int* nn) {
00174     assert (etat != ETATNONDEF) ;
00175     if ( (etat == ETATZERO) || (etat == ETATUN) )
00176     return ;
00177     
00178     del_deriv() ;
00179     
00180     va.coef() ;
00181     va.set_etat_cf_qcq() ;
00182     int nz = mp->get_mg()->get_nzone() ;
00183     int* nr = new int[nz];
00184     int* nt = new int[nz];
00185     int* np = new int[nz];
00186     for (int l=0; l<=nz-1; l++) {
00187     nr[l] = mp->get_mg()->get_nr(l) ; 
00188     nt[l] = mp->get_mg()->get_nt(l) ; 
00189     np[l] = mp->get_mg()->get_np(l) ; 
00190     }
00191  
00192     for (int l=0; l<=nz-1; l++)
00193     for (int k=0 ; k<np[l]+1 ; k++)
00194         if (k!=1)
00195         for (int j=0 ; j<nt[l] ; j++)
00196             for (int i=nr[l]-1; i>nr[l]-1-nn[l] ; i--)
00197             va.c_cf->set(l, k, j, i) = 0 ;
00198 
00199     if (va.c != 0x0) {
00200         delete va.c ;
00201         va.c = 0x0 ;
00202     }
00203     
00204 }
00205 
00206 
00207 /*
00208  * Annule les n derniers coefficients en r dans zone nz
00209  */
00210 
00211 void Scalar::filtre_r (int n, int nz) {
00212     assert (etat != ETATNONDEF) ;
00213     if ( (etat == ETATZERO) || (etat == ETATUN) )
00214     return ;
00215     
00216     del_deriv() ;
00217     
00218     va.coef() ;
00219     va.set_etat_cf_qcq() ;
00220     int nr = mp->get_mg()->get_nr(nz) ; 
00221     int nt = mp->get_mg()->get_nt(nz) ; 
00222     int np = mp->get_mg()->get_np(nz) ; 
00223 
00224     for (int k=0 ; k<np+1 ; k++)
00225         if (k!=1)
00226         for (int j=0 ; j<nt ; j++)
00227             for (int i=nr-1; i>nr-1-n ; i--)
00228             va.c_cf->set(nz, k, j, i) = 0 ;
00229 
00230     if (va.c != 0x0) {
00231         delete va.c ;
00232     va.c = 0x0 ;
00233     }
00234 
00235 }
00236 
00237 
00238 /*
00239  * Annule les n derniers coefficients en phi dans zone nz
00240  */
00241  
00242 void Scalar::filtre_phi (int n, int nz) {
00243     assert (etat != ETATNONDEF) ;
00244     if ( (etat == ETATZERO) || (etat == ETATUN) )
00245     return ;
00246     
00247     del_deriv() ;
00248     
00249     va.coef() ;
00250     va.set_etat_cf_qcq() ;
00251     int np = mp->get_mg()->get_np(nz) ;
00252     int nt = mp->get_mg()->get_nt(nz) ;
00253     int nr = mp->get_mg()->get_nr(nz) ;
00254     
00255     for (int k=np+1-n ; k<np+1 ; k++)
00256     for (int j=0 ; j<nt ; j++)
00257         for (int i=0 ; i<nr ; i++)
00258         va.c_cf->set(nz, k, j, i) = 0 ;
00259 
00260 }
00261 
00262 
00263 void Scalar::filtre_tp(int nn, int nz1, int nz2) {
00264 
00265     va.filtre_tp(nn, nz1, nz2) ;
00266 
00267 }
00268 
00269 
00270   /* Sets the value of the {\tt Scalar} at the inner boundary of a given 
00271    * domain. 
00272    * @param l [input] domain index
00273    * @param x [input] (constant) value at the inner boundary of domain no. {\tt l}
00274    */
00275 
00276 void Scalar::set_inner_boundary(int l0, double x0) {
00277     
00278     assert (etat != ETATNONDEF) ;
00279     if (etat == ETATZERO) {
00280         if (x0 == double(0)) return ;
00281         else annule_hard() ;
00282     }
00283     
00284     if (etat == ETATUN) {
00285         if (x0 == double(1)) return ;
00286         else etat = ETATQCQ ;
00287     }
00288     
00289     del_deriv() ;
00290     
00291     int nt = mp->get_mg()->get_nt(l0) ;
00292     int np = mp->get_mg()->get_np(l0) ;
00293     
00294     va.coef_i() ;
00295     va.set_etat_c_qcq() ;
00296     
00297     for (int k=0 ; k<np ; k++)
00298     for (int j=0 ; j<nt ; j++)
00299         va.set(l0, k, j, 0) = x0 ;
00300 }
00301 
00302   /* Sets the value of the {\tt Scalar} at the outer boundary of a given 
00303    * domain. 
00304    * @param l [input] domain index
00305    * @param x [input] (constant) value at the outer boundary of domain no. {\tt l}
00306    */
00307 
00308 void Scalar::set_outer_boundary(int l0, double x0) {
00309     
00310     assert (etat != ETATNONDEF) ;
00311     if (etat == ETATZERO) {
00312         if (x0 == double(0)) return ;
00313         else annule_hard() ;
00314     }
00315     
00316     if (etat == ETATUN) {
00317         if (x0 == double(1)) return ;
00318         else etat = ETATQCQ ;
00319     }
00320     
00321     del_deriv() ;
00322     
00323     int nrm1 = mp->get_mg()->get_nr(l0) - 1 ;
00324     int nt = mp->get_mg()->get_nt(l0) ;
00325     int np = mp->get_mg()->get_np(l0) ;
00326     
00327     va.coef_i() ;
00328     va.set_etat_c_qcq() ;
00329     
00330     for (int k=0 ; k<np ; k++)
00331     for (int j=0 ; j<nt ; j++)
00332         va.set(l0, k, j, nrm1) = x0 ;
00333 }
00334 
00335 /*
00336  * Permet de fixer la decroissance du cmp a l infini en viurant les 
00337  * termes en 1/r^n
00338  */
00339 void Scalar::fixe_decroissance (int puis) {
00340     
00341     if (puis<dzpuis)
00342     return ;
00343     else {
00344     
00345     int nbre = puis-dzpuis ;
00346     
00347     // le confort avant tout ! (c'est bien le confort ...)
00348     int nz = mp->get_mg()->get_nzone() ;
00349     int np = mp->get_mg()->get_np(nz-1) ;
00350     int nt = mp->get_mg()->get_nt(nz-1) ;
00351     int nr = mp->get_mg()->get_nr(nz-1) ;
00352     
00353     const Map_af* map  = dynamic_cast<const Map_af*>(mp) ;
00354     if (map == 0x0) {
00355         cout << "Le mapping doit etre affine" << endl ;
00356         abort() ;
00357     }
00358     
00359     double alpha = map->get_alpha()[nz-1] ;
00360     
00361     Scalar courant (*this) ;
00362     
00363     va.coef() ;
00364     va.set_etat_cf_qcq() ;
00365     
00366     for (int conte=0 ; conte<nbre ; conte++) {
00367         
00368         int base_r = courant.va.base.get_base_r(nz-1) ;
00369         
00370         courant.va.coef() ;
00371         
00372         // On calcul les coefficients de 1/r^conte
00373         double* coloc = new double [nr] ;
00374         int * deg = new int[3] ;
00375         deg[0] = 1 ; 
00376         deg[1] = 1 ;
00377         deg[2] = nr ;
00378             
00379         for (int i=0 ; i<nr ; i++)
00380         coloc[i] =pow(alpha, double(conte))*
00381             pow(-1-cos(M_PI*i/(nr-1)), double(conte)) ;
00382             
00383         cfrcheb(deg, deg, coloc, deg, coloc) ;
00384         
00385         for (int k=0 ; k<np+1 ; k++)
00386         if (k != 1)
00387         for (int j=0 ; j<nt ; j++) {
00388             
00389             // On doit determiner le coefficient du truc courant :
00390             double* coef = new double [nr] ;
00391             double* auxi = new double[1] ;
00392             for (int i=0 ; i<nr ; i++)
00393             coef[i] = (*courant.va.c_cf)(nz-1, k, j, i) ;
00394             switch (base_r) {
00395             case R_CHEBU :
00396             som_r_chebu (coef, nr, 1, 1, 1, auxi) ;
00397             break ;
00398             default :
00399             som_r_pas_prevu (coef, nr, 1, 1, 1, auxi) ;
00400             break ;
00401             }
00402             
00403             // On modifie le cmp courant :
00404             courant.va.coef() ;
00405             courant.va.set_etat_cf_qcq() ;
00406             courant.va.c_cf->set(nz-1, k, j, 0) -= *auxi ;  
00407             
00408             for (int i=0 ; i<nr ; i++)
00409                 this->va.c_cf->set(nz-1, k, j, i) -= *auxi * coloc[i] ;
00410 
00411               
00412             delete [] coef ;
00413             delete [] auxi ;
00414         }
00415         delete [] coloc ;
00416         delete [] deg ;
00417         
00418         courant.mult_r_ced() ;
00419     }
00420     }
00421 }
00422 
00423 Tbl Scalar::tbl_out_bound(int l_zone, bool output_ylm) {
00424     va.coef() ;
00425     if (output_ylm) va.ylm() ;
00426 
00427     int np = mp->get_mg()->get_np(l_zone) ;
00428     int nt = mp->get_mg()->get_nt(l_zone) ;
00429     Tbl resu(np+2, nt) ;
00430     if (etat == ETATZERO) resu.set_etat_zero() ;
00431     else {
00432     assert(etat == ETATQCQ) ;
00433     resu.set_etat_qcq() ;
00434     for (int k=0; k<np+2; k++)
00435         for (int j=0; j<nt; j++)
00436         resu.set(k, j) = va.c_cf->val_out_bound_jk(l_zone, j, k) ;
00437     }
00438     return resu ;
00439 }
00440 
00441 Tbl Scalar::tbl_in_bound(int l_zone, bool output_ylm) {
00442     assert(mp->get_mg()->get_type_r(l_zone) != RARE) ;
00443     va.coef() ;
00444     if (output_ylm) va.ylm() ;
00445 
00446     int np = mp->get_mg()->get_np(l_zone) ;
00447     int nt = mp->get_mg()->get_nt(l_zone) ;
00448     Tbl resu(np+2, nt) ;
00449     if (etat == ETATZERO) resu.set_etat_zero() ;
00450     else {
00451     assert(etat == ETATQCQ) ;
00452     resu.set_etat_qcq() ;
00453     for (int k=0; k<np+2; k++)
00454         for (int j=0; j<nt; j++)
00455         resu.set(k, j) = va.c_cf->val_in_bound_jk(l_zone, j, k) ;
00456     }
00457     return resu ;
00458 }
00459 
00460 Scalar Scalar::scalar_out_bound(int l_zone, bool output_ylm) {
00461     va.coef() ;
00462     if (output_ylm) va.ylm() ;
00463 
00464     Scalar resu(mp->mp_angu(l_zone)) ;
00465     resu.std_spectral_base() ;
00466     Base_val base = resu.get_spectral_base() ;
00467     base.set_base_t(va.base.get_base_t(l_zone)) ;
00468     resu.set_spectral_base(base) ;
00469     if (etat == ETATZERO) resu.set_etat_zero() ;
00470     else {
00471     assert(etat == ETATQCQ) ;
00472     resu.annule_hard() ;
00473     int np = mp->get_mg()->get_np(l_zone) ;
00474     int nt = mp->get_mg()->get_nt(l_zone) ;
00475     for (int k=0; k<np+2; k++)
00476         for (int j=0; j<nt; j++)
00477         resu.set_spectral_va().c_cf->set(0, k, j, 0) 
00478             = va.c_cf->val_out_bound_jk(l_zone, j, k) ;
00479     delete resu.set_spectral_va().c ;
00480     resu.set_spectral_va().c = 0x0 ;
00481     }
00482     return resu ;
00483 }

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