MFIX  2016-1
dbg_mod.f
Go to the documentation of this file.
1 !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvv!
2 ! !
3 ! Module: DBG !
4 ! Author: J.Musser Date: 11-Nov-13 !
5 ! !
6 ! Purpose: Contains routines used for extracting the Am matrix, Bm !
7 ! vector in .csv format. Additionally, routines for extracting array !
8 ! subsets are included. !
9 ! !
10 ! Comments: !
11 ! !
12 ! > initExtract must be invoked prior to invoking any other routine !
13 ! in this file. The initialization routine needs to set several !
14 ! values required by other functions. !
15 ! !
16 ! !
17 !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^!
18  MODULE dbg
19 
20  use exit, only: mfix_exit
21  use mpi_utility
22  use param1, only: undefined_c
23  use param1, only: undefined_i
24 
25  IMPLICIT NONE
26 
27  PRIVATE
28 
29 ! Subroutine Access:
30 !---------------------------------------------------------------------//
31  PUBLIC :: initextract, matrixextract, arrayextract
32 
33 ! Interface
34 !---------------------------------------------------------------------//
35  interface arrayextract
36  module procedure arrayextract_int
37  module procedure arrayextract_dbl
38  module procedure arrayextract_log
39  end interface
40 
41 ! Module variables (private):
42 !---------------------------------------------------------------------//
43 
44 ! Coeffs for the space filling curve.
45  INTEGER :: am_c0, am_c1, am_c2
46 
47 ! The extent to extract the Am/Bm matrix/vector.
48  INTEGER :: iub, ilb
49  INTEGER :: jub, jlb
50  INTEGER :: kub, klb
51 
52 ! Logical flag to inflate A_M matrix to from sparce to full NxN.
53  LOGICAL, parameter :: inflateam = .false.
54 ! Debug Mode:
55  LOGICAL, parameter :: dbgmode = .false.
56 
57 ! Logical to verify that the initialization routine was called.
58  LOGICAL :: initnotcalled = .true.
59 
60 ! Data storage buffer size.
61  INTEGER :: dbgdimn
62 
63 ! Buffers for I/J/K indices
64  INTEGER, allocatable :: i_ofbuff(:)
65  INTEGER, allocatable :: j_ofbuff(:)
66  INTEGER, allocatable :: k_ofbuff(:)
67  INTEGER, allocatable :: ijk_buff(:)
68 
69 ! Buffers for various data types.
70 ! Buffer for generic integer array.
71  INTEGER, allocatable :: outbuff_i(:)
72 ! Buffer for generic double precision array.
73  DOUBLE PRECISION, allocatable :: outbuff_dp(:)
74 ! Buffer for Am matrix.
75  DOUBLE PRECISION, allocatable :: outam(:)
76 
77 ! Flag stating that the dgb routines are in-use
78  LOGICAL :: dbglock
79 
80 ! File Unit
81  INTEGER, parameter :: dbgunit = 9659
82 
83 ! File structure information: (optional)
84 ! Pass count if dumping at multiple locations or during iterations.
85  INTEGER :: fpass
86 ! Flag to append to previous output.
87  LOGICAL :: fapnd
88 ! Flag to write IJK values with output.
89  LOGICAL :: fwijk, pwijk
90 
91  contains
92 
93 !``````````````````````````````````````````````````````````````````````!
94 ! !
95 ! !
96 !``````````````````````````````````````````````````````````````````````!
97  INTEGER FUNCTION dbg_funijk(xxi, xxj, xxk)
98  INTEGER, intent(in) :: xxi, xxj, xxk
99 
100  dbg_funijk = xxj + am_c0 + xxi*am_c1 + xxk*am_c2
101 
102  END FUNCTION dbg_funijk
103 
104 !----------------------------------------------------------------------!
105 ! !
106 ! !
107 ! !
108 !----------------------------------------------------------------------!
109  SUBROUTINE initextract(iLow, iHgh, jLow, jHgh, kLow, kHgh)
111  use compar
112  use funits, only: dmp_log
113  use geometry, only: imin3, imax3
114  use geometry, only: jmin3, jmax3
115  use geometry, only: kmin3, kmax3, do_k
116 
117  INTEGER, optional, intent(in) :: iLow, iHgh
118  INTEGER, optional, intent(in) :: jLow, jHgh
119  INTEGER, optional, intent(in) :: kLow, kHgh
120 
121  INTEGER :: PROC
122 
123 ! Lock the dbg routines.
124  dbglock = .true.
125 
126  ilb = merge(ilow, imin3, present(ilow))
127  iub = merge(ihgh, imax3, present(ihgh))
128 
129  jlb = merge(jlow, jmin3, present(jlow))
130  jub = merge(jhgh, jmax3, present(jhgh))
131 
132  klb = merge(klow, kmin3, present(klow))
133  kub = merge(khgh, kmax3, present(khgh))
134 
135 ! Some basic logical checks:
136  IF(ilb > iub .OR. ilb<imin3 .OR. iub>imax3) THEN
137  IF(dmp_log) WRITE(*,1000)'I',ilb,imin3,iub,imax3
138  CALL mfix_exit(mype)
139  ELSEIF(jlb > jub .OR. jlb<jmin3 .OR. jub>jmax3) THEN
140  IF(dmp_log) WRITE(*,1000)'J',jlb,jmin3,jub,jmax3
141  CALL mfix_exit(mype)
142  ELSEIF(klb > kub .OR. klb<kmin3 .OR. kub>kmax3) THEN
143  IF(dmp_log) WRITE(*,1000)'K',klb,kmin3,kub,kmax3
144  CALL mfix_exit(mype)
145  ENDIF
146 
147 
148 ! Calculate the coeffs for the space filling curve.
149  am_c0 = 1 - jlb
150  am_c1 = (jub-jlb+1)
151  am_c2 = (jub-jlb+1)*(iub-ilb+1)
152  am_c0 = am_c0 - am_c1*ilb - am_c2*klb
153 
154  dbgdimn = (1+iub-ilb) * (1+jub-jlb) * (1+kub-klb)
155 
156 ! Allocate the output vector for Am
157  if(inflateam) then
158  Allocate( outam( dbgdimn) )
159  else
160  if(do_k) then
161  Allocate( outam(-3:3) )
162  else
163  Allocate( outam(-2:2) )
164  endif
165  endif
166 
167 ! Report some basic data to the screen.
168  IF(dmp_log) THEN
169  WRITE(*,"(4/,3x,'Matrix Map:')")
170  write(*,"(/,5X,'Domain Limits >')")
171  write(*,"(7X,'I: ',I4,2x,I4)") imin3, imax3
172  write(*,"(7X,'J: ',I4,2x,I4)") jmin3, jmax3
173  if(do_k) write(*,"(7X,'K: ',I4,2x,I4)") kmin3, kmax3
174 
175  write(*,"(/5x,'Local IJK Coeffs >')")
176  write(*,"(7x,'C0: ',I9)") c0
177  write(*,"(7x,'C1: ',I9)") c1
178  write(*,"(7x,'C2: ',I9)") c2
179 
180  write(*,"(/5X,'Extraction Region: >')")
181  write(*,"(7X,'I: ',I4,2x,I4)") ilb, iub
182  write(*,"(7X,'J: ',I4,2x,I4)") jlb, jub
183  if(do_k) write(*,"(7X,'K: ',I4,2x,I4)") klb, kub
184 
185  write(*,"(/5x,'dbgIJK Coeffs >')")
186  write(*,"(7x,'Am_C0: ',I9)") am_c0
187  write(*,"(7x,'Am_C1: ',I9)") am_c1
188  write(*,"(7x,'Am_C2: ',I9)") am_c2
189 
190  if(inflateam) then
191  write(*,"(/5x,'A_M is going to be inflated.')")
192  else
193  write(*,"(/5x,'A_M is NOT going to be inflated.')")
194  endif
195  ENDIF
196 
197 ! Set the flag indicating that the initialization routine was called.
198  initnotcalled = .false.
199 
200  do proc = 0, numpes-1
201  if(mype == proc) then
202  write(*,"(/3x,'Proc: ',I3)")proc
203  write(*,"(5x,'I start/end 1:',2(2x,I3))") istart1, iend1
204  write(*,"(5x,'J start/end 1:',2(2x,I3))") jstart1, jend1
205  if(do_k)write(*,"(5x,'K start/end 1:',2(2x,I3))") kstart1, kend1
206  endif
207 #ifdef MPI
208  CALL mpi_barrier(mpi_comm_world,mpierr)
209 #endif
210  enddo
211 
212 ! Unlock dbg routines.
213  dbglock = .false.
214 
215  RETURN
216 
217  1000 FORMAT(/1x,70('*')/' From: initExtract',/' Error 1000:', &
218  ' Invalid parameters. Axis: ',a1,/8x,'LB =',i4,3x,'Min2 =', &
219  i12,/8x,'UB =',i4,3x,'Max2 =',i12,/1x,70('*'),2/)
220 
221  END SUBROUTINE initextract
222 
223 !----------------------------------------------------------------------!
224 ! !
225 ! !
226 ! !
227 !----------------------------------------------------------------------!
228  SUBROUTINE matrixextract(A_m, B_m, M, VAR, PASS)
230  use compar
231  use fldvar
232 ! use geometry
233 ! use indices
234 ! use run
235 ! use sendrecv
236  use functions
237  use param, only: dimension_3, dimension_m
238  use param1, only: zero
239 
240  IMPLICIT NONE
241 
242 ! Septadiagonal matrix A_m
243  DOUBLE PRECISION, intent(in) :: A_m(dimension_3, -3:3, 0:dimension_m)
244 ! Vector b_m
245  DOUBLE PRECISION, intent(in) :: B_m(dimension_3, 0:dimension_m)
246 
247  INTEGER, intent(in) :: M
248 
249  CHARACTER(len=*), intent(in), optional :: VAR
250 
251  INTEGER, intent(in), optional :: PASS
252 
253 ! Septadiagonal matrix A_m
254  DOUBLE PRECISION :: lA_m(-3:3)
255 ! Vector b_m
256  DOUBLE PRECISION :: lB_m
257 ! Neighbor info for debugging.
258  INTEGER :: NBGHS(-3:3)
259 
260  LOGICAL :: lexist
261  CHARACTER(len=64) :: AmFName
262  CHARACTER(len=64) :: BmFName
263 
264  INTEGER, parameter :: AmUnit = 9657
265  INTEGER, parameter :: BmUnit = 9658
266 
267  INTEGER :: IJK, I, J, K, OWNER
268 
269 ! If the initialization routine was not called, flag the error and exit.
270  IF(initnotcalled)THEN
271  IF(dmp_log) THEN
272  WRITE(*,"(3/,' Fatal Error in matrixExtract!')")
273  WRITE(*,"(' The initialization routine was never called.')")
274  WRITE(*,"(' USR0 should contain a call to initExtract.')")
275  WRITE(*,"(' Forcing a hard stop.')")
276  ENDIF
277  CALL mfix_exit(mype)
278  ENDIF
279 
280  amfname=''
281  bmfname=''
282  IF(present(var) .AND. present(pass)) THEN
283  write(amfname,"('Am_',A,'_',I6.6,'.csv')")trim(adjustl(var)),pass
284  write(bmfname,"('Bm_',A,'_',I6.6,'.csv')")trim(adjustl(var)),pass
285  ELSEIF(present(var)) THEN
286  write(amfname,"('Am_',A,'.csv')")trim(adjustl(var))
287  write(bmfname,"('Bm_',A,'.csv')")trim(adjustl(var))
288  ELSEIF(present(pass)) THEN
289  write(amfname,"('Am_',I6.6,'.csv')") pass
290  write(bmfname,"('Bm_',I6.6,'.csv')") pass
291  ELSE
292  amfname='Am.csv'
293  bmfname='Bm.csv'
294  ENDIF
295 
296 ! Rank 0 opens the file.
297  IF(mype == pe_io) THEN
298  inquire(file=trim(amfname),exist=lexist)
299  IF(lexist) THEN
300  open(amunit,file=trim(amfname), status='replace',convert='big_endian')
301  ELSE
302  open(amunit,file=trim(amfname), status='new',convert='big_endian')
303  ENDIF
304 
305  inquire(file=trim(bmfname),exist=lexist)
306  IF(lexist) THEN
307  open(bmunit,file=trim(bmfname), status='replace',convert='big_endian')
308  ELSE
309  open(bmunit,file=trim(bmfname), status='new',convert='big_endian')
310  ENDIF
311  ENDIF
312 
313 ! Everyone waits for the file to open.
314 #ifdef MPI
315  CALL mpi_barrier(mpi_comm_world,mpierr)
316 #endif
317 
318  do k = klb, kub
319  do i = ilb, iub
320  do j = jlb, jub
321 
322 ! All ranks set lAm and lBm to zero.
323  owner = 0
324  nbghs = 0
325  lb_m = zero
326  la_m(-3:3) = zero
327 
328 ! Only the rank that owns IJK populates lAm and lBm.
329  if(is_on_mype_owns(i,j,k)) then
330  ijk = funijk(i,j,k)
331  la_m(-3:3) = a_m(ijk,-3:3,m)
332  lb_m = b_m(ijk,m)
333 ! Incorporate neighbor information for debugging.
334  if(dbgmode) then
335  owner = mype
336  if(do_k) nbghs(-3) = bottom_of(ijk)
337  nbghs(-2) = south_of(ijk)
338  nbghs(-1) = west_of(ijk)
339  nbghs( 0) = ijk
340  nbghs( 1) = east_of(ijk)
341  nbghs( 2) = north_of(ijk)
342  if(do_k) nbghs( 3) = top_of(ijk)
343  endif
344  endif
345 
346 ! Global sum updates lAm and lBm for all ranks.
347  CALL global_all_sum(la_m)
348  CALL global_all_sum(lb_m)
349  if(dbgmode) then
350  CALL global_all_sum(owner)
351  CALL global_all_sum(nbghs)
352  endif
353 
354  if(mype==pe_io) then
355  CALL am_to_aout(i, j, k, nbghs, owner, la_m)
356  CALL write_aoutbm(lb_m)
357  endif
358 
359 ! Everyone waits for the owner to complete this cycle.
360 #ifdef MPI
361  CALL mpi_barrier(mpi_comm_world, mpierr)
362 #endif
363 
364  enddo
365  enddo
366  enddo
367 
368 ! Rank 0 closes the file.
369  IF(mype == pe_io) then
370  close(amunit)
371  close(bmunit)
372  ENDIF
373 
374  return
375 
376  contains
377 
378 !``````````````````````````````````````````````````````````````````````!
379 ! !
380 ! !
381 !``````````````````````````````````````````````````````````````````````!
382  SUBROUTINE am_to_aout(I, J, K, lNBGHS, lOWNER, lA_m)
384  INTEGER, intent(in) :: I, J, K
385  INTEGER, intent(in) :: lNBGHS(-3:3), lOWNER
386  DOUBLE PRECISION, intent(in) :: lA_m(-3:3)
387 
388 ! Local process IJK for neighbor cell mapped to global Am output matrix.
389  INTEGER :: nIJK
390 ! Increments from I,J,K for neighbor calculations.
391  INTEGER :: ii, jj, kk
392 
393 ! Bounds for extracted section or walls.
394  INTEGER :: sMin, wMin, bMin
395  INTEGER :: eMax, nMax, tMax
396 
397 ! Initialize:
398  outam = 0.0d0
399 
400 ! Set up bounds check values.
401  wmin = ilb; emax = (1+iub-ilb)
402  smin = jlb; nmax = (1+jub-jlb)
403  bmin = klb; tmax = (1+kub-klb)
404 
405 ! Print mapping info.
406  if(dbgmode) write(*,9003) lnbghs(0), dbg_funijk(i,j,k), &
407  i, j, k, lowner
408 
409 ! Bottom Neighboring Fluid Cell
410 !---------------------------------------------------------------------//
411  if(do_k) then
412  if(k > klb) then
413  ii = i; jj = j; kk = k-1
414  nijk = merge(dbg_funijk(ii,jj,kk), -3, inflateam)
415  outam(nijk) = la_m(-3)
416 
417  if(dbgmode) write(*,9000)'Bottom of ', lnbghs(-3), &
418  nijk, ii, jj, kk, la_m(-3)
419  else
420  if(dbgmode) write(*,9001)'Bottom of ', la_m(-3)
421  endif
422  endif
423 
424 ! South Neighboring Fluid Cell
425 !---------------------------------------------------------------------//
426  if(j > jlb) then
427 
428  ii = i; jj = j-1; kk = k
429  nijk = merge(dbg_funijk(ii,jj,kk), -2, inflateam)
430  outam(nijk) = la_m(-2)
431 
432  if(dbgmode) write(*,9000)'South of ', lnbghs(-2), &
433  nijk, ii, jj, kk, la_m(-2)
434  else
435  if(dbgmode) write(*,9001)'South of ', la_m(-2)
436  endif
437 
438 ! West Neighboring Fluid Cell
439 !---------------------------------------------------------------------//
440  if(i > ilb) then
441 
442  ii = i-1; jj = j; kk = k
443  nijk = merge(dbg_funijk(ii,jj,kk), -1, inflateam)
444  outam(nijk) = la_m(-1)
445 
446  if(dbgmode) write(*,9000)'West of ', lnbghs(-1), &
447  nijk, ii, jj, kk, la_m(-1)
448  else
449  if(dbgmode) write(*,9001)'West of ', la_m(-1)
450  endif
451 
452 ! Center Coefficient
453 !---------------------------------------------------------------------//
454 
455  ii = i; jj = j; kk = k
456  nijk = merge(dbg_funijk(ii,jj,kk), 0, inflateam)
457  outam(nijk) = la_m(0)
458 
459  if(dbgmode) write(*,9000)'Cntr Coef ', lnbghs(0), &
460  nijk, ii, jj, kk, la_m(0)
461 
462 ! East Neighboring Fluid Cell
463 !---------------------------------------------------------------------//
464  if(i < iub) then
465 
466  ii = i+1; jj = j; kk = k
467  nijk = merge(dbg_funijk(ii,jj,kk), 1, inflateam)
468  outam(nijk) = la_m( 1)
469 
470  if(dbgmode) write(*,9000)'East of ', lnbghs( 1), &
471  nijk, ii, jj, kk, la_m( 1)
472  else
473  if(dbgmode) write(*,9001)'East of ', la_m( 1)
474  endif
475 
476 ! North Neighboring Fluid Cell
477 !---------------------------------------------------------------------//
478  if(j < jub) then
479 
480  ii = i; jj = j+1; kk = k
481  nijk = merge(dbg_funijk(ii,jj,kk), 2, inflateam)
482  outam(nijk) = la_m( 2)
483 
484  if(dbgmode) write(*,9000)'North of ', lnbghs(2), &
485  nijk, ii, jj, kk, la_m( 2)
486  else
487  if(dbgmode) write(*,9001)'North of ', la_m( 2)
488  endif
489 
490 ! Top Neighboring Fluid Cell
491 !---------------------------------------------------------------------//
492  if(do_k) then
493  if(k < kub) then
494 
495  ii = i; jj = j; kk = k+1
496  nijk = merge(dbg_funijk(ii,jj,kk), 3, inflateam)
497  outam(nijk) = la_m( 3)
498 
499  if(dbgmode) write(*,9000)'Top of ', lnbghs( 3), &
500  nijk, ii, jj, kk, la_m( 3)
501  else
502  if(dbgmode) write(*,9001)'Top of ', la_m( 3)
503  endif
504  endif
505 
506  return
507 
508  9000 Format(5x,a,':: ',i4,' --> ',i4,3x,'(',i3,',',i3,',',i3, &
509  ') = ',f12.4)
510 
511  9001 Format(5x,a,':: ............ OoR ............ = ',f12.4)
512 
513 
514  9003 Format(//3x,'Mapping: ',i4,' --> ',i4,3x,'(',i3,',',i3,',',i3, &
515  ')',4x,'Rank: ',i5)
516 
517  END SUBROUTINE am_to_aout
518 
519 
520 !``````````````````````````````````````````````````````````````````````!
521 ! !
522 ! !
523 !``````````````````````````````````````````````````````````````````````!
524  SUBROUTINE write_aoutbm(lB_m)
525 
526  DOUBLE PRECISION, intent(in) :: lB_m
527 
528  INTEGER :: lStart
529  INTEGER :: lEnd
530  INTEGER :: IJK
531 
532  lstart = lbound(outam,1)
533  lend = ubound(outam,1)
534 
535  if(dbgmode) then
536  do ijk = lstart, lend-1
537  write(amunit,"(F12.4,',')",advance='no')outam(ijk)
538  enddo
539  write(amunit,"(F12.4)")outam(lend)
540  write(bmunit,"(F12.4)")lb_m
541 
542  else
543  do ijk = lstart, lend-1
544  write(amunit,"(e14.6,',')",advance='no')outam(ijk)
545  enddo
546  write(amunit,"(e14.6)")outam(lend)
547  write(bmunit,"(e14.6)")lb_m
548  endif
549 
550  END SUBROUTINE write_aoutbm
551 
552  END SUBROUTINE matrixextract
553 
554 !----------------------------------------------------------------------!
555 ! !
556 ! !
557 ! !
558 !----------------------------------------------------------------------!
559  SUBROUTINE dbg_write(lMsg, Flush)
560 
561  use compar
562 ! use sendrecv
563 
564  implicit none
565 
566  CHARACTER(len=*), intent(in) :: lMsg
567 
568  LOGICAL, optional, intent(in) :: FLUSH
569  LOGICAL :: lFlush
570 
571  if(dbgmode)then
572  lflush = merge(FLUSH, .false., present(flush))
573 
574  if(mype == pe_io) then
575  if(lflush) write(*,"(' ')")
576  write(*,"(3x,A)") trim(lmsg)
577  endif
578 #ifdef MPI
579  CALL mpi_barrier(mpi_comm_world,mpierr)
580 #endif
581  endif
582 
583  RETURN
584  END SUBROUTINE dbg_write
585 
586 !----------------------------------------------------------------------!
587 ! !
588 ! !
589 ! !
590 !----------------------------------------------------------------------!
591  SUBROUTINE arrayextract_init(vName)
593  use compar, only: mype
594  use funits, only: dmp_log
595 
596  implicit none
597 
598 ! Variable named used to create file name.
599  CHARACTER(len=*), intent(in) :: vName
600 
601 ! Variables for opening files.
602  LOGICAL :: lexist
603  CHARACTER(len=64) :: VarFName
604 ! Integer Error Flag.
605  INTEGER :: iErr
606 
607 ! If the initialization routine was not called, flag the error and exit.
608  IF(initnotcalled)THEN
609  IF(dmp_log) WRITE(*,1000)
610  CALL mfix_exit(mype)
611  ENDIF
612 
613 ! Verify that the files are unlocked.
614  IF(dbglock)THEN
615  IF(dmp_log) WRITE(*,1001)
616  CALL mfix_exit(mype)
617  ELSE
618  dbglock = .true.
619  ENDIF
620 
621 ! Allocate indice arrays.
622  IF(allocated(i_ofbuff))THEN
623  IF(dmp_log) WRITE(*,1002) 'i_ofBuff'
624  CALL mfix_exit(mype)
625  ELSE
626  allocate( i_ofbuff(dbgdimn) )
627  i_ofbuff = 0
628  ENDIF
629 
630  IF(allocated(j_ofbuff))THEN
631  IF(dmp_log) WRITE(*,1002) 'j_ofBuff'
632  CALL mfix_exit(mype)
633  ELSE
634  allocate( j_ofbuff(dbgdimn) )
635  j_ofbuff = 0
636  ENDIF
637 
638  IF(allocated(k_ofbuff))THEN
639  IF(dmp_log) WRITE(*,1002) 'k_ofBuff'
640  CALL mfix_exit(mype)
641  ELSE
642  allocate( k_ofbuff(dbgdimn) )
643  k_ofbuff = 0
644  ENDIF
645 
646  IF(allocated(ijk_buff))THEN
647  IF(dmp_log) WRITE(*,1002) 'ijk_Buff'
648  CALL mfix_exit(mype)
649  ELSE
650  allocate( ijk_buff(dbgdimn) )
651  ijk_buff = 0
652  ENDIF
653 
654 ! Construct the file name.
655  varfname=''
656  IF(fapnd) THEN
657  write(varfname,"(A,'.csv')") &
658  trim(adjustl(vname))
659  ELSEIF(fpass /= undefined_i) THEN
660  write(varfname,"(A,'_',I6.6,'.csv')") &
661  trim(adjustl(vname)), fpass
662  ELSE
663  write(varfname,"(A,'.csv')") &
664  trim(adjustl(vname))
665  ENDIF
666 
667 ! Open the file -- Rank 0 ONLY.
668  inquire(file=trim(varfname), exist=lexist)
669  IF(mype == pe_io) THEN
670  IF(lexist) THEN
671  IF(fapnd) THEN
672  open(dbgunit,file=trim(varfname), &
673  status='old', position='append', iostat=ierr,convert='big_endian')
674  ELSE
675  open(dbgunit,file=trim(varfname), &
676  status='replace', iostat=ierr,convert='big_endian')
677  ENDIF
678  ELSE
679  open(dbgunit,file=trim(varfname), &
680  status='new', iostat=ierr,convert='big_endian')
681  ENDIF
682  ENDIF
683  CALL bcast(ierr, pe_io)
684  IF(ierr /= 0)THEN
685  IF(mype == pe_io) write(*,1003) trim(varfname)
686  CALL mfix_exit(mype)
687  ENDIF
688 
689 ! Set printIJK Flag. Only print the IJK values if this is the first time
690 ! for an append file.
691  pwijk = merge(.false., fwijk, fwijk.AND.lexist)
692 
693  RETURN
694 
695  1000 FORMAT(//1x,70('*')/' From: dbg_mod -> arrayExtract_init',/, &
696  ' Error 1000: The initialization routine was never called.', &
697  ' Include the',/' following in USR0: CALL initExtract.',2/, &
698  ' These arguments are used to specify a domain subset to', &
699  ' extract. If',/' not defined, the entire domain (MIN3/MAX3)',&
700  ' is extracted.',2/,' Optional arguments:',/, &
701  3x,'iLow - lower I index; iHgh - Upper I index (X-axis)',/, &
702  3x,'jLow - lower J index; jHgh - Upper J index (Y-axis)',/, &
703  3x,'kLow - lower K index; kHgh - Upper K index (Z-axis)',/ &
704  1x,70('*'),2/)
705 
706  1001 FORMAT(//1x,70('*')/' From: dbg_mod -> arrayExtract_init',/, &
707  ' Error 1001: dgbLock is set. Something must have failed.',/ &
708  1x,70('*'),2/)
709 
710  1002 FORMAT(//1x,70('*')/' From: dbg_mod -> arrayExtract_init',/, &
711  ' Error 1002: Buffer already allocated: ',a,/1x,70('*'),2/)
712 
713  1003 FORMAT(//1x,70('*')/' From: dbg_mod -> arrayExtract_init',/, &
714  ' Error 1002: Error opening file: ',a,/1x,70('*'),2/)
715 
716  END SUBROUTINE arrayextract_init
717 
718 !----------------------------------------------------------------------!
719 ! !
720 ! !
721 ! !
722 !----------------------------------------------------------------------!
723  SUBROUTINE arrayextract_finl
725  implicit none
726 
727  IF(allocated(i_ofbuff)) deallocate(i_ofbuff)
728  IF(allocated(j_ofbuff)) deallocate(j_ofbuff)
729  IF(allocated(k_ofbuff)) deallocate(k_ofbuff)
730  IF(allocated(ijk_buff)) deallocate(ijk_buff)
731 
732  IF(mype == pe_io) close(dbgunit)
733 
734  dbglock = .false.
735 
736  RETURN
737  END SUBROUTINE arrayextract_finl
738 
739 !----------------------------------------------------------------------!
740 ! !
741 ! !
742 ! !
743 !----------------------------------------------------------------------!
744  SUBROUTINE arrayextract_prnt(dType)
745 
746  implicit none
747 
748 ! Data type
749  CHARACTER(len=3), intent(in) :: dType
750 
751 ! Loop counter
752  INTEGER :: IJK
753 
754 ! Rank 0 writes and closes the file.
755  IF(mype /= pe_io) RETURN
756 
757  IF(fapnd) THEN
758 ! Print header info.
759  IF(pwijk) THEN
760  CALL write_index(i_ofbuff)
761  CALL write_index(j_ofbuff)
762  if(do_k) CALL write_index(k_ofbuff)
763  CALL write_index(ijk_buff)
764  ENDIF
765 
766 ! Store the pass count.
767  IF(fpass /= undefined_i) &
768  WRITE(dbgunit,"(I14,',')",advance='no') fpass
769 ! Output the data.
770  SELECT CASE(dtype)
771  CASE('INT'); CALL write_int
772  CASE('DBL'); CALL write_dbl
773  CASE('LOG'); CALL write_log
774  END SELECT
775 
776  ELSE
777 
778  IF(fpass /= undefined_i) &
779  WRITE(dbgunit,"(2x,'Pass: ',I8)") fpass
780 
781  DO ijk=1, dbgdimn
782 ! Output the IJK info.
783  IF(pwijk)THEN
784  WRITE(dbgunit,"(4(I14,','))",advance='no') ijk, &
785  ijk_buff(ijk), i_ofbuff(ijk), j_ofbuff(ijk)
786  if(do_k)WRITE(dbgunit,"(I14,',')",advance='no') &
787  k_ofbuff(ijk)
788  ENDIF
789 ! Write the formatted output.
790  SELECT CASE(dtype)
791  CASE('INT'); CALL write_int(ijk)
792  CASE('DBL'); CALL write_dbl(ijk)
793  CASE('LOG'); CALL write_log(ijk)
794  END SELECT
795  ENDDO
796  ENDIF
797 
798  RETURN
799  END SUBROUTINE arrayextract_prnt
800 
801 !----------------------------------------------------------------------!
802 ! !
803 ! !
804 ! !
805 !----------------------------------------------------------------------!
806  SUBROUTINE write_index(intArray)
807 
808  implicit none
809 
810  INTEGER, intent(in) :: intArray(dbgdimn)
811 ! Looped IJK
812  INTEGER :: IJK
813 
814 ! Loop through array entries writing them in one long continuous line.
815  IF(fpass /= undefined_i) WRITE(dbgunit,"(14X,',')",advance='no')
816  DO ijk = 1, dbgdimn-1
817  WRITE(dbgunit,"(I14,',')",advance='no')intarray(ijk)
818  ENDDO
819  WRITE(dbgunit,"(I14)")intarray(dbgdimn)
820 
821  RETURN
822  END SUBROUTINE write_index
823 
824 !----------------------------------------------------------------------!
825 ! !
826 ! !
827 ! !
828 !----------------------------------------------------------------------!
829  SUBROUTINE write_int(pIJK)
830 
831  implicit none
832 
833 ! Dummy arguments.
834  INTEGER, optional, intent(in) :: pIJK
835 ! Looped IJK
836  INTEGER :: IJK
837 
838  IF(present(pijk)) THEN
839 ! Write the entry and return.
840  WRITE(dbgunit,"(I14)") outbuff_i(pijk)
841 
842  ELSE
843 ! Loop through array entries writing them in one long continuous line.
844  DO ijk = 1, dbgdimn-1
845  WRITE(dbgunit,"(I14,',')",advance='no')outbuff_i(ijk)
846  ENDDO
847  WRITE(dbgunit,"(I14)")outbuff_i(dbgdimn)
848  ENDIF
849 
850  RETURN
851  END SUBROUTINE write_int
852 
853 !----------------------------------------------------------------------!
854 ! !
855 ! !
856 ! !
857 !----------------------------------------------------------------------!
858  SUBROUTINE write_dbl(pIJK)
859 
860  implicit none
861 
862 ! Dummy arguments.
863  INTEGER, optional, intent(in) :: pIJK
864 ! Looped IJK
865  INTEGER :: IJK
866 
867  IF(present(pijk)) THEN
868 ! Write the entry and return.
869  WRITE(dbgunit,"(E14.6)") outbuff_dp(pijk)
870  ELSE
871 ! Loop through array entries writing them in one long continuous line.
872  DO ijk = 1, dbgdimn-1
873  WRITE(dbgunit,"(E14.6,',')",advance='no')outbuff_dp(ijk)
874  ENDDO
875  WRITE(dbgunit,"(E14.6)")outbuff_dp(dbgdimn)
876  ENDIF
877 
878  RETURN
879  END SUBROUTINE write_dbl
880 
881 !----------------------------------------------------------------------!
882 ! !
883 ! !
884 ! !
885 !----------------------------------------------------------------------!
886  SUBROUTINE write_log(pIJK)
887 
888  implicit none
889 
890 ! Dummy arguments.
891  INTEGER, optional, intent(in) :: pIJK
892 ! Looped IJK
893  INTEGER :: IJK
894  LOGICAL :: INTtoLOG
895 
896  IF(present(pijk)) THEN
897 ! Write the entry and return.
898  inttolog= (outbuff_i(pijk) .eq. 1)
899  WRITE(dbgunit,"(L14)") inttolog
900 
901  ELSE
902 ! Loop through array entries writing them in one long continuous line.
903  DO ijk = 1, dbgdimn-1
904  inttolog=(outbuff_i(ijk) .eq. 1)
905  WRITE(dbgunit,"(L14,',')",advance='no') inttolog
906  ENDDO
907  inttolog = (outbuff_i(dbgdimn) .eq. 1)
908  WRITE(dbgunit,"(L14)") inttolog
909 
910  ENDIF
911 
912  RETURN
913  END SUBROUTINE write_log
914 
915 !----------------------------------------------------------------------!
916 ! !
917 ! !
918 ! !
919 !----------------------------------------------------------------------!
920  SUBROUTINE arrayextract_int(Array, VAR, PASS, APND, withIJK)
921 
922  use compar
923  use fldvar
924  use geometry
925  use indices
926  use mflux
927  use parallel
928  use param
929  use param1
930  use physprop
931  use run
932  use sendrecv
933  USE mpi_utility
934  USE functions
935 
936  implicit none
937 
938 ! Array to be extracted.
939  INTEGER, intent(in) :: Array(dimension_3)
940 ! Variable named used to create file name.
941  CHARACTER(len=*), intent(in) :: VAR
942 ! Pass count if dumping at multiple locations or during iterations.
943  INTEGER, intent(in), optional :: PASS
944 ! Flag to append to previous output.
945  LOGICAL, intent(in), optional :: APND
946 ! Flag to write IJK values with output.
947  LOGICAL, intent(in), optional :: withIJK
948 ! Loop counters:
949  INTEGER :: I, J, K, IJK, dbgIJK
950 ! Debugging message.
951  CHARACTER(len=64) :: MSG
952 
953  msg='Entered arrayExtract_int'
954  if(dbgmode) CALL dbg_write(trim(msg),flush=.true.)
955 
956 ! Set local flags based on optional arguments.
957  fpass = merge(pass, undefined_i, present(pass))
958  fapnd = merge(apnd, .false., present(apnd))
959  fwijk = merge(withijk, .false., present(withijk))
960 
961  msg=' > Calling arrayExtract_INIT'
962  if(dbgmode) CALL dbg_write(trim(msg))
963  CALL arrayextract_init(var)
964 
965  msg=' > Allocating outBuff_i'
966  if(dbgmode) CALL dbg_write(trim(msg))
967  allocate( outbuff_i(dbgdimn) ); outbuff_i = 0
968 
969  msg=' > Extracting array data.'
970  if(dbgmode) CALL dbg_write(trim(msg))
971 
972  do k = klb, kub
973  do i = ilb, iub
974  do j = jlb, jub
975  if(is_on_mype_owns(i,j,k)) then
976  IF (dead_cell_at(i,j,k)) cycle ! skip dead cells
977  ijk = funijk(i,j,k)
978  dbgijk = dbg_funijk(i,j,k)
979 ! Only the rank that owns IJK populates lAm and lBm.
980  outbuff_i(dbgijk) = array(ijk)
981 ! Store indicies.
982  i_ofbuff(dbgijk) = i
983  j_ofbuff(dbgijk) = j
984  k_ofbuff(dbgijk) = k
985  ijk_buff(dbgijk) = ijk
986  endif
987  enddo
988  enddo
989  enddo
990 
991  msg=' > Collecting array data.'
992  if(dbgmode) CALL dbg_write(trim(msg))
993 ! Global sum updates lAm and lBm for all ranks.
994  CALL global_all_sum(outbuff_i)
995  if(pwijk) then
996  CALL global_all_sum(i_ofbuff)
997  CALL global_all_sum(j_ofbuff)
998  CALL global_all_sum(ijk_buff)
999  if(do_k)CALL global_all_sum(k_ofbuff)
1000  endif
1001 
1002 ! Write the data.
1003  msg=' > Calling arrayExtract_prnt.'
1004  if(dbgmode) CALL dbg_write(trim(msg))
1005  CALL arrayextract_prnt('INT')
1006 
1007 ! Clean up the used memory.
1008  msg=' > Calling arrayExtract_finl.'
1009  if(dbgmode) CALL dbg_write(trim(msg))
1010 
1011  if(allocated(outbuff_i)) deallocate(outbuff_i)
1012  CALL arrayextract_finl
1013 
1014  msg='Leaving arrayExtract_int.'
1015  if(dbgmode) CALL dbg_write(trim(msg))
1016 
1017  RETURN
1018  END SUBROUTINE arrayextract_int
1019 
1020 !----------------------------------------------------------------------!
1021 ! !
1022 ! !
1023 ! !
1024 !----------------------------------------------------------------------!
1025  SUBROUTINE arrayextract_dbl(Array, VAR, PASS, APND, withIJK)
1027  use compar
1028  use fldvar
1029  use geometry
1030  use indices
1031  use mflux
1032  use parallel
1033  use param
1034  use param1
1035  use physprop
1036  use run
1037  use sendrecv
1038  USE mpi_utility
1039  USE functions
1040 
1041  implicit none
1042 
1043 ! Array to be extracted.
1044  DOUBLE PRECISION, intent(in) :: Array(dimension_3)
1045 ! Variable named used to create file name.
1046  CHARACTER(len=*), intent(in) :: VAR
1047 ! Pass count if dumping at multiple locations or during iterations.
1048  INTEGER, intent(in), optional :: PASS
1049 ! Flag to append to previous output.
1050  LOGICAL, intent(in), optional :: APND
1051 ! Flag to write IJK values with output.
1052  LOGICAL, intent(in), optional :: withIJK
1053 
1054 ! Loop indices:
1055  INTEGER :: I, J, K, IJK, dbgIJK
1056 ! Debugging message.
1057  CHARACTER(len=64) :: MSG
1058 
1059  msg='Entered arrayExtract_dbl'
1060  CALL dbg_write(trim(msg), flush=.true.)
1061 
1062 ! Set local flags based on optional arguments.
1063  fpass = merge(pass, undefined_i, present(pass))
1064  fapnd = merge(apnd, .false., present(apnd))
1065  fwijk = merge(withijk, .false., present(withijk))
1066 
1067  msg=' > Calling arrayExtract_INIT'
1068  CALL dbg_write(trim(msg))
1069  CALL arrayextract_init(var)
1070 
1071 ! Allocate the global storage array.
1072  allocate( outbuff_dp(dbgdimn) ); outbuff_dp = 0
1073 
1074  msg=' > Extracting array data.'
1075  CALL dbg_write(trim(msg))
1076 
1077 ! Extract the data... one cell at a time.
1078  do k = klb, kub
1079  do i = ilb, iub
1080  do j = jlb, jub
1081  if(is_on_mype_owns(i,j,k)) then
1082  IF (dead_cell_at(i,j,k)) cycle ! skip dead cells
1083  ijk = funijk(i,j,k)
1084  dbgijk = dbg_funijk(i,j,k)
1085 ! Only the rank that owns IJK populates lAm and lBm.
1086  outbuff_dp(dbgijk) = array(ijk)
1087 ! Store indicies.
1088  i_ofbuff(dbgijk) = i
1089  j_ofbuff(dbgijk) = j
1090  k_ofbuff(dbgijk) = k
1091  ijk_buff(dbgijk) = ijk
1092  endif
1093  enddo
1094  enddo
1095  enddo
1096 
1097 ! Global sum updates lAm and lBm for all ranks.
1098  msg=' > Collecting array data.'
1099  CALL dbg_write(trim(msg))
1100 
1101  CALL global_all_sum(outbuff_dp)
1102  if(pwijk) then
1103  CALL global_all_sum(i_ofbuff)
1104  CALL global_all_sum(j_ofbuff)
1105  CALL global_all_sum(ijk_buff)
1106  if(do_k)CALL global_all_sum(k_ofbuff)
1107  endif
1108 
1109 ! Write the data.
1110  msg=' > Calling arrayExtract_prnt.'
1111  CALL dbg_write(trim(msg))
1112  CALL arrayextract_prnt('DBL')
1113 
1114 ! Clean up the used memory.
1115  msg=' > Calling arrayExtract_finl.'
1116  CALL dbg_write(trim(msg))
1117  if(allocated(outbuff_dp)) deallocate(outbuff_dp)
1118  CALL arrayextract_finl
1119 
1120  msg='Leaving arrayExtract_dbl.'
1121  CALL dbg_write(trim(msg))
1122 
1123  RETURN
1124  END SUBROUTINE arrayextract_dbl
1125 
1126 !----------------------------------------------------------------------!
1127 ! !
1128 ! !
1129 ! !
1130 !----------------------------------------------------------------------!
1131  SUBROUTINE arrayextract_log(Array, VAR, PASS, APND, withIJK)
1133  use compar
1134  use fldvar
1135  use geometry
1136  use indices
1137  use mflux
1138  use parallel
1139  use param
1140  use param1
1141  use physprop
1142  use run
1143  use sendrecv
1144  USE mpi_utility
1145  USE functions
1146 
1147  implicit none
1148 
1149 ! Array to be extracted.
1150  LOGICAL, intent(in) :: Array(dimension_3)
1151 ! Variable named used to create file name.
1152  CHARACTER(len=*), intent(in) :: VAR
1153 ! Pass count if dumping at multiple locations or during iterations.
1154  INTEGER, intent(in), optional :: PASS
1155 ! Flag to append to previous output.
1156  LOGICAL, intent(in), optional :: APND
1157 ! Flag to write IJK values with output.
1158  LOGICAL, intent(in), optional :: withIJK
1159 ! Loop counters:
1160  INTEGER :: I, J, K, IJK, dbgIJK
1161 ! Debugging message.
1162  CHARACTER(len=64) :: MSG
1163 
1164  msg='Entered arrayExtract_log'
1165  if(dbgmode) CALL dbg_write(trim(msg),flush=.true.)
1166 
1167 ! Set local flags based on optional arguments.
1168  fpass = merge(pass, undefined_i, present(pass))
1169  fapnd = merge(apnd, .false., present(apnd))
1170  fwijk = merge(withijk, .false., present(withijk))
1171 
1172  msg=' > Calling arrayExtract_INIT'
1173  if(dbgmode) CALL dbg_write(trim(msg))
1174  CALL arrayextract_init(var)
1175 
1176  msg=' > Allocating outBuff_i'
1177  if(dbgmode) CALL dbg_write(trim(msg))
1178  allocate( outbuff_i(dbgdimn) ); outbuff_i = 0
1179 
1180  msg=' > Extracting array data.'
1181  if(dbgmode) CALL dbg_write(trim(msg))
1182 
1183  do k = klb, kub
1184  do i = ilb, iub
1185  do j = jlb, jub
1186  if(is_on_mype_owns(i,j,k)) then
1187  IF (dead_cell_at(i,j,k)) cycle ! skip dead cells
1188  ijk = funijk(i,j,k)
1189  dbgijk = dbg_funijk(i,j,k)
1190 ! Only the rank that owns IJK populates lAm and lBm.
1191 ! Convert the logical to interger.
1192  outbuff_i(dbgijk) = merge(1,0,array(ijk))
1193 ! Store indicies.
1194  i_ofbuff(dbgijk) = i
1195  j_ofbuff(dbgijk) = j
1196  k_ofbuff(dbgijk) = k
1197  ijk_buff(dbgijk) = ijk
1198  endif
1199  enddo
1200  enddo
1201  enddo
1202 
1203  msg=' > Collecting array data.'
1204  if(dbgmode) CALL dbg_write(trim(msg))
1205 ! Global sum updates lAm and lBm for all ranks.
1206  CALL global_all_sum(outbuff_i)
1207  if(pwijk) then
1208  CALL global_all_sum(i_ofbuff)
1209  CALL global_all_sum(j_ofbuff)
1210  CALL global_all_sum(ijk_buff)
1211  if(do_k)CALL global_all_sum(k_ofbuff)
1212  endif
1213 
1214 ! Write the data.
1215  msg=' > Calling arrayExtract_prnt.'
1216  if(dbgmode) CALL dbg_write(trim(msg))
1217  CALL arrayextract_prnt('LOG')
1218 
1219 ! Clean up the used memory.
1220  msg=' > Calling arrayExtract_finl.'
1221  if(dbgmode) CALL dbg_write(trim(msg))
1222 
1223  if(allocated(outbuff_i)) deallocate(outbuff_i)
1224  CALL arrayextract_finl
1225 
1226  msg='Leaving arrayExtract_log.'
1227  if(dbgmode) CALL dbg_write(trim(msg))
1228 
1229  RETURN
1230  END SUBROUTINE arrayextract_log
1231 
1232  END MODULE dbg
logical dmp_log
Definition: funits_mod.f:6
integer c0
Definition: compar_mod.f:104
integer kend1
Definition: compar_mod.f:80
integer istart1
Definition: compar_mod.f:80
integer iend1
Definition: compar_mod.f:80
integer imax3
Definition: geometry_mod.f:91
integer dimension_3
Definition: param_mod.f:11
Definition: dbg_mod.f:18
subroutine arrayextract_init(vName)
Definition: dbg_mod.f:592
integer mpierr
Definition: compar_mod.f:27
subroutine, public matrixextract(A_m, B_m, M, VAR, PASS)
Definition: dbg_mod.f:229
integer imin3
Definition: geometry_mod.f:90
subroutine, public initextract(iLow, iHgh, jLow, jHgh, kLow, kHgh)
Definition: dbg_mod.f:110
integer kstart1
Definition: compar_mod.f:80
integer numpes
Definition: compar_mod.f:24
subroutine am_to_aout(I, J, K, lNBGHS, lOWNER, lA_m)
Definition: dbg_mod.f:383
integer pe_io
Definition: compar_mod.f:30
integer c2
Definition: compar_mod.f:104
subroutine mfix_exit(myID, normal_termination)
Definition: exit.f:5
logical, dimension(:,:,:), allocatable dead_cell_at
Definition: compar_mod.f:127
Definition: exit.f:2
integer jmax3
Definition: geometry_mod.f:91
subroutine arrayextract_log(Array, VAR, PASS, APND, withIJK)
Definition: dbg_mod.f:1132
Definition: run_mod.f:13
subroutine arrayextract_finl
Definition: dbg_mod.f:724
Definition: param_mod.f:2
integer jmin3
Definition: geometry_mod.f:90
integer kmax3
Definition: geometry_mod.f:91
logical do_k
Definition: geometry_mod.f:30
integer mype
Definition: compar_mod.f:24
integer, parameter undefined_i
Definition: param1_mod.f:19
integer kmin3
Definition: geometry_mod.f:90
integer dimension_m
Definition: param_mod.f:18
double precision, parameter zero
Definition: param1_mod.f:27
integer c1
Definition: compar_mod.f:104
integer jend1
Definition: compar_mod.f:80
integer jstart1
Definition: compar_mod.f:80
subroutine arrayextract_dbl(Array, VAR, PASS, APND, withIJK)
Definition: dbg_mod.f:1026
character, parameter undefined_c
Definition: param1_mod.f:20