18 use scale_prc,
only: &
40 real(rp),
allocatable :: send_buf(:,:)
41 real(rp),
allocatable :: recv_buf(:,:)
42 integer :: nnode_lcmeshface
47 integer :: s_tilelocalid
50 procedure,
public :: init => localmeshcommdata_init
51 procedure,
public :: final => localmeshcommdata_final
52 procedure,
public :: sendrecv => localmeshcommdata_sendrecv
53 procedure,
public :: pc_init_recv => localmeshcommdata_pc_init_recv
54 procedure,
public :: pc_init_send => localmeshcommdata_pc_init_send
60 integer :: hvfield_num
61 integer :: htensorfield_num
62 integer :: field_num_tot
63 integer :: nfaces_comm
66 real(rp),
allocatable :: send_buf(:,:,:)
67 real(rp),
allocatable :: recv_buf(:,:,:)
68 integer,
allocatable :: request_send(:)
69 integer,
allocatable :: request_recv(:)
72 integer,
allocatable :: is_f(:,:)
74 logical :: mpi_pc_flag
75 logical :: use_mpi_pc_fujitsu_ext
76 integer,
allocatable :: request_pc(:)
78 integer :: req_counter
79 logical :: call_wait_flag_sub_get
84 procedure(meshfieldcommbase_get),
public,
deferred :: get
85 procedure(meshfieldcommbase_exchange),
public,
deferred :: exchange
86 procedure,
public :: prepare_pc => meshfieldcommbase_prepare_pc
109 integer,
intent(in) :: varid_s
112 subroutine meshfieldcommbase_get(this, field_list, varid_s)
117 integer,
intent(in) :: varid_s
118 end subroutine meshfieldcommbase_get
120 subroutine meshfieldcommbase_exchange(this, do_wait)
123 logical,
intent(in),
optional :: do_wait
124 end subroutine meshfieldcommbase_exchange
126 subroutine meshfieldcomm_final(this)
129 end subroutine meshfieldcomm_final
147 integer :: obj_ind = 0
148 integer,
parameter :: OBJ_INDEX_MAX = 8192
162 sfield_num, hvfield_num, htensorfield_num, bufsize_per_field, comm_face_num, &
163 Nnode_LCMeshFace, mesh )
168 class(
meshbase),
intent(in),
target :: mesh
169 integer,
intent(in) :: sfield_num
170 integer,
intent(in) :: hvfield_num
171 integer,
intent(in) :: htensorfield_num
172 integer,
intent(in) :: bufsize_per_field
173 integer,
intent(in) :: comm_face_num
174 integer,
intent(in) :: nnode_lcmeshface(comm_face_num,mesh%local_mesh_num)
181 this%sfield_num = sfield_num
182 this%hvfield_num = hvfield_num
183 this%htensorfield_num = htensorfield_num
184 this%field_num_tot = sfield_num + hvfield_num*2 + htensorfield_num*4
185 this%nfaces_comm = comm_face_num
187 if (this%field_num_tot > 0)
then
188 allocate( this%send_buf(bufsize_per_field, this%field_num_tot, mesh%LOCAL_MESH_NUM) )
189 allocate( this%recv_buf(bufsize_per_field, this%field_num_tot, mesh%LOCAL_MESH_NUM) )
190 allocate( this%request_send(comm_face_num*mesh%LOCAL_MESH_NUM) )
191 allocate( this%request_recv(comm_face_num*mesh%LOCAL_MESH_NUM) )
193 allocate( this%commdata_list(comm_face_num,mesh%LOCAL_MESH_NUM) )
194 allocate( this%is_f(comm_face_num,mesh%LOCAL_MESH_NUM) )
196 do n=1, mesh%LOCAL_MESH_NUM
198 do f=2, this%nfaces_comm
199 this%is_f(f,n) = this%is_f(f-1,n) + nnode_lcmeshface(f-1,n)
202 call mesh%GetLocalMesh(n, lcmesh)
203 do f=1, this%nfaces_comm
204 call this%commdata_list(f,n)%Init( this, lcmesh, f, nnode_lcmeshface(f,n) )
209 this%MPI_pc_flag = .false.
211 obj_ind = obj_ind + 1
212 if ( obj_ind > obj_index_max )
then
213 log_error(
"MeshFieldCommBase_Init",*)
'obj_ind > OBJ_INDEX_MAX. Check!'
216 this%obj_ind = obj_ind
234 if (this%field_num_tot > 0)
then
235 deallocate( this%send_buf, this%recv_buf )
236 deallocate( this%request_send, this%request_recv )
238 do n=1, this%mesh%LOCAL_MESH_NUM
239 do f=1, this%nfaces_comm
240 call this%commdata_list(f,n)%Final()
243 deallocate( this%commdata_list, this%is_f )
245 if ( this%MPI_pc_flag )
then
246 do ireq=1, this%req_counter
247 call mpi_request_free( this%request_pc(ireq), ierr )
249 deallocate( this%request_pc )
261 subroutine meshfieldcommbase_prepare_pc( this, &
262 use_mpi_pc_fujitsu_ext )
265 logical,
intent(in),
optional :: use_mpi_pc_fujitsu_ext
270 this%MPI_pc_flag = .true.
273 this%use_mpi_pc_fujitsu_ext = .true.
275 this%use_mpi_pc_fujitsu_ext = .false.
277 if (
present(use_mpi_pc_fujitsu_ext) )
then
278 this%use_mpi_pc_fujitsu_ext = use_mpi_pc_fujitsu_ext
280 if ( use_mpi_pc_fujitsu_ext )
then
281 log_error(
"MeshFieldCommBase_prepare_PC",*)
'use_mpi_pc_fujitsu_ext=.true., but Fujitsu MPI is unavailable. Check!'
287 if (this%field_num_tot > 0)
then
288 allocate( this%request_pc(2*this%nfaces_comm*this%mesh%LOCAL_MESH_NUM) )
292 do n=1, this%mesh%LOCAL_MESH_NUM
293 do f=1, this%nfaces_comm
294 call this%commdata_list(f,n)%PC_Init_recv( &
295 this%req_counter, this%request_pc(:), &
296 this%obj_ind, this%use_mpi_pc_fujitsu_ext )
299 do n=1, this%mesh%LOCAL_MESH_NUM
300 do f=1, this%nfaces_comm
301 call this%commdata_list(f,n)%PC_Init_send( &
302 this%req_counter, this%request_pc(:), &
303 this%obj_ind, this%use_mpi_pc_fujitsu_ext )
308 end subroutine meshfieldcommbase_prepare_pc
318 fjmpi_prequest_startall
326 type(
localmeshcommdata),
intent(inout),
target :: commdata_list(this%nfaces_comm, this%mesh%local_mesh_num)
327 logical,
intent(in),
optional :: do_wait
337 if (
present(do_wait) )
then
345 if ( this%MPI_pc_flag )
then
348 if ( this%use_mpi_pc_fujitsu_ext )
then
351 call fjmpi_prequest_startall( this%req_counter, this%request_pc(1:this%req_counter), ierr )
354 call mpi_startall( this%req_counter, this%request_pc(1:this%req_counter), ierr )
358 do n=1, this%mesh%LOCAL_MESH_NUM
359 do f=1, this%nfaces_comm
360 lcommdata => commdata_list(f,n)
361 if ( lcommdata%s_rank == lcommdata%lcmesh%PRC_myrank )
then
362 commdata_list(abs(lcommdata%s_faceID), lcommdata%s_tilelocalID)%recv_buf(:,:) &
363 = lcommdata%send_buf(:,:)
370 do n=1, this%mesh%LOCAL_MESH_NUM
371 do f=1, this%nfaces_comm
372 call commdata_list(f,n)%SendRecv( &
373 this%req_counter, this%request_send(:), this%request_recv(:), &
383 this%call_wait_flag_sub_get = .false.
385 this%call_wait_flag_sub_get = .true.
396 field_list, dim, varid_s, lcmesh_list )
404 type(
localmeshcommdata),
intent(inout),
target :: commdata_list(this%nfaces_comm, this%mesh%local_mesh_num)
406 integer,
intent(in),
optional :: dim
407 integer,
intent(in),
optional :: varid_s
408 class(
localmeshbase),
intent(in),
optional,
target :: lcmesh_list(:)
411 integer :: stat_send(mpi_status_size, this%req_counter)
412 integer :: stat_recv(mpi_status_size, this%req_counter)
413 integer :: stat_pc(mpi_status_size, this%req_counter)
417 integer :: irs(this%nfaces_comm,this%mesh%local_mesh_num), ire(this%nfaces_comm,this%mesh%local_mesh_num)
420 integer :: val_size(this%mesh%local_mesh_num)
424 if ( this%MPI_pc_flag )
then
425 if (this%req_counter > 0)
then
426 call mpi_waitall( this%req_counter, this%request_pc(1:this%req_counter), stat_pc, ierr )
429 if ( this%req_counter > 0 )
then
430 call mpi_waitall( this%req_counter, this%request_recv(1:this%req_counter), stat_recv, ierr )
431 call mpi_waitall( this%req_counter, this%request_send(1:this%req_counter), stat_send, ierr )
440 if (
present(field_list) )
then
441 do n=1, this%mesh%LOCAL_MESH_NUM
442 lcmesh => lcmesh_list(n)
443 val_size(n) = lcmesh%refElem%Np * lcmesh%NeA
444 irs(1,n) = lcmesh%refElem%Np * lcmesh%Ne + 1
445 do f=1, this%nfaces_comm
446 ire(f,n) = irs(f,n) + commdata_list(f,n)%Nnode_LCMeshFace - 1
447 if (f<this%nfaces_comm) irs(f+1,n) = ire(f,n) + 1
452 do n=1, this%mesh%LOCAL_MESH_NUM
453 do i=1,
size(field_list)
454 do f=1, this%nfaces_comm
455 var_id = varid_s + i - 1
457 call set_bounddata( field_list(var_id)%field1d%local(n)%val, val_size(n), irs(f,n), ire(f,n), commdata_list(f,n)%recv_buf(:,var_id) )
458 else if (dim==2)
then
459 call set_bounddata( field_list(var_id)%field2d%local(n)%val, val_size(n), irs(f,n), ire(f,n), commdata_list(f,n)%recv_buf(:,var_id) )
460 else if (dim==3)
then
461 call set_bounddata( field_list(var_id)%field3d%local(n)%val, val_size(n), irs(f,n), ire(f,n), commdata_list(f,n)%recv_buf(:,var_id) )
467 do n=1, this%mesh%LOCAL_MESH_NUM
469 do f=1, this%nfaces_comm
470 ire(f,n) = irs(f,n) + commdata_list(f,n)%Nnode_LCMeshFace - 1
471 if (f<this%nfaces_comm) irs(f+1,n) = ire(f,n) + 1
476 do n=1, this%mesh%LOCAL_MESH_NUM
477 do var_id=1, this%field_num_tot
478 do f=1, this%nfaces_comm
479 this%recv_buf(irs(f,n):ire(f,n),var_id,n) = commdata_list(f,n)%recv_buf(:,var_id)
490 integer,
intent(in) :: IA
491 real(RP),
intent(inout) :: var(IA)
492 integer,
intent(in) :: irs_, ire_
493 real(RP),
intent(in) :: recv_buf(ire_-irs_+1)
495 var(irs_:ire_) = recv_buf(:)
507 real(rp),
intent(in) :: var(refelem%np * mesh%nea)
508 real(rp),
intent(out) :: buf(size(mesh%vmapb))
515 buf(i) = var(mesh%VmapB(i))
525 integer,
intent(in) :: varid_s
526 integer,
intent(in) :: dim
528 real(rp),
intent(out) :: buf(size(lcmesh_list(1)%vmapb),size(field_list),size(lcmesh_list))
535 do n=1,
size(lcmesh_list)
536 lcmesh => lcmesh_list(n)
538 do while( i <=
size(field_list) )
539 varid = varid_s + i - 1
542 call extract_bounddata_var2( buf(:,varid,n), buf(:,varid+1,n), field_list(varid)%field1d%local(n)%val, field_list(varid+1)%field1d%local(n)%val, lcmesh, lcmesh%refElem )
544 call extract_bounddata_var2( buf(:,varid,n), buf(:,varid+1,n), field_list(varid)%field2d%local(n)%val, field_list(varid+1)%field2d%local(n)%val, lcmesh, lcmesh%refElem )
546 call extract_bounddata_var2( buf(:,varid,n), buf(:,varid+1,n), field_list(varid)%field3d%local(n)%val, field_list(varid+1)%field3d%local(n)%val, lcmesh, lcmesh%refElem )
564 subroutine extract_bounddata_var2( buf1_, buf2_, var1, var2, lmesh, elem )
568 real(rp),
intent(out) :: buf1_(size(lmesh%vmapb))
569 real(rp),
intent(out) :: buf2_(size(lmesh%vmapb))
570 real(rp),
intent(inout) :: var1(elem%np*lmesh%nea)
571 real(rp),
intent(inout) :: var2(elem%np*lmesh%nea)
578 buf1_(ii) = var1(lmesh%vmapB(ii))
579 buf2_(ii) = var2(lmesh%vmapB(ii))
582 end subroutine extract_bounddata_var2
591 real(rp),
intent(in) :: buf(size(mesh%vmapb))
592 real(rp),
intent(inout) :: var(refelem%np * mesh%nea)
595 var(refelem%Np*mesh%NeE+1:refelem%Np*mesh%NeE+
size(buf)) = buf(:)
601 subroutine localmeshcommdata_init( this, comm, lcmesh, faceID, Nnode_LCMeshFace )
607 integer,
intent(in) :: faceid
608 integer,
intent(in) :: nnode_lcmeshface
611 this%lcmesh => lcmesh
612 this%Nnode_LCMeshFace = nnode_lcmeshface
614 allocate( this%send_buf(nnode_lcmeshface, comm%field_num_tot))
615 allocate( this%recv_buf(nnode_lcmeshface, comm%field_num_tot))
617 this%s_tileID = comm%mesh%tileID_globalMap(faceid, lcmesh%tileID)
618 this%s_faceID = comm%mesh%tileFaceID_globalMap(faceid, lcmesh%tileID)
619 this%s_panelID = comm%mesh%tilePanelID_globalMap(faceid, lcmesh%tileID)
620 this%s_rank = comm%mesh%PRCrank_globalMap(this%s_tileID)
621 this%s_tilelocalID = comm%mesh%tileID_global2localMap(this%s_tileID)
625 end subroutine localmeshcommdata_init
627 subroutine localmeshcommdata_sendrecv( this, &
628 req_counter, req_send, req_recv, &
631 use scale_prc,
only: &
639 integer,
intent(inout) :: req_counter
640 integer,
intent(inout) :: req_send(:)
641 integer,
intent(inout) :: req_recv(:)
649 if ( this%s_rank /= this%lcmesh%PRC_myrank )
then
651 req_counter = req_counter + 1
653 tag = 10 * this%lcmesh%tileID + this%faceID
654 bufsize =
size(this%recv_buf)
655 call mpi_irecv( this%recv_buf(1,1), bufsize, mpi_double_precision, &
656 this%s_rank, tag, prc_local_comm_world, &
657 req_recv(req_counter), ierr)
659 bufsize =
size(this%send_buf)
660 tag = 10 * this%s_tileID + abs(this%s_faceID)
661 call mpi_isend( this%send_buf(1,1), bufsize, mpi_double_precision, &
662 this%s_rank, tag, prc_local_comm_world, &
663 req_send(req_counter), ierr )
665 else if ( this%s_rank == this%lcmesh%PRC_myrank )
then
666 lccommdat_list(abs(this%s_faceID), this%s_tilelocalID)%recv_buf(:,:) &
671 end subroutine localmeshcommdata_sendrecv
673 subroutine localmeshcommdata_pc_init_send( this, &
674 req_counter, req, obj_ind_ , &
675 use_mpi_pc_fujisu_ext )
677 use scale_prc,
only: &
684 fjmpi_prequest_send_init
689 integer,
intent(inout) :: req_counter
690 integer,
intent(inout) :: req(:)
691 integer,
intent(in) :: obj_ind_
692 logical,
intent(in) :: use_mpi_pc_fujisu_ext
699 if ( this%s_rank /= this%lcmesh%PRC_myrank )
then
700 req_counter = req_counter + 1
702 bufsize =
size(this%send_buf)
706 if ( use_mpi_pc_fujisu_ext )
then
708 call fjmpi_prequest_send_init( this%send_buf(1,1), bufsize, mpi_double_precision, &
709 this%s_rank, tag, prc_local_comm_world, &
710 req(req_counter), ierr )
713 call mpi_send_init( this%send_buf(1,1), bufsize, mpi_double_precision, &
714 this%s_rank, tag, prc_local_comm_world, &
715 req(req_counter), ierr )
720 end subroutine localmeshcommdata_pc_init_send
721 subroutine localmeshcommdata_pc_init_recv( this, &
722 req_counter, req, obj_ind_, &
723 use_mpi_pc_fujisu_ext )
725 use scale_prc,
only: &
732 fjmpi_prequest_recv_init
737 integer,
intent(inout) :: req_counter
738 integer,
intent(inout) :: req(:)
739 integer,
intent(in) :: obj_ind_
740 logical,
intent(in) :: use_mpi_pc_fujisu_ext
747 if ( this%s_rank /= this%lcmesh%PRC_myrank )
then
748 req_counter = req_counter + 1
750 bufsize =
size(this%recv_buf)
754 if ( use_mpi_pc_fujisu_ext )
then
756 call fjmpi_prequest_recv_init( this%recv_buf(1,1), bufsize, mpi_double_precision, &
757 this%s_rank, tag, prc_local_comm_world, &
758 req(req_counter), ierr )
761 call mpi_recv_init( this%recv_buf(1,1), bufsize, mpi_double_precision, &
762 this%s_rank, tag, prc_local_comm_world, &
763 req(req_counter), ierr )
768 end subroutine localmeshcommdata_pc_init_recv
770 subroutine localmeshcommdata_final( this )
776 deallocate( this%send_buf )
777 deallocate( this%recv_buf )
780 end subroutine localmeshcommdata_final
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)
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_extract_bounddata_2(field_list, dim, varid_s, lcmesh_list, buf)
Extract halo data from data array with MeshField object and set it to the recieving 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.
subroutine set_bounddata(var, ia, irs_, ire_, recv_buf)
Derived type representing an arbitrary finite element.
Derived type to manage a local computational domain (base type)
Derived type representing a field with 1D mesh.
Derived type representing a field with 2D mesh.
Derived type representing a field with 3D mesh.
Derived type representing a field (base type)
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.