MFIX  2016-1
stiff_chem_stats_mod.f
Go to the documentation of this file.
1 !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvv!
2 ! !
3 ! Module name: STIFF_CHEM_DEBUG !
4 ! !
5 ! Purpose: !
6 ! !
7 ! Author: J.Musser Date: !
8 ! !
9 ! Comments: !
10 ! !
11 !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^!
13 
14  PRIVATE
15 
16 ! Variable Access:
17 !---------------------------------------------------------------------//
18  PUBLIC :: failedcount
19 
20 
21 ! Subroutine Access:
22 !---------------------------------------------------------------------//
24 
25  PUBLIC :: init_stiff_chem_stats
26  PUBLIC :: update_stiff_chem_stats
27  PUBLIC :: write_stiff_chem_stats
28 
29 
30 ! Routine used to compare to values.
31  LOGICAL, external :: compare
32 
33 
34 ! Static variables/parameters.
35 !---------------------------------------------------------------------//
36 ! Frequency to report the number of steps distribution.
37  INTEGER, parameter :: reportnst_freq = 10
38 
39 ! Variables updated once each call to the stiff solver.
40 !---------------------------------------------------------------------//
41 ! Frequency to report the number of steps distribution.
42  INTEGER :: reportnst
43 
44  INTEGER :: failedcount_total
45  INTEGER :: countincpt_total
46 
47 
48 ! Variables updated every IJK loop cycle.
49 !---------------------------------------------------------------------//
50 
51 ! The minimum number of integrations needed (over all IJK)
52  INTEGER, allocatable :: minnst(:) ! local
53  INTEGER, allocatable :: minnst_all(:) ! global
54 
55 ! The maximum number of integrations needed (over all IJK)
56  INTEGER, allocatable :: maxnst(:) ! local
57  INTEGER, allocatable :: maxnst_all(:) ! global
58 
59 ! An array that stores the distrubtion of the number of steps needed
60 ! to integrate ODES.
61  INTEGER, allocatable :: countnst(:) ! local
62  INTEGER, allocatable :: countnst_all(:) ! global
63 
64 ! Number of cells that only have homogeneous chemical reactions.
65  INTEGER, allocatable :: homogns(:) ! local
66  INTEGER, allocatable :: homogns_all(:) ! global
67 
68 ! Number of cells that only have homogeneous and/or heterogeneous
69 ! chemical reactions.
70  INTEGER, allocatable :: hetrgns(:) ! local
71  INTEGER, allocatable :: hetrgns_all(:) ! global
72 
73 ! Number of cells that failed to successfully integration ODEs.
74  INTEGER, allocatable :: failedcount(:) ! local
75  INTEGER, allocatable :: failedcount_all(:) ! global
76 
77 ! Maximum number of attempts to integrate.
78  INTEGER, allocatable :: maxattempts(:) ! local
79  INTEGER, allocatable :: maxattempts_all(:) ! global
80 
81 ! Maximum number of incomplete integrations.
82  INTEGER, allocatable :: countincpt(:) ! local
83  INTEGER, allocatable :: countincpt_all(:) ! global
84 
85  DOUBLE PRECISION :: ode_starttime
86 
87  contains
88 
89 
90 !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvv!
91 ! !
92 ! Module name: ALLOCATE_STIFF_CHEM_STATS !
93 ! !
94 ! Purpose: !
95 ! !
96 ! Author: J.Musser Date: !
97 ! !
98 ! Comments: !
99 ! !
100 !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^!
101  SUBROUTINE allocate_stiff_chem_stats
103  use compar, only: mype
104  use compar, only: numpes
105  use output, only: full_log
106 
107  implicit none
108 
109  reportnst = 1
110  failedcount_total = 0
111  countincpt_total = 0
112 
113 ! Number of cells that failed to successfully integration ODEs.
114  allocate( failedcount(0:numpes-1) ); failedcount = 0 ! local
115  allocate( failedcount_all(0:numpes-1) ) ! global
116 
117  if(.NOT.full_log) return
118 
119 ! The minimum number of integrations needed (over all IJK)
120  allocate( minnst(0:numpes-1)); minnst = 0 ! local
121  allocate( minnst_all(0:numpes-1) ) ! global
122  minnst(mype) = 5000
123 
124 ! The maximum number of integrations needed (over all IJK)
125  allocate( maxnst(0:numpes-1) ); maxnst = 0 ! local
126  allocate( maxnst_all(0:numpes-1) ) ! global
127 
128 ! An array that stores the distrubtion of the number of steps needed
129 ! to integrate ODES.
130  allocate( countnst(5) ); IF(reportnst==1) countnst = 0 ! local
131  allocate( countnst_all(5) ) ! global
132 
133 ! Number of cells that only have homogeneous chemical reactions.
134  allocate( homogns(0:numpes-1) ); homogns = 0 ! local
135  allocate( homogns_all(0:numpes-1) ) ! global
136 
137 ! Number of cells that only have homogeneous and/or heterogeneous
138 ! chemical reactions.
139  allocate( hetrgns(0:numpes-1) ); hetrgns = 0; ! local
140  allocate( hetrgns_all(0:numpes-1) ) ! global
141 
142 ! Maximum number of attempts to integrate.
143  allocate( maxattempts(0:numpes-1) ); maxattempts = 0 ! local
144  allocate( maxattempts_all(0:numpes-1) ) ! global
145 
146 ! Number of cells that fail to completely integrate the time step
147 ! given the maximum number of steps.
148  allocate( countincpt(0:numpes-1) ); countincpt = 0 ! local
149  allocate( countincpt_all(0:numpes-1) ) ! global
150 
151 
152  RETURN
153  END SUBROUTINE allocate_stiff_chem_stats
154 
155 
156 !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvv!
157 ! !
158 ! Module name: INIT_ODE_STATS0 !
159 ! !
160 ! Purpose: !
161 ! !
162 ! Author: J.Musser Date: !
163 ! !
164 ! Comments: !
165 ! !
166 !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^!
167  SUBROUTINE init_stiff_chem_stats
169  use compar, only: mype
170  use compar, only: pe_io
171  use output, only: full_log
172 
173  implicit none
174 
175  if(.NOT.full_log) return
176 
177  CALL cpu_time(ode_starttime)
178 
179  hetrgns = 0
180  homogns = 0
181  failedcount = 0
182  countincpt = 0
183 
184  if(mype == pe_io) &
185  write(*,"(/3x,'Integrating stiff chemistry...')",advance="NO")
186 
187  RETURN
188  END SUBROUTINE init_stiff_chem_stats
189 
190 
191 !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvv!
192 ! !
193 ! Module name: UPDATE_STIFF_CHEM_STATS !
194 ! !
195 ! Purpose: !
196 ! !
197 ! Author: J.Musser Date: !
198 ! !
199 ! Comments: !
200 ! !
201 !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^!
202  SUBROUTINE update_stiff_chem_stats(lNEQ, lNEQ_DIMN, lNST, &
203  lode_dimn, latps, lincpt)
205  use compar, only: mype
206 
207  implicit none
208 
209 ! The number of steps needed to integrate.
210  INTEGER, intent(in) :: lNEQ_DIMN
211 
212 ! (1) :: Number of ODEs
213 ! (2) :: Fluid cell index (IJK) passed into ODEPACK
214  INTEGER, dimension(lNEQ_DIMN), intent(in) :: lNEQ
215 ! The number of steps needed to integrate.
216  INTEGER, intent(in) :: lNST
217 ! The number ODEs (maximum).
218  INTEGER, intent(in) :: lODE_DIMN
219 
220 ! The number of attempts.
221  INTEGER, intent(in) :: lAtps
222 
223 ! Flag that the integration is incomplete
224  LOGICAL, intent(in) :: lIncpt
225 
226 
227  IF(lneq(1) == lode_dimn) THEN
228  hetrgns(mype) = hetrgns(mype) + 1
229  ELSE
230  homogns(mype) = homogns(mype) + 1
231  ENDIF
232 
233  maxattempts(mype) = max(latps, maxattempts(mype))
234 
235  minnst(mype) = min(minnst(mype), lnst)
236  maxnst(mype) = max(maxnst(mype), lnst)
237 
238  IF (lnst < 10) THEN
239  countnst(1) = countnst(1) + 1
240  ELSE IF (lnst < 100) THEN
241  countnst(2) = countnst(2) + 1
242  ELSE IF (lnst < 1000) THEN
243  countnst(3) = countnst(3) + 1
244  ELSE IF (lnst < 10000) THEN
245  countnst(4) = countnst(4) + 1
246  ELSE
247  countnst(5) = countnst(5) + 1
248  ENDIF
249 
250  IF(lincpt) countincpt(mype) = countincpt(mype) + 1
251 
252  RETURN
253  END SUBROUTINE update_stiff_chem_stats
254 
255 
256 !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvv!
257 ! Module name: WRITE_ODE_STATS !
258 ! !
259 ! Purpose: !
260 ! !
261 ! Author: J.Musser Date: !
262 ! !
263 ! Comments: !
264 ! !
265 !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^!
266  SUBROUTINE write_stiff_chem_stats
268  use compar, only: mype
269  use compar, only: pe_io
270  use output, only: full_log
271 
272  use mpi_utility
273 
274  implicit none
275 
276 ! Message buffer.
277  CHARACTER(LEN=64) :: lMsg0, lMsg1
278 
279  DOUBLE PRECISION :: lODE_EndTime, lODE_RunTime
280 
281  IF(.NOT.full_log) return
282 
283 
284 ! Update screen message.
285  IF(mype == pe_io) WRITE(*,"(2x,'DONE.',/)")
286 
287 
288  CALL cpu_time(lode_endtime)
289  lode_runtime = lode_endtime - ode_starttime
290 
291 
292 ! Collect stats on min/max number of steps.
293  minnst_all = 0; CALL global_sum(minnst, minnst_all)
294  maxnst_all = 0; CALL global_sum(maxnst, maxnst_all)
295 
296 ! Collect stats on the number of cells with pure homogeneous reactions.
297  homogns_all = 0;
298  CALL global_sum(homogns, homogns_all)
299 
300 ! Collect stats on the number of cells with heterogeneous and
301 ! homogeneous reactions.
302  hetrgns_all = 0
303  CALL global_sum(hetrgns, hetrgns_all)
304 
305 ! Collect stats on the maximum number of integration attempts.
306  maxattempts_all = 0
307  CALL global_sum(maxattempts, maxattempts_all)
308 
309 ! Collect stats on the maximum number of incomplete integrations.
310  countincpt_all = 0
311  CALL global_sum(countincpt, countincpt_all)
312 
313 ! Collect stats on the number of failed integrations.
314  failedcount_all = 0
315  CALL global_sum(failedcount, failedcount_all)
316 
317 
318 ! Display stiff solver summary.
319  IF(mype == pe_io) THEN
320 
321 ! Report Min/Max steps:
322  lmsg0=''; write(lmsg0,*) minval(minnst_all)
323  lmsg1=''; write(lmsg1,*) maxval(maxnst_all)
324  write(*,1000) trim(adjustl(lmsg0)), trim(adjustl(lmsg1))
325 
326 ! Report Homogeneous/Heterogeneous reactions:
327  lmsg0=''; write(lmsg0,*) sum(homogns_all)
328  lmsg1=''; write(lmsg1,*) sum(hetrgns_all)
329  write(*,1001) trim(adjustl(lmsg0)), trim(adjustl(lmsg1))
330 
331 ! Report Max attempts:
332  lmsg0=''; write(lmsg0,*) maxval(maxattempts_all)
333  write(*,1004) trim(adjustl(lmsg0))
334 
335 ! Report incomplete integrations:
336  countincpt_total = countincpt_total + sum(countincpt_all)
337 
338  IF(countincpt_total > 0) THEN
339  lmsg0=''; write(lmsg0,*) sum(countincpt_all)
340  lmsg1=''; write(lmsg1,*) countincpt_total
341  write(*,1002) 'incomplete', trim(adjustl(lmsg0)), trim(adjustl(lmsg1))
342  ENDIF
343 
344 ! Report failed integrations:
345  failedcount_total = failedcount_total + sum(failedcount_all)
346 
347  IF(failedcount_total > 0) THEN
348  lmsg0=''; write(lmsg0,*) sum(failedcount_all)
349  lmsg1=''; write(lmsg1,*) failedcount_total
350  write(*,1002) 'failed', trim(adjustl(lmsg0)), trim(adjustl(lmsg1))
351  ENDIF
352 
353  IF(lode_runtime > 3.6d3) THEN
354  lmsg0=''; write(lmsg0,"(f8.4)") lode_runtime/3.6d3
355  lmsg1='hrs'
356  ELSEIF(lode_runtime > 6.0d1) THEN
357  lmsg0=''; write(lmsg0,"(f8.4)") lode_runtime/6.0d1
358  lmsg1='min'
359  ELSE
360  lmsg0=''; write(lmsg0,"(f8.4)") lode_runtime
361  lmsg1='sec'
362  ENDIF
363  write(*,1003) trim(adjustl(lmsg0)), trim(adjustl(lmsg1))
364 
365  ENDIF
366 
367 
368  if(reportnst == reportnst_freq) then
369 ! Collect the number of steps distributions.
370  countnst_all = 0;
371  CALL global_sum(countnst, countnst_all)
372 
373  countnst_all = int(countnst_all/reportnst_freq)
374 
375  if(mype == pe_io) then
376  write(*,"(/5x,'Average Integration Distribution:')")
377  write(*,"(7x,'NST < 10^1: ', I6)")countnst_all(1)
378  write(*,"(7x,'NST < 10^2: ', I6)")countnst_all(2)
379  write(*,"(7x,'NST < 10^3: ', I6)")countnst_all(3)
380  write(*,"(7x,'NST < 10^4: ', I6)")countnst_all(4)
381  write(*,"(7x,'NST > 10^5: ', I6)")countnst_all(5)
382  endif
383 ! Reset the reporting counter.
384  reportnst = 1
385 ! Clear out old data.
386  countnst = 0
387  countnst_all = 0
388  else
389 ! Increment the reporting counter.
390  reportnst = reportnst + 1
391  endif
392 
393  if(mype == pe_io)write(*,"(/' ')")
394 
395  RETURN
396 
397  1000 Format(5x,'Minimum/Maximum number of steps over all cells: ',a,'/',a)
398  1001 Format(5x,'Number of cells with Homogeneous/Heterogeneous reactions: ',a,'/',a)
399  1002 Format(5x,'Number of Current/Cumulative ',a,' integrations: ',a,'/',a)
400  1003 Format(5x,'CPU Time Used: ',a,' ',a)
401  1004 Format(5x,'Maximum number of integration attempts: ',a)
402 
403  END SUBROUTINE write_stiff_chem_stats
404 
405  END MODULE stiff_chem_stats
logical full_log
Definition: output_mod.f:31
subroutine, public write_stiff_chem_stats
subroutine, public update_stiff_chem_stats(lNEQ, lNEQ_DIMN, lNST, lODE_DIMN, lAtps, lIncpt)
integer numpes
Definition: compar_mod.f:24
integer, dimension(:), allocatable, public failedcount
integer pe_io
Definition: compar_mod.f:30
subroutine, public init_stiff_chem_stats
integer mype
Definition: compar_mod.f:24
subroutine, public allocate_stiff_chem_stats