30char map_C[] =
"$Header: /cvsroot/Lorene/C++/Source/Map/map.C,v 1.11 2014/10/13 08:53:02 j_novak Exp $" ;
130#include "utilitaires.h"
139Map::Map(
const Mg3d& mgi) : mg(&mgi),
140 ori_x(0), ori_y(0), ori_z(0), rot_phi(0),
142 "Mapping orthonormal spherical basis"),
143 bvect_cart(rot_phi,
"Mapping Cartesian basis"),
157Map::Map(
const Map& mp) : mg(mp.mg),
158 ori_x(mp.ori_x), ori_y(mp.ori_y), ori_z(mp.ori_z),
161 "Mapping orthonormal spherical basis"),
162 bvect_cart(rot_phi,
"Mapping Cartesian basis"),
176Map::Map(
const Mg3d& mgi, FILE* fd) : mg(&mgi),
178 "Mapping orthonormal spherical basis"),
185 if (*mg != *mg_tmp) {
186 cout <<
"Map::Map(const Mg3d&, FILE*): grid not consistent !"
192 fread_be(&ori_x,
sizeof(
double), 1, fd) ;
193 fread_be(&ori_y,
sizeof(
double), 1, fd) ;
194 fread_be(&ori_z,
sizeof(
double), 1, fd) ;
195 fread_be(&rot_phi,
sizeof(
double), 1, fd) ;
224void Map::sauve(FILE* fd)
const {
228 fwrite_be(&ori_x,
sizeof(
double), 1, fd) ;
229 fwrite_be(&ori_y,
sizeof(
double), 1, fd) ;
230 fwrite_be(&ori_z,
sizeof(
double), 1, fd) ;
231 fwrite_be(&rot_phi,
sizeof(
double), 1, fd) ;
239ostream& operator<<(ostream& o,
const Map & cv) {
240 o <<
"Absolute coordinates of the mapping origin: " << endl ;
241 o <<
" X_0, Y_0, Z_0 : " << cv.get_ori_x() <<
" "
242 << cv.get_ori_y() <<
" " << cv.get_ori_z() << endl ;
243 o <<
"Rotation angle between the x-axis and X-axis : "
244 << cv.get_rot_phi() << endl ;
253void Map::set_ori(
double xo,
double yo,
double zo) {
258 bvect_spher.
set_ori(ori_x, ori_y, ori_z) ;
263void Map::set_rot_phi(
double newphi) {
276void Map::reset_coord() {
302void Map::convert_absolute(
double xx,
double yy,
double zz,
303 double& rr,
double& theta,
double& pphi)
const {
307 double x1 = xx - ori_x ;
308 double y1 = yy - ori_y ;
309 double z1 = zz - ori_z ;
312 double rho2 = x1*x1 + y1*y1 ;
313 double rho =
sqrt( rho2 ) ;
314 rr =
sqrt(rho2 + z1*z1) ;
315 theta = atan2(rho, z1) ;
316 pphi = atan2(y1, x1) - rot_phi ;
317 if (pphi < 0) pphi += 2*M_PI ;
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.
void set_ori(double xa0, double ya0, double za0)
Sets a new origin.
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.
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
void del_t() const
Logical destructor (deletes the Mtbl member *c ).
Flat metric for tensor calculation.
Cmp sqrt(const Cmp &)
Square root.
Base_vect_spher bvect_spher
Base class for coordinate mappings.
int fread_be(int *aa, int size, int nb, FILE *fich)
Reads integer(s) from a binary file according to the big endian convention.
Map_af * p_mp_angu
Pointer on the "angular" mapping.
virtual void reset_coord()
Resets all the member Coords.
Coord z
z coordinate centered on the grid
Metric_flat * p_flat_met_spher
Pointer onto the flat metric associated with the spherical coordinates and with components expressed ...
Base_vect_cart bvect_cart
Cartesian basis associated with the coordinates (x,y,z) of the mapping, i.e.
Metric_flat * p_flat_met_cart
Pointer onto the flat metric associated with the Cartesian coordinates and with components expressed ...
Cmp * p_cmp_zero
The null Cmp.
Coord y
y coordinate centered on the grid
Map(const Mg3d &)
Constructor from a multi-domain 3D grid.
Coord phi
coordinate centered on the grid
Coord tet
coordinate centered on the grid
Coord x
x coordinate centered on the grid
Coord xa
Absolute x coordinate.
Coord za
Absolute z coordinate.
Coord r
r coordinate centered on the grid
Coord ya
Absolute y coordinate.