scalar.h

00001 /*
00002  *  Definition of Lorene class Scalar
00003  *
00004  */
00005 
00006 /*
00007  *   Copyright (c) 2003-2004 Eric Gourgoulhon & Jerome Novak
00008  *   
00009  *   Copyright (c) 1999-2000 Jean-Alain Marck (for previous class Cmp)
00010  *   Copyright (c) 1999-2002 Eric Gourgoulhon (for previous class Cmp)
00011  *   Copyright (c) 1999-2001 Philippe Grandclement (for previous class Cmp)
00012  *   Copyright (c) 2000-2002 Jerome Novak (for previous class Cmp)
00013  *   Copyright (c) 2000-2001 Keisuke Taniguchi (for previous class Cmp)
00014  *
00015  *   This file is part of LORENE.
00016  *
00017  *   LORENE is free software; you can redistribute it and/or modify
00018  *   it under the terms of the GNU General Public License as published by
00019  *   the Free Software Foundation; either version 2 of the License, or
00020  *   (at your option) any later version.
00021  *
00022  *   LORENE is distributed in the hope that it will be useful,
00023  *   but WITHOUT ANY WARRANTY; without even the implied warranty of
00024  *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00025  *   GNU General Public License for more details.
00026  *
00027  *   You should have received a copy of the GNU General Public License
00028  *   along with LORENE; if not, write to the Free Software
00029  *   Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00030  *
00031  */
00032 
00033 
00034 #ifndef __SCALAR_H_ 
00035 #define __SCALAR_H_ 
00036 
00037 
00038 /*
00039  * $Id: scalar.h,v 1.88 2012/01/17 15:05:46 j_penner Exp $
00040  * $Log: scalar.h,v $
00041  * Revision 1.88  2012/01/17 15:05:46  j_penner
00042  * *** empty log message ***
00043  *
00044  * Revision 1.87  2012/01/17 10:16:27  j_penner
00045  * functions added: sarra_filter_r, sarra_filter_r_all_domains, Heaviside
00046  *
00047  * Revision 1.86  2011/04/08 13:13:09  e_gourgoulhon
00048  * Changed the comment of function val_point to indicate specifically the
00049  * division by r^dzpuis in the compactified external domain.
00050  *
00051  * Revision 1.85  2008/09/29 13:23:51  j_novak
00052  * Implementation of the angular mapping associated with an affine
00053  * mapping. Things must be improved to take into account the domain index.
00054  *
00055  * Revision 1.84  2008/09/22 19:08:01  j_novak
00056  * New methods to deal with boundary conditions
00057  *
00058  * Revision 1.83  2008/05/24 15:05:22  j_novak
00059  * New method Scalar::match_tau to match the output of an explicit time-marching scheme with the tau method.
00060  *
00061  * Revision 1.82  2007/12/21 16:06:16  j_novak
00062  * Methods to filter Tensor, Vector and Sym_tensor objects.
00063  *
00064  * Revision 1.81  2007/10/31 10:33:11  j_novak
00065  * Added exponential filters to smooth Gibbs-type phenomena.
00066  *
00067  * Revision 1.80  2007/06/21 19:56:36  k_taniguchi
00068  * Introduction of another filtre_r.
00069  *
00070  * Revision 1.79  2007/05/06 10:48:08  p_grandclement
00071  * Modification of a few operators for the vorton project
00072  *
00073  * Revision 1.78  2007/01/16 15:05:59  n_vasset
00074  * New constructor (taking a Scalar in mono-domain angular grid for
00075  * boundary) for function sol_elliptic_boundary
00076  *
00077  * Revision 1.77  2006/05/26 09:00:09  j_novak
00078  * New members for multiplication or division by cos(theta).
00079  *
00080  * Revision 1.76  2005/11/30 13:48:06  e_gourgoulhon
00081  * Replaced M_PI/2 by 1.57... in argument list of sol_elliptic_sin_zec
00082  * (in order not to require the definition of M_PI).
00083  *
00084  * Revision 1.75  2005/11/30 11:09:03  p_grandclement
00085  * Changes for the Bin_ns_bh project
00086  *
00087  * Revision 1.74  2005/11/17 15:29:46  e_gourgoulhon
00088  * Added arithmetics with Mtbl.
00089  *
00090  * Revision 1.73  2005/10/25 08:56:34  p_grandclement
00091  * addition of std_spectral_base in the case of odd functions near the origin
00092  *
00093  * Revision 1.72  2005/09/07 13:10:47  j_novak
00094  * Added a filter setting to zero all mulitpoles in a given range.
00095  *
00096  * Revision 1.71  2005/08/26 14:02:38  p_grandclement
00097  * Modification of the elliptic solver that matches with an oscillatory exterior solution
00098  * small correction in Poisson tau also...
00099  *
00100  * Revision 1.70  2005/08/25 12:14:07  p_grandclement
00101  * Addition of a new method to solve the scalar Poisson equation, based on a multi-domain Tau-method
00102  *
00103  * Revision 1.69  2005/06/09 07:56:25  f_limousin
00104  * Implement a new function sol_elliptic_boundary() and
00105  * Vector::poisson_boundary(...) which solve the vectorial poisson
00106  * equation (method 6) with an inner boundary condition.
00107  *
00108  * Revision 1.68  2005/06/08 12:35:18  j_novak
00109  * New method for solving divergence-like ODEs.
00110  *
00111  * Revision 1.67  2005/05/20 14:42:30  j_novak
00112  * Added the method Scalar::get_spectral_base().
00113  *
00114  * Revision 1.66  2005/04/04 21:28:57  e_gourgoulhon
00115  * Added argument lambda to method poisson_angu
00116  * to deal with the generalized angular Poisson equation:
00117  *    Lap_ang u + lambda u = source.
00118  *
00119  * Revision 1.65  2004/12/14 09:09:39  f_limousin
00120  * Modif. comments.
00121  *
00122  * Revision 1.64  2004/11/23 12:41:53  f_limousin
00123  * Intoduce function poisson_dir_neu(...) to solve a scalar poisson
00124  * equation with a mixed boundary condition (Dirichlet + Neumann).
00125  *
00126  * Revision 1.63  2004/10/11 15:09:00  j_novak
00127  * The radial manipulation functions take Scalar as arguments, instead of Cmp.
00128  * Added a conversion operator from Scalar to Cmp.
00129  * The Cmp radial manipulation function make conversion to Scalar, call to the
00130  * Map_radial version with a Scalar argument and back.
00131  *
00132  * Revision 1.62  2004/08/24 09:14:40  p_grandclement
00133  * Addition of some new operators, like Poisson in 2d... It now requieres the
00134  * GSL library to work.
00135  *
00136  * Also, the way a variable change is stored by a Param_elliptic is changed and
00137  * no longer uses Change_var but rather 2 Scalars. The codes using that feature
00138  * will requiere some modification. (It should concern only the ones about monopoles)
00139  *
00140  * Revision 1.61  2004/07/27 08:24:26  j_novak
00141  * Modif. comments
00142  *
00143  * Revision 1.60  2004/07/26 16:02:21  j_novak
00144  * Added a flag to specify whether the primitive should be zero either at r=0
00145  * or at r going to infinity.
00146  *
00147  * Revision 1.59  2004/07/06 13:36:27  j_novak
00148  * Added methods for desaliased product (operator |) only in r direction.
00149  *
00150  * Revision 1.58  2004/06/22 08:49:57  p_grandclement
00151  * Addition of everything needed for using the logarithmic mapping
00152  *
00153  * Revision 1.57  2004/06/14 15:24:23  e_gourgoulhon
00154  * Added method primr (radial primitive).
00155  *
00156  * Revision 1.56  2004/06/11 14:29:56  j_novak
00157  * Scalar::multipole_spectrum() is now a const method.
00158  *
00159  * Revision 1.55  2004/05/24 14:07:31  e_gourgoulhon
00160  * Method set_domain now includes a call to del_deriv() for safety.
00161  *
00162  * Revision 1.54  2004/05/07 11:26:10  f_limousin
00163  * New method filtre_r(int*)
00164  *
00165  * Revision 1.53  2004/03/22 13:12:43  j_novak
00166  * Modification of comments to use doxygen instead of doc++
00167  *
00168  * Revision 1.52  2004/03/17 15:58:47  p_grandclement
00169  * Slight modification of sol_elliptic_no_zec
00170  *
00171  * Revision 1.51  2004/03/11 12:07:30  e_gourgoulhon
00172  * Added method visu_section_anim.
00173  *
00174  * Revision 1.50  2004/03/08 15:45:38  j_novak
00175  * Modif. comment
00176  *
00177  * Revision 1.49  2004/03/05 15:09:40  e_gourgoulhon
00178  * Added method smooth_decay.
00179  *
00180  * Revision 1.48  2004/03/01 09:57:02  j_novak
00181  * the wave equation is solved with Scalars. It now accepts a grid with a
00182  * compactified external domain, which the solver ignores and where it copies
00183  * the values of the field from one time-step to the next.
00184  *
00185  * Revision 1.47  2004/02/27 09:43:58  f_limousin
00186  * New methods filtre_phi(int) and filtre_theta(int).
00187  *
00188  * Revision 1.46  2004/02/26 22:46:26  e_gourgoulhon
00189  * Added methods derive_cov, derive_con and derive_lie.
00190  *
00191  * Revision 1.45  2004/02/21 17:03:49  e_gourgoulhon
00192  * -- Method "point" renamed "val_grid_point".
00193  * -- Method "set_point" renamed "set_grid_point".
00194  *
00195  * Revision 1.44  2004/02/19 22:07:35  e_gourgoulhon
00196  * Added argument "comment" in method spectral_display.
00197  *
00198  * Revision 1.43  2004/02/11 09:47:44  p_grandclement
00199  * Addition of a new elliptic solver, matching with the homogeneous solution
00200  * at the outer shell and not solving in the external domain (more details
00201  * coming soon ; check your local Lorene dealer...)
00202  *
00203  * Revision 1.42  2004/01/28 16:46:22  p_grandclement
00204  * Addition of the sol_elliptic_fixe_der_zero stuff
00205  *
00206  * Revision 1.41  2004/01/28 13:25:40  j_novak
00207  * The ced_mult_r arguments have been suppressed from the Scalar::*dsd* methods.
00208  * In the div/mult _r_dzpuis, there is no more default value.
00209  *
00210  * Revision 1.40  2004/01/28 10:39:17  j_novak
00211  * Comments modified.
00212  *
00213  * Revision 1.39  2004/01/27 15:10:01  j_novak
00214  * New methods Scalar::div_r_dzpuis(int) and Scalar_mult_r_dzpuis(int)
00215  * which replace div_r_inc*. Tried to clean the dzpuis handling.
00216  * WARNING: no testing at this point!!
00217  *
00218  * Revision 1.38  2004/01/23 13:25:44  e_gourgoulhon
00219  * Added methods set_inner_boundary and set_outer_boundary.
00220  * Methods set_val_inf and set_val_hor, which are particular cases of
00221  * the above, have been suppressed.
00222  *
00223  * Revision 1.37  2004/01/22 16:10:09  e_gourgoulhon
00224  * Added (provisory) method div_r_inc().
00225  *
00226  * Revision 1.36  2003/12/16 06:32:20  e_gourgoulhon
00227  * Added method visu_box.
00228  *
00229  * Revision 1.35  2003/12/14 21:46:35  e_gourgoulhon
00230  * Added argument start_dx in visu_section.
00231  *
00232  * Revision 1.34  2003/12/11 16:19:38  e_gourgoulhon
00233  * Added method visu_section for visualization with OpenDX.
00234  *
00235  * Revision 1.33  2003/12/11 14:48:47  p_grandclement
00236  * Addition of ALL (and that is a lot !) the files needed for the general elliptic solver ... UNDER DEVELOPEMENT...
00237  *
00238  * Revision 1.32  2003/11/13 13:43:53  p_grandclement
00239  * Addition of things needed for Bhole::update_metric (const Etoile_bin&, double, double)
00240  *
00241  * Revision 1.31  2003/11/06 14:43:37  e_gourgoulhon
00242  * Gave a name to const arguments in certain method prototypes (e.g.
00243  * constructors) to correct a bug of DOC++.
00244  *
00245  * Revision 1.30  2003/11/04 22:55:50  e_gourgoulhon
00246  * Added new methods mult_cost(), mult_sint() and div_sint().
00247  *
00248  * Revision 1.29  2003/10/29 13:09:11  e_gourgoulhon
00249  * -- Added integer argument to derivative functions dsdr, etc...
00250  *    so that one can choose the dzpuis of the result (default=2).
00251  * -- Change of method name: laplacien --> laplacian.
00252  *
00253  * Revision 1.28  2003/10/29 11:00:42  e_gourgoulhon
00254  * Virtual functions dec_dzpuis and inc_dzpuis have now an integer argument to
00255  *  specify by which amount dzpuis is to be increased.
00256  * Accordingly virtual methods dec2_dzpuis and inc2_dzpuis have been suppressed.
00257  *
00258  * Revision 1.27  2003/10/20 14:26:02  j_novak
00259  * New assignement operators.
00260  *
00261  * Revision 1.26  2003/10/19 19:46:33  e_gourgoulhon
00262  * -- Method spectral_display now virtual (from Tensor), list of argument
00263  *    changed.
00264  *
00265  * Revision 1.25  2003/10/17 13:46:14  j_novak
00266  * The argument is now between 1 and 3 (instead of 0->2)
00267  *
00268  * Revision 1.24  2003/10/16 15:23:41  e_gourgoulhon
00269  * Name of method div_r_ced() changed to div_r_inc2().
00270  * Name of method div_rsint_ced() changed to div_rsint_inc2().
00271  *
00272  * Revision 1.23  2003/10/15 21:10:11  e_gourgoulhon
00273  * Added method poisson_angu().
00274  *
00275  * Revision 1.22  2003/10/15 16:03:35  j_novak
00276  * Added the angular Laplace operator for Scalar.
00277  *
00278  * Revision 1.21  2003/10/15 10:29:05  e_gourgoulhon
00279  * Added new members p_dsdt and p_stdsdp.
00280  * Added new methods dsdt(), stdsdp() and div_tant().
00281  *
00282  * Revision 1.20  2003/10/13 13:52:39  j_novak
00283  * Better managment of derived quantities.
00284  *
00285  * Revision 1.19  2003/10/10 15:57:27  j_novak
00286  * Added the state one (ETATUN) to the class Scalar
00287  *
00288  * Revision 1.18  2003/10/08 14:24:08  j_novak
00289  * replaced mult_r_zec with mult_r_ced
00290  *
00291  * Revision 1.17  2003/10/06 16:16:02  j_novak
00292  * New constructor from a Tensor.
00293  *
00294  * Revision 1.16  2003/10/06 13:58:45  j_novak
00295  * The memory management has been improved.
00296  * Implementation of the covariant derivative with respect to the exact Tensor
00297  * type.
00298  *
00299  * Revision 1.15  2003/10/05 21:06:31  e_gourgoulhon
00300  * - Added new methods div_r_ced() and div_rsint_ced().
00301  * - Added new virtual method std_spectral_base()
00302  * - Removed method std_spectral_base_scal()
00303  *
00304  * Revision 1.14  2003/10/01 13:02:58  e_gourgoulhon
00305  * Suppressed the constructor from Map* .
00306  *
00307  * Revision 1.13  2003/09/29 12:52:56  j_novak
00308  * Methods for changing the triad are implemented.
00309  *
00310  * Revision 1.12  2003/09/25 09:33:36  j_novak
00311  * Added methods for integral calculation and various manipulations
00312  *
00313  * Revision 1.11  2003/09/25 09:11:21  e_gourgoulhon
00314  * Added functions for radial operations (divr, etc...)
00315  *
00316  * Revision 1.10  2003/09/25 08:55:23  e_gourgoulhon
00317  * Added members raccord*.
00318  *
00319  * Revision 1.9  2003/09/25 08:50:11  j_novak
00320  * Added the members import
00321  *
00322  * Revision 1.8  2003/09/25 08:13:51  j_novak
00323  * Added method for calculating derivatives
00324  *
00325  * Revision 1.7  2003/09/25 07:59:26  e_gourgoulhon
00326  * Added prototypes for PDE resolutions.
00327  *
00328  * Revision 1.6  2003/09/25 07:17:58  j_novak
00329  * Method asymptot implemented.
00330  *
00331  * Revision 1.5  2003/09/24 20:53:38  e_gourgoulhon
00332  * Added  -- constructor by conversion from a Cmp
00333  *        -- assignment from Cmp
00334  *
00335  * Revision 1.4  2003/09/24 15:10:54  j_novak
00336  * Suppression of the etat flag in class Tensor (still present in Scalar)
00337  *
00338  * Revision 1.3  2003/09/24 12:01:44  j_novak
00339  * Added friend functions for math.
00340  *
00341  * Revision 1.2  2003/09/24 10:22:01  e_gourgoulhon
00342  * still in progress...
00343  *
00344  * Revision 1.1  2003/09/22 12:50:47  e_gourgoulhon
00345  * First version: not ready yet!
00346  *
00347  *
00348  * $Header: /cvsroot/Lorene/C++/Include/scalar.h,v 1.88 2012/01/17 15:05:46 j_penner Exp $
00349  *
00350  */
00351 
00352 // Headers Lorene 
00353 #include "valeur.h"
00354 #include "tensor.h"
00355 
00356 class Param ; 
00357 class Cmp ;
00358 class Param_elliptic ;
00359 
00367 class Scalar : public Tensor {
00368   
00369   // Data : 
00370   // -----
00371  protected:
00372   
00376   int etat ; 
00377   
00383   int dzpuis ;  
00384   
00385   Valeur va ;       
00386   
00387   // Derived data : 
00388   // ------------
00389  protected:
00391   mutable Scalar* p_dsdr ;  
00392 
00396   mutable Scalar* p_srdsdt ;    
00397 
00401   mutable Scalar* p_srstdsdp ;
00402 
00404   mutable Scalar* p_dsdt ;  
00405 
00409   mutable Scalar* p_stdsdp ;    
00410   
00414   mutable Scalar* p_dsdx ;  
00415   
00419   mutable Scalar* p_dsdy ;  
00420 
00424   mutable Scalar* p_dsdz ;  
00425   
00428   mutable Scalar* p_lap ;   
00429   
00432   mutable Scalar* p_lapang ;    
00433    
00435   mutable Scalar* p_dsdradial ; 
00436 
00438   mutable Scalar* p_dsdrho ;    
00439 
00443   mutable int ind_lap ; 
00444 
00448   mutable Tbl* p_integ ; 
00449   
00450   // Constructors - Destructor
00451   // -------------------------
00452   
00453  public:
00454   
00455   explicit Scalar(const Map& mpi) ; 
00456 
00458   Scalar(const Tensor& a) ;     
00459 
00460   Scalar(const Scalar& a) ;     
00461   
00462   Scalar(const Cmp& a) ;    
00463   
00465   Scalar(const Map&, const Mg3d&, FILE* ) ;         
00466   
00467   virtual ~Scalar() ;           
00468   
00469   
00470   // Memory management
00471   // -----------------
00472  protected:
00473   void del_t() ;            
00474   virtual void del_deriv() const;       
00475   void set_der_0x0() const;     
00476   
00477  public:
00478   
00483   virtual void set_etat_nondef() ;   
00484   
00490   virtual void set_etat_zero() ;        
00491   
00498   virtual void set_etat_qcq() ;     
00499   
00505   void set_etat_one() ;     
00506   
00515   virtual void allocate_all() ; 
00516   
00525   void annule_hard() ;
00526   
00527   // Extraction of information
00528   // -------------------------
00529     public:
00533   int get_etat() const {return etat;} ; 
00534   
00536   int get_dzpuis() const {return dzpuis;} ; 
00537   
00541   bool dz_nonzero() const ; 
00542     
00547   bool check_dzpuis(int dzi) const ; 
00548   
00549   // Assignment
00550   // -----------
00551  public: 
00553   void operator=(const Scalar& a) ; 
00554   
00556   virtual void operator=(const Tensor& a) ; 
00557 
00558   void operator=(const Cmp& a) ;    
00559   void operator=(const Valeur& a) ; 
00560   void operator=(const Mtbl& a) ;    
00561   void operator=(double ) ;  
00562   void operator=(int ) ;         
00563   
00564   // Conversion oprators
00565   //---------------------
00566   operator Cmp() const ; 
00567 
00568   // Access to individual elements
00569   // -----------------------------
00570     public:
00571   
00573   const Valeur& get_spectral_va() const {return va;} ; 
00574   
00576   Valeur& set_spectral_va() {return va;} ; 
00577   
00587   Tbl& set_domain(int l) {
00588     assert(etat == ETATQCQ) ;
00589     del_deriv() ; 
00590     return va.set(l) ;
00591   };
00592   
00597   const Tbl& domain(int l) const {
00598     assert( (etat == ETATQCQ) || (etat == ETATUN) ) ;
00599     return va(l) ;
00600   };
00601   
00602   
00609   double val_grid_point(int l, int k, int j, int i) const {
00610     assert(etat != ETATNONDEF) ;
00611     if (etat == ETATZERO) {
00612       double zero = 0. ;
00613       return zero ; 
00614     }
00615     else {
00616       if (etat == ETATUN) {
00617       double one = 1. ;
00618       return one ;
00619       }
00620       else{         
00621     return va(l, k, j, i) ;
00622       }
00623     }
00624   };
00625   
00640   double val_point(double r, double theta, double phi) const ; 
00641   
00642   
00656   double& set_grid_point(int l, int k, int j, int i) {
00657     assert(etat == ETATQCQ) ;
00658     return va.set(l, k, j, i) ;
00659   };
00660   
00661     
00672   virtual void annule(int l_min, int l_max) ; 
00673   
00679   void set_inner_boundary(int l, double x) ;
00680     
00686   void set_outer_boundary(int l, double x) ;
00687 
00696   Tbl multipole_spectrum () const ;
00697 
00704   Tbl tbl_out_bound(int l_dom, bool leave_ylm = false) ;
00705   
00712   Tbl tbl_in_bound(int n, bool leave_ylm = false) ;
00713   
00720   Scalar scalar_out_bound(int n, bool leave_ylm = false) ;
00721   
00722   // Differential operators and others
00723   // ---------------------------------
00724  public:
00729   const Scalar& dsdr() const ; 
00730  
00735   const Scalar& srdsdt() const ; 
00736   
00741   const Scalar& srstdsdp() const ; 
00742   
00745   const Scalar& dsdt() const ; 
00746   
00754   const Scalar& dsdradial() const ; 
00755   
00760   const Scalar& dsdrho() const ;
00761  
00764   const Scalar& stdsdp() const ; 
00765   
00771   const Scalar& dsdx() const ;  
00772   
00778   const Scalar& dsdy() const ;  
00779   
00785   const Scalar& dsdz() const ;  
00786   
00793   const Scalar& deriv(int i) const ;    
00794   
00799     const Vector& derive_cov(const Metric& gam) const ; 
00800 
00801 
00807     const Vector& derive_con(const Metric& gam) const ; 
00808 
00810     Scalar derive_lie(const Vector& v) const ; 
00811 
00812 
00821   const Scalar& laplacian(int ced_mult_r = 4) const ; 
00822   
00830   const Scalar& lapang() const ; 
00831   
00833   void div_r() ;    
00834  
00839   void div_r_dzpuis(int ced_mult_r) ; 
00840   
00844   void div_r_ced() ;
00845 
00847   void mult_r() ;  
00848   
00853   void mult_r_dzpuis(int ced_mult_r) ; 
00854   
00858   void mult_r_ced() ;
00859   
00861   void mult_rsint() ;   
00862   
00867   void mult_rsint_dzpuis(int ced_mult_r) ; 
00868 
00870   void div_rsint() ;    
00871   
00876   void div_rsint_dzpuis(int ced_mult_r) ; 
00877   
00878   void mult_cost() ;   
00879 
00880   void div_cost() ;    
00881 
00882   void mult_sint() ;   
00883 
00884   void div_sint() ;    
00885 
00886   void div_tant() ;    
00887   
00899   Scalar primr(bool null_infty = true) const  ;         
00900 
00908   double integrale() const ; 
00909     
00919   const Tbl& integrale_domains() const ; 
00920     
00925   virtual void dec_dzpuis(int dec = 1) ; 
00926 
00931   virtual void inc_dzpuis(int inc = 1) ; 
00932     
00936   virtual void change_triad(const Base_vect& new_triad) ; 
00937     
00942   void filtre (int n) ;
00943     
00948   void filtre_r (int* nn) ;
00949 
00953   void filtre_r (int n, int nzone) ;
00954 
00965 //  virtual void exponential_filter_r(int lzmin, int lzmax, int p, 
00966 //              double alpha= -16.) ;
00967   virtual void exponential_filter_r(int lzmin, int lzmax, int p, 
00968                 double alpha= -16.) ;
00969 
00980   void sarra_filter_r(int lzmin, int lzmax, double p, 
00981                 double alpha= -1E-16) ;
00982 
00988   void exp_filter_r_all_domains(Scalar &ss, int p, double alpha=-16.) ;
00989 
00996   void sarra_filter_r_all_domains(double p, double alpha=1E-16) ;
00997 
01005   virtual void exponential_filter_ylm(int lzmin, int lzmax, int p, 
01006                 double alpha= -16.) ;
01007 
01015   void annule_l (int l_min, int l_max, bool ylm_output= false ) ;
01016 
01021   void filtre_phi (int n, int zone) ;
01022  
01027   void filtre_tp(int nn, int nz1, int nz2) ;
01028   
01029    
01035   void fixe_decroissance (int puis) ;
01036 
01042     void smooth_decay(int k, int n) ; 
01043 
01048   void raccord(int n) ;
01049     
01056   void raccord_c1_zec(int puis, int nbre, int lmax) ;
01057 
01061   void raccord_externe(int puis, int nbre, int lmax) ;
01062 
01071   void match_tau(Param& par_bc, Param* par_mat=0x0) ;
01072 
01073   // Outputs
01074   // -------
01075  public:
01076   virtual void sauve(FILE *) const ;        
01077     
01088     virtual void spectral_display(const char* comment = 0x0, 
01089                             double threshold = 1.e-7, int precision = 4, 
01090                 ostream& ostr = cout) const ;
01091 
01093   friend ostream& operator<<(ostream& , const Scalar & ) ;  
01094   
01118     void visu_section(const char section_type, double aa, double umin, double umax, double vmin,
01119         double vmax, const char* title = 0x0, const char* filename = 0x0,
01120         bool start_dx = true, int nu = 200, int nv = 200) const ;   
01121 
01150     void visu_section(const Tbl& plane, double umin, double umax, double vmin,
01151         double vmax, const char* title = 0x0, const char* filename = 0x0,
01152         bool start_dx = true, int nu = 200, int nv = 200) const ;   
01153 
01182    void visu_section_anim(const char section_type, double aa, double umin, 
01183         double umax, double vmin, double vmax, int jtime, double ttime, 
01184         int jgraph = 1, const char* title = 0x0, const char* filename_root = 0x0, 
01185         bool start_dx = false, int nu = 200, int nv = 200) const ;   
01186 
01208     void visu_box(double xmin, double xmax, double ymin, double ymax,
01209         double zmin, double zmax, const char* title0 = 0x0, 
01210         const char* filename0 = 0x0, bool start_dx = true, int nx = 40, int ny = 40, 
01211         int nz = 40) const ;      
01212         
01213 
01214 
01215   // Member arithmetics
01216   // ------------------
01217  public:
01218   void operator+=(const Scalar &) ;         
01219   void operator-=(const Scalar &) ;         
01220   void operator*=(const Scalar &) ;         
01221 
01222   // Manipulation of spectral bases
01223   // ------------------------------    
01227   virtual void std_spectral_base() ;     
01228   
01232   virtual void std_spectral_base_odd() ;     
01233   
01236   void set_spectral_base(const Base_val& ) ;     
01237 
01239   const Base_val& get_spectral_base( ) const {return va.base ;} ;    
01240 
01246   void set_dzpuis(int ) ; 
01247 
01263   Valeur** asymptot(int n, const int flag = 0) const ; 
01264     
01265 
01266   // PDE resolution 
01267   // --------------
01268  public:
01278   Scalar poisson() const ;
01279 
01291   void poisson(Param& par, Scalar& uu) const ;
01292       
01302   Scalar poisson_tau() const ;
01303    
01314   void poisson_tau(Param& par, Scalar& uu) const ;
01315 
01330   Scalar poisson_dirichlet (const Valeur& limite, int num) const ;
01331     
01336   Scalar poisson_neumann   (const Valeur&, int) const ;
01337 
01338 
01357   Scalar poisson_dir_neu  (const Valeur& limite , int num, 
01358                double fact_dir, double fact_neu) const ;
01359 
01366   Scalar poisson_frontiere_double   (const Valeur&, const Valeur&, int) const ;
01367 
01392   void poisson_regular(int k_div, int nzet, double unsgam1, Param& par,
01393                Scalar& uu, Scalar& uu_regu, Scalar& uu_div,
01394                Tensor& duu_div,
01395                Scalar& source_regu, Scalar& source_div) const ;
01396 
01431   Tbl test_poisson(const Scalar& uu, ostream& ostr, 
01432            bool detail = false) const ;  
01433 
01447     Scalar poisson_angu(double lambda =0) const ;
01448 
01481   Scalar avance_dalembert(Param& par, const Scalar& fJm1, const Scalar& source) 
01482     const ;
01483 
01488   Scalar sol_elliptic(Param_elliptic& params) const ;
01489  
01499   Scalar sol_elliptic_boundary(Param_elliptic& params, const Mtbl_cf& bound,
01500                    double fact_dir, double fact_neu) const ;
01501  
01506   Scalar sol_elliptic_boundary(Param_elliptic& params, const Scalar& bound,
01507                    double fact_dir, double fact_neu) const ;
01508  
01509 
01516   Scalar sol_elliptic_2d(Param_elliptic&) const ;
01517  
01522   Scalar sol_elliptic_pseudo_1d(Param_elliptic&) const ;
01523 
01531   Scalar sol_elliptic_no_zec(Param_elliptic& params, double val = 0) const ;
01532   
01539   Scalar sol_elliptic_only_zec(Param_elliptic& params, double val) const ;
01540 
01550   Scalar sol_elliptic_sin_zec(Param_elliptic& params, double* coefs, double* phases) const ;
01551 
01552    
01560   Scalar sol_elliptic_fixe_der_zero(double val, 
01561                     Param_elliptic& params) const ;
01562   
01563 
01574   Scalar sol_divergence(int n) const ;
01575     
01576   // Import from other mapping 
01577   // -------------------------
01578 
01584   void import(const Scalar& ci) ;    
01585 
01592   void import_symy(const Scalar& ci) ;   
01593 
01601   void import_asymy(const Scalar& ci) ;  
01602 
01614   void import(int nzet, const Scalar& ci) ;  
01615 
01628   void import_symy(int nzet, const Scalar& ci) ;     
01629 
01643   void import_asymy(int nzet, const Scalar& ci) ;    
01644 
01645  protected:
01659   void import_gal(int nzet, const Scalar& ci) ;  
01660 
01674   void import_align(int nzet, const Scalar& ci) ;    
01675 
01690   void import_anti(int nzet, const Scalar& ci) ;     
01691 
01706   void import_align_symy(int nzet, const Scalar& ci) ;   
01707 
01723   void import_anti_symy(int nzet, const Scalar& ci) ;    
01724 
01740   void import_align_asymy(int nzet, const Scalar& ci) ;  
01741 
01758   void import_anti_asymy(int nzet, const Scalar& ci) ;   
01759     
01760 
01761   friend Scalar operator-(const Scalar& ) ;         
01762   friend Scalar operator+(const Scalar&, const Scalar &) ;  
01763   friend Scalar operator+(const Scalar&, const Mtbl&) ; 
01764   friend Scalar operator+(const Scalar&, double ) ;     
01765   friend Scalar operator-(const Scalar &, const Scalar &) ;
01766   friend Scalar operator-(const Scalar&, const Mtbl&) ; 
01767   friend Scalar operator-(const Scalar&, double ) ;     
01768   friend Scalar operator*(const Scalar &, const Scalar &) ;
01769   friend Scalar operator%(const Scalar &, const Scalar &) ;
01770   friend Scalar operator|(const Scalar &, const Scalar &) ;
01771   friend Scalar operator*(const Mtbl&, const Scalar &) ;        
01772   friend Scalar operator*(double, const Scalar &) ;     
01773   friend Scalar operator/(const Scalar &, const Scalar &) ;
01774   friend Scalar operator/(const Scalar &, const Mtbl &) ;
01775   friend Scalar operator/(const Mtbl &, const Scalar &) ;
01776   friend Scalar operator/(const Scalar&, double ) ;        
01777   friend Scalar operator/(double, const Scalar &) ;
01778 
01779   friend Scalar sin(const Scalar& ) ;
01780   friend Scalar cos(const Scalar& ) ;
01781   friend Scalar tan(const Scalar& ) ;
01782   friend Scalar asin(const Scalar& ) ;
01783   friend Scalar acos(const Scalar& ) ;
01784   friend Scalar atan(const Scalar& ) ;
01785   friend Scalar exp(const Scalar& ) ;   
01786   friend Scalar Heaviside(const Scalar& ) ; 
01787   friend Scalar log(const Scalar& ) ;   
01788   friend Scalar log10(const Scalar& ) ; 
01789   friend Scalar sqrt(const Scalar& ) ;  
01790   friend Scalar racine_cubique (const Scalar& ) ;
01791   friend Scalar pow(const Scalar& , int ) ; 
01792   friend Scalar pow(const Scalar& , double ) ; 
01793   friend Scalar abs(const Scalar& ) ;   
01794 
01795   friend double totalmax(const Scalar& ) ;   
01796   friend double totalmin(const Scalar& ) ;   
01797   friend Tbl max(const Scalar& ) ;   
01798   friend Tbl min(const Scalar& ) ;   
01799   friend Tbl norme(const Scalar& ) ;   
01800   friend Tbl diffrel(const Scalar& a, const Scalar& b) ; 
01801   friend Tbl diffrelmax(const Scalar& a, const Scalar& b) ; 
01802 
01803 };
01804 
01805 ostream& operator<<(ostream& , const Scalar & ) ;   
01806 
01807 // Prototypage de l'arithmetique
01814 Scalar operator+(const Scalar& ) ;          
01815 Scalar operator-(const Scalar& ) ;          
01816 Scalar operator+(const Scalar&, const Scalar &) ;   
01817 Scalar operator+(const Scalar&, const Mtbl&) ;  
01818 Scalar operator+(const Mtbl&, const Scalar&) ;  
01819 Scalar operator+(const Scalar&, double ) ;      
01820 Scalar operator+(double, const Scalar& ) ;      
01821 Scalar operator+(const Scalar&, int ) ;     
01822 Scalar operator+(int, const Scalar& ) ;     
01823 Scalar operator-(const Scalar &, const Scalar &) ;  
01824 Scalar operator-(const Scalar&, const Mtbl&) ;  
01825 Scalar operator-(const Mtbl&, const Scalar&) ;  
01826 Scalar operator-(const Scalar&, double ) ;      
01827 Scalar operator-(double, const Scalar& ) ;      
01828 Scalar operator-(const Scalar&, int ) ;     
01829 Scalar operator-(int, const Scalar& ) ;     
01830 Scalar operator*(const Scalar &, const Scalar &) ;  
01831 
01833 Scalar operator%(const Scalar &, const Scalar &) ;  
01834 
01836 Scalar operator|(const Scalar &, const Scalar &) ;
01837 
01838 Scalar operator*(const Mtbl&, const Scalar&) ; 
01839 Scalar operator*(const Scalar&, const Mtbl&) ; 
01840     
01841 Scalar operator*(const Scalar&, double ) ;      
01842 Scalar operator*(double, const Scalar &) ;      
01843 Scalar operator*(const Scalar&, int ) ;     
01844 Scalar operator*(int, const Scalar& ) ;     
01845 Scalar operator/(const Scalar &, const Scalar &) ;  
01846 Scalar operator/(const Scalar&, double ) ;      
01847 Scalar operator/(double, const Scalar &) ;      
01848 Scalar operator/(const Scalar&, int ) ;     
01849 Scalar operator/(int, const Scalar &) ;     
01850 Scalar operator/(const Scalar &, const Mtbl&) ; 
01851 Scalar operator/(const Mtbl&, const Scalar &) ; 
01852     
01853 
01854 Scalar sin(const Scalar& ) ;        
01855 Scalar cos(const Scalar& ) ;        
01856 Scalar tan(const Scalar& ) ;        
01857 Scalar asin(const Scalar& ) ;       
01858 Scalar acos(const Scalar& ) ;       
01859 Scalar atan(const Scalar& ) ;       
01860 Scalar exp(const Scalar& ) ;        
01861 Scalar Heaviside(const Scalar& ) ;      
01862 Scalar log(const Scalar& ) ;        
01863 Scalar log10(const Scalar& ) ;  
01864 Scalar sqrt(const Scalar& ) ;       
01865 Scalar racine_cubique (const Scalar& ) ;        
01866 Scalar pow(const Scalar& , int ) ;  
01867 Scalar pow(const Scalar& , double ) ; 
01868 Scalar abs(const Scalar& ) ;        
01869 
01875 double totalmax(const Scalar& ) ;   
01876 
01882 double totalmin(const Scalar& ) ;   
01883 
01889 Tbl max(const Scalar& ) ;   
01890 
01896 Tbl min(const Scalar& ) ;   
01897 
01904 Tbl norme(const Scalar& ) ;   
01905 
01914 Tbl diffrel(const Scalar& a, const Scalar& b) ; 
01915 
01924 Tbl diffrelmax(const Scalar& a, const Scalar& b) ; 
01925 
01930 void exp_filter_ylm_all_domains(Scalar& ss, int p, double alpha=-16.) ;
01931 
01933 #endif

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