FE-Project
Loading...
Searching...
No Matches
scale_localmesh_base.F90
Go to the documentation of this file.
1!-------------------------------------------------------------------------------
9#include "scaleFElib.h"
11
12 !-----------------------------------------------------------------------------
13 !
14 !++ used modules
15 !
16 use scale_precision
17 use scale_io
18
20
21 !-----------------------------------------------------------------------------
22 implicit none
23 private
24
25 !-----------------------------------------------------------------------------
26 !
27 !++ Public type & procedure
28 !
29
30 type, public :: localmeshbase
31 integer :: ne
32 integer :: nes
33 integer :: nee
34 integer :: nea
35 integer :: nv
36
37 class(elementbase), pointer :: refelem
38
39 real(rp), allocatable :: pos_ev(:,:)
40 real(rp), allocatable :: pos_en(:,:,:)
41 real(rp), allocatable :: normal_fn(:,:,:)
42
43 real(rp), allocatable :: sj(:,:)
44 real(rp), allocatable :: j(:,:)
45
46 real(rp), allocatable :: escale(:,:,:,:)
47 real(rp), allocatable :: fscale(:,:)
48
49 integer, allocatable :: etov(:,:)
50 integer, allocatable :: etoe(:,:)
51 integer, allocatable :: etof(:,:)
52 integer, allocatable :: vmapm(:,:)
53 integer, allocatable :: vmapp(:,:)
54 integer, allocatable :: mapm(:,:)
55 integer, allocatable :: mapp(:,:)
56
57 integer, allocatable :: bctype(:,:)
58 integer, allocatable :: mapb(:)
59 integer, allocatable :: vmapb(:)
60 ! integer, allocatable :: MapD(:)
61 ! integer, allocatable :: VMapD(:)
62 ! integer, allocatable :: MapN(:)
63 ! integer, allocatable :: VMapN(:)
64
65 integer :: tileid
66 integer :: panelid
67 integer :: prc_myrank
68 integer :: lcdomid
69
70 real(rp), allocatable :: g_ij(:,:,:,:)
71 real(rp), allocatable :: gij(:,:,:,:)
72 real(rp), allocatable :: gsqrt(:,:)
73 end type localmeshbase
74
75 public :: localmeshbase_init
76 public :: localmeshbase_final
77
78 !-----------------------------------------------------------------------------
79 !
80 !++ Public parameters & variables
81 !
82 integer, public, parameter :: bctype_interior = 0
83 integer, public, parameter :: bctype_undefbc = 1
84 integer, public, parameter :: bctype_periodic = 2
85 integer, public, parameter :: bctype_dirchlet = 3
86 integer, public, parameter :: bctype_neuman = 4
87 integer, public, parameter :: bctype_shoreline = 5
88
89 !-----------------------------------------------------------------------------
90 !
91 !++ Private procedure
92 !
93
94 !-----------------------------------------------------------------------------
95 !
96 !++ Private parameters & variables
97 !
98
99contains
100 subroutine localmeshbase_init( this, lcdomID, refElem, ndim, myrank )
101
102 use scale_prc, only: prc_myrank
103 implicit none
104
105 class(localmeshbase), intent(inout) :: this
106 integer, intent(in) :: lcdomid
107 class(elementbase), intent(in), target :: refelem
108 integer, intent(in) :: ndim
109 integer, intent(in), optional :: myrank
110 !-----------------------------------------------------------------------------
111
112 this%lcdomID = lcdomid
113 this%refElem => refelem
114
115 if (present(myrank)) then
116 this%PRC_myrank = myrank
117 else
118 this%PRC_myrank = prc_myrank
119 end if
120
121 return
122 end subroutine localmeshbase_init
123
124 subroutine localmeshbase_final( this, is_generated )
125 implicit none
126
127 class(localmeshbase), intent(inout) :: this
128 logical, intent(in) :: is_generated
129 !-----------------------------------------------------------------------------
130
131 if ( is_generated ) then
132 deallocate( this%pos_ev, this%pos_en, this%normal_fn )
133 deallocate( this%sJ, this%J )
134
135 deallocate( this%Escale, this%Fscale )
136
137 deallocate( this%EToV, this%EToE, this%EToF )
138
139 if ( allocated(this%VMapM) ) then
140 deallocate( this%VMapM, this%VMapP, this%MapM, this%MapP )
141 end if
142 if ( allocated(this%VMapB) ) then
143 deallocate( this%BCType )
144 deallocate( this%VMapB, this%MapB )
145 end if
146 if ( allocated(this%G_ij) ) then
147 deallocate( this%G_ij, this%GIJ )
148 deallocate( this%Gsqrt )
149 end if
150 end if
151
152 return
153 end subroutine localmeshbase_final
154end module scale_localmesh_base
module FElib / Element / Base
module FElib / Mesh / Local, Base
integer, parameter, public bctype_interior
integer, parameter, public bctype_shoreline
integer, parameter, public bctype_neuman
subroutine, public localmeshbase_final(this, is_generated)
integer, parameter, public bctype_periodic
integer, parameter, public bctype_dirchlet
subroutine, public localmeshbase_init(this, lcdomid, refelem, ndim, myrank)
integer, parameter, public bctype_undefbc