FE-Project
Loading...
Searching...
No Matches
scale_localmesh_3d.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_localmesh_base, only: &
21
22 !-----------------------------------------------------------------------------
23 implicit none
24 private
25
26 !-----------------------------------------------------------------------------
27 !
28 !++ Public type & procedure
29 !
30 type, extends(localmeshbase), public :: localmesh3d
31
32 type(elementbase3d), pointer :: refelem3d
33 real(rp) :: xmin, xmax
34 real(rp) :: ymin, ymax
35 real(rp) :: zmin, zmax
36 integer :: nex
37 integer :: ney
38 integer :: nez
39 integer :: ne2d
40 integer :: ne2da
41
42 real(rp), allocatable :: sz(:,:)
43 real(rp), allocatable :: zs(:,:)
44 real(rp), allocatable :: gi3(:,:,:)
45 real(rp), allocatable :: gsqrth(:,:)
46 real(rp), allocatable :: zlev(:,:)
47 real(rp), allocatable :: gam(:,:)
48
49 class(localmesh2d), pointer :: lcmesh2d
50 real(rp), allocatable :: lon2d(:,:)
51 real(rp), allocatable :: lat2d(:,:)
52
53 integer, allocatable :: emap3dto2d(:)
54 contains
55 procedure :: setlocalmesh2d => localmesh3d_setlocalmesh2d
56 procedure :: getvmapz1d => localmesh3d_getvmapz1d
57 end type localmesh3d
58
60
61 !-----------------------------------------------------------------------------
62 !
63 !++ Public parameters & variables
64 !
65
66 !-----------------------------------------------------------------------------
67 !
68 !++ Private procedure
69 !
70
71 !-----------------------------------------------------------------------------
72 !
73 !++ Private parameters & variables
74 !
75
76contains
77!OCL SERIAL
78 subroutine localmesh3d_init( this, &
79 lcdomID, refElem, myrank )
80 implicit none
81
82 class(localmesh3d), intent(inout) :: this
83 integer, intent(in) :: lcdomid
84 class(elementbase3d), intent(in), target :: refelem
85 integer, intent(in), optional :: myrank
86 !-------------------------------------------------
87
88 this%refElem3D => refelem
89 nullify( this%lcmesh2D )
90
91 call localmeshbase_init( this, lcdomid, refelem, 3, myrank )
92
93 return
94 end subroutine localmesh3d_init
95
96!OCL SERIAL
97 subroutine localmesh3d_final( this, is_generated )
98 implicit none
99 type(localmesh3d), intent(inout) :: this
100 logical, intent(in) :: is_generated
101 !-------------------------------------------------
102
103 call localmeshbase_final( this, is_generated )
104 if (is_generated) then
105 deallocate( this%zS, this%Sz )
106 deallocate( this%GI3, this%GsqrtH )
107 deallocate( this%zlev )
108 deallocate( this%gam )
109 deallocate( this%lon2D, this%lat2D )
110 deallocate( this%EMap3Dto2D )
111 end if
112
113 return
114 end subroutine localmesh3d_final
115
116!OCL SERIAL
117 subroutine localmesh3d_setlocalmesh2d( this, lcmesh2D )
118 implicit none
119 class(localmesh3d), intent(inout) :: this
120 class(localmesh2d), intent(in), target :: lcmesh2d
121 !-------------------------------------------------
122
123 this%lcmesh2D => lcmesh2d
124
125 return
126 end subroutine localmesh3d_setlocalmesh2d
127
128!OCL SERIAL
129 subroutine localmesh3d_getvmapz1d( this, vmapM, vmapP )
130
131 implicit none
132 class(localmesh3d), intent(in), target :: this
133 integer, intent(out) :: vmapm(this%refelem3d%nfptot,this%nez)
134 integer, intent(out) :: vmapp(this%refelem3d%nfptot,this%nez)
135
136 integer :: ke_z
137 integer :: f
138 integer :: vs, ve
139
140 class(elementbase3d), pointer :: elem
141 !------------------------------
142
143 elem => this%refElem3D
144
145 do ke_z=1, this%NeZ
146 do f=1, elem%Nfaces_h
147 vs = 1 + (f-1)*elem%Nfp_h
148 ve = vs + elem%Nfp_h - 1
149 vmapm(vs:ve,ke_z) = elem%Fmask_h(:,f) + (ke_z-1)*elem%Np
150 end do
151 do f=1, elem%Nfaces_v
152 vs = elem%Nfp_h*elem%Nfaces_h + 1 + (f-1)*elem%Nfp_v
153 ve = vs + elem%Nfp_v - 1
154 vmapm(vs:ve,ke_z) = elem%Fmask_v(:,f) + (ke_z-1)*elem%Np
155 end do
156 vmapp(:,ke_z) = vmapm(:,ke_z)
157 end do
158
159 do ke_z=1, this%NeZ
160 vs = elem%Nfp_h*elem%Nfaces_h + 1
161 ve = vs + elem%Nfp_v - 1
162 if (ke_z > 1) &
163 vmapp(vs:ve,ke_z) = elem%Fmask_v(:,2) + (ke_z-2)*elem%Np
164
165 vs = elem%Nfp_h*elem%Nfaces_h + elem%Nfp_v + 1
166 ve = vs + elem%Nfp_v - 1
167 if (ke_z < this%NeZ) &
168 vmapp(vs:ve,ke_z) = elem%Fmask_v(:,1) + ke_z*elem%Np
169 end do
170
171 return
172 end subroutine localmesh3d_getvmapz1d
173
174end 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)
subroutine, public localmesh3d_init(this, lcdomid, refelem, myrank)
module FElib / Mesh / Local, Base
subroutine, public localmeshbase_final(this, is_generated)
subroutine, public localmeshbase_init(this, lcdomid, refelem, ndim, myrank)