Skip navigation links
(Wind-US Documentation Home Page) (Wind-US User's Guide) (GMAN User's Guide) (MADCAP User's Guide) (CFPOST User's Guide) (Wind-US Utilities) (Common File User's Guide) (Wind-US Installation Guide) (Wind-US Developer's Reference) (Guidelines Documents)

(Common File Versions) (Basic Concepts) (Reference and Scaling Data) (Subroutine Descriptions) (Status Codes) (Default Library Values and Limits) (Definitions for CFD Applications) (Conversion Factors) (Transferring Common Files Between Computer Systems) (CFCRAY - Convert Version 1 or 2 Files to and from Cray Format) (Using Common Files Directly in PLOT3D)

Definitions for CFD Applications

This section describes the conventions to be used when using a common file for CFD applications. Consistent application of these conventions will assure that data will be accessible to a developing base of pre- and post- processing tools. For a picture representation of the following definitions see the figure representing the common file layout for .cfl files.

Node Identifiers

Zone Node Identifiers

Zone nodes are subnodes of the root node. A zone node is accessed and the header read through the following Fortran code fragment:

   CHARACTER*32 ZNAME
   INTEGER ZONE
   INTEGER STATUS, RSTATE, ZSTATE

   ! Note that the default sizes of IPARZ, FPARZ and CPARZ are used.

   INTEGER IPARZ(64)
   REAL    FPARZ(64)
   CHARACTER*80 CPARZ(2)

   ! RSTATE is the state of the root node returned by CFOPEN
   ! All zone nodes under the root node are named 'ZONE nnn'

   ZONE = zone number (1, 2, 3, etc.)
   WRITE (ZNAME, '(''ZONE '', I3)') ZONE
   CALL CFSNOD (STATUS, RSTATE, ZNAME, ZSTATE)
   CALL CFRNOD (STATUS, ZSTATE, IPARZ, FPARZ, CPARZ)
After executing the above code, you may access the data for the zone through the state variable ZSTATE.

Boundary Node Identifiers

Boundary nodes are subnodes of a zone node. A boundary node is accessed through the following Fortran construct:

   CHARACTER*32 BNAME
   INTEGER BNDRY
   INTEGER STATUS, BSTATE, ZSTATE

   ! All boundary nodes for a given zone are stored under the zone
   ! node with a name of 'BNDRY n' where n is a number defined below.

   BNDRY = boundary (1=I1, 2=IMAX, 3=J1, 4=JMAX, 5=K1, 6=KMAX, 7=Chimera, 
                     >= 8=Unstructured)
   WRITE (BNAME, '(''BNDRY  '', I3)') BNDRY
   CALL CFSNOD (STATUS, ZSTATE, BNAME, BSTATE)

Interior Node Identifiers

Interior nodes are subnodes of a zone node (for 3-D unstructured grids only). An interior node is accessed through the following Fortran construct:

   CHARACTER*32 INAME
   INTEGER STATUS, ISTATE, ZSTATE

   ! An interior node for a given zone is stored under the zone
   ! node with the name 'INTERIOR'.

   INAME = 'INTERIOR'
   CALL CFSNOD (STATUS, ZSTATE, INAME, ISTATE)

Node Header Definitions

Note: All dimensional values stored in the common sections of the node headers must be in SI units. Data stored in the application-specific sections can be in any units the application desires.

Root Node Common Definitions

For unstructured grids, IPAR(1), IPAR(2), and IPAR(3) are not used. For a file containing both structured and unstructured zones, these parameters will have the correct values for the structured zones.

    Parameter   Definition

IPAR(1)   Maximum I in all zones
IPAR(2) Maximum J in all zones
IPAR(3) Maximum K in all zones
IPAR(4) Number of zones
IPAR(5) Maximum number of points in any one zone
IPAR(6) Symmetry flag:
0 = no symmetry assumed,
1 = x-y axi-symmetric,
2 = y-z axi-symmetric,
3 = x-z axi-symmetric
IPAR(7) Mxset
IPAR(8) Rotation system flag:
0 = none,
1 = rotating system,
2 = gravity system
IPAR(9-39) Reserved for future use
IPAR(40-64) Reserved for application-specific data

