gauss.F90       coverage:  50.00 %func     16.38 %block


     1) module Gauss_module
     2) 
     3)   use Grid_Unstructured_Cell_module
     4)   use PFLOTRAN_Constants_module
     5) 
     6)   implicit none
     7) 
     8) #include "petsc/finclude/petscsys.h"
     9) 
    10)   private
    11)   
    12)   PetscInt, parameter, public :: LINE_TYPE          = 7
    13)   
    14)   type, public :: gauss_type
    15)     PetscInt :: dim                 ! dimension
    16)     PetscInt :: EleType             ! Element type
    17)     PetscInt :: NGPTS               ! Number of gauss points
    18)     PetscReal, pointer :: r(:,:)    ! location of points
    19)     PetscReal, pointer :: w(:)      ! weights
    20)   end type gauss_type
    21)     
    22)   public :: GaussCalculatePoints, GaussDestroy, GaussInitialize
    23)   
    24)   contains
    25) 
    26) ! ************************************************************************** !
    27) 
    28) subroutine GaussInitialize(gauss)
    29)   ! 
    30)   ! Initializes Gauss type
    31)   ! 
    32)   ! Author: Satish Karra, LANL
    33)   ! Date: 6/19/2013
    34)   ! 
    35) 
    36)   type(gauss_type) :: gauss
    37)   
    38)   gauss%dim = 0 
    39)   gauss%EleType = 0
    40)   gauss%NGPTS = 0
    41)   nullify(gauss%r)
    42)   nullify(gauss%w)
    43) 
    44) end subroutine GaussInitialize   
    45) 
    46) ! ************************************************************************** !
    47) 
    48) subroutine GaussCalculatePoints(gauss)
    49)   ! 
    50)   ! Calculates Gauss points
    51)   ! 
    52)   ! Author: Satish Karra, LANL
    53)   ! Date: 5/17/2013
    54)   ! 
    55) 
    56)   use Utility_module, only: DeallocateArray
    57) 
    58)   type(gauss_type) :: gauss
    59)   PetscReal, pointer :: r(:,:)
    60)   PetscReal, pointer :: w(:)
    61)   
    62)   select case(gauss%dim)
    63)     case(ONE_DIM_GRID)
    64)       call Gauss1D(gauss%EleType,gauss%NGPTS,r,w)
    65)     case(TWO_DIM_GRID)
    66)       call Gauss2D(gauss%EleType,gauss%NGPTS,r,w)
    67)     case(THREE_DIM_GRID)
    68)       call Gauss3D(gauss%EleType,gauss%NGPTS,r,w)
    69)     case default
    70)       print *, 'Error: Invalid dimension for Gauss point calculation'
    71)       stop
    72)   end select  
    73)   
    74)   allocate(gauss%r(size(r,1),size(r,2)))
    75)   allocate(gauss%w(size(w)))
    76)   
    77)   gauss%r = r
    78)   gauss%w = w
    79)   
    80)   deallocate(r)
    81)   deallocate(w)
    82)     
    83) end subroutine GaussCalculatePoints  
    84) 
    85) ! ************************************************************************** !
    86) 
    87) subroutine Gauss1D(EleType,NGPTS,r,w)
    88)   ! 
    89)   ! Calculates Gauss points for 1D elements
    90)   ! 
    91)   ! Author: Satish Karra, LANL
    92)   ! Date: 5/17/2013
    93)   ! 
    94) 
    95)   PetscInt :: EleType
    96)   PetscInt :: NGPTS
    97)   PetscReal, pointer :: r(:,:)
    98)   PetscReal, pointer :: w(:)
    99)   
   100)   allocate(r(NGPTS,1))
   101)   allocate(w(NGPTS))
   102)   
   103)   if (EleType /= LINE_TYPE) then
   104)     print *, 'Error: in Element type. Only L2 ' // &
   105)              '(line type) can be used for 1D Gauss quadrature.'
   106)   endif
   107)     
   108)   select case(NGPTS)
   109)     case(1)
   110) 
   111)     !------------------------------
   112)     ! No of Gauss Points = 1
   113)     !------------------------------
   114)         
   115)       r(1,1) = 0.0
   116)       !
   117)       w(1) = 2.0
   118) 
   119)     !------------------------------
   120)     ! No of Gauss Points = 2
   121)     !------------------------------
   122) 
   123)     case(2)
   124)         
   125)       r(1,1) = -1.0/sqrt(3.0)
   126)       r(2,1) = -r(1,1) 
   127)       !
   128)       w(1) = 1.0
   129)       w(2) = 1.0
   130) 
   131)     !-------------------------------
   132)     ! No of Gauss Points = 3
   133)     !-------------------------------
   134) 
   135)     case(3)
   136)         
   137)       r(1,1) = -sqrt(0.6)
   138)       r(2,1) = 0.0
   139)       r(3,1) = -r(1,1)
   140)       !
   141)       w(1) = 5.0/9.0
   142)       w(2) = 8.0/9.0
   143)       w(3) = w(1)
   144) 
   145)     !--------------------------------
   146)     ! No of Gauss Points = 4
   147)     !--------------------------------
   148) 
   149)     case(4)
   150)         
   151)       r(1,1) = -0.861136311594053
   152)       r(2,1) = -0.339981043584856
   153)       r(3,1) =  0.339981043584856
   154)       r(4,1) =  0.861136311594053
   155)       ! 
   156)       w(1) = 0.347854845137454
   157)       w(2) = 0.652145154862546
   158)       w(3) = 0.652145154862546
   159)       w(4) = 0.347854845137454
   160) 
   161)     !----------------------------------
   162)     ! No of Gauss Points = 5
   163)     !----------------------------------
   164) 
   165)     case(5)
   166)         
   167)       r(1,1) = -0.906179845938664
   168)       r(2,1) = -0.538469310105683
   169)       r(3,1) =  0.000000000000000
   170)       r(4,1) =  0.538469310105683
   171)       r(5,1) =  0.906179845938664
   172)       !
   173)       w(1) =  0.236926885056189
   174)       w(2) =  0.478628670499366 
   175)       w(3) =  0.568888888888889
   176)       w(4) =  0.478628670499366
   177)       w(5) =  0.236926885056189
   178) 
   179)     !----------------------------------
   180)     ! No of Gauss Points = 6
   181)     !----------------------------------
   182)    
   183)     case(6)
   184)         
   185)       r(1,1) = -0.932469514203152
   186)       r(2,1) = -0.661209386466265
   187)       r(3,1) = -0.238619186083197
   188)       r(4,1) =  0.238619186083197
   189)       r(5,1) =  0.661209386466265
   190)       r(6,1) =  0.932469514203152
   191)       !
   192)       w(1) =  0.171324492379170
   193)       w(2) =  0.360761573048139
   194)       w(3) =  0.467913934572691
   195)       w(4) =  0.467913934572691
   196)       w(5) =  0.360761573048139
   197)       w(6) =  0.171324492379170
   198) 
   199)     !------------------------------------
   200)     ! No of Gauss Points = 7
   201)     !------------------------------------
   202) 
   203)     case(7)
   204)         
   205)       r(1,1) = -0.949107912342759
   206)       r(2,1) = -0.741531185599394
   207)       r(3,1) = -0.405845151377397
   208)       r(4,1) =  0.000000000000000
   209)       r(5,1) =  0.405845151377397
   210)       r(6,1) =  0.741531185599394
   211)       r(7,1) =  0.949107912342759
   212)       !
   213)       w(1) =  0.129484966168870
   214)       w(2) =  0.279705391489277
   215)       w(3) =  0.381830050505119
   216)       w(4) =  0.417959183673469
   217)       w(5) =  0.381830050505119
   218)       w(6) =  0.279705391489277
   219)       w(7) =  0.129484966168870
   220) 
   221)     !------------------------------------
   222)     ! No of Gauss Points = 8
   223)     !------------------------------------
   224)     
   225)     case(8)
   226)         
   227)       r(1,1) = -0.960289856497536
   228)       r(2,1) = -0.796666477413627
   229)       r(3,1) = -0.525532409916329
   230)       r(4,1) = -0.183434642495650
   231)       r(5,1) =  0.183434642495650
   232)       r(6,1) =  0.525532409916329
   233)       r(7,1) =  0.796666477413627
   234)       r(8,1) =  0.960289856497536
   235)       !
   236)       w(1) =  0.101228536290376
   237)       w(2) =  0.222381034453374
   238)       w(3) =  0.313706645877887
   239)       w(4) =  0.362683783378362
   240)       w(5) =  0.362683783378362
   241)       w(6) =  0.313706645877887
   242)       w(7) =  0.222381034453374
   243)       w(8) =  0.101228536290376
   244) 
   245)     !------------------------------------
   246)     ! No of Gauss Points = 9
   247)     !------------------------------------
   248)  
   249)     case(9)
   250)         
   251)       r(1,1) = -0.968160239507626
   252)       r(2,1) = -0.836031170326636
   253)       r(3,1) = -0.613371432700590
   254)       r(4,1) = -0.324253423403809
   255)       r(5,1) =  0.000000000000000
   256)       r(6,1) =  0.324253423403809
   257)       r(7,1) =  0.613371432700590
   258)       r(8,1) =  0.836031107326636
   259)       r(9,1) =  0.968160239507626
   260) 
   261)       w(1) =  0.081274388361574
   262)       w(2) =  0.180648160694857
   263)       w(3) =  0.260610696402935
   264)       w(4) =  0.312347077040003
   265)       w(5) =  0.330239355001260
   266)       w(6) =  0.312347077040003
   267)       w(7) =  0.260610696402935
   268)       w(8) =  0.180648160694857
   269)       w(9) =  0.081274388361574
   270) 
   271)    case default
   272)      print *, 'Error in NGPTS for 1D Gauss quadrature'
   273)      stop
   274)    end select
   275) 
   276) end subroutine Gauss1D  
   277) 
   278) ! ************************************************************************** !
   279) 
   280) subroutine Gauss2D(EleType,NGPTS,r,w)
   281)   ! 
   282)   ! Calculates Gauss points for 2D elements
   283)   ! 
   284)   ! Author: Satish Karra, LANL
   285)   ! Date: 5/17/2013
   286)   ! 
   287) 
   288)   PetscInt :: EleType
   289)   PetscInt :: NGPTS
   290)   PetscReal, pointer :: r(:,:)
   291)   PetscReal, pointer :: w(:)
   292)     
   293)   select case(EleType)
   294)     case(QUAD_TYPE)
   295)       call GaussSquare(NGPTS,r,w)
   296)     case(TRI_TYPE)
   297)       call GaussTriangle(NGPTS,r,w)
   298)     case default
   299)       print *, 'Error: Only T3 and Q4 elements available for 2D.'
   300)       stop
   301)   end select
   302) 
   303) end subroutine Gauss2D
   304) 
   305) ! ************************************************************************** !
   306) 
   307) subroutine GaussSquare(NGPTS,r,w)
   308)   ! 
   309)   ! Calculates Gauss points for Q4 element
   310)   ! 
   311)   ! Author: Satish Karra, LANL
   312)   ! Date: 5/17/2013
   313)   ! 
   314) 
   315)   PetscInt :: NGPTS
   316)   PetscReal, pointer :: r(:,:)
   317)   PetscReal, pointer :: w(:)
   318)   PetscReal, pointer :: l(:,:)
   319)   PetscReal, pointer :: m(:)
   320)   PetscInt :: counter,i,j 
   321)   
   322)   allocate(r(NGPTS*NGPTS,2))
   323)   allocate(w(NGPTS*NGPTS))
   324)   allocate(l(NGPTS,1))
   325)   allocate(m(NGPTS))
   326) 
   327)   ! 1D Gauss points are stored in l vector and weights are stored in m vector
   328) 
   329)   call Gauss1D(LINE_TYPE,NGPTS,l,m)
   330)   
   331)   ! Generate the Q4 Gauss points and weights using for loops
   332)   counter = 1
   333)   do i = 1, NGPTS
   334)     do j = 1, NGPTS
   335)       r(counter,1) = l(i,1)
   336)       r(counter,2) = l(j,1)
   337)       w(counter) = m(i)*m(j)
   338)       counter = counter + 1
   339)     enddo
   340)   enddo
   341) 
   342)   deallocate(l)
   343)   deallocate(m)
   344) 
   345) end subroutine GaussSquare
   346) 
   347) ! ************************************************************************** !
   348) 
   349) subroutine GaussTriangle(NGPTS,r,w)
   350)   ! 
   351)   ! Calculates Gauss points for T3 elements
   352)   ! 
   353)   ! Author: Satish Karra, LANL
   354)   ! Date: 5/17/2013
   355)   ! 
   356) 
   357)   PetscInt :: EleType
   358)   PetscInt :: NGPTS
   359)   PetscReal, pointer :: r(:,:)
   360)   PetscReal, pointer :: w(:)
   361)   
   362)   allocate(r(NGPTS,2))
   363)   allocate(w(NGPTS))
   364)     
   365)   select case(NGPTS)
   366)     case(1)
   367) 
   368)     !------------------------------
   369)     ! No of Gauss Points = 1
   370)     !------------------------------
   371)       r(1,:) = 1.0/3.0*(/1.0,1.0/)
   372)       !
   373)       w(1) = 0.5
   374) 
   375)     !-------------------------------
   376)     ! No of Gauss Points = 3
   377)     !-------------------------------
   378) 
   379)     case(3)
   380)         
   381)       r(1,:) = 0.5*(/1.0,1.0/)
   382)       r(2,:) = 0.5*(/1.0,0.0/)
   383)       r(3,:) = 0.5*(/0.0,1.0/) 
   384)       !
   385)       w = 1.0/6.0*(/1.0,1.0,1.0/)
   386) 
   387)     !--------------------------------
   388)     ! No of Gauss Points = 4
   389)     !--------------------------------
   390) 
   391)     case(4)
   392)     
   393)       r(1,:) = (/0.333333333333333,0.333333333333333/)
   394)       r(2,:) = (/0.600000000000000,0.200000000000000/)
   395)       r(3,:) = (/0.200000000000000,0.600000000000000/)
   396)       r(4,:) = (/0.200000000000000,0.200000000000000/)
   397)       ! 
   398)       w = 0.5*(/-0.562500000000000, &
   399)                  0.520833333333333, &
   400)                  0.520833333333333, &
   401)                  0.520833333333333/)
   402) 
   403)     !------------------------------------
   404)     ! No of Gauss Points = 7
   405)     !------------------------------------
   406) 
   407)     case(7)
   408)     
   409)       r(1,:) = (/0.333333333333333,0.333333333333333/)
   410)       r(2,:) = (/0.797426985353087,0.101286507323456/)
   411)       r(3,:) = (/0.101286507323456,0.797426985353087/)
   412)       r(4,:) = (/0.101286507323456,0.101286507323456/)
   413)       r(5,:) = (/0.470142064105115,0.059715871789770/)
   414)       r(6,:) = (/0.059715871789770,0.470142064105115/)
   415)       r(7,:) = (/0.470142064105115,0.470142064105115/)
   416)       !
   417)       w = 0.5*(/0.225000000000000, &
   418)                 0.125939180544827, &
   419)                 0.125939180544827, &
   420)                 0.125939180544827, &
   421)                 0.132394152788506, &
   422)                 0.132394152788506, &
   423)                 0.132394152788506/)
   424) 
   425)     !------------------------------------
   426)     ! No of Gauss Points = 9
   427)     !------------------------------------
   428)  
   429)     case(9)
   430)        
   431)       r(1,:) = (/0.124949503233232,0.437525248383384/)
   432)       r(2,:) = (/0.437525248383384,0.124949503233232/)
   433)       r(3,:) = (/0.437525248383384,0.437525248383384/)
   434)       r(4,:) = (/0.797112651860071,0.165409927389841/)
   435)       r(5,:) = (/0.797112651860071,0.037477420750088/)
   436)       r(6,:) = (/0.165409927389841,0.797112651860071/)
   437)       r(7,:) = (/0.165409927389841,0.037477420750088/)
   438)       r(8,:) = (/0.037477420750088,0.797112651860071/)
   439)       r(9,:) = (/0.037477420750088,0.165409927389841/)
   440) 
   441)       w = 0.5*(/0.205950504760887, &
   442)                 0.205950504760887, &
   443)                 0.205950504760887, &
   444)                 0.063691414286223, &
   445)                 0.063691414286223, &
   446)                 0.063691414286223, &
   447)                 0.063691414286223, &
   448)                 0.063691414286223, &
   449)                 0.063691414286223/)
   450)                 
   451)     !------------------------------------
   452)     ! No of Gauss Points = 12
   453)     !------------------------------------
   454)  
   455)     case(12)
   456)        
   457)       r(1,:) = (/0.873821971016996,0.063089014491502/)
   458)       r(2,:) = (/0.063089014491502,0.873821971016996/) 
   459)       r(3,:) = (/0.063089014491502,0.063089014491502/)
   460)       r(4,:) = (/0.501426509658179,0.249286745170910/)
   461)       r(5,:) = (/0.249286745170910,0.501426509658179/) 
   462)       r(6,:) = (/0.249286745170910,0.249286745170910/)
   463)       r(7,:) = (/0.636502499121399,0.310352451033785/)
   464)       r(8,:) = (/0.636502499121399,0.053145049844816/)
   465)       r(9,:) = (/0.310352451033785,0.636502499121399/)
   466)       r(10,:) = (/0.310352451033785,0.053145049844816/)
   467)       r(11,:) = (/0.053145049844816,0.636502499121399/) 
   468)       r(12,:) = (/0.053145049844816,0.310352451033785/)
   469) 
   470)       w = 0.5*(/0.050844906370207, &
   471)                 0.050844906370207, &
   472)                 0.050844906370207, &
   473)                 0.116786275726379, &
   474)                 0.116786275726379, &
   475)                 0.116786275726379, &
   476)                 0.082851075618374, &
   477)                 0.082851075618374, &
   478)                 0.082851075618374, &
   479)                 0.082851075618374, &
   480)                 0.082851075618374, &
   481)                 0.082851075618374/)
   482)                 
   483)     !------------------------------------
   484)     ! No of Gauss Points = 13
   485)     !------------------------------------
   486)  
   487)     case(13)
   488)        
   489)       r(1,:) = (/0.333333333333333,0.333333333333333/)
   490)       r(2,:) = (/0.479308067841923,0.260345966079038/)
   491)       r(3,:) = (/0.260345966079038,0.479308067841923/)
   492)       r(4,:) = (/0.260345966079038,0.260345966079038/)
   493)       r(5,:) = (/0.869739794195568,0.065130102902216/)
   494)       r(6,:) = (/0.065130102902216,0.869739794195568/)
   495)       r(7,:) = (/0.065130102902216,0.065130102902216/)
   496)       r(8,:) = (/0.638444188569809,0.312865496004875/)
   497)       r(9,:) = (/0.638444188569809,0.086903154253160/)
   498)       r(10,:) = (/0.312865496004875,0.638444188569809/)
   499)       r(11,:) = (/0.312865496004875,0.086903154253160/)
   500)       r(12,:) = (/0.086903154253160,0.638444188569809/)
   501)       r(13,:) = (/0.086903154253160,0.312865496004875/) 
   502)       
   503)       w = 0.5*(/-0.149570044467670, &
   504)                 +0.175615257433204, &
   505)                 +0.175615257433204, &
   506)                 +0.175615257433204, &
   507)                 +0.053347235608839, &
   508)                 +0.053347235608839, &
   509)                 +0.053347235608839, &
   510)                 +0.077113760890257, &
   511)                 +0.077113760890257, &
   512)                 +0.077113760890257, &
   513)                 +0.077113760890257, &
   514)                 +0.077113760890257, &
   515)                 +0.077113760890257/)
   516)               
   517)    case default
   518)      print *, 'Invalid NGPTS for T3 Gauss quadrature'
   519)      stop
   520)    end select
   521)   
   522)   
   523) end subroutine GaussTriangle
   524) 
   525) ! ************************************************************************** !
   526) 
   527) subroutine GaussTetrahedra(NGPTS,r,w)
   528)   ! 
   529)   ! Calculates Gauss points for tetrahedra elements
   530)   ! 
   531)   ! Author: Satish Karra, LANL
   532)   ! Date: 7/11/2013
   533)   ! 
   534) 
   535)   PetscInt :: EleType
   536)   PetscInt :: NGPTS
   537)   PetscReal, pointer :: r(:,:)
   538)   PetscReal, pointer :: w(:)
   539)   PetscReal :: p1,p2,p3,p4,p5,p6,p7
   540)   PetscReal :: q1,q2,q3,q4
   541)   
   542)   allocate(r(NGPTS,3))
   543)   allocate(w(NGPTS))
   544)     
   545)   select case(NGPTS)
   546)     case(1)
   547) 
   548)     !------------------------------
   549)     ! No of Gauss Points = 1
   550)     !------------------------------
   551)       r(1,:) = 1.0/4.0*(/1.0,1.0,1.0/)
   552)       !
   553)       w(1) = 1.0/6.0
   554) 
   555)     !--------------------------------
   556)     ! No of Gauss Points = 4
   557)     !--------------------------------
   558) 
   559)     case(4)
   560)     
   561)       p1 = 0.5854101966249638
   562)       p2 = 0.1381966011250105
   563)       r(1,:) = (/p1,p2,p2/)
   564)       r(2,:) = (/p2,p1,p2/)
   565)       r(3,:) = (/p2,p2,p1/)
   566)       r(4,:) = (/p1,p1,p1/)
   567)       ! 
   568)       w = 1.0/(6.0*4.0)*(/1.0,1.0,1.0,1.0/)
   569)       
   570)     !------------------------------------
   571)     ! No of Gauss Points = 5
   572)     !------------------------------------     
   573)     
   574)     case(5)
   575)     
   576)       r(1,:) = 1.0/4.0*(/1.0,1.0,1.0/)
   577)       r(2,:) = (/1.0/2.0,1.0/6.0,1.0/6.0/)
   578)       r(3,:) = (/1.0/6.0,1.0/2.0,1.0/6.0/)
   579)       r(4,:) = (/1.0/6.0,1.0/6.0,1.0/2.0/)
   580)       r(5,:) = (/1.0/6.0,1.0/6.0,1.0/6.0/)
   581)       !
   582)       w = 1.0/6.0*(/-4.0/5.0,9.0/20.0,9.0/20.0,9.0/20.0,9.0/20.0/)
   583) 
   584)     !------------------------------------
   585)     ! No of Gauss Points = 11
   586)     !------------------------------------
   587) 
   588)     case(11)
   589)     
   590)       p1 = 0.250000000000000  
   591)       p2 = 0.785714285714286
   592)       p3 = 0.071428571428571
   593)       p4 = 0.399403576166799
   594)       p5 = 0.100596423833201
   595) 
   596)       r(1,:) = (/p1,p1,p1/)
   597)       r(2,:) = (/p2,p3,p3/)
   598)       r(3,:) = (/p3,p2,p3/)
   599)       r(4,:) = (/p3,p3,p2/)
   600)       r(5,:) = (/p3,p3,p3/)
   601)       r(6,:)  = (/p4,p5,p5/)
   602)       r(7,:)  = (/p5,p4,p5/)
   603)       r(8,:)  = (/p5,p5,p4/)
   604)       r(9,:)  = (/p5,p4,p4/)
   605)       r(10,:) = (/p4,p5,p4/)
   606)       r(11,:) = (/p4,p4,p5/)
   607)       !
   608)       q1 = -0.013155555555556
   609)       q2 =  0.007622222222222
   610)       q3 =  0.024888888888889
   611) 
   612)       w = (/q1,q2,q2,q2,q2,q3,q3,q3,q3,q3,q3/)
   613)       
   614) 
   615)     !------------------------------------
   616)     ! No of Gauss Points = 15
   617)     !------------------------------------
   618)  
   619)     case(15)
   620)  
   621)       p1 = 0.250000000000000
   622)       p2 = 0.000000000000000
   623)       p3 = 0.333333333333333
   624)       p4 = 0.727272727272727
   625)       p5 = 0.090909090909091
   626)       p6 = 0.066550153573664
   627)       p7 = 0.433449846426336
   628) 
   629)       r(1,:) = (/p1,p1,p1/)
   630)       r(2,:) = (/p2,p3,p3/)
   631)       r(3,:) = (/p3,p2,p3/)
   632)       r(4,:) = (/p3,p3,p2/)
   633)       r(5,:) = (/p3,p3,p3/)
   634)       r(6,:) = (/p4,p5,p5/)
   635)       r(7,:) = (/p5,p4,p5/)
   636)       r(8,:) = (/p5,p5,p4/)
   637)       r(9,:) = (/p5,p5,p5/)
   638)       r(10,:) = (/p6,p7,p7/)
   639)       r(11,:) = (/p7,p6,p7/)
   640)       r(12,:) = (/p7,p7,p6/)
   641)       r(13,:) = (/p7,p6,p6/)
   642)       r(14,:) = (/p6,p7,p6/)
   643)       r(15,:) = (/p6,p6,p7/)
   644)       !
   645)       q1 = 0.030283678097089
   646)       q2 = 0.006026785714286
   647)       q3 = 0.011645249086029
   648)       q4 = 0.010949141561386
   649) 
   650)       w = (/q1,q2,q2,q2, q2,q3,q3,q3,q3,q4,q4,q4,q4,q4,q4/)
   651) 
   652)     case default
   653)      print *, 'Invalid NGPTS for Tetrahedra Gauss quadrature'
   654)      stop
   655)    end select
   656)   
   657)   
   658) end subroutine GaussTetrahedra
   659) 
   660) ! ************************************************************************** !
   661) 
   662) subroutine GaussPyramid(NGPTS,r,w)
   663)   ! 
   664)   ! Calculates Gauss points for tetrahedra elements
   665)   ! Reference:
   666)   ! http://people.sc.fsu.edu/~jburkardt/datasets/quadrature_rules_pyramid/quadrature_rules_pyramid.html
   667)   ! 
   668)   ! Author: Satish Karra, LANL
   669)   ! Date: 7/11/2013
   670)   ! 
   671) 
   672)   PetscInt :: EleType
   673)   PetscInt :: NGPTS
   674)   PetscReal, pointer :: r(:,:)
   675)   PetscReal, pointer :: w(:)
   676)   
   677)   allocate(r(NGPTS,3))
   678)   allocate(w(NGPTS))
   679)     
   680)   select case(NGPTS)
   681)     case(1)
   682) 
   683)     !------------------------------
   684)     ! No of Gauss Points = 1
   685)     !------------------------------
   686)       r(1,:) = (/0.0,0.0,0.25/)
   687)       !
   688)       w(1) = 1.0
   689) 
   690)     !--------------------------------
   691)     ! No of Gauss Points = 5
   692)     !--------------------------------
   693) 
   694)     case(5)
   695)     
   696)     
   697)       r(1,:) = (/-0.48686449556014765641,-0.48686449556014765641,0.16666666666666666667/) 
   698)       r(2,:) = (/ 0.48686449556014765641,-0.48686449556014765641,0.16666666666666666667/) 
   699)       r(3,:) = (/ 0.48686449556014765641, 0.48686449556014765641,0.16666666666666666667/) 
   700)       r(4,:) = (/-0.48686449556014765641, 0.48686449556014765641,0.16666666666666666667/)
   701)       r(5,:) = (/ 0.00000000000000000000, 0.00000000000000000000,0.70000000000000000000/)
   702)       !
   703)       w = (/0.2109375, &
   704)             0.2109375, &
   705)             0.2109375, &
   706)             0.2109375, &
   707)             0.15625/)
   708)       
   709)     !------------------------------------
   710)     ! No of Gauss Points = 6
   711)     !------------------------------------     
   712)     
   713)     case(6)
   714)     
   715)       r(1,:) = (/-0.48795003647426658968,-0.48795003647426658968,0.16666666666666666667/)
   716)       r(2,:) = (/ 0.48795003647426658968,-0.48795003647426658968,0.16666666666666666667/)
   717)       r(3,:) = (/ 0.48795003647426658968, 0.48795003647426658968,0.16666666666666666667/)
   718)       r(4,:) = (/-0.48795003647426658968, 0.48795003647426658968,0.16666666666666666667/)
   719)       r(5,:) = (/ 0.00000000000000000000, 0.00000000000000000000,0.58333333333333333333/)
   720)       r(6,:) = (/ 0.00000000000000000000, 0.00000000000000000000,0.75000000000000000000/)
   721)       !
   722)       w = (/0.21000000000000000000, &
   723)             0.21000000000000000000, &
   724)             0.21000000000000000000, &
   725)             0.21000000000000000000, &
   726)             0.06000000000000000000, &
   727)             0.10000000000000000000/)
   728) 
   729)      case default
   730)      print *, 'Invalid NGPTS for Tetrahedra Gauss quadrature'
   731)      stop
   732)    end select
   733)   
   734)   
   735) end subroutine GaussPyramid
   736) 
   737) ! ************************************************************************** !
   738) 
   739) subroutine Gauss3D(EleType,NGPTS,r,w)
   740)   ! 
   741)   ! Calculates Gauss points for 3D element
   742)   ! 
   743)   ! Author: Satish Karra, LANL
   744)   ! Date: 5/17/2013
   745)   ! 
   746) 
   747)   PetscInt :: EleType
   748)   PetscInt :: NGPTS
   749)   PetscReal, pointer :: r(:,:)
   750)   PetscReal, pointer :: w(:)
   751)   
   752)   select case(EleType)
   753)     case(HEX_TYPE)
   754)       call GaussBrick(NGPTS,r,w)
   755)     case(WEDGE_TYPE)
   756)       call GaussWedge(NGPTS,r,w)
   757)     case(TET_TYPE)
   758)       call GaussTetrahedra(NGPTS,r,w)
   759)     case(PYR_TYPE)
   760)       call GaussPyramid(NGPTS,r,w)
   761)     case default
   762)       print *, 'Error: Only B8, W6, P5 and TET4 elements available for 3D.'
   763)       stop
   764)   end select
   765)   
   766) end subroutine Gauss3D
   767) 
   768) ! ************************************************************************** !
   769) 
   770) subroutine GaussBrick(NGPTS,r,w)
   771)   ! 
   772)   ! Calculates Gauss points for B8 element
   773)   ! 
   774)   ! Author: Satish Karra, LANL
   775)   ! Date: 5/17/2013
   776)   ! 
   777) 
   778)   PetscInt :: NGPTS
   779)   PetscReal, pointer :: r(:,:)
   780)   PetscReal, pointer :: w(:)
   781)   PetscReal, pointer :: l(:,:)
   782)   PetscReal, pointer :: m(:)
   783)   PetscInt :: counter, i, j, k
   784)   
   785)   allocate(r(NGPTS*NGPTS*NGPTS,3))
   786)   allocate(w(NGPTS*NGPTS*NGPTS))
   787)   
   788)   call Gauss1D(LINE_TYPE,NGPTS,l,m)
   789)   
   790)   ! Generate the B8 Gauss points and weights using for loops
   791) 
   792)   counter = 1
   793)   do i = 1, NGPTS
   794)     do j = 1, NGPTS
   795)       do k = 1, NGPTS
   796)         r(counter,1) = l(i,1)
   797)         r(counter,2) = l(j,1)
   798)         r(counter,3) = l(k,1)
   799)         w(counter) = m(i)*m(j)*m(k)
   800)         counter = counter + 1
   801)       enddo
   802)     enddo
   803)   enddo
   804)   
   805)   deallocate(l)
   806)   deallocate(m)
   807)   
   808) end subroutine GaussBrick
   809) 
   810) ! ************************************************************************** !
   811) 
   812) subroutine GaussWedge(NGPTS,r,w)
   813)   ! 
   814)   ! Calculates Gauss points for wedge element
   815)   ! 
   816)   ! Author: Satish Karra, LANL
   817)   ! Date: 7/10//2013
   818)   ! 
   819) 
   820)   PetscInt :: NGPTS
   821)   PetscReal, pointer :: r(:,:)
   822)   PetscReal, pointer :: w(:)
   823)   PetscReal, pointer :: rT3(:,:),rL2(:,:)  
   824)   PetscReal, pointer :: wT3(:),wL2(:)
   825)   PetscInt :: counter, i, j, k
   826)   
   827)   allocate(r(NGPTS*NGPTS,3))
   828)   allocate(w(NGPTS*NGPTS))
   829)   
   830)   call Gauss1D(LINE_TYPE,NGPTS,rL2,wL2)
   831)   call Gauss2D(TRI_TYPE,NGPTS,rT3,wT3)
   832)   
   833)   ! Generate the wedge Gauss points and weights using for loops
   834)   do i = 1, NGPTS
   835)     do j = 1, NGPTS
   836)       r((i-1)*NGPTS+j,1) = rT3(i,1)
   837)       r((i-1)*NGPTS+j,2) = rT3(i,2)
   838)       r((i-1)*NGPTS+j,3) = rL2(j,1)
   839)       w((i-1)*NGPTS+j) = wT3(i)*wL2(j)
   840)     enddo
   841)   enddo
   842)   
   843)   deallocate(rL2)
   844)   deallocate(rT3)
   845)   deallocate(wL2)
   846)   deallocate(wT3)
   847)   
   848) end subroutine GaussWedge
   849) 
   850) ! ************************************************************************** !
   851) 
   852) subroutine GaussDestroy(gauss)
   853)   ! 
   854)   ! Deallocate gauss type
   855)   ! 
   856)   ! Author: Satish Karra, LANL
   857)   ! Date: 5/17/2013
   858)   ! 
   859) 
   860)   type(gauss_type) :: gauss
   861)   
   862)   deallocate(gauss%r)
   863)   nullify(gauss%r)
   864)   deallocate(gauss%w)
   865)   nullify(gauss%w)
   866) 
   867) end subroutine GaussDestroy
   868)      
   869) end module Gauss_module

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