forked from AMReX-Astro/Microphysics
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathreact_util.H
293 lines (242 loc) · 6.23 KB
/
react_util.H
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
#ifndef REACT_UTIL_H
#define REACT_UTIL_H
#include <extern_parameters.H>
using namespace amrex::literals;
using namespace unit_test_rp;
struct init_t {
int nprim;
int is1;
int is2;
int is3;
int n1;
int n2;
int n3;
amrex::Real Xp_min;
amrex::Real Xp_max;
};
AMREX_INLINE
init_t setup_composition(const int nz) {
// get the primary species indices
init_t comp_data;
comp_data.nprim = 0;
// we absolutely require primary_species_1 to be defined
comp_data.is1 = network_spec_index(primary_species_1);
if (comp_data.is1 >= 0) {
comp_data.nprim++;
} else {
amrex::Error("Error: invalid primary_species_1");
}
// we'll check the next species, but if it is not valid,
// we'll just move on
comp_data.is2 = network_spec_index(primary_species_2);
if (comp_data.is2 >= 0) {
comp_data.nprim++;
// only consider primary_species_3 if primary_species_2
// was defined
comp_data.is3 = network_spec_index(primary_species_3);
if (comp_data.is3 >= 0) {
comp_data.nprim++;
}
}
if (comp_data.nprim == 0) {
amrex::Error("ERROR: no primary species set");
}
// figure out how many zones to allocate to the each of the primary
// species and the extrema for the primary species
if (comp_data.nprim == 1) {
comp_data.n1 = nz;
comp_data.n2 = 0;
comp_data.n3 = 0;
comp_data.Xp_min = 0.2_rt;
comp_data.Xp_max = 0.9_rt;
} else if (comp_data.nprim == 2) {
comp_data.n1 = nz/2;
comp_data.n2 = nz - comp_data.n1;
comp_data.n3 = 0;
comp_data.Xp_min = 0.2_rt;
comp_data.Xp_max = 0.8_rt;
} else if (comp_data.nprim == 3) {
comp_data.n1 = nz/3;
comp_data.n2 = nz/3;
comp_data.n3 = nz - comp_data.n1 - comp_data.n2;
comp_data.Xp_min = 0.2_rt;
comp_data.Xp_max = 0.7_rt;
}
return comp_data;
}
AMREX_INLINE AMREX_GPU_HOST_DEVICE
void get_xn(const int k, const init_t cd, amrex::Real *xn_zone, bool uniform_composition=false) {
for (int n = 0; n < NumSpec; n++) {
xn_zone[n] = 0.0_rt;
}
if (uniform_composition) {
for (int n = 0; n < NumSpec; ++n) {
xn_zone[n] = 1.0_rt / NumSpec;
}
}
else {
if (k < cd.n1) {
if (cd.nprim >= 2) xn_zone[cd.is2] = cd.Xp_min/2;
if (cd.nprim >= 3) xn_zone[cd.is3] = cd.Xp_min/2;
amrex::Real dX = (cd.Xp_max - cd.Xp_min) /
amrex::max((cd.n1 - 1), 1);
xn_zone[cd.is1] = cd.Xp_min + k * dX;
} else if (cd.nprim >= 2 && k < cd.n1 + cd.n2) {
xn_zone[cd.is1] = cd.Xp_min/2;
if (cd.nprim >= 3) {
xn_zone[cd.is3] = cd.Xp_min/2;
}
amrex::Real dX = (cd.Xp_max - cd.Xp_min) /
amrex::max((cd.n2 - 1), 1);
xn_zone[cd.is2] = cd.Xp_min + (k - cd.n1) * dX;
} else {
xn_zone[cd.is1] = cd.Xp_min/2;
xn_zone[cd.is2] = cd.Xp_min/2;
amrex::Real dX = (cd.Xp_max - cd.Xp_min) /
amrex::max((cd.n3 - 1), 1);
xn_zone[cd.is3] = cd.Xp_min + (k - (cd.n1 + cd.n2)) * dX;
}
amrex::Real excess = 0.0_rt;
for (int n = 0; n < NumSpec; n++) {
excess += xn_zone[n];
}
excess = 1.0_rt - excess;
for (int n = 0; n < NumSpec; n++) {
if (n == cd.is1 ||
(cd.nprim >= 2 && n == cd.is2) ||
(cd.nprim >= 3 && n == cd.is3)) {
continue;
}
xn_zone[n] = excess / (NumSpec - cd.nprim);
}
}
// normalize -- just in case
amrex::Real sum_X = 0.0_rt;
for (int n = 0; n < NumSpec; n++) {
sum_X += xn_zone[n];
}
for (int n = 0; n < NumSpec; n++) {
xn_zone[n] = xn_zone[n] / sum_X;
}
}
///
/// return the value of the runtime parameter XN, where N is
/// the integer index passed in
///
AMREX_INLINE AMREX_GPU_HOST_DEVICE
amrex::Real get_xn(const int index, bool uniform_composition=false) {
amrex::Real mass_fraction{};
if (uniform_composition) {
mass_fraction = 1.0_rt / NumSpec;
return mass_fraction;
}
switch (index) {
case 1:
mass_fraction = X1;
break;
case 2:
mass_fraction = X2;
break;
case 3:
mass_fraction = X3;
break;
case 4:
mass_fraction = X4;
break;
case 5:
mass_fraction = X5;
break;
case 6:
mass_fraction = X6;
break;
case 7:
mass_fraction = X7;
break;
case 8:
mass_fraction = X8;
break;
case 9:
mass_fraction = X9;
break;
case 10:
mass_fraction = X10;
break;
case 11:
mass_fraction = X11;
break;
case 12:
mass_fraction = X12;
break;
case 13:
mass_fraction = X13;
break;
case 14:
mass_fraction = X14;
break;
case 15:
mass_fraction = X15;
break;
case 16:
mass_fraction = X16;
break;
case 17:
mass_fraction = X17;
break;
case 18:
mass_fraction = X18;
break;
case 19:
mass_fraction = X19;
break;
case 20:
mass_fraction = X20;
break;
case 21:
mass_fraction = X21;
break;
case 22:
mass_fraction = X22;
break;
case 23:
mass_fraction = X23;
break;
case 24:
mass_fraction = X24;
break;
case 25:
mass_fraction = X25;
break;
case 26:
mass_fraction = X26;
break;
case 27:
mass_fraction = X27;
break;
case 28:
mass_fraction = X28;
break;
case 29:
mass_fraction = X29;
break;
case 30:
mass_fraction = X30;
break;
case 31:
mass_fraction = X31;
break;
case 32:
mass_fraction = X32;
break;
case 33:
mass_fraction = X33;
break;
case 34:
mass_fraction = X34;
break;
case 35:
mass_fraction = X35;
break;
}
return mass_fraction;
}
#endif