PLUTO
entropy_solve.c
Go to the documentation of this file.
1 /* ///////////////////////////////////////////////////////////////////// */
2 /*!
3  \file
4  \brief Inversion scheme for RMHD using entropy.
5 
6  Convert the conservative variables u=[D, m, B, sigma_c]
7  (where sigma_c = D*sigma is the conserved entropy) to
8  primitive variable v=[rho,v,B,p] using a Newton-Raphson/Bisection scheme.
9 
10  \authors C. Zanni \n
11  A. Mignone
12 
13  \date June 25, 2015
14 */
15 /* ///////////////////////////////////////////////////////////////////// */
16 #include "pluto.h"
17 
18 #define MAX_ITER 50
19 /* ********************************************************************* */
21 /*!
22  *
23  *********************************************************************** */
24 {
25  int k, done;
26  double v2, v2max, v2min, acc = 1.e-13;
27  double h, W, W2, W3, S2_W2;
28  double scrh, x, temp;
29  double lor, lor2, rho, sigma, th;
30  double dv2_dW, dW_dlor;
31  double fv2, dfv2, dv2, dv2old;
32  double m2, B2, S, S2, D;
33 
34  D = par->D;
35  m2 = par->m2;
36  B2 = par->B2;
37  S = par->S;
38  S2 = S*S;
39 
40  sigma = par->sigma_c/D;
41  v2max = 1.0 - 1.e-8;
42  v2min = 0.0;
43 
44  done = 0;
45 
46  v2 = 0.5*(v2max+v2min);
47  dv2old = v2max-v2min;
48  dv2 = dv2old;
49 
50  lor2 = 1.0/(1.0 - v2);
51  lor = sqrt(lor2);
52  rho = D/lor;
53 
54  #if EOS == IDEAL
55  x = pow(rho,g_gamma - 1.0);
56  h = 1.0 + g_gamma/(g_gamma - 1.0)*sigma*x;
57  W = D*h*lor;
58  dW_dlor = (2 - g_gamma)*W/lor + (g_gamma - 1.0)*D;
59  #elif EOS == TAUB
60  x = pow(rho,2.0/3.0);
61  th = sigma*x/sqrt(1.0 + 3.0*sigma*x);
62  h = 2.5*th + sqrt(2.25*th*th + 1.0);
63  W = D*h*lor;
64 
65  scrh = -2.0*x/(3.0*lor)*sigma*(2.0*sigma*x - 3.0*th*th);
66  scrh /= 2.0*th*(1.0 + 3.0*sigma*x);
67  dW_dlor = W/lor + D*lor*(5.0*h - 8.0*th)/(2.0*h - 5.0*th)*scrh;
68  #endif
69 
70  W2 = W*W;
71  W3 = W2*W;
72  S2_W2 = S2/W2;
73  scrh = W + B2;
74 
75  fv2 = v2 - (m2 + S2_W2*(scrh + W))/(scrh*scrh);
76 
77  dv2_dW = S2*(3.0*W*scrh + B2*B2) + m2*W3;
78  dv2_dW *= -2.0/(W3*scrh*scrh*scrh);
79 
80  dfv2 = 1.0 - 0.5*dv2_dW*dW_dlor*lor2*lor;
81 
82  for (k = 1; k < MAX_ITER; k++) {
83  if ((((v2-v2max)*dfv2-fv2)*((v2-v2min)*dfv2-fv2) > 0.0)
84  || (fabs(2.*fv2) > fabs(dv2old*dfv2))) {
85  dv2old = dv2;
86  dv2 = 0.5*(v2max-v2min);
87  v2 = v2min+dv2;
88  if (v2min == v2) done = 1;
89  } else {
90  dv2old = dv2;
91  dv2 = fv2/dfv2;
92  temp = v2;
93  v2 -= dv2;
94  if (temp == v2) done = 1;
95  }
96 
97  if (fabs(dv2) < acc*v2 || fabs(fv2) < acc) done = 1;
98 
99  lor2 = 1.0/(1.0 - v2);
100  lor = sqrt(lor2);
101 
102  rho = D/lor;
103 
104  #if EOS == IDEAL
105  x = pow(rho,g_gamma - 1.0);
106  h = 1.0 + g_gamma/(g_gamma - 1.0)*sigma*x;
107  W = D*h*lor;
108  dW_dlor = (2 - g_gamma)*W/lor + (g_gamma - 1.0)*D;
109  #elif EOS == TAUB
110  x = pow(rho,2.0/3.0);
111  th = sigma*x/sqrt(1.0 + 3.0*sigma*x);
112  h = 2.5*th + sqrt(2.25*th*th + 1.0);
113  W = D*h*lor;
114 
115  scrh = -2.0*x/(3.0*lor)*sigma*(2.0*sigma*x - 3.0*th*th);
116  scrh /= 2.0*th*(1.0 + 3.0*sigma*x);
117  dW_dlor = W/lor + D*lor*(5.0*h-8.0*th)/(2.0*h-5.0*th)*scrh;
118  #endif
119 
120  W2 = W*W;
121  S2_W2 = S2/W2;
122 
123  if (done) break;
124 
125  W3 = W2*W;
126  scrh = W + B2;
127 
128  fv2 = v2 - (m2 + S2_W2*(scrh + W))/(scrh*scrh);
129 
130  if (fv2 < 0.0) v2min = v2;
131  else v2max = v2;
132 
133  dv2_dW = S2*(3.0*W*scrh + B2*B2) + m2*W3;
134  dv2_dW *= -2.0/(W3*scrh*scrh*scrh);
135 
136  dfv2 = 1.-0.5*dv2_dW*dW_dlor*lor2*lor;
137  }
138 
139 /* -- Compure pressure */
140 
141  #if EOS == IDEAL
142  par->prs = sigma*x*rho;
143  #elif EOS == TAUB
144  par->prs = th*rho;
145  #endif
146 
147  if (k == MAX_ITER) {
148  WARNING(
149  print ("! EntropySolve: too many iterations,%d, ",k);
150  )
151  return(1);
152  }
153 
154  if (par->prs < 0.0 || sigma < 0.0 || W < 0.0) {
155  WARNING(
156  print ("! EntropySolve: negative pressure, p = %12.6e\n",par->prs);
157  )
158  return(1);
159  }
160 
161  par->W = W;
162  par->lor = lor;
163  par->rho = rho;
164 
165 /* -- Redefine energy -- */
166 
167  #if RMHD_REDUCED_ENERGY == YES
168  #if EOS == IDEAL
169  par->E = par->prs*(g_gamma/(g_gamma-1.0)*lor2 - 1.0) + D*lor2*v2/(1.0 + lor)
170  + 0.5*(1.0 + v2)*B2 - 0.5*S2_W2;
171  #elif EOS == TAUB
172  scrh = 2.25*th*th + v2;
173  par->E = par->prs*(2.5*lor2 - 1.0)
174  + D*lor2/(lor*sqrt(2.25*th*th + 1.0) + 1.0)*scrh
175  + 0.5*(1.0 + v2)*B2 - 0.5*S2_W2;
176  #endif
177  #else
178  par->E = W - par->prs + 0.5*(1.0 + v2)*B2 - 0.5*S2_W2;
179  #endif
180  return(0); /* -- success -- */
181 }
182 #undef MAX_ITER
double g_gamma
Definition: globals.h:112
double m2
Square of total momentum (input).
Definition: mod_defs.h:89
double D
Lab density (input).
Definition: mod_defs.h:86
tuple scrh
Definition: configure.py:200
#define WARNING(a)
Definition: macros.h:217
#define MAX_ITER
Definition: entropy_solve.c:18
#define B2
double rho
proper density (output)
Definition: mod_defs.h:91
int EntropySolve(Map_param *par)
Definition: entropy_solve.c:21
double W
D*h*lor (output).
Definition: mod_defs.h:92
double S
mB (input).
Definition: mod_defs.h:94
int k
Definition: analysis.c:2
void print(const char *fmt,...)
Definition: amrPluto.cpp:497
PLUTO main header file.
double lor
Lorentz factor (output).
Definition: mod_defs.h:93
double prs
Thermal pressure (output).
Definition: mod_defs.h:94
double sigma_c
Conserved entropy (input).
Definition: mod_defs.h:87
double B2
Square of magnetic field (input).
Definition: mod_defs.h:96
double E
Total energy (input).
Definition: mod_defs.h:88