vector_poisson_block.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_block_C[] = "$Header: /cvsroot/Lorene/C++/Source/Tensor/vector_poisson_block.C,v 1.6 2007/12/21 10:45:06 j_novak Exp $" ;
00029 
00030 /*
00031  * $Id: vector_poisson_block.C,v 1.6 2007/12/21 10:45:06 j_novak Exp $
00032  * $Log: vector_poisson_block.C,v $
00033  * Revision 1.6  2007/12/21 10:45:06  j_novak
00034  * Better treatment of spherical symmetry
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_block.C,v 1.6 2007/12/21 10:45:06 j_novak 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 
00069 void Vector::poisson_block(double lam, Vector& resu) const {
00070 
00071     const Map_af* mpaff = dynamic_cast<const Map_af*>(mp) ;
00072 #ifndef NDEBUG 
00073     for (int i=0; i<3; i++)
00074     assert(cmp[i]->check_dzpuis(4)) ;
00075     // All this has a meaning only for spherical components:
00076     const Base_vect_spher* bvs = dynamic_cast<const Base_vect_spher*>(triad) ;
00077     assert(bvs != 0x0) ; 
00078     //## ... and affine mapping, for the moment!
00079     assert( mpaff != 0x0) ;
00080 #endif
00081 
00082      if (fabs(lam + 1.) < 1.e-8) {
00083     const Metric_flat& mets = mp->flat_met_spher() ;
00084     Vector_divfree sou(*mp, *triad, mets) ;
00085     for (int i=1; i<=3; i++) sou.set(i) = *cmp[i-1] ;
00086     resu = sou.poisson() ;
00087     return ;
00088      }
00089 
00090     // Some working objects
00091     //---------------------
00092   const Mg3d& mg = *mpaff->get_mg() ;
00093   int nz = mg.get_nzone() ; int nzm1 = nz - 1;
00094   int nt = mg.get_nt(0) ;
00095   int np = mg.get_np(0) ;
00096   assert(mg.get_type_r(0) == RARE) ;
00097   assert(mg.get_type_r(nzm1) == UNSURR) ;
00098   Scalar S_r = *cmp[0] ;
00099   Scalar S_eta = eta() ;
00100   Scalar het(*mpaff) ; 
00101   Scalar vr(*mpaff) ; 
00102  bool all_zero = false ;
00103  if (S_r.get_etat() == ETATZERO) {
00104      if (S_eta.get_etat() == ETATZERO) {
00105      vr.set_etat_zero() ;
00106      het.set_etat_zero() ;
00107      all_zero = true ;
00108      }
00109      else {
00110      S_r.annule_hard() ;
00111      S_r.set_spectral_base(S_eta.get_spectral_base()) ;
00112      }
00113  }
00114  if ((S_eta.get_etat() == ETATZERO)&&(!all_zero)) {
00115      S_eta.annule_hard() ;
00116      S_eta.set_spectral_base(S_r.get_spectral_base()) ;
00117  }
00118  if (!all_zero) {
00119   S_r.set_spectral_va().ylm() ;
00120   S_eta.set_spectral_va().ylm() ;
00121   const Base_val& base = S_eta.get_spectral_va().base ;
00122   Mtbl_cf sol_part_eta(mg, base) ; sol_part_eta.annule_hard() ;
00123   Mtbl_cf sol_part_vr(mg, base) ; sol_part_vr.annule_hard() ;
00124   Mtbl_cf sol_hom_un_eta(mg, base) ; sol_hom_un_eta.annule_hard() ;
00125   Mtbl_cf sol_hom_deux_eta(mg, base) ; sol_hom_deux_eta.annule_hard() ;
00126   Mtbl_cf sol_hom_trois_eta(mg, base) ; sol_hom_trois_eta.annule_hard() ;
00127   Mtbl_cf sol_hom_quatre_eta(mg, base) ; sol_hom_quatre_eta.annule_hard() ;
00128   Mtbl_cf sol_hom_un_vr(mg, base) ; sol_hom_un_vr.annule_hard() ;
00129   Mtbl_cf sol_hom_deux_vr(mg, base) ; sol_hom_deux_vr.annule_hard() ;
00130   Mtbl_cf sol_hom_trois_vr(mg, base) ; sol_hom_trois_vr.annule_hard() ;
00131   Mtbl_cf sol_hom_quatre_vr(mg, base) ; sol_hom_quatre_vr.annule_hard() ;
00132  
00133   // Solution of the l=0 part for Vr
00134   //---------------------------------
00135   Scalar sou_l0 = (*cmp[0]) / (1. + lam) ;
00136   sou_l0.set_spectral_va().ylm() ;
00137   Scalar vrl0 = pois_vect_r0(sou_l0) ;
00138 
00139   // Build-up & inversion of the system for (eta, V^r) in each domain
00140   //-----------------------------------------------------------------
00141 
00142   // Nucleus
00143   //--------
00144   int nr = mg.get_nr(0) ;
00145   double alpha = mpaff->get_alpha()[0] ; double alp2 = alpha*alpha ;
00146   double beta = mpaff->get_beta()[0] ;
00147   int l_q = 0 ; int m_q = 0; int base_r = 0 ;
00148 
00149   // Loop on l and m
00150   //----------------
00151   for (int k=0 ; k<np+1 ; k++) {
00152       for (int j=0 ; j<nt ; j++) { 
00153       base.give_quant_numbers(0, k, j, m_q, l_q, base_r) ;
00154       if ( (nullite_plm(j, nt, k, np, base) == 1) && (l_q > 0) ) {
00155           int aa = 0 ; int bb = 0 ;  int nr0 = 0 ;
00156           if (base_r == R_CHEBP) {
00157           nr0 = nr - 1 ;
00158           aa = 0 ; bb = 1 ;
00159           }
00160           else {
00161           assert (base_r == R_CHEBI) ;
00162           nr0 = nr - 2 ;
00163           aa = 2 ; bb = 1 ;
00164           }
00165           int d0 = nr - nr0 ;
00166           int nrtot = 2*nr0 ;
00167           Matrice oper(2*nr, 2*nr) ; oper.set_etat_qcq() ;
00168           Tbl sec_membre(nrtot) ; sec_membre.set_etat_qcq() ;
00169           Diff_dsdx2 d2(base_r, nr) ; const Matrice& md2 = d2.get_matrice() ;
00170           Diff_sxdsdx xd(base_r, nr) ; const Matrice& mxd = xd.get_matrice() ;
00171           Diff_sx2 s2(base_r, nr) ; const Matrice& ms2 = s2.get_matrice() ;
00172 
00173           // Building the operator
00174           //----------------------
00175           for (int lin=0; lin<nr; lin++) { //eq.1 
00176           for (int col=0; col<nr; col++) 
00177               oper.set(lin,col) = 
00178               (md2(lin,col) + 2*mxd(lin,col) 
00179                -(lam+1)*l_q*(l_q+1)*ms2(lin,col)) / alp2 ;
00180           for (int col=0; col<nr; col++) 
00181               oper.set(lin,col+nr) = 
00182               (lam*mxd(lin,col) + 2*(1+lam)*ms2(lin,col)) / alp2 ;
00183           }
00184           for (int lin=0; lin<nr; lin++) { //eq.2
00185           for (int col=0; col<nr; col++)
00186               oper.set(lin+nr,col) = 
00187               (-lam*l_q*(l_q+1)*mxd(lin,col) 
00188                +(lam+2)*l_q*(l_q+1)*ms2(lin,col)) / alp2 ;
00189           for (int col=0; col<nr; col++)
00190               oper.set(lin+nr, col+nr) = 
00191               ((lam+1)*(md2(lin,col) + 2*mxd(lin,col)) 
00192                - (2*(lam+1) + l_q*(l_q+1))*ms2(lin,col)) / alp2 ;
00193           }
00194           bool pb_eta = ( ( fabs( lam*double(l_q+3) + 2 ) < 0.01) && (l_q <=2) ) ;
00195           if (!pb_eta) {
00196           for (int col=0; col<nr; col++) 
00197               oper.set(nr0-1, col) = 1 ;
00198           for (int col=0; col<nr; col++) 
00199               oper.set(nr0-1,nr+col) = 0 ;        
00200           }
00201           if ((l_q > 2)||pb_eta) {
00202           for (int col=0; col<nr; col++) 
00203               oper.set(nr+nr0-1, col) = 0 ;
00204           for (int col=0; col<nr; col++) 
00205               oper.set(nr+nr0-1,nr+col) = 1 ;         
00206           }
00207 
00208           Matrice op2(nrtot, nrtot) ;
00209           op2.set_etat_qcq() ;
00210           for (int i=0; i<nr0; i++) {
00211           for (int col=0; col<nr0;col++) 
00212               op2.set(i,col) = (aa*col+bb)*oper(i,col+1) 
00213               + (aa*(col+1)+bb)*oper(i,col) ;
00214           for (int col=0; col<nr0;col++) 
00215               op2.set(i,col+nr0) = (aa*col+bb)*oper(i,col+nr+1)
00216               + (aa*(col+1)+bb)*oper(i,col+nr) ;
00217           }
00218           
00219           for (int i=nr0; i<nrtot; i++) {
00220           for (int col=0; col<nr0;col++) 
00221               op2.set(i,col) = (aa*col+bb)*oper(i+d0,col+1) 
00222               + (aa*(col+1)+bb)*oper(i+d0,col) ;
00223           for (int col=0; col<nr0;col++) 
00224               op2.set(i,col+nr0) = (aa*col+bb)*oper(i+d0,col+nr+1)
00225               + (aa*(col+1)+bb)*oper(i+d0,col+nr) ;
00226           }
00227           op2.set_lu() ;
00228 
00229           // Filling the r.h.s
00230           //------------------
00231           for (int i=0; i<nr0; i++)  //eq.1
00232           sec_membre.set(i) = (*S_eta.get_spectral_va().c_cf)(0, k, j, i) ;
00233           if (!pb_eta) sec_membre.set(nr0-1) = 0 ;
00234           for (int i=0; i<nr0; i++) //eq.2
00235           sec_membre.set(i+nr0) 
00236               = (*S_r.get_spectral_va().c_cf)(0, k, j, i) ;
00237           if ((l_q > 2)||pb_eta) sec_membre.set(nrtot-1) = 0 ;
00238 
00239           // Inversion of the "big" operator
00240           //--------------------------------
00241           Tbl big_res = op2.inverse(sec_membre) ;
00242           
00243           // Putting coefficients of eta and Vr to individual arrays
00244           //--------------------------------------------------------
00245           sol_part_eta.set(0, k, j, 0) = (aa+bb)*big_res(0) ;
00246           for (int i=1; i<nr0; i++)
00247           sol_part_eta.set(0, k, j, i) = (aa*(i-1)+bb)*big_res(i-1) 
00248               + (aa*(i+1)+bb)*big_res(i);
00249           sol_part_eta.set(0, k, j, nr0) = big_res(nr0-1)*(aa*(nr0-1) + bb) ;
00250           sol_part_vr.set(0, k, j, 0) = (aa+bb)*big_res(nr0) ;
00251           for (int i=1; i<nr0; i++)
00252           sol_part_vr.set(0, k, j, i) = (aa*(i-1)+bb)*big_res(nr0+i-1) 
00253               + (aa*(i+1)+bb)*big_res(nr0+i);
00254           sol_part_vr.set(0, k, j, nr0) = big_res(nrtot-1)*(aa*(nr0-1)+bb) ;
00255           if (base_r == R_CHEBI) {
00256           sol_part_eta.set(0, k, j, nr-1) = 0 ;
00257           sol_part_vr.set(0, k, j, nr-1) = 0 ;
00258           }
00259 
00260           // Homogeneous solutions (only r^(l-1) and r^(l+1) in the nucleus)
00261           if (l_q<=2) {
00262           double fac_eta = lam*double(l_q+3) + 2. ;
00263           double fac_vr = double(l_q+1)*(lam*l_q - 2.) ;
00264           Tbl sol_hom1 = solh(nr, l_q-1, 0., base_r) ;
00265           Tbl sol_hom2 = solh(nr, l_q+1, 0., base_r) ;
00266           for (int i=0 ; i<nr ; i++) {
00267               sol_hom_un_eta.set(0, k, j, i) = sol_hom1(i) ;
00268               sol_hom_un_vr.set(0, k, j, i) = l_q*sol_hom1(i) ;
00269               sol_hom_trois_eta.set(0, k, j, i) = fac_eta*sol_hom2(i) ; 
00270               sol_hom_trois_vr.set(0, k, j, i) = fac_vr*sol_hom2(i) ; 
00271           }
00272           }
00273           else {
00274           for (int i=0; i<nrtot; i++) sec_membre.set(i) = 0 ;
00275           // First homogeneous sol.
00276           sec_membre.set(nr0-1) = 1 ;
00277           big_res = op2.inverse(sec_membre) ;
00278 
00279           sol_hom_un_eta.set(0, k, j, 0) = (aa+bb)*big_res(0) ;
00280           for (int i=1; i<nr0; i++)
00281               sol_hom_un_eta.set(0, k, j, i) = (aa*(i-1)+bb)*big_res(i-1) 
00282               + (aa*(i+1)+bb)*big_res(i);
00283           sol_hom_un_eta.set(0, k, j, nr0) = big_res(nr0-1)*(aa*(nr0-1) + bb) ;
00284           sol_hom_un_vr.set(0, k, j, 0) = (aa+bb)*big_res(nr0) ;
00285           for (int i=1; i<nr0; i++)
00286               sol_hom_un_vr.set(0, k, j, i) = (aa*(i-1)+bb)*big_res(nr0+i-1) 
00287               + (aa*(i+1)+bb)*big_res(nr0+i);
00288           sol_hom_un_vr.set(0, k, j, nr0) = big_res(nrtot-1)*(aa*(nr0-1)+bb) ;
00289           if (base_r == R_CHEBI) {
00290               sol_hom_un_eta.set(0, k, j, nr-1) = 0 ;
00291               sol_hom_un_vr.set(0, k, j, nr-1) = 0 ;
00292           }
00293 
00294           sec_membre.set(nr0-1) = 0 ;
00295           sec_membre.set(nrtot-1) = 1 ;
00296           big_res = op2.inverse(sec_membre) ;
00297 
00298           sol_hom_trois_eta.set(0, k, j, 0) = (aa+bb)*big_res(0) ;
00299           for (int i=1; i<nr0; i++)
00300               sol_hom_trois_eta.set(0, k, j, i) = (aa*(i-1)+bb)*big_res(i-1) 
00301               + (aa*(i+1)+bb)*big_res(i);
00302           sol_hom_trois_eta.set(0, k, j, nr0) = big_res(nr0-1)*(aa*(nr0-1)+bb) ;
00303           sol_hom_trois_vr.set(0, k, j, 0) = (aa+bb)*big_res(nr0) ;
00304           for (int i=1; i<nr0; i++)
00305               sol_hom_trois_vr.set(0, k, j, i) = (aa*(i-1)+bb)*big_res(nr0+i-1) 
00306               + (aa*(i+1)+bb)*big_res(nr0+i);
00307           sol_hom_trois_vr.set(0, k, j, nr0) = big_res(nrtot-1)*(aa*(nr0-1)+bb) ;
00308           if (base_r == R_CHEBI) {
00309               sol_hom_trois_eta.set(0, k, j, nr-1) = 0 ;
00310               sol_hom_trois_vr.set(0, k, j, nr-1) = 0 ;
00311           }
00312           
00313           }
00314       }
00315       }
00316   }    
00317 
00318 
00319   // Shells
00320   //-------
00321   for (int zone=1 ; zone<nzm1 ; zone++) {
00322       nr = mg.get_nr(zone) ; 
00323       assert (nr > 5) ;
00324       assert(nt == mg.get_nt(zone)) ;
00325       assert(np == mg.get_np(zone)) ;
00326       alpha = mpaff->get_alpha()[zone] ;
00327       beta = mpaff->get_beta()[zone] ;
00328       double ech = beta / alpha ;
00329 
00330   // Loop on l and m
00331   //----------------
00332       for (int k=0 ; k<np+1 ; k++) {
00333       for (int j=0 ; j<nt ; j++) {
00334       base.give_quant_numbers(zone, k, j, m_q, l_q, base_r) ;
00335       if ( (nullite_plm(j, nt, k, np, base) == 1) && (l_q != 0) ) {
00336           int dege1 = 2 ; //degeneracy of eq.1
00337           int dege2 = 2 ;  //degeneracy of eq.2
00338           int nr_eq1 = nr - dege1 ; //Eq.1 is for eta
00339           int nr_eq2 = nr - dege2 ; //Eq.2 is for V^r
00340           int nrtot = 2*nr ;
00341           Matrice oper(nrtot, nrtot) ; oper.set_etat_qcq() ;
00342           Tbl sec_membre(nrtot) ; sec_membre.set_etat_qcq() ;
00343           Diff_x2dsdx2 x2d2(base_r, nr); const Matrice& m2d2 = x2d2.get_matrice() ;
00344           Diff_xdsdx2 xd2(base_r, nr) ; const Matrice& mxd2 = xd2.get_matrice() ;
00345           Diff_dsdx2 d2(base_r, nr) ; const Matrice& md2 = d2.get_matrice() ;
00346           Diff_xdsdx xd(base_r, nr) ; const Matrice& mxd = xd.get_matrice() ;
00347           Diff_dsdx d1(base_r, nr) ; const Matrice& md = d1.get_matrice() ;
00348           Diff_id id(base_r, nr) ; const Matrice& mid = id.get_matrice() ;
00349 
00350           // Building the operator
00351           //----------------------
00352           for (int lin=0; lin<nr_eq1; lin++) { 
00353           for (int col=0; col<nr; col++) 
00354               oper.set(lin,col) 
00355               = m2d2(lin,col) + 2*ech*mxd2(lin,col) + ech*ech*md2(lin,col) 
00356               + 2*(mxd(lin,col) + ech*md(lin,col)) 
00357               - (lam+1)*l_q*(l_q+1)*mid(lin,col) ;
00358           for (int col=0; col<nr; col++) 
00359               oper.set(lin,col+nr) 
00360               = lam*(mxd(lin,col) + ech*md(lin,col)) 
00361               + 2*(1+lam)*mid(lin,col) ; 
00362           }
00363           oper.set(nr_eq1, 0) = 1 ;
00364           for (int col=1; col<2*nr; col++) 
00365           oper.set(nr_eq1, col) = 0 ;
00366           oper.set(nr_eq1+1, 0) = 0 ;
00367           oper.set(nr_eq1+1, 1) = 1 ;
00368           for (int col=2; col<2*nr; col++) 
00369           oper.set(nr_eq1+1, 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) 
00373               = -l_q*(l_q+1)*(lam*(mxd(lin,col) + ech*md(lin,col))
00374                       - (lam+2)*mid(lin,col)) ;
00375           for (int col=0; col<nr; col++)
00376               oper.set(lin+nr, col+nr) 
00377               = (lam+1)*(m2d2(lin,col) + 2*ech*mxd2(lin,col) 
00378                      + ech*ech*md2(lin,col) 
00379                      + 2*(mxd(lin,col) + ech*md(lin,col)))
00380               -(2*(lam+1)+l_q*(l_q+1))*mid(lin,col) ;         
00381           }
00382           for (int col=0; col<nr; col++) 
00383           oper.set(nr+nr_eq2, col) = 0 ;
00384           oper.set(nr+nr_eq2, nr) = 1 ;
00385           for (int col=nr+1; col<2*nr; col++) 
00386           oper.set(nr+nr_eq2, col) = 0 ;
00387           for (int col=0; col<nr+1; col++) 
00388           oper.set(nr+nr_eq2+1, col) = 0 ;
00389           oper.set(nr+nr_eq2+1, nr+1) = 1 ;
00390           for (int col=nr+2; col<2*nr; col++) 
00391           oper.set(nr+nr_eq2+1, col) = 0 ;
00392 
00393           oper.set_lu() ;
00394           
00395           // Filling the r.h.s
00396           //------------------
00397           Tbl sr(nr) ; sr.set_etat_qcq() ; 
00398           Tbl seta(nr) ; seta.set_etat_qcq() ;
00399           for (int i=0; i<nr; i++) {
00400           sr.set(i) = (*S_r.get_spectral_va().c_cf)(zone, k, j, i);
00401           seta.set(i) = (*S_eta.get_spectral_va().c_cf)(zone, k, j, i) ;
00402           }
00403           Tbl xsr= sr ;  Tbl x2sr= sr ;
00404           Tbl xseta= seta ; Tbl x2seta = seta ;
00405           multx2_1d(nr, &x2sr.t, base_r) ; multx_1d(nr, &xsr.t, base_r) ;
00406           multx2_1d(nr, &x2seta.t, base_r) ; multx_1d(nr, &xseta.t, base_r) ;
00407           
00408           for (int i=0; i<nr_eq1; i++) 
00409           sec_membre.set(i) = alpha*(alpha*x2seta(i) + 2*beta*xseta(i))
00410               + beta*beta*seta(i);
00411           sec_membre.set(nr_eq1) = 0 ;
00412           sec_membre.set(nr_eq1+1) = 0 ;
00413           for (int i=0; i<nr_eq2; i++)
00414           sec_membre.set(i+nr) = beta*beta*sr(i) 
00415               + alpha*(alpha*x2sr(i) + 2*beta*xsr(i)) ;
00416           sec_membre.set(nr+nr_eq2) = 0 ;
00417           sec_membre.set(nr+nr_eq2+1) = 0 ;
00418 
00419           // Inversion of the "big" operator
00420           //--------------------------------
00421           Tbl big_res = oper.inverse(sec_membre) ;
00422           for (int i=0; i<nr; i++) {
00423           sol_part_eta.set(zone, k, j, i) = big_res(i) ;
00424           sol_part_vr.set(zone, k, j, i) = big_res(nr+i) ;        
00425           }
00426           
00427           // Getting the homogeneous solutions
00428           //----------------------------------
00429           sec_membre.annule_hard() ;
00430           sec_membre.set(nr_eq1) = 1 ;
00431           big_res = oper.inverse(sec_membre) ;
00432           for (int i=0 ; i<nr ; i++) {
00433           sol_hom_un_eta.set(zone, k, j, i) = big_res(i) ;
00434           sol_hom_un_vr.set(zone, k, j, i) = big_res(nr+i) ;
00435           }
00436           sec_membre.set(nr_eq1) = 0 ;
00437           sec_membre.set(nr_eq1+1) = 1 ;
00438           big_res = oper.inverse(sec_membre) ;
00439           for (int i=0 ; i<nr ; i++) {
00440           sol_hom_deux_eta.set(zone, k, j, i) = big_res(i) ;
00441           sol_hom_deux_vr.set(zone, k, j, i) = big_res(nr+i) ;
00442           }
00443           sec_membre.set(nr_eq1+1) = 0 ;
00444           sec_membre.set(nr+nr_eq2) = 1 ;
00445           big_res = oper.inverse(sec_membre) ;
00446           for (int i=0 ; i<nr ; i++) {
00447           sol_hom_trois_eta.set(zone, k, j, i) = big_res(i) ;
00448           sol_hom_trois_vr.set(zone, k, j, i) = big_res(nr+i) ;
00449           }
00450           sec_membre.set(nr+nr_eq2) = 0 ;
00451           sec_membre.set(nr+nr_eq2+1) = 1 ;
00452           big_res = oper.inverse(sec_membre) ;
00453           for (int i=0 ; i<nr ; i++) {
00454           sol_hom_quatre_eta.set(zone, k, j, i) = big_res(i) ;
00455           sol_hom_quatre_vr.set(zone, k, j, i) = big_res(nr+i) ;
00456           }
00457 
00458       } 
00459       }
00460       }
00461   }
00462 
00463   // Compactified external domain
00464   //-----------------------------
00465   nr = mg.get_nr(nzm1) ;
00466   assert(nt == mg.get_nt(nzm1)) ;
00467   assert(np == mg.get_np(nzm1)) ;
00468   alpha = mpaff->get_alpha()[nzm1] ; alp2 = alpha*alpha ;
00469   assert (nr > 4) ;
00470   
00471   // Loop on l and m
00472   //----------------
00473   for (int k=0 ; k<np+1 ; k++) {
00474       for (int j=0 ; j<nt ; j++) { 
00475       base.give_quant_numbers(nzm1, k, j, m_q, l_q, base_r) ;
00476       if ( (nullite_plm(j, nt, k, np, base) == 1) && (l_q != 0) ) {
00477           int dege1 = 3; //degeneracy of eq.1
00478           int dege2 = (l_q == 1 ? 2 : 3); //degeneracy of eq.2
00479           int nr_eq1 = nr - dege1 ; //Eq.1 is for eta
00480           int nr_eq2 = nr - dege2 ; //Eq.2 is the div-free condition
00481           int nrtot = 2*nr ;
00482           Matrice oper(nrtot, nrtot) ; oper.set_etat_qcq() ;
00483           Tbl sec_membre(nrtot) ; sec_membre.set_etat_qcq() ;
00484           Diff_dsdx2 d2(base_r, nr) ; const Matrice& md2 = d2.get_matrice() ;
00485           Diff_sxdsdx sxd(base_r, nr) ; const Matrice& mxd = sxd.get_matrice() ;
00486           Diff_sx2 sx2(base_r, nr) ; const Matrice& ms2 = sx2.get_matrice() ;
00487 
00488           // Building the operator
00489           //----------------------
00490           for (int lin=0; lin<nr_eq1; lin++) { 
00491           for (int col=0; col<nr; col++) 
00492               oper.set(lin,col) 
00493               = (md2(lin,col) - (lam+1)*l_q*(l_q+1)*ms2(lin,col))/alp2 ;
00494               for (int col=0; col<nr; col++) 
00495                   oper.set(lin,col+nr) = 
00496                   (-lam*mxd(lin,col) + 2*(1+lam)*ms2(lin,col)) / alp2 ;
00497           }
00498           for (int col=0; col<nr; col++)
00499           oper.set(nr_eq1, col) = 1 ;
00500           for (int col=nr; col<nrtot; col++)
00501           oper.set(nr_eq1, col) = 0 ;
00502           int pari = -1 ;
00503           for (int col=0; col<nr; col++) {
00504           oper.set(nr_eq1+1, col) = pari*col*col ;
00505           pari = pari ;
00506           }
00507           for (int col=nr; col<nrtot; col++)
00508           oper.set(nr_eq1+1, col) = 0 ;
00509           oper.set(nr_eq1+2, 0) = 1 ;
00510           for (int col=1; col<nrtot; col++)
00511           oper.set(nr_eq1+2, col) = 0 ;
00512           for (int lin=0; lin<nr_eq2; lin++) {
00513           for (int col=0; col<nr; col++)
00514               oper.set(lin+nr,col) = (l_q*(l_q+1)*(lam*mxd(lin,col) 
00515                    + (lam+2)*ms2(lin,col))) / alp2 ;
00516           for (int col=0; col<nr; col++)
00517               oper.set(lin+nr, col+nr) = ((lam+1)*md2(lin,col) 
00518                  - (2*(lam+1) + l_q*(l_q+1))*ms2(lin,col)) / alp2 ;
00519           }
00520           for (int col=0; col<nr; col++)
00521           oper.set(nr+nr_eq2, col) = 0 ;
00522           for (int col=nr; col<nrtot; col++)
00523           oper.set(nr+nr_eq2, col) = 1 ;
00524           for (int col=0; col<nr; col++)
00525           oper.set(nr+nr_eq2+1, col) = 0 ;
00526           pari = -1 ;
00527           for (int col=0; col<nr; col++) {
00528           oper.set(nr+nr_eq2+1, nr+col) = pari*col*col ;
00529           pari = pari ;
00530           }
00531           if (l_q > 1) {
00532           for (int col=0; col<nr; col++)
00533               oper.set(nr+nr_eq2+2, col) = 0 ;
00534           oper.set(nr+nr_eq2+2, nr) = 1 ;
00535           for (int col=nr+1; col<nrtot; col++)
00536               oper.set(nr+nr_eq2+2, col) = 0 ;
00537           }
00538           oper.set_lu() ;
00539 
00540           // Filling the r.h.s
00541           //------------------
00542           for (int i=0; i<nr_eq1; i++) 
00543           sec_membre.set(i) = (*S_eta.get_spectral_va().c_cf)(nzm1, k, j, i) ;
00544           for (int i=nr_eq1; i<nr; i++)
00545           sec_membre.set(i) = 0 ;
00546           for (int i=0; i<nr_eq2; i++)
00547           sec_membre.set(i+nr) =(*S_r.get_spectral_va().c_cf)(nzm1, k, j, i);
00548           for (int i=nr_eq2; i<nr; i++)
00549           sec_membre.set(nr+i) = 0 ;
00550           Tbl big_res = oper.inverse(sec_membre) ;
00551           
00552           // Putting coefficients of h and v to individual arrays
00553           //-----------------------------------------------------
00554           for (int i=0; i<nr; i++) {
00555           sol_part_eta.set(nzm1, k, j, i) = big_res(i) ;
00556           sol_part_vr.set(nzm1, k, j, i) = big_res(i+nr) ;
00557           }
00558 
00559           // Homogeneous solutions (only 1/r^(l+2) and 1/r^l in the CED)
00560           //------------------------------------------------------------
00561           if (l_q == 1) {         
00562           Tbl sol_hom1 = solh(nr, 0, 0., base_r) ;
00563           Tbl sol_hom2 = solh(nr, 2, 0., base_r) ;
00564           double fac_eta = lam + 2. ;
00565           double fac_vr = 2*lam + 2. ;
00566           for (int i=0 ; i<nr ; i++) {
00567               sol_hom_deux_eta.set(nzm1, k, j, i) = sol_hom2(i) ;
00568               sol_hom_quatre_eta.set(nzm1, k, j, i) = fac_eta*sol_hom1(i) ;
00569               sol_hom_deux_vr.set(nzm1, k, j, i) = -2*sol_hom2(i) ;
00570               sol_hom_quatre_vr.set(nzm1, k, j, i) = fac_vr*sol_hom1(i) ;
00571           }
00572           }
00573           else {
00574           sec_membre.annule_hard() ;
00575           sec_membre.set(nr-1) = 1 ;
00576           big_res = oper.inverse(sec_membre) ;
00577 
00578           for (int i=0; i<nr; i++) {
00579               sol_hom_deux_eta.set(nzm1, k, j, i) = big_res(i) ;
00580               sol_hom_deux_vr.set(nzm1, k, j, i) = big_res(nr+i) ;
00581           }
00582           sec_membre.set(nr-1) = 0 ;
00583           sec_membre.set(2*nr-1) = 1 ;
00584           big_res = oper.inverse(sec_membre) ;
00585           for (int i=0; i<nr; i++) {
00586               sol_hom_quatre_eta.set(nzm1, k, j, i) = big_res(i) ;
00587               sol_hom_quatre_vr.set(nzm1, k, j, i) = big_res(nr+i) ;
00588           }
00589           }
00590       } 
00591       }
00592   }
00593 
00594   // Now let's match everything ...
00595   //-------------------------------
00596 
00597   // Resulting V^r & eta
00598   vr.set_etat_qcq() ;
00599   vr.set_spectral_base(base) ;
00600   vr.set_spectral_va().set_etat_cf_qcq() ;
00601   Mtbl_cf& cf_vr = *vr.set_spectral_va().c_cf ;
00602   cf_vr.annule_hard() ;
00603   het.set_etat_qcq() ;
00604   het.set_spectral_base(base) ;
00605   het.set_spectral_va().set_etat_cf_qcq() ;
00606   Mtbl_cf& cf_eta = *het.set_spectral_va().c_cf ;
00607   cf_eta.annule_hard() ;
00608   int taille = 4*nzm1 ;
00609   Tbl sec_membre(taille) ; 
00610   Matrice systeme(taille, taille) ; 
00611   systeme.set_etat_qcq() ;
00612   int ligne ;  int colonne ;
00613   
00614   // Loop on l and m
00615   //----------------
00616   for (int k=0 ; k<np+1 ; k++)
00617       for (int j=0 ; j<nt ; j++) {
00618       base.give_quant_numbers(0, k, j, m_q, l_q, base_r) ;
00619       if ((nullite_plm(j, nt, k, np, base) == 1)&&(l_q != 0)) {
00620         
00621           ligne = 0 ;
00622           colonne = 0 ;
00623           systeme.annule_hard() ;
00624           sec_membre.annule_hard() ;
00625 
00626           //Nucleus 
00627           nr = mg.get_nr(0) ;
00628           alpha = mpaff->get_alpha()[0] ;
00629           // value of at x=1 of eta ...
00630           systeme.set(ligne, colonne) = sol_hom_un_eta.val_out_bound_jk(0, j, k) ;
00631           systeme.set(ligne, colonne+1) = sol_hom_trois_eta.val_out_bound_jk(0, j, k) ;
00632           sec_membre.set(ligne) = -sol_part_eta.val_out_bound_jk(0, j, k) ;
00633           ligne++ ;
00634           // ... and of its couterpart for V^r
00635           systeme.set(ligne, colonne) = sol_hom_un_vr.val_out_bound_jk(0, j, k) ;
00636           systeme.set(ligne, colonne+1) = sol_hom_trois_vr.val_out_bound_jk(0, j, k) ;
00637           sec_membre.set(ligne) = -sol_part_vr.val_out_bound_jk(0,j,k) ; 
00638           ligne++ ; 
00639 
00640               //derivatives
00641           int pari = (base_r == R_CHEBP ? 0 : 1) ;
00642           for (int i=0; i<nr; i++) {
00643           systeme.set(ligne, colonne) 
00644               += (2*i+pari)*(2*i+pari)*sol_hom_un_eta(0, k, j, i)/alpha ;
00645           systeme.set(ligne, colonne+1) 
00646               += (2*i+pari)*(2*i+pari)*sol_hom_trois_eta(0, k, j, i)/alpha ;
00647           sec_membre.set(ligne) 
00648               -= (2*i+pari)*(2*i+pari)* sol_part_eta(0, k, j, i)/alpha ;
00649           }
00650           ligne++ ;
00651           // ... and of its couterpart for V^r
00652           for (int i=0; i<nr; i++) {
00653           systeme.set(ligne, colonne) 
00654               += (2*i+pari)*(2*i+pari)*sol_hom_un_vr(0, k, j, i)/alpha ;
00655           systeme.set(ligne, colonne+1) 
00656               += (2*i+pari)*(2*i+pari)*sol_hom_trois_vr(0, k, j, i)/alpha ;
00657           sec_membre.set(ligne) 
00658               -= (2*i+pari)*(2*i+pari)* sol_part_vr(0, k, j, i)/alpha ;
00659           }
00660           colonne += 2 ; 
00661 
00662               //shells
00663           for (int zone=1 ; zone<nzm1 ; zone++) {
00664           nr = mg.get_nr(zone) ;
00665           alpha = mpaff->get_alpha()[zone] ;
00666           ligne -= 3 ;
00667           systeme.set(ligne, colonne) 
00668               = -sol_hom_un_eta.val_in_bound_jk(zone, j, k)  ;
00669           systeme.set(ligne, colonne+1) 
00670               = -sol_hom_deux_eta.val_in_bound_jk(zone, j, k)  ;
00671           systeme.set(ligne, colonne+2) 
00672               = -sol_hom_trois_eta.val_in_bound_jk(zone, j, k)  ;
00673           systeme.set(ligne, colonne+3) 
00674               = -sol_hom_quatre_eta.val_in_bound_jk(zone, j, k)  ;
00675           sec_membre.set(ligne) += sol_part_eta.val_in_bound_jk(zone, j, k) ;
00676           ligne++ ;
00677           // ... and their couterparts for V^r
00678           systeme.set(ligne, colonne) 
00679               = -sol_hom_un_vr.val_in_bound_jk(zone, j, k)  ;
00680           systeme.set(ligne, colonne+1) 
00681               = -sol_hom_deux_vr.val_in_bound_jk(zone, j, k)  ;
00682           systeme.set(ligne, colonne+2) 
00683               = -sol_hom_trois_vr.val_in_bound_jk(zone, j, k)  ;
00684           systeme.set(ligne, colonne+3) 
00685               = -sol_hom_quatre_vr.val_in_bound_jk(zone, j, k)  ;
00686           sec_membre.set(ligne) += sol_part_vr.val_in_bound_jk(zone, j, k) ;
00687           ligne++ ;
00688 
00689           //derivative of (x+echelle)^(l-1) at -1 
00690           pari = -1 ;
00691           for (int i=0; i<nr; i++) {
00692               systeme.set(ligne, colonne) 
00693               -= pari*i*i*sol_hom_un_eta(zone, k, j, i)/alpha ;
00694               systeme.set(ligne, colonne+1) 
00695               -= pari*i*i*sol_hom_deux_eta(zone, k, j, i)/alpha ;
00696               systeme.set(ligne, colonne+2) 
00697               -= pari*i*i*sol_hom_trois_eta(zone, k, j, i)/alpha ;
00698               systeme.set(ligne, colonne+3) 
00699               -= pari*i*i*sol_hom_quatre_eta(zone, k, j, i)/alpha ;
00700               sec_membre.set(ligne) 
00701               += pari*i*i* sol_part_eta(zone, k, j, i)/alpha ;
00702               pari = -pari ;
00703           }
00704           ligne++ ;
00705           // ... and their couterparts for V^r
00706           pari = -1 ;
00707           for (int i=0; i<nr; i++) {
00708               systeme.set(ligne, colonne) 
00709               -= pari*i*i*sol_hom_un_vr(zone, k, j, i)/alpha ;
00710               systeme.set(ligne, colonne+1) 
00711               -= pari*i*i*sol_hom_deux_vr(zone, k, j, i)/alpha ;
00712               systeme.set(ligne, colonne+2) 
00713               -= pari*i*i*sol_hom_trois_vr(zone, k, j, i)/alpha ;
00714               systeme.set(ligne, colonne+3) 
00715               -= pari*i*i*sol_hom_quatre_vr(zone, k, j, i)/alpha ;
00716               sec_membre.set(ligne) 
00717               += pari*i*i* sol_part_vr(zone, k, j, i)/alpha ;
00718               pari = -pari ;
00719           }
00720           ligne++ ;
00721             
00722           systeme.set(ligne, colonne) 
00723               += sol_hom_un_eta.val_out_bound_jk(zone, j, k)  ;
00724           systeme.set(ligne, colonne+1) 
00725               += sol_hom_deux_eta.val_out_bound_jk(zone, j, k)  ;
00726           systeme.set(ligne, colonne+2) 
00727               += sol_hom_trois_eta.val_out_bound_jk(zone, j, k)  ;
00728           systeme.set(ligne, colonne+3) 
00729               += sol_hom_quatre_eta.val_out_bound_jk(zone, j, k)  ;
00730           sec_membre.set(ligne) -= sol_part_eta.val_out_bound_jk(zone, j, k) ;
00731           ligne++ ;
00732           // ... and their couterparts for V^r
00733           systeme.set(ligne, colonne) 
00734               += sol_hom_un_vr.val_out_bound_jk(zone, j, k)  ;
00735           systeme.set(ligne, colonne+1) 
00736               += sol_hom_deux_vr.val_out_bound_jk(zone, j, k)  ;
00737           systeme.set(ligne, colonne+2) 
00738               += sol_hom_trois_vr.val_out_bound_jk(zone, j, k)  ;
00739           systeme.set(ligne, colonne+3) 
00740               += sol_hom_quatre_vr.val_out_bound_jk(zone, j, k)  ;
00741           sec_membre.set(ligne) -= sol_part_vr.val_out_bound_jk(zone, j, k) ;
00742           ligne++ ;
00743 
00744           //derivative at 1 
00745           pari = 1 ;
00746           for (int i=0; i<nr; i++) {
00747               systeme.set(ligne, colonne) 
00748               += pari*i*i*sol_hom_un_eta(zone, k, j, i)/alpha ;
00749               systeme.set(ligne, colonne+1) 
00750               += pari*i*i*sol_hom_deux_eta(zone, k, j, i)/alpha ;
00751               systeme.set(ligne, colonne+2) 
00752               += pari*i*i*sol_hom_trois_eta(zone, k, j, i)/alpha ;
00753               systeme.set(ligne, colonne+3) 
00754               += pari*i*i*sol_hom_quatre_eta(zone, k, j, i)/alpha ;
00755               sec_membre.set(ligne) 
00756               -= pari*i*i* sol_part_eta(zone, k, j, i)/alpha ;
00757               pari = pari ;
00758           }
00759           ligne++ ;
00760           // ... and their couterparts for V^r
00761           pari = 1 ;
00762           for (int i=0; i<nr; i++) {
00763               systeme.set(ligne, colonne) 
00764               += pari*i*i*sol_hom_un_vr(zone, k, j, i)/alpha ;
00765               systeme.set(ligne, colonne+1) 
00766               += pari*i*i*sol_hom_deux_vr(zone, k, j, i)/alpha ;
00767               systeme.set(ligne, colonne+2) 
00768               += pari*i*i*sol_hom_trois_vr(zone, k, j, i)/alpha ;
00769               systeme.set(ligne, colonne+3) 
00770               += pari*i*i*sol_hom_quatre_vr(zone, k, j, i)/alpha ;
00771               sec_membre.set(ligne) 
00772               -= pari*i*i* sol_part_vr(zone, k, j, i)/alpha ;
00773               pari = pari ;
00774           }
00775           colonne += 4 ;
00776           }         
00777           //Compactified external domain
00778           nr = mg.get_nr(nzm1) ;
00779 
00780           alpha = mpaff->get_alpha()[nzm1] ;
00781           ligne -= 3 ;
00782           //value of (x-1)^(l+2) at -1 :
00783           systeme.set(ligne, colonne) 
00784           -= sol_hom_deux_eta.val_in_bound_jk(nzm1, j, k) ;
00785           systeme.set(ligne, colonne+1) 
00786           -= sol_hom_quatre_eta.val_in_bound_jk(nzm1, j, k) ;
00787           sec_membre.set(ligne) += sol_part_eta.val_in_bound_jk(nzm1, j, k) ;
00788           //... and of its couterpart for V^r
00789           systeme.set(ligne+1, colonne) 
00790           -= sol_hom_deux_vr.val_in_bound_jk(nzm1, j, k) ;
00791           systeme.set(ligne+1, colonne+1) 
00792           -= sol_hom_quatre_vr.val_in_bound_jk(nzm1, j, k) ;
00793           sec_membre.set(ligne+1) += sol_part_vr.val_in_bound_jk(nzm1, j, k) ;
00794             
00795           ligne += 2 ;
00796           //derivative of (x-1)^(l+2) at -1 :
00797           pari = 1 ;
00798           for (int i=0; i<nr; i++) {
00799               systeme.set(ligne, colonne) 
00800               -= 4*alpha*pari*i*i*sol_hom_deux_eta(nzm1, k, j, i) ;
00801               systeme.set(ligne, colonne+1) 
00802               -= 4*alpha*pari*i*i*sol_hom_quatre_eta(nzm1, k, j, i) ;
00803               sec_membre.set(ligne) 
00804               += 4*alpha*pari*i*i* sol_part_eta(nzm1, k, j, i) ;
00805               pari = -pari ;
00806           }
00807           ligne++ ;
00808           // ... and their couterparts for V^r
00809           pari = 1 ;
00810           for (int i=0; i<nr; i++) {
00811               systeme.set(ligne, colonne) 
00812               -= 4*alpha*pari*i*i*sol_hom_deux_vr(nzm1, k, j, i) ;
00813               systeme.set(ligne, colonne+1) 
00814               -= 4*alpha*pari*i*i*sol_hom_quatre_vr(nzm1, k, j, i) ;
00815               sec_membre.set(ligne) 
00816               += 4*alpha*pari*i*i* sol_part_vr(nzm1, k, j, i) ;
00817               pari = -pari ;
00818           }
00819             
00820           // Solution of the system giving the coefficients for the homogeneous 
00821           // solutions
00822           //-------------------------------------------------------------------
00823           systeme.set_lu() ;
00824           Tbl facteurs(systeme.inverse(sec_membre)) ;
00825           int conte = 0 ;
00826 
00827           // everything is put to the right place, the same combination of hom.
00828           // solutions (with some l or -(l+1) factors) must be used for V^r
00829           //-------------------------------------------------------------------
00830           nr = mg.get_nr(0) ; //nucleus
00831           for (int i=0 ; i<nr ; i++) {
00832           cf_eta.set(0, k, j, i) = sol_part_eta(0, k, j, i)
00833               +facteurs(conte)*sol_hom_un_eta(0, k, j, i) 
00834               +facteurs(conte+1)*sol_hom_trois_eta(0, k, j, i) ;
00835           cf_vr.set(0, k, j, i) = sol_part_vr(0, k, j, i)
00836               +facteurs(conte)*sol_hom_un_vr(0, k, j, i) 
00837               +facteurs(conte+1)*sol_hom_trois_vr(0, k, j, i) ;
00838           }
00839           conte += 2 ;
00840           for (int zone=1 ; zone<nzm1 ; zone++) { //shells
00841           nr = mg.get_nr(zone) ;
00842           for (int i=0 ; i<nr ; i++) {
00843               cf_eta.set(zone, k, j, i) = 
00844               sol_part_eta(zone, k, j, i)
00845               +facteurs(conte)*sol_hom_un_eta(zone, k, j, i) 
00846               +facteurs(conte+1)*sol_hom_deux_eta(zone, k, j, i) 
00847               +facteurs(conte+2)*sol_hom_trois_eta(zone, k, j, i) 
00848               +facteurs(conte+3)*sol_hom_quatre_eta(zone, k, j, i) ;
00849               cf_vr.set(zone, k, j, i) = sol_part_vr(zone, k, j, i)
00850               +facteurs(conte)*sol_hom_un_vr(zone, k, j, i) 
00851               +facteurs(conte+1)*sol_hom_deux_vr(zone, k, j, i) 
00852               +facteurs(conte+2)*sol_hom_trois_vr(zone, k, j, i) 
00853               +facteurs(conte+3)*sol_hom_quatre_vr(zone, k, j, i) ;
00854           }
00855           conte+=4 ;
00856           }
00857           nr = mg.get_nr(nz-1) ; //compactified external domain
00858           for (int i=0 ; i<nr ; i++) {
00859           cf_eta.set(nzm1, k, j, i) = sol_part_eta(nzm1, k, j, i)
00860               +facteurs(conte)*sol_hom_deux_eta(nzm1, k, j, i) 
00861               +facteurs(conte+1)*sol_hom_quatre_eta(nzm1, k, j, i) ;
00862           cf_vr.set(nzm1, k, j, i) = sol_part_vr(nzm1, k, j, i)
00863               +facteurs(conte)*sol_hom_deux_vr(nzm1, k, j, i) 
00864               +facteurs(conte+1)*sol_hom_quatre_vr(nzm1, k, j, i) ;
00865 
00866           }
00867       } // End of nullite_plm  
00868       } //End of loop on theta
00869   vr.set_spectral_va().ylm_i() ;
00870   vr += vrl0 ;
00871   het.set_spectral_va().ylm_i() ;
00872  }
00873  resu.set_vr_eta_mu(vr, het, mu().poisson()) ;
00874  if ((nt==1)&&(np==1)) {
00875      resu.set(2).set_etat_zero() ;
00876      resu.set(3).set_etat_zero() ;
00877  }
00878  
00879  return ;
00880   
00881 }

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