-
Notifications
You must be signed in to change notification settings - Fork 17
/
readinput.c
390 lines (334 loc) · 13.9 KB
/
readinput.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
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
/* ------- file: -------------------------- readinput.c -------------
Version: rh2.0, 1-D plane-parallel
Author: Han Uitenbroek ([email protected])
Last modified: Wed Jul 24 13:00:16 2013 --
-------------------------- ----------RH-- */
/* --- Reads input data for and defines keywords for 1-D
plane-parallel version -- -------------- */
#include <string.h>
#include "rh.h"
#include "atom.h"
#include "atmos.h"
#include "spectrum.h"
#include "constant.h"
#include "statistics.h"
#include "inputs.h"
#include "error.h"
/* --- Global variables -- -------------- */
extern enum Topology topology;
extern Atmosphere atmos;
extern Spectrum spectrum;
extern InputData input;
extern CommandLine commandline;
extern ProgramStats stats;
extern char messageStr[];
/* --- Function prototypes -- -------------- */
/* ------- begin -------------------------- readInput.c ------------- */
void readInput(char *input_string)
{
const char routineName[] = "readInput";
static char atom_input[MAX_VALUE_LENGTH], molecule_input[MAX_VALUE_LENGTH];
int Nkeyword;
char *keyword_string;
Keyword theKeywords[] = {
{"ATMOS_FILE", "", FALSE, KEYWORD_REQUIRED, input.atmos_input,
setcharValue},
{"ABUND_FILE", "", FALSE, KEYWORD_REQUIRED, input.abund_input,
setcharValue},
{"NRAYS", "0", FALSE, KEYWORD_OPTIONAL, &atmos.Nrays, setintValue},
{"ANGLE_SET", "NO_SET", FALSE, KEYWORD_OPTIONAL, &atmos.angleSet,
setAngleSet},
{"EDDINGTON", "FALSE", FALSE, KEYWORD_OPTIONAL, &input.Eddington,
setboolValue},
{"ATMOS_ITOP", "none", FALSE, KEYWORD_OPTIONAL, input.Itop, setcharValue},
{"WAVETABLE", "none", FALSE, KEYWORD_OPTIONAL, input.wavetable_input,
setcharValue},
{"ATOMS_FILE", "", FALSE, KEYWORD_REQUIRED, input.atoms_input,
setcharValue},
{"MOLECULES_FILE", "", FALSE, KEYWORD_REQUIRED, input.molecules_input,
setcharValue},
{"NON_ICE", "FALSE", FALSE, KEYWORD_OPTIONAL, &input.NonICE,
setboolValue},
{"N_MAX_SCATTER", "5", FALSE, KEYWORD_OPTIONAL, &input.NmaxScatter,
setintValue},
{"I_SUM", "0", FALSE, KEYWORD_REQUIRED, &input.isum, setintValue},
{"N_MAX_ITER", "", FALSE, KEYWORD_REQUIRED, &input.NmaxIter,
setintValue},
{"ITER_LIMIT", "", FALSE, KEYWORD_REQUIRED, &input.iterLimit,
setdoubleValue},
{"NG_DELAY", "0", FALSE, KEYWORD_OPTIONAL, &input.Ngdelay, setintValue},
{"NG_ORDER", "0", FALSE, KEYWORD_OPTIONAL, &input.Ngorder, setintValue},
{"NG_PERIOD", "1", FALSE, KEYWORD_OPTIONAL, &input.Ngperiod,
setintValue},
{"NG_MOLECULES", "FALSE", FALSE, KEYWORD_DEFAULT, &input.accelerate_mols,
setboolValue},
{"PRD_N_MAX_ITER", "3", FALSE, KEYWORD_OPTIONAL, &input.PRD_NmaxIter,
setintValue},
{"PRD_ITER_LIMIT", "1.0E-2", FALSE, KEYWORD_OPTIONAL, &input.PRDiterLimit,
setdoubleValue},
{"PRD_NG_DELAY", "0", FALSE, KEYWORD_OPTIONAL, &input.PRD_Ngdelay,
setintValue},
{"PRD_NG_ORDER", "0", FALSE, KEYWORD_OPTIONAL, &input.PRD_Ngorder,
setintValue},
{"PRD_NG_PERIOD", "0", FALSE, KEYWORD_OPTIONAL, &input.PRD_Ngperiod,
setintValue},
{"PRD_ANGLE_DEP", "PRD_ANGLE_INDEP", FALSE, KEYWORD_DEFAULT, &input.PRD_angle_dep,
setPRDangle},
{"PRDH_LIMIT_MEM", "FALSE", FALSE, KEYWORD_OPTIONAL, &input.prdh_limit_mem,
setboolValue},
{"BACKGR_IN_MEM", "FALSE", FALSE, KEYWORD_OPTIONAL, &input.backgr_in_mem,
setboolValue},
{"XRD", "FALSE", FALSE, KEYWORD_DEFAULT, &input.XRD, setboolValue},
{"J_FILE", "", FALSE, KEYWORD_REQUIRED, input.JFile, setcharValue},
{"BACKGROUND_FILE", "", FALSE, KEYWORD_REQUIRED, input.background_File,
setcharValue},
{"BACKGROUND_RAY_FILE", "background.ray", FALSE,
KEYWORD_OPTIONAL, input.background_ray_File,
setcharValue},
{"OLD_BACKGROUND", "FALSE", FALSE, KEYWORD_OPTIONAL,
&input.old_background, setboolValue},
{"STARTING_J", "", FALSE, KEYWORD_REQUIRED, &input.startJ,
setstartValue},
{"HYDROGEN_LTE", "FALSE", FALSE, KEYWORD_DEFAULT, &atmos.H_LTE,
setboolValue},
{"HYDROSTATIC", "FALSE", FALSE, KEYWORD_DEFAULT, &atmos.hydrostatic,
setboolValue},
{"KURUCZ_DATA", "none", FALSE, KEYWORD_OPTIONAL, &input.KuruczData,
setcharValue},
{"BARKLEM_DATA_DIR", "../../Atoms", FALSE, KEYWORD_OPTIONAL,
&input.BarklemDir, setcharValue},
{"RLK_SCATTER", "FALSE", FALSE, KEYWORD_DEFAULT, &input.rlkscatter,
setboolValue},
{"KURUCZ_PF_DATA", "../../Atoms/pf_Kurucz.input", FALSE,
KEYWORD_REQUIRED, &input.pfData, setcharValue},
{"SOLVE_NE", "NONE", FALSE, KEYWORD_DEFAULT, &input.solve_ne,
setnesolution},
{"OPACITY_FUDGE", "none", FALSE, KEYWORD_OPTIONAL, &input.fudgeData,
setcharValue},
{"METALLICITY", "0.0", FALSE, KEYWORD_DEFAULT, &input.metallicity,
setdoubleValue},
{"ATMOS_OUTPUT", "atmos.out", FALSE, KEYWORD_DEFAULT, input.atmos_output,
setcharValue},
{"GEOMETRY_OUTPUT", "geometry.out", FALSE, KEYWORD_OPTIONAL,
input.geometry_output, setcharValue},
{"SPECTRUM_OUTPUT", "spectrum.out", FALSE, KEYWORD_DEFAULT,
input.spectrum_output, setcharValue},
{"OPACITY_OUTPUT", "none", FALSE, KEYWORD_OPTIONAL, input.opac_output,
setcharValue},
{"RADRATE_OUTPUT", "none", FALSE, KEYWORD_OPTIONAL, input.radrateFile,
setcharValue},
{"COLLRATE_OUTPUT", "none", FALSE, KEYWORD_OPTIONAL, input.collrateFile,
setcharValue},
{"DAMPING_OUTPUT", "none", FALSE, KEYWORD_OPTIONAL, input.dampingFile,
setcharValue},
{"COOLING_OUTPUT", "none", FALSE, KEYWORD_OPTIONAL, input.coolingFile,
setcharValue},
{"VMICRO_CHAR", "", FALSE, KEYWORD_REQUIRED, &atmos.vmicro_char,
setdoubleValue},
{"VMACRO_TRESH", "0.1", FALSE, KEYWORD_OPTIONAL, &atmos.vmacro_tresh,
setdoubleValue},
{"LAMBDA_REF", "500.0", FALSE, KEYWORD_DEFAULT, &atmos.lambda_ref,
setdoubleValue},
{"VACUUM_TO_AIR", "0", FALSE, KEYWORD_OPTIONAL, &spectrum.vacuum_to_air,
setboolValue},
/* --- Magnetic field related inputs go here -- ------------- */
{"STOKES_INPUT", "none", FALSE, KEYWORD_OPTIONAL, input.Stokes_input,
setcharValue},
{"B_STRENGTH_CHAR", "0.0", FALSE, KEYWORD_DEFAULT, &atmos.B_char,
setdoubleValue},
{"STOKES_MODE", "NO_STOKES", FALSE, KEYWORD_OPTIONAL,
&input.StokesMode, setStokesMode},
{"MAGNETO_OPTICAL", "TRUE", FALSE, KEYWORD_DEFAULT,
&input.magneto_optical, setboolValue},
{"BACKGROUND_POLARIZATION", "FALSE", FALSE, KEYWORD_DEFAULT,
&input.backgr_pol, setboolValue},
{"XDR_ENDIAN", "TRUE", FALSE, KEYWORD_OPTIONAL,
&input.xdr_endian, setboolValue},
{"S_INTERPOLATION", "LINEAR", FALSE, KEYWORD_DEFAULT,
&input.S_interpolation, set_S_interpolation},
{"S_INTERPOLATION_STOKES", "DELO_BEZIER3", FALSE, KEYWORD_DEFAULT,
&input.S_interpolation_stokes, set_S_interpolation_stokes},
{"INTERPOLATE_3D", "LINEAR_3D", FALSE, KEYWORD_DEFAULT,
&input.interpolate_3D, setInterpolate_3D},
{"PRINT_CPU", "0", FALSE, KEYWORD_OPTIONAL, &stats.printCPU,
setboolValue},
{"N_THREADS", "0", FALSE, KEYWORD_OPTIONAL, &input.Nthreads,
setThreadValue},
{"LIMIT_MEMORY", "FALSE", FALSE, KEYWORD_DEFAULT, &input.limit_memory,
setboolValue},
{"ALLOW_PASSIVE_BB", "TRUE", FALSE, KEYWORD_DEFAULT,
&input.allow_passive_bb, setboolValue},
/* --- 1.5D version related inputs go here -- --------------- */
{"SNAPSHOT", "0", FALSE, KEYWORD_OPTIONAL, &input.p15d_nt,
setintValue},
{"X_START", "0", FALSE, KEYWORD_OPTIONAL, &input.p15d_x0,
setintValue},
{"X_END", "-1", FALSE, KEYWORD_OPTIONAL, &input.p15d_x1,
setintValue},
{"X_STEP", "1", FALSE, KEYWORD_OPTIONAL, &input.p15d_xst,
setintValue},
{"Y_START", "0", FALSE, KEYWORD_OPTIONAL, &input.p15d_y0,
setintValue},
{"Y_END", "-1", FALSE, KEYWORD_OPTIONAL, &input.p15d_y1,
setintValue},
{"Y_STEP", "1", FALSE, KEYWORD_OPTIONAL, &input.p15d_yst,
setintValue},
{"COLLRAD_SWITCH", "0.0", FALSE, KEYWORD_OPTIONAL, &input.crsw,
setdoubleValue},
{"COLLRAD_SWITCH_INI", "1.0", FALSE, KEYWORD_OPTIONAL, &input.crsw_ini,
setdoubleValue},
{"PRD_SWITCH", "0.0", FALSE, KEYWORD_OPTIONAL, &input.prdsw,
setdoubleValue},
{"N_PESC_ITER", "3", FALSE, KEYWORD_OPTIONAL, &input.NpescIter,
setintValue},
{"VTURB_MULTIPLIER", "1.0", FALSE, KEYWORD_OPTIONAL, &input.vturb_mult,
setdoubleValue},
{"VTURB_ADD", "0.0", FALSE, KEYWORD_OPTIONAL, &input.vturb_add,
setdoubleValue},
{"15D_DEPTH_REFINE", "FALSE", FALSE, KEYWORD_OPTIONAL, &input.p15d_refine,
setboolValue},
{"15D_DEPTH_ZCUT", "TRUE", FALSE, KEYWORD_OPTIONAL, &input.p15d_zcut,
setboolValue},
{"15D_RERUN", "FALSE", FALSE, KEYWORD_OPTIONAL, &input.p15d_rerun,
setboolValue},
{"15D_TMAX_CUT", "-1.0", FALSE, KEYWORD_OPTIONAL, &input.p15d_tmax,
setdoubleValue},
{"15D_WRITE_POPS", "FALSE", FALSE, KEYWORD_OPTIONAL, &input.p15d_wpop,
setboolValue},
{"15D_WRITE_RRATES", "FALSE", FALSE, KEYWORD_OPTIONAL, &input.p15d_wrates,
setboolValue},
{"15D_WRITE_CRATES", "FALSE", FALSE, KEYWORD_OPTIONAL, &input.p15d_wcrates,
setboolValue},
{"15D_WRITE_TAU1", "FALSE", FALSE, KEYWORD_OPTIONAL, &input.p15d_wtau,
setboolValue},
{"15D_WRITE_EXTRA", "TRUE", FALSE, KEYWORD_OPTIONAL, &input.p15d_wxtra,
setboolValue}
};
Nkeyword = sizeof(theKeywords) / sizeof(Keyword);
if (input_string == NULL) { /* Read file from disk */
input.keyword_file_contents = readWholeFile(commandline.keyword_input);
} else {
input.keyword_file_contents = input_string;
}
keyword_string = input.keyword_file_contents;
readValues(keyword_string, Nkeyword, theKeywords);
/* --- Perform some sanity checks -- -------------- */
switch (topology) {
case ONE_D_PLANE:
if (atmos.Nrays == 0) {
Error(ERROR_LEVEL_2, routineName,
"Must set keyword NRAYS in 1-D plane parallel geometry");
}
if (atmos.angleSet.set != NO_SET) {
Error(WARNING, routineName,
"Ignoring value of keyword ANGLE_SET in 1-D plane geometry");
}
break;
case SPHERICAL_SYMMETRIC:
if (atmos.Nrays > 0 || atmos.angleSet.set != NO_SET) {
Error(WARNING, routineName,
"Ignoring value of keywords ANGLE_SET and NRAYS in "
"spherical geometry");
}
break;
case TWO_D_PLANE:
if (input.Eddington && atmos.angleSet.set != NO_SET) {
Error(WARNING, routineName,
"Ignoring value of keywords ANGLE_SET > NO_VERTICAL when\n "
" EDDINGTON is set to TRUE\n"
" Using SET_VERTICAL with muz = 1/sqrt(3)");
atmos.angleSet.set = SET_VERTICAL;
}
break;
case THREE_D_PLANE:
if (atmos.angleSet.set == NO_SET)
Error(ERROR_LEVEL_2, routineName,
"Must set keyword ANGLE_SET in multi-D plane geometry");
if (atmos.Nrays > 0)
Error(WARNING, routineName,
"Ignoring value of keyword NRAYS in multi-D plane geometry");
break;
}
/* --- Lambda_ref should be positive since it is used for the conversion
of depth scales in routine convertScales -- ------------- */
switch (topology) {
case ONE_D_PLANE:
case SPHERICAL_SYMMETRIC:
if (atmos.lambda_ref <= 0.0)
Error(ERROR_LEVEL_2, routineName,
"Value of LAMBDA_REF should be larger than 0.0");
break;
default:
if (atmos.lambda_ref < 0.0)
Error(ERROR_LEVEL_2, routineName,
"Value of LAMBDA_REF should be larger than or equal 0.0");
break;
}
/* --- Stokes for the moment only in 1D plane -- -------------- */
if (strcmp(input.Stokes_input, "none")) {
switch (topology) {
case ONE_D_PLANE:
case TWO_D_PLANE:
if (atmos.B_char == 0.0) {
Error(WARNING, routineName,
"Parameter atmos.B_char not set or set to zero\n"
" Wavelength grids of line profiles do not take account "
" of Zeeman splitting");
}
if (input.StokesMode == NO_STOKES) {
sprintf(messageStr, "%s",
"Keyword STOKES_MODE == NO_STOKES.\n"
" Set to FIELD_FREE, POLARIZATION_FREE, or FULL_STOKES\n"
" when doing polarization calculations");
Error(ERROR_LEVEL_1, routineName, messageStr);
}
break;
case THREE_D_PLANE:
if (atmos.B_char == 0.0) {
Error(WARNING, routineName,
"Parameter atmos.B_char not set or set to zero\n"
" Wavelength grids of line profiles do not take account "
" of Zeeman splitting");
}
if (input.StokesMode == NO_STOKES) {
sprintf(messageStr, "%s",
"Keyword STOKES_MODE == NO_STOKES.\n"
" Set to FIELD_FREE, POLARIZATION_FREE, or FULL_STOKES\n"
" when doing polarization calculations");
Error(ERROR_LEVEL_1, routineName, messageStr);
}
break;
case SPHERICAL_SYMMETRIC:
Error(ERROR_LEVEL_2, routineName,
"Cannot accomodate magnetic fields in this topology");
}
} else {
if (atmos.B_char != 0.0) {
Error(WARNING, routineName,
"Ignoring value of keyword B_STRENGTH_CHAR when no "
"magnetic field is read");
}
if (input.StokesMode > NO_STOKES) {
Error(WARNING, routineName,
"Ignoring value of keyword STOKES_MODE when no "
"magnetic field is read");
}
}
/* --- Hydrostatic equilibrium only in 1-D plane parallel -- ------ */
switch (topology) {
case ONE_D_PLANE:
break;
default:
if (atmos.hydrostatic)
Error(ERROR_LEVEL_2, routineName,
"Can only perform hydrostatic equilibrium calculation"
" in 1-D Cartesian geometry");
break;
}
/* --- If called with -showkeywords commandline option -- --------- */
if (commandline.showkeywords) showValues(Nkeyword, theKeywords);
/* --- Convert to MKSA units where necessary -- ------------- */
atmos.vmicro_char *= KM_TO_M;
atmos.vmacro_tresh *= KM_TO_M;
}
/* ------- end ---------------------------- readInput.c ------------- */