FE-Project
Loading...
Searching...
No Matches
scale_meshfieldcomm_base.F90
Go to the documentation of this file.
1!-------------------------------------------------------------------------------
9#include "scaleFElib.h"
11
12 !-----------------------------------------------------------------------------
13 !
14 !++ used modules
15 !
16 use scale_precision
17 use scale_io
18
20 use scale_mesh_base, only: meshbase
21 use scale_meshfield_base, only: &
23 use scale_localmesh_base, only: &
25
26 !-----------------------------------------------------------------------------
27 implicit none
28 private
29
30 !-----------------------------------------------------------------------------
31 !
32 !++ Public type & procedure
33 !
34
36 type, public :: localmeshcommdata
37 class(localmeshbase), pointer :: lcmesh
38 real(rp), allocatable :: send_buf(:,:)
39 real(rp), allocatable :: recv_buf(:,:)
40 integer :: nnode_lcmeshface
41 integer :: s_tileid
42 integer :: s_faceid
43 integer :: s_panelid
44 integer :: s_rank
45 integer :: s_tilelocalid
46 integer :: faceid
47 contains
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
52 end type
53
55 type, abstract, public :: meshfieldcommbase
56 integer :: sfield_num
57 integer :: hvfield_num
58 integer :: htensorfield_num
59 integer :: field_num_tot
60 integer :: nfaces_comm
61
62 class(meshbase), pointer :: mesh
63 real(rp), allocatable :: send_buf(:,:,:)
64 real(rp), allocatable :: recv_buf(:,:,:)
65 integer, allocatable :: request_send(:)
66 integer, allocatable :: request_recv(:)
67
68 type(localmeshcommdata), allocatable :: commdata_list(:,:)
69 integer, allocatable :: is_f(:,:)
70
71 logical :: mpi_pc_flag
72 integer, allocatable :: request_pc(:)
73
74 integer :: req_counter
75 logical :: call_wait_flag_sub_get
76 contains
77 procedure(meshfieldcommbase_put), public, deferred :: put
78 procedure(meshfieldcommbase_get), public, deferred :: get
79 procedure(meshfieldcommbase_exchange), public, deferred :: exchange
80 procedure, public :: prepare_pc => meshfieldcommbase_prepare_pc
81 end type meshfieldcommbase
82
88
89 type, public :: meshfieldcontainer
90 class(meshfield1d), pointer :: field1d
91 class(meshfield2d), pointer :: field2d
92 class(meshfield3d), pointer :: field3d
93 end type
94
95 interface
96 subroutine meshfieldcommbase_put(this, field_list, varid_s)
99 class(meshfieldcommbase), intent(inout) :: this
100 type(meshfieldcontainer), intent(in) :: field_list(:)
101 integer, intent(in) :: varid_s
102 end subroutine meshfieldcommbase_put
103
104 subroutine meshfieldcommbase_get(this, field_list, varid_s)
105 import meshfieldcommbase
106 import meshfieldcontainer
107 class(meshfieldcommbase), intent(inout) :: this
108 type(meshfieldcontainer), intent(inout) :: field_list(:)
109 integer, intent(in) :: varid_s
110 end subroutine meshfieldcommbase_get
111
112 subroutine meshfieldcommbase_exchange(this, do_wait)
113 import meshfieldcommbase
114 class(meshfieldcommbase), intent(inout), target :: this
115 logical, intent(in), optional :: do_wait
116 end subroutine meshfieldcommbase_exchange
117
118 subroutine meshfieldcomm_final(this)
119 import meshfieldcommbase
120 class(meshfieldcommbase), intent(inout) :: this
121 end subroutine meshfieldcomm_final
122 end interface
123
124 !-----------------------------------------------------------------------------
125 !
126 !++ Public parameters & variables
127 !
128
129 !-----------------------------------------------------------------------------
130 !
131 !++ Private procedure
132 !
133
134 !-----------------------------------------------------------------------------
135 !
136 !++ Private parameters & variables
137 !
138
139contains
140
149!OCL SERIAL
150 subroutine meshfieldcommbase_init( this, &
151 sfield_num, hvfield_num, htensorfield_num, bufsize_per_field, comm_face_num, &
152 Nnode_LCMeshFace, mesh )
153
154 implicit none
155
156 class(meshfieldcommbase), intent(inout) :: this
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)
164
165 integer :: n, f
166 class(localmeshbase), pointer :: lcmesh
167 !-----------------------------------------------------------------------------
168
169 this%mesh => mesh
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
175
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) )
181
182 allocate( this%commdata_list(comm_face_num,mesh%LOCAL_MESH_NUM) )
183 allocate( this%is_f(comm_face_num,mesh%LOCAL_MESH_NUM) )
184
185 do n=1, mesh%LOCAL_MESH_NUM
186 this%is_f(1,n) = 1
187 do f=2, this%nfaces_comm
188 this%is_f(f,n) = this%is_f(f-1,n) + nnode_lcmeshface(f-1,n)
189 end do
190
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) )
194 end do
195 end do
196 end if
197
198 this%MPI_pc_flag = .false.
199
200 return
201 end subroutine meshfieldcommbase_init
202
205!OCL SERIAL
206 subroutine meshfieldcommbase_final( this )
207 implicit none
208
209 class(meshfieldcommbase), intent(inout) :: this
210
211 integer :: n, f
212 integer :: ireq
213 integer :: ierr
214 !-----------------------------------------------------------------------------
215
216 if (this%field_num_tot > 0) then
217 deallocate( this%send_buf, this%recv_buf )
218 deallocate( this%request_send, this%request_recv )
219
220 do n=1, this%mesh%LOCAL_MESH_NUM
221 do f=1, this%nfaces_comm
222 call this%commdata_list(f,n)%Final()
223 end do
224 end do
225 deallocate( this%commdata_list, this%is_f )
226
227 if (this%MPI_pc_flag) then
228 do ireq=1, this%req_counter
229 call mpi_request_free( this%request_pc(ireq), ierr )
230 end do
231 deallocate( this%request_pc )
232 end if
233 end if
234
235 return
236 end subroutine meshfieldcommbase_final
237
240!OCL SERIAL
241 subroutine meshfieldcommbase_prepare_pc( this )
242 implicit none
243 class(meshfieldcommbase), intent(inout) :: this
244
245 integer :: n, f
246 !-----------------------------------------------------------------------------
247
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) )
251 end if
252
253 this%req_counter = 0
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(:) ) ! (inout)
258 end do
259 end do
260
261 return
262 end subroutine meshfieldcommbase_prepare_pc
263
268!OCL SERIAL
269 subroutine meshfieldcommbase_exchange_core( this, commdata_list, do_wait )
270! use scale_prof
271 implicit none
272
273 class(meshfieldcommbase), intent(inout) :: this
274 type(localmeshcommdata), intent(inout), target :: commdata_list(this%nfaces_comm, this%mesh%local_mesh_num)
275 logical, intent(in), optional :: do_wait
276
277 integer :: n, f
278 integer :: ierr
279 logical :: do_wait_
280 !-----------------------------------------------------------------------------
281
282 if ( present(do_wait) ) then
283 do_wait_ = do_wait
284 else
285 do_wait_ = .true.
286 end if
287
288! call PROF_rapstart( 'meshfiled_comm_ex_core', 3)
289 if ( this%MPI_pc_flag ) then
290 call mpi_startall( this%req_counter, this%request_pc(1:this%req_counter), ierr )
291 else
292 this%req_counter = 0
293 end if
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(:), & ! (inout)
298 commdata_list(:,:), & ! (inout)
299 this%MPI_pc_flag ) ! (in)
300 end do
301 end do
302! call PROF_rapend( 'meshfiled_comm_ex_core', 3)
303
304 if ( do_wait_ ) then
305 call meshfieldcommbase_wait_core( this, commdata_list )
306 this%call_wait_flag_sub_get = .false.
307 else
308 this%call_wait_flag_sub_get = .true.
309 end if
310
311 return
313
317!OCL SERIAL
318 subroutine meshfieldcommbase_wait_core( this, commdata_list )
319 use mpi, only: &
320 !MPI_waitall, &
321 mpi_status_size
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
328 integer :: ierr
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)
332
333 integer :: n, f
334 integer :: var_id
335 integer :: irs, ire
336 !----------------------------
337
338! call PROF_rapstart( 'meshfiled_comm_wait_core', 3)
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 )
342 end if
343 else
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 )
347 end if
348 end if
349! call PROF_rapend( 'meshfiled_comm_wait_core', 3)
350
351 !---------------------
352
353! call PROF_rapstart( 'meshfiled_comm_wait_post', 3)
354 do n=1, this%mesh%LOCAL_MESH_NUM
355 !$omp parallel do private(var_id,f,irs,ire)
356 do var_id=1, this%field_num_tot
357 irs = 1
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)
361 irs = ire + 1
362 end do ! end loop for face
363 end do
364 end do
365! call PROF_rapend( 'meshfiled_comm_wait_post', 3)
366
367 return
368 end subroutine meshfieldcommbase_wait_core
369
370!OCL SERIAL
371 subroutine meshfieldcommbase_extract_bounddata(var, refElem, mesh, buf)
372 implicit none
373
374 class(elementbase), intent(in) :: refelem
375 class(localmeshbase), intent(in) :: mesh
376 real(rp), intent(in) :: var(refelem%np * mesh%nea)
377 real(rp), intent(out) :: buf(size(mesh%vmapb))
378
379 integer :: i
380 !-----------------------------------------------------------------------------
381 !$omp parallel do
382!OCL PREFETCH
383 do i=1, size(buf)
384 buf(i) = var(mesh%VmapB(i))
385 end do
386 return
388
389 subroutine meshfieldcommbase_set_bounddata(buf, refElem, mesh, var)
390 implicit none
391
392 class(elementbase), intent(in) :: refelem
393 class(localmeshbase), intent(in) :: mesh
394 real(rp), intent(in) :: buf(size(mesh%vmapb))
395 real(rp), intent(inout) :: var(refelem%np * mesh%nea)
396 !-----------------------------------------------------------------------------
397
398 var(refelem%Np*mesh%NeE+1:refelem%Np*mesh%NeE+size(buf)) = buf(:)
399 return
400 end subroutine meshfieldcommbase_set_bounddata
401
402 !-------------------------------------------
403
404 subroutine localmeshcommdata_init( this, comm, lcmesh, faceID, Nnode_LCMeshFace )
405 implicit none
406
407 class(localmeshcommdata), intent(inout) :: this
408 class(meshfieldcommbase), intent(in) :: comm
409 class(localmeshbase), intent(in), target :: lcmesh
410 integer, intent(in) :: faceid
411 integer, intent(in) :: nnode_lcmeshface
412 !-----------------------------------------------------------------------------
413
414 this%lcmesh => lcmesh
415 this%Nnode_LCMeshFace = nnode_lcmeshface
416
417 allocate( this%send_buf(nnode_lcmeshface, comm%field_num_tot))
418 allocate( this%recv_buf(nnode_lcmeshface, comm%field_num_tot))
419
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)
425 this%faceID = faceid
426
427 return
428 end subroutine localmeshcommdata_init
429
430 subroutine localmeshcommdata_sendrecv( this, &
431 req_counter, req_send, req_recv, &
432 lccommdat_list, MPI_pc_flag )
433
434 use scale_prc, only: &
435 prc_local_comm_world
436 use mpi, only: &
437 !MPI_isend, MPI_irecv, &
438 mpi_double_precision
439 implicit none
440
441 class(localmeshcommdata), intent(inout) :: this
442 integer, intent(inout) :: req_counter
443 integer, intent(inout) :: req_send(:)
444 integer, intent(inout) :: req_recv(:)
445 type(localmeshcommdata), intent(inout) :: lccommdat_list(:,:)
446 logical, intent(in) :: mpi_pc_flag
447
448 integer :: tag
449 integer :: bufsize
450 integer :: ierr
451 !-------------------------------------------
452
453 if ( this%s_rank /= this%lcmesh%PRC_myrank &
454 .and. (.not. mpi_pc_flag ) ) then
455
456 req_counter = req_counter + 1
457
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)
463
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 )
469
470 else if ( this%s_rank == this%lcmesh%PRC_myrank ) then
471 lccommdat_list(abs(this%s_faceID), this%s_tilelocalID)%recv_buf(:,:) &
472 = this%send_buf(:,:)
473 end if
474
475 return
476 end subroutine localmeshcommdata_sendrecv
477
478 subroutine localmeshcommdata_pc_init( this, &
479 req_counter, req )
480
481 use scale_prc, only: &
482 prc_local_comm_world
483 use mpi, only: &
484 !MPI_isend, MPI_irecv, &
485 mpi_double_precision
486 implicit none
487
488 class(localmeshcommdata), intent(inout) :: this
489 integer, intent(inout) :: req_counter
490 integer, intent(inout) :: req(:)
491
492 integer :: tag
493 integer :: bufsize
494 integer :: ierr
495 !-------------------------------------------
496
497 if ( this%s_rank /= this%lcmesh%PRC_myrank ) then
498
499 req_counter = req_counter + 1
500
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)
506
507 !--
508 req_counter = req_counter + 1
509
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 )
515 end if
516
517 return
518 end subroutine localmeshcommdata_pc_init
519
520 subroutine localmeshcommdata_final( this )
521 implicit none
522
523 class(localmeshcommdata), intent(inout) :: this
524 !-----------------------------------------------------------------------------
525
526 deallocate( this%send_buf )
527 deallocate( this%recv_buf )
528
529 return
530 end subroutine localmeshcommdata_final
531
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.