FPAR(1) p0, freestream stagnation pressure
FPAR(2) T0, freestream stagnation temperature
FPAR(3) a0, freestream stagnation speed of sound
FPAR(4) ρ0, freestream stagnation density
FPAR(5) M, freestream Mach number
FPAR(6) p, freestream static pressure
FPAR(7) T, freestream static temperature
FPAR(8) a, freestream speed of sound
FPAR(9) ρ, freestream density
FPAR(10) k, freestream k of k-ε or k-ω turbulence models
FPAR(11) ε or ω, freestream ε or ω of k-ε or k-ω turbulence model
FPAR(12) μ, freestream viscosity
FPAR(13) Re, Reynolds number
FPAR(14) α, angle of attack
FPAR(15) β, yaw angle
FPAR(16-23) Reserved for future use
FPAR(24) β, effective γ
FPAR(25) R, gas constant
FPAR(26) γ, ratio of specific heats (1.4)
FPAR(27) Pr, Prandtl number (0.72)
FPAR(28) Prt, turbulent Prandtl number (0.9)
FPAR(29) x, x base point of axi-symmetric line
FPAR(30) y, y base point of axi-symmetric line
FPAR(31) z, z base point of axi-symmetric line
FPAR(32) M, slope of axi-symmetric line
FPAR(33) A, angle of rotation about axi-symmetric line
FPAR(34-36) Rotation system xyz center (IPAR(8) = 1)
FPAR(37-39) Rotation system xyz rotation rate (IPAR(8) = 1), or gravity system xyz terms (IPAR(8) = 2)
FPAR(40-64) Reserved for application-specific data

CPAR(1) Grid title
CPAR(2) Flowfield title

Wind-US Root Node Application Specific Data

    Parameter   Definition

IPAR(59)   Last time level completed (Global Newton)
IPAR(60) −1 for new format global header
IPAR(61) Last I plane completed (marching)
IPAR(62) Last zone completed
IPAR(63) Number of k-ε cycles
IPAR(64) Number of cycles

FPAR(60) Global Newton big norm
FPAR(61) Global Newton L2 norm
FPAR(62) Global Newton max convergence variable ever

Zone Node Header Common Definitions

Zone nodes can contain either structured or unstructured grids. These two types are distinguished based on IPAR(9) in the Zone Node Header. IPAR(9) will be 0 for structured grids and 1 or 2 for unstructured grids.

Structured Grid
    Parameter   Definition

IPAR(1)   I dimension
IPAR(2) J dimension
IPAR(3) K dimension
IPAR(4) Number of fringe points
IPAR(5) Number of overlapping tracking definitions
IPAR(6) Size of overlapping definition work array
IPAR(7) Reserved for future use
IPAR(8) Rotation flag: 0 = Non-rotating; 1 = Rotating (see FPAR(11-16))
IPAR(9) Grid type: 0 = Structured; 1,2 = Unstructured
IPAR(10) Gas model:
0 = Ideal gas
1 = Thermally perfect (frozen chemistry)
2 = Equilibrium air
3 = Finite rate
IPAR(11) Turbulence model:
0 = Euler (Inviscid)
1 = Laminar
2 = Baldwin-Lomax
3 = Cebeci-Smith
4 = k-ε (obsolete)
5 = Baldwin-Lomax and PDT
6 = Baldwin-Barth (1 equation)
7 = Spalart-Allmaras (1 equation)
8 = SST-Menter (2 equation, k-ω)
10 = Chien k-ε
IPAR(12) Cell/node/variable relationships
0 = Variable values are at node points
1 = Variable values are at cell centers
IPAR(13) Wall function mode
0 = No wall function
1 = White/Christoph model
IPAR(20) Boundary type for I = 1 boundary
0 = Undefined
1 = Reflection/symmetry
2 = Adiabatic wall
3 = Freestream
4 = Viscous wall
5 = Unused
6 = Unused
7 = Arbitrary inflow
8 = Outflow
9 = Inviscid wall
10 = Self-closing
11 = Singular axis
12 = Inviscid axis and wall (not used in Wind-US)
13 = Coupled/point by point
14 = Unused
15 = Bleed
16 = Pinwheel axis
17 = Frozen
18 = Chimera
IPAR(21) Boundary type for I = IMAX boundary
IPAR(22) Boundary type for J = 1 boundary
IPAR(23) Boundary type for J = JMAX boundary
IPAR(24) Boundary type for K = 1 boundary
IPAR(25) Boundary type for K = KMAX boundary
IPAR(26) Grid velocity flag; 0 = none, 1 = global, 13 = point by point
IPAR(27-39) Reserved for future use
IPAR(40-64) Reserved for application-specific data

