43 type :: veccovariantcomp
44 type(MeshField3D),
pointer :: u1 => null()
45 type(MeshField3D),
pointer :: u2 => null()
51 type(veccovariantcomp),
allocatable :: vec_covariant_comp_ptrlist(:)
52 integer,
allocatable :: nnode_lcmeshallface(:)
54 procedure,
public :: init => meshfieldcommcubedspheredom3d_init
55 procedure,
public :: put => meshfieldcommcubedspheredom3d_put
56 procedure,
public :: get => meshfieldcommcubedspheredom3d_get
57 procedure,
public :: exchange => meshfieldcommcubedspheredom3d_exchange
58 procedure,
public :: setcovariantvec => meshfieldcommcubedspheredom3d_set_covariantvec
59 procedure,
public :: final => meshfieldcommcubedspheredom3d_final
71 private :: post_exchange_core
72 private :: push_localsendbuf
79 integer,
parameter :: comm_face_num = 6
83 subroutine meshfieldcommcubedspheredom3d_init( this, &
84 sfield_num, hvfield_num, htensorfield_num, mesh3d )
89 integer,
intent(in) :: sfield_num
90 integer,
intent(in) :: hvfield_num
91 integer,
intent(in) :: htensorfield_num
97 integer :: nnode_lcmeshface(comm_face_num,mesh3d%local_mesh_num)
100 this%mesh3d => mesh3d
101 lcmesh => mesh3d%lcmesh_list(1)
102 elem => lcmesh%refElem3D
105 + 2*lcmesh%NeX*lcmesh%NeY*elem%Nfp_v
107 allocate( this%Nnode_LCMeshAllFace(mesh3d%LOCAL_MESH_NUM) )
108 do n=1, this%mesh3d%LOCAL_MESH_NUM
109 lcmesh => this%mesh3d%lcmesh_list(n)
110 nnode_lcmeshface(:,n) = &
111 (/ lcmesh%NeX, lcmesh%NeY, lcmesh%NeX, lcmesh%NeY, 0, 0 /) * lcmesh%NeZ * lcmesh%refElem3D%Nfp_h &
112 + (/ 0, 0, 0, 0, 1, 1 /) * lcmesh%NeX*lcmesh%NeY * lcmesh%refElem3D%Nfp_v
113 this%Nnode_LCMeshAllFace(n) = sum(nnode_lcmeshface(:,n))
118 if (hvfield_num > 0)
then
119 allocate( this%vec_covariant_comp_ptrlist(hvfield_num) )
123 end subroutine meshfieldcommcubedspheredom3d_init
126 subroutine meshfieldcommcubedspheredom3d_final( this )
133 if ( this%hvfield_num > 0 )
then
134 deallocate( this%vec_covariant_comp_ptrlist )
140 end subroutine meshfieldcommcubedspheredom3d_final
142 subroutine meshfieldcommcubedspheredom3d_set_covariantvec( &
143 this, hvfield_ID, u1, u2 )
146 integer,
intent(in) :: hvfield_id
151 this%vec_covariant_comp_ptrlist(hvfield_id)%u1 => u1
152 this%vec_covariant_comp_ptrlist(hvfield_id)%u2 => u2
155 end subroutine meshfieldcommcubedspheredom3d_set_covariantvec
158 subroutine meshfieldcommcubedspheredom3d_put(this, field_list, varid_s)
162 integer,
intent(in) :: varid_s
169 do i=1,
size(field_list)
170 do n=1, this%mesh%LOCAL_MESH_NUM
171 lcmesh => this%mesh3d%lcmesh_list(n)
173 this%send_buf(:,varid_s+i-1,n) )
178 end subroutine meshfieldcommcubedspheredom3d_put
181 subroutine meshfieldcommcubedspheredom3d_get(this, field_list, varid_s)
188 integer,
intent(in) :: varid_s
199 integer :: varid_vec_s
201 real(rp),
allocatable :: g_ij(:,:,:,:)
204 varnum =
size(field_list)
207 if ( this%call_wait_flag_sub_get )
then
209 call post_exchange_core( this )
214 do n=1, this%mesh3d%LOCAL_MESH_NUM
215 lcmesh => this%mesh3d%lcmesh_list(n)
217 field_list(i)%field3d%local(n)%val )
221 varid_e = varid_s + varnum - 1
222 if ( varid_e > this%sfield_num )
then
223 do i=1, this%hvfield_num
225 varid_vec_s = this%sfield_num + 2*i - 1
226 if ( varid_vec_s > varid_e )
exit
228 if (
associated(this%vec_covariant_comp_ptrlist(i)%u1 ) &
229 .and.
associated(this%vec_covariant_comp_ptrlist(i)%u2 ) )
then
231 do n=1, this%mesh3d%LOCAL_MESH_NUM
232 lcmesh => this%mesh3d%lcmesh_list(n)
233 elem => lcmesh%refElem3D
235 allocate( g_ij(elem%Np,lcmesh%Ne,2,2) )
237 do ke=lcmesh%NeS, lcmesh%NeE
238 ke2d = lcmesh%EMap3Dto2D(ke)
239 g_ij(:,ke,1,1) = lcmesh%G_ij(elem%IndexH2Dto3D(:),ke2d,1,1)
240 g_ij(:,ke,2,1) = lcmesh%G_ij(elem%IndexH2Dto3D(:),ke2d,2,1)
241 g_ij(:,ke,1,2) = lcmesh%G_ij(elem%IndexH2Dto3D(:),ke2d,1,2)
242 g_ij(:,ke,2,2) = lcmesh%G_ij(elem%IndexH2Dto3D(:),ke2d,2,2)
245 call set_boundary_data3d_u1u2( &
246 this%recv_buf(:,varid_vec_s,n), this%recv_buf(:,varid_vec_s+1,n), &
247 lcmesh%refElem3D, lcmesh, g_ij(:,:,:,:), &
248 this%vec_covariant_comp_ptrlist(i)%u1%local(n)%val, &
249 this%vec_covariant_comp_ptrlist(i)%u2%local(n)%val )
256 end subroutine meshfieldcommcubedspheredom3d_get
262 subroutine meshfieldcommcubedspheredom3d_exchange( this, do_wait )
274 logical,
intent(in),
optional :: do_wait
279 real(rp),
allocatable :: fpos3d(:,:)
280 real(rp),
allocatable :: lcfpos3d(:,:)
281 real(rp),
allocatable :: unity_fac(:)
282 real(rp),
allocatable :: tmp_svec3d(:,:)
283 real(rp),
allocatable :: tmp1_htensor3d(:,:,:)
284 real(rp),
allocatable :: tmp2_htensor3d(:,:,:)
293 do n=1, this%mesh%LOCAL_MESH_NUM
294 lcmesh => this%mesh3d%lcmesh_list(n)
295 elem => lcmesh%refElem3D
297 allocate( fpos3d(this%Nnode_LCMeshAllFace(n),2) )
298 call extract_boundary_data3d( lcmesh%pos_en(:,:,1), elem, lcmesh, fpos3d(:,1) )
299 call extract_boundary_data3d( lcmesh%pos_en(:,:,2), elem, lcmesh, fpos3d(:,2) )
301 do f=1, this%nfaces_comm
302 commdata => this%commdata_list(f,n)
303 call push_localsendbuf( commdata%send_buf(:,:), &
304 this%send_buf(:,:,n), commdata%s_faceID, this%is_f(f,n), &
305 commdata%Nnode_LCMeshFace, this%field_num_tot, &
308 if ( commdata%s_panelID /= lcmesh%panelID )
then
309 if ( this%hvfield_num > 0 )
then
311 allocate( lcfpos3d(commdata%Nnode_LCMeshFace,2), unity_fac(commdata%Nnode_LCMeshFace) )
312 allocate( tmp_svec3d(commdata%Nnode_LCMeshFace,2) )
314 call push_localsendbuf( lcfpos3d, &
315 fpos3d, commdata%s_faceID, this%is_f(f,n), commdata%Nnode_LCMeshFace, 2, &
317 unity_fac(:) = 1.0_rp
319 do varid=this%sfield_num+1, this%sfield_num+2*this%hvfield_num-1, 2
320 tmp_svec3d(:,1) = commdata%send_buf(:,varid )
321 tmp_svec3d(:,2) = commdata%send_buf(:,varid+1)
324 lcmesh%panelID, lcfpos3d(:,1), lcfpos3d(:,2), unity_fac(:), &
325 commdata%Nnode_LCMeshFace, &
326 tmp_svec3d(:,1), tmp_svec3d(:,2), &
327 commdata%send_buf(:,varid), commdata%send_buf(:,varid+1) )
330 deallocate( lcfpos3d, unity_fac, tmp_svec3d )
333 if ( this%htensorfield_num > 0 )
then
334 allocate( lcfpos3d(commdata%Nnode_LCMeshFace,2), unity_fac(commdata%Nnode_LCMeshFace) )
335 allocate( tmp1_htensor3d(commdata%Nnode_LCMeshFace,2,2) )
336 allocate( tmp2_htensor3d(commdata%Nnode_LCMeshFace,2,2) )
338 call push_localsendbuf( lcfpos3d, &
339 fpos3d, commdata%s_faceID, this%is_f(f,n), commdata%Nnode_LCMeshFace, 2, &
341 unity_fac(:) = 1.0_rp
343 do varid=this%sfield_num+2*this%hvfield_num+1, this%field_num_tot-3, 4
344 tmp1_htensor3d(:,1,1) = commdata%send_buf(:,varid )
345 tmp1_htensor3d(:,2,1) = commdata%send_buf(:,varid+1)
346 tmp1_htensor3d(:,1,2) = commdata%send_buf(:,varid+2)
347 tmp1_htensor3d(:,2,2) = commdata%send_buf(:,varid+3)
350 lcmesh%panelID, lcfpos3d(:,1), lcfpos3d(:,2), unity_fac(:), &
351 commdata%Nnode_LCMeshFace, &
352 tmp1_htensor3d(:,1,1), tmp1_htensor3d(:,2,1), &
353 tmp2_htensor3d(:,1,1), tmp2_htensor3d(:,2,1) )
355 lcmesh%panelID, lcfpos3d(:,1), lcfpos3d(:,2), unity_fac(:), &
356 commdata%Nnode_LCMeshFace, &
357 tmp1_htensor3d(:,1,2), tmp1_htensor3d(:,2,2), &
358 tmp2_htensor3d(:,1,2), tmp2_htensor3d(:,2,2) )
360 lcmesh%panelID, lcfpos3d(:,1), lcfpos3d(:,2), unity_fac(:), &
361 commdata%Nnode_LCMeshFace, &
362 tmp2_htensor3d(:,1,1), tmp2_htensor3d(:,1,2), &
363 commdata%send_buf(:,varid), commdata%send_buf(:,varid+2) )
365 lcmesh%panelID, lcfpos3d(:,1), lcfpos3d(:,2), unity_fac(:), &
366 commdata%Nnode_LCMeshFace, &
367 tmp2_htensor3d(:,2,1), tmp2_htensor3d(:,2,2), &
368 commdata%send_buf(:,varid+1), commdata%send_buf(:,varid+3) )
370 deallocate( lcfpos3d, unity_fac, tmp1_htensor3d, tmp2_htensor3d )
382 if ( .not. this%call_wait_flag_sub_get ) &
383 call post_exchange_core( this )
386 end subroutine meshfieldcommcubedspheredom3d_exchange
390 subroutine post_exchange_core( this )
401 real(rp),
allocatable :: fpos3d(:,:)
402 real(rp),
allocatable :: lcfpos3d(:,:)
403 real(rp),
allocatable :: unity_fac(:)
404 real(rp),
allocatable :: tmp1_htensor3d(:,:,:)
413 do n=1, this%mesh%LOCAL_MESH_NUM
414 lcmesh => this%mesh3d%lcmesh_list(n)
415 elem => lcmesh%refElem3D
417 allocate( fpos3d(this%Nnode_LCMeshAllFace(n),2) )
418 call extract_boundary_data3d( lcmesh%pos_en(:,:,1), elem, lcmesh, fpos3d(:,1) )
419 call extract_boundary_data3d( lcmesh%pos_en(:,:,2), elem, lcmesh, fpos3d(:,2) )
422 do f=1, this%nfaces_comm
423 commdata => this%commdata_list(f,n)
424 ire = irs + commdata%Nnode_LCMeshFace - 1
426 if ( commdata%s_panelID /= lcmesh%panelID )
then
427 if ( this%hvfield_num > 0 )
then
429 allocate( lcfpos3d(commdata%Nnode_LCMeshFace,2), unity_fac(commdata%Nnode_LCMeshFace) )
430 call push_localsendbuf( lcfpos3d, &
431 fpos3d, f, this%is_f(f,n), commdata%Nnode_LCMeshFace, 2, &
433 unity_fac(:) = 1.0_rp
435 do varid=this%sfield_num+1, this%sfield_num+2*this%hvfield_num-1, 2
437 lcmesh%panelID, lcfpos3d(:,1), lcfpos3d(:,2), unity_fac(:), &
438 commdata%Nnode_LCMeshFace, &
439 commdata%recv_buf(:,varid), commdata%recv_buf(:,varid+1), &
440 this%recv_buf(irs:ire,varid,n), this%recv_buf(irs:ire,varid+1,n) )
442 deallocate( lcfpos3d, unity_fac )
445 if ( this%htensorfield_num > 0 )
then
446 allocate( lcfpos3d(commdata%Nnode_LCMeshFace,2), unity_fac(commdata%Nnode_LCMeshFace) )
447 allocate( tmp1_htensor3d(commdata%Nnode_LCMeshFace,2,2) )
448 unity_fac(:) = 1.0_rp
450 call push_localsendbuf( lcfpos3d, &
451 fpos3d, f, this%is_f(f,n), commdata%Nnode_LCMeshFace, 2, &
454 do varid=this%sfield_num+2*this%hvfield_num+1, this%field_num_tot-3, 4
456 lcmesh%panelID, lcfpos3d(:,1), lcfpos3d(:,2), unity_fac(:), &
457 commdata%Nnode_LCMeshFace, &
458 commdata%recv_buf(:,varid), commdata%recv_buf(:,varid+1), &
459 tmp1_htensor3d(:,1,1), tmp1_htensor3d(:,2,1) )
461 lcmesh%panelID, lcfpos3d(:,1), lcfpos3d(:,2), unity_fac(:), &
462 commdata%Nnode_LCMeshFace, &
463 commdata%recv_buf(:,varid+2), commdata%recv_buf(:,varid+3), &
464 tmp1_htensor3d(:,1,2), tmp1_htensor3d(:,2,2) )
466 lcmesh%panelID, lcfpos3d(:,1), lcfpos3d(:,2), unity_fac(:), &
467 commdata%Nnode_LCMeshFace, &
468 tmp1_htensor3d(:,1,1), tmp1_htensor3d(:,1,2), &
469 this%recv_buf(irs:ire,varid,n), this%recv_buf(irs:ire,varid+2,n) )
471 lcmesh%panelID, lcfpos3d(:,1), lcfpos3d(:,2), unity_fac(:), &
472 commdata%Nnode_LCMeshFace, &
473 tmp1_htensor3d(:,2,1), tmp1_htensor3d(:,2,2), &
474 this%recv_buf(irs:ire,varid+1,n), this%recv_buf(irs:ire,varid+3,n) )
476 deallocate( lcfpos3d, unity_fac, tmp1_htensor3d )
487 end subroutine post_exchange_core
490 subroutine push_localsendbuf( lc_send_buf, send_buf, s_faceID, is, Nnode_LCMeshFace, var_num, &
494 integer,
intent(in) :: var_num
495 integer,
intent(in) :: nnode_lcmeshface
496 real(rp),
intent(inout) :: lc_send_buf(nnode_lcmeshface,var_num)
498 integer,
intent(in) :: s_faceid, is
505 ie = is + nnode_lcmeshface - 1
506 if ( s_faceid > 0 )
then
507 lc_send_buf(:,:) = send_buf(is:ie,:)
509 call revert_hori( lc_send_buf, send_buf(is:ie,:), mesh3d, elem3d, mesh3d%lcmesh2D )
514 subroutine revert_hori( revert, ori, mesh, e3D, mesh2D )
519 real(rp),
intent(out) :: revert(e3d%nnode_h1d,e3d%nnode_v,mesh2d%nex,mesh%nez, var_num)
520 real(rp),
intent(in) :: ori(e3d%nnode_h1d,e3d%nnode_v,mesh2d%nex,mesh%nez, var_num)
522 integer :: p1, p3, i, k, n
529 i_ = mesh2d%NeX - i + 1
531 do p1=1, e3d%Nnode_h1D
532 p1_ = e3d%Nnode_h1D - p1 + 1
533 revert(p1,p3,i,k,n) = ori(p1_,p3,i_,k,n)
541 end subroutine revert_hori
542 end subroutine push_localsendbuf
545 subroutine extract_boundary_data3d( var, elem, mesh, buf )
550 real(dp),
intent(in) :: var(elem%np * mesh%ne)
551 real(dp),
intent(inout) :: buf( 2*(mesh%nex + mesh%ney)*mesh%nez*elem%nfp_h &
552 + 2*mesh%nex*mesh%ney*elem%nfp_v )
555 buf(:) = var(mesh%VmapB(:))
558 end subroutine extract_boundary_data3d
561 subroutine set_boundary_data3d_u1u2( buf_U, buf_V, &
569 real(dp),
intent(in) :: buf_u(2*(mesh%nex + mesh%ney)*mesh%nez*elem%nfp_h &
570 + 2*mesh%nex*mesh%ney*elem%nfp_v)
571 real(dp),
intent(in) :: buf_v(2*(mesh%nex + mesh%ney)*mesh%nez*elem%nfp_h &
572 + 2*mesh%nex*mesh%ney*elem%nfp_v)
573 real(dp),
intent(in) :: g_ij(elem%np * mesh%ne,2,2)
574 real(dp),
intent(inout) :: u1(elem%np * mesh%nea)
575 real(dp),
intent(inout) :: u2(elem%np * mesh%nea)
578 u1(elem%Np*mesh%NeE+1:elem%Np*mesh%NeE+
size(buf_u)) &
579 = g_ij(mesh%VmapB,1,1) * buf_u(:) + g_ij(mesh%VmapB,1,2) * buf_v(:)
580 u2(elem%Np*mesh%NeE+1:elem%Np*mesh%NeE+
size(buf_u)) &
581 = g_ij(mesh%VmapB,2,1) * buf_u(:) + g_ij(mesh%VmapB,2,2) * buf_v(:)
584 end subroutine set_boundary_data3d_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 / 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)
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 3D cubed-sphere domain
integer bufsize_per_field
Buffer size per a field.
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 representing a field with 3D 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 3D cubed-sphere domain.