/
NRPyEOS_readtable_set_EOS_params.c
218 lines (181 loc) · 9.37 KB
/
NRPyEOS_readtable_set_EOS_params.c
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
#include "NRPyEOS.h"
// mini NoMPI
#ifdef HAVE_CAPABILITY_MPI
#include <mpi.h>
#define BCAST(buffer, size) MPI_Bcast(buffer, size, MPI_BYTE, my_reader_process, MPI_COMM_WORLD)
#else
#define BCAST(buffer, size) do { /* do nothing */ } while(0)
#endif
// If on the IO proc (doIO == True) actually perform HDF5 IO, catch possible
// HDF5 errors
#define HDF5_DO_IO(fn_call) \
{ \
int _error_code = fn_call; \
if (_error_code < 0) { \
fprintf(stderr,"(NRPyEOS) HDF5 call '%s' returned error code %d", \
#fn_call, _error_code); \
} \
}
void NRPyEOS_readtable_set_EOS_params(const char *nuceos_table_name, NRPyEOS_params_tabulated *restrict eos_params) {
printf("(NRPyEOS) *******************************\n");
printf("(NRPyEOS) Reading EOS table from file:\n");
printf("(NRPyEOS) %s\n",nuceos_table_name);
printf("(NRPyEOS) *******************************\n");
hid_t file;
HDF5_DO_IO(file = H5Fopen(nuceos_table_name, H5F_ACC_RDONLY, H5P_DEFAULT));
// Use these two defines to easily read in a lot of variables in the same way
// The first reads in one variable of a given type completely
#define READ_BCAST_EOS_HDF5(NAME,VAR,TYPE,MEM,NELEMS) \
do { \
hid_t dataset; \
HDF5_DO_IO(dataset = H5Dopen(file, NAME, H5P_DEFAULT)); \
HDF5_DO_IO(H5Dread(dataset, TYPE, MEM, H5S_ALL, H5P_DEFAULT, VAR)); \
BCAST (VAR, sizeof(*(VAR))*(NELEMS)); \
HDF5_DO_IO(H5Dclose(dataset)); \
} while (0)
// The second reads a given variable into a hyperslab of the alltables_temp array
#define READ_BCAST_EOSTABLE_HDF5(NAME,OFF,DIMS) \
do { \
READ_BCAST_EOS_HDF5(NAME,&alltables_temp[(OFF)*(DIMS)[1]],H5T_NATIVE_DOUBLE,H5S_ALL,(DIMS)[1]); \
} while (0)
// Read size of tables
READ_BCAST_EOS_HDF5("pointsrho", &eos_params->nrho, H5T_NATIVE_INT, H5S_ALL, 1);
READ_BCAST_EOS_HDF5("pointstemp", &eos_params->ntemp, H5T_NATIVE_INT, H5S_ALL, 1);
READ_BCAST_EOS_HDF5("pointsye", &eos_params->nye, H5T_NATIVE_INT, H5S_ALL, 1);
// Allocate memory for tables
double* alltables_temp;
if (!(alltables_temp = (double*)malloc(eos_params->nrho * eos_params->ntemp * eos_params->nye * NRPyEOS_ntablekeys * sizeof(double)))) {
fprintf(stderr,"(NRPyEOS) Cannot allocate memory for EOS table");
}
if (!(eos_params->logrho = (double*)malloc(eos_params->nrho * sizeof(double)))) {
fprintf(stderr,"(NRPyEOS) Cannot allocate memory for EOS table");
}
if (!(eos_params->logtemp = (double*)malloc(eos_params->ntemp * sizeof(double)))) {
fprintf(stderr,"(NRPyEOS) Cannot allocate memory for EOS table");
}
if (!(eos_params->yes = (double*)malloc(eos_params->nye * sizeof(double)))) {
fprintf(stderr,"(NRPyEOS) Cannot allocate memory for EOS table");
}
// Prepare HDF5 to read hyperslabs into alltables_temp
hsize_t table_dims[2] = {NRPyEOS_ntablekeys, (hsize_t)eos_params->nrho * eos_params->ntemp * eos_params->nye};
hid_t mem3 = H5Screate_simple(2, table_dims, NULL);
// Read alltables_temp
READ_BCAST_EOSTABLE_HDF5("logpress", 0, table_dims);
READ_BCAST_EOSTABLE_HDF5("logenergy", 1, table_dims);
READ_BCAST_EOSTABLE_HDF5("entropy", 2, table_dims);
READ_BCAST_EOSTABLE_HDF5("munu", 3, table_dims);
READ_BCAST_EOSTABLE_HDF5("cs2", 4, table_dims);
READ_BCAST_EOSTABLE_HDF5("dedt", 5, table_dims);
READ_BCAST_EOSTABLE_HDF5("dpdrhoe", 6, table_dims);
READ_BCAST_EOSTABLE_HDF5("dpderho", 7, table_dims);
// chemical potentials
READ_BCAST_EOSTABLE_HDF5("muhat", 8, table_dims);
READ_BCAST_EOSTABLE_HDF5("mu_e", 9, table_dims);
READ_BCAST_EOSTABLE_HDF5("mu_p", 10, table_dims);
READ_BCAST_EOSTABLE_HDF5("mu_n", 11, table_dims);
// compositions
READ_BCAST_EOSTABLE_HDF5("Xa", 12, table_dims);
READ_BCAST_EOSTABLE_HDF5("Xh", 13, table_dims);
READ_BCAST_EOSTABLE_HDF5("Xn", 14, table_dims);
READ_BCAST_EOSTABLE_HDF5("Xp", 15, table_dims);
// average nucleus
READ_BCAST_EOSTABLE_HDF5("Abar", 16, table_dims);
READ_BCAST_EOSTABLE_HDF5("Zbar", 17, table_dims);
// Gamma
READ_BCAST_EOSTABLE_HDF5("gamma", 18, table_dims);
// Read additional tables and variables
READ_BCAST_EOS_HDF5("logrho", eos_params->logrho, H5T_NATIVE_DOUBLE, H5S_ALL, eos_params->nrho);
READ_BCAST_EOS_HDF5("logtemp", eos_params->logtemp, H5T_NATIVE_DOUBLE, H5S_ALL, eos_params->ntemp);
READ_BCAST_EOS_HDF5("ye", eos_params->yes, H5T_NATIVE_DOUBLE, H5S_ALL, eos_params->nye);
READ_BCAST_EOS_HDF5("energy_shift", &eos_params->energy_shift, H5T_NATIVE_DOUBLE, H5S_ALL, 1);
HDF5_DO_IO(H5Sclose(mem3));
HDF5_DO_IO(H5Fclose(file));
// change ordering of alltables array so that
// the table kind is the fastest changing index
if (!(eos_params->alltables = (double*)malloc(eos_params->nrho * eos_params->ntemp * eos_params->nye * NRPyEOS_ntablekeys
* sizeof(double)))) {
fprintf(stderr,"(NRPyEOS) Cannot allocate memory for EOS table");
}
for(int iv = 0;iv<NRPyEOS_ntablekeys;iv++)
for(int k = 0; k<eos_params->nye;k++)
for(int j = 0; j<eos_params->ntemp; j++)
for(int i = 0; i<eos_params->nrho; i++) {
int indold = i + eos_params->nrho*(j + eos_params->ntemp*(k + eos_params->nye*iv));
int indnew = iv + NRPyEOS_ntablekeys*(i + eos_params->nrho*(j + eos_params->ntemp*k));
eos_params->alltables[indnew] = alltables_temp[indold];
}
// free memory of temporary array
free(alltables_temp);
// convert units, convert logs to natural log
// The latter is great, because exp() is way faster than pow()
// pressure
eos_params->energy_shift = eos_params->energy_shift * EPSGF;
for(int i=0;i<eos_params->nrho;i++) {
// rewrite:
//logrho[i] = log(pow(10.0,logrho[i]) * RHOGF);
// by using log(a^b*c) = b*log(a)+log(c)
eos_params->logrho[i] = eos_params->logrho[i] * log(10.) + log(RHOGF);
}
for(int i=0;i<eos_params->ntemp;i++) {
//logtemp[i] = log(pow(10.0,logtemp[i]));
eos_params->logtemp[i] = eos_params->logtemp[i]*log(10.0);
}
// allocate epstable; a linear-scale eps table
// that allows us to extrapolate to negative eps
if (!(eos_params->epstable = (double*)malloc(eos_params->nrho * eos_params->ntemp * eos_params->nye
* sizeof(double)))) {
fprintf(stderr,"(NRPyEOS) Cannot allocate memory for eps table\n");
}
// convert units
for(int i=0;i<eos_params->nrho*eos_params->ntemp*eos_params->nye;i++) {
{ // pressure
int idx = 0 + NRPyEOS_ntablekeys*i;
eos_params->alltables[idx] = eos_params->alltables[idx] * log(10.0) + log(PRESSGF);
}
{ // eps
int idx = 1 + NRPyEOS_ntablekeys*i;
eos_params->alltables[idx] = eos_params->alltables[idx] * log(10.0) + log(EPSGF);
eos_params->epstable[i] = exp(eos_params->alltables[idx]);
}
{ // cs2
int idx = 4 + NRPyEOS_ntablekeys*i;
eos_params->alltables[idx] *= LENGTHGF*LENGTHGF/TIMEGF/TIMEGF;
}
{ // dedT
int idx = 5 + NRPyEOS_ntablekeys*i;
eos_params->alltables[idx] *= EPSGF;
}
{ // dpdrhoe
int idx = 6 + NRPyEOS_ntablekeys*i;
eos_params->alltables[idx] *= PRESSGF/RHOGF;
}
{ // dpderho
int idx = 7 + NRPyEOS_ntablekeys*i;
eos_params->alltables[idx] *= PRESSGF/EPSGF;
}
}
eos_params->temp0 = exp(eos_params->logtemp[0]);
eos_params->temp1 = exp(eos_params->logtemp[1]);
// set up some vars
eos_params->dtemp = (eos_params->logtemp[eos_params->ntemp-1] - eos_params->logtemp[0]) / (1.0*(eos_params->ntemp-1));
eos_params->dtempi = 1.0/eos_params->dtemp;
eos_params->dlintemp = eos_params->temp1-eos_params->temp0;
eos_params->dlintempi = 1.0/eos_params->dlintemp;
eos_params->drho = (eos_params->logrho[eos_params->nrho-1] - eos_params->logrho[0]) / (1.0*(eos_params->nrho-1));
eos_params->drhoi = 1.0/eos_params->drho;
eos_params->dye = (eos_params->yes[eos_params->nye-1] - eos_params->yes[0]) / (1.0*(eos_params->nye-1));
eos_params->dyei = 1.0/eos_params->dye;
eos_params->drhotempi = eos_params->drhoi * eos_params->dtempi;
eos_params->drholintempi = eos_params->drhoi * eos_params->dlintempi;
eos_params->drhoyei = eos_params->drhoi * eos_params->dyei;
eos_params->dtempyei = eos_params->dtempi * eos_params->dyei;
eos_params->dlintempyei = eos_params->dlintempi * eos_params->dyei;
eos_params->drhotempyei = eos_params->drhoi * eos_params->dtempi * eos_params->dyei;
eos_params->drholintempyei = eos_params->drhoi * eos_params->dlintempi * eos_params->dyei;
eos_params->eos_rhomax = exp(eos_params->logrho[eos_params->nrho-1]);
eos_params->eos_rhomin = exp(eos_params->logrho[0]);
eos_params->eos_tempmax = exp(eos_params->logtemp[eos_params->ntemp-1]);
eos_params->eos_tempmin = exp(eos_params->logtemp[0]);
eos_params->eos_yemax = eos_params->yes[eos_params->nye-1];
eos_params->eos_yemin = eos_params->yes[0];
}