Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
RFLU_ModFaceList.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: Suite of routines to build face list.
26 !
27 ! Description: None.
28 !
29 ! Notes: None.
30 !
31 ! ******************************************************************************
32 !
33 ! $Id: RFLU_ModFaceList.F90,v 1.43 2008/12/06 08:44:21 mtcampbe Exp $
34 !
35 ! Copyright: (c) 2002-2006 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 modbndpatch, ONLY: t_patch
46  USE moddatastruct, ONLY: t_region
47  USE modgrid, ONLY: t_grid
48  USE modborder, ONLY: t_border
49  USE modmpi
50 
51  IMPLICIT NONE
52 
53  PRIVATE
54  PUBLIC :: rflu_buildavface2borderlist, &
69 
70  SAVE
71 
72 ! ******************************************************************************
73 ! Declarations and definitions
74 ! ******************************************************************************
75 
76  INTEGER, PARAMETER, PUBLIC :: DESTROY_FACE_INT = 1, &
77  DESTROY_FACE_EXT = 2, &
78  DESTROY_FACE_ALL = 4
79 
80  CHARACTER(CHRLEN) :: &
81  RCSIdentString = '$RCSfile: RFLU_ModFaceList.F90,v $ $Revision: 1.43 $'
82 
83 ! ******************************************************************************
84 ! Routines
85 ! ******************************************************************************
86 
87  CONTAINS
88 
89 
90 
91 
92 
93 
94 
95 ! ******************************************************************************
96 !
97 ! Purpose: Build actual-virtual-face-to-border list.
98 !
99 ! Description: None.
100 !
101 ! Input:
102 ! pRegion Pointer to region
103 !
104 ! Output: None.
105 !
106 ! Notes: None.
107 !
108 ! ******************************************************************************
109 
110  SUBROUTINE rflu_buildavface2borderlist(pRegion)
111 
112  USE modsortsearch
113 
114  IMPLICIT NONE
115 
116 ! ******************************************************************************
117 ! Declarations and definitions
118 ! ******************************************************************************
119 
120 ! ==============================================================================
121 ! Arguments
122 ! ==============================================================================
123 
124  TYPE(t_region), POINTER :: pregion
125 
126 ! ==============================================================================
127 ! Locals
128 ! ==============================================================================
129 
130  INTEGER :: c1,c2,errorflag,iborder,icg,icgmax,icl,ifg,ifl,iloc, &
131  ncellsrecvtot
132  INTEGER, DIMENSION(:,:), ALLOCATABLE :: vc,vc2b
133  TYPE(t_global), POINTER :: global
134  TYPE(t_border), POINTER :: pborder
135  TYPE(t_grid), POINTER :: pgrid
136  TYPE(t_patch), POINTER :: ppatch
137 
138 ! ******************************************************************************
139 ! Start
140 ! ******************************************************************************
141 
142  global => pregion%global
143 
144  CALL registerfunction(global,'RFLU_BuildAVFace2BorderList',&
145  'RFLU_ModFaceList.F90')
146 
147  IF ( global%myProcid == masterproc .AND. &
148  global%verbLevel >= verbose_high ) THEN
149  WRITE(stdout,'(A,1X,A)') solver_name, &
150  'Building av-face-to-border list...'
151  END IF ! global%verbLevel
152 
153 ! ******************************************************************************
154 ! Set grid pointer
155 ! ******************************************************************************
156 
157  pgrid => pregion%grid
158 
159 ! ******************************************************************************
160 ! Build temporary lists
161 ! ******************************************************************************
162 
163 ! ==============================================================================
164 ! Determine maximum cell index of actual-virtual faces so can limit searching
165 ! ==============================================================================
166 
167  icgmax = 0
168 
169  DO ifg = pgrid%nFaces-pgrid%nFacesAV+1,pgrid%nFaces
170  c1 = pgrid%f2c(1,ifg)
171  c2 = pgrid%f2c(2,ifg)
172 
173  IF ( c1 > pgrid%nCells ) THEN ! c1 is virtual cell
174  icgmax = max(c1,icgmax)
175  ELSE IF ( c2 > pgrid%nCells ) THEN ! c2 is virtual cell
176  icgmax = max(c2,icgmax)
177  ELSE ! Error
178 ! TEMPORARY
179  WRITE(*,*) 'ERROR 1!'
180  stop
181 ! END TEMPORARY
182  END IF ! c1
183  END DO ! ifg
184 
185 ! ==============================================================================
186 ! Determine number of border cells whose global cell indices are smaller than
187 ! maximum cell index determined above and allocate temporary memory
188 ! ==============================================================================
189 
190  ncellsrecvtot = 0
191 
192  DO iborder = 1,pgrid%nBorders
193  pborder => pgrid%borders(iborder)
194 
195  DO icl = 1,pborder%nCellsRecv
196  IF ( pborder%icgRecv(icl) <= icgmax ) THEN
197  ncellsrecvtot = ncellsrecvtot + 1
198  END IF ! pBorder%icgRecv
199  END DO ! icl
200  END DO ! iBorder
201 
202  ALLOCATE(vc(2,ncellsrecvtot),stat=errorflag)
203  global%error = errorflag
204  IF ( global%error /= err_none ) THEN
205  CALL errorstop(global,err_allocate,__line__,'vc')
206  END IF ! global%error
207 
208  ALLOCATE(vc2b(2,ncellsrecvtot),stat=errorflag)
209  global%error = errorflag
210  IF ( global%error /= err_none ) THEN
211  CALL errorstop(global,err_allocate,__line__,'vc2b')
212  END IF ! global%error
213 
214 ! ==============================================================================
215 ! Build temporary lists and sort to allow searching
216 ! ==============================================================================
217 
218  ncellsrecvtot = 0 ! Reinitialize
219 
220  DO iborder = 1,pgrid%nBorders
221  pborder => pgrid%borders(iborder)
222 
223  DO icl = 1,pborder%nCellsRecv
224  IF ( pborder%icgRecv(icl) <= icgmax ) THEN
225  ncellsrecvtot = ncellsrecvtot + 1
226 
227  vc(1,ncellsrecvtot) = pborder%icgRecv(icl)
228  vc(2,ncellsrecvtot) = ncellsrecvtot
229 
230  vc2b(1,ncellsrecvtot) = iborder
231  vc2b(2,ncellsrecvtot) = icl
232  END IF ! pBorder%icgRecv
233  END DO ! icl
234  END DO ! iBorder
235 
236  CALL quicksortintegerinteger(vc(1:1,1:ncellsrecvtot), &
237  vc(2:2,1:ncellsrecvtot), &
238  ncellsrecvtot)
239 
240 ! ******************************************************************************
241 ! Build
242 ! ******************************************************************************
243 
244  DO ifg = pgrid%nFaces-pgrid%nFacesAV+1,pgrid%nFaces
245  c1 = pgrid%f2c(1,ifg)
246  c2 = pgrid%f2c(2,ifg)
247 
248  IF ( c1 > pgrid%nCells ) THEN ! c1 is virtual cell
249  icg = c1
250  ELSE IF ( c2 > pgrid%nCells ) THEN ! c2 is virtual cell
251  icg = c2
252  ELSE ! Error
253 ! TEMPORARY
254  WRITE(*,*) 'ERROR 2!'
255  stop
256 ! END TEMPORARY
257  END IF ! c1
258 
259  CALL binarysearchinteger(vc(1:1,1:ncellsrecvtot),ncellsrecvtot,icg,iloc)
260 
261  IF ( iloc /= element_not_found ) THEN
262  icl = vc(2,iloc)
263 
264  ifl = ifg - (pgrid%nFaces-pgrid%nFacesAV)
265 
266  pgrid%avf2b(1,ifl) = vc2b(1,icl)
267  pgrid%avf2b(2,ifl) = vc2b(2,icl)
268 
269  ELSE
270 ! TEMPORARY
271  WRITE(*,*) 'ERROR 3!'
272  stop
273 ! END TEMPORARY
274  END IF ! iLoc
275  END DO ! ifg
276 
277 ! ==============================================================================
278 ! Deallocate temporary memory
279 ! ==============================================================================
280 
281  DEALLOCATE(vc,stat=errorflag)
282  global%error = errorflag
283  IF ( global%error /= err_none ) THEN
284  CALL errorstop(global,err_deallocate,__line__,'vc')
285  END IF ! global%error
286 
287  DEALLOCATE(vc2b,stat=errorflag)
288  global%error = errorflag
289  IF ( global%error /= err_none ) THEN
290  CALL errorstop(global,err_deallocate,__line__,'vc2b')
291  END IF ! global%error
292 
293 ! ******************************************************************************
294 ! End
295 ! ******************************************************************************
296 
297  IF ( global%myProcid == masterproc .AND. &
298  global%verbLevel >= verbose_high ) THEN
299  WRITE(stdout,'(A,1X,A)') solver_name, &
300  'Building av-face-to-border list done.'
301  END IF ! global%verbLevel
302 
303  CALL deregisterfunction(global)
304 
305  END SUBROUTINE rflu_buildavface2borderlist
306 
307 
308 
309 
310 
311 
312 
313 
314 ! ******************************************************************************
315 !
316 ! Purpose: Build actual-virtual-face-to-patch list.
317 !
318 ! Description: None.
319 !
320 ! Input:
321 ! pRegion Pointer to region
322 !
323 ! Output: None.
324 !
325 ! Notes:
326 ! 1. Acess avf2b lists, so these must be built first.
327 !
328 ! ******************************************************************************
329 
330  SUBROUTINE rflu_buildavface2patchlist(pRegion)
331 
332  USE modsortsearch
333 
334  IMPLICIT NONE
335 
336 ! ******************************************************************************
337 ! Declarations and definitions
338 ! ******************************************************************************
339 
340 ! ==============================================================================
341 ! Arguments
342 ! ==============================================================================
343 
344  TYPE(t_region), POINTER :: pregion
345 
346 ! ==============================================================================
347 ! Locals
348 ! ==============================================================================
349 
350  INTEGER :: errorflag,iborder,icg,icl,ifg,ifl,iloc,ipatch
351  TYPE(t_global), POINTER :: global
352  TYPE(t_border), POINTER :: pborder
353  TYPE(t_grid), POINTER :: pgrid
354  TYPE(t_patch), POINTER :: ppatch
355 
356 ! ******************************************************************************
357 ! Start
358 ! ******************************************************************************
359 
360  global => pregion%global
361 
362  CALL registerfunction(global,'RFLU_BuildAVFace2PatchList',&
363  'RFLU_ModFaceList.F90')
364 
365  IF ( global%myProcid == masterproc .AND. &
366  global%verbLevel >= verbose_high ) THEN
367  WRITE(stdout,'(A,1X,A)') solver_name, &
368  'Building av-face-to-patch list...'
369  END IF ! global%verbLevel
370 
371 ! ******************************************************************************
372 ! Set grid pointer
373 ! ******************************************************************************
374 
375  pgrid => pregion%grid
376 
377 ! ******************************************************************************
378 ! Build temporary lists of sorted virtual cells adjacent to patches
379 ! ******************************************************************************
380 
381  DO ipatch = 1,pgrid%nPatches
382  ppatch => pregion%patches(ipatch)
383 
384  IF ( ppatch%nBCellsVirt > 0 ) THEN
385  ALLOCATE(ppatch%bvcSorted(ppatch%nBCellsVirt),stat=errorflag)
386  global%error = errorflag
387  IF ( global%error /= err_none ) THEN
388  CALL errorstop(global,err_allocate,__line__,'pPatch%bvcSorted')
389  END IF ! global%error
390 
391  DO icl = 1,ppatch%nBCellsVirt
392  ppatch%bvcSorted(icl) = ppatch%bvc(icl)
393  END DO ! icl
394 
395  CALL quicksortinteger(ppatch%bvcSorted,ppatch%nBCellsVirt)
396  ELSE
397  nullify(ppatch%bvcSorted)
398  END IF ! pPatch%nBCellsVirt
399  END DO ! iPatch
400 
401 ! ******************************************************************************
402 ! Build lists of patch indices for each av face
403 ! ******************************************************************************
404 
405  DO ifl = 1,pgrid%nFacesAV
406  iborder = pgrid%avf2b(1,ifl)
407  pborder => pgrid%borders(iborder)
408 
409  icl = pgrid%avf2b(2,ifl)
410 
411  IF ( icl <= pborder%nCellsRecv ) THEN
412  icg = pborder%icgRecv(icl)
413  ELSE
414 ! TEMPORARY
415  WRITE(*,*) 'ERROR! Exceeding dims of pBorder%icgRecv'
416  stop
417 ! END TEMPORARY
418  END IF ! icl
419 
420  pgrid%avf2p(ifl) = crazy_value_int
421 
422  patchloop: DO ipatch = 1,pgrid%nPatches
423  ppatch => pregion%patches(ipatch)
424 
425  IF ( ppatch%nBCellsVirt > 0 ) THEN
426  CALL binarysearchinteger(ppatch%bvcSorted,ppatch%nBCellsVirt, &
427  icg,iloc)
428 
429  IF ( iloc /= element_not_found ) THEN
430  pgrid%avf2p(ifl) = ipatch
431 
432  EXIT patchloop
433  END IF ! iLoc
434  END IF ! pPatch%nBCellsVirt
435  END DO patchloop
436  END DO ! ifl
437 
438 ! ******************************************************************************
439 ! Destroy temporary lists
440 ! ******************************************************************************
441 
442  DO ipatch = 1,pgrid%nPatches
443  ppatch => pregion%patches(ipatch)
444 
445  IF ( ppatch%nBCellsVirt > 0 ) THEN
446  DEALLOCATE(ppatch%bvcSorted,stat=errorflag)
447  global%error = errorflag
448  IF ( global%error /= err_none ) THEN
449  CALL errorstop(global,err_allocate,__line__,'pPatch%bvcSorted')
450  END IF ! global%error
451  END IF ! pPatch%nBCellsVirt
452  END DO ! iPatch
453 
454 ! ******************************************************************************
455 ! End
456 ! ******************************************************************************
457 
458  IF ( global%myProcid == masterproc .AND. &
459  global%verbLevel >= verbose_high ) THEN
460  WRITE(stdout,'(A,1X,A)') solver_name, &
461  'Building av-face-to-patch list done.'
462  END IF ! global%verbLevel
463 
464  CALL deregisterfunction(global)
465 
466  END SUBROUTINE rflu_buildavface2patchlist
467 
468 
469 
470 
471 
472 
473 
474 
475 ! ******************************************************************************
476 !
477 ! Purpose: Build cell-to-face list.
478 !
479 ! Description: None.
480 !
481 ! Input:
482 ! pRegion Pointer to region
483 !
484 ! Output: None.
485 !
486 ! Notes: None.
487 !
488 ! ******************************************************************************
489 
490  SUBROUTINE rflu_buildcell2facelist(pRegion)
491 
493 
494  IMPLICIT NONE
495 
496 ! ******************************************************************************
497 ! Declarations and definitions
498 ! ******************************************************************************
499 
500 ! ==============================================================================
501 ! Arguments
502 ! ==============================================================================
503 
504  TYPE(t_region), POINTER :: pregion
505 
506 ! ==============================================================================
507 ! Locals
508 ! ==============================================================================
509 
510  INTEGER :: errorflag,icl,icg,ick,ifg,ifl,ipatch
511  TYPE(t_global), POINTER :: global
512  TYPE(t_grid), POINTER :: pgrid
513  TYPE(t_patch), POINTER :: ppatch
514 
515 ! ******************************************************************************
516 ! Start
517 ! ******************************************************************************
518 
519  global => pregion%global
520 
521  CALL registerfunction(global,'RFLU_BuildCell2FaceList',&
522  'RFLU_ModFaceList.F90')
523 
524  IF ( global%myProcid == masterproc .AND. &
525  global%verbLevel >= verbose_high ) THEN
526  WRITE(stdout,'(A,1X,A)') solver_name,'Building cell-to-face list...'
527  END IF ! global%verbLevel
528 
529 ! ******************************************************************************
530 ! Set grid pointer
531 ! ******************************************************************************
532 
533  pgrid => pregion%grid
534 
535 ! ******************************************************************************
536 ! Initialize
537 ! ******************************************************************************
538 
539  IF ( pgrid%nTetsTot > 0 ) THEN
540  DO icl = 1,pgrid%nTetsTot ! Explicit loops because of ASCI White
541  DO ifl = 1,4
542  pgrid%tet2f(1,ifl,icl) = c2f_init
543  pgrid%tet2f(2,ifl,icl) = c2f_init
544  END DO ! ifl
545  END DO ! icl
546  END IF ! pGrid%nTetsTot
547 
548  IF ( pgrid%nHexsTot > 0 ) THEN
549  DO icl = 1,pgrid%nHexsTot ! Explicit loops because of ASCI White
550  DO ifl = 1,6
551  pgrid%hex2f(1,ifl,icl) = c2f_init
552  pgrid%hex2f(2,ifl,icl) = c2f_init
553  END DO ! ifl
554  END DO ! icl
555  END IF ! pGrid%nHexsTot
556 
557  IF ( pgrid%nPrisTot > 0 ) THEN
558  DO icl = 1,pgrid%nPrisTot ! Explicit loops because of ASCI White
559  DO ifl = 1,5
560  pgrid%pri2f(1,ifl,icl) = c2f_init
561  pgrid%pri2f(2,ifl,icl) = c2f_init
562  END DO ! ifl
563  END DO ! icl
564  END IF ! pGrid%nPrisTot
565 
566  IF ( pgrid%nPyrsTot > 0 ) THEN
567  DO icl = 1,pgrid%nPyrsTot ! Explicit loops because of ASCI White
568  DO ifl = 1,5
569  pgrid%pyr2f(1,ifl,icl) = c2f_init
570  pgrid%pyr2f(2,ifl,icl) = c2f_init
571  END DO ! ifl
572  END DO ! icl
573  END IF ! pGrid%nPyrsTot
574 
575 ! ******************************************************************************
576 ! Build cell-to-face list
577 ! ******************************************************************************
578 
579 ! ==============================================================================
580 ! Interior faces
581 ! ==============================================================================
582 
583  ipatch = 0 ! Indicates interior face
584 
585  DO ifg = 1,pgrid%nFacesTot
586  DO icl = 1,2
587  icg = pgrid%f2c(icl,ifg)
588  ick = rflu_getglobalcellkind(global,pgrid,icg)
589 
590  IF ( ick /= cell_kind_ext .AND. ick /= cell_kind_bnd ) THEN
591  CALL rflu_insertintocell2facelist(pregion,ipatch,icg,ifg)
592  END IF ! ick
593  END DO ! icl
594  END DO ! ifg
595 
596 ! ==============================================================================
597 ! Boundary faces. NOTE skip periodic and symmetry patches because otherwise
598 ! get duplicated entries and hence error message from algorithm in routine
599 ! RFLU_InsertIntoCell2FaceList
600 ! ==============================================================================
601 
602  DO ipatch = 1,pgrid%nPatches
603  ppatch => pregion%patches(ipatch)
604 
605  IF ( ppatch%bcType /= bc_periodic .AND. &
606  ppatch%bcType /= bc_symmetry ) THEN
607  DO ifl = 1,ppatch%nBFacesTot
608  icg = ppatch%bf2c(ifl)
609  ick = rflu_getglobalcellkind(global,pgrid,icg)
610 
611  CALL rflu_insertintocell2facelist(pregion,ipatch,icg,ifl)
612  END DO ! ifl
613  END IF ! pPatch%bcType
614  END DO ! iPatch
615 
616 ! ******************************************************************************
617 ! Check that all cells have all faces defined
618 ! ******************************************************************************
619 
620  IF ( global%checkLevel > check_none ) THEN
621  IF ( pgrid%nTetsTot > 0 ) THEN
622  DO icl = 1,pgrid%nTetsTot ! Explicit loops because of ASCI White
623  DO ifl = 1,4
624  IF ( pgrid%tet2f(1,ifl,icl) == c2f_init .OR. &
625  pgrid%tet2f(2,ifl,icl) == c2f_init ) THEN
626  CALL errorstop(global,err_c2flist_invalid,__line__)
627  END IF ! pGrid%tet2f
628  END DO ! ifl
629  END DO ! icl
630  END IF ! pGrid%nTetsTot
631 
632  IF ( pgrid%nHexsTot > 0 ) THEN
633  DO icl = 1,pgrid%nHexsTot ! Explicit loops because of ASCI White
634  DO ifl = 1,6
635  IF ( pgrid%hex2f(1,ifl,icl) == c2f_init .OR. &
636  pgrid%hex2f(2,ifl,icl) == c2f_init ) THEN
637  CALL errorstop(global,err_c2flist_invalid,__line__)
638  END IF ! pGrid%hex2f
639  END DO ! ifl
640  END DO ! icl
641  END IF ! pGrid%nHexsTot
642 
643  IF ( pgrid%nPrisTot > 0 ) THEN
644  DO icl = 1,pgrid%nPrisTot ! Explicit loops because of ASCI White
645  DO ifl = 1,5
646  IF ( pgrid%pri2f(1,ifl,icl) == c2f_init .OR. &
647  pgrid%pri2f(2,ifl,icl) == c2f_init ) THEN
648  CALL errorstop(global,err_c2flist_invalid,__line__)
649  END IF ! pGrid%pri2f
650  END DO ! ifl
651  END DO ! icl
652  END IF ! pGrid%nPrisTot
653 
654  IF ( pgrid%nPyrsTot > 0 ) THEN
655  DO icl = 1,pgrid%nPyrsTot ! Explicit loops because of ASCI White
656  DO ifl = 1,5
657  IF ( pgrid%pyr2f(1,ifl,icl) == c2f_init .OR. &
658  pgrid%pyr2f(2,ifl,icl) == c2f_init ) THEN
659  CALL errorstop(global,err_c2flist_invalid,__line__)
660  END IF ! pGrid%pyr2f
661  END DO ! ifl
662  END DO ! icl
663  END IF ! pGrid%nPyrsTot
664  END IF ! global%checkLevel
665 
666 #ifdef CHECK_DATASTRUCT
667 ! ******************************************************************************
668 ! Data structure output for checking
669 ! ******************************************************************************
670 
671  WRITE(stdout,'(A)') solver_name
672  WRITE(stdout,'(A,1X,A)') solver_name,'### START CHECK OUTPUT ###'
673  WRITE(stdout,'(A,1X,A)') solver_name,'Face-to-cell list'
674 
675  IF ( pgrid%nTetsTot > 0 ) THEN
676  WRITE(stdout,'(A,1X,A)') solver_name,'Tetrahedra:'
677  DO icl = 1,pgrid%nTetsTot
678  DO ifl = 1,4
679  WRITE(stdout,'(A,4(1X,I6))') solver_name,icl,ifl, &
680  pgrid%tet2f(1:2,ifl,icl)
681  END DO ! ifc
682  END DO ! icl
683  END IF ! pGrid%nTetsTot
684 
685  IF ( pgrid%nHexsTot > 0 ) THEN
686  WRITE(stdout,'(A,1X,A)') solver_name,'Hexahedra:'
687  DO icl = 1,pgrid%nHexsTot
688  DO ifl = 1,6
689  WRITE(stdout,'(A,4(1X,I6))') solver_name,icl,ifl, &
690  pgrid%hex2f(1:2,ifl,icl)
691  END DO ! ifl
692  END DO ! icl
693  END IF ! pGrid%nHexsTot
694 
695  IF ( pgrid%nPrisTot > 0 ) THEN
696  WRITE(stdout,'(A,1X,A)') solver_name,'Prisms:'
697  DO icl = 1,pgrid%nPrisTot
698  DO ifl = 1,5
699  WRITE(stdout,'(A,4(1X,I6))') solver_name,icl,ifl, &
700  pgrid%pri2f(1:2,ifl,icl)
701  END DO ! ifc
702  END DO ! ic
703  END IF ! pGrid%nPrisTot
704 
705  IF ( pgrid%nPyrsTot > 0 ) THEN
706  WRITE(stdout,'(A,1X,A)') solver_name,'Pyramids:'
707  DO icl = 1,pgrid%nPyrsTot
708  DO ifl = 1,5
709  WRITE(stdout,'(A,4(1X,I6))') solver_name,icl,ifl, &
710  pgrid%pyr2f(1:2,ifl,icl)
711  END DO ! ifl
712  END DO ! icl
713  END IF ! pGrid%nPyrsTot
714 
715  WRITE(stdout,'(A,1X,A)') solver_name,'### END CHECK OUTPUT ###'
716  WRITE(stdout,'(A)') solver_name
717 #endif
718 
719 ! ******************************************************************************
720 ! End
721 ! ******************************************************************************
722 
723  IF ( global%myProcid == masterproc .AND. &
724  global%verbLevel >= verbose_high ) THEN
725  WRITE(stdout,'(A,1X,A)') solver_name,'Building cell-to-face list done.'
726  END IF ! global%verbLevel
727 
728  CALL deregisterfunction(global)
729 
730  END SUBROUTINE rflu_buildcell2facelist
731 
732 
733 
734 
735 
736 
737 
738 
739 ! ******************************************************************************
740 !
741 ! Purpose: Build face list.
742 !
743 ! Description: None.
744 !
745 ! Input:
746 ! pRegion Pointer to region
747 !
748 ! Output: None.
749 !
750 ! Notes: None.
751 !
752 ! ******************************************************************************
753 
754  SUBROUTINE rflu_buildfacelist(pRegion)
755 
756  USE modsortsearch
757 
759  USE rflu_modgrid
761 
762  IMPLICIT NONE
763 
764 ! ******************************************************************************
765 ! Declarations and definitions
766 ! ******************************************************************************
767 
768 ! ==============================================================================
769 ! Arguments
770 ! ==============================================================================
771 
772  TYPE(t_region), POINTER :: pregion
773 
774 ! ==============================================================================
775 ! Locals
776 ! ==============================================================================
777 
778  INTEGER :: c1,c1k,c1t,c2,c2k,c2t,errorflag,fcntr,ftype,facetype,fksum, &
779  fvsize,icg,icl,ict,ifc,ifg,ifk,ifl,ipatch,iq,it,key, &
780  nbfacesquad,nbfacestri,nfacesint,nfacesquad,nfacestri, &
781  nhexslow,nhexsupp,nprislow,nprisupp,npyrslow,npyrsupp, &
782  ntetslow,ntetsupp,v1g,v2g,v3g,v4g
783  INTEGER :: fv(4)
784  INTEGER :: fkcntr(face_kind_aa:face_kind_ab,face_type_tri:face_type_quad), &
785  fkoffs(face_kind_aa:face_kind_ab)
786  INTEGER, DIMENSION(:), ALLOCATABLE :: fkind
787  TYPE(t_grid), POINTER :: pgrid
788  TYPE(t_patch), POINTER :: ppatch
789  TYPE(t_global), POINTER :: global
790 
791 ! TEMPORARY
792  INTEGER :: nbfaces
793 ! END TEMPORARY
794 
795 ! ******************************************************************************
796 ! Start
797 ! ******************************************************************************
798 
799  global => pregion%global
800 
801  CALL registerfunction(global,'RFLU_BuildFaceList',&
802  'RFLU_ModFaceList.F90')
803 
804  IF ( global%myProcid == masterproc .AND. &
805  global%verbLevel >= verbose_high ) THEN
806  WRITE(stdout,'(A,1X,A)') solver_name,'Building face list...'
807  WRITE(stdout,'(A,3X,A,1X,I5.5)') solver_name,'Global region:', &
808  pregion%iRegionGlobal
809  WRITE(stdout,'(A,3X,A)') solver_name,'Building hash table...'
810  END IF ! global%verbLevel
811 
812 ! ******************************************************************************
813 ! Set grid pointer and initialize
814 ! ******************************************************************************
815 
816  pgrid => pregion%grid
817 
818  pgrid%nFacesTot = 0
819 
820  nbfacesquad = 0
821  nbfacestri = 0
822  nfacestri = 0
823  nfacesquad = 0
824 
825 ! ******************************************************************************
826 ! Create hash table
827 ! ******************************************************************************
828 
829  CALL rflu_createhashtable(global,pgrid%nFacesEst)
830 
831  IF ( global%myProcid == masterproc .AND. &
832  global%verbLevel >= verbose_high ) THEN
833  WRITE(stdout,'(A,5X,A,12X,I9)') solver_name,'Hash table size:', &
834  hashtablesize
835  END IF ! global%verbLevel
836 
837 ! ******************************************************************************
838 ! Loop over cell types, construct hash table of faces
839 ! ******************************************************************************
840 
841  IF ( global%myProcid == masterproc .AND. &
842  global%verbLevel >= verbose_high ) THEN
843  WRITE(stdout,'(A,5X,A)') solver_name,'Looping over cell types...'
844  END IF ! global%verbLevel
845 
846 ! ==============================================================================
847 ! Tetrahedra
848 ! ==============================================================================
849 
850  IF ( pgrid%nTetsTot /= 0 ) THEN
851  IF ( global%myProcid == masterproc .AND. &
852  global%verbLevel >= verbose_high ) THEN
853  WRITE(stdout,'(A,7X,A)') solver_name,'Tetrahedra...'
854  END IF ! global%verbLevel
855  END IF ! pGrid%nTetsTot
856 
857  fvsize = 3 ! Tetrahedra only have 3 vertices per face
858 
859  DO icl = 1,pgrid%nTetsTot
860  icg = pgrid%tet2CellGlob(icl)
861 
862  DO ifl = 1,4
863  fv(1) = pgrid%tet2v(f2vtet(1,ifl),icl)
864  fv(2) = pgrid%tet2v(f2vtet(2,ifl),icl)
865  fv(3) = pgrid%tet2v(f2vtet(3,ifl),icl)
866 
867  CALL quicksortinteger(fv(1:fvsize),fvsize)
868  CALL rflu_hashbuildkey(fv(1:3),3,key)
869  CALL rflu_hashface(global,key,pgrid,icg,ifl,fv(1:3),fvsize,facetype)
870 
871  IF ( facetype == face_type_new ) THEN
872  nfacestri = nfacestri + 1
873  END IF ! faceType
874  END DO ! ifl
875  END DO ! icl
876 
877 ! ==============================================================================
878 ! Hexahedra
879 ! ==============================================================================
880 
881  IF ( pgrid%nHexsTot /= 0 ) THEN
882  IF ( global%myProcid == masterproc .AND. &
883  global%verbLevel >= verbose_high ) THEN
884  WRITE(stdout,'(A,7X,A)') solver_name,'Hexahedra...'
885  END IF ! global%verbLevel
886  END IF ! pGrid%nHexsTot
887 
888  fvsize = 4 ! Hexahedra only have 4 vertices per face
889 
890  DO icl = 1,pgrid%nHexsTot
891  icg = pgrid%hex2CellGlob(icl)
892 
893  DO ifl = 1,6
894  fv(1) = pgrid%hex2v(f2vhex(1,ifl),icl)
895  fv(2) = pgrid%hex2v(f2vhex(2,ifl),icl)
896  fv(3) = pgrid%hex2v(f2vhex(3,ifl),icl)
897  fv(4) = pgrid%hex2v(f2vhex(4,ifl),icl)
898 
899  CALL quicksortinteger(fv(1:fvsize),fvsize)
900  CALL rflu_hashbuildkey(fv(1:3),3,key)
901  CALL rflu_hashface(global,key,pgrid,icg,ifl,fv(1:3),fvsize,facetype)
902 
903  IF ( facetype == face_type_new ) THEN
904  nfacesquad = nfacesquad + 1
905  END IF ! faceType
906  END DO ! ifl
907  END DO ! icl
908 
909 ! ==============================================================================
910 ! Prisms
911 ! ==============================================================================
912 
913  IF ( pgrid%nPrisTot /= 0 ) THEN
914  IF ( global%myProcid == masterproc .AND. &
915  global%verbLevel >= verbose_high ) THEN
916  WRITE(stdout,'(A,7X,A)') solver_name,'Prisms...'
917  END IF ! global%verbLevel
918  END IF ! pGrid%nPrisTot
919 
920  DO icl = 1,pgrid%nPrisTot
921  icg = pgrid%pri2CellGlob(icl)
922 
923  DO ifl = 1,5
924  fvsize = 3
925  fv(1) = pgrid%pri2v(f2vpri(1,ifl),icl)
926  fv(2) = pgrid%pri2v(f2vpri(2,ifl),icl)
927  fv(3) = pgrid%pri2v(f2vpri(3,ifl),icl)
928 
929  IF ( f2vpri(4,ifl) /= vert_none ) THEN
930  fvsize = 4
931  fv(4) = pgrid%pri2v(f2vpri(4,ifl),icl)
932  END IF ! f2v%pri
933 
934  CALL quicksortinteger(fv(1:fvsize),fvsize)
935  CALL rflu_hashbuildkey(fv(1:3),3,key)
936  CALL rflu_hashface(global,key,pgrid,icg,ifl,fv(1:3),fvsize,facetype)
937 
938  IF ( facetype == face_type_new ) THEN
939  IF ( fvsize == 4 ) THEN
940  nfacesquad = nfacesquad + 1
941  ELSE
942  nfacestri = nfacestri + 1
943  END IF ! fvSize
944  END IF ! faceType
945  END DO ! ifl
946  END DO ! icl
947 
948 ! ==============================================================================
949 ! Pyramids
950 ! ==============================================================================
951 
952  IF ( pgrid%nPyrsTot /= 0 ) THEN
953  IF ( global%myProcid == masterproc .AND. &
954  global%verbLevel >= verbose_high ) THEN
955  WRITE(stdout,'(A,7X,A)') solver_name,'Pyramids...'
956  END IF ! global%verbLevel
957  END IF ! pGrid%nPyrsTot
958 
959  DO icl = 1,pgrid%nPyrsTot
960  icg = pgrid%pyr2CellGlob(icl)
961 
962  DO ifl = 1,5
963  fvsize = 3
964  fv(1) = pgrid%pyr2v(f2vpyr(1,ifl),icl)
965  fv(2) = pgrid%pyr2v(f2vpyr(2,ifl),icl)
966  fv(3) = pgrid%pyr2v(f2vpyr(3,ifl),icl)
967 
968  IF ( f2vpyr(4,ifl) /= vert_none ) THEN
969  fvsize = 4
970  fv(4) = pgrid%pyr2v(f2vpyr(4,ifl),icl)
971  END IF ! f2v%pyr
972 
973  CALL quicksortinteger(fv(1:fvsize),fvsize)
974  CALL rflu_hashbuildkey(fv(1:3),3,key)
975  CALL rflu_hashface(global,key,pgrid,icg,ifl,fv(1:3),fvsize,facetype)
976 
977  IF ( facetype == face_type_new ) THEN
978  IF ( fvsize == 3 ) THEN
979  nfacestri = nfacestri + 1
980  ELSE
981  nfacesquad = nfacesquad + 1
982  END IF ! fvSize
983  END IF ! faceType
984  END DO ! ifl
985  END DO ! icl
986 
987  IF ( global%myProcid == masterproc .AND. &
988  global%verbLevel >= verbose_high ) THEN
989  WRITE(stdout,'(A,5X,A,5X,I9)') solver_name,'Hash table collisions: ', &
990  hashtablecollisions
991  END IF ! global%myProcid
992 
993 ! ******************************************************************************
994 ! Number of internal faces (not on boundaries). NOTE nFacesTot is computed
995 ! in RFLU_HashFace. At this point it is the total number of faces. Further
996 ! below, it will be redefined to include only those faces which are not on
997 ! boundaries.
998 ! ******************************************************************************
999 
1000  nfacesint = pgrid%nFacesTot - pgrid%nBFacesTot
1001 
1002 ! ******************************************************************************
1003 ! Build actual facelist from hash table
1004 ! ******************************************************************************
1005 
1006 ! ==============================================================================
1007 ! Actual facelist - boundary faces. NOTE combined face list is ordered such
1008 ! that actual triangular faces appear first, then quadrilateral actual
1009 ! faces, virtual triangular faces, and quadrilateral virtual faces last.
1010 ! ==============================================================================
1011 
1012  IF ( global%myProcid == masterproc .AND. &
1013  global%verbLevel >= verbose_high ) THEN
1014  WRITE(stdout,'(A,3X,A)') solver_name,'Building boundary face lists...'
1015  END IF ! global%verbLevel
1016 
1017  DO ipatch = 1,pgrid%nPatches
1018  ppatch => pregion%patches(ipatch)
1019 
1020  nbfacestri = nbfacestri + ppatch%nBTrisTot
1021  nbfacesquad = nbfacesquad + ppatch%nBQuadsTot
1022 
1023  IF ( global%myProcid == masterproc .AND. &
1024  global%verbLevel >= verbose_high ) THEN
1025  WRITE(stdout,'(A,5X,A,I4,1X,A,1X,I4,A)') solver_name,'Patch: ',ipatch, &
1026  '(Global patch:',ppatch%iPatchGlobal,')'
1027 
1028  WRITE(stdout,'(A,7X,A,8X,2(1X,I9))') solver_name, &
1029  'Triangular faces:', &
1030  ppatch%nBTris, &
1031  ppatch%nBTrisTot-ppatch%nBTris
1032  WRITE(stdout,'(A,7X,A,5X,2(1X,I9))') solver_name, &
1033  'Quadrilateral faces:', &
1034  ppatch%nBQuads, &
1035  ppatch%nBQuadsTot-ppatch%nBQuads
1036  WRITE(stdout,'(A,7X,A,3X,2(1X,I9))') solver_name, &
1037  'Total number of faces:', &
1038  ppatch%nBFaces, &
1039  ppatch%nBFacesTot-ppatch%nBFaces
1040  END IF ! global%verbLevel
1041 
1042 ! ------------------------------------------------------------------------------
1043 ! Triangular faces - NOTE sorting of actual and virtual faces...
1044 ! ------------------------------------------------------------------------------
1045 
1046  fvsize = 3
1047 
1048  DO it = 1,ppatch%nBTrisTot
1049  fv(1) = ppatch%bTri2v(1,it)
1050  fv(2) = ppatch%bTri2v(2,it)
1051  fv(3) = ppatch%bTri2v(3,it)
1052 
1053  CALL quicksortinteger(fv(1:fvsize),fvsize)
1054  CALL rflu_hashbuildkey(fv(1:3),3,key)
1055  CALL rflu_unhashbface(global,key,pgrid,fv(1:3),fvsize,ppatch%bcType, &
1056  icg,ifg)
1057 
1058  fcntr = it
1059 
1060  IF ( it > ppatch%nBTris ) THEN
1061  fcntr = it + ppatch%nBQuads
1062  END IF ! it
1063 
1064  ppatch%bf2c( fcntr) = icg
1065  ppatch%bf2v(1,fcntr) = ppatch%bTri2v(1,it)
1066  ppatch%bf2v(2,fcntr) = ppatch%bTri2v(2,it)
1067  ppatch%bf2v(3,fcntr) = ppatch%bTri2v(3,it)
1068  END DO ! it
1069 
1070 ! ------------------------------------------------------------------------------
1071 ! Quadrilateral faces - NOTE sorting of actual and virtual faces...
1072 ! ------------------------------------------------------------------------------
1073 
1074  fvsize = 4
1075 
1076  DO iq = 1,ppatch%nBQuadsTot
1077  fv(1) = ppatch%bQuad2v(1,iq)
1078  fv(2) = ppatch%bQuad2v(2,iq)
1079  fv(3) = ppatch%bQuad2v(3,iq)
1080  fv(4) = ppatch%bQuad2v(4,iq)
1081 
1082  CALL quicksortinteger(fv(1:fvsize),fvsize)
1083  CALL rflu_hashbuildkey(fv(1:3),3,key)
1084  CALL rflu_unhashbface(global,key,pgrid,fv(1:3),fvsize,ppatch%bcType, &
1085  icg,ifg)
1086 
1087  IF ( iq <= ppatch%nBQuads ) THEN
1088  fcntr = iq + ppatch%nBTris
1089  ELSE
1090  fcntr = iq + ppatch%nBTrisTot
1091  END IF ! iq
1092 
1093  ppatch%bf2c( fcntr) = icg
1094  ppatch%bf2v(1,fcntr) = ppatch%bQuad2v(1,iq)
1095  ppatch%bf2v(2,fcntr) = ppatch%bQuad2v(2,iq)
1096  ppatch%bf2v(3,fcntr) = ppatch%bQuad2v(3,iq)
1097  ppatch%bf2v(4,fcntr) = ppatch%bQuad2v(4,iq)
1098  END DO ! iq
1099  END DO ! iPatch
1100 
1101 ! ==============================================================================
1102 ! Destroy hash table
1103 ! ==============================================================================
1104 
1105  CALL rflu_destroyhashtable(global)
1106 
1107 ! ==============================================================================
1108 ! Deallocate f2v array, not needed, regenerated below automatically
1109 ! ==============================================================================
1110 
1111  DEALLOCATE(pgrid%f2v,stat=errorflag)
1112  global%error = errorflag
1113  IF ( global%error /= err_none ) THEN
1114  CALL errorstop(global,err_deallocate,__line__,'pGrid%f2v')
1115  END IF ! global%error
1116 
1117 ! ******************************************************************************
1118 ! Determine face kinds. NOTE do not modify fkCntr array, it is used below.
1119 ! NOTE also: These statistics must be computed AFTER the boundary face lists
1120 ! are built, because before cannot distinguish between faces of kinds
1121 ! FACE_KIND_VB and FACE_KIND_VX. This is because f2c array is initialized
1122 ! with CELL_TYPE_EXT, and boundary-face lists, which are given explicitly
1123 ! in grid file, have not yet been used to mark those faces.
1124 ! ******************************************************************************
1125 
1126  IF ( global%myProcid == masterproc .AND. &
1127  global%verbLevel >= verbose_high ) THEN
1128  WRITE(stdout,'(A,3X,A)') solver_name,'Determining face statistics...'
1129  END IF ! global%verbLevel
1130 
1131  fkcntr(face_kind_aa,face_type_tri:face_type_quad) = 0
1132  fkcntr(face_kind_av,face_type_tri:face_type_quad) = 0
1133  fkcntr(face_kind_vv,face_type_tri:face_type_quad) = 0
1134  fkcntr(face_kind_vb,face_type_tri:face_type_quad) = 0
1135  fkcntr(face_kind_vx,face_type_tri:face_type_quad) = 0
1136  fkcntr(face_kind_ab,face_type_tri:face_type_quad) = 0
1137 
1138  DO ifc = 1,pgrid%nFacesTot
1139  c1 = pgrid%f2c(1,ifc)
1140  c2 = pgrid%f2c(2,ifc)
1141 
1142  ifl = pgrid%f2c(3,ifc)
1143  c1t = rflu_getglobalcelltype(global,pgrid,c1)
1144 
1145  SELECT CASE ( c1t )
1146  CASE ( cell_type_tet ) ! Tetrahedral cell
1147  ftype = face_type_tri
1148  CASE ( cell_type_hex ) ! Hexahedral cell
1149  ftype = face_type_quad
1150  CASE ( cell_type_pri ) ! Prismatic cell
1151  IF ( f2vpri(4,ifl) /= vert_none ) THEN
1152  ftype = face_type_quad
1153  ELSE
1154  ftype = face_type_tri
1155  END IF ! f2vPri
1156  CASE ( cell_type_pyr ) ! Pyramidal cell
1157  IF ( f2vpyr(4,ifl) == vert_none ) THEN
1158  ftype = face_type_tri
1159  ELSE
1160  ftype = face_type_quad
1161  END IF ! f2vPyr
1162  CASE default
1163  CALL errorstop(global,err_reached_default,__line__)
1164  END SELECT ! cellType
1165 
1166  c1k = rflu_getglobalcellkind(global,pgrid,c1)
1167  c2k = rflu_getglobalcellkind(global,pgrid,c2)
1168 
1169  SELECT CASE ( rflu_getfacekind(global,c1k,c2k) )
1170  CASE ( face_kind_aa )
1171  fkcntr(face_kind_aa,ftype) = fkcntr(face_kind_aa,ftype) + 1
1172  CASE ( face_kind_av )
1173  fkcntr(face_kind_av,ftype) = fkcntr(face_kind_av,ftype) + 1
1174  CASE ( face_kind_vv )
1175  fkcntr(face_kind_vv,ftype) = fkcntr(face_kind_vv,ftype) + 1
1176  CASE ( face_kind_vb )
1177  fkcntr(face_kind_vb,ftype) = fkcntr(face_kind_vb,ftype) + 1
1178  CASE ( face_kind_ab )
1179  fkcntr(face_kind_ab,ftype) = fkcntr(face_kind_ab,ftype) + 1
1180  CASE ( face_kind_vx )
1181  fkcntr(face_kind_vx,ftype) = fkcntr(face_kind_vx,ftype) + 1
1182  CASE default
1183  CALL errorstop(global,err_reached_default,__line__)
1184  END SELECT ! fKind
1185  END DO ! ifc
1186 
1187 ! ==============================================================================
1188 ! Write info
1189 ! ==============================================================================
1190 
1191  IF ( global%myProcid == masterproc .AND. &
1192  global%verbLevel >= verbose_high ) THEN
1193  WRITE(stdout,'(A,5X,A)') solver_name,'Face-type statistics:'
1194  WRITE(stdout,'(A,7X,A,5X,I9)') solver_name,'Total faces: ', &
1195  pgrid%nFacesTot
1196  WRITE(stdout,'(A,9X,A,3X,I9)') solver_name,'Triangular faces: ', &
1197  nfacestri
1198  WRITE(stdout,'(A,9X,A,3X,I9)') solver_name,'Quadrilateral faces: ', &
1199  nfacesquad
1200  WRITE(stdout,'(A,7X,A,5X,I9)') solver_name,'Total boundary faces: ', &
1201  pgrid%nBFacesTot
1202  WRITE(stdout,'(A,9X,A,3X,I9)') solver_name,'Triangular faces: ', &
1203  nbfacestri
1204  WRITE(stdout,'(A,9X,A,3X,I9)') solver_name,'Quadrilateral faces: ', &
1205  nbfacesquad
1206  WRITE(stdout,'(A,7X,A,5X,I9)') solver_name,'Non-boundary faces: ', &
1207  nfacesint
1208  WRITE(stdout,'(A,9X,A,3X,I9)') solver_name,'Triangular faces: ', &
1209  nfacestri - nbfacestri
1210  WRITE(stdout,'(A,9X,A,3X,I9)') solver_name,'Quadrilateral faces: ', &
1211  nfacesquad - nbfacesquad
1212  WRITE(stdout,'(A,5X,A)') solver_name,'Face-kind statistics:'
1213  WRITE(stdout,'(A,7X,A,5X,I9)') solver_name,'Actual-actual faces: ', &
1214  fkcntr(face_kind_aa,face_type_tri ) + &
1215  fkcntr(face_kind_aa,face_type_quad)
1216  WRITE(stdout,'(A,9X,A,3X,I9)') solver_name,'Triangular faces: ', &
1217  fkcntr(face_kind_aa,face_type_tri)
1218  WRITE(stdout,'(A,9X,A,3X,I9)') solver_name,'Quadrilateral faces: ', &
1219  fkcntr(face_kind_aa,face_type_quad)
1220  WRITE(stdout,'(A,7X,A,5X,I9)') solver_name,'Actual-virtual faces: ', &
1221  fkcntr(face_kind_av,face_type_tri ) + &
1222  fkcntr(face_kind_av,face_type_quad)
1223  WRITE(stdout,'(A,9X,A,3X,I9)') solver_name,'Triangular faces: ', &
1224  fkcntr(face_kind_av,face_type_tri)
1225  WRITE(stdout,'(A,9X,A,3X,I9)') solver_name,'Quadrilateral faces: ', &
1226  fkcntr(face_kind_av,face_type_quad)
1227  WRITE(stdout,'(A,7X,A,5X,I9)') solver_name,'Virtual-virtual faces: ', &
1228  fkcntr(face_kind_vv,face_type_tri ) + &
1229  fkcntr(face_kind_vv,face_type_quad)
1230  WRITE(stdout,'(A,9X,A,3X,I9)') solver_name,'Triangular faces: ', &
1231  fkcntr(face_kind_vv,face_type_tri)
1232  WRITE(stdout,'(A,9X,A,3X,I9)') solver_name,'Quadrilateral faces: ', &
1233  fkcntr(face_kind_vv,face_type_quad)
1234  WRITE(stdout,'(A,7X,A,5X,I9)') solver_name,'Virtual-boundary faces:', &
1235  fkcntr(face_kind_vb,face_type_tri ) + &
1236  fkcntr(face_kind_vb,face_type_quad)
1237  WRITE(stdout,'(A,9X,A,3X,I9)') solver_name,'Triangular faces: ', &
1238  fkcntr(face_kind_vb,face_type_tri)
1239  WRITE(stdout,'(A,9X,A,3X,I9)') solver_name,'Quadrilateral faces: ', &
1240  fkcntr(face_kind_vb,face_type_quad)
1241  WRITE(stdout,'(A,7X,A,5X,I9)') solver_name,'Actual-boundary faces: ', &
1242  fkcntr(face_kind_ab,face_type_tri ) + &
1243  fkcntr(face_kind_ab,face_type_quad)
1244  WRITE(stdout,'(A,9X,A,3X,I9)') solver_name,'Triangular faces: ', &
1245  fkcntr(face_kind_ab,face_type_tri)
1246  WRITE(stdout,'(A,9X,A,3X,I9)') solver_name,'Quadrilateral faces: ', &
1247  fkcntr(face_kind_ab,face_type_quad)
1248  WRITE(stdout,'(A,7X,A,5X,I9)') solver_name,'Virtual-external faces:', &
1249  fkcntr(face_kind_vx,face_type_tri ) + &
1250  fkcntr(face_kind_vx,face_type_quad)
1251  WRITE(stdout,'(A,9X,A,3X,I9)') solver_name,'Triangular faces: ', &
1252  fkcntr(face_kind_vx,face_type_tri)
1253  WRITE(stdout,'(A,9X,A,3X,I9)') solver_name,'Quadrilateral faces: ', &
1254  fkcntr(face_kind_vx,face_type_quad)
1255  END IF ! global%myProcid
1256 
1257  IF ( global%myProcid == masterproc .AND. &
1258  global%verbLevel >= verbose_high ) THEN
1259  WRITE(stdout,'(A,3X,A)') solver_name,'Determining face statistics done.'
1260  END IF ! global%verbLevel
1261 
1262 ! ******************************************************************************
1263 ! Actual facelist - interior faces
1264 ! ******************************************************************************
1265 
1266  IF ( global%myProcid == masterproc .AND. &
1267  global%verbLevel >= verbose_high ) THEN
1268  WRITE(stdout,'(A,3X,A)') solver_name,'Building non-boundary face '// &
1269  'lists...'
1270  END IF ! global%verbLevel
1271 
1272 ! ==============================================================================
1273 ! Determine face kind offsets for sorting of faces
1274 ! ==============================================================================
1275 
1276  fkoffs(face_kind_aa) = 0
1277  fkoffs(face_kind_av) = fkcntr(face_kind_aa,face_type_tri ) &
1278  + fkcntr(face_kind_aa,face_type_quad)
1279  fkoffs(face_kind_vv) = fkoffs(face_kind_av) &
1280  + fkcntr(face_kind_av,face_type_tri ) &
1281  + fkcntr(face_kind_av,face_type_quad)
1282  fkoffs(face_kind_vx) = fkoffs(face_kind_vv) &
1283  + fkcntr(face_kind_vv,face_type_tri ) &
1284  + fkcntr(face_kind_vv,face_type_quad)
1285 
1286 ! ==============================================================================
1287 ! Determine number of flux faces (nFaces) and total internal faces. NOTE need
1288 ! to set nBFaces and nBFacesTot here (in addition to RFLU_CreateGrid) because
1289 ! when adding virtual cells for symmetry or periodic patches, face-list is
1290 ! regenerated because of new cells, but RFLU_CreateGrid is not called again.
1291 ! ==============================================================================
1292 
1293  pgrid%nFaces = fkcntr(face_kind_aa,face_type_tri ) &
1294  + fkcntr(face_kind_aa,face_type_quad) &
1295  + fkcntr(face_kind_av,face_type_tri ) &
1296  + fkcntr(face_kind_av,face_type_quad)
1297 
1298  pgrid%nFacesAV = fkcntr(face_kind_av,face_type_tri ) &
1299  + fkcntr(face_kind_av,face_type_quad)
1300  pgrid%nFacesVV = fkcntr(face_kind_vv,face_type_tri ) &
1301  + fkcntr(face_kind_vv,face_type_quad)
1302 
1303  nfacesint = fkcntr(face_kind_aa,face_type_tri ) &
1304  + fkcntr(face_kind_aa,face_type_quad) &
1305  + fkcntr(face_kind_av,face_type_tri ) &
1306  + fkcntr(face_kind_av,face_type_quad) &
1307  + fkcntr(face_kind_vv,face_type_tri ) &
1308  + fkcntr(face_kind_vv,face_type_quad) &
1309  + fkcntr(face_kind_vx,face_type_tri ) &
1310  + fkcntr(face_kind_vx,face_type_quad)
1311 
1312  pgrid%nBFaces = fkcntr(face_kind_ab,face_type_tri ) &
1313  + fkcntr(face_kind_ab,face_type_quad)
1314 
1315  pgrid%nBFacesTot = fkcntr(face_kind_ab,face_type_tri ) &
1316  + fkcntr(face_kind_ab,face_type_quad) &
1317  + fkcntr(face_kind_vb,face_type_tri ) &
1318  + fkcntr(face_kind_vb,face_type_quad)
1319 
1320 ! ==============================================================================
1321 ! Allocate temporary face arrays and initialize arrays and counters. NOTE
1322 ! in the following, write counter information into fkCntr(ifk,FACE_TYPE_TRI)
1323 ! for convenience even if face is not a triangle
1324 ! ==============================================================================
1325 
1326  ALLOCATE(pgrid%f2cTemp(2,nfacesint),stat=errorflag) ! NOTE nFacesInt
1327  IF ( global%error /= err_none ) THEN
1328  CALL errorstop(global,err_allocate,__line__,'pGrid%f2cTemp')
1329  END IF ! global%error
1330 
1331  DO ifc = 1,nfacesint ! Explicit loop because of Frost
1332  pgrid%f2cTemp(1,ifc) = 0 ! Initial value immaterial
1333  pgrid%f2cTemp(2,ifc) = 0 ! Initial value immaterial
1334  END DO ! ifc
1335 
1336  ALLOCATE(pgrid%f2vTemp(4,nfacesint),stat=errorflag) ! NOTE nFacesInt
1337  global%error = errorflag
1338  IF ( global%error /= err_none ) THEN
1339  CALL errorstop(global,err_allocate,__line__,'pGrid%f2vTemp')
1340  END IF ! global%error
1341 
1342  DO ifc = 1,nfacesint ! Explicit loop because of Frost
1343  pgrid%f2vTemp(1,ifc) = 0 ! Initial value immaterial
1344  pgrid%f2vTemp(2,ifc) = 0 ! Initial value immaterial
1345  pgrid%f2vTemp(3,ifc) = 0 ! Initial value immaterial
1346  pgrid%f2vTemp(4,ifc) = 0 ! Initial value immaterial
1347  END DO ! ifc
1348 
1349  fcntr = 0
1350 
1351  fkcntr(face_kind_aa,face_type_tri:face_type_quad) = 0
1352  fkcntr(face_kind_av,face_type_tri:face_type_quad) = 0
1353  fkcntr(face_kind_vv,face_type_tri:face_type_quad) = 0
1354  fkcntr(face_kind_vx,face_type_tri:face_type_quad) = 0
1355 
1356 ! ==============================================================================
1357 ! Loop over all faces
1358 ! ==============================================================================
1359 
1360  DO ifg = 1,pgrid%nFacesTot
1361 
1362 ! ------------------------------------------------------------------------------
1363 ! Determine face kind
1364 ! ------------------------------------------------------------------------------
1365 
1366  c1 = pgrid%f2c(1,ifg)
1367  c2 = pgrid%f2c(2,ifg)
1368 
1369  c1k = rflu_getglobalcellkind(global,pgrid,c1)
1370  c2k = rflu_getglobalcellkind(global,pgrid,c2)
1371  ifk = rflu_getfacekind(global,c1k,c2k)
1372 
1373 ! ------------------------------------------------------------------------------
1374 ! If not on boundary, add to face list
1375 ! ------------------------------------------------------------------------------
1376 
1377  IF ( (ifk /= face_kind_ab) .AND. (ifk /= face_kind_vb) ) THEN
1378  ict = pgrid%cellGlob2Loc(1,c1) ! cell type
1379  icl = pgrid%cellGlob2Loc(2,c1) ! local cell index
1380  ifl = pgrid%f2c(3,ifg) ! local face index
1381 
1382  SELECT CASE ( ict )
1383  CASE ( cell_type_tet ) ! Tetrahedral cell
1384  v1g = pgrid%tet2v(f2vtet(1,ifl),icl)
1385  v2g = pgrid%tet2v(f2vtet(2,ifl),icl)
1386  v3g = pgrid%tet2v(f2vtet(3,ifl),icl)
1387  v4g = vert_none
1388  CASE ( cell_type_hex ) ! Hexahedral cell
1389  v1g = pgrid%hex2v(f2vhex(1,ifl),icl)
1390  v2g = pgrid%hex2v(f2vhex(2,ifl),icl)
1391  v3g = pgrid%hex2v(f2vhex(3,ifl),icl)
1392  v4g = pgrid%hex2v(f2vhex(4,ifl),icl)
1393  CASE ( cell_type_pri ) ! Prismatic cell
1394  v1g = pgrid%pri2v(f2vpri(1,ifl),icl)
1395  v2g = pgrid%pri2v(f2vpri(2,ifl),icl)
1396  v3g = pgrid%pri2v(f2vpri(3,ifl),icl)
1397  v4g = vert_none
1398 
1399  IF ( f2vpri(4,ifl) /= vert_none ) THEN
1400  v4g = pgrid%pri2v(f2vpri(4,ifl),icl)
1401  END IF ! f2vPri
1402  CASE ( cell_type_pyr ) ! Pyramidal cell
1403  v1g = pgrid%pyr2v(f2vpyr(1,ifl),icl)
1404  v2g = pgrid%pyr2v(f2vpyr(2,ifl),icl)
1405  v3g = pgrid%pyr2v(f2vpyr(3,ifl),icl)
1406  v4g = vert_none
1407 
1408  IF ( f2vpyr(4,ifl) /= vert_none ) THEN
1409  v4g = pgrid%pyr2v(f2vpyr(4,ifl),icl)
1410  END IF ! f2vPyr
1411  CASE default
1412  CALL errorstop(global,err_reached_default,__line__)
1413  END SELECT ! cellType
1414 
1415 ! ------------------------------------------------------------------------------
1416 ! Insert into face list
1417 ! ------------------------------------------------------------------------------
1418 
1419  fkcntr(ifk,face_type_tri) = fkcntr(ifk,face_type_tri) + 1
1420  fcntr = fkcntr(ifk,face_type_tri) + fkoffs(ifk)
1421 
1422  pgrid%f2cTemp(1,fcntr) = pgrid%f2c(1,ifg)
1423  pgrid%f2cTemp(2,fcntr) = pgrid%f2c(2,ifg)
1424 
1425  pgrid%f2vTemp(1,fcntr) = v1g
1426  pgrid%f2vTemp(2,fcntr) = v2g
1427  pgrid%f2vTemp(3,fcntr) = v3g
1428  pgrid%f2vTemp(4,fcntr) = v4g
1429  END IF ! ifk
1430  END DO ! ifg
1431 
1432 ! ------------------------------------------------------------------------------
1433 ! Check for consistency
1434 ! ------------------------------------------------------------------------------
1435 
1436  fksum = fkcntr(face_kind_aa,face_type_tri) &
1437  + fkcntr(face_kind_av,face_type_tri) &
1438  + fkcntr(face_kind_vv,face_type_tri) &
1439  + fkcntr(face_kind_vx,face_type_tri)
1440 
1441  IF ( fksum /= nfacesint ) THEN ! Incorrect number of internal faces
1442  CALL errorstop(global,err_nfaces_wrong,__line__)
1443  END IF ! fkSum
1444 
1445  pgrid%nFacesTot = nfacesint
1446 
1447 ! ==============================================================================
1448 ! Deallocate arrays
1449 ! ==============================================================================
1450 
1451  DEALLOCATE(pgrid%f2c,stat=errorflag)
1452  global%error = errorflag
1453  IF ( global%error /= err_none ) THEN
1454  CALL errorstop(global,err_deallocate,__line__,'pGrid%f2c')
1455  END IF ! global%error
1456 
1457 ! ==============================================================================
1458 ! Reallocate arrays, copy from temporary arrays, and deallocate those
1459 ! ==============================================================================
1460 
1461  ALLOCATE(pgrid%f2c(2,pgrid%nFacesTot),stat=errorflag) ! NOTE nFacesInt
1462  IF ( global%error /= err_none ) THEN
1463  CALL errorstop(global,err_allocate,__line__,'pGrid%f2c')
1464  END IF ! global%error
1465 
1466  DO ifc = 1,pgrid%nFacesTot ! Explicit loop because of Frost
1467  pgrid%f2c(1,ifc) = pgrid%f2cTemp(1,ifc)
1468  pgrid%f2c(2,ifc) = pgrid%f2cTemp(2,ifc)
1469  END DO ! ifc
1470 
1471  DEALLOCATE(pgrid%f2cTemp,stat=errorflag)
1472  global%error = errorflag
1473  IF ( global%error /= err_none ) THEN
1474  CALL errorstop(global,err_deallocate,__line__,'pGrid%f2cTemp')
1475  END IF ! global%error
1476 
1477  ALLOCATE(pgrid%f2v(4,pgrid%nFacesTot),stat=errorflag) ! NOTE nFacesInt
1478  global%error = errorflag
1479  IF ( global%error /= err_none ) THEN
1480  CALL errorstop(global,err_allocate,__line__,'pGrid%f2v')
1481  END IF ! global%error
1482 
1483  DO ifc = 1,pgrid%nFacesTot ! Explicit loop because of Frost
1484  pgrid%f2v(1,ifc) = pgrid%f2vTemp(1,ifc)
1485  pgrid%f2v(2,ifc) = pgrid%f2vTemp(2,ifc)
1486  pgrid%f2v(3,ifc) = pgrid%f2vTemp(3,ifc)
1487  pgrid%f2v(4,ifc) = pgrid%f2vTemp(4,ifc)
1488  END DO ! ifc
1489 
1490  DEALLOCATE(pgrid%f2vTemp,stat=errorflag)
1491  global%error = errorflag
1492  IF ( global%error /= err_none ) THEN
1493  CALL errorstop(global,err_deallocate,__line__,'pGrid%f2vTemp')
1494  END IF ! global%error
1495 
1496 #ifdef CHECK_DATASTRUCT
1497 ! ******************************************************************************
1498 ! Data structure output for checking
1499 ! ******************************************************************************
1500 
1501  WRITE(stdout,'(A)') solver_name
1502  WRITE(stdout,'(A,1X,A)') solver_name,'### START CHECK OUTPUT ###'
1503 
1504  WRITE(stdout,'(A,1X,A)') solver_name,'Face neighbours and vertices:'
1505 
1506  DO ifc = 1,pgrid%nFacesTot
1507  WRITE(stdout,'(A,7(1X,I6))') solver_name,ifc,pgrid%f2c(1:2,ifc), &
1508  pgrid%f2v(1:4,ifc)
1509  END DO ! ifc
1510 
1511  WRITE(stdout,'(A,1X,A)') solver_name,'### END CHECK OUTPUT ###'
1512  WRITE(stdout,'(A)') solver_name
1513 
1514  IF ( pgrid%nPatches > 0 ) THEN
1515  WRITE(stdout,'(A,1X,A)') solver_name,'### START CHECK OUTPUT ###'
1516 
1517  DO ipatch = 1,pgrid%nPatches
1518  ppatch => pregion%patches(ipatch)
1519  WRITE(stdout,'(A,1X,A,1X,I6)') solver_name,'Patch:',ipatch
1520 
1521  WRITE(stdout,'(A,3X,A)') solver_name,'Face neighbour and vertices:'
1522 
1523  DO ifc = 1,ppatch%nBFacesTot
1524  WRITE(stdout,'(A,6(1X,I6))') solver_name,ifc,ppatch%bf2c(ifc), &
1525  ppatch%bf2v(1:4,ifc)
1526  END DO ! ifc
1527  END DO ! iPatch
1528 
1529  WRITE(stdout,'(A,1X,A)') solver_name,'### END CHECK OUTPUT ###'
1530  WRITE(stdout,'(A)') solver_name
1531  END IF ! pGrid%nPatches
1532 #endif
1533 
1534 ! ******************************************************************************
1535 ! End
1536 ! ******************************************************************************
1537 
1538  IF ( global%myProcid == masterproc .AND. &
1539  global%verbLevel >= verbose_high ) THEN
1540  WRITE(stdout,'(A,1X,A)') solver_name,'Building face list done.'
1541  END IF ! global%verbLevel
1542 
1543  CALL deregisterfunction(global)
1544 
1545  END SUBROUTINE rflu_buildfacelist
1546 
1547 
1548 
1549 
1550 
1551 
1552 
1553 
1554 ! ******************************************************************************
1555 !
1556 ! Purpose: Create actual-virtual-face-to-border list.
1557 !
1558 ! Description: None.
1559 !
1560 ! Input:
1561 ! pRegion Pointer to region
1562 !
1563 ! Output: None.
1564 !
1565 ! Notes: None.
1566 !
1567 ! ******************************************************************************
1568 
1569  SUBROUTINE rflu_createavface2borderlist(pRegion)
1570 
1571  IMPLICIT NONE
1572 
1573 ! ******************************************************************************
1574 ! Declarations and definitions
1575 ! ******************************************************************************
1576 
1577 ! ==============================================================================
1578 ! Arguments
1579 ! ==============================================================================
1580 
1581  TYPE(t_region), POINTER :: pregion
1582 
1583 ! ==============================================================================
1584 ! Locals
1585 ! ==============================================================================
1586 
1587  INTEGER :: errorflag
1588  TYPE(t_grid), POINTER :: pgrid
1589  TYPE(t_global), POINTER :: global
1590 
1591 ! ******************************************************************************
1592 ! Start
1593 ! ******************************************************************************
1594 
1595  global => pregion%global
1596 
1597  CALL registerfunction(global,'RFLU_CreateAVFace2BorderList',&
1598  'RFLU_ModFaceList.F90')
1599 
1600  IF ( global%myProcid == masterproc .AND. &
1601  global%verbLevel >= verbose_high ) THEN
1602  WRITE(stdout,'(A,1X,A)') solver_name,'Creating av-face-to-border list...'
1603  END IF ! global%verbLevel
1604 
1605 ! ******************************************************************************
1606 ! Nullify memory
1607 ! ******************************************************************************
1608 
1609  CALL rflu_nullifyavface2borderlist(pregion)
1610 
1611 ! ******************************************************************************
1612 ! Set grid pointer
1613 ! ******************************************************************************
1614 
1615  pgrid => pregion%grid
1616 
1617 ! ******************************************************************************
1618 ! Allocate memory
1619 ! ******************************************************************************
1620 
1621  ALLOCATE(pgrid%avf2b(2,pgrid%nFacesAV),stat=errorflag)
1622  global%error = errorflag
1623  IF ( global%error /= err_none ) THEN
1624  CALL errorstop(global,err_allocate,__line__,'pGrid%avf2b')
1625  END IF ! global%error
1626 
1627 ! ******************************************************************************
1628 ! End
1629 ! ******************************************************************************
1630 
1631  IF ( global%myProcid == masterproc .AND. &
1632  global%verbLevel >= verbose_high ) THEN
1633  WRITE(stdout,'(A,1X,A)') solver_name, &
1634  'Creating av-face-to-border list done.'
1635  END IF ! global%verbLevel
1636 
1637  CALL deregisterfunction(global)
1638 
1639  END SUBROUTINE rflu_createavface2borderlist
1640 
1641 
1642 
1643 
1644 
1645 
1646 
1647 ! ******************************************************************************
1648 !
1649 ! Purpose: Create actual-virtual-face-to-patch list.
1650 !
1651 ! Description: None.
1652 !
1653 ! Input:
1654 ! pRegion Pointer to region
1655 !
1656 ! Output: None.
1657 !
1658 ! Notes: None.
1659 !
1660 ! ******************************************************************************
1661 
1662  SUBROUTINE rflu_createavface2patchlist(pRegion)
1663 
1664  IMPLICIT NONE
1665 
1666 ! ******************************************************************************
1667 ! Declarations and definitions
1668 ! ******************************************************************************
1669 
1670 ! ==============================================================================
1671 ! Arguments
1672 ! ==============================================================================
1673 
1674  TYPE(t_region), POINTER :: pregion
1675 
1676 ! ==============================================================================
1677 ! Locals
1678 ! ==============================================================================
1679 
1680  INTEGER :: errorflag
1681  TYPE(t_grid), POINTER :: pgrid
1682  TYPE(t_global), POINTER :: global
1683 
1684 ! ******************************************************************************
1685 ! Start
1686 ! ******************************************************************************
1687 
1688  global => pregion%global
1689 
1690  CALL registerfunction(global,'RFLU_CreateAVFace2PatchList',&
1691  'RFLU_ModFaceList.F90')
1692 
1693  IF ( global%myProcid == masterproc .AND. &
1694  global%verbLevel >= verbose_high ) THEN
1695  WRITE(stdout,'(A,1X,A)') solver_name,'Creating av-face-to-patch list...'
1696  END IF ! global%verbLevel
1697 
1698 ! ******************************************************************************
1699 ! Nullify memory
1700 ! ******************************************************************************
1701 
1702  CALL rflu_nullifyavface2patchlist(pregion)
1703 
1704 ! ******************************************************************************
1705 ! Set grid pointer
1706 ! ******************************************************************************
1707 
1708  pgrid => pregion%grid
1709 
1710 ! ******************************************************************************
1711 ! Allocate memory
1712 ! ******************************************************************************
1713 
1714  ALLOCATE(pgrid%avf2p(pgrid%nFacesAV),stat=errorflag)
1715  global%error = errorflag
1716  IF ( global%error /= err_none ) THEN
1717  CALL errorstop(global,err_allocate,__line__,'pGrid%avf2p')
1718  END IF ! global%error
1719 
1720 ! ******************************************************************************
1721 ! End
1722 ! ******************************************************************************
1723 
1724  IF ( global%myProcid == masterproc .AND. &
1725  global%verbLevel >= verbose_high ) THEN
1726  WRITE(stdout,'(A,1X,A)') solver_name, &
1727  'Creating av-face-to-patch list done.'
1728  END IF ! global%verbLevel
1729 
1730  CALL deregisterfunction(global)
1731 
1732  END SUBROUTINE rflu_createavface2patchlist
1733 
1734 
1735 
1736 
1737 
1738 
1739 
1740 ! ******************************************************************************
1741 !
1742 ! Purpose: Create cell-to-face lists.
1743 !
1744 ! Description: None.
1745 !
1746 ! Input:
1747 ! pRegion Pointer to region
1748 !
1749 ! Output: None.
1750 !
1751 ! Notes: None.
1752 !
1753 ! ******************************************************************************
1754 
1755  SUBROUTINE rflu_createcell2facelist(pRegion)
1756 
1757  IMPLICIT NONE
1758 
1759 ! ******************************************************************************
1760 ! Declarations and definitions
1761 ! ******************************************************************************
1762 
1763 ! ==============================================================================
1764 ! Arguments
1765 ! ==============================================================================
1766 
1767  TYPE(t_region), POINTER :: pregion
1768 
1769 ! ==============================================================================
1770 ! Locals
1771 ! ==============================================================================
1772 
1773  INTEGER :: errorflag
1774  TYPE(t_grid), POINTER :: pgrid
1775  TYPE(t_global), POINTER :: global
1776 
1777 ! ******************************************************************************
1778 ! Start
1779 ! ******************************************************************************
1780 
1781  global => pregion%global
1782 
1783  CALL registerfunction(global,'RFLU_CreateCell2FaceList',&
1784  'RFLU_ModFaceList.F90')
1785 
1786  IF ( global%myProcid == masterproc .AND. &
1787  global%verbLevel >= verbose_high ) THEN
1788  WRITE(stdout,'(A,1X,A)') solver_name,'Creating cell-to-face list...'
1789  END IF ! global%verbLevel
1790 
1791 ! ******************************************************************************
1792 ! Nullify memory
1793 ! ******************************************************************************
1794 
1795  CALL rflu_nullifycell2facelist(pregion)
1796 
1797 ! ******************************************************************************
1798 ! Set grid pointer
1799 ! ******************************************************************************
1800 
1801  pgrid => pregion%grid
1802 
1803 ! ******************************************************************************
1804 ! Allocate memory
1805 ! ******************************************************************************
1806 
1807  IF ( pgrid%nTetsTot > 0 ) THEN
1808  ALLOCATE(pgrid%tet2f(2,4,pgrid%nTetsTot),stat=errorflag)
1809  global%error = errorflag
1810  IF ( global%error /= err_none ) THEN
1811  CALL errorstop(global,err_allocate,__line__,'pGrid%tet2f')
1812  END IF ! global%error
1813  END IF ! pGrid%nTetsTot
1814 
1815  IF ( pgrid%nHexsTot > 0 ) THEN
1816  ALLOCATE(pgrid%hex2f(2,6,pgrid%nHexsTot),stat=errorflag)
1817  global%error = errorflag
1818  IF ( global%error /= err_none ) THEN
1819  CALL errorstop(global,err_allocate,__line__,'pGrid%hex2f')
1820  END IF ! global%error
1821  END IF ! pGrid%nHexsTot
1822 
1823  IF ( pgrid%nPrisTot > 0 ) THEN
1824  ALLOCATE(pgrid%pri2f(2,5,pgrid%nPrisTot),stat=errorflag)
1825  global%error = errorflag
1826  IF ( global%error /= err_none ) THEN
1827  CALL errorstop(global,err_allocate,__line__,'pGrid%pri2f')
1828  END IF ! global%error
1829  END IF ! pGrid%nPrisTot
1830 
1831  IF ( pgrid%nPyrsTot > 0 ) THEN
1832  ALLOCATE(pgrid%pyr2f(2,5,pgrid%nPyrsTot),stat=errorflag)
1833  global%error = errorflag
1834  IF ( global%error /= err_none ) THEN
1835  CALL errorstop(global,err_allocate,__line__,'pGrid%pyr2f')
1836  END IF ! global%error
1837  END IF ! pGrid%nPyrsTot
1838 
1839 ! ******************************************************************************
1840 ! End
1841 ! ******************************************************************************
1842 
1843  IF ( global%myProcid == masterproc .AND. &
1844  global%verbLevel >= verbose_high ) THEN
1845  WRITE(stdout,'(A,1X,A)') solver_name,'Creating cell-to-face list done.'
1846  END IF ! global%verbLevel
1847 
1848  CALL deregisterfunction(global)
1849 
1850  END SUBROUTINE rflu_createcell2facelist
1851 
1852 
1853 
1854 
1855 
1856 
1857 
1858 
1859 
1860 
1861 ! ******************************************************************************
1862 !
1863 ! Purpose: Create face list.
1864 !
1865 ! Description: None.
1866 !
1867 ! Input:
1868 ! pRegion Pointer to region
1869 !
1870 ! Output: None.
1871 !
1872 ! Notes:
1873 ! 1. The bf2c and bf2v arrays are not created here, because it is more
1874 ! convenient to do so when creating the grid.
1875 !
1876 ! ******************************************************************************
1877 
1878  SUBROUTINE rflu_createfacelist(pRegion)
1879 
1880  IMPLICIT NONE
1881 
1882 ! ******************************************************************************
1883 ! Declarations and definitions
1884 ! ******************************************************************************
1885 
1886 ! ==============================================================================
1887 ! Arguments
1888 ! ==============================================================================
1889 
1890  TYPE(t_region), POINTER :: pregion
1891 
1892 ! ==============================================================================
1893 ! Locals
1894 ! ==============================================================================
1895 
1896  INTEGER :: errorflag,ifc,ipatch,nbfaces
1897  TYPE(t_grid), POINTER :: pgrid
1898  TYPE(t_global), POINTER :: global
1899  TYPE(t_patch), POINTER :: ppatch
1900 
1901 ! ******************************************************************************
1902 ! Start
1903 ! ******************************************************************************
1904 
1905  global => pregion%global
1906 
1907  CALL registerfunction(global,'RFLU_CreateFaceList',&
1908  'RFLU_ModFaceList.F90')
1909 
1910  IF ( global%myProcid == masterproc .AND. &
1911  global%verbLevel >= verbose_high ) THEN
1912  WRITE(stdout,'(A,1X,A)') solver_name,'Creating face list...'
1913  END IF ! global%verbLevel
1914 
1915 ! ******************************************************************************
1916 ! Nullify memory
1917 ! ******************************************************************************
1918 
1919  CALL rflu_nullifyfacelist(pregion)
1920 
1921 ! ******************************************************************************
1922 ! Set grid pointer
1923 ! ******************************************************************************
1924 
1925  pgrid => pregion%grid
1926 
1927 ! ******************************************************************************
1928 ! Estimate number of faces for allocation of hash table. The estimation
1929 ! formula shown below assumes that boundary effects are negligible, so for
1930 ! very small grids, where boundary faces dominate the total number of faces ,
1931 ! need a kludge. Kludge is also needed if running in parallel and boundary
1932 ! faces become a significant fraction of total number of faces, and when
1933 ! running code with periodic hack, where the missing boundaries mean that
1934 ! number of faces is underestimated.
1935 ! ******************************************************************************
1936 
1937  nbfaces = 0
1938 
1939  DO ipatch = 1,pgrid%nPatches
1940  ppatch => pregion%patches(ipatch)
1941 
1942  nbfaces = nbfaces + ppatch%nBTrisTot + ppatch%nBQuadsTot
1943  END DO ! iPatch
1944 
1945 ! ==============================================================================
1946 ! Estimate total number of faces
1947 ! ==============================================================================
1948 
1949  pgrid%nFacesEst = nbfaces + 2*pgrid%nTetsTot + 3*pgrid%nHexsTot &
1950  + 5*pgrid%nPrisTot/2
1951 
1952  IF ( nbfaces/REAL(pGrid%nFacesEst,KIND=RFREAL) > 0.8_rfreal .OR. &
1953  nbfaces/REAL(pGrid%nFacesEst,KIND=RFREAL) < 0.3_rfreal ) then
1954  pgrid%nFacesEst = 2*pgrid%nFacesEst
1955 
1956  IF ( global%myProcid == masterproc .AND. &
1957  global%verbLevel >= verbose_high ) THEN
1958  WRITE(stdout,'(A,3X,A,1X,A)') solver_name,'Corrected estimate', &
1959  'of number of faces.'
1960  END IF ! global%verbLevel
1961  END IF ! nBFaces
1962 
1963  IF ( global%myProcid == masterproc .AND. &
1964  global%verbLevel >= verbose_high ) THEN
1965  WRITE(stdout,'(A,3X,A,3X,I9)') solver_name,'Estimated number of '// &
1966  'faces: ',pgrid%nFacesEst
1967  END IF ! global%verbLevel
1968 
1969 ! ******************************************************************************
1970 ! Allocate memory for interior faces
1971 ! ******************************************************************************
1972 
1973  ALLOCATE(pgrid%f2c(4,pgrid%nFacesEst),stat=errorflag)
1974  global%error = errorflag
1975  IF ( global%error /= err_none ) THEN
1976  CALL errorstop(global,err_allocate,__line__,'pGrid%f2c')
1977  END IF ! global%error
1978 
1979  DO ifc = 1,pgrid%nFacesEst
1980  pgrid%f2c(1,ifc) = cell_type_ext ! Initial value NOT immaterial
1981  pgrid%f2c(2,ifc) = cell_type_ext ! Initial value NOT immaterial
1982  pgrid%f2c(3,ifc) = 0 ! Initial value immaterial
1983  pgrid%f2c(4,ifc) = 0 ! Initial value immaterial
1984  END DO ! ifc
1985 
1986  ALLOCATE(pgrid%f2v(3,pgrid%nFacesEst),stat=errorflag)
1987  global%error = errorflag
1988  IF ( global%error /= err_none ) THEN
1989  CALL errorstop(global,err_allocate,__line__,'pGrid%f2v')
1990  END IF ! global%error
1991 
1992  DO ifc = 1,pgrid%nFacesEst
1993  pgrid%f2v(1,ifc) = vert_none ! Initial value NOT immaterial
1994  pgrid%f2v(2,ifc) = vert_none ! Initial value NOT immaterial
1995  pgrid%f2v(3,ifc) = vert_none ! Initial value NOT immaterial
1996  END DO ! ifc
1997 
1998 ! ******************************************************************************
1999 ! Allocate memory for boundary faces
2000 ! ******************************************************************************
2001 
2002  DO ipatch = 1,pgrid%nPatches
2003  ppatch => pregion%patches(ipatch)
2004 
2005  ALLOCATE(ppatch%bf2c(ppatch%nBFacesTot),stat=errorflag)
2006  global%error = errorflag
2007  IF ( global%error /= err_none ) THEN
2008  CALL errorstop(global,err_allocate,__line__,'pPatch%bf2c')
2009  END IF ! global%error
2010 
2011  ALLOCATE(ppatch%bf2v(4,ppatch%nBFacesTot),stat=errorflag)
2012  global%error = errorflag
2013  IF ( global%error /= err_none ) THEN
2014  CALL errorstop(global,err_allocate,__line__,'pPatch%bf2v')
2015  END IF ! global%error
2016 
2017  DO ifc = 1,ppatch%nBFacesTot
2018  ppatch%bf2c(ifc) = 0 ! Initial value immaterial, not used
2019  ppatch%bf2v(1,ifc) = vert_none ! Initial value NOT immaterial
2020  ppatch%bf2v(2,ifc) = vert_none ! Initial value NOT immaterial
2021  ppatch%bf2v(3,ifc) = vert_none ! Initial value NOT immaterial
2022  ppatch%bf2v(4,ifc) = vert_none ! Initial value NOT immaterial
2023  END DO ! ifc
2024  END DO ! iPatch
2025 
2026 ! ******************************************************************************
2027 ! End
2028 ! ******************************************************************************
2029 
2030  IF ( global%myProcid == masterproc .AND. &
2031  global%verbLevel >= verbose_high ) THEN
2032  WRITE(stdout,'(A,1X,A)') solver_name,'Creating face list done.'
2033  END IF ! global%verbLevel
2034 
2035  CALL deregisterfunction(global)
2036 
2037  END SUBROUTINE rflu_createfacelist
2038 
2039 
2040 
2041 
2042 
2043 
2044 ! ******************************************************************************
2045 !
2046 ! Purpose: Destroy actual-virtual-face-to-border list.
2047 !
2048 ! Description: None.
2049 !
2050 ! Input:
2051 ! pRegion Pointer to region
2052 !
2053 ! Output: None.
2054 !
2055 ! Notes: None.
2056 !
2057 ! ******************************************************************************
2058 
2059  SUBROUTINE rflu_destroyavface2borderlist(pRegion)
2060 
2061 ! ******************************************************************************
2062 ! Declarations and definitions
2063 ! ******************************************************************************
2064 
2065 ! ==============================================================================
2066 ! Arguments
2067 ! ==============================================================================
2068 
2069  TYPE(t_region), POINTER :: pregion
2070 
2071 ! ==============================================================================
2072 ! Locals
2073 ! ==============================================================================
2074 
2075  INTEGER :: errorflag
2076  TYPE(t_grid), POINTER :: pgrid
2077  TYPE(t_global), POINTER :: global
2078 
2079 ! ******************************************************************************
2080 ! Start
2081 ! ******************************************************************************
2082 
2083  global => pregion%global
2084 
2085  CALL registerfunction(global,'RFLU_DestroyAVFace2BorderList',&
2086  'RFLU_ModFaceList.F90')
2087 
2088  IF ( global%myProcid == masterproc .AND. &
2089  global%verbLevel >= verbose_high ) THEN
2090  WRITE(stdout,'(A,1X,A)') solver_name, &
2091  'Destroying av-face-to-border list...'
2092  END IF ! global%verbLevel
2093 
2094  pgrid => pregion%grid
2095 
2096 ! ******************************************************************************
2097 ! Deallocate memory
2098 ! ******************************************************************************
2099 
2100  DEALLOCATE(pgrid%avf2b,stat=errorflag)
2101  global%error = errorflag
2102  IF ( global%error /= err_none ) THEN
2103  CALL errorstop(global,err_deallocate,__line__,'pGrid%avf2b')
2104  END IF ! global%error
2105 
2106 ! ******************************************************************************
2107 ! Nullify memory
2108 ! ******************************************************************************
2109 
2110  CALL rflu_nullifyavface2borderlist(pregion)
2111 
2112 ! ******************************************************************************
2113 ! End
2114 ! ******************************************************************************
2115 
2116  IF ( global%myProcid == masterproc .AND. &
2117  global%verbLevel >= verbose_high ) THEN
2118  WRITE(stdout,'(A,1X,A)') solver_name, &
2119  'Destroying av-face-to-border list done.'
2120  END IF ! global%verbLevel
2121 
2122  CALL deregisterfunction(global)
2123 
2124  END SUBROUTINE rflu_destroyavface2borderlist
2125 
2126 
2127 
2128 
2129 
2130 
2131 
2132 ! ******************************************************************************
2133 !
2134 ! Purpose: Destroy actual-virtual-face-to-patch list.
2135 !
2136 ! Description: None.
2137 !
2138 ! Input:
2139 ! pRegion Pointer to region
2140 !
2141 ! Output: None.
2142 !
2143 ! Notes: None.
2144 !
2145 ! ******************************************************************************
2146 
2147  SUBROUTINE rflu_destroyavface2patchlist(pRegion)
2148 
2149 ! ******************************************************************************
2150 ! Declarations and definitions
2151 ! ******************************************************************************
2152 
2153 ! ==============================================================================
2154 ! Arguments
2155 ! ==============================================================================
2156 
2157  TYPE(t_region), POINTER :: pregion
2158 
2159 ! ==============================================================================
2160 ! Locals
2161 ! ==============================================================================
2162 
2163  INTEGER :: errorflag
2164  TYPE(t_grid), POINTER :: pgrid
2165  TYPE(t_global), POINTER :: global
2166 
2167 ! ******************************************************************************
2168 ! Start
2169 ! ******************************************************************************
2170 
2171  global => pregion%global
2172 
2173  CALL registerfunction(global,'RFLU_DestroyAVFace2PatchList',&
2174  'RFLU_ModFaceList.F90')
2175 
2176  IF ( global%myProcid == masterproc .AND. &
2177  global%verbLevel >= verbose_high ) THEN
2178  WRITE(stdout,'(A,1X,A)') solver_name, &
2179  'Destroying av-face-to-patch list...'
2180  END IF ! global%verbLevel
2181 
2182  pgrid => pregion%grid
2183 
2184 ! ******************************************************************************
2185 ! Deallocate memory
2186 ! ******************************************************************************
2187 
2188  DEALLOCATE(pgrid%avf2p,stat=errorflag)
2189  global%error = errorflag
2190  IF ( global%error /= err_none ) THEN
2191  CALL errorstop(global,err_deallocate,__line__,'pGrid%avf2p')
2192  END IF ! global%error
2193 
2194 ! ******************************************************************************
2195 ! Nullify memory
2196 ! ******************************************************************************
2197 
2198  CALL rflu_nullifyavface2patchlist(pregion)
2199 
2200 ! ******************************************************************************
2201 ! End
2202 ! ******************************************************************************
2203 
2204  IF ( global%myProcid == masterproc .AND. &
2205  global%verbLevel >= verbose_high ) THEN
2206  WRITE(stdout,'(A,1X,A)') solver_name, &
2207  'Destroying av-face-to-patch list done.'
2208  END IF ! global%verbLevel
2209 
2210  CALL deregisterfunction(global)
2211 
2212  END SUBROUTINE rflu_destroyavface2patchlist
2213 
2214 
2215 
2216 
2217 
2218 ! ******************************************************************************
2219 !
2220 ! Purpose: Destroy cell-to-face lists.
2221 !
2222 ! Description: None.
2223 !
2224 ! Input:
2225 ! pRegion Pointer to region
2226 !
2227 ! Output: None.
2228 !
2229 ! Notes: None.
2230 !
2231 ! ******************************************************************************
2232 
2233  SUBROUTINE rflu_destroycell2facelist(pRegion)
2234 
2235 ! ******************************************************************************
2236 ! Declarations and definitions
2237 ! ******************************************************************************
2238 
2239 ! ==============================================================================
2240 ! Arguments
2241 ! ==============================================================================
2242 
2243  TYPE(t_region), POINTER :: pregion
2244 
2245 ! ==============================================================================
2246 ! Locals
2247 ! ==============================================================================
2248 
2249  INTEGER :: errorflag
2250  TYPE(t_grid), POINTER :: pgrid
2251  TYPE(t_global), POINTER :: global
2252 
2253 ! ******************************************************************************
2254 ! Start
2255 ! ******************************************************************************
2256 
2257  global => pregion%global
2258 
2259  CALL registerfunction(global,'RFLU_DestroyCell2FaceList',&
2260  'RFLU_ModFaceList.F90')
2261 
2262  IF ( global%myProcid == masterproc .AND. &
2263  global%verbLevel >= verbose_high ) THEN
2264  WRITE(stdout,'(A,1X,A)') solver_name,'Destroying cell-to-face list...'
2265  END IF ! global%verbLevel
2266 
2267  pgrid => pregion%grid
2268 
2269 ! ******************************************************************************
2270 ! Deallocate memory
2271 ! ******************************************************************************
2272 
2273  IF ( pgrid%nTetsTot > 0 ) THEN
2274  DEALLOCATE(pgrid%tet2f,stat=errorflag)
2275  global%error = errorflag
2276  IF ( global%error /= err_none ) THEN
2277  CALL errorstop(global,err_deallocate,__line__,'pGrid%tet2f')
2278  END IF ! global%error
2279  END IF ! pGrid%nTetsTot
2280 
2281  IF ( pgrid%nHexsTot > 0 ) THEN
2282  DEALLOCATE(pgrid%hex2f,stat=errorflag)
2283  global%error = errorflag
2284  IF ( global%error /= err_none ) THEN
2285  CALL errorstop(global,err_deallocate,__line__,'pGrid%hex2f')
2286  END IF ! global%error
2287  END IF ! pGrid%nHexsTot
2288 
2289  IF ( pgrid%nPrisTot > 0 ) THEN
2290  DEALLOCATE(pgrid%pri2f,stat=errorflag)
2291  global%error = errorflag
2292  IF ( global%error /= err_none ) THEN
2293  CALL errorstop(global,err_deallocate,__line__,'pGrid%pri2f')
2294  END IF ! global%error
2295  END IF ! pGrid%nPrisTot
2296 
2297  IF ( pgrid%nPyrsTot > 0 ) THEN
2298  DEALLOCATE(pgrid%pyr2f,stat=errorflag)
2299  global%error = errorflag
2300  IF ( global%error /= err_none ) THEN
2301  CALL errorstop(global,err_deallocate,__line__,'pGrid%pyr2f')
2302  END IF ! global%error
2303  END IF ! pGrid%nPyrsTot
2304 
2305 ! ******************************************************************************
2306 ! Nullify memory
2307 ! ******************************************************************************
2308 
2309  CALL rflu_nullifycell2facelist(pregion)
2310 
2311 ! ******************************************************************************
2312 ! End
2313 ! ******************************************************************************
2314 
2315  IF ( global%myProcid == masterproc .AND. &
2316  global%verbLevel >= verbose_high ) THEN
2317  WRITE(stdout,'(A,1X,A)') solver_name, &
2318  'Destroying cell-to-face list done.'
2319  END IF ! global%verbLevel
2320 
2321  CALL deregisterfunction(global)
2322 
2323  END SUBROUTINE rflu_destroycell2facelist
2324 
2325 
2326 
2327 
2328 
2329 
2330 
2331 
2332 
2333 
2334 
2335 ! ******************************************************************************
2336 !
2337 ! Purpose: Destroy face list.
2338 !
2339 ! Description: None.
2340 !
2341 ! Input:
2342 ! pRegion Pointer to region
2343 !
2344 ! Output: None.
2345 !
2346 ! Notes: None.
2347 !
2348 ! ******************************************************************************
2349 
2350  SUBROUTINE rflu_destroyfacelist(pRegion)
2351 
2352  IMPLICIT NONE
2353 
2354 ! ******************************************************************************
2355 ! Declarations and definitions
2356 ! ******************************************************************************
2357 
2358 ! ==============================================================================
2359 ! Arguments
2360 ! ==============================================================================
2361 
2362  TYPE(t_region), POINTER :: pregion
2363 
2364 ! ==============================================================================
2365 ! Locals
2366 ! ==============================================================================
2367 
2368  INTEGER :: errorflag,ipatch
2369  TYPE(t_grid), POINTER :: pgrid
2370  TYPE(t_patch), POINTER :: ppatch
2371  TYPE(t_global), POINTER :: global
2372 
2373 ! ******************************************************************************
2374 ! Start
2375 ! ******************************************************************************
2376 
2377  global => pregion%global
2378 
2379  CALL registerfunction(global,'RFLU_DestroyFaceList',&
2380  'RFLU_ModFaceList.F90')
2381 
2382  IF ( global%myProcid == masterproc .AND. &
2383  global%verbLevel >= verbose_high ) THEN
2384  WRITE(stdout,'(A,1X,A)') solver_name,'Destroying face lists...'
2385  END IF ! global%verbLevel
2386 
2387 ! ******************************************************************************
2388 ! Set grid pointer
2389 ! ******************************************************************************
2390 
2391  pgrid => pregion%grid
2392 
2393 ! ******************************************************************************
2394 ! Deallocate and nullify memory for interior face lists
2395 ! ******************************************************************************
2396 
2397  DEALLOCATE(pgrid%f2c,stat=errorflag)
2398  global%error = errorflag
2399  IF ( global%error /= err_none ) THEN
2400  CALL errorstop(global,err_deallocate,__line__,'pGrid%f2c')
2401  END IF ! global%error
2402 
2403  DEALLOCATE(pgrid%f2v,stat=errorflag)
2404  global%error = errorflag
2405  IF ( global%error /= err_none ) THEN
2406  CALL errorstop(global,err_deallocate,__line__,'pGrid%f2v')
2407  END IF ! global%error
2408 
2409 ! ******************************************************************************
2410 ! Deallocate and nullify memory for boundary face lists. NOTE need to set
2411 ! renumbering flag to false because when have symmetry or periodic patches
2412 ! need to recreate boundary face lists after adding virtual cells on patches
2413 ! and renumber them. If not reset here, lists do not get renumbered because
2414 ! will get skipped on account of flag already indicating renumbering having
2415 ! been done.
2416 ! ******************************************************************************
2417 
2418  DO ipatch = 1,pgrid%nPatches
2419  ppatch => pregion%patches(ipatch)
2420 
2421  DEALLOCATE(ppatch%bf2c,stat=errorflag)
2422  global%error = errorflag
2423  IF ( global%error /= err_none ) THEN
2424  CALL errorstop(global,err_deallocate,__line__,'pPatch%bf2c')
2425  END IF ! global%error
2426 
2427  DEALLOCATE(ppatch%bf2v,stat=errorflag)
2428  global%error = errorflag
2429  IF ( global%error /= err_none ) THEN
2430  CALL errorstop(global,err_deallocate,__line__,'pPatch%bf2v')
2431  END IF ! global%error
2432 
2433  ppatch%renumFlag = .false.
2434  END DO ! iPatch
2435 
2436 ! ******************************************************************************
2437 ! Nullify memory
2438 ! ******************************************************************************
2439 
2440  CALL rflu_nullifyfacelist(pregion)
2441 
2442 ! ******************************************************************************
2443 ! End
2444 ! ******************************************************************************
2445 
2446  IF ( global%myProcid == masterproc .AND. &
2447  global%verbLevel >= verbose_high ) THEN
2448  WRITE(stdout,'(A,1X,A)') solver_name,'Destroying face lists done.'
2449  END IF ! global%verbLevel
2450 
2451  CALL deregisterfunction(global)
2452 
2453  END SUBROUTINE rflu_destroyfacelist
2454 
2455 
2456 
2457 
2458 
2459 
2460 
2461 
2462 ! ******************************************************************************
2463 !
2464 ! Purpose: Get faces opposing given face.
2465 !
2466 ! Description: Compare vertices of given face with vertices of the faces of the
2467 ! two cells adjacent to given face.
2468 !
2469 ! Input:
2470 ! pRegion Pointer to region
2471 ! iPatch Patch index
2472 ! ifg Face index
2473 !
2474 ! Output:
2475 ! nFacesOpp Number of opposing faces
2476 ! oppFaceInfo Information about opposing faces
2477 !
2478 ! Notes:
2479 ! 1. There are in general two faces opposing a given face because two cells
2480 ! share a face. For faces of a boundary cell, there will be at most one
2481 ! opposing face.
2482 ! 2. The term opposing is used in the sense that two faces of a cell are
2483 ! opposing if they do not share any vertices. Hence tetrahedra or pyramids
2484 ! do not have any opposing faces.
2485 ! 3. The array oppFaceInfo contains two entries for each opposing face. The
2486 ! first entry is the patch number of the opposing face (0 if the opposing
2487 ! face is an internal face), and the second entry is the index of the face
2488 ! in the respective face list (internal face list or patch face list). If
2489 ! a given face has no or only one opposing face, the corresponding entries
2490 ! are set to OPP_FACE_NONE.
2491 ! 4. Looping down when determining whether global and local face are identical
2492 ! because missing vertex (which could have entry 0 depending on parameter
2493 ! VERT_NONE) would mean that could only check for future unnecessary
2494 ! comparisons for local face 2, so by starting from end can save one
2495 ! comparison.
2496 !
2497 ! ******************************************************************************
2498 
2499  SUBROUTINE rflu_getopposingfaces(pRegion,iPatch,ifg,nFacesOpp,faceOppInfo)
2500 
2501  USE modsortsearch
2502 
2503  USE rflu_modgrid
2504 
2505 ! ******************************************************************************
2506 ! Declarations and definitions
2507 ! ******************************************************************************
2508 
2509 ! ==============================================================================
2510 ! Arguments
2511 ! ==============================================================================
2512 
2513  INTEGER, INTENT(IN) :: ifg,ipatch
2514  INTEGER, INTENT(OUT) :: nfacesopp
2515  INTEGER, INTENT(OUT) :: faceoppinfo(2,2)
2516  TYPE(t_region), POINTER :: pregion
2517 
2518 ! ==============================================================================
2519 ! Locals
2520 ! ==============================================================================
2521 
2522  INTEGER :: ic,icg,icl,ict,ifgopp,ifl,iflopp,ipatchopp,ivl,ncells,term
2523  INTEGER :: c(2),f2vsort(4),f2vsort2(4)
2524  TYPE(t_global), POINTER :: global
2525  TYPE(t_grid), POINTER :: pgrid
2526  TYPE(t_patch), POINTER :: ppatch
2527 
2528 ! ******************************************************************************
2529 ! Start
2530 ! ******************************************************************************
2531 
2532  global => pregion%global
2533 
2534  CALL registerfunction(global,'RFLU_GetOpposingFaces',&
2535  'RFLU_ModFaceList.F90')
2536 
2537  pgrid => pregion%grid
2538 
2539 ! ******************************************************************************
2540 ! Set patch pointer if face is located on a patch
2541 ! ******************************************************************************
2542 
2543  IF ( ipatch > 0 ) THEN
2544  ppatch => pregion%patches(ipatch)
2545  ELSE IF ( ipatch == 0 ) THEN
2546  nullify(ppatch)
2547  ELSE
2548  CALL errorstop(global,err_reached_default,__line__)
2549  END IF ! iPatch
2550 
2551 ! ******************************************************************************
2552 ! Get cells sharing face and the vertices of face
2553 ! ******************************************************************************
2554 
2555  IF ( ipatch == 0 ) THEN ! Interior face
2556  ncells = 2
2557 
2558  c(1) = pgrid%f2c(1,ifg)
2559  c(2) = pgrid%f2c(2,ifg)
2560 
2561  DO ivl = 1,4
2562  f2vsort(ivl) = pgrid%f2v(ivl,ifg)
2563  END DO ! ivl
2564  ELSE IF ( ipatch > 0 ) THEN ! Boundary face
2565  ncells = 1
2566 
2567  c(1) = ppatch%bf2c(ifg)
2568 
2569  DO ivl = 1,4
2570  IF ( ppatch%bf2v(ivl,ifg) /= vert_none ) THEN
2571  f2vsort(ivl) = ppatch%bv(ppatch%bf2v(ivl,ifg))
2572  ELSE
2573  f2vsort(ivl) = vert_none
2574  END IF ! pPatch%bf2v
2575  END DO ! ivl
2576  END IF ! iPatch
2577 
2578  CALL quicksortinteger(f2vsort,4)
2579 
2580 ! ******************************************************************************
2581 ! Find local face of cells which matches vertices of face in question
2582 ! ******************************************************************************
2583 
2584  nfacesopp = 0
2585 
2586  DO ic = 1,ncells
2587  icg = c(ic)
2588 
2589  ict = pgrid%cellGlob2Loc(1,icg) ! cell type
2590  icl = pgrid%cellGlob2Loc(2,icg) ! local cell index
2591 
2592  ipatchopp = opp_face_none
2593  ifgopp = opp_face_none
2594 
2595 ! ==============================================================================
2596 ! Select cell type
2597 ! ==============================================================================
2598 
2599  SELECT CASE ( ict )
2600 
2601 ! ------------------------------------------------------------------------------
2602 ! Tetrahedron and pyramid (do not have opposing faces)
2603 ! ------------------------------------------------------------------------------
2604 
2605  CASE ( cell_type_tet, cell_type_pyr )
2606 
2607 ! ------------------------------------------------------------------------------
2608 ! Hexahedron
2609 ! ------------------------------------------------------------------------------
2610 
2611  CASE ( cell_type_hex )
2612  hexfaceloop: DO ifl = 1,6
2613  DO ivl = 1,4 ! Get vertices of local face
2614  f2vsort2(ivl) = pgrid%hex2v(f2vhex(ivl,ifl),icl)
2615  END DO ! ivl
2616 
2617  CALL quicksortinteger(f2vsort2,4)
2618 
2619  term = 0
2620 
2621  DO ivl = 4,1,-1 ! See note about looping down
2622  term = term + abs(f2vsort(ivl) - f2vsort2(ivl))
2623 
2624  IF ( term /= 0 ) THEN
2625  cycle hexfaceloop
2626  END IF ! term
2627  END DO ! ivl
2628 
2629  IF ( term == 0 ) THEN ! Found opposing face
2630  iflopp = f2fopphex(ifl)
2631 
2632  nfacesopp = nfacesopp + 1
2633 
2634  faceoppinfo(1,nfacesopp) = pgrid%hex2f(1,iflopp,icl)
2635  faceoppinfo(2,nfacesopp) = pgrid%hex2f(2,iflopp,icl)
2636 
2637  EXIT hexfaceloop
2638  END IF ! term
2639  END DO hexfaceloop
2640 
2641 ! ------------------------------------------------------------------------------
2642 ! Prism
2643 ! ------------------------------------------------------------------------------
2644 
2645  CASE ( cell_type_pri )
2646  prifaceloop: DO ifl = 1,5 ! Get vertices of local face
2647  DO ivl = 1,4
2648  IF ( f2vpri(ivl,ifl) /= vert_none ) THEN
2649  f2vsort2(ivl) = pgrid%pri2v(f2vpri(ivl,ifl),icl)
2650  ELSE
2651  f2vsort2(ivl) = vert_none
2652  END IF ! f2vPri
2653  END DO ! ivl
2654 
2655  CALL quicksortinteger(f2vsort2,4)
2656 
2657  term = 0
2658 
2659  DO ivl = 4,1,-1 ! See note about looping down
2660  term = term + abs(f2vsort(ivl) - f2vsort2(ivl))
2661 
2662  IF ( term /= 0 ) THEN
2663  cycle prifaceloop
2664  END IF ! term
2665  END DO ! ivl
2666 
2667  IF ( term == 0 ) THEN ! Found opposing face
2668  iflopp = f2fopppri(ifl)
2669 
2670  nfacesopp = nfacesopp + 1
2671 
2672  faceoppinfo(1,nfacesopp) = pgrid%pri2f(1,iflopp,icl)
2673  faceoppinfo(2,nfacesopp) = pgrid%pri2f(2,iflopp,icl)
2674 
2675  EXIT prifaceloop
2676  END IF ! term
2677  END DO prifaceloop
2678 
2679 ! ------------------------------------------------------------------------------
2680 ! Default
2681 ! ------------------------------------------------------------------------------
2682 
2683  CASE default
2684  CALL errorstop(global,err_reached_default,__line__)
2685  END SELECT ! ict
2686  END DO ! ic
2687 
2688 ! ******************************************************************************
2689 ! End
2690 ! ******************************************************************************
2691 
2692  CALL deregisterfunction(global)
2693 
2694  END SUBROUTINE rflu_getopposingfaces
2695 
2696 
2697 
2698 
2699 
2700 
2701 
2702 ! ******************************************************************************
2703 !
2704 ! Purpose: Insert face into cell-to-face list.
2705 !
2706 ! Description: None.
2707 !
2708 ! Input:
2709 ! pRegion Pointer to region data
2710 ! iPatch Patch index
2711 ! icg Cell index
2712 ! ifg Face index
2713 !
2714 ! Output: None.
2715 !
2716 ! Notes:
2717 ! 1. Looping down when determining whether global and local face are identical
2718 ! because missing vertex (which could have entry 0 depending on parameter
2719 ! VERT_NONE) would mean that could only check for future unnecessary
2720 ! comparisons for local face 2, so by starting from end can save one
2721 ! comparison.
2722 !
2723 ! ******************************************************************************
2724 
2725  SUBROUTINE rflu_insertintocell2facelist(pRegion,iPatch,icg,ifg)
2726 
2727  USE modsortsearch, ONLY: quicksortinteger
2728 
2729  USE rflu_modgrid
2730 
2731 ! ******************************************************************************
2732 ! Declarations and definitions
2733 ! ******************************************************************************
2734 
2735 ! ==============================================================================
2736 ! Arguments
2737 ! ==============================================================================
2738 
2739  INTEGER, INTENT(IN) :: icg,ifg,ipatch
2740  TYPE(t_region), POINTER :: pregion
2741 
2742 ! ==============================================================================
2743 ! Locals
2744 ! ==============================================================================
2745 
2746  INTEGER :: icl,ict,ifl,ivl,nfaces,term
2747  INTEGER, DIMENSION(4) :: f2vsort,f2vsort2
2748  TYPE(t_global), POINTER :: global
2749  TYPE(t_grid), POINTER :: pgrid
2750  TYPE(t_patch), POINTER :: ppatch
2751 
2752 ! ******************************************************************************
2753 ! Start
2754 ! ******************************************************************************
2755 
2756  global => pregion%global
2757 
2758  CALL registerfunction(global,'RFLU_InsertIntoCell2FaceList',&
2759  'RFLU_ModFaceList.F90')
2760 
2761  pgrid => pregion%grid
2762 
2763 ! ******************************************************************************
2764 ! Get cell type and local cell index
2765 ! ******************************************************************************
2766 
2767  ict = pgrid%cellGlob2Loc(1,icg) ! cell type
2768  icl = pgrid%cellGlob2Loc(2,icg) ! local cell index
2769 
2770 ! ******************************************************************************
2771 ! Get vertices of face
2772 ! ******************************************************************************
2773 
2774  IF ( ipatch == 0 ) THEN
2775  nullify(ppatch)
2776 
2777  DO ivl = 1,4
2778  f2vsort(ivl) = pgrid%f2v(ivl,ifg)
2779  END DO ! ivl
2780  ELSE IF ( ipatch > 0 ) THEN
2781  ppatch => pregion%patches(ipatch)
2782 
2783  DO ivl = 1,4
2784  IF ( ppatch%bf2v(ivl,ifg) /= vert_none ) THEN
2785  f2vsort(ivl) = ppatch%bv(ppatch%bf2v(ivl,ifg))
2786  ELSE
2787  f2vsort(ivl) = vert_none
2788  END IF ! pPatch%bf2v
2789  END DO ! ivl
2790  ELSE
2791  CALL errorstop(global,err_reached_default,__line__)
2792  END IF ! iPatch
2793 
2794  CALL quicksortinteger(f2vsort,4)
2795 
2796 ! ******************************************************************************
2797 ! Insert face into cell-to-face list
2798 ! ******************************************************************************
2799 
2800  SELECT CASE ( ict )
2801 
2802 ! ==============================================================================
2803 ! Tetrahedra
2804 ! ==============================================================================
2805 
2806  CASE ( cell_type_tet )
2807  nfaces = 4
2808 
2809  tetfaceloop: DO ifl = 1,nfaces
2810  DO ivl = 1,4
2811  IF ( f2vtet(ivl,ifl) /= vert_none ) THEN
2812  f2vsort2(ivl) = pgrid%tet2v(f2vtet(ivl,ifl),icl)
2813  ELSE
2814  f2vsort2(ivl) = vert_none
2815  END IF ! f2vTet
2816  END DO ! ivl
2817 
2818  CALL quicksortinteger(f2vsort2,4)
2819 
2820  term = 0
2821 
2822  tetfacevertloop: DO ivl = 4,1,-1 ! See note about looping down
2823  term = term + abs(f2vsort(ivl) - f2vsort2(ivl))
2824 
2825  IF ( term /= 0 ) THEN
2826  EXIT tetfacevertloop
2827  END IF ! term
2828  END DO tetfacevertloop
2829 
2830  IF ( term == 0 ) THEN
2831  IF ( pgrid%tet2f(1,ifl,icl) == c2f_init .OR. &
2832  pgrid%tet2f(2,ifl,icl) == c2f_init ) THEN
2833  pgrid%tet2f(1,ifl,icl) = ipatch
2834  pgrid%tet2f(2,ifl,icl) = ifg
2835 
2836  EXIT tetfaceloop
2837  ELSE
2838  CALL errorstop(global,err_reached_default,__line__)
2839  END IF ! pGrid%tet2f
2840  ELSE
2841  IF ( ifl == nfaces ) THEN
2842  CALL errorstop(global,err_reached_default,__line__)
2843  END IF ! ifl
2844  END IF ! term
2845  END DO tetfaceloop
2846 
2847 ! ==============================================================================
2848 ! Hexahedra
2849 ! ==============================================================================
2850 
2851  CASE ( cell_type_hex )
2852  nfaces = 6
2853 
2854  hexfaceloop: DO ifl = 1,nfaces
2855  DO ivl = 1,4
2856  f2vsort2(ivl) = pgrid%hex2v(f2vhex(ivl,ifl),icl)
2857  END DO ! ivl
2858 
2859  CALL quicksortinteger(f2vsort2,4)
2860 
2861  term = 0
2862 
2863  hexfacevertloop: DO ivl = 4,1,-1 ! See note about looping down
2864  term = term + abs(f2vsort(ivl) - f2vsort2(ivl))
2865 
2866  IF ( term /= 0 ) THEN
2867  EXIT hexfacevertloop
2868  END IF ! term
2869  END DO hexfacevertloop
2870 
2871  IF ( term == 0 ) THEN
2872  IF ( pgrid%hex2f(1,ifl,icl) == c2f_init .OR. &
2873  pgrid%hex2f(2,ifl,icl) == c2f_init ) THEN
2874  pgrid%hex2f(1,ifl,icl) = ipatch
2875  pgrid%hex2f(2,ifl,icl) = ifg
2876 
2877  EXIT hexfaceloop
2878  ELSE
2879  CALL errorstop(global,err_reached_default,__line__)
2880  END IF ! pGrid%hex2f
2881  ELSE
2882  IF ( ifl == nfaces ) THEN
2883  CALL errorstop(global,err_reached_default,__line__)
2884  END IF ! ifl
2885  END IF ! term
2886  END DO hexfaceloop
2887 
2888 ! ==============================================================================
2889 ! Prisms
2890 ! ==============================================================================
2891 
2892  CASE ( cell_type_pri )
2893  nfaces = 5
2894 
2895  prifaceloop: DO ifl = 1,nfaces
2896  DO ivl = 1,4
2897  IF ( f2vpri(ivl,ifl) /= vert_none ) THEN
2898  f2vsort2(ivl) = pgrid%pri2v(f2vpri(ivl,ifl),icl)
2899  ELSE
2900  f2vsort2(ivl) = vert_none
2901  END IF ! f2vPri
2902  END DO ! ivl
2903 
2904  CALL quicksortinteger(f2vsort2,4)
2905 
2906  term = 0
2907 
2908  prifacevertloop: DO ivl = 4,1,-1 ! See note about looping down
2909  term = term + abs(f2vsort(ivl) - f2vsort2(ivl))
2910 
2911  IF ( term /= 0 ) THEN
2912  EXIT prifacevertloop
2913  END IF ! term
2914  END DO prifacevertloop
2915 
2916  IF ( term == 0 ) THEN
2917  IF ( pgrid%pri2f(1,ifl,icl) == c2f_init .OR. &
2918  pgrid%pri2f(2,ifl,icl) == c2f_init ) THEN
2919  pgrid%pri2f(1,ifl,icl) = ipatch
2920  pgrid%pri2f(2,ifl,icl) = ifg
2921 
2922  EXIT prifaceloop
2923  ELSE
2924  CALL errorstop(global,err_reached_default,__line__)
2925  END IF ! pGrid%pri2f
2926  ELSE
2927  IF ( ifl == nfaces ) THEN
2928  CALL errorstop(global,err_reached_default,__line__)
2929  END IF ! ifl
2930  END IF ! term
2931  END DO prifaceloop
2932 
2933 ! ==============================================================================
2934 ! Pyramids
2935 ! ==============================================================================
2936 
2937  CASE ( cell_type_pyr )
2938  nfaces = 5
2939 
2940  pyrfaceloop: DO ifl = 1,nfaces
2941  DO ivl = 1,4
2942  IF ( f2vpyr(ivl,ifl) /= vert_none ) THEN
2943  f2vsort2(ivl) = pgrid%pyr2v(f2vpyr(ivl,ifl),icl)
2944  ELSE
2945  f2vsort2(ivl) = vert_none
2946  END IF ! f2vPyr
2947  END DO ! ivl
2948 
2949  CALL quicksortinteger(f2vsort2,4)
2950 
2951  term = 0
2952 
2953  pyrfacevertloop: DO ivl = 4,1,-1 ! See note about looping down
2954  term = term + abs(f2vsort(ivl) - f2vsort2(ivl))
2955 
2956  IF ( term /= 0 ) THEN
2957  EXIT pyrfacevertloop
2958  END IF ! term
2959  END DO pyrfacevertloop
2960 
2961  IF ( term == 0 ) THEN
2962  IF ( pgrid%pyr2f(1,ifl,icl) == c2f_init .OR. &
2963  pgrid%pyr2f(2,ifl,icl) == c2f_init ) THEN
2964  pgrid%pyr2f(1,ifl,icl) = ipatch
2965  pgrid%pyr2f(2,ifl,icl) = ifg
2966 
2967  EXIT pyrfaceloop
2968  ELSE
2969  CALL errorstop(global,err_reached_default,__line__)
2970  END IF ! pGrid%pyr2f
2971  ELSE
2972  IF ( ifl == nfaces ) THEN
2973  CALL errorstop(global,err_reached_default,__line__)
2974  END IF ! ifl
2975  END IF ! term
2976  END DO pyrfaceloop
2977 
2978 ! ==============================================================================
2979 ! Default
2980 ! ==============================================================================
2981 
2982  CASE default
2983  CALL errorstop(global,err_reached_default,__line__)
2984  END SELECT ! ict
2985 
2986 ! ******************************************************************************
2987 ! End
2988 ! ******************************************************************************
2989 
2990  CALL deregisterfunction(global)
2991 
2992  END SUBROUTINE rflu_insertintocell2facelist
2993 
2994 
2995 
2996 
2997 
2998 
2999 
3000 ! ******************************************************************************
3001 !
3002 ! Purpose: Nullify actual-virtual-face-to-border lists.
3003 !
3004 ! Description: None.
3005 !
3006 ! Input:
3007 ! pRegion Pointer to region
3008 !
3009 ! Output: None.
3010 !
3011 ! Notes: None.
3012 !
3013 ! ******************************************************************************
3014 
3015  SUBROUTINE rflu_nullifyavface2borderlist(pRegion)
3016 
3017  IMPLICIT NONE
3018 
3019 ! ******************************************************************************
3020 ! Declarations and definitions
3021 ! ******************************************************************************
3022 
3023 ! ==============================================================================
3024 ! Arguments
3025 ! ==============================================================================
3026 
3027  TYPE(t_region), POINTER :: pregion
3028 
3029 ! ==============================================================================
3030 ! Locals
3031 ! ==============================================================================
3032 
3033  INTEGER :: errorflag
3034  TYPE(t_grid), POINTER :: pgrid
3035  TYPE(t_global), POINTER :: global
3036 
3037 ! ******************************************************************************
3038 ! Start
3039 ! ******************************************************************************
3040 
3041  global => pregion%global
3042 
3043  CALL registerfunction(global,'RFLU_NullifyAVFace2BorderList',&
3044  'RFLU_ModFaceList.F90')
3045 
3046  pgrid => pregion%grid
3047 
3048 ! ******************************************************************************
3049 ! Nullify memory
3050 ! ******************************************************************************
3051 
3052  nullify(pgrid%avf2b)
3053 
3054 ! ******************************************************************************
3055 ! End
3056 ! ******************************************************************************
3057 
3058  CALL deregisterfunction(global)
3059 
3060  END SUBROUTINE rflu_nullifyavface2borderlist
3061 
3062 
3063 
3064 
3065 
3066 
3067 ! ******************************************************************************
3068 !
3069 ! Purpose: Nullify actual-virtual-face-to-patch lists.
3070 !
3071 ! Description: None.
3072 !
3073 ! Input:
3074 ! pRegion Pointer to region
3075 !
3076 ! Output: None.
3077 !
3078 ! Notes: None.
3079 !
3080 ! ******************************************************************************
3081 
3082  SUBROUTINE rflu_nullifyavface2patchlist(pRegion)
3083 
3084  IMPLICIT NONE
3085 
3086 ! ******************************************************************************
3087 ! Declarations and definitions
3088 ! ******************************************************************************
3089 
3090 ! ==============================================================================
3091 ! Arguments
3092 ! ==============================================================================
3093 
3094  TYPE(t_region), POINTER :: pregion
3095 
3096 ! ==============================================================================
3097 ! Locals
3098 ! ==============================================================================
3099 
3100  INTEGER :: errorflag
3101  TYPE(t_grid), POINTER :: pgrid
3102  TYPE(t_global), POINTER :: global
3103 
3104 ! ******************************************************************************
3105 ! Start
3106 ! ******************************************************************************
3107 
3108  global => pregion%global
3109 
3110  CALL registerfunction(global,'RFLU_NullifyAVFace2PatchList',&
3111  'RFLU_ModFaceList.F90')
3112 
3113  pgrid => pregion%grid
3114 
3115 ! ******************************************************************************
3116 ! Nullify memory
3117 ! ******************************************************************************
3118 
3119  nullify(pgrid%avf2p)
3120 
3121 ! ******************************************************************************
3122 ! End
3123 ! ******************************************************************************
3124 
3125  CALL deregisterfunction(global)
3126 
3127  END SUBROUTINE rflu_nullifyavface2patchlist
3128 
3129 
3130 
3131 
3132 
3133 
3134 ! ******************************************************************************
3135 !
3136 ! Purpose: Nullify cell-to-face lists.
3137 !
3138 ! Description: None.
3139 !
3140 ! Input:
3141 ! pRegion Pointer to region
3142 !
3143 ! Output: None.
3144 !
3145 ! Notes: None.
3146 !
3147 ! ******************************************************************************
3148 
3149  SUBROUTINE rflu_nullifycell2facelist(pRegion)
3150 
3151  IMPLICIT NONE
3152 
3153 ! ******************************************************************************
3154 ! Declarations and definitions
3155 ! ******************************************************************************
3156 
3157 ! ==============================================================================
3158 ! Arguments
3159 ! ==============================================================================
3160 
3161  TYPE(t_region), POINTER :: pregion
3162 
3163 ! ==============================================================================
3164 ! Locals
3165 ! ==============================================================================
3166 
3167  INTEGER :: errorflag
3168  TYPE(t_grid), POINTER :: pgrid
3169  TYPE(t_global), POINTER :: global
3170 
3171 ! ******************************************************************************
3172 ! Start
3173 ! ******************************************************************************
3174 
3175  global => pregion%global
3176 
3177  CALL registerfunction(global,'RFLU_NullifyCell2FaceList',&
3178  'RFLU_ModFaceList.F90')
3179 
3180  pgrid => pregion%grid
3181 
3182 ! ******************************************************************************
3183 ! Nullify memory
3184 ! ******************************************************************************
3185 
3186  nullify(pgrid%tet2f)
3187  nullify(pgrid%hex2f)
3188  nullify(pgrid%pri2f)
3189  nullify(pgrid%pyr2f)
3190 
3191 ! ******************************************************************************
3192 ! End
3193 ! ******************************************************************************
3194 
3195  CALL deregisterfunction(global)
3196 
3197  END SUBROUTINE rflu_nullifycell2facelist
3198 
3199 
3200 
3201 
3202 
3203 
3204 
3205 
3206 
3207 ! ******************************************************************************
3208 !
3209 ! Purpose: Nullify face list.
3210 !
3211 ! Description: None.
3212 !
3213 ! Input:
3214 ! pRegion Pointer to region
3215 !
3216 ! Output: None.
3217 !
3218 ! Notes: None.
3219 !
3220 ! ******************************************************************************
3221 
3222  SUBROUTINE rflu_nullifyfacelist(pRegion)
3223 
3224  IMPLICIT NONE
3225 
3226 ! ******************************************************************************
3227 ! Declarations and definitions
3228 ! ******************************************************************************
3229 
3230 ! ==============================================================================
3231 ! Arguments
3232 ! ==============================================================================
3233 
3234  TYPE(t_region), POINTER :: pregion
3235 
3236 ! ==============================================================================
3237 ! Locals
3238 ! ==============================================================================
3239 
3240  INTEGER :: errorflag,ipatch
3241  TYPE(t_grid), POINTER :: pgrid
3242  TYPE(t_global), POINTER :: global
3243  TYPE(t_patch), POINTER :: ppatch
3244 
3245 ! ******************************************************************************
3246 ! Start
3247 ! ******************************************************************************
3248 
3249  global => pregion%global
3250 
3251  CALL registerfunction(global,'RFLU_NullifyFaceList',&
3252  'RFLU_ModFaceList.F90')
3253 
3254 ! ******************************************************************************
3255 ! Set grid pointer
3256 ! ******************************************************************************
3257 
3258  pgrid => pregion%grid
3259 
3260 ! ******************************************************************************
3261 ! Nullify memory
3262 ! ******************************************************************************
3263 
3264  nullify(pgrid%f2c)
3265  nullify(pgrid%f2v)
3266 
3267 ! TEMPORARY
3268 ! DO iPatch = 1,pGrid%nPatches
3269 ! pPatch => pRegion%patches(iPatch)
3270 !
3271 ! NULLIFY(pPatch%bf2c)
3272 ! NULLIFY(pPatch%bf2v)
3273 ! END DO ! iPatch
3274 ! END TEMPORARY
3275 
3276 ! ******************************************************************************
3277 ! End
3278 ! ******************************************************************************
3279 
3280  CALL deregisterfunction(global)
3281 
3282  END SUBROUTINE rflu_nullifyfacelist
3283 
3284 
3285 
3286 
3287 
3288 
3289 
3290 
3291 ! ******************************************************************************
3292 !
3293 ! Purpose: Reorient faces.
3294 !
3295 ! Description: None.
3296 !
3297 ! Input:
3298 ! pRegion Pointer to region
3299 !
3300 ! Output: None.
3301 !
3302 ! Notes:
3303 ! 1. Need to make sure that normals of interpartition faces are identical for
3304 ! the two domains to which such faces belong. This is because, for
3305 ! example, the Roe dissipative flux depends on ABS(qh-ah), and this will
3306 ! give different results for the two domains if qh has different signs in
3307 ! the two domains, which it will if the normal vector does not point in
3308 ! the same direction. Here make normal vector point in the same direction
3309 ! as on serial grid.
3310 !
3311 ! ******************************************************************************
3312 
3313  SUBROUTINE rflu_reorientfaces(pRegion)
3314 
3316 
3317  IMPLICIT NONE
3318 
3319 ! ******************************************************************************
3320 ! Declarations and definitions
3321 ! ******************************************************************************
3322 
3323 ! ==============================================================================
3324 ! Arguments
3325 ! ==============================================================================
3326 
3327  TYPE(t_region), POINTER :: pregion
3328 
3329 ! ==============================================================================
3330 ! Locals
3331 ! ==============================================================================
3332 
3333  INTEGER :: cntr,c1,c1s,c2,c2s,ifg,v1g,v2g,v3g,v4g
3334  TYPE(t_grid), POINTER :: pgrid
3335  TYPE(t_global), POINTER :: global
3336 
3337 ! ******************************************************************************
3338 ! Start
3339 ! ******************************************************************************
3340 
3341  global => pregion%global
3342 
3343  CALL registerfunction(global,'RFLU_ReorientFaces',&
3344  'RFLU_ModFaceList.F90')
3345 
3346  IF ( global%myProcid == masterproc .AND. &
3347  global%verbLevel >= verbose_high ) THEN
3348  WRITE(stdout,'(A,1X,A)') solver_name, &
3349  'Reorienting faces...'
3350  END IF ! global%verbLevel
3351 
3352 ! ******************************************************************************
3353 ! Set grid pointer and initialize counter
3354 ! ******************************************************************************
3355 
3356  pgrid => pregion%grid
3357 
3358  cntr = 0
3359 
3360 ! ******************************************************************************
3361 ! Loop over faces and reorient depending on indices in serial grid
3362 ! ******************************************************************************
3363 
3364  DO ifg = 1,pgrid%nFaces
3365  c1 = pgrid%f2c(1,ifg)
3366  c2 = pgrid%f2c(2,ifg)
3367 
3368  c1s = pgrid%pc2sc(c1)
3369  c2s = pgrid%pc2sc(c2)
3370 
3371  IF ( c1s > c2s ) THEN
3372  pgrid%f2c(1,ifg) = c2
3373  pgrid%f2c(2,ifg) = c1
3374 
3375  v1g = pgrid%f2v(1,ifg)
3376  v2g = pgrid%f2v(2,ifg)
3377  v3g = pgrid%f2v(3,ifg)
3378  v4g = pgrid%f2v(4,ifg)
3379 
3380  pgrid%f2v(1,ifg) = v3g
3381  pgrid%f2v(2,ifg) = v2g
3382  pgrid%f2v(3,ifg) = v1g
3383  pgrid%f2v(4,ifg) = v4g
3384 
3385  cntr = cntr + 1
3386  END IF ! c1s
3387  END DO ! ifg
3388 
3389 ! ******************************************************************************
3390 ! Write info
3391 ! ******************************************************************************
3392 
3393  IF ( global%myProcid == masterproc .AND. &
3394  global%verbLevel >= verbose_high ) THEN
3395  WRITE(stdout,'(A,3X,A,1X,I6)') solver_name, &
3396  'Number of reoriented faces:',cntr
3397  END IF ! global%verbLevel
3398 
3399 ! ******************************************************************************
3400 ! End
3401 ! ******************************************************************************
3402 
3403  IF ( global%myProcid == masterproc .AND. &
3404  global%verbLevel >= verbose_high ) THEN
3405  WRITE(stdout,'(A,1X,A)') solver_name, &
3406  'Reorienting faces done.'
3407  END IF ! global%verbLevel
3408 
3409  CALL deregisterfunction(global)
3410 
3411  END SUBROUTINE rflu_reorientfaces
3412 
3413 
3414 
3415 
3416 
3417 
3418 ! ******************************************************************************
3419 ! End
3420 ! ******************************************************************************
3421 
3422 END MODULE rflu_modfacelist
3423 
3424 
3425 ! ******************************************************************************
3426 !
3427 ! RCS Revision history:
3428 !
3429 ! $Log: RFLU_ModFaceList.F90,v $
3430 ! Revision 1.43 2008/12/06 08:44:21 mtcampbe
3431 ! Updated license.
3432 !
3433 ! Revision 1.42 2008/11/19 22:17:32 mtcampbe
3434 ! Added Illinois Open Source License/Copyright
3435 !
3436 ! Revision 1.41 2006/08/18 14:01:48 haselbac
3437 ! Added routines for AVFace2Patch list, removed Nullify routines from IF
3438 !
3439 ! Revision 1.40 2006/03/25 21:52:12 haselbac
3440 ! Substantial changes because of sype patches
3441 !
3442 ! Revision 1.39 2006/03/20 13:54:04 haselbac
3443 ! Added output of region index
3444 !
3445 ! Revision 1.38 2005/05/18 22:11:10 fnajjar
3446 ! ACH: Fixed bug in building of av-face list, added setting of nFacesAV
3447 !
3448 ! Revision 1.37 2005/04/29 23:00:11 haselbac
3449 ! Added routines for AVFace2Border list
3450 !
3451 ! Revision 1.36 2005/04/15 15:06:52 haselbac
3452 ! Removed Charm/FEM stuff and changed routine to reorient faces
3453 !
3454 ! Revision 1.35 2004/12/29 21:07:33 haselbac
3455 ! Removed setting of pGrid%nBFaces, now done when reading dims
3456 !
3457 ! Revision 1.34 2004/12/04 03:29:48 haselbac
3458 ! Cosmetics only
3459 !
3460 ! Revision 1.33 2004/11/09 00:28:37 haselbac
3461 ! Cosmetics only
3462 !
3463 ! Revision 1.32 2004/11/03 15:04:53 haselbac
3464 ! Changed correction of estimation of boundary faces, cosmetics
3465 !
3466 ! Revision 1.31 2004/10/19 19:40:37 haselbac
3467 ! Adapted to changes in boundary connectivity, GEN3
3468 !
3469 ! Revision 1.30 2004/10/07 14:55:55 haselbac
3470 ! Bug fix: Prism limit for pyramid loop
3471 !
3472 ! Revision 1.29 2004/10/06 15:35:19 haselbac
3473 ! Bug fixes: VERT_NONE was used to access face and cell conn
3474 !
3475 ! Revision 1.28 2004/09/27 02:43:22 haselbac
3476 ! Bug fix in RFLU_InsertIntoCell2FaceList
3477 !
3478 ! Revision 1.27 2004/09/27 01:39:03 haselbac
3479 ! Added proper headers and routine for getting opposing faces
3480 !
3481 ! Revision 1.26 2004/07/06 15:14:38 haselbac
3482 ! Adapted to changes in libflu and modflu, cosmetics
3483 !
3484 ! Revision 1.25 2004/06/16 20:01:01 haselbac
3485 ! Now store nBFaces directly in grid data type
3486 !
3487 ! Revision 1.24 2004/01/22 16:03:59 haselbac
3488 ! Made contents of modules PRIVATE, only procs PUBLIC, to avoid errors on ALC
3489 ! and titan
3490 !
3491 ! Revision 1.23 2003/12/04 03:28:47 haselbac
3492 ! Cleaned up
3493 !
3494 ! Revision 1.22 2003/11/03 03:49:51 haselbac
3495 ! Removed building and deallocation of bf2bg list
3496 !
3497 ! Revision 1.21 2003/08/20 02:09:58 haselbac
3498 ! Changed verbosity conditions to reduce solver output in GENx runs
3499 !
3500 ! Revision 1.20 2003/08/19 22:48:43 haselbac
3501 ! Added explicit initialization (to avoid Frost problems)
3502 !
3503 ! Revision 1.19 2003/06/04 22:08:30 haselbac
3504 ! Added Nullify routines, some cosmetics
3505 !
3506 ! Revision 1.18 2003/04/16 19:14:24 mtcampbe
3507 ! ACH: Changed implied to explicit loops bcos of Frost problems
3508 !
3509 ! Revision 1.17 2003/04/07 14:24:45 haselbac
3510 ! Added cell-to-face lists
3511 !
3512 ! Revision 1.16 2003/04/02 17:27:22 haselbac
3513 ! Changed limits of modified estimate of no of faces
3514 !
3515 ! Revision 1.15 2003/03/25 19:14:46 haselbac
3516 ! Added new subroutine for reorienting AV faces
3517 !
3518 ! Revision 1.14 2003/03/15 18:07:29 haselbac
3519 ! Added creation, removed face splitting, numerous other fixes
3520 !
3521 ! Revision 1.13 2003/01/28 16:29:53 haselbac
3522 ! Added creation, removed renumbering (bcos of RFLU_InitFlowSolver changes),
3523 ! and more
3524 !
3525 ! Revision 1.12 2002/10/27 19:04:50 haselbac
3526 ! Bug fix, added CHECK_DATASTRUCT output, and use explicit copies (ASCI White)
3527 !
3528 ! Revision 1.11 2002/10/08 15:49:21 haselbac
3529 ! {IO}STAT=global%error replaced by {IO}STAT=errorFlag - SGI problem
3530 !
3531 ! Revision 1.10 2002/09/10 20:26:24 haselbac
3532 ! Corrected bug in RFLU_BuildExtraBFaceList
3533 !
3534 ! Revision 1.9 2002/09/09 15:05:09 haselbac
3535 ! global now under regions, compute nBFaces, added construction of bf2bg access
3536 ! array
3537 !
3538 ! Revision 1.8 2002/07/27 18:10:39 haselbac
3539 ! More conservative estimate of nFacesEst when no boundaries present, got
3540 ! overflow
3541 !
3542 ! Revision 1.7 2002/07/25 15:00:19 haselbac
3543 ! Only write out for MASTERPROC, more output for CHECK_DATASTRUCT
3544 !
3545 ! Revision 1.6 2002/06/27 15:49:47 haselbac
3546 ! Modifications for parallelization, flagging of boundary faces, added
3547 ! CHECK_DATASTRUCT flag
3548 !
3549 ! Revision 1.5 2002/06/17 13:39:45 haselbac
3550 ! Prefixed SOLVER_NAME to all screen output
3551 !
3552 ! Revision 1.4 2002/06/14 20:13:18 haselbac
3553 ! Added destroy flag, routines for parallelization (for Charm)
3554 !
3555 ! Revision 1.3 2002/06/05 18:51:58 haselbac
3556 ! Cosmetic change only, added empty line after destroying face list
3557 !
3558 ! Revision 1.2 2002/05/04 16:39:51 haselbac
3559 ! Cosmetic changes
3560 !
3561 ! Revision 1.1 2002/04/11 18:48:48 haselbac
3562 ! Initial revision
3563 !
3564 ! ******************************************************************************
3565 
3566 
3567 
3568 
3569 
3570 
3571 
3572 
3573 
3574 
3575 
3576 
3577 
3578 
3579 
3580 
3581 
3582 
3583 
3584 
3585 
3586 
3587 
3588 
3589 
subroutine, public rflu_buildcell2facelist(pRegion)
subroutine, public rflu_destroycell2facelist(pRegion)
subroutine rflu_nullifyfacelist(pRegion)
subroutine, public rflu_createhashtable(global, size)
subroutine, public rflu_destroyfacelist(pRegion)
INTEGER function, public rflu_getglobalcellkind(global, pGrid, icg)
subroutine, public rflu_buildavface2borderlist(pRegion)
subroutine, public rflu_destroyhashtable(global)
subroutine, public rflu_destroyavface2patchlist(pRegion)
subroutine, public rflu_createcell2facelist(pRegion)
Vector_n max(const Array_n_const &v1, const Array_n_const &v2)
Definition: Vector_n.h:354
subroutine, public rflu_createavface2borderlist(pRegion)
subroutine registerfunction(global, funName, fileName)
Definition: ModError.F90:449
subroutine, public rflu_hashbuildkey(a, aSize, key)
**********************************************************************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 ic
subroutine, public rflu_unhashbface(global, key, pGrid, fv, nVert, bcType, icg, ifg)
subroutine quicksortinteger(a, n)
RT c() const
Definition: Line_2.h:150
subroutine binarysearchinteger(a, n, v, i, j)
subroutine, public rflu_createfacelist(pRegion)
subroutine rflu_nullifyavface2borderlist(pRegion)
IndexType nfaces() const
Definition: Mesh.H:641
subroutine quicksortintegerinteger(a, b, n)
subroutine, public rflu_buildfacelist(pRegion)
INTEGER function, public rflu_getfacekind(global, c1k, c2k)
subroutine, public rflu_destroyavface2borderlist(pRegion)
subroutine, public rflu_buildavface2patchlist(pRegion)
subroutine, public rflu_getopposingfaces(pRegion, iPatch, ifg, nFacesOpp, faceOppInfo)
subroutine rflu_nullifyavface2patchlist(pRegion)
subroutine errorstop(global, errorCode, errorLine, addMessage)
Definition: ModError.F90:483
subroutine rflu_nullifycell2facelist(pRegion)
subroutine, public rflu_createavface2patchlist(pRegion)
static T_Key key
Definition: vinci_lass.c:76
subroutine deregisterfunction(global)
Definition: ModError.F90:469
subroutine, public rflu_hashface(global, key, pGrid, icg, ifl, fv, nVert, faceType)
subroutine, public rflu_reorientfaces(pRegion)
INTEGER function, public rflu_getglobalcelltype(global, pGrid, icg)
subroutine, public rflu_insertintocell2facelist(pRegion, iPatch, icg, ifg)