LORENE
map.C
1 /*
2  * Methods of class Map
3  *
4  * (see file map.h for documentation)
5  *
6  */
7 
8 /*
9  * Copyright (c) 1999-2003 Eric Gourgoulhon
10  *
11  * This file is part of LORENE.
12  *
13  * LORENE is free software; you can redistribute it and/or modify
14  * it under the terms of the GNU General Public License as published by
15  * the Free Software Foundation; either version 2 of the License, or
16  * (at your option) any later version.
17  *
18  * LORENE is distributed in the hope that it will be useful,
19  * but WITHOUT ANY WARRANTY; without even the implied warranty of
20  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21  * GNU General Public License for more details.
22  *
23  * You should have received a copy of the GNU General Public License
24  * along with LORENE; if not, write to the Free Software
25  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
26  *
27  */
28 
29 
30 
31 
32 
33 /*
34  * $Id: map.C,v 1.13 2025/03/28 09:12:06 j_novak Exp $
35  * $Log: map.C,v $
36  * Revision 1.13 2025/03/28 09:12:06 j_novak
37  * New method to change the multi-grid, to be used for de-aliasing.
38  *
39  * Revision 1.12 2016/12/05 16:17:56 j_novak
40  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
41  *
42  * Revision 1.11 2014/10/13 08:53:02 j_novak
43  * Lorene classes and functions now belong to the namespace Lorene.
44  *
45  * Revision 1.10 2014/10/06 15:13:12 j_novak
46  * Modified #include directives to use c++ syntax.
47  *
48  * Revision 1.9 2008/09/29 13:23:51 j_novak
49  * Implementation of the angular mapping associated with an affine
50  * mapping. Things must be improved to take into account the domain index.
51  *
52  * Revision 1.8 2004/01/29 08:50:03 p_grandclement
53  * Modification of Map::operator==(const Map&) and addition of the surface
54  * integrales using Scalar.
55  *
56  * Revision 1.7 2003/12/30 22:53:23 e_gourgoulhon
57  * Added methods flat_met_spher() and flat_met_cart() to get
58  * flat metric associated with the coordinates described by the mapping.
59  *
60  * Revision 1.6 2003/10/15 10:30:46 e_gourgoulhon
61  * Map::set_ori: changed x,y,z to xo,yo,zo not to shadow Coord's x,y,z
62  * Map::set_rot_phi: changed phi to newphi not to shadow Coord phi.
63  *
64  * Revision 1.5 2003/06/20 14:45:06 f_limousin
65  * Add operator==
66  *
67  * Revision 1.4 2002/10/16 14:36:40 j_novak
68  * Reorganization of #include instructions of standard C++, in order to
69  * use experimental version 3 of gcc.
70  *
71  * Revision 1.3 2002/05/22 12:44:04 f_limousin
72  * Added print of ori_x, ori_y, ori_z and rot_phi in operator<<.
73  *
74  * Revision 1.2 2001/12/04 21:27:53 e_gourgoulhon
75  *
76  * All writing/reading to a binary file are now performed according to
77  * the big endian convention, whatever the system is big endian or
78  * small endian, thanks to the functions fwrite_be and fread_be
79  *
80  * Revision 1.1.1.1 2001/11/20 15:19:27 e_gourgoulhon
81  * LORENE
82  *
83  * Revision 2.12 2000/02/11 14:26:54 eric
84  * *** empty log message ***
85  *
86  * Revision 2.11 2000/02/11 13:38:23 eric
87  * Ajout de la fonction convert_absolute.
88  *
89  * Revision 2.10 2000/02/09 13:25:29 eric
90  * Mise en conformite avec le nouveau constructeur de Base_vect_spher
91  * (passage des arguments ori_x, ori_y et ori_z).
92  *
93  * Revision 2.9 2000/01/12 12:54:53 eric
94  * Ajout du Cmp null, *p_cmp_zero, et de la methode associee cmp_zero().
95  *
96  * Revision 2.8 2000/01/10 13:28:47 eric
97  * Ajout des bases vectorielles associees aux coordonnees :
98  * membres bvect_spher et bvect_cart.
99  *
100  * Revision 2.7 1999/11/24 14:39:34 eric
101  * Appel de reset_coord() dans les fonctions set_ori et set_rot_phi.
102  *
103  * Revision 2.6 1999/11/22 10:35:36 eric
104  * Fonction del_coord() rebaptisee reset_coord().
105  *
106  * Revision 2.5 1999/10/15 14:12:29 eric
107  * Suppression de l'appel a del_coord() dans le destructeur de Map.
108  *
109  * Revision 2.4 1999/10/15 09:23:04 eric
110  * *** empty log message ***
111  *
112  * Revision 2.3 1999/10/14 14:26:51 eric
113  * Depoussierage.
114  * Documentation.
115  *
116  * Revision 2.2 1999/03/01 16:57:13 eric
117  * Operateur <<
118  *
119  * Revision 2.1 1999/03/01 14:57:52 eric
120  * *** empty log message ***
121  *
122  * Revision 2.0 1999/02/15 10:42:45 hyc
123  * *** empty log message ***
124  *
125  * $Header: /cvsroot/Lorene/C++/Source/Map/map.C,v 1.13 2025/03/28 09:12:06 j_novak Exp $
126  *
127  */
128 
129 // headers C
130 #include <cmath>
131 
132 // headers Lorene
133 #include "map.h"
134 #include "cmp.h"
135 #include "metric.h"
136 #include "utilitaires.h"
137 
138  //---------------//
139  // Constructeurs //
140  //---------------//
141 
142 // Constructor from a grid
143 // -----------------------
144 namespace Lorene {
145 Map::Map(const Mg3d& mgi) : mg(&mgi),
146  ori_x(0), ori_y(0), ori_z(0), rot_phi(0),
147  bvect_spher(ori_x, ori_y, ori_z, rot_phi,
148  "Mapping orthonormal spherical basis"),
149  bvect_cart(rot_phi, "Mapping Cartesian basis"),
150  p_flat_met_spher(0x0),
151  p_flat_met_cart(0x0),
152  p_mp_angu(0x0)
153 {
154  // The Coord's are constructed by the default Coord constructor
155 
156  // The null Cmp :
157  p_cmp_zero = new Cmp(this) ;
159 }
160 
161 // Copy constructor
162 // ----------------
163 Map::Map(const Map& mp) : mg(mp.mg),
164  ori_x(mp.ori_x), ori_y(mp.ori_y), ori_z(mp.ori_z),
165  rot_phi(mp.rot_phi),
166  bvect_spher(ori_x, ori_y, ori_z, rot_phi,
167  "Mapping orthonormal spherical basis"),
168  bvect_cart(rot_phi, "Mapping Cartesian basis"),
169  p_flat_met_spher(0x0),
170  p_flat_met_cart(0x0),
171  p_mp_angu(0x0)
172 {
173  // The Coord's are constructed by the default Coord constructor
174 
175  // The null Cmp :
176  p_cmp_zero = new Cmp(this) ;
178 }
179 
180 // Constructor from file
181 // ---------------------
182 Map::Map(const Mg3d& mgi, FILE* fd) : mg(&mgi),
183  bvect_spher(0., 0., 0., 0.,
184  "Mapping orthonormal spherical basis"),
185  bvect_cart(0., "Mapping Cartesian basis"),
186  p_flat_met_spher(0x0),
187  p_flat_met_cart(0x0),
188  p_mp_angu(0x0)
189 {
190  Mg3d* mg_tmp = new Mg3d(fd) ; // la multi-grille d'origine
191  if (*mg != *mg_tmp) {
192  cout << "Map::Map(const Mg3d&, FILE*): grid not consistent !"
193  << endl ;
194  abort() ;
195  }
196  delete mg_tmp ;
197 
198  fread_be(&ori_x, sizeof(double), 1, fd) ; // ori_x
199  fread_be(&ori_y, sizeof(double), 1, fd) ; // ori_y
200  fread_be(&ori_z, sizeof(double), 1, fd) ; // ori_z
201  fread_be(&rot_phi, sizeof(double), 1, fd) ; // rot_phi
202 
206 
207  // The Coord's are constructed by the default Coord constructor
208 
209  // The null Cmp :
210  p_cmp_zero = new Cmp(this) ;
212 }
213 
214  //--------------//
215  // Destructeurs //
216  //--------------//
217 
218 // Destructeur
220  if (p_flat_met_spher != 0x0) delete p_flat_met_spher ;
221  if (p_flat_met_cart != 0x0) delete p_flat_met_cart ;
222  if (p_mp_angu != 0x0) delete p_mp_angu ;
223  delete p_cmp_zero ;
224 }
225 
226  //------------//
227  // Sauvegarde //
228  //------------//
229 
230 void Map::sauve(FILE* fd) const {
231 
232  mg->sauve(fd) ; // la multi-grille
233 
234  fwrite_be(&ori_x, sizeof(double), 1, fd) ; // ori_x
235  fwrite_be(&ori_y, sizeof(double), 1, fd) ; // ori_y
236  fwrite_be(&ori_z, sizeof(double), 1, fd) ; // ori_z
237  fwrite_be(&rot_phi, sizeof(double), 1, fd) ; // rot_phi
238 }
239 
240  //------------//
241  // Impression //
242  //------------//
243 
244 // Operateurs <<
245 ostream& operator<<(ostream& o, const Map & cv) {
246  o << "Absolute coordinates of the mapping origin: " << endl ;
247  o << " X_0, Y_0, Z_0 : " << cv.get_ori_x() << " "
248  << cv.get_ori_y() << " " << cv.get_ori_z() << endl ;
249  o << "Rotation angle between the x-axis and X-axis : "
250  << cv.get_rot_phi() << endl ;
251  cv >> o ;
252  return o ;
253 }
254 
255  //-------------------//
256  // Methodes diverses //
257  //-------------------//
258 
259 void Map::set_ori(double xo, double yo, double zo) {
260  ori_x = xo ;
261  ori_y = yo ;
262  ori_z = zo ;
263 
265 
266  reset_coord() ; // Mise a jour des Coords
267 }
268 
269 void Map::set_rot_phi(double newphi) {
270 
271  rot_phi = newphi ;
272 
275 
276  reset_coord() ; // Mise a jour des Coords
277 
278 }
279 
280 // Gestion de la memoire
281 // ---------------------
283 
284  // Coordonnees spheriques centrees sur la grille :
285  r.del_t() ;
286  tet.del_t() ;
287  phi.del_t() ;
288  sint.del_t() ;
289  cost.del_t() ;
290  sinp.del_t() ;
291  cosp.del_t() ;
292 
293  // Coordonnees cartesiennes centrees sur la grille :
294  x.del_t() ;
295  y.del_t() ;
296  z.del_t() ;
297 
298  // Coordonnees cartesiennes "absolues" :
299  xa.del_t() ;
300  ya.del_t() ;
301  za.del_t() ;
302 }
303 
304 
305 // Conversion coordonnees (X,Y,Z) --> (r, theta, phi)
306 // --------------------------------------------------
307 
308 void Map::convert_absolute(double xx, double yy, double zz,
309  double& rr, double& theta, double& pphi) const {
310 
311  // Cartesian coordinates aligned with the absolute ones but centered
312  // on the mapping origin :
313  double x1 = xx - ori_x ;
314  double y1 = yy - ori_y ;
315  double z1 = zz - ori_z ;
316 
317  // Spherical coordinates :
318  double rho2 = x1*x1 + y1*y1 ;
319  double rho = sqrt( rho2 ) ;
320  rr = sqrt(rho2 + z1*z1) ;
321  theta = atan2(rho, z1) ;
322  pphi = atan2(y1, x1) - rot_phi ; // (rotation)
323  if (pphi < 0) pphi += 2*M_PI ;
324 
325 }
326 
328 
329  if (p_flat_met_spher == 0x0) {
330  p_flat_met_spher = new Metric_flat(*this, bvect_spher) ;
331  }
332 
333  return *p_flat_met_spher ;
334 
335 }
336 
338 
339  if (p_flat_met_cart == 0x0) {
340  p_flat_met_cart = new Metric_flat(*this, bvect_cart) ;
341  }
342 
343  return *p_flat_met_cart ;
344 
345 }
346 
347 
348  void Map::set_new_grid(const Mg3d& mg2) {
349 
350 #ifndef NDEBUG
351  int nz = mg->get_nzone( );
352  assert(mg2.get_nzone() == nz) ;
353  for (int l=0; l<nz; l++) {
354  assert(mg2.get_type_r(l) == mg->get_type_r(l)) ;
355  }
356  assert(mg2.get_type_t() == mg->get_type_t()) ;
357  assert(mg2.get_type_p() == mg->get_type_p()) ;
358 #endif
359  mg = &mg2 ;
360  }
361 
362 
363 
364 
365 
366 
367 
368 
369 }
Coord xa
Absolute x coordinate.
Definition: map.h:757
int get_type_p() const
Returns the type of sampling in the direction: SYM : : symmetry with respect to the transformatio...
Definition: grilles.h:514
void set_ori(double xa0, double ya0, double za0)
Sets a new origin.
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Definition: cmp.h:446
virtual void sauve(FILE *) const
Save in a file.
Definition: map.C:230
Cmp sqrt(const Cmp &)
Square root.
Definition: cmp_math.C:223
void set_rot_phi(double phi0)
Sets a new rotation angle.
Definition: map.C:269
void set_rot_phi(double rot_phi_i)
Sets a new value to the angle rot_phi between the x –axis and the absolute frame X –axis...
double get_ori_y() const
Returns the y coordinate of the origin.
Definition: map.h:797
Lorene prototypes.
Definition: app_hor.h:67
Cmp * p_cmp_zero
The null Cmp.
Definition: map.h:740
Flat metric for tensor calculation.
Definition: metric.h:261
Base class for coordinate mappings.
Definition: map.h:697
double get_ori_x() const
Returns the x coordinate of the origin.
Definition: map.h:795
int get_type_t() const
Returns the type of sampling in the direction: SYM : : symmetry with respect to the equatorial pl...
Definition: grilles.h:504
void set_rot_phi(double rot_phi_i)
Sets a new value to the angle rot_phi between the x –axis and the absolute frame X –axis...
Map_af * p_mp_angu
Pointer on the "angular" mapping.
Definition: map.h:742
double get_rot_phi() const
Returns the angle between the x –axis and X –axis.
Definition: map.h:802
void set_new_grid(const Mg3d &new_mg)
Sets a new grid to a Map.
Definition: map.C:348
const Metric_flat & flat_met_cart() const
Returns the flat metric associated with the Cartesian coordinates and with components expressed in th...
Definition: map.C:337
Coord tet
coordinate centered on the grid
Definition: map.h:746
void set_ori(double xa0, double ya0, double za0)
Sets a new origin.
Definition: map.C:259
Coord phi
coordinate centered on the grid
Definition: map.h:747
Coord sint
Definition: map.h:748
void set_etat_zero()
Sets the logical state to ETATZERO (zero).
Definition: cmp.C:292
Metric_flat * p_flat_met_spher
Pointer onto the flat metric associated with the spherical coordinates and with components expressed ...
Definition: map.h:729
void sauve(FILE *fd, bool save_base=false) const
Saves into a file.
Definition: mg3d.C:503
double ori_y
Absolute coordinate y of the origin.
Definition: map.h:706
double ori_z
Absolute coordinate z of the origin.
Definition: map.h:707
Coord sinp
Definition: map.h:750
int get_nzone() const
Returns the number of domains.
Definition: grilles.h:467
Metric_flat * p_flat_met_cart
Pointer onto the flat metric associated with the Cartesian coordinates and with components expressed ...
Definition: map.h:734
int fwrite_be(const int *aa, int size, int nb, FILE *fich)
Writes integer(s) into a binary file according to the big endian convention.
Definition: fwrite_be.C:73
virtual void reset_coord()
Resets all the member Coords.
Definition: map.C:282
int fread_be(int *aa, int size, int nb, FILE *fich)
Reads integer(s) from a binary file according to the big endian convention.
Definition: fread_be.C:72
Coord ya
Absolute y coordinate.
Definition: map.h:758
const Mg3d * mg
Pointer on the multi-grid Mgd3 on which this is defined.
Definition: map.h:703
Multi-domain grid.
Definition: grilles.h:276
double ori_x
Absolute coordinate x of the origin.
Definition: map.h:705
Coord y
y coordinate centered on the grid
Definition: map.h:754
Coord za
Absolute z coordinate.
Definition: map.h:759
Coord cosp
Definition: map.h:751
Coord x
x coordinate centered on the grid
Definition: map.h:753
Base_vect_spher bvect_spher
Orthonormal vectorial basis associated with the coordinates of the mapping.
Definition: map.h:716
double get_ori_z() const
Returns the z coordinate of the origin.
Definition: map.h:799
int get_type_r(int l) const
Returns the type of sampling in the radial direction in domain no.
Definition: grilles.h:493
Map(const Mg3d &)
Constructor from a multi-domain 3D grid.
Definition: map.C:145
void convert_absolute(double xx, double yy, double zz, double &rr, double &theta, double &pphi) const
Determines the coordinates corresponding to given absolute Cartesian coordinates (X...
Definition: map.C:308
virtual ~Map()
Destructor.
Definition: map.C:219
Coord z
z coordinate centered on the grid
Definition: map.h:755
double rot_phi
Angle between the x –axis and X –axis.
Definition: map.h:708
const Metric_flat & flat_met_spher() const
Returns the flat metric associated with the spherical coordinates and with components expressed in th...
Definition: map.C:327
Base_vect_cart bvect_cart
Cartesian basis associated with the coordinates (x,y,z) of the mapping, i.e.
Definition: map.h:724
void del_t() const
Logical destructor (deletes the Mtbl member *c ).
Definition: coord.C:128
Coord r
r coordinate centered on the grid
Definition: map.h:745
Coord cost
Definition: map.h:749