aboutsummaryrefslogtreecommitdiff
path: root/src/nuc_eos/nuc_eos.F90
blob: e84ac82b2c41a233a6d3049697aba81f6aa52ebe (plain)
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
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
! #########################################################
!
! Copyright C. D. Ott and Evan O'Connor, July 2009
!
! UNITS: density              g/cm^3
!        temperature          MeV
!        ye                   number fraction per baryon
!        energy               erg/g
!        pressure             dyn/cm^2
!        chemical potentials  MeV
!        entropy              k_B / baryon
!        cs2                  cm^2/s^2 (not relativistic)
!
! keyerr --> error output; should be 0
! rfeps --> root finding relative accuracy, set around 1.0d-10
! keytemp: 0 -> coming in with eps
!          1 -> coming in with temperature
!          2 -> coming in with entropy
!
subroutine nuc_eos_full(xrho,xtemp,xye,xenr,xprs,xent,xcs2,xdedt,&
     xdpderho,xdpdrhoe,xxa,xxh,xxn,xxp,xabar,xzbar,xmu_e,xmu_n,xmu_p, &
     xmuhat,keytemp,keyerr,rfeps)

  use eosmodule
  implicit none

  real*8, intent(in)    :: xrho,xye
  real*8, intent(inout) :: xtemp,xenr,xent
  real*8, intent(out)   :: xprs,xcs2,xdedt
  real*8, intent(out)   :: xdpderho,xdpdrhoe,xxa,xxh,xxn,xxp
  real*8, intent(out)   :: xabar,xzbar,xmu_e,xmu_n,xmu_p,xmuhat
  real*8, intent(in)    :: rfeps
  integer, intent(in)   :: keytemp
  integer, intent(out)  :: keyerr


  ! local variables
  real*8 :: lr,lt,y,xx,xeps,leps,xs
  real*8 :: d1,d2,d3
  real*8 :: ff(nvars)
  integer :: keyerrt = 0

  if(xrho.gt.eos_rhomax) then
     stop "nuc_eos: rho > rhomax"
  endif

  if(xrho.lt.eos_rhomin) then
     stop "nuc_eos: rho < rhomin"
  endif

  if(xye.gt.eos_yemax) then
     stop "nuc_eos: ye > yemax"
  endif

  if(xye.lt.eos_yemin) then
     stop "nuc_eos: ye < yemin"
  endif

  if(keytemp.eq.1) then
     if(xtemp.gt.eos_tempmax) then
        stop "nuc_eos: temp > tempmax"
     endif
     
     if(xtemp.lt.eos_tempmin) then
        stop "nuc_eos: temp < tempmin"
     endif
  endif

  lr = log10(xrho)
  lt = log10(xtemp)
  y = xye
  xeps = xenr + energy_shift
  leps = log10(max(xeps,1.0d0))

  keyerr = 0

  if(keytemp.eq.0) then
     !need to find temperature based on xeps
     call findtemp(lr,lt,y,leps,keyerrt,rfeps)
     if(keyerrt.ne.0) then
        stop "Did not find temperature"
     endif
     xtemp = 10.0d0**lt

  elseif(keytemp.eq.2) then
     !need to find temperature based on xent
     xs = xent
     call findtemp_entropy(lr,lt,y,xs,keyerrt,rfeps)
     xtemp = 10.0d0**lt
  endif

  ! have temperature, proceed:
  call findall(lr,lt,y,ff)

  xprs = 10.0d0**ff(1)

  xeps = 10.0d0**ff(2) - energy_shift
  if(keytemp.eq.1.or.keytemp.eq.2) then
     xenr = xeps
  endif

  if(keytemp.eq.1.or.keytemp.eq.0) then
     xent = ff(3)
  endif

  xcs2 = ff(5)

! derivatives
  xdedt = ff(6)

  xdpdrhoe = ff(7)

  xdpderho = ff(8)

! chemical potentials
  xmuhat = ff(9)

  xmu_e = ff(10)

  xmu_p = ff(11)

  xmu_n = ff(12)

! compositions
  xxa = ff(13)

  xxh = ff(14)

  xxn = ff(15)

  xxp = ff(16)

  xabar = ff(17)

  xzbar = ff(18)


end subroutine nuc_eos_full

