FE-Project
All Classes Namespaces Files Functions Variables Pages
scale_meshfieldcomm_cubedom3d.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
22 use scale_meshfieldcomm_base, only: &
29 use scale_prof
30 !-----------------------------------------------------------------------------
31 implicit none
32 private
33
34 !-----------------------------------------------------------------------------
35 !
36 !++ Public type & procedure
37 !
38
39
40 type, public, extends(meshfieldcommbase) :: meshfieldcommcubedom3d
41 class(meshcubedom3d), pointer :: mesh3d
42 contains
43 procedure, public :: init => meshfieldcommcubedom3d_init
44 procedure, public :: put => meshfieldcommcubedom3d_put
45 procedure, public :: get => meshfieldcommcubedom3d_get
46 procedure, public :: exchange => meshfieldcommcubedom3d_exchange
47 procedure, public :: final => meshfieldcommcubedom3d_final
49
50 !-----------------------------------------------------------------------------
51 !
52 !++ Public parameters & variables
53 !
54
55 !-----------------------------------------------------------------------------
56 !
57 !++ Private procedure
58 !
59 private :: push_localsendbuf
60
61 !-----------------------------------------------------------------------------
62 !
63 !++ Private parameters & variables
64 !
66 integer, parameter :: comm_face_num = 6
67
68contains
69 subroutine meshfieldcommcubedom3d_init( this, &
70 sfield_num, hvfield_num, htensorfield_num, mesh3d )
71
72 implicit none
73
74 class(meshfieldcommcubedom3d), intent(inout) :: this
75 integer, intent(in) :: sfield_num
76 integer, intent(in) :: hvfield_num
77 integer, intent(in) :: htensorfield_num
78 class(meshcubedom3d), intent(in), target :: mesh3d
79
80 type(localmesh3d), pointer :: lcmesh
81 type(elementbase3d), pointer :: elem
82 integer :: n
83 integer :: nnode_lcmeshface(comm_face_num,mesh3d%local_mesh_num)
84 !-----------------------------------------------------------------------------
85
86 this%mesh3d => mesh3d
87 lcmesh => mesh3d%lcmesh_list(1)
88 elem => lcmesh%refElem3D
89 bufsize_per_field = 2*(lcmesh%NeX + lcmesh%NeY)*lcmesh%NeZ*elem%Nfp_h &
90 + 2*lcmesh%NeX*lcmesh%NeY*elem%Nfp_v
91
92 do n=1, this%mesh3d%LOCAL_MESH_NUM
93 lcmesh => this%mesh3d%lcmesh_list(n)
94 nnode_lcmeshface(:,n) = &
95 (/ lcmesh%NeX, lcmesh%NeY, lcmesh%NeX, lcmesh%NeY, 0, 0 /) * lcmesh%NeZ * lcmesh%refElem3D%Nfp_h &
96 + (/ 0, 0, 0, 0, 1, 1 /) * lcmesh%NeX*lcmesh%NeY * lcmesh%refElem3D%Nfp_v
97 end do
98
99 call meshfieldcommbase_init( this, sfield_num, hvfield_num, htensorfield_num, bufsize_per_field, comm_face_num, nnode_lcmeshface, mesh3d )
100
101 return
102 end subroutine meshfieldcommcubedom3d_init
103
104 subroutine meshfieldcommcubedom3d_final( this )
105 implicit none
106 class(meshfieldcommcubedom3d), intent(inout) :: this
107 !-----------------------------------------------------------------------------
108
109 call meshfieldcommbase_final( this )
110
111 return
112 end subroutine meshfieldcommcubedom3d_final
113
114 subroutine meshfieldcommcubedom3d_put(this, field_list, varid_s)
115 implicit none
116 class(meshfieldcommcubedom3d), intent(inout) :: this
117 type(meshfieldcontainer), intent(in) :: field_list(:)
118 integer, intent(in) :: varid_s
119
120 integer :: i
121 integer :: n
122 type(localmesh3d), pointer :: lcmesh
123 integer :: field_num
124 !-----------------------------------------------------------------------------
125
126! call PROF_rapstart( 'meshfiled_comm_put', 3)
127 field_num = size(field_list)
128 do n=1, this%mesh%LOCAL_MESH_NUM
129 lcmesh => this%mesh3d%lcmesh_list(n)
130 do i=1, field_num
131 call meshfieldcommbase_extract_bounddata( field_list(i)%field3d%local(n)%val, lcmesh%refElem, lcmesh, & ! (in)
132 this%send_buf(:,varid_s+i-1,n) ) ! (out)
133 end do
134 end do
135! call PROF_rapend( 'meshfiled_comm_put', 3)
136
137 return
138 end subroutine meshfieldcommcubedom3d_put
139
140 subroutine meshfieldcommcubedom3d_get(this, field_list, varid_s)
141 use scale_meshfieldcomm_base, only: &
143 implicit none
144
145 class(meshfieldcommcubedom3d), intent(inout) :: this
146 type(meshfieldcontainer), intent(inout) :: field_list(:)
147 integer, intent(in) :: varid_s
148
149 integer :: i
150 integer :: n
151 type(localmesh3d), pointer :: lcmesh
152 !-----------------------------------------------------------------------------
153
154 if ( this%call_wait_flag_sub_get ) then
155 call meshfieldcommbase_wait_core( this, this%commdata_list )
156 end if
157
158! call PROF_rapstart( 'meshfiled_comm_get', 3)
159 do i=1, size(field_list)
160 do n=1, this%mesh3d%LOCAL_MESH_NUM
161 lcmesh => this%mesh3d%lcmesh_list(n)
162 call meshfieldcommbase_set_bounddata( this%recv_buf(:,varid_s+i-1,n), lcmesh%refElem, lcmesh, & !(in)
163 field_list(i)%field3d%local(n)%val ) !(out)
164 end do
165 end do
166! call PROF_rapend( 'meshfiled_comm_get', 3)
167
168 return
169 end subroutine meshfieldcommcubedom3d_get
170
171!OCL SERIAL
172 subroutine meshfieldcommcubedom3d_exchange( this, do_wait )
173 use scale_meshfieldcomm_base, only: &
176 implicit none
177
178 class(meshfieldcommcubedom3d), intent(inout), target :: this
179 logical, intent(in), optional :: do_wait
180
181 integer :: n, f
182 type(localmesh3d), pointer :: lcmesh
183 type(localmeshcommdata), pointer :: commdata
184 !-----------------------------------------------------------------------------
185
186! call PROF_rapstart( 'meshfiled_comm_ex_push_buf', 3)
187 do n=1, this%mesh%LOCAL_MESH_NUM
188 lcmesh => this%mesh3d%lcmesh_list(n)
189 do f=1, this%nfaces_comm
190 commdata => this%commdata_list(f,n)
191 call push_localsendbuf( commdata%send_buf(:,:), & ! (inout)
192 this%send_buf(:,:,n), commdata%s_faceID, this%is_f(f,n), & ! (in)
193 commdata%Nnode_LCMeshFace, this%field_num_tot, lcmesh ) ! (in)
194 end do
195 end do
196! call PROF_rapend( 'meshfiled_comm_ex_push_buf', 3)
197
198 !-----------------------
199
200 call meshfieldcommbase_exchange_core(this, this%commdata_list(:,:), do_wait )
201
202 return
203 end subroutine meshfieldcommcubedom3d_exchange
204
205!----------------------------
206
207!OCL SERIAL
208 subroutine push_localsendbuf( lc_send_buf, send_buf, s_faceID, is, Nnode_LCMeshFace, var_num, lcmesh )
209 implicit none
210
211 integer, intent(in) :: nnode_lcmeshface
212 integer, intent(in) :: var_num
213 type(localmesh3d), intent(in) :: lcmesh
214 real(rp), intent(out) :: lc_send_buf(nnode_lcmeshface,var_num)
215 real(rp), intent(in) :: send_buf(bufsize_per_field,var_num)
216 integer, intent(in) :: s_faceid, is
217
218 integer :: ie
219 integer :: vid, i
220 class(elementbase3d), pointer :: e3d
221 !-----------------------------------------------------------------------------
222
223 if ( s_faceid > 0 ) then
224 !$omp parallel do
225 do vid=1, var_num
226 do i=1, nnode_lcmeshface
227 lc_send_buf(i,vid) = send_buf(is+i-1,vid)
228 end do
229 end do
230 else if ( -5 < s_faceid .and. s_faceid < 0) then
231 e3d => lcmesh%refElem3D
232 ie = is + nnode_lcmeshface - 1
233 call revert_hori( lc_send_buf(:,:), send_buf(is:ie,:), nnode_lcmeshface/lcmesh%NeZ, lcmesh%NeZ )
234 end if
235
236 return
237 contains
238 subroutine revert_hori(revert, ori, Ne_h1D, NeZ)
239 integer, intent(in) :: ne_h1d, nez
240 real(rp), intent(out) :: revert(e3d%nnode_h1d, e3d%nnode_v, ne_h1d, nez, var_num)
241 real(rp), intent(in) :: ori(e3d%nnode_h1d, e3d%nnode_v, ne_h1d, nez, var_num)
242
243 integer :: p1, p3, i, k, n
244 integer :: i_, p1_
245 !-----------------------------------------------------------------------------
246
247 do n=1, var_num
248 do k=1, nez
249 do i=1, ne_h1d
250 i_ = ne_h1d - i + 1
251 do p3=1, e3d%Nnode_v
252 do p1=1, e3d%Nnode_h1D
253 p1_ = e3d%Nnode_h1D - p1 + 1
254 revert(p1,p3,i,k,n) = ori(p1_,p3,i_,k,n)
255 end do
256 end do
257 end do
258 end do
259 end do
260 end subroutine revert_hori
261 end subroutine push_localsendbuf
262end module scale_meshfieldcomm_cubedom3d
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)
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.
module FElib / Data / Communication 3D cubic domain
Derived type to manage data communication between a face of local mesh.
Base derived type to manage data communication.