PLUTO
pressure_fix.c File Reference

Inversion scheme for RMHD using a pressure fix. More...

#include "pluto.h"
Include dependency graph for pressure_fix.c:

Go to the source code of this file.

Macros

#define MAX_ITER   50
 

Functions

static double VelocitySquareFunc (double, void *)
 
int PressureFix (Map_param *par)
 

Detailed Description

Inversion scheme for RMHD using a pressure fix.

Fix p to a small value, solve for the square of velocity by using secant or bisection algorithm applied to Eq (A3). This step involved re-computing W at each step of the iteration. Once the root has been found, we recompute total energy E. Return 0 if succesful, 1 otherwise.

Authors
A. Mignone
C. Zanni
Date
June 25, 2015

Definition in file pressure_fix.c.

Macro Definition Documentation

#define MAX_ITER   50

Definition at line 20 of file pressure_fix.c.

Function Documentation

int PressureFix ( Map_param par)

Fix p to a small value, solve for the square of velocity by using secant algorithm applied to Eq (9) of Mignone, Plewa & Bodo (2005). This step involved re-computing W at each step of the iteration. Once the root has been found, we recompute total energy E. Return 0 if succesful, 1 otherwise.

Definition at line 23 of file pressure_fix.c.

27 {
28  int k, done=0;
29  double v2, f, dW, S2_W2;
30  double fmin, fmax, v2min, v2max;
31 
32  par->prs = g_smallPressure;
33 
34 /* ------------------------------------------------
35  Use secant algorithm
36  ------------------------------------------------ */
37 
38 /*
39  v2max = 1.0-1.e-8;
40  v2c = 0.95;
41  fc = VelocitySquareFunc(v2c, par);
42  v2 = 0.96;
43  for (k = 1; k < MAX_ITER; k++){
44  f = VelocitySquareFunc(v2, par);
45 print ("%d, v2 = %12.6e f = %12.6e\n",k,v2,f);
46  if (done == 1) break;
47  dW = (v2 - v2c)/(f - fc)*f;
48  v2c = v2; fc = f;
49  v2 -= dW;
50  v2 = MIN(v2max,v2);
51  v2 = MAX(v2, 0.0);
52  if (fabs(f) < 1.e-9) done = 1;
53  }
54 
55  if (v2 >= 1.0 || k >= MAX_ITER) {
56  FILE *fp;
57  fp = fopen("pressure_fix.dat","w");
58  for (v2 = 0.99; v2 < 1.0; v2 += 1.e-5){
59  f = VelocitySquareFunc (v2, par);
60  fprintf (fp, "%12.6e %12.6e\n", v2, f);
61  }
62  fclose(fp);
63  print ("! PressureFix: too many iter while fixing p , v^2 = %f\n", v2);
64 
65  return (1);
66  }
67 */
68 
69 /* ------------------------------------------------
70  Use bisection algorithm
71  ------------------------------------------------ */
72 
73  v2min = 0.0; fmin = VelocitySquareFunc(v2min, par);
74  v2max = 1.0-1.e-9; fmax = VelocitySquareFunc(v2max, par);
75 
76  for (k = 1; k < MAX_ITER; k++){
77 
78  v2 = 0.5*(v2min + v2max);
79  f = VelocitySquareFunc(v2, par);
80 
81  if (f*fmin > 0.0){
82  v2min = v2; fmin = f;
83  }else{
84  v2max = v2; fmax = f;
85  }
86  if (fabs(f) < 1.e-9) break;
87 
88  if (v2 >= 1.0 || k >= MAX_ITER) {
89 /*
90  FILE *fp;
91  fp = fopen("pressure_fix.dat","w");
92  for (v2 = 0.99; v2 < 1.0; v2 += 1.e-5){
93  f = VelocitySquareFunc (v2, par);
94  fprintf (fp, "%12.6e %12.6e\n", v2, f);
95  }
96  fclose(fp);
97 */
98  print ("! PressureFix(): too many iterations while fixing p , v^2 = %f\n", v2);
99  return (1);
100  }
101  }
102 
103 /* -- Redefine energy, proper density and entropy -- */
104 
105  S2_W2 = par->S2/(par->W*par->W);
106 #if RMHD_REDUCED_ENERGY == YES
107  par->E = par->W - par->D - par->prs + 0.5*(1.0 + v2)*par->B2 - 0.5*S2_W2;
108 #else
109  par->E = par->W - par->prs + 0.5*(1.0 + v2)*par->B2 - 0.5*S2_W2;
110 #endif
111  par->rho = par->D/par->lor;
112 
113 /* -- Recompute entropy consistently -- */
114 
115 #if ENTROPY_SWITCH
116 {
117  double rho = par->rho;
118  double th = par->prs/rho;
119  #if EOS == IDEAL
120  par->sigma_c = par->prs*par->lor/pow(rho,g_gamma-1);
121  #elif EOS == TAUB
122  par->sigma_c = par->prs*par->lor/pow(rho,2.0/3.0)*(1.5*th + sqrt(2.25*th*th + 1.0));
123  #endif
124 }
125 #endif
126 
127  return(0); /* -- success -- */
128 }
double g_gamma
Definition: globals.h:112
double D
Lab density (input).
Definition: mod_defs.h:86
static double VelocitySquareFunc(double, void *)
Definition: pressure_fix.c:131
double g_smallPressure
Small value for pressure fix.
Definition: globals.h:110
#define MAX_ITER
Definition: pressure_fix.c:20
double S2
Square of S (input).
Definition: mod_defs.h:95
double rho
proper density (output)
Definition: mod_defs.h:91
double W
D*h*lor (output).
Definition: mod_defs.h:92
int k
Definition: analysis.c:2
void print(const char *fmt,...)
Definition: amrPluto.cpp:497
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

Here is the call graph for this function:

Here is the caller graph for this function:

double VelocitySquareFunc ( double  v2,
void *  par 
)
static

Implement Eq (A3) of Mignone & McKinney (2007).

Definition at line 131 of file pressure_fix.c.

137 {
138  double lor2, W2, f, pg;
139  Map_param *p = (Map_param *) par;
140 
141  lor2 = 1.0/(1.0 - v2);
142  p->lor = sqrt(lor2);
143 
144  pg = p->prs*p->lor;
145  #if EOS == IDEAL
146  p->W = (p->D + pg*g_gamma/(g_gamma - 1.0))*p->lor;
147  #elif EOS == TAUB
148  p->W = (2.5*pg + sqrt(2.25*pg*pg + p->D*p->D))*p->lor;
149  #endif
150 
151  W2 = p->W*p->W;
152 
153  f = p->S2*(2.0*p->W + p->B2) + p->m2*W2;
154  f /= (p->W + p->B2)*(p->W + p->B2)*W2;
155  f -= v2;
156  return f;
157 }
#define EOS
Definition: pluto.h:341
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
#define TAUB
Definition: pluto.h:47
double S2
Square of S (input).
Definition: mod_defs.h:95
double W
D*h*lor (output).
Definition: mod_defs.h:92
double lor
Lorentz factor (output).
Definition: mod_defs.h:93
double prs
Thermal pressure (output).
Definition: mod_defs.h:94
double B2
Square of magnetic field (input).
Definition: mod_defs.h:96

Here is the caller graph for this function: