ADflow  v1.0
ADflow is a finite volume RANS solver tailored for gradient-based aerodynamic design optimization.
surfaceFamilies.F90
Go to the documentation of this file.
2 
3  use constants
4 #ifndef USE_TAPENADE
5 #include <petsc/finclude/petsc.h>
6  use petsc
7  implicit none
8 
10  ! Vectors for global traction calc
11 
12  ! Parallel vector of un-uniqufied values concatenated for all
13  ! included surfaces
14  vec nodevallocal
15 
16  ! Parallel vector of uniqueifed values.
17  vec nodevalglobal
18 
19  ! Sum global. Same size as nodeValGlobal
20  vec sumglobal
21 
22  ! Scatter from nodeValLocal to NodeValGlobal
23  vecscatter scatter
24 
25  ! Flag for allocated petsc variables
26  logical :: allocated = .false.
27  integer(kind=intType), dimension(:), allocatable :: famlist
28 
29  integer(Kind=intType) :: sps
30  integer(kind=intType) :: nnodes
31  end type familyexchange
32 
34  integer(kind=intType), pointer, dimension(:) :: famlist
35  end type bcgrouptype
36 
37  ! Generic PETSc scatters
38  is is1, is2
39 
40  ! The list of exchanges based on boundary condition type
41  type(familyexchange), dimension(:, :), allocatable, target :: bcfamexchange
42 
43  ! List of familis grouped by BC. See constants.F90 for the indices
44  ! to use for this array.
45  type(bcgrouptype), dimension(nFamExchange) :: bcfamgroups
46 
47  ! The full list of the family names
48  character(len=maxCGNSNameLen), dimension(:), allocatable :: famnames
49 
50  ! List of all families. This is just 1,2,3,4...nFam. It is just used
51  ! in fortran when a specific family is not required.
52  integer(kind=intType), dimension(:), allocatable :: fullfamlist
53 #endif
54 
55  ! Special BC array's that are sometime required for reducitons.
56  real(kind=realtype), dimension(:, :), allocatable, target :: zerocellval
57  real(kind=realtype), dimension(:, :), allocatable, target :: onecellval
58  real(kind=realtype), dimension(:, :), allocatable, target :: zeronodeval
59 
60 #ifndef USE_TAPENADE
61 contains
62 
63  subroutine getnfam(nfam)
64 
65  implicit none
66  integer(kind=inttype), intent(out) :: nfam
67  if (allocated(famnames)) then
68  nfam = size(famnames)
69  else
70  nfam = 0
71  end if
72 
73  end subroutine getnfam
74 
75  subroutine getfam(i, fam)
76  implicit none
77  character(len=maxCGNSNameLen), intent(out) :: fam
78  integer(kind=intType), intent(in) :: i
79 
80  if (allocated(famnames)) then
81  fam = famnames(i)
82  end if
83 
84  end subroutine getfam
85 
86  subroutine destroyfamilyexchange(exch)
87 
88  use constants
89  type(familyexchange) :: exch
90  integer(kind=intType) :: ierr
91  if (exch%allocated) then
92 
93  call vecdestroy(exch%nodeValLocal, ierr)
94  call vecdestroy(exch%nodeValGlobal, ierr)
95  call vecdestroy(exch%sumGlobal, ierr)
96  call vecscatterdestroy(exch%scatter, ierr)
97 
98  end if
99 
100  exch%allocated = .false.
101 
102  end subroutine destroyfamilyexchange
103 #endif
104 end module surfacefamilies
type(bcgrouptype), dimension(nfamexchange) bcfamgroups
subroutine destroyfamilyexchange(exch)
character(len=maxcgnsnamelen), dimension(:), allocatable famnames
real(kind=realtype), dimension(:, :), allocatable, target zeronodeval
real(kind=realtype), dimension(:, :), allocatable, target onecellval
subroutine getfam(i, fam)
real(kind=realtype), dimension(:, :), allocatable, target zerocellval
integer(kind=inttype), dimension(:), allocatable fullfamlist
subroutine getnfam(nfam)
type(familyexchange), dimension(:, :), allocatable, target bcfamexchange