FPAR(1-7) Zone min/max (xmin, xmax, ymin, ymax, zmin, zmax, checksum)
FPAR(8-10) Global translation velocity (u, v, w)
FPAR(11-13) Global rotation velocity (vr, vθ, vψ)
FPAR(14-16) Global rotation center (x, y, z)
FPAR(17-39) Reserved for future use
FPAR(40-64) Reserved for application-specific data

CPAR(1) Zone title
CPAR(2) Characters 1:40, grid generator identification
Characters 41:80, flow solver identification
Unstructured Grid
    Parameter   Definition

IPAR(1)   Number of points
IPAR(2) Total number of edges, including internal as well as boundary edges. (Optional; only needed if edge data structure is stored.)
IPAR(3) Total number of faces, including internal as well as boundary faces.
IPAR(4) Number of cells (0 --> surface)
IPAR(5) Max number of nodes per face
IPAR(6) Max number of faces per cell
IPAR(7) Number of surface offset data entries
IPAR(8) Max number of nodes per cell
IPAR(9) Grid type:
0 = Structured
1 = Unstructured (tetrahedral)
2 = Hybrid (mixed element types)
IPAR(10) Gas model:
0 = Ideal gas
1 = Thermally perfect (frozen chemistry)
2 = Equilibrium air
3 = Finite rate
IPAR(11) Turbulence model:
0 = Euler (Inviscid)
1 = Laminar
2 = Baldwin-Lomax
3 = Cebeci-Smith
4 = k-ε
5 = Baldwin-Lomax and PDT
IPAR(12) Cell/node/variable relationships:
0 = Variable values are at node points
1 = Variable values are at cell centers
IPAR(13) Number of fringe points
IPAR(14) Number of sequenced cells
IPAR(15-18) Reserved for flow code use
IPAR(19) Total number of elements (CGNS)
IPAR(20) Number of boundary points
IPAR(21) Number of boundary edges
IPAR(22) Number of boundary faces
IPAR(23) Boundary condition type (also used to distinguish between 2-D and 3-D grids):
0 = Face (3D)
1 = Node (3D)
2 = Edge (2D)
3 = Node (2D)
IPAR(24) Number of closed boundaries (Curves-2D, Surfaces-3D) up to a maximum of 20
IPAR(25-58) Array of boundary size pairs:
(25,27,...) - No. of points (2D or 3D surface)
(26,28,...) - No. of faces (Full 3D boundary)
IPAR(60) Number of overlapping tracking definitions
IPAR(61) Size of overlapping definition work array

FPAR(1-39) Reserved for future use
FPAR(40-64) Reserved for application-specific data

CPAR(1) Zone title
CPAR(2) Characters 1:40, grid generator identification
Characters 41:80, flow solver identification

Wind-US Zonal Node Application Specific Data

    Parameter   Definition

IPAR(55)   Compressor face mode
IPAR(56) I location of downstream pressure
IPAR(57) J location of downstream pressure
IPAR(58) K location of downstream pressure
IPAR(59) Downstream pressure variable flag (I, J, or K)
IPAR(63) Unstructured grid type:
0 = Mixed cells
1 = Tets only
5 = Pyramid only
11 = Prism only
111 = Hex only

