Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
RFLU_ModProbes.F90
Go to the documentation of this file.
1 ! *********************************************************************
2 ! * Rocstar Simulation Suite *
3 ! * Copyright@2015, Illinois Rocstar LLC. All rights reserved. *
4 ! * *
5 ! * Illinois Rocstar LLC *
6 ! * Champaign, IL *
7 ! * www.illinoisrocstar.com *
8 ! * sales@illinoisrocstar.com *
9 ! * *
10 ! * License: See LICENSE file in top level of distribution package or *
11 ! * http://opensource.org/licenses/NCSA *
12 ! *********************************************************************
13 ! *********************************************************************
14 ! * THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, *
15 ! * EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES *
16 ! * OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND *
17 ! * NONINFRINGEMENT. IN NO EVENT SHALL THE CONTRIBUTORS OR *
18 ! * COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER *
19 ! * LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, *
20 ! * Arising FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE *
21 ! * USE OR OTHER DEALINGS WITH THE SOFTWARE. *
22 ! *********************************************************************
23 ! ******************************************************************************
24 !
25 ! Purpose: Collection of routines related to probes.
26 !
27 ! Description: None.
28 !
29 ! Notes: None.
30 !
31 ! ******************************************************************************
32 !
33 ! $Id: RFLU_ModProbes.F90,v 1.6 2008/12/06 08:44:23 mtcampbe Exp $
34 !
35 ! Copyright: (c) 2005 by the University of Illinois
36 !
37 ! ******************************************************************************
38 
40 
41  USE modglobal, ONLY: t_global
42  USE moddatatypes
43  USE modparameters
44  USE moderror
45  USE moddatastruct, ONLY: t_region
46  USE modgrid, ONLY: t_grid
47  USE modmpi
48 
49  IMPLICIT NONE
50 
51  PRIVATE
52  PUBLIC :: rflu_closeprobefiles, &
57 
58 ! ******************************************************************************
59 ! Declarations and definitions
60 ! ******************************************************************************
61 
62  CHARACTER(CHRLEN) :: &
63  RCSIdentString = '$RCSfile: RFLU_ModProbes.F90,v $ $Revision: 1.6 $'
64 
65 ! *****************************************************************************
66 ! Routines
67 ! *****************************************************************************
68 
69  CONTAINS
70 
71 
72 
73 
74 
75 
76 
77 ! ******************************************************************************
78 !
79 ! Purpose: Close probe files.
80 !
81 ! Description: None.
82 !
83 ! Input:
84 ! pRegion Pointer to region
85 !
86 ! Output: None.
87 !
88 ! Notes: None.
89 !
90 ! ******************************************************************************
91 
92 SUBROUTINE rflu_closeprobefiles(pRegion)
93 
94  IMPLICIT NONE
95 
96 ! ******************************************************************************
97 ! Declarations
98 ! ******************************************************************************
99 
100 ! ==============================================================================
101 ! Arguments
102 ! ==============================================================================
103 
104  TYPE (t_region), POINTER :: pregion
105 
106 ! ******************************************************************************
107 ! Locals
108 ! ******************************************************************************
109 
110  INTEGER :: errorflag,iprobe
111  TYPE(t_global), POINTER :: global
112 
113 ! ******************************************************************************
114 ! Start
115 ! ******************************************************************************
116 
117  global => pregion%global
118 
119  CALL registerfunction(global,'RFLU_CloseProbeFiles',&
120  'RFLU_ModProbes.F90')
121 
122 ! ******************************************************************************
123 ! Loop over probes
124 ! ******************************************************************************
125 
126  DO iprobe = 1,global%nProbes
127 
128 ! ==============================================================================
129 ! Check if probe region on current processor
130 ! ==============================================================================
131 
132  IF ( global%probePos(iprobe,probe_region) == pregion%iRegionGlobal ) THEN
133  CLOSE(if_probe+iprobe-1,iostat=errorflag)
134  global%error = errorflag
135  IF ( global%error /= err_none ) THEN
136  CALL errorstop(global,err_file_close,__line__)
137  END IF ! global%error
138  END IF ! global
139  END DO ! iProbe
140 
141 ! ******************************************************************************
142 ! End
143 ! ******************************************************************************
144 
145  CALL deregisterfunction(global)
146 
147 END SUBROUTINE rflu_closeprobefiles
148 
149 
150 
151 
152 
153 
154 
155 
156 
157 ! ******************************************************************************
158 !
159 ! Purpose: Determine whether to write probe data to files.
160 !
161 ! Description: None.
162 !
163 ! Input:
164 ! global Pointer to global data
165 !
166 ! Output:
167 ! RFLU_DecideWrite = .TRUE. If should write probe data to files
168 ! RFLU_DecideWrite = .FALSE. If should not write probe data to files
169 !
170 ! Notes: None.
171 !
172 ! ******************************************************************************
173 
174 LOGICAL FUNCTION rflu_decidewriteprobes(global)
175 
176  IMPLICIT NONE
177 
178 ! ******************************************************************************
179 ! Declarations and definitions
180 ! ******************************************************************************
181 
182 ! ==============================================================================
183 ! Arguments
184 ! ==============================================================================
185 
186  TYPE(t_global), POINTER :: global
187 
188 ! ==============================================================================
189 ! Locals
190 ! ==============================================================================
191 
192  LOGICAL :: logical1,logical2,logical3
193 
194 ! ******************************************************************************
195 ! Start
196 ! ******************************************************************************
197 
198  CALL registerfunction(global,'RFLU_DecideWriteProbes',&
199  'RFLU_ModProbes.F90')
200 
201 ! ******************************************************************************
202 ! Initialize
203 ! ******************************************************************************
204 
205  rflu_decidewriteprobes = .false.
206 
207 ! ******************************************************************************
208 ! Determine whether should print to screen
209 ! ******************************************************************************
210 
211 ! ==============================================================================
212 ! Unsteady flow
213 ! ==============================================================================
214 
215  IF ( global%flowType == flow_unsteady ) THEN
216  logical1 = abs(global%timeSinceProbe-global%probeSaveTime) &
217  < 0.1_rfreal*global%dtMin
218  logical2 = (global%timeSinceProbe > global%probeSaveTime)
219  logical3 = (global%iterSinceRestart == 1)
220 
221  IF ( logical1 .OR. logical2 .OR. logical3 ) THEN
222  rflu_decidewriteprobes = .true.
223  END IF ! logical1
224 
225 ! ==============================================================================
226 ! Steady flow
227 ! ==============================================================================
228 
229  ELSE
230  rflu_decidewriteprobes = (mod(global%currentIter,global%probeSaveIter) == 0)
231 
232  IF ( global%currentIter == 1 ) THEN
233  rflu_decidewriteprobes = .true.
234  END IF ! global%currentIter
235  END IF ! global%flowType
236 
237 ! ******************************************************************************
238 ! End
239 ! ******************************************************************************
240 
241  CALL deregisterfunction(global)
242 
243 END FUNCTION rflu_decidewriteprobes
244 
245 
246 
247 
248 
249 
250 
251 
252 ! ******************************************************************************
253 !
254 ! Purpose: Determine cell centroids which are closest to probe coordinates.
255 !
256 ! Description: None.
257 !
258 ! Input:
259 ! pRegion Pointer to region
260 !
261 ! Output: None.
262 !
263 ! Notes:
264 ! 1. It is not enough to simply get closest cell from Octree, because
265 ! proximity to cell centroid does not guarantee that the specified
266 ! location is actually inside the cell.
267 ! 2. If the probe location is very close to a vertex, a very large number
268 ! of cells may have to be searched for tetrahedral grids...
269 !
270 ! ******************************************************************************
271 
272 SUBROUTINE rflu_findprobecells(pRegion)
273 
276  USE rflu_modoctree
277 
278  IMPLICIT NONE
279 
280 ! ******************************************************************************
281 ! Declarations
282 ! ******************************************************************************
283 
284 ! ==============================================================================
285 ! Arguments
286 ! ===============================================================================
287 
288  TYPE(t_region), POINTER :: pregion
289 
290 ! ==============================================================================
291 ! Locals
292 ! ==============================================================================
293 
294  INTEGER :: ccsize,errorflag,icg,iprobe
295  INTEGER, DIMENSION(:), ALLOCATABLE :: cc
296  REAL(RFREAL) :: delfrac,xdel,xp,xmax,xmin,ydel,yp,ymax,ymin,zdel,zp,zmax,zmin
297  TYPE(t_global), POINTER :: global
298  TYPE(t_grid), POINTER :: pgrid
299 
300 ! ******************************************************************************
301 ! Start
302 ! ******************************************************************************
303 
304  global => pregion%global
305 
306  CALL registerfunction(global,'RFLU_FindProbeCells',&
307  'RFLU_ModProbes.F90')
308 
309  IF ( global%myProcid == masterproc .AND. &
310  global%verbLevel >= verbose_high ) THEN
311  WRITE(stdout,'(A,1X,A)') solver_name,'Finding cells containing probes...'
312  END IF ! global%verbLevel
313 
314 ! ******************************************************************************
315 ! Set variables and pointers
316 ! ******************************************************************************
317 
318  pgrid => pregion%grid
319 
320  delfrac = 0.01_rfreal
321 
322  ccsize = min(100,pgrid%nCells) ! Must be larger than unity
323 
324  ALLOCATE(cc(ccsize),stat=errorflag)
325  global%error = errorflag
326  IF ( global%error /= err_none ) THEN
327  CALL errorstop(global,err_allocate,__line__,'cc')
328  END IF ! global%error
329 
330 ! ******************************************************************************
331 ! Build bounding box
332 ! ******************************************************************************
333 
334  xmin = minval(pgrid%xyz(xcoord,1:pgrid%nVert))
335  xmax = maxval(pgrid%xyz(xcoord,1:pgrid%nVert))
336  ymin = minval(pgrid%xyz(ycoord,1:pgrid%nVert))
337  ymax = maxval(pgrid%xyz(ycoord,1:pgrid%nVert))
338  zmin = minval(pgrid%xyz(zcoord,1:pgrid%nVert))
339  zmax = maxval(pgrid%xyz(zcoord,1:pgrid%nVert))
340 
341  xdel = xmax - xmin
342  ydel = ymax - ymin
343  zdel = zmax - zmin
344 
345  xmin = xmin - delfrac*xdel
346  xmax = xmax + delfrac*xdel
347  ymin = ymin - delfrac*ydel
348  ymax = ymax + delfrac*ydel
349  zmin = zmin - delfrac*zdel
350  zmax = zmax + delfrac*zdel
351 
352 ! ******************************************************************************
353 ! Build Octree with cell centroids: NOTE only for actual cells
354 ! ******************************************************************************
355 
356  CALL rflu_createoctree(global,pgrid%nCells)
357 
358  CALL rflu_buildoctree(pgrid%cofg(xcoord,1:pgrid%nCells), &
359  pgrid%cofg(ycoord,1:pgrid%nCells), &
360  pgrid%cofg(zcoord,1:pgrid%nCells), &
362 
363 ! ******************************************************************************
364 ! Loop over probes
365 ! ******************************************************************************
366 
367  DO iprobe = 1,global%nProbes
368 
369 ! ==============================================================================
370 ! Test probe location against bounding box of partition
371 ! ==============================================================================
372 
373  xp = global%probeXyz(iprobe,1)
374  yp = global%probeXyz(iprobe,2)
375  zp = global%probeXyz(iprobe,3)
376 
377  IF ( rflu_testinboundbox(global,xp,yp,zp,xmin,xmax,ymin,ymax, &
378  zmin,zmax) .EQV. .true. ) THEN
379 
380 ! ------------------------------------------------------------------------------
381 ! Query octree to get closest cells
382 ! ------------------------------------------------------------------------------
383 
384  CALL rflu_queryoctree(xp,yp,zp,ccsize,cc)
385 
386 ! ------------------------------------------------------------------------------
387 ! Test cells obtained from Octree if they contain specified location
388 ! ------------------------------------------------------------------------------
389 
390  cellloop: DO icg = 1,ccsize
391  IF ( rflu_ict_testincell(pregion,xp,yp,zp,cc(icg)) .EQV. .true. ) THEN
392  global%probePos(iprobe,probe_region) = pregion%iRegionGlobal
393  global%probePos(iprobe,probe_cell) = cc(icg)
394 
395  EXIT cellloop
396  END IF ! RFLU_ICT_TestInCell
397  END DO cellloop
398  END IF ! RFLU_TestInBoundNox
399  END DO ! iProbe
400 
401 ! ******************************************************************************
402 ! Destroy Octree and deallocate temporary memory
403 ! ******************************************************************************
404 
405  CALL rflu_destroyoctree(global)
406 
407  DEALLOCATE(cc,stat=errorflag)
408  global%error = errorflag
409  IF ( global%error /= err_none ) THEN
410  CALL errorstop(global,err_deallocate,__line__,'cc')
411  END IF ! global%error
412 
413 ! ******************************************************************************
414 ! End
415 ! ******************************************************************************
416 
417  IF ( global%myProcid == masterproc .AND. &
418  global%verbLevel >= verbose_high ) THEN
419  WRITE(stdout,'(A,1X,A)') solver_name,'Finding cells containing probes done.'
420  END IF ! global%verbLevel
421 
422  CALL deregisterfunction(global)
423 
424 END SUBROUTINE rflu_findprobecells
425 
426 
427 
428 
429 
430 
431 
432 
433 
434 
435 ! ******************************************************************************
436 !
437 ! Purpose: Open probe files and position files at right line when restarting.
438 !
439 ! Description: None.
440 !
441 ! Input:
442 ! pRegion Pointer to region
443 !
444 ! Output: None.
445 !
446 ! Notes: None.
447 !
448 ! ******************************************************************************
449 
450 SUBROUTINE rflu_openprobefiles(pRegion)
451 
453  USE modtools, ONLY: floatless
454 
455  IMPLICIT NONE
456 
457 ! ******************************************************************************
458 ! Declarations
459 ! ******************************************************************************
460 
461 ! ==============================================================================
462 ! Arguments
463 ! ==============================================================================
464 
465  TYPE (t_region), POINTER :: pregion
466 
467 ! ******************************************************************************
468 ! Locals
469 ! ******************************************************************************
470 
471  LOGICAL :: fileappend,fileexists
472  CHARACTER(CHRLEN+9) :: fname
473  INTEGER :: errorflag,iprobe,probeiter
474  REAL(RFREAL) :: probetime
475  TYPE(t_global), POINTER :: global
476 
477 ! ******************************************************************************
478 ! Start
479 ! ******************************************************************************
480 
481  global => pregion%global
482 
483  CALL registerfunction(global,'RFLU_OpenProbeFiles',&
484  'RFLU_ModProbes.F90')
485 
486 ! ******************************************************************************
487 ! Initialize
488 ! ******************************************************************************
489 
490  probetime = huge(1.0_rfreal)
491 
492 ! ******************************************************************************
493 ! Loop over probes
494 ! ******************************************************************************
495 
496  DO iprobe = 1,global%nProbes
497 
498 ! ==============================================================================
499 ! Check if probe region on current processor
500 ! ==============================================================================
501 
502  IF ( global%probePos(iprobe,probe_region) == pregion%iRegionGlobal ) THEN
503 
504 ! ------------------------------------------------------------------------------
505 ! Generate file name
506 ! ------------------------------------------------------------------------------
507 
508  WRITE(fname,'(A,I4.4)') trim(global%outDir)// &
509  trim(global%casename)//'.prb_',iprobe
510 
511 ! ------------------------------------------------------------------------------
512 ! Open file with appropriate options
513 ! ------------------------------------------------------------------------------
514 
515  IF ( (global%flowType == flow_unsteady .AND. &
516  global%currentTime > 0.0_rfreal) .OR. &
517  (global%flowType == flow_steady .AND. &
518  global%currentIter > 1) ) THEN
519  INQUIRE(file=fname,exist=fileexists)
520 
521  IF ( fileexists .EQV. .true. ) THEN
522  fileappend = .true.
523 
524  OPEN(if_probe+iprobe-1,file=fname,form='FORMATTED',status='OLD', &
525  position='APPEND',iostat=errorflag)
526  ELSE
527  fileappend = .false.
528 
529  OPEN(if_probe+iprobe-1,file=fname,form='FORMATTED',status='UNKNOWN', &
530  iostat=errorflag)
531  END IF ! fileExists
532  ELSE
533  fileappend = .false.
534 
535  OPEN(if_probe+iprobe-1,file=fname,form='FORMATTED',status='UNKNOWN', &
536  iostat=errorflag)
537  END IF ! global
538 
539  global%error = errorflag
540  IF (global%error /= err_none ) THEN
541  CALL errorstop(global,err_file_open,__line__,'File: '//trim(fname))
542  END IF ! global%error
543 
544 ! ------------------------------------------------------------------------------
545 ! Back up to right time when appending to existing file. NOTE if encounter
546 ! errors when backspacing or reading file (for example because it may be
547 ! empty), just quit loop instead of treating error properly.
548 ! ------------------------------------------------------------------------------
549 
550  IF ( fileappend .EQV. .true. ) THEN
551 
552 ! ----- Unsteady flow ----------------------------------------------------------
553 
554  IF ( global%flowType == flow_unsteady ) THEN
555  emptyloopunsteady: DO
556  backspace(if_probe+iprobe-1,iostat=errorflag)
557  IF ( errorflag /= err_none ) THEN
558  EXIT emptyloopunsteady
559  END IF ! errorFlag
560 
561  READ(if_probe+iprobe-1,*,iostat=errorflag) probetime
562  IF ( errorflag /= err_none ) THEN
563  EXIT emptyloopunsteady
564  END IF ! errorFlag
565 
566  IF ( floatless(probetime,global%currentTime) .EQV. .true. ) THEN
567  EXIT emptyloopunsteady
568  ELSE
569  backspace(if_probe+iprobe-1,iostat=errorflag)
570  IF ( errorflag /= err_none ) THEN
571  EXIT emptyloopunsteady
572  END IF ! errorFlag
573  END IF ! probeTime
574  END DO emptyloopunsteady
575 
576 ! ----- Steady flow ------------------------------------------------------------
577 
578  ELSE
579  emptyloopsteady: DO
580  backspace(if_probe+iprobe-1,iostat=errorflag)
581  IF ( errorflag /= err_none ) THEN
582  EXIT emptyloopsteady
583  END IF ! errorFlag
584 
585  READ(if_probe+iprobe-1,*,iostat=errorflag) probeiter
586  IF ( errorflag /= err_none ) THEN
587  EXIT emptyloopsteady
588  END IF ! errorFlag
589 
590  IF ( probeiter < global%currentIter ) THEN
591  EXIT emptyloopsteady
592  ELSE
593  backspace(if_probe+iprobe-1,iostat=errorflag)
594  IF ( errorflag /= err_none ) THEN
595  EXIT emptyloopsteady
596  END IF ! errorFlag
597  END IF ! probeTime
598  END DO emptyloopsteady
599  END IF ! global%flowType
600  END IF ! fileAppend
601  END IF ! global
602  END DO ! iProbe
603 
604 ! ******************************************************************************
605 ! End
606 ! ******************************************************************************
607 
608  CALL deregisterfunction(global)
609 
610 END SUBROUTINE rflu_openprobefiles
611 
612 
613 
614 
615 
616 
617 
618 ! ******************************************************************************
619 !
620 ! Purpose: Print information on probe cells.
621 !
622 ! Description: None.
623 !
624 ! Input:
625 ! global Pointer to global data
626 !
627 ! Output: None.
628 !
629 ! Notes:
630 ! 1. It is important to note that the use of the max reduction requires the
631 ! global%probePos array to be initialized with negative integers (see
632 ! readProbeSection.F90). Hence if the entries of a given probe are still
633 ! equal to the initial values, know that probe was not located in any of
634 ! the regions.
635 !
636 ! ******************************************************************************
637 
638 SUBROUTINE rflu_printprobeinfo(global)
639 
640  IMPLICIT NONE
641 
642 ! ******************************************************************************
643 ! Declarations
644 ! ******************************************************************************
645 
646 ! ==============================================================================
647 ! Arguments
648 ! ==============================================================================
649 
650  TYPE(t_global), POINTER :: global
651 
652 ! ==============================================================================
653 ! Locals
654 ! ==============================================================================
655 
656  CHARACTER(CHRLEN) :: errorstring
657  INTEGER :: errorflag,iprobe
658  INTEGER, DIMENSION(:), ALLOCATABLE :: globalvals,localvals
659 
660 ! ******************************************************************************
661 ! Start
662 ! ******************************************************************************
663 
664  CALL registerfunction(global,'RFLU_PrintProbeInfo',&
665  'RFLU_ModProbes.F90')
666 
667  IF ( global%myProcid == masterproc .AND. &
668  global%verbLevel >= verbose_high ) THEN
669  WRITE(stdout,'(A,1X,A)') solver_name,'Printing probe information...'
670  END IF ! global%verbLevel
671 
672 ! ******************************************************************************
673 ! Gather information
674 ! ******************************************************************************
675 
676 ! ==============================================================================
677 ! Allocate temporary memory
678 ! ==============================================================================
679 
680  ALLOCATE(localvals(global%nProbes),stat=errorflag)
681  global%error = errorflag
682  IF ( global%error /= err_none ) THEN
683  CALL errorstop(global,err_allocate,__line__,'localVals')
684  END IF ! global%error
685 
686  ALLOCATE(globalvals(global%nProbes),stat=errorflag)
687  global%error = errorflag
688  IF ( global%error /= err_none ) THEN
689  CALL errorstop(global,err_allocate,__line__,'globalVals')
690  END IF ! global%error
691 
692 ! ==============================================================================
693 ! Peform reduction operations
694 ! ==============================================================================
695 
696  DO iprobe = 1,global%nProbes
697  localvals(iprobe) = global%probePos(iprobe,probe_region)
698  END DO ! iProbe
699 
700  CALL mpi_allreduce(localvals,globalvals,global%nProbes,mpi_integer,mpi_max, &
701  global%mpiComm,errorflag)
702  global%error = errorflag
703  IF ( global%error /= err_none ) THEN
704  CALL errorstop(global,err_mpi_trouble,__line__)
705  END IF ! global%errorFlag
706 
707  DO iprobe = 1,global%nProbes
708  global%probePos(iprobe,probe_region) = globalvals(iprobe)
709  END DO ! iProbe
710 
711  DO iprobe = 1,global%nProbes
712  localvals(iprobe) = global%probePos(iprobe,probe_cell)
713  END DO ! iProbe
714 
715  CALL mpi_allreduce(localvals,globalvals,global%nProbes,mpi_integer,mpi_max, &
716  global%mpiComm,errorflag)
717  global%error = errorflag
718  IF ( global%error /= err_none ) THEN
719  CALL errorstop(global,err_mpi_trouble,__line__)
720  END IF ! global%errorFlag
721 
722  DO iprobe = 1,global%nProbes
723  global%probePos(iprobe,probe_cell) = globalvals(iprobe)
724  END DO ! iProbe
725 
726 ! ==============================================================================
727 ! Deallocate temporary memory
728 ! ==============================================================================
729 
730  DEALLOCATE(localvals,stat=errorflag)
731  global%error = errorflag
732  IF ( global%error /= err_none ) THEN
733  CALL errorstop(global,err_deallocate,__line__,'localVals')
734  END IF ! global%error
735 
736  DEALLOCATE(globalvals,stat=errorflag)
737  global%error = errorflag
738  IF ( global%error /= err_none ) THEN
739  CALL errorstop(global,err_deallocate,__line__,'globalVals')
740  END IF ! global%error
741 
742 ! ******************************************************************************
743 ! If probe position values are still equal to initial values, could not find
744 ! a cell containing probe position
745 ! ******************************************************************************
746 
747  IF ( global%myProcid == masterproc ) THEN
748  DO iprobe = 1,global%nProbes
749  IF ( (global%probePos(iprobe,probe_region) == crazy_value_int) .OR. &
750  (global%probePos(iprobe,probe_cell ) == crazy_value_int) ) THEN
751  WRITE(errorstring,'(A,1X,I3)') 'Probe:',iprobe
752  CALL errorstop(global,err_probe_location,__line__,trim(errorstring))
753  END IF ! global%probePos
754  END DO ! iProbe
755  END IF ! global%myProcid
756 
757 ! ******************************************************************************
758 ! Print information
759 ! ******************************************************************************
760 
761  IF ( global%myProcid == masterproc .AND. &
762  global%verbLevel >= verbose_high ) THEN
763  WRITE(stdout,'(A,5X,A,3X,A,8X,A)') solver_name,'#','Region','Cell'
764 
765  DO iprobe = 1,global%nProbes
766  WRITE(stdout,'(A,3X,I3,3X,I6,3X,I9)') solver_name,iprobe, &
767  global%probePos(iprobe,probe_region),global%probePos(iprobe,probe_cell)
768  END DO ! iProbe
769  END IF ! global%myProcid
770 
771 ! ******************************************************************************
772 ! End
773 ! ******************************************************************************
774 
775  IF ( global%myProcid == masterproc .AND. &
776  global%verbLevel >= verbose_high ) THEN
777  WRITE(stdout,'(A,1X,A)') solver_name,'Printing probe information done.'
778  END IF ! global%verbLevel
779 
780  CALL deregisterfunction(global)
781 
782 END SUBROUTINE rflu_printprobeinfo
783 
784 
785 
786 
787 
788 
789 ! ******************************************************************************
790 ! End
791 ! ******************************************************************************
792 
793 END MODULE rflu_modprobes
794 
795 
796 ! ******************************************************************************
797 !
798 ! RCS Revision history:
799 !
800 ! $Log: RFLU_ModProbes.F90,v $
801 ! Revision 1.6 2008/12/06 08:44:23 mtcampbe
802 ! Updated license.
803 !
804 ! Revision 1.5 2008/11/19 22:17:34 mtcampbe
805 ! Added Illinois Open Source License/Copyright
806 !
807 ! Revision 1.4 2006/12/15 13:25:08 haselbac
808 ! Fixed bug in format statement, found by ifort
809 !
810 ! Revision 1.3 2006/04/07 15:19:20 haselbac
811 ! Removed tabs
812 !
813 ! Revision 1.2 2005/12/24 21:30:45 haselbac
814 ! Adapted to changes in ICT
815 !
816 ! Revision 1.1 2005/04/29 12:41:04 haselbac
817 ! Initial revision
818 !
819 ! ******************************************************************************
820 
821 
822 
823 
824 
825 
826 
827 
828 
829 
830 
831 
double ymin() const
double xmax() const
subroutine, public rflu_queryoctree(XPT, YPT, ZPT, NUMP, NEIGHP)
LOGICAL function, public rflu_testinboundbox(global, xLoc, yLoc, zLoc, xMin, xMax, yMin, yMax, zMin, zMax)
double xmin() const
subroutine registerfunction(global, funName, fileName)
Definition: ModError.F90:449
int status() const
Obtain the status of the attribute.
Definition: Attribute.h:240
double zmin() const
subroutine, public rflu_buildoctree(XI, YI, ZI, XLOW, XUPP, YLOW, YUPP, ZLOW, ZUPP)
subroutine, public rflu_findprobecells(pRegion)
LOGICAL function, public rflu_decidewriteprobes(global)
subroutine buildfilenameplain(global, dest, ext, fileName)
LOGICAL function floatless(a, b)
Definition: ModTools.F90:140
subroutine, public rflu_openprobefiles(pRegion)
subroutine, public rflu_printprobeinfo(global)
subroutine, public rflu_destroyoctree(global)
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE SUBROUTINE form
double zmax() const
double ymax() const
LOGICAL function, public rflu_ict_testincell(pRegion, xLoc, yLoc, zLoc, icg)
Vector_n min(const Array_n_const &v1, const Array_n_const &v2)
Definition: Vector_n.h:346
subroutine errorstop(global, errorCode, errorLine, addMessage)
Definition: ModError.F90:483
subroutine deregisterfunction(global)
Definition: ModError.F90:469
subroutine, public rflu_createoctree(global, nPoints)
subroutine, public rflu_closeprobefiles(pRegion)