Source

dengo / templates / rates_and_rate_tables.c.template

The default branch has multiple heads

Full commit
  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
{% block solver_comment_header %}
/* THIS FILE HAS BEEN AUTO-GENERATED.  DO NOT EDIT. */

/* This is C++ code to read HDF5 files for
   reaction rates, cooling rates, and initial
   conditions for the chemical network defined
   by the user.  In addition, this contains
   code for calculating temperature from the
   gas energy and computing the RHS and the
   Jacobian of the system of equations which
   will be fed into the solver.
*/
{% endblock %}

{% block includes %}

/* stdlib, hdf5, local includes */

#include "time.h"
#include "sys/time.h"
#include "stdlib.h"
#include "math.h"
#include "hdf5.h"
#include "hdf5_hl.h"
#include "stdio.h"

#define MAX_NCELLS 1024

typedef int(*rhs_f)(double *, double *, int, int, void *);
typedef int(*jac_f)(double *, double *, int, int, void *);
int BE_chem_solve(rhs_f f, jac_f J, double *u, double dt, double *rtol, 
                  double *atol, int nstrip, int nchem, void *sdata);

{% endblock %} {# includes #}

{% block data_structs_for_rates %}

typedef struct {{ solver_name }}_data_ {

    /* All of the network bins will be the same width */
    double dbin;
    double idbin;
    double bounds[2];
    int nbins;

    double Ts[MAX_NCELLS];
    double Tdef[MAX_NCELLS]; /* t1, t2, tdef */
    double dT[MAX_NCELLS]; /* t1, t2, tdef */
    double logTs[MAX_NCELLS];
    double dTs_{{ network.energy_term.name }}[MAX_NCELLS];
    /* Now we do all of our cooling and chemical tables */
    {%- for name, rate in network.reactions | dictsort %}
    double r_{{name}}[{{ network.T | length }}];
    double rs_{{name}}[MAX_NCELLS];
    double drs_{{name}}[MAX_NCELLS];
    {% endfor %}
    {%- for name, rate in network.cooling_actions | dictsort %}
    {%- for name2 in rate.tables | sort %}
    double c_{{name}}_{{name2}}[{{ network.T | length }}];
    double cs_{{name}}_{{name2}}[MAX_NCELLS];
    double dcs_{{name}}_{{name2}}[MAX_NCELLS];
    {% endfor %}
    {% endfor %}
    int bin_id[MAX_NCELLS];
    int ncells;
} {{ solver_name }}_data;

void {{ solver_name }}_read_rate_tables({{solver_name}}_data*);
void {{ solver_name }}_read_cooling_tables({{solver_name}}_data*);
{% endblock %} {# data_struct_for_rates #}

{% block function_declarations %}

int calculate_jacobian_{{solver_name}}(double *input, double *Joutput,
        int nstrip, int nchem, void *sdata);
int calculate_rhs_{{solver_name}}(double *input, double *rhs, int nstrip,
                  int nchem, void *sdata);

{% endblock %} {# function_declarations #}

{% block entry_point %}
int main(int argc, char** argv)
{
    {{ solver_name }}_data *data;

    data = ({{solver_name}}_data *) malloc(sizeof({{solver_name}}_data));

    data->bounds[0] = {{ network.T_bounds[0] }};
    data->bounds[1] = {{ network.T_bounds[1] }};
    data->nbins = {{ network.T | length }};
    data->dbin = (log(data->bounds[1]) - log(data->bounds[0])) / data->nbins;
    data->idbin = 1.0L / data->dbin;
    
    {{ solver_name }}_read_rate_tables(data);
    fprintf(stderr, "Successfully read in rate tables.\n");

    {{ solver_name }}_read_cooling_tables(data);
    fprintf(stderr, "Successfully read in cooling rate tables.\n");

    /* Initial conditions */

    hid_t file_id = H5Fopen("{{ solver_name }}_initial_conditions.h5", H5F_ACC_RDONLY, H5P_DEFAULT);
    if (file_id < 0) {fprintf(stderr, "Failed to open "
        "{{ solver_name }}_initial_conditions.h5 so dying.\n");
        return(1);}

    /* Allocate the correct number of cells */
    hsize_t dims; /* We have flat versus number of species */

    /* Check gas energy to get the number of cells */
    fprintf(stderr, "Getting dimensionality from {{ network.energy_term.name }}:\n");
    herr_t status = H5LTget_dataset_info(file_id, "/{{ network.energy_term.name }}", &dims, NULL, NULL);
    fprintf(stderr, "  ncells = % 3i\n", (int) dims);
    data->ncells = dims;

    int N = {{network.required_species | length}};

    double *atol, *rtol;
    atol = (double *) alloca(N * dims * sizeof(double));
    rtol = (double *) alloca(N * dims * sizeof(double));

    double *tics = (double *) alloca(dims * sizeof(double));
    double *ics = (double *) alloca(dims * N * sizeof(double));
    int i = 0, j;
    {% for s in network.required_species | sort %}
    fprintf(stderr, "Reading I.C. for /{{ s.name }}\n");
    H5LTread_dataset_double(file_id, "/{{ s.name }}", tics);
    for (j = 0; j < dims; j++) {
        ics[j * N + i] = tics[j];
        atol[j] = tics[j] * 1e-7;
        rtol[j] = 1e-7;
        if(j==0) {
            fprintf(stderr, "{{s.name}}[0] = %0.3g, atol => %0.5g\n",
                    tics[j], atol[j]);
        }
    }
    i++;
    {% endfor %}

    H5Fclose(file_id);

    float dt = 1e8;
    rhs_f f = calculate_rhs_{{solver_name}};
    jac_f jf = calculate_jacobian_{{solver_name}};
    int rv = BE_chem_solve(f, jf, ics, dt, rtol, atol, dims, N, (void *) data);
    fprintf(stderr, "Return value: %i\n", rv);
}
{% endblock %} {# entry_point #}

{% block read_tables %}
void {{ solver_name }}_read_rate_tables({{solver_name}}_data *data)
{
    herr_t status;
    hid_t file_id = H5Fopen("{{ solver_name }}_rate_tables.h5", H5F_ACC_RDONLY, H5P_DEFAULT);
    /* Allocate the correct number of rate tables */
    hsize_t dims;

    {%- for name, rate in network.reactions | dictsort %}
    H5LTread_dataset_double(file_id, "/{{ name }}", data->r_{{name}});
    {%- endfor %}

    H5Fclose(file_id);
}

void {{ solver_name }}_read_cooling_tables({{solver_name}}_data *data)
{

    herr_t status;
    hid_t file_id = H5Fopen("{{ solver_name }}_cooling_tables.h5", H5F_ACC_RDONLY, H5P_DEFAULT);
    /* Allocate the correct number of rate tables */
    hsize_t dims;

    {%- for name, rate in network.cooling_actions | dictsort %}
    {%- for name2 in rate.tables | sort %}
    H5LTread_dataset_double(file_id, "/{{name}}_{{name2}}",
                            data->c_{{name}}_{{name2}});
    {%- endfor %}
    {%- endfor %}

    H5Fclose(file_id);
}

{% endblock %} {# read_tables #}

{% block calculate_temperature %}
void {{ solver_name }}_calculate_temperature({{ solver_name }}_data *data,
                        double *input, int nstrip, int nchem)
{
    int i, j;
    double density;
    double kb = 1.3806504e-16; // Boltzmann constant [erg/K]
    double mh = 1.67e-24;
    double gamma = 5.e0/3.e0;
    {% if 'H2I' in network.species_list() %}
    double gammaH2 = 7.e0/5.e0; // Should be a function of temperature
    	   	     		// this is a temporary solution
    {% endif %}
    double temperature;
    /* Calculate total density */

    {%- for species in network.required_species | sort %}
    double {{species.name}};
    {%- endfor %}

    for (i = 0; i<nstrip; i++) {
        j = i * nchem;
    {%- for species in network.required_species | sort %}
        {{species.name}} = input[j];
        j++;
    {% endfor %}
        density = {{network.print_mass_density()}};
        data->Ts[i] = {{ network.temperature_calculation() }};
        if (data->Ts[i] < data->bounds[0]) data->Ts[i] = data->bounds[0];
        if (data->Ts[i] > data->bounds[1]) data->Ts[i] = data->bounds[1];
        data->logTs[i] = log(data->Ts[i]);
	data->dTs_{{ network.energy_term.name }}[i] = 
        {{ network.temperature_calculation(derivative=True) }};
        fprintf(stderr, "T[%d] = %0.5g, LogT[%d] = %0.5g, dTdee = %0.5g\n", i,
        data->Ts[i], i, data->logTs[i], data->dTs_{{network.energy_term.name}}[i]);
    }
         
}
{% endblock %} {# calculate_temperature #}

{% block interpolate_rates %}
/*
   This setup may be different than the user may anticipate, as a result
   of the lockstep timestep we use for a pencil beam through the grid.
   As such, it accepts the number of things to interpolate and makes
   assumptions about the sizes of the rates.
*/

/* This also requires no templating other than for the solver name...*/
void {{ solver_name }}_interpolate_rates({{ solver_name }}_data *data,
                    int nstrip)
{
    int i, ri, bin_id;
    static int ncalls = 0;
    double dd, bv, dy, lb, ub;
    double t1, t2, tdef;
    lb = log(data->bounds[0]);
    ub = log(data->bounds[1]);
    fprintf(stderr, "lb = %0.5g, ub = %0.5g\n", lb, ub);
    for (i = 0; i < nstrip; i++) {
        data->bin_id[i] = bin_id = (int) (data->idbin * (data->logTs[i] - lb));
        t1 = (lb + (bin_id - 1) * data->dbin);
        t2 = (lb + (bin_id    ) * data->dbin);
        data->Tdef[i] = (data->logTs[i] - t1)/(t2 - t1);
        data->dT[i] = (t2 - t1);
        fprintf(stderr, "INTERP: %d, bin_id = %d, dT = %0.5g, T = %0.5g, logT = %0.5g\n",
                i, data->bin_id[i], data->dT[i], data->Ts[i], data->logTs[i]);
    }

    {%- for name, rate in network.reactions | dictsort %}
    for (i = 0; i < nstrip; i++) {
        bin_id = data->bin_id[i];
        data->rs_{{name}}[i] = data->r_{{name}}[bin_id] +
            data->Tdef[i] * (data->r_{{name}}[bin_id+1] - data->r_{{name}}[bin_id]);
        data->drs_{{name}}[i] = (data->r_{{name}}[bin_id+1] - data->r_{{name}}[bin_id]);
        data->drs_{{name}}[i] /= data->dT[i];
    }
    {% endfor %}
    {%- for name, rate in network.cooling_actions | dictsort %}
    {%- for name2 in rate.tables | sort %}
    for (i = 0; i < nstrip; i++) {
        bin_id = data->bin_id[i];
        data->cs_{{name}}_{{name2}}[i] = data->c_{{name}}_{{name2}}[bin_id] +
            data->Tdef[i] * (data->c_{{name}}_{{name2}}[bin_id+1] - data->c_{{name}}_{{name2}}[bin_id]);
        data->dcs_{{name}}_{{name2}}[i] = (data->c_{{name}}_{{name2}}[bin_id+1] - data->c_{{name}}_{{name2}}[bin_id]);;
        data->dcs_{{name}}_{{name2}}[i] /= data->dT[i];
    }
    {% endfor %}
    {% endfor %}

}
{% endblock %} {# interpolate_rates #}

{% block calculate_rhs %}

int calculate_rhs_{{solver_name}}(double *input, double *rhs, int nstrip,
                  int nchem, void *sdata)
{
    /* We iterate over all of the rates */
    /* Calculate temperature first */
    {{solver_name}}_data *data = ({{solver_name}}_data*)sdata;
    int i, j;
    {{solver_name}}_calculate_temperature(data, input, nstrip, nchem);

    {{solver_name}}_interpolate_rates(data, nstrip);

    /* Now we set up some temporaries */
    {%- for name, rate in network.reactions | dictsort %}
    double *{{name}} = data->rs_{{name}};
    {%- endfor %}

    {%- for name, rate in network.cooling_actions | dictsort %}
    {%- for name2 in rate.tables | sort %}
    double *{{name}}_{{name2}} = data->cs_{{name}}_{{name2}};
    {%- endfor %}
    {%- endfor %}
    {%- for species in network.required_species | sort %}
    double {{species.name}};
    {%- endfor %}

    double total;
    for (i = 0; i<nstrip; i++) {
        j = i * nchem;
        total = 0.0;
    {%- for species in network.required_species | sort %}
        {{species.name}} = input[j];
        j++;
    {% endfor %}
        j = i * nchem;
    {%- for species in network.required_species | sort %}
        // 
        // Species: {{species.name}}
        // 
        {{ network.print_ccode(species, assign_to = "rhs[j]") }}
        fprintf(stderr, "rhs[{{species.name}}][%d] = %0.5g\n", i, rhs[j]);
        {% if species.name != "ge" %}
        total += rhs[j];
        {% else %}
        rhs[j] = 0.0;
        {% endif %}
        j++;
    {% endfor %}
        fprintf(stderr, "rhs_total[%d] = %0.5g\n", i, total);
    }
    return 0;
}

{% endblock %}

{% block calculate_jacobian %}
int calculate_jacobian_{{solver_name}}(double *input, double *Joutput,
        int nstrip, int nchem, void *sdata)
{
    /* We iterate over all of the rates */
    /* Calculate temperature first */
    {{solver_name}}_data *data = ({{solver_name}}_data*)sdata;

    int i, j;
    {{solver_name}}_calculate_temperature(data, input, nstrip, nchem);

    {{solver_name}}_interpolate_rates(data, nstrip);

    /* Now we set up some temporaries */
    double *T{{ network.energy_term.name }} = data->dTs_{{ network.energy_term.name }};
    {%- for name, rate in network.reactions | dictsort %}
    double *{{name}} = data->rs_{{name}};
    double *r{{name}} = data->drs_{{name}};
    {%- endfor %}

    {%- for name, rate in network.cooling_actions | dictsort %}
    {%- for name2 in rate.tables | sort %}
    double *{{name}}_{{name2}} = data->cs_{{name}}_{{name2}};
    double *r{{name}}_{{name2}} = data->dcs_{{name}}_{{name2}};
    {%- endfor %}
    {%- endfor %}
    {%- for species in network.required_species | sort %}
    double {{species.name}};
    {%- endfor %}

    for (i = 0; i<nstrip; i++) {
        j = i * nchem;
    {%- for species in network.required_species | sort %}
        {{species.name}} = input[j];
        j++;
    {% endfor %}
        j = i * nchem * nchem;
    {%- for s1 in network.required_species | sort %}
        // 
        // Species: {{s1.name }}
        // 
        {%- for s2 in network.required_species | sort %}
            // {{s1.name}} by {{s2.name}}
            {{ network.print_jacobian_component(s1, s2, assign_to="Joutput[j]") }}
	        {% if s2.name == 'ge' %}
            Joutput[j] *= T{{ network.energy_term.name }}[i];
            {% endif %}
            if(Joutput[j] != -100.0) {
	    	{% if s2.name != 'geaaaa' %}
                fprintf(stderr, "Joutput[{{s1.name}}][{{s2.name}}][%d] = %0.5g\n",
                    i, Joutput[j]);
		{% else %}
		fprintf(stderr, "Joutput[{{s1.name}}][{{s2.name}}][%d] = %0.5g, T{{ network.energy_term.name }}[%d] = %0.5g\n",
                    i, Joutput[j], i, Tge[i]);
		{% endif %}
            }
            j++;
        {%- endfor %}
    {% endfor %}
    }
    int nzero = 0;
    double mi = 1e30;
    double ma = -1e30;
    double avg = 0.0;
    for (i = 0; i < nchem * nchem * nstrip; i++) {
        if (Joutput[i] == 0.0) nzero++;
        if (Joutput[i] < mi) mi = Joutput[i];
        if (Joutput[i] > ma) ma = Joutput[i];
        avg += Joutput[i];
    }
    fprintf(stderr, "JOUTPUTSTATS %d %0.5g %0.5g %d %0.5g\n",
        nchem*nchem*nstrip, mi, ma, nzero, avg / (nchem*nchem*nstrip));
    return 0;
    
}
{% endblock %}