FPAR(37-39) .dat file SMO1-3 smoothing parameters
FPAR(40-42) Current SMO1-3 levels in solution file
FPAR(43-54) Load convergence FPxyz, FVxyz, MPxyz, MVxyz
FPAR(55) Compressor face T0
FPAR(56) Compressor face mass flow
FPAR(57) Compressor face Mach
FPAR(58) Downstream pressure for IJK mode
FPAR(59) Downstream Mach for IJK mode
FPAR(60)
FPAR(61) Downstream pressure
FPAR(62) Mass flow ratio
FPAR(63)
FPAR(64) Capture area used with specified mass flow ratio and bleed rate boundary conditions

Boundary Node Application Specific Data

    Parameter   Definition

IPAR(1)   BC format: 0 = old, 1 = new
IPAR(8) Unstructured offset (CGNS)
IPAR(20) BC code (CGNS)
IPAR(40) Wall function mode
0 = No wall function
1 = White/Christoph model
IPAR(58) Boundary rotation mode
IPAR(59-64) Rotation specific data depending on mode

FPAR(40) Capture area (see zonal FPAR(64))
FPAR(41) Bleed area (see BLAREA variable)

CPAR(1) Capture area region name
CPAR(2) Bleed area region name

Variable Identifiers

Geometry

    Parameter   Definition

x   x, x coordinate
y y, y coordinate
z z, z coordinate
IBLANK Blanking data asscociated with overlapping grids
ui x direction grid velocity
vi y direction grid velocity
wi z direction grid velocity

The variables x, y, and z are located under the zone node for structured grids. For unstructured grids, x, y, and z for the boundary points are located under the zone node and x, y, and z for the interior points (3D) are located under the interior node. Note that IBLANK has labels encoded in it as follows: 1 = a normal point, −(label + 1) = a hole point, and +(label + 1) = a fringe point, where label is a positive number.

The following variables are applicable to unstructured grids only, and are located under the zone node for boundary surfaces and under the interior node for volumes. These variables describe the connectivity of the faces and/or edges.

    Variable Definition
facpi   Index for point i of a face, i = 1 to the number of points per face. facp4 = 0 for triangular faces in a mixed quad/tri face grid. The facp values for a particular face should be oriented such that the face normal points into the domain.
edgpi Index for point i of an edge segment
facei Index for edge i of a face, i = 1 to the number of edges per face

The edgpi and facei variable definitions are provided for future extensions and will not be initially supported by common applications (such as post-processing).

The following variables are applicable to unstructured grids only, and are located under the the interior node. These variables describe the connectivity of the cells.

    Variable Definition
celpi   Index for point i of a cell, i = 1 to the number of points per cell. celpi = 0 for nodes of dimensionality greater than the current cell. (I.e., celp5 = 0 for a tetrahedral cell.)
celfi Index for face i of a cell, i = 1 to the number of faces per cell
celtyp The type of cell (see definition of IPAR(63) under "Wind-US Zonal Node Application Specific Data"). Used for hybrid unstructured grids only.

The celfi variable definitions are provided for future extensions and will not be initially supported by common applications (such as post-processing)

The following variables are applicable to unstructured grids only, and are located under the zone node. These variables describe collections of faces or surfaces. Currently, each surface must contain only a single face type (i.e., all triangles or all quads). Each of these variables are of size nsurfs = ipar(7).

    Variable Definition
srfoff   Starting index in the face list of a surface segment
srfsiz Number of faces on a surface segment
srfid ID number for a surface segment
srfbc Boundary condition for a surface segment, BC codes same as structured grid
srfvrt For hybrid grids, the number of points per face for all faces of that boundary surface
srfpts Number of points referenced by all of the faces on that boundary surface

General Flow Variables

    Variable Definition
rho   ρ, density
rho*u ρu, x component of momentum per unit volume
rho*v ρv, y component of momentum per unit volume
rho*w ρw, z component of momentum per unit volume
rho*e0 ρe0, stagnation energy per unit volume
P p, pressure
T T, temperature
u u, x component of velocity
v v, y component of velocity
w w, z component of velocity
e0 e0, stagnation energy per unit mass
M M, Mach number
s s, entropy
h h, enthalpy
omegax Ωx, x component of vorticity
omegay Ωy, y component of vorticity
omegaz Ωz, z component of vorticity

