38 real(rp),
allocatable :: send_buf(:,:)
39 real(rp),
allocatable :: recv_buf(:,:)
40 integer :: nnode_lcmeshface
45 integer :: s_tilelocalid
48 procedure,
public :: init => localmeshcommdata_init
49 procedure,
public :: final => localmeshcommdata_final
50 procedure,
public :: sendrecv => localmeshcommdata_sendrecv
51 procedure,
public :: pc_init => localmeshcommdata_pc_init
57 integer :: hvfield_num
58 integer :: htensorfield_num
59 integer :: field_num_tot
60 integer :: nfaces_comm
63 real(rp),
allocatable :: send_buf(:,:,:)
64 real(rp),
allocatable :: recv_buf(:,:,:)
65 integer,
allocatable :: request_send(:)
66 integer,
allocatable :: request_recv(:)
69 integer,
allocatable :: is_f(:,:)
71 logical :: mpi_pc_flag
72 integer,
allocatable :: request_pc(:)
74 integer :: req_counter
75 logical :: call_wait_flag_sub_get
78 procedure(meshfieldcommbase_get),
public,
deferred :: get
79 procedure(meshfieldcommbase_exchange),
public,
deferred :: exchange
80 procedure,
public :: prepare_pc => meshfieldcommbase_prepare_pc
101 integer,
intent(in) :: varid_s
104 subroutine meshfieldcommbase_get(this, field_list, varid_s)
109 integer,
intent(in) :: varid_s
110 end subroutine meshfieldcommbase_get
112 subroutine meshfieldcommbase_exchange(this, do_wait)
115 logical,
intent(in),
optional :: do_wait
116 end subroutine meshfieldcommbase_exchange
118 subroutine meshfieldcomm_final(this)
121 end subroutine meshfieldcomm_final
151 sfield_num, hvfield_num, htensorfield_num, bufsize_per_field, comm_face_num, &
152 Nnode_LCMeshFace, mesh )
157 class(
meshbase),
intent(in),
target :: mesh
158 integer,
intent(in) :: sfield_num
159 integer,
intent(in) :: hvfield_num
160 integer,
intent(in) :: htensorfield_num
161 integer,
intent(in) :: bufsize_per_field
162 integer,
intent(in) :: comm_face_num
163 integer,
intent(in) :: nnode_lcmeshface(comm_face_num,mesh%local_mesh_num)
170 this%sfield_num = sfield_num
171 this%hvfield_num = hvfield_num
172 this%htensorfield_num = htensorfield_num
173 this%field_num_tot = sfield_num + hvfield_num*2 + htensorfield_num*4
174 this%nfaces_comm = comm_face_num
176 if (this%field_num_tot > 0)
then
177 allocate( this%send_buf(bufsize_per_field, this%field_num_tot, mesh%LOCAL_MESH_NUM) )
178 allocate( this%recv_buf(bufsize_per_field, this%field_num_tot, mesh%LOCAL_MESH_NUM) )
179 allocate( this%request_send(comm_face_num*mesh%LOCAL_MESH_NUM) )
180 allocate( this%request_recv(comm_face_num*mesh%LOCAL_MESH_NUM) )
182 allocate( this%commdata_list(comm_face_num,mesh%LOCAL_MESH_NUM) )
183 allocate( this%is_f(comm_face_num,mesh%LOCAL_MESH_NUM) )
185 do n=1, mesh%LOCAL_MESH_NUM
187 do f=2, this%nfaces_comm
188 this%is_f(f,n) = this%is_f(f-1,n) + nnode_lcmeshface(f-1,n)
191 call mesh%GetLocalMesh(n, lcmesh)
192 do f=1, this%nfaces_comm
193 call this%commdata_list(f,n)%Init( this, lcmesh, f, nnode_lcmeshface(f,n) )
198 this%MPI_pc_flag = .false.
216 if (this%field_num_tot > 0)
then
217 deallocate( this%send_buf, this%recv_buf )
218 deallocate( this%request_send, this%request_recv )
220 do n=1, this%mesh%LOCAL_MESH_NUM
221 do f=1, this%nfaces_comm
222 call this%commdata_list(f,n)%Final()
225 deallocate( this%commdata_list, this%is_f )
227 if (this%MPI_pc_flag)
then
228 do ireq=1, this%req_counter
229 call mpi_request_free( this%request_pc(ireq), ierr )
231 deallocate( this%request_pc )
241 subroutine meshfieldcommbase_prepare_pc( this )
248 this%MPI_pc_flag = .true.
249 if (this%field_num_tot > 0)
then
250 allocate( this%request_pc(2*this%nfaces_comm*this%mesh%LOCAL_MESH_NUM) )
254 do n=1, this%mesh%LOCAL_MESH_NUM
255 do f=1, this%nfaces_comm
256 call this%commdata_list(f,n)%PC_Init( &
257 this%req_counter, this%request_pc(:) )
262 end subroutine meshfieldcommbase_prepare_pc
274 type(
localmeshcommdata),
intent(inout),
target :: commdata_list(this%nfaces_comm, this%mesh%local_mesh_num)
275 logical,
intent(in),
optional :: do_wait
282 if (
present(do_wait) )
then
289 if ( this%MPI_pc_flag )
then
290 call mpi_startall( this%req_counter, this%request_pc(1:this%req_counter), ierr )
294 do n=1, this%mesh%LOCAL_MESH_NUM
295 do f=1, this%nfaces_comm
296 call commdata_list(f,n)%SendRecv( &
297 this%req_counter, this%request_send(:), this%request_recv(:), &
298 commdata_list(:,:), &
306 this%call_wait_flag_sub_get = .false.
308 this%call_wait_flag_sub_get = .true.
326 type(
localmeshcommdata),
intent(inout),
target :: commdata_list(this%nfaces_comm, this%mesh%local_mesh_num)
329 integer :: stat_send(mpi_status_size, this%nfaces_comm * this%mesh%local_mesh_num)
330 integer :: stat_recv(mpi_status_size, this%nfaces_comm * this%mesh%local_mesh_num)
331 integer :: stat_pc(mpi_status_size, 2*this%nfaces_comm * this%mesh%local_mesh_num)
339 if ( this%MPI_pc_flag )
then
340 if (this%req_counter > 0)
then
341 call mpi_waitall( this%req_counter, this%request_pc(1:this%req_counter), stat_pc(:,1:this%req_counter), ierr )
344 if ( this%req_counter > 0 )
then
345 call mpi_waitall( this%req_counter, this%request_recv(1:this%req_counter), stat_recv(:,1:this%req_counter), ierr )
346 call mpi_waitall( this%req_counter, this%request_send(1:this%req_counter), stat_send(:,1:this%req_counter), ierr )
354 do n=1, this%mesh%LOCAL_MESH_NUM
356 do var_id=1, this%field_num_tot
358 do f=1, this%nfaces_comm
359 ire = irs + commdata_list(f,n)%Nnode_LCMeshFace - 1
360 this%recv_buf(irs:ire,var_id,n) = commdata_list(f,n)%recv_buf(:,var_id)
376 real(rp),
intent(in) :: var(refelem%np * mesh%nea)
377 real(rp),
intent(out) :: buf(size(mesh%vmapb))
384 buf(i) = var(mesh%VmapB(i))
394 real(rp),
intent(in) :: buf(size(mesh%vmapb))
395 real(rp),
intent(inout) :: var(refelem%np * mesh%nea)
398 var(refelem%Np*mesh%NeE+1:refelem%Np*mesh%NeE+
size(buf)) = buf(:)
404 subroutine localmeshcommdata_init( this, comm, lcmesh, faceID, Nnode_LCMeshFace )
410 integer,
intent(in) :: faceid
411 integer,
intent(in) :: nnode_lcmeshface
414 this%lcmesh => lcmesh
415 this%Nnode_LCMeshFace = nnode_lcmeshface
417 allocate( this%send_buf(nnode_lcmeshface, comm%field_num_tot))
418 allocate( this%recv_buf(nnode_lcmeshface, comm%field_num_tot))
420 this%s_tileID = comm%mesh%tileID_globalMap(faceid, lcmesh%tileID)
421 this%s_faceID = comm%mesh%tileFaceID_globalMap(faceid, lcmesh%tileID)
422 this%s_panelID = comm%mesh%tilePanelID_globalMap(faceid, lcmesh%tileID)
423 this%s_rank = comm%mesh%PRCrank_globalMap(this%s_tileID)
424 this%s_tilelocalID = comm%mesh%tileID_global2localMap(this%s_tileID)
428 end subroutine localmeshcommdata_init
430 subroutine localmeshcommdata_sendrecv( this, &
431 req_counter, req_send, req_recv, &
432 lccommdat_list, MPI_pc_flag )
434 use scale_prc,
only: &
442 integer,
intent(inout) :: req_counter
443 integer,
intent(inout) :: req_send(:)
444 integer,
intent(inout) :: req_recv(:)
446 logical,
intent(in) :: mpi_pc_flag
453 if ( this%s_rank /= this%lcmesh%PRC_myrank &
454 .and. (.not. mpi_pc_flag ) )
then
456 req_counter = req_counter + 1
458 tag = 10 * this%lcmesh%tileID + this%faceID
459 bufsize =
size(this%recv_buf)
460 call mpi_irecv( this%recv_buf(1,1), bufsize, mpi_double_precision, &
461 this%s_rank, tag, prc_local_comm_world, &
462 req_recv(req_counter), ierr)
464 bufsize =
size(this%send_buf)
465 tag = 10 * this%s_tileID + abs(this%s_faceID)
466 call mpi_isend( this%send_buf(1,1), bufsize, mpi_double_precision, &
467 this%s_rank, tag, prc_local_comm_world, &
468 req_send(req_counter), ierr )
470 else if ( this%s_rank == this%lcmesh%PRC_myrank )
then
471 lccommdat_list(abs(this%s_faceID), this%s_tilelocalID)%recv_buf(:,:) &
476 end subroutine localmeshcommdata_sendrecv
478 subroutine localmeshcommdata_pc_init( this, &
481 use scale_prc,
only: &
489 integer,
intent(inout) :: req_counter
490 integer,
intent(inout) :: req(:)
497 if ( this%s_rank /= this%lcmesh%PRC_myrank )
then
499 req_counter = req_counter + 1
501 tag = 100 * this%lcmesh%tileID + this%faceID
502 bufsize =
size(this%recv_buf)
503 call mpi_recv_init( this%recv_buf(1,1), bufsize, mpi_double_precision, &
504 this%s_rank, tag, prc_local_comm_world, &
505 req(req_counter), ierr)
508 req_counter = req_counter + 1
510 bufsize =
size(this%send_buf)
511 tag = 100 * this%s_tileID + abs(this%s_faceID)
512 call mpi_send_init( this%send_buf(1,1), bufsize, mpi_double_precision, &
513 this%s_rank, tag, prc_local_comm_world, &
514 req(req_counter), ierr )
518 end subroutine localmeshcommdata_pc_init
520 subroutine localmeshcommdata_final( this )
526 deallocate( this%send_buf )
527 deallocate( this%recv_buf )
530 end subroutine localmeshcommdata_final
532end module scale_meshfieldcomm_base
module FElib / Element / Base
module FElib / Mesh / Local, Base
module FElib / Mesh / Base
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.
Derived type to manage data communication between a face of local mesh.
Base derived type to manage data communication.