Line data Source code
1 : ! ***********************************************************************
2 : !
3 : ! Copyright (C) 2010 The MESA Team
4 : !
5 : ! This program is free software: you can redistribute it and/or modify
6 : ! it under the terms of the GNU Lesser General Public License
7 : ! as published by the Free Software Foundation,
8 : ! either version 3 of the License, or (at your option) any later version.
9 : !
10 : ! This program is distributed in the hope that it will be useful,
11 : ! but WITHOUT ANY WARRANTY; without even the implied warranty of
12 : ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
13 : ! See the GNU Lesser General Public License for more details.
14 : !
15 : ! You should have received a copy of the GNU Lesser General Public License
16 : ! along with this program. If not, see <https://www.gnu.org/licenses/>.
17 : !
18 : ! ***********************************************************************
19 :
20 : module history_specs
21 :
22 : use star_private_def
23 : use star_history_def
24 : use const_def, only: dp
25 : use chem_def
26 : use num_lib, only: linear_interp, find0
27 :
28 : implicit none
29 :
30 : public ! history.f90 uses most of this
31 :
32 : ! spacing between these must be larger than max number of nuclides
33 : integer, parameter :: idel = 10000
34 :
35 : integer, parameter :: center_xa_offset = idel
36 : integer, parameter :: surface_xa_offset = center_xa_offset + idel
37 : integer, parameter :: average_xa_offset = surface_xa_offset + idel
38 : integer, parameter :: category_offset = average_xa_offset + idel
39 : integer, parameter :: total_mass_offset = category_offset + idel
40 : integer, parameter :: log_total_mass_offset = total_mass_offset + idel
41 : integer, parameter :: log_average_xa_offset = log_total_mass_offset + idel
42 : integer, parameter :: log_center_xa_offset = log_average_xa_offset + idel
43 : integer, parameter :: log_surface_xa_offset = log_center_xa_offset + idel
44 : integer, parameter :: cz_max_offset = log_surface_xa_offset + idel
45 : integer, parameter :: cz_top_max_offset = cz_max_offset + idel
46 : integer, parameter :: max_eps_nuc_offset = cz_top_max_offset + idel
47 : integer, parameter :: c_log_eps_burn_offset = max_eps_nuc_offset + idel
48 :
49 : integer, parameter :: bc_offset = c_log_eps_burn_offset + idel
50 : integer, parameter :: abs_mag_offset = bc_offset + idel
51 : integer, parameter :: lum_band_offset = abs_mag_offset + idel
52 : integer, parameter :: log_lum_band_offset = lum_band_offset + idel
53 :
54 : integer, parameter :: raw_rate_offset = log_lum_band_offset + idel
55 : integer, parameter :: screened_rate_offset = raw_rate_offset + idel
56 : integer, parameter :: eps_nuc_rate_offset = screened_rate_offset + idel
57 : integer, parameter :: eps_neu_rate_offset = eps_nuc_rate_offset + idel
58 :
59 : integer, parameter :: start_of_special_cases = eps_neu_rate_offset + idel
60 : ! mixing and burning regions must be given the largest offsets
61 : ! so they can be distinguished from the other ones
62 :
63 : integer, parameter :: mixing_offset = start_of_special_cases
64 : integer, parameter :: mix_relr_offset = mixing_offset + idel
65 : integer, parameter :: burning_offset = mix_relr_offset + idel
66 : integer, parameter :: burn_relr_offset = burning_offset + idel
67 :
68 : !integer, parameter :: next_available_offset = burning_offset + idel
69 :
70 : contains
71 :
72 1 : recursive subroutine add_history_columns( &
73 : s, level, capacity, spec, history_columns_file, report, ierr)
74 : use utils_lib
75 : use utils_def
76 : use chem_def
77 : use chem_lib
78 : use const_def, only: mesa_dir
79 : use colors_def, only: bc_total_num_colors
80 : use colors_lib
81 : type (star_info), pointer :: s
82 : integer, intent(in) :: level
83 : integer, intent(inout) :: capacity
84 : integer, pointer :: spec(:)
85 : character (len=*), intent(in) :: history_columns_file
86 : logical, intent(in) :: report
87 : integer, intent(out) :: ierr
88 :
89 : integer :: iounit, n, i, t, j, cnt, ii, nxt_spec, spec_err
90 : character (len=strlen) :: buffer, string, filename
91 : integer, parameter :: max_level = 20
92 : logical :: special_case, exists
93 : logical, parameter :: dbg = .false.
94 :
95 : include 'formats'
96 :
97 1 : if (level > max_level) then
98 0 : write(*,*) 'too many levels of nesting for log column files', level
99 0 : ierr = -1
100 0 : return
101 : end if
102 :
103 1 : ierr = 0
104 1 : spec_err = 0
105 :
106 :
107 : ! first try local directory
108 1 : filename = history_columns_file
109 :
110 1 : if(level==1) then ! First pass either the user set the file or we load the defaults
111 1 : if (len_trim(filename) == 0) filename = 'history_columns.list'
112 :
113 1 : exists=.false.
114 1 : inquire(file=filename,exist=exists)
115 :
116 1 : if(.not.exists) filename = trim(mesa_dir) // '/star/defaults/history_columns.list'
117 : else
118 : ! User had include '' in their history_columns.list file, so dont try to load the local one, jump to the defaults
119 0 : if (len_trim(filename) == 0) filename =trim(mesa_dir) // '/star/defaults/history_columns.list'
120 : end if
121 :
122 :
123 :
124 1 : open(newunit=iounit, file=trim(filename), action='read', status='old', iostat=ierr)
125 1 : if (ierr /= 0) then
126 0 : write(*,*) 'failed to open ' // trim(history_columns_file)
127 0 : return
128 : end if
129 :
130 : if (dbg) then
131 : write(*,'(A)')
132 : write(*,*) 'history_columns_file <' // trim(filename) // '>'
133 : write(*,'(A)')
134 : end if
135 :
136 1 : call count_specs
137 :
138 1 : n = 0
139 1 : i = 0
140 :
141 : do
142 :
143 61 : t = token(iounit, n, i, buffer, string)
144 61 : if (t == eof_token) exit
145 60 : if (t /= name_token) then
146 0 : call error; return
147 : end if
148 :
149 60 : if (string == 'include') then
150 0 : t = token(iounit, n, i, buffer, string)
151 0 : if (t /= string_token) then
152 0 : call error; return
153 : end if
154 0 : call add_history_columns(s, level+1, capacity, spec, string, report, spec_err)
155 0 : if (spec_err /= 0) then
156 0 : write(*,*) 'failed for included log columns list ' // trim(string)
157 0 : call error; return
158 : end if
159 0 : call count_specs
160 0 : cycle
161 : end if
162 :
163 : spec_err = 0
164 :
165 60 : if (string == 'add_center_abundances') then
166 0 : call do_abundances(center_xa_offset, spec_err)
167 0 : if (spec_err /= 0) then
168 0 : call error; return
169 : end if
170 0 : call count_specs
171 0 : cycle
172 : end if
173 :
174 60 : if (string == 'add_surface_abundances') then
175 0 : call do_abundances(surface_xa_offset, spec_err)
176 0 : if (spec_err /= 0) then
177 0 : call error; return
178 : end if
179 0 : call count_specs
180 0 : cycle
181 : end if
182 :
183 60 : if (string == 'add_average_abundances') then
184 0 : call do_abundances(average_xa_offset, spec_err)
185 0 : if (spec_err /= 0) then
186 0 : call error; return
187 : end if
188 0 : call count_specs
189 0 : cycle
190 : end if
191 :
192 60 : if (string == 'add_log_center_abundances') then
193 0 : call do_abundances(log_center_xa_offset, spec_err)
194 0 : if (spec_err /= 0) then
195 0 : call error; return
196 : end if
197 0 : call count_specs
198 0 : cycle
199 : end if
200 :
201 60 : if (string == 'add_log_surface_abundances') then
202 0 : call do_abundances(log_surface_xa_offset, spec_err)
203 0 : if (spec_err /= 0) then
204 0 : call error; return
205 : end if
206 0 : call count_specs
207 0 : cycle
208 : end if
209 :
210 60 : if (string == 'add_log_average_abundances') then
211 0 : call do_abundances(log_average_xa_offset, spec_err)
212 0 : if (spec_err /= 0) then
213 0 : call error; return
214 : end if
215 0 : call count_specs
216 0 : cycle
217 : end if
218 :
219 60 : if (string == 'add_total_mass') then
220 0 : call do_abundances(total_mass_offset, spec_err)
221 0 : if (spec_err /= 0) then
222 0 : call error; return
223 : end if
224 0 : call count_specs
225 0 : cycle
226 : end if
227 :
228 60 : if (string == 'add_log_total_mass') then
229 0 : call do_abundances(log_total_mass_offset, spec_err)
230 0 : if (spec_err /= 0) then
231 0 : call error; return
232 : end if
233 0 : call count_specs
234 0 : cycle
235 : end if
236 :
237 60 : if (string == 'add_bc') then
238 0 : call do_colors(bc_offset,'bc_', spec_err)
239 0 : if (spec_err /= 0) then
240 0 : call error; return
241 : end if
242 0 : call count_specs
243 0 : cycle
244 : end if
245 :
246 60 : if (string == 'add_abs_mag') then
247 0 : call do_colors(abs_mag_offset,'abs_mag_', spec_err)
248 0 : if (spec_err /= 0) then
249 0 : call error; return
250 : end if
251 0 : call count_specs
252 0 : cycle
253 : end if
254 :
255 60 : if (string == 'add_lum_band') then
256 0 : call do_colors(lum_band_offset,'lum_band_', spec_err)
257 0 : if (spec_err /= 0) then
258 0 : call error; return
259 : end if
260 0 : call count_specs
261 0 : cycle
262 : end if
263 :
264 60 : if (string == 'add_log_lum_band') then
265 0 : call do_colors(log_lum_band_offset,'log_lum_band_', spec_err)
266 0 : if (spec_err /= 0) then
267 0 : call error; return
268 : end if
269 0 : call count_specs
270 0 : cycle
271 : end if
272 :
273 60 : if (string == 'add_raw_rates') then
274 0 : call do_rate(raw_rate_offset,'raw_rate_', spec_err)
275 0 : if (spec_err /= 0) then
276 0 : call error; return
277 : end if
278 0 : call count_specs
279 0 : cycle
280 : end if
281 :
282 60 : if (string == 'add_screened_rates') then
283 0 : call do_rate(screened_rate_offset,'screened_rate_', spec_err)
284 0 : if (spec_err /= 0) then
285 0 : call error; return
286 : end if
287 0 : call count_specs
288 0 : cycle
289 : end if
290 :
291 60 : if (string == 'add_eps_nuc_rates') then
292 0 : call do_rate(eps_nuc_rate_offset,'eps_nuc_rate_', spec_err)
293 0 : if (spec_err /= 0) then
294 0 : call error; return
295 : end if
296 0 : call count_specs
297 0 : cycle
298 : end if
299 :
300 60 : if (string == 'add_eps_neu_rates') then
301 0 : call do_rate(eps_neu_rate_offset,'eps_neu_rate_', spec_err)
302 0 : if (spec_err /= 0) then
303 0 : call error; return
304 : end if
305 0 : call count_specs
306 0 : cycle
307 : end if
308 :
309 : nxt_spec = do1_history_spec( &
310 60 : s, iounit, t, n, i, string, buffer, special_case, report, spec_err)
311 61 : if (spec_err /= 0) then
312 0 : ierr = spec_err
313 : else
314 60 : if (.not. special_case) then
315 60 : call insert_spec(nxt_spec, string, spec_err)
316 60 : if (spec_err /= 0) then
317 0 : write(*,*) 'failed for history item ' // trim(string)
318 0 : ierr = -1; cycle
319 : end if
320 0 : else if (nxt_spec == h_mixing_regions) then
321 0 : t = token(iounit, n, i, buffer, string)
322 0 : if (t /= name_token) then
323 0 : ierr = -1; cycle
324 : end if
325 0 : read(string,fmt=*,iostat=spec_err) cnt
326 0 : if (spec_err /= 0 .or. cnt <= 0 .or. cnt > 1000) then
327 0 : write(*,*) 'bad integer count for mixing regions: ' // trim(string)
328 0 : ierr = -1; cycle
329 : end if
330 0 : do ii=1,2*cnt
331 0 : call insert_spec(mixing_offset + ii, string, spec_err)
332 0 : if (spec_err /= 0) then
333 0 : call error; return
334 : end if
335 : end do
336 0 : else if (nxt_spec == h_mix_relr_regions) then
337 0 : t = token(iounit, n, i, buffer, string)
338 0 : if (t /= name_token) then
339 0 : ierr = -1; cycle
340 : end if
341 0 : read(string,fmt=*,iostat=spec_err) cnt
342 0 : if (spec_err /= 0 .or. cnt <= 0 .or. cnt > 1000) then
343 0 : write(*,*) 'bad integer count for mix_relr regions: ' // trim(string)
344 0 : ierr = -1; cycle
345 : end if
346 0 : do ii=1,2*cnt
347 0 : call insert_spec(mix_relr_offset + ii, string, spec_err)
348 0 : if (spec_err /= 0) then
349 0 : call error; return
350 : end if
351 : end do
352 0 : else if (nxt_spec == h_burning_regions) then
353 0 : t = token(iounit, n, i, buffer, string)
354 0 : if (t /= name_token) then
355 0 : ierr = -1; cycle
356 : end if
357 0 : read(string,fmt=*,iostat=spec_err) cnt
358 0 : if (spec_err /= 0 .or. cnt <= 0 .or. cnt > 1000) then
359 0 : write(*,*) 'bad integer count following burning regions: ' // trim(string)
360 0 : ierr = -1
361 : end if
362 0 : do ii=1,2*cnt
363 0 : call insert_spec(burning_offset + ii, string, spec_err)
364 0 : if (spec_err /= 0) then
365 0 : call error; return
366 : end if
367 : end do
368 0 : else if (nxt_spec == h_burn_relr_regions) then
369 0 : t = token(iounit, n, i, buffer, string)
370 0 : if (t /= name_token) then
371 0 : ierr = -1; cycle
372 : end if
373 0 : read(string,fmt=*,iostat=spec_err) cnt
374 0 : if (spec_err /= 0 .or. cnt <= 0 .or. cnt > 1000) then
375 0 : write(*,*) 'bad integer count following burn_relr regions: ' // trim(string)
376 0 : ierr = -1
377 : end if
378 0 : do ii=1,2*cnt
379 0 : call insert_spec(burn_relr_offset + ii, string, spec_err)
380 0 : if (spec_err /= 0) then
381 0 : call error; return
382 : end if
383 : end do
384 : else
385 0 : write(*,*) 'confusion in history specs'
386 0 : ierr = -1
387 : end if
388 : end if
389 :
390 : end do
391 :
392 : if (dbg) write(*,*) 'finished ' // trim(filename)
393 :
394 1 : close(iounit)
395 :
396 2 : if (dbg) then
397 : write(*,'(A)')
398 : write(*,*) 'done add_history_columns ' // trim(filename)
399 : write(*,'(A)')
400 : end if
401 :
402 :
403 : contains
404 :
405 :
406 1 : subroutine count_specs
407 : integer :: i
408 1 : j = 1
409 1 : do i=1, capacity
410 1 : if (spec(i) == 0) then
411 1 : j = i; exit
412 : end if
413 : end do
414 1 : end subroutine count_specs
415 :
416 :
417 60 : subroutine make_room(ierr)
418 : integer, intent(out) :: ierr
419 60 : if (j < capacity) return
420 0 : capacity = 50 + (3*capacity)/2
421 0 : call realloc_integer(spec,capacity,ierr)
422 0 : spec(j+1:capacity) = 0
423 : end subroutine make_room
424 :
425 :
426 120 : subroutine insert_spec(c, name, ierr)
427 : integer, intent(in) :: c
428 : character (len=*) :: name
429 : integer, intent(out) :: ierr
430 : integer :: i
431 : include 'formats'
432 1830 : do i=1,j-1
433 1830 : if (spec(i) == c) return
434 : end do
435 60 : call make_room(ierr)
436 60 : if (ierr /= 0) return
437 60 : spec(j) = c
438 : if (dbg) write(*,2) trim(name), spec(j)
439 60 : j = j+1
440 : end subroutine insert_spec
441 :
442 :
443 0 : subroutine do_abundances(offset, ierr)
444 : integer, intent(in) :: offset
445 : integer, intent(out) :: ierr
446 : integer :: k
447 0 : ierr = 0
448 0 : do k=1,s% species
449 : call insert_spec( &
450 : offset + s% chem_id(k), &
451 0 : chem_isos% name(s% chem_id(k)), ierr)
452 0 : if (ierr /= 0) return
453 : end do
454 : end subroutine do_abundances
455 :
456 0 : subroutine do_colors(offset,prefix,ierr)
457 : integer, intent(in) :: offset
458 : character(len=*) :: prefix
459 : integer, intent(out) :: ierr
460 : integer :: k
461 0 : ierr = 0
462 0 : do k=1,bc_total_num_colors
463 : call insert_spec( &
464 0 : offset + k,trim(prefix)//trim(get_bc_name_by_id(k,ierr)), ierr)
465 0 : if (ierr /= 0) return
466 : end do
467 : end subroutine do_colors
468 :
469 0 : subroutine do_rate(offset,prefix,ierr) ! raw_rate, screened_rate, eps_nuc_rate, eps_neu_rate
470 : use rates_def, only: reaction_name
471 : use net_def, only: get_net_ptr
472 : integer, intent(in) :: offset
473 : character(len=*) :: prefix
474 : integer, intent(out) :: ierr
475 : integer :: k, ir
476 : type(net_general_info), pointer :: g=>null()
477 : ierr = 0
478 0 : call get_net_ptr(s% net_handle, g, ierr)
479 0 : if(ierr/=0) return
480 0 : do k=1,s% num_reactions
481 0 : ir = g% reaction_id(k)
482 : call insert_spec( &
483 0 : offset + k,trim(prefix)//trim(reaction_name(ir)), ierr)
484 0 : if (ierr /= 0) return
485 : end do
486 0 : end subroutine do_rate
487 :
488 :
489 0 : subroutine error
490 0 : ierr = -1
491 0 : close(iounit)
492 0 : end subroutine error
493 :
494 :
495 : end subroutine add_history_columns
496 :
497 :
498 60 : integer function do1_history_spec( &
499 : s, iounit, t, n, i, string, buffer, special_case, report, ierr) result(spec)
500 :
501 : use utils_lib
502 : use utils_def
503 : use chem_def
504 : use chem_lib
505 : use net_def
506 : integer :: iounit, t, n, i
507 :
508 : type(star_info), pointer :: s
509 : character (len=*) :: string, buffer
510 : logical, intent(out) :: special_case
511 : logical, intent(in) :: report
512 : integer, intent(out) :: ierr
513 :
514 : integer :: id
515 : type(Net_General_Info), pointer :: g
516 :
517 60 : ierr = 0
518 60 : spec = -1
519 60 : special_case = .false.
520 :
521 60 : call get_net_ptr(s% net_handle, g, ierr)
522 60 : if(ierr/=0) return
523 :
524 : ! These must come first otherwise things like center_mu will be caught by the
525 : ! center abundances check
526 60 : id = do_get_history_id(string)
527 60 : if (id > 0) then
528 49 : spec = id
529 : if (id == h_mixing_regions .or. &
530 : id == h_mix_relr_regions .or. &
531 49 : id == h_burning_regions .or. &
532 : id == h_burn_relr_regions) then
533 0 : special_case = .true.
534 : end if
535 49 : return
536 : end if
537 11 : id = rates_category_id(string)
538 11 : if (id > 0) then
539 3 : spec = category_offset + id
540 3 : return
541 : end if
542 :
543 68 : if (do1(string, 'raw_rate', raw_rate_offset, do1_rate)) then
544 0 : return
545 :
546 8 : else if (do1(string, 'screened_rate', screened_rate_offset, do1_rate)) then
547 0 : return
548 :
549 8 : else if (do1(string, 'eps_nuc_rate', eps_nuc_rate_offset, do1_rate)) then
550 0 : return
551 :
552 8 : else if (do1(string, 'eps_neu_rate', eps_neu_rate_offset, do1_rate)) then
553 0 : return
554 :
555 8 : else if (do1(string, 'abs_mag', abs_mag_offset, do1_color)) then
556 0 : return
557 :
558 8 : else if (do1(string, 'bc', bc_offset, do1_color)) then
559 0 : return
560 :
561 8 : else if (do1(string, 'lum_band', lum_band_offset, do1_color)) then
562 0 : return
563 :
564 8 : else if (do1(string, 'log_lum_band', log_lum_band_offset, do1_color)) then
565 0 : return
566 :
567 8 : else if (do1(string, 'center', center_xa_offset, do1_nuclide)) then
568 4 : return
569 :
570 4 : else if (do1(string, 'surface', surface_xa_offset, do1_nuclide)) then
571 2 : return
572 :
573 2 : else if (do1(string, 'average', average_xa_offset, do1_nuclide)) then
574 0 : return
575 :
576 2 : else if (do1(string, 'total_mass', total_mass_offset, do1_nuclide)) then
577 2 : return
578 :
579 0 : else if (do1(string, 'log_total_mass', log_total_mass_offset, do1_nuclide)) then
580 0 : return
581 :
582 0 : else if (do1(string, 'log_average', log_average_xa_offset, do1_nuclide)) then
583 0 : return
584 :
585 0 : else if (do1(string, 'log_center', log_center_xa_offset, do1_nuclide)) then
586 0 : return
587 :
588 0 : else if (do1(string, 'log_surface', log_surface_xa_offset, do1_nuclide)) then
589 0 : return
590 :
591 0 : else if (do1(string, 'max_eps_nuc_log_xa', max_eps_nuc_offset, do1_nuclide)) then
592 0 : return
593 :
594 0 : else if (do1(string, 'cz_top_log_xa', cz_top_max_offset, do1_nuclide)) then
595 0 : return
596 :
597 0 : else if (do1(string, 'cz_log_xa', cz_max_offset, do1_nuclide)) then
598 0 : return
599 :
600 0 : else if (do1(string, 'c_log_eps_burn', c_log_eps_burn_offset, do1_rates_category)) then
601 0 : return
602 :
603 : else
604 0 : if (report) write(*,*) 'bad history list name: ' // trim(string)
605 0 : ierr = -1
606 : end if
607 :
608 :
609 : contains
610 :
611 80 : logical function do1(string, name, offset, func)
612 : character(len=*) :: string, name
613 : integer :: offset, k
614 : interface
615 : subroutine func(offset)
616 : implicit none
617 : integer, intent(in) :: offset
618 : end subroutine func
619 : end interface
620 :
621 80 : if(string == name) then
622 : ! We have string value (i.e total_mass c12)
623 8 : call func(offset)
624 8 : do1 = .true.
625 72 : else if(string(1:len_trim(name)+1) == trim(name)//'_') then
626 : ! We have string_value (i.e total_mass_c12)
627 : ! Rewrite string so its in the form string value (i.e total_mass c12)
628 : ! By finding the last _ and replacing with a space
629 0 : k = index(string,'_',.true.)
630 0 : string(k:) = ' '
631 0 : buffer(k:k) = ' '
632 0 : i = len_trim(name)
633 0 : call func(offset)
634 0 : do1 = .true.
635 : else
636 : do1 = .false.
637 : end if
638 60 : end function do1
639 :
640 :
641 :
642 8 : subroutine do1_nuclide(offset)
643 : integer, intent(in) :: offset
644 8 : t = token(iounit, n, i, buffer, string)
645 8 : if (t /= name_token) then
646 0 : ierr = -1; return
647 : end if
648 8 : id = chem_get_iso_id(string)
649 8 : if (id > 0) then
650 8 : spec = offset + id
651 8 : return
652 : end if
653 0 : write(*,*) 'bad iso name: ' // trim(string)
654 0 : ierr = -1
655 : end subroutine do1_nuclide
656 :
657 :
658 0 : subroutine do1_rates_category(offset)
659 : integer, intent(in) :: offset
660 0 : t = token(iounit, n, i, buffer, string)
661 0 : if (t /= name_token) then
662 0 : ierr = -1; return
663 : end if
664 0 : id = rates_category_id(string)
665 0 : if (id > 0) then
666 0 : spec = offset + id
667 0 : return
668 : end if
669 0 : write(*,*) 'bad rates category name: ' // trim(string)
670 0 : ierr = -1
671 : end subroutine do1_rates_category
672 :
673 0 : subroutine do1_color(offset)
674 : use colors_lib, only : get_bc_id_by_name
675 : integer, intent(in) :: offset
676 0 : t = token(iounit, n, i, buffer, string)
677 0 : if (t /= name_token) then
678 0 : ierr = -1; return
679 : end if
680 0 : id = get_bc_id_by_name(string,ierr)
681 0 : if (ierr/=0) return
682 0 : if (id > 0) then
683 0 : spec = offset + id
684 0 : return
685 : end if
686 0 : write(*,*) 'bad filter name: ' // trim(string)
687 0 : ierr = -1
688 0 : end subroutine do1_color
689 :
690 0 : subroutine do1_rate(offset)
691 0 : use rates_lib, only: rates_reaction_id
692 : integer, intent(in) :: offset
693 0 : t = token(iounit, n, i, buffer, string)
694 0 : if (t /= name_token) then
695 0 : ierr = -1; return
696 : end if
697 0 : id = rates_reaction_id(string)
698 0 : id = g% net_reaction(id) ! Convert to net id not the global rate id
699 0 : if (ierr/=0) return
700 0 : if (id > 0) then
701 0 : spec = offset + id
702 0 : return
703 : end if
704 0 : write(*,*) 'bad reaction name: ' // trim(string)
705 0 : ierr = -1
706 : end subroutine do1_rate
707 :
708 :
709 : end function do1_history_spec
710 :
711 :
712 1 : subroutine set_history_columns(id, history_columns_file, report, ierr)
713 : use utils_lib, only: realloc_integer
714 : integer, intent(in) :: id
715 : character (len=*), intent(in) :: history_columns_file
716 : logical, intent(in) :: report
717 : integer, intent(out) :: ierr
718 :
719 : type (star_info), pointer :: s
720 : integer :: capacity, cnt, i
721 : logical, parameter :: dbg = .false.
722 : integer, pointer :: old_history_column_spec(:) => null()
723 : character (len=strlen) :: fname
724 : logical :: history_file_exists
725 : if (dbg) write(*,*) 'set_history_columns'
726 : ierr = 0
727 1 : call get_star_ptr(id, s, ierr)
728 1 : if (ierr /= 0) return
729 1 : old_history_column_spec => null()
730 1 : if (associated(s% history_column_spec)) old_history_column_spec => s% history_column_spec
731 : nullify(s% history_column_spec)
732 1 : capacity = 100 ! will increase if needed
733 1 : allocate(s% history_column_spec(capacity), stat=ierr)
734 1 : if (ierr /= 0) return
735 101 : s% history_column_spec(:) = 0
736 : call add_history_columns( &
737 1 : s, 1, capacity, s% history_column_spec, history_columns_file, report, ierr)
738 1 : if (ierr /= 0) then
739 0 : if (associated(old_history_column_spec)) deallocate(old_history_column_spec)
740 0 : return
741 : end if
742 : ! delete trailing 0's
743 1 : cnt = capacity+1
744 61 : do i=1, capacity
745 61 : if (s% history_column_spec(i) == 0) then
746 : cnt = i; exit
747 : end if
748 : end do
749 1 : capacity = cnt-1
750 1 : call realloc_integer(s% history_column_spec, capacity, ierr)
751 1 : if (ierr /= 0) return
752 1 : if (associated(old_history_column_spec)) then
753 : ! check that haven't changed the cols specs for an existing log file
754 0 : fname = trim(s% log_directory) // '/' // trim(s% star_history_name)
755 0 : inquire(file=trim(fname), exist=history_file_exists)
756 0 : if (history_file_exists) then
757 0 : if (capacity /= size(old_history_column_spec)) then
758 0 : ierr = -1
759 0 : write(*,*) 'new size of history col specs', capacity
760 0 : write(*,*) 'old size of history col specs', size(old_history_column_spec)
761 : else
762 0 : do i=1,capacity
763 0 : if (old_history_column_spec(i) /= s% history_column_spec(i)) then
764 0 : write(*,*) 'change in history col spec', &
765 0 : i, old_history_column_spec(i), s% history_column_spec(i)
766 0 : ierr = -1
767 : end if
768 : end do
769 : end if
770 0 : if (ierr /= 0) then
771 0 : write(*,*) 'ERROR: cannot change history columns when have an existing history file'
772 0 : write(*,*) 'please delete the history file or go back to previous history columns list'
773 : end if
774 : end if
775 0 : deallocate(old_history_column_spec)
776 0 : if (ierr /= 0) return
777 : end if
778 :
779 2 : end subroutine set_history_columns
780 :
781 : end module history_specs
|