For rotating systems (IPAR(8) = 1), an "r" is appended to the name of the following variables to indicate the variable is in a rotating reference system: rho*u, rho*v, rho*w, rho*e0, u, v, w, e0, M.

Turbulence Variables

    Variable Definition
mul   μl, laminar viscosity
mut μt, turbulent viscosity
k k, kinetic energy (k-ε model)
rho*k ρk (k-ε model)
epsilon ε, rate of dissipation (k-ε model)
rho*epsi ρε (k-ε model)
K k (k-ε and SST models)
omega ω (SST model)
anut νt (Baldwin-Barth and Spalart-Allmaras models)

Chemistry Variables

    Variable Definition
a   a, local speed of sound
beta β, effective gamma
Z Z, compressibility
kappa κ, thermal conductivity
etake νK.E., kinetic energy efficiency
PHI Mass fraction of non-reacting species
H Mass fraction of H
N Mass fraction of N
O Mass fraction of O
H2 Mass fraction of H2
N2 Mass fraction of N2
O2 Mass fraction of O2
OH Mass fraction of OH
NO Mass fraction of NO
H2O Mass fraction of H2O
CO2 Mass fraction of CO2
rho*H ρH, ρ times the mass fraction of H; similarly for N, O, etc.

Miscellaneous Variables

    Variable Definition
Cp   Cp, pressure coefficient
delta* δ*, displacement thickness of boundary layer
THETA θ, momentum thickness of boundary layer
Redelta* Reδ*, Reynolds number based on δ*
Cf1 Cf1, skin friction in I direction
Cf2 Cf2, skin friction in J direction
Cf3 Cf3, skin friction in K direction

Wind-US Application Specific Variables

The following variables are used by Wind-US. For the new BC format, the "_n" variables are for node-centered grids, and the "_c" variables are for cell-centered grids. For cell vertex grids the size of the boundary is the number of points, and for cell-centered grids it is the number of faces.

Variable   Description   Location   Type

maxr Max residual for all zones for each equation Root Real
BLAREA Bleed area Root Real
SETDATA(3,*) For each set; (1,*) = # iterations, (2,*) = max residual, and 3 = integrated time (located under zone node) Zone Real
LABLST Label for overlapping tracking Zone Integer
LABTYP Hole/fringe type for overlapping tracking Zone Integer
LABIND Index into LABDEF for overlapping tracking Zone Integer
LABDEF Generation data for overlapping tracking Zone Real
NZN If > 0, zone coupled to; if ≤ 0, boundary condition Boundary Integer
NBD Boundary coupled to Boundary Integer
IFRG The I value of the fringe point (FRNGBND) Boundary Integer
JFRG The J value of the fringe point (FRNGBND) Boundary Integer
KFRG The K value of the fringe point (FRNGBND) Boundary Integer
I1 Node coupled to (first coordinate) Boundary Integer
I2 Node coupled to (second coordinate) Boundary Integer
I3 Node coupled to (third coordinate) Boundary Integer
F1 First coordinate tri/bilinear interpolating factor Boundary Real
F2 Second coordinate tri/bilinear interpolating factor Boundary Real
F3 Third coordinate trilinear interpolating factor Boundary Real
Zone_[nc] If > 0, zone coupled to; if ≤ 0, boundary condition Boundary Integer
BndryNo_[nc] The number of the boundary coupled to Boundary Integer
Cell_n,Cell_c The cell number containing the coupled point/face. For structured grids it is the IJK product for the point. Boundary Integer
CCell_n,CCell_c The cell number in the coupled zone containing the coupled point/face. Boundary Integer
CellLoc_n1,2,3 For structured grids, the IJK location plus the bilinear factor in the coupled-to cell. For unstructured grids, the 1,2,3 weighting factor of the 1,2,3 boundary face point. For quads the fourth weighting factor is 1.0 − ∑(1,2,3). Boundary Double
CellLoc_c1,2,3 For structured grids, the IJK location plus the bilinear factor in the coupled-to cell. For unstructured grids, the X,Y,Z value of the cell center of the coupled face. Boundary Double
CFace_c Face used to parameterize XYZ in CellLoc_c Boundary Integer
Trans Turbulent transition specification array Boundary Real
Temp Temperature specification array, K Boundary Real

