vector_df_poisson.C

00001 /*
00002  *  Resolution of the divergence-free vector Poisson equation
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 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_df_poisson_C[] = "$Header: /cvsroot/Lorene/C++/Source/Tensor/vector_df_poisson.C,v 1.13 2005/02/15 09:45:22 j_novak Exp $" ;
00029 
00030 /*
00031  * $Id: vector_df_poisson.C,v 1.13 2005/02/15 09:45:22 j_novak Exp $
00032  * $Log: vector_df_poisson.C,v $
00033  * Revision 1.13  2005/02/15 09:45:22  j_novak
00034  * Correction of an error in the matching.
00035  *
00036  * Revision 1.12  2005/02/09 16:53:11  j_novak
00037  * Now V^r and eta are matched across domains, but not any of their derivatives.
00038  *
00039  * Revision 1.11  2005/02/09 14:52:01  j_novak
00040  * Better solution in the shells.
00041  *
00042  * Revision 1.10  2005/02/09 13:20:27  j_novak
00043  * Completely new way of solving the vector Poisson equation in the Div_free
00044  * case: the system is inverted "as a whole" for V^r and eta. This works only
00045  * with Map_af...
00046  *
00047  *
00048  * $Header: /cvsroot/Lorene/C++/Source/Tensor/vector_df_poisson.C,v 1.13 2005/02/15 09:45:22 j_novak Exp $
00049  *
00050  */
00051 
00052 // C headers
00053 #include <assert.h>
00054 #include <math.h>
00055 
00056 // Lorene headers
00057 #include "tensor.h"
00058 #include "diff.h"
00059 #include "proto.h"
00060 #include "param.h"
00061 
00062 Vector_divfree Vector_divfree::poisson(Param& par ) const {
00063 
00064    // All this has a meaning only for spherical components:
00065 #ifndef NDEBUG 
00066     const Base_vect_spher* bvs = dynamic_cast<const Base_vect_spher*>(triad) ;
00067     assert(bvs != 0x0) ; 
00068 #endif
00069     
00070     int nitermax = par.get_int() ; 
00071     int& niter = par.get_int_mod() ; 
00072     double relax = par.get_double() ; 
00073     double precis = par.get_double(1) ;     
00074     Cmp& ss_khi = par.get_cmp_mod(0) ;
00075     Cmp& ss_mu = par.get_cmp_mod(1) ;
00076       
00077     // Solution for the r-component
00078     // ----------------------------
00079     
00080     Scalar source_r = *(cmp[0]) ; 
00081     source_r.mult_r() ; 
00082     
00083     Param par_khi ;
00084     par_khi.add_int(nitermax, 0) ;
00085     par_khi.add_int_mod(niter, 0) ;
00086     par_khi.add_double(relax, 0) ;
00087     par_khi.add_double(precis, 1) ;
00088     par_khi.add_cmp_mod(ss_khi, 0) ;
00089 
00090     Scalar khi (*mp) ;
00091     khi.set_etat_zero() ;
00092 
00093     source_r.poisson(par_khi, khi) ; 
00094     khi.div_r() ;   // khi now contains V^r
00095     
00096     // Solution for mu
00097     // ---------------
00098     
00099     Param par_mu ;
00100     par_mu.add_int(nitermax, 0) ;
00101     par_mu.add_int_mod(niter, 0) ;
00102     par_mu.add_double(relax, 0) ;
00103     par_mu.add_double(precis, 1) ;
00104     par_mu.add_cmp_mod(ss_mu, 0) ;
00105 
00106     Scalar mu_resu (*mp) ;
00107     mu_resu.set_etat_zero() ;
00108  
00109     mu().poisson(par_mu, mu_resu) ;
00110 
00111     // Final result
00112     // ------------
00113     
00114     Vector_divfree resu(*mp, *triad, *met_div) ; 
00115     
00116     resu.set_vr_mu(khi, mu_resu) ; 
00117     
00118     return resu ;
00119 
00120 }
00121 
00122 /*
00123  * In the case without parameters, first is solved the equation for mu and then
00124  * the system of equations for (eta, V^r) is inverted as a whole:
00125  * d2 eta / dr2 + 2/r d eta / dr - 1/r d V^r / dr = S(eta)
00126  * d V^r / dr + 2/r V^r - l(l+1)/r eta = 0 (div free condition)
00127  * 
00128  * There is no l=0 contribution (divergence free in all space!).
00129  * In the nucleus and the CED the system is inverted for h(r) and v(r) , 
00130  * such that eta = r^2 h and V^r = r^2 v in the nucleus,
00131  * in the compactified domain one has eta = u^2 h and V^r = u^2 v (where u=1/r); 
00132  * In the shells, both equations are  multiplied by r.
00133  * These methods are used only to get particular solutions.
00134  * 
00135  * Homogeneous solutions are known analitically: r^(l-1) and/or 1/r^(l+2)
00136  * It is then only possible to match eta and V^r, but not their derivatives,
00137  * due to the div-free condition.
00138  */
00139 Vector_divfree Vector_divfree::poisson() const {
00140 
00141     const Map_af* mpaff = dynamic_cast<const Map_af*>(mp) ;
00142 #ifndef NDEBUG 
00143     for (int i=0; i<3; i++)
00144     assert(cmp[i]->check_dzpuis(4)) ;
00145     // All this has a meaning only for spherical components:
00146     const Base_vect_spher* bvs = dynamic_cast<const Base_vect_spher*>(triad) ;
00147     assert(bvs != 0x0) ; 
00148     //## ... and affine mapping, for the moment!
00149     assert( mpaff != 0x0) ;
00150 #endif
00151 
00152   // Final result
00153   // ------------
00154     Vector_divfree resu(*mpaff, *triad, *met_div) ; 
00155 
00156   // Solution for mu
00157   // ---------------
00158     Scalar mu_resu = mu().poisson() ;
00159     
00160     Scalar f_r(*mpaff) ;
00161     if (cmp[0]->get_etat() == ETATZERO) { // no work needed ...
00162     f_r.set_etat_zero() ;
00163     resu.set_vr_mu(f_r, mu_resu) ; 
00164     return resu ;
00165     }
00166   
00167     // Some working objects
00168     //---------------------
00169   const Mg3d& mg = *mpaff->get_mg() ;
00170   int nz = mg.get_nzone() ; int nzm1 = nz - 1;
00171   assert(mg.get_type_r(0) == RARE) ;
00172   assert(mg.get_type_r(nzm1) == UNSURR) ;
00173   Scalar S_r = *cmp[0] ;
00174   Scalar S_eta = eta() ;
00175   S_r.set_spectral_va().ylm() ;
00176   S_eta.set_spectral_va().ylm() ;
00177   const Base_val& base = S_eta.get_spectral_va().base ;
00178   Mtbl_cf sol_part_eta(mg, base) ; sol_part_eta.annule_hard() ;
00179   Mtbl_cf sol_part_vr(mg, base) ; sol_part_vr.annule_hard() ;
00180   Mtbl_cf solution_hom_un(mg, base) ; solution_hom_un.annule_hard() ;
00181   Mtbl_cf solution_hom_deux(mg, base) ; solution_hom_deux.annule_hard() ;
00182 
00183   // Build-up & inversion of the system for (eta, V^r) in each domain
00184   //-----------------------------------------------------------------
00185 
00186   // Nucleus
00187   //--------
00188   int nr = mg.get_nr(0) ;
00189   int nt = mg.get_nt(0) ;
00190   int np = mg.get_np(0) ;
00191   double alpha = mpaff->get_alpha()[0] ;
00192   double beta = mpaff->get_beta()[0] ;
00193   int l_q = 0 ; int m_q = 0; int base_r = 0 ;
00194   int nr0 = nr - 1 ; //one degree of freedom less because of division by r^2
00195 
00196   // Loop on l and m
00197   //----------------
00198   for (int k=0 ; k<np+1 ; k++) {
00199       for (int j=0 ; j<nt ; j++) { 
00200       base.give_quant_numbers(0, k, j, m_q, l_q, base_r) ;
00201       if ( (nullite_plm(j, nt, k, np, base) == 1) && (l_q != 0) ) {
00202           int dege1 = (l_q <3 ? 0 : 1) ; //degeneracy of eq.1
00203           int dege2 = 0 ;  //degeneracy of eq.2
00204           int nr_eq1 = nr0 - dege1 ; //Eq.1 is for h (eta/r^2)
00205           int nr_eq2 = nr0 - dege2 ; //Eq.2 is the div-free condition
00206           int nrtot = nr_eq1 + nr_eq2 ;
00207           Matrice oper(nrtot, nrtot) ; oper.set_etat_qcq() ;
00208           Tbl sec_membre(nrtot) ; sec_membre.set_etat_qcq() ;
00209           Diff_x2dsdx2 d2(base_r, nr) ; const Matrice& md2 = d2.get_matrice() ;
00210           Diff_xdsdx xd(base_r, nr) ; const Matrice& mxd = xd.get_matrice() ;
00211           Diff_id id(base_r, nr) ; const Matrice& mid = id.get_matrice() ;
00212 
00213           // Building the operator
00214           //----------------------
00215           for (int lin=0; lin<nr_eq1; lin++) { //eq.1 
00216           for (int col=dege1; col<nr0; col++) 
00217               oper.set(lin,col-dege1) 
00218               = md2(lin,col) + 6*mxd(lin,col) + 6*mid(lin,col) ;
00219           for (int col=dege2; col<nr0; col++) 
00220               oper.set(lin,col-dege2+nr_eq1) = -mxd(lin,col) -2*mid(lin,col) ;
00221           }
00222           for (int lin=0; lin<nr0; lin++) { //eq.2
00223           for (int col=dege1; col<nr0; col++)
00224               oper.set(lin+nr_eq1,col-dege1) = -l_q*(l_q+1)*mid(lin,col) ;
00225           for (int col=dege2; col<nr0; col++)
00226               oper.set(lin+nr_eq1, col-dege2+nr_eq1) = mxd(lin,col) + 4*mid(lin,col) ;
00227           }
00228           oper.set_lu() ;
00229 
00230           // Filling the r.h.s
00231           //------------------
00232           for (int i=0; i<nr_eq1; i++)  //eq.1
00233           sec_membre.set(i) = (*S_eta.get_spectral_va().c_cf)(0, k, j, i) ;
00234           for (int i=0; i<nr0; i++) //eq.2
00235           sec_membre.set(i+nr_eq1) = 0 ;
00236 
00237           // Inversion of the "big" operator
00238           //--------------------------------
00239           Tbl big_res = oper.inverse(sec_membre) ;
00240           Tbl res_eta(nr) ;  res_eta.set_etat_qcq() ;
00241           Tbl res_vr(nr) ;  res_vr.set_etat_qcq() ;
00242           
00243           // Putting coefficients of h and v to individual arrays
00244           //-----------------------------------------------------
00245           for (int i=0; i<dege1; i++)
00246           res_eta.set(i) = 0 ;
00247           for (int i=dege1; i<nr0; i++)
00248           res_eta.set(i) = big_res(i-dege1) ;
00249           res_eta.set(nr0) = 0 ;
00250           for (int i=0; i<dege2; i++)
00251           res_vr.set(i) = 0 ;
00252           for (int i=dege2; i<nr0; i++)
00253           res_vr.set(i) = big_res(i-dege2+nr_eq1) ;
00254           res_vr.set(nr0) = 0 ;
00255 
00256           // Multiplication by xi^2 (the alpha^2 factor comes soon)
00257           //-------------------------------------------------------
00258           multx2_1d(nr, &res_eta.t, base_r) ;
00259           multx2_1d(nr, &res_vr.t, base_r) ;
00260 
00261           // Homogeneous solution (only r^(l-1) in the nucleus)
00262           Tbl sol_hom = solh(nr, l_q-1, 0., base_r) ;
00263           for (int i=0 ; i<nr ; i++) {
00264           sol_part_eta.set(0, k, j, i) = alpha*alpha*res_eta(i) ;
00265           sol_part_vr.set(0, k, j, i) = alpha*alpha*res_vr(i) ;
00266           solution_hom_un.set(0, k, j, i) = sol_hom(i) ;
00267           solution_hom_deux.set(0, k, j, i) = 0. ; 
00268           }
00269       }
00270       }
00271   }     
00272 
00273 
00274   // Shells
00275   //-------
00276   for (int zone=1 ; zone<nzm1 ; zone++) {
00277       nr = mg.get_nr(zone) ; 
00278       assert (nr > 5) ;
00279       assert(nt == mg.get_nt(zone)) ;
00280       assert(np == mg.get_np(zone)) ;
00281       alpha = mpaff->get_alpha()[zone] ;
00282       beta = mpaff->get_beta()[zone] ;
00283       double ech = beta / alpha ;
00284 
00285   // Loop on l and m
00286   //----------------
00287       for (int k=0 ; k<np+1 ; k++) {
00288       for (int j=0 ; j<nt ; j++) {
00289       base.give_quant_numbers(zone, k, j, m_q, l_q, base_r) ;
00290       if ( (nullite_plm(j, nt, k, np, base) == 1) && (l_q != 0) ) {
00291           int dege1 = 2 ; //degeneracy of eq.1
00292           int dege2 = 0 ;  //degeneracy of eq.2
00293           int nr_eq1 = nr - dege1 ; //Eq.1 is for eta
00294           int nr_eq2 = nr - dege2 ; //Eq.2 is the div-free condition
00295           int nrtot = nr_eq1 + nr_eq2 + 1;
00296           Matrice oper(nrtot, nrtot) ; oper.set_etat_qcq() ;
00297           Tbl sec_membre(nrtot) ; sec_membre.set_etat_qcq() ;
00298           Diff_x2dsdx2 x2d2(base_r, nr+1); const Matrice& m2d2 = x2d2.get_matrice() ;
00299           Diff_xdsdx2 xd2(base_r, nr+1) ; const Matrice& mxd2 = xd2.get_matrice() ;
00300           Diff_dsdx2 d2(base_r, nr+1) ; const Matrice& md2 = d2.get_matrice() ;
00301           Diff_xdsdx xd(base_r, nr+1) ; const Matrice& mxd = xd.get_matrice() ;
00302           Diff_dsdx d1(base_r, nr+1) ; const Matrice& md = d1.get_matrice() ;
00303           Diff_id id(base_r, nr+1) ; const Matrice& mid = id.get_matrice() ;
00304 
00305           // Building the operator
00306           //----------------------
00307           for (int lin=0; lin<nr_eq1; lin++) { 
00308           for (int col=dege1; col<nr; col++) 
00309               oper.set(lin,col-dege1) 
00310               = mxd2(lin,col) + ech*md2(lin,col) + 2*md(lin,col) ;
00311               for (int col=dege2; col<nr+1; col++) 
00312                   oper.set(lin,col-dege2+nr_eq1) = -md(lin,col) ; 
00313           }
00314           for (int lin=0; lin<nr_eq2; lin++) {
00315           for (int col=dege1; col<nr; col++)
00316               oper.set(lin+nr_eq1,col-dege1) = -l_q*(l_q+1)*mid(lin,col) ;
00317           for (int col=dege2; col<nr+1; col++)
00318               oper.set(lin+nr_eq1, col-dege2+nr_eq1) = 
00319               mxd(lin,col) + ech*md(lin,col) + 2*mid(lin,col) ;
00320           }
00321           //Additional line to avoid spurious homogeneous solutions
00322           //this line is the first one of the V^r eq.
00323           for (int col=dege1; col<nr; col++)
00324           oper.set(nrtot-1, col-dege1) = 0 ;
00325           for (int col=dege2; col<nr+1; col++)
00326           oper.set(nrtot-1, col-dege2+nr_eq1) = 
00327               m2d2(0,col) + ech*(2*mxd2(0,col) + ech*md2(0,col))
00328               +4*(mxd(0,col) +ech*md(0,col)) 
00329               +(2 - l_q*(l_q+1))*mid(0,col) ;
00330           oper.set_lu() ;
00331           
00332           // Filling the r.h.s
00333           //------------------
00334           Tbl sr(5) ; sr.set_etat_qcq() ; 
00335           Tbl seta(nr) ; seta.set_etat_qcq() ;
00336           for (int i=0; i<5; i++) {
00337           sr.set(i) = (*S_r.get_spectral_va().c_cf)(zone, k, j, i);
00338           seta.set(i) = (*S_eta.get_spectral_va().c_cf)(zone, k, j, i) ;
00339           }
00340           for (int i=5; i<nr; i++) 
00341           seta.set(i) = (*S_eta.get_spectral_va().c_cf)(zone, k, j, i) ;
00342           Tbl xsr= sr ;  Tbl x2sr= sr ;
00343           Tbl xseta= seta ;
00344           multx2_1d(5, &x2sr.t, base_r) ; multx_1d(5, &xsr.t, base_r) ;
00345           multx_1d(nr, &xseta.t, base_r) ;
00346           
00347           for (int i=0; i<nr_eq1; i++) 
00348           sec_membre.set(i) = alpha*(alpha*xseta(i) + beta*seta(i)) ;
00349           for (int i=0; i<nr_eq2; i++)
00350           sec_membre.set(i+nr_eq1) = 0 ;
00351           sec_membre.set(nr_eq1+nr_eq2) = alpha*alpha*x2sr(0) + 2*alpha*beta*xsr(0)
00352           + beta*beta*sr(0) ;
00353 
00354           // Inversion of the "big" operator
00355           //--------------------------------
00356           Tbl big_res = oper.inverse(sec_membre) ;
00357           Tbl res_eta(nr) ;  res_eta.set_etat_qcq() ;
00358           Tbl res_vr(nr) ;  res_vr.set_etat_qcq() ;
00359           
00360           // Putting coefficients of h and v to individual arrays
00361           //-----------------------------------------------------
00362           for (int i=0; i<dege1; i++)
00363           res_eta.set(i) = 0 ;
00364           for (int i=dege1; i<nr; i++)
00365           res_eta.set(i) = big_res(i-dege1) ;
00366           for (int i=0; i<dege2; i++)
00367           res_vr.set(i) = 0 ;
00368           for (int i=dege2; i<nr; i++)
00369           res_vr.set(i) = big_res(i-dege2+nr_eq1) ;
00370 
00371           //homogeneous solutions
00372           Tbl sol_hom1 = solh(nr, l_q-1, ech, base_r) ;
00373           Tbl sol_hom2 = solh(nr, l_q+1, ech, base_r) ;
00374           for (int i=0 ; i<nr ; i++) {
00375           sol_part_eta.set(zone, k, j, i) = res_eta(i) ;
00376           sol_part_vr.set(zone, k, j, i) = res_vr(i) ;
00377           solution_hom_un.set(zone, k, j, i) = sol_hom1(0,i) ;
00378           solution_hom_deux.set(zone, k, j, i) = sol_hom2(1,i) ;
00379           }
00380       } 
00381       }
00382       }
00383   }
00384 
00385   // Compactified external domain
00386   //-----------------------------
00387   nr = mg.get_nr(nzm1) ;
00388   assert(nt == mg.get_nt(nzm1)) ;
00389   assert(np == mg.get_np(nzm1)) ;
00390   alpha = mpaff->get_alpha()[nzm1] ;
00391   assert (nr > 4) ;
00392   nr0 = nr - 2 ;  //two degrees of freedom less because of division by r^2
00393 
00394   // Loop on l and m
00395   //----------------
00396   for (int k=0 ; k<np+1 ; k++) {
00397       for (int j=0 ; j<nt ; j++) { 
00398       base.give_quant_numbers(nzm1, k, j, m_q, l_q, base_r) ;
00399       if ( (nullite_plm(j, nt, k, np, base) == 1) && (l_q != 0) ) {
00400           int dege1 = 0; //degeneracy of eq.1
00401           int dege2 = 1; //degeneracy of eq.2
00402           int nr_eq1 = nr0 - dege1 ; //Eq.1 is for eta
00403           int nr_eq2 = nr0 - dege2 ; //Eq.2 is the div-free condition
00404           int nrtot = nr_eq1 + nr_eq2 ;
00405           Matrice oper(nrtot, nrtot) ; oper.set_etat_qcq() ;
00406           Tbl sec_membre(nrtot) ; sec_membre.set_etat_qcq() ;
00407           Diff_x2dsdx2 x2d2(base_r, nr) ; const Matrice& m2d2 = x2d2.get_matrice() ;
00408           Diff_xdsdx xd(base_r, nr) ; const Matrice& mxd = xd.get_matrice() ;
00409           Diff_id id(base_r, nr) ; const Matrice& mid = id.get_matrice() ;
00410 
00411           // Building the operator
00412           //----------------------
00413           for (int lin=0; lin<nr_eq1; lin++) { 
00414           for (int col=dege1; col<nr0; col++) 
00415               oper.set(lin,col-dege1) 
00416               = m2d2(lin,col) + 4*mxd(lin,col) + 2*mid(lin,col) ;
00417               for (int col=dege2; col<nr0; col++) 
00418                   oper.set(lin,col-dege2+nr_eq1) = 
00419                   mxd(lin,col) + 2*mid(lin,col) ;
00420           }
00421           for (int lin=0; lin<nr_eq2; lin++) {
00422           for (int col=dege1; col<nr0; col++)
00423               oper.set(lin+nr_eq1,col-dege1) = l_q*(l_q+1)*mid(lin,col) ;
00424           for (int col=dege2; col<nr0; col++)
00425               oper.set(lin+nr_eq1, col-dege2+nr_eq1) = mxd(lin,col) ;
00426           }
00427           oper.set_lu() ;
00428           
00429           // Filling the r.h.s
00430           //------------------
00431           for (int i=0; i<nr_eq1; i++) 
00432           sec_membre.set(i) = (*S_eta.get_spectral_va().c_cf)(nzm1, k, j, i) ;
00433           for (int i=0; i<nr_eq2; i++)
00434           sec_membre.set(i+nr_eq1) = 0 ;
00435           Tbl big_res = oper.inverse(sec_membre) ;
00436           Tbl res_eta(nr) ;  res_eta.set_etat_qcq() ;
00437           Tbl res_vr(nr) ;  res_vr.set_etat_qcq() ;
00438           
00439           // Putting coefficients of h and v to individual arrays
00440           //-----------------------------------------------------
00441           for (int i=0; i<dege1; i++)
00442           res_eta.set(i) = 0 ;
00443           for (int i=dege1; i<nr0; i++)
00444           res_eta.set(i) = big_res(i-dege1) ;
00445           res_eta.set(nr0) = 0 ;
00446           res_eta.set(nr0+1) = 0 ;
00447           for (int i=0; i<dege2; i++)
00448           res_vr.set(i) = 0 ;
00449           for (int i=dege2; i<nr0; i++)
00450           res_vr.set(i) = big_res(i-dege2+nr_eq1) ;
00451           res_vr.set(nr0) = 0 ;
00452           res_vr.set(nr0+1) = 0 ;
00453 
00454           // Multiplication by r^2 
00455           //-----------------------
00456           Tbl res_v2(nr) ;  res_v2.set_etat_qcq() ;
00457           Tbl res_e2(nr) ; res_e2.set_etat_qcq() ;
00458           mult2_xm1_1d_cheb(nr, res_eta.t, res_e2.t) ; 
00459           mult2_xm1_1d_cheb(nr, res_vr.t, res_v2.t) ; 
00460 
00461           // Homogeneous solution (only 1/r^(l+2) in the CED)
00462           Tbl sol_hom = solh(nr, l_q+1, 0., base_r) ;
00463           for (int i=0 ; i<nr ; i++) {
00464           sol_part_eta.set(nzm1, k, j, i) = alpha*alpha*res_e2(i) ;
00465           sol_part_vr.set(nzm1, k, j, i) = alpha*alpha*res_v2(i) ;
00466           solution_hom_un.set(nzm1, k, j, i) = sol_hom(i) ;
00467           solution_hom_deux.set(nzm1, k, j, i) = 0. ; 
00468           }
00469       }     
00470       }
00471   }
00472 
00473   // Now let's match everything ...
00474   //-------------------------------
00475 
00476   // Resulting V^r & eta
00477   Scalar vr(*mpaff) ; vr.set_etat_qcq() ;
00478   vr.set_spectral_base(base) ;
00479   vr.set_spectral_va().set_etat_cf_qcq() ;
00480   Mtbl_cf& cf_vr = *vr.set_spectral_va().c_cf ;
00481   cf_vr.annule_hard() ;
00482   Scalar het(*mpaff) ; het.set_etat_qcq() ;
00483   het.set_spectral_base(base) ;
00484   het.set_spectral_va().set_etat_cf_qcq() ;
00485   Mtbl_cf& cf_eta = *het.set_spectral_va().c_cf ;
00486   cf_eta.annule_hard() ;
00487   int taille = 2*nzm1 ;
00488   Tbl sec_membre(taille) ; 
00489   Matrice systeme(taille, taille) ; 
00490   systeme.set_etat_qcq() ;
00491   int ligne ;  int colonne ;
00492   
00493   // Loop on l and m
00494   //----------------
00495   for (int k=0 ; k<np+1 ; k++)
00496       for (int j=0 ; j<nt ; j++) {
00497       base.give_quant_numbers(0, k, j, m_q, l_q, base_r) ;
00498       if ((nullite_plm(j, nt, k, np, base) == 1)&&(l_q != 0)) {
00499         
00500           ligne = 0 ;
00501           colonne = 0 ;
00502           sec_membre.annule_hard() ;
00503           for (int l=0; l<taille; l++) 
00504           for (int c=0; c<taille; c++)
00505               systeme.set(l,c) = 0 ;
00506           //Nucleus 
00507           nr = mg.get_nr(0) ;
00508           alpha = mpaff->get_alpha()[0] ;
00509           // value of x^(l-1) at 1 ...
00510           systeme.set(ligne, colonne) = 1. ;
00511           for (int i=0 ; i<nr ; i++)
00512           sec_membre.set(ligne) -= sol_part_eta(0, k, j, i) ;
00513           ligne++ ;
00514           // ... and of its couterpart for V^r
00515           systeme.set(ligne, colonne) = l_q;
00516           for (int i=0; i<nr; i++)
00517           sec_membre.set(ligne) -= sol_part_vr(0,k,j,i) ; 
00518           colonne++ ; 
00519               //shells
00520           for (int zone=1 ; zone<nzm1 ; zone++) {
00521           nr = mg.get_nr(zone) ;
00522           alpha = mpaff->get_alpha()[zone] ;
00523           double echelle = mpaff->get_beta()[zone]/alpha ;
00524           ligne -- ;
00525           //value of (x+echelle)^(l-1) at -1 
00526           systeme.set(ligne, colonne) = -pow(echelle-1., double(l_q-1)) ;
00527           // value of 1/(x+echelle) ^(l+2) at -1 
00528           systeme.set(ligne, colonne+1) = -1/pow(echelle-1., double(l_q+2)) ;
00529           for (int i=0 ; i<nr ; i++)
00530               if (i%2 == 0)
00531               sec_membre.set(ligne) += sol_part_eta(zone, k, j, i) ;
00532               else sec_membre.set(ligne) -= sol_part_eta(zone, k, j, i) ;
00533           ligne++ ;
00534           // ... and their couterparts for V^r
00535           systeme.set(ligne, colonne) = -l_q*pow(echelle-1., double(l_q-1)) ;
00536           systeme.set(ligne, colonne+1) = (l_q+1)/pow(echelle-1., double(l_q+2));
00537           for (int i=0 ; i<nr ; i++)
00538               if (i%2 == 0)
00539               sec_membre.set(ligne) += sol_part_vr(zone, k, j, i) ;
00540               else sec_membre.set(ligne) -= sol_part_vr(zone, k, j, i) ;
00541           ligne++ ;
00542             
00543           //value of (x+echelle)^(l-1) at 1 :
00544           systeme.set(ligne, colonne) = pow(echelle+1., double(l_q-1)) ;
00545           // value of 1/(x+echelle)^(l+2) at 1 :
00546           systeme.set(ligne, colonne+1) = 1./pow(echelle+1., double(l_q+2)) ;
00547           for (int i=0 ; i<nr ; i++)
00548               sec_membre.set(ligne) -= sol_part_eta(zone, k, j, i) ;
00549           ligne ++ ;
00550           //... and their couterparts for V^r
00551           systeme.set(ligne, colonne) = l_q*pow(echelle+1., double(l_q-1)) ;
00552           systeme.set(ligne, colonne+1) = -double(l_q+1)
00553               / pow(echelle+1., double(l_q+2)) ;
00554           for (int i=0 ; i<nr ; i++)
00555               sec_membre.set(ligne) -= sol_part_vr(zone, k, j, i); 
00556           colonne += 2 ;
00557           }         
00558           //Compactified external domain
00559           nr = mg.get_nr(nzm1) ;
00560 
00561           alpha = mpaff->get_alpha()[nzm1] ;
00562           ligne -- ;
00563           //value of (x-1)^(l+2) at -1 :
00564           systeme.set(ligne, colonne) = -pow(-2, double(l_q+2)) ;
00565           for (int i=0 ; i<nr ; i++)
00566           if (i%2 == 0) sec_membre.set(ligne) += sol_part_eta(nzm1, k, j, i) ;
00567           else sec_membre.set(ligne) -= sol_part_eta(nzm1, k, j, i) ;
00568           //... and of its couterpart for V^r
00569           systeme.set(ligne+1, colonne) = double(l_q+1)*pow(-2, double(l_q+2)) ;
00570           for (int i=0 ; i<nr ; i++)
00571           if (i%2 == 0) sec_membre.set(ligne+1) += sol_part_vr(nzm1, k, j, i) ;
00572           else sec_membre.set(ligne+1) -= sol_part_vr(nzm1, k, j, i) ;
00573             
00574           // Solution of the system giving the coefficients for the homogeneous 
00575           // solutions
00576           //-------------------------------------------------------------------
00577           if (taille > 2) systeme.set_band(2,2) ;
00578           else systeme.set_band(1,1) ;
00579           systeme.set_lu() ;
00580           Tbl facteurs(systeme.inverse(sec_membre)) ;
00581           int conte = 0 ;
00582 
00583           // everything is put to the right place, the same combination of hom.
00584           // solutions (with some l or -(l+1) factors) must be used for V^r
00585           //-------------------------------------------------------------------
00586           nr = mg.get_nr(0) ; //nucleus
00587           for (int i=0 ; i<nr ; i++) {
00588           cf_eta.set(0, k, j, i) = sol_part_eta(0, k, j, i)
00589               +facteurs(conte)*solution_hom_un(0, k, j, i) ;
00590           cf_vr.set(0, k, j, i) = sol_part_vr(0, k, j, i)
00591               +double(l_q)*facteurs(conte)*solution_hom_un(0, k, j, i) ;
00592           }
00593           conte++ ;
00594           for (int zone=1 ; zone<nzm1 ; zone++) { //shells
00595           nr = mg.get_nr(zone) ;
00596           for (int i=0 ; i<nr ; i++) {
00597               cf_eta.set(zone, k, j, i) = 
00598               sol_part_eta(zone, k, j, i)
00599               +facteurs(conte)*solution_hom_un(zone, k, j, i) 
00600               +facteurs(conte+1)*solution_hom_deux(zone, k, j, i) ;
00601               cf_vr.set(zone, k, j, i) = sol_part_vr(zone, k, j, i)
00602               +double(l_q)*facteurs(conte)*solution_hom_un(zone, k, j, i) 
00603               -double(l_q+1)*facteurs(conte+1)*solution_hom_deux(zone, k, j, i) ;
00604           }
00605           conte+=2 ;
00606           }
00607           nr = mg.get_nr(nz-1) ; //compactified external domain
00608           for (int i=0 ; i<nr ; i++) {
00609           cf_eta.set(nzm1, k, j, i) = sol_part_eta(nzm1, k, j, i)
00610               +facteurs(conte)*solution_hom_un(nzm1, k, j, i) ;
00611           cf_vr.set(nzm1, k, j, i) = sol_part_vr(nzm1, k, j, i)
00612               -double(l_q+1)*facteurs(conte)*solution_hom_un(nzm1, k, j, i) ;
00613           }
00614       } // End of nullite_plm  
00615       } //End of loop on theta
00616 
00617   vr.set_spectral_va().ylm_i() ;
00618   het.set_spectral_va().ylm_i() ;
00619 
00620   resu.set_vr_eta_mu(vr, het, mu_resu) ;
00621 
00622   return resu ;
00623 
00624 }
00625 

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