subroutine nuc_low_eos(xrho,xenr,xprs,xcs2,xdpderho,xdpdrhoe,keytemp)

  implicit none
  real*8 xrho,xenr,xprs,xcs2,xdpderho,xdpdrhoe
  real*8,parameter :: idealK1 =  1.2435d15 * (0.5d0**(4.d0/3.d0))
  real*8,parameter :: idealgamma = 1.41d0
  integer keytemp

  if(keytemp.eq.1) then
!	energy wanted
     xprs=idealK1*xrho**(idealgamma)
     xenr=xprs/xrho/(idealgamma-1.d0)
     xcs2=idealgamma*xprs/xrho
  endif

  xprs = (idealgamma - 1.0d0) * xrho * xenr

  xdpderho = (idealgamma - 1.0d0 ) * xrho
  xdpdrhoe = (idealgamma - 1.0d0 ) * xenr

  xcs2= xdpdrhoe+xdpderho*xprs/xrho**2

end subroutine nuc_low_eos

subroutine nuc_eos_short(xrho,xtemp,xye,xenr,xprs,xent,xcs2,xdedt,&
     xdpderho,xdpdrhoe,xmunu,keytemp,keyerr,rfeps)

  use eosmodule
  implicit none

  real*8, intent(in)    :: xrho,xye
  real*8, intent(inout) :: xtemp,xenr,xent
  real*8, intent(in)    :: rfeps
  real*8, intent(out)   :: xprs,xmunu,xcs2,xdedt
  real*8, intent(out)   :: xdpderho,xdpdrhoe
  integer, intent(in)   :: keytemp
  integer, intent(out)  :: keyerr

  ! local variables
  real*8 :: lr,lt,y,xx,xeps,leps,xs
  real*8 :: d1,d2,d3,ff(8)
  integer :: keyerrt = 0

  if(xrho.gt.eos_rhomax) then
     stop "nuc_eos: rho > rhomax"
  endif

  if(xrho.lt.eos_rhomin*1.2d0) then
     call nuc_low_eos(xrho,xenr,xprs,xcs2,xdpderho,xdpdrhoe,keytemp)
     xent = 4.0d0
     return
  endif

  if(xye.gt.eos_yemax) then
     stop "nuc_eos: ye > yemax"
  endif

  if(xye.lt.eos_yemin) then
     keyerr = 102
     write(6,*) "ye: ", xye
     stop "nuc_eos ye < yemin"
  endif

  if(keytemp.eq.1) then
     if(xtemp.gt.eos_tempmax) then
        stop "nuc_eos: temp > tempmax"
     endif
     
     if(xtemp.lt.eos_tempmin) then
        call nuc_low_eos(xrho,xenr,xprs,xcs2,xdpderho,xdpdrhoe,keytemp)
        xent = 4.0d0
        return
     endif
  endif

  lr = log10(xrho)
  lt = log10(xtemp)
  y = xye
  xeps = xenr + energy_shift
  leps = log10(max(xeps,1.0d0))

  keyerr = 0

  if(keytemp.eq.0) then
     !need to find temperature based on xeps
     call findtemp(lr,lt,y,leps,keyerrt,rfeps)
     if(keyerrt.ne.0) then
        keyerr = keyerrt
        return
     endif
     xtemp = 10.0d0**lt

  elseif(keytemp.eq.2) then
     !need to find temperature based on xent
     xs = xent
     call findtemp_entropy(lr,lt,y,xs,keyerrt,rfeps)
     if(keyerrt.ne.0) then
        keyerr = keyerrt
        return
     endif
     xtemp = 10.0d0**lt

  endif

  ! have temperature, proceed:
  call findall_short(lr,lt,y,ff)
  xprs = 10.0d0**ff(1)

  xeps = 10.0d0**ff(2) - energy_shift
  if((keytemp.eq.1).or.(keytemp.eq.2)) then
     xenr = xeps
  endif

  if(keytemp.eq.1.or.keytemp.eq.0) then
     xent = ff(3)
  endif

  xmunu = ff(4)

  xcs2 = ff(5)

  xdedt = ff(6)

  xdpdrhoe = ff(7)

  xdpderho = ff(8)

end subroutine nuc_eos_short


