FE-Project
Loading...
Searching...
No Matches
scale_meshfieldcomm_base.F90
Go to the documentation of this file.
1!-------------------------------------------------------------------------------
2!> module FElib / Data / Communication base
3!!
4!! @par Description
5!! Base module to manage data communication for element-based methods
6!!
7!! @author Yuta Kawai, Team SCALE
8!<
9#include "scaleFElib.h"
11
12 !-----------------------------------------------------------------------------
13 !
14 !++ used modules
15 !
16 use scale_precision
17 use scale_io
18 use scale_prc, only: &
19 prc_abort
20
22 use scale_mesh_base, only: meshbase
23 use scale_meshfield_base, only: &
25 use scale_localmesh_base, only: &
27
28 !-----------------------------------------------------------------------------
29 implicit none
30 private
31
32 !-----------------------------------------------------------------------------
33 !
34 !++ Public type & procedure
35 !
36
37 !> Derived type to manage data communication at a face between adjacent local meshes
38 type, public :: localmeshcommdata
39 class(localmeshbase), pointer :: lcmesh
40 real(rp), allocatable :: send_buf(:,:) !< Buffer for sending data
41 real(rp), allocatable :: recv_buf(:,:) !< Buffer for receiving data
42 integer :: nnode_lcmeshface !< Number of nodes at a face
43 integer :: s_tileid !< Destination tile ID
44 integer :: s_faceid !< Destination face ID
45 integer :: s_panelid !< Destination panel ID
46 integer :: s_rank !< Destination MPI rank
47 integer :: s_tilelocalid !< Destination local tile ID
48 integer :: faceid !< Own face ID
49 contains
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
55 end type
56
57 !> Base derived type to manage data communication
58 type, abstract, public :: meshfieldcommbase
59 integer :: sfield_num !< Number of scalar fields
60 integer :: hvfield_num !< Number of horizontal vector fields
61 integer :: htensorfield_num !< Number of horizontal tensor fields
62 integer :: field_num_tot !< Total number of fields
63 integer :: nfaces_comm !< Number of faces where halo data is communicated
64
65 class(meshbase), pointer :: mesh
66 real(rp), allocatable :: send_buf(:,:,:) !< Buffer for sending data
67 real(rp), allocatable :: recv_buf(:,:,:) !< Buffer for receiving data
68 integer, allocatable :: request_send(:)
69 integer, allocatable :: request_recv(:)
70
71 type(localmeshcommdata), allocatable :: commdata_list(:,:)
72 integer, allocatable :: is_f(:,:)
73
74 logical :: mpi_pc_flag !< Flag whether persistent communication is used
75 logical :: use_mpi_pc_fujitsu_ext !< Flag whether Fujitsu extension routines are used for persistent communication
76 integer, allocatable :: request_pc(:)
77
78 integer :: req_counter
79 logical :: call_wait_flag_sub_get !< Flag whether MPI_wait need to be called before getting halo data from recv_buf
80
81 integer :: obj_ind
82 contains
83 procedure(meshfieldcommbase_put), public, deferred :: put
84 procedure(meshfieldcommbase_get), public, deferred :: get
85 procedure(meshfieldcommbase_exchange), public, deferred :: exchange
86 procedure, public :: prepare_pc => meshfieldcommbase_prepare_pc
87 end type meshfieldcommbase
88
95
96 !> Container to save a pointer of MeshField(1D, 2D, 3D) object
97 type, public :: meshfieldcontainer
98 class(meshfield1d), pointer :: field1d
99 class(meshfield2d), pointer :: field2d
100 class(meshfield3d), pointer :: field3d
101 end type
102
103 interface
104 subroutine meshfieldcommbase_put(this, field_list, varid_s)
105 import meshfieldcommbase
106 import meshfieldcontainer
107 class(meshfieldcommbase), intent(inout) :: this
108 type(meshfieldcontainer), intent(in) :: field_list(:)
109 integer, intent(in) :: varid_s
110 end subroutine meshfieldcommbase_put
111
112 subroutine meshfieldcommbase_get(this, field_list, varid_s)
113 import meshfieldcommbase
114 import meshfieldcontainer
115 class(meshfieldcommbase), intent(inout) :: this
116 type(meshfieldcontainer), intent(inout) :: field_list(:)
117 integer, intent(in) :: varid_s
118 end subroutine meshfieldcommbase_get
119
120 subroutine meshfieldcommbase_exchange(this, do_wait)
121 import meshfieldcommbase
122 class(meshfieldcommbase), intent(inout), target :: this
123 logical, intent(in), optional :: do_wait
124 end subroutine meshfieldcommbase_exchange
125
126 subroutine meshfieldcomm_final(this)
127 import meshfieldcommbase
128 class(meshfieldcommbase), intent(inout) :: this
129 end subroutine meshfieldcomm_final
130 end interface
131
132 !-----------------------------------------------------------------------------
133 !
134 !++ Public parameters & variables
135 !
136
137 !-----------------------------------------------------------------------------
138 !
139 !++ Private procedure
140 !
141
142 !-----------------------------------------------------------------------------
143 !
144 !++ Private parameters & variables
145 !
146
147 integer :: obj_ind = 0
148 integer, parameter :: OBJ_INDEX_MAX = 8192
149
150contains
151
152!> Initialize a base object to manage data communication of fields
153!!
154!! @param sfield_num Number of scalar fields
155!! @param hvfield_num Number of horizontal vector fields
156!! @param htensorfield_num Number of horizontal tensor fields
157!! @param bufsize_per_field Buffer size per a field
158!! @param comm_face_num Number of faces on a local mesh which perform data communication
159!! @param Nnode_LCMeshFace Array to store the number of nodes
160!OCL SERIAL
161 subroutine meshfieldcommbase_init( this, &
162 sfield_num, hvfield_num, htensorfield_num, bufsize_per_field, comm_face_num, &
163 Nnode_LCMeshFace, mesh )
164
165 implicit none
166
167 class(meshfieldcommbase), intent(inout) :: this
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)
175
176 integer :: n, f
177 class(localmeshbase), pointer :: lcmesh
178 !-----------------------------------------------------------------------------
179
180 this%mesh => mesh
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
186
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) )
192
193 allocate( this%commdata_list(comm_face_num,mesh%LOCAL_MESH_NUM) )
194 allocate( this%is_f(comm_face_num,mesh%LOCAL_MESH_NUM) )
195
196 do n=1, mesh%LOCAL_MESH_NUM
197 this%is_f(1,n) = 1
198 do f=2, this%nfaces_comm
199 this%is_f(f,n) = this%is_f(f-1,n) + nnode_lcmeshface(f-1,n)
200 end do
201
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) )
205 end do
206 end do
207 end if
208
209 this%MPI_pc_flag = .false.
210
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!'
214 call prc_abort
215 end if
216 this%obj_ind = obj_ind
217
218 return
219 end subroutine meshfieldcommbase_init
220
221!> Finalize a base object to manage data communication of fields
222!!
223!OCL SERIAL
224 subroutine meshfieldcommbase_final( this )
225 implicit none
226
227 class(meshfieldcommbase), intent(inout) :: this
228
229 integer :: n, f
230 integer :: ireq
231 integer :: ierr
232 !-----------------------------------------------------------------------------
233
234 if (this%field_num_tot > 0) then
235 deallocate( this%send_buf, this%recv_buf )
236 deallocate( this%request_send, this%request_recv )
237
238 do n=1, this%mesh%LOCAL_MESH_NUM
239 do f=1, this%nfaces_comm
240 call this%commdata_list(f,n)%Final()
241 end do
242 end do
243 deallocate( this%commdata_list, this%is_f )
244
245 if ( this%MPI_pc_flag ) then
246 do ireq=1, this%req_counter
247 call mpi_request_free( this%request_pc(ireq), ierr )
248 end do
249 deallocate( this%request_pc )
250 end if
251 end if
252
253 return
254 end subroutine meshfieldcommbase_final
255
256!> Prepare persistent communication
257!!
258!! @param use_mpi_pc_fujitsu_ext Flag whether the extension routines of Fujitsu MPI are used
259!!
260!OCL SERIAL
261 subroutine meshfieldcommbase_prepare_pc( this, &
262 use_mpi_pc_fujitsu_ext )
263 implicit none
264 class(meshfieldcommbase), intent(inout) :: this
265 logical, intent(in), optional :: use_mpi_pc_fujitsu_ext
266
267 integer :: n, f
268 !-----------------------------------------------------------------------------
269
270 this%MPI_pc_flag = .true.
271
272#ifdef __FUJITSU
273 this%use_mpi_pc_fujitsu_ext = .true.
274#else
275 this%use_mpi_pc_fujitsu_ext = .false.
276#endif
277 if ( present(use_mpi_pc_fujitsu_ext) ) then
278 this%use_mpi_pc_fujitsu_ext = use_mpi_pc_fujitsu_ext
279#ifndef __FUJITSU
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!'
282 call prc_abort
283 end if
284#endif
285 end if
286
287 if (this%field_num_tot > 0) then
288 allocate( this%request_pc(2*this%nfaces_comm*this%mesh%LOCAL_MESH_NUM) )
289 end if
290
291 this%req_counter = 0
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(:), & ! (inout)
296 this%obj_ind, this%use_mpi_pc_fujitsu_ext ) ! (in)
297 end do
298 end do
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(:), & ! (inout)
303 this%obj_ind, this%use_mpi_pc_fujitsu_ext ) ! (in)
304 end do
305 end do
306
307 return
308 end subroutine meshfieldcommbase_prepare_pc
309
310!> Exchange halo data
311!!
312!! @param commdata_list Array of LocalMeshCommData objects which manage information and halo data
313!! @param do_wait Flag whether MPI_waitall is called and move tmp data of LocalMeshCommData object to a recv buffer
314!OCL SERIAL
315 subroutine meshfieldcommbase_exchange_core( this, commdata_list, do_wait )
316#ifdef __FUJITSU
317 use mpi_ext, only: &
318 fjmpi_prequest_startall
319#endif
320! use mpi, only: &
321! MPI_startall
322 use scale_prof
323 implicit none
324
325 class(meshfieldcommbase), intent(inout) :: this
326 type(localmeshcommdata), intent(inout), target :: commdata_list(this%nfaces_comm, this%mesh%local_mesh_num)
327 logical, intent(in), optional :: do_wait
328
329 integer :: n, f
330 integer :: ierr
331 logical :: do_wait_
332
333 integer :: i
334 type(localmeshcommdata), pointer :: lcommdata
335 !-----------------------------------------------------------------------------
336
337 if ( present(do_wait) ) then
338 do_wait_ = do_wait
339 else
340 do_wait_ = .true.
341 end if
342
343! call PROF_rapstart( 'meshfiled_comm_ex_core', 3)
344 !
345 if ( this%MPI_pc_flag ) then
346 !$omp parallel
347 !$omp master
348 if ( this%use_mpi_pc_fujitsu_ext ) then
349#ifdef __FUJITSU
350 ! Use Fujitsu MPI extension for persistent communication
351 call fjmpi_prequest_startall( this%req_counter, this%request_pc(1:this%req_counter), ierr )
352#endif
353 else
354 call mpi_startall( this%req_counter, this%request_pc(1:this%req_counter), ierr )
355 end if
356 !$omp end master
357 !$omp do collapse(2) private(lcommdata)
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(:,:)
364 end if
365 end do
366 end do
367 !$omp end parallel
368 else
369 this%req_counter = 0
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(:), & ! (inout)
374 commdata_list(:,:) ) ! (inout)
375 end do
376 end do
377 end if
378
379! call PROF_rapend( 'meshfiled_comm_ex_core', 3)
380
381 if ( do_wait_ ) then
382 call meshfieldcommbase_wait_core( this, commdata_list )
383 this%call_wait_flag_sub_get = .false.
384 else
385 this%call_wait_flag_sub_get = .true.
386 end if
387
388 return
390
391!> Wait data communication and move tmp data of LocalMeshCommData object to a recv buffer
392!!
393!! @param commdata_list Array of LocalMeshCommData objects which manage information and halo data
394!OCL SERIAL
395 subroutine meshfieldcommbase_wait_core( this, commdata_list, &
396 field_list, dim, varid_s, lcmesh_list )
397 use mpi, only: &
398! MPI_waitall, &
399 mpi_status_size
400 use scale_prof
401 implicit none
402
403 class(meshfieldcommbase), intent(inout) :: this
404 type(localmeshcommdata), intent(inout), target :: commdata_list(this%nfaces_comm, this%mesh%local_mesh_num)
405 type(meshfieldcontainer), intent(inout), optional :: field_list(:)
406 integer, intent(in), optional :: dim
407 integer, intent(in), optional :: varid_s
408 class(localmeshbase), intent(in), optional, target :: lcmesh_list(:)
409
410 integer :: ierr
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)
414
415 integer :: n, f
416 integer :: i, var_id
417 integer :: irs(this%nfaces_comm,this%mesh%local_mesh_num), ire(this%nfaces_comm,this%mesh%local_mesh_num)
418
419 class(localmeshbase), pointer :: lcmesh
420 integer :: val_size(this%mesh%local_mesh_num)
421 !----------------------------
422
423! call PROF_rapstart( 'meshfiled_comm_wait_core', 2)
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 )
427 end if
428 else
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 )
432 end if
433 end if
434! call PROF_rapend( 'meshfiled_comm_wait_core', 2)
435
436 !---------------------
437
438! call PROF_rapstart( 'meshfiled_comm_wait_post', 2)
439
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
448 end do
449 end do
450
451 !$omp parallel do private(var_id,n,i,f) collapse(3)
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
456 if (dim==1) then
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) )
462 end if
463 end do ! end loop for face
464 end do
465 end do
466 else
467 do n=1, this%mesh%LOCAL_MESH_NUM
468 irs(1,n) = 1
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
472 end do
473 end do
474
475 !$omp parallel do private(n,var_id,f) collapse(3)
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)
480 end do ! end loop for face
481 end do
482 end do
483 end if
484! call PROF_rapend( 'meshfiled_comm_wait_post', 2)
485 return
486 contains
487!OCL SERIAL
488 subroutine set_bounddata( var, IA, irs_, ire_, recv_buf )
489 implicit none
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)
494 !-----------------------------
495 var(irs_:ire_) = recv_buf(:)
496 return
497 end subroutine set_bounddata
498 end subroutine meshfieldcommbase_wait_core
499
500!> Extract halo data from data array with MeshField object and set it to the recieving buffer
501!OCL SERIAL
502 subroutine meshfieldcommbase_extract_bounddata(var, refElem, mesh, buf)
503 implicit none
504
505 class(elementbase), intent(in) :: refelem
506 class(localmeshbase), intent(in) :: mesh
507 real(rp), intent(in) :: var(refelem%np * mesh%nea)
508 real(rp), intent(out) :: buf(size(mesh%vmapb))
509
510 integer :: i
511 !-----------------------------------------------------------------------------
512 !$omp parallel do
513!OCL PREFETCH
514 do i=1, size(buf)
515 buf(i) = var(mesh%VmapB(i))
516 end do
517 return
519
520!> Extract halo data from data array with MeshField object and set it to the recieving buffer
521!OCL SERIAL
522 subroutine meshfieldcommbase_extract_bounddata_2(field_list, dim, varid_s, lcmesh_list, buf)
523 implicit none
524 class(meshfieldcontainer), intent(in), target :: field_list(:)
525 integer, intent(in) :: varid_s
526 integer, intent(in) :: dim
527 class(localmeshbase), intent(in), target :: lcmesh_list(:)
528 real(rp), intent(out) :: buf(size(lcmesh_list(1)%vmapb),size(field_list),size(lcmesh_list))
529
530 class(localmeshbase), pointer :: lcmesh
531 integer :: varid
532 integer :: i, n
533 !-----------------------------------------------------------------------------
534
535 do n=1, size(lcmesh_list)
536 lcmesh => lcmesh_list(n)
537 i = 1
538 do while( i <= size(field_list) )
539 varid = varid_s + i - 1
540 if ( i+1 <= n ) then
541 if (dim==1) then
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 )
543 else if(dim==2) then
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 )
545 else if(dim==3) then
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 )
547 end if
548 i = i + 2
549 else
550 if (dim==1) then
551 call meshfieldcommbase_extract_bounddata( field_list(varid)%field1d%local(n)%val, lcmesh%refElem, lcmesh, buf(:,varid,n) )
552 else if(dim==2) then
553 call meshfieldcommbase_extract_bounddata( field_list(varid)%field2d%local(n)%val, lcmesh%refElem, lcmesh, buf(:,varid,n) )
554 else if(dim==3) then
555 call meshfieldcommbase_extract_bounddata( field_list(varid)%field3d%local(n)%val, lcmesh%refElem, lcmesh, buf(:,varid,n) )
556 end if
557 i = i + 1
558 end if
559 end do
560 end do
561 return
562 contains
563!OCL SERIAL
564 subroutine extract_bounddata_var2( buf1_, buf2_, var1, var2, lmesh, elem )
565 implicit none
566 class(localmeshbase), intent(in) :: lmesh
567 class(elementbase), intent(in) :: 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)
572
573 integer :: ii
574 !-----------------------------
575 !$omp parallel do
576!OCL PREFETCH
577 do ii=1, size(buf)
578 buf1_(ii) = var1(lmesh%vmapB(ii))
579 buf2_(ii) = var2(lmesh%vmapB(ii))
580 end do
581 return
582 end subroutine extract_bounddata_var2
584
585!> Extract halo data from the recieving buffer and set it to data array with MeshField object
586 subroutine meshfieldcommbase_set_bounddata(buf, refElem, mesh, var)
587 implicit none
588
589 class(elementbase), intent(in) :: refelem
590 class(localmeshbase), intent(in) :: mesh
591 real(rp), intent(in) :: buf(size(mesh%vmapb))
592 real(rp), intent(inout) :: var(refelem%np * mesh%nea)
593 !-----------------------------------------------------------------------------
594
595 var(refelem%Np*mesh%NeE+1:refelem%Np*mesh%NeE+size(buf)) = buf(:)
596 return
597 end subroutine meshfieldcommbase_set_bounddata
598
599 !-------------------------------------------
600
601 subroutine localmeshcommdata_init( this, comm, lcmesh, faceID, Nnode_LCMeshFace )
602 implicit none
603
604 class(localmeshcommdata), intent(inout) :: this
605 class(meshfieldcommbase), intent(in) :: comm
606 class(localmeshbase), intent(in), target :: lcmesh
607 integer, intent(in) :: faceid
608 integer, intent(in) :: nnode_lcmeshface
609 !-----------------------------------------------------------------------------
610
611 this%lcmesh => lcmesh
612 this%Nnode_LCMeshFace = nnode_lcmeshface
613
614 allocate( this%send_buf(nnode_lcmeshface, comm%field_num_tot))
615 allocate( this%recv_buf(nnode_lcmeshface, comm%field_num_tot))
616
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)
622 this%faceID = faceid
623
624 return
625 end subroutine localmeshcommdata_init
626
627 subroutine localmeshcommdata_sendrecv( this, &
628 req_counter, req_send, req_recv, &
629 lccommdat_list )
630
631 use scale_prc, only: &
632 prc_local_comm_world
633 use mpi, only: &
634 !MPI_isend, MPI_irecv, &
635 mpi_double_precision
636 implicit none
637
638 class(localmeshcommdata), intent(inout) :: this
639 integer, intent(inout) :: req_counter
640 integer, intent(inout) :: req_send(:)
641 integer, intent(inout) :: req_recv(:)
642 type(localmeshcommdata), intent(inout) :: lccommdat_list(:,:)
643
644 integer :: tag
645 integer :: bufsize
646 integer :: ierr
647 !-------------------------------------------
648
649 if ( this%s_rank /= this%lcmesh%PRC_myrank ) then
650
651 req_counter = req_counter + 1
652
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)
658
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 )
664
665 else if ( this%s_rank == this%lcmesh%PRC_myrank ) then
666 lccommdat_list(abs(this%s_faceID), this%s_tilelocalID)%recv_buf(:,:) &
667 = this%send_buf(:,:)
668 end if
669
670 return
671 end subroutine localmeshcommdata_sendrecv
672
673 subroutine localmeshcommdata_pc_init_send( this, &
674 req_counter, req, obj_ind_ , &
675 use_mpi_pc_fujisu_ext )
676
677 use scale_prc, only: &
678 prc_local_comm_world
679 use mpi, only: &
680 mpi_double_precision!, &
681! MPI_send_init
682#ifdef __FUJITSU
683 use mpi_ext, only: &
684 fjmpi_prequest_send_init
685#endif
686 implicit none
687
688 class(localmeshcommdata), intent(inout) :: this
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
693
694 integer :: tag
695 integer :: bufsize
696 integer :: ierr
697 !-------------------------------------------
698
699 if ( this%s_rank /= this%lcmesh%PRC_myrank ) then
700 req_counter = req_counter + 1
701
702 bufsize = size(this%send_buf)
703! tag = 1000 * this%s_tileID + 10 * obj_ind_ + abs(this%s_faceID)
704 tag = obj_ind_
705
706 if ( use_mpi_pc_fujisu_ext ) then
707#ifdef __FUJITSU
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 )
711#endif
712 else
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 )
716 end if
717 end if
718
719 return
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 )
724
725 use scale_prc, only: &
726 prc_local_comm_world
727 use mpi, only: &
728 mpi_double_precision!, &
729! MPI_recv_init
730#ifdef __FUJITSU
731 use mpi_ext, only: &
732 fjmpi_prequest_recv_init
733#endif
734 implicit none
735
736 class(localmeshcommdata), intent(inout) :: this
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
741
742 integer :: tag
743 integer :: bufsize
744 integer :: ierr
745 !-------------------------------------------
746
747 if ( this%s_rank /= this%lcmesh%PRC_myrank ) then
748 req_counter = req_counter + 1
749
750 bufsize = size(this%recv_buf)
751! tag = 1000 * this%lcmesh%tileID + 10 * obj_ind_ + this%faceID
752 tag = obj_ind_
753
754 if ( use_mpi_pc_fujisu_ext ) then
755#ifdef __FUJITSU
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 )
759#endif
760 else
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 )
764 end if
765 end if
766
767 return
768 end subroutine localmeshcommdata_pc_init_recv
769
770 subroutine localmeshcommdata_final( this )
771 implicit none
772
773 class(localmeshcommdata), intent(inout) :: this
774 !-----------------------------------------------------------------------------
775
776 deallocate( this%send_buf )
777 deallocate( this%recv_buf )
778
779 return
780 end subroutine localmeshcommdata_final
781
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.