LORENE
et_rot_diff.h
1/*
2 * Definition of Lorene class Et_rot_diff
3 *
4 */
5
6/*
7 * Copyright (c) 2001 Eric Gourgoulhon
8 *
9 * This file is part of LORENE.
10 *
11 * LORENE is free software; you can redistribute it and/or modify
12 * it under the terms of the GNU General Public License as published by
13 * the Free Software Foundation; either version 2 of the License, or
14 * (at your option) any later version.
15 *
16 * LORENE is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
19 * GNU General Public License for more details.
20 *
21 * You should have received a copy of the GNU General Public License
22 * along with LORENE; if not, write to the Free Software
23 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
24 *
25 */
26
27
28#ifndef __ET_ROT_DIFF_H_
29#define __ET_ROT_DIFF_H_
30
31/*
32 * $Id: et_rot_diff.h,v 1.6 2021/04/13 11:24:53 j_novak Exp $
33 * $Log: et_rot_diff.h,v $
34 * Revision 1.6 2021/04/13 11:24:53 j_novak
35 * Corrected a bug in Etoile_rot::fait_shift() which was missing the np=1 case.
36 *
37 * Revision 1.5 2014/10/13 08:52:34 j_novak
38 * Lorene classes and functions now belong to the namespace Lorene.
39 *
40 * Revision 1.4 2005/10/05 15:14:47 j_novak
41 * Added a Param* as parameter of Etoile_rot::equilibrium
42 *
43 * Revision 1.3 2004/03/22 13:12:41 j_novak
44 * Modification of comments to use doxygen instead of doc++
45 *
46 * Revision 1.2 2002/09/13 09:17:33 j_novak
47 * Modif. commentaires
48 *
49 * Revision 1.1.1.1 2001/11/20 15:19:27 e_gourgoulhon
50 * LORENE
51 *
52 * Revision 1.2 2001/10/25 09:20:35 eric
53 * Ajout de la fonction virtuelle display_poly.
54 *
55 * Revision 1.1 2001/10/19 08:17:59 eric
56 * Initial revision
57 *
58 *
59 * $Header: /cvsroot/Lorene/C++/Include/et_rot_diff.h,v 1.6 2021/04/13 11:24:53 j_novak Exp $
60 *
61 */
62
63// Headers Lorene
64#include "etoile.h"
65
66namespace Lorene {
72
73class Et_rot_diff : public Etoile_rot {
74
75 // Data :
76 // -----
77 protected:
88 double (*frot)(double, const Tbl&) ;
89
96 double (*primfrot)(double, const Tbl&) ;
97
106
109
110 double omega_min ;
111 double omega_max ;
112
115
116 // Constructors - Destructor
117 // -------------------------
118 public:
134 Et_rot_diff(Map& mp_i, int nzet_i, bool relat, const Eos& eos_i,
135 double (*frot_i)(double, const Tbl&),
136 double (*primfrot_i)(double, const Tbl&),
137 const Tbl& par_frot_i) ;
138
139 Et_rot_diff(const Et_rot_diff& ) ;
140
151 Et_rot_diff(Map& mp_i, const Eos& eos_i, FILE* fich,
152 double (*frot_i)(double, const Tbl&),
153 double (*primfrot_i)(double, const Tbl&) ) ;
154
155 virtual ~Et_rot_diff() ;
156
157
158 // Memory management
159 // -----------------
160
161 // Everything is inherited from Etoile_rot
162
163 // Mutators / assignment
164 // ---------------------
165 public:
167 void operator=(const Et_rot_diff& ) ;
168
169 // Accessors
170 // ---------
171 public:
173 const Tenseur& get_omega_field() const {return omega_field;} ;
174
178 virtual double get_omega_c() const ;
179
180 // Outputs
181 // -------
182 public:
183 virtual void sauve(FILE *) const ;
184
186 virtual void display_poly(ostream& ) const ;
187
188 protected:
190 virtual ostream& operator>>(ostream& ) const ;
191
192
193 // Computational routines
194 // ----------------------
195 public:
196
197 virtual double tsw() const ;
198
210 virtual void hydro_euler() ;
211
230 void fait_omega_field(double omeg_min, double omeg_max,
231 double precis, int nitermax) ;
232
234 void fait_prim_field() ;
235
251 double funct_omega(double omeg) const ;
252
261 double prim_funct_omega(double omeg) const ;
262
342 virtual void equilibrium(double ent_c, double omega0, double fact_omega,
343 int nzadapt, const Tbl& ent_limit,
344 const Itbl& icontrol, const Tbl& control,
345 double mbar_wanted, double aexp_mass,
346 Tbl& diff, Param* = 0x0) ;
347
348 };
349
350}
351#endif
Equation of state base class.
Definition eos.h:206
double omega_min
Minimum value of .
void fait_omega_field(double omeg_min, double omeg_max, double precis, int nitermax)
Computes (member omega_field ).
Tenseur prim_field
Field .
virtual void hydro_euler()
Computes the hydrodynamical quantities relative to the Eulerian observer from those in the fluid fram...
virtual void display_poly(ostream &) const
Display in polytropic units.
double omega_max
Maximum value of .
virtual double tsw() const
Ratio T/W.
Tenseur omega_field
Field .
void fait_prim_field()
Computes the member prim_field from omega_field .
Et_rot_diff(Map &mp_i, int nzet_i, bool relat, const Eos &eos_i, double(*frot_i)(double, const Tbl &), double(*primfrot_i)(double, const Tbl &), const Tbl &par_frot_i)
Standard constructor.
Definition et_rot_diff.C:85
virtual void sauve(FILE *) const
Save in a file.
virtual ~Et_rot_diff()
Destructor.
Tbl par_frot
Parameters of the function .
double(* primfrot)(double, const Tbl &)
Primitive of the function , which vanishes at the stellar center.
Definition et_rot_diff.h:96
double prim_funct_omega(double omeg) const
Evaluates the primitive of , where F is the function defining the rotation profile.
void operator=(const Et_rot_diff &)
Assignment to another Et_rot_diff.
virtual double get_omega_c() const
Returns the central value of the rotation angular velocity ([f_unit] ).
const Tenseur & get_omega_field() const
Returns the angular velocity field .
double funct_omega(double omeg) const
Evaluates , where F is the function defining the rotation profile.
virtual ostream & operator>>(ostream &) const
Operator >> (virtual function called by the operator <<).
virtual void equilibrium(double ent_c, double omega0, double fact_omega, int nzadapt, const Tbl &ent_limit, const Itbl &icontrol, const Tbl &control, double mbar_wanted, double aexp_mass, Tbl &diff, Param *=0x0)
Computes an equilibrium configuration.
double(* frot)(double, const Tbl &)
Function defining the rotation profile.
Definition et_rot_diff.h:88
Etoile_rot(Map &mp_i, int nzet_i, bool relat, const Eos &eos_i)
Standard constructor.
Definition etoile_rot.C:146
Basic integer array class.
Definition itbl.h:122
Parameter storage.
Definition param.h:125
Basic array class.
Definition tbl.h:161
Tensor handling *** DEPRECATED : use class Tensor instead ***.
Definition tenseur.h:304
Lorene prototypes.
Definition app_hor.h:67
Map(const Mg3d &)
Constructor from a multi-domain 3D grid.
Definition map.C:142