8.1
general documentation
Mesh Fortran structure, pointers to the C structure

Functions/Subroutines

elemental pure integer function ifacel (iside, ifac)
  Index-numbers of the two (only) neighboring cells for each internal face More...
 
elemental pure integer function ifabor (ifac)
  index-number of the (unique) neighboring cell for each boundary face More...
 
elemental pure integer function isolid (iporos, iel)
  integer to mark out the "solid" cells (where the fluid volume is 0). Only available when iposros > 0. When iporos = 0, this array has a unique value (isolid_0(1:1)=0). More...
 

Variables

integer ndim
  spatial dimension (3) More...
 
integer, save ncelet = 0
  number of extended (real + ghost of the 'halo') cells. See Note 1: ghost cells - (halos) More...
 
integer, save ncel = 0
  number of real cells in the mesh More...
 
integer, save nfac = 0
  number of internal faces (see Note 2: internal faces) More...
 
integer, save nfabor = 0
  number of boundary faces (see Note 2: internal faces) More...
 
integer, dimension(:,:), pointer ifacel_0
 
integer, dimension(:), pointer ifabor_0
 
double precision, dimension(:,:), pointer xyzcen
  coordinate of the cell centers More...
 
double precision, dimension(:,:), pointer surfac
  surface vector of the internal faces. Its norm is the surface of the face and it is oriented from ifacel(1,.) to ifacel(2,.) More...
 
double precision, dimension(:,:), pointer surfbo
  surface vector of the boundary faces. Its norm is the surface of the face and it is oriented outwards More...
 
double precision, dimension(:,:), pointer suffac
  fluid surface vector of the internal faces. Its norm is the surface of the face and it is oriented from ifacel(1,.) to ifacel(2,.) More...
 
double precision, dimension(:,:), pointer suffbo
  fluid surface vector of the boundary faces. Its norm is the surface of the face and it is oriented outwards More...
 
integer, dimension(:), pointer isolid_0
 
double precision, dimension(:,:), pointer cdgfac
  coordinates of the centers of the internal faces More...
 
double precision, dimension(:,:), pointer cdgfbo
  coordinates of the centers of the boundary faces More...
 
double precision, dimension(:), pointer volume
  volume of each cell More...
 
double precision, dimension(:), pointer cell_f_vol
  fluid volume of each cell More...
 
double precision, dimension(:), pointer surfan
  norm of the surface vector of the internal faces More...
 
double precision, dimension(:), pointer surfbn
  norm of the surface of the boundary faces More...
 
double precision, dimension(:), pointer suffan
  norm of the fluid surface vector of the internal faces More...
 
double precision, dimension(:), pointer suffbn
  norm of the fluid surface of the boundary faces More...
 
double precision, dimension(:), pointer dist
  for every internal face, dot product of the vectors $ \vect{IJ}$ and $\vect{n}$. I and J are respectively the centers of the first and the second neighboring cell. The vector $\vect{n}$ is the unit vector normal to the face and oriented from the first to the second cell More...
 
double precision, dimension(:), pointer distb
  For every boundary face, dot product between the vectors $\vect{IF}$ and $\vect{n}$. I is the center of the neighboring cell. F is the face center. The vector $\vect{n}$ is the unit vector normal to the face and oriented to the exterior of the domain More...
 
double precision, dimension(:), pointer pond
  weighting (Aij=pond Ai+(1-pond)Aj) for every internal face, $\displaystyle\frac{\vect{FJ}.\vect{n}}{\vect{IJ}.\vect{n}}$. With regard to the mesh quality, its ideal value is 0.5 More...
 
double precision, dimension(:,:), pointer dijpf
  vector I'J' for interior faces for every internal face, the three components of the vector $\vect{I'J'}$, where I' and J' are respectively the orthogonal projections of the neighboring cell centers I and J on a straight line orthogonal to the face and passing through its center More...
 
double precision, dimension(:,:), pointer diipb
  vector II' for interior faces for every boundary face, the three components of the vector $\vect{II'}$. I' is the orthogonal projection of I, center of the neighboring cell, on the straight line perpendicular to the face and passing through its center More...
 
double precision, dimension(:,:), pointer dofij
  vector OF for interior faces for every internal face, the three components of the vector $\vect{OF}$. O is the intersection point between the face and the straight line joining the centers of the two neighboring cells. F is the face center More...
 

Detailed Description

Function/Subroutine Documentation

◆ ifabor()

elemental pure integer function mesh::ifabor ( integer, intent(in)  ifac)

index-number of the (unique) neighboring cell for each boundary face