subroutine nuc_eos_press_eps(xrho,xtemp,xye,xenr,xprs,&
     keytemp,keyerr,rfeps)

  use eosmodule
  implicit none

  real*8, intent(in)    :: xrho,xye
  real*8, intent(inout) :: xtemp,xenr
  real*8, intent(in)    :: rfeps
  real*8, intent(out)   :: xprs
  integer, intent(in)   :: keytemp
  integer, intent(out)  :: keyerr

  ! local variables
  real*8 :: xcs2,xdpderho,xdpdrhoe
  real*8 :: lr,lt,y,xx,xeps,leps,xs
  real*8 :: d1,d2,d3,ff(8)
  integer :: keyerrt = 0

  if(xrho.gt.eos_rhomax) then
     stop "nuc_eos: rho > rhomax"
  endif

  if(xrho.lt.eos_rhomin*1.2d0) then
     call nuc_low_eos(xrho,xenr,xprs,xcs2,xdpderho,xdpdrhoe,keytemp)
     return
  endif

  if(xye.gt.eos_yemax) then
     stop "nuc_eos: ye > yemax"
  endif

  if(xye.lt.eos_yemin) then
     keyerr = 102
     write(6,*) "ye: ", xye
     stop "nuc_eos ye < yemin"
  endif

  if(keytemp.eq.1) then
     if(xtemp.gt.eos_tempmax) then
        stop "nuc_eos: temp > tempmax"
     endif
     
     if(xtemp.lt.eos_tempmin) then
        call nuc_low_eos(xrho,xenr,xprs,xcs2,xdpderho,xdpdrhoe,keytemp)
        return
     endif
  endif

  lr = log10(xrho)
  lt = log10(xtemp)
  y = xye
  xeps = xenr + energy_shift
  leps = log10(max(xeps,1.0d0))

  keyerr = 0

  if(keytemp.eq.0) then
     !need to find temperature based on xeps
     call findtemp(lr,lt,y,leps,keyerrt,rfeps)
     if(keyerrt.ne.0) then
        keyerr = keyerrt
        return
     endif
     xtemp = 10.0d0**lt

  elseif(keytemp.eq.2) then
     stop "eos_nuc_press does not support keytemp.eq.2"
  endif

  ! have temperature, proceed:
  call findall_press_eps(lr,lt,y,ff)
  xprs = 10.0d0**ff(1)

  xeps = 10.0d0**ff(2) - energy_shift
  if((keytemp.eq.1).or.(keytemp.eq.2)) then
     xenr = xeps
  endif

end subroutine nuc_eos_press_eps


subroutine findthis(lr,lt,y,value,array,d1,d2,d3)

  use eosmodule

  implicit none

  integer rip,rim
  integer tip,tim
  integer yip,yim

  real*8 lr,lt,y,value,d1,d2,d3
  real*8 array(*)

! Ewald's interpolator           
  call intp3d(lr,lt,y,value,1,array,nrho,ntemp,nye,logrho,logtemp,ye,d1,d2,d3)


end subroutine findthis


subroutine findall(lr,lt,y,ff)

  use eosmodule
  implicit none

  real*8 ff(nvars)
  real*8 ffx(nvars,1)
  real*8 lr,lt,y
  integer i

! Ewald's interpolator           
  call intp3d_many(lr,lt,y,ffx,1,alltables,&
       nrho,ntemp,nye,nvars,logrho,logtemp,ye)
  ff(:) = ffx(:,1)


end subroutine findall


subroutine findall_short(lr,lt,y,ff)

  use eosmodule
  implicit none

  real*8 ffx(8,1)
  real*8 ff(8)
  real*8 lr,lt,y
  integer i
  integer :: nvarsx = 8


! Ewald's interpolator           
  call intp3d_many(lr,lt,y,ffx,1,alltables(:,:,:,1:8), &
       nrho,ntemp,nye,nvarsx,logrho,logtemp,ye)
  ff(:) = ffx(:,1)

end subroutine findall_short


subroutine findall_press_eps(lr,lt,y,ff)

  use eosmodule
  implicit none

  real*8 ffx(8,1)
  real*8 ff(8)
  real*8 lr,lt,y
  integer i
  integer :: nvarsx = 2


! Ewald's interpolator           
  call intp3d_many(lr,lt,y,ffx,1,alltables(:,:,:,1:2), &
       nrho,ntemp,nye,nvarsx,logrho,logtemp,ye)
  ff(:) = ffx(:,1)

end subroutine findall_press_eps