vector.C

00001 /*
00002  *  Methods of class Vector
00003  *
00004  *   (see file vector.h for documentation)
00005  *
00006  */
00007 
00008 /*
00009  *   Copyright (c) 2003 Eric Gourgoulhon & Jerome Novak
00010  *
00011  *   This file is part of LORENE.
00012  *
00013  *   LORENE is free software; you can redistribute it and/or modify
00014  *   it under the terms of the GNU General Public License as published by
00015  *   the Free Software Foundation; either version 2 of the License, or
00016  *   (at your option) any later version.
00017  *
00018  *   LORENE is distributed in the hope that it will be useful,
00019  *   but WITHOUT ANY WARRANTY; without even the implied warranty of
00020  *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00021  *   GNU General Public License for more details.
00022  *
00023  *   You should have received a copy of the GNU General Public License
00024  *   along with LORENE; if not, write to the Free Software
00025  *   Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00026  *
00027  */
00028 
00029 
00030 char vector_C[] = "$Header: /cvsroot/Lorene/C++/Source/Tensor/vector.C,v 1.28 2008/10/29 14:09:14 jl_cornou Exp $" ;
00031 
00032 /*
00033  * $Id: vector.C,v 1.28 2008/10/29 14:09:14 jl_cornou Exp $
00034  * $Log: vector.C,v $
00035  * Revision 1.28  2008/10/29 14:09:14  jl_cornou
00036  * Spectral bases for pseudo vectors and curl added
00037  *
00038  * Revision 1.27  2008/08/27 08:52:23  jl_cornou
00039  * Added fonctions for angular potential A
00040  *
00041  * Revision 1.26  2007/12/21 16:07:08  j_novak
00042  * Methods to filter Tensor, Vector and Sym_tensor objects.
00043  *
00044  * Revision 1.25  2005/02/14 13:01:50  j_novak
00045  * p_eta and p_mu are members of the class Vector. Most of associated functions
00046  * have been moved from the class Vector_divfree to the class Vector.
00047  *
00048  * Revision 1.24  2005/01/25 15:37:35  j_novak
00049  * Solved some dzpuis problem...
00050  *
00051  * Revision 1.23  2005/01/12 16:48:23  j_novak
00052  * Better treatment of the case where all vector components are null in
00053  * decompose_div .
00054  *
00055  * Revision 1.22  2004/10/12 09:58:25  j_novak
00056  * Better memory management.
00057  *
00058  * Revision 1.21  2004/10/11 09:46:31  j_novak
00059  * Speed improvements.
00060  *
00061  * Revision 1.20  2004/05/09 20:55:05  e_gourgoulhon
00062  * Added method flux.
00063  *
00064  * Revision 1.19  2004/03/29 11:57:45  e_gourgoulhon
00065  * Added methods ope_killing and ope_killing_conf.
00066  *
00067  * Revision 1.18  2004/02/26 22:48:50  e_gourgoulhon
00068  * -- Method divergence: call to Tensor::divergence and cast of the
00069  *    result.
00070  * -- Added method derive_lie.
00071  *
00072  * Revision 1.17  2004/02/24 09:46:20  j_novak
00073  * Correction to cope with SGI compiler's warnings.
00074  *
00075  * Revision 1.16  2004/02/20 10:53:04  j_novak
00076  * Added the matching of the potential adapted to the case where the
00077  * vector is the source of a Poisson equation (dzpuis = 4).
00078  *
00079  * Revision 1.15  2004/01/30 10:30:17  j_novak
00080  * Changed dzpuis handling in Vector::decompose_div (this may be temporary).
00081  *
00082  * Revision 1.14  2003/12/30 23:09:47  e_gourgoulhon
00083  * Change in methods derive_cov() and divergence() to take into account
00084  *  the change of name: Metric::get_connect() --> Metric::connect().
00085  *
00086  * Revision 1.13  2003/12/19 15:18:16  j_novak
00087  * Shadow variables hunt
00088  *
00089  * Revision 1.12  2003/10/29 11:04:34  e_gourgoulhon
00090  * dec2_dpzuis() replaced by dec_dzpuis(2).
00091  * inc2_dpzuis() replaced by inc_dzpuis(2).
00092  *
00093  * Revision 1.11  2003/10/22 14:24:19  j_novak
00094  * *** empty log message ***
00095  *
00096  * Revision 1.9  2003/10/20 13:00:38  j_novak
00097  * Memory error corrected
00098  *
00099  * Revision 1.8  2003/10/20 09:32:12  j_novak
00100  * Members p_potential and p_div_free of the Helmholtz decomposition
00101  * + the method decompose_div(Metric).
00102  *
00103  * Revision 1.7  2003/10/16 14:21:37  j_novak
00104  * The calculation of the divergence of a Tensor is now possible.
00105  *
00106  * Revision 1.6  2003/10/13 13:52:40  j_novak
00107  * Better managment of derived quantities.
00108  *
00109  * Revision 1.5  2003/10/06 13:58:48  j_novak
00110  * The memory management has been improved.
00111  * Implementation of the covariant derivative with respect to the exact Tensor
00112  * type.
00113  *
00114  * Revision 1.4  2003/10/05 21:14:20  e_gourgoulhon
00115  * Added method std_spectral_base().
00116  *
00117  * Revision 1.3  2003/10/03 14:10:32  e_gourgoulhon
00118  * Added constructor from Tensor.
00119  *
00120  * Revision 1.2  2003/10/03 14:08:46  j_novak
00121  * Removed old change_trid...
00122  *
00123  * Revision 1.1  2003/09/26 08:05:31  j_novak
00124  * New class Vector.
00125  *
00126  *
00127  * $Header: /cvsroot/Lorene/C++/Source/Tensor/vector.C,v 1.28 2008/10/29 14:09:14 jl_cornou Exp $
00128  *
00129  */
00130 
00131 // Headers C
00132 #include <stdlib.h>
00133 #include <assert.h>
00134 #include <math.h>
00135 
00136 // Headers Lorene
00137 #include "metric.h"
00138 #include "proto.h"
00139 #include "matrice.h"
00140 #include "nbr_spx.h"
00141 
00142 
00143             //--------------//
00144             // Constructors //
00145             //--------------//
00146 
00147 // Standard constructor 
00148 // --------------------
00149 Vector::Vector(const Map& map, int tipe, const Base_vect& triad_i) 
00150         : Tensor(map, 1, tipe, triad_i) {
00151         
00152   set_der_0x0() ;
00153 
00154 }
00155 
00156 // Standard constructor with the triad passed as a pointer
00157 // -------------------------------------------------------
00158 Vector::Vector(const Map& map, int tipe, const Base_vect* triad_i) 
00159         : Tensor(map, 1, tipe, *triad_i) {
00160         
00161   set_der_0x0() ;
00162 }
00163     
00164 // Copy constructor
00165 // ----------------
00166 Vector::Vector (const Vector& source) : 
00167     Tensor(source) {
00168   
00169   assert(valence == 1) ;
00170   set_der_0x0() ;
00171 
00172 }   
00173 
00174 
00175 // Constructor from a {\tt Tensor}.
00176 //--------------------------------
00177 Vector::Vector(const Tensor& uu) : Tensor(uu) {
00178 
00179   assert(valence == 1) ;
00180   set_der_0x0() ;
00181 
00182 }
00183 
00184 
00185 // Constructor from a file
00186 // -----------------------
00187 Vector::Vector(const Map& mapping, const Base_vect& triad_i, FILE* fd) : 
00188   Tensor(mapping, triad_i, fd) {
00189    
00190   assert ( (valence == 1) && (n_comp == 3) ) ;
00191   set_der_0x0() ;
00192 
00193 }
00194 
00195 
00196             //--------------//
00197             //  Destructor  //
00198             //--------------//
00199 
00200 
00201 Vector::~Vector () {
00202 
00203   Vector::del_deriv() ;
00204 
00205 }
00206 
00207 
00208             //-------------------//
00209             // Memory managment  //
00210             //-------------------//
00211 
00212 void Vector::del_deriv() const {
00213 
00214   for (int i=0; i<N_MET_MAX; i++) 
00215     del_derive_met(i) ;
00216   
00217   if (p_A != 0x0) delete p_A ;
00218   if (p_eta != 0x0) delete p_eta ; 
00219   if (p_mu != 0x0) delete p_mu ; 
00220   set_der_0x0() ;
00221   Tensor::del_deriv() ;
00222 
00223 }
00224 
00225 void Vector::set_der_0x0() const {
00226 
00227   for (int i=0; i<N_MET_MAX; i++) 
00228     set_der_met_0x0(i) ;
00229   p_A = 0x0 ;
00230   p_eta = 0x0 ; 
00231   p_mu = 0x0 ; 
00232 
00233 }
00234 
00235 void Vector::del_derive_met(int j) const {
00236 
00237   assert( (j>=0) && (j<N_MET_MAX) ) ;
00238 
00239   if (met_depend[j] != 0x0) {
00240     if (p_potential[j] != 0x0)
00241       delete p_potential[j] ;
00242     if (p_div_free[j] != 0x0)
00243       delete p_div_free[j] ;
00244     
00245     set_der_met_0x0(j) ;
00246     
00247     Tensor::del_derive_met(j) ;
00248   }
00249 }
00250 
00251 void Vector::set_der_met_0x0(int i) const {
00252 
00253   assert( (i>=0) && (i<N_MET_MAX) ) ;
00254 
00255   p_potential[i] = 0x0 ;
00256   p_div_free[i] = 0x0 ;
00257 
00258 }
00259 
00260 void Vector::operator=(const Vector& t) {
00261     
00262     triad = t.triad ; 
00263 
00264     assert(t.type_indice(0) == type_indice(0)) ;
00265 
00266     for (int i=0 ; i<3 ; i++) {
00267       *cmp[i] = *t.cmp[i] ;
00268     }
00269     del_deriv() ;
00270 }
00271 
00272 void Vector::operator=(const Tensor& t) {
00273     
00274     assert (t.valence == 1) ;
00275 
00276     triad = t.triad ; 
00277 
00278     assert(t.type_indice(0) == type_indice(0)) ;
00279 
00280     for (int i=0 ; i<3 ; i++) {
00281       *cmp[i] = *t.cmp[i] ;
00282     }
00283     del_deriv() ;
00284 }
00285 
00286 
00287 
00288 // Affectation d'une composante :
00289 Scalar& Vector::set(int index) {
00290     
00291   assert ( (index>=1) && (index<=3) ) ;
00292 
00293   del_deriv() ;
00294   
00295   return *cmp[index - 1] ;
00296 }
00297 
00298 const Scalar& Vector::operator()(int index) const {
00299     
00300   assert ( (index>=1) && (index<=3) ) ;
00301 
00302   return *cmp[index - 1] ;
00303 
00304 }
00305 
00306 
00307 // Sets the standard spectal bases of decomposition for each component
00308 
00309 void Vector::std_spectral_base() {
00310 
00311     Base_val** bases = 0x0 ;
00312 
00313     if ( triad->identify() == (mp->get_bvect_cart()).identify() ) {
00314 
00315         // Cartesian case
00316         bases = mp->get_mg()->std_base_vect_cart() ;
00317 
00318     }
00319     else {
00320         // Spherical case
00321         assert( triad->identify() == (mp->get_bvect_spher()).identify()) ;
00322         bases = mp->get_mg()->std_base_vect_spher() ;
00323     }
00324         
00325     for (int i=0 ; i<3 ; i++) {
00326         cmp[i]->set_spectral_base( *bases[i] ) ; 
00327     }
00328         
00329     for (int i=0 ; i<3 ; i++) {
00330         delete bases[i] ;
00331     }
00332     delete [] bases ;
00333 
00334 
00335 }
00336 
00337 // Sets the standard spectral bases of decomposition for each component for a pseudo vector
00338 
00339 void Vector::pseudo_spectral_base() {
00340 
00341     Base_val** bases = 0x0 ;
00342 
00343     if ( triad->identify() == (mp->get_bvect_cart()).identify() ) {
00344 
00345         // Cartesian case
00346         bases = mp->get_mg()->pseudo_base_vect_cart() ;
00347 
00348     }
00349     else {
00350         // Spherical case
00351         assert( triad->identify() == (mp->get_bvect_spher()).identify()) ;
00352         bases = mp->get_mg()->pseudo_base_vect_spher() ;
00353     }
00354         
00355     for (int i=0 ; i<3 ; i++) {
00356         cmp[i]->set_spectral_base( *bases[i] ) ; 
00357     }
00358         
00359     for (int i=0 ; i<3 ; i++) {
00360         delete bases[i] ;
00361     }
00362     delete [] bases ;
00363 
00364 
00365 }
00366 
00367 
00368 
00369                     //-------------------------------//
00370                     //    Computational methods      //
00371                     //-------------------------------//
00372                     
00373 
00374 const Scalar& Vector::divergence(const Metric& metre) const {
00375   
00376     const Scalar* pscal = 
00377         dynamic_cast<const Scalar*>( &(Tensor::divergence(metre)) ) ;
00378 
00379     assert(pscal != 0x0) ;
00380 
00381     return *pscal ;
00382 }
00383 
00384 
00385 Vector Vector::derive_lie(const Vector& vv) const {
00386 
00387     Vector resu(*mp, type_indice(0), triad) ; 
00388     
00389     compute_derive_lie(vv, resu) ;
00390     
00391     return resu ; 
00392     
00393 }
00394 
00395 const Vector_divfree Vector::curl() const {
00396 
00397     if ( triad->identify() == (mp->get_bvect_cart()).identify() ) {
00398     const Metric_flat& metc = mp->flat_met_cart() ;
00399     Vector_divfree resu(*mp, mp->get_bvect_cart(), metc) ;
00400     resu.set(1)= cmp[3]->dsdy() - cmp[2]->dsdz();
00401     resu.set(2)= cmp[1]->dsdz() - cmp[3]->dsdx();
00402     resu.set(3)= cmp[2]->dsdx() - cmp[1]->dsdy();
00403     resu.pseudo_spectral_base();
00404     return resu ;
00405     }
00406     else    {
00407         assert( triad->identify() == (mp->get_bvect_spher()).identify()) ;
00408         const Metric_flat& mets = mp->flat_met_spher() ;
00409         Vector_divfree resu(*mp, mp->get_bvect_spher(), mets);
00410         Scalar tmp = *cmp[3] ;
00411         tmp.div_tant();
00412         tmp += cmp[3]->dsdt();
00413         tmp.div_r();
00414         resu.set(1) = tmp - cmp[2]->srstdsdp() ;
00415         tmp = *cmp[3] ;
00416         tmp.mult_r();
00417         tmp = tmp.dsdr();
00418         tmp.div_r();
00419         resu.set(2) = cmp[1]->srstdsdp() - tmp ;
00420         tmp = *cmp[2];
00421         tmp.mult_r();   
00422         resu.set(3) = tmp.dsdr() - cmp[1]->dsdt() ;
00423         resu.set(3).div_r(); 
00424         resu.pseudo_spectral_base();
00425         return resu ;
00426         }
00427 
00428 }
00429 
00430 
00431 Sym_tensor Vector::ope_killing(const Metric& gam) const {
00432 
00433     Sym_tensor resu(*mp, type_indice(0), *triad) ; 
00434 
00435     if (type_indice(0) == CON ) {
00436         for (int i=1; i<=3; i++) {
00437             for (int j=i; j<=3; j++) {
00438                 resu.set(i,j) = derive_con(gam)(i,j) + derive_con(gam)(j,i)  ;
00439             }
00440         }
00441     }
00442     else {
00443         for (int i=1; i<=3; i++) {
00444             for (int j=i; j<=3; j++) {
00445                 resu.set(i,j) = derive_cov(gam)(i,j) + derive_cov(gam)(j,i)  ;
00446             }
00447         }
00448     }
00449     
00450     return resu ; 
00451 
00452 } 
00453 
00454 
00455 Sym_tensor Vector::ope_killing_conf(const Metric& gam) const {
00456 
00457     Sym_tensor resu(*mp, type_indice(0), *triad) ; 
00458 
00459     if (type_indice(0) == CON ) {
00460         for (int i=1; i<=3; i++) {
00461             for (int j=i; j<=3; j++) {
00462                 resu.set(i,j) = derive_con(gam)(i,j) + derive_con(gam)(j,i)  
00463                                 - 0.6666666666666666* divergence(gam) 
00464                                             * gam.con()(i,j) ;
00465             }
00466         }
00467     }
00468     else {
00469         for (int i=1; i<=3; i++) {
00470             for (int j=i; j<=3; j++) {
00471                 resu.set(i,j) = derive_cov(gam)(i,j) + derive_cov(gam)(j,i)  
00472                                 - 0.6666666666666666* derive_con(gam).trace() 
00473                                             * gam.cov()(i,j) ;
00474             }
00475         }
00476     }
00477 
00478     return resu ; 
00479 
00480 } 
00481 
00482 
00483 
00484 
00485 const Scalar& Vector::potential(const Metric& metre) const {
00486 
00487   set_dependance(metre) ;
00488   int j = get_place_met(metre) ;
00489   assert ((j>=0) && (j<N_MET_MAX)) ;
00490   if (p_potential[j] == 0x0) {
00491     decompose_div(metre) ;
00492   }
00493 
00494   return *p_potential[j] ;
00495 }
00496 
00497 const Vector_divfree& Vector::div_free(const Metric& metre) const {
00498 
00499   set_dependance(metre) ;
00500   int j = get_place_met(metre) ;
00501   assert ((j>=0) && (j<N_MET_MAX)) ;
00502   if (p_div_free[j] == 0x0) {
00503     decompose_div(metre) ;
00504   }
00505 
00506   return *p_div_free[j] ;
00507 }
00508 
00509 void Vector::decompose_div(const Metric& metre) const {
00510 
00511   assert( type_indice(0) == CON ) ; //Only for contravariant vectors...
00512 
00513   set_dependance(metre) ;
00514   int ind =  get_place_met(metre) ;
00515   assert ((ind>=0) && (ind<N_MET_MAX)) ;
00516 
00517   if ( (p_potential[ind] != 0x0) && (p_div_free[ind] != 0x0) ) 
00518     return ; // Nothing to do ...
00519 
00520   else {
00521     if (p_div_free[ind] != 0x0)
00522       delete p_div_free[ind] ;  
00523     if (p_potential[ind] != 0x0) 
00524       delete p_potential[ind] ;
00525 
00526     const Mg3d* mmg = mp->get_mg() ;
00527     int nz = mmg->get_nzone() ;
00528 
00529     int dzp = cmp[0]->get_dzpuis() ;
00530 #ifndef NDEBUG
00531     bool dz_zero = cmp[0]->check_dzpuis(0) ;
00532     bool dz_four = cmp[0]->check_dzpuis(4) ;
00533 #endif
00534     assert( dz_zero || dz_four) ;
00535     assert (cmp[1]->check_dzpuis(dzp)) ;
00536     assert (cmp[2]->check_dzpuis(dzp)) ;
00537 
00538     Scalar dive = divergence(metre) ;
00539 
00540     if (dive.get_etat() == ETATZERO) {
00541     p_potential[ind] = new Scalar(*mp) ;
00542     p_potential[ind]->set_etat_zero() ;
00543     p_potential[ind]->set_dzpuis(dzp) ;
00544     }
00545     else {
00546         //No matching of the solution at this point
00547     p_potential[ind] = new Scalar(dive.poisson()) ; 
00548 
00549     if (dzp == 4) {
00550         assert (mmg->get_type_r(nz-1) == UNSURR) ;
00551         // Let's now do the matching ... in the case dzpuis = 4
00552         const Map_af* mapping = dynamic_cast<const Map_af*>(mp) ;
00553         assert (mapping != 0x0) ;
00554         Valeur& val = p_potential[ind]->set_spectral_va() ;
00555         val.ylm() ;
00556         Mtbl_cf& mtc = *val.c_cf ;
00557         if (nz > 1) {
00558         int np = mmg->get_np(0) ;
00559         int nt = mmg->get_nt(0) ;
00560 #ifndef NDEBUG
00561         for (int lz=0; lz<nz; lz++) {
00562             assert (mmg->get_np(lz) == np) ;
00563             assert (mmg->get_nt(lz) == nt) ;
00564         }
00565 #endif
00566         int lmax = 0 ;
00567         for (int k=0; k<np+1; k++) 
00568             for (int j=0; j<nt; j++) 
00569             if ( nullite_plm(j, nt, k, np, val.base)) {
00570                 int m_quant, l_quant, base_r ;
00571                 donne_lm (nz, 0, j, k, val.base, m_quant, 
00572                       l_quant, base_r) ; 
00573                 lmax = (l_quant > lmax ? l_quant : lmax) ;
00574             }
00575         Scalar** ri = new Scalar*[lmax+1] ;
00576         Scalar** rmi = new Scalar*[lmax+1] ;
00577         Scalar erre(*mp) ;
00578         erre = mp->r ;
00579         for (int l=0; l<=lmax; l++) {
00580             ri[l] = new Scalar(*mp) ;
00581             rmi[l] = new Scalar(*mp) ;
00582             if (l == 0) *(ri[l]) = 1. ;
00583             else *(ri[l]) = pow(erre, l) ;
00584             ri[l]->annule_domain(nz - 1) ;
00585             ri[l]->std_spectral_base() ; //exact base in r will be set later
00586             if (l==2) *(rmi[l]) = 1 ;
00587             else *(rmi[l]) = pow(erre, 2-l) ;
00588             rmi[l]->annule(0,nz-2) ;
00589             Scalar tmp = pow(erre, -l-1) ;
00590             tmp.annule_domain(nz-1) ;
00591             tmp.annule_domain(0) ;
00592             *(rmi[l]) += tmp ;
00593             rmi[l]->set_dzpuis(3) ;
00594             rmi[l]->std_spectral_base() ;//exact base in r will be set later
00595         }
00596 
00597         for (int k=0; k<np+1; k++) {
00598             for (int j=0; j<nt; j++) {
00599             if ( nullite_plm(j, nt, k, np, val.base)) {
00600                 int m_quant, l_quant, base_r, l, c ;
00601                 donne_lm (nz, 0, j, k, val.base, m_quant, l_quant, base_r) ; 
00602                 bool lzero = (l_quant == 0) || (l_quant == 1) ;
00603                 int nb_hom_ced  = (l_quant < 2 ? 0 : 1) ; 
00604           
00605                 int taille = 2*(nz-1) - 1 + nb_hom_ced ;
00606                 Tbl deuz(taille) ;
00607                 deuz.set_etat_qcq() ;
00608                 Matrice systeme(taille,taille) ;
00609                 systeme.set_etat_qcq() ;
00610                 for (l=0; l<taille; l++) 
00611                 for (c=0; c<taille; c++) systeme.set(l,c) = 0 ;
00612                 for (l=0; l<taille; l++) deuz.set(l) = 0 ;
00613           
00614                 //---------
00615                 // Nucleus
00616                 //---------
00617                 assert(mmg->get_type_r(0) == RARE) ;
00618                 assert ((base_r == R_CHEBP)||(base_r == R_CHEBI)) ;
00619                 ri[l_quant]->set_spectral_va().set_base_r(0, base_r) ;
00620       
00621                 double alpha = mapping->get_alpha()[0] ;
00622                 int nr = mmg->get_nr(0) ;
00623                 Tbl partn(nr) ;
00624                 partn.set_etat_qcq() ;
00625                 for (int i=0; i<nr; i++)
00626                 partn.set(i) = mtc(0, k, j, i) ;
00627                 l=0 ; c=0 ;
00628 
00629                 systeme.set(l,c) += pow(alpha, double(l_quant)) ;
00630 
00631                 deuz.set(l) -= val1_dern_1d(0, partn, base_r) ;
00632                 l++ ;
00633 
00634                 if (!lzero) {
00635                 systeme.set(l,c) += l_quant * pow(alpha, double(l_quant - 1)) ;
00636    
00637                 deuz.set(l) -= val1_dern_1d(1, partn, base_r) / alpha ;
00638                 }
00639 
00640                 //----------
00641                 //  Shells
00642                 //----------
00643 
00644                 for (int lz=1; lz<nz-1; lz++) {
00645                 alpha = mapping->get_alpha()[lz] ;
00646                 double beta = mapping->get_beta()[lz] ;
00647                 double xm1 = beta - alpha ;
00648                 double xp1 = alpha + beta ;
00649                 nr = mmg->get_nr(lz) ;
00650                 Tbl parts(nr) ;
00651                 parts.set_etat_qcq() ;
00652                 for (int i=0; i<nr; i++)
00653                     parts.set(i) = mtc(lz, k, j, i) ;
00654   
00655                 //Function at x = -1
00656                 l-- ; 
00657                 c = 2*lz - 1 ;
00658                 systeme.set(l,c) -= pow(xm1, double(l_quant)) ;
00659                 c++;
00660                 systeme.set(l,c) -= pow(xm1, double(-l_quant-1)) ;
00661           
00662                 deuz.set(l) += valm1_dern_1d(0, parts, R_CHEB) ;
00663         
00664                 if ((lz>1) || (!lzero)) {
00665                     //First derivative at x=-1
00666                     l++ ;
00667                     c-- ;
00668                     systeme.set(l,c) -=  l_quant*pow(xm1, double(l_quant - 1)) ;
00669                     c++;
00670                     systeme.set(l,c) -= (-l_quant - 1)*
00671                     pow(xm1, double(-l_quant - 2)) ;
00672 
00673                     deuz.set(l) += valm1_dern_1d(1, parts, R_CHEB) / alpha ;
00674                 }
00675 
00676                 //Function at x = 1
00677                 l++ ; 
00678                 c-- ;
00679                 systeme.set(l,c) += pow(xp1, double(l_quant)) ;
00680                 c++;
00681                 systeme.set(l,c) += pow(xp1, double(-l_quant-1)) ;
00682 
00683                 deuz.set(l) -= val1_dern_1d(0, parts, R_CHEB) ;
00684         
00685                 //First derivative at x = 1
00686                 l++ ; 
00687                 c-- ;
00688                 systeme.set(l,c) +=  l_quant*pow(xp1, double(l_quant - 1)) ;
00689                 c++;
00690                 systeme.set(l,c) += (-l_quant - 1)*
00691                     pow(xp1, double(-l_quant - 2)) ;
00692         
00693                 deuz.set(l) -= val1_dern_1d(1, parts, R_CHEB) / alpha ;
00694         
00695                 } //End of the loop on different shells
00696 
00697                 //-------------------------------
00698                 //  Compactified external domain
00699                 //-------------------------------
00700                 assert(mmg->get_type_r(nz-1) == UNSURR) ;
00701                 nr = mmg->get_nr(nz-1) ;
00702                 Tbl partc(nr) ;
00703                 partc.set_etat_qcq() ;
00704                 for (int i=0; i<nr; i++)
00705                 partc.set(i) = mtc(nz-1, k, j, i) ;
00706 
00707                 alpha = mapping->get_alpha()[nz-1] ;
00708                 double beta =  mapping->get_beta()[nz-1] ;
00709                 double xm1 = beta - alpha ; // 1 / r_left
00710                 double part0, part1 ;
00711                 part0 = valm1_dern_1d(0, partc, R_CHEB) ;
00712                 part1 = xm1*valm1_dern_1d(1, partc, R_CHEB) / alpha ;
00713                 assert (p_potential[ind]->get_dzpuis() == 3) ;
00714         
00715                 //Function at x = -1
00716                 l--;
00717                 if (!lzero) {
00718                 c++;
00719                 systeme.set(l,c) -= pow(xm1, double(l_quant+1)) ;
00720                 }
00721                 deuz.set(l) += part0*xm1*xm1*xm1 ;
00722 
00723                 // First derivative at x = -1
00724                 l++ ;
00725                 if (!lzero) {
00726                 systeme.set(l,c) -= (-l_quant - 1)*
00727                     pow(xm1, double(l_quant + 2)) ;
00728                 }
00729                 if ( (nz > 2) || (!lzero))
00730                 deuz.set(l) += -xm1*xm1*xm1*xm1*(3*part0 + part1) ;
00731 
00732                 //--------------------------------------
00733                 //   Solution of the linear system
00734                 //--------------------------------------
00735       
00736                 int inf = 1 + nb_hom_ced;
00737                 int sup = 3 - nb_hom_ced;
00738                 systeme.set_band(sup, inf) ;
00739                 systeme.set_lu() ;
00740                 Tbl facteur(systeme.inverse(deuz)) ;
00741                 ri[l_quant]->set_spectral_va().coef() ;
00742                 rmi[l_quant]->set_spectral_va().coef() ;
00743 
00744                 //Linear combination in the nucleus
00745                 nr = mmg->get_nr(0) ;
00746                 for (int i=0; i<nr; i++) 
00747                 mtc.set(0, k, j, i) += 
00748          facteur(0)*(*(ri[l_quant]->get_spectral_va().c_cf))(0, 0, 0, i) ;
00749         
00750                 //Linear combination in the shells
00751                 for (int lz=1; lz<nz-1; lz++) {
00752                 nr = mmg->get_nr(lz) ;
00753                 for (int i=0; i<nr; i++) 
00754                     mtc.set(lz, k, j, i) += facteur(2*lz - 1)*
00755            (*(ri[l_quant]->get_spectral_va().c_cf))(lz, 0, 0, i) ;
00756                 for (int i=0; i<nr; i++) 
00757                     mtc.set(lz, k, j, i) += facteur(2*lz)*
00758            (*(rmi[l_quant]->get_spectral_va().c_cf))(lz, 0, 0, i) ;
00759                 }
00760 
00761                 //Linear combination in the CED
00762                 nr = mmg->get_nr(nz-1) ;
00763                 if (!lzero) {
00764                 for (int i=0; i<nr; i++) 
00765                     mtc.set(nz-1, k, j, i) += 
00766                     facteur(taille - 1)*
00767           (*(rmi[l_quant]->get_spectral_va().c_cf))(nz-1, 0, 0, i) ;
00768                 }
00769             } //End of nullite_plm ...
00770 
00771             } //End of j/theta loop   
00772         } //End of k/phi loop 
00773 
00774         for (int l=0; l<=lmax; l++) {
00775             delete ri[l] ;
00776             delete rmi[l] ;
00777         }
00778         delete [] ri ;
00779         delete [] rmi ;
00780 
00781         } //End of the case of more than one domain
00782 
00783         val.ylm_i() ;
00784 
00785     } //End of the case dzp = 4
00786     }
00787     p_div_free[ind] = new Vector_divfree(*mp, *triad, metre) ;
00788 
00789     Vector gradient = p_potential[ind]->derive_con(metre) ;
00790     if (dzp != 4)  gradient.dec_dzpuis(2) ;
00791 
00792     *p_div_free[ind] = ( *this - gradient ) ;
00793 
00794   }
00795   
00796 }
00797 
00798 
00799 
00800 double Vector::flux(double radius, const Metric& met) const {
00801 
00802     assert(type_indice(0) == CON) ; 
00803     
00804     const Map_af* mp_af = dynamic_cast<const Map_af*>(mp) ; 
00805     if (mp_af == 0x0) {
00806         cerr << 
00807         "Vector::flux : the case of a mapping outside the class Map_af\n"
00808             << " is not implemented yet !" << endl ; 
00809         abort() ; 
00810     } 
00811     
00812     const Metric_flat* ff = dynamic_cast<const Metric_flat*>(&met) ;
00813     if (ff == 0x0) {
00814         cerr << 
00815         "Vector::flux : the case of a non flat metric is not implemented yet !"
00816              << endl ; 
00817         abort() ; 
00818     }
00819     
00820     const Base_vect_cart* bcart = dynamic_cast<const Base_vect_cart*>(triad) ; 
00821     Vector* vspher = 0x0 ; 
00822     if (bcart != 0x0) { // switch to spherical components:
00823         vspher = new Vector(*this) ; 
00824         vspher->change_triad(mp->get_bvect_spher()) ;
00825     }
00826     else assert( dynamic_cast<const Base_vect_spher*>(triad) != 0x0 ) ; 
00827 
00828     const Vector* vv = (bcart != 0x0) ? vspher : this ; 
00829 
00830     const Scalar& vr = vv->operator()(1) ; 
00831     
00832     double resu ; 
00833     if (radius == __infinity) {
00834         resu = mp_af->integrale_surface_infini(vr) ; 
00835     }
00836     else {
00837         resu = mp_af->integrale_surface(vr, radius) ;
00838     }
00839     
00840     return resu ; 
00841 }
00842 
00843 void Vector::exponential_filter_r(int lzmin, int lzmax, int p, 
00844               double alpha) {
00845     if( triad->identify() == (mp->get_bvect_cart()).identify() ) 
00846     for (int i=0; i<n_comp; i++)
00847         cmp[i]->exponential_filter_r(lzmin, lzmax, p, alpha) ;
00848     else {
00849     assert( triad->identify() == (mp->get_bvect_spher()).identify()) ;
00850     Scalar vr_tmp = operator()(1) ; 
00851     vr_tmp.exponential_filter_r(lzmin, lzmax, p, alpha) ;
00852     Scalar eta_tmp = eta() ;
00853     eta_tmp.exponential_filter_r(lzmin, lzmax, p, alpha) ;
00854     Scalar mu_tmp = mu() ;
00855     mu_tmp.exponential_filter_r(lzmin, lzmax, p, alpha) ;
00856     set_vr_eta_mu(vr_tmp, eta_tmp, mu_tmp) ;
00857     }
00858 }
00859 
00860 void Vector::exponential_filter_ylm(int lzmin, int lzmax, int p, 
00861               double alpha) {
00862     if( triad->identify() == (mp->get_bvect_cart()).identify() ) 
00863     for (int i=0; i<n_comp; i++)
00864         cmp[i]->exponential_filter_ylm(lzmin, lzmax, p, alpha) ;
00865     else {
00866     assert( triad->identify() == (mp->get_bvect_spher()).identify()) ;
00867     Scalar vr_tmp = operator()(1) ; 
00868     vr_tmp.exponential_filter_ylm(lzmin, lzmax, p, alpha) ;
00869     Scalar eta_tmp = eta() ;
00870     eta_tmp.exponential_filter_ylm(lzmin, lzmax, p, alpha) ;
00871     Scalar mu_tmp = mu() ;
00872     mu_tmp.exponential_filter_ylm(lzmin, lzmax, p, alpha) ;
00873     set_vr_eta_mu(vr_tmp, eta_tmp, mu_tmp) ;
00874     }
00875 }

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