FE-Project
Loading...
Searching...
No Matches
scale_localmesh_3d.F90
Go to the documentation of this file.
1!-------------------------------------------------------------------------------
2!> module FElib / Mesh / Local 3D
3!!
4!! @par Description
5!! Module to manage 3D local mesh for element-based methods
6!!
7!! @author Yuta Kawai, Team SCALE
8!<
9#include "scaleFElib.h"
11
12 !-----------------------------------------------------------------------------
13 !
14 !++ used modules
15 !
16 use scale_precision
17 use scale_localmesh_base, only: &
21
22 !-----------------------------------------------------------------------------
23 implicit none
24 private
25
26 !-----------------------------------------------------------------------------
27 !
28 !++ Public type & procedure
29 !
30
31 !> Derived type to manage a local 3D computational domain
32 type, extends(localmeshbase), public :: localmesh3d
33
34 type(elementbase3d), pointer :: refelem3d !< Pointer to a object of 3D reference element
35 real(rp) :: xmin !< Minimum x-coordinate value of the local computational domain
36 real(rp) :: xmax !< Maximum x-coordinate value of the local computational domain
37 real(rp) :: ymin !< Minimum y-coordinate value of the local computational domain
38 real(rp) :: ymax !< Maximum y-coordinate value of the local computational domain
39 real(rp) :: zmin !< Minimum z-coordinate value of the local computational domain
40 real(rp) :: zmax !< Maximum z-coordinate value of the local computational domain
41
42 integer :: nex !< Number of finite elements for the x-direction in the local computational domain
43 integer :: ney !< Number of finite elements for the y-direction in the local computational domain
44 integer :: nez !< Number of finite elements for the z-direction in the local computational domain
45 integer :: ne2d !< NeX * NeY
46 integer :: ne2da !< NeX * NeY + halo size
47
48 real(rp), allocatable :: sz(:,:)
49 real(rp), allocatable :: zs(:,:)
50 real(rp), allocatable :: gi3(:,:,:) !< The contravariant component of metric tensor with vertical general coordinate
51 real(rp), allocatable :: gsqrth(:,:) !< The Jacobian of horizontal transformation in the computational coordinate
52 real(rp), allocatable :: zlev(:,:) !< z-coordinates (actual level)
53 real(rp), allocatable :: gam(:,:) !< Factor for approximation with spherical shell domain (= r/a)
54
55 class(localmesh2d), pointer :: lcmesh2d
56 real(rp), allocatable :: lon2d(:,:) !< Longitude coordinate with 2D mesh
57 real(rp), allocatable :: lat2d(:,:) !< Latitude coordinate with 2D mesh
58
59 integer, allocatable :: emap3dto2d(:) !< Array to map 3D element ID to 2D element ID
60 contains
61 procedure :: setlocalmesh2d => localmesh3d_setlocalmesh2d
62 procedure :: getvmapz1d => localmesh3d_getvmapz1d
63 procedure :: getvmapz3d => localmesh3d_getvmapz3d
64 end type localmesh3d
65
67
68 !-----------------------------------------------------------------------------
69 !
70 !++ Public parameters & variables
71 !
72
73 !-----------------------------------------------------------------------------
74 !
75 !++ Private procedure
76 !
77
78 !-----------------------------------------------------------------------------
79 !
80 !++ Private parameters & variables
81 !
82
83contains
84!> Initialize an object to manage a 3D local computational domain
85!OCL SERIAL
86 subroutine localmesh3d_init( this, &
87 lcdomID, refElem, myrank )
88 implicit none
89
90 class(localmesh3d), intent(inout) :: this
91 integer, intent(in) :: lcdomid
92 class(elementbase3d), intent(in), target :: refelem
93 integer, intent(in), optional :: myrank
94 !-------------------------------------------------
95
96 this%refElem3D => refelem
97 nullify( this%lcmesh2D )
98
99 call localmeshbase_init( this, lcdomid, refelem, 3, myrank )
100
101 return
102 end subroutine localmesh3d_init
103
104!> Finalize an object to manage a 3D local computational domain
105!OCL SERIAL
106 subroutine localmesh3d_final( this, is_generated )
107 implicit none
108 type(localmesh3d), intent(inout) :: this
109 logical, intent(in) :: is_generated
110 !-------------------------------------------------
111
112 call localmeshbase_final( this, is_generated )
113 if (is_generated) then
114 deallocate( this%zS, this%Sz )
115 deallocate( this%GI3, this%GsqrtH )
116 deallocate( this%zlev )
117 deallocate( this%gam )
118 deallocate( this%lon2D, this%lat2D )
119 deallocate( this%EMap3Dto2D )
120 end if
121
122 return
123 end subroutine localmesh3d_final
124
125!> Set a pointer to an object to mange 2D local computational domain
126!OCL SERIAL
127 subroutine localmesh3d_setlocalmesh2d( this, lcmesh2D )
128 implicit none
129 class(localmesh3d), intent(inout) :: this
130 class(localmesh2d), intent(in), target :: lcmesh2d
131 !-------------------------------------------------
132
133 this%lcmesh2D => lcmesh2d
134
135 return
136 end subroutine localmesh3d_setlocalmesh2d
137
138!> Get arrays for vertical mapping node index with vertical element boundary to that with DG data
139!OCL SERIAL
140 subroutine localmesh3d_getvmapz1d( this, vmapM, vmapP )
141
142 implicit none
143 class(localmesh3d), intent(in), target :: this
144 integer, intent(out) :: vmapm(this%refelem3d%nfptot,this%nez)
145 integer, intent(out) :: vmapp(this%refelem3d%nfptot,this%nez)
146
147 integer :: ke_z
148 integer :: f
149 integer :: vs, ve
150
151 class(elementbase3d), pointer :: elem
152 !------------------------------
153
154 elem => this%refElem3D
155
156 do ke_z=1, this%NeZ
157 do f=1, elem%Nfaces_h
158 vs = 1 + (f-1)*elem%Nfp_h
159 ve = vs + elem%Nfp_h - 1
160 vmapm(vs:ve,ke_z) = elem%Fmask_h(:,f) + (ke_z-1)*elem%Np
161 end do
162 do f=1, elem%Nfaces_v
163 vs = elem%Nfp_h*elem%Nfaces_h + 1 + (f-1)*elem%Nfp_v
164 ve = vs + elem%Nfp_v - 1
165 vmapm(vs:ve,ke_z) = elem%Fmask_v(:,f) + (ke_z-1)*elem%Np
166 end do
167 vmapp(:,ke_z) = vmapm(:,ke_z)
168 end do
169
170 do ke_z=1, this%NeZ
171 vs = elem%Nfp_h*elem%Nfaces_h + 1
172 ve = vs + elem%Nfp_v - 1
173 if (ke_z > 1) &
174 vmapp(vs:ve,ke_z) = elem%Fmask_v(:,2) + (ke_z-2)*elem%Np
175
176 vs = elem%Nfp_h*elem%Nfaces_h + elem%Nfp_v + 1
177 ve = vs + elem%Nfp_v - 1
178 if (ke_z < this%NeZ) &
179 vmapp(vs:ve,ke_z) = elem%Fmask_v(:,1) + ke_z*elem%Np
180 end do
181
182 return
183 end subroutine localmesh3d_getvmapz1d
184
185!> Get arrays for vertical mapping node index with vertical element boundary to that with DG data
186!OCL SERIAL
187 subroutine localmesh3d_getvmapz3d( this, vmapM, vmapP )
188
189 implicit none
190 class(localmesh3d), intent(in), target :: this
191 integer, intent(out) :: vmapm(this%refelem3d%nfptot,this%ne)
192 integer, intent(out) :: vmapp(this%refelem3d%nfptot,this%ne)
193
194 integer :: ke, ke_xy, ke_z
195 integer :: ke_neigh
196 integer :: f
197 integer :: vs, ve
198
199 class(elementbase3d), pointer :: elem
200 !------------------------------
201
202 elem => this%refElem3D
203
204 !$omp parallel do collapse(2) private(ke,ke_xy,ke_z,ke_neigh,f,vs,ve)
205 do ke_z=1, this%NeZ
206 do ke_xy=1, this%NeX * this%NeY
207 ke = ke_xy + (ke_z-1) * this%NeX * this%NeY
208
209 do f=1, elem%Nfaces_h
210 vs = 1 + (f-1)*elem%Nfp_h
211 ve = vs + elem%Nfp_h - 1
212 vmapm(vs:ve,ke) = elem%Fmask_h(:,f) + (ke-1)*elem%Np
213 end do
214 do f=1, elem%Nfaces_v
215 vs = elem%Nfp_h*elem%Nfaces_h + 1 + (f-1)*elem%Nfp_v
216 ve = vs + elem%Nfp_v - 1
217 vmapm(vs:ve,ke) = elem%Fmask_v(:,f) + (ke-1)*elem%Np
218 end do
219 vmapp(:,ke) = vmapm(:,ke)
220
221 !--
222 vs = elem%Nfp_h*elem%Nfaces_h + 1
223 ve = vs + elem%Nfp_v - 1
224 if (ke_z > 1) then
225 ke_neigh = ke_xy + (ke_z-2) * this%NeX * this%NeY
226 vmapp(vs:ve,ke) = elem%Fmask_v(:,2) + (ke_neigh-1)*elem%Np
227 end if
228 vs = elem%Nfp_h*elem%Nfaces_h + elem%Nfp_v + 1
229 ve = vs + elem%Nfp_v - 1
230 if (ke_z < this%NeZ) then
231 ke_neigh = ke_xy + ke_z * this%NeX * this%NeY
232 vmapp(vs:ve,ke) = elem%Fmask_v(:,1) + (ke_neigh-1)*elem%Np
233 end if
234 end do
235 end do
236
237 return
238 end subroutine localmesh3d_getvmapz3d
239
240end module scale_localmesh_3d
module FElib / Element / Base
module FElib / Mesh / Local 2D
module FElib / Mesh / Local 3D
subroutine, public localmesh3d_final(this, is_generated)
Finalize an object to manage a 3D local computational domain.
subroutine, public localmesh3d_init(this, lcdomid, refelem, myrank)
Initialize an object to manage a 3D local computational domain.
module FElib / Mesh / Local, Base
subroutine, public localmeshbase_final(this, is_generated)
Finalize an object to manage a local computational mesh.
subroutine, public localmeshbase_init(this, lcdomid, refelem, ndim, myrank)
Setup an object to manage a local computational mesh.
Derived type representing a 3D reference element.
Derived type representing an arbitrary finite element.
Derived type to manage a local 3D computational domain.
Derived type to manage a local computational domain (base type)