FE-Project
Loading...
Searching...
No Matches
scale_meshfieldcomm_cubedom3d.F90
Go to the documentation of this file.
1!-------------------------------------------------------------------------------
2!> module FElib / Data / Communication 3D cubic domain
3!!
4!! @par Description
5!! A module to manage data communication with 3D cubic domain 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
22 use scale_meshfieldcomm_base, only: &
30 use scale_prof
31 !-----------------------------------------------------------------------------
32 implicit none
33 private
34
35 !-----------------------------------------------------------------------------
36 !
37 !++ Public type & procedure
38 !
39
40 !> Base derived type to manage data communication with 3D cubic domain
41 type, public, extends(meshfieldcommbase) :: meshfieldcommcubedom3d
42 class(meshcubedom3d), pointer :: mesh3d !< Pointer to an object representing 3D cubic computational mesh
43 contains
44 procedure, public :: init => meshfieldcommcubedom3d_init
45 procedure, public :: put => meshfieldcommcubedom3d_put
46 procedure, public :: get => meshfieldcommcubedom3d_get
47 procedure, public :: exchange => meshfieldcommcubedom3d_exchange
48 procedure, public :: final => meshfieldcommcubedom3d_final
50
51 !-----------------------------------------------------------------------------
52 !
53 !++ Public parameters & variables
54 !
55
56 !-----------------------------------------------------------------------------
57 !
58 !++ Private procedure
59 !
60 private :: push_localsendbuf
61
62 !-----------------------------------------------------------------------------
63 !
64 !++ Private parameters & variables
65 !
66 integer :: bufsize_per_field !< Buffer size per a field
67 integer, parameter :: comm_face_num = 6 !< Number of faces with data communication
68
69contains
70!> Initialize an object to manage data communication with 3D cubic domain
71 subroutine meshfieldcommcubedom3d_init( this, &
72 sfield_num, hvfield_num, htensorfield_num, mesh3d )
73
74 implicit none
75
76 class(meshfieldcommcubedom3d), intent(inout) :: this
77 integer, intent(in) :: sfield_num !< Number of scalar fields
78 integer, intent(in) :: hvfield_num !< Number of horizontal vector fields
79 integer, intent(in) :: htensorfield_num !< Number of horizontal vector fields
80 class(meshcubedom3d), intent(in), target :: mesh3d !< Object to manage a 3D cubic computational mesh
81
82 type(localmesh3d), pointer :: lcmesh
83 type(elementbase3d), pointer :: elem
84 integer :: n
85 integer :: nnode_lcmeshface(comm_face_num,mesh3d%local_mesh_num)
86 !-----------------------------------------------------------------------------
87
88 this%mesh3d => mesh3d
89 lcmesh => mesh3d%lcmesh_list(1)
90 elem => lcmesh%refElem3D
91 bufsize_per_field = 2*(lcmesh%NeX + lcmesh%NeY)*lcmesh%NeZ*elem%Nfp_h &
92 + 2*lcmesh%NeX*lcmesh%NeY*elem%Nfp_v
93
94 do n=1, this%mesh3d%LOCAL_MESH_NUM
95 lcmesh => this%mesh3d%lcmesh_list(n)
96 nnode_lcmeshface(:,n) = &
97 (/ lcmesh%NeX, lcmesh%NeY, lcmesh%NeX, lcmesh%NeY, 0, 0 /) * lcmesh%NeZ * lcmesh%refElem3D%Nfp_h &
98 + (/ 0, 0, 0, 0, 1, 1 /) * lcmesh%NeX*lcmesh%NeY * lcmesh%refElem3D%Nfp_v
99 end do
100
101 call meshfieldcommbase_init( this, sfield_num, hvfield_num, htensorfield_num, bufsize_per_field, comm_face_num, nnode_lcmeshface, mesh3d )
102
103 return
104 end subroutine meshfieldcommcubedom3d_init
105
106!> Finalize an object to manage data communication with 3D cubic domain
107 subroutine meshfieldcommcubedom3d_final( this )
108 implicit none
109 class(meshfieldcommcubedom3d), intent(inout) :: this
110 !-----------------------------------------------------------------------------
111
112 call meshfieldcommbase_final( this )
113
114 return
115 end subroutine meshfieldcommcubedom3d_final
116
117!> Put field data into temporary buffers
118 subroutine meshfieldcommcubedom3d_put(this, field_list, varid_s)
119 implicit none
120 class(meshfieldcommcubedom3d), intent(inout) :: this
121 type(meshfieldcontainer), intent(in) :: field_list(:) !< Array of objects with 3D mesh field
122 integer, intent(in) :: varid_s !< Start index with variables when field_list(1) is written to buffers for data communication
123
124 ! integer :: i
125 ! integer :: n
126 ! type(LocalMesh3d), pointer :: lcmesh
127 ! integer :: field_num
128 !-----------------------------------------------------------------------------
129
130! call PROF_rapstart( 'meshfiled_comm_put', 3)
131 ! field_num = size(field_list)
132 ! do n=1, this%mesh%LOCAL_MESH_NUM
133 ! lcmesh => this%mesh3d%lcmesh_list(n)
134 ! do i=1, field_num
135 ! call MeshFieldCommBase_extract_bounddata( field_list(i)%field3d%local(n)%val, lcmesh%refElem, lcmesh, & ! (in)
136 ! this%send_buf(:,varid_s+i-1,n) ) ! (out)
137 ! end do
138 ! end do
139 call meshfieldcommbase_extract_bounddata_2( field_list, 3, varid_s, this%mesh3d%lcmesh_list, this%send_buf )
140! call PROF_rapend( 'meshfiled_comm_put', 3)
141
142 return
143 end subroutine meshfieldcommcubedom3d_put
144
145!> Extract field data from temporary buffers
146 subroutine meshfieldcommcubedom3d_get(this, field_list, varid_s)
147 use scale_meshfieldcomm_base, only: &
149 implicit none
150
151 class(meshfieldcommcubedom3d), intent(inout) :: this
152 type(meshfieldcontainer), intent(inout) :: field_list(:) !< Array of objects with 3D mesh field
153 integer, intent(in) :: varid_s !< Start index with variables when field_list(1) is written to buffers for data communication
154
155 integer :: i
156 integer :: n
157 type(localmesh3d), pointer :: lcmesh
158 !-----------------------------------------------------------------------------
159
160 if ( this%call_wait_flag_sub_get ) then
161 ! call PROF_rapstart( 'meshfiled_comm_wait_get', 2)
162 call meshfieldcommbase_wait_core( this, this%commdata_list, &
163 field_list, 3, varid_s, this%mesh3d%lcmesh_list )
164 ! call PROF_rapend( 'meshfiled_comm_wait_get', 2)
165 else
166 ! call PROF_rapstart( 'meshfiled_comm_get', 2)
167 do i=1, size(field_list)
168 do n=1, this%mesh3d%LOCAL_MESH_NUM
169 lcmesh => this%mesh3d%lcmesh_list(n)
170 call meshfieldcommbase_set_bounddata( this%recv_buf(:,varid_s+i-1,n), lcmesh%refElem, lcmesh, & !(in)
171 field_list(i)%field3d%local(n)%val ) !(out)
172 end do
173 end do
174 ! call PROF_rapend( 'meshfiled_comm_get', 2)
175 end if
176
177 return
178 end subroutine meshfieldcommcubedom3d_get
179
180!> Exchange field data between neighboring MPI processes
181!!
182!! @param do_wait Flag whether MPI_waitall is called and move tmp data of LocalMeshCommData object to a recv buffer
183!OCL SERIAL
184 subroutine meshfieldcommcubedom3d_exchange( this, do_wait )
185 use scale_meshfieldcomm_base, only: &
188 use scale_prof
189 implicit none
190
191 class(meshfieldcommcubedom3d), intent(inout), target :: this
192 logical, intent(in), optional :: do_wait
193
194 integer :: n, f
195 type(localmesh3d), pointer :: lcmesh
196 type(localmeshcommdata), pointer :: commdata
197 !-----------------------------------------------------------------------------
198
199! call PROF_rapstart( 'meshfiled_comm_ex_push_buf', 3)
200 do n=1, this%mesh%LOCAL_MESH_NUM
201 lcmesh => this%mesh3d%lcmesh_list(n)
202 do f=1, this%nfaces_comm
203 commdata => this%commdata_list(f,n)
204 call push_localsendbuf( commdata%send_buf(:,:), & ! (inout)
205 this%send_buf(:,:,n), commdata%s_faceID, this%is_f(f,n), & ! (in)
206 commdata%Nnode_LCMeshFace, this%field_num_tot, lcmesh ) ! (in)
207 end do
208 end do
209! call PROF_rapend( 'meshfiled_comm_ex_push_buf', 3)
210
211 !-----------------------
212
213 call meshfieldcommbase_exchange_core(this, this%commdata_list(:,:), do_wait )
214
215 return
216 end subroutine meshfieldcommcubedom3d_exchange
217
218!----------------------------
219
220!OCL SERIAL
221 subroutine push_localsendbuf( lc_send_buf, send_buf, s_faceID, is, Nnode_LCMeshFace, var_num, lcmesh )
222 implicit none
223
224 integer, intent(in) :: nnode_lcmeshface
225 integer, intent(in) :: var_num
226 type(localmesh3d), intent(in) :: lcmesh
227 real(rp), intent(out) :: lc_send_buf(nnode_lcmeshface,var_num)
228 real(rp), intent(in) :: send_buf(bufsize_per_field,var_num)
229 integer, intent(in) :: s_faceid, is
230
231 integer :: ie
232 integer :: vid, i
233 class(elementbase3d), pointer :: e3d
234 !-----------------------------------------------------------------------------
235
236 if ( s_faceid > 0 ) then
237 !$omp parallel do
238 do vid=1, var_num
239 do i=1, nnode_lcmeshface
240 lc_send_buf(i,vid) = send_buf(is+i-1,vid)
241 end do
242 end do
243 else if ( -5 < s_faceid .and. s_faceid < 0) then
244 e3d => lcmesh%refElem3D
245 ie = is + nnode_lcmeshface - 1
246 call revert_hori( lc_send_buf(:,:), send_buf(is:ie,:), nnode_lcmeshface/lcmesh%NeZ, lcmesh%NeZ )
247 end if
248
249 return
250 contains
251 subroutine revert_hori(revert, ori, Ne_h1D, NeZ)
252 integer, intent(in) :: ne_h1d, nez
253 real(rp), intent(out) :: revert(e3d%nnode_h1d, e3d%nnode_v, ne_h1d, nez, var_num)
254 real(rp), intent(in) :: ori(e3d%nnode_h1d, e3d%nnode_v, ne_h1d, nez, var_num)
255
256 integer :: p1, p3, i, k, n
257 integer :: i_, p1_
258 !-----------------------------------------------------------------------------
259
260 do n=1, var_num
261 do k=1, nez
262 do i=1, ne_h1d
263 i_ = ne_h1d - i + 1
264 do p3=1, e3d%Nnode_v
265 do p1=1, e3d%Nnode_h1D
266 p1_ = e3d%Nnode_h1D - p1 + 1
267 revert(p1,p3,i,k,n) = ori(p1_,p3,i_,k,n)
268 end do
269 end do
270 end do
271 end do
272 end do
273 end subroutine revert_hori
274 end subroutine push_localsendbuf
module FElib / Element / Base
module FElib / Mesh / Local 3D
module FElib / Mesh / Cubic 3D domain
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.
module FElib / Data / Communication 3D cubic domain
integer bufsize_per_field
Buffer size per a field.
Derived type representing a 3D reference element.
Derived type representing an arbitrary finite element.
Derived type to manage a local 3D computational domain.
Derived type representing a field with 3D mesh.
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.
Base derived type to manage data communication with 3D cubic domain.