/usr/share/doc/lorene-doc/refguide/star__bhns__vel__pot_8C_source.html is in lorene-doc 0.0.0~cvs20161116+dfsg-1.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 | <!DOCTYPE html PUBLIC "-//W3C//DTD XHTML 1.0 Transitional//EN" "http://www.w3.org/TR/xhtml1/DTD/xhtml1-transitional.dtd">
<html xmlns="http://www.w3.org/1999/xhtml">
<head>
<meta http-equiv="Content-Type" content="text/xhtml;charset=UTF-8"/>
<meta http-equiv="X-UA-Compatible" content="IE=9"/>
<meta name="generator" content="Doxygen 1.8.12"/>
<meta name="viewport" content="width=device-width, initial-scale=1"/>
<title>LORENE: star_bhns_vel_pot.C Source File</title>
<link href="tabs.css" rel="stylesheet" type="text/css"/>
<script type="text/javascript" src="jquery.js"></script>
<script type="text/javascript" src="dynsections.js"></script>
<link href="navtree.css" rel="stylesheet" type="text/css"/>
<script type="text/javascript" src="resize.js"></script>
<script type="text/javascript" src="navtreedata.js"></script>
<script type="text/javascript" src="navtree.js"></script>
<script type="text/javascript">
$(document).ready(initResizable);
</script>
<link href="doxygen.css" rel="stylesheet" type="text/css" />
</head>
<body>
<div id="top"><!-- do not remove this div, it is closed by doxygen! -->
<div id="titlearea">
<table cellspacing="0" cellpadding="0">
<tbody>
<tr style="height: 56px;">
<td id="projectalign" style="padding-left: 0.5em;">
<div id="projectname">LORENE
</div>
</td>
</tr>
</tbody>
</table>
</div>
<!-- end header part -->
<!-- Generated by Doxygen 1.8.12 -->
<script type="text/javascript" src="menudata.js"></script>
<script type="text/javascript" src="menu.js"></script>
<script type="text/javascript">
$(function() {
initMenu('',false,false,'search.php','Search');
});
</script>
<div id="main-nav"></div>
</div><!-- top -->
<div id="side-nav" class="ui-resizable side-nav-resizable">
<div id="nav-tree">
<div id="nav-tree-contents">
<div id="nav-sync" class="sync"></div>
</div>
</div>
<div id="splitbar" style="-moz-user-select:none;"
class="ui-resizable-handle">
</div>
</div>
<script type="text/javascript">
$(document).ready(function(){initNavTree('star__bhns__vel__pot_8C_source.html','');});
</script>
<div id="doc-content">
<div class="header">
<div class="headertitle">
<div class="title">star_bhns_vel_pot.C</div> </div>
</div><!--header-->
<div class="contents">
<div class="fragment"><div class="line"><a name="l00001"></a><span class="lineno"> 1</span> <span class="comment">/*</span></div><div class="line"><a name="l00002"></a><span class="lineno"> 2</span> <span class="comment"> * Method of class Star_bhns to compute the velocity scalar potential $\psi$</span></div><div class="line"><a name="l00003"></a><span class="lineno"> 3</span> <span class="comment"> * by solving the continuity equation.</span></div><div class="line"><a name="l00004"></a><span class="lineno"> 4</span> <span class="comment"> *</span></div><div class="line"><a name="l00005"></a><span class="lineno"> 5</span> <span class="comment"> * (see file star_bhns.h for documentation).</span></div><div class="line"><a name="l00006"></a><span class="lineno"> 6</span> <span class="comment"> *</span></div><div class="line"><a name="l00007"></a><span class="lineno"> 7</span> <span class="comment"> */</span></div><div class="line"><a name="l00008"></a><span class="lineno"> 8</span> </div><div class="line"><a name="l00009"></a><span class="lineno"> 9</span> <span class="comment">/*</span></div><div class="line"><a name="l00010"></a><span class="lineno"> 10</span> <span class="comment"> * Copyright (c) 2005-2007 Keisuke Taniguchi</span></div><div class="line"><a name="l00011"></a><span class="lineno"> 11</span> <span class="comment"> *</span></div><div class="line"><a name="l00012"></a><span class="lineno"> 12</span> <span class="comment"> * This file is part of LORENE.</span></div><div class="line"><a name="l00013"></a><span class="lineno"> 13</span> <span class="comment"> *</span></div><div class="line"><a name="l00014"></a><span class="lineno"> 14</span> <span class="comment"> * LORENE is free software; you can redistribute it and/or modify</span></div><div class="line"><a name="l00015"></a><span class="lineno"> 15</span> <span class="comment"> * it under the terms of the GNU General Public License version 2</span></div><div class="line"><a name="l00016"></a><span class="lineno"> 16</span> <span class="comment"> * as published by the Free Software Foundation.</span></div><div class="line"><a name="l00017"></a><span class="lineno"> 17</span> <span class="comment"> *</span></div><div class="line"><a name="l00018"></a><span class="lineno"> 18</span> <span class="comment"> * LORENE is distributed in the hope that it will be useful,</span></div><div class="line"><a name="l00019"></a><span class="lineno"> 19</span> <span class="comment"> * but WITHOUT ANY WARRANTY; without even the implied warranty of</span></div><div class="line"><a name="l00020"></a><span class="lineno"> 20</span> <span class="comment"> * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the</span></div><div class="line"><a name="l00021"></a><span class="lineno"> 21</span> <span class="comment"> * GNU General Public License for more details.</span></div><div class="line"><a name="l00022"></a><span class="lineno"> 22</span> <span class="comment"> *</span></div><div class="line"><a name="l00023"></a><span class="lineno"> 23</span> <span class="comment"> * You should have received a copy of the GNU General Public License</span></div><div class="line"><a name="l00024"></a><span class="lineno"> 24</span> <span class="comment"> * along with LORENE; if not, write to the Free Software</span></div><div class="line"><a name="l00025"></a><span class="lineno"> 25</span> <span class="comment"> * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA</span></div><div class="line"><a name="l00026"></a><span class="lineno"> 26</span> <span class="comment"> *</span></div><div class="line"><a name="l00027"></a><span class="lineno"> 27</span> <span class="comment"> */</span></div><div class="line"><a name="l00028"></a><span class="lineno"> 28</span> </div><div class="line"><a name="l00029"></a><span class="lineno"> 29</span> <span class="keywordtype">char</span> star_bhns_vel_pot_C[] = <span class="stringliteral">"$Header: /cvsroot/Lorene/C++/Source/Star_bhns/star_bhns_vel_pot.C,v 1.3 2014/10/13 08:53:41 j_novak Exp $"</span> ;</div><div class="line"><a name="l00030"></a><span class="lineno"> 30</span> </div><div class="line"><a name="l00031"></a><span class="lineno"> 31</span> <span class="comment">/*</span></div><div class="line"><a name="l00032"></a><span class="lineno"> 32</span> <span class="comment"> * $Id: star_bhns_vel_pot.C,v 1.3 2014/10/13 08:53:41 j_novak Exp $</span></div><div class="line"><a name="l00033"></a><span class="lineno"> 33</span> <span class="comment"> * $Log: star_bhns_vel_pot.C,v $</span></div><div class="line"><a name="l00034"></a><span class="lineno"> 34</span> <span class="comment"> * Revision 1.3 2014/10/13 08:53:41 j_novak</span></div><div class="line"><a name="l00035"></a><span class="lineno"> 35</span> <span class="comment"> * Lorene classes and functions now belong to the namespace Lorene.</span></div><div class="line"><a name="l00036"></a><span class="lineno"> 36</span> <span class="comment"> *</span></div><div class="line"><a name="l00037"></a><span class="lineno"> 37</span> <span class="comment"> * Revision 1.2 2008/05/15 19:20:29 k_taniguchi</span></div><div class="line"><a name="l00038"></a><span class="lineno"> 38</span> <span class="comment"> * Change of some parameters.</span></div><div class="line"><a name="l00039"></a><span class="lineno"> 39</span> <span class="comment"> *</span></div><div class="line"><a name="l00040"></a><span class="lineno"> 40</span> <span class="comment"> * Revision 1.1 2007/06/22 01:33:14 k_taniguchi</span></div><div class="line"><a name="l00041"></a><span class="lineno"> 41</span> <span class="comment"> * *** empty log message ***</span></div><div class="line"><a name="l00042"></a><span class="lineno"> 42</span> <span class="comment"> *</span></div><div class="line"><a name="l00043"></a><span class="lineno"> 43</span> <span class="comment"> *</span></div><div class="line"><a name="l00044"></a><span class="lineno"> 44</span> <span class="comment"> * $Header: /cvsroot/Lorene/C++/Source/Star_bhns/star_bhns_vel_pot.C,v 1.3 2014/10/13 08:53:41 j_novak Exp $</span></div><div class="line"><a name="l00045"></a><span class="lineno"> 45</span> <span class="comment"> *</span></div><div class="line"><a name="l00046"></a><span class="lineno"> 46</span> <span class="comment"> */</span></div><div class="line"><a name="l00047"></a><span class="lineno"> 47</span> </div><div class="line"><a name="l00048"></a><span class="lineno"> 48</span> <span class="comment">// C++ headers</span></div><div class="line"><a name="l00049"></a><span class="lineno"> 49</span> <span class="comment">//#include <></span></div><div class="line"><a name="l00050"></a><span class="lineno"> 50</span> </div><div class="line"><a name="l00051"></a><span class="lineno"> 51</span> <span class="comment">// C headers</span></div><div class="line"><a name="l00052"></a><span class="lineno"> 52</span> <span class="comment">//#include <></span></div><div class="line"><a name="l00053"></a><span class="lineno"> 53</span> </div><div class="line"><a name="l00054"></a><span class="lineno"> 54</span> <span class="comment">// Lorene headers</span></div><div class="line"><a name="l00055"></a><span class="lineno"> 55</span> <span class="preprocessor">#include "star_bhns.h"</span></div><div class="line"><a name="l00056"></a><span class="lineno"> 56</span> <span class="preprocessor">#include "eos.h"</span></div><div class="line"><a name="l00057"></a><span class="lineno"> 57</span> <span class="preprocessor">#include "param.h"</span></div><div class="line"><a name="l00058"></a><span class="lineno"> 58</span> <span class="preprocessor">#include "cmp.h"</span></div><div class="line"><a name="l00059"></a><span class="lineno"> 59</span> <span class="preprocessor">#include "tenseur.h"</span></div><div class="line"><a name="l00060"></a><span class="lineno"> 60</span> <span class="preprocessor">#include "utilitaires.h"</span></div><div class="line"><a name="l00061"></a><span class="lineno"> 61</span> <span class="preprocessor">#include "unites.h"</span></div><div class="line"><a name="l00062"></a><span class="lineno"> 62</span> </div><div class="line"><a name="l00063"></a><span class="lineno"> 63</span> <span class="comment">// Local prototype</span></div><div class="line"><a name="l00064"></a><span class="lineno"> 64</span> <span class="keyword">namespace </span><a class="code" href="namespaceLorene.html">Lorene</a> {</div><div class="line"><a name="l00065"></a><span class="lineno"> 65</span> Cmp raccord_c1(<span class="keyword">const</span> Cmp& uu, <span class="keywordtype">int</span> l1) ;</div><div class="line"><a name="l00066"></a><span class="lineno"> 66</span> </div><div class="line"><a name="l00067"></a><span class="lineno"><a class="line" href="classLorene_1_1Star__bhns.html#a26f0c79991a8a68ed50e4425271b762d"> 67</a></span> <span class="keywordtype">double</span> <a class="code" href="classLorene_1_1Star__bhns.html#a26f0c79991a8a68ed50e4425271b762d">Star_bhns::velo_pot_bhns</a>(<span class="keyword">const</span> <span class="keywordtype">double</span>& mass_bh, <span class="keyword">const</span> <span class="keywordtype">double</span>& sepa,</div><div class="line"><a name="l00068"></a><span class="lineno"> 68</span>  <span class="keywordtype">bool</span> kerrschild, <span class="keywordtype">int</span> mermax, <span class="keywordtype">double</span> precis,</div><div class="line"><a name="l00069"></a><span class="lineno"> 69</span>  <span class="keywordtype">double</span> relax) {</div><div class="line"><a name="l00070"></a><span class="lineno"> 70</span> </div><div class="line"><a name="l00071"></a><span class="lineno"> 71</span>  <span class="comment">// Fundamental constants and units</span></div><div class="line"><a name="l00072"></a><span class="lineno"> 72</span>  <span class="comment">// -------------------------------</span></div><div class="line"><a name="l00073"></a><span class="lineno"> 73</span>  <span class="keyword">using namespace </span><a class="code" href="namespaceUnites.html">Unites</a> ;</div><div class="line"><a name="l00074"></a><span class="lineno"> 74</span> </div><div class="line"><a name="l00075"></a><span class="lineno"> 75</span>  <span class="keywordtype">int</span> nzm1 = <a class="code" href="classLorene_1_1Star.html#ac8b026b43b98b1f8f25cc57d96cb3243">mp</a>.<a class="code" href="classLorene_1_1Map.html#a4df8c6a943a85ef314ea636240ef31a7">get_mg</a>()-><a class="code" href="classLorene_1_1Mg3d.html#a99fb0e7a3f1f0f3e405624855200055e">get_nzone</a>() - 1 ;</div><div class="line"><a name="l00076"></a><span class="lineno"> 76</span> </div><div class="line"><a name="l00077"></a><span class="lineno"> 77</span>  <span class="comment">//--------------------------------</span></div><div class="line"><a name="l00078"></a><span class="lineno"> 78</span>  <span class="comment">// Specific relativistic enthalpy ---> hhh</span></div><div class="line"><a name="l00079"></a><span class="lineno"> 79</span>  <span class="comment">//--------------------------------</span></div><div class="line"><a name="l00080"></a><span class="lineno"> 80</span> </div><div class="line"><a name="l00081"></a><span class="lineno"> 81</span>  <a class="code" href="classLorene_1_1Scalar.html">Scalar</a> hhh = <a class="code" href="group__cmp__m.html#ga66a0d18b41ceae25049869b5ccdca44d">exp</a>(<a class="code" href="classLorene_1_1Star.html#abd29e4e86f0f9f65b941948e7baa286f">ent</a>) ;</div><div class="line"><a name="l00082"></a><span class="lineno"> 82</span>  hhh.<a class="code" href="classLorene_1_1Scalar.html#a4eae068572578d6fdd0ab74b90c17b76">std_spectral_base</a>() ;</div><div class="line"><a name="l00083"></a><span class="lineno"> 83</span> </div><div class="line"><a name="l00084"></a><span class="lineno"> 84</span>  <span class="comment">//---------------------------------------------------</span></div><div class="line"><a name="l00085"></a><span class="lineno"> 85</span>  <span class="comment">// Computation of W^i = psi4 h Gamma_n B^i/lapse_tot</span></div><div class="line"><a name="l00086"></a><span class="lineno"> 86</span>  <span class="comment">// See Eq (62) from Gourgoulhon et al. (2001)</span></div><div class="line"><a name="l00087"></a><span class="lineno"> 87</span>  <span class="comment">//---------------------------------------------------</span></div><div class="line"><a name="l00088"></a><span class="lineno"> 88</span> </div><div class="line"><a name="l00089"></a><span class="lineno"> 89</span>  <a class="code" href="classLorene_1_1Vector.html">Vector</a> www = hhh * <a class="code" href="classLorene_1_1Star.html#a8a0796aa1b273029e32cdfb104f330cb">gam_euler</a> * <a class="code" href="classLorene_1_1Star__bhns.html#a285e3e765784ecabd7b5043df2ce5c9c">bsn</a> * <a class="code" href="classLorene_1_1Star__bhns.html#a6d677eed1bff187093e1fb8c19d11538">psi4</a> ;</div><div class="line"><a name="l00090"></a><span class="lineno"> 90</span> </div><div class="line"><a name="l00091"></a><span class="lineno"> 91</span>  www.<a class="code" href="classLorene_1_1Vector.html#ad916ce65d0f862b7ea905af6e0ab9785">change_triad</a>( <a class="code" href="classLorene_1_1Star.html#ac8b026b43b98b1f8f25cc57d96cb3243">mp</a>.<a class="code" href="classLorene_1_1Map.html#adfc2b4aa1ef4f22a109bce235e7eb6b7">get_bvect_cart</a>() ) ; <span class="comment">// components on the mapping</span></div><div class="line"><a name="l00092"></a><span class="lineno"> 92</span>  <span class="comment">// Cartesian basis</span></div><div class="line"><a name="l00093"></a><span class="lineno"> 93</span> </div><div class="line"><a name="l00094"></a><span class="lineno"> 94</span>  <span class="comment">//-------------------------------------------------</span></div><div class="line"><a name="l00095"></a><span class="lineno"> 95</span>  <span class="comment">// Constant value of W^i at the center of the star</span></div><div class="line"><a name="l00096"></a><span class="lineno"> 96</span>  <span class="comment">//-------------------------------------------------</span></div><div class="line"><a name="l00097"></a><span class="lineno"> 97</span> </div><div class="line"><a name="l00098"></a><span class="lineno"> 98</span>  <a class="code" href="classLorene_1_1Vector.html">Vector</a> v_orb(<a class="code" href="classLorene_1_1Star.html#ac8b026b43b98b1f8f25cc57d96cb3243">mp</a>, CON, <a class="code" href="classLorene_1_1Star.html#ac8b026b43b98b1f8f25cc57d96cb3243">mp</a>.<a class="code" href="classLorene_1_1Map.html#adfc2b4aa1ef4f22a109bce235e7eb6b7">get_bvect_cart</a>()) ;</div><div class="line"><a name="l00099"></a><span class="lineno"> 99</span>  v_orb.<a class="code" href="classLorene_1_1Tensor.html#ad0c0b10eae3f7f145600569be438ba80">set_etat_qcq</a>() ;</div><div class="line"><a name="l00100"></a><span class="lineno"> 100</span> </div><div class="line"><a name="l00101"></a><span class="lineno"> 101</span>  <span class="keywordflow">for</span> (<span class="keywordtype">int</span> i=1; i<=3; i++) {</div><div class="line"><a name="l00102"></a><span class="lineno"> 102</span>  v_orb.set(i) = www(i).val_grid_point(0, 0, 0, 0) ;</div><div class="line"><a name="l00103"></a><span class="lineno"> 103</span>  }</div><div class="line"><a name="l00104"></a><span class="lineno"> 104</span> </div><div class="line"><a name="l00105"></a><span class="lineno"> 105</span>  v_orb.<a class="code" href="classLorene_1_1Tensor.html#a230bf6b9bfad003434a5384e5c559217">set_triad</a>( *(www.<a class="code" href="classLorene_1_1Tensor.html#afc084e117c3f5b816afc88a6ffd01c56">get_triad</a>()) ) ;</div><div class="line"><a name="l00106"></a><span class="lineno"> 106</span>  v_orb.std_spectral_base() ;</div><div class="line"><a name="l00107"></a><span class="lineno"> 107</span> </div><div class="line"><a name="l00108"></a><span class="lineno"> 108</span>  <span class="comment">//-------------------------------------------------</span></div><div class="line"><a name="l00109"></a><span class="lineno"> 109</span>  <span class="comment">// Source and coefficients a,b for poisson_compact (independent from psi0)</span></div><div class="line"><a name="l00110"></a><span class="lineno"> 110</span>  <span class="comment">//-------------------------------------------------</span></div><div class="line"><a name="l00111"></a><span class="lineno"> 111</span> </div><div class="line"><a name="l00112"></a><span class="lineno"> 112</span>  <a class="code" href="classLorene_1_1Scalar.html">Scalar</a> dndh_log = <a class="code" href="classLorene_1_1Star.html#a50a6ced39ae7a7d5602ef5953414be79">eos</a>.<a class="code" href="classLorene_1_1Eos.html#afa6ed337a87a4025e5c0e84ad6b83ed6">der_nbar_ent</a>(<a class="code" href="classLorene_1_1Star.html#abd29e4e86f0f9f65b941948e7baa286f">ent</a>, <a class="code" href="classLorene_1_1Star.html#ae5b312dbec065adce16692d1fe834252">nzet</a>) ;</div><div class="line"><a name="l00113"></a><span class="lineno"> 113</span> </div><div class="line"><a name="l00114"></a><span class="lineno"> 114</span>  <span class="comment">// In order to avoid any division by zero in the computation of zeta_h</span></div><div class="line"><a name="l00115"></a><span class="lineno"> 115</span>  <span class="comment">// the value of dndh_log is set to 1 in the external domains:</span></div><div class="line"><a name="l00116"></a><span class="lineno"> 116</span>  <span class="keywordflow">for</span> (<span class="keywordtype">int</span> l=<a class="code" href="classLorene_1_1Star.html#ae5b312dbec065adce16692d1fe834252">nzet</a>; l<=nzm1; l++) {</div><div class="line"><a name="l00117"></a><span class="lineno"> 117</span>  dndh_log.<a class="code" href="classLorene_1_1Scalar.html#a44e1e1f7fc197b7229932aa80c98c479">set_domain</a>(l) = 1. ;</div><div class="line"><a name="l00118"></a><span class="lineno"> 118</span>  }</div><div class="line"><a name="l00119"></a><span class="lineno"> 119</span> </div><div class="line"><a name="l00120"></a><span class="lineno"> 120</span>  <a class="code" href="classLorene_1_1Scalar.html">Scalar</a> zeta_h( <a class="code" href="classLorene_1_1Star.html#abd29e4e86f0f9f65b941948e7baa286f">ent</a> / dndh_log ) ;</div><div class="line"><a name="l00121"></a><span class="lineno"> 121</span>  zeta_h.<a class="code" href="classLorene_1_1Scalar.html#a4eae068572578d6fdd0ab74b90c17b76">std_spectral_base</a>() ;</div><div class="line"><a name="l00122"></a><span class="lineno"> 122</span> </div><div class="line"><a name="l00123"></a><span class="lineno"> 123</span>  <span class="keywordtype">double</span> erreur ;</div><div class="line"><a name="l00124"></a><span class="lineno"> 124</span> </div><div class="line"><a name="l00125"></a><span class="lineno"> 125</span>  <a class="code" href="classLorene_1_1Scalar.html">Scalar</a> source(<a class="code" href="classLorene_1_1Star.html#ac8b026b43b98b1f8f25cc57d96cb3243">mp</a>) ;</div><div class="line"><a name="l00126"></a><span class="lineno"> 126</span>  <a class="code" href="classLorene_1_1Vector.html">Vector</a> bb(<a class="code" href="classLorene_1_1Star.html#ac8b026b43b98b1f8f25cc57d96cb3243">mp</a>, CON, <a class="code" href="classLorene_1_1Star.html#ac8b026b43b98b1f8f25cc57d96cb3243">mp</a>.<a class="code" href="classLorene_1_1Map.html#a6066565316b9953e46b1326b10eff8d4">get_bvect_spher</a>()) ;</div><div class="line"><a name="l00127"></a><span class="lineno"> 127</span>  <a class="code" href="classLorene_1_1Metric__flat.html">Metric_flat</a> flat_spher( <a class="code" href="classLorene_1_1Star.html#ac8b026b43b98b1f8f25cc57d96cb3243">mp</a>.<a class="code" href="classLorene_1_1Map.html#af7ea501bf6b9706de056032bfb807194">flat_met_spher</a>() ) ;</div><div class="line"><a name="l00128"></a><span class="lineno"> 128</span> </div><div class="line"><a name="l00129"></a><span class="lineno"> 129</span>  <span class="keywordflow">if</span> (kerrschild) {</div><div class="line"><a name="l00130"></a><span class="lineno"> 130</span> </div><div class="line"><a name="l00131"></a><span class="lineno"> 131</span>  cout << <span class="stringliteral">"!!!!! WARNING: Not yet available !!!!!"</span> << endl ;</div><div class="line"><a name="l00132"></a><span class="lineno"> 132</span>  abort() ;</div><div class="line"><a name="l00133"></a><span class="lineno"> 133</span> </div><div class="line"><a name="l00134"></a><span class="lineno"> 134</span>  } <span class="comment">// End of Kerr-Schild</span></div><div class="line"><a name="l00135"></a><span class="lineno"> 135</span>  <span class="keywordflow">else</span> { <span class="comment">// Isotropic coordinates with the maximal slicing</span></div><div class="line"><a name="l00136"></a><span class="lineno"> 136</span> </div><div class="line"><a name="l00137"></a><span class="lineno"> 137</span>  <a class="code" href="classLorene_1_1Scalar.html">Scalar</a> lnlappsi(<a class="code" href="classLorene_1_1Star.html#ac8b026b43b98b1f8f25cc57d96cb3243">mp</a>) ;</div><div class="line"><a name="l00138"></a><span class="lineno"> 138</span>  lnlappsi = <a class="code" href="group__cmp__m.html#gafdf39b6dcd55c040de3866406aefe648">log</a>( <a class="code" href="classLorene_1_1Star__bhns.html#a9887ba973118c540900e92b2c29b631b">lapconf_tot</a> * <a class="code" href="classLorene_1_1Star__bhns.html#a6d49d03aab650292b00d2f1b90b65e0e">confo_tot</a> ) ;</div><div class="line"><a name="l00139"></a><span class="lineno"> 139</span>  lnlappsi.<a class="code" href="classLorene_1_1Scalar.html#a4eae068572578d6fdd0ab74b90c17b76">std_spectral_base</a>() ;</div><div class="line"><a name="l00140"></a><span class="lineno"> 140</span> </div><div class="line"><a name="l00141"></a><span class="lineno"> 141</span>  bb = (1. - zeta_h) * <a class="code" href="classLorene_1_1Star.html#abd29e4e86f0f9f65b941948e7baa286f">ent</a>.<a class="code" href="classLorene_1_1Scalar.html#ac7adde08612d8ecd42e4ef4d26224f7c">derive_con</a>(flat_spher)</div><div class="line"><a name="l00142"></a><span class="lineno"> 142</span>  + zeta_h * lnlappsi.<a class="code" href="classLorene_1_1Scalar.html#ac7adde08612d8ecd42e4ef4d26224f7c">derive_con</a>(flat_spher) ;</div><div class="line"><a name="l00143"></a><span class="lineno"> 143</span> </div><div class="line"><a name="l00144"></a><span class="lineno"> 144</span>  <a class="code" href="classLorene_1_1Vector.html">Vector</a> dentmb(<a class="code" href="classLorene_1_1Star.html#ac8b026b43b98b1f8f25cc57d96cb3243">mp</a>, COV, <a class="code" href="classLorene_1_1Star.html#ac8b026b43b98b1f8f25cc57d96cb3243">mp</a>.<a class="code" href="classLorene_1_1Map.html#adfc2b4aa1ef4f22a109bce235e7eb6b7">get_bvect_cart</a>()) ;</div><div class="line"><a name="l00145"></a><span class="lineno"> 145</span>  dentmb.<a class="code" href="classLorene_1_1Tensor.html#ad0c0b10eae3f7f145600569be438ba80">set_etat_qcq</a>() ;</div><div class="line"><a name="l00146"></a><span class="lineno"> 146</span>  <span class="keywordflow">for</span> (<span class="keywordtype">int</span> i=1; i<=3; i++) {</div><div class="line"><a name="l00147"></a><span class="lineno"> 147</span>  dentmb.set(i) = <a class="code" href="classLorene_1_1Star.html#abd29e4e86f0f9f65b941948e7baa286f">ent</a>.<a class="code" href="classLorene_1_1Scalar.html#a1c9daa4fb7f429eb7dc07cfe0c63b759">deriv</a>(i)</div><div class="line"><a name="l00148"></a><span class="lineno"> 148</span>  - (<a class="code" href="classLorene_1_1Star__bhns.html#a71e58dd05a6fa4a759e199ba5e9d568b">d_lapconf_auto</a>(i) + <a class="code" href="classLorene_1_1Star__bhns.html#a9f755b2c65abe93653cb426d6437cec1">d_lapconf_comp</a>(i)) / <a class="code" href="classLorene_1_1Star__bhns.html#a9887ba973118c540900e92b2c29b631b">lapconf_tot</a></div><div class="line"><a name="l00149"></a><span class="lineno"> 149</span>  - (<a class="code" href="classLorene_1_1Star__bhns.html#a7452dd044e729aa703af306e1313769b">d_confo_auto</a>(i) + <a class="code" href="classLorene_1_1Star__bhns.html#a6bca9b5dc3a93a15b7d0afe2c7f1deeb">d_confo_comp</a>(i)) / confo_tot ;</div><div class="line"><a name="l00150"></a><span class="lineno"> 150</span>  }</div><div class="line"><a name="l00151"></a><span class="lineno"> 151</span>  dentmb.<a class="code" href="classLorene_1_1Scalar.html#a4eae068572578d6fdd0ab74b90c17b76">std_spectral_base</a>() ;</div><div class="line"><a name="l00152"></a><span class="lineno"> 152</span> </div><div class="line"><a name="l00153"></a><span class="lineno"> 153</span>  source = <a class="code" href="group__tens__cal.html#ga39112c7788de4bfd90424c7d552bae21">contract</a>(www - v_orb, 0, <a class="code" href="classLorene_1_1Star.html#abd29e4e86f0f9f65b941948e7baa286f">ent</a>.<a class="code" href="classLorene_1_1Scalar.html#a174b781161f0c232c80d8cd1f1ba5615">derive_cov</a>(<a class="code" href="classLorene_1_1Star__bhns.html#a0067ae1d842b4c67441b5463590ef51e">flat</a>), 0)</div><div class="line"><a name="l00154"></a><span class="lineno"> 154</span>  +zeta_h*(<a class="code" href="group__tens__cal.html#ga39112c7788de4bfd90424c7d552bae21">contract</a>(v_orb, 0, dentmb, 0)</div><div class="line"><a name="l00155"></a><span class="lineno"> 155</span>  + <a class="code" href="group__tens__cal.html#ga39112c7788de4bfd90424c7d552bae21">contract</a>(www/<a class="code" href="classLorene_1_1Star.html#a8a0796aa1b273029e32cdfb104f330cb">gam_euler</a>, 0, <a class="code" href="classLorene_1_1Star.html#a8a0796aa1b273029e32cdfb104f330cb">gam_euler</a>.<a class="code" href="classLorene_1_1Scalar.html#a174b781161f0c232c80d8cd1f1ba5615">derive_cov</a>(<a class="code" href="classLorene_1_1Star__bhns.html#a0067ae1d842b4c67441b5463590ef51e">flat</a>), 0)</div><div class="line"><a name="l00156"></a><span class="lineno"> 156</span>  ) ;</div><div class="line"><a name="l00157"></a><span class="lineno"> 157</span> </div><div class="line"><a name="l00158"></a><span class="lineno"> 158</span>  source.<a class="code" href="classLorene_1_1Scalar.html#a7d65bc7f080d1974a8429de810e0eff3">annule</a>(<a class="code" href="classLorene_1_1Star.html#ae5b312dbec065adce16692d1fe834252">nzet</a>, nzm1) ;</div><div class="line"><a name="l00159"></a><span class="lineno"> 159</span> </div><div class="line"><a name="l00160"></a><span class="lineno"> 160</span>  } <span class="comment">// End of Isotropic coordinates</span></div><div class="line"><a name="l00161"></a><span class="lineno"> 161</span> </div><div class="line"><a name="l00162"></a><span class="lineno"> 162</span> </div><div class="line"><a name="l00163"></a><span class="lineno"> 163</span>  <span class="comment">//----------------------------------------------------</span></div><div class="line"><a name="l00164"></a><span class="lineno"> 164</span>  <span class="comment">// Resolution by means of Map_radial::poisson_compact</span></div><div class="line"><a name="l00165"></a><span class="lineno"> 165</span>  <span class="comment">//----------------------------------------------------</span></div><div class="line"><a name="l00166"></a><span class="lineno"> 166</span> </div><div class="line"><a name="l00167"></a><span class="lineno"> 167</span>  <a class="code" href="classLorene_1_1Param.html">Param</a> par ;</div><div class="line"><a name="l00168"></a><span class="lineno"> 168</span>  <span class="keywordtype">int</span> niter ;</div><div class="line"><a name="l00169"></a><span class="lineno"> 169</span>  par.<a class="code" href="classLorene_1_1Param.html#aefd8de56123122b81d919e877a7a5cb7">add_int</a>(mermax) ;</div><div class="line"><a name="l00170"></a><span class="lineno"> 170</span>  par.<a class="code" href="classLorene_1_1Param.html#a049f47c98bc2f7162a94b5a8c9215fde">add_double</a>(precis, 0) ;</div><div class="line"><a name="l00171"></a><span class="lineno"> 171</span>  par.<a class="code" href="classLorene_1_1Param.html#a049f47c98bc2f7162a94b5a8c9215fde">add_double</a>(relax, 1) ;</div><div class="line"><a name="l00172"></a><span class="lineno"> 172</span>  par.<a class="code" href="classLorene_1_1Param.html#a8ea6e8bc0a7166c622764641bbc3744e">add_int_mod</a>(niter) ;</div><div class="line"><a name="l00173"></a><span class="lineno"> 173</span> </div><div class="line"><a name="l00174"></a><span class="lineno"> 174</span>  <span class="keywordflow">if</span> (<a class="code" href="classLorene_1_1Star__bhns.html#abfc578637808e1dd52261e2d8030e9e7">psi0</a>.<a class="code" href="classLorene_1_1Scalar.html#a957e9163ea5a022bef63bdd295cc473f">get_etat</a>() == ETATZERO) {</div><div class="line"><a name="l00175"></a><span class="lineno"> 175</span>  <a class="code" href="classLorene_1_1Star__bhns.html#abfc578637808e1dd52261e2d8030e9e7">psi0</a> = 0. ;</div><div class="line"><a name="l00176"></a><span class="lineno"> 176</span>  }</div><div class="line"><a name="l00177"></a><span class="lineno"> 177</span> </div><div class="line"><a name="l00178"></a><span class="lineno"> 178</span>  <a class="code" href="classLorene_1_1Cmp.html">Cmp</a> source_cmp(source) ;</div><div class="line"><a name="l00179"></a><span class="lineno"> 179</span>  <a class="code" href="classLorene_1_1Cmp.html">Cmp</a> zeta_h_cmp(zeta_h) ;</div><div class="line"><a name="l00180"></a><span class="lineno"> 180</span>  <a class="code" href="classLorene_1_1Cmp.html">Cmp</a> psi0_cmp(<a class="code" href="classLorene_1_1Star__bhns.html#abfc578637808e1dd52261e2d8030e9e7">psi0</a>) ;</div><div class="line"><a name="l00181"></a><span class="lineno"> 181</span> </div><div class="line"><a name="l00182"></a><span class="lineno"> 182</span>  <a class="code" href="classLorene_1_1Tenseur.html">Tenseur</a> bb_cmp(<a class="code" href="classLorene_1_1Star.html#ac8b026b43b98b1f8f25cc57d96cb3243">mp</a>, 1, CON, <a class="code" href="classLorene_1_1Star.html#ac8b026b43b98b1f8f25cc57d96cb3243">mp</a>.<a class="code" href="classLorene_1_1Map.html#a6066565316b9953e46b1326b10eff8d4">get_bvect_spher</a>()) ;</div><div class="line"><a name="l00183"></a><span class="lineno"> 183</span>  bb_cmp.<a class="code" href="classLorene_1_1Tenseur.html#a48a4d68fcd77f7c3e3a8d074b6f0648a">set_etat_qcq</a>() ;</div><div class="line"><a name="l00184"></a><span class="lineno"> 184</span>  <a class="code" href="classLorene_1_1Cmp.html">Cmp</a> bb_cmp1(bb(1)) ;</div><div class="line"><a name="l00185"></a><span class="lineno"> 185</span>  <a class="code" href="classLorene_1_1Cmp.html">Cmp</a> bb_cmp2(bb(2)) ;</div><div class="line"><a name="l00186"></a><span class="lineno"> 186</span>  <a class="code" href="classLorene_1_1Cmp.html">Cmp</a> bb_cmp3(bb(3)) ;</div><div class="line"><a name="l00187"></a><span class="lineno"> 187</span>  bb_cmp.set(0) = bb_cmp1 ;</div><div class="line"><a name="l00188"></a><span class="lineno"> 188</span>  bb_cmp.<a class="code" href="classLorene_1_1Cmp.html#aff912841b755c3dfc7b4ac020d0fa208">set</a>(1) = bb_cmp2 ;</div><div class="line"><a name="l00189"></a><span class="lineno"> 189</span>  bb_cmp.<a class="code" href="classLorene_1_1Cmp.html#aff912841b755c3dfc7b4ac020d0fa208">set</a>(2) = bb_cmp3 ;</div><div class="line"><a name="l00190"></a><span class="lineno"> 190</span> </div><div class="line"><a name="l00191"></a><span class="lineno"> 191</span>  source_cmp.<a class="code" href="classLorene_1_1Cmp.html#a353a6fdc2e07cdc85716ac4f833da9aa">va</a>.<a class="code" href="classLorene_1_1Valeur.html#a611a9887c6903ff084043da54195a760">ylm</a>() ;</div><div class="line"><a name="l00192"></a><span class="lineno"> 192</span> </div><div class="line"><a name="l00193"></a><span class="lineno"> 193</span>  <a class="code" href="classLorene_1_1Star.html#ac8b026b43b98b1f8f25cc57d96cb3243">mp</a>.<a class="code" href="classLorene_1_1Map.html#a12b8a0a54a3b7a30f8cbc55d1c68711a">poisson_compact</a>(source_cmp, zeta_h_cmp, bb_cmp, par, psi0_cmp) ;</div><div class="line"><a name="l00194"></a><span class="lineno"> 194</span> </div><div class="line"><a name="l00195"></a><span class="lineno"> 195</span>  <a class="code" href="classLorene_1_1Star__bhns.html#abfc578637808e1dd52261e2d8030e9e7">psi0</a> = psi0_cmp ;</div><div class="line"><a name="l00196"></a><span class="lineno"> 196</span> </div><div class="line"><a name="l00197"></a><span class="lineno"> 197</span>  <span class="comment">//-----------------------</span></div><div class="line"><a name="l00198"></a><span class="lineno"> 198</span>  <span class="comment">// Check of the solution</span></div><div class="line"><a name="l00199"></a><span class="lineno"> 199</span>  <span class="comment">//-----------------------</span></div><div class="line"><a name="l00200"></a><span class="lineno"> 200</span> </div><div class="line"><a name="l00201"></a><span class="lineno"> 201</span>  <a class="code" href="classLorene_1_1Scalar.html">Scalar</a> bb_dpsi0 = <a class="code" href="group__tens__cal.html#ga39112c7788de4bfd90424c7d552bae21">contract</a>(bb, 0, <a class="code" href="classLorene_1_1Star__bhns.html#abfc578637808e1dd52261e2d8030e9e7">psi0</a>.<a class="code" href="classLorene_1_1Scalar.html#a174b781161f0c232c80d8cd1f1ba5615">derive_cov</a>(flat_spher), 0) ;</div><div class="line"><a name="l00202"></a><span class="lineno"> 202</span> </div><div class="line"><a name="l00203"></a><span class="lineno"> 203</span>  <a class="code" href="classLorene_1_1Scalar.html">Scalar</a> oper = zeta_h * <a class="code" href="classLorene_1_1Star__bhns.html#abfc578637808e1dd52261e2d8030e9e7">psi0</a>.<a class="code" href="classLorene_1_1Scalar.html#a311eda1e80ff4e136043530022ec2f45">laplacian</a>() + bb_dpsi0 ;</div><div class="line"><a name="l00204"></a><span class="lineno"> 204</span> </div><div class="line"><a name="l00205"></a><span class="lineno"> 205</span>  source.<a class="code" href="classLorene_1_1Scalar.html#aa9c1f146171d343ba77204581f3db771">set_spectral_va</a>().<a class="code" href="classLorene_1_1Valeur.html#ad5d93b89e41be4b7af1ab55df0e4a50c">ylm_i</a>() ;</div><div class="line"><a name="l00206"></a><span class="lineno"> 206</span> </div><div class="line"><a name="l00207"></a><span class="lineno"> 207</span>  erreur = <a class="code" href="group__cmp__m.html#ga8796af0ce7ead5666ec0a06ad1829b8a">diffrel</a>(oper, source)(0) ;</div><div class="line"><a name="l00208"></a><span class="lineno"> 208</span> </div><div class="line"><a name="l00209"></a><span class="lineno"> 209</span>  cout << <span class="stringliteral">"Check of the resolution of the continuity equation : "</span></div><div class="line"><a name="l00210"></a><span class="lineno"> 210</span>  << endl ;</div><div class="line"><a name="l00211"></a><span class="lineno"> 211</span>  cout << <span class="stringliteral">" norme(source) : "</span> << <a class="code" href="group__cmp__m.html#ga87b26d4a1417bc5eedd182ec2daa16e8">norme</a>(source)(0)</div><div class="line"><a name="l00212"></a><span class="lineno"> 212</span>  << <span class="stringliteral">" diff oper/source : "</span> << erreur << endl ;</div><div class="line"><a name="l00213"></a><span class="lineno"> 213</span> </div><div class="line"><a name="l00214"></a><span class="lineno"> 214</span>  <span class="comment">//--------------------------</span></div><div class="line"><a name="l00215"></a><span class="lineno"> 215</span>  <span class="comment">// Computation of grad(psi)</span></div><div class="line"><a name="l00216"></a><span class="lineno"> 216</span>  <span class="comment">//--------------------------</span></div><div class="line"><a name="l00217"></a><span class="lineno"> 217</span> </div><div class="line"><a name="l00218"></a><span class="lineno"> 218</span>  <span class="keywordflow">for</span> (<span class="keywordtype">int</span> i=1; i<=3; i++)</div><div class="line"><a name="l00219"></a><span class="lineno"> 219</span>  <a class="code" href="classLorene_1_1Star__bhns.html#a7337c680b6bd8de2367ba78248570f81">d_psi</a>.<a class="code" href="classLorene_1_1Vector.html#ae6ed605a9527b3b5ae0b810fbef713a6">set</a>(i) = (<a class="code" href="classLorene_1_1Star__bhns.html#abfc578637808e1dd52261e2d8030e9e7">psi0</a>.<a class="code" href="classLorene_1_1Scalar.html#a174b781161f0c232c80d8cd1f1ba5615">derive_cov</a>(<a class="code" href="classLorene_1_1Star__bhns.html#a0067ae1d842b4c67441b5463590ef51e">flat</a>))(i) + v_orb(i) ;</div><div class="line"><a name="l00220"></a><span class="lineno"> 220</span> </div><div class="line"><a name="l00221"></a><span class="lineno"> 221</span> </div><div class="line"><a name="l00222"></a><span class="lineno"> 222</span>  <span class="comment">// C^1 continuation of d_psi outside the star</span></div><div class="line"><a name="l00223"></a><span class="lineno"> 223</span>  <span class="comment">// (to ensure a smooth enthalpy field accross the stellar surface)</span></div><div class="line"><a name="l00224"></a><span class="lineno"> 224</span>  <span class="comment">// ----------------------------------------------------------------</span></div><div class="line"><a name="l00225"></a><span class="lineno"> 225</span> </div><div class="line"><a name="l00226"></a><span class="lineno"> 226</span>  <a class="code" href="classLorene_1_1Star__bhns.html#a7337c680b6bd8de2367ba78248570f81">d_psi</a>.<a class="code" href="classLorene_1_1Tensor.html#a814e10dd902c7f7a56b9da0c0caf2b53">annule</a>(<a class="code" href="classLorene_1_1Star.html#ae5b312dbec065adce16692d1fe834252">nzet</a>, nzm1) ;</div><div class="line"><a name="l00227"></a><span class="lineno"> 227</span> </div><div class="line"><a name="l00228"></a><span class="lineno"> 228</span>  <span class="keywordflow">for</span> (<span class="keywordtype">int</span> i=1; i<=3; i++) {</div><div class="line"><a name="l00229"></a><span class="lineno"> 229</span>  <a class="code" href="classLorene_1_1Cmp.html">Cmp</a> d_psi_i( <a class="code" href="classLorene_1_1Star__bhns.html#a7337c680b6bd8de2367ba78248570f81">d_psi</a>(i) ) ;</div><div class="line"><a name="l00230"></a><span class="lineno"> 230</span>  d_psi_i = raccord_c1(d_psi_i, <a class="code" href="classLorene_1_1Star.html#ae5b312dbec065adce16692d1fe834252">nzet</a>) ;</div><div class="line"><a name="l00231"></a><span class="lineno"> 231</span>  <a class="code" href="classLorene_1_1Star__bhns.html#a7337c680b6bd8de2367ba78248570f81">d_psi</a>.<a class="code" href="classLorene_1_1Vector.html#ae6ed605a9527b3b5ae0b810fbef713a6">set</a>(i) = d_psi_i ;</div><div class="line"><a name="l00232"></a><span class="lineno"> 232</span>  }</div><div class="line"><a name="l00233"></a><span class="lineno"> 233</span> </div><div class="line"><a name="l00234"></a><span class="lineno"> 234</span>  <span class="keywordflow">return</span> erreur ;</div><div class="line"><a name="l00235"></a><span class="lineno"> 235</span> </div><div class="line"><a name="l00236"></a><span class="lineno"> 236</span> }</div><div class="line"><a name="l00237"></a><span class="lineno"> 237</span> }</div><div class="ttc" id="classLorene_1_1Star__bhns_html_a6d677eed1bff187093e1fb8c19d11538"><div class="ttname"><a href="classLorene_1_1Star__bhns.html#a6d677eed1bff187093e1fb8c19d11538">Lorene::Star_bhns::psi4</a></div><div class="ttdeci">Scalar psi4</div><div class="ttdoc">Fourth power of the total conformal factor. </div><div class="ttdef"><b>Definition:</b> <a href="star__bhns_8h_source.html#l00176">star_bhns.h:176</a></div></div>
<div class="ttc" id="classLorene_1_1Star__bhns_html_a71e58dd05a6fa4a759e199ba5e9d568b"><div class="ttname"><a href="classLorene_1_1Star__bhns.html#a71e58dd05a6fa4a759e199ba5e9d568b">Lorene::Star_bhns::d_lapconf_auto</a></div><div class="ttdeci">Vector d_lapconf_auto</div><div class="ttdoc">Derivative of the lapconf function generated by the star . </div><div class="ttdef"><b>Definition:</b> <a href="star__bhns_8h_source.html#l00130">star_bhns.h:130</a></div></div>
<div class="ttc" id="classLorene_1_1Star__bhns_html_a9887ba973118c540900e92b2c29b631b"><div class="ttname"><a href="classLorene_1_1Star__bhns.html#a9887ba973118c540900e92b2c29b631b">Lorene::Star_bhns::lapconf_tot</a></div><div class="ttdeci">Scalar lapconf_tot</div><div class="ttdoc">Total lapconf function. </div><div class="ttdef"><b>Definition:</b> <a href="star__bhns_8h_source.html#l00119">star_bhns.h:119</a></div></div>
<div class="ttc" id="group__cmp__m_html_gafdf39b6dcd55c040de3866406aefe648"><div class="ttname"><a href="group__cmp__m.html#gafdf39b6dcd55c040de3866406aefe648">Lorene::log</a></div><div class="ttdeci">Cmp log(const Cmp &)</div><div class="ttdoc">Neperian logarithm. </div><div class="ttdef"><b>Definition:</b> <a href="cmp__math_8C_source.html#l00296">cmp_math.C:296</a></div></div>
<div class="ttc" id="classLorene_1_1Tensor_html_ad0c0b10eae3f7f145600569be438ba80"><div class="ttname"><a href="classLorene_1_1Tensor.html#ad0c0b10eae3f7f145600569be438ba80">Lorene::Tensor::set_etat_qcq</a></div><div class="ttdeci">virtual void set_etat_qcq()</div><div class="ttdoc">Sets the logical state of all components to ETATQCQ (ordinary state). </div><div class="ttdef"><b>Definition:</b> <a href="tensor_8C_source.html#l00481">tensor.C:481</a></div></div>
<div class="ttc" id="classLorene_1_1Cmp_html"><div class="ttname"><a href="classLorene_1_1Cmp.html">Lorene::Cmp</a></div><div class="ttdoc">Component of a tensorial field *** DEPRECATED : use class Scalar instead ***. </div><div class="ttdef"><b>Definition:</b> <a href="cmp_8h_source.html#l00446">cmp.h:446</a></div></div>
<div class="ttc" id="group__cmp__m_html_ga66a0d18b41ceae25049869b5ccdca44d"><div class="ttname"><a href="group__cmp__m.html#ga66a0d18b41ceae25049869b5ccdca44d">Lorene::exp</a></div><div class="ttdeci">Cmp exp(const Cmp &)</div><div class="ttdoc">Exponential. </div><div class="ttdef"><b>Definition:</b> <a href="cmp__math_8C_source.html#l00270">cmp_math.C:270</a></div></div>
<div class="ttc" id="classLorene_1_1Param_html_aefd8de56123122b81d919e877a7a5cb7"><div class="ttname"><a href="classLorene_1_1Param.html#aefd8de56123122b81d919e877a7a5cb7">Lorene::Param::add_int</a></div><div class="ttdeci">void add_int(const int &n, int position=0)</div><div class="ttdoc">Adds the address of a new int to the list. </div><div class="ttdef"><b>Definition:</b> <a href="param_8C_source.html#l00246">param.C:246</a></div></div>
<div class="ttc" id="classLorene_1_1Star_html_ac8b026b43b98b1f8f25cc57d96cb3243"><div class="ttname"><a href="classLorene_1_1Star.html#ac8b026b43b98b1f8f25cc57d96cb3243">Lorene::Star::mp</a></div><div class="ttdeci">Map & mp</div><div class="ttdoc">Mapping associated with the star. </div><div class="ttdef"><b>Definition:</b> <a href="star_8h_source.html#l00180">star.h:180</a></div></div>
<div class="ttc" id="classLorene_1_1Valeur_html_ad5d93b89e41be4b7af1ab55df0e4a50c"><div class="ttname"><a href="classLorene_1_1Valeur.html#ad5d93b89e41be4b7af1ab55df0e4a50c">Lorene::Valeur::ylm_i</a></div><div class="ttdeci">void ylm_i()</div><div class="ttdoc">Inverse of ylm() </div><div class="ttdef"><b>Definition:</b> <a href="valeur__ylm__i_8C_source.html#l00131">valeur_ylm_i.C:131</a></div></div>
<div class="ttc" id="classLorene_1_1Scalar_html_a7d65bc7f080d1974a8429de810e0eff3"><div class="ttname"><a href="classLorene_1_1Scalar.html#a7d65bc7f080d1974a8429de810e0eff3">Lorene::Scalar::annule</a></div><div class="ttdeci">virtual void annule(int l_min, int l_max)</div><div class="ttdoc">Sets the Scalar to zero in several domains. </div><div class="ttdef"><b>Definition:</b> <a href="scalar_8C_source.html#l00391">scalar.C:391</a></div></div>
<div class="ttc" id="classLorene_1_1Star__bhns_html_a285e3e765784ecabd7b5043df2ce5c9c"><div class="ttname"><a href="classLorene_1_1Star__bhns.html#a285e3e765784ecabd7b5043df2ce5c9c">Lorene::Star_bhns::bsn</a></div><div class="ttdeci">Vector bsn</div><div class="ttdoc">3-vector shift, divided by N , of the rotating coordinates, . </div><div class="ttdef"><b>Definition:</b> <a href="star__bhns_8h_source.html#l00099">star_bhns.h:99</a></div></div>
<div class="ttc" id="classLorene_1_1Map_html_a6066565316b9953e46b1326b10eff8d4"><div class="ttname"><a href="classLorene_1_1Map.html#a6066565316b9953e46b1326b10eff8d4">Lorene::Map::get_bvect_spher</a></div><div class="ttdeci">const Base_vect_spher & get_bvect_spher() const</div><div class="ttdoc">Returns the orthonormal vectorial basis associated with the coordinates of the mapping. </div><div class="ttdef"><b>Definition:</b> <a href="map_8h_source.html#l00783">map.h:783</a></div></div>
<div class="ttc" id="namespaceLorene_html"><div class="ttname"><a href="namespaceLorene.html">Lorene</a></div><div class="ttdoc">Lorene prototypes. </div><div class="ttdef"><b>Definition:</b> <a href="app__hor_8h_source.html#l00064">app_hor.h:64</a></div></div>
<div class="ttc" id="classLorene_1_1Star_html_a50a6ced39ae7a7d5602ef5953414be79"><div class="ttname"><a href="classLorene_1_1Star.html#a50a6ced39ae7a7d5602ef5953414be79">Lorene::Star::eos</a></div><div class="ttdeci">const Eos & eos</div><div class="ttdoc">Equation of state of the stellar matter. </div><div class="ttdef"><b>Definition:</b> <a href="star_8h_source.html#l00185">star.h:185</a></div></div>
<div class="ttc" id="namespaceUnites_html"><div class="ttname"><a href="namespaceUnites.html">Unites</a></div><div class="ttdoc">Standard units of space, time and mass. </div></div>
<div class="ttc" id="classLorene_1_1Valeur_html_a611a9887c6903ff084043da54195a760"><div class="ttname"><a href="classLorene_1_1Valeur.html#a611a9887c6903ff084043da54195a760">Lorene::Valeur::ylm</a></div><div class="ttdeci">void ylm()</div><div class="ttdoc">Computes the coefficients of *this. </div><div class="ttdef"><b>Definition:</b> <a href="valeur__ylm_8C_source.html#l00138">valeur_ylm.C:138</a></div></div>
<div class="ttc" id="classLorene_1_1Map_html_a4df8c6a943a85ef314ea636240ef31a7"><div class="ttname"><a href="classLorene_1_1Map.html#a4df8c6a943a85ef314ea636240ef31a7">Lorene::Map::get_mg</a></div><div class="ttdeci">const Mg3d * get_mg() const</div><div class="ttdoc">Gives the Mg3d on which the mapping is defined. </div><div class="ttdef"><b>Definition:</b> <a href="map_8h_source.html#l00765">map.h:765</a></div></div>
<div class="ttc" id="classLorene_1_1Scalar_html_a311eda1e80ff4e136043530022ec2f45"><div class="ttname"><a href="classLorene_1_1Scalar.html#a311eda1e80ff4e136043530022ec2f45">Lorene::Scalar::laplacian</a></div><div class="ttdeci">const Scalar & laplacian(int ced_mult_r=4) const</div><div class="ttdoc">Returns the Laplacian of *this. </div><div class="ttdef"><b>Definition:</b> <a href="scalar__deriv_8C_source.html#l00436">scalar_deriv.C:436</a></div></div>
<div class="ttc" id="classLorene_1_1Metric__flat_html"><div class="ttname"><a href="classLorene_1_1Metric__flat.html">Lorene::Metric_flat</a></div><div class="ttdoc">Flat metric for tensor calculation. </div><div class="ttdef"><b>Definition:</b> <a href="metric_8h_source.html#l00261">metric.h:261</a></div></div>
<div class="ttc" id="classLorene_1_1Star__bhns_html_a9f755b2c65abe93653cb426d6437cec1"><div class="ttname"><a href="classLorene_1_1Star__bhns.html#a9f755b2c65abe93653cb426d6437cec1">Lorene::Star_bhns::d_lapconf_comp</a></div><div class="ttdeci">Vector d_lapconf_comp</div><div class="ttdoc">Derivative of the lapconf function generated by the companion black hole. </div><div class="ttdef"><b>Definition:</b> <a href="star__bhns_8h_source.html#l00135">star_bhns.h:135</a></div></div>
<div class="ttc" id="classLorene_1_1Scalar_html"><div class="ttname"><a href="classLorene_1_1Scalar.html">Lorene::Scalar</a></div><div class="ttdoc">Tensor field of valence 0 (or component of a tensorial field). </div><div class="ttdef"><b>Definition:</b> <a href="scalar_8h_source.html#l00387">scalar.h:387</a></div></div>
<div class="ttc" id="classLorene_1_1Map_html_a12b8a0a54a3b7a30f8cbc55d1c68711a"><div class="ttname"><a href="classLorene_1_1Map.html#a12b8a0a54a3b7a30f8cbc55d1c68711a">Lorene::Map::poisson_compact</a></div><div class="ttdeci">virtual void poisson_compact(const Cmp &source, const Cmp &aa, const Tenseur &bb, const Param &par, Cmp &psi) const =0</div><div class="ttdoc">Resolution of the elliptic equation in the case where the stellar interior is covered by a single do...</div></div>
<div class="ttc" id="classLorene_1_1Scalar_html_a4eae068572578d6fdd0ab74b90c17b76"><div class="ttname"><a href="classLorene_1_1Scalar.html#a4eae068572578d6fdd0ab74b90c17b76">Lorene::Scalar::std_spectral_base</a></div><div class="ttdeci">virtual void std_spectral_base()</div><div class="ttdoc">Sets the spectral bases of the Valeur va to the standard ones for a scalar field. ...</div><div class="ttdef"><b>Definition:</b> <a href="scalar_8C_source.html#l00784">scalar.C:784</a></div></div>
<div class="ttc" id="classLorene_1_1Vector_html_ad916ce65d0f862b7ea905af6e0ab9785"><div class="ttname"><a href="classLorene_1_1Vector.html#ad916ce65d0f862b7ea905af6e0ab9785">Lorene::Vector::change_triad</a></div><div class="ttdeci">virtual void change_triad(const Base_vect &)</div><div class="ttdoc">Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly. </div><div class="ttdef"><b>Definition:</b> <a href="vector__change__triad_8C_source.html#l00075">vector_change_triad.C:75</a></div></div>
<div class="ttc" id="classLorene_1_1Scalar_html_a957e9163ea5a022bef63bdd295cc473f"><div class="ttname"><a href="classLorene_1_1Scalar.html#a957e9163ea5a022bef63bdd295cc473f">Lorene::Scalar::get_etat</a></div><div class="ttdeci">int get_etat() const</div><div class="ttdoc">Returns the logical state ETATNONDEF (undefined), ETATZERO (null) or ETATQCQ (ordinary). </div><div class="ttdef"><b>Definition:</b> <a href="scalar_8h_source.html#l00554">scalar.h:554</a></div></div>
<div class="ttc" id="classLorene_1_1Star_html_abd29e4e86f0f9f65b941948e7baa286f"><div class="ttname"><a href="classLorene_1_1Star.html#abd29e4e86f0f9f65b941948e7baa286f">Lorene::Star::ent</a></div><div class="ttdeci">Scalar ent</div><div class="ttdoc">Log-enthalpy. </div><div class="ttdef"><b>Definition:</b> <a href="star_8h_source.html#l00190">star.h:190</a></div></div>
<div class="ttc" id="classLorene_1_1Scalar_html_ac7adde08612d8ecd42e4ef4d26224f7c"><div class="ttname"><a href="classLorene_1_1Scalar.html#ac7adde08612d8ecd42e4ef4d26224f7c">Lorene::Scalar::derive_con</a></div><div class="ttdeci">const Vector & derive_con(const Metric &gam) const</div><div class="ttdoc">Returns the "contravariant" derivative of *this with respect to some metric , by raising the index of...</div><div class="ttdef"><b>Definition:</b> <a href="scalar__deriv_8C_source.html#l00402">scalar_deriv.C:402</a></div></div>
<div class="ttc" id="classLorene_1_1Vector_html"><div class="ttname"><a href="classLorene_1_1Vector.html">Lorene::Vector</a></div><div class="ttdoc">Tensor field of valence 1. </div><div class="ttdef"><b>Definition:</b> <a href="vector_8h_source.html#l00188">vector.h:188</a></div></div>
<div class="ttc" id="classLorene_1_1Scalar_html_a44e1e1f7fc197b7229932aa80c98c479"><div class="ttname"><a href="classLorene_1_1Scalar.html#a44e1e1f7fc197b7229932aa80c98c479">Lorene::Scalar::set_domain</a></div><div class="ttdeci">Tbl & set_domain(int l)</div><div class="ttdoc">Read/write of the value in a given domain. </div><div class="ttdef"><b>Definition:</b> <a href="scalar_8h_source.html#l00615">scalar.h:615</a></div></div>
<div class="ttc" id="group__cmp__m_html_ga8796af0ce7ead5666ec0a06ad1829b8a"><div class="ttname"><a href="group__cmp__m.html#ga8796af0ce7ead5666ec0a06ad1829b8a">Lorene::diffrel</a></div><div class="ttdeci">Tbl diffrel(const Cmp &a, const Cmp &b)</div><div class="ttdoc">Relative difference between two Cmp (norme version). </div><div class="ttdef"><b>Definition:</b> <a href="cmp__math_8C_source.html#l00504">cmp_math.C:504</a></div></div>
<div class="ttc" id="classLorene_1_1Star__bhns_html_a7337c680b6bd8de2367ba78248570f81"><div class="ttname"><a href="classLorene_1_1Star__bhns.html#a7337c680b6bd8de2367ba78248570f81">Lorene::Star_bhns::d_psi</a></div><div class="ttdeci">Vector d_psi</div><div class="ttdoc">Gradient of (in the irrotational case) (Spherical components with respect to the mapping of the star...</div><div class="ttdef"><b>Definition:</b> <a href="star__bhns_8h_source.html#l00082">star_bhns.h:82</a></div></div>
<div class="ttc" id="classLorene_1_1Star_html_ae5b312dbec065adce16692d1fe834252"><div class="ttname"><a href="classLorene_1_1Star.html#ae5b312dbec065adce16692d1fe834252">Lorene::Star::nzet</a></div><div class="ttdeci">int nzet</div><div class="ttdoc">Number of domains of *mp occupied by the star. </div><div class="ttdef"><b>Definition:</b> <a href="star_8h_source.html#l00183">star.h:183</a></div></div>
<div class="ttc" id="classLorene_1_1Star__bhns_html_a0067ae1d842b4c67441b5463590ef51e"><div class="ttname"><a href="classLorene_1_1Star__bhns.html#a0067ae1d842b4c67441b5463590ef51e">Lorene::Star_bhns::flat</a></div><div class="ttdeci">Metric_flat flat</div><div class="ttdoc">Flat metric defined on the mapping (Spherical components with respect to the mapping of the star )...</div><div class="ttdef"><b>Definition:</b> <a href="star__bhns_8h_source.html#l00192">star_bhns.h:192</a></div></div>
<div class="ttc" id="classLorene_1_1Star__bhns_html_a6d49d03aab650292b00d2f1b90b65e0e"><div class="ttname"><a href="classLorene_1_1Star__bhns.html#a6d49d03aab650292b00d2f1b90b65e0e">Lorene::Star_bhns::confo_tot</a></div><div class="ttdeci">Scalar confo_tot</div><div class="ttdoc">Total conformal factor. </div><div class="ttdef"><b>Definition:</b> <a href="star__bhns_8h_source.html#l00163">star_bhns.h:163</a></div></div>
<div class="ttc" id="classLorene_1_1Star_html_a8a0796aa1b273029e32cdfb104f330cb"><div class="ttname"><a href="classLorene_1_1Star.html#a8a0796aa1b273029e32cdfb104f330cb">Lorene::Star::gam_euler</a></div><div class="ttdeci">Scalar gam_euler</div><div class="ttdoc">Lorentz factor between the fluid and Eulerian observers. </div><div class="ttdef"><b>Definition:</b> <a href="star_8h_source.html#l00204">star.h:204</a></div></div>
<div class="ttc" id="group__cmp__m_html_ga87b26d4a1417bc5eedd182ec2daa16e8"><div class="ttname"><a href="group__cmp__m.html#ga87b26d4a1417bc5eedd182ec2daa16e8">Lorene::norme</a></div><div class="ttdeci">Tbl norme(const Cmp &)</div><div class="ttdoc">Sums of the absolute values of all the values of the Cmp in each domain. </div><div class="ttdef"><b>Definition:</b> <a href="cmp__math_8C_source.html#l00481">cmp_math.C:481</a></div></div>
<div class="ttc" id="classLorene_1_1Tensor_html_afc084e117c3f5b816afc88a6ffd01c56"><div class="ttname"><a href="classLorene_1_1Tensor.html#afc084e117c3f5b816afc88a6ffd01c56">Lorene::Tensor::get_triad</a></div><div class="ttdeci">const Base_vect * get_triad() const</div><div class="ttdoc">Returns the vectorial basis (triad) on which the components are defined. </div><div class="ttdef"><b>Definition:</b> <a href="tensor_8h_source.html#l00866">tensor.h:866</a></div></div>
<div class="ttc" id="classLorene_1_1Star__bhns_html_a6bca9b5dc3a93a15b7d0afe2c7f1deeb"><div class="ttname"><a href="classLorene_1_1Star__bhns.html#a6bca9b5dc3a93a15b7d0afe2c7f1deeb">Lorene::Star_bhns::d_confo_comp</a></div><div class="ttdeci">Vector d_confo_comp</div><div class="ttdoc">Derivative of the conformal factor generated by the companion black hole. </div><div class="ttdef"><b>Definition:</b> <a href="star__bhns_8h_source.html#l00173">star_bhns.h:173</a></div></div>
<div class="ttc" id="classLorene_1_1Param_html"><div class="ttname"><a href="classLorene_1_1Param.html">Lorene::Param</a></div><div class="ttdoc">Parameter storage. </div><div class="ttdef"><b>Definition:</b> <a href="param_8h_source.html#l00125">param.h:125</a></div></div>
<div class="ttc" id="classLorene_1_1Mg3d_html_a99fb0e7a3f1f0f3e405624855200055e"><div class="ttname"><a href="classLorene_1_1Mg3d.html#a99fb0e7a3f1f0f3e405624855200055e">Lorene::Mg3d::get_nzone</a></div><div class="ttdeci">int get_nzone() const</div><div class="ttdoc">Returns the number of domains. </div><div class="ttdef"><b>Definition:</b> <a href="grilles_8h_source.html#l00448">grilles.h:448</a></div></div>
<div class="ttc" id="classLorene_1_1Scalar_html_a1c9daa4fb7f429eb7dc07cfe0c63b759"><div class="ttname"><a href="classLorene_1_1Scalar.html#a1c9daa4fb7f429eb7dc07cfe0c63b759">Lorene::Scalar::deriv</a></div><div class="ttdeci">const Scalar & deriv(int i) const</div><div class="ttdoc">Returns of *this , where . </div><div class="ttdef"><b>Definition:</b> <a href="scalar__deriv_8C_source.html#l00359">scalar_deriv.C:359</a></div></div>
<div class="ttc" id="group__tens__cal_html_ga39112c7788de4bfd90424c7d552bae21"><div class="ttname"><a href="group__tens__cal.html#ga39112c7788de4bfd90424c7d552bae21">Lorene::contract</a></div><div class="ttdeci">Tenseur contract(const Tenseur &, int id1, int id2)</div><div class="ttdoc">Self contraction of two indices of a Tenseur . </div><div class="ttdef"><b>Definition:</b> <a href="tenseur__operateur_8C_source.html#l00279">tenseur_operateur.C:279</a></div></div>
<div class="ttc" id="classLorene_1_1Cmp_html_aff912841b755c3dfc7b4ac020d0fa208"><div class="ttname"><a href="classLorene_1_1Cmp.html#aff912841b755c3dfc7b4ac020d0fa208">Lorene::Cmp::set</a></div><div class="ttdeci">Tbl & set(int l)</div><div class="ttdoc">Read/write of the value in a given domain. </div><div class="ttdef"><b>Definition:</b> <a href="cmp_8h_source.html#l00724">cmp.h:724</a></div></div>
<div class="ttc" id="classLorene_1_1Map_html_adfc2b4aa1ef4f22a109bce235e7eb6b7"><div class="ttname"><a href="classLorene_1_1Map.html#adfc2b4aa1ef4f22a109bce235e7eb6b7">Lorene::Map::get_bvect_cart</a></div><div class="ttdeci">const Base_vect_cart & get_bvect_cart() const</div><div class="ttdoc">Returns the Cartesian basis associated with the coordinates (x,y,z) of the mapping, i.e. </div><div class="ttdef"><b>Definition:</b> <a href="map_8h_source.html#l00791">map.h:791</a></div></div>
<div class="ttc" id="classLorene_1_1Star__bhns_html_a7452dd044e729aa703af306e1313769b"><div class="ttname"><a href="classLorene_1_1Star__bhns.html#a7452dd044e729aa703af306e1313769b">Lorene::Star_bhns::d_confo_auto</a></div><div class="ttdeci">Vector d_confo_auto</div><div class="ttdoc">Derivative of the conformal factor generated by the star . </div><div class="ttdef"><b>Definition:</b> <a href="star__bhns_8h_source.html#l00168">star_bhns.h:168</a></div></div>
<div class="ttc" id="classLorene_1_1Param_html_a049f47c98bc2f7162a94b5a8c9215fde"><div class="ttname"><a href="classLorene_1_1Param.html#a049f47c98bc2f7162a94b5a8c9215fde">Lorene::Param::add_double</a></div><div class="ttdeci">void add_double(const double &x, int position=0)</div><div class="ttdoc">Adds the the address of a new double to the list. </div><div class="ttdef"><b>Definition:</b> <a href="param_8C_source.html#l00315">param.C:315</a></div></div>
<div class="ttc" id="classLorene_1_1Tensor_html_a230bf6b9bfad003434a5384e5c559217"><div class="ttname"><a href="classLorene_1_1Tensor.html#a230bf6b9bfad003434a5384e5c559217">Lorene::Tensor::set_triad</a></div><div class="ttdeci">void set_triad(const Base_vect &new_triad)</div><div class="ttdoc">Assigns a new vectorial basis (triad) of decomposition. </div><div class="ttdef"><b>Definition:</b> <a href="tensor_8C_source.html#l00519">tensor.C:519</a></div></div>
<div class="ttc" id="classLorene_1_1Scalar_html_aa9c1f146171d343ba77204581f3db771"><div class="ttname"><a href="classLorene_1_1Scalar.html#aa9c1f146171d343ba77204581f3db771">Lorene::Scalar::set_spectral_va</a></div><div class="ttdeci">Valeur & set_spectral_va()</div><div class="ttdoc">Returns va (read/write version) </div><div class="ttdef"><b>Definition:</b> <a href="scalar_8h_source.html#l00604">scalar.h:604</a></div></div>
<div class="ttc" id="classLorene_1_1Vector_html_ae6ed605a9527b3b5ae0b810fbef713a6"><div class="ttname"><a href="classLorene_1_1Vector.html#ae6ed605a9527b3b5ae0b810fbef713a6">Lorene::Vector::set</a></div><div class="ttdeci">Scalar & set(int)</div><div class="ttdoc">Read/write access to a component. </div><div class="ttdef"><b>Definition:</b> <a href="vector_8C_source.html#l00296">vector.C:296</a></div></div>
<div class="ttc" id="classLorene_1_1Tenseur_html_a48a4d68fcd77f7c3e3a8d074b6f0648a"><div class="ttname"><a href="classLorene_1_1Tenseur.html#a48a4d68fcd77f7c3e3a8d074b6f0648a">Lorene::Tenseur::set_etat_qcq</a></div><div class="ttdeci">void set_etat_qcq()</div><div class="ttdoc">Sets the logical state to ETATQCQ (ordinary state). </div><div class="ttdef"><b>Definition:</b> <a href="tenseur_8C_source.html#l00636">tenseur.C:636</a></div></div>
<div class="ttc" id="classLorene_1_1Eos_html_afa6ed337a87a4025e5c0e84ad6b83ed6"><div class="ttname"><a href="classLorene_1_1Eos.html#afa6ed337a87a4025e5c0e84ad6b83ed6">Lorene::Eos::der_nbar_ent</a></div><div class="ttdeci">Cmp der_nbar_ent(const Cmp &ent, int nzet, int l_min=0, const Param *par=0x0) const</div><div class="ttdoc">Computes the logarithmic derivative from the log-enthalpy and extra parameters. </div><div class="ttdef"><b>Definition:</b> <a href="eos_8C_source.html#l00407">eos.C:407</a></div></div>
<div class="ttc" id="classLorene_1_1Scalar_html_a174b781161f0c232c80d8cd1f1ba5615"><div class="ttname"><a href="classLorene_1_1Scalar.html#a174b781161f0c232c80d8cd1f1ba5615">Lorene::Scalar::derive_cov</a></div><div class="ttdeci">const Vector & derive_cov(const Metric &gam) const</div><div class="ttdoc">Returns the gradient (1-form = covariant vector) of *this. </div><div class="ttdef"><b>Definition:</b> <a href="scalar__deriv_8C_source.html#l00390">scalar_deriv.C:390</a></div></div>
<div class="ttc" id="classLorene_1_1Tensor_html_a814e10dd902c7f7a56b9da0c0caf2b53"><div class="ttname"><a href="classLorene_1_1Tensor.html#a814e10dd902c7f7a56b9da0c0caf2b53">Lorene::Tensor::annule</a></div><div class="ttdeci">virtual void annule(int l_min, int l_max)</div><div class="ttdoc">Sets the Tensor to zero in several domains. </div><div class="ttdef"><b>Definition:</b> <a href="tensor_8C_source.html#l00671">tensor.C:671</a></div></div>
<div class="ttc" id="classLorene_1_1Star__bhns_html_abfc578637808e1dd52261e2d8030e9e7"><div class="ttname"><a href="classLorene_1_1Star__bhns.html#abfc578637808e1dd52261e2d8030e9e7">Lorene::Star_bhns::psi0</a></div><div class="ttdeci">Scalar psi0</div><div class="ttdoc">Scalar potential of the non-translational part of fluid 4-velocity (in the irrotational case) ...</div><div class="ttdef"><b>Definition:</b> <a href="star__bhns_8h_source.html#l00077">star_bhns.h:77</a></div></div>
<div class="ttc" id="classLorene_1_1Cmp_html_a353a6fdc2e07cdc85716ac4f833da9aa"><div class="ttname"><a href="classLorene_1_1Cmp.html#a353a6fdc2e07cdc85716ac4f833da9aa">Lorene::Cmp::va</a></div><div class="ttdeci">Valeur va</div><div class="ttdoc">The numerical value of the Cmp. </div><div class="ttdef"><b>Definition:</b> <a href="cmp_8h_source.html#l00464">cmp.h:464</a></div></div>
<div class="ttc" id="classLorene_1_1Map_html_af7ea501bf6b9706de056032bfb807194"><div class="ttname"><a href="classLorene_1_1Map.html#af7ea501bf6b9706de056032bfb807194">Lorene::Map::flat_met_spher</a></div><div class="ttdeci">const Metric_flat & flat_met_spher() const</div><div class="ttdoc">Returns the flat metric associated with the spherical coordinates and with components expressed in th...</div><div class="ttdef"><b>Definition:</b> <a href="map_8C_source.html#l00321">map.C:321</a></div></div>
<div class="ttc" id="classLorene_1_1Tenseur_html"><div class="ttname"><a href="classLorene_1_1Tenseur.html">Lorene::Tenseur</a></div><div class="ttdoc">Tensor handling *** DEPRECATED : use class Tensor instead ***. </div><div class="ttdef"><b>Definition:</b> <a href="tenseur_8h_source.html#l00301">tenseur.h:301</a></div></div>
<div class="ttc" id="classLorene_1_1Param_html_a8ea6e8bc0a7166c622764641bbc3744e"><div class="ttname"><a href="classLorene_1_1Param.html#a8ea6e8bc0a7166c622764641bbc3744e">Lorene::Param::add_int_mod</a></div><div class="ttdeci">void add_int_mod(int &n, int position=0)</div><div class="ttdoc">Adds the address of a new modifiable int to the list. </div><div class="ttdef"><b>Definition:</b> <a href="param_8C_source.html#l00385">param.C:385</a></div></div>
<div class="ttc" id="classLorene_1_1Star__bhns_html_a26f0c79991a8a68ed50e4425271b762d"><div class="ttname"><a href="classLorene_1_1Star__bhns.html#a26f0c79991a8a68ed50e4425271b762d">Lorene::Star_bhns::velo_pot_bhns</a></div><div class="ttdeci">double velo_pot_bhns(const double &mass_bh, const double &sepa, bool kerrschild, int mermax, double precis, double relax)</div><div class="ttdoc">Computes the non-translational part of the velocity scalar potential by solving the continuity equat...</div><div class="ttdef"><b>Definition:</b> <a href="star__bhns__vel__pot_8C_source.html#l00067">star_bhns_vel_pot.C:67</a></div></div>
</div><!-- fragment --></div><!-- contents -->
</div><!-- doc-content -->
<!-- start footer part -->
<div id="nav-path" class="navpath"><!-- id is needed for treeview function! -->
<ul>
<li class="navelem"><a class="el" href="dir_2916f9fbcacdf2e5f0583aa159c707b6.html">C++</a></li><li class="navelem"><a class="el" href="dir_be81ae66a82603ee109e545ee9e7219c.html">Source</a></li><li class="navelem"><a class="el" href="dir_eb75cb0379775c1decf9799397229ad0.html">Star_bhns</a></li><li class="navelem"><b>star_bhns_vel_pot.C</b></li>
<li class="footer">Generated by
<a href="http://www.doxygen.org/index.html">
<img class="footer" src="doxygen.png" alt="doxygen"/></a> 1.8.12 </li>
</ul>
</div>
</body>
</html>
|