40 type :: veccovariantcomp
41 type(MeshField2D),
pointer :: u1 => null()
42 type(MeshField2D),
pointer :: u2 => null()
48 type(veccovariantcomp),
allocatable :: vec_covariant_comp_ptrlist(:)
49 integer,
allocatable :: nnode_lcmeshallface(:)
51 procedure,
public :: init => meshfieldcommcubedspheredom2d_init
52 procedure,
public :: put => meshfieldcommcubedspheredom2d_put
53 procedure,
public :: get => meshfieldcommcubedspheredom2d_get
54 procedure,
public :: exchange => meshfieldcommcubedspheredom2d_exchange
55 procedure,
public :: setcovariantvec => meshfieldcommcubedspheredom2d_set_covariantvec
56 procedure,
public :: final => meshfieldcommcubedspheredom2d_final
68 private :: post_exchange_core
69 private :: push_localsendbuf
76 integer,
parameter :: comm_face_num = 4
79 subroutine meshfieldcommcubedspheredom2d_init( this, &
80 sfield_num, hvfield_num, htensorfield_num, mesh2d )
85 integer,
intent(in) :: sfield_num
86 integer,
intent(in) :: hvfield_num
87 integer,
intent(in) :: htensorfield_num
92 integer :: nnode_lcmeshface(comm_face_num,mesh2d%local_mesh_num)
96 lcmesh => mesh2d%lcmesh_list(1)
99 allocate( this%Nnode_LCMeshAllFace(mesh2d%LOCAL_MESH_NUM) )
100 do n=1, this%mesh2d%LOCAL_MESH_NUM
101 lcmesh => this%mesh2d%lcmesh_list(n)
102 nnode_lcmeshface(:,n) = (/ lcmesh%NeX, lcmesh%NeY, lcmesh%NeX, lcmesh%NeY /) * lcmesh%refElem2D%Nfp
103 this%Nnode_LCMeshAllFace(n) = sum(nnode_lcmeshface(:,n))
108 if (hvfield_num > 0)
then
109 allocate( this%vec_covariant_comp_ptrlist(hvfield_num) )
113 end subroutine meshfieldcommcubedspheredom2d_init
115 subroutine meshfieldcommcubedspheredom2d_final( this )
121 deallocate( this%Nnode_LCMeshAllFace )
123 if ( this%hvfield_num > 0 )
then
124 deallocate( this%vec_covariant_comp_ptrlist )
130 end subroutine meshfieldcommcubedspheredom2d_final
132 subroutine meshfieldcommcubedspheredom2d_set_covariantvec( &
133 this, hvfield_ID, u1, u2 )
136 integer,
intent(in) :: hvfield_id
141 this%vec_covariant_comp_ptrlist(hvfield_id)%u1 => u1
142 this%vec_covariant_comp_ptrlist(hvfield_id)%u2 => u2
145 end subroutine meshfieldcommcubedspheredom2d_set_covariantvec
147 subroutine meshfieldcommcubedspheredom2d_put(this, field_list, varid_s)
151 integer,
intent(in) :: varid_s
158 do i=1,
size(field_list)
159 do n=1, this%mesh%LOCAL_MESH_NUM
160 lcmesh => this%mesh2d%lcmesh_list(n)
162 this%send_buf(:,varid_s+i-1,n) )
167 end subroutine meshfieldcommcubedspheredom2d_put
169 subroutine meshfieldcommcubedspheredom2d_get(this, field_list, varid_s)
176 integer,
intent(in) :: varid_s
184 integer :: varid_vec_s
188 varnum =
size(field_list)
191 if ( this%call_wait_flag_sub_get )
then
193 call post_exchange_core( this )
198 do n=1, this%mesh2d%LOCAL_MESH_NUM
199 lcmesh => this%mesh2d%lcmesh_list(n)
201 field_list(i)%field2d%local(n)%val )
205 varid_e = varid_s + varnum - 1
206 if ( varid_e > this%sfield_num )
then
207 do i=1, this%hvfield_num
209 varid_vec_s = this%sfield_num + 2*i - 1
210 if ( varid_vec_s > varid_e )
exit
212 if (
associated(this%vec_covariant_comp_ptrlist(i)%u1 ) &
213 .and.
associated(this%vec_covariant_comp_ptrlist(i)%u2 ) )
then
215 do n=1, this%mesh2d%LOCAL_MESH_NUM
216 call set_boundary_data2d_u1u2( &
217 this%recv_buf(:,varid_vec_s,n), this%recv_buf(:,varid_vec_s+1,n), &
218 lcmesh%refElem2D, lcmesh, lcmesh%G_ij, &
219 this%vec_covariant_comp_ptrlist(i)%u1%local(n)%val, &
220 this%vec_covariant_comp_ptrlist(i)%u2%local(n)%val )
227 end subroutine meshfieldcommcubedspheredom2d_get
230 subroutine meshfieldcommcubedspheredom2d_exchange( this, do_wait )
242 logical,
intent(in),
optional :: do_wait
247 real(rp),
allocatable :: fpos2d(:,:)
248 real(rp),
allocatable :: lcfpos2d(:,:)
249 real(rp),
allocatable :: unity_fac(:)
250 real(rp),
allocatable :: tmp_svec2d(:,:)
259 do n=1, this%mesh%LOCAL_MESH_NUM
260 lcmesh => this%mesh2d%lcmesh_list(n)
261 elem => lcmesh%refElem2D
263 allocate( fpos2d(this%Nnode_LCMeshAllFace(n),2) )
264 call extract_boundary_data2d( lcmesh%pos_en(:,:,1), elem, lcmesh, fpos2d(:,1) )
265 call extract_boundary_data2d( lcmesh%pos_en(:,:,2), elem, lcmesh, fpos2d(:,2) )
268 do f=1, this%nfaces_comm
269 commdata => this%commdata_list(f,n)
270 call push_localsendbuf( commdata%send_buf(:,:), &
271 this%send_buf(:,:,n), commdata%s_faceID, this%is_f(f,n), &
272 commdata%Nnode_LCMeshFace, this%field_num_tot)
274 if ( commdata%s_panelID /= lcmesh%panelID )
then
275 if ( this%hvfield_num > 0 )
then
277 allocate( lcfpos2d(commdata%Nnode_LCMeshFace,2), unity_fac(commdata%Nnode_LCMeshFace) )
278 allocate( tmp_svec2d(commdata%Nnode_LCMeshFace,2) )
280 call push_localsendbuf( lcfpos2d, &
281 fpos2d, commdata%s_faceID, this%is_f(f,n), commdata%Nnode_LCMeshFace, 2 )
282 unity_fac(:) = 1.0_rp
284 ire = irs + commdata%Nnode_LCMeshFace - 1
286 do varid=this%sfield_num+1, this%field_num_tot-1,2
287 tmp_svec2d(:,1) = commdata%send_buf(:,varid )
288 tmp_svec2d(:,2) = commdata%send_buf(:,varid+1)
290 lcmesh%panelID, lcfpos2d(:,1), lcfpos2d(:,2), unity_fac, &
291 commdata%Nnode_LCMeshFace, &
292 tmp_svec2d(:,1), tmp_svec2d(:,2), &
293 commdata%send_buf(:,varid), commdata%send_buf(:,varid+1) )
296 deallocate( lcfpos2d, unity_fac, tmp_svec2d )
308 if ( .not. this%call_wait_flag_sub_get ) &
309 call post_exchange_core( this )
312 end subroutine meshfieldcommcubedspheredom2d_exchange
316 subroutine post_exchange_core( this )
327 real(rp),
allocatable :: fpos2d(:,:)
328 real(rp),
allocatable :: lcfpos2d(:,:)
329 real(rp),
allocatable :: unity_fac(:)
338 do n=1, this%mesh%LOCAL_MESH_NUM
339 lcmesh => this%mesh2d%lcmesh_list(n)
340 elem => lcmesh%refElem2D
342 allocate( fpos2d(this%Nnode_LCMeshAllFace(n),2) )
343 call extract_boundary_data2d( lcmesh%pos_en(:,:,1), elem, lcmesh, fpos2d(:,1) )
344 call extract_boundary_data2d( lcmesh%pos_en(:,:,2), elem, lcmesh, fpos2d(:,2) )
347 do f=1, this%nfaces_comm
348 commdata => this%commdata_list(f,n)
349 ire = irs + commdata%Nnode_LCMeshFace - 1
351 if ( commdata%s_panelID /= lcmesh%panelID )
then
352 if ( this%hvfield_num > 0 )
then
354 allocate( lcfpos2d(commdata%Nnode_LCMeshFace,2), unity_fac(commdata%Nnode_LCMeshFace) )
356 call push_localsendbuf( lcfpos2d, &
357 fpos2d, f, this%is_f(f,n), commdata%Nnode_LCMeshFace, 2 )
358 unity_fac(:) = 1.0_rp
360 do varid=this%sfield_num+1, this%field_num_tot-1, 2
362 lcmesh%panelID, lcfpos2d(:,1), lcfpos2d(:,2), unity_fac(:), &
363 commdata%Nnode_LCMeshFace, &
364 commdata%recv_buf(:,varid), commdata%recv_buf(:,varid+1), &
365 this%recv_buf(irs:ire,varid,n), this%recv_buf(irs:ire,varid+1,n) )
367 deallocate( lcfpos2d, unity_fac )
378 end subroutine post_exchange_core
380 subroutine push_localsendbuf( lc_send_buf, send_buf, s_faceID, is, Nnode_LCMeshFace, var_num )
383 integer,
intent(in) :: var_num
384 integer,
intent(in) :: nnode_lcmeshface
385 real(rp),
intent(inout) :: lc_send_buf(nnode_lcmeshface,var_num)
387 integer,
intent(in) :: s_faceid, is
392 ie = is + nnode_lcmeshface - 1
393 if ( s_faceid > 0 )
then
394 lc_send_buf(:,:) = send_buf(is:ie,:)
396 lc_send_buf(:,:) = send_buf(ie:is:-1,:)
400 end subroutine push_localsendbuf
402 subroutine extract_boundary_data2d( var, refElem, mesh, buf )
407 real(dp),
intent(in) :: var(refelem%np * mesh%ne)
408 real(dp),
intent(inout) :: buf(refelem%nfp * mesh%nex * 4)
411 buf(:) = var(mesh%VmapB(:))
414 end subroutine extract_boundary_data2d
416 subroutine set_boundary_data2d_u1u2( buf_U, buf_V, &
424 real(dp),
intent(in) :: buf_u(elem%nfp * mesh%nex * 4)
425 real(dp),
intent(in) :: buf_v(elem%nfp * mesh%nex * 4)
426 real(dp),
intent(in) :: g_ij(elem%np * mesh%ne,2,2)
427 real(dp),
intent(inout) :: u1(elem%np * mesh%nea)
428 real(dp),
intent(inout) :: u2(elem%np * mesh%nea)
431 u1(elem%Np*mesh%NeE+1:elem%Np*mesh%NeE+
size(buf_u)) &
432 = g_ij(mesh%VmapB,1,1) * buf_u(:) + g_ij(mesh%VmapB,1,2) * buf_v(:)
433 u2(elem%Np*mesh%NeE+1:elem%Np*mesh%NeE+
size(buf_u)) &
434 = g_ij(mesh%VmapB,2,1) * buf_u(:) + g_ij(mesh%VmapB,2,2) * buf_v(:)
437 end subroutine set_boundary_data2d_u1u2
Module common / Coordinate conversion with cubed-sphere projection.
subroutine, public cubedspherecoordcnv_cs2lonlatvec(panelid, alpha, beta, gam, np, vecalpha, vecbeta, veclon, veclat, lat)
Covert the components of a vector in local coordinates with an equiangular gnomonic cubed-sphere proj...
subroutine, public cubedspherecoordcnv_lonlat2csvec(panelid, alpha, beta, gam, np, veclon, veclat, vecalpha, vecbeta, lat)
Covert the components of a vector in longitude and latitude coordinates to those in local coordinates...
module FElib / Element / Base
module FElib / Mesh / Local 2D
module FElib / Mesh / Cubed-sphere 2D domain
module FElib / Data / base
module FElib / Data / Communication base
subroutine, public meshfieldcommbase_extract_bounddata(var, refelem, mesh, buf)
Extract halo data from data array with MeshField object and set it to the recieving buffer.
subroutine, public meshfieldcommbase_final(this)
Finalize a base object to manage data communication of fields.
subroutine, public meshfieldcommbase_wait_core(this, commdata_list, field_list, dim, varid_s, lcmesh_list)
Wait data communication and move tmp data of LocalMeshCommData object to a recv buffer.
subroutine, public meshfieldcommbase_set_bounddata(buf, refelem, mesh, var)
Extract halo data from the recieving buffer and set it to data array with MeshField object.
subroutine, public meshfieldcommbase_init(this, sfield_num, hvfield_num, htensorfield_num, bufsize_per_field, comm_face_num, nnode_lcmeshface, mesh)
Initialize a base object to manage data communication of fields.
subroutine, public meshfieldcommbase_exchange_core(this, commdata_list, do_wait)
Exchange halo data.
module FElib / Data / Communication in 2D cubed-sphere domain
integer bufsize_per_field
Derived type representing a 2D reference element.
Derived type representing an arbitrary finite element.
Derived type representing a field with 2D mesh.
Derived type to manage data communication at a face between adjacent local meshes.
Base derived type to manage data communication.
Container to save a pointer of MeshField(1D, 2D, 3D) object.
Base derived type to manage data communication with 2D cubed-sphere domain.