The variables Trans and Temp are stored in boundary nodes in the .cfl file since they are solution-specific data.

For hybrid grids we have the following. When an unstructured cell-centered grid is coupled to a structured vertex-centered grid, the CellLoc_c data contains the structured bilinear factors used by the structured zone to interpolate the data to the unstructured face-centered location. For the other way, CellLoc_n contains the parameterized value of the structured point being coupled. The parameterization is based on the unstructured zone's face that is coupled to. This face number is stored in the CFace_n array. This is done since this face is not unique. Note this is not a problem for cell-centered unstructured-to-unstructured, since the parameterization is based on the unstructured face being coupled, and so there is no ambiguity. The follow code segments define the parameterization of the XYZ values into CellLoc_x.

The following takes the parameterization back to XYZ:

   !-----  Copy out the parameterization
           e123(1:3) = paramxyz(n,1:3)

   !-----  Calculate the coordinate vectors
           if ( VrtcsPerBndryFace >= 4 .and. facep(4,n) > 0 ) then
   !-----    Quad face
             v123(1,1:3) = VrtxXyz(facep(3,n),1:3) - VrtxXyz(facep(1,n),1:3)
             v123(2,1:3) = VrtxXyz(facep(4,n),1:3) - VrtxXyz(facep(2,n),1:3)
           else
   !-----    Tri face
             v123(1,1:3) = VrtxXyz(facep(2,n),1:3) - VrtxXyz(facep(1,n),1:3)
             v123(2,1:3) = VrtxXyz(facep(3,n),1:3) - VrtxXyz(facep(1,n),1:3)
           end if
           norm(1)   = (v123(1,2)*v123(2,3) - v123(2,2)*v123(1,3))
           norm(2)   = (v123(1,3)*v123(2,1) - v123(2,3)*v123(1,1))
           norm(3)   = (v123(1,1)*v123(2,2) - v123(2,1)*v123(1,2))
           v123(3,:) = norm / sqrt(sum(norm**2)) * sqrt(sum(v123(1,:)**2))

   !-----  Calculate the xyz cell center
           xyz(n,1) = VrtxXyz(facep(1,n),1) + sum(e123*v123(:,1))
           xyz(n,2) = VrtxXyz(facep(1,n),2) + sum(e123*v123(:,2))
           xyz(n,3) = VrtxXyz(facep(1,n),3) + sum(e123*v123(:,3))

The following parameterizes XYZ based on the input face:

   !-----  Setup the equations
           if ( VrtcsPerBndryFace >= 4 .and. facep(4) > 0 ) then
   !-----    Quad face
             v123(1,1:3) = VrtxXyz(1:3,facep(3)) - VrtxXyz(1:3,facep(1))
             v123(2,1:3) = VrtxXyz(1:3,facep(4)) - VrtxXyz(1:3,facep(2))
           else
   !-----    Tri face
             v123(1,1:3) = VrtxXyz(1:3,facep(2)) - VrtxXyz(1:3,facep(1))
             v123(2,1:3) = VrtxXyz(1:3,facep(3)) - VrtxXyz(1:3,facep(1))
           end if
           norm(1)   = (v123(1,2)*v123(2,3) - v123(2,2)*v123(1,3))
           norm(2)   = (v123(1,3)*v123(2,1) - v123(2,3)*v123(1,1))
           norm(3)   = (v123(1,1)*v123(2,2) - v123(2,1)*v123(1,2))
           v123(3,:) = norm / sqrt(sum(norm**2)) * sqrt(sum(v123(1,:)**2))
           b123 = xyzc - xyz(:,facep(1))

   !-----  Use Cramer's rule to solve the equations Ax = B
           call bg_ggutil_linsolve3x3_f ( status, v123, b123, e123 )
           if (status /= 0 ) e123 = 0

CFD Application Figures

Common File Layout for .cfl Files

Figure showing common file layout for .cfl files

Common File Layout for .cgd Files (Old BC Format)

Figure showing common file layout for .cgd files


Last updated 14 July 2005