sym_tensor_trans_pde.C

00001 /*
00002  *  Functions to solve various PDEs for a divergence-free symmetric tensor.
00003  *
00004  *    (see file sym_tensor.h for documentation).
00005  *
00006  */
00007 
00008 /*
00009  *   Copyright (c) 2005-2006  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 sym_tensor_trans_pde_C[] = "$Header: /cvsroot/Lorene/C++/Source/Tensor/sym_tensor_trans_pde.C,v 1.14 2010/10/11 10:38:34 j_novak Exp $" ;
00029 
00030 /*
00031  * $Id: sym_tensor_trans_pde.C,v 1.14 2010/10/11 10:38:34 j_novak Exp $
00032  * $Log: sym_tensor_trans_pde.C,v $
00033  * Revision 1.14  2010/10/11 10:38:34  j_novak
00034  * *** empty log message ***
00035  *
00036  * Revision 1.13  2010/10/11 10:23:03  j_novak
00037  * Removed methods Sym_tensor_trans::solve_hrr() and Sym_tensor_trans::set_WX_det_one(), as they are no longer relevant.
00038  *
00039  * Revision 1.12  2006/09/05 15:38:45  j_novak
00040  * The fuctions sol_Dirac... are in a seperate file, with new parameters to
00041  * control the boundary conditions.
00042  *
00043  * Revision 1.11  2006/08/31 12:13:22  j_novak
00044  * Added an argument of type Param to Sym_tensor_trans::sol_Dirac_A().
00045  *
00046  * Revision 1.10  2006/06/28 07:48:26  j_novak
00047  * Better treatment of some null cases.
00048  *
00049  * Revision 1.9  2006/06/21 15:42:47  j_novak
00050  * Minor changes.
00051  *
00052  * Revision 1.8  2006/06/20 12:07:15  j_novak
00053  * Improved execution speed for sol_Dirac_tildeB...
00054  *
00055  * Revision 1.7  2006/06/14 10:04:21  j_novak
00056  * New methods sol_Dirac_l01, set_AtB_det_one and set_AtB_trace_zero.
00057  *
00058  * Revision 1.6  2006/06/13 13:30:12  j_novak
00059  * New members sol_Dirac_A and sol_Dirac_tildeB (see documentation).
00060  *
00061  * Revision 1.5  2006/06/12 13:37:23  j_novak
00062  * Added bounds in l (multipolar momentum) for Sym_tensor_trans::solve_hrr.
00063  *
00064  * Revision 1.4  2005/11/28 14:45:17  j_novak
00065  * Improved solution of the Poisson tensor equation in the case of a transverse
00066  * tensor.
00067  *
00068  * Revision 1.3  2005/11/24 14:07:54  j_novak
00069  * Use of Matrice::annule_hard()
00070  *
00071  * Revision 1.2  2005/11/24 09:24:25  j_novak
00072  * Corrected some missing references.
00073  *
00074  * Revision 1.1  2005/09/16 13:58:11  j_novak
00075  * New Poisson solver for a Sym_tensor_trans.
00076  *
00077  *
00078  * $Header: /cvsroot/Lorene/C++/Source/Tensor/sym_tensor_trans_pde.C,v 1.14 2010/10/11 10:38:34 j_novak Exp $
00079  *
00080  */
00081 
00082 // C headers
00083 #include <assert.h>
00084 #include <math.h>
00085 
00086 // Lorene headers
00087 #include "tensor.h"
00088 #include "diff.h"
00089 #include "proto.h"
00090 #include "param.h"
00091 
00092 Sym_tensor_trans Sym_tensor_trans::poisson(const Scalar* h_guess) const {
00093 
00094     // All this has a meaning only for spherical components...
00095     assert(dynamic_cast<const Base_vect_spher*>(triad) != 0x0) ; 
00096     //## ... and affine mapping, for the moment!
00097     const Map_af* mpaff = dynamic_cast<const Map_af*>(mp) ;
00098     assert( mpaff!= 0x0) ;
00099     Sym_tensor_trans resu(*mp, *triad, *met_div) ;
00100 
00101     const Mg3d& gri = *mp->get_mg() ;
00102     int np = gri.get_np(0) ;
00103     int nt = gri.get_nt(0) ;
00104     assert (nt > 4) ;
00105     if (np == 1) {
00106     int nz = gri.get_nzone() ;
00107     double* bornes = new double[nz+1] ;
00108     const double* alp = mpaff->get_alpha() ;
00109     const double* bet = mpaff->get_beta() ; 
00110     for (int lz=0; lz<nz; lz++) {
00111         assert (gri.get_np(lz) == np) ;
00112         assert (gri.get_nt(lz) == nt) ;
00113         switch (gri.get_type_r(lz)) {
00114         case RARE: {
00115             bornes[lz] = bet[lz] ;
00116             break ;
00117         }
00118         case FIN: {
00119             bornes[lz] = bet[lz] - alp[lz] ;
00120             break ;
00121         }
00122         case UNSURR: {
00123             bornes[lz] = double(1) / ( bet[lz] - alp[lz] ) ;
00124             break ;
00125         }
00126         default: {
00127             cout << "Sym_tensor_trans::poisson() : problem with the grid!" 
00128              << endl ;
00129             abort() ;
00130             break ;
00131         }
00132         }
00133     }
00134     if (gri.get_type_r(nz-1) == UNSURR) 
00135         bornes[nz] = 1./(alp[nz-1] + bet[nz-1]) ;
00136     else
00137         bornes[nz] = alp[nz-1] + bet[nz-1] ;
00138     
00139     const Mg3d& gr2 = *gri.get_non_axi() ;
00140     Map_af mp2(gr2, bornes) ;
00141     int np2 = ( np > 3 ? np : 4 ) ;
00142     
00143     Sym_tensor sou_cart(mp2, CON, mp2.get_bvect_spher()) ;
00144     for (int l=1; l<=3; l++)
00145         for (int c=l; c<=3; c++) {
00146         switch (this->operator()(l,c).get_etat() ) {
00147             case ETATZERO: {
00148             sou_cart.set(l,c).set_etat_zero() ;
00149             break ;
00150             }
00151             case ETATUN: {
00152             sou_cart.set(l,c).set_etat_one() ;
00153             break ;
00154             }
00155             case ETATQCQ : {
00156             sou_cart.set(l,c).allocate_all() ;
00157             for (int lz=0; lz<nz; lz++)             
00158                 for (int k=0; k<np2; k++)
00159                 for (int j=0; j<nt; j++)
00160                     for(int i=0; i<gr2.get_nr(lz); i++)
00161                     sou_cart.set(l,c).set_grid_point(lz, k, j, i)
00162                    = this->operator()(l,c).val_grid_point(lz, 0, j, i) ;
00163             break ;
00164             }
00165             default: {
00166             cout << 
00167                 "Sym_tensor_trans::poisson() : source in undefined state!" 
00168                  << endl ;
00169             abort() ;
00170             break ; 
00171             }
00172         }
00173         sou_cart.set(l,c).set_dzpuis(this->operator()(l,c).get_dzpuis()) ;
00174         }
00175     sou_cart.std_spectral_base() ;
00176     sou_cart.change_triad(mp2.get_bvect_cart()) ;
00177     Sym_tensor res_cart(mp2, CON, mp2.get_bvect_cart()) ;
00178     for (int i=1; i<=3; i++)
00179         for(int j=i; j<=3; j++) 
00180         res_cart.set(i,j) = sou_cart(i,j).poisson() ;
00181     res_cart.change_triad(mp2.get_bvect_spher()) ;
00182     Scalar res_A(*mp) ; Scalar big_A = res_cart.compute_A() ;
00183     Scalar res_B(*mp) ; Scalar big_B = res_cart.compute_tilde_B_tt() ;
00184 
00185     switch (big_A.get_etat() ) {
00186         case ETATZERO: {
00187         res_A.set_etat_zero() ;
00188         break ;
00189         }
00190         case ETATUN : {
00191         res_A.set_etat_one() ;
00192         break ;
00193         }
00194         case ETATQCQ : {
00195         res_A.allocate_all() ;
00196         for (int lz=0; lz<nz; lz++)             
00197             for (int k=0; k<np; k++)
00198             for (int j=0; j<nt; j++)
00199                 for(int i=0; i<gri.get_nr(lz); i++)
00200                 res_A.set_grid_point(lz, k, j, i)
00201                     = big_A.val_grid_point(lz, k, j, i) ;
00202         break ;
00203         }
00204         default: {
00205         cout << 
00206             "Sym_tensor_trans::poisson() : res_A in undefined state!" 
00207              << endl ;
00208         abort() ;
00209         break ; 
00210         }
00211     }
00212     res_A.set_spectral_base(big_A.get_spectral_base()) ;
00213     int dzA = big_A.get_dzpuis() ;
00214     res_A.set_dzpuis(dzA) ;
00215 
00216     switch (big_B.get_etat() ) {
00217         case ETATZERO: {
00218         res_B.set_etat_zero() ;
00219         break ;
00220         }
00221         case ETATUN : {
00222         res_B.set_etat_one() ;
00223         break ;
00224         }
00225         case ETATQCQ : {
00226         res_B.allocate_all() ;
00227         for (int lz=0; lz<nz; lz++)             
00228             for (int k=0; k<np; k++)
00229             for (int j=0; j<nt; j++)
00230                 for(int i=0; i<gri.get_nr(lz); i++)
00231                 res_B.set_grid_point(lz, k, j, i)
00232                     = big_B.val_grid_point(lz, k, j, i) ;
00233         break ;
00234         }
00235         default: {
00236         cout << 
00237             "Sym_tensor_trans::poisson() : res_B in undefined state!" 
00238              << endl ;
00239         abort() ;
00240         break ; 
00241         }
00242     }
00243     res_B.set_spectral_base(big_B.get_spectral_base()) ;
00244     int dzB = big_B.get_dzpuis() ;
00245     res_B.set_dzpuis(dzB) ;
00246     
00247     resu.set_AtBtt_det_one(res_A, res_B, h_guess) ;
00248     
00249     delete [] bornes ;
00250     }
00251     else {
00252     assert (np >=4) ;
00253     Sym_tensor_trans sou_cart = *this ;
00254     sou_cart.change_triad(mp->get_bvect_cart()) ;
00255     
00256     Sym_tensor res_cart(*mp, CON, mp->get_bvect_cart()) ;
00257     for (int i=1; i<=3; i++)
00258         for(int j=i; j<=3; j++) 
00259         res_cart.set(i,j) = sou_cart(i,j).poisson() ;
00260 
00261     res_cart.change_triad(*triad) ;
00262     
00263     resu.set_AtBtt_det_one(res_cart.compute_A(), res_cart.compute_tilde_B_tt(), h_guess) ;
00264     
00265     }
00266 #ifndef NDEBUG
00267     Vector dive = resu.divergence(*met_div) ;
00268     dive.dec_dzpuis(2) ;
00269     maxabs(dive, "Sym_tensor_trans::poisson : divergence of the solution") ;
00270 #endif    
00271     return resu ;   
00272 }
00273 
00274 

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