PLUTO
jacobian.c
Go to the documentation of this file.
1 #include "pluto.h"
2 #include "cooling_defs.h"
3 
4 /* ********************************************************************* */
5 void Jacobian (double *v, double *rhs, double **dfdy)
6 /*!
7  * Compute the jacobian J(k,l) = dfdy
8  *
9  * k = row index
10  * l = column index
11  *
12  * J(0,0) J(0,1) ... J(0, n-1)
13  * J(1,0) J(1,1) .... J(1, n-1)
14  * . . .
15  * . . .
16  * . . .
17  * J(n-1,0) .... J(n-1, n-1)
18  *
19  *
20  * or,
21  *
22  * +-----------------------+
23  * + | |
24  * + | |
25  * + | |
26  * + dX'/dX | dX'/dp |
27  * + (JXX) | (JXp) |
28  * + | |
29  * + | |
30  * +--------------+--------+
31  * + dp'/dX | dp'/dp |
32  * + (JpX) | Jpp |
33  * +-----------------------+
34  *
35  *
36  *
37  *********************************************************************** */
38 {
39  int k, l, nv, n;
40  double dmu_dX[NIONS], N, mu;
41  double T, dLdX, dCdX, dRdX, eps, scrh;
42  double vp[NVAR], vm[NVAR], rhs_m[NVAR], rhs_p[NVAR];
43  double dfdp[NIONS];
44  double *L, *La, *Lb, *Lc;
45  double *C, *Ca, *Cb, *Cc;
46  double *R, *Ra, *Rb, *Rc;
47  double *X, *dnel_dX, *de;
48  double Unit_Time, E_cost;
49  static double **delta;
50 
51  n = NIONS + 1;
52 
53  Unit_Time = UNIT_LENGTH/UNIT_VELOCITY;
54 
55  E_cost = UNIT_LENGTH/UNIT_DENSITY/
57 
58  if (delta == NULL) {
59  delta = ARRAY_2D(NIONS, NIONS, double);
60  for (k = 0; k < NIONS; k++){
61  for (l = 0; l < NIONS; l++){
62  delta[k][l] = (k == l ? 1.0:0.0);
63  }}
64  }
65 
66  Radiat (v, rhs);
67 
68  L = CoolCoeffs.Lrate; La = CoolCoeffs.La; Lb = CoolCoeffs.Lb; Lc = CoolCoeffs.Lc;
69  C = CoolCoeffs.Crate; Ca = CoolCoeffs.Ca; Cb = CoolCoeffs.Cb; Cc = CoolCoeffs.Cc;
70  R = CoolCoeffs.Rrate; Ra = CoolCoeffs.Ra; Rb = CoolCoeffs.Rb; Rc = CoolCoeffs.Rc;
71 
72  de = CoolCoeffs.de;
73 
74  X = v + NFLX;
75  dnel_dX = CoolCoeffs.dnel_dX;
76 
77  N = v[RHO]*find_N_rho(); /* -- Total number density -- */
78  mu = MeanMolecularWeight(v);
79  T = v[PRS]/v[RHO]*KELVIN*mu;
80 
81 /* -- compute the vector grad_X (mu) -- */
82 
83  scrh = 1.0/CoolCoeffs.muD;
84  for (k = 0; k < n - 1; k++){
85  dmu_dX[k] = (CoolCoeffs.dmuN_dX[k]
86  - CoolCoeffs.dmuD_dX[k]*mu)*scrh;
87  }
88 
89 /* -------------------------------------------------
90  Compute dX'/dX, dp'/dX with
91 
92  k = 0....n - 2
93  l = 0....n - 2
94  ------------------------------------------------- */
95 
96  /* -- compute d\dot{X}_0/dXl (H right hand side) separately -- */
97 
98  for (l = 0; l < n - 1; l++) { /* -- first row -- */
99  dfdy[0][l] = - (R[0] + C[0])*delta[0][l] +
100  (1.0 - X[0])*CoolCoeffs.fRH*dnel_dX[l]
101  - X[0] *CoolCoeffs.fCH*dnel_dX[l];
102  }
103 
104  for (k = 1; k < n - 1; k++) { /* -- C's and L's contribution -- */
105  for (l = 0; l < n - 1; l++) {
106  dLdX = La[k]*dnel_dX[l] + Lb[k]*delta[0][l] + Lc[k]*delta[1][l];
107  dCdX = Ca[k]*dnel_dX[l] + Cb[k]*delta[0][l] + Cc[k]*delta[1][l];
108  dfdy[k][l] = L[k]*delta[k - 1][l] + dLdX*X[k - 1] - C[k]*delta[k][l] - dCdX*X[k];
109  }}
110 
111  for (k = 1; k < n - 2; k++) { /* -- R's contribution -- */
112  for (l = 0; l < n - 1; l++) {
113  dRdX = Ra[k]*dnel_dX[l] + Rb[k]*delta[0][l] + Rc[k]*delta[1][l];
114  dfdy[k][l] += R[k]*delta[k + 1][l] + dRdX*X[k + 1];
115  }}
116 
117  for (k = 0; k < n - 1; k++) { /* -- Get correct dimensions -- */
118  for (l = 0; l < n - 1; l++) {
119  dfdy[k][l] *= Unit_Time;
120  }}
121 
122 /* -------------------------------------------------------
123  Compute PDEs of cooling function with respect to X's
124  ------------------------------------------------------- */
125 
126  for (l = 0; l < n - 1; l++){
127 
128  dfdy[n - 1][l] = dnel_dX[l]*rhs[PRS]/CoolCoeffs.Ne;
129 
130  scrh = 0.0;
131  for (nv = 0; nv < NIONS; nv++){
132  scrh += X[nv]*CoolCoeffs.de_dne[nv]*elem_ab[elem_part[nv]];
133  }
134  scrh *= dnel_dX[l];
135  scrh += de[l]*elem_ab[elem_part[l]];
136  scrh += CoolCoeffs.dLIR_dX[0]*delta[0][l];
137  scrh += CoolCoeffs.dLIR_dX[2]*delta[2][l];
138  scrh *= N*CoolCoeffs.Ne*E_cost;
139 
140  dfdy[n - 1][l] -= (g_gamma-1.0)*scrh;
141  }
142 
143 /* --------------------------------------------------------
144  Compute partial derivatives with respect to T
145  using numerical differentiation.
146  -------------------------------------------------------- */
147 
148  eps = 1.e-4;
149 
150  for (nv = 0; nv < NVAR; nv++){
151  vm[nv] = vp[nv] = v[nv];
152  }
153  vp[PRS] = v[PRS]*(1.0 + eps);
154  vm[PRS] = v[PRS]*(1.0 - eps);
155 
156  Radiat (vp, rhs_p);
157  Radiat (vm, rhs_m);
158 
159 /* -- Compute last column (Jxp and Jpp) drhs/dp -- */
160 
161  for (k = 0; k < n - 1; k++){
162  dfdp[k] = (rhs_p[k + NFLX] - rhs_m[k + NFLX])/(2.0*eps*v[PRS]);
163  dfdy[k][n - 1] = dfdp[k];
164  }
165  dfdp[n - 1] = (rhs_p[PRS] - rhs_m[PRS])/(2.0*eps*v[PRS]);
166  dfdy[n - 1][n - 1] = dfdp[n - 1];
167 
168 /* -- Add df/dT*dT/dX term to all columns except last one -- */
169 
170  scrh = v[PRS]/mu;
171  for (k = 0; k < n ; k++) {
172  for (l = 0; l < n - 1; l++) {
173  dfdy[k][l] += dfdp[k]*scrh*dmu_dX[l];
174  }}
175 
176 }
177 
178 
double dnel_dX[NIONS]
Definition: cooling_defs.h:51
tuple T
Definition: Sph_disk.py:33
double dLIR_dX[NIONS]
Definition: cooling_defs.h:55
double Crate[NIONS]
Definition: cooling_defs.h:31
double g_gamma
Definition: globals.h:112
double de_dne[NIONS]
Definition: cooling_defs.h:47
#define UNIT_DENSITY
Unit density in gr/cm^3.
Definition: pluto.h:369
double Cc[NIONS]
Definition: cooling_defs.h:41
double Rb[NIONS]
Definition: cooling_defs.h:43
static int n
Definition: analysis.c:3
void Radiat(double *, double *)
Definition: radiat.c:94
const double elem_ab[8]
Definition: radiat.c:9
#define RHO
Definition: mod_defs.h:19
tuple scrh
Definition: configure.py:200
COOL_COEFF CoolCoeffs
Definition: radiat.c:38
double de[NIONS]
Definition: cooling_defs.h:46
double Lrate[NIONS]
Definition: cooling_defs.h:30
#define KELVIN
Definition: pluto.h:401
double Lb[NIONS]
Definition: cooling_defs.h:37
double La[NIONS]
Definition: cooling_defs.h:36
#define UNIT_VELOCITY
Unit velocity in cm/sec.
Definition: pluto.h:377
#define eps
#define X
#define NFLX
Definition: mod_defs.h:32
#define UNIT_LENGTH
Unit Length in cm.
Definition: pluto.h:373
double dmuN_dX[NIONS]
Definition: cooling_defs.h:53
#define NIONS
Definition: cooling.h:28
void Jacobian(real *v, real *rhs, real **dfdy)
Definition: jacobian.c:4
int k
Definition: analysis.c:2
const int elem_part[31]
Definition: radiat.c:16
double Rrate[NIONS]
Definition: cooling_defs.h:32
double find_N_rho()
Definition: radiat.c:694
double Rc[NIONS]
Definition: cooling_defs.h:44
double dmuD_dX[NIONS]
Definition: cooling_defs.h:54
PLUTO main header file.
double MeanMolecularWeight(double *V)
double Ra[NIONS]
Definition: cooling_defs.h:42
double Lc[NIONS]
Definition: cooling_defs.h:38
double Ca[NIONS]
Definition: cooling_defs.h:39
#define ARRAY_2D(nx, ny, type)
Definition: prototypes.h:171
#define NVAR
Definition: pluto.h:609
double Cb[NIONS]
Definition: cooling_defs.h:40