Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
RFLU_ModTECPLOTUtils.F90
Go to the documentation of this file.
1 ! *********************************************************************
2 ! * Rocstar Simulation Suite *
3 ! * Copyright@2015, Illinois Rocstar LLC. All rights reserved. *
4 ! * *
5 ! * Illinois Rocstar LLC *
6 ! * Champaign, IL *
7 ! * www.illinoisrocstar.com *
8 ! * sales@illinoisrocstar.com *
9 ! * *
10 ! * License: See LICENSE file in top level of distribution package or *
11 ! * http://opensource.org/licenses/NCSA *
12 ! *********************************************************************
13 ! *********************************************************************
14 ! * THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, *
15 ! * EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES *
16 ! * OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND *
17 ! * NONINFRINGEMENT. IN NO EVENT SHALL THE CONTRIBUTORS OR *
18 ! * COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER *
19 ! * LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, *
20 ! * Arising FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE *
21 ! * USE OR OTHER DEALINGS WITH THE SOFTWARE. *
22 ! *********************************************************************
23 
24 !
25 ! Purpose: Collection of utility routines for writing TECPLOT files.
26 !
27 ! Description: None.
28 !
29 ! Input: None.
30 !
31 ! Output: None
32 !
33 ! Notes: None.
34 !
35 ! ******************************************************************************
36 !
37 ! $Id: RFLU_ModTECPLOTUtils.F90,v 1.14 2008/12/06 08:45:06 mtcampbe Exp $
38 !
39 ! Copyright: (c) 2004 by the University of Illinois
40 !
41 ! ******************************************************************************
42 
44 
45  USE moddatatypes
46  USE modparameters
47  USE moderror
48  USE modglobal, ONLY: t_global
49  USE modgrid, ONLY: t_grid
50  USE modbndpatch, ONLY: t_patch
51  USE moddatastruct, ONLY: t_region
52 
53  IMPLICIT NONE
54 
55 ! ******************************************************************************
56 ! Definitions and declarations
57 ! ******************************************************************************
58 
59  PRIVATE
60 
61 ! ==============================================================================
62 ! Private data
63 ! ==============================================================================
64 
65  CHARACTER(CHRLEN), PARAMETER, PRIVATE :: &
66  RCSIdentString = '$RCSfile: RFLU_ModTECPLOTUtils.F90,v $ $Revision: 1.14 $'
67 
68 ! ==============================================================================
69 ! Public data
70 ! ==============================================================================
71 
72  CHARACTER(1), PARAMETER, PUBLIC :: NULLCHAR = CHAR(0)
73  CHARACTER(2*CHRLEN), PUBLIC :: headerTEC
74  INTEGER, PARAMETER, PUBLIC :: ZONE_TYPE_TRI = 2, &
75  ZONE_TYPE_QUAD = 3, &
76  ZONE_TYPE_TET = 4, &
77  ZONE_TYPE_HEX = 5
78  INTEGER, PARAMETER, PUBLIC :: ZONE_FORM_BLOCK = 1
79  INTEGER, PARAMETER, PUBLIC :: VAR_POS_CELL = 0, &
80  VAR_POS_FACE = 0, &
81  VAR_POS_VERT = 1
82  INTEGER, PARAMETER, PUBLIC :: FILE_TYPE_FIELD = 1, &
83  FILE_TYPE_PATCH = 2, &
84  FILE_TYPE_PATCH_STATS = 3
85  INTEGER, PARAMETER, PUBLIC :: FILE_CNTR_TEC_MAX = 10
86  INTEGER, PUBLIC :: debugFlagTEC,doubleFlagTEC,fileCntrTEC,nVarsTEC, &
87  nVarsTotTEC,nVarsCellTEC,nVarsFaceTEC,nVarsVertTEC
88  INTEGER, PUBLIC :: fileType2CntrTEC(FILE_CNTR_TEC_MAX)
89  INTEGER, DIMENSION(8), PARAMETER, PUBLIC :: tet2vTEC = (/1,2,3,3,4,4,4,4/)
90  INTEGER, DIMENSION(8), PARAMETER, PUBLIC :: hex2vTEC = (/1,2,3,4,5,6,7,8/)
91  INTEGER, DIMENSION(8), PARAMETER, PUBLIC :: pri2vTEC = (/1,2,3,1,4,5,6,4/)
92  INTEGER, DIMENSION(8), PARAMETER, PUBLIC :: pyr2vTEC = (/1,2,3,4,5,5,5,5/)
93  INTEGER, ALLOCATABLE, DIMENSION(:), PUBLIC :: posTEC
94  REAL(RFREAL), ALLOCATABLE, DIMENSION(:,:), PUBLIC :: varCellTEC,varFaceTEC, &
95  varVertTEC
96 
97 ! ==============================================================================
98 ! Public procedures
99 ! ==============================================================================
100 
101  PUBLIC :: rflu_tec_closefile, &
109 
110 
111 ! ******************************************************************************
112 ! Subroutines and functions
113 ! ******************************************************************************
114 
115  CONTAINS
116 
117 
118 
119 
120 
121 
122 ! ******************************************************************************
123 !
124 ! Purpose: Close TECPLOT file
125 !
126 ! Description: None.
127 !
128 ! Input:
129 ! global Pointer to global type
130 !
131 ! Output: None.
132 !
133 ! Notes:
134 ! 1. Isolated into separate routine so that can write several zones to same
135 ! TECPLOT file.
136 !
137 ! ******************************************************************************
138 
139 SUBROUTINE rflu_tec_closefile(global)
140 
141 ! ******************************************************************************
142 ! Declarations and definitions
143 ! ******************************************************************************
144 
145 ! ==============================================================================
146 ! Locals
147 ! ==============================================================================
148 
149  INTEGER :: errorflag
150 
151 ! ==============================================================================
152 ! Arguments
153 ! ==============================================================================
154 
155  TYPE(t_global), POINTER :: global
156 
157 ! ==============================================================================
158 ! Externals: TECPLOT functions
159 ! ==============================================================================
160 
161  INTEGER, EXTERNAL :: tecend100
162 
163 ! ******************************************************************************
164 ! Start
165 ! ******************************************************************************
166 
167  CALL registerfunction(global,'RFLU_TEC_CloseFile', &
168  'RFLU_ModTECPLOTUtils.F90')
169 
170 ! ******************************************************************************
171 ! Close TECPLOT file
172 ! ******************************************************************************
173 
174  errorflag = tecend100()
175  global%error = errorflag
176  IF ( global%error /= err_none ) THEN
177  CALL errorstop(global,err_tecplot_output,__line__)
178  END IF ! global%error
179 
180 ! ******************************************************************************
181 ! End
182 ! ******************************************************************************
183 
184  CALL deregisterfunction(global)
185 
186 END SUBROUTINE rflu_tec_closefile
187 
188 
189 
190 
191 
192 
193 
194 ! ******************************************************************************
195 !
196 ! Purpose: Open TECPLOT file
197 !
198 ! Description: None.
199 !
200 ! Input:
201 ! global Pointer to global data
202 ! title Title of TECPLOT dataset
203 ! fileName Name of TECPLOT file
204 !
205 ! Output: None.
206 !
207 ! Notes:
208 ! 1. Isolated into separate routine so that can write several zones to same
209 ! TECPLOT file.
210 !
211 ! ******************************************************************************
212 
213 SUBROUTINE rflu_tec_openfile(global,title,fileName)
214 
215  IMPLICIT NONE
216 
217 ! ******************************************************************************
218 ! Declarations and definitions
219 ! ******************************************************************************
220 
221 ! ==============================================================================
222 ! Arguments
223 ! ==============================================================================
224 
225  CHARACTER(CHRLEN) :: filename,title
226  TYPE(t_global), POINTER :: global
227 
228 ! ==============================================================================
229 ! Locals
230 ! ==============================================================================
231 
232  CHARACTER(CHRLEN) :: scratchdir
233  INTEGER :: errorflag
234 
235 ! ==============================================================================
236 ! Externals: TECPLOT functions
237 ! ==============================================================================
238 
239  INTEGER, EXTERNAL :: tecini100
240 
241 ! ******************************************************************************
242 ! Start
243 ! ******************************************************************************
244 
245  CALL registerfunction(global,'RFLU_TEC_OpenFile', &
246  'RFLU_ModTECPLOTUtils.F90')
247 
248 ! ==============================================================================
249 ! Set some TECPLOT variables
250 ! ==============================================================================
251 
252  doubleflagtec = 1 ! Flag for single/double precision data
253  debugflagtec = 0 ! Flag for debugging
254 
255  scratchdir = '.'
256 
257 ! ==============================================================================
258 ! Open TECPLOT file
259 ! ==============================================================================
260 
261  errorflag = tecini100(trim(title)//nullchar, &
262  trim(headertec)//nullchar, &
263  trim(filename)//nullchar, &
264  trim(scratchdir)//nullchar, &
265  debugflagtec, &
266  doubleflagtec)
267  global%error = errorflag
268  IF ( global%error /= err_none ) THEN
269  CALL errorstop(global,err_tecplot_output,__line__)
270  END IF ! global%error
271 
272 ! ******************************************************************************
273 ! End
274 ! ******************************************************************************
275 
276  CALL deregisterfunction(global)
277 
278 END SUBROUTINE rflu_tec_openfile
279 
280 
281 
282 
283 
284 
285 
286 ! ******************************************************************************
287 !
288 ! Purpose: Write zone for special cells.
289 !
290 ! Description: None.
291 !
292 ! Input:
293 ! pRegion Pointer to region
294 !
295 ! Output: None.
296 !
297 ! Notes:
298 ! 1. Special cells are written as degenerate hexahedra.
299 !
300 ! ******************************************************************************
301 
303 
304  USE modsortsearch
305 
306  USE rflu_modgrid
307 
308  IMPLICIT NONE
309 
310 ! ******************************************************************************
311 ! Declarations and definitions
312 ! ******************************************************************************
313 
314 ! ==============================================================================
315 ! Arguments
316 ! ==============================================================================
317 
318  TYPE(t_region), POINTER :: pregion
319 
320 ! ==============================================================================
321 ! Locals
322 ! ==============================================================================
323 
324  CHARACTER(CHRLEN) :: zoneprefix,zonetitle
325  INTEGER :: cs2vcard,errorflag,icg,icl,ics,ict,iloc,ivar,ivarcell,ivarvert, &
326  ivcntr,ivl,nvert,zonetype
327  INTEGER :: cs2v(8),cs2vcopy(8),cs2vflag(8),cs2vtec(8)
328  INTEGER, DIMENSION(:), ALLOCATABLE :: nullarray
329  REAL(RFREAL), DIMENSION(:), ALLOCATABLE :: varcelltmp
330  REAL(RFREAL), DIMENSION(:,:), ALLOCATABLE :: varverttmp
331  TYPE(t_global), POINTER :: global
332  TYPE(t_grid), POINTER :: pgrid
333 
334 ! ==============================================================================
335 ! Externals: TECPLOT functions
336 ! ==============================================================================
337 
338  INTEGER, EXTERNAL :: tecdat100,tecnod100,teczne100
339 
340 ! ******************************************************************************
341 ! Start
342 ! ******************************************************************************
343 
344  global => pregion%global
345 
346  CALL registerfunction(global,'RFLU_TEC_WriteZoneCellsSpecial', &
347  'RFLU_ModTECPLOTUtils.F90')
348 
349 ! ******************************************************************************
350 ! Set grid pointer
351 ! ******************************************************************************
352 
353  pgrid => pregion%grid
354 
355 ! ******************************************************************************
356 ! Create null array
357 ! ******************************************************************************
358 
359  ALLOCATE(nullarray(nvarstec),stat=errorflag)
360  global%error = errorflag
361  IF ( global%error /= err_none ) THEN
362  CALL errorstop(global,err_allocate,__line__,'nullArray')
363  END IF ! global%error
364 
365  DO ivar = 1,nvarstec
366  nullarray(ivar) = 0
367  END DO ! iVar
368 
369 ! ******************************************************************************
370 ! Loop over special cells
371 ! ******************************************************************************
372 
373  DO ics = 1,pgrid%nCellsSpecial
374  icg = pgrid%cellsSpecial(ics)
375 
376 ! ==============================================================================
377 ! Determine cell type and local cell index
378 ! ==============================================================================
379 
380  ict = pgrid%cellGlob2Loc(1,icg)
381  icl = pgrid%cellGlob2Loc(2,icg)
382 
383  IF ( global%verbLevel > verbose_none ) THEN
384  WRITE(stdout,'(A,5X,A,1X,I3)') solver_name,'Writing special cell:',ics
385  WRITE(stdout,'(A,7X,A,1X,I9)') solver_name,'Global cell index:',icg
386  WRITE(stdout,'(A,7X,A,1X,I9)') solver_name,'Local cell index: ',icl
387  END IF ! global%verbLevel
388 
389 ! ==============================================================================
390 ! Set local cell connectivity. Note the local cell connectivity for TECPLOT
391 ! contains duplicated vertices...
392 ! ==============================================================================
393 
394  SELECT CASE ( ict )
395  CASE ( cell_type_tet )
396  IF ( global%verbLevel > verbose_none ) THEN
397  WRITE(stdout,'(A,7X,A)') solver_name,'Local cell type: Tetrahedron'
398  END IF ! global%verbLevel
399 
400  nvert = 4
401 
402  DO ivl = 1,8
403  cs2v(ivl) = pgrid%tet2v(tet2vtec(ivl),icl)
404  cs2vtec(ivl) = tet2vtec(ivl)
405  END DO ! ivl
406  CASE ( cell_type_hex )
407  IF ( global%verbLevel > verbose_none ) THEN
408  WRITE(stdout,'(A,7X,A)') solver_name,'Local cell type: Hexahedron'
409  END IF ! global%verbLevel
410 
411  nvert = 8
412 
413  DO ivl = 1,8
414  cs2v(ivl) = pgrid%hex2v(hex2vtec(ivl),icl)
415  cs2vtec(ivl) = hex2vtec(ivl)
416  END DO ! ivl
417  CASE ( cell_type_pri )
418  IF ( global%verbLevel > verbose_none ) THEN
419  WRITE(stdout,'(A,7X,A)') solver_name,'Local cell type: Prism'
420  END IF ! global%verbLevel
421 
422  nvert = 6
423 
424  DO ivl = 1,8
425  cs2v(ivl) = pgrid%pri2v(pri2vtec(ivl),icl)
426  cs2vtec(ivl) = pri2vtec(ivl)
427  END DO ! ivl
428  CASE ( cell_type_pyr )
429  IF ( global%verbLevel > verbose_none ) THEN
430  WRITE(stdout,'(A,7X,A)') solver_name,'Local cell type: Pyramid'
431  END IF ! global%verbLevel
432 
433  nvert = 5
434 
435  DO ivl = 1,8
436  cs2v(ivl) = pgrid%pyr2v(pyr2vtec(ivl),icl)
437  cs2vtec(ivl) = pyr2vtec(ivl)
438  END DO ! ivl
439  CASE default
440  CALL errorstop(global,err_reached_default,__line__)
441  END SELECT ! ict
442 
443  IF ( global%verbLevel > verbose_low ) THEN
444  WRITE(stdout,'(A,7X,A,8(1X,I9))') solver_name,'Vertices:',cs2v(1:8)
445  END IF ! global%verbLevel
446 
447 ! ==============================================================================
448 ! Allocate temporary memory
449 ! ==============================================================================
450 
451  ALLOCATE(varverttmp(nvert,nvarsverttec),stat=errorflag)
452  global%error = errorflag
453  IF ( global%error /= err_none ) THEN
454  CALL errorstop(global,err_allocate,__line__,'varVertTmp')
455  END IF ! global%error
456 
457  IF ( nvarscelltec > 0 ) THEN
458  ALLOCATE(varcelltmp(nvarscelltec),stat=errorflag)
459  global%error = errorflag
460  IF ( global%error /= err_none ) THEN
461  CALL errorstop(global,err_allocate,__line__,'varCellTmp')
462  END IF ! global%error
463  END IF ! nVarsCellTEC
464 
465 ! ==============================================================================
466 ! Build additional data structures: Sort cs2v array and delete duplicated
467 ! vertices. This is required because solution is not written in duplicated
468 ! manner. Also check that vertex list without duplicated vertices has correct
469 ! length (same as nVertPerCell defined above).
470 ! ==============================================================================
471 
472  cs2vcopy(1:8) = cs2v(1:8)
473 
474  CALL quicksortinteger(cs2vcopy,8)
475  CALL simplifysortedintegers(cs2vcopy,8,cs2vcard)
476 
477  IF ( cs2vcard /= nvert ) THEN
478  CALL errorstop(global,err_svert_list_invalid,__line__)
479  END IF ! cs2vCard
480 
481 ! ==============================================================================
482 ! Extract vertex variables
483 ! ==============================================================================
484 
485  IF ( global%verbLevel > verbose_low ) THEN
486  WRITE(stdout,'(A,7X,A)') solver_name,'Cell vertex coordinates:'
487  WRITE(stdout,'(A,9X,A,1X,A,4X,A,2X,6X,A,2(12X,A))') solver_name,'#', &
488  'Local','Global','x-coordinate','y-coordinate','z-coordinate'
489  END IF ! global%verbLevel
490 
491  cs2vflag(1:8) = 0
492  ivcntr = 0
493 
494  DO ivl = 1,8
495  CALL binarysearchinteger(cs2vcopy(1:nvert),nvert,cs2v(ivl),iloc)
496 
497  IF ( cs2vflag(iloc) == 0 ) THEN
498  ivcntr = ivcntr + 1
499 
500  cs2vflag(iloc) = ivcntr
501 
502  ivarvert = 0
503 
504  DO ivar = 1,nvarstec
505  IF ( postec(ivar) == var_pos_vert ) THEN
506  ivarvert = ivarvert + 1
507  varverttmp(ivcntr,ivarvert) = pregion%varVertTEC(cs2v(ivl),ivarvert)
508  END IF ! posTEC
509  END DO ! iVar
510 
511  IF ( global%verbLevel > verbose_low ) THEN
512  WRITE(stdout,'(A,9X,I1,3X,I1,3X,I9,1X,3(1X,E23.16))') solver_name, &
513  ivcntr,ivl,cs2v(ivl),pregion%varVertTEC(cs2v(ivl),1:3)
514  END IF ! global%verbLevel
515  END IF ! cs2vFlag
516  END DO ! ivl
517 
518 ! ==============================================================================
519 ! Extract cell variables
520 ! ==============================================================================
521 
522  IF ( nvarscelltec > 0 ) THEN
523  ivarcell = 0
524 
525  DO ivar = 1,nvarstec
526  IF ( postec(ivar) == var_pos_cell ) THEN
527  ivarcell = ivarcell + 1
528  varcelltmp(ivarcell) = pregion%varCellTEC(icg,ivarcell)
529  END IF ! posTEC
530  END DO ! iVar
531  END IF ! nVarsCellTEC
532 
533 ! ==============================================================================
534 ! Write to TECPLOT file
535 ! ==============================================================================
536 
537 ! ------------------------------------------------------------------------------
538 ! Zone information
539 ! ------------------------------------------------------------------------------
540 
541  WRITE(zonetitle,'(A,I8.8,A,I5.5)') 'C_',icg,'_',pregion%iRegionGlobal
542  errorflag = teczne100(trim(zonetitle)//nullchar,zone_type_hex,nvert,1, &
543  0,0,0,0,zone_form_block,0,0,postec,nullarray,0)
544  global%error = errorflag
545  IF ( global%error /= err_none ) THEN
546  CALL errorstop(global,err_tecplot_output,__line__)
547  END IF ! global%error
548 
549 ! ------------------------------------------------------------------------------
550 ! Variables
551 ! ------------------------------------------------------------------------------
552 
553  ivarcell = 0
554  ivarvert = 0
555 
556  DO ivar = 1,nvarstec
557  IF ( postec(ivar) == var_pos_cell ) THEN
558  ivarcell = ivarcell + 1
559  errorflag = tecdat100(1,varcelltmp(ivarcell),doubleflagtec)
560  ELSE IF ( postec(ivar) == var_pos_vert ) THEN
561  ivarvert = ivarvert + 1
562  errorflag = tecdat100(nvert,varverttmp(1,ivarvert),doubleflagtec)
563  ELSE ! defensive coding
564  CALL errorstop(global,err_reached_default,__line__)
565  END IF ! posTEC
566 
567  global%error = errorflag
568  IF ( global%error /= err_none ) THEN
569  CALL errorstop(global,err_tecplot_output,__line__)
570  END IF ! global%error
571  END DO ! iVar
572 
573 ! ------------------------------------------------------------------------------
574 ! Connectivity
575 ! ------------------------------------------------------------------------------
576 
577  errorflag = tecnod100(cs2vtec)
578  global%error = errorflag
579  IF ( global%error /= err_none ) THEN
580  CALL errorstop(global,err_tecplot_output,__line__)
581  END IF ! global%error
582 
583 ! ==============================================================================
584 ! Deallocate temporary memory
585 ! ==============================================================================
586 
587  DEALLOCATE(varverttmp,stat=errorflag)
588  global%error = errorflag
589  IF ( global%error /= err_none ) THEN
590  CALL errorstop(global,err_deallocate,__line__,'varVertTmp')
591  END IF ! global%error
592 
593  IF ( nvarscelltec > 0 ) THEN
594  DEALLOCATE(varcelltmp,stat=errorflag)
595  global%error = errorflag
596  IF ( global%error /= err_none ) THEN
597  CALL errorstop(global,err_deallocate,__line__,'varCellTmp')
598  END IF ! global%error
599  END IF ! nVarsCellTEC
600  END DO ! ics
601 
602 ! ******************************************************************************
603 ! Destroy null array
604 ! ******************************************************************************
605 
606  DEALLOCATE(nullarray,stat=errorflag)
607  global%error = errorflag
608  IF ( global%error /= err_none ) THEN
609  CALL errorstop(global,err_deallocate,__line__,'nullArray')
610  END IF ! global%error
611 
612 ! ******************************************************************************
613 ! End
614 ! ******************************************************************************
615 
616  CALL deregisterfunction(global)
617 
618 END SUBROUTINE rflu_tec_writezonecellsspecial
619 
620 
621 
622 
623 
624 
625 
626 
627 
628 ! ******************************************************************************
629 !
630 ! Purpose: Write zone for special faces.
631 !
632 ! Description: None.
633 !
634 ! Input:
635 ! pRegion Pointer to region
636 !
637 ! Output: None.
638 !
639 ! Notes:
640 ! 1. Special faces are written as degenerate quadrilaterals.
641 ! 2. Cannot write special faces when have cell-centered data because do not
642 ! have any data for interior faces.
643 !
644 ! ******************************************************************************
645 
647 
648  USE modsortsearch
649 
650  USE rflu_modgrid
651 
652  IMPLICIT NONE
653 
654 ! ******************************************************************************
655 ! Declarations and definitions
656 ! ******************************************************************************
657 
658 ! ==============================================================================
659 ! Arguments
660 ! ==============================================================================
661 
662  TYPE(t_region), POINTER :: pregion
663 
664 ! ==============================================================================
665 ! Locals
666 ! ==============================================================================
667 
668  CHARACTER(CHRLEN) :: zoneprefix,zonetitle
669  INTEGER :: fs2vcard,errorflag,ifl,ifs,iloc,ipatch,ivar,ivarcell,ivarvert, &
670  ivcntr,ivl,nvert,zonetype
671  INTEGER :: fs2v(4),fs2vcopy(4),fs2vflag(4),fs2vtec(4)
672  INTEGER, DIMENSION(:), ALLOCATABLE :: nullarray
673  REAL(RFREAL), DIMENSION(:,:), ALLOCATABLE :: varverttmp
674  TYPE(t_global), POINTER :: global
675  TYPE(t_grid), POINTER :: pgrid
676  TYPE(t_patch), POINTER :: ppatch
677 
678 ! ==============================================================================
679 ! Externals: TECPLOT functions
680 ! ==============================================================================
681 
682  INTEGER, EXTERNAL :: tecdat100,tecnod100,teczne100
683 
684 ! ******************************************************************************
685 ! Start
686 ! ******************************************************************************
687 
688  global => pregion%global
689 
690  CALL registerfunction(global,'RFLU_TEC_WriteZoneFacesSpecial', &
691  'RFLU_ModTECPLOTUtils.F90')
692 
693 ! ******************************************************************************
694 ! Set grid pointer
695 ! ******************************************************************************
696 
697  pgrid => pregion%grid
698 
699 ! ******************************************************************************
700 ! Create null array
701 ! ******************************************************************************
702 
703  ALLOCATE(nullarray(nvarstec),stat=errorflag)
704  global%error = errorflag
705  IF ( global%error /= err_none ) THEN
706  CALL errorstop(global,err_allocate,__line__,'nullArray')
707  END IF ! global%error
708 
709  DO ivar = 1,nvarstec
710  nullarray(ivar) = 0
711  END DO ! iVar
712 
713 ! ******************************************************************************
714 ! Loop over special cells
715 ! ******************************************************************************
716 
717  DO ifs = 1,pgrid%nFacesSpecial
718  ipatch = pgrid%facesSpecial(1,ifs)
719  ifl = pgrid%facesSpecial(2,ifs)
720 
721  IF ( global%verbLevel > verbose_none ) THEN
722  WRITE(stdout,'(A,5X,A,1X,I3)') solver_name,'Writing special face:',ifs
723  WRITE(stdout,'(A,7X,A,1X,I9)') solver_name,'Patch index:',ipatch
724  WRITE(stdout,'(A,7X,A,1X,I9)') solver_name,'Face index: ',ifl
725  END IF ! global%verbLevel
726 
727 ! ==============================================================================
728 ! Set local face connectivity. Note the local face connectivity for TECPLOT
729 ! contains duplicated vertices for triangles
730 ! ==============================================================================
731 
732  IF ( ipatch == 0 ) THEN
733  DO ivl = 1,4
734  fs2v(ivl) = pgrid%f2v(ivl,ifl)
735  fs2vtec(ivl) = ivl
736  END DO ! ivl
737  ELSE IF ( ipatch > 0 .AND. ipatch <= pgrid%nPatches ) THEN
738  ppatch => pregion%patches(ipatch)
739 
740  DO ivl = 1,4
741  IF ( ppatch%bf2v(ivl,ifl) /= vert_none ) THEN
742  fs2v(ivl) = ppatch%bv(ppatch%bf2v(ivl,ifl))
743  ELSE
744  fs2v(ivl) = vert_none
745  END IF ! pPatch%bf2v
746 
747  fs2vtec(ivl) = ivl
748  END DO ! ivl
749  ELSE
750  CALL errorstop(global,err_reached_default,__line__)
751  END IF ! iPatch
752 
753  IF ( fs2v(4) == vert_none ) THEN ! Triangle
754  nvert = 3
755 
756  fs2v(4) = fs2v(1)
757  fs2vtec(4) = 1
758  ELSE ! Quadrilateral
759  nvert = 4
760  END IF ! fs2v
761 
762  IF ( global%verbLevel > verbose_low ) THEN
763  WRITE(stdout,'(A,7X,A,4(1X,I9))') solver_name,'Vertices:',fs2v(1:4)
764  END IF ! global%verbLevel
765 
766 ! ==============================================================================
767 ! Allocate temporary memory
768 ! ==============================================================================
769 
770  ALLOCATE(varverttmp(nvert,nvarsverttec),stat=errorflag)
771  global%error = errorflag
772  IF ( global%error /= err_none ) THEN
773  CALL errorstop(global,err_allocate,__line__,'varVertTmp')
774  END IF ! global%error
775 
776 ! ==============================================================================
777 ! Build additional data structures: Sort cs2v array and delete duplicated
778 ! vertices. This is required because solution is not written in duplicated
779 ! manner. Also check that vertex list without duplicated vertices has correct
780 ! length (same as nVertPerCell defined above).
781 ! ==============================================================================
782 
783  fs2vcopy(1:4) = fs2v(1:4)
784 
785  CALL quicksortinteger(fs2vcopy,4)
786  CALL simplifysortedintegers(fs2vcopy,4,fs2vcard)
787 
788  IF ( fs2vcard /= nvert ) THEN
789  CALL errorstop(global,err_svert_list_invalid,__line__)
790  END IF ! fs2vCard
791 
792 ! ==============================================================================
793 ! Extract vertex variables
794 ! ==============================================================================
795 
796  IF ( global%verbLevel > verbose_low ) THEN
797  WRITE(stdout,'(A,7X,A)') solver_name,'Cell vertex coordinates:'
798  WRITE(stdout,'(A,9X,A,1X,A,4X,A,2X,6X,A,2(12X,A))') solver_name,'#', &
799  'Local','Global','x-coordinate','y-coordinate','z-coordinate'
800  END IF ! global%verbLevel
801 
802  fs2vflag(1:4) = 0
803  ivcntr = 0
804 
805  DO ivl = 1,4
806  CALL binarysearchinteger(fs2vcopy(1:nvert),nvert,fs2v(ivl),iloc)
807 
808  IF ( fs2vflag(iloc) == 0 ) THEN
809  ivcntr = ivcntr + 1
810 
811  fs2vflag(iloc) = ivcntr
812 
813  ivarvert = 0
814 
815  DO ivar = 1,nvarstec
816  IF ( postec(ivar) == var_pos_vert ) THEN
817  ivarvert = ivarvert + 1
818  varverttmp(ivcntr,ivarvert) = pregion%varVertTEC(fs2v(ivl),ivarvert)
819  END IF ! posTEC
820  END DO ! iVar
821 
822  IF ( global%verbLevel > verbose_low ) THEN
823  WRITE(stdout,'(A,9X,I1,3X,I1,3X,I9,1X,3(1X,E23.16))') solver_name, &
824  ivcntr,ivl,fs2v(ivl),pregion%varVertTEC(fs2v(ivl),1:3)
825  END IF ! global%verbLevel
826  END IF ! cs2vFlag
827  END DO ! ivl
828 
829 ! ==============================================================================
830 ! Write to TECPLOT file
831 ! ==============================================================================
832 
833 ! ------------------------------------------------------------------------------
834 ! Zone information
835 ! ------------------------------------------------------------------------------
836 
837  WRITE(zonetitle,'(A,I8.8,A,I3.3,A,I5.5)') 'F_',ifl,'_PAT_',ipatch,'_', &
838  pregion%iRegionGlobal
839  errorflag = teczne100(trim(zonetitle)//nullchar,zone_type_quad,nvert,1, &
840  0,0,0,0,zone_form_block,0,0,postec,nullarray,0)
841  global%error = errorflag
842  IF ( global%error /= err_none ) THEN
843  CALL errorstop(global,err_tecplot_output,__line__)
844  END IF ! global%error
845 
846 ! ------------------------------------------------------------------------------
847 ! Variables
848 ! ------------------------------------------------------------------------------
849 
850  ivarvert = 0
851 
852  DO ivar = 1,nvarstec
853  IF ( postec(ivar) == var_pos_vert ) THEN
854  ivarvert = ivarvert + 1
855  errorflag = tecdat100(nvert,varverttmp(1,ivarvert),doubleflagtec)
856  ELSE ! defensive coding
857  CALL errorstop(global,err_reached_default,__line__)
858  END IF ! posTEC
859 
860  global%error = errorflag
861  IF ( global%error /= err_none ) THEN
862  CALL errorstop(global,err_tecplot_output,__line__)
863  END IF ! global%error
864  END DO ! iVar
865 
866 ! ------------------------------------------------------------------------------
867 ! Connectivity
868 ! ------------------------------------------------------------------------------
869 
870  errorflag = tecnod100(fs2vtec)
871  global%error = errorflag
872  IF ( global%error /= err_none ) THEN
873  CALL errorstop(global,err_tecplot_output,__line__)
874  END IF ! global%error
875 
876 ! ==============================================================================
877 ! Deallocate temporary memory
878 ! ==============================================================================
879 
880  DEALLOCATE(varverttmp,stat=errorflag)
881  global%error = errorflag
882  IF ( global%error /= err_none ) THEN
883  CALL errorstop(global,err_deallocate,__line__,'varVertTmp')
884  END IF ! global%error
885  END DO ! ifs
886 
887 ! ******************************************************************************
888 ! Destroy null array
889 ! ******************************************************************************
890 
891  DEALLOCATE(nullarray,stat=errorflag)
892  global%error = errorflag
893  IF ( global%error /= err_none ) THEN
894  CALL errorstop(global,err_deallocate,__line__,'nullArray')
895  END IF ! global%error
896 
897 ! ******************************************************************************
898 ! End
899 ! ******************************************************************************
900 
901  CALL deregisterfunction(global)
902 
903 END SUBROUTINE rflu_tec_writezonefacesspecial
904 
905 
906 
907 
908 
909 
910 
911 
912 
913 ! ******************************************************************************
914 !
915 ! Purpose: Write zone for surface data.
916 !
917 ! Description: None.
918 !
919 ! Input:
920 ! pRegion Pointer to region
921 !
922 ! Output: None.
923 !
924 ! Notes: None.
925 !
926 ! ******************************************************************************
927 
928 SUBROUTINE rflu_tec_writezoneinterf(pRegion)
929 
931 
932  IMPLICIT NONE
933 
934 ! ******************************************************************************
935 ! Declarations and definitions
936 ! ******************************************************************************
937 
938 ! ==============================================================================
939 ! Arguments
940 ! ==============================================================================
941 
942  TYPE(t_region), POINTER :: pregion
943 
944 ! ==============================================================================
945 ! Locals
946 ! ==============================================================================
947 
948  CHARACTER(CHRLEN) :: zoneprefix,zonetitle
949  INTEGER :: c1,c1k,c2,c2k,errorflag,ifg,ifl,iquad,itri,ivar, &
950  ivarface,ivarvert,ivl,nfaces,nfacesest,nquads,ntris,nvert, &
951  nvertest,v1,v2,v3,v4,zonetype
952  INTEGER, DIMENSION(:), ALLOCATABLE :: bvflag,nullarray
953  INTEGER, DIMENSION(:,:), ALLOCATABLE :: ibface,f2v,f2vtmp
954  REAL(RFREAL), DIMENSION(:,:), ALLOCATABLE :: varfacetmp,varverttmp
955  TYPE(t_global), POINTER :: global
956  TYPE(t_grid), POINTER :: pgrid
957 
958 ! ==============================================================================
959 ! Externals: TECPLOT functions
960 ! ==============================================================================
961 
962  INTEGER, EXTERNAL :: tecdat100,tecnod100,teczne100
963 
964 ! ******************************************************************************
965 ! Start
966 ! ******************************************************************************
967 
968  global => pregion%global
969 
970  CALL registerfunction(global,'RFLU_TEC_WriteZoneInterf', &
971  'RFLU_ModTECPLOTUtils.F90')
972 
973 ! ******************************************************************************
974 ! Set grid pointer
975 ! ******************************************************************************
976 
977  pgrid => pregion%grid
978 
979 ! ******************************************************************************
980 ! Create null array
981 ! ******************************************************************************
982 
983  ALLOCATE(nullarray(nvarstec),stat=errorflag)
984  global%error = errorflag
985  IF ( global%error /= err_none ) THEN
986  CALL errorstop(global,err_allocate,__line__,'nullArray')
987  END IF ! global%error
988 
989  DO ivar = 1,nvarstec
990  nullarray(ivar) = 0
991  END DO ! iVar
992 
993 ! ******************************************************************************
994 ! Allocate temporary memory. NOTE include pretty large safety factor, because
995 ! the difference between the total and actual number of vertices can lead to
996 ! underestimation of interface vertices if have small number of virtual cells,
997 ! (as for first-order scheme with no grid motion), and hence small number of
998 ! virtual vertices.
999 ! ******************************************************************************
1000 
1001  nvertest = int(10.0_rfreal*(pgrid%nVertTot - pgrid%nVert))
1002 
1003  ALLOCATE(bvflag(pgrid%nVertTot),stat=errorflag)
1004  global%error = errorflag
1005  IF ( global%error /= err_none ) THEN
1006  CALL errorstop(global,err_allocate,__line__,'bvFlag')
1007  END IF ! global%error
1008 
1009  ALLOCATE(varverttmp(nvertest,nvarsverttec),stat=errorflag)
1010  global%error = errorflag
1011  IF ( global%error /= err_none ) THEN
1012  CALL errorstop(global,err_allocate,__line__,'varVertTmp')
1013  END IF ! global%error
1014 
1015 ! ******************************************************************************
1016 ! Allocate memory for interpartition face list
1017 ! ******************************************************************************
1018 
1019  IF ( pgrid%nFacesTot /= pgrid%nFaces ) THEN
1020  nfacesest = pgrid%nFacesTot - pgrid%nFaces
1021  ELSE ! May have just one face
1022  nfacesest = 1
1023  END IF ! pGrid%nFacesTot
1024 
1025  ALLOCATE(f2v(4,nfacesest),stat=errorflag)
1026  global%error = errorflag
1027  IF ( global%error /= err_none ) THEN
1028  CALL errorstop(global,err_allocate,__line__,'f2v')
1029  END IF ! global%error
1030 
1031 ! ******************************************************************************
1032 ! Build temporary interpartition face list
1033 ! ******************************************************************************
1034 
1035  nfaces = 0
1036  ntris = 0
1037  nquads = 0
1038 
1039  DO ifg = 1,pgrid%nFaces ! NOTE upper limit
1040  c1 = pgrid%f2c(1,ifg)
1041  c2 = pgrid%f2c(2,ifg)
1042 
1043  c1k = rflu_getglobalcellkind(global,pgrid,c1)
1044  c2k = rflu_getglobalcellkind(global,pgrid,c2)
1045 
1046  IF ( rflu_getfacekind(global,c1k,c2k) == face_kind_av ) THEN
1047  nfaces = nfaces + 1
1048 
1049  IF ( nfaces > nfacesest ) THEN
1050  CALL errorstop(global,err_nfaces_estimate,__line__)
1051  END IF ! nFaces
1052 
1053  f2v(1:4,nfaces) = pgrid%f2v(1:4,ifg)
1054 
1055  IF ( f2v(4,nfaces) == vert_none ) THEN
1056  ntris = ntris + 1
1057  ELSE
1058  nquads = nquads + 1
1059  END IF ! f2v
1060  END IF ! RFLU_GetFaceKind
1061  END DO ! ifg
1062 
1063 ! ******************************************************************************
1064 ! Build final interpartition face list (note sorting)
1065 ! ******************************************************************************
1066 
1067  ALLOCATE(f2vtmp(4,nfaces),stat=errorflag)
1068  global%error = errorflag
1069  IF ( global%error /= err_none ) THEN
1070  CALL errorstop(global,err_allocate,__line__,'f2vTmp')
1071  END IF ! global%error
1072 
1073  f2vtmp(1:4,1:nfaces) = f2v(1:4,1:nfaces)
1074 
1075  DEALLOCATE(f2v,stat=errorflag)
1076  global%error = errorflag
1077  IF ( global%error /= err_none ) THEN
1078  CALL errorstop(global,err_deallocate,__line__,'f2v')
1079  END IF ! global%error
1080 
1081  ALLOCATE(f2v(4,nfaces),stat=errorflag)
1082  global%error = errorflag
1083  IF ( global%error /= err_none ) THEN
1084  CALL errorstop(global,err_allocate,__line__,'f2v')
1085  END IF ! global%error
1086 
1087  itri = 0
1088  iquad = ntris
1089 
1090  DO ifg = 1,nfaces
1091  IF ( f2vtmp(4,ifg) == vert_none ) THEN
1092  itri = itri + 1
1093  f2v(1:4,itri) = f2vtmp(1:4,ifg)
1094  ELSE
1095  iquad = iquad + 1
1096  f2v(1:4,iquad) = f2vtmp(1:4,ifg)
1097  END IF ! f2vTmp
1098  END DO ! ifg
1099 
1100  DEALLOCATE(f2vtmp,stat=errorflag)
1101  global%error = errorflag
1102  IF ( global%error /= err_none ) THEN
1103  CALL errorstop(global,err_deallocate,__line__,'f2vTmp')
1104  END IF ! global%error
1105 
1106 ! ******************************************************************************
1107 ! Allocate memory for face variables
1108 ! ******************************************************************************
1109 
1110  IF ( nvarsfacetec > 0 ) THEN
1111  ALLOCATE(varfacetmp(nfaces,nvarsfacetec),stat=errorflag)
1112  global%error = errorflag
1113  IF ( global%error /= err_none ) THEN
1114  CALL errorstop(global,err_allocate,__line__,'varFaceTmp')
1115  END IF ! global%error
1116  END IF ! nVarsFaceTEC
1117 
1118 ! ******************************************************************************
1119 ! Triangular interpartition faces
1120 ! ******************************************************************************
1121 
1122  IF ( ntris > 0 ) THEN
1123  ALLOCATE(f2vtmp(3,ntris),stat=errorflag)
1124  global%error = errorflag
1125  IF ( global%error /= err_none ) THEN
1126  CALL errorstop(global,err_allocate,__line__,'f2vTmp')
1127  END IF ! global%error
1128 
1129  DO ivl = 1,SIZE(bvflag,1)
1130  bvflag(ivl) = 0
1131  END DO ! ivl
1132 
1133  nvert = 0
1134 
1135  DO ifl = 1,ntris
1136  v1 = f2v(1,ifl)
1137  v2 = f2v(2,ifl)
1138  v3 = f2v(3,ifl)
1139 
1140  IF ( bvflag(v1) == 0 ) THEN
1141  nvert = nvert + 1
1142  bvflag(v1) = nvert
1143 
1144  ivarvert = 0
1145 
1146  DO ivar = 1,nvarstec
1147  IF ( postec(ivar) == var_pos_vert ) THEN
1148  ivarvert = ivarvert + 1
1149 
1150  IF ( nvert > nvertest ) THEN
1151  CALL errorstop(global,err_nbvert_estimate,__line__)
1152  END IF ! nVert
1153 
1154  varverttmp(nvert,ivarvert) = pregion%varVertTEC(v1,ivarvert)
1155  END IF ! posTEC
1156  END DO ! iVar
1157  END IF ! bvFlag
1158 
1159  IF ( bvflag(v2) == 0 ) THEN
1160  nvert = nvert + 1
1161  bvflag(v2) = nvert
1162 
1163  ivarvert = 0
1164 
1165  DO ivar = 1,nvarstec
1166  IF ( postec(ivar) == var_pos_vert ) THEN
1167  ivarvert = ivarvert + 1
1168 
1169  IF ( nvert > nvertest ) THEN
1170  CALL errorstop(global,err_nbvert_estimate,__line__)
1171  END IF ! nVert
1172 
1173  varverttmp(nvert,ivarvert) = pregion%varVertTEC(v2,ivarvert)
1174  END IF ! posTEC
1175  END DO ! iVar
1176  END IF ! bvFlag
1177 
1178  IF ( bvflag(v3) == 0 ) THEN
1179  nvert = nvert + 1
1180  bvflag(v3) = nvert
1181 
1182  ivarvert = 0
1183 
1184  DO ivar = 1,nvarstec
1185  IF ( postec(ivar) == var_pos_vert ) THEN
1186  ivarvert = ivarvert + 1
1187 
1188  IF ( nvert > nvertest ) THEN
1189  CALL errorstop(global,err_nbvert_estimate,__line__)
1190  END IF ! nVert
1191 
1192  varverttmp(nvert,ivarvert) = pregion%varVertTEC(v3,ivarvert)
1193  END IF ! posTEC
1194  END DO ! iVar
1195  END IF ! bvFlag
1196 
1197  f2vtmp(1,ifl) = bvflag(v1)
1198  f2vtmp(2,ifl) = bvflag(v2)
1199  f2vtmp(3,ifl) = bvflag(v3)
1200  END DO ! ifl
1201 
1202 ! ******************************************************************************
1203 ! Write data
1204 ! ******************************************************************************
1205 
1206 ! ==============================================================================
1207 ! Zone information
1208 ! ==============================================================================
1209 
1210  WRITE(zonetitle,'(A,I5.5)') 'INT_TRI_',pregion%iRegionGlobal
1211  errorflag = teczne100(trim(zonetitle)//nullchar,zone_type_tri,nvert, &
1212  ntris,0,0,0,0,zone_form_block,0,0,postec,nullarray,0)
1213  global%error = errorflag
1214  IF ( global%error /= err_none ) THEN
1215  CALL errorstop(global,err_tecplot_output,__line__)
1216  END IF ! global%error
1217 
1218 ! ==============================================================================
1219 ! Variables
1220 ! ==============================================================================
1221 
1222  ivarface = 0
1223  ivarvert = 0
1224 
1225  DO ivar = 1,nvarstec
1226  IF ( postec(ivar) == var_pos_vert ) THEN
1227  ivarvert = ivarvert + 1
1228  errorflag = tecdat100(nvert,varverttmp(1,ivarvert),doubleflagtec)
1229  ELSE IF ( postec(ivar) == var_pos_face ) THEN
1230  ivarface = ivarface + 1
1231  errorflag = tecdat100(ntris,varfacetmp(1,ivarface),doubleflagtec)
1232  ELSE
1233  CALL errorstop(global,err_reached_default,__line__)
1234  END IF ! posTEC
1235 
1236  global%error = errorflag
1237  IF ( global%error /= err_none ) THEN
1238  CALL errorstop(global,err_tecplot_output,__line__)
1239  END IF ! global%error
1240  END DO ! iVar
1241 
1242 ! ==============================================================================
1243 ! Connectivity
1244 ! ==============================================================================
1245 
1246  errorflag = tecnod100(f2vtmp)
1247  global%error = errorflag
1248  IF ( global%error /= err_none ) THEN
1249  CALL errorstop(global,err_tecplot_output,__line__)
1250  END IF ! global%error
1251 
1252  DEALLOCATE(f2vtmp,stat=errorflag)
1253  global%error = errorflag
1254  IF ( global%error /= err_none ) THEN
1255  CALL errorstop(global,err_deallocate,__line__,'f2vTmp')
1256  END IF ! global%error
1257  END IF ! nTris
1258 
1259 ! ******************************************************************************
1260 ! Quadrilateral interpartition faces
1261 ! ******************************************************************************
1262 
1263  IF ( nquads > 0 ) THEN
1264  ALLOCATE(f2vtmp(4,nquads),stat=errorflag)
1265  global%error = errorflag
1266  IF ( global%error /= err_none ) THEN
1267  CALL errorstop(global,err_allocate,__line__,'f2vTmp')
1268  END IF ! global%error
1269 
1270  DO ivl = 1,SIZE(bvflag,1)
1271  bvflag(ivl) = 0
1272  END DO ! ivl
1273 
1274  nvert = 0
1275 
1276  DO ifl = 1,nquads
1277  v1 = f2v(1,ntris+ifl)
1278  v2 = f2v(2,ntris+ifl)
1279  v3 = f2v(3,ntris+ifl)
1280  v4 = f2v(4,ntris+ifl)
1281 
1282  IF ( bvflag(v1) == 0 ) THEN
1283  nvert = nvert + 1
1284  bvflag(v1) = nvert
1285 
1286  ivarvert = 0
1287 
1288  DO ivar = 1,nvarstec
1289  IF ( postec(ivar) == var_pos_vert ) THEN
1290  ivarvert = ivarvert + 1
1291 
1292  IF ( nvert > nvertest ) THEN
1293  CALL errorstop(global,err_nbvert_estimate,__line__)
1294  END IF ! nVert
1295 
1296  varverttmp(nvert,ivarvert) = pregion%varVertTEC(v1,ivarvert)
1297  END IF ! posTEC
1298  END DO ! iVar
1299  END IF ! bvFlag
1300 
1301  IF ( bvflag(v2) == 0 ) THEN
1302  nvert = nvert + 1
1303  bvflag(v2) = nvert
1304 
1305  ivarvert = 0
1306 
1307  DO ivar = 1,nvarstec
1308  IF ( postec(ivar) == var_pos_vert ) THEN
1309  ivarvert = ivarvert + 1
1310 
1311  IF ( nvert > nvertest ) THEN
1312  CALL errorstop(global,err_nbvert_estimate,__line__)
1313  END IF ! nVert
1314 
1315  varverttmp(nvert,ivarvert) = pregion%varVertTEC(v2,ivarvert)
1316  END IF ! posTEC
1317  END DO ! iVar
1318  END IF ! bvFlag
1319 
1320  IF ( bvflag(v3) == 0 ) THEN
1321  nvert = nvert + 1
1322  bvflag(v3) = nvert
1323 
1324  ivarvert = 0
1325 
1326  DO ivar = 1,nvarstec
1327  IF ( postec(ivar) == var_pos_vert ) THEN
1328  ivarvert = ivarvert + 1
1329  varverttmp(nvert,ivarvert) = pregion%varVertTEC(v3,ivarvert)
1330  END IF ! posTEC
1331  END DO ! iVar
1332  END IF ! bvFlag
1333 
1334  IF ( bvflag(v4) == 0 ) THEN
1335  nvert = nvert + 1
1336  bvflag(v4) = nvert
1337 
1338  ivarvert = 0
1339 
1340  DO ivar = 1,nvarstec
1341  IF ( postec(ivar) == var_pos_vert ) THEN
1342  ivarvert = ivarvert + 1
1343 
1344  IF ( nvert > nvertest ) THEN
1345  CALL errorstop(global,err_nbvert_estimate,__line__)
1346  END IF ! nVert
1347 
1348  varverttmp(nvert,ivarvert) = pregion%varVertTEC(v4,ivarvert)
1349  END IF ! posTEC
1350  END DO ! iVar
1351  END IF ! bvFlag
1352 
1353  f2vtmp(1,ifl) = bvflag(v1)
1354  f2vtmp(2,ifl) = bvflag(v2)
1355  f2vtmp(3,ifl) = bvflag(v3)
1356  f2vtmp(4,ifl) = bvflag(v4)
1357  END DO ! ifl
1358 
1359 ! ******************************************************************************
1360 ! Write data
1361 ! ******************************************************************************
1362 
1363 ! ==============================================================================
1364 ! Zone information
1365 ! ==============================================================================
1366 
1367  WRITE(zonetitle,'(A,I5.5)') 'INT_QUAD_',pregion%iRegionGlobal
1368  errorflag = teczne100(trim(zonetitle)//nullchar,zone_type_quad,nvert, &
1369  nquads,0,0,0,0,zone_form_block,0,0,postec,nullarray,0)
1370  global%error = errorflag
1371  IF ( global%error /= err_none ) THEN
1372  CALL errorstop(global,err_tecplot_output,__line__)
1373  END IF ! global%error
1374 
1375 ! ==============================================================================
1376 ! Variables
1377 ! ==============================================================================
1378 
1379  ivarface = 0
1380  ivarvert = 0
1381 
1382  DO ivar = 1,nvarstec
1383  IF ( postec(ivar) == var_pos_vert ) THEN
1384  ivarvert = ivarvert + 1
1385  errorflag = tecdat100(nvert,varverttmp(1,ivarvert),doubleflagtec)
1386  ELSE IF ( postec(ivar) == var_pos_face ) THEN
1387  ivarface = ivarface + 1
1388  errorflag = tecdat100(nquads,varfacetmp(1,ivarface),doubleflagtec)
1389  ELSE
1390  CALL errorstop(global,err_reached_default,__line__)
1391  END IF ! posTEC
1392 
1393  global%error = errorflag
1394  IF ( global%error /= err_none ) THEN
1395  CALL errorstop(global,err_tecplot_output,__line__)
1396  END IF ! global%error
1397  END DO ! iVar
1398 
1399 ! ==============================================================================
1400 ! Connectivity
1401 ! ==============================================================================
1402 
1403  errorflag = tecnod100(f2vtmp)
1404  global%error = errorflag
1405  IF ( global%error /= err_none ) THEN
1406  CALL errorstop(global,err_tecplot_output,__line__)
1407  END IF ! global%error
1408 
1409  DEALLOCATE(f2vtmp,stat=errorflag)
1410  global%error = errorflag
1411  IF ( global%error /= err_none ) THEN
1412  CALL errorstop(global,err_deallocate,__line__,'f2vTmp')
1413  END IF ! global%error
1414  END IF ! nQuads
1415 
1416  DEALLOCATE(f2v,stat=errorflag)
1417  global%error = errorflag
1418  IF ( global%error /= err_none ) THEN
1419  CALL errorstop(global,err_deallocate,__line__,'f2v')
1420  END IF ! global%error
1421 
1422 ! ******************************************************************************
1423 ! Deallocate temporary memory
1424 ! ******************************************************************************
1425 
1426  DEALLOCATE(bvflag,stat=errorflag)
1427  global%error = errorflag
1428  IF ( global%error /= err_none ) THEN
1429  CALL errorstop(global,err_deallocate,__line__,'bvFlag')
1430  END IF ! global%error
1431 
1432  DEALLOCATE(varverttmp,stat=errorflag)
1433  global%error = errorflag
1434  IF ( global%error /= err_none ) THEN
1435  CALL errorstop(global,err_deallocate,__line__,'varVertTmp')
1436  END IF ! global%error
1437 
1438  IF ( nvarsfacetec > 0 ) THEN
1439  DEALLOCATE(varfacetmp,stat=errorflag)
1440  global%error = errorflag
1441  IF ( global%error /= err_none ) THEN
1442  CALL errorstop(global,err_deallocate,__line__,'varFaceTmp')
1443  END IF ! global%error
1444  END IF ! nVarsFaceTEC
1445 
1446 ! ******************************************************************************
1447 ! Destroy null array
1448 ! ******************************************************************************
1449 
1450  DEALLOCATE(nullarray,stat=errorflag)
1451  global%error = errorflag
1452  IF ( global%error /= err_none ) THEN
1453  CALL errorstop(global,err_deallocate,__line__,'nullArray')
1454  END IF ! global%error
1455 
1456 ! ******************************************************************************
1457 ! End
1458 ! ******************************************************************************
1459 
1460  CALL deregisterfunction(global)
1461 
1462 END SUBROUTINE rflu_tec_writezoneinterf
1463 
1464 
1465 
1466 
1467 
1468 
1469 
1470 
1471 
1472 ! ******************************************************************************
1473 !
1474 ! Purpose: Write zone for surface data.
1475 !
1476 ! Description: None.
1477 !
1478 ! Input:
1479 ! pRegion Pointer to region
1480 ! pPatch Pointer to patch
1481 ! faceType Face type
1482 ! faceKind Face kind
1483 !
1484 ! Output: None.
1485 !
1486 ! Notes:
1487 ! 1. This routine combines writing of surface data irrespective of whether
1488 ! volume data is written.
1489 !
1490 ! ******************************************************************************
1491 
1492 SUBROUTINE rflu_tec_writezonesurf(pRegion,pPatch,faceType,faceKind)
1493 
1494  USE modsortsearch
1495 
1496  IMPLICIT NONE
1497 
1498 ! ******************************************************************************
1499 ! Declarations and definitions
1500 ! ******************************************************************************
1501 
1502 ! ==============================================================================
1503 ! Arguments
1504 ! ==============================================================================
1505 
1506  INTEGER, INTENT(IN) :: facekind,facetype
1507  TYPE(t_patch), POINTER :: ppatch
1508  TYPE(t_region), POINTER :: pregion
1509 
1510 ! ==============================================================================
1511 ! Locals
1512 ! ==============================================================================
1513 
1514  CHARACTER(CHRLEN) :: zoneprefix,zonetitle
1515  INTEGER :: errorflag,ifg,ifgbeg,ifgend,ifgoff,ifg2,ifl,iloc,ivar,ivarface, &
1516  ivarvert,ivg,ivlf,ivlp,nbfaces,nbvert,zonetype
1517  INTEGER, DIMENSION(:), ALLOCATABLE :: bvflag,nullarray
1518  INTEGER, DIMENSION(:,:), ALLOCATABLE :: bf2vtmp
1519  REAL(RFREAL), DIMENSION(:,:), ALLOCATABLE :: varfacetmp,varverttmp
1520  TYPE(t_global), POINTER :: global
1521  TYPE(t_grid), POINTER :: pgrid
1522 
1523 ! ==============================================================================
1524 ! Externals: TECPLOT functions
1525 ! ==============================================================================
1526 
1527  INTEGER, EXTERNAL :: tecdat100,tecnod100,teczne100
1528 
1529 ! ******************************************************************************
1530 ! Start
1531 ! ******************************************************************************
1532 
1533  global => pregion%global
1534 
1535  CALL registerfunction(global,'RFLU_TEC_WriteZoneSurf', &
1536  'RFLU_ModTECPLOTUtils.F90')
1537 
1538 ! ******************************************************************************
1539 ! Set grid pointer
1540 ! ******************************************************************************
1541 
1542  pgrid => pregion%grid
1543 
1544 ! ******************************************************************************
1545 ! Create null array
1546 ! ******************************************************************************
1547 
1548  ALLOCATE(nullarray(nvarstec),stat=errorflag)
1549  global%error = errorflag
1550  IF ( global%error /= err_none ) THEN
1551  CALL errorstop(global,err_allocate,__line__,'nullArray')
1552  END IF ! global%error
1553 
1554  DO ivar = 1,nvarstec
1555  nullarray(ivar) = 0
1556  END DO ! iVar
1557 
1558 ! ******************************************************************************
1559 ! Set zone type, number of faces, and zone prefix. Face offset must be defined
1560 ! so that boundary face data can be accessed properly. The difference arises
1561 ! because the connectivity is stored for each face type, whereas the data is
1562 ! stored for the combined face list (actual triangles, actual quads, virtual
1563 ! triangles, virtual quads).
1564 ! ******************************************************************************
1565 
1566  SELECT CASE ( facekind )
1567  CASE ( face_kind_ab )
1568  SELECT CASE ( facetype )
1569  CASE ( face_type_tri )
1570  ifgbeg = 1
1571  ifgend = ppatch%nBTris
1572  ifgoff = 0
1573  zonetype = zone_type_tri
1574  zoneprefix = 'TRI-A_'
1575  CASE ( face_type_quad )
1576  ifgbeg = 1
1577  ifgend = ppatch%nBQuads
1578  ifgoff = ppatch%nBTris
1579  zonetype = zone_type_quad
1580  zoneprefix = 'QUAD-A_'
1581  CASE default
1582  CALL errorstop(global,err_reached_default,__line__)
1583  END SELECT ! faceType
1584  CASE ( face_kind_vb )
1585  SELECT CASE ( facetype )
1586  CASE ( face_type_tri )
1587  ifgbeg = ppatch%nBTris+1
1588  ifgend = ppatch%nBTrisTot
1589  ifgoff = ppatch%nBTris+ppatch%nBQuads
1590  zonetype = zone_type_tri
1591  zoneprefix = 'TRI-V_'
1592  CASE ( face_type_quad )
1593  ifgbeg = ppatch%nBQuads+1
1594  ifgend = ppatch%nBQuadsTot
1595  ifgoff = ppatch%nBTrisTot+ppatch%nBQuads
1596  zonetype = zone_type_quad
1597  zoneprefix = 'QUAD-V_'
1598  CASE default
1599  CALL errorstop(global,err_reached_default,__line__)
1600  END SELECT ! cellType
1601  CASE default
1602  CALL errorstop(global,err_reached_default,__line__)
1603  END SELECT ! cellKind
1604 
1605  nbfaces = ifgend - ifgbeg + 1
1606 
1607 ! ******************************************************************************
1608 ! Allocate temporary memory
1609 ! ******************************************************************************
1610 
1611  ALLOCATE(bvflag(ppatch%nBVertTot),stat=errorflag)
1612  global%error = errorflag
1613  IF ( global%error /= err_none ) THEN
1614  CALL errorstop(global,err_allocate,__line__,'bvFlag')
1615  END IF ! global%error
1616 
1617  IF ( nvarsfacetec > 0 ) THEN
1618  ALLOCATE(varfacetmp(nbfaces,nvarsfacetec),stat=errorflag)
1619  global%error = errorflag
1620  IF ( global%error /= err_none ) THEN
1621  CALL errorstop(global,err_allocate,__line__,'varFaceTmp')
1622  END IF ! global%error
1623  END IF ! nVarsFaceTEC
1624 
1625 ! ******************************************************************************
1626 ! Construct local connectivity array. NOTE this loop computes nBVert, which will
1627 ! differ from pPatch%nBVert if have triangular and quadrilateral faces on same
1628 ! patch or actual and virtual faces on same patch.
1629 ! ******************************************************************************
1630 
1631 ! ==============================================================================
1632 ! Triangles
1633 ! ==============================================================================
1634 
1635  IF ( facetype == face_type_tri ) THEN
1636  ALLOCATE(bf2vtmp(3,nbfaces),stat=errorflag)
1637  global%error = errorflag
1638  IF ( global%error /= err_none ) THEN
1639  CALL errorstop(global,err_allocate,__line__,'bf2vTmp')
1640  END IF ! global%error
1641 
1642  DO ivlp = 1,ppatch%nBVertTot
1643  bvflag(ivlp) = 0
1644  END DO ! ivlp
1645 
1646  ifl = 0
1647  nbvert = 0
1648 
1649 ! ------------------------------------------------------------------------------
1650 ! Loop over faces
1651 ! ------------------------------------------------------------------------------
1652 
1653  DO ifg = ifgbeg,ifgend
1654  ifl = ifl + 1
1655  ifg2 = ifg - ifgbeg + ifgoff + 1
1656 
1657  DO ivlf = 1,3
1658  ivg = ppatch%bTri2v(ivlf,ifg)
1659 
1660  CALL binarysearchinteger(ppatch%bv,ppatch%nBVertTot,ivg,ivlp)
1661 
1662  IF ( ivlp == element_not_found ) THEN
1663  CALL errorstop(global,err_binary_search,__line__)
1664  END IF ! ivlp
1665 
1666  IF ( bvflag(ivlp) == 0 ) THEN
1667  nbvert = nbvert + 1
1668  bvflag(ivlp) = nbvert
1669  END IF ! bvFlag
1670 
1671 ! ----- Define connectivity ----------------------------------------------------
1672 
1673  bf2vtmp(ivlf,ifl) = bvflag(ivlp)
1674 
1675 ! ----- Copy face data ---------------------------------------------------------
1676 
1677  ivarface = 0
1678 
1679  DO ivar = 1,nvarstec
1680  IF ( postec(ivar) == var_pos_face ) THEN
1681  ivarface = ivarface + 1
1682  varfacetmp(ifl,ivarface) = ppatch%varFaceTEC(ifg2,ivarface)
1683  END IF ! posTEC
1684  END DO ! iVar
1685  END DO ! ivlf
1686  END DO ! ifg
1687 
1688 ! ==============================================================================
1689 ! Quadrilaterals
1690 ! ==============================================================================
1691 
1692  ELSE
1693  ALLOCATE(bf2vtmp(4,nbfaces),stat=errorflag)
1694  global%error = errorflag
1695  IF ( global%error /= err_none ) THEN
1696  CALL errorstop(global,err_allocate,__line__,'bf2vTmp')
1697  END IF ! global%error
1698 
1699  DO ivlp = 1,ppatch%nBVertTot
1700  bvflag(ivlp) = 0
1701  END DO ! ivlp
1702 
1703  ifl = 0
1704  nbvert = 0
1705 
1706 ! ------------------------------------------------------------------------------
1707 ! Loop over faces
1708 ! ------------------------------------------------------------------------------
1709 
1710  DO ifg = ifgbeg,ifgend
1711  ifl = ifl + 1
1712  ifg2 = ifg - ifgbeg + ifgoff + 1
1713 
1714  DO ivlf = 1,4
1715  ivg = ppatch%bQuad2v(ivlf,ifg)
1716 
1717  CALL binarysearchinteger(ppatch%bv,ppatch%nBVertTot,ivg,ivlp)
1718 
1719  IF ( ivlp == element_not_found ) THEN
1720  CALL errorstop(global,err_binary_search,__line__)
1721  END IF ! ivlp
1722 
1723  IF ( bvflag(ivlp) == 0 ) THEN
1724  nbvert = nbvert + 1
1725  bvflag(ivlp) = nbvert
1726  END IF ! bvFlag
1727 
1728 ! ----- Define connectivity ----------------------------------------------------
1729 
1730  bf2vtmp(ivlf,ifl) = bvflag(ivlp)
1731 
1732 ! ----- Copy face data ---------------------------------------------------------
1733 
1734  ivarface = 0
1735 
1736  DO ivar = 1,nvarstec
1737  IF ( postec(ivar) == var_pos_face ) THEN
1738  ivarface = ivarface + 1
1739  varfacetmp(ifl,ivarface) = ppatch%varFaceTEC(ifg2,ivarface)
1740  END IF ! posTEC
1741  END DO ! iVar
1742  END DO ! ivlf
1743  END DO ! ifg
1744  END IF ! faceType
1745 
1746 ! ******************************************************************************
1747 ! Allocate memory for and copy vertex data
1748 ! ******************************************************************************
1749 
1750  ALLOCATE(varverttmp(nbvert,nvarsverttec),stat=errorflag)
1751  global%error = errorflag
1752  IF ( global%error /= err_none ) THEN
1753  CALL errorstop(global,err_allocate,__line__,'varVertTmp')
1754  END IF ! global%error
1755 
1756  DO ivlp = 1,ppatch%nBVertTot
1757  IF ( bvflag(ivlp) /= 0 ) THEN
1758  ivarvert = 0
1759 
1760  DO ivar = 1,nvarstec
1761  IF ( postec(ivar) == var_pos_vert ) THEN
1762  ivarvert = ivarvert + 1
1763  varverttmp(bvflag(ivlp),ivarvert) = ppatch%varVertTEC(ivlp,ivarvert)
1764  END IF ! posTEC
1765  END DO ! iVar
1766  END IF ! bvFlag
1767  END DO ! ivlp
1768 
1769 ! ******************************************************************************
1770 ! Write to TECPLOT file
1771 ! ******************************************************************************
1772 
1773 ! ==============================================================================
1774 ! Zone information
1775 ! ==============================================================================
1776 
1777  WRITE(zonetitle,'(A,I3.3,A,I5.5)') 'PAT_',ppatch%iPatchGlobal, &
1778  '_'//trim(zoneprefix), &
1779  pregion%iRegionGlobal
1780 
1781  errorflag = teczne100(trim(zonetitle)//nullchar,zonetype,nbvert,nbfaces, &
1782  0,0,0,0,zone_form_block,0,0,postec,nullarray,0)
1783  global%error = errorflag
1784  IF ( global%error /= err_none ) THEN
1785  CALL errorstop(global,err_tecplot_output,__line__)
1786  END IF ! global%error
1787 
1788 ! ==============================================================================
1789 ! Variables
1790 ! ==============================================================================
1791 
1792  ivarface = 0
1793  ivarvert = 0
1794 
1795  DO ivar = 1,nvarstec
1796  IF ( postec(ivar) == var_pos_face ) THEN
1797  ivarface = ivarface + 1
1798  errorflag = tecdat100(nbfaces,varfacetmp(1,ivarface),doubleflagtec)
1799  ELSE IF ( postec(ivar) == var_pos_vert ) THEN
1800  ivarvert = ivarvert + 1
1801  errorflag = tecdat100(nbvert,varverttmp(1,ivarvert),doubleflagtec)
1802  ELSE ! defensive coding
1803  CALL errorstop(global,err_reached_default,__line__)
1804  END IF ! posTEC
1805 
1806  global%error = errorflag
1807  IF ( global%error /= err_none ) THEN
1808  CALL errorstop(global,err_tecplot_output,__line__)
1809  END IF ! global%error
1810  END DO ! iVar
1811 
1812 ! ==============================================================================
1813 ! Connectivity
1814 ! ==============================================================================
1815 
1816  errorflag = tecnod100(bf2vtmp)
1817  global%error = errorflag
1818  IF ( global%error /= err_none ) THEN
1819  CALL errorstop(global,err_tecplot_output,__line__)
1820  END IF ! global%error
1821 
1822 ! ******************************************************************************
1823 ! Deallocate temporary memory
1824 ! ******************************************************************************
1825 
1826  DEALLOCATE(bvflag,stat=errorflag)
1827  global%error = errorflag
1828  IF ( global%error /= err_none ) THEN
1829  CALL errorstop(global,err_deallocate,__line__,'bvFlag')
1830  END IF ! global%error
1831 
1832  IF ( ALLOCATED(varfacetmp) .EQV. .true. ) THEN
1833  DEALLOCATE(varfacetmp,stat=errorflag)
1834  global%error = errorflag
1835  IF ( global%error /= err_none ) THEN
1836  CALL errorstop(global,err_deallocate,__line__,'varFaceTmp')
1837  END IF ! global%error
1838  END IF ! ALLOCATED
1839 
1840  DEALLOCATE(varverttmp,stat=errorflag)
1841  global%error = errorflag
1842  IF ( global%error /= err_none ) THEN
1843  CALL errorstop(global,err_deallocate,__line__,'varVertTmp')
1844  END IF ! global%error
1845 
1846  DEALLOCATE(bf2vtmp,stat=errorflag)
1847  global%error = errorflag
1848  IF ( global%error /= err_none ) THEN
1849  CALL errorstop(global,err_deallocate,__line__,'bf2vTmp')
1850  END IF ! global%error
1851 
1852 ! ******************************************************************************
1853 ! Destroy null array
1854 ! ******************************************************************************
1855 
1856  DEALLOCATE(nullarray,stat=errorflag)
1857  global%error = errorflag
1858  IF ( global%error /= err_none ) THEN
1859  CALL errorstop(global,err_deallocate,__line__,'nullArray')
1860  END IF ! global%error
1861 
1862 ! ******************************************************************************
1863 ! End
1864 ! ******************************************************************************
1865 
1866  CALL deregisterfunction(global)
1867 
1868 END SUBROUTINE rflu_tec_writezonesurf
1869 
1870 
1871 
1872 
1873 
1874 
1875 
1876 
1877 
1878 ! ******************************************************************************
1879 !
1880 ! Purpose: Write zone for surface data in mixed format.
1881 !
1882 ! Description: None.
1883 !
1884 ! Input:
1885 ! pRegion Pointer to region
1886 ! pPatch Pointer to patch
1887 !
1888 ! Output: None.
1889 !
1890 ! Notes:
1891 ! 1. The term mixed format refers to the list of faces being written out in
1892 ! a single block, without distinguishing between triangular and
1893 ! quadrilateral faces.
1894 ! 2. Triangular faces are written as degenerate quadrilateral faces.
1895 ! 3. Only actual faces are written because face data is typically not defined
1896 ! for virtual faces.
1897 !
1898 ! ******************************************************************************
1899 
1900 SUBROUTINE rflu_tec_writezonesurfmixed(pRegion,pPatch)
1901 
1902  USE modsortsearch
1903 
1904  IMPLICIT NONE
1905 
1906 ! ******************************************************************************
1907 ! Declarations and definitions
1908 ! ******************************************************************************
1909 
1910 ! ==============================================================================
1911 ! Arguments
1912 ! ==============================================================================
1913 
1914  TYPE(t_patch), POINTER :: ppatch
1915  TYPE(t_region), POINTER :: pregion
1916 
1917 ! ==============================================================================
1918 ! Locals
1919 ! ==============================================================================
1920 
1921  CHARACTER(CHRLEN) :: zonetitle
1922  INTEGER :: errorflag,ifl,ivar,ivarface,ivarvert,ivl
1923  INTEGER, DIMENSION(:), ALLOCATABLE :: nullarray
1924  INTEGER, DIMENSION(:,:), ALLOCATABLE :: bf2vtmp
1925  TYPE(t_global), POINTER :: global
1926  TYPE(t_grid), POINTER :: pgrid
1927 
1928 ! ==============================================================================
1929 ! Externals: TECPLOT functions
1930 ! ==============================================================================
1931 
1932  INTEGER, EXTERNAL :: tecdat100,tecnod100,teczne100
1933 
1934 ! ******************************************************************************
1935 ! Start
1936 ! ******************************************************************************
1937 
1938  global => pregion%global
1939 
1940  CALL registerfunction(global,'RFLU_TEC_WriteZoneSurfMixed', &
1941  'RFLU_ModTECPLOTUtils.F90')
1942 
1943 ! ******************************************************************************
1944 ! Set grid pointer
1945 ! ******************************************************************************
1946 
1947  pgrid => pregion%grid
1948 
1949 ! ******************************************************************************
1950 ! Create null array
1951 ! ******************************************************************************
1952 
1953  ALLOCATE(nullarray(nvarstec),stat=errorflag)
1954  global%error = errorflag
1955  IF ( global%error /= err_none ) THEN
1956  CALL errorstop(global,err_allocate,__line__,'nullArray')
1957  END IF ! global%error
1958 
1959  DO ivar = 1,nvarstec
1960  nullarray(ivar) = 0
1961  END DO ! iVar
1962 
1963 ! ******************************************************************************
1964 ! Allocate temporary memory
1965 ! ******************************************************************************
1966 
1967  ALLOCATE(bf2vtmp(4,ppatch%nBFaces),stat=errorflag)
1968  global%error = errorflag
1969  IF ( global%error /= err_none ) THEN
1970  CALL errorstop(global,err_allocate,__line__,'bf2vTmp')
1971  END IF ! global%error
1972 
1973 ! ******************************************************************************
1974 ! Convert connectivity array
1975 ! ******************************************************************************
1976 
1977  DO ifl = 1,ppatch%nBFaces
1978  bf2vtmp(1,ifl) = ppatch%bf2v(1,ifl)
1979  bf2vtmp(2,ifl) = ppatch%bf2v(2,ifl)
1980  bf2vtmp(3,ifl) = ppatch%bf2v(3,ifl)
1981 
1982  IF ( ppatch%bf2v(4,ifl) /= vert_none ) THEN
1983  bf2vtmp(4,ifl) = ppatch%bf2v(4,ifl)
1984  ELSE
1985  bf2vtmp(4,ifl) = bf2vtmp(1,ifl)
1986  END IF ! v4l
1987  END DO ! ifl
1988 
1989 ! ******************************************************************************
1990 ! Write to TECPLOT file
1991 ! ******************************************************************************
1992 
1993 ! ==============================================================================
1994 ! Zone information
1995 ! ==============================================================================
1996 
1997  WRITE(zonetitle,'(A,I3.3,A,I5.5)') 'PAT_',ppatch%iPatchGlobal, &
1998  '_MIX-A_',pregion%iRegionGlobal
1999 
2000  errorflag = teczne100(trim(zonetitle)//nullchar,zone_type_quad, &
2001  ppatch%nBVert,ppatch%nBFaces,0,0,0,0,zone_form_block, &
2002  0,0,postec,nullarray,0)
2003  global%error = errorflag
2004  IF ( global%error /= err_none ) THEN
2005  CALL errorstop(global,err_tecplot_output,__line__)
2006  END IF ! global%error
2007 
2008 ! ==============================================================================
2009 ! Variables
2010 ! ==============================================================================
2011 
2012  ivarface = 0
2013  ivarvert = 0
2014 
2015  DO ivar = 1,nvarstec
2016  IF ( postec(ivar) == var_pos_face ) THEN
2017  ivarface = ivarface + 1
2018  errorflag = tecdat100(ppatch%nBFaces,ppatch%varFaceTEC(1,ivarface), &
2019  doubleflagtec)
2020  ELSE IF ( postec(ivar) == var_pos_vert ) THEN
2021  ivarvert = ivarvert + 1
2022  errorflag = tecdat100(ppatch%nBVert,ppatch%varVertTEC(1,ivarvert), &
2023  doubleflagtec)
2024  ELSE ! defensive coding
2025  CALL errorstop(global,err_reached_default,__line__)
2026  END IF ! posTEC
2027 
2028  global%error = errorflag
2029  IF ( global%error /= err_none ) THEN
2030  CALL errorstop(global,err_tecplot_output,__line__)
2031  END IF ! global%error
2032  END DO ! iVar
2033 
2034 ! ==============================================================================
2035 ! Connectivity
2036 ! ==============================================================================
2037 
2038  errorflag = tecnod100(bf2vtmp)
2039  global%error = errorflag
2040  IF ( global%error /= err_none ) THEN
2041  CALL errorstop(global,err_tecplot_output,__line__)
2042  END IF ! global%error
2043 
2044 ! ******************************************************************************
2045 ! Deallocate temporary memory
2046 ! ******************************************************************************
2047 
2048  DEALLOCATE(bf2vtmp,stat=errorflag)
2049  global%error = errorflag
2050  IF ( global%error /= err_none ) THEN
2051  CALL errorstop(global,err_deallocate,__line__,'bf2vTmp')
2052  END IF ! global%error
2053 
2054 ! ******************************************************************************
2055 ! Destroy null array
2056 ! ******************************************************************************
2057 
2058  DEALLOCATE(nullarray,stat=errorflag)
2059  global%error = errorflag
2060  IF ( global%error /= err_none ) THEN
2061  CALL errorstop(global,err_deallocate,__line__,'nullArray')
2062  END IF ! global%error
2063 
2064 ! ******************************************************************************
2065 ! End
2066 ! ******************************************************************************
2067 
2068  CALL deregisterfunction(global)
2069 
2070 END SUBROUTINE rflu_tec_writezonesurfmixed
2071 
2072 
2073 
2074 
2075 
2076 
2077 ! ******************************************************************************
2078 !
2079 ! Purpose: Write zone for volume data.
2080 !
2081 ! Description: None.
2082 !
2083 ! Input:
2084 ! pRegion Pointer to region
2085 ! cellType Cell type
2086 ! cellKind Cell kind
2087 !
2088 ! Output: None.
2089 !
2090 ! Notes:
2091 ! 1. Prisms and pyramids are treated as degenerate hexahedra.
2092 !
2093 ! ******************************************************************************
2094 
2095 SUBROUTINE rflu_tec_writezonevol(pRegion,cellType,cellKind)
2096 
2097  IMPLICIT NONE
2098 
2099 ! ******************************************************************************
2100 ! Declarations and definitions
2101 ! ******************************************************************************
2102 
2103 ! ==============================================================================
2104 ! Arguments
2105 ! ==============================================================================
2106 
2107  INTEGER, INTENT(IN) :: cellkind,celltype
2108  TYPE(t_region), POINTER :: pregion
2109 
2110 ! ==============================================================================
2111 ! Locals
2112 ! ==============================================================================
2113 
2114  CHARACTER(CHRLEN) :: zoneprefix,zonetitle
2115  INTEGER :: errorflag,icg,icgbeg,icgend,icl,ivar,ivarcell,ivarvert,ncells, &
2116  nvert,zonetype
2117  INTEGER, DIMENSION(:), ALLOCATABLE :: nullarray
2118  INTEGER, DIMENSION(:), POINTER :: pxyz2cellglob
2119  INTEGER, DIMENSION(:,:), ALLOCATABLE :: c2v
2120  REAL(RFREAL), DIMENSION(:,:), ALLOCATABLE :: varcelltmp
2121  TYPE(t_global), POINTER :: global
2122  TYPE(t_grid), POINTER :: pgrid
2123 
2124 ! ==============================================================================
2125 ! Externals: TECPLOT functions
2126 ! ==============================================================================
2127 
2128  INTEGER, EXTERNAL :: tecdat100,tecnod100,teczne100
2129 
2130 ! ******************************************************************************
2131 ! Start
2132 ! ******************************************************************************
2133 
2134  global => pregion%global
2135 
2136  CALL registerfunction(global,'RFLU_TEC_WriteZoneVol', &
2137  'RFLU_ModTECPLOTUtils.F90')
2138 
2139 ! ******************************************************************************
2140 ! Set grid pointer
2141 ! ******************************************************************************
2142 
2143  pgrid => pregion%grid
2144 
2145 ! ******************************************************************************
2146 ! Create null array
2147 ! ******************************************************************************
2148 
2149  ALLOCATE(nullarray(nvarstec),stat=errorflag)
2150  global%error = errorflag
2151  IF ( global%error /= err_none ) THEN
2152  CALL errorstop(global,err_allocate,__line__,'nullArray')
2153  END IF ! global%error
2154 
2155  DO ivar = 1,nvarstec
2156  nullarray(ivar) = 0
2157  END DO ! iVar
2158 
2159 ! ******************************************************************************
2160 ! Set zone type, number of cells, and zone title
2161 ! ******************************************************************************
2162 
2163  SELECT CASE ( cellkind )
2164  CASE ( cell_kind_actual )
2165  SELECT CASE ( celltype )
2166  CASE ( cell_type_tet )
2167  icgbeg = 1
2168  icgend = pgrid%nTets
2169  nvert = pgrid%nVertTot
2170  zonetype = zone_type_tet
2171  zoneprefix = 'TET-A_'
2172  CASE ( cell_type_hex )
2173  icgbeg = 1
2174  icgend = pgrid%nHexs
2175  nvert = pgrid%nVertTot
2176  zonetype = zone_type_hex
2177  zoneprefix = 'HEX-A_'
2178  CASE ( cell_type_pri )
2179  icgbeg = 1
2180  icgend = pgrid%nPris
2181  nvert = pgrid%nVertTot
2182  zonetype = zone_type_hex
2183  zoneprefix = 'PRI-A_'
2184  CASE ( cell_type_pyr )
2185  icgbeg = 1
2186  icgend = pgrid%nPyrs
2187  nvert = pgrid%nVertTot
2188  zonetype = zone_type_hex
2189  zoneprefix = 'PYR-A_'
2190  CASE default
2191  CALL errorstop(global,err_reached_default,__line__)
2192  END SELECT ! cellType
2193  CASE ( cell_kind_virtual )
2194  SELECT CASE ( celltype )
2195  CASE ( cell_type_tet )
2196  icgbeg = pgrid%nTets+1
2197  icgend = pgrid%nTetsTot
2198  nvert = pgrid%nVertTot
2199  zonetype = zone_type_tet
2200  zoneprefix = 'TET-V_'
2201  CASE ( cell_type_hex )
2202  icgbeg = pgrid%nHexs+1
2203  icgend = pgrid%nHexsTot
2204  nvert = pgrid%nVertTot
2205  zonetype = zone_type_hex
2206  zoneprefix = 'HEX-V_'
2207  CASE ( cell_type_pri )
2208  icgbeg = pgrid%nPris+1
2209  icgend = pgrid%nPrisTot
2210  nvert = pgrid%nVertTot
2211  zonetype = zone_type_hex
2212  zoneprefix = 'PRI-V_'
2213  CASE ( cell_type_pyr )
2214  icgbeg = pgrid%nPyrs+1
2215  icgend = pgrid%nPyrsTot
2216  nvert = pgrid%nVertTot
2217  zonetype = zone_type_hex
2218  zoneprefix = 'PYR-V_'
2219  CASE default
2220  CALL errorstop(global,err_reached_default,__line__)
2221  END SELECT ! cellType
2222  CASE default
2223  CALL errorstop(global,err_reached_default,__line__)
2224  END SELECT ! cellKind
2225 
2226  ncells = icgend - icgbeg + 1
2227 
2228 ! ******************************************************************************
2229 ! Convert cell array. This is required whenever have mixed grid because nCells
2230 ! would then not match pGrid%nCellsTot, which is the dimension of the array
2231 ! pRegion%varCellTEC.
2232 ! ******************************************************************************
2233 
2234  IF ( nvarscelltec > 0 ) THEN
2235  ALLOCATE(varcelltmp(ncells,nvarscelltec),stat=errorflag)
2236  global%error = errorflag
2237  IF ( global%error /= err_none ) THEN
2238  CALL errorstop(global,err_allocate,__line__,'varCellTmp')
2239  END IF ! global%error
2240 
2241  SELECT CASE ( celltype )
2242  CASE ( cell_type_tet )
2243  pxyz2cellglob => pgrid%tet2CellGlob
2244  CASE ( cell_type_hex )
2245  pxyz2cellglob => pgrid%hex2CellGlob
2246  CASE ( cell_type_pri )
2247  pxyz2cellglob => pgrid%pri2CellGlob
2248  CASE ( cell_type_pyr )
2249  pxyz2cellglob => pgrid%pyr2CellGlob
2250  CASE default
2251  CALL errorstop(global,err_reached_default,__line__)
2252  END SELECT ! cellType
2253 
2254  DO icg = icgbeg,icgend
2255  icl = icg - icgbeg + 1
2256 
2257  DO ivarcell = 1,nvarscelltec
2258  varcelltmp(icl,ivarcell) = &
2259  pregion%varCellTEC(pxyz2cellglob(icg),ivarcell)
2260  END DO ! iVarCell
2261  END DO ! icl
2262  END IF ! nVarsCellTEC
2263 
2264 ! ******************************************************************************
2265 ! Convert connectivity for prisms and pyramids
2266 ! ******************************************************************************
2267 
2268  IF ( celltype == cell_type_pri .OR. celltype == cell_type_pyr ) THEN
2269  ALLOCATE(c2v(8,ncells),stat=errorflag)
2270  global%error = errorflag
2271  IF ( global%error /= err_none ) THEN
2272  CALL errorstop(global,err_allocate,__line__,'c2v')
2273  END IF ! global%error
2274 
2275  IF ( celltype == cell_type_pri ) THEN ! Prisms
2276  DO icg = icgbeg,icgend
2277  icl = icg - icgbeg + 1
2278  c2v(1:3,icl) = pgrid%pri2v(1:3,icg)
2279  c2v(4 ,icl) = pgrid%pri2v(1 ,icg)
2280  c2v(5:7,icl) = pgrid%pri2v(4:6,icg)
2281  c2v(8 ,icl) = pgrid%pri2v(4 ,icg)
2282  END DO ! icg
2283  ELSE ! Pyramids
2284  DO icg = icgbeg,icgend
2285  icl = icg - icgbeg + 1
2286  c2v(1:5,icl) = pgrid%pyr2v(1:5,icg)
2287  c2v(6 ,icl) = pgrid%pyr2v(5 ,icg)
2288  c2v(7 ,icl) = pgrid%pyr2v(5 ,icg)
2289  c2v(8 ,icl) = pgrid%pyr2v(5 ,icg)
2290  END DO ! icg
2291  END IF ! cellType
2292  END IF ! cellType
2293 
2294 ! ******************************************************************************
2295 ! Write to TECPLOT file
2296 ! ******************************************************************************
2297 
2298 ! ==============================================================================
2299 ! Zone information
2300 ! ==============================================================================
2301 
2302  WRITE(zonetitle,'(A,I5.5)') trim(zoneprefix),pregion%iRegionGlobal
2303 
2304  errorflag = teczne100(trim(zonetitle)//nullchar,zonetype,nvert,ncells, &
2305  0,0,0,0,zone_form_block,0,0,postec,nullarray,0)
2306  global%error = errorflag
2307  IF ( global%error /= err_none ) THEN
2308  CALL errorstop(global,err_tecplot_output,__line__)
2309  END IF ! global%error
2310 
2311 ! ==============================================================================
2312 ! Variables
2313 ! ==============================================================================
2314 
2315  ivarcell = 0
2316  ivarvert = 0
2317 
2318  DO ivar = 1,nvarstec
2319  IF ( postec(ivar) == var_pos_cell ) THEN
2320  ivarcell = ivarcell + 1
2321  errorflag = tecdat100(ncells,varcelltmp(1,ivarcell),doubleflagtec)
2322  ELSE IF ( postec(ivar) == var_pos_vert ) THEN
2323  ivarvert = ivarvert + 1
2324  errorflag = tecdat100(pgrid%nVertTot,pregion%varVertTEC(1,ivarvert), &
2325  doubleflagtec)
2326  ELSE ! defensive coding
2327  CALL errorstop(global,err_reached_default,__line__)
2328  END IF ! posTEC
2329 
2330  global%error = errorflag
2331  IF ( global%error /= err_none ) THEN
2332  CALL errorstop(global,err_tecplot_output,__line__)
2333  END IF ! global%error
2334  END DO ! iVar
2335 
2336 ! ==============================================================================
2337 ! Connectivity
2338 ! ==============================================================================
2339 
2340  SELECT CASE ( celltype )
2341  CASE ( cell_type_tet )
2342  errorflag = tecnod100(pgrid%tet2v(1:4,icgbeg:icgend))
2343  CASE ( cell_type_hex )
2344  errorflag = tecnod100(pgrid%hex2v(1:8,icgbeg:icgend))
2345  CASE ( cell_type_pri )
2346  errorflag = tecnod100(c2v)
2347  CASE ( cell_type_pyr )
2348  errorflag = tecnod100(c2v)
2349  CASE default
2350  CALL errorstop(global,err_reached_default,__line__)
2351  END SELECT ! cellType
2352 
2353  global%error = errorflag
2354  IF ( global%error /= err_none ) THEN
2355  CALL errorstop(global,err_tecplot_output,__line__)
2356  END IF ! global%error
2357 
2358 ! ******************************************************************************
2359 ! Deallocate temporary memory
2360 ! ******************************************************************************
2361 
2362  IF ( nvarscelltec > 0 ) THEN
2363  DEALLOCATE(varcelltmp,stat=errorflag)
2364  global%error = errorflag
2365  IF ( global%error /= err_none ) THEN
2366  CALL errorstop(global,err_deallocate,__line__,'varCellTmp')
2367  END IF ! global%error
2368  END IF ! nVarsCellTEC
2369 
2370  IF ( celltype == cell_type_pri .OR. celltype == cell_type_pyr ) THEN
2371  DEALLOCATE(c2v,stat=errorflag)
2372  global%error = errorflag
2373  IF ( global%error /= err_none ) THEN
2374  CALL errorstop(global,err_deallocate,__line__,'c2v')
2375  END IF ! global%error
2376  END IF ! cellType
2377 
2378 ! ******************************************************************************
2379 ! Destroy null array
2380 ! ******************************************************************************
2381 
2382  DEALLOCATE(nullarray,stat=errorflag)
2383  global%error = errorflag
2384  IF ( global%error /= err_none ) THEN
2385  CALL errorstop(global,err_deallocate,__line__,'nullArray')
2386  END IF ! global%error
2387 
2388 ! ******************************************************************************
2389 ! End
2390 ! ******************************************************************************
2391 
2392  CALL deregisterfunction(global)
2393 
2394 END SUBROUTINE rflu_tec_writezonevol
2395 
2396 
2397 
2398 
2399 
2400 
2401 END MODULE rflu_modtecplotutils
2402 
2403 ! ******************************************************************************
2404 !
2405 ! RCS Revision history:
2406 !
2407 ! $Log: RFLU_ModTECPLOTUtils.F90,v $
2408 ! Revision 1.14 2008/12/06 08:45:06 mtcampbe
2409 ! Updated license.
2410 !
2411 ! Revision 1.13 2008/11/19 22:18:17 mtcampbe
2412 ! Added Illinois Open Source License/Copyright
2413 !
2414 ! Revision 1.12 2006/04/07 15:19:26 haselbac
2415 ! Removed tabs
2416 !
2417 ! Revision 1.11 2005/09/23 19:02:21 haselbac
2418 ! Added FILE_TYPE_PATCH_STATS
2419 !
2420 ! Revision 1.10 2005/03/16 20:31:44 haselbac
2421 ! Changed naming of special face zones to be more consistent
2422 !
2423 ! Revision 1.9 2005/01/06 04:43:17 haselbac
2424 ! Two bug fixes to allow cell data to be written correctly for mixed grids
2425 !
2426 ! Revision 1.8 2004/12/27 23:34:42 haselbac
2427 ! Added writing of field cell and face data
2428 !
2429 ! Revision 1.7 2004/12/21 15:11:04 fnajjar
2430 ! Increased Tecplot header size to accommodate for PLAG surface statistics
2431 !
2432 ! Revision 1.6 2004/10/19 19:30:18 haselbac
2433 ! Modified RFLU_TEC_WriteZoneSurf, improved estimate of interface vertices
2434 !
2435 ! Revision 1.5 2004/09/27 01:42:54 haselbac
2436 ! Added procedure to write zone with special faces
2437 !
2438 ! Revision 1.4 2004/08/09 16:22:28 haselbac
2439 ! Bug fix: Increased safety margin for alloc
2440 !
2441 ! Revision 1.3 2004/07/06 15:15:09 haselbac
2442 ! Added xyz2vTEC mapping arrays
2443 !
2444 ! Revision 1.2 2004/07/02 03:05:39 haselbac
2445 ! Bug fix: Changed pNull to nullArray
2446 !
2447 ! Revision 1.1 2004/06/16 20:01:40 haselbac
2448 ! Initial revision
2449 !
2450 ! ******************************************************************************
2451 
2452 
2453 
2454 
2455 
2456 
2457 
2458 
2459 
2460 
2461 
2462 
2463 
2464 
subroutine, public rflu_tec_writezonesurf(pRegion, pPatch, faceType, faceKind)
INTEGER function, public rflu_getglobalcellkind(global, pGrid, icg)
subroutine registerfunction(global, funName, fileName)
Definition: ModError.F90:449
subroutine simplifysortedintegers(a, na, nb)
subroutine quicksortinteger(a, n)
subroutine binarysearchinteger(a, n, v, i, j)
subroutine, public rflu_tec_writezonefacesspecial(pRegion)
IndexType nfaces() const
Definition: Mesh.H:641
subroutine, public rflu_tec_writezoneinterf(pRegion)
INTEGER function, public rflu_getfacekind(global, c1k, c2k)
subroutine, public rflu_tec_closefile(global)
subroutine, public rflu_tec_openfile(global, title, fileName)
subroutine, public rflu_tec_writezonesurfmixed(pRegion, pPatch)
subroutine errorstop(global, errorCode, errorLine, addMessage)
Definition: ModError.F90:483
subroutine deregisterfunction(global)
Definition: ModError.F90:469
subroutine, public rflu_tec_writezonecellsspecial(pRegion)
IndexType nvert() const
Definition: Mesh.H:565
subroutine, public rflu_tec_writezonevol(pRegion, cellType, cellKind)