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 reaclib_print
21 : use rates_def
22 :
23 : implicit none
24 :
25 : contains
26 :
27 :
28 0 : subroutine write_reaction_data(unitno, rates, ierr)
29 : integer, intent(in) :: unitno
30 : type(reaction_data), intent(in) :: rates
31 : integer, intent(out) :: ierr
32 : integer :: i
33 : character(len=*), parameter :: fmt = &
34 : '(a36, i3, tr1, 6(i5, tr1), a1, tr1, f5.1, tr1, /, 4es14.6, /, 3es14.6, /, 2es14.6, tr1, i2, tr1, 3(/, 8(f12.9, tr1)))'
35 0 : write(unitno, *) rates% nreactions
36 0 : write(unitno, *) rates% nchapters_included
37 0 : write(unitno, '(i6)') rates% chapters_present(:)
38 0 : write(unitno, '(2i6, tr1)') (rates% bookmarks(:, i), i=1, nchapters)
39 0 : do i = 1, rates% nreactions
40 : write(unitno, fmt=fmt, iostat=ierr) &
41 0 : rates% reaction_handle(i), rates% chapter(i), rates% pspecies(:, i), rates% reaction_flag(i), &
42 0 : 1.0d0/rates% weight(i), rates% coefficients(:, i), rates% inverse_coefficients(:, i), rates% inverse_exp(i), &
43 0 : rates% inverse_part(:, i)
44 : end do
45 0 : end subroutine write_reaction_data
46 :
47 :
48 0 : subroutine pretty_print_reactions(unitno, rates, nuclides, ierr)
49 : integer, intent(in) :: unitno
50 : type(reaction_data), intent(in) :: rates
51 : type(nuclide_data), intent(in) :: nuclides
52 : integer, intent(out) :: ierr
53 : integer :: i, pass
54 : logical :: reverse
55 : character (len=100) :: str
56 0 : write(unitno,'(/,a,/)') 'forward reactions'
57 0 : do pass = 1, 2
58 0 : str = ''
59 0 : reverse = (pass == 2)
60 0 : do i = 1, rates% nreactions
61 0 : call do_pretty_print_reaction(unitno, i, rates, nuclides, reverse, str, ierr)
62 0 : if (ierr /= 0) exit
63 : end do
64 0 : if (pass == 1) write(unitno,'(/,a,/)') 'reverse reactions'
65 : end do
66 0 : end subroutine pretty_print_reactions
67 :
68 :
69 0 : subroutine do_pretty_print_reaction(unitno, i, rates, nuclides, reverse, str, ierr)
70 : integer, intent(in) :: unitno, i
71 : type(reaction_data), intent(in) :: rates
72 : type(nuclide_data), intent(in) :: nuclides
73 : logical, intent(in) :: reverse
74 : character (len=100), intent(inout) :: str
75 : integer, intent(out) :: ierr
76 : character (len=100) :: str_nxt
77 0 : ierr = 0
78 0 : str_nxt = ''
79 0 : if (reverse .and. i <= rates% num_from_weaklib) return
80 0 : select case(rates% chapter(i))
81 : case(r_one_one)
82 0 : call write_n_to_m(1,1)
83 : case(r_one_two)
84 0 : call write_n_to_m(1,2)
85 : case(r_one_three)
86 0 : call write_n_to_m(1,3)
87 : case(r_two_one)
88 0 : call write_n_to_m(2,1)
89 : case(r_two_two)
90 0 : call write_n_to_m(2,2)
91 : case(r_two_three)
92 0 : call write_n_to_m(2,3)
93 : case(r_two_four)
94 0 : call write_n_to_m(2,4)
95 : case(r_three_one)
96 0 : call write_n_to_m(3,1)
97 : case(r_three_two)
98 0 : call write_n_to_m(3,2)
99 : case(r_four_two)
100 0 : call write_n_to_m(4,2)
101 : case(r_one_four)
102 0 : call write_n_to_m(1,4)
103 : end select
104 0 : if (trim(str_nxt) /= trim(str) .and. len_trim(str_nxt) > 0) then
105 0 : if (i <= rates% num_from_weaklib) str_nxt = trim(str_nxt) // ' weaklib'
106 0 : write(unitno,fmt='(a)') trim(str_nxt)
107 0 : str = str_nxt
108 : end if
109 :
110 : contains
111 :
112 0 : subroutine write_n_to_m(n,m)
113 : integer, intent(in) :: n, m
114 : integer :: j
115 0 : if (.not. reverse) then
116 0 : do j=1,n
117 0 : str_nxt = trim(str_nxt) // ' ' // nuclides% name(rates% pspecies(j, i))
118 : end do
119 0 : str_nxt = trim(str_nxt) // " -> "
120 0 : do j=1,m
121 0 : str_nxt = trim(str_nxt) // ' ' // nuclides% name(rates% pspecies(n+j, i))
122 : end do
123 0 : else if (rates% reaction_flag(i) /= 'w' .and. rates% reaction_flag(i) /= 'e') then
124 0 : do j=1,m
125 0 : str_nxt = trim(str_nxt) // ' ' // nuclides% name(rates% pspecies(n+j, i))
126 : end do
127 0 : str_nxt = trim(str_nxt) // " -> "
128 0 : do j=1,n
129 0 : str_nxt = trim(str_nxt) // ' ' // nuclides% name(rates% pspecies(j, i))
130 : end do
131 : end if
132 0 : end subroutine write_n_to_m
133 :
134 :
135 : end subroutine do_pretty_print_reaction
136 :
137 :
138 0 : subroutine print_short_format_reactions(unitno, rates, nuclides, ierr)
139 : integer, intent(in) :: unitno
140 : type(reaction_data), intent(in) :: rates
141 : type(nuclide_data), intent(in) :: nuclides
142 : integer, intent(out) :: ierr
143 : integer :: i, pass
144 : logical :: reverse
145 : character (len=100) :: str
146 0 : write(unitno,'(/,a,/)') 'forward reactions'
147 0 : do pass = 1, 2
148 0 : str = ''
149 0 : reverse = (pass == 2)
150 0 : do i = 1, rates% nreactions
151 0 : call do_print_short_format_reaction(unitno, i, rates, nuclides, reverse, str, ierr)
152 0 : if (ierr /= 0) exit
153 : end do
154 0 : if (pass == 1) write(unitno,'(/,a,/)') 'reverse reactions'
155 : end do
156 0 : end subroutine print_short_format_reactions
157 :
158 :
159 0 : subroutine do_print_short_format_reaction(unitno, i, rates, nuclides, reverse, str, ierr)
160 : use reaclib_support, only: get1_reaction_handle
161 : integer, intent(in) :: unitno, i
162 : type(reaction_data), intent(in) :: rates
163 : type(nuclide_data), intent(in) :: nuclides
164 : logical, intent(in) :: reverse
165 : character (len=100), intent(inout) :: str
166 : integer, intent(out) :: ierr
167 :
168 : character (len=100) :: str_nxt
169 : integer :: chapter, num_in, num_out
170 :
171 0 : ierr = 0
172 0 : str = ''
173 :
174 0 : if (reverse .and. &
175 : (rates% reaction_flag(i) == 'w' .or. rates% reaction_flag(i) == 'e')) return
176 :
177 0 : chapter = rates% chapter(i)
178 0 : num_in = Nin(chapter)
179 0 : num_out = Nout(chapter)
180 :
181 : call get1_reaction_handle( &
182 : num_in, num_out, rates% pspecies(:,i), nuclides, reverse, &
183 0 : rates% reaction_flag(i), str_nxt)
184 0 : if (trim(str_nxt) /= trim(str) .and. len_trim(str_nxt) > 0) then
185 0 : write(unitno,fmt='(a)') trim(str_nxt)
186 0 : str = str_nxt
187 : end if
188 :
189 0 : end subroutine do_print_short_format_reaction
190 :
191 :
192 : end module reaclib_print
|