42 type :: veccovariantcomp
43 type(MeshField3D),
pointer :: u1 => null()
44 type(MeshField3D),
pointer :: u2 => null()
49 type(veccovariantcomp),
allocatable :: vec_covariant_comp_ptrlist(:)
50 integer,
allocatable :: nnode_lcmeshallface(:)
52 procedure,
public :: init => meshfieldcommcubedspheredom3d_init
53 procedure,
public :: put => meshfieldcommcubedspheredom3d_put
54 procedure,
public :: get => meshfieldcommcubedspheredom3d_get
55 procedure,
public :: exchange => meshfieldcommcubedspheredom3d_exchange
56 procedure,
public :: setcovariantvec => meshfieldcommcubedspheredom3d_set_covariantvec
57 procedure,
public :: final => meshfieldcommcubedspheredom3d_final
69 private :: post_exchange_core
70 private :: push_localsendbuf
77 integer,
parameter :: comm_face_num = 6
80 subroutine meshfieldcommcubedspheredom3d_init( this, &
81 sfield_num, hvfield_num, htensorfield_num, mesh3d )
86 integer,
intent(in) :: sfield_num
87 integer,
intent(in) :: hvfield_num
88 integer,
intent(in) :: htensorfield_num
94 integer :: nnode_lcmeshface(comm_face_num,mesh3d%local_mesh_num)
98 lcmesh => mesh3d%lcmesh_list(1)
99 elem => lcmesh%refElem3D
102 + 2*lcmesh%NeX*lcmesh%NeY*elem%Nfp_v
104 allocate( this%Nnode_LCMeshAllFace(mesh3d%LOCAL_MESH_NUM) )
105 do n=1, this%mesh3d%LOCAL_MESH_NUM
106 lcmesh => this%mesh3d%lcmesh_list(n)
107 nnode_lcmeshface(:,n) = &
108 (/ lcmesh%NeX, lcmesh%NeY, lcmesh%NeX, lcmesh%NeY, 0, 0 /) * lcmesh%NeZ * lcmesh%refElem3D%Nfp_h &
109 + (/ 0, 0, 0, 0, 1, 1 /) * lcmesh%NeX*lcmesh%NeY * lcmesh%refElem3D%Nfp_v
110 this%Nnode_LCMeshAllFace(n) = sum(nnode_lcmeshface(:,n))
115 if (hvfield_num > 0)
then
116 allocate( this%vec_covariant_comp_ptrlist(hvfield_num) )
120 end subroutine meshfieldcommcubedspheredom3d_init
122 subroutine meshfieldcommcubedspheredom3d_final( this )
129 if ( this%hvfield_num > 0 )
then
130 deallocate( this%vec_covariant_comp_ptrlist )
136 end subroutine meshfieldcommcubedspheredom3d_final
138 subroutine meshfieldcommcubedspheredom3d_set_covariantvec( &
139 this, hvfield_ID, u1, u2 )
142 integer,
intent(in) :: hvfield_id
147 this%vec_covariant_comp_ptrlist(hvfield_id)%u1 => u1
148 this%vec_covariant_comp_ptrlist(hvfield_id)%u2 => u2
151 end subroutine meshfieldcommcubedspheredom3d_set_covariantvec
153 subroutine meshfieldcommcubedspheredom3d_put(this, field_list, varid_s)
157 integer,
intent(in) :: varid_s
164 do i=1,
size(field_list)
165 do n=1, this%mesh%LOCAL_MESH_NUM
166 lcmesh => this%mesh3d%lcmesh_list(n)
168 this%send_buf(:,varid_s+i-1,n) )
173 end subroutine meshfieldcommcubedspheredom3d_put
175 subroutine meshfieldcommcubedspheredom3d_get(this, field_list, varid_s)
182 integer,
intent(in) :: varid_s
193 integer :: varid_vec_s
195 real(rp),
allocatable :: g_ij(:,:,:,:)
198 varnum =
size(field_list)
201 if ( this%call_wait_flag_sub_get )
then
202 call post_exchange_core( this )
208 do n=1, this%mesh3d%LOCAL_MESH_NUM
209 lcmesh => this%mesh3d%lcmesh_list(n)
211 field_list(i)%field3d%local(n)%val )
215 varid_e = varid_s + varnum - 1
216 if ( varid_e > this%sfield_num )
then
217 do i=1, this%hvfield_num
219 varid_vec_s = this%sfield_num + 2*i - 1
220 if ( varid_vec_s > varid_e )
exit
222 if (
associated(this%vec_covariant_comp_ptrlist(i)%u1 ) &
223 .and.
associated(this%vec_covariant_comp_ptrlist(i)%u2 ) )
then
225 do n=1, this%mesh3d%LOCAL_MESH_NUM
226 lcmesh => this%mesh3d%lcmesh_list(n)
227 elem => lcmesh%refElem3D
229 allocate( g_ij(elem%Np,lcmesh%Ne,2,2) )
231 do ke=lcmesh%NeS, lcmesh%NeE
232 ke2d = lcmesh%EMap3Dto2D(ke)
233 g_ij(:,ke,1,1) = lcmesh%G_ij(elem%IndexH2Dto3D(:),ke2d,1,1)
234 g_ij(:,ke,2,1) = lcmesh%G_ij(elem%IndexH2Dto3D(:),ke2d,2,1)
235 g_ij(:,ke,1,2) = lcmesh%G_ij(elem%IndexH2Dto3D(:),ke2d,1,2)
236 g_ij(:,ke,2,2) = lcmesh%G_ij(elem%IndexH2Dto3D(:),ke2d,2,2)
239 call set_boundary_data3d_u1u2( &
240 this%recv_buf(:,varid_vec_s,n), this%recv_buf(:,varid_vec_s+1,n), &
241 lcmesh%refElem3D, lcmesh, g_ij(:,:,:,:), &
242 this%vec_covariant_comp_ptrlist(i)%u1%local(n)%val, &
243 this%vec_covariant_comp_ptrlist(i)%u2%local(n)%val )
250 end subroutine meshfieldcommcubedspheredom3d_get
253 subroutine meshfieldcommcubedspheredom3d_exchange( this, do_wait )
265 logical,
intent(in),
optional :: do_wait
270 real(rp),
allocatable :: fpos3d(:,:)
271 real(rp),
allocatable :: lcfpos3d(:,:)
272 real(rp),
allocatable :: unity_fac(:)
273 real(rp),
allocatable :: tmp_svec3d(:,:)
274 real(rp),
allocatable :: tmp1_htensor3d(:,:,:)
275 real(rp),
allocatable :: tmp2_htensor3d(:,:,:)
284 do n=1, this%mesh%LOCAL_MESH_NUM
285 lcmesh => this%mesh3d%lcmesh_list(n)
286 elem => lcmesh%refElem3D
288 allocate( fpos3d(this%Nnode_LCMeshAllFace(n),2) )
289 call extract_boundary_data3d( lcmesh%pos_en(:,:,1), elem, lcmesh, fpos3d(:,1) )
290 call extract_boundary_data3d( lcmesh%pos_en(:,:,2), elem, lcmesh, fpos3d(:,2) )
292 do f=1, this%nfaces_comm
293 commdata => this%commdata_list(f,n)
294 call push_localsendbuf( commdata%send_buf(:,:), &
295 this%send_buf(:,:,n), commdata%s_faceID, this%is_f(f,n), &
296 commdata%Nnode_LCMeshFace, this%field_num_tot, &
299 if ( commdata%s_panelID /= lcmesh%panelID )
then
300 if ( this%hvfield_num > 0 )
then
302 allocate( lcfpos3d(commdata%Nnode_LCMeshFace,2), unity_fac(commdata%Nnode_LCMeshFace) )
303 allocate( tmp_svec3d(commdata%Nnode_LCMeshFace,2) )
305 call push_localsendbuf( lcfpos3d, &
306 fpos3d, commdata%s_faceID, this%is_f(f,n), commdata%Nnode_LCMeshFace, 2, &
308 unity_fac(:) = 1.0_rp
310 do varid=this%sfield_num+1, this%sfield_num+2*this%hvfield_num-1, 2
311 tmp_svec3d(:,1) = commdata%send_buf(:,varid )
312 tmp_svec3d(:,2) = commdata%send_buf(:,varid+1)
315 lcmesh%panelID, lcfpos3d(:,1), lcfpos3d(:,2), unity_fac(:), &
316 commdata%Nnode_LCMeshFace, &
317 tmp_svec3d(:,1), tmp_svec3d(:,2), &
318 commdata%send_buf(:,varid), commdata%send_buf(:,varid+1) )
321 deallocate( lcfpos3d, unity_fac, tmp_svec3d )
324 if ( this%htensorfield_num > 0 )
then
325 allocate( lcfpos3d(commdata%Nnode_LCMeshFace,2), unity_fac(commdata%Nnode_LCMeshFace) )
326 allocate( tmp1_htensor3d(commdata%Nnode_LCMeshFace,2,2) )
327 allocate( tmp2_htensor3d(commdata%Nnode_LCMeshFace,2,2) )
329 call push_localsendbuf( lcfpos3d, &
330 fpos3d, commdata%s_faceID, this%is_f(f,n), commdata%Nnode_LCMeshFace, 2, &
332 unity_fac(:) = 1.0_rp
334 do varid=this%sfield_num+2*this%hvfield_num+1, this%field_num_tot-3, 4
335 tmp1_htensor3d(:,1,1) = commdata%send_buf(:,varid )
336 tmp1_htensor3d(:,2,1) = commdata%send_buf(:,varid+1)
337 tmp1_htensor3d(:,1,2) = commdata%send_buf(:,varid+2)
338 tmp1_htensor3d(:,2,2) = commdata%send_buf(:,varid+3)
341 lcmesh%panelID, lcfpos3d(:,1), lcfpos3d(:,2), unity_fac(:), &
342 commdata%Nnode_LCMeshFace, &
343 tmp1_htensor3d(:,1,1), tmp1_htensor3d(:,2,1), &
344 tmp2_htensor3d(:,1,1), tmp2_htensor3d(:,2,1) )
346 lcmesh%panelID, lcfpos3d(:,1), lcfpos3d(:,2), unity_fac(:), &
347 commdata%Nnode_LCMeshFace, &
348 tmp1_htensor3d(:,1,2), tmp1_htensor3d(:,2,2), &
349 tmp2_htensor3d(:,1,2), tmp2_htensor3d(:,2,2) )
351 lcmesh%panelID, lcfpos3d(:,1), lcfpos3d(:,2), unity_fac(:), &
352 commdata%Nnode_LCMeshFace, &
353 tmp2_htensor3d(:,1,1), tmp2_htensor3d(:,1,2), &
354 commdata%send_buf(:,varid), commdata%send_buf(:,varid+2) )
356 lcmesh%panelID, lcfpos3d(:,1), lcfpos3d(:,2), unity_fac(:), &
357 commdata%Nnode_LCMeshFace, &
358 tmp2_htensor3d(:,2,1), tmp2_htensor3d(:,2,2), &
359 commdata%send_buf(:,varid+1), commdata%send_buf(:,varid+3) )
361 deallocate( lcfpos3d, unity_fac, tmp1_htensor3d, tmp2_htensor3d )
373 if ( .not. this%call_wait_flag_sub_get ) &
374 call post_exchange_core( this )
377 end subroutine meshfieldcommcubedspheredom3d_exchange
381 subroutine post_exchange_core( this )
392 real(rp),
allocatable :: fpos3d(:,:)
393 real(rp),
allocatable :: lcfpos3d(:,:)
394 real(rp),
allocatable :: unity_fac(:)
395 real(rp),
allocatable :: tmp1_htensor3d(:,:,:)
404 do n=1, this%mesh%LOCAL_MESH_NUM
405 lcmesh => this%mesh3d%lcmesh_list(n)
406 elem => lcmesh%refElem3D
408 allocate( fpos3d(this%Nnode_LCMeshAllFace(n),2) )
409 call extract_boundary_data3d( lcmesh%pos_en(:,:,1), elem, lcmesh, fpos3d(:,1) )
410 call extract_boundary_data3d( lcmesh%pos_en(:,:,2), elem, lcmesh, fpos3d(:,2) )
413 do f=1, this%nfaces_comm
414 commdata => this%commdata_list(f,n)
415 ire = irs + commdata%Nnode_LCMeshFace - 1
417 if ( commdata%s_panelID /= lcmesh%panelID )
then
418 if ( this%hvfield_num > 0 )
then
420 allocate( lcfpos3d(commdata%Nnode_LCMeshFace,2), unity_fac(commdata%Nnode_LCMeshFace) )
421 call push_localsendbuf( lcfpos3d, &
422 fpos3d, f, this%is_f(f,n), commdata%Nnode_LCMeshFace, 2, &
424 unity_fac(:) = 1.0_rp
426 do varid=this%sfield_num+1, this%sfield_num+2*this%hvfield_num-1, 2
428 lcmesh%panelID, lcfpos3d(:,1), lcfpos3d(:,2), unity_fac(:), &
429 commdata%Nnode_LCMeshFace, &
430 commdata%recv_buf(:,varid), commdata%recv_buf(:,varid+1), &
431 this%recv_buf(irs:ire,varid,n), this%recv_buf(irs:ire,varid+1,n) )
433 deallocate( lcfpos3d, unity_fac )
436 if ( this%htensorfield_num > 0 )
then
437 allocate( lcfpos3d(commdata%Nnode_LCMeshFace,2), unity_fac(commdata%Nnode_LCMeshFace) )
438 allocate( tmp1_htensor3d(commdata%Nnode_LCMeshFace,2,2) )
439 unity_fac(:) = 1.0_rp
441 call push_localsendbuf( lcfpos3d, &
442 fpos3d, f, this%is_f(f,n), commdata%Nnode_LCMeshFace, 2, &
445 do varid=this%sfield_num+2*this%hvfield_num+1, this%field_num_tot-3, 4
447 lcmesh%panelID, lcfpos3d(:,1), lcfpos3d(:,2), unity_fac(:), &
448 commdata%Nnode_LCMeshFace, &
449 commdata%recv_buf(:,varid), commdata%recv_buf(:,varid+1), &
450 tmp1_htensor3d(:,1,1), tmp1_htensor3d(:,2,1) )
452 lcmesh%panelID, lcfpos3d(:,1), lcfpos3d(:,2), unity_fac(:), &
453 commdata%Nnode_LCMeshFace, &
454 commdata%recv_buf(:,varid+2), commdata%recv_buf(:,varid+3), &
455 tmp1_htensor3d(:,1,2), tmp1_htensor3d(:,2,2) )
457 lcmesh%panelID, lcfpos3d(:,1), lcfpos3d(:,2), unity_fac(:), &
458 commdata%Nnode_LCMeshFace, &
459 tmp1_htensor3d(:,1,1), tmp1_htensor3d(:,1,2), &
460 this%recv_buf(irs:ire,varid,n), this%recv_buf(irs:ire,varid+2,n) )
462 lcmesh%panelID, lcfpos3d(:,1), lcfpos3d(:,2), unity_fac(:), &
463 commdata%Nnode_LCMeshFace, &
464 tmp1_htensor3d(:,2,1), tmp1_htensor3d(:,2,2), &
465 this%recv_buf(irs:ire,varid+1,n), this%recv_buf(irs:ire,varid+3,n) )
467 deallocate( lcfpos3d, unity_fac, tmp1_htensor3d )
478 end subroutine post_exchange_core
481 subroutine push_localsendbuf( lc_send_buf, send_buf, s_faceID, is, Nnode_LCMeshFace, var_num, &
485 integer,
intent(in) :: var_num
486 integer,
intent(in) :: nnode_lcmeshface
487 real(rp),
intent(inout) :: lc_send_buf(nnode_lcmeshface,var_num)
489 integer,
intent(in) :: s_faceid, is
496 ie = is + nnode_lcmeshface - 1
497 if ( s_faceid > 0 )
then
498 lc_send_buf(:,:) = send_buf(is:ie,:)
500 call revert_hori( lc_send_buf, send_buf(is:ie,:), mesh3d, elem3d, mesh3d%lcmesh2D )
505 subroutine revert_hori( revert, ori, mesh, e3D, mesh2D )
510 real(rp),
intent(out) :: revert(e3d%nnode_h1d,e3d%nnode_v,mesh2d%nex,mesh%nez, var_num)
511 real(rp),
intent(in) :: ori(e3d%nnode_h1d,e3d%nnode_v,mesh2d%nex,mesh%nez, var_num)
513 integer :: p1, p3, i, k, n
520 i_ = mesh2d%NeX - i + 1
522 do p1=1, e3d%Nnode_h1D
523 p1_ = e3d%Nnode_h1D - p1 + 1
524 revert(p1,p3,i,k,n) = ori(p1_,p3,i_,k,n)
532 end subroutine revert_hori
533 end subroutine push_localsendbuf
536 subroutine extract_boundary_data3d( var, elem, mesh, buf )
541 real(dp),
intent(in) :: var(elem%np * mesh%ne)
542 real(dp),
intent(inout) :: buf( 2*(mesh%nex + mesh%ney)*mesh%nez*elem%nfp_h &
543 + 2*mesh%nex*mesh%ney*elem%nfp_v )
546 buf(:) = var(mesh%VmapB(:))
549 end subroutine extract_boundary_data3d
552 subroutine set_boundary_data3d_u1u2( buf_U, buf_V, &
560 real(dp),
intent(in) :: buf_u(2*(mesh%nex + mesh%ney)*mesh%nez*elem%nfp_h &
561 + 2*mesh%nex*mesh%ney*elem%nfp_v)
562 real(dp),
intent(in) :: buf_v(2*(mesh%nex + mesh%ney)*mesh%nez*elem%nfp_h &
563 + 2*mesh%nex*mesh%ney*elem%nfp_v)
564 real(dp),
intent(in) :: g_ij(elem%np * mesh%ne,2,2)
565 real(dp),
intent(inout) :: u1(elem%np * mesh%nea)
566 real(dp),
intent(inout) :: u2(elem%np * mesh%nea)
569 u1(elem%Np*mesh%NeE+1:elem%Np*mesh%NeE+
size(buf_u)) &
570 = g_ij(mesh%VmapB,1,1) * buf_u(:) + g_ij(mesh%VmapB,1,2) * buf_v(:)
571 u2(elem%Np*mesh%NeE+1:elem%Np*mesh%NeE+
size(buf_u)) &
572 = g_ij(mesh%VmapB,2,1) * buf_u(:) + g_ij(mesh%VmapB,2,2) * buf_v(:)
575 end subroutine set_boundary_data3d_u1u2
577end module scale_meshfieldcomm_cubedspheredom3d
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 / Local 3D
module FElib / Mesh / Cubed-sphere 3D 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 3D 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.