KADATH
domain_spheric_time_nucleus.cpp
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 #include "headcpp.hpp"
21 #include "utilities.hpp"
22 #include "spheric_time.hpp"
23 #include "point.hpp"
24 #include "array_math.hpp"
25 #include "val_domain.hpp"
26 
27 namespace Kadath {
28 void coef_1d (int, Array<double>&) ;
29 void coef_i_1d (int, Array<double>&) ;
30 int der_1d (int, Array<double>&) ;
31 
32 // Standard constructor
33 Domain_spheric_time_nucleus::Domain_spheric_time_nucleus (int numdom, int ttype, double tmmin, double tmmax, double r, const Dim_array& nbr) :
34  Domain(numdom, ttype, nbr), alpha(r), tmin(tmmin), tmax(tmmax) {
35  assert (nbr.get_ndim()==2) ;
36  do_coloc() ;
37 }
38 
39 // Constructor by copy
40 Domain_spheric_time_nucleus::Domain_spheric_time_nucleus (const Domain_spheric_time_nucleus& so) : Domain(so), alpha(so.alpha), tmin(so.tmin), tmax(so.tmax) {
41 }
42 
44  fread_be (&alpha, sizeof(double), 1, fd) ;
45  fread_be (&tmin, sizeof(double), 1, fd) ;
46  fread_be (&tmax, sizeof(double), 1, fd) ;
47  do_coloc() ;
48 }
49 
50 // Destructor
51 Domain_spheric_time_nucleus::~Domain_spheric_time_nucleus() {}
52 
53 void Domain_spheric_time_nucleus::save (FILE* fd) const {
54  nbr_points.save(fd) ;
55  nbr_coefs.save(fd) ;
56  fwrite_be (&ndim, sizeof(int), 1, fd) ;
57  fwrite_be (&type_base, sizeof(int), 1, fd) ;
58  fwrite_be (&alpha, sizeof(double), 1, fd) ;
59  fwrite_be (&tmin, sizeof(double), 1, fd) ;
60  fwrite_be (&tmax, sizeof(double), 1, fd) ;
61 }
62 
63 ostream& Domain_spheric_time_nucleus::print (ostream& o) const {
64  o << "Spheric-time nucleus" << endl ;
65  o << "time goes from " << tmin << " to " << tmax << endl ;
66  o << "Rmax = " << alpha << endl ;
67  o << "Nbr pts = " << nbr_points << endl ;
68  o << endl ;
69  return o ;
70 }
71 
72 
74 
75  switch (bound) {
76  case OUTER_BC : {
77  Val_domain res (so.der_var(1)/alpha) ;
78  return res;
79  }
80  case TIME_INIT : {
81  Val_domain res (so.der_var(2)*2./(tmax-tmin)) ;
82  return res ;
83  }
84  default:
85  cerr << "Unknown boundary case in Domain_spheric_time_nucleus::der_normal" << endl ;
86  abort() ;
87  }
88 }
89 
90 // Computes the cartesian coordinates
92  for (int i=0 ; i<2 ; i++)
93  assert (coloc[i] != 0x0) ;
94  for (int i=0 ; i<2 ; i++)
95  assert (absol[i] == 0x0) ;
96  for (int i=0 ; i<2 ; i++) {
97  absol[i] = new Val_domain(this) ;
98  absol[i]->allocate_conf() ;
99  }
100  Index index (nbr_points) ;
101  do {
102  absol[0]->set(index) = alpha* ((*coloc[0])(index(0))) ;
103  absol[1]->set(index) = (*coloc[1])(index(1)) * (tmax-tmin)/2. + (tmax+tmin)/2. ;
104  }
105  while (index.inc()) ;
106 
107 }
108 
109 // Computes the radius
111 
112  for (int i=0 ; i<2 ; i++)
113  assert (coloc[i] != 0x0) ;
114  assert (radius == 0x0) ;
115  radius = new Val_domain(this) ;
116  radius->allocate_conf() ;
117  Index index (nbr_points) ;
118  do
119  radius->set(index) = alpha* ((*coloc[0])(index(0))) ;
120  while (index.inc()) ;
121 }
122 
123 
124 // Is a point inside this domain ?
125 bool Domain_spheric_time_nucleus::is_in (const Point& xx, double prec) const {
126 
127  assert (xx.get_ndim()==2) ;
128 
129  bool res = (xx(1) <= alpha+prec) ? true : false ;
130  if ((xx(2)<tmin-prec) || (xx(2)>tmax + prec))
131  res = false ;
132  return res ;
133 }
134 
135 // Convert absolute coordinates to numerical ones
137 
138  assert (is_in(abs)) ;
139  Point num(2) ;
140 
141  num.set(1) = abs(1)/alpha ;
142  num.set(2) = 2./(tmax-tmin)*(abs(2)- (tmax+tmin)/2.) ;
143 
144  return num ;
145 }
146 
147 double coloc_leg(int, int) ;
148 double coloc_leg_parity(int, int) ;
150 
151  switch (type_base) {
152  case CHEB_TYPE:
154  del_deriv() ;
155  for (int i=0 ; i<ndim ; i++)
156  coloc[i] = new Array<double> (nbr_points(i)) ;
157  for (int i=0 ; i<nbr_points(0) ; i++)
158  coloc[0]->set(i) = sin(M_PI/2.*i/(nbr_points(0)-1)) ;
159  for (int j=0 ; j<nbr_points(1) ; j++)
160  coloc[1]->set(j) = -cos(M_PI*j/(nbr_points(1)-1)) ;
161  break ;
162  case LEG_TYPE:
164  del_deriv() ;
165  for (int i=0 ; i<ndim ; i++)
166  coloc[i] = new Array<double> (nbr_points(i)) ;
167  for (int i=0 ; i<nbr_points(0) ; i++)
168  coloc[0]->set(i) = coloc_leg_parity(i, nbr_points(0)) ;
169  for (int j=0 ; j<nbr_points(1) ; j++)
170  coloc[1]->set(j) = coloc_leg(j, nbr_points(1)) ;
171  break ;
172  default :
173  cerr << "Unknown type of basis in Domain_spheric_time_nucleus::do_coloc" << endl ;
174  abort() ;
175  }
176 }
177 
178 // Base for a function symetric in z, using Chebyshev
180 
181  assert (type_base == CHEB_TYPE) ;
182  base.allocate(nbr_coefs) ;
183 
184  Index index(base.bases_1d[0]->get_dimensions()) ;
185 
186  base.def=true ;
187  base.bases_1d[1]->set(0) = CHEB ;
188  for (int j=0 ; j<nbr_coefs(1) ; j++)
189  base.bases_1d[0]->set(j) = CHEB_EVEN ;
190 }
191 
192 // Base for a function symetric in z, using Legendre
194 
195  assert (type_base == LEG_TYPE) ;
196  base.allocate(nbr_coefs) ;
197 
198  Index index(base.bases_1d[0]->get_dimensions()) ;
199  base.def=true ;
200  base.bases_1d[1]->set(0) = LEG ;
201  for (int j=0 ; j<nbr_coefs(1) ; j++)
202  base.bases_1d[0]->set(j) = LEG_EVEN ;
203  }
204 
205 // Computes the derivatives with respect to rho,Z as a function of the numerical ones.
206 void Domain_spheric_time_nucleus::do_der_abs_from_der_var(const Val_domain *const *const der_var, Val_domain **const der_abs) const {
207  // d/dr
208  der_abs[0] = new Val_domain (*der_var[0]/alpha) ;
209  // d/dt :
210  der_abs[1] = new Val_domain (*der_var[1]*2./(tmax-tmin)) ;
211 }
212 
213 // Rules for the multiplication of two basis.
215 
216  assert (a.ndim==2) ;
217  assert (b.ndim==2) ;
218 
219  Base_spectral res(2) ;
220  bool res_def = true ;
221 
222  if (!a.def)
223  res_def=false ;
224  if (!b.def)
225  res_def=false ;
226 
227  if (res_def) {
228 
229 
230  // Bases in time :
231  res.bases_1d[1] = new Array<int> (a.bases_1d[1]->get_dimensions()) ;
232  switch ((*a.bases_1d[1])(0)) {
233  case CHEB:
234  switch ((*b.bases_1d[1])(0)) {
235  case CHEB:
236  res.bases_1d[1]->set(0) = CHEB ;
237  break ;
238  default:
239  res_def = false ;
240  break ;
241  }
242  break ;
243  case LEG:
244  switch ((*b.bases_1d[1])(0)) {
245  case LEG:
246  res.bases_1d[1]->set(0) = LEG ;
247  break ;
248  default:
249  res_def = false ;
250  break ;
251  }
252  break ;
253 
254  default:
255  res_def = false ;
256  break ;
257  }
258 
259  // Base in r :
260  Index index_0 (a.bases_1d[0]->get_dimensions()) ;
261  res.bases_1d[0] = new Array<int> (a.bases_1d[0]->get_dimensions()) ;
262  do {
263  switch ((*a.bases_1d[0])(index_0)) {
264  case CHEB_EVEN:
265  switch ((*b.bases_1d[0])(index_0)) {
266  case CHEB_EVEN:
267  res.bases_1d[0]->set(index_0) = CHEB_EVEN ;
268  break ;
269  case CHEB_ODD:
270  res.bases_1d[0]->set(index_0) = CHEB_ODD ;
271  break ;
272  default:
273  res_def = false ;
274  break ;
275  }
276  break ;
277  case CHEB_ODD:
278  switch ((*b.bases_1d[0])(index_0)) {
279  case CHEB_EVEN:
280  res.bases_1d[0]->set(index_0) = CHEB_ODD ;
281  break ;
282  case CHEB_ODD:
283  res.bases_1d[0]->set(index_0) = CHEB_EVEN ;
284  break ;
285  default:
286  res_def = false ;
287  break ;
288  }
289  break ;
290  case LEG_EVEN:
291  switch ((*b.bases_1d[0])(index_0)) {
292  case LEG_EVEN:
293  res.bases_1d[0]->set(index_0) = LEG_EVEN ;
294  break ;
295  case LEG_ODD:
296  res.bases_1d[0]->set(index_0) = LEG_ODD ;
297  break ;
298  default:
299  res_def = false ;
300  break ;
301  }
302  break ;
303  case LEG_ODD:
304  switch ((*b.bases_1d[0])(index_0)) {
305  case LEG_EVEN:
306  res.bases_1d[0]->set(index_0) = LEG_ODD ;
307  break ;
308  case LEG_ODD:
309  res.bases_1d[0]->set(index_0) = LEG_EVEN ;
310  break ;
311  default:
312  res_def = false ;
313  break ;
314  }
315  break ;
316  default:
317  res_def = false ;
318  break ;
319  }
320  }
321  while (index_0.inc()) ;
322  }
323  if (!res_def)
324  for (int dim=0 ; dim<a.ndim ; dim++)
325  if (res.bases_1d[dim]!= 0x0) {
326  delete res.bases_1d[dim] ;
327  res.bases_1d[dim] = 0x0 ;
328  }
329  res.def = res_def ;
330  return res ;
331 }
332 }
Class for storing the basis of decompositions of a field.
Bases_container bases_1d
Arrays containing the various basis of decomposition.
void allocate(const Dim_array &nbr_coefs)
Allocates the various arrays, for a given number of coefficients.
bool def
true if the Base_spectral is defined and false otherwise.
int ndim
Number of dimensions.
Class for storing the dimensions of an array.
Definition: dim_array.hpp:34
int get_ndim() const
Returns the number of dimensions.
Definition: dim_array.hpp:63
void save(FILE *) const
Save function.
Definition: dim_array.cpp:32
Class for a 2-dimensional spherical domain containing the origin and a symetry with respect to the pl...
virtual ostream & print(ostream &o) const
Delegate function to virtualize the << operator.
virtual const Point absol_to_num(const Point &) const
Computes the numerical coordinates from the physical ones.
Domain_spheric_time_nucleus(int num, int ttype, double tmmin, double tmmax, double radius, const Dim_array &nbr)
Standard constructor :
virtual bool is_in(const Point &xx, double prec=1e-13) const
Check whether a point lies inside Domain.
virtual void do_absol() const
Computes the absolute coordinates.
virtual void set_legendre_base(Base_spectral &) const
Gives the standard base for Legendre polynomials.
virtual void set_cheb_base(Base_spectral &) const
Gives the standard base for Chebyshev polynomials.
virtual Base_spectral mult(const Base_spectral &, const Base_spectral &) const
Method for the multiplication of two Base_spectral.
double alpha
Relates the numerical to the physical radii.
virtual void do_coloc()
Computes the colocation points.
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 save(FILE *) const
Saving function.
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.
Abstract class that implements the fonctionnalities common to all the type of domains.
Definition: space.hpp:60
virtual void del_deriv()
Destroys the derivated members (like coloc, cart and radius), when changing the type of colocation po...
Definition: domain.cpp:77
Val_domain * radius
The generalized radius.
Definition: space.hpp:78
Memory_mapped_array< Val_domain * > absol
Asbolute coordinates (if defined ; usually Cartesian-like)
Definition: space.hpp:76
int ndim
Number of dimensions.
Definition: space.hpp:64
Dim_array nbr_coefs
Number of coefficients.
Definition: space.hpp:66
Dim_array nbr_points
Number of colocation points.
Definition: space.hpp:65
int type_base
Type of colocation point :
Definition: space.hpp:73
Memory_mapped_array< Array< double > * > coloc
Colocation points in each dimension (stored in ndim 1d- arrays)
Definition: space.hpp:75
Class that gives the position inside a multi-dimensional Array.
Definition: index.hpp:38
bool inc(int increm, int var=0)
Increments the position of the Index.
Definition: index.hpp:99
The class Point is used to store the coordinates of a point.
Definition: point.hpp:30
const int & get_ndim() const
Returns the number of dimensions.
Definition: point.hpp:51
double & set(int i)
Read/write of a coordinate.
Definition: point.hpp:47
Class for storing the basis of decompositions of a field and its values on both the configuration and...
Definition: val_domain.hpp:69
double & set(const Index &pos)
Read/write the value of the field in the configuration space.
Definition: val_domain.cpp:171
Val_domain der_var(int i) const
Computes the derivative with respect to a numerical coordinate.
Definition: val_domain.cpp:670
void allocate_conf()
Allocates the values in the configuration space and destroys the values in the coefficients space.
Definition: val_domain.cpp:209