grid_unstructured_cell.F90       coverage:  85.71 %func     73.27 %block


     1) module Grid_Unstructured_Cell_module
     2)  
     3)   use PFLOTRAN_Constants_module
     4) 
     5)   implicit none
     6) 
     7)   private
     8) 
     9) #include "petsc/finclude/petscsys.h"
    10) 
    11)   ! grid cell type
    12)   PetscInt, parameter, public :: HEX_TYPE          = 1
    13)   PetscInt, parameter, public :: TET_TYPE          = 2
    14)   PetscInt, parameter, public :: WEDGE_TYPE        = 3
    15)   PetscInt, parameter, public :: PYR_TYPE          = 4
    16)   PetscInt, parameter, public :: POLY_TYPE         = 8 ! move LINE_TYPE from gauss.F90
    17)   ! 2D cell types:
    18)   PetscInt, parameter, public :: TRI_TYPE          = 5
    19)   PetscInt, parameter, public :: QUAD_TYPE         = 6
    20) 
    21)   ! grid cell properties
    22)   PetscInt, parameter, public :: LINE_FACE_TYPE    = 1
    23)   PetscInt, parameter, public :: TRI_FACE_TYPE     = 2
    24)   PetscInt, parameter, public :: QUAD_FACE_TYPE    = 3
    25)   PetscInt, parameter, public :: MAX_VERT_PER_FACE = 4
    26)   PetscInt, parameter, public :: MAX_FACE_PER_CELL = 6
    27)   PetscInt, parameter, public :: MAX_FACE_PER_CELL_SURF = 4
    28) 
    29)   type, public :: point_type
    30)     PetscInt :: id
    31)     PetscReal :: x
    32)     PetscReal :: y
    33)     PetscReal :: z
    34)   end type point_type
    35)   
    36)   type, public :: plane_type
    37)     PetscReal :: A
    38)     PetscReal :: B
    39)     PetscReal :: C
    40)     PetscReal :: D
    41)   end type plane_type
    42)   
    43)   public :: UCellComputeCentroid, &
    44)             UCellComputeVolume, &
    45)             UCellComputeArea, &
    46)             UCellComputePlane, &
    47)             UCellGetPlaneIntercept, &
    48)             UCellProjectPointOntoPlane, &
    49)             UCellComputeDistanceFromPlane, &
    50)             UCellGetLineIntercept, &
    51)             UCellGetNVertices, &
    52)             UCellGetNFaces, &
    53)             UCellGetNFaceVertices, &
    54)             UCellGetFaceType, &
    55)             UCellGetFaceVertices, &
    56)             UCellGetNFaceVertsandVerts, &
    57)             UCellTypeToWord, &
    58)             UCellFaceTypeToWord, &
    59)             UCellQuality
    60)             
    61) contains
    62) 
    63) ! ************************************************************************** !
    64) 
    65) function UCellComputeCentroid(cell_type,vertices,option)
    66)   ! 
    67)   ! Computes the centroid a grid cell
    68)   ! 
    69)   ! Author: Glenn Hammond
    70)   ! Date: 10/30/09
    71)   ! 
    72) 
    73)   use Option_module
    74) 
    75)   implicit none
    76)   
    77)   PetscInt :: cell_type
    78)   type(point_type) :: vertices(*)
    79)   type(option_type) :: option
    80)   
    81)   PetscReal :: UCellComputeCentroid(3)
    82)   PetscInt :: ivertex
    83)   
    84)   UCellComputeCentroid = 0.d0
    85)   select case(cell_type)
    86)     case(HEX_TYPE)
    87)       ! need something more sophisticated, but for now, just use average
    88)       do ivertex = 1, 8
    89)         UCellComputeCentroid(1) = UCellComputeCentroid(1) + vertices(ivertex)%x
    90)         UCellComputeCentroid(2) = UCellComputeCentroid(2) + vertices(ivertex)%y
    91)         UCellComputeCentroid(3) = UCellComputeCentroid(3) + vertices(ivertex)%z
    92)       enddo
    93)       UCellComputeCentroid = UCellComputeCentroid / 8.d0
    94)     case(WEDGE_TYPE)
    95)       ! need something more sophisticated, but for now, just use average
    96)       do ivertex = 1, 6
    97)         UCellComputeCentroid(1) = UCellComputeCentroid(1) + vertices(ivertex)%x
    98)         UCellComputeCentroid(2) = UCellComputeCentroid(2) + vertices(ivertex)%y
    99)         UCellComputeCentroid(3) = UCellComputeCentroid(3) + vertices(ivertex)%z
   100)       enddo
   101)       UCellComputeCentroid = UCellComputeCentroid / 6.d0
   102)     case(PYR_TYPE)
   103)       ! need something more sophisticated, but for now, just use average
   104)       do ivertex = 1, 5
   105)         UCellComputeCentroid(1) = UCellComputeCentroid(1) + vertices(ivertex)%x
   106)         UCellComputeCentroid(2) = UCellComputeCentroid(2) + vertices(ivertex)%y
   107)         UCellComputeCentroid(3) = UCellComputeCentroid(3) + vertices(ivertex)%z
   108)       enddo
   109)       UCellComputeCentroid = UCellComputeCentroid / 5.d0
   110)     case(TET_TYPE, QUAD_TYPE)
   111)       do ivertex = 1, 4
   112)         UCellComputeCentroid(1) = UCellComputeCentroid(1) + vertices(ivertex)%x
   113)         UCellComputeCentroid(2) = UCellComputeCentroid(2) + vertices(ivertex)%y
   114)         UCellComputeCentroid(3) = UCellComputeCentroid(3) + vertices(ivertex)%z
   115)       enddo
   116)       UCellComputeCentroid = UCellComputeCentroid / 4.d0
   117)     case(TRI_TYPE)
   118)       do ivertex = 1, 3
   119)         UCellComputeCentroid(1) = UCellComputeCentroid(1) + vertices(ivertex)%x
   120)         UCellComputeCentroid(2) = UCellComputeCentroid(2) + vertices(ivertex)%y
   121)         UCellComputeCentroid(3) = UCellComputeCentroid(3) + vertices(ivertex)%z
   122)       enddo
   123)       UCellComputeCentroid = UCellComputeCentroid / 3.d0
   124)     case default
   125)       option%io_buffer = 'Cell type not recognized'
   126)       call printErrMsg(option)
   127)   end select
   128) 
   129) end function UCellComputeCentroid
   130) 
   131) ! ************************************************************************** !
   132) 
   133) function UCellComputeVolume(cell_type,vertices,option)
   134)   ! 
   135)   ! Computes the volume a grid cell
   136)   ! 
   137)   ! Author: Glenn Hammond
   138)   ! Date: 11/06/09
   139)   ! 
   140) 
   141)   use Utility_module, only : DotProduct, CrossProduct
   142)   use Option_module
   143) 
   144)   implicit none
   145)   
   146)   PetscInt :: cell_type
   147)   type(point_type) :: vertices(*)
   148)   type(option_type) :: option
   149)   
   150)   PetscReal :: UCellComputeVolume
   151)   PetscReal :: v(3)
   152)   PetscReal :: l1, l2, l3
   153)   PetscReal :: n1(3), n2(3), v1(3), v2(3)
   154)   PetscReal :: area1, area2, dz
   155)   PetscReal :: vv(3,8)
   156)   PetscInt :: i, j
   157)   
   158)   UCellComputeVolume = 0.d0
   159)   select case(cell_type)
   160)     case(HEX_TYPE)
   161)       ! split into 5 tetrahedron
   162)       UCellComputeVolume = &
   163)         UCellComputeVolumeOfTetrahedron(vertices(1),vertices(2), &
   164)                                         vertices(3),vertices(6))
   165)       UCellComputeVolume = &
   166)         UCellComputeVolume + &
   167)         UCellComputeVolumeOfTetrahedron(vertices(1),vertices(6), &
   168)                                         vertices(8),vertices(5))
   169)       UCellComputeVolume = &
   170)         UCellComputeVolume + &
   171)         UCellComputeVolumeOfTetrahedron(vertices(1),vertices(6), &
   172)                                         vertices(3),vertices(8))
   173)       UCellComputeVolume = &
   174)         UCellComputeVolume + &
   175)         UCellComputeVolumeOfTetrahedron(vertices(8),vertices(6), &
   176)                                         vertices(3),vertices(7))
   177)       UCellComputeVolume = &
   178)         UCellComputeVolume + &
   179)         UCellComputeVolumeOfTetrahedron(vertices(1),vertices(3), &
   180)                                         vertices(4),vertices(8))
   181)     case(WEDGE_TYPE)
   182)       ! split into 3 tetrahedrons: v1v2v3v4, v3v4v5v6, v2v3v4v5
   183)       UCellComputeVolume = &
   184)         UCellComputeVolumeOfTetrahedron(vertices(1),vertices(2),vertices(3), &
   185)                                         vertices(4))
   186)       UCellComputeVolume = &
   187)         UCellComputeVolume + &
   188)         UCellComputeVolumeOfTetrahedron(vertices(3),vertices(4),vertices(5), &
   189)                                         vertices(6))
   190)       UCellComputeVolume = &
   191)         UCellComputeVolume + &
   192)         UCellComputeVolumeOfTetrahedron(vertices(2),vertices(3),vertices(4), &
   193)                                         vertices(5))
   194)     case(PYR_TYPE)
   195)       ! split pyramid into two tets and compute
   196)       UCellComputeVolume = &
   197)         UCellComputeVolumeOfTetrahedron(vertices(1),vertices(2), &
   198)                                         vertices(3),vertices(5))
   199)       UCellComputeVolume = &
   200)         UCellComputeVolume + &
   201)         UCellComputeVolumeOfTetrahedron(vertices(3),vertices(4), &
   202)                                         vertices(1),vertices(5))
   203)     case(TET_TYPE)
   204)       UCellComputeVolume = &
   205)         UCellComputeVolumeOfTetrahedron(vertices(1),vertices(2),vertices(3), &
   206)                                         vertices(4))
   207)     case(QUAD_TYPE, TRI_TYPE)
   208)       option%io_buffer = 'Cell type is QUAD or TRI, thus one should call '// &
   209)         'UCellComputeArea instead of UCellComputeVolume'
   210)       call printErrMsg(option)
   211)     case default
   212)       option%io_buffer = 'Cell type not recognized'
   213)       call printErrMsg(option)
   214)   end select
   215) 
   216) end function UCellComputeVolume
   217) 
   218) ! ************************************************************************** !
   219) 
   220) function UCellComputeArea(cell_type,vertices,option)
   221)   ! 
   222)   ! Computes the area a 2D grid cell
   223)   ! 
   224)   ! Author: Gautam Bisht
   225)   ! Date: 03/17/12
   226)   ! 
   227) 
   228)   use Utility_module, only : DotProduct, CrossProduct
   229)   use Option_module
   230) 
   231)   implicit none
   232)   
   233)   PetscInt :: cell_type
   234)   type(point_type) :: vertices(*)
   235)   type(option_type) :: option
   236)   
   237)   PetscReal :: UCellComputeArea
   238)   PetscReal :: v(3)
   239)   PetscReal :: l1, l2, l3
   240)   PetscReal :: n1(3), n2(3), v1(3), v2(3)
   241)   PetscReal :: area1, area2, dz
   242)   PetscReal :: vv(3,8)
   243)   PetscInt :: i, j
   244)   
   245)   UCellComputeArea = 0.d0
   246)   select case(cell_type)
   247)     case(QUAD_TYPE)
   248)       v1(1) = vertices(3)%x-vertices(2)%x
   249)       v1(2) = vertices(3)%y-vertices(2)%y
   250)       v1(3) = vertices(3)%z-vertices(2)%z
   251)       v2(1) = vertices(1)%x-vertices(2)%x
   252)       v2(2) = vertices(1)%y-vertices(2)%y
   253)       v2(3) = vertices(1)%z-vertices(2)%z
   254)       n1 = CrossProduct(v1,v2)
   255)       area1 = 0.5d0*sqrt(DotProduct(n1,n1))
   256)       
   257)       v1(1) = vertices(3)%x-vertices(4)%x
   258)       v1(2) = vertices(3)%y-vertices(4)%y
   259)       v1(3) = vertices(3)%z-vertices(4)%z
   260)       v2(1) = vertices(1)%x-vertices(4)%x
   261)       v2(2) = vertices(1)%y-vertices(4)%y
   262)       v2(3) = vertices(1)%z-vertices(4)%z
   263)       n2 = CrossProduct(v1,v2)
   264)       area2 = 0.5d0*sqrt(DotProduct(n2,n2))
   265)       
   266)       UCellComputeArea = (area1 + area2)
   267)     case(TRI_TYPE)
   268)       v1(1) = vertices(3)%x-vertices(2)%x
   269)       v1(2) = vertices(3)%y-vertices(2)%y
   270)       v1(3) = vertices(3)%z-vertices(2)%z
   271)       v2(1) = vertices(1)%x-vertices(2)%x
   272)       v2(2) = vertices(1)%y-vertices(2)%y
   273)       v2(3) = vertices(1)%z-vertices(2)%z
   274)       n1 = CrossProduct(v1,v2)
   275)       area1 = 0.5d0*sqrt(DotProduct(n1,n1))
   276) 
   277)       UCellComputeArea = area1
   278)     case default
   279)       option%io_buffer = 'Cell type not recognized'
   280)       call printErrMsg(option)
   281)   end select
   282) 
   283) end function UCellComputeArea
   284) 
   285) ! ************************************************************************** !
   286) 
   287) function UCellComputeVolumeOfTetrahedron(point1,point2,point3,point4)
   288)   ! 
   289)   ! Computes the voluem of a tetrahedron
   290)   ! given four points
   291)   ! 
   292)   ! Author: Glenn Hammond
   293)   ! Date: 12/06/11
   294)   ! 
   295) 
   296)   use Utility_module, only : DotProduct, CrossProduct
   297) 
   298)   implicit none
   299)   
   300)   type(point_type) :: point1, point2, point3, point4
   301)   
   302)   PetscReal :: vv(3,4)
   303)   PetscReal :: vv1_minus_vv4(3)
   304)   PetscReal :: vv2_minus_vv4(3)
   305)   PetscReal :: vv3_minus_vv4(3)
   306)   PetscReal :: cross_2_minus_4_X_3_minus_4(3)
   307)   PetscReal :: UCellComputeVolumeOfTetrahedron
   308)   PetscInt :: i
   309) 
   310)   vv(1,1) = point1%x
   311)   vv(2,1) = point1%y
   312)   vv(3,1) = point1%z
   313)   vv(1,2) = point2%x
   314)   vv(2,2) = point2%y
   315)   vv(3,2) = point2%z
   316)   vv(1,3) = point3%x
   317)   vv(2,3) = point3%y
   318)   vv(3,3) = point3%z
   319)   vv(1,4) = point4%x
   320)   vv(2,4) = point4%y
   321)   vv(3,4) = point4%z
   322) 
   323)   ! V = |(a-d).((b-d)x(c-d))| / 6
   324) 
   325)   !geh: Intel Visual Fortran creates temporary arrays and reports warnings
   326)   !     to the screen.  Therefore, I will use temporary variables below
   327) 
   328)   !UCellComputeVolumeOfTetrahedron = dabs(DotProduct(vv(:,1)-vv(:,4), &
   329)   !                                       CrossProduct(vv(:,2)-vv(:,4), &
   330)   !                                                    vv(:,3)-vv(:,4)))) / &
   331)   !                                  6.d0
   332) 
   333)   vv1_minus_vv4 = vv(:,1)-vv(:,4)
   334)   vv2_minus_vv4 = vv(:,2)-vv(:,4)
   335)   vv3_minus_vv4 = vv(:,3)-vv(:,4)
   336)   cross_2_minus_4_X_3_minus_4 = CrossProduct(vv2_minus_vv4, &
   337)                                              vv3_minus_vv4)
   338)   UCellComputeVolumeOfTetrahedron = dabs(DotProduct(vv1_minus_vv4, &
   339)                                          cross_2_minus_4_X_3_minus_4)) / &
   340)                                     6.d0
   341) 
   342) end function UCellComputeVolumeOfTetrahedron
   343) 
   344) ! ************************************************************************** !
   345) 
   346) subroutine UCellComputePlane(plane,point1,point2,point3)
   347)   ! 
   348)   ! Calculates the plane intersected by 3 points
   349)   ! 
   350)   ! Author: Glenn Hammond
   351)   ! Date: 10/30/09
   352)   ! 
   353) 
   354)   implicit none
   355)   
   356)   type(plane_type) :: plane
   357)   type(point_type) :: point1, point2, point3
   358)   
   359)   PetscReal :: x1,y1,z1
   360)   PetscReal :: x2,y2,z2
   361)   PetscReal :: x3,y3,z3
   362)   x1 = point1%x
   363)   y1 = point1%y
   364)   z1 = point1%z
   365)   x2 = point2%x
   366)   y2 = point2%y
   367)   z2 = point2%z
   368)   x3 = point3%x
   369)   y3 = point3%y
   370)   z3 = point3%z
   371)   
   372)   ! this grabbed from python script
   373)   plane%A = y1*(z2-z3)+y2*(z3-z1)+y3*(z1-z2)
   374)   plane%B = z1*(x2-x3)+z2*(x3-x1)+z3*(x1-x2)
   375)   plane%C = x1*(y2-y3)+x2*(y3-y1)+x3*(y1-y2)
   376)   plane%D = -1.*(x1*(y2*z3-y3*z2)+x2*(y3*z1-y1*z3)+x3*(y1*z2-y2*z1))
   377) 
   378) end subroutine UCellComputePlane
   379) 
   380) ! ************************************************************************** !
   381) 
   382) subroutine UCellProjectPointOntoPlane(plane,point,intercept)
   383)   ! 
   384)   ! Calculates the intercept of a point with a plane
   385)   ! 
   386)   ! Author: Glenn Hammond
   387)   ! Date: 11/22/11
   388)   ! 
   389) 
   390)   implicit none
   391)   
   392)   type(plane_type) :: plane
   393)   type(point_type) :: point
   394)   type(point_type) :: intercept
   395) 
   396)   PetscReal :: scalar
   397)   
   398)   ! plane equation:
   399)   !   A*x + B*y + C*z + D = 0
   400) 
   401)   scalar = (plane%A*point%x + plane%B*point%y + plane%C*point%z + plane%D) / &
   402)            (plane%A*plane%A + plane%B*plane%B + plane%C*plane%C)
   403)   
   404)   intercept%x = point%x - plane%A * scalar
   405)   intercept%y = point%y - plane%B * scalar
   406)   intercept%z = point%z - plane%C * scalar
   407) 
   408) end subroutine UCellProjectPointOntoPlane
   409) 
   410) ! ************************************************************************** !
   411) 
   412) subroutine UCellGetPlaneIntercept(plane,point1,point2,intercept)
   413)   ! 
   414)   ! Calculates the intercept of a line with a plane
   415)   ! 
   416)   ! Author: Glenn Hammond
   417)   ! Date: 10/30/09
   418)   ! 
   419) 
   420)   implicit none
   421)   
   422)   type(plane_type) :: plane
   423)   type(point_type) :: point1, point2
   424)   type(point_type) :: intercept
   425) 
   426)   PetscReal :: x1,y1,z1
   427)   PetscReal :: x2,y2,z2
   428)   PetscReal :: u
   429)     
   430)   x1 = point1%x
   431)   y1 = point1%y
   432)   z1 = point1%z
   433)   x2 = point2%x
   434)   y2 = point2%y
   435)   z2 = point2%z
   436)  
   437)  
   438)   u = (plane%A*x1 + plane%B*y1 + plane%C*z1 + plane%D) / &
   439)       (plane%A*(x1-x2) + plane%B*(y1-y2) + plane%C*(z1-z2))
   440) 
   441)   intercept%x = point1%x + u*(point2%x-point1%x)
   442)   intercept%y = point1%y + u*(point2%y-point1%y)
   443)   intercept%z = point1%z + u*(point2%z-point1%z)
   444) 
   445) end subroutine UCellGetPlaneIntercept
   446) 
   447) ! ************************************************************************** !
   448) 
   449) function UCellComputeDistanceFromPlane(plane,point)
   450)   ! 
   451)   ! Calculates the distance of a point from a plane
   452)   ! 
   453)   ! Author: Glenn Hammond
   454)   ! Date: 10/24/11
   455)   ! 
   456) 
   457)   implicit none
   458)   
   459)   type(plane_type) :: plane
   460)   type(point_type) :: point
   461)   PetscReal :: UCellComputeDistanceFromPlane
   462) 
   463)   UCellComputeDistanceFromPlane = &
   464)     (plane%A*point%x + plane%B*point%y + plane%C*point%z + plane%D) / &
   465)     sqrt(plane%A*plane%A+plane%B*plane%B+plane%C*plane%C)
   466)   
   467) end function UCellComputeDistanceFromPlane
   468) 
   469) ! ************************************************************************** !
   470) 
   471) function UCellGetNVertices(cell_type,option)
   472)   ! 
   473)   ! Returns number of vertices in a cell
   474)   ! 
   475)   ! Author: Glenn Hammond
   476)   ! Date: 10/24/11
   477)   ! 
   478) 
   479)   use Option_module
   480)   implicit none
   481)   
   482)   PetscInt :: cell_type
   483)   PetscInt :: UCellGetNVertices
   484)   type(option_type) :: option
   485)   
   486)   select case(cell_type)
   487)     case(HEX_TYPE)
   488)       UCellGetNVertices = 8
   489)     case(WEDGE_TYPE)
   490)       UCellGetNVertices = 6
   491)     case(PYR_TYPE)
   492)       UCellGetNVertices = 5
   493)     case(TET_TYPE, QUAD_TYPE)
   494)       UCellGetNVertices = 4
   495)     case(TRI_TYPE)
   496)       UCellGetNVertices = 3
   497)     case default
   498)       option%io_buffer = 'Cell type not recognized'
   499)       call printErrMsg(option)
   500)   end select  
   501)   
   502) end function UCellGetNVertices
   503) 
   504) ! ************************************************************************** !
   505) 
   506) function UCellGetNEdges(cell_type)
   507)   ! 
   508)   ! Returns number of edges in a cell
   509)   ! 
   510)   ! Author: Glenn Hammond
   511)   ! Date: 01/17/12
   512)   ! 
   513) 
   514)   implicit none
   515)   
   516)   PetscInt :: cell_type
   517)   PetscInt :: UCellGetNEdges
   518)   
   519)   select case(cell_type)
   520)     case(HEX_TYPE)
   521)       UCellGetNEdges = 12
   522)     case(WEDGE_TYPE)
   523)       UCellGetNEdges = 9
   524)     case(PYR_TYPE)
   525)       UCellGetNEdges = 8
   526)     case(TET_TYPE)
   527)       UCellGetNEdges = 6
   528)     case(QUAD_TYPE)
   529)       UCellGetNEdges = 4
   530)     case(TRI_TYPE)
   531)       UCellGetNEdges = 3
   532)   end select  
   533)   
   534) end function UCellGetNEdges
   535) 
   536) ! ************************************************************************** !
   537) 
   538) function UCellGetNFaces(cell_type,option)
   539)   ! 
   540)   ! Returns number of faces in a cell
   541)   ! 
   542)   ! Author: Glenn Hammond
   543)   ! Date: 10/24/11
   544)   ! 
   545) 
   546)   use Option_module
   547)   implicit none
   548)   
   549)   PetscInt :: cell_type
   550)   PetscInt :: UCellGetNFaces
   551)   type(option_type) :: option
   552)   
   553)   select case(cell_type)
   554)     case(HEX_TYPE)
   555)       UCellGetNFaces = 6
   556)     case(WEDGE_TYPE,PYR_TYPE)
   557)       UCellGetNFaces = 5
   558)     case(TET_TYPE)
   559)       UCellGetNFaces = 4
   560)     case(QUAD_TYPE)
   561)       UCellGetNFaces = 4
   562)     case(TRI_TYPE)
   563)       UCellGetNFaces = 3
   564)     case default
   565)       option%io_buffer = 'Cell type not recognized'
   566)       call printErrMsg(option)
   567)   end select  
   568)   
   569) end function UCellGetNFaces
   570) 
   571) ! ************************************************************************** !
   572) 
   573) function UCellGetNFaceVertices(cell_type,iface,option)
   574)   ! 
   575)   ! Returns number of vertices in a cell face
   576)   ! 
   577)   ! Author: Glenn Hammond
   578)   ! Date: 10/24/11
   579)   ! 
   580) 
   581)   use Option_module
   582)   implicit none
   583)   
   584)   PetscInt :: cell_type
   585)   PetscInt :: iface
   586)   PetscInt :: UCellGetNFaceVertices
   587)   type(option_type) :: option
   588)   
   589)   select case(cell_type)
   590)     case(HEX_TYPE)
   591)       UCellGetNFaceVertices = 4
   592)     case(WEDGE_TYPE)
   593)       if (iface > 3) then
   594)         UCellGetNFaceVertices = 3
   595)       else 
   596)         UCellGetNFaceVertices = 4
   597)       endif
   598)     case(PYR_TYPE)
   599)       if (iface > 4) then
   600)         UCellGetNFaceVertices = 4
   601)       else 
   602)         UCellGetNFaceVertices = 3
   603)       endif
   604)     case(TET_TYPE)
   605)       UCellGetNFaceVertices = 3
   606)     case(QUAD_TYPE)
   607)       UCellGetNFaceVertices = 2
   608)     case(TRI_TYPE)
   609)       UCellGetNFaceVertices = 2
   610)     case default
   611)       option%io_buffer = 'Cell type not recognized'
   612)       call printErrMsg(option)
   613)   end select
   614)       
   615) end function UCellGetNFaceVertices
   616) 
   617) ! ************************************************************************** !
   618) 
   619) function UCellGetFaceType(cell_type,iface,option)
   620)   ! 
   621)   ! Returns type of cell face
   622)   ! 
   623)   ! Author: Glenn Hammond
   624)   ! Date: 10/24/11
   625)   ! 
   626) 
   627)   use Option_module
   628)   implicit none
   629)   
   630)   PetscInt :: cell_type
   631)   PetscInt :: iface
   632)   PetscInt :: UCellGetFaceType
   633)   type(option_type) :: option
   634)   
   635)   select case(cell_type)
   636)     case(HEX_TYPE)
   637)       UCellGetFaceType = QUAD_FACE_TYPE
   638)     case(WEDGE_TYPE)
   639)       if (iface > 3) then
   640)         UCellGetFaceType = TRI_FACE_TYPE
   641)       else 
   642)         UCellGetFaceType = QUAD_FACE_TYPE
   643)       endif
   644)     case(PYR_TYPE)
   645)       if (iface > 4) then
   646)         UCellGetFaceType = QUAD_FACE_TYPE
   647)       else 
   648)         UCellGetFaceType = TRI_FACE_TYPE
   649)       endif
   650)     case(TET_TYPE)
   651)       UCellGetFaceType = TRI_FACE_TYPE
   652)     case(QUAD_TYPE)
   653)       UCellGetFaceType = LINE_FACE_TYPE
   654)     case(TRI_TYPE)
   655)       UCellGetFaceType = LINE_FACE_TYPE
   656)     case default
   657)       option%io_buffer = 'Cell type not recognized'
   658)       call printErrMsg(option)
   659)   end select
   660)   
   661) end function UCellGetFaceType
   662) 
   663) ! ************************************************************************** !
   664) 
   665) function UCellTypeToWord(cell_type,option)
   666)   ! 
   667)   ! Returns type of cell as a string
   668)   ! 
   669)   ! Author: Glenn Hammond
   670)   ! Date: 12/09/11
   671)   ! 
   672) 
   673)   use Option_module
   674)   implicit none
   675)   
   676)   PetscInt :: cell_type
   677)   type(option_type) :: option
   678) 
   679)   character(len=MAXWORDLENGTH) :: UCellTypeToWord
   680)   
   681)   select case(cell_type)
   682)     case(HEX_TYPE)
   683)       UCellTypeToWord = 'hexahedron'
   684)     case(WEDGE_TYPE)
   685)       UCellTypeToWord = 'wedge'
   686)     case(PYR_TYPE)
   687)       UCellTypeToWord = 'pyramid'
   688)     case(TET_TYPE)
   689)       UCellTypeToWord = 'tetrahedron'
   690)     case(QUAD_TYPE)
   691)       UCellTypeToWord = 'quadrilateral'
   692)     case(TRI_TYPE)
   693)       UCellTypeToWord = 'triangle'
   694)     case default
   695)       option%io_buffer = 'Cell type not recognized'
   696)       call printErrMsg(option)
   697)   end select
   698)   
   699) end function UCellTypeToWord
   700) 
   701) ! ************************************************************************** !
   702) 
   703) function UCellFaceTypeToWord(face_type,option)
   704)   ! 
   705)   ! Returns type of cell face as a string
   706)   ! 
   707)   ! Author: Glenn Hammond
   708)   ! Date: 12/09/11
   709)   ! 
   710) 
   711)   use Option_module
   712)   implicit none
   713)   
   714)   PetscInt :: face_type
   715)   type(option_type) :: option
   716) 
   717)   character(len=MAXWORDLENGTH) :: UCellFaceTypeToWord
   718)   
   719)   select case(face_type)
   720)     case(TRI_FACE_TYPE)
   721)       UCellFaceTypeToWord = 'triangle'
   722)     case(QUAD_FACE_TYPE)
   723)       UCellFaceTypeToWord = 'quadrilateral'
   724)     case(LINE_FACE_TYPE)
   725)       UCellFaceTypeToWord = 'line'
   726)     case default
   727)       option%io_buffer = 'Face type not recognized'
   728)       call printErrMsg(option)
   729)   end select
   730)   
   731) end function UCellFaceTypeToWord
   732) 
   733) ! ************************************************************************** !
   734) 
   735) subroutine UCellGetNFaceVertsandVerts(option,cell_type,iface,nvertices, &
   736)                                       vertex_ids)
   737)   ! 
   738)   ! returns the numbber of vertices for a face and
   739)   ! the vertices
   740)   ! 
   741)   ! Author: Glenn Hammond
   742)   ! Date: 12/06/11
   743)   ! 
   744)   use Option_module
   745)   
   746)   implicit none
   747)   
   748)   type(option_type) :: option
   749)   PetscInt :: cell_type
   750)   PetscInt :: iface
   751)   PetscInt :: nvertices
   752)   PetscInt :: vertex_ids(*)
   753)   
   754)   nvertices = UCellGetNFaceVertices(cell_type,iface,option)
   755)   call UCellGetFaceVertices(option,cell_type,iface,vertex_ids)
   756) 
   757) end subroutine UCellGetNFaceVertsandVerts
   758) 
   759) ! ************************************************************************** !
   760) 
   761) subroutine UCellGetFaceVertices(option,cell_type,iface,vertex_ids)
   762)   ! 
   763)   ! returns vertex ids of a face
   764)   ! 
   765)   ! Author: Glenn Hammond
   766)   ! Date: 11/24/11
   767)   ! 
   768) 
   769)   use Option_module
   770)   
   771)   implicit none
   772)   
   773)   type(option_type) :: option
   774)   PetscInt :: cell_type
   775)   PetscInt :: iface
   776)   PetscInt :: vertex_ids(*)
   777)   
   778)   select case(cell_type)
   779)     case(HEX_TYPE)
   780)       select case(iface)
   781)         case(1)
   782)           vertex_ids(1) = 1
   783)           vertex_ids(2) = 2
   784)           vertex_ids(3) = 6
   785)           vertex_ids(4) = 5
   786)         case(2)
   787)           vertex_ids(1) = 2
   788)           vertex_ids(2) = 3
   789)           vertex_ids(3) = 7
   790)           vertex_ids(4) = 6
   791)         case(3)
   792)           vertex_ids(1) = 3
   793)           vertex_ids(2) = 4
   794)           vertex_ids(3) = 8
   795)           vertex_ids(4) = 7
   796)         case(4)
   797)           vertex_ids(1) = 4
   798)           vertex_ids(2) = 1
   799)           vertex_ids(3) = 5
   800)           vertex_ids(4) = 8
   801)         case(5)
   802)           vertex_ids(1) = 1
   803)           vertex_ids(2) = 4
   804)           vertex_ids(3) = 3
   805)           vertex_ids(4) = 2
   806)         case(6)
   807)           vertex_ids(1) = 5
   808)           vertex_ids(2) = 6
   809)           vertex_ids(3) = 7
   810)           vertex_ids(4) = 8
   811)       end select
   812)     case(WEDGE_TYPE)
   813)       select case(iface)
   814)         case(1)
   815)           vertex_ids(1) = 1
   816)           vertex_ids(2) = 2
   817)           vertex_ids(3) = 5
   818)           vertex_ids(4) = 4
   819)        case(2)
   820)           vertex_ids(1) = 2
   821)           vertex_ids(2) = 3
   822)           vertex_ids(3) = 6
   823)           vertex_ids(4) = 5
   824)        case(3)
   825)           vertex_ids(1) = 3
   826)           vertex_ids(2) = 1
   827)           vertex_ids(3) = 4
   828)           vertex_ids(4) = 6
   829)        case(4)
   830)           vertex_ids(1) = 1
   831)           vertex_ids(2) = 3
   832)           vertex_ids(3) = 2
   833)        case(5)
   834)           vertex_ids(1) = 4
   835)           vertex_ids(2) = 5
   836)           vertex_ids(3) = 6
   837)        case default
   838)           option%io_buffer='Cell WEDGE_TYPE has only 5 faces'
   839)           call printErrMsg(option)
   840)        end select
   841)     case(PYR_TYPE)
   842)       select case(iface)
   843)         case(1)
   844)           vertex_ids(1) = 1
   845)           vertex_ids(2) = 2
   846)           vertex_ids(3) = 5
   847)        case(2)
   848)           vertex_ids(1) = 2
   849)           vertex_ids(2) = 3
   850)           vertex_ids(3) = 5
   851)        case(3)
   852)           vertex_ids(1) = 3
   853)           vertex_ids(2) = 4
   854)           vertex_ids(3) = 5
   855)        case(4)
   856)           vertex_ids(1) = 4
   857)           vertex_ids(2) = 1
   858)           vertex_ids(3) = 5
   859)        case(5)
   860)           vertex_ids(1) = 1
   861)           vertex_ids(2) = 4
   862)           vertex_ids(3) = 3
   863)           vertex_ids(4) = 2
   864)        case default
   865)           option%io_buffer='Cell PYR_TYPE has only 5 faces'
   866)           call printErrMsg(option)
   867)        end select
   868)     case(TET_TYPE)
   869)       select case(iface)
   870)         case(1)
   871)           vertex_ids(1) = 1
   872)           vertex_ids(2) = 2
   873)           vertex_ids(3) = 4
   874)         case(2)
   875)           vertex_ids(1) = 2
   876)           vertex_ids(2) = 3
   877)           vertex_ids(3) = 4
   878)         case(3)
   879)           vertex_ids(1) = 1
   880)           vertex_ids(2) = 4
   881)           vertex_ids(3) = 3
   882)         case(4)
   883)           vertex_ids(1) = 1
   884)           vertex_ids(2) = 3
   885)           vertex_ids(3) = 2
   886)         case default
   887)           option%io_buffer='Cell TET_TYPE has only 4 faces'
   888)           call printErrMsg(option)
   889)       end select       
   890)     case(QUAD_TYPE)
   891)       select case(iface)
   892)         case(1)
   893)           vertex_ids(1) = 1
   894)           vertex_ids(2) = 2
   895)         case(2)
   896)           vertex_ids(1) = 2
   897)           vertex_ids(2) = 3
   898)         case(3)
   899)           vertex_ids(1) = 3
   900)           vertex_ids(2) = 4
   901)         case(4)
   902)           vertex_ids(1) = 4
   903)           vertex_ids(2) = 1
   904)         case default
   905)           option%io_buffer='Cell QUAD_TYPE has only 4 faces'
   906)           call printErrMsg(option)
   907)       end select
   908)     case(TRI_TYPE)
   909)       select case(iface)
   910)         case(1)
   911)           vertex_ids(1) = 1
   912)           vertex_ids(2) = 2
   913)         case(2)
   914)           vertex_ids(1) = 2
   915)           vertex_ids(2) = 3
   916)         case(3)
   917)           vertex_ids(1) = 3
   918)           vertex_ids(2) = 1
   919)         case default
   920)           option%io_buffer='Cell TRI_TYPE has only 3 faces'
   921)           call printErrMsg(option)
   922)       end select
   923)     case default
   924)       option%io_buffer = 'Cell type not recognized'
   925)       call printErrMsg(option)
   926)   end select
   927) 
   928) end subroutine UCellGetFaceVertices
   929) 
   930) ! ************************************************************************** !
   931) 
   932) subroutine UCellGetEdgeVertices(cell_type,iedge,vertex_ids)
   933)   ! 
   934)   ! returns vertex ids of an edge
   935)   ! 
   936)   ! Author: Glenn Hammond
   937)   ! Date: 01/17/12
   938)   ! 
   939) 
   940)   implicit none
   941)   
   942)   PetscInt :: cell_type
   943)   PetscInt :: iedge
   944)   PetscInt :: vertex_ids(2)
   945)   
   946)   vertex_ids = UNINITIALIZED_INTEGER
   947)   
   948)   select case(cell_type)
   949)     case(HEX_TYPE)
   950)       select case(iedge)
   951)         case(1)
   952)           vertex_ids(1) = 1
   953)           vertex_ids(2) = 2
   954)         case(2)
   955)           vertex_ids(1) = 2
   956)           vertex_ids(2) = 3
   957)         case(3)
   958)           vertex_ids(1) = 3
   959)           vertex_ids(2) = 4
   960)         case(4)
   961)           vertex_ids(1) = 4
   962)           vertex_ids(2) = 1
   963)         case(5)
   964)           vertex_ids(1) = 5
   965)           vertex_ids(2) = 6
   966)         case(6)
   967)           vertex_ids(1) = 6
   968)           vertex_ids(2) = 7
   969)         case(7)
   970)           vertex_ids(1) = 7
   971)           vertex_ids(2) = 8
   972)         case(8)
   973)           vertex_ids(1) = 8
   974)           vertex_ids(2) = 5
   975)         case(9)
   976)           vertex_ids(1) = 1
   977)           vertex_ids(2) = 5
   978)         case(10)
   979)           vertex_ids(1) = 2
   980)           vertex_ids(2) = 6
   981)         case(11)
   982)           vertex_ids(1) = 3
   983)           vertex_ids(2) = 7
   984)         case(12)
   985)           vertex_ids(1) = 4
   986)           vertex_ids(2) = 8
   987)       end select
   988)     case(WEDGE_TYPE)
   989)       select case(iedge)
   990)         case(1)
   991)           vertex_ids(1) = 1
   992)           vertex_ids(2) = 2
   993)         case(2)
   994)           vertex_ids(1) = 2
   995)           vertex_ids(2) = 3
   996)         case(3)
   997)           vertex_ids(1) = 3
   998)           vertex_ids(2) = 1
   999)         case(4)
  1000)           vertex_ids(1) = 4
  1001)           vertex_ids(2) = 5
  1002)         case(5)
  1003)           vertex_ids(1) = 5
  1004)           vertex_ids(2) = 6
  1005)         case(6)
  1006)           vertex_ids(1) = 6
  1007)           vertex_ids(2) = 4
  1008)         case(7)
  1009)           vertex_ids(1) = 1
  1010)           vertex_ids(2) = 4
  1011)         case(8)
  1012)           vertex_ids(1) = 2
  1013)           vertex_ids(2) = 5
  1014)         case(9)
  1015)           vertex_ids(1) = 3
  1016)           vertex_ids(2) = 6
  1017)        end select
  1018)     case(PYR_TYPE)
  1019)       select case(iedge)
  1020)         case(1)
  1021)           vertex_ids(1) = 1
  1022)           vertex_ids(2) = 2
  1023)         case(2)
  1024)           vertex_ids(1) = 2
  1025)           vertex_ids(2) = 3
  1026)         case(3)
  1027)           vertex_ids(1) = 3
  1028)           vertex_ids(2) = 4
  1029)         case(4)
  1030)           vertex_ids(1) = 4
  1031)           vertex_ids(2) = 1
  1032)         case(5)
  1033)           vertex_ids(1) = 1
  1034)           vertex_ids(2) = 5
  1035)         case(6)
  1036)           vertex_ids(1) = 2
  1037)           vertex_ids(2) = 5
  1038)         case(7)
  1039)           vertex_ids(1) = 3
  1040)           vertex_ids(2) = 5
  1041)         case(8)
  1042)           vertex_ids(1) = 4
  1043)           vertex_ids(2) = 5
  1044)        end select
  1045)     case(TET_TYPE)
  1046)       select case(iedge)
  1047)         case(1)
  1048)           vertex_ids(1) = 1
  1049)           vertex_ids(2) = 2
  1050)         case(2)
  1051)           vertex_ids(1) = 2
  1052)           vertex_ids(2) = 3
  1053)         case(3)
  1054)           vertex_ids(1) = 3
  1055)           vertex_ids(2) = 1
  1056)         case(4)
  1057)           vertex_ids(1) = 1
  1058)           vertex_ids(2) = 4
  1059)         case(5)
  1060)           vertex_ids(1) = 2
  1061)           vertex_ids(2) = 4
  1062)         case(6)
  1063)           vertex_ids(1) = 3
  1064)           vertex_ids(2) = 4
  1065)       end select       
  1066)     case(QUAD_TYPE)
  1067)       select case(iedge)
  1068)         case(1)
  1069)           vertex_ids(1) = 1
  1070)           vertex_ids(2) = 2
  1071)         case(2)
  1072)           vertex_ids(1) = 2
  1073)           vertex_ids(2) = 3
  1074)         case(3)
  1075)           vertex_ids(1) = 3
  1076)           vertex_ids(2) = 4
  1077)         case(4)
  1078)           vertex_ids(1) = 4
  1079)           vertex_ids(2) = 1
  1080)       end select
  1081)     case(TRI_TYPE)
  1082)       select case(iedge)
  1083)         case(1)
  1084)           vertex_ids(1) = 1
  1085)           vertex_ids(2) = 2
  1086)         case(2)
  1087)           vertex_ids(1) = 2
  1088)           vertex_ids(2) = 3
  1089)         case(3)
  1090)           vertex_ids(1) = 3
  1091)           vertex_ids(2) = 1
  1092)       end select
  1093)   end select
  1094) 
  1095) end subroutine UCellGetEdgeVertices
  1096) 
  1097) ! ************************************************************************** !
  1098) 
  1099) function UCellGetEdgeLength(cell_type,iedge,vertices)
  1100)   ! 
  1101)   ! returns length of an edge
  1102)   ! 
  1103)   ! Author: Glenn Hammond
  1104)   ! Date: 01/17/12
  1105)   ! 
  1106) 
  1107)   implicit none
  1108)   
  1109)   PetscInt :: cell_type
  1110)   PetscInt :: iedge
  1111)   type(point_type) :: vertices(*)
  1112)   
  1113)   PetscReal :: UCellGetEdgeLength
  1114)   
  1115)   PetscReal :: len_x, len_y, len_z
  1116)   PetscInt :: vertex_ids(2)
  1117)   
  1118)   call UCellGetEdgeVertices(cell_type,iedge,vertex_ids)
  1119)   
  1120)   len_x = vertices(vertex_ids(1))%x-vertices(vertex_ids(2))%x
  1121)   len_y = vertices(vertex_ids(1))%y-vertices(vertex_ids(2))%y
  1122)   len_z = vertices(vertex_ids(1))%z-vertices(vertex_ids(2))%z
  1123)   UCellGetEdgeLength = sqrt(len_x*len_x+len_y*len_y+len_z*len_z)
  1124) 
  1125) end function UCellGetEdgeLength
  1126) 
  1127) ! ************************************************************************** !
  1128) 
  1129) function UCellQuality(cell_type,vertices,option)
  1130)   ! 
  1131)   ! returns vertex ids of a face
  1132)   ! 
  1133)   ! Author: Glenn Hammond
  1134)   ! Date: 01/17/12
  1135)   ! 
  1136) 
  1137)   use Option_module
  1138) 
  1139)   implicit none
  1140)   
  1141)   PetscInt :: cell_type
  1142)   type(point_type) :: vertices(*)
  1143)   type(option_type) :: option
  1144)   
  1145)   PetscReal :: UCellQuality
  1146)   
  1147)   PetscInt :: iedge
  1148)   PetscReal :: max_side, min_side
  1149)   PetscReal :: len
  1150) 
  1151)   max_side = -1.d20
  1152)   min_side = 1.d20
  1153)   do iedge = 1, UCellGetNVertices(cell_type,option)
  1154)     len = UCellGetEdgeLength(cell_type,iedge,vertices)
  1155)     if (len > max_side) max_side = len
  1156)     if (len < min_side) min_side = len
  1157)   enddo
  1158)   
  1159)   UCellQuality = max_side / min_side
  1160) 
  1161) end function UCellQuality
  1162) 
  1163) ! ************************************************************************** !
  1164) 
  1165) subroutine UCellGetLineIntercept(line_start,line_end,point,intercept)
  1166)   ! 
  1167)   ! Computes the intercept of a point with a line
  1168)   ! 
  1169)   ! Author: Gautam Bisht
  1170)   ! Date: 02/26/12
  1171)   ! 
  1172) 
  1173)   implicit none
  1174)   type(point_type) :: line_start
  1175)   type(point_type) :: line_end
  1176)   type(point_type) :: point
  1177)   type(point_type) :: intercept
  1178) 
  1179)   PetscReal :: dx,dy,dz
  1180)   PetscReal :: u, line_mag
  1181)   
  1182)   dx = (line_end%x - line_start%x)
  1183)   dy = (line_end%y - line_start%y)
  1184)   dz = (line_end%z - line_start%z)
  1185)   
  1186)   line_mag = (dx*dx + dy*dy + dz*dz)
  1187)   
  1188)   u = ((point%x - line_start%x)*(line_end%x - line_start%x) + &
  1189)        (point%y - line_start%y)*(line_end%y - line_start%y) + &
  1190)        (point%z - line_start%z)*(line_end%z - line_start%z))/ line_mag
  1191) 
  1192)   intercept%x = line_start%x + u*(line_end%x - line_start%x)
  1193)   intercept%y = line_start%y + u*(line_end%y - line_start%y)
  1194)   intercept%z = line_start%z + u*(line_end%z - line_start%z)
  1195) 
  1196) end subroutine
  1197) 
  1198) end module Grid_Unstructured_Cell_module

generated by
Intel(R) C++/Fortran Compiler code-coverage tool
Web-Page Owner: Nobody