PLUTO
entropy_solve.c
Go to the documentation of this file.
1 /* ///////////////////////////////////////////////////////////////////// */
2 /*!
3  \file
4  \brief Inversion scheme for RHD using entropy.
5 
6  Convert the conservative variables u=[D, m, sigma_c]
7  (where sigma_c = D*sigma is the conserved entropy) to
8  primitive variable using a Newton-Raphson/Bisection scheme.
9 
10  \authors C. Zanni \n
11  A. Mignone
12 
13  \date Oct 5, 2012
14 */
15 /* ///////////////////////////////////////////////////////////////////// */
16 
17 #include "pluto.h"
18 
19 #define MAX_ITER 20
20 /* ********************************************************************* */
22 /*!
23  * Convert the conservative variables {D, m, sigma_c}
24  * (where sigma_c = D*sigma is the conserved entropy) to
25  * primitive variable using a Newton-Raphson/Bisection scheme.
26  *
27  * LAST MODIFIED
28  *
29  * February 5th (2012) by C. Zanni & A. Mignone
30  * (zanni@oato.inaf.it, mignone@oato.inaf.it)
31  *
32  *********************************************************************** */
33 {
34  int k, done;
35  double v2, v2max, v2min, acc = 1.e-13;
36  double h, W, W2, W3;
37  double scrh, x, temp;
38  double lor, lor2, rho, sigma, th;
39  double dW_dlor;
40  double fv2, dfv2, dv2, dv2old;
41  double m2, D;
42 
43  D = par->D;
44  m2 = par->m2;
45 
46  sigma = par->sigma_c/D;
47  v2max = 1.0 - 1.e-8;
48  v2min = 0.0;
49 
50  done = 0;
51 
52  v2 = 0.5*(v2max+v2min);
53  dv2old = v2max-v2min;
54  dv2 = dv2old;
55 
56  lor2 = 1.0/(1.0 - v2);
57  lor = sqrt(lor2);
58  rho = D/lor;
59 
60  #if EOS == IDEAL
61  x = pow(rho,g_gamma - 1.0);
62  h = 1.0 + g_gamma/(g_gamma - 1.0)*sigma*x;
63  W = D*h*lor;
64  dW_dlor = (2 - g_gamma)*W/lor + (g_gamma - 1.0)*D;
65  #elif EOS == TAUB
66  x = pow(rho,2.0/3.0);
67  th = sigma*x/sqrt(1.0 + 3.0*sigma*x);
68  h = 2.5*th + sqrt(2.25*th*th + 1.0);
69  W = D*h*lor;
70 
71  scrh = -2.0*x/(3.0*lor)*sigma*(2.0*sigma*x - 3.0*th*th);
72  scrh /= 2.0*th*(1.0 + 3.0*sigma*x);
73  dW_dlor = W/lor + D*lor*(5.0*h - 8.0*th)/(2.0*h - 5.0*th)*scrh;
74  #endif
75 
76  W2 = W*W;
77  W3 = W2*W;
78 
79  fv2 = v2 - m2/W2;
80  dfv2 = 1.0 + m2/W3*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 
122  if (done) break;
123 
124  W3 = W2*W;
125  fv2 = v2 - m2/W2;
126 
127  if (fv2 < 0.0) v2min = v2;
128  else v2max = v2;
129 
130  dfv2 = 1.0 + m2/W3*dW_dlor*lor2*lor;
131  }
132 
133  #if EOS == IDEAL
134  par->prs = sigma*x*rho;
135  #elif EOS == TAUB
136  par->prs = th*rho;
137  #endif
138 
139  if (k == MAX_ITER) {
140  WARNING(
141  print ("! EntropySolve: too many iterations,%d, ",k);
142  )
143  par->prs = g_smallPressure;
144  par->E = W - par->prs;
145  return(1);
146  }
147 
148  if (par->prs < 0.0 || sigma < 0.0 || W < 0.0) {
149  WARNING(
150  print ("! EntropySolve: negative pressure, p = %12.6e\n",par->prs);
151  )
152  par->prs = g_smallPressure;
153  par->E = W - par->prs;
154  return(1);
155  }
156 
157  par->rho = par->D/lor;
158  par->lor = lor;
159  par->E = W - par->prs; /* redefine energy */
160  return(0); /* -- success -- */
161 }
162 #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
double g_smallPressure
Small value for pressure fix.
Definition: globals.h:110
double rho
proper density (output)
Definition: mod_defs.h:91
int EntropySolve(Map_param *par)
Definition: entropy_solve.c:21
int k
Definition: analysis.c:2
void print(const char *fmt,...)
Definition: amrPluto.cpp:497
#define MAX_ITER
Definition: entropy_solve.c:19
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 E
Total energy (input).
Definition: mod_defs.h:88