FE-Project
Loading...
Searching...
No Matches
Data Types | Functions/Subroutines
scale_meshfieldcomm_base Module Reference

module FElib / Data / Communication base More...

Data Types

type  localmeshcommdata
 Derived type to manage data communication between a face of local mesh. More...
 
type  meshfieldcommbase
 Base derived type to manage data communication. More...
 
interface  meshfieldcommbase_put
 
type  meshfieldcontainer
 

Functions/Subroutines

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_final (this)
 Finalize a base object to manage data communication of fields.
 
subroutine, public meshfieldcommbase_exchange_core (this, commdata_list, do_wait)
 Exchange halo data.
 
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_extract_bounddata (var, refelem, mesh, buf)
 
subroutine, public meshfieldcommbase_set_bounddata (buf, refelem, mesh, var)
 

Detailed Description

module FElib / Data / Communication base

Description
Base module to manage data communication for element-based methods
Author
Yuta Kawai, Team SCALE

Function/Subroutine Documentation

◆ meshfieldcommbase_init()

subroutine, public scale_meshfieldcomm_base::meshfieldcommbase_init ( class(meshfieldcommbase), intent(inout) this,
integer, intent(in) sfield_num,
integer, intent(in) hvfield_num,
integer, intent(in) htensorfield_num,
integer, intent(in) bufsize_per_field,
integer, intent(in) comm_face_num,
integer, dimension(comm_face_num,mesh%local_mesh_num), intent(in) nnode_lcmeshface,
class(meshbase), intent(in), target mesh )

Initialize a base object to manage data communication of fields.

Parameters
sfield_numNumber of scalar fields
hvfield_numNumber of horizontal vector fields
htensorfield_numNumber of horizontal tensor fields
bufsize_per_fieldBuffer size per a field
comm_face_numNumber of faces on a local mesh which perform data communication
Nnode_LCMeshFaceArray to store the number of nodes

Definition at line 150 of file scale_meshfieldcomm_base.F90.

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

Referenced by scale_meshfieldcomm_base::meshfieldcommbase::prepare_pc().

◆ meshfieldcommbase_final()

subroutine, public scale_meshfieldcomm_base::meshfieldcommbase_final ( class(meshfieldcommbase), intent(inout) this)

Finalize a base object to manage data communication of fields.

Definition at line 206 of file scale_meshfieldcomm_base.F90.

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

Referenced by scale_meshfieldcomm_base::meshfieldcommbase::prepare_pc().

◆ meshfieldcommbase_exchange_core()

subroutine, public scale_meshfieldcomm_base::meshfieldcommbase_exchange_core ( class(meshfieldcommbase), intent(inout) this,
type(localmeshcommdata), dimension(this%nfaces_comm, this%mesh%local_mesh_num), intent(inout), target commdata_list,
logical, intent(in), optional do_wait )

Exchange halo data.

Parameters
commdata_listArray of LocalMeshCommData objects which manage information and halo data
do_waitFlag whether MPI_waitall is called and move tmp data of LocalMeshCommData object to a recv buffer

Definition at line 269 of file scale_meshfieldcomm_base.F90.

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

References meshfieldcommbase_wait_core().

Referenced by scale_meshfieldcomm_base::meshfieldcommbase::prepare_pc().

◆ meshfieldcommbase_wait_core()

subroutine, public scale_meshfieldcomm_base::meshfieldcommbase_wait_core ( class(meshfieldcommbase), intent(inout) this,
type(localmeshcommdata), dimension(this%nfaces_comm, this%mesh%local_mesh_num), intent(inout), target commdata_list )

Wait data communication and move tmp data of LocalMeshCommData object to a recv buffer.

Parameters
commdata_listArray of LocalMeshCommData objects which manage information and halo data

Definition at line 318 of file scale_meshfieldcomm_base.F90.

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

Referenced by meshfieldcommbase_exchange_core(), and scale_meshfieldcomm_base::meshfieldcommbase::prepare_pc().

◆ meshfieldcommbase_extract_bounddata()

subroutine, public scale_meshfieldcomm_base::meshfieldcommbase_extract_bounddata ( real(rp), dimension(refelem%np * mesh%nea), intent(in) var,
class(elementbase), intent(in) refelem,
class(localmeshbase), intent(in) mesh,
real(rp), dimension(size(mesh%vmapb)), intent(out) buf )

Definition at line 371 of file scale_meshfieldcomm_base.F90.

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

Referenced by scale_meshfieldcomm_base::meshfieldcommbase::prepare_pc().

◆ meshfieldcommbase_set_bounddata()

subroutine, public scale_meshfieldcomm_base::meshfieldcommbase_set_bounddata ( real(rp), dimension(size(mesh%vmapb)), intent(in) buf,
class(elementbase), intent(in) refelem,
class(localmeshbase), intent(in) mesh,
real(rp), dimension(refelem%np * mesh%nea), intent(inout) var )

Definition at line 389 of file scale_meshfieldcomm_base.F90.

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

Referenced by scale_meshfieldcomm_base::meshfieldcommbase::prepare_pc().