vector_poisson_boundary2.C

00001 /*
00002  *  Method for vector Poisson equation inverting eqs. for V^r and eta as a block.
00003  *
00004  *    (see file vector.h for documentation).
00005  *
00006  */
00007 
00008 /*
00009  *   Copyright (c) 2005  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 version 2
00015  *   as published by the Free Software Foundation.
00016  *
00017  *   LORENE is distributed in the hope that it will be useful,
00018  *   but WITHOUT ANY WARRANTY; without even the implied warranty of
00019  *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00020  *   GNU General Public License for more details.
00021  *
00022  *   You should have received a copy of the GNU General Public License
00023  *   along with LORENE; if not, write to the Free Software
00024  *   Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00025  *
00026  */
00027 
00028 char vector_poisson_boundary2_C[] = "$Header: /cvsroot/Lorene/C++/Source/Tensor/vector_poisson_boundary2.C,v 1.7 2008/08/20 15:07:36 n_vasset Exp $" ;
00029 
00030 /*
00031  * $Id: vector_poisson_boundary2.C,v 1.7 2008/08/20 15:07:36 n_vasset Exp $
00032  * $Log: vector_poisson_boundary2.C,v $
00033  * Revision 1.7  2008/08/20 15:07:36  n_vasset
00034  * Cleaning up the code...
00035  *
00036  * Revision 1.5  2007/09/05 12:35:18  j_novak
00037  * Homogeneous solutions are no longer obtained through the analytic formula, but
00038  * by solving (again) the operator system with almost zero as r.h.s. This is to
00039  * lower the condition number of the matching system.
00040  *
00041  * Revision 1.4  2007/01/23 17:08:46  j_novak
00042  * New function pois_vect_r0.C to solve the l=0 part of the vector Poisson
00043  * equation, which involves only the r-component.
00044  *
00045  * Revision 1.3  2005/12/30 13:39:38  j_novak
00046  * Changed the Galerkin base in the nucleus to (hopefully) stabilise the solver
00047  * when used in an iteration. Similar changes in the CED too.
00048  *
00049  * Revision 1.2  2005/02/15 15:43:18  j_novak
00050  * First version of the block inversion for the vector Poisson equation (method 6).
00051  *
00052  *
00053  * $Header: /cvsroot/Lorene/C++/Source/Tensor/vector_poisson_boundary2.C,v 1.7 2008/08/20 15:07:36 n_vasset Exp $
00054  *
00055  */
00056 
00057 // C headers
00058 #include <assert.h>
00059 #include <stdlib.h>
00060 #include <math.h>
00061 
00062 // Lorene headers
00063 #include "metric.h"
00064 #include "diff.h"
00065 #include "param_elliptic.h"
00066 #include "proto.h"
00067 #include "graphique.h"
00068 #include "map.h"
00069 #include "utilitaires.h"
00070 
00071 
00072 
00073 void Vector::poisson_boundary2(double lam, Vector& resu, Scalar boundvr, Scalar boundeta, Scalar boundmu, double dir_vr, double neum_vr, double dir_eta, double neum_eta, double dir_mu, double neum_mu) const {
00074 
00075     const Map_af* mpaff = dynamic_cast<const Map_af*>(mp) ;
00076 #ifndef NDEBUG 
00077     for (int i=0; i<3; i++)
00078     assert(cmp[i]->check_dzpuis(4)) ;
00079     // All this has a meaning only for spherical components:
00080     const Base_vect_spher* bvs = dynamic_cast<const Base_vect_spher*>(triad) ;
00081     assert(bvs != 0x0) ; 
00082     //## ... and affine mapping, for the moment!
00083     assert( mpaff != 0x0) ;
00084 #endif
00085 
00086      if (fabs(lam + 1.) < 1.e-8) {
00087 
00088        cout <<" lambda = -1 is not supported by this code!" << endl;
00089        abort(); 
00090      }
00091 
00092 
00093 
00094      // Extraction of Mtbl_cf objects from boundary informations; 
00095      Scalar bound_vr2 = boundvr; 
00096      bound_vr2.set_spectral_va().ylm(); 
00097      Scalar bound_eta2 = boundeta; 
00098      bound_eta2.set_spectral_va().ylm();
00099      Scalar bound_mu2 = boundmu; bound_mu2.set_spectral_va().ylm();
00100      
00101      const  Mtbl_cf *bound_vr = bound_vr2.get_spectral_va().c_cf; 
00102      const  Mtbl_cf *bound_mu = bound_mu2.get_spectral_va().c_cf;
00103      const Mtbl_cf *bound_eta = bound_eta2.get_spectral_va().c_cf;
00104 
00105 
00106 
00107 
00108     // Some working objects
00109     //---------------------
00110   const Mg3d& mg = *mpaff->get_mg() ;
00111   int nz = mg.get_nzone() ; int nzm1 = nz - 1;
00112   assert(mg.get_type_r(0) == RARE) ;
00113   assert(mg.get_type_r(nzm1) == UNSURR) ;
00114   Scalar S_r = *cmp[0] ;
00115   Scalar S_eta = eta() ;
00116   Scalar het(*mpaff) ; 
00117   Scalar vr(*mpaff) ; 
00118  
00119   if (S_r.get_etat() == ETATZERO) {
00120       if (S_eta.get_etat() == ETATZERO) {
00121 
00122       S_r.annule_hard() ;
00123       S_r.set_spectral_base(S_eta.get_spectral_base()) ;
00124       S_eta.annule_hard() ;
00125       S_eta.set_spectral_base(S_r.get_spectral_base()) ;
00126       }
00127       else {
00128       S_r.annule_hard() ;
00129       S_r.set_spectral_base(S_eta.get_spectral_base()) ;
00130       }
00131   }
00132   if (S_eta.get_etat() == ETATZERO) {
00133       S_eta.annule_hard() ;
00134       S_eta.set_spectral_base(S_r.get_spectral_base()) ;
00135   }
00136 
00137   S_r.set_spectral_va().ylm() ;
00138   S_eta.set_spectral_va().ylm() ;
00139   const Base_val& base = S_eta.get_spectral_va().base ;
00140   Mtbl_cf sol_part_eta(mg, base) ; sol_part_eta.annule_hard() ;
00141   Mtbl_cf sol_part_vr(mg, base) ; sol_part_vr.annule_hard() ;
00142   Mtbl_cf sol_hom_un_eta(mg, base) ; sol_hom_un_eta.annule_hard() ;
00143   Mtbl_cf sol_hom_deux_eta(mg, base) ; sol_hom_deux_eta.annule_hard() ;
00144   Mtbl_cf sol_hom_trois_eta(mg, base) ; sol_hom_trois_eta.annule_hard() ;
00145   Mtbl_cf sol_hom_quatre_eta(mg, base) ; sol_hom_quatre_eta.annule_hard() ;
00146   Mtbl_cf sol_hom_un_vr(mg, base) ; sol_hom_un_vr.annule_hard() ;
00147   Mtbl_cf sol_hom_deux_vr(mg, base) ; sol_hom_deux_vr.annule_hard() ;
00148   Mtbl_cf sol_hom_trois_vr(mg, base) ; sol_hom_trois_vr.annule_hard() ;
00149   Mtbl_cf sol_hom_quatre_vr(mg, base) ; sol_hom_quatre_vr.annule_hard() ;
00150  
00151   // Solution of the l=0 part for Vr
00152   //---------------------------------
00153   Scalar sou_l0 = (*cmp[0]) / (1. + lam) ;
00154   sou_l0.set_spectral_va().ylm() ;
00155 
00156  Param_elliptic param_l0(sou_l0) ;
00157   for (int l=0; l<nz; l++){
00158       param_l0.set_poisson_vect_r(l, true) ;
00159   }
00160 
00161 
00162    Scalar vrl0 = sou_l0.sol_elliptic_boundary(param_l0, *bound_vr, dir_vr, neum_vr) ;
00163 
00164 
00165   // Build-up & inversion of the system for (eta, V^r) in each domain (except for the nucleus!)
00166   //-----------------------------------------------------------------
00167 
00168 
00169   // Shells
00170   //-------
00171 
00172   int nr = mg.get_nr(1) ;
00173   int nt = mg.get_nt(1) ;
00174   int np = mg.get_np(1) ;
00175   double alpha = mpaff->get_alpha()[1] ; double alp2 = alpha*alpha ;
00176   double beta = mpaff->get_beta()[1] ;
00177  
00178   int l_q = 0 ; int m_q = 0; int base_r = 0 ;
00179   for (int zone=1 ; zone<nzm1 ; zone++) {
00180       nr = mg.get_nr(zone) ; 
00181       assert (nr > 5) ;
00182       assert(nt == mg.get_nt(zone)) ;
00183       assert(np == mg.get_np(zone)) ;
00184       alpha = mpaff->get_alpha()[zone] ;
00185       beta = mpaff->get_beta()[zone] ;
00186       double ech = beta / alpha ;
00187 
00188   // Loop on l and m
00189   //----------------
00190       for (int k=0 ; k<np+1 ; k++) {
00191       for (int j=0 ; j<nt ; j++) {
00192       base.give_quant_numbers(zone, k, j, m_q, l_q, base_r) ;
00193       if ( (nullite_plm(j, nt, k, np, base) == 1) && (l_q != 0) ) {
00194           int dege1 = 2 ; //degeneracy of eq.1
00195           int dege2 = 2 ;  //degeneracy of eq.2
00196           int nr_eq1 = nr - dege1 ; //Eq.1 is for eta
00197           int nr_eq2 = nr - dege2 ; //Eq.2 is for V^r
00198           int nrtot = 2*nr ;
00199           Matrice oper(nrtot, nrtot) ; oper.set_etat_qcq() ;
00200           Tbl sec_membre(nrtot) ; sec_membre.set_etat_qcq() ;
00201           Diff_x2dsdx2 x2d2(base_r, nr); const Matrice& m2d2 = x2d2.get_matrice() ;
00202           Diff_xdsdx2 xd2(base_r, nr) ; const Matrice& mxd2 = xd2.get_matrice() ;
00203           Diff_dsdx2 d2(base_r, nr) ; const Matrice& md2 = d2.get_matrice() ;
00204           Diff_xdsdx xd(base_r, nr) ; const Matrice& mxd = xd.get_matrice() ;
00205           Diff_dsdx d1(base_r, nr) ; const Matrice& md = d1.get_matrice() ;
00206           Diff_id id(base_r, nr) ; const Matrice& mid = id.get_matrice() ;
00207 
00208           // Building the operator
00209           //----------------------
00210           for (int lin=0; lin<nr_eq1; lin++) { 
00211           for (int col=0; col<nr; col++) 
00212               oper.set(lin,col) 
00213               = m2d2(lin,col) + 2*ech*mxd2(lin,col) + ech*ech*md2(lin,col) 
00214               + 2*(mxd(lin,col) + ech*md(lin,col)) 
00215               - (lam+1)*l_q*(l_q+1)*mid(lin,col) ;
00216           for (int col=0; col<nr; col++) 
00217               oper.set(lin,col+nr) 
00218               = lam*(mxd(lin,col) + ech*md(lin,col)) 
00219               + 2*(1+lam)*mid(lin,col) ; 
00220           }
00221           oper.set(nr_eq1, 0) = 1 ;
00222           for (int col=1; col<2*nr; col++) 
00223           oper.set(nr_eq1, col) = 0 ;
00224           oper.set(nr_eq1+1, 0) = 0 ;
00225           oper.set(nr_eq1+1, 1) = 1 ;
00226           for (int col=2; col<2*nr; col++) 
00227           oper.set(nr_eq1+1, col) = 0 ;
00228           for (int lin=0; lin<nr_eq2; lin++) {
00229           for (int col=0; col<nr; col++)
00230               oper.set(lin+nr,col) 
00231               = -l_q*(l_q+1)*(lam*(mxd(lin,col) + ech*md(lin,col))
00232                       - (lam+2)*mid(lin,col)) ;
00233           for (int col=0; col<nr; col++)
00234               oper.set(lin+nr, col+nr) 
00235               = (lam+1)*(m2d2(lin,col) + 2*ech*mxd2(lin,col) 
00236                      + ech*ech*md2(lin,col) 
00237                      + 2*(mxd(lin,col) + ech*md(lin,col)))
00238               -(2*(lam+1)+l_q*(l_q+1))*mid(lin,col) ;         
00239           }
00240           for (int col=0; col<nr; col++) 
00241           oper.set(nr+nr_eq2, col) = 0 ;
00242           oper.set(nr+nr_eq2, nr) = 1 ;
00243           for (int col=nr+1; col<2*nr; col++) 
00244           oper.set(nr+nr_eq2, col) = 0 ;
00245           for (int col=0; col<nr+1; col++) 
00246           oper.set(nr+nr_eq2+1, col) = 0 ;
00247           oper.set(nr+nr_eq2+1, nr+1) = 1 ;
00248           for (int col=nr+2; col<2*nr; col++) 
00249           oper.set(nr+nr_eq2+1, col) = 0 ;
00250 
00251           oper.set_lu() ;
00252           
00253           // Filling the r.h.s
00254           //------------------
00255           Tbl sr(nr) ; sr.set_etat_qcq() ; 
00256           Tbl seta(nr) ; seta.set_etat_qcq() ;
00257           for (int i=0; i<nr; i++) {
00258           sr.set(i) = (*S_r.get_spectral_va().c_cf)(zone, k, j, i);
00259           seta.set(i) = (*S_eta.get_spectral_va().c_cf)(zone, k, j, i) ;
00260           }
00261           Tbl xsr= sr ;  Tbl x2sr= sr ;
00262           Tbl xseta= seta ; Tbl x2seta = seta ;
00263           multx2_1d(nr, &x2sr.t, base_r) ; multx_1d(nr, &xsr.t, base_r) ;
00264           multx2_1d(nr, &x2seta.t, base_r) ; multx_1d(nr, &xseta.t, base_r) ;
00265           
00266           for (int i=0; i<nr_eq1; i++) 
00267           sec_membre.set(i) = alpha*(alpha*x2seta(i) + 2*beta*xseta(i))
00268               + beta*beta*seta(i);
00269           sec_membre.set(nr_eq1) = 0 ;
00270           sec_membre.set(nr_eq1+1) = 0 ;
00271           for (int i=0; i<nr_eq2; i++)
00272           sec_membre.set(i+nr) = beta*beta*sr(i) 
00273               + alpha*(alpha*x2sr(i) + 2*beta*xsr(i)) ;
00274           sec_membre.set(nr+nr_eq2) = 0 ;
00275           sec_membre.set(nr+nr_eq2+1) = 0 ;
00276 
00277           // Inversion of the "big" operator
00278           //--------------------------------
00279           Tbl big_res = oper.inverse(sec_membre) ;
00280           for (int i=0; i<nr; i++) {
00281           sol_part_eta.set(zone, k, j, i) = big_res(i) ;
00282           sol_part_vr.set(zone, k, j, i) = big_res(nr+i) ;        
00283           }
00284           
00285           // Getting the homogeneous solutions
00286           //----------------------------------
00287           sec_membre.annule_hard() ;
00288           sec_membre.set(nr_eq1) = 1 ;
00289           big_res = oper.inverse(sec_membre) ;
00290           for (int i=0 ; i<nr ; i++) {
00291           sol_hom_un_eta.set(zone, k, j, i) = big_res(i) ;
00292           sol_hom_un_vr.set(zone, k, j, i) = big_res(nr+i) ;
00293           }
00294           sec_membre.set(nr_eq1) = 0 ;
00295           sec_membre.set(nr_eq1+1) = 1 ;
00296           big_res = oper.inverse(sec_membre) ;
00297           for (int i=0 ; i<nr ; i++) {
00298           sol_hom_deux_eta.set(zone, k, j, i) = big_res(i) ;
00299           sol_hom_deux_vr.set(zone, k, j, i) = big_res(nr+i) ;
00300           }
00301           sec_membre.set(nr_eq1+1) = 0 ;
00302           sec_membre.set(nr+nr_eq2) = 1 ;
00303           big_res = oper.inverse(sec_membre) ;
00304           for (int i=0 ; i<nr ; i++) {
00305           sol_hom_trois_eta.set(zone, k, j, i) = big_res(i) ;
00306           sol_hom_trois_vr.set(zone, k, j, i) = big_res(nr+i) ;
00307           }
00308           sec_membre.set(nr+nr_eq2) = 0 ;
00309           sec_membre.set(nr+nr_eq2+1) = 1 ;
00310           big_res = oper.inverse(sec_membre) ;
00311           for (int i=0 ; i<nr ; i++) {
00312           sol_hom_quatre_eta.set(zone, k, j, i) = big_res(i) ;
00313           sol_hom_quatre_vr.set(zone, k, j, i) = big_res(nr+i) ;
00314           }
00315 
00316       } 
00317       }
00318       }
00319   }
00320 
00321   // Compactified external domain
00322   //-----------------------------
00323   nr = mg.get_nr(nzm1) ;
00324   assert(nt == mg.get_nt(nzm1)) ;
00325   assert(np == mg.get_np(nzm1)) ;
00326   alpha = mpaff->get_alpha()[nzm1] ; alp2 = alpha*alpha ;
00327   assert (nr > 4) ;
00328   
00329   // Loop on l and m
00330   //----------------
00331   for (int k=0 ; k<np+1 ; k++) {
00332       for (int j=0 ; j<nt ; j++) { 
00333       base.give_quant_numbers(nzm1, k, j, m_q, l_q, base_r) ;
00334       if ( (nullite_plm(j, nt, k, np, base) == 1) && (l_q != 0) ) {
00335           int dege1 = 3; //degeneracy of eq.1
00336           int dege2 = (l_q == 1 ? 2 : 3); //degeneracy of eq.2
00337           int nr_eq1 = nr - dege1 ; //Eq.1 is for eta
00338           int nr_eq2 = nr - dege2 ; //Eq.2 is the div-free condition
00339           int nrtot = 2*nr ;
00340           Matrice oper(nrtot, nrtot) ; oper.set_etat_qcq() ;
00341           Tbl sec_membre(nrtot) ; sec_membre.set_etat_qcq() ;
00342           Diff_dsdx2 d2(base_r, nr) ; const Matrice& md2 = d2.get_matrice() ;
00343           Diff_sxdsdx sxd(base_r, nr) ; const Matrice& mxd = sxd.get_matrice() ;
00344           Diff_sx2 sx2(base_r, nr) ; const Matrice& ms2 = sx2.get_matrice() ;
00345 
00346           // Building the operator
00347           //----------------------
00348           for (int lin=0; lin<nr_eq1; lin++) { 
00349           for (int col=0; col<nr; col++) 
00350               oper.set(lin,col) 
00351               = (md2(lin,col) - (lam+1)*l_q*(l_q+1)*ms2(lin,col))/alp2 ;
00352               for (int col=0; col<nr; col++) 
00353                   oper.set(lin,col+nr) = 
00354                   (-lam*mxd(lin,col) + 2*(1+lam)*ms2(lin,col)) / alp2 ;
00355           }
00356           for (int col=0; col<nr; col++)
00357           oper.set(nr_eq1, col) = 1 ;
00358           for (int col=nr; col<nrtot; col++)
00359           oper.set(nr_eq1, col) = 0 ;
00360           int pari = -1 ;
00361           for (int col=0; col<nr; col++) {
00362           oper.set(nr_eq1+1, col) = pari*col*col ;
00363           pari = pari ;
00364           }
00365           for (int col=nr; col<nrtot; col++)
00366           oper.set(nr_eq1+1, col) = 0 ;
00367           oper.set(nr_eq1+2, 0) = 1 ;
00368           for (int col=1; col<nrtot; col++)
00369           oper.set(nr_eq1+2, col) = 0 ;
00370           for (int lin=0; lin<nr_eq2; lin++) {
00371           for (int col=0; col<nr; col++)
00372               oper.set(lin+nr,col) = (l_q*(l_q+1)*(lam*mxd(lin,col) 
00373                    + (lam+2)*ms2(lin,col))) / alp2 ;
00374           for (int col=0; col<nr; col++)
00375               oper.set(lin+nr, col+nr) = ((lam+1)*md2(lin,col) 
00376                  - (2*(lam+1) + l_q*(l_q+1))*ms2(lin,col)) / alp2 ;
00377           }
00378           for (int col=0; col<nr; col++)
00379           oper.set(nr+nr_eq2, col) = 0 ;
00380           for (int col=nr; col<nrtot; col++)
00381           oper.set(nr+nr_eq2, col) = 1 ;
00382           for (int col=0; col<nr; col++)
00383           oper.set(nr+nr_eq2+1, col) = 0 ;
00384           pari = -1 ;
00385           for (int col=0; col<nr; col++) {
00386           oper.set(nr+nr_eq2+1, nr+col) = pari*col*col ;
00387           pari = pari ;
00388           }
00389           if (l_q > 1) {
00390           for (int col=0; col<nr; col++)
00391               oper.set(nr+nr_eq2+2, col) = 0 ;
00392           oper.set(nr+nr_eq2+2, nr) = 1 ;
00393           for (int col=nr+1; col<nrtot; col++)
00394               oper.set(nr+nr_eq2+2, col) = 0 ;
00395           }
00396           oper.set_lu() ;
00397 
00398           // Filling the r.h.s
00399           //------------------
00400           for (int i=0; i<nr_eq1; i++) 
00401           sec_membre.set(i) = (*S_eta.get_spectral_va().c_cf)(nzm1, k, j, i) ;
00402           for (int i=nr_eq1; i<nr; i++)
00403           sec_membre.set(i) = 0 ;
00404           for (int i=0; i<nr_eq2; i++)
00405           sec_membre.set(i+nr) =(*S_r.get_spectral_va().c_cf)(nzm1, k, j, i);
00406           for (int i=nr_eq2; i<nr; i++)
00407           sec_membre.set(nr+i) = 0 ;
00408           Tbl big_res = oper.inverse(sec_membre) ;
00409           
00410           // Putting coefficients of h and v to individual arrays
00411           //-----------------------------------------------------
00412           for (int i=0; i<nr; i++) {
00413           sol_part_eta.set(nzm1, k, j, i) = big_res(i) ;
00414           sol_part_vr.set(nzm1, k, j, i) = big_res(i+nr) ;
00415           }
00416 
00417           // Homogeneous solutions (only 1/r^(l+2) and 1/r^l in the CED)
00418           //------------------------------------------------------------
00419           if (l_q == 1) {         
00420           Tbl sol_hom1 = solh(nr, 0, 0., base_r) ;
00421           Tbl sol_hom2 = solh(nr, 2, 0., base_r) ;
00422           double fac_eta = lam + 2. ;
00423           double fac_vr = 2*lam + 2. ;
00424           for (int i=0 ; i<nr ; i++) {
00425               sol_hom_deux_eta.set(nzm1, k, j, i) = sol_hom2(i) ;
00426               sol_hom_quatre_eta.set(nzm1, k, j, i) = fac_eta*sol_hom1(i) ;
00427               sol_hom_deux_vr.set(nzm1, k, j, i) = -2*sol_hom2(i) ;
00428               sol_hom_quatre_vr.set(nzm1, k, j, i) = fac_vr*sol_hom1(i) ;
00429           }
00430           }
00431           else {
00432           sec_membre.annule_hard() ;
00433           sec_membre.set(nr-1) = 1 ;
00434           big_res = oper.inverse(sec_membre) ;
00435 
00436           for (int i=0; i<nr; i++) {
00437               sol_hom_deux_eta.set(nzm1, k, j, i) = big_res(i) ;
00438               sol_hom_deux_vr.set(nzm1, k, j, i) = big_res(nr+i) ;
00439           }
00440           sec_membre.set(nr-1) = 0 ;
00441           sec_membre.set(2*nr-1) = 1 ;
00442           big_res = oper.inverse(sec_membre) ;
00443           for (int i=0; i<nr; i++) {
00444               sol_hom_quatre_eta.set(nzm1, k, j, i) = big_res(i) ;
00445               sol_hom_quatre_vr.set(nzm1, k, j, i) = big_res(nr+i) ;
00446           }
00447           }
00448       }
00449       }
00450   }
00451 
00452   // Now let's match everything ...
00453   //-------------------------------
00454 
00455   // Resulting V^r & eta
00456   vr.set_etat_qcq() ;
00457   vr.set_spectral_base(base) ;
00458   vr.set_spectral_va().set_etat_cf_qcq() ;
00459   Mtbl_cf& cf_vr = *vr.set_spectral_va().c_cf ;
00460   cf_vr.annule_hard() ;
00461   het.set_etat_qcq() ;
00462   het.set_spectral_base(base) ;
00463   het.set_spectral_va().set_etat_cf_qcq() ;
00464   Mtbl_cf& cf_eta = *het.set_spectral_va().c_cf ;
00465   cf_eta.annule_hard() ;
00466   int taille = 4*nzm1 -2 ; // Two less than R3 case 
00467   Tbl sec_membre(taille) ; 
00468   Matrice systeme(taille, taille) ; 
00469   systeme.set_etat_qcq() ;
00470   int ligne ;  int colonne ;
00471   
00472   // Loop on l and m
00473   //----------------
00474   for (int k=0 ; k<np+1 ; k++)
00475       for (int j=0 ; j<nt ; j++) {
00476       base.give_quant_numbers(0, k, j, m_q, l_q, base_r) ;
00477       if ((nullite_plm(j, nt, k, np, base) == 1)&&(l_q != 0)) {
00478         
00479           ligne = 0 ;
00480           colonne = 0 ;
00481           systeme.annule_hard() ;
00482           sec_membre.annule_hard() ;
00483 
00484 
00485               //shell 1 (boundary condition)
00486           int zone1 = 1;
00487           nr = mg.get_nr(zone1) ;
00488           alpha = mpaff->get_alpha()[zone1] ;
00489 
00490           // Here we prepare for a Robyn-type boundary condition on eta. 
00491           // Parameters are dir_eta and neum_eta
00492           systeme.set(ligne, colonne) 
00493               = -dir_eta*sol_hom_un_eta.val_in_bound_jk(zone1, j, k)  ;
00494           systeme.set(ligne, colonne+1) 
00495               = -dir_eta * sol_hom_deux_eta.val_in_bound_jk(zone1, j, k)  ;
00496           systeme.set(ligne, colonne+2) 
00497               = -dir_eta*sol_hom_trois_eta.val_in_bound_jk(zone1, j, k)  ;
00498           systeme.set(ligne, colonne+3) 
00499               = -dir_eta*sol_hom_quatre_eta.val_in_bound_jk(zone1, j, k)  ;
00500  
00501  
00502 
00503 
00504           sec_membre.set(ligne) += dir_eta* sol_part_eta.val_in_bound_jk(zone1, j, k) ;
00505 
00506 
00507 
00508           int pari = -1 ;
00509           for (int i=0; i<nr; i++) {
00510               systeme.set(ligne, colonne) 
00511               -= neum_eta*pari*i*i*sol_hom_un_eta(zone1, k, j, i)/alpha ;
00512               systeme.set(ligne, colonne+1) 
00513               -= neum_eta*pari*i*i*sol_hom_deux_eta(zone1, k, j, i)/alpha ;
00514               systeme.set(ligne, colonne+2) 
00515               -= neum_eta*pari*i*i*sol_hom_trois_eta(zone1, k, j, i)/alpha ;
00516               systeme.set(ligne, colonne+3) 
00517               -= neum_eta*pari*i*i*sol_hom_quatre_eta(zone1, k, j, i)/alpha ;
00518               sec_membre.set(ligne) 
00519               += neum_eta*pari*i*i* sol_part_eta(zone1, k, j, i)/alpha ;
00520               pari = -pari ;
00521           }
00522 
00523           sec_membre.set(ligne) -= (*bound_eta).val_in_bound_jk(1,j,k);
00524           ligne++ ;
00525  
00526                  
00527 
00528           // ... and their couterparts for V^r
00529           systeme.set(ligne, colonne) 
00530               = - dir_vr*sol_hom_un_vr.val_in_bound_jk(zone1, j, k)  ;
00531           systeme.set(ligne, colonne+1) 
00532               = - dir_vr*sol_hom_deux_vr.val_in_bound_jk(zone1, j, k)  ;
00533           systeme.set(ligne, colonne+2) 
00534               = -dir_vr*sol_hom_trois_vr.val_in_bound_jk(zone1, j, k)  ;
00535           systeme.set(ligne, colonne+3) 
00536               = -dir_vr*sol_hom_quatre_vr.val_in_bound_jk(zone1, j, k)  ;
00537           sec_membre.set(ligne) += dir_vr*sol_part_vr.val_in_bound_jk(zone1, j, k) ;
00538     
00539 
00540     
00541     
00542           pari = -1 ;
00543           for (int i=0; i<nr; i++) {
00544               systeme.set(ligne, colonne) 
00545               -= neum_vr*pari*i*i*sol_hom_un_vr(zone1, k, j, i)/alpha ;
00546               systeme.set(ligne, colonne+1) 
00547               -= neum_vr*pari*i*i*sol_hom_deux_vr(zone1, k, j, i)/alpha ;
00548               systeme.set(ligne, colonne+2) 
00549               -= neum_vr*pari*i*i*sol_hom_trois_vr(zone1, k, j, i)/alpha ;
00550               systeme.set(ligne, colonne+3) 
00551               -= neum_vr*pari*i*i*sol_hom_quatre_vr(zone1, k, j, i)/alpha ;
00552               sec_membre.set(ligne) 
00553               += neum_vr*pari*i*i* sol_part_vr(zone1, k, j, i)/alpha ;
00554               pari = -pari ;
00555           }
00556 
00557                   sec_membre.set(ligne) -= (*bound_vr).val_in_bound_jk(1,j,k);
00558           ligne++ ;
00559 
00560 
00561 
00562     
00563           // Values in 1: beginning with solution in eta...
00564             
00565           systeme.set(ligne, colonne) 
00566               += sol_hom_un_eta.val_out_bound_jk(zone1, j, k)  ;
00567           systeme.set(ligne, colonne+1) 
00568               += sol_hom_deux_eta.val_out_bound_jk(zone1, j, k)  ;
00569           systeme.set(ligne, colonne+2) 
00570               += sol_hom_trois_eta.val_out_bound_jk(zone1, j, k)  ;
00571           systeme.set(ligne, colonne+3) 
00572               += sol_hom_quatre_eta.val_out_bound_jk(zone1, j, k)  ;
00573           sec_membre.set(ligne) -= sol_part_eta.val_out_bound_jk(zone1, j, k) ;
00574           ligne++ ;
00575           // ... and their counterparts for V^r
00576           systeme.set(ligne, colonne) 
00577               += sol_hom_un_vr.val_out_bound_jk(zone1, j, k)  ;
00578           systeme.set(ligne, colonne+1) 
00579               += sol_hom_deux_vr.val_out_bound_jk(zone1, j, k)  ;
00580           systeme.set(ligne, colonne+2) 
00581               += sol_hom_trois_vr.val_out_bound_jk(zone1, j, k)  ;
00582           systeme.set(ligne, colonne+3) 
00583               += sol_hom_quatre_vr.val_out_bound_jk(zone1, j, k)  ;
00584           sec_membre.set(ligne) -= sol_part_vr.val_out_bound_jk(zone1, j, k) ;
00585           ligne++ ;
00586 
00587           //derivative at 1 
00588           pari = 1 ;
00589           for (int i=0; i<nr; i++) {
00590               systeme.set(ligne, colonne) 
00591               += pari*i*i*sol_hom_un_eta(zone1, k, j, i)/alpha ;
00592               systeme.set(ligne, colonne+1) 
00593               += pari*i*i*sol_hom_deux_eta(zone1, k, j, i)/alpha ;
00594               systeme.set(ligne, colonne+2) 
00595               += pari*i*i*sol_hom_trois_eta(zone1, k, j, i)/alpha ;
00596               systeme.set(ligne, colonne+3) 
00597               += pari*i*i*sol_hom_quatre_eta(zone1, k, j, i)/alpha ;
00598               sec_membre.set(ligne) 
00599               -= pari*i*i* sol_part_eta(zone1, k, j, i)/alpha ;
00600               pari = pari ;
00601           }
00602           ligne++ ;
00603           // ... and their counterparts for V^r
00604           pari = 1 ;
00605           for (int i=0; i<nr; i++) {
00606               systeme.set(ligne, colonne) 
00607               += pari*i*i*sol_hom_un_vr(zone1, k, j, i)/alpha ;
00608               systeme.set(ligne, colonne+1) 
00609               += pari*i*i*sol_hom_deux_vr(zone1, k, j, i)/alpha ;
00610               systeme.set(ligne, colonne+2) 
00611               += pari*i*i*sol_hom_trois_vr(zone1, k, j, i)/alpha ;
00612               systeme.set(ligne, colonne+3) 
00613               += pari*i*i*sol_hom_quatre_vr(zone1, k, j, i)/alpha ;
00614               sec_membre.set(ligne) 
00615               -= pari*i*i* sol_part_vr(zone1, k, j, i)/alpha ;
00616               pari = pari ;
00617           }
00618           colonne += 4 ;
00619 
00620           // Other shells
00621           if ( nz <= 2){
00622             cout <<"WARNING!! Mapping must contain at least 2 shells!!" << endl;}
00623           for (int zone=2 ; zone<nzm1 ; zone++) {
00624 
00625       nr = mg.get_nr(zone) ;
00626           alpha = mpaff->get_alpha()[zone] ;
00627           ligne -= 3 ;
00628           systeme.set(ligne, colonne) 
00629               = -sol_hom_un_eta.val_in_bound_jk(zone, j, k)  ;
00630           systeme.set(ligne, colonne+1) 
00631               = -sol_hom_deux_eta.val_in_bound_jk(zone, j, k)  ;
00632           systeme.set(ligne, colonne+2) 
00633               = -sol_hom_trois_eta.val_in_bound_jk(zone, j, k)  ;
00634           systeme.set(ligne, colonne+3) 
00635               = -sol_hom_quatre_eta.val_in_bound_jk(zone, j, k)  ;
00636           sec_membre.set(ligne) += sol_part_eta.val_in_bound_jk(zone, j, k) ;
00637           ligne++ ;
00638           // ... and their counterparts for V^r
00639           systeme.set(ligne, colonne) 
00640               = -sol_hom_un_vr.val_in_bound_jk(zone, j, k)  ;
00641           systeme.set(ligne, colonne+1) 
00642               = -sol_hom_deux_vr.val_in_bound_jk(zone, j, k)  ;
00643           systeme.set(ligne, colonne+2) 
00644               = -sol_hom_trois_vr.val_in_bound_jk(zone, j, k)  ;
00645           systeme.set(ligne, colonne+3) 
00646               = -sol_hom_quatre_vr.val_in_bound_jk(zone, j, k)  ;
00647           sec_membre.set(ligne) += sol_part_vr.val_in_bound_jk(zone, j, k) ;
00648           ligne++ ;
00649 
00650           //derivative of (x+echelle)^(l-1) at -1 
00651           pari = -1 ;
00652           for (int i=0; i<nr; i++) {
00653               systeme.set(ligne, colonne) 
00654               -= pari*i*i*sol_hom_un_eta(zone, k, j, i)/alpha ;
00655               systeme.set(ligne, colonne+1) 
00656               -= pari*i*i*sol_hom_deux_eta(zone, k, j, i)/alpha ;
00657               systeme.set(ligne, colonne+2) 
00658               -= pari*i*i*sol_hom_trois_eta(zone, k, j, i)/alpha ;
00659               systeme.set(ligne, colonne+3) 
00660               -= pari*i*i*sol_hom_quatre_eta(zone, k, j, i)/alpha ;
00661               sec_membre.set(ligne) 
00662               += pari*i*i* sol_part_eta(zone, k, j, i)/alpha ;
00663               pari = -pari ;
00664           }
00665           ligne++ ;
00666           // ... and their counterparts for V^r
00667           pari = -1 ;
00668           for (int i=0; i<nr; i++) {
00669               systeme.set(ligne, colonne) 
00670               -= pari*i*i*sol_hom_un_vr(zone, k, j, i)/alpha ;
00671               systeme.set(ligne, colonne+1) 
00672               -= pari*i*i*sol_hom_deux_vr(zone, k, j, i)/alpha ;
00673               systeme.set(ligne, colonne+2) 
00674               -= pari*i*i*sol_hom_trois_vr(zone, k, j, i)/alpha ;
00675               systeme.set(ligne, colonne+3) 
00676               -= pari*i*i*sol_hom_quatre_vr(zone, k, j, i)/alpha ;
00677               sec_membre.set(ligne) 
00678               += pari*i*i* sol_part_vr(zone, k, j, i)/alpha ;
00679               pari = -pari ;
00680           }
00681           ligne++ ;
00682             
00683           systeme.set(ligne, colonne) 
00684               += sol_hom_un_eta.val_out_bound_jk(zone, j, k)  ;
00685           systeme.set(ligne, colonne+1) 
00686               += sol_hom_deux_eta.val_out_bound_jk(zone, j, k)  ;
00687           systeme.set(ligne, colonne+2) 
00688               += sol_hom_trois_eta.val_out_bound_jk(zone, j, k)  ;
00689           systeme.set(ligne, colonne+3) 
00690               += sol_hom_quatre_eta.val_out_bound_jk(zone, j, k)  ;
00691           sec_membre.set(ligne) -= sol_part_eta.val_out_bound_jk(zone, j, k) ;
00692           ligne++ ;
00693           // ... and their counterparts for V^r
00694           systeme.set(ligne, colonne) 
00695               += sol_hom_un_vr.val_out_bound_jk(zone, j, k)  ;
00696           systeme.set(ligne, colonne+1) 
00697               += sol_hom_deux_vr.val_out_bound_jk(zone, j, k)  ;
00698           systeme.set(ligne, colonne+2) 
00699               += sol_hom_trois_vr.val_out_bound_jk(zone, j, k)  ;
00700           systeme.set(ligne, colonne+3) 
00701               += sol_hom_quatre_vr.val_out_bound_jk(zone, j, k)  ;
00702           sec_membre.set(ligne) -= sol_part_vr.val_out_bound_jk(zone, j, k) ;
00703           ligne++ ;
00704 
00705           //derivative at 1 
00706           pari = 1 ;
00707           for (int i=0; i<nr; i++) {
00708               systeme.set(ligne, colonne) 
00709               += pari*i*i*sol_hom_un_eta(zone, k, j, i)/alpha ;
00710               systeme.set(ligne, colonne+1) 
00711               += pari*i*i*sol_hom_deux_eta(zone, k, j, i)/alpha ;
00712               systeme.set(ligne, colonne+2) 
00713               += pari*i*i*sol_hom_trois_eta(zone, k, j, i)/alpha ;
00714               systeme.set(ligne, colonne+3) 
00715               += pari*i*i*sol_hom_quatre_eta(zone, k, j, i)/alpha ;
00716               sec_membre.set(ligne) 
00717               -= pari*i*i* sol_part_eta(zone, k, j, i)/alpha ;
00718               pari = pari ;
00719           }
00720           ligne++ ;
00721           // ... and their couterparts for V^r
00722           pari = 1 ;
00723           for (int i=0; i<nr; i++) {
00724               systeme.set(ligne, colonne) 
00725               += pari*i*i*sol_hom_un_vr(zone, k, j, i)/alpha ;
00726               systeme.set(ligne, colonne+1) 
00727               += pari*i*i*sol_hom_deux_vr(zone, k, j, i)/alpha ;
00728               systeme.set(ligne, colonne+2) 
00729               += pari*i*i*sol_hom_trois_vr(zone, k, j, i)/alpha ;
00730               systeme.set(ligne, colonne+3) 
00731               += pari*i*i*sol_hom_quatre_vr(zone, k, j, i)/alpha ;
00732               sec_membre.set(ligne) 
00733               -= pari*i*i* sol_part_vr(zone, k, j, i)/alpha ;
00734               pari = pari ;
00735           }
00736           colonne += 4 ;
00737           }         
00738                     
00739           //Compactified external domain
00740           nr = mg.get_nr(nzm1) ;
00741 
00742           alpha = mpaff->get_alpha()[nzm1] ;
00743           ligne -= 3 ;
00744           //value of (x-1)^(l+2) at -1 :
00745           systeme.set(ligne, colonne) 
00746           -= sol_hom_deux_eta.val_in_bound_jk(nzm1, j, k) ;
00747           systeme.set(ligne, colonne+1) 
00748           -= sol_hom_quatre_eta.val_in_bound_jk(nzm1, j, k) ;
00749           sec_membre.set(ligne) += sol_part_eta.val_in_bound_jk(nzm1, j, k) ;
00750           //... and of its counterpart for V^r
00751           systeme.set(ligne+1, colonne) 
00752           -= sol_hom_deux_vr.val_in_bound_jk(nzm1, j, k) ;
00753           systeme.set(ligne+1, colonne+1) 
00754           -= sol_hom_quatre_vr.val_in_bound_jk(nzm1, j, k) ;
00755           sec_membre.set(ligne+1) += sol_part_vr.val_in_bound_jk(nzm1, j, k) ;
00756             
00757           ligne += 2 ;
00758           //derivative of (x-1)^(l+2) at -1 :
00759           pari = 1 ;
00760           for (int i=0; i<nr; i++) {
00761               systeme.set(ligne, colonne) 
00762               -= 4*alpha*pari*i*i*sol_hom_deux_eta(nzm1, k, j, i) ;
00763               systeme.set(ligne, colonne+1) 
00764               -= 4*alpha*pari*i*i*sol_hom_quatre_eta(nzm1, k, j, i) ;
00765               sec_membre.set(ligne) 
00766               += 4*alpha*pari*i*i* sol_part_eta(nzm1, k, j, i) ;
00767               pari = -pari ;
00768           }
00769           ligne++ ;
00770           // ... and their counterparts for V^r
00771           pari = 1 ;
00772           for (int i=0; i<nr; i++) {
00773               systeme.set(ligne, colonne) 
00774               -= 4*alpha*pari*i*i*sol_hom_deux_vr(nzm1, k, j, i) ;
00775               systeme.set(ligne, colonne+1) 
00776               -= 4*alpha*pari*i*i*sol_hom_quatre_vr(nzm1, k, j, i) ;
00777               sec_membre.set(ligne) 
00778               += 4*alpha*pari*i*i* sol_part_vr(nzm1, k, j, i) ;
00779               pari = -pari ;
00780           }
00781             
00782           // Solution of the system giving the coefficients for the homogeneous 
00783           // solutions
00784           //-------------------------------------------------------------------
00785           systeme.set_lu() ;
00786           Tbl facteurs(systeme.inverse(sec_membre)) ;
00787           int conte = 0 ;
00788 
00789           // everything is put to the right place, the same combination of hom.
00790           // solutions (with some l or -(l+1) factors) must be used for V^r
00791           //-------------------------------------------------------------------
00792           nr = mg.get_nr(0) ; //nucleus
00793           for (int i=0 ; i<nr ; i++) {
00794         cf_eta.set(0, k, j, i) = 0.;
00795         cf_vr.set(0, k, j, i) = 0.;
00796           }
00797 
00798           for (int zone=1 ; zone<nzm1 ; zone++) { //shells
00799           nr = mg.get_nr(zone) ;
00800           for (int i=0 ; i<nr ; i++) {
00801               cf_eta.set(zone, k, j, i) = 
00802               sol_part_eta(zone, k, j, i)
00803               +facteurs(conte)*sol_hom_un_eta(zone, k, j, i) 
00804               +facteurs(conte+1)*sol_hom_deux_eta(zone, k, j, i) 
00805               +facteurs(conte+2)*sol_hom_trois_eta(zone, k, j, i) 
00806               +facteurs(conte+3)*sol_hom_quatre_eta(zone, k, j, i) ;
00807               cf_vr.set(zone, k, j, i) = sol_part_vr(zone, k, j, i)
00808               +facteurs(conte)*sol_hom_un_vr(zone, k, j, i) 
00809               +facteurs(conte+1)*sol_hom_deux_vr(zone, k, j, i) 
00810               +facteurs(conte+2)*sol_hom_trois_vr(zone, k, j, i) 
00811               +facteurs(conte+3)*sol_hom_quatre_vr(zone, k, j, i) ;
00812           }
00813           conte+=4 ;
00814           }
00815           nr = mg.get_nr(nz-1) ; //compactified external domain
00816           for (int i=0 ; i<nr ; i++) {
00817           cf_eta.set(nzm1, k, j, i) = sol_part_eta(nzm1, k, j, i)
00818               +facteurs(conte)*sol_hom_deux_eta(nzm1, k, j, i) 
00819               +facteurs(conte+1)*sol_hom_quatre_eta(nzm1, k, j, i) ;
00820           cf_vr.set(nzm1, k, j, i) = sol_part_vr(nzm1, k, j, i)
00821               +facteurs(conte)*sol_hom_deux_vr(nzm1, k, j, i) 
00822               +facteurs(conte+1)*sol_hom_quatre_vr(nzm1, k, j, i) ;
00823 
00824           }
00825       } // End of nullite_plm  
00826       } //End of loop on theta
00827   vr.set_spectral_va().ylm_i() ;
00828   vr += vrl0 ;
00829   het.set_spectral_va().ylm_i() ;
00830 
00831 
00832   
00833   Scalar pipo(*mpaff); 
00834   pipo.annule_hard(); 
00835   pipo.std_spectral_base(); 
00836   pipo.set_dzpuis(3);
00837   Param_elliptic param_mu(mu());
00838 
00839 
00840   
00841 
00842     resu.set_vr_eta_mu(vr, het, mu().sol_elliptic_boundary(param_mu, (*bound_mu), dir_mu , neum_mu)) ;
00843 
00844   return ;
00845   
00846 }

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