KADATH
spheric.hpp
1 /*
2  Copyright 2017 Philippe Grandclement
3 
4  This file is part of Kadath.
5 
6  Kadath is free software: you can redistribute it and/or modify
7  it under the terms of the GNU General Public License as published by
8  the Free Software Foundation, either version 3 of the License, or
9  (at your option) any later version.
10 
11  Kadath is distributed in the hope that it will be useful,
12  but WITHOUT ANY WARRANTY; without even the implied warranty of
13  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  GNU General Public License for more details.
15 
16  You should have received a copy of the GNU General Public License
17  along with Kadath. If not, see <http://www.gnu.org/licenses/>.
18 */
19 
20 #ifndef __SPHERIC_HPP_
21 #define __SPHERIC_HPP_
22 
23 #define STD_TYPE 0
24 #define LOG_TYPE 1
25 #define SURR_TYPE 2
26 
27 #include "point.hpp"
28 #include "space.hpp"
29 
30 namespace Kadath {
31 
32 class Domain_shell_log ;
33 class Domain_shell_surr ;
34 class Term_eq ;
35 class Metric ;
36 
66 class Domain_nucleus : public Domain {
67 
68  private:
69  double alpha ;
71 
72  public:
81  Domain_nucleus (int num, int ttype, double radius, const Point& cr, const Dim_array& nbr) ;
82  Domain_nucleus (const Domain_nucleus& so) ;
88  Domain_nucleus (int num, FILE* ff) ;
89 
90  virtual ~Domain_nucleus() ;
91  virtual void save (FILE*) const ;
92 
93  private:
94  virtual void do_absol () const ;
95  virtual void do_cart () const ;
96  virtual void do_cart_surr () const ;
97  virtual void do_radius () const ;
98  private:
110  virtual void set_cheb_base(Base_spectral& so) const ;
111  virtual void set_cheb_base_with_m(Base_spectral& so, int m) const ;
123  virtual void set_legendre_base(Base_spectral& so) const ;
135  virtual void set_anti_cheb_base(Base_spectral& so) const ;
147  virtual void set_anti_legendre_base(Base_spectral& so) const ;
148 
149  virtual void set_cheb_base_r_spher(Base_spectral&) const ;
150  virtual void set_cheb_base_t_spher(Base_spectral&) const ;
151  virtual void set_cheb_base_p_spher(Base_spectral&) const ;
152  virtual void set_cheb_base_r_mtz(Base_spectral&) const ;
153  virtual void set_cheb_base_t_mtz(Base_spectral&) const ;
154  virtual void set_cheb_base_p_mtz(Base_spectral&) const ;
155 
156  virtual void set_legendre_base_r_spher(Base_spectral&) const ;
157  virtual void set_legendre_base_t_spher(Base_spectral&) const ;
158  virtual void set_legendre_base_p_spher(Base_spectral&) const ;
159  virtual void set_legendre_base_r_mtz(Base_spectral&) const ;
160  virtual void set_legendre_base_t_mtz(Base_spectral&) const ;
161  virtual void set_legendre_base_p_mtz(Base_spectral&) const ;
162 
163  virtual void set_cheb_r_base(Base_spectral&) const ;
164  virtual void set_legendre_r_base(Base_spectral&) const ;
165 
166  virtual void do_coloc () ;
167  virtual int give_place_var (char*) const ;
168  public:
169  virtual double get_rmax() const {return alpha ; } ;
170  virtual Point get_center () const {return center ;} ;
171 
172 
173  virtual bool is_in(const Point&xx, double prec=1e-13) const ;
182  virtual const Point absol_to_num(const Point& xxx) const;
183 
184  virtual const Point absol_to_num_bound(const Point&, int) const;
185 
203  virtual void do_der_abs_from_der_var(const Val_domain *const *const der_var, Val_domain **const der_abs) const ;
204 
205  virtual Base_spectral mult (const Base_spectral&, const Base_spectral&) const ;
206 
207  public:
208 
209  virtual Val_domain mult_cos_phi (const Val_domain&) const ;
210  virtual Val_domain mult_sin_phi (const Val_domain&) const ;
211  virtual Val_domain mult_cos_theta (const Val_domain&) const ;
212  virtual Val_domain mult_sin_theta (const Val_domain&) const ;
213  virtual Val_domain div_sin_theta (const Val_domain&) const ;
214  virtual Val_domain div_cos_theta (const Val_domain&) const ;
215 
216  virtual Tensor change_basis_cart_to_spher (int dd, const Tensor&) const ;
217  virtual Tensor change_basis_spher_to_cart (int dd, const Tensor&) const ;
218 
219 
220  virtual Val_domain div_x (const Val_domain&) const ;
221  virtual Val_domain mult_r (const Val_domain&) const ;
222  virtual Val_domain div_r (const Val_domain&) const ;
223  virtual Val_domain der_r (const Val_domain&) const ;
224  virtual Val_domain ddp (const Val_domain&) const ;
225  virtual Val_domain srdr (const Val_domain&) const ;
226  virtual Val_domain div_1mx2 (const Val_domain&) const ;
227  virtual Val_domain laplacian2 (const Val_domain&, int) const ;
228 
229  virtual double val_boundary (int, const Val_domain&, const Index&) const ;
230  virtual void find_other_dom (int, int, int&, int&) const ;
231  virtual Val_domain der_normal (const Val_domain&, int) const ;
232 
233  virtual int nbr_unknowns (const Tensor&, int) const ;
242  int nbr_unknowns_val_domain (const Val_domain& so, int mlim, int llim) const ;
249  int nbr_unknowns_val_domain_vr (const Val_domain& so) const ;
256  int nbr_unknowns_val_domain_vt (const Val_domain& so) const ;
263  int nbr_unknowns_val_domain_vp (const Val_domain& so) const ;
264 
265  virtual Array<int> nbr_conditions (const Tensor&, int, int, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
275  int nbr_conditions_val_domain (const Val_domain& so, int mlim, int llim, int order) const ;
284  int nbr_conditions_val_domain_vr (const Val_domain& so, int order) const ;
293  int nbr_conditions_val_domain_vp (const Val_domain& so, int order) const ;
302  int nbr_conditions_val_domain_vt (const Val_domain& so, int order) const ;
303 
304  virtual Array<int> nbr_conditions_boundary (const Tensor&, int, int, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
313  int nbr_conditions_val_domain_boundary (const Val_domain& eq, int mlim) const ;
322  int nbr_conditions_val_domain_boundary_vr (const Val_domain& eq) const ;
331  int nbr_conditions_val_domain_boundary_vt (const Val_domain& eq) const ;
340  int nbr_conditions_val_domain_boundary_vp (const Val_domain& eq) const ;
341 
342  virtual void export_tau (const Tensor&, int, int, Array<double>&, int&, const Array<int>&, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
354  void export_tau_val_domain (const Val_domain& eq, int mlim, int llim, int order, Array<double>& res, int& pos_res, int ncond) const ;
365  void export_tau_val_domain_vr (const Val_domain& eq, int order, Array<double>& res, int& pos_res, int ncond) const ;
376  void export_tau_val_domain_vt (const Val_domain& eq, int order, Array<double>& res, int& pos_res, int ncond) const ;
387  void export_tau_val_domain_vp (const Val_domain& eq, int order, Array<double>& res, int& pos_res, int ncond) const ;
388 
389  virtual void export_tau_boundary (const Tensor&, int, int, Array<double>&, int&, const Array<int>&, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
400  void export_tau_val_domain_boundary (const Val_domain& eq, int mlim, int bound, Array<double>& res, int& pos_res, int ncond) const ;
411  void export_tau_val_domain_boundary_vr (const Val_domain& eq, int bound, Array<double>& res, int& pos_res, int ncond) const ;
422  void export_tau_val_domain_boundary_vt (const Val_domain& eq, int bound, Array<double>& res, int& pos_res, int ncond) const ;
433  void export_tau_val_domain_boundary_vp (const Val_domain& eq, int bound, Array<double>& res, int& pos_res, int ncond) const ;
434 
435  virtual void affecte_tau (Tensor&, int, const Array<double>&, int&) const ;
445  void affecte_tau_val_domain (Val_domain& so, int mlim, int llim, const Array<double>& cf, int& pos_cf) const ;
454  void affecte_tau_val_domain_vr (Val_domain& so, const Array<double>& cf, int& pos_cf) const ;
463  void affecte_tau_val_domain_vt (Val_domain& so, const Array<double>& cf, int& pos_cf) const ;
472  void affecte_tau_val_domain_vp (Val_domain& so, const Array<double>& cf, int& pos_cf) const ;
473 
474  virtual void affecte_tau_one_coef (Tensor&, int, int, int&) const ;
484  void affecte_tau_one_coef_val_domain (Val_domain& so, int mlim, int llim, int cc, int& pos_cf) const ;
493  void affecte_tau_one_coef_val_domain_vr (Val_domain& so, int cc, int& pos_cf) const ;
502  void affecte_tau_one_coef_val_domain_vt (Val_domain& so, int cc, int& pos_cf) const ;
511  void affecte_tau_one_coef_val_domain_vp (Val_domain& so, int cc, int& pos_cf) const ;
512 
513  virtual int nbr_points_boundary (int, const Base_spectral&) const ;
514  virtual void do_which_points_boundary (int, const Base_spectral&, Index**, int) const ;
515 
516  virtual double integ_volume (const Val_domain& so) const ;
517 
518  virtual Term_eq derive_flat_spher (int, char, const Term_eq&, const Metric*) const ;
519  virtual Term_eq derive_flat_cart (int, char, const Term_eq&, const Metric*) const ;
520  virtual double integ(const Val_domain& so, int bound) const ;
521  virtual Tensor import (int, int, int, const Array<int>&, Tensor**) const ;
522 
523 public:
524  virtual ostream& print (ostream& o) const ;
525 } ;
555 class Domain_shell : public Domain {
556 
557  private:
558  double alpha ;
559  double beta ;
561 
562  public:
572  Domain_shell (int num, int ttype, double r_int, double r_ext, const Point& cr, const Dim_array& nbr) ;
573  Domain_shell (const Domain_shell& so) ;
579  Domain_shell (int num, FILE* ff) ;
580 
581  virtual ~Domain_shell() ;
582  virtual void save (FILE*) const ;
583 
584  private:
585  virtual void do_absol () const ;
586  virtual void do_radius () const ;
587  virtual void do_cart () const ;
588  virtual void do_cart_surr () const ;
589  private :
601  virtual void set_cheb_base(Base_spectral& so) const ;
602  virtual void set_cheb_base_with_m(Base_spectral& so, int m) const ;
614  virtual void set_legendre_base(Base_spectral& so) const ;
626  virtual void set_anti_cheb_base(Base_spectral& so) const ;
638  virtual void set_anti_legendre_base(Base_spectral& so) const ;
639 
640  virtual void set_cheb_base_r_spher(Base_spectral&) const ;
641  virtual void set_cheb_base_t_spher(Base_spectral&) const ;
642  virtual void set_cheb_base_p_spher(Base_spectral&) const ;
643  virtual void set_cheb_base_r_mtz(Base_spectral&) const ;
644  virtual void set_cheb_base_t_mtz(Base_spectral&) const ;
645  virtual void set_cheb_base_p_mtz(Base_spectral&) const ;
646 
647  virtual void set_legendre_base_r_spher(Base_spectral&) const ;
648  virtual void set_legendre_base_t_spher(Base_spectral&) const ;
649  virtual void set_legendre_base_p_spher(Base_spectral&) const ;
650  virtual void set_legendre_base_r_mtz(Base_spectral&) const ;
651  virtual void set_legendre_base_t_mtz(Base_spectral&) const ;
652  virtual void set_legendre_base_p_mtz(Base_spectral&) const ;
653 
654  virtual void set_cheb_r_base(Base_spectral&) const ;
655  virtual void set_legendre_r_base(Base_spectral&) const ;
656 
657  virtual void do_coloc () ;
658 
659  virtual int give_place_var (char*) const ;
660  public:
661  virtual double get_rmin() const {return beta - alpha ;} ;
662  virtual double get_rmax() const {return beta+alpha ; } ;
663 
664  virtual Point get_center () const {return center ;} ;
665  virtual bool is_in(const Point& xx, double prec=1e-13) const ;
674  virtual const Point absol_to_num(const Point& xxx) const;
675 
676  virtual const Point absol_to_num_bound(const Point&, int) const;
694  virtual void do_der_abs_from_der_var(const Val_domain *const *const der_var, Val_domain **const der_abs) const ;
695 
696  virtual Base_spectral mult (const Base_spectral&, const Base_spectral&) const ;
697 
698  public:
699 
700  virtual Val_domain mult_cos_phi (const Val_domain&) const ;
701  virtual Val_domain mult_sin_phi (const Val_domain&) const ;
702  virtual Val_domain mult_cos_theta (const Val_domain&) const ;
703  virtual Val_domain mult_sin_theta (const Val_domain&) const ;
704  virtual Val_domain div_sin_theta (const Val_domain&) const ;
705  virtual Val_domain div_cos_theta (const Val_domain&) const ;
706 
707  virtual Tensor change_basis_cart_to_spher (int dd, const Tensor&) const ;
708  virtual Tensor change_basis_spher_to_cart (int dd, const Tensor&) const ;
709 
710  virtual Val_domain div_xp1 (const Val_domain&) const ;
711  virtual Val_domain div_1mrsL (const Val_domain&) const ;
712  virtual Val_domain mult_1mrsL (const Val_domain&) const ;
713 
714  virtual Val_domain mult_r (const Val_domain&) const ;
715  virtual Val_domain div_r (const Val_domain&) const ;
716  virtual Val_domain der_r (const Val_domain&) const ;
717  virtual Val_domain der_partial_var (const Val_domain&, int) const ;
718  virtual Val_domain ddp (const Val_domain&) const ;
719  virtual Val_domain div_xm1 (const Val_domain&) const ;
720  virtual Val_domain mult_xm1 (const Val_domain&) const ;
721  virtual Val_domain div_1mx2 (const Val_domain&) const ;
722  virtual Val_domain laplacian2 (const Val_domain&, int) const ;
723  virtual Val_domain dt (const Val_domain&) const ;
724 
725  virtual double multipoles_sym (int, int, int, const Val_domain&, const Array<double>&) const ;
726  virtual double multipoles_asym (int, int, int, const Val_domain&, const Array<double>&) const ;
727 
728  virtual Term_eq multipoles_sym (int, int, int, const Term_eq&, const Array<double>&) const ;
729  virtual Term_eq multipoles_asym (int, int, int, const Term_eq&, const Array<double>&) const ;
730 
731  virtual Term_eq radial_part_sym (const Space&, int, int, const Term_eq&, Term_eq (*f) (const Space&, int, int, const Term_eq&, const Param&), const Param&) const ;
732  virtual Term_eq radial_part_asym (const Space&, int, int, const Term_eq&, Term_eq (*f) (const Space&, int, int, const Term_eq&, const Param&), const Param&) const ;
733  virtual Term_eq harmonics_sym (const Term_eq&, const Term_eq&, int, Term_eq (*f) (const Space&, int, int, const Term_eq&, const Param&), const Param&, const Array<double>&) const ;
734  virtual Term_eq harmonics_asym (const Term_eq&, const Term_eq&, int, Term_eq (*f) (const Space&, int, int, const Term_eq&, const Param&), const Param&, const Array<double>&) const ;
735 
736 
737  virtual Term_eq der_radial_part_sym (const Space&, int, int, const Term_eq&, Term_eq (*f) (const Space&, int, int, const Term_eq&, const Param& param), const Param&) const ;
738  virtual Term_eq der_radial_part_asym (const Space&, int, int, const Term_eq&, Term_eq (*f) (const Space&, int, int, const Term_eq&, const Param&), const Param&) const ;
739  virtual Term_eq der_harmonics_sym (const Term_eq&, const Term_eq&, int, Term_eq (*f) (const Space&, int, int, const Term_eq&, const Param&), const Param&, const Array<double>&) const ;
740  virtual Term_eq der_harmonics_asym (const Term_eq&, const Term_eq&, int, Term_eq (*f) (const Space&, int, int, const Term_eq&, const Param&), const Param&, const Array<double>&) const ;
741  virtual Term_eq der_multipoles_sym (int, int, int, const Term_eq&, const Array<double>&) const ;
742  virtual Term_eq der_multipoles_asym (int, int, int, const Term_eq&, const Array<double>&) const ;
743 
750  Term_eq bc_waves (const Term_eq& gamma, const Term_eq& omega) const ;
759  Tensor bc_waves (int dom, const Tensor& gamma, const double omega, bool toinf = true) const ;
760 
767  Val_domain fitschwarz (const Val_domain&, int) const ;
768 
775  Term_eq fitschwarz (const Term_eq&, int) const ;
776 
777  virtual void find_other_dom (int, int, int&, int&) const ;
778  virtual Val_domain der_normal (const Val_domain&, int) const ;
779  virtual double val_boundary (int, const Val_domain&, const Index&) const ;
780 
781  virtual int nbr_unknowns (const Tensor&, int) const ;
789  int nbr_unknowns_val_domain_mquant (const Val_domain& so, int mquant) const ;
797  int nbr_unknowns_val_domain (const Val_domain& so, int mlim) const ;
798 
799  virtual Array<int> nbr_conditions (const Tensor&, int, int, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
808  int nbr_conditions_val_domain_mquant (const Val_domain& so, int mquant, int order) const ;
817  int nbr_conditions_val_domain (const Val_domain& so, int mlim, int order) const ;
818 
819  virtual Array<int> nbr_conditions_boundary (const Tensor&, int, int, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
828  int nbr_conditions_val_domain_boundary_mquant (const Val_domain& so, int mquant) const ;
837  int nbr_conditions_val_domain_boundary (const Val_domain& so, int mlim) const ;
838 
839  virtual void export_tau (const Tensor&, int, int, Array<double>&, int&, const Array<int>&, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
850  void export_tau_val_domain_mquant (const Val_domain& eq, int mquant, int order, Array<double>& res, int& pos_res, int ncond) const ;
861  void export_tau_val_domain (const Val_domain& eq, int mlim, int order, Array<double>& res, int& pos_res, int ncond) const ;
862 
863  virtual void export_tau_boundary (const Tensor&, int, int, Array<double>&, int&, const Array<int>&, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
874  void export_tau_val_domain_boundary_mquant (const Val_domain& eq, int mquant, int bound, Array<double>& res, int& pos_res, int ncond) const ;
885  void export_tau_val_domain_boundary (const Val_domain& eq, int mlim, int bound, Array<double>& res, int& pos_res, int ncond) const ;
886 
887  virtual void affecte_tau (Tensor&, int, const Array<double>&, int&) const ;
896  void affecte_tau_val_domain_mquant (Val_domain& so, int mquant, const Array<double>& cf, int& pos_cf) const ;
905  void affecte_tau_val_domain (Val_domain& so, int mlim, const Array<double>& cf, int& pos_cf) const ;
906 
907  virtual void affecte_tau_one_coef (Tensor&, int, int, int&) const ;
916  void affecte_tau_one_coef_val_domain_mquant (Val_domain& so, int mquant, int cc, int& pos_cf) const ;
925  void affecte_tau_one_coef_val_domain (Val_domain& so, int mlim, int cc, int& pos_cf) const ;
926 
927  virtual void export_tau_boundary_exception (const Tensor&, int, int, Array<double>&, int&, const Array<int>&, const Param&, int,
928  const Tensor&, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
943  void export_tau_val_domain_boundary_exception_mquant (const Val_domain& eq, int mquant, int bound, Array<double>& res, int& pos_res, int ncond, const Param& param,
944  int type_exception, const Val_domain& exception) const ;
958  void export_tau_val_domain_boundary_exception (const Val_domain& eq, int mlim, int bound, Array<double>& res, int& pos_res, int ncond, const Param& param,
959  int type_exception, const Val_domain& exception) const ;
960 
961  virtual int nbr_points_boundary (int, const Base_spectral&) const ;
962  virtual void do_which_points_boundary (int, const Base_spectral&, Index**, int) const ;
963  virtual Tensor import (int, int, int, const Array<int>&, Tensor**) const ;
964 
965  virtual double integ_volume (const Val_domain& so) const ;
966  virtual Term_eq derive_flat_spher (int, char, const Term_eq&, const Metric*) const ;
967  virtual Term_eq derive_flat_cart (int, char, const Term_eq&, const Metric*) const ;
968  virtual Term_eq derive_flat_mtz (int, char, const Term_eq&, const Metric*) const ;
969  virtual double integ(const Val_domain& so, int bound) const ;
970 
971 public:
972  virtual ostream& print (ostream& o) const ;
973 
974  friend class Domain_shell_log ;
975  friend class Domain_shell_surr ;
976 } ;
977 
1007 class Domain_compact : public Domain {
1008 
1009  private:
1010  double alpha ;
1012 
1013  public:
1022  Domain_compact (int num, int ttype, double r_int, const Point& cr, const Dim_array& nbr) ;
1023  Domain_compact (const Domain_compact& so) ;
1029  Domain_compact (int num, FILE* ff) ;
1030 
1031  virtual ~Domain_compact() ;
1032  virtual void save(FILE*) const ;
1033 
1034  private:
1035  virtual void do_absol () const ;
1036  virtual void do_radius () const ;
1037  virtual void do_cart () const ;
1038  virtual void do_cart_surr () const ;
1039 
1040  private:
1052  virtual void set_cheb_base(Base_spectral& so) const ;
1053  virtual void set_cheb_base_with_m(Base_spectral& so, int m) const ;
1065  virtual void set_legendre_base(Base_spectral& so) const ;
1077  virtual void set_anti_cheb_base(Base_spectral& so) const ;
1089  virtual void set_anti_legendre_base(Base_spectral& so) const ;
1090 
1091  virtual void set_cheb_base_r_spher(Base_spectral&) const ;
1092  virtual void set_cheb_base_t_spher(Base_spectral&) const ;
1093  virtual void set_cheb_base_p_spher(Base_spectral&) const ;
1094  virtual void set_legendre_base_r_spher(Base_spectral&) const ;
1095  virtual void set_legendre_base_t_spher(Base_spectral&) const ;
1096  virtual void set_legendre_base_p_spher(Base_spectral&) const ;
1097  virtual void set_legendre_base_r_mtz(Base_spectral&) const ;
1098  virtual void set_legendre_base_t_mtz(Base_spectral&) const ;
1099  virtual void set_legendre_base_p_mtz(Base_spectral&) const ;
1100  virtual void set_cheb_base_r_mtz(Base_spectral&) const ;
1101  virtual void set_cheb_base_t_mtz(Base_spectral&) const ;
1102  virtual void set_cheb_base_p_mtz(Base_spectral&) const ;
1103 
1104  virtual void set_cheb_r_base(Base_spectral&) const ;
1105  virtual void set_legendre_r_base(Base_spectral&) const ;
1106 
1107  virtual void do_coloc() ;
1108  virtual int give_place_var (char*) const ;
1109  public:
1110  virtual Point get_center () const {return center ;} ;
1111  double get_alpha() const
1112  {return alpha ;} ;
1113 
1120  virtual bool is_in(const Point& xx, double prec=1e-13) const ;
1129  virtual const Point absol_to_num(const Point& xxx) const;
1130 
1131  virtual const Point absol_to_num_bound(const Point&, int) const;
1149  virtual void do_der_abs_from_der_var(const Val_domain *const *const der_var, Val_domain **const der_abs) const ;
1155  virtual Base_spectral mult (const Base_spectral&, const Base_spectral&) const ;
1156 
1157  public:
1158 
1159  virtual Tensor change_basis_cart_to_spher (int dd, const Tensor&) const ;
1160  virtual Tensor change_basis_spher_to_cart (int dd, const Tensor&) const ;
1161 
1162 
1163 
1164  virtual Val_domain mult_cos_phi (const Val_domain&) const ;
1165  virtual Val_domain mult_sin_phi (const Val_domain&) const ;
1166  virtual Val_domain mult_cos_theta (const Val_domain&) const ;
1167  virtual Val_domain mult_sin_theta (const Val_domain&) const ;
1168  virtual Val_domain div_sin_theta (const Val_domain&) const ;
1169  virtual Val_domain div_cos_theta (const Val_domain&) const ;
1170  virtual Val_domain div_xm1 (const Val_domain&) const ;
1171  virtual Val_domain mult_xm1 (const Val_domain&) const ;
1172  virtual Val_domain mult_r (const Val_domain&) const ;
1173  virtual Val_domain div_r (const Val_domain&) const ;
1174  virtual Val_domain der_r (const Val_domain&) const ;
1175  virtual Val_domain der_r_rtwo (const Val_domain&) const ;
1176  virtual Val_domain der_partial_var (const Val_domain&, int) const ;
1177  virtual Val_domain laplacian2 (const Val_domain&, int) const ;
1178  virtual Val_domain dt (const Val_domain&) const ;
1179 
1180  virtual void set_val_inf (Val_domain&, double) const ;
1181 
1182  virtual void find_other_dom (int, int, int&, int&) const ;
1183  virtual Val_domain der_normal (const Val_domain&, int) const ;
1184  virtual double integ (const Val_domain&, int) const ;
1185 
1186  virtual double val_boundary (int, const Val_domain&, const Index&) const ;
1187  virtual int nbr_points_boundary (int, const Base_spectral&) const ;
1188  virtual void do_which_points_boundary (int, const Base_spectral&, Index**, int) const ;
1189 
1190  virtual int nbr_unknowns (const Tensor&, int) const ;
1198  int nbr_unknowns_val_domain_mquant (const Val_domain& so , int mquant) const ;
1206  int nbr_unknowns_val_domain (const Val_domain& so , int mlim) const ;
1207 
1208  virtual Array<int> nbr_conditions (const Tensor&, int, int, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
1217  int nbr_conditions_val_domain_mquant (const Val_domain& eq, int mquant, int order) const ;
1226  int nbr_conditions_val_domain (const Val_domain& eq, int mlim, int order) const ;
1227 
1228  virtual Array<int> nbr_conditions_boundary (const Tensor&, int, int, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
1237  int nbr_conditions_val_domain_boundary_mquant (const Val_domain& eq, int mquant) const ;
1246  int nbr_conditions_val_domain_boundary (const Val_domain& eq, int mlim) const ;
1247 
1248  virtual void export_tau (const Tensor&, int, int, Array<double>&, int&, const Array<int>&, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
1259  void export_tau_val_domain_mquant (const Val_domain& eq, int mquant, int order, Array<double>& res, int& pos_res, int ncond) const ;
1270  void export_tau_val_domain (const Val_domain& eq, int mlim, int order, Array<double>& res, int& pos_res, int ncond) const ;
1271 
1272  virtual void export_tau_boundary (const Tensor&, int, int, Array<double>&, int&, const Array<int>&, int n_cmp=-1, Array<int>** p_cmp=0x0) const ;
1283  void export_tau_val_domain_boundary_mquant (const Val_domain& eq, int mquant, int bound, Array<double>& res, int& pos_res, int ncond) const ;
1294  void export_tau_val_domain_boundary (const Val_domain& eq, int mlim, int bound, Array<double>& res, int& pos_res, int ncond) const ;
1295 
1296  virtual void affecte_tau (Tensor&, int, const Array<double>&, int&) const ;
1305  void affecte_tau_val_domain_mquant (Val_domain& so, int mquant, const Array<double>& cf, int& pos_cf) const ;
1314  void affecte_tau_val_domain (Val_domain& so, int mlim, const Array<double>& cf, int& pos_cf) const ;
1315 
1316  virtual void affecte_tau_one_coef (Tensor&, int, int, int&) const ;
1325  void affecte_tau_one_coef_val_domain_mquant (Val_domain& so, int mquant, int cc, int& pos_cf) const ;
1334  void affecte_tau_one_coef_val_domain (Val_domain& so, int mlim, int cc, int& pos_cf) const ;
1335 
1336  virtual Tensor import (int, int, int, const Array<int>&, Tensor**) const ;
1337 
1338  virtual double integ_volume (const Val_domain& so) const ;
1339  virtual Term_eq derive_flat_spher (int, char, const Term_eq&, const Metric*) const ;
1340  virtual Term_eq derive_flat_cart (int, char, const Term_eq&, const Metric*) const ;
1341 
1342 public:
1343  virtual ostream& print (ostream& o) const ;
1344 } ;
1345 
1350 class Space_spheric : public Space {
1351  public:
1360  Space_spheric (int ttype, const Point& cr, const Dim_array& nbr, const Array<double>& bounds, bool withzec=true) ;
1361 
1369  Space_spheric (int ttype, const Point& cr, Dim_array** nbr, const Array<double>& bounds) ;
1370 
1371 
1380  Space_spheric (int ttype, const Point& cr, const Dim_array& nbr, const Array<double>& bounds, const Array<int>& type_shells) ;
1381 
1382  Space_spheric (FILE*, bool withzec=true) ;
1383  virtual ~Space_spheric() ;
1384  virtual void save(FILE*) const ;
1385 
1394  void add_inner_bc (System_of_eqs& syst, const char* eq, int nused=-1, Array<int>** pused=0x0) ;
1395 
1403  void add_outer_bc (System_of_eqs& syst, const char* eq, int nused=-1, Array<int>** pused=0x0) ;
1404 
1412  void add_eq (System_of_eqs& syst, const char* eq, int nused=-1, Array<int>** pused=0x0) ;
1413 
1421  void add_matching (System_of_eqs& syst, const char* eq, int nused=-1, Array<int>** pused=0x0) ;
1422 
1432  void add_eq (System_of_eqs& syst, const char* eq, const char* rac, const char* rac_der, int nused=-1, Array<int>** pused=0x0) ;
1433 
1441  void add_eq_one_side (System_of_eqs& syst, const char* eq, int nused=-1, Array<int>** pused=0x0) ;
1442 
1450  void add_eq_full (System_of_eqs& syst, const char* eq, int nused=-1, Array<int>** pused=0x0) ;
1451 
1461  void add_eq_one_side (System_of_eqs& syst, const char* eq, const char* rac, int nused=-1, Array<int>** pused=0x0) ;
1467  void add_eq_int_volume (System_of_eqs& syst, const char* eq) ;
1468 
1474  void add_eq_ori (System_of_eqs& syst, const char* eq) ;
1483  void add_eq_mode_mid (System_of_eqs& syst, const char* f, int itarget, int jtarget, int ktarget) ;
1484 
1490  double int_inf (const Scalar&) const ;
1497  void add_eq_point (System_of_eqs& syst, const Point& pp, const char* eq) ;
1498 
1499  virtual Array<int> get_indices_matching_non_std(int, int) const ;
1500 } ;
1501 
1508 
1509 
1510  public:
1520  Domain_shell_log (int num, int ttype, double r_int, double r_ext, const Point& cr, const Dim_array& nbr) ;
1521  Domain_shell_log (const Domain_shell_log& so) ;
1527  Domain_shell_log (int num, FILE* ff) ;
1528 
1529  virtual ~Domain_shell_log() ;
1530  virtual void save (FILE*) const ;
1531 
1532  private:
1533  virtual void do_absol () const ;
1534  virtual void do_radius () const ;
1535  virtual void do_cart () const ;
1536 
1537  public:
1538  virtual double get_rmin() const {return exp(beta - alpha) ;} ;
1539  virtual double get_rmax() const {return exp(beta+alpha) ; } ;
1540 
1541 
1548  virtual bool is_in(const Point& xx, double prec=1e-13) const ;
1549 
1558  virtual const Point absol_to_num(const Point& xxx) const;
1559 
1560  virtual const Point absol_to_num_bound(const Point&, int) const;
1561 
1579  virtual void do_der_abs_from_der_var(const Val_domain *const *const der_var, Val_domain **const der_abs) const ;
1580 
1581  virtual Val_domain mult_r (const Val_domain&) const ;
1582  virtual Val_domain div_r (const Val_domain&) const ;
1583  virtual Val_domain der_r (const Val_domain&) const ;
1584 
1585  virtual Val_domain der_normal (const Val_domain&, int) const ;
1586 
1587 
1588 public:
1589  virtual ostream& print (ostream& o) const ;
1590 } ;
1591 
1598 
1599 
1600  public:
1610  Domain_shell_surr (int num, int ttype, double r_int, double r_ext, const Point& cr, const Dim_array& nbr) ;
1617  Domain_shell_surr (int num , FILE* ff) ;
1618 
1619  virtual ~Domain_shell_surr() ;
1620  virtual void save (FILE*) const ;
1621 
1622  private:
1623  virtual void do_absol () const ;
1624  virtual void do_radius () const ;
1625  virtual void do_cart () const ;
1626 
1627  public:
1628  virtual double get_rmin() const {return 1./(beta - alpha) ;} ;
1629  virtual double get_rmax() const {return 1./(beta+alpha) ; } ;
1630 
1631 
1638  virtual bool is_in(const Point& xx, double prec=1e-13) const ;
1639 
1648  virtual const Point absol_to_num(const Point& xxx) const;
1649 
1650  virtual const Point absol_to_num_bound(const Point&, int) const;
1651 
1669  virtual void do_der_abs_from_der_var(const Val_domain *const *const der_var, Val_domain **const der_abs) const ;
1670 
1671  virtual Val_domain mult_r (const Val_domain&) const ;
1672  virtual Val_domain div_r (const Val_domain&) const ;
1673  virtual Val_domain der_r (const Val_domain&) const ;
1674 
1675  virtual Val_domain der_normal (const Val_domain&, int) const ;
1676 
1677 
1678 public:
1679  virtual ostream& print (ostream& o) const ;
1680 
1681 } ;
1682 
1683 }
1684 #endif
Class for storing the basis of decompositions of a field.
Class for storing the dimensions of an array.
Definition: dim_array.hpp:34
Class for a spherical compactified domain and a symmetry with respect to the plane .
Definition: spheric.hpp:1007
virtual Tensor change_basis_spher_to_cart(int dd, const Tensor &) const
Changes the tensorial basis from spherical to Cartesian in a given domain.
int nbr_unknowns_val_domain(const Val_domain &so, int mlim) const
Computes the number of true unknowns of a Val_domain.
virtual double val_boundary(int, const Val_domain &, const Index &) const
Computes the value of a field at a boundary.
virtual void set_legendre_r_base(Base_spectral &) const
Gives the base using odd Legendre polynomials$ for the radius.
void affecte_tau_one_coef_val_domain(Val_domain &so, int mlim, int cc, int &pos_cf) const
Sets at most one coefficient of a Val_domain to 1.
void affecte_tau_val_domain_mquant(Val_domain &so, int mquant, const Array< double > &cf, int &pos_cf) const
Affects some coefficients to a Val_domain.
virtual ostream & print(ostream &o) const
Delegate function to virtualize the << operator.
virtual Term_eq derive_flat_cart(int, char, const Term_eq &, const Metric *) const
Computes the flat derivative of a Term_eq, in Cartesian coordinates.
virtual Val_domain div_xm1(const Val_domain &) const
Division by .
virtual void set_cheb_base(Base_spectral &so) const
Sets the base to the standard one for Chebyshev polynomials.
virtual void do_coloc()
Computes the colocation points.
virtual Val_domain der_normal(const Val_domain &, int) const
Normal derivative with respect to a given surface.
virtual Val_domain mult_sin_phi(const Val_domain &) const
Multiplication by .
int nbr_unknowns_val_domain_mquant(const Val_domain &so, int mquant) const
Computes the number of true unknowns of a Val_domain.
virtual void set_val_inf(Val_domain &, double) const
Sets the value at infinity of a Val_domain : not implemented for this type of Domain.
virtual Val_domain der_r_rtwo(const Val_domain &) const
Compute the radial derivative multiplied by of a scalar field.
virtual const Point absol_to_num_bound(const Point &, int) const
Computes the numerical coordinates from the physical ones for a point lying on a boundary.
virtual void save(FILE *) const
Saving function.
virtual void do_cart() const
Computes the Cartesian coordinates.
virtual void set_cheb_base_r_mtz(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the radial component of a vector in the MTZ setting.
virtual void set_legendre_base_t_mtz(Base_spectral &) const
Gives the base using Legendre polynomials, for the component of a vector in the MTZ context.
virtual Val_domain mult_cos_theta(const Val_domain &) const
Multiplication by .
double alpha
Relates the numerical to the physical radii.
Definition: spheric.hpp:1010
virtual Val_domain div_cos_theta(const Val_domain &) const
Division by .
virtual void find_other_dom(int, int, int &, int &) const
Gives the informations corresponding the a touching neighboring domain.
virtual Val_domain div_sin_theta(const Val_domain &) const
Division by .
virtual Val_domain mult_r(const Val_domain &) const
Multiplication by .
Point center
Absolute coordinates of the center.
Definition: spheric.hpp:1011
int nbr_conditions_val_domain(const Val_domain &eq, int mlim, int order) const
Computes number of discretized equations associated with a given tensorial equation in the bulk.
virtual Val_domain dt(const Val_domain &) const
Compute the derivative with respect to of a scalar field.
Domain_compact(int num, int ttype, double r_int, const Point &cr, const Dim_array &nbr)
Standard constructor :
virtual void do_cart_surr() const
Computes the Cartesian coordinates over the radius.
virtual const Point absol_to_num(const Point &xxx) const
Computes the numerical coordinates from the physical ones.
virtual Val_domain mult_sin_theta(const Val_domain &) const
Multiplication by .
virtual Val_domain der_partial_var(const Val_domain &, int) const
Partial derivative with respect to a coordinate.
virtual void set_cheb_base_p_spher(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the component of a vector.
double get_alpha() const
Returns the center ;.
Definition: spheric.hpp:1111
int nbr_conditions_val_domain_boundary(const Val_domain &eq, int mlim) const
Computes number of discretized equations associated with a given equation on a boundary.
virtual void set_legendre_base_r_spher(Base_spectral &) const
Gives the base using Legendre polynomials, for the radial component of a vector.
int nbr_conditions_val_domain_boundary_mquant(const Val_domain &eq, int mquant) const
Computes number of discretized equations associated with a given equation on a boundary.
virtual void set_legendre_base_r_mtz(Base_spectral &) const
Gives the base using Legendre polynomials, for the radial component of a vector in the MTZ context.
void export_tau_val_domain(const Val_domain &eq, int mlim, int order, Array< double > &res, int &pos_res, int ncond) const
Exports a residual equation in the bulk.
virtual Array< int > nbr_conditions(const Tensor &, int, int, int n_cmp=-1, Array< int > **p_cmp=0x0) const
Computes number of discretized equations associated with a given tensorial equation in the bulk.
virtual void do_absol() const
Computes the absolute coordinates.
virtual Val_domain laplacian2(const Val_domain &, int) const
Computes the ordinary flat 2dè- Laplacian for a scalar field with an harmonic index m.
virtual void set_legendre_base_p_mtz(Base_spectral &) const
Gives the base using Legendre polynomials, for the component of a vector in the MTZ context.
virtual Val_domain mult_cos_phi(const Val_domain &) const
Multiplication by .
void export_tau_val_domain_mquant(const Val_domain &eq, int mquant, int order, Array< double > &res, int &pos_res, int ncond) const
Exports a residual equation in the bulk.
virtual void affecte_tau(Tensor &, int, const Array< double > &, int &) const
Affects some coefficients to a Tensor.
void affecte_tau_val_domain(Val_domain &so, int mlim, const Array< double > &cf, int &pos_cf) const
Affects some coefficients to a Val_domain.
virtual Tensor change_basis_cart_to_spher(int dd, const Tensor &) const
Changes the tensorial basis from Cartsian to spherical in a given domain.
void export_tau_val_domain_boundary_mquant(const Val_domain &eq, int mquant, int bound, Array< double > &res, int &pos_res, int ncond) const
Exports all the residual equations corresponding to a tensorial one on a given boundary It makes use ...
virtual void export_tau(const Tensor &, int, int, Array< double > &, int &, const Array< int > &, int n_cmp=-1, Array< int > **p_cmp=0x0) const
Exports all the residual equations corresponding to a tensorial one in the bulk.
virtual void set_cheb_base_p_mtz(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the component of a vector in the MTZ setting.
int nbr_conditions_val_domain_mquant(const Val_domain &eq, int mquant, int order) const
Computes number of discretized equations associated with a given tensorial equation in the bulk.
virtual void set_cheb_base_t_mtz(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the component of a vector in the MTZ setting.
virtual void set_cheb_base_r_spher(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the radial component of a vector.
virtual void set_anti_legendre_base(Base_spectral &so) const
Sets the base to the standard one for Legendre polynomials for functions antisymetric in The bases a...
virtual Base_spectral mult(const Base_spectral &, const Base_spectral &) const
Method for the multiplication of two Base_spectral.
virtual void set_cheb_base_t_spher(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the component of a vector.
virtual int give_place_var(char *) const
Translates a name of a coordinate into its corresponding numerical name.
virtual void set_anti_cheb_base(Base_spectral &so) const
Sets the base to the standard one for Chebyshev polynomials for functions antisymetric in The bases ...
virtual int nbr_unknowns(const Tensor &, int) const
Computes the number of true unknowns of a Tensor, in a given domain.
virtual double integ_volume(const Val_domain &so) const
Volume integral.
virtual void set_legendre_base(Base_spectral &so) const
Sets the base to the standard one for Legendre polynomials.
virtual Val_domain mult_xm1(const Val_domain &) const
Multiplication by .
virtual bool is_in(const Point &xx, double prec=1e-13) const
Check whether a point lies inside Domain.
virtual Array< int > nbr_conditions_boundary(const Tensor &, int, int, int n_cmp=-1, Array< int > **p_cmp=0x0) const
Computes number of discretized equations associated with a given tensorial equation on a boundary.
virtual void export_tau_boundary(const Tensor &, int, int, Array< double > &, int &, const Array< int > &, int n_cmp=-1, Array< int > **p_cmp=0x0) const
Exports all the residual equations corresponding to a tensorial one on a given boundary It makes use ...
virtual Term_eq derive_flat_spher(int, char, const Term_eq &, const Metric *) const
Computes the flat derivative of a Term_eq, in spherical orthonormal coordinates.
virtual Point get_center() const
Returns the center.
Definition: spheric.hpp:1110
virtual void do_der_abs_from_der_var(const Val_domain *const *const der_var, Val_domain **const der_abs) const
Computes the derivative with respect to the absolute Cartesian coordinates from the derivative with r...
virtual void set_legendre_base_t_spher(Base_spectral &) const
Gives the base using Legendre polynomials, for the component of a vector.
void affecte_tau_one_coef_val_domain_mquant(Val_domain &so, int mquant, int cc, int &pos_cf) const
Sets at most one coefficient of a Val_domain to 1.
virtual void do_which_points_boundary(int, const Base_spectral &, Index **, int) const
Lists all the indices corresponding to true collocation points on a boundary.
virtual void affecte_tau_one_coef(Tensor &, int, int, int &) const
Sets at most one coefficient of a Tensor to 1.
virtual void set_cheb_r_base(Base_spectral &) const
Gives the base using odd Chebyshev polynomials$ for the radius.
virtual double integ(const Val_domain &, int) const
Surface integral on a given boundary.
virtual void do_radius() const
Computes the generalized radius.
virtual Val_domain div_r(const Val_domain &) const
Division by .
void export_tau_val_domain_boundary(const Val_domain &eq, int mlim, int bound, Array< double > &res, int &pos_res, int ncond) const
Exports all the residual equations corresponding to a tensorial one on a given boundary It makes use ...
virtual Val_domain der_r(const Val_domain &) const
Compute the radial derivative of a scalar field.
virtual int nbr_points_boundary(int, const Base_spectral &) const
Computes the number of relevant collocation points on a boundary.
virtual void set_cheb_base_with_m(Base_spectral &so, int m) const
Gives the standard base using Chebyshev polynomials.
virtual void set_legendre_base_p_spher(Base_spectral &) const
Gives the base using Legendre polynomials, for the component of a vector.
Class for a spherical domain containing the origin and a symmetry with respect to the plane .
Definition: spheric.hpp:66
void export_tau_val_domain_boundary_vr(const Val_domain &eq, int bound, Array< double > &res, int &pos_res, int ncond) const
Exports all the residual equations corresponding to a tensorial one on a given boundary It makes use ...
virtual void do_radius() const
Computes the generalized radius.
void export_tau_val_domain(const Val_domain &eq, int mlim, int llim, int order, Array< double > &res, int &pos_res, int ncond) const
Exports a residual equation in the bulk.
void export_tau_val_domain_boundary(const Val_domain &eq, int mlim, int bound, Array< double > &res, int &pos_res, int ncond) const
Exports all the residual equations corresponding to a tensorial one on a given boundary It makes use ...
Point center
Absolute coordinates of the center.
Definition: spheric.hpp:70
virtual Val_domain der_normal(const Val_domain &, int) const
Normal derivative with respect to a given surface.
virtual bool is_in(const Point &xx, double prec=1e-13) const
Check whether a point lies inside Domain.
virtual void affecte_tau(Tensor &, int, const Array< double > &, int &) const
Affects some coefficients to a Tensor.
virtual void do_cart() const
Computes the Cartesian coordinates.
virtual void set_cheb_r_base(Base_spectral &) const
Gives the base using odd Chebyshev polynomials$ for the radius.
virtual Val_domain div_x(const Val_domain &) const
Division by .
virtual Val_domain der_r(const Val_domain &) const
Compute the radial derivative of a scalar field.
virtual ostream & print(ostream &o) const
Delegate function to virtualize the << operator.
virtual void set_legendre_base_t_spher(Base_spectral &) const
Gives the base using Legendre polynomials, for the component of a vector.
void export_tau_val_domain_vr(const Val_domain &eq, int order, Array< double > &res, int &pos_res, int ncond) const
Exports a residual equation in the bulk.
virtual int nbr_unknowns(const Tensor &, int) const
Computes the number of true unknowns of a Tensor, in a given domain.
int nbr_conditions_val_domain_boundary_vt(const Val_domain &eq) const
Computes number of discretized equations associated with a given equation on a boundary.
virtual void set_legendre_base_p_spher(Base_spectral &) const
Gives the base using Legendre polynomials, for the component of a vector.
int nbr_unknowns_val_domain_vp(const Val_domain &so) const
Computes the number of true unknowns of a Val_domain.
void export_tau_val_domain_vt(const Val_domain &eq, int order, Array< double > &res, int &pos_res, int ncond) const
Exports a residual equation in the bulk.
void affecte_tau_val_domain_vp(Val_domain &so, const Array< double > &cf, int &pos_cf) const
Affects some coefficients to a Val_domain.
virtual Val_domain div_r(const Val_domain &) const
Division by .
int nbr_unknowns_val_domain_vt(const Val_domain &so) const
Computes the number of true unknowns of a Val_domain.
void affecte_tau_val_domain_vt(Val_domain &so, const Array< double > &cf, int &pos_cf) const
Affects some coefficients to a Val_domain.
virtual double integ(const Val_domain &so, int bound) const
Surface integral on a given boundary.
virtual Tensor change_basis_cart_to_spher(int dd, const Tensor &) const
Changes the tensorial basis from Cartsian to spherical in a given domain.
virtual Tensor change_basis_spher_to_cart(int dd, const Tensor &) const
Changes the tensorial basis from spherical to Cartesian in a given domain.
virtual void set_legendre_base(Base_spectral &so) const
Sets the base to the standard one for Legendre polynomials.
virtual Val_domain srdr(const Val_domain &) const
Compute the of a scalar field .
virtual const Point absol_to_num_bound(const Point &, int) const
Computes the numerical coordinates from the physical ones for a point lying on a boundary.
virtual void do_cart_surr() const
Computes the Cartesian coordinates over the radius.
virtual void set_cheb_base_t_spher(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the component of a vector.
virtual void set_anti_cheb_base(Base_spectral &so) const
Sets the base to the standard one for Chebyshev polynomials for functions antisymetric in The bases ...
void affecte_tau_one_coef_val_domain(Val_domain &so, int mlim, int llim, int cc, int &pos_cf) const
Sets at most one coefficient of a Val_domain to 1.
int nbr_conditions_val_domain_boundary(const Val_domain &eq, int mlim) const
Computes number of discretized equations associated with a given equation on a boundary.
virtual void set_cheb_base(Base_spectral &so) const
Sets the base to the standard one for Chebyshev polynomials.
virtual Base_spectral mult(const Base_spectral &, const Base_spectral &) const
Method for the multiplication of two Base_spectral.
void affecte_tau_one_coef_val_domain_vp(Val_domain &so, int cc, int &pos_cf) const
Sets at most one coefficient of a Val_domain to 1.
virtual Array< int > nbr_conditions_boundary(const Tensor &, int, int, int n_cmp=-1, Array< int > **p_cmp=0x0) const
Computes number of discretized equations associated with a given tensorial equation on a boundary.
int nbr_conditions_val_domain_vt(const Val_domain &so, int order) const
Computes number of discretized equations associated with a given tensorial equation in the bulk.
void affecte_tau_val_domain_vr(Val_domain &so, const Array< double > &cf, int &pos_cf) const
Affects some coefficients to a Val_domain.
virtual void set_legendre_base_t_mtz(Base_spectral &) const
Gives the base using Legendre polynomials, for the component of a vector in the MTZ context.
int nbr_conditions_val_domain_boundary_vr(const Val_domain &eq) const
Computes number of discretized equations associated with a given equation on a boundary.
virtual void set_cheb_base_with_m(Base_spectral &so, int m) const
Gives the standard base using Chebyshev polynomials.
virtual void set_cheb_base_r_spher(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the radial component of a vector.
virtual Val_domain div_sin_theta(const Val_domain &) const
Division by .
virtual Val_domain ddp(const Val_domain &) const
Compute the second derivative with respect to of a scalar field.
virtual double val_boundary(int, const Val_domain &, const Index &) const
Computes the value of a field at a boundary.
void export_tau_val_domain_vp(const Val_domain &eq, int order, Array< double > &res, int &pos_res, int ncond) const
Exports a residual equation in the bulk.
virtual void set_legendre_base_r_spher(Base_spectral &) const
Gives the base using Legendre polynomials, for the radial component of a vector.
int nbr_conditions_val_domain(const Val_domain &so, int mlim, int llim, int order) const
Computes number of discretized equations associated with a given tensorial equation in the bulk.
virtual void affecte_tau_one_coef(Tensor &, int, int, int &) const
Sets at most one coefficient of a Tensor to 1.
virtual Val_domain mult_sin_theta(const Val_domain &) const
Multiplication by .
virtual void find_other_dom(int, int, int &, int &) const
Gives the informations corresponding the a touching neighboring domain.
virtual int give_place_var(char *) const
Translates a name of a coordinate into its corresponding numerical name.
Domain_nucleus(int num, int ttype, double radius, const Point &cr, const Dim_array &nbr)
Standard constructor :
virtual void save(FILE *) const
Saving function.
virtual void set_cheb_base_p_spher(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the component of a vector.
void export_tau_val_domain_boundary_vp(const Val_domain &eq, int bound, Array< double > &res, int &pos_res, int ncond) const
Exports all the residual equations corresponding to a tensorial one on a given boundary It makes use ...
virtual Point get_center() const
Returns the center.
Definition: spheric.hpp:170
virtual Val_domain mult_cos_theta(const Val_domain &) const
Multiplication by .
int nbr_conditions_val_domain_boundary_vp(const Val_domain &eq) const
Computes number of discretized equations associated with a given equation on a boundary.
virtual double integ_volume(const Val_domain &so) const
Volume integral.
void export_tau_val_domain_boundary_vt(const Val_domain &eq, int bound, Array< double > &res, int &pos_res, int ncond) const
Exports all the residual equations corresponding to a tensorial one on a given boundary It makes use ...
virtual void set_legendre_base_r_mtz(Base_spectral &) const
Gives the base using Legendre polynomials, for the radial component of a vector in the MTZ context.
int nbr_conditions_val_domain_vr(const Val_domain &so, int order) const
Computes number of discretized equations associated with a given tensorial equation in the bulk.
virtual Term_eq derive_flat_cart(int, char, const Term_eq &, const Metric *) const
Computes the flat derivative of a Term_eq, in Cartesian coordinates.
virtual Val_domain div_cos_theta(const Val_domain &) const
Division by .
virtual const Point absol_to_num(const Point &xxx) const
Computes the numerical coordinates from the physical ones.
virtual Val_domain mult_sin_phi(const Val_domain &) const
Multiplication by .
virtual void set_legendre_r_base(Base_spectral &) const
Gives the base using odd Legendre polynomials$ for the radius.
virtual void do_absol() const
Computes the absolute coordinates.
virtual Val_domain mult_r(const Val_domain &) const
Multiplication by .
int nbr_unknowns_val_domain_vr(const Val_domain &so) const
Computes the number of true unknowns of a Val_domain.
double alpha
Relates the numerical to the physical radii.
Definition: spheric.hpp:69
virtual void do_coloc()
Computes the colocation points.
void affecte_tau_one_coef_val_domain_vr(Val_domain &so, int cc, int &pos_cf) const
Sets at most one coefficient of a Val_domain to 1.
virtual int nbr_points_boundary(int, const Base_spectral &) const
Computes the number of relevant collocation points on a boundary.
int nbr_unknowns_val_domain(const Val_domain &so, int mlim, int llim) const
Computes the number of true unknowns of a Val_domain.
virtual double get_rmax() const
Returns the maximum radius.
Definition: spheric.hpp:169
virtual Array< int > nbr_conditions(const Tensor &, int, int, int n_cmp=-1, Array< int > **p_cmp=0x0) const
Computes number of discretized equations associated with a given tensorial equation in the bulk.
virtual void export_tau(const Tensor &, int, int, Array< double > &, int &, const Array< int > &, int n_cmp=-1, Array< int > **p_cmp=0x0) const
Exports all the residual equations corresponding to a tensorial one in the bulk.
virtual void set_cheb_base_t_mtz(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the component of a vector in the MTZ setting.
virtual Val_domain laplacian2(const Val_domain &, int) const
Computes the ordinary flat 2dè- Laplacian for a scalar field with an harmonic index m.
void affecte_tau_one_coef_val_domain_vt(Val_domain &so, int cc, int &pos_cf) const
Sets at most one coefficient of a Val_domain to 1.
virtual Val_domain mult_cos_phi(const Val_domain &) const
Multiplication by .
virtual void export_tau_boundary(const Tensor &, int, int, Array< double > &, int &, const Array< int > &, int n_cmp=-1, Array< int > **p_cmp=0x0) const
Exports all the residual equations corresponding to a tensorial one on a given boundary It makes use ...
virtual void set_legendre_base_p_mtz(Base_spectral &) const
Gives the base using Legendre polynomials, for the component of a vector in the MTZ context.
void affecte_tau_val_domain(Val_domain &so, int mlim, int llim, const Array< double > &cf, int &pos_cf) const
Affects some coefficients to a Val_domain.
int nbr_conditions_val_domain_vp(const Val_domain &so, int order) const
Computes number of discretized equations associated with a given tensorial equation in the bulk.
virtual void do_der_abs_from_der_var(const Val_domain *const *const der_var, Val_domain **const der_abs) const
Computes the derivative with respect to the absolute Cartesian coordinates from the derivative with r...
virtual Val_domain div_1mx2(const Val_domain &) const
Division by .
virtual void do_which_points_boundary(int, const Base_spectral &, Index **, int) const
Lists all the indices corresponding to true collocation points on a boundary.
virtual void set_cheb_base_p_mtz(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the component of a vector in the MTZ setting.
virtual Term_eq derive_flat_spher(int, char, const Term_eq &, const Metric *) const
Computes the flat derivative of a Term_eq, in spherical orthonormal coordinates.
virtual void set_cheb_base_r_mtz(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the radial component of a vector in the MTZ setting.
virtual void set_anti_legendre_base(Base_spectral &so) const
Sets the base to the standard one for Legendre polynomials for functions antisymetric in The bases a...
Class for a spherical shell and a symmetry with respect to the plane .
Definition: spheric.hpp:1507
virtual Val_domain div_r(const Val_domain &) const
Division by .
virtual const Point absol_to_num_bound(const Point &, int) const
Computes the numerical coordinates from the physical ones for a point lying on a boundary.
virtual void do_absol() const
Computes the absolute coordinates.
virtual void do_der_abs_from_der_var(const Val_domain *const *const der_var, Val_domain **const der_abs) const
Computes the derivative with respect to the absolute Cartesian coordinates from the derivative with r...
virtual Val_domain der_normal(const Val_domain &, int) const
Normal derivative with respect to a given surface.
virtual void do_radius() const
Computes the generalized radius.
virtual Val_domain mult_r(const Val_domain &) const
Multiplication by .
virtual Val_domain der_r(const Val_domain &) const
Compute the radial derivative of a scalar field.
virtual const Point absol_to_num(const Point &xxx) const
Computes the numerical coordinates from the physical ones.
virtual double get_rmin() const
Returns the minimum radius.
Definition: spheric.hpp:1538
virtual double get_rmax() const
Returns the maximum radius.
Definition: spheric.hpp:1539
virtual bool is_in(const Point &xx, double prec=1e-13) const
Check whether a point lies inside Domain.
virtual ostream & print(ostream &o) const
Delegate function to virtualize the << operator.
virtual void save(FILE *) const
Saving function.
virtual void do_cart() const
Computes the Cartesian coordinates.
Class for a spherical shell and a symmetry with respect to the plane .
Definition: spheric.hpp:1597
virtual void do_cart() const
Computes the Cartesian coordinates.
virtual const Point absol_to_num(const Point &xxx) const
Computes the numerical coordinates from the physical ones.
virtual const Point absol_to_num_bound(const Point &, int) const
Computes the numerical coordinates from the physical ones for a point lying on a boundary.
virtual Val_domain der_r(const Val_domain &) const
Compute the radial derivative of a scalar field.
virtual void do_absol() const
Computes the absolute coordinates.
virtual ostream & print(ostream &o) const
Delegate function to virtualize the << operator.
virtual Val_domain mult_r(const Val_domain &) const
Multiplication by .
virtual double get_rmin() const
Returns the minimum radius.
Definition: spheric.hpp:1628
virtual Val_domain div_r(const Val_domain &) const
Division by .
virtual void save(FILE *) const
Saving function.
virtual double get_rmax() const
Returns the maximum radius.
Definition: spheric.hpp:1629
virtual void do_der_abs_from_der_var(const Val_domain *const *const der_var, Val_domain **const der_abs) const
Computes the derivative with respect to the absolute Cartesian coordinates from the derivative with r...
virtual Val_domain der_normal(const Val_domain &, int) const
Normal derivative with respect to a given surface.
virtual bool is_in(const Point &xx, double prec=1e-13) const
Check whether a point lies inside Domain.
virtual void do_radius() const
Computes the generalized radius.
Class for a spherical shell and a symmetry with respect to the plane .
Definition: spheric.hpp:555
virtual Term_eq der_harmonics_sym(const Term_eq &, const Term_eq &, int, Term_eq(*f)(const Space &, int, int, const Term_eq &, const Param &), const Param &, const Array< double > &) const
Fit, spherical harmonic by spherical harmonic, for the radial derivative of a symmetric function.
virtual Val_domain div_xp1(const Val_domain &) const
Division by .
virtual Val_domain div_r(const Val_domain &) const
Division by .
void export_tau_val_domain_boundary_mquant(const Val_domain &eq, int mquant, int bound, Array< double > &res, int &pos_res, int ncond) const
Exports all the residual equations corresponding to a tensorial one on a given boundary It makes use ...
virtual Term_eq derive_flat_cart(int, char, const Term_eq &, const Metric *) const
Computes the flat derivative of a Term_eq, in Cartesian coordinates.
virtual void affecte_tau(Tensor &, int, const Array< double > &, int &) const
Affects some coefficients to a Tensor.
virtual Val_domain mult_cos_theta(const Val_domain &) const
Multiplication by .
int nbr_unknowns_val_domain_mquant(const Val_domain &so, int mquant) const
Computes the number of true unknowns of a Val_domain.
virtual void set_legendre_base_t_spher(Base_spectral &) const
Gives the base using Legendre polynomials, for the component of a vector.
void affecte_tau_one_coef_val_domain(Val_domain &so, int mlim, int cc, int &pos_cf) const
Sets at most one coefficient of a Val_domain to 1.
virtual Array< int > nbr_conditions_boundary(const Tensor &, int, int, int n_cmp=-1, Array< int > **p_cmp=0x0) const
Computes number of discretized equations associated with a given tensorial equation on a boundary.
virtual void set_cheb_base_with_m(Base_spectral &so, int m) const
Gives the standard base using Chebyshev polynomials.
virtual bool is_in(const Point &xx, double prec=1e-13) const
Check whether a point lies inside Domain.
virtual Val_domain div_1mrsL(const Val_domain &) const
Division by .
virtual double integ_volume(const Val_domain &so) const
Volume integral.
virtual void affecte_tau_one_coef(Tensor &, int, int, int &) const
Sets at most one coefficient of a Tensor to 1.
int nbr_conditions_val_domain_boundary_mquant(const Val_domain &so, int mquant) const
Computes number of discretized equations associated with a given equation on a boundary.
virtual double integ(const Val_domain &so, int bound) const
Surface integral on a given boundary.
double beta
Relates the numerical to the physical radii.
Definition: spheric.hpp:559
void export_tau_val_domain_boundary_exception(const Val_domain &eq, int mlim, int bound, Array< double > &res, int &pos_res, int ncond, const Param &param, int type_exception, const Val_domain &exception) const
Exports all the residual equations corresponding to one tensorial one on a given boundary,...
virtual void set_cheb_base_r_mtz(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the radial component of a vector in the MTZ setting.
virtual void set_cheb_base_p_mtz(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the component of a vector in the MTZ setting.
virtual Val_domain div_sin_theta(const Val_domain &) const
Division by .
virtual Base_spectral mult(const Base_spectral &, const Base_spectral &) const
Method for the multiplication of two Base_spectral.
void export_tau_val_domain(const Val_domain &eq, int mlim, int order, Array< double > &res, int &pos_res, int ncond) const
Exports a residual equation in the bulk.
virtual void export_tau_boundary_exception(const Tensor &, int, int, Array< double > &, int &, const Array< int > &, const Param &, int, const Tensor &, int n_cmp=-1, Array< int > **p_cmp=0x0) const
Exports all the residual equations corresponding to one tensorial one on a given boundary,...
virtual Val_domain ddp(const Val_domain &) const
Compute the second derivative with respect to of a scalar field.
virtual void set_cheb_base_t_mtz(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the component of a vector in the MTZ setting.
virtual Val_domain der_partial_var(const Val_domain &, int) const
Partial derivative with respect to a coordinate.
virtual Val_domain der_r(const Val_domain &) const
Compute the radial derivative of a scalar field.
double alpha
Relates the numerical to the physical radii.
Definition: spheric.hpp:558
virtual double multipoles_asym(int, int, int, const Val_domain &, const Array< double > &) const
Extraction of a given multipole, at some boundary, for a anti-symmetric scalar function.
virtual void set_cheb_r_base(Base_spectral &) const
Gives the base using odd Chebyshev polynomials$ for the radius.
virtual Val_domain div_cos_theta(const Val_domain &) const
Division by .
virtual void do_coloc()
Computes the colocation points.
virtual void save(FILE *) const
Saving function.
virtual int give_place_var(char *) const
Translates a name of a coordinate into its corresponding numerical name.
virtual Term_eq der_multipoles_sym(int, int, int, const Term_eq &, const Array< double > &) const
Extraction of a given multipole, at some boundary, for the radial derivative a symmetric scalar funct...
virtual Term_eq harmonics_sym(const Term_eq &, const Term_eq &, int, Term_eq(*f)(const Space &, int, int, const Term_eq &, const Param &), const Param &, const Array< double > &) const
Fit, spherical harmonic by spherical harmonic, for a symmetric function.
Definition: multipoles.cpp:539
Domain_shell(int num, int ttype, double r_int, double r_ext, const Point &cr, const Dim_array &nbr)
Standard constructor :
virtual void set_legendre_base_r_spher(Base_spectral &) const
Gives the base using Legendre polynomials, for the radial component of a vector.
virtual Val_domain div_xm1(const Val_domain &) const
Division by .
virtual void do_which_points_boundary(int, const Base_spectral &, Index **, int) const
Lists all the indices corresponding to true collocation points on a boundary.
void export_tau_val_domain_mquant(const Val_domain &eq, int mquant, int order, Array< double > &res, int &pos_res, int ncond) const
Exports a residual equation in the bulk.
virtual void set_cheb_base_p_spher(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the component of a vector.
virtual void set_legendre_base_t_mtz(Base_spectral &) const
Gives the base using Legendre polynomials, for the component of a vector in the MTZ context.
Val_domain fitschwarz(const Val_domain &, int) const
Fit some field with a decay (Val_domain version).
virtual Val_domain mult_cos_phi(const Val_domain &) const
Multiplication by .
virtual Term_eq der_harmonics_asym(const Term_eq &, const Term_eq &, int, Term_eq(*f)(const Space &, int, int, const Term_eq &, const Param &), const Param &, const Array< double > &) const
Fit, spherical harmonic by spherical harmonic, for the radial derivative of an anti-symmetric functio...
virtual Term_eq der_multipoles_asym(int, int, int, const Term_eq &, const Array< double > &) const
Extraction of a given multipole, at some boundary, for the radial derivative of an anti-symmetric sca...
virtual void find_other_dom(int, int, int &, int &) const
Gives the informations corresponding the a touching neighboring domain.
virtual Val_domain der_normal(const Val_domain &, int) const
Normal derivative with respect to a given surface.
virtual void set_legendre_base(Base_spectral &so) const
Sets the base to the standard one for Legendre polynomials.
virtual Point get_center() const
Returns the center.
Definition: spheric.hpp:664
virtual void set_cheb_base(Base_spectral &so) const
Sets the base to the standard one for Chebyshev polynomials.
virtual Term_eq der_radial_part_asym(const Space &, int, int, const Term_eq &, Term_eq(*f)(const Space &, int, int, const Term_eq &, const Param &), const Param &) const
Gives some radial fit for a given multipole, intended for the radial derivative of an anti-symmetric ...
int nbr_conditions_val_domain_mquant(const Val_domain &so, int mquant, int order) const
Computes number of discretized equations associated with a given tensorial equation in the bulk.
virtual void set_anti_legendre_base(Base_spectral &so) const
Sets the base to the standard one for Legendre polynomials for functions antisymetric in The bases a...
virtual double multipoles_sym(int, int, int, const Val_domain &, const Array< double > &) const
Extraction of a given multipole, at some boundary, for a symmetric scalar function.
Term_eq bc_waves(const Term_eq &gamma, const Term_eq &omega) const
Gives an matching of the spatial metric, based on homogeneous solutions of outgoing waves.
Definition: bc_waves.cpp:33
virtual void set_cheb_base_r_spher(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the radial component of a vector.
int nbr_conditions_val_domain_boundary(const Val_domain &so, int mlim) const
Computes number of discretized equations associated with a given equation on a boundary.
virtual double get_rmin() const
Returns the minimum radius.
Definition: spheric.hpp:661
virtual Term_eq der_radial_part_sym(const Space &, int, int, const Term_eq &, Term_eq(*f)(const Space &, int, int, const Term_eq &, const Param &param), const Param &) const
Gives some radial fit for a given multipole, intended for the radial derivative of a symmetric scalar...
void affecte_tau_val_domain(Val_domain &so, int mlim, const Array< double > &cf, int &pos_cf) const
Affects some coefficients to a Val_domain.
virtual Val_domain dt(const Val_domain &) const
Compute the derivative with respect to of a scalar field.
virtual Term_eq radial_part_asym(const Space &, int, int, const Term_eq &, Term_eq(*f)(const Space &, int, int, const Term_eq &, const Param &), const Param &) const
Gives some radial fit for a given multipole, intended for anti-symmetric scalar function.
Definition: multipoles.cpp:524
virtual void set_cheb_base_t_spher(Base_spectral &) const
Gives the base using Chebyshev polynomials, for the component of a vector.
virtual const Point absol_to_num(const Point &xxx) const
Computes the numerical coordinates from the physical ones.
int nbr_unknowns_val_domain(const Val_domain &so, int mlim) const
Computes the number of true unknowns of a Val_domain.
virtual Val_domain mult_sin_phi(const Val_domain &) const
Multiplication by .
virtual Term_eq harmonics_asym(const Term_eq &, const Term_eq &, int, Term_eq(*f)(const Space &, int, int, const Term_eq &, const Param &), const Param &, const Array< double > &) const
Fit, spherical harmonic by spherical harmonic, for an anti-symmetric function.
Definition: multipoles.cpp:579
virtual void set_anti_cheb_base(Base_spectral &so) const
Sets the base to the standard one for Chebyshev polynomials for functions antisymetric in The bases ...
void export_tau_val_domain_boundary_exception_mquant(const Val_domain &eq, int mquant, int bound, Array< double > &res, int &pos_res, int ncond, const Param &param, int type_exception, const Val_domain &exception) const
Exports all the residual equations corresponding to one tensorial one on a given boundary,...
virtual void do_absol() const
Computes the absolute coordinates.
virtual double val_boundary(int, const Val_domain &, const Index &) const
Computes the value of a field at a boundary.
virtual Val_domain laplacian2(const Val_domain &, int) const
Computes the ordinary flat 2dè- Laplacian for a scalar field with an harmonic index m.
virtual Tensor change_basis_cart_to_spher(int dd, const Tensor &) const
Changes the tensorial basis from Cartsian to spherical in a given domain.
virtual ostream & print(ostream &o) const
Delegate function to virtualize the << operator.
int nbr_conditions_val_domain(const Val_domain &so, int mlim, int order) const
Computes number of discretized equations associated with a given tensorial equation in the bulk.
virtual Val_domain mult_r(const Val_domain &) const
Multiplication by .
Point center
Absolute coordinates of the center.
Definition: spheric.hpp:560
virtual int nbr_points_boundary(int, const Base_spectral &) const
Computes the number of relevant collocation points on a boundary.
virtual void set_legendre_base_p_spher(Base_spectral &) const
Gives the base using Legendre polynomials, for the component of a vector.
virtual void do_der_abs_from_der_var(const Val_domain *const *const der_var, Val_domain **const der_abs) const
Computes the derivative with respect to the absolute Cartesian coordinates from the derivative with r...
void export_tau_val_domain_boundary(const Val_domain &eq, int mlim, int bound, Array< double > &res, int &pos_res, int ncond) const
Exports all the residual equations corresponding to a tensorial one on a given boundary It makes use ...
virtual void set_legendre_base_r_mtz(Base_spectral &) const
Gives the base using Legendre polynomials, for the radial component of a vector in the MTZ context.
virtual void set_legendre_r_base(Base_spectral &) const
Gives the base using odd Legendre polynomials$ for the radius.
virtual Val_domain mult_1mrsL(const Val_domain &) const
Multiplication by .
virtual Term_eq radial_part_sym(const Space &, int, int, const Term_eq &, Term_eq(*f)(const Space &, int, int, const Term_eq &, const Param &), const Param &) const
Gives some radial fit for a given multipole, intended for symmetric scalar function.
Definition: multipoles.cpp:510
void affecte_tau_val_domain_mquant(Val_domain &so, int mquant, const Array< double > &cf, int &pos_cf) const
Affects some coefficients to a Val_domain.
virtual double get_rmax() const
Returns the maximum radius.
Definition: spheric.hpp:662
virtual void export_tau_boundary(const Tensor &, int, int, Array< double > &, int &, const Array< int > &, int n_cmp=-1, Array< int > **p_cmp=0x0) const
Exports all the residual equations corresponding to a tensorial one on a given boundary It makes use ...
virtual const Point absol_to_num_bound(const Point &, int) const
Computes the numerical coordinates from the physical ones for a point lying on a boundary.
virtual void set_legendre_base_p_mtz(Base_spectral &) const
Gives the base using Legendre polynomials, for the component of a vector in the MTZ context.
virtual void export_tau(const Tensor &, int, int, Array< double > &, int &, const Array< int > &, int n_cmp=-1, Array< int > **p_cmp=0x0) const
Exports all the residual equations corresponding to a tensorial one in the bulk.
virtual void do_cart() const
Computes the Cartesian coordinates.
virtual Val_domain div_1mx2(const Val_domain &) const
Division by .
void affecte_tau_one_coef_val_domain_mquant(Val_domain &so, int mquant, int cc, int &pos_cf) const
Sets at most one coefficient of a Val_domain to 1.
virtual Term_eq derive_flat_spher(int, char, const Term_eq &, const Metric *) const
Computes the flat derivative of a Term_eq, in spherical orthonormal coordinates.
virtual Tensor change_basis_spher_to_cart(int dd, const Tensor &) const
Changes the tensorial basis from spherical to Cartesian in a given domain.
virtual Array< int > nbr_conditions(const Tensor &, int, int, int n_cmp=-1, Array< int > **p_cmp=0x0) const
Computes number of discretized equations associated with a given tensorial equation in the bulk.
virtual Term_eq derive_flat_mtz(int, char, const Term_eq &, const Metric *) const
Computes the flat derivative of a Term_eq, in spherical coordinates where the constant radii sections...
virtual Val_domain mult_sin_theta(const Val_domain &) const
Multiplication by .
virtual int nbr_unknowns(const Tensor &, int) const
Computes the number of true unknowns of a Tensor, in a given domain.
virtual void do_radius() const
Computes the generalized radius.
virtual Val_domain mult_xm1(const Val_domain &) const
Multiplication by .
virtual void do_cart_surr() const
Computes the Cartesian coordinates over the radius.
Abstract class that implements the fonctionnalities common to all the type of domains.
Definition: space.hpp:60
Val_domain * radius
The generalized radius.
Definition: space.hpp:78
Class that gives the position inside a multi-dimensional Array.
Definition: index.hpp:38
Purely abstract class for metric handling.
Definition: metric.hpp:39
Parameter storage.
Definition: param.hpp:30
The class Point is used to store the coordinates of a point.
Definition: point.hpp:30
The class Scalar does not really implements scalars in the mathematical sense but rather tensorial co...
Definition: scalar.hpp:67
The Space_spheric class fills the space with one nucleus, several shells and a compactified domain,...
Definition: spheric.hpp:1350
virtual Array< int > get_indices_matching_non_std(int, int) const
Gives the number of the other domains, touching a given boundary.
void add_eq(System_of_eqs &syst, const char *eq, int nused=-1, Array< int > **pused=0x0)
Adds a bulk equation in all the domains of a given system (equation is assumed to be second order)
void add_eq_point(System_of_eqs &syst, const Point &pp, const char *eq)
Adds an equation being the value of some field at a given point.
double int_inf(const Scalar &) const
Computes the surface integral at infinity.
void add_eq_ori(System_of_eqs &syst, const char *eq)
Adds an equation being the value of some field at the origin.
Space_spheric(int ttype, const Point &cr, const Dim_array &nbr, const Array< double > &bounds, bool withzec=true)
Standard constructor.
void add_eq_mode_mid(System_of_eqs &syst, const char *f, int itarget, int jtarget, int ktarget)
Adds an equation saying that one coefficient of a field is zero (in the first shell)
void add_eq_int_volume(System_of_eqs &syst, const char *eq)
Adds an equation being a volume integral in the whole computational space.
void add_eq_full(System_of_eqs &syst, const char *eq, int nused=-1, Array< int > **pused=0x0)
Adds a bulk equation in all the domains of a given system (equation is assumed to be 0th order)
void add_outer_bc(System_of_eqs &syst, const char *eq, int nused=-1, Array< int > **pused=0x0)
Sets a boundary condition at the outer radius of the compactified domain.
virtual void save(FILE *) const
Saving function.
void add_matching(System_of_eqs &syst, const char *eq, int nused=-1, Array< int > **pused=0x0)
Adds a matching condition, at all the interface present in a given system.
virtual ~Space_spheric()
Destructor
void add_inner_bc(System_of_eqs &syst, const char *eq, int nused=-1, Array< int > **pused=0x0)
Sets a boundary condition at the inner radius of the innermost shell.
void add_eq_one_side(System_of_eqs &syst, const char *eq, int nused=-1, Array< int > **pused=0x0)
Adds a bulk equation in all the domains of a given system (equation is assumed to be first order)
The Space class is an ensemble of domains describing the whole space of the computation.
Definition: space.hpp:1362
Class used to describe and solve a system of equations.
Tensor handling.
Definition: tensor.hpp:149
This class is intended to describe the manage objects appearing in the equations.
Definition: term_eq.hpp:62
Class for storing the basis of decompositions of a field and its values on both the configuration and...
Definition: val_domain.hpp:69