map_af_elliptic.C

00001 /*
00002  *   Copyright (c) 1999-2001 Eric Gourgoulhon
00003  *   Copyright (c) 1999-2001 Philippe Grandclement
00004  *
00005  *   This file is part of LORENE.
00006  *
00007  *   LORENE is free software; you can redistribute it and/or modify
00008  *   it under the terms of the GNU General Public License as published by
00009  *   the Free Software Foundation; either version 2 of the License, or
00010  *   (at your option) any later version.
00011  *
00012  *   LORENE is distributed in the hope that it will be useful,
00013  *   but WITHOUT ANY WARRANTY; without even the implied warranty of
00014  *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00015  *   GNU General Public License for more details.
00016  *
00017  *   You should have received a copy of the GNU General Public License
00018  *   along with LORENE; if not, write to the Free Software
00019  *   Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00020  *
00021  */
00022 
00023 
00024 char map_af_elliptic_C[] = "$Header: /cvsroot/Lorene/C++/Source/Map/map_af_elliptic.C,v 1.11 2007/05/06 10:48:11 p_grandclement Exp $" ;
00025 
00026 /*
00027  * $Id: map_af_elliptic.C,v 1.11 2007/05/06 10:48:11 p_grandclement Exp $
00028  * $Log: map_af_elliptic.C,v $
00029  * Revision 1.11  2007/05/06 10:48:11  p_grandclement
00030  * Modification of a few operators for the vorton project
00031  *
00032  * Revision 1.10  2007/01/16 15:08:07  n_vasset
00033  * New constructor, usn Scalar on mono-domain angular grid for boundary,
00034  * for function sol_elliptic_boundary()
00035  *
00036  * Revision 1.9  2005/11/30 11:09:07  p_grandclement
00037  * Changes for the Bin_ns_bh project
00038  *
00039  * Revision 1.8  2005/08/26 14:02:40  p_grandclement
00040  * Modification of the elliptic solver that matches with an oscillatory exterior solution
00041  * small correction in Poisson tau also...
00042  *
00043  * Revision 1.7  2005/06/09 07:57:31  f_limousin
00044  * Implement a new function sol_elliptic_boundary() and
00045  * Vector::poisson_boundary(...) which solve the vectorial poisson
00046  * equation (method 6) with an inner boundary condition.
00047  *
00048  * Revision 1.6  2004/08/24 09:14:42  p_grandclement
00049  * Addition of some new operators, like Poisson in 2d... It now requieres the
00050  * GSL library to work.
00051  *
00052  * Also, the way a variable change is stored by a Param_elliptic is changed and
00053  * no longer uses Change_var but rather 2 Scalars. The codes using that feature
00054  * will requiere some modification. (It should concern only the ones about monopoles)
00055  *
00056  * Revision 1.5  2004/06/22 08:49:58  p_grandclement
00057  * Addition of everything needed for using the logarithmic mapping
00058  *
00059  * Revision 1.4  2004/03/17 15:58:48  p_grandclement
00060  * Slight modification of sol_elliptic_no_zec
00061  *
00062  * Revision 1.3  2004/02/11 09:47:46  p_grandclement
00063  * Addition of a new elliptic solver, matching with the homogeneous solution
00064  * at the outer shell and not solving in the external domain (more details
00065  * coming soon ; check your local Lorene dealer...)
00066  *
00067  * Revision 1.2  2004/01/28 16:46:23  p_grandclement
00068  * Addition of the sol_elliptic_fixe_der_zero stuff
00069  *
00070  * Revision 1.1  2003/12/11 14:48:48  p_grandclement
00071  * Addition of ALL (and that is a lot !) the files needed for the general elliptic solver ... UNDER DEVELOPEMENT...
00072  *
00073  * 
00074  * $Header: /cvsroot/Lorene/C++/Source/Map/map_af_elliptic.C,v 1.11 2007/05/06 10:48:11 p_grandclement Exp $
00075  *
00076  */
00077 
00078 // Header C : 
00079 #include <stdlib.h>
00080 #include <math.h>
00081 
00082 // Headers Lorene :
00083 #include "tbl.h"
00084 #include "mtbl_cf.h"
00085 #include "map.h"
00086 #include "param_elliptic.h"
00087 #include "proto.h"         
00088 
00089             //----------------------------------------------
00090        //       General elliptic solver
00091       //----------------------------------------------
00092 
00093 void Map_af::sol_elliptic(Param_elliptic& ope_var, const Scalar& source, 
00094               Scalar& pot) const {
00095     
00096   assert(source.get_etat() != ETATNONDEF) ; 
00097   assert(source.get_mp().get_mg() == mg) ; 
00098   assert(pot.get_mp().get_mg() == mg) ; 
00099   
00100   assert(source.check_dzpuis(2) || source.check_dzpuis(3) || 
00101      source.check_dzpuis(4)) ; 
00102   // Spherical harmonic expansion of the source
00103   // ------------------------------------------
00104   
00105   const Valeur& sourva = source.get_spectral_va() ; 
00106   
00107   if (sourva.get_etat() == ETATZERO) {
00108     pot.set_etat_zero() ;
00109     return ;  
00110     }
00111   
00112   // Spectral coefficients of the source
00113   assert(sourva.get_etat() == ETATQCQ) ; 
00114   
00115   Valeur rho(sourva.get_mg()) ; 
00116   sourva.coef() ; 
00117   rho = *(sourva.c_cf) ;    // copy of the coefficients of the source
00118   
00119   rho.ylm() ;           // spherical harmonic transforms 
00120   
00121   // On met les bonnes bases dans le changement de variable 
00122   // de ope_var :
00123   ope_var.var_F.set_spectral_va().coef() ;
00124   ope_var.var_F.set_spectral_va().ylm() ;
00125   ope_var.var_G.set_spectral_va().coef() ;
00126   ope_var.var_G.set_spectral_va().ylm() ;
00127 
00128   // Call to the Mtbl_cf version
00129   // ---------------------------
00130   Mtbl_cf resu = elliptic_solver (ope_var, *(rho.c_cf)) ;
00131   
00132   // Final result returned as a Scalar
00133   // ---------------------------------
00134   
00135   pot.set_etat_zero() ;  // to call Scalar::del_t().
00136   
00137   pot.set_etat_qcq() ; 
00138   
00139   pot.set_spectral_va() = resu ;
00140   pot.set_spectral_va().ylm_i() ; // On repasse en base standard.       
00141   
00142   pot.set_dzpuis(0) ; 
00143 }
00144 
00145 
00146  
00147             //-----------------------------------------------
00148        //       General elliptic solver with boundary as Mtbl-cf
00149       //-------------------------------------------------
00150 
00151 
00152 void Map_af::sol_elliptic_boundary(Param_elliptic& ope_var, const Scalar& source, 
00153                    Scalar& pot,  const Mtbl_cf& bound, 
00154                    double fact_dir, double fact_neu ) const {
00155     
00156   assert(source.get_etat() != ETATNONDEF) ; 
00157   assert(source.get_mp().get_mg() == mg) ; 
00158   assert(pot.get_mp().get_mg() == mg) ; 
00159   
00160   assert(source.check_dzpuis(2) || source.check_dzpuis(3) || 
00161      source.check_dzpuis(4)) ; 
00162   // Spherical harmonic expansion of the source
00163   // ------------------------------------------
00164   
00165   const Valeur& sourva = source.get_spectral_va() ; 
00166   
00167   if (sourva.get_etat() == ETATZERO) {
00168     pot.set_etat_zero() ;
00169     return ;  
00170     }
00171   
00172   // Spectral coefficients of the source
00173   assert(sourva.get_etat() == ETATQCQ) ; 
00174   
00175   Valeur rho(sourva.get_mg()) ; 
00176   sourva.coef() ; 
00177   rho = *(sourva.c_cf) ;    // copy of the coefficients of the source
00178   
00179   rho.ylm() ;           // spherical harmonic transforms 
00180   
00181   // On met les bonnes bases dans le changement de variable 
00182   // de ope_var :
00183   ope_var.var_F.set_spectral_va().coef() ;
00184   ope_var.var_F.set_spectral_va().ylm() ;
00185   ope_var.var_G.set_spectral_va().coef() ;
00186   ope_var.var_G.set_spectral_va().ylm() ;
00187 
00188   // Call to the Mtbl_cf version
00189   // ---------------------------
00190   Mtbl_cf resu = elliptic_solver_boundary (ope_var, *(rho.c_cf),  bound, 
00191                        fact_dir, fact_neu) ;
00192   
00193   // Final result returned as a Scalar
00194   // ---------------------------------
00195   
00196   pot.set_etat_zero() ;  // to call Scalar::del_t().
00197   
00198   pot.set_etat_qcq() ; 
00199   
00200   pot.set_spectral_va() = resu ;
00201   pot.set_spectral_va().ylm_i() ; // On repasse en base standard.       
00202   
00203   pot.set_dzpuis(0) ; 
00204 }
00205 
00206 
00207 
00208 
00209 
00210             //-----------------------------------------------
00211        //       General elliptic solver with boundary as Scalar
00212       //-------------------------------------------------
00213 
00214 
00215 void Map_af::sol_elliptic_boundary(Param_elliptic& ope_var, const Scalar& source, 
00216                    Scalar& pot,  const Scalar& bound, 
00217                    double fact_dir, double fact_neu ) const {
00218     
00219   assert(source.get_etat() != ETATNONDEF) ; 
00220   assert(source.get_mp().get_mg() == mg) ; 
00221   assert(pot.get_mp().get_mg() == mg) ; 
00222   
00223   assert(source.check_dzpuis(2) || source.check_dzpuis(3) || 
00224      source.check_dzpuis(4)) ; 
00225   // Spherical harmonic expansion of the source
00226   // ------------------------------------------
00227   
00228   const Valeur& sourva = source.get_spectral_va() ; 
00229   
00230   if (sourva.get_etat() == ETATZERO) {
00231     pot.set_etat_zero() ;
00232     return ;  
00233     }
00234   
00235   // Spectral coefficients of the source
00236   assert(sourva.get_etat() == ETATQCQ) ; 
00237   
00238   Valeur rho(sourva.get_mg()) ; 
00239   sourva.coef() ; 
00240   rho = *(sourva.c_cf) ;    // copy of the coefficients of the source
00241   
00242   rho.ylm() ;           // spherical harmonic transforms 
00243   
00244   // On met les bonnes bases dans le changement de variable 
00245   // de ope_var :
00246   ope_var.var_F.set_spectral_va().coef() ;
00247   ope_var.var_F.set_spectral_va().ylm() ;
00248   ope_var.var_G.set_spectral_va().coef() ;
00249   ope_var.var_G.set_spectral_va().ylm() ;
00250 
00251   // Call to the Mtbl_cf version
00252   // ---------------------------
00253  
00254 // REMINDER : The scalar bound must be defined on a mono-domain angular grid corresponding with the full grid of the scalar source (following assert()) 
00255  
00256   Scalar bbound = bound;
00257   bbound.set_spectral_va().ylm() ;
00258   const Map& mapp = bbound.get_mp();
00259  
00260   const Mg3d& gri2d = *mapp.get_mg();
00261 
00262   assert(&gri2d == source.get_mp().get_mg()->get_angu_1dom()) ;  // Attention à cet assert !! 
00263   
00264   Mtbl_cf bound2 (gri2d , bbound.get_spectral_base()) ;
00265   bound2.annule_hard() ;  
00266 
00267   if (bbound.get_etat() != ETATZERO){ 
00268  
00269       int nr = gri2d.get_nr(0) ;
00270       int nt = gri2d.get_nt(0) ; 
00271       int np = gri2d.get_np(0) ; 
00272        
00273       for(int k=0; k<np+2; k++)
00274           for (int j=0; j<=nt-1; j++)
00275           for(int xi=0; xi<= nr-1; xi++)
00276           {
00277 
00278   bound2.set(0, k , j , xi) = (*bbound.get_spectral_va().c_cf)(0, k, j, xi) ;   
00279   }
00280   }
00281   Mtbl_cf resu = elliptic_solver_boundary (ope_var, *(rho.c_cf),  bound2, 
00282                        fact_dir, fact_neu) ;
00283   
00284   // Final result returned as a Scalar
00285   // ---------------------------------
00286   
00287   pot.set_etat_zero() ;  // to call Scalar::del_t().
00288   
00289   pot.set_etat_qcq() ; 
00290   
00291   pot.set_spectral_va() = resu ;
00292   pot.set_spectral_va().ylm_i() ; // On repasse en base standard.       
00293   
00294   pot.set_dzpuis(0) ; 
00295 }
00296 
00297 
00298 
00299 
00300 
00301             //----------------------------------------------
00302        //      General elliptic solver with no ZEC
00303       //----------------------------------------------
00304 
00305 void Map_af::sol_elliptic_no_zec(Param_elliptic& ope_var, const Scalar& source, 
00306               Scalar& pot, double val) const {
00307     
00308   assert(source.get_etat() != ETATNONDEF) ; 
00309   assert(source.get_mp().get_mg() == mg) ; 
00310   assert(pot.get_mp().get_mg() == mg) ; 
00311   
00312   assert(source.check_dzpuis(2) || source.check_dzpuis(3) || 
00313      source.check_dzpuis(4)) ; 
00314   // Spherical harmonic expansion of the source
00315   // ------------------------------------------
00316   
00317   const Valeur& sourva = source.get_spectral_va() ; 
00318   
00319   if (sourva.get_etat() == ETATZERO) {
00320     pot.set_etat_zero() ;
00321     return ;  
00322     }
00323   
00324   // Spectral coefficients of the source
00325   assert(sourva.get_etat() == ETATQCQ) ; 
00326   
00327   Valeur rho(sourva.get_mg()) ; 
00328   sourva.coef() ; 
00329   rho = *(sourva.c_cf) ;    // copy of the coefficients of the source
00330   
00331   rho.ylm() ;           // spherical harmonic transforms 
00332     
00333   // On met les bonnes bases dans le changement de variable 
00334   // de ope_var :
00335   ope_var.var_F.set_spectral_va().coef() ;
00336   ope_var.var_F.set_spectral_va().ylm() ;
00337   ope_var.var_G.set_spectral_va().coef() ;
00338   ope_var.var_G.set_spectral_va().ylm() ;
00339 
00340   // Call to the Mtbl_cf version
00341   // ---------------------------
00342   Mtbl_cf resu = elliptic_solver_no_zec (ope_var, *(rho.c_cf), val) ;
00343   
00344   // Final result returned as a Scalar
00345   // ---------------------------------
00346   
00347   pot.set_etat_zero() ;  // to call Scalar::del_t().
00348   
00349   pot.set_etat_qcq() ; 
00350   
00351   pot.set_spectral_va() = resu ;
00352   pot.set_spectral_va().ylm_i() ; // On repasse en base standard.       
00353   
00354   pot.set_dzpuis(0) ; 
00355 }
00356 
00357           //----------------------------------------------
00358        //      General elliptic solver only in the ZEC
00359       //----------------------------------------------
00360 
00361 void Map_af::sol_elliptic_only_zec(Param_elliptic& ope_var, 
00362                    const Scalar& source, 
00363                    Scalar& pot, double val) const {
00364     
00365   assert(source.get_etat() != ETATNONDEF) ; 
00366   assert(source.get_mp().get_mg() == mg) ; 
00367   assert(pot.get_mp().get_mg() == mg) ; 
00368   
00369   assert(source.check_dzpuis(2) || source.check_dzpuis(3) || 
00370      source.check_dzpuis(4)) ; 
00371   // Spherical harmonic expansion of the source
00372   // ------------------------------------------
00373   
00374   const Valeur& sourva = source.get_spectral_va() ; 
00375   
00376   if (sourva.get_etat() == ETATZERO) {
00377     pot.set_etat_zero() ;
00378     return ;  
00379     }
00380   
00381   // Spectral coefficients of the source
00382   assert(sourva.get_etat() == ETATQCQ) ; 
00383   
00384   Valeur rho(sourva.get_mg()) ; 
00385   sourva.coef() ; 
00386   rho = *(sourva.c_cf) ;    // copy of the coefficients of the source
00387   
00388   rho.ylm() ;           // spherical harmonic transforms 
00389     
00390   // On met les bonnes bases dans le changement de variable 
00391   // de ope_var :
00392   ope_var.var_F.set_spectral_va().coef() ;
00393   ope_var.var_F.set_spectral_va().ylm() ;
00394   ope_var.var_G.set_spectral_va().coef() ;
00395   ope_var.var_G.set_spectral_va().ylm() ;
00396 
00397   // Call to the Mtbl_cf version
00398   // ---------------------------
00399   Mtbl_cf resu = elliptic_solver_only_zec (ope_var, *(rho.c_cf), val) ;
00400   
00401   // Final result returned as a Scalar
00402   // ---------------------------------
00403   
00404   pot.set_etat_zero() ;  // to call Scalar::del_t().
00405   
00406   pot.set_etat_qcq() ; 
00407   
00408   pot.set_spectral_va() = resu ;
00409   pot.set_spectral_va().ylm_i() ; // On repasse en base standard.       
00410   
00411   pot.set_dzpuis(0) ; 
00412 }
00413 
00414             //----------------------------------------------
00415        //      General elliptic solver with no ZEC
00416            //  and a mtaching with sin(r)/r
00417       //----------------------------------------------
00418 
00419 void Map_af::sol_elliptic_sin_zec(Param_elliptic& ope_var, 
00420                   const Scalar& source, Scalar& pot, double* amplis, double* phases) const {
00421     
00422   assert(source.get_etat() != ETATNONDEF) ; 
00423   assert(source.get_mp().get_mg() == mg) ; 
00424   assert(pot.get_mp().get_mg() == mg) ; 
00425   
00426   assert(source.check_dzpuis(2) || source.check_dzpuis(3) || 
00427      source.check_dzpuis(4)) ; 
00428   // Spherical harmonic expansion of the source
00429   // ------------------------------------------
00430   
00431   const Valeur& sourva = source.get_spectral_va() ; 
00432   
00433   if (sourva.get_etat() == ETATZERO) {
00434     pot.set_etat_zero() ;
00435     return ;  
00436     }
00437   
00438   // Spectral coefficients of the source
00439   assert(sourva.get_etat() == ETATQCQ) ; 
00440   
00441   Valeur rho(sourva.get_mg()) ; 
00442   sourva.coef() ; 
00443   rho = *(sourva.c_cf) ;    // copy of the coefficients of the source
00444   
00445   rho.ylm() ;           // spherical harmonic transforms 
00446     
00447   // On met les bonnes bases dans le changement de variable 
00448   // de ope_var :
00449   ope_var.var_F.set_spectral_va().coef() ;
00450   ope_var.var_F.set_spectral_va().ylm() ;
00451   ope_var.var_G.set_spectral_va().coef() ;
00452   ope_var.var_G.set_spectral_va().ylm() ;
00453 
00454   // Call to the Mtbl_cf version
00455   // ---------------------------
00456   Mtbl_cf resu = elliptic_solver_sin_zec (ope_var, *(rho.c_cf), amplis, phases) ;
00457   
00458   // Final result returned as a Scalar
00459   // ---------------------------------
00460   
00461   pot.set_etat_zero() ;  // to call Scalar::del_t().
00462   
00463   pot.set_etat_qcq() ; 
00464   
00465   pot.set_spectral_va() = resu ;
00466   pot.set_spectral_va().ylm_i() ; // On repasse en base standard.       
00467   
00468   pot.set_dzpuis(0) ; 
00469 }
00470 
00471 
00472             //----------------------------------------------
00473        //      General elliptic solver with no ZEC
00474       //----------------------------------------------
00475 
00476 void Map_af::sol_elliptic_fixe_der_zero (double valeur, 
00477                      Param_elliptic& ope_var, 
00478                      const Scalar& source, 
00479                      Scalar& pot) const {
00480     
00481   assert(source.get_etat() != ETATNONDEF) ; 
00482   assert(source.get_mp().get_mg() == mg) ; 
00483   assert(pot.get_mp().get_mg() == mg) ; 
00484   
00485   assert(source.check_dzpuis(2) || source.check_dzpuis(3) || 
00486      source.check_dzpuis(4)) ; 
00487   // Spherical harmonic expansion of the source
00488   // ------------------------------------------
00489   
00490   const Valeur& sourva = source.get_spectral_va() ; 
00491   
00492   if (sourva.get_etat() == ETATZERO) {
00493     pot.set_etat_zero() ;
00494     return ;  
00495     }
00496   
00497   // Spectral coefficients of the source
00498   assert(sourva.get_etat() == ETATQCQ) ; 
00499   
00500   Valeur rho(sourva.get_mg()) ; 
00501   sourva.coef() ; 
00502   rho = *(sourva.c_cf) ;    // copy of the coefficients of the source
00503   
00504   rho.ylm() ;           // spherical harmonic transforms 
00505     
00506   // On met les bonnes bases dans le changement de variable 
00507   // de ope_var :
00508   ope_var.var_F.set_spectral_va().coef() ;
00509   ope_var.var_F.set_spectral_va().ylm() ;
00510   ope_var.var_G.set_spectral_va().coef() ;
00511   ope_var.var_G.set_spectral_va().ylm() ;
00512 
00513   // Call to the Mtbl_cf version
00514   // ---------------------------
00515   valeur *= alpha[0] ;
00516   Mtbl_cf resu = elliptic_solver_fixe_der_zero (valeur, ope_var, *(rho.c_cf)) ;
00517   
00518   // Final result returned as a Scalar
00519   // ---------------------------------
00520   
00521   pot.set_etat_zero() ;  // to call Scalar::del_t().
00522   
00523   pot.set_etat_qcq() ; 
00524   
00525   pot.set_spectral_va() = resu ;
00526   pot.set_spectral_va().ylm_i() ; // On repasse en base standard.       
00527   
00528   pot.set_dzpuis(0) ; 
00529 }
00530 

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