48 real(rp),
allocatable :: sz(:,:)
49 real(rp),
allocatable :: zs(:,:)
50 real(rp),
allocatable :: gi3(:,:,:)
51 real(rp),
allocatable :: gsqrth(:,:)
52 real(rp),
allocatable :: zlev(:,:)
53 real(rp),
allocatable :: gam(:,:)
56 real(rp),
allocatable :: lon2d(:,:)
57 real(rp),
allocatable :: lat2d(:,:)
59 integer,
allocatable :: emap3dto2d(:)
61 procedure :: setlocalmesh2d => localmesh3d_setlocalmesh2d
62 procedure :: getvmapz1d => localmesh3d_getvmapz1d
63 procedure :: getvmapz3d => localmesh3d_getvmapz3d
87 lcdomID, refElem, myrank )
91 integer,
intent(in) :: lcdomid
93 integer,
intent(in),
optional :: myrank
96 this%refElem3D => refelem
97 nullify( this%lcmesh2D )
109 logical,
intent(in) :: 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 )
127 subroutine localmesh3d_setlocalmesh2d( this, lcmesh2D )
133 this%lcmesh2D => lcmesh2d
136 end subroutine localmesh3d_setlocalmesh2d
140 subroutine localmesh3d_getvmapz1d( this, vmapM, vmapP )
144 integer,
intent(out) :: vmapm(this%refelem3d%nfptot,this%nez)
145 integer,
intent(out) :: vmapp(this%refelem3d%nfptot,this%nez)
154 elem => this%refElem3D
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
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
167 vmapp(:,ke_z) = vmapm(:,ke_z)
171 vs = elem%Nfp_h*elem%Nfaces_h + 1
172 ve = vs + elem%Nfp_v - 1
174 vmapp(vs:ve,ke_z) = elem%Fmask_v(:,2) + (ke_z-2)*elem%Np
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
183 end subroutine localmesh3d_getvmapz1d
187 subroutine localmesh3d_getvmapz3d( this, vmapM, vmapP )
191 integer,
intent(out) :: vmapm(this%refelem3d%nfptot,this%ne)
192 integer,
intent(out) :: vmapp(this%refelem3d%nfptot,this%ne)
194 integer :: ke, ke_xy, ke_z
202 elem => this%refElem3D
206 do ke_xy=1, this%NeX * this%NeY
207 ke = ke_xy + (ke_z-1) * this%NeX * this%NeY
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
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
219 vmapp(:,ke) = vmapm(:,ke)
222 vs = elem%Nfp_h*elem%Nfaces_h + 1
223 ve = vs + elem%Nfp_v - 1
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
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
238 end subroutine localmesh3d_getvmapz3d
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)