39 type :: veccovariantcomp
40 type(MeshField2D),
pointer :: u1 => null()
41 type(MeshField2D),
pointer :: u2 => null()
46 type(veccovariantcomp),
allocatable :: vec_covariant_comp_ptrlist(:)
47 integer,
allocatable :: nnode_lcmeshallface(:)
49 procedure,
public :: init => meshfieldcommcubedspheredom2d_init
50 procedure,
public :: put => meshfieldcommcubedspheredom2d_put
51 procedure,
public :: get => meshfieldcommcubedspheredom2d_get
52 procedure,
public :: exchange => meshfieldcommcubedspheredom2d_exchange
53 procedure,
public :: setcovariantvec => meshfieldcommcubedspheredom2d_set_covariantvec
54 procedure,
public :: final => meshfieldcommcubedspheredom2d_final
66 private :: post_exchange_core
67 private :: push_localsendbuf
74 integer,
parameter :: comm_face_num = 4
77 subroutine meshfieldcommcubedspheredom2d_init( this, &
78 sfield_num, hvfield_num, htensorfield_num, mesh2d )
83 integer,
intent(in) :: sfield_num
84 integer,
intent(in) :: hvfield_num
85 integer,
intent(in) :: htensorfield_num
90 integer :: nnode_lcmeshface(comm_face_num,mesh2d%local_mesh_num)
94 lcmesh => mesh2d%lcmesh_list(1)
97 allocate( this%Nnode_LCMeshAllFace(mesh2d%LOCAL_MESH_NUM) )
98 do n=1, this%mesh2d%LOCAL_MESH_NUM
99 lcmesh => this%mesh2d%lcmesh_list(n)
100 nnode_lcmeshface(:,n) = (/ lcmesh%NeX, lcmesh%NeY, lcmesh%NeX, lcmesh%NeY /) * lcmesh%refElem2D%Nfp
101 this%Nnode_LCMeshAllFace(n) = sum(nnode_lcmeshface(:,n))
106 if (hvfield_num > 0)
then
107 allocate( this%vec_covariant_comp_ptrlist(hvfield_num) )
111 end subroutine meshfieldcommcubedspheredom2d_init
113 subroutine meshfieldcommcubedspheredom2d_final( this )
119 deallocate( this%Nnode_LCMeshAllFace )
121 if ( this%hvfield_num > 0 )
then
122 deallocate( this%vec_covariant_comp_ptrlist )
128 end subroutine meshfieldcommcubedspheredom2d_final
130 subroutine meshfieldcommcubedspheredom2d_set_covariantvec( &
131 this, hvfield_ID, u1, u2 )
134 integer,
intent(in) :: hvfield_id
139 this%vec_covariant_comp_ptrlist(hvfield_id)%u1 => u1
140 this%vec_covariant_comp_ptrlist(hvfield_id)%u2 => u2
143 end subroutine meshfieldcommcubedspheredom2d_set_covariantvec
145 subroutine meshfieldcommcubedspheredom2d_put(this, field_list, varid_s)
149 integer,
intent(in) :: varid_s
156 do i=1,
size(field_list)
157 do n=1, this%mesh%LOCAL_MESH_NUM
158 lcmesh => this%mesh2d%lcmesh_list(n)
160 this%send_buf(:,varid_s+i-1,n) )
165 end subroutine meshfieldcommcubedspheredom2d_put
167 subroutine meshfieldcommcubedspheredom2d_get(this, field_list, varid_s)
174 integer,
intent(in) :: varid_s
182 integer :: varid_vec_s
186 varnum =
size(field_list)
189 if ( this%call_wait_flag_sub_get )
then
190 call post_exchange_core( this )
196 do n=1, this%mesh2d%LOCAL_MESH_NUM
197 lcmesh => this%mesh2d%lcmesh_list(n)
199 field_list(i)%field2d%local(n)%val )
203 varid_e = varid_s + varnum - 1
204 if ( varid_e > this%sfield_num )
then
205 do i=1, this%hvfield_num
207 varid_vec_s = this%sfield_num + 2*i - 1
208 if ( varid_vec_s > varid_e )
exit
210 if (
associated(this%vec_covariant_comp_ptrlist(i)%u1 ) &
211 .and.
associated(this%vec_covariant_comp_ptrlist(i)%u2 ) )
then
213 do n=1, this%mesh2d%LOCAL_MESH_NUM
214 call set_boundary_data2d_u1u2( &
215 this%recv_buf(:,varid_vec_s,n), this%recv_buf(:,varid_vec_s+1,n), &
216 lcmesh%refElem2D, lcmesh, lcmesh%G_ij, &
217 this%vec_covariant_comp_ptrlist(i)%u1%local(n)%val, &
218 this%vec_covariant_comp_ptrlist(i)%u2%local(n)%val )
225 end subroutine meshfieldcommcubedspheredom2d_get
228 subroutine meshfieldcommcubedspheredom2d_exchange( this, do_wait )
240 logical,
intent(in),
optional :: do_wait
245 real(rp),
allocatable :: fpos2d(:,:)
246 real(rp),
allocatable :: lcfpos2d(:,:)
247 real(rp),
allocatable :: unity_fac(:)
248 real(rp),
allocatable :: tmp_svec2d(:,:)
257 do n=1, this%mesh%LOCAL_MESH_NUM
258 lcmesh => this%mesh2d%lcmesh_list(n)
259 elem => lcmesh%refElem2D
261 allocate( fpos2d(this%Nnode_LCMeshAllFace(n),2) )
262 call extract_boundary_data2d( lcmesh%pos_en(:,:,1), elem, lcmesh, fpos2d(:,1) )
263 call extract_boundary_data2d( lcmesh%pos_en(:,:,2), elem, lcmesh, fpos2d(:,2) )
266 do f=1, this%nfaces_comm
267 commdata => this%commdata_list(f,n)
268 call push_localsendbuf( commdata%send_buf(:,:), &
269 this%send_buf(:,:,n), commdata%s_faceID, this%is_f(f,n), &
270 commdata%Nnode_LCMeshFace, this%field_num_tot)
272 if ( commdata%s_panelID /= lcmesh%panelID )
then
273 if ( this%hvfield_num > 0 )
then
275 allocate( lcfpos2d(commdata%Nnode_LCMeshFace,2), unity_fac(commdata%Nnode_LCMeshFace) )
276 allocate( tmp_svec2d(commdata%Nnode_LCMeshFace,2) )
278 call push_localsendbuf( lcfpos2d, &
279 fpos2d, commdata%s_faceID, this%is_f(f,n), commdata%Nnode_LCMeshFace, 2 )
280 unity_fac(:) = 1.0_rp
282 ire = irs + commdata%Nnode_LCMeshFace - 1
284 do varid=this%sfield_num+1, this%field_num_tot-1,2
285 tmp_svec2d(:,1) = commdata%send_buf(:,varid )
286 tmp_svec2d(:,2) = commdata%send_buf(:,varid+1)
288 lcmesh%panelID, lcfpos2d(:,1), lcfpos2d(:,2), unity_fac, &
289 commdata%Nnode_LCMeshFace, &
290 tmp_svec2d(:,1), tmp_svec2d(:,2), &
291 commdata%send_buf(:,varid), commdata%send_buf(:,varid+1) )
294 deallocate( lcfpos2d, unity_fac, tmp_svec2d )
306 if ( .not. this%call_wait_flag_sub_get ) &
307 call post_exchange_core( this )
310 end subroutine meshfieldcommcubedspheredom2d_exchange
314 subroutine post_exchange_core( this )
325 real(rp),
allocatable :: fpos2d(:,:)
326 real(rp),
allocatable :: lcfpos2d(:,:)
327 real(rp),
allocatable :: unity_fac(:)
336 do n=1, this%mesh%LOCAL_MESH_NUM
337 lcmesh => this%mesh2d%lcmesh_list(n)
338 elem => lcmesh%refElem2D
340 allocate( fpos2d(this%Nnode_LCMeshAllFace(n),2) )
341 call extract_boundary_data2d( lcmesh%pos_en(:,:,1), elem, lcmesh, fpos2d(:,1) )
342 call extract_boundary_data2d( lcmesh%pos_en(:,:,2), elem, lcmesh, fpos2d(:,2) )
345 do f=1, this%nfaces_comm
346 commdata => this%commdata_list(f,n)
347 ire = irs + commdata%Nnode_LCMeshFace - 1
349 if ( commdata%s_panelID /= lcmesh%panelID )
then
350 if ( this%hvfield_num > 0 )
then
352 allocate( lcfpos2d(commdata%Nnode_LCMeshFace,2), unity_fac(commdata%Nnode_LCMeshFace) )
354 call push_localsendbuf( lcfpos2d, &
355 fpos2d, f, this%is_f(f,n), commdata%Nnode_LCMeshFace, 2 )
356 unity_fac(:) = 1.0_rp
358 do varid=this%sfield_num+1, this%field_num_tot-1, 2
360 lcmesh%panelID, lcfpos2d(:,1), lcfpos2d(:,2), unity_fac(:), &
361 commdata%Nnode_LCMeshFace, &
362 commdata%recv_buf(:,varid), commdata%recv_buf(:,varid+1), &
363 this%recv_buf(irs:ire,varid,n), this%recv_buf(irs:ire,varid+1,n) )
365 deallocate( lcfpos2d, unity_fac )
376 end subroutine post_exchange_core
378 subroutine push_localsendbuf( lc_send_buf, send_buf, s_faceID, is, Nnode_LCMeshFace, var_num )
381 integer,
intent(in) :: var_num
382 integer,
intent(in) :: nnode_lcmeshface
383 real(rp),
intent(inout) :: lc_send_buf(nnode_lcmeshface,var_num)
385 integer,
intent(in) :: s_faceid, is
390 ie = is + nnode_lcmeshface - 1
391 if ( s_faceid > 0 )
then
392 lc_send_buf(:,:) = send_buf(is:ie,:)
394 lc_send_buf(:,:) = send_buf(ie:is:-1,:)
398 end subroutine push_localsendbuf
400 subroutine extract_boundary_data2d( var, refElem, mesh, buf )
405 real(dp),
intent(in) :: var(refelem%np * mesh%ne)
406 real(dp),
intent(inout) :: buf(refelem%nfp * mesh%nex * 4)
409 buf(:) = var(mesh%VmapB(:))
412 end subroutine extract_boundary_data2d
414 subroutine set_boundary_data2d_u1u2( buf_U, buf_V, &
422 real(dp),
intent(in) :: buf_u(elem%nfp * mesh%nex * 4)
423 real(dp),
intent(in) :: buf_v(elem%nfp * mesh%nex * 4)
424 real(dp),
intent(in) :: g_ij(elem%np * mesh%ne,2,2)
425 real(dp),
intent(inout) :: u1(elem%np * mesh%nea)
426 real(dp),
intent(inout) :: u2(elem%np * mesh%nea)
429 u1(elem%Np*mesh%NeE+1:elem%Np*mesh%NeE+
size(buf_u)) &
430 = g_ij(mesh%VmapB,1,1) * buf_u(:) + g_ij(mesh%VmapB,1,2) * buf_v(:)
431 u2(elem%Np*mesh%NeE+1:elem%Np*mesh%NeE+
size(buf_u)) &
432 = g_ij(mesh%VmapB,2,1) * buf_u(:) + g_ij(mesh%VmapB,2,2) * buf_v(:)
435 end subroutine set_boundary_data2d_u1u2
437end module scale_meshfieldcomm_cubedspheredom2d
module common / Coordinate conversion with a cubed-sphere
subroutine, public cubedspherecoordcnv_cs2lonlatvec(panelid, alpha, beta, gam, np, vecalpha, vecbeta, veclon, veclat, lat)
subroutine, public cubedspherecoordcnv_lonlat2csvec(panelid, alpha, beta, gam, np, veclon, veclat, vecalpha, vecbeta, lat)
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)
subroutine, public meshfieldcommbase_final(this)
Finalize a base object to manage data communication of fields.
subroutine, public meshfieldcommbase_wait_core(this, commdata_list)
Wait data communication and move tmp data of LocalMeshCommData object to a recv buffer.
subroutine, public meshfieldcommbase_set_bounddata(buf, refelem, mesh, var)
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 to manage data communication between a face of local mesh.
Base derived type to manage data communication.