◆ ifacel()

elemental pure integer function mesh::ifacel ( integer, intent(in)  iside,
integer, intent(in)  ifac 
)

Index-numbers of the two (only) neighboring cells for each internal face

◆ isolid()

elemental pure integer function mesh::isolid ( integer, intent(in)  iporos,
integer, intent(in)  iel 
)

integer to mark out the "solid" cells (where the fluid volume is 0). Only available when iposros > 0. When iporos = 0, this array has a unique value (isolid_0(1:1)=0).

Variable Documentation

◆ cdgfac

double precision, dimension(:,:), pointer cdgfac

coordinates of the centers of the internal faces

◆ cdgfbo

double precision, dimension(:,:), pointer cdgfbo

coordinates of the centers of the boundary faces

◆ cell_f_vol

double precision, dimension(:), pointer cell_f_vol

fluid volume of each cell

◆ diipb

double precision, dimension(:,:), pointer diipb

vector II' for interior faces for every boundary face, the three components of the vector $\vect{II'}$. I' is the orthogonal projection of I, center of the neighboring cell, on the straight line perpendicular to the face and passing through its center

◆ dijpf

double precision, dimension(:,:), pointer dijpf

vector I'J' for interior faces for every internal face, the three components of the vector $\vect{I'J'}$, where I' and J' are respectively the orthogonal projections of the neighboring cell centers I and J on a straight line orthogonal to the face and passing through its center

◆ dist

double precision, dimension(:), pointer dist

for every internal face, dot product of the vectors $ \vect{IJ}$ and $\vect{n}$. I and J are respectively the centers of the first and the second neighboring cell. The vector $\vect{n}$ is the unit vector normal to the face and oriented from the first to the second cell

◆ distb

double precision, dimension(:), pointer distb

For every boundary face, dot product between the vectors $\vect{IF}$ and $\vect{n}$. I is the center of the neighboring cell. F is the face center. The vector $\vect{n}$ is the unit vector normal to the face and oriented to the exterior of the domain

◆ dofij

double precision, dimension(:,:), pointer dofij

vector OF for interior faces for every internal face, the three components of the vector $\vect{OF}$. O is the intersection point between the face and the straight line joining the centers of the two neighboring cells. F is the face center

◆ ifabor_0

integer, dimension(:), pointer ifabor_0

◆ ifacel_0

integer, dimension(:,:), pointer ifacel_0

◆ isolid_0

integer, dimension(:), pointer isolid_0

◆ ncel

integer, save ncel = 0

number of real cells in the mesh

◆ ncelet

integer, save ncelet = 0

number of extended (real + ghost of the 'halo') cells. See Note 1: ghost cells - (halos)

◆ ndim

integer ndim

spatial dimension (3)

◆ nfabor

integer, save nfabor = 0

number of boundary faces (see Note 2: internal faces)

◆ nfac

integer, save nfac = 0

number of internal faces (see Note 2: internal faces)

◆ pond

double precision, dimension(:), pointer pond

weighting (Aij=pond Ai+(1-pond)Aj) for every internal face, $\displaystyle\frac{\vect{FJ}.\vect{n}}{\vect{IJ}.\vect{n}}$. With regard to the mesh quality, its ideal value is 0.5

◆ suffac

double precision, dimension(:,:), pointer suffac

fluid surface vector of the internal faces. Its norm is the surface of the face and it is oriented from ifacel(1,.) to ifacel(2,.)

◆ suffan

double precision, dimension(:), pointer suffan

norm of the fluid surface vector of the internal faces

◆ suffbn

double precision, dimension(:), pointer suffbn

norm of the fluid surface of the boundary faces

◆ suffbo

double precision, dimension(:,:), pointer suffbo

fluid surface vector of the boundary faces. Its norm is the surface of the face and it is oriented outwards

◆ surfac

double precision, dimension(:,:), pointer surfac

surface vector of the internal faces. Its norm is the surface of the face and it is oriented from ifacel(1,.) to ifacel(2,.)

◆ surfan

double precision, dimension(:), pointer surfan

norm of the surface vector of the internal faces

◆ surfbn

double precision, dimension(:), pointer surfbn

norm of the surface of the boundary faces

◆ surfbo

double precision, dimension(:,:), pointer surfbo

surface vector of the boundary faces. Its norm is the surface of the face and it is oriented outwards

◆ volume

double precision, dimension(:), pointer volume

volume of each cell

◆ xyzcen

double precision, dimension(:,:), pointer xyzcen

coordinate of the cell centers