LORENE
scalar.h
1 /*
2  * Definition of Lorene class Scalar
3  *
4  */
5 
6 /*
7  * Copyright (c) 2003-2004 Eric Gourgoulhon & Jerome Novak
8  *
9  * Copyright (c) 1999-2000 Jean-Alain Marck (for previous class Cmp)
10  * Copyright (c) 1999-2002 Eric Gourgoulhon (for previous class Cmp)
11  * Copyright (c) 1999-2001 Philippe Grandclement (for previous class Cmp)
12  * Copyright (c) 2000-2002 Jerome Novak (for previous class Cmp)
13  * Copyright (c) 2000-2001 Keisuke Taniguchi (for previous class Cmp)
14  *
15  * This file is part of LORENE.
16  *
17  * LORENE is free software; you can redistribute it and/or modify
18  * it under the terms of the GNU General Public License as published by
19  * the Free Software Foundation; either version 2 of the License, or
20  * (at your option) any later version.
21  *
22  * LORENE is distributed in the hope that it will be useful,
23  * but WITHOUT ANY WARRANTY; without even the implied warranty of
24  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
25  * GNU General Public License for more details.
26  *
27  * You should have received a copy of the GNU General Public License
28  * along with LORENE; if not, write to the Free Software
29  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
30  *
31  */
32 
33 
34 #ifndef __SCALAR_H_
35 #define __SCALAR_H_
36 
37 
38 /*
39  * $Id: scalar.h,v 1.98 2025/03/31 08:29:50 j_novak Exp $
40  * $Log: scalar.h,v $
41  * Revision 1.98 2025/03/31 08:29:50 j_novak
42  * Methods to copyt Tensors and Scalars to smaller grids...
43  *
44  * Revision 1.97 2025/03/28 09:12:50 j_novak
45  * News methods to copy Tensors and Scalars to larger grids. To be used for de-aliasing.
46  *
47  * Revision 1.96 2023/08/31 08:27:26 g_servignat
48  * Added the possibility to filter in the r direction within the ylm filter. An order filtering of 0 means no filtering (for all 3 directions).
49  *
50  * Revision 1.95 2023/08/28 09:53:33 g_servignat
51  * Added ylm filter for Tensor and Scalar in theta + phi directions
52  *
53  * Revision 1.94 2015/12/18 15:52:51 j_novak
54  * New method is_nan() for class Scalar.
55  *
56  * Revision 1.93 2015/09/10 13:28:05 j_novak
57  * New methods for the class Hot_Eos
58  *
59  * Revision 1.92 2014/10/13 08:52:36 j_novak
60  * Lorene classes and functions now belong to the namespace Lorene.
61  *
62  * Revision 1.91 2013/06/05 15:06:10 j_novak
63  * Legendre bases are treated as standard bases, when the multi-grid
64  * (Mg3d) is built with BASE_LEG.
65  *
66  * Revision 1.90 2013/01/11 15:44:54 j_novak
67  * Addition of Legendre bases (part 2).
68  *
69  * Revision 1.89 2012/12/19 13:59:56 j_penner
70  * added a few lines to the documentation for scalar_match_tau function
71  *
72  * Revision 1.88 2012/01/17 15:05:46 j_penner
73  * *** empty log message ***
74  *
75  * Revision 1.87 2012/01/17 10:16:27 j_penner
76  * functions added: sarra_filter_r, sarra_filter_r_all_domains, Heaviside
77  *
78  * Revision 1.86 2011/04/08 13:13:09 e_gourgoulhon
79  * Changed the comment of function val_point to indicate specifically the
80  * division by r^dzpuis in the compactified external domain.
81  *
82  * Revision 1.85 2008/09/29 13:23:51 j_novak
83  * Implementation of the angular mapping associated with an affine
84  * mapping. Things must be improved to take into account the domain index.
85  *
86  * Revision 1.84 2008/09/22 19:08:01 j_novak
87  * New methods to deal with boundary conditions
88  *
89  * Revision 1.83 2008/05/24 15:05:22 j_novak
90  * New method Scalar::match_tau to match the output of an explicit time-marching scheme with the tau method.
91  *
92  * Revision 1.82 2007/12/21 16:06:16 j_novak
93  * Methods to filter Tensor, Vector and Sym_tensor objects.
94  *
95  * Revision 1.81 2007/10/31 10:33:11 j_novak
96  * Added exponential filters to smooth Gibbs-type phenomena.
97  *
98  * Revision 1.80 2007/06/21 19:56:36 k_taniguchi
99  * Introduction of another filtre_r.
100  *
101  * Revision 1.79 2007/05/06 10:48:08 p_grandclement
102  * Modification of a few operators for the vorton project
103  *
104  * Revision 1.78 2007/01/16 15:05:59 n_vasset
105  * New constructor (taking a Scalar in mono-domain angular grid for
106  * boundary) for function sol_elliptic_boundary
107  *
108  * Revision 1.77 2006/05/26 09:00:09 j_novak
109  * New members for multiplication or division by cos(theta).
110  *
111  * Revision 1.76 2005/11/30 13:48:06 e_gourgoulhon
112  * Replaced M_PI/2 by 1.57... in argument list of sol_elliptic_sin_zec
113  * (in order not to require the definition of M_PI).
114  *
115  * Revision 1.75 2005/11/30 11:09:03 p_grandclement
116  * Changes for the Bin_ns_bh project
117  *
118  * Revision 1.74 2005/11/17 15:29:46 e_gourgoulhon
119  * Added arithmetics with Mtbl.
120  *
121  * Revision 1.73 2005/10/25 08:56:34 p_grandclement
122  * addition of std_spectral_base in the case of odd functions near the origin
123  *
124  * Revision 1.72 2005/09/07 13:10:47 j_novak
125  * Added a filter setting to zero all mulitpoles in a given range.
126  *
127  * Revision 1.71 2005/08/26 14:02:38 p_grandclement
128  * Modification of the elliptic solver that matches with an oscillatory exterior solution
129  * small correction in Poisson tau also...
130  *
131  * Revision 1.70 2005/08/25 12:14:07 p_grandclement
132  * Addition of a new method to solve the scalar Poisson equation, based on a multi-domain Tau-method
133  *
134  * Revision 1.69 2005/06/09 07:56:25 f_limousin
135  * Implement a new function sol_elliptic_boundary() and
136  * Vector::poisson_boundary(...) which solve the vectorial poisson
137  * equation (method 6) with an inner boundary condition.
138  *
139  * Revision 1.68 2005/06/08 12:35:18 j_novak
140  * New method for solving divergence-like ODEs.
141  *
142  * Revision 1.67 2005/05/20 14:42:30 j_novak
143  * Added the method Scalar::get_spectral_base().
144  *
145  * Revision 1.66 2005/04/04 21:28:57 e_gourgoulhon
146  * Added argument lambda to method poisson_angu
147  * to deal with the generalized angular Poisson equation:
148  * Lap_ang u + lambda u = source.
149  *
150  * Revision 1.65 2004/12/14 09:09:39 f_limousin
151  * Modif. comments.
152  *
153  * Revision 1.64 2004/11/23 12:41:53 f_limousin
154  * Intoduce function poisson_dir_neu(...) to solve a scalar poisson
155  * equation with a mixed boundary condition (Dirichlet + Neumann).
156  *
157  * Revision 1.63 2004/10/11 15:09:00 j_novak
158  * The radial manipulation functions take Scalar as arguments, instead of Cmp.
159  * Added a conversion operator from Scalar to Cmp.
160  * The Cmp radial manipulation function make conversion to Scalar, call to the
161  * Map_radial version with a Scalar argument and back.
162  *
163  * Revision 1.62 2004/08/24 09:14:40 p_grandclement
164  * Addition of some new operators, like Poisson in 2d... It now requieres the
165  * GSL library to work.
166  *
167  * Also, the way a variable change is stored by a Param_elliptic is changed and
168  * no longer uses Change_var but rather 2 Scalars. The codes using that feature
169  * will requiere some modification. (It should concern only the ones about monopoles)
170  *
171  * Revision 1.61 2004/07/27 08:24:26 j_novak
172  * Modif. comments
173  *
174  * Revision 1.60 2004/07/26 16:02:21 j_novak
175  * Added a flag to specify whether the primitive should be zero either at r=0
176  * or at r going to infinity.
177  *
178  * Revision 1.59 2004/07/06 13:36:27 j_novak
179  * Added methods for desaliased product (operator |) only in r direction.
180  *
181  * Revision 1.58 2004/06/22 08:49:57 p_grandclement
182  * Addition of everything needed for using the logarithmic mapping
183  *
184  * Revision 1.57 2004/06/14 15:24:23 e_gourgoulhon
185  * Added method primr (radial primitive).
186  *
187  * Revision 1.56 2004/06/11 14:29:56 j_novak
188  * Scalar::multipole_spectrum() is now a const method.
189  *
190  * Revision 1.55 2004/05/24 14:07:31 e_gourgoulhon
191  * Method set_domain now includes a call to del_deriv() for safety.
192  *
193  * Revision 1.54 2004/05/07 11:26:10 f_limousin
194  * New method filtre_r(int*)
195  *
196  * Revision 1.53 2004/03/22 13:12:43 j_novak
197  * Modification of comments to use doxygen instead of doc++
198  *
199  * Revision 1.52 2004/03/17 15:58:47 p_grandclement
200  * Slight modification of sol_elliptic_no_zec
201  *
202  * Revision 1.51 2004/03/11 12:07:30 e_gourgoulhon
203  * Added method visu_section_anim.
204  *
205  * Revision 1.50 2004/03/08 15:45:38 j_novak
206  * Modif. comment
207  *
208  * Revision 1.49 2004/03/05 15:09:40 e_gourgoulhon
209  * Added method smooth_decay.
210  *
211  * Revision 1.48 2004/03/01 09:57:02 j_novak
212  * the wave equation is solved with Scalars. It now accepts a grid with a
213  * compactified external domain, which the solver ignores and where it copies
214  * the values of the field from one time-step to the next.
215  *
216  * Revision 1.47 2004/02/27 09:43:58 f_limousin
217  * New methods filtre_phi(int) and filtre_theta(int).
218  *
219  * Revision 1.46 2004/02/26 22:46:26 e_gourgoulhon
220  * Added methods derive_cov, derive_con and derive_lie.
221  *
222  * Revision 1.45 2004/02/21 17:03:49 e_gourgoulhon
223  * -- Method "point" renamed "val_grid_point".
224  * -- Method "set_point" renamed "set_grid_point".
225  *
226  * Revision 1.44 2004/02/19 22:07:35 e_gourgoulhon
227  * Added argument "comment" in method spectral_display.
228  *
229  * Revision 1.43 2004/02/11 09:47:44 p_grandclement
230  * Addition of a new elliptic solver, matching with the homogeneous solution
231  * at the outer shell and not solving in the external domain (more details
232  * coming soon ; check your local Lorene dealer...)
233  *
234  * Revision 1.42 2004/01/28 16:46:22 p_grandclement
235  * Addition of the sol_elliptic_fixe_der_zero stuff
236  *
237  * Revision 1.41 2004/01/28 13:25:40 j_novak
238  * The ced_mult_r arguments have been suppressed from the Scalar::*dsd* methods.
239  * In the div/mult _r_dzpuis, there is no more default value.
240  *
241  * Revision 1.40 2004/01/28 10:39:17 j_novak
242  * Comments modified.
243  *
244  * Revision 1.39 2004/01/27 15:10:01 j_novak
245  * New methods Scalar::div_r_dzpuis(int) and Scalar_mult_r_dzpuis(int)
246  * which replace div_r_inc*. Tried to clean the dzpuis handling.
247  * WARNING: no testing at this point!!
248  *
249  * Revision 1.38 2004/01/23 13:25:44 e_gourgoulhon
250  * Added methods set_inner_boundary and set_outer_boundary.
251  * Methods set_val_inf and set_val_hor, which are particular cases of
252  * the above, have been suppressed.
253  *
254  * Revision 1.37 2004/01/22 16:10:09 e_gourgoulhon
255  * Added (provisory) method div_r_inc().
256  *
257  * Revision 1.36 2003/12/16 06:32:20 e_gourgoulhon
258  * Added method visu_box.
259  *
260  * Revision 1.35 2003/12/14 21:46:35 e_gourgoulhon
261  * Added argument start_dx in visu_section.
262  *
263  * Revision 1.34 2003/12/11 16:19:38 e_gourgoulhon
264  * Added method visu_section for visualization with OpenDX.
265  *
266  * Revision 1.33 2003/12/11 14:48:47 p_grandclement
267  * Addition of ALL (and that is a lot !) the files needed for the general elliptic solver ... UNDER DEVELOPEMENT...
268  *
269  * Revision 1.32 2003/11/13 13:43:53 p_grandclement
270  * Addition of things needed for Bhole::update_metric (const Etoile_bin&, double, double)
271  *
272  * Revision 1.31 2003/11/06 14:43:37 e_gourgoulhon
273  * Gave a name to const arguments in certain method prototypes (e.g.
274  * constructors) to correct a bug of DOC++.
275  *
276  * Revision 1.30 2003/11/04 22:55:50 e_gourgoulhon
277  * Added new methods mult_cost(), mult_sint() and div_sint().
278  *
279  * Revision 1.29 2003/10/29 13:09:11 e_gourgoulhon
280  * -- Added integer argument to derivative functions dsdr, etc...
281  * so that one can choose the dzpuis of the result (default=2).
282  * -- Change of method name: laplacien --> laplacian.
283  *
284  * Revision 1.28 2003/10/29 11:00:42 e_gourgoulhon
285  * Virtual functions dec_dzpuis and inc_dzpuis have now an integer argument to
286  * specify by which amount dzpuis is to be increased.
287  * Accordingly virtual methods dec2_dzpuis and inc2_dzpuis have been suppressed.
288  *
289  * Revision 1.27 2003/10/20 14:26:02 j_novak
290  * New assignement operators.
291  *
292  * Revision 1.26 2003/10/19 19:46:33 e_gourgoulhon
293  * -- Method spectral_display now virtual (from Tensor), list of argument
294  * changed.
295  *
296  * Revision 1.25 2003/10/17 13:46:14 j_novak
297  * The argument is now between 1 and 3 (instead of 0->2)
298  *
299  * Revision 1.24 2003/10/16 15:23:41 e_gourgoulhon
300  * Name of method div_r_ced() changed to div_r_inc2().
301  * Name of method div_rsint_ced() changed to div_rsint_inc2().
302  *
303  * Revision 1.23 2003/10/15 21:10:11 e_gourgoulhon
304  * Added method poisson_angu().
305  *
306  * Revision 1.22 2003/10/15 16:03:35 j_novak
307  * Added the angular Laplace operator for Scalar.
308  *
309  * Revision 1.21 2003/10/15 10:29:05 e_gourgoulhon
310  * Added new members p_dsdt and p_stdsdp.
311  * Added new methods dsdt(), stdsdp() and div_tant().
312  *
313  * Revision 1.20 2003/10/13 13:52:39 j_novak
314  * Better managment of derived quantities.
315  *
316  * Revision 1.19 2003/10/10 15:57:27 j_novak
317  * Added the state one (ETATUN) to the class Scalar
318  *
319  * Revision 1.18 2003/10/08 14:24:08 j_novak
320  * replaced mult_r_zec with mult_r_ced
321  *
322  * Revision 1.17 2003/10/06 16:16:02 j_novak
323  * New constructor from a Tensor.
324  *
325  * Revision 1.16 2003/10/06 13:58:45 j_novak
326  * The memory management has been improved.
327  * Implementation of the covariant derivative with respect to the exact Tensor
328  * type.
329  *
330  * Revision 1.15 2003/10/05 21:06:31 e_gourgoulhon
331  * - Added new methods div_r_ced() and div_rsint_ced().
332  * - Added new virtual method std_spectral_base()
333  * - Removed method std_spectral_base_scal()
334  *
335  * Revision 1.14 2003/10/01 13:02:58 e_gourgoulhon
336  * Suppressed the constructor from Map* .
337  *
338  * Revision 1.13 2003/09/29 12:52:56 j_novak
339  * Methods for changing the triad are implemented.
340  *
341  * Revision 1.12 2003/09/25 09:33:36 j_novak
342  * Added methods for integral calculation and various manipulations
343  *
344  * Revision 1.11 2003/09/25 09:11:21 e_gourgoulhon
345  * Added functions for radial operations (divr, etc...)
346  *
347  * Revision 1.10 2003/09/25 08:55:23 e_gourgoulhon
348  * Added members raccord*.
349  *
350  * Revision 1.9 2003/09/25 08:50:11 j_novak
351  * Added the members import
352  *
353  * Revision 1.8 2003/09/25 08:13:51 j_novak
354  * Added method for calculating derivatives
355  *
356  * Revision 1.7 2003/09/25 07:59:26 e_gourgoulhon
357  * Added prototypes for PDE resolutions.
358  *
359  * Revision 1.6 2003/09/25 07:17:58 j_novak
360  * Method asymptot implemented.
361  *
362  * Revision 1.5 2003/09/24 20:53:38 e_gourgoulhon
363  * Added -- constructor by conversion from a Cmp
364  * -- assignment from Cmp
365  *
366  * Revision 1.4 2003/09/24 15:10:54 j_novak
367  * Suppression of the etat flag in class Tensor (still present in Scalar)
368  *
369  * Revision 1.3 2003/09/24 12:01:44 j_novak
370  * Added friend functions for math.
371  *
372  * Revision 1.2 2003/09/24 10:22:01 e_gourgoulhon
373  * still in progress...
374  *
375  * Revision 1.1 2003/09/22 12:50:47 e_gourgoulhon
376  * First version: not ready yet!
377  *
378  *
379  * $Header: /cvsroot/Lorene/C++/Include/scalar.h,v 1.98 2025/03/31 08:29:50 j_novak Exp $
380  *
381  */
382 
383 // Headers Lorene
384 #include "valeur.h"
385 #include "tensor.h"
386 
387 namespace Lorene {
388 class Param ;
389 class Cmp ;
390 class Param_elliptic ;
391 
399 class Scalar : public Tensor {
400 
401  // Data :
402  // -----
403  protected:
404 
408  int etat ;
409 
415  int dzpuis ;
416 
418 
419  // Derived data :
420  // ------------
421  protected:
423  mutable Scalar* p_dsdr ;
424 
428  mutable Scalar* p_srdsdt ;
429 
433  mutable Scalar* p_srstdsdp ;
434 
436  mutable Scalar* p_dsdt ;
437 
441  mutable Scalar* p_stdsdp ;
442 
446  mutable Scalar* p_dsdx ;
447 
451  mutable Scalar* p_dsdy ;
452 
456  mutable Scalar* p_dsdz ;
457 
460  mutable Scalar* p_lap ;
461 
464  mutable Scalar* p_lapang ;
465 
467  mutable Scalar* p_dsdradial ;
468 
470  mutable Scalar* p_dsdrho ;
471 
475  mutable int ind_lap ;
476 
480  mutable Tbl* p_integ ;
481 
482  // Constructors - Destructor
483  // -------------------------
484 
485  public:
486 
487  explicit Scalar(const Map& mpi) ;
488 
490  Scalar(const Tensor& a) ;
491 
492  Scalar(const Scalar& a) ;
493 
494  Scalar(const Cmp& a) ;
495 
497  Scalar(const Map&, const Mg3d&, FILE* ) ;
498 
499  virtual ~Scalar() ;
500 
501 
502  // Memory management
503  // -----------------
504  protected:
505  void del_t() ;
506  virtual void del_deriv() const;
507  void set_der_0x0() const;
508 
509  public:
510 
516  virtual void set_etat_nondef() ;
517 
523  virtual void set_etat_zero() ;
524 
531  virtual void set_etat_qcq() ;
532 
538  void set_etat_one() ;
539 
548  virtual void allocate_all() ;
549 
558  void annule_hard() ;
559 
560  // Extraction of information
561  // -------------------------
562  public:
566  int get_etat() const {return etat;} ;
567 
569  int get_dzpuis() const {return dzpuis;} ;
570 
574  bool dz_nonzero() const ;
575 
580  bool check_dzpuis(int dzi) const ;
581 
587  bool is_nan(bool verbose=false) const ;
588 
589  // Assignment
590  // -----------
591  public:
593  void operator=(const Scalar& a) ;
594 
596  virtual void operator=(const Tensor& a) ;
597 
598  void operator=(const Cmp& a) ;
599  void operator=(const Valeur& a) ;
600  void operator=(const Mtbl& a) ;
601  void operator=(double ) ;
602  void operator=(int ) ;
603 
604  // Conversion oprators
605  //---------------------
606  operator Cmp() const ;
607 
608  // Access to individual elements
609  // -----------------------------
610  public:
611 
613  const Valeur& get_spectral_va() const {return va;} ;
614 
616  Valeur& set_spectral_va() {return va;} ;
617 
627  Tbl& set_domain(int l) {
628  assert(etat == ETATQCQ) ;
629  del_deriv() ;
630  return va.set(l) ;
631  };
632 
637  const Tbl& domain(int l) const {
638  assert( (etat == ETATQCQ) || (etat == ETATUN) ) ;
639  return va(l) ;
640  };
641 
642 
649  double val_grid_point(int l, int k, int j, int i) const {
650  assert(etat != ETATNONDEF) ;
651  if (etat == ETATZERO) {
652  double zero = 0. ;
653  return zero ;
654  }
655  else {
656  if (etat == ETATUN) {
657  double one = 1. ;
658  return one ;
659  }
660  else{
661  return va(l, k, j, i) ;
662  }
663  }
664  };
665 
680  double val_point(double r, double theta, double phi) const ;
681 
682 
696  double& set_grid_point(int l, int k, int j, int i) {
697  assert(etat == ETATQCQ) ;
698  return va.set(l, k, j, i) ;
699  };
700 
701 
712  virtual void annule(int l_min, int l_max) ;
713 
719  void set_inner_boundary(int l, double x) ;
720 
726  void set_outer_boundary(int l, double x) ;
727 
736  Tbl multipole_spectrum () const ;
737 
744  Tbl tbl_out_bound(int l_dom, bool leave_ylm = false) ;
745 
752  Tbl tbl_in_bound(int n, bool leave_ylm = false) ;
753 
760  Scalar scalar_out_bound(int n, bool leave_ylm = false) ;
761 
762  // Differential operators and others
763  // ---------------------------------
764  public:
769  const Scalar& dsdr() const ;
770 
775  const Scalar& srdsdt() const ;
776 
781  const Scalar& srstdsdp() const ;
782 
785  const Scalar& dsdt() const ;
786 
794  const Scalar& dsdradial() const ;
795 
800  const Scalar& dsdrho() const ;
801 
804  const Scalar& stdsdp() const ;
805 
811  const Scalar& dsdx() const ;
812 
818  const Scalar& dsdy() const ;
819 
825  const Scalar& dsdz() const ;
826 
833  const Scalar& deriv(int i) const ;
834 
839  const Vector& derive_cov(const Metric& gam) const ;
840 
841 
847  const Vector& derive_con(const Metric& gam) const ;
848 
850  Scalar derive_lie(const Vector& v) const ;
851 
852 
861  const Scalar& laplacian(int ced_mult_r = 4) const ;
862 
870  const Scalar& lapang() const ;
871 
873  void div_r() ;
874 
879  void div_r_dzpuis(int ced_mult_r) ;
880 
884  void div_r_ced() ;
885 
887  void mult_r() ;
888 
893  void mult_r_dzpuis(int ced_mult_r) ;
894 
898  void mult_r_ced() ;
899 
901  void mult_rsint() ;
902 
907  void mult_rsint_dzpuis(int ced_mult_r) ;
908 
910  void div_rsint() ;
911 
916  void div_rsint_dzpuis(int ced_mult_r) ;
917 
918  void mult_cost() ;
919 
920  void div_cost() ;
921 
922  void mult_sint() ;
923 
924  void div_sint() ;
925 
926  void div_tant() ;
927 
938  Scalar primr(bool null_infty = true) const ;
939 
947  double integrale() const ;
948 
958  const Tbl& integrale_domains() const ;
959 
964  virtual void dec_dzpuis(int dec = 1) ;
965 
970  virtual void inc_dzpuis(int inc = 1) ;
971 
975  virtual void change_triad(const Base_vect& new_triad) ;
976 
977 
978  // Filters and de-aliasing tools
979  //---------------------------------
984  void filtre (int n) ;
985 
990  void filtre_r (int* nn) ;
991 
995  void filtre_r (int n, int nzone) ;
996 
1007 // virtual void exponential_filter_r(int lzmin, int lzmax, int p,
1008 // double alpha= -16.) ;
1009  virtual void exponential_filter_r(int lzmin, int lzmax, int p,
1010  double alpha= -16.) ;
1011 
1022  void sarra_filter_r(int lzmin, int lzmax, double p,
1023  double alpha= -1E-16) ;
1024 
1030  void exp_filter_r_all_domains(Scalar &ss, int p, double alpha=-16.) ;
1031 
1038  void sarra_filter_r_all_domains(double p, double alpha=1E-16) ;
1039 
1047  virtual void exponential_filter_ylm(int lzmin, int lzmax, int p,
1048  double alpha= -16.) ;
1049 
1057  virtual void exponential_filter_ylm_phi(int lzmin, int lzmax, int p_r, int p_tet, int p_phi,
1058  double alpha= -16.) ;
1059 
1067  void annule_l (int l_min, int l_max, bool ylm_output= false ) ;
1068 
1073  void filtre_phi (int n, int zone) ;
1074 
1079  void filtre_tp(int nn, int nz1, int nz2) ;
1080 
1089  void copy_coefs_from_small_grid(const Scalar& scal_small_grid) ;
1090 
1099  void copy_coefs_from_large_grid(const Scalar& scal_large_grid) ;
1100 
1106  void fixe_decroissance (int puis) ;
1107 
1113  void smooth_decay(int k, int n) ;
1114 
1119  void raccord(int n) ;
1120 
1127  void raccord_c1_zec(int puis, int nbre, int lmax) ;
1128 
1132  void raccord_externe(int puis, int nbre, int lmax) ;
1133 
1150  void match_tau(Param& par_bc, Param* par_mat=0x0) ;
1151 
1168  void match_tau_star(Param& par_bc, Param* par_mat=0x0) ;
1169 
1178  void match_tau_shell(Param& par_bc, Param* par_mat=0x0) ;
1179 
1188  void match_collocation(Param& par_bc, Param* par_mat=0x0) ;
1189 
1190  // Outputs
1191  // -------
1192  public:
1193  virtual void sauve(FILE *) const ;
1194 
1205  virtual void spectral_display(const char* comment = 0x0,
1206  double threshold = 1.e-7, int precision = 4,
1207  ostream& ostr = cout) const ;
1208 
1210  friend ostream& operator<<(ostream& , const Scalar & ) ;
1211 
1235  void visu_section(const char section_type, double aa, double umin, double umax, double vmin,
1236  double vmax, const char* title = 0x0, const char* filename = 0x0,
1237  bool start_dx = true, int nu = 200, int nv = 200) const ;
1238 
1267  void visu_section(const Tbl& plane, double umin, double umax, double vmin,
1268  double vmax, const char* title = 0x0, const char* filename = 0x0,
1269  bool start_dx = true, int nu = 200, int nv = 200) const ;
1270 
1299  void visu_section_anim(const char section_type, double aa, double umin,
1300  double umax, double vmin, double vmax, int jtime, double ttime,
1301  int jgraph = 1, const char* title = 0x0, const char* filename_root = 0x0,
1302  bool start_dx = false, int nu = 200, int nv = 200) const ;
1303 
1325  void visu_box(double xmin, double xmax, double ymin, double ymax,
1326  double zmin, double zmax, const char* title0 = 0x0,
1327  const char* filename0 = 0x0, bool start_dx = true, int nx = 40, int ny = 40,
1328  int nz = 40) const ;
1329 
1330 
1331 
1332  // Member arithmetics
1333  // ------------------
1334  public:
1335  void operator+=(const Scalar &) ;
1336  void operator-=(const Scalar &) ;
1337  void operator*=(const Scalar &) ;
1338 
1339  // Manipulation of spectral bases
1340  // ------------------------------
1344  virtual void std_spectral_base() ;
1345 
1349  virtual void std_spectral_base_odd() ;
1350 
1353  void set_spectral_base(const Base_val& ) ;
1354 
1356  const Base_val& get_spectral_base( ) const {return va.base ;} ;
1357 
1363  void set_dzpuis(int ) ;
1364 
1380  Valeur** asymptot(int n, const int flag = 0) const ;
1381 
1382 
1383  // PDE resolution
1384  // --------------
1385  public:
1395  Scalar poisson() const ;
1396 
1408  void poisson(Param& par, Scalar& uu) const ;
1409 
1419  Scalar poisson_tau() const ;
1420 
1431  void poisson_tau(Param& par, Scalar& uu) const ;
1432 
1447  Scalar poisson_dirichlet (const Valeur& limite, int num) const ;
1448 
1453  Scalar poisson_neumann (const Valeur&, int) const ;
1454 
1455 
1474  Scalar poisson_dir_neu (const Valeur& limite , int num,
1475  double fact_dir, double fact_neu) const ;
1476 
1483  Scalar poisson_frontiere_double (const Valeur&, const Valeur&, int) const ;
1484 
1509  void poisson_regular(int k_div, int nzet, double unsgam1, Param& par,
1510  Scalar& uu, Scalar& uu_regu, Scalar& uu_div,
1511  Tensor& duu_div,
1512  Scalar& source_regu, Scalar& source_div) const ;
1513 
1548  Tbl test_poisson(const Scalar& uu, ostream& ostr,
1549  bool detail = false) const ;
1550 
1564  Scalar poisson_angu(double lambda =0) const ;
1565 
1598  Scalar avance_dalembert(Param& par, const Scalar& fJm1, const Scalar& source)
1599  const ;
1600 
1605  Scalar sol_elliptic(Param_elliptic& params) const ;
1606 
1616  Scalar sol_elliptic_boundary(Param_elliptic& params, const Mtbl_cf& bound,
1617  double fact_dir, double fact_neu) const ;
1618 
1623  Scalar sol_elliptic_boundary(Param_elliptic& params, const Scalar& bound,
1624  double fact_dir, double fact_neu) const ;
1625 
1626 
1634 
1640 
1648  Scalar sol_elliptic_no_zec(Param_elliptic& params, double val = 0) const ;
1649 
1656  Scalar sol_elliptic_only_zec(Param_elliptic& params, double val) const ;
1657 
1667  Scalar sol_elliptic_sin_zec(Param_elliptic& params, double* coefs, double* phases) const ;
1668 
1669 
1677  Scalar sol_elliptic_fixe_der_zero(double val,
1678  Param_elliptic& params) const ;
1679 
1680 
1691  Scalar sol_divergence(int n) const ;
1692 
1693  // Import from other mapping
1694  // -------------------------
1695 
1701  void import(const Scalar& ci) ;
1702 
1709  void import_symy(const Scalar& ci) ;
1710 
1718  void import_asymy(const Scalar& ci) ;
1719 
1731  void import(int nzet, const Scalar& ci) ;
1732 
1745  void import_symy(int nzet, const Scalar& ci) ;
1746 
1760  void import_asymy(int nzet, const Scalar& ci) ;
1761 
1762  protected:
1776  void import_gal(int nzet, const Scalar& ci) ;
1777 
1791  void import_align(int nzet, const Scalar& ci) ;
1792 
1807  void import_anti(int nzet, const Scalar& ci) ;
1808 
1823  void import_align_symy(int nzet, const Scalar& ci) ;
1824 
1840  void import_anti_symy(int nzet, const Scalar& ci) ;
1841 
1857  void import_align_asymy(int nzet, const Scalar& ci) ;
1858 
1875  void import_anti_asymy(int nzet, const Scalar& ci) ;
1876 
1877 
1878  friend Scalar operator-(const Scalar& ) ;
1879  friend Scalar operator+(const Scalar&, const Scalar &) ;
1880  friend Scalar operator+(const Scalar&, const Mtbl&) ;
1881  friend Scalar operator+(const Scalar&, double ) ;
1882  friend Scalar operator-(const Scalar &, const Scalar &) ;
1883  friend Scalar operator-(const Scalar&, const Mtbl&) ;
1884  friend Scalar operator-(const Scalar&, double ) ;
1885  friend Scalar operator*(const Scalar &, const Scalar &) ;
1886  friend Scalar operator%(const Scalar &, const Scalar &) ;
1887  friend Scalar operator|(const Scalar &, const Scalar &) ;
1888  friend Scalar operator*(const Mtbl&, const Scalar &) ;
1889  friend Scalar operator*(double, const Scalar &) ;
1890  friend Scalar operator/(const Scalar &, const Scalar &) ;
1891  friend Scalar operator/(const Scalar &, const Mtbl &) ;
1892  friend Scalar operator/(const Mtbl &, const Scalar &) ;
1893  friend Scalar operator/(const Scalar&, double ) ;
1894  friend Scalar operator/(double, const Scalar &) ;
1895 
1896  friend Scalar sin(const Scalar& ) ;
1897  friend Scalar cos(const Scalar& ) ;
1898  friend Scalar tan(const Scalar& ) ;
1899  friend Scalar asin(const Scalar& ) ;
1900  friend Scalar acos(const Scalar& ) ;
1901  friend Scalar atan(const Scalar& ) ;
1902  friend Scalar exp(const Scalar& ) ;
1903  friend Scalar Heaviside(const Scalar& ) ;
1904  friend Scalar log(const Scalar& ) ;
1905  friend Scalar log10(const Scalar& ) ;
1906  friend Scalar sqrt(const Scalar& ) ;
1907  friend Scalar racine_cubique (const Scalar& ) ;
1908  friend Scalar pow(const Scalar& , int ) ;
1909  friend Scalar pow(const Scalar& , double ) ;
1910  friend Scalar abs(const Scalar& ) ;
1911 
1912  friend double totalmax(const Scalar& ) ;
1913  friend double totalmin(const Scalar& ) ;
1914  friend Tbl max(const Scalar& ) ;
1915  friend Tbl min(const Scalar& ) ;
1916  friend Tbl norme(const Scalar& ) ;
1917  friend Tbl diffrel(const Scalar& a, const Scalar& b) ;
1918  friend Tbl diffrelmax(const Scalar& a, const Scalar& b) ;
1919 
1920 };
1921 
1922 ostream& operator<<(ostream& , const Scalar & ) ;
1923 
1924 // Prototypage de l'arithmetique
1931 Scalar operator+(const Scalar& ) ;
1932 Scalar operator-(const Scalar& ) ;
1933 Scalar operator+(const Scalar&, const Scalar &) ;
1934 Scalar operator+(const Scalar&, const Mtbl&) ;
1935 Scalar operator+(const Mtbl&, const Scalar&) ;
1936 Scalar operator+(const Scalar&, double ) ;
1937 Scalar operator+(double, const Scalar& ) ;
1938 Scalar operator+(const Scalar&, int ) ;
1939 Scalar operator+(int, const Scalar& ) ;
1940 Scalar operator-(const Scalar &, const Scalar &) ;
1941 Scalar operator-(const Scalar&, const Mtbl&) ;
1942 Scalar operator-(const Mtbl&, const Scalar&) ;
1943 Scalar operator-(const Scalar&, double ) ;
1944 Scalar operator-(double, const Scalar& ) ;
1945 Scalar operator-(const Scalar&, int ) ;
1946 Scalar operator-(int, const Scalar& ) ;
1947 Scalar operator*(const Scalar &, const Scalar &) ;
1948 
1950 Scalar operator%(const Scalar &, const Scalar &) ;
1951 
1953 Scalar operator|(const Scalar &, const Scalar &) ;
1954 
1955 Scalar operator*(const Mtbl&, const Scalar&) ;
1956 Scalar operator*(const Scalar&, const Mtbl&) ;
1957 
1958 Scalar operator*(const Scalar&, double ) ;
1959 Scalar operator*(double, const Scalar &) ;
1960 Scalar operator*(const Scalar&, int ) ;
1961 Scalar operator*(int, const Scalar& ) ;
1962 Scalar operator/(const Scalar &, const Scalar &) ;
1963 Scalar operator/(const Scalar&, double ) ;
1964 Scalar operator/(double, const Scalar &) ;
1965 Scalar operator/(const Scalar&, int ) ;
1966 Scalar operator/(int, const Scalar &) ;
1967 Scalar operator/(const Scalar &, const Mtbl&) ;
1968 Scalar operator/(const Mtbl&, const Scalar &) ;
1969 
1970 
1971 Scalar sin(const Scalar& ) ;
1972 Scalar cos(const Scalar& ) ;
1973 Scalar tan(const Scalar& ) ;
1974 Scalar asin(const Scalar& ) ;
1975 Scalar acos(const Scalar& ) ;
1976 Scalar atan(const Scalar& ) ;
1977 Scalar exp(const Scalar& ) ;
1978 Scalar Heaviside(const Scalar& ) ;
1979 Scalar log(const Scalar& ) ;
1980 Scalar log10(const Scalar& ) ;
1981 Scalar sqrt(const Scalar& ) ;
1982 Scalar racine_cubique (const Scalar& ) ;
1983 Scalar pow(const Scalar& , int ) ;
1984 Scalar pow(const Scalar& , double ) ;
1985 Scalar abs(const Scalar& ) ;
1986 
1992 double totalmax(const Scalar& ) ;
1993 
1999 double totalmin(const Scalar& ) ;
2000 
2006 Tbl max(const Scalar& ) ;
2007 
2013 Tbl min(const Scalar& ) ;
2014 
2021 Tbl norme(const Scalar& ) ;
2022 
2031 Tbl diffrel(const Scalar& a, const Scalar& b) ;
2032 
2041 Tbl diffrelmax(const Scalar& a, const Scalar& b) ;
2042 
2047 void exp_filter_ylm_all_domains(Scalar& ss, int p, double alpha=-16.) ;
2048 
2053 void exp_filter_ylm_all_domains(Scalar& ss, int p_tet, int p_phi, double alpha=-16.) ;
2054 
2056 }
2057 #endif
const Base_val & get_spectral_base() const
Returns the spectral bases of the Valeur va.
Definition: scalar.h:1356
Scalar * p_srdsdt
Pointer on of *this (0x0 if not up to date)
Definition: scalar.h:428
void annule_l(int l_min, int l_max, bool ylm_output=false)
Sets all the multipolar components between l_min and l_max to zero.
Definition: scalar_manip.C:123
Metric for tensor calculation.
Definition: metric.h:90
Scalar * p_srstdsdp
Pointer on of *this (0x0 if not up to date)
Definition: scalar.h:433
Cmp log(const Cmp &)
Neperian logarithm.
Definition: cmp_math.C:299
Scalar poisson_dir_neu(const Valeur &limite, int num, double fact_dir, double fact_neu) const
Is identicall to Scalar::poisson() .
Cmp asin(const Cmp &)
Arcsine.
Definition: cmp_math.C:147
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Definition: cmp.h:446
Cmp exp(const Cmp &)
Exponential.
Definition: cmp_math.C:273
void filtre_r(int *nn)
Sets the n lasts coefficients in r to 0 in all domains.
Definition: scalar_manip.C:192
Scalar poisson_frontiere_double(const Valeur &, const Valeur &, int) const
Idem as Scalar::poisson_dirichlet , the boundary condition being on both the function and its radial ...
const Scalar & dsdradial() const
Returns of *this if the mapping is affine (class Map_af) and of *this if the mapping is logarithmic...
Definition: scalar_deriv.C:490
void mult_r()
Multiplication by r everywhere; dzpuis is not changed.
void copy_coefs_from_small_grid(const Scalar &scal_small_grid)
Copies the content of the argument Scalar to this defined on a larger grid, with similar mappings...
Definition: scalar_manip.C:505
const Scalar & lapang() const
Returns the angular Laplacian of *this , where .
Definition: scalar_deriv.C:460
Cmp sqrt(const Cmp &)
Square root.
Definition: cmp_math.C:223
virtual void set_etat_zero()
Sets the logical state to ETATZERO (zero).
Definition: scalar.C:330
void set_der_0x0() const
Sets the pointers for derivatives to 0x0.
Definition: scalar.C:312
virtual void annule(int l_min, int l_max)
Sets the Scalar to zero in several domains.
Definition: scalar.C:397
Multi-domain array.
Definition: mtbl.h:118
const Tbl & domain(int l) const
Read-only of the value in a given domain.
Definition: scalar.h:637
friend Scalar exp(const Scalar &)
Exponential.
Definition: scalar_math.C:326
friend double totalmax(const Scalar &)
Maximum values of a Scalar in each domain.
Definition: scalar_math.C:556
Scalar(const Map &mpi)
Constructor from mapping.
Definition: scalar.C:210
Lorene prototypes.
Definition: app_hor.h:67
void div_cost()
Division by .
void visu_section_anim(const char section_type, double aa, double umin, double umax, double vmin, double vmax, int jtime, double ttime, int jgraph=1, const char *title=0x0, const char *filename_root=0x0, bool start_dx=false, int nu=200, int nv=200) const
3D visualization via time evolving plane section (animation).
Definition: scalar_visu.C:538
const Scalar & dsdt() const
Returns of *this .
Definition: scalar_deriv.C:207
const Scalar & laplacian(int ced_mult_r=4) const
Returns the Laplacian of *this.
Definition: scalar_deriv.C:435
Cmp racine_cubique(const Cmp &)
Cube root.
Definition: cmp_math.C:248
void visu_box(double xmin, double xmax, double ymin, double ymax, double zmin, double zmax, const char *title0=0x0, const char *filename0=0x0, bool start_dx=true, int nx=40, int ny=40, int nz=40) const
3D visualization (volume rendering) via OpenDX.
Definition: scalar_visu.C:345
Tbl multipole_spectrum() const
Gives the spectrum in terms of multipolar modes l .
Definition: scalar.C:962
const Scalar & dsdz() const
Returns of *this , where .
Definition: scalar_deriv.C:327
Tensor field of valence 0 (or component of a tensorial field).
Definition: scalar.h:399
virtual void change_triad(const Base_vect &new_triad)
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
Definition: scalar.C:1003
Scalar poisson() const
Solves the scalar Poisson equation with *this as a source.
Definition: scalar_pde.C:139
Base_val operator*(const Base_val &, const Base_val &)
This operator is used when calling multiplication or division of Valeur .
Base class for coordinate mappings.
Definition: map.h:697
void import_gal(int nzet, const Scalar &ci)
Assignment to another Scalar defined on a different mapping, when the two mappings do not have a part...
void operator-=(const Scalar &)
-= Scalar
double integrale() const
Computes the integral over all space of *this .
Definition: scalar_integ.C:67
friend Tbl max(const Scalar &)
Maximum values of a Scalar in each domain.
Definition: scalar_math.C:614
friend ostream & operator<<(ostream &, const Scalar &)
Display.
Definition: scalar.C:708
Scalar * p_dsdt
Pointer on of *this (0x0 if not up to date)
Definition: scalar.h:436
Scalar primr(bool null_infty=true) const
Computes the radial primitive which vanishes for .
Definition: scalar_integ.C:106
friend Scalar atan(const Scalar &)
Arctangent.
Definition: scalar_math.C:234
void match_tau(Param &par_bc, Param *par_mat=0x0)
Method for matching accross domains and imposing boundary condition.
void mult_sint()
Multiplication by .
virtual void exponential_filter_ylm_phi(int lzmin, int lzmax, int p_r, int p_tet, int p_phi, double alpha=-16.)
Applies an exponential filter to the spectral coefficients in the angular directions.
void poisson_regular(int k_div, int nzet, double unsgam1, Param &par, Scalar &uu, Scalar &uu_regu, Scalar &uu_div, Tensor &duu_div, Scalar &source_regu, Scalar &source_div) const
Solves the scalar Poisson equation with *this as a source (version with parameters to control the res...
int ind_lap
Power of r by which the last computed Laplacian has been multiplied in the compactified external doma...
Definition: scalar.h:475
const Scalar & dsdx() const
Returns of *this , where .
Definition: scalar_deriv.C:265
virtual void exponential_filter_r(int lzmin, int lzmax, int p, double alpha=-16.)
Applies an exponential filter to the spectral coefficients in the radial direction.
virtual void std_spectral_base()
Sets the spectral bases of the Valeur va to the standard ones for a scalar field. ...
Definition: scalar.C:790
void import_anti_symy(int nzet, const Scalar &ci)
Assignment to another Scalar defined on a different mapping, when the two mappings have anti-aligned ...
friend Tbl min(const Scalar &)
Minimum values of a Scalar in each domain.
Definition: scalar_math.C:642
Values and coefficients of a (real-value) function.
Definition: valeur.h:297
Cmp operator%(const Cmp &, const Cmp &)
Cmp * Cmp with desaliasing.
Definition: cmp_arithm.C:367
Scalar scalar_out_bound(int n, bool leave_ylm=false)
Returns the Scalar containing the values of angular coefficients at the outer boundary.
Definition: scalar_manip.C:479
virtual void allocate_all()
Sets the logical state to ETATQCQ (ordinary state) and performs the memory allocation of all the elem...
Definition: scalar.C:371
virtual void set_etat_nondef()
Sets the logical state to ETATNONDEF (undefined).
Definition: scalar.C:350
int dzpuis
Power of r by which the quantity represented by this must be divided in the compactified external do...
Definition: scalar.h:415
int get_etat() const
Returns the logical state ETATNONDEF (undefined), ETATZERO (null) or ETATQCQ (ordinary).
Definition: scalar.h:566
void operator*=(const Scalar &)
*= Scalar
Scalar sol_elliptic_fixe_der_zero(double val, Param_elliptic &params) const
Resolution of a general elliptic equation fixing the dericative at the origin and relaxing one contin...
Definition: scalar_pde.C:375
Cmp operator/(const Cmp &, const Cmp &)
Cmp / Cmp.
Definition: cmp_arithm.C:460
void operator+=(const Scalar &)
+= Scalar
Tbl min(const Cmp &)
Minimum values of a Cmp in each domain.
Definition: cmp_math.C:461
Scalar sol_divergence(int n) const
Resolution of a divergence-like equation.
virtual void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: scalar.C:359
friend Scalar sqrt(const Scalar &)
Square root.
Definition: scalar_math.C:266
Tbl tbl_in_bound(int n, bool leave_ylm=false)
Returns the Tbl containing the values of angular coefficients at the inner boundary.
Definition: scalar_manip.C:460
Scalar operator|(const Scalar &, const Scalar &)
Scalar * Scalar with desaliasing only in r.
Tbl test_poisson(const Scalar &uu, ostream &ostr, bool detail=false) const
Checks if a Poisson equation with *this as a source has been correctly solved.
const Vector & derive_con(const Metric &gam) const
Returns the "contravariant" derivative of *this with respect to some metric , by raising the index of...
Definition: scalar_deriv.C:401
Tensor field of valence 1.
Definition: vector.h:188
friend Scalar operator*(const Scalar &, const Scalar &)
Scalar * Scalar.
Scalar poisson_angu(double lambda=0) const
Solves the (generalized) angular Poisson equation with *this as source.
Definition: scalar_pde.C:189
Vectorial bases (triads) with respect to which the tensorial components are defined.
Definition: base_vect.h:105
Cmp cos(const Cmp &)
Cosine.
Definition: cmp_math.C:97
Scalar sol_elliptic_only_zec(Param_elliptic &params, double val) const
Resolution of a general elliptic equation solving in the compactified domain and putting a given valu...
Definition: scalar_pde.C:330
Tbl & set_domain(int l)
Read/write of the value in a given domain.
Definition: scalar.h:627
Tbl diffrel(const Cmp &a, const Cmp &b)
Relative difference between two Cmp (norme version).
Definition: cmp_math.C:507
friend Scalar acos(const Scalar &)
Arccosine.
Definition: scalar_math.C:202
void set_dzpuis(int)
Modifies the dzpuis flag.
Definition: scalar.C:814
double val_grid_point(int l, int k, int j, int i) const
Returns the value of the field at a specified grid point.
Definition: scalar.h:649
void annule_hard()
Sets the Scalar to zero in a hard way.
Definition: scalar.C:386
Scalar * p_lap
Pointer on the Laplacian of *this (0x0 if not up to date)
Definition: scalar.h:460
void exp_filter_ylm_all_domains(Scalar &ss, int p, double alpha=-16.)
Applies an exponential filter in angular directions in all domains.
void raccord(int n)
Performs the matching of the nucleus with respect to the first shell.
Cmp tan(const Cmp &)
Tangent.
Definition: cmp_math.C:123
Scalar avance_dalembert(Param &par, const Scalar &fJm1, const Scalar &source) const
Performs one time-step integration (from ) of the scalar d&#39;Alembert equation with *this being the val...
Definition: scalar_pde.C:206
friend Scalar pow(const Scalar &, int)
Power .
Definition: scalar_math.C:454
Tbl * p_integ
Pointer on the space integral of *this (values in each domain) (0x0 if not up to date) ...
Definition: scalar.h:480
void import_anti(int nzet, const Scalar &ci)
Assignment to another Scalar defined on a different mapping, when the two mappings have anti-aligned ...
void smooth_decay(int k, int n)
Performs a matching of the last non-compactified shell with a decaying function where is the spher...
friend Scalar cos(const Scalar &)
Cosine.
Definition: scalar_math.C:107
void div_r()
Division by r everywhere; dzpuis is not changed.
const Scalar & dsdrho() const
Returns of *this .
Definition: scalar_deriv.C:521
virtual void sauve(FILE *) const
Save in a file.
Definition: scalar.C:692
const Tbl & integrale_domains() const
Computes the integral in each domain of *this .
Definition: scalar_integ.C:85
virtual void del_deriv() const
Logical destructor of the derivatives.
Definition: scalar.C:293
Scalar * p_dsdradial
Pointer on of *this.
Definition: scalar.h:467
Tbl norme(const Cmp &)
Sums of the absolute values of all the values of the Cmp in each domain.
Definition: cmp_math.C:484
Cmp atan(const Cmp &)
Arctangent.
Definition: cmp_math.C:198
Scalar sol_elliptic_no_zec(Param_elliptic &params, double val=0) const
Resolution of a general elliptic equation, putting a given value at the outermost shell and not solvi...
Definition: scalar_pde.C:303
virtual void inc_dzpuis(int inc=1)
Increases by inc units the value of dzpuis and changes accordingly the values of the Scalar in the co...
Scalar * p_dsdrho
Pointer on of *this.
Definition: scalar.h:470
const Scalar & dsdy() const
Returns of *this , where .
Definition: scalar_deriv.C:296
Scalar * p_dsdy
Pointer on of *this , where (0x0 if not up to date)
Definition: scalar.h:451
void operator=(const Scalar &a)
Assignment to another Scalar defined on the same mapping.
Definition: scalar.C:452
Cmp operator+(const Cmp &)
Definition: cmp_arithm.C:107
Scalar * p_lapang
Pointer on the Laplacian of *this (0x0 if not up to date)
Definition: scalar.h:464
Scalar poisson_neumann(const Valeur &, int) const
Idem as Scalar::poisson_dirichlet , the boundary condition being on the radial derivative of the solu...
friend Scalar operator+(const Scalar &, const Scalar &)
Scalar + Scalar.
Base_val base
Bases on which the spectral expansion is performed.
Definition: valeur.h:315
int get_dzpuis() const
Returns dzpuis.
Definition: scalar.h:569
Scalar sol_elliptic(Param_elliptic &params) const
Resolution of a general elliptic equation, putting zero at infinity.
Definition: scalar_pde.C:223
virtual void spectral_display(const char *comment=0x0, double threshold=1.e-7, int precision=4, ostream &ostr=cout) const
Displays the spectral coefficients and the associated basis functions.
Definition: scalar.C:747
void set_outer_boundary(int l, double x)
Sets the value of the Scalar at the outer boundary of a given domain.
Definition: scalar_manip.C:327
friend Scalar sin(const Scalar &)
Sine.
Definition: scalar_math.C:75
void mult_r_ced()
Multiplication by r in the compactified external domain (CED), the dzpuis flag is not changed...
Scalar sol_elliptic_boundary(Param_elliptic &params, const Mtbl_cf &bound, double fact_dir, double fact_neu) const
Resolution of a general elliptic equation, putting zero at infinity and with inner boundary condition...
Definition: scalar_pde.C:245
Parameter storage.
Definition: param.h:125
void set_etat_one()
Sets the logical state to ETATUN (one).
Definition: scalar.C:340
friend Tbl diffrelmax(const Scalar &a, const Scalar &b)
Relative difference between two Scalar (max version).
Definition: scalar_math.C:733
Valeur va
The numerical value of the Scalar.
Definition: scalar.h:417
void del_t()
Logical destructor.
Definition: scalar.C:285
Tbl max(const Cmp &)
Maximum values of a Cmp in each domain.
Definition: cmp_math.C:438
virtual void std_spectral_base_odd()
Sets the spectral bases of the Valeur va to the standard odd ones for a scalar field.
Definition: scalar.C:797
const Scalar & deriv(int i) const
Returns of *this , where .
Definition: scalar_deriv.C:358
double totalmax(const Mtbl &)
Maximum value of the Mtbl elements in all domains.
Definition: mtbl_math.C:497
friend Tbl diffrel(const Scalar &a, const Scalar &b)
Relative difference between two Scalar (norme version).
Definition: scalar_math.C:698
friend Scalar operator|(const Scalar &, const Scalar &)
Scalar * Scalar with desaliasing only in r.
const Scalar & stdsdp() const
Returns of *this .
Definition: scalar_deriv.C:237
void mult_cost()
Multiplication by .
void div_rsint_dzpuis(int ced_mult_r)
Division by but with the output flag dzpuis set to ced_mult_r .
friend double totalmin(const Scalar &)
Minimum values of a Scalar in each domain.
Definition: scalar_math.C:585
void visu_section(const char section_type, double aa, double umin, double umax, double vmin, double vmax, const char *title=0x0, const char *filename=0x0, bool start_dx=true, int nu=200, int nv=200) const
3D visualization via a plane section.
Definition: scalar_visu.C:81
void mult_rsint()
Multiplication by everywhere; dzpuis is not changed.
Cmp pow(const Cmp &, int)
Power .
Definition: cmp_math.C:351
Tensor handling.
Definition: tensor.h:300
void sarra_filter_r_all_domains(double p, double alpha=1E-16)
Applies an exponential filter in radial direction in all domains for the case where p is a double (se...
void filtre(int n)
Sets the n lasts coefficients in r to 0 in the external domain.
Definition: scalar_manip.C:163
This class contains the parameters needed to call the general elliptic solver.
friend Scalar operator-(const Scalar &)
- Scalar
Definition: scalar_arithm.C:93
friend Scalar log(const Scalar &)
Neperian logarithm.
Definition: scalar_math.C:388
double totalmin(const Mtbl &)
Minimum value of the Mtbl elements in all domain.
Definition: mtbl_math.C:525
friend Scalar racine_cubique(const Scalar &)
Cube root.
Definition: scalar_math.C:296
Scalar derive_lie(const Vector &v) const
Computes the derivative of this along a vector field v.
Definition: scalar_deriv.C:418
void import_align_symy(int nzet, const Scalar &ci)
Assignment to another Scalar defined on a different mapping, when the two mappings have aligned Carte...
double val_point(double r, double theta, double phi) const
Computes the value of the field at an arbitrary point , by means of the spectral expansion.
Definition: scalar.C:896
void import_align(int nzet, const Scalar &ci)
Assignment to another Scalar defined on a different mapping, when the two mappings have aligned Carte...
void set_inner_boundary(int l, double x)
Sets the value of the Scalar at the inner boundary of a given domain.
Definition: scalar_manip.C:295
void filtre_tp(int nn, int nz1, int nz2)
Sets the n lasts coefficients in to 0 in the domain nz1 to nz2 when expressed in spherical harmonics...
Definition: scalar_manip.C:282
void fixe_decroissance(int puis)
Substracts all the components behaving like in the external domain, with n strictly lower than puis ...
Definition: scalar_manip.C:358
void set_spectral_base(const Base_val &)
Sets the spectral bases of the Valeur va
Definition: scalar.C:803
void copy_coefs_from_large_grid(const Scalar &scal_large_grid)
Copies the content of the argument Scalar to this defined on a smaller grid, with similar mappings...
Definition: scalar_manip.C:590
Multi-domain grid.
Definition: grilles.h:276
Bases of the spectral expansions.
Definition: base_val.h:325
void import_anti_asymy(int nzet, const Scalar &ci)
Assignment to another Scalar defined on a different mapping, when the two mappings have anti-aligned ...
Scalar * p_dsdz
Pointer on of *this , where (0x0 if not up to date)
Definition: scalar.h:456
double & set_grid_point(int l, int k, int j, int i)
Setting the value of the field at a given grid point.
Definition: scalar.h:696
friend Scalar operator/(const Scalar &, const Scalar &)
Scalar / Scalar.
void filtre_phi(int n, int zone)
Sets the n lasts coefficients in to 0 in the domain zone .
Definition: scalar_manip.C:261
bool dz_nonzero() const
Returns true if the last domain is compactified and *this is not zero in this domain.
Definition: scalar.C:820
void sarra_filter_r(int lzmin, int lzmax, double p, double alpha=-1E-16)
Applies an exponential filter to the spectral coefficients in the radial direction.
const Scalar & srdsdt() const
Returns of *this .
Definition: scalar_deriv.C:144
virtual ~Scalar()
Destructor.
Definition: scalar.C:273
void raccord_c1_zec(int puis, int nbre, int lmax)
Performs the matching of the external domain with respect to the last shell using function like wit...
Cmp log10(const Cmp &)
Basis 10 logarithm.
Definition: cmp_math.C:325
void exp_filter_r_all_domains(Scalar &ss, int p, double alpha=-16.)
Applies an exponential filter in radial direction in all domains.
Scalar sol_elliptic_pseudo_1d(Param_elliptic &) const
Solves a pseudo-1d elliptic equation with *this as a source.
Definition: scalar_pde.C:419
Cmp acos(const Cmp &)
Arccosine.
Definition: cmp_math.C:172
Cmp abs(const Cmp &)
Absolute value.
Definition: cmp_math.C:413
friend Scalar log10(const Scalar &)
Basis 10 logarithm.
Definition: scalar_math.C:421
int etat
The logical state ETATNONDEF (undefined), ETATZERO (null), ETATUN (one), or ETATQCQ (ordinary)...
Definition: scalar.h:408
Scalar * p_stdsdp
Pointer on of *this (0x0 if not up to date)
Definition: scalar.h:441
Coefficients storage for the multi-domain spectral method.
Definition: mtbl_cf.h:196
const Scalar & dsdr() const
Returns of *this .
Definition: scalar_deriv.C:112
void import_symy(const Scalar &ci)
Assignment to another Scalar defined on a different mapping.
friend Tbl norme(const Scalar &)
Sums of the absolute values of all the values of the Scalar in each domain.
Definition: scalar_math.C:670
Scalar * p_dsdr
Pointer on of *this (0x0 if not up to date)
Definition: scalar.h:423
void div_sint()
Division by .
void div_rsint()
Division by everywhere; dzpuis is not changed.
Valeur ** asymptot(int n, const int flag=0) const
Asymptotic expansion at r = infinity.
Tbl tbl_out_bound(int l_dom, bool leave_ylm=false)
Returns the Tbl containing the values of angular coefficients at the outer boundary.
Definition: scalar_manip.C:442
Cmp operator-(const Cmp &)
- Cmp
Definition: cmp_arithm.C:111
void mult_rsint_dzpuis(int ced_mult_r)
Multiplication by but with the output flag dzpuis set to ced_mult_r .
Cmp sin(const Cmp &)
Sine.
Definition: cmp_math.C:72
Mtbl Heaviside(const Mtbl &)
Heaviside function.
Definition: mtbl_math.C:320
Basic array class.
Definition: tbl.h:164
Valeur & set_spectral_va()
Returns va (read/write version)
Definition: scalar.h:616
bool check_dzpuis(int dzi) const
Returns false if the last domain is compactified and *this is not zero in this domain and dzpuis is n...
Definition: scalar.C:879
bool is_nan(bool verbose=false) const
Looks for NaNs (not a number) in the scalar field.
Definition: scalar.C:1014
void match_tau_star(Param &par_bc, Param *par_mat=0x0)
Method for matching accross domains and imposing boundary condition.
void match_collocation(Param &par_bc, Param *par_mat=0x0)
Method for matching accross domains and imposing boundary condition.
friend Scalar tan(const Scalar &)
Tangent.
Definition: scalar_math.C:139
friend Scalar abs(const Scalar &)
Absolute value.
Definition: scalar_math.C:526
Scalar poisson_dirichlet(const Valeur &limite, int num) const
Is identicall to Scalar::poisson() .
virtual void exponential_filter_ylm(int lzmin, int lzmax, int p, double alpha=-16.)
Applies an exponential filter to the spectral coefficients in the angular directions.
void match_tau_shell(Param &par_bc, Param *par_mat=0x0)
Method for matching accross domains and imposing boundary condition.
Scalar sol_elliptic_sin_zec(Param_elliptic &params, double *coefs, double *phases) const
General elliptic solver.
Definition: scalar_pde.C:352
void import_align_asymy(int nzet, const Scalar &ci)
Assignment to another Scalar defined on a different mapping, when the two mappings have aligned Carte...
void div_r_dzpuis(int ced_mult_r)
Division by r everywhere but with the output flag dzpuis set to ced_mult_r .
const Vector & derive_cov(const Metric &gam) const
Returns the gradient (1-form = covariant vector) of *this
Definition: scalar_deriv.C:389
void div_tant()
Division by .
void raccord_externe(int puis, int nbre, int lmax)
Matching of the external domain with the outermost shell.
void div_r_ced()
Division by r in the compactified external domain (CED), the dzpuis flag is not changed.
friend Scalar operator%(const Scalar &, const Scalar &)
Scalar * Scalar with desaliasing.
virtual void dec_dzpuis(int dec=1)
Decreases by dec units the value of dzpuis and changes accordingly the values of the Scalar in the co...
Scalar sol_elliptic_2d(Param_elliptic &) const
Solves the scalar 2-dimensional elliptic equation with *this as a source.
Definition: scalar_pde.C:398
void import_asymy(const Scalar &ci)
Assignment to another Scalar defined on a different mapping.
friend Scalar asin(const Scalar &)
Arcsine.
Definition: scalar_math.C:170
friend Scalar Heaviside(const Scalar &)
Heaviside function.
Definition: scalar_math.C:358
void mult_r_dzpuis(int ced_mult_r)
Multiplication by r everywhere but with the output flag dzpuis set to ced_mult_r ...
Tbl diffrelmax(const Cmp &a, const Cmp &b)
Relative difference between two Cmp (max version).
Definition: cmp_math.C:542
const Valeur & get_spectral_va() const
Returns va (read only version)
Definition: scalar.h:613
const Scalar & srstdsdp() const
Returns of *this .
Definition: scalar_deriv.C:176
Tbl & set(int l)
Read/write of the value in a given domain (configuration space).
Definition: valeur.h:373
Scalar poisson_tau() const
Solves the scalar Poisson equation with *this as a source using a real Tau method The source of the ...
Definition: scalar_pde.C:165
Scalar * p_dsdx
Pointer on of *this , where (0x0 if not up to date)
Definition: scalar.h:446