aboutsummaryrefslogtreecommitdiff
path: root/src/nuc_eos/readtable.F90
blob: b865b39d9f290dbef45065d2d42ccaf56b9f3172 (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
#include "cctk.h"
#include "cctk_Parameters.h"


subroutine nuc_eos_readtable(eos_filename)
! This routine reads the table and initializes
! all variables in the module. 

  use eosmodule
  use hdf5 

  implicit none

  DECLARE_CCTK_PARAMETERS

  character(*) eos_filename

  character(len=100) message

! HDF5 vars
  integer(HID_T) file_id,dset_id,dspace_id
  integer(HSIZE_T) dims1(1), dims3(3)
  integer error,rank,accerr
  integer i,j,k

  real*8 amu_cgs_andi
  real*8 buffer1,buffer2,buffer3,buffer4
  accerr=0

!  write(*,*) "Reading Nuclear EOS Table"

  call h5open_f(error)

  call h5fopen_f (trim(adjustl(eos_filename)), H5F_ACC_RDONLY_F, file_id, error)

!  write(6,*) trim(adjustl(eos_filename))

! read scalars
  dims1(1)=1
  call h5dopen_f(file_id, "pointsrho", dset_id, error)
  call h5dread_f(dset_id, H5T_NATIVE_INTEGER, nrho, dims1, error)
  call h5dclose_f(dset_id,error)

  if(error.ne.0) then
     call CCTK_WARN (0, "Could not read EOS table file")
  endif

  dims1(1)=1
  call h5dopen_f(file_id, "pointstemp", dset_id, error)
  call h5dread_f(dset_id, H5T_NATIVE_INTEGER, ntemp, dims1, error)
  call h5dclose_f(dset_id,error)

  if(error.ne.0) then
     call CCTK_WARN (0, "Could not read EOS table file")
  endif

  dims1(1)=1
  call h5dopen_f(file_id, "pointsye", dset_id, error)
  call h5dread_f(dset_id, H5T_NATIVE_INTEGER, nye, dims1, error)
  call h5dclose_f(dset_id,error)

  if(error.ne.0) then
     call CCTK_WARN (0, "Could not read EOS table file")
  endif

!  write(message,"(a25,1P3i5)") "We have nrho ntemp nye: ", nrho,ntemp,nye
!  write(*,*) message

  allocate(alltables(nrho,ntemp,nye,nvars))

  ! index variable mapping:
  !  1 -> logpress
  !  2 -> logenergy
  !  3 -> entropy
  !  4 -> munu
  !  5 -> cs2
  !  6 -> dedT
  !  7 -> dpdrhoe
  !  8 -> dpderho
  !  9 -> muhat
  ! 10 -> mu_e
  ! 11 -> mu_p
  ! 12 -> mu_n
  ! 13 -> xa
  ! 14 -> xh
  ! 15 -> xn
  ! 16 -> xp
  ! 17 -> abar
  ! 18 -> zbar


  dims3(1)=nrho
  dims3(2)=ntemp
  dims3(3)=nye
  call h5dopen_f(file_id, "logpress", dset_id, error)
  call h5dread_f(dset_id, H5T_NATIVE_DOUBLE, alltables(:,:,:,1), dims3, error)
  call h5dclose_f(dset_id,error)
  accerr=accerr+error
  call h5dopen_f(file_id, "logenergy", dset_id, error)
  call h5dread_f(dset_id, H5T_NATIVE_DOUBLE, alltables(:,:,:,2), dims3, error)
  call h5dclose_f(dset_id,error)
  accerr=accerr+error
  call h5dopen_f(file_id, "entropy", dset_id, error)
  call h5dread_f(dset_id, H5T_NATIVE_DOUBLE, alltables(:,:,:,3), dims3, error)
  call h5dclose_f(dset_id,error)
  accerr=accerr+error
  call h5dopen_f(file_id, "munu", dset_id, error)
  call h5dread_f(dset_id, H5T_NATIVE_DOUBLE, alltables(:,:,:,4), dims3, error)
  call h5dclose_f(dset_id,error)
  accerr=accerr+error
  call h5dopen_f(file_id, "cs2", dset_id, error)
  call h5dread_f(dset_id, H5T_NATIVE_DOUBLE, alltables(:,:,:,5), dims3, error)
  call h5dclose_f(dset_id,error)
  accerr=accerr+error
  call h5dopen_f(file_id, "dedt", dset_id, error)
  call h5dread_f(dset_id, H5T_NATIVE_DOUBLE, alltables(:,:,:,6), dims3, error)
  call h5dclose_f(dset_id,error)
  accerr=accerr+error
  call h5dopen_f(file_id, "dpdrhoe", dset_id, error)
  call h5dread_f(dset_id, H5T_NATIVE_DOUBLE, alltables(:,:,:,7), dims3, error)
  call h5dclose_f(dset_id,error)
  accerr=accerr+error
  call h5dopen_f(file_id, "dpderho", dset_id, error)
  call h5dread_f(dset_id, H5T_NATIVE_DOUBLE, alltables(:,:,:,8), dims3, error)
  call h5dclose_f(dset_id,error)
  accerr=accerr+error

! chemical potentials
  call h5dopen_f(file_id, "muhat", dset_id, error)
  call h5dread_f(dset_id, H5T_NATIVE_DOUBLE, alltables(:,:,:,9), dims3, error)
  call h5dclose_f(dset_id,error)
  accerr=accerr+error

  call h5dopen_f(file_id, "mu_e", dset_id, error)
  call h5dread_f(dset_id, H5T_NATIVE_DOUBLE, alltables(:,:,:,10), dims3, error)
  call h5dclose_f(dset_id,error)
  accerr=accerr+error

  call h5dopen_f(file_id, "mu_p", dset_id, error)
  call h5dread_f(dset_id, H5T_NATIVE_DOUBLE, alltables(:,:,:,11), dims3, error)
  call h5dclose_f(dset_id,error)
  accerr=accerr+error

  call h5dopen_f(file_id, "mu_n", dset_id, error)
  call h5dread_f(dset_id, H5T_NATIVE_DOUBLE, alltables(:,:,:,12), dims3, error)
  call h5dclose_f(dset_id,error)
  accerr=accerr+error

! compositions
  call h5dopen_f(file_id, "Xa", dset_id, error)
  call h5dread_f(dset_id, H5T_NATIVE_DOUBLE, alltables(:,:,:,13), dims3, error)
  call h5dclose_f(dset_id,error)
  accerr=accerr+error

  call h5dopen_f(file_id, "Xh", dset_id, error)
  call h5dread_f(dset_id, H5T_NATIVE_DOUBLE, alltables(:,:,:,14), dims3, error)
  call h5dclose_f(dset_id,error)
  accerr=accerr+error

  call h5dopen_f(file_id, "Xn", dset_id, error)
  call h5dread_f(dset_id, H5T_NATIVE_DOUBLE, alltables(:,:,:,15), dims3, error)
  call h5dclose_f(dset_id,error)
  accerr=accerr+error

  call h5dopen_f(file_id, "Xp", dset_id, error)
  call h5dread_f(dset_id, H5T_NATIVE_DOUBLE, alltables(:,:,:,16), dims3, error)
  call h5dclose_f(dset_id,error)
  accerr=accerr+error


! average nucleus
  call h5dopen_f(file_id, "Abar", dset_id, error)
  call h5dread_f(dset_id, H5T_NATIVE_DOUBLE, alltables(:,:,:,17), dims3, error)
  call h5dclose_f(dset_id,error)
  accerr=accerr+error

  call h5dopen_f(file_id, "Zbar", dset_id, error)
  call h5dread_f(dset_id, H5T_NATIVE_DOUBLE, alltables(:,:,:,18), dims3, error)
  call h5dclose_f(dset_id,error)
  accerr=accerr+error

! Gamma
  call h5dopen_f(file_id, "gamma", dset_id, error)
  call h5dread_f(dset_id, H5T_NATIVE_DOUBLE, alltables(:,:,:,19), dims3, error)
  call h5dclose_f(dset_id,error)
  accerr=accerr+error

  allocate(logrho(nrho))
  dims1(1)=nrho
  call h5dopen_f(file_id, "logrho", dset_id, error)
  call h5dread_f(dset_id, H5T_NATIVE_DOUBLE, logrho, dims1, error)
  call h5dclose_f(dset_id,error)
  accerr=accerr+error

  allocate(logtemp(ntemp))
  dims1(1)=ntemp
  call h5dopen_f(file_id, "logtemp", dset_id, error)
  call h5dread_f(dset_id, H5T_NATIVE_DOUBLE, logtemp, dims1, error)
  call h5dclose_f(dset_id,error)
  accerr=accerr+error

  allocate(ye(nye))
  dims1(1)=nye
  call h5dopen_f(file_id, "ye", dset_id, error)
  call h5dread_f(dset_id, H5T_NATIVE_DOUBLE, ye, dims1, error)
  call h5dclose_f(dset_id,error)
  accerr=accerr+error


  call h5dopen_f(file_id, "energy_shift", dset_id, error)
  call h5dread_f(dset_id, H5T_NATIVE_DOUBLE, energy_shift, dims1, error)
  call h5dclose_f(dset_id,error)
  accerr=accerr+error
  
  if(accerr.ne.0) then
    call CCTK_WARN (0, "Problem reading EOS table file")
  endif

  call h5fclose_f (file_id,error)

  if(do_energy_shift.ne.1) then
     energy_shift = 0.0d0
  endif

! not a good idea to call this function as it will screw recovery
!  call h5close_f (error)

  ! set min-max values:

  eos_rhomin = 10.0d0**logrho(1)
  eos_rhomax = 10.0d0**logrho(nrho)

  eos_yemin = ye(1)
  eos_yemax = ye(nye)

  eos_tempmin = 10.0d0**logtemp(1)
  eos_tempmax = 10.0d0**logtemp(ntemp)

!  write(6,*) "Done reading eos tables"


end subroutine nuc_eos_readtable