FE-Project
Loading...
Searching...
No Matches
scale_meshfieldcomm_cubedspheredom2d.F90
Go to the documentation of this file.
1!-------------------------------------------------------------------------------
2!> module FElib / Data / Communication in 2D cubed-sphere domain
3!!
4!! @par Description
5!! A module to manage data communication with 2D cubed-sphere 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: &
29
30 !-----------------------------------------------------------------------------
31 implicit none
32 private
33
34 !-----------------------------------------------------------------------------
35 !
36 !++ Public type & procedure
37 !
38
39 !> Derived type to represent covariant vector components
40 type :: veccovariantcomp
41 type(MeshField2D), pointer :: u1 => null()
42 type(MeshField2D), pointer :: u2 => null()
43 end type
44
45 !> Base derived type to manage data communication with 2D cubed-sphere domain
47 class(meshcubedspheredom2d), pointer :: mesh2d !< Pointer to an object representing 2D cubed-sphere computational mesh
48 type(veccovariantcomp), allocatable :: vec_covariant_comp_ptrlist(:)
49 integer, allocatable :: nnode_lcmeshallface(:)
50 contains
51 procedure, public :: init => meshfieldcommcubedspheredom2d_init
52 procedure, public :: put => meshfieldcommcubedspheredom2d_put
53 procedure, public :: get => meshfieldcommcubedspheredom2d_get
54 procedure, public :: exchange => meshfieldcommcubedspheredom2d_exchange
55 procedure, public :: setcovariantvec => meshfieldcommcubedspheredom2d_set_covariantvec
56 procedure, public :: final => meshfieldcommcubedspheredom2d_final
58
59 !-----------------------------------------------------------------------------
60 !
61 !++ Public parameters & variables
62 !
63
64 !-----------------------------------------------------------------------------
65 !
66 !++ Private type & procedure
67 !
68 private :: post_exchange_core
69 private :: push_localsendbuf
70
71 !-----------------------------------------------------------------------------
72 !
73 !++ Private parameters & variables
74 !
76 integer, parameter :: comm_face_num = 4
77
78contains
79 subroutine meshfieldcommcubedspheredom2d_init( this, &
80 sfield_num, hvfield_num, htensorfield_num, mesh2d )
81
82 implicit none
83
84 class(meshfieldcommcubedspheredom2d), intent(inout) :: this
85 integer, intent(in) :: sfield_num
86 integer, intent(in) :: hvfield_num
87 integer, intent(in) :: htensorfield_num
88 class(meshcubedspheredom2d), intent(in), target :: mesh2d
89
90 type(localmesh2d), pointer :: lcmesh
91 integer :: n
92 integer :: nnode_lcmeshface(comm_face_num,mesh2d%local_mesh_num)
93 !-----------------------------------------------------------------------------
94
95 this%mesh2d => mesh2d
96 lcmesh => mesh2d%lcmesh_list(1)
97 bufsize_per_field = 2*(lcmesh%NeX + lcmesh%NeY)*lcmesh%refElem2D%Nfp
98
99 allocate( this%Nnode_LCMeshAllFace(mesh2d%LOCAL_MESH_NUM) )
100 do n=1, this%mesh2d%LOCAL_MESH_NUM
101 lcmesh => this%mesh2d%lcmesh_list(n)
102 nnode_lcmeshface(:,n) = (/ lcmesh%NeX, lcmesh%NeY, lcmesh%NeX, lcmesh%NeY /) * lcmesh%refElem2D%Nfp
103 this%Nnode_LCMeshAllFace(n) = sum(nnode_lcmeshface(:,n))
104 end do
105
106 call meshfieldcommbase_init( this, sfield_num, hvfield_num, htensorfield_num, bufsize_per_field, comm_face_num, nnode_lcmeshface, mesh2d)
107
108 if (hvfield_num > 0) then
109 allocate( this%vec_covariant_comp_ptrlist(hvfield_num) )
110 end if
111
112 return
113 end subroutine meshfieldcommcubedspheredom2d_init
114
115 subroutine meshfieldcommcubedspheredom2d_final( this )
116 implicit none
117
118 class(meshfieldcommcubedspheredom2d), intent(inout) :: this
119 !-----------------------------------------------------------------------------
120
121 deallocate( this%Nnode_LCMeshAllFace )
122
123 if ( this%hvfield_num > 0 ) then
124 deallocate( this%vec_covariant_comp_ptrlist )
125 end if
126
127 call meshfieldcommbase_final( this )
128
129 return
130 end subroutine meshfieldcommcubedspheredom2d_final
131
132 subroutine meshfieldcommcubedspheredom2d_set_covariantvec( &
133 this, hvfield_ID, u1, u2 )
134 implicit none
135 class(meshfieldcommcubedspheredom2d), intent(inout) :: this
136 integer, intent(in) :: hvfield_id
137 type(meshfield2d), intent(in), target :: u1
138 type(meshfield2d), intent(in), target :: u2
139 !--------------------------------------------------------------
140
141 this%vec_covariant_comp_ptrlist(hvfield_id)%u1 => u1
142 this%vec_covariant_comp_ptrlist(hvfield_id)%u2 => u2
143
144 return
145 end subroutine meshfieldcommcubedspheredom2d_set_covariantvec
146
147 subroutine meshfieldcommcubedspheredom2d_put(this, field_list, varid_s)
148 implicit none
149 class(meshfieldcommcubedspheredom2d), intent(inout) :: this
150 type(meshfieldcontainer), intent(in) :: field_list(:)
151 integer, intent(in) :: varid_s
152
153 integer :: i
154 integer :: n
155 type(localmesh2d), pointer :: lcmesh
156 !-----------------------------------------------------------------------------
157
158 do i=1, size(field_list)
159 do n=1, this%mesh%LOCAL_MESH_NUM
160 lcmesh => this%mesh2d%lcmesh_list(n)
161 call meshfieldcommbase_extract_bounddata( field_list(i)%field2d%local(n)%val, lcmesh%refElem, lcmesh, & ! (in)
162 this%send_buf(:,varid_s+i-1,n) ) ! (out)
163 end do
164 end do
165
166 return
167 end subroutine meshfieldcommcubedspheredom2d_put
168
169 subroutine meshfieldcommcubedspheredom2d_get(this, field_list, varid_s)
170 use scale_meshfieldcomm_base, only: &
172 implicit none
173
174 class(meshfieldcommcubedspheredom2d), intent(inout) :: this
175 type(meshfieldcontainer), intent(inout) :: field_list(:)
176 integer, intent(in) :: varid_s
177
178 integer :: i
179 integer :: n
180 type(localmesh2d), pointer :: lcmesh
181
182 integer :: varnum
183 integer :: varid_e
184 integer :: varid_vec_s
185 type(meshfield2d), pointer :: u1, u2
186 !-----------------------------------------------------------------------------
187
188 varnum = size(field_list)
189
190 !--
191 if ( this%call_wait_flag_sub_get ) then
192 call meshfieldcommbase_wait_core( this, this%commdata_list )
193 call post_exchange_core( this )
194 end if
195
196 !--
197 do i=1, varnum
198 do n=1, this%mesh2d%LOCAL_MESH_NUM
199 lcmesh => this%mesh2d%lcmesh_list(n)
200 call meshfieldcommbase_set_bounddata( this%recv_buf(:,varid_s+i-1,n), lcmesh%refElem, lcmesh, & !(in)
201 field_list(i)%field2d%local(n)%val ) !(out)
202 end do
203 end do
204
205 varid_e = varid_s + varnum - 1
206 if ( varid_e > this%sfield_num ) then
207 do i=1, this%hvfield_num
208
209 varid_vec_s = this%sfield_num + 2*i - 1
210 if ( varid_vec_s > varid_e ) exit
211
212 if ( associated(this%vec_covariant_comp_ptrlist(i)%u1 ) &
213 .and. associated(this%vec_covariant_comp_ptrlist(i)%u2 ) ) then
214
215 do n=1, this%mesh2d%LOCAL_MESH_NUM
216 call set_boundary_data2d_u1u2( &
217 this%recv_buf(:,varid_vec_s,n), this%recv_buf(:,varid_vec_s+1,n), & ! (in)
218 lcmesh%refElem2D, lcmesh, lcmesh%G_ij, & ! (in)
219 this%vec_covariant_comp_ptrlist(i)%u1%local(n)%val, & ! (out)
220 this%vec_covariant_comp_ptrlist(i)%u2%local(n)%val ) ! (out)
221 end do
222 end if
223 end do
224 end if
225
226 return
227 end subroutine meshfieldcommcubedspheredom2d_get
228
229!OCL SERIAL
230 subroutine meshfieldcommcubedspheredom2d_exchange( this, do_wait )
231 use scale_meshfieldcomm_base, only: &
234
235 use scale_cubedsphere_coord_cnv, only: &
238
239 implicit none
240
241 class(meshfieldcommcubedspheredom2d), intent(inout), target :: this
242 logical, intent(in), optional :: do_wait
243
244 integer :: n, f
245 integer :: varid
246
247 real(rp), allocatable :: fpos2d(:,:)
248 real(rp), allocatable :: lcfpos2d(:,:)
249 real(rp), allocatable :: unity_fac(:)
250 real(rp), allocatable :: tmp_svec2d(:,:)
251
252 class(elementbase2d), pointer :: elem
253 type(localmesh2d), pointer :: lcmesh
254 type(localmeshcommdata), pointer :: commdata
255
256 integer :: irs, ire
257 !-----------------------------------------------------------------------------
258
259 do n=1, this%mesh%LOCAL_MESH_NUM
260 lcmesh => this%mesh2d%lcmesh_list(n)
261 elem => lcmesh%refElem2D
262
263 allocate( fpos2d(this%Nnode_LCMeshAllFace(n),2) )
264 call extract_boundary_data2d( lcmesh%pos_en(:,:,1), elem, lcmesh, fpos2d(:,1) )
265 call extract_boundary_data2d( lcmesh%pos_en(:,:,2), elem, lcmesh, fpos2d(:,2) )
266
267 irs = 1
268 do f=1, this%nfaces_comm
269 commdata => this%commdata_list(f,n)
270 call push_localsendbuf( commdata%send_buf(:,:), & ! (inout)
271 this%send_buf(:,:,n), commdata%s_faceID, this%is_f(f,n), & ! (in)
272 commdata%Nnode_LCMeshFace, this%field_num_tot) ! (in)
273
274 if ( commdata%s_panelID /= lcmesh%panelID ) then
275 if ( this%hvfield_num > 0 ) then
276
277 allocate( lcfpos2d(commdata%Nnode_LCMeshFace,2), unity_fac(commdata%Nnode_LCMeshFace) )
278 allocate( tmp_svec2d(commdata%Nnode_LCMeshFace,2) )
279
280 call push_localsendbuf( lcfpos2d, &
281 fpos2d, commdata%s_faceID, this%is_f(f,n), commdata%Nnode_LCMeshFace, 2 )
282 unity_fac(:) = 1.0_rp
283
284 ire = irs + commdata%Nnode_LCMeshFace - 1
285
286 do varid=this%sfield_num+1, this%field_num_tot-1,2
287 tmp_svec2d(:,1) = commdata%send_buf(:,varid )
288 tmp_svec2d(:,2) = commdata%send_buf(:,varid+1)
290 lcmesh%panelID, lcfpos2d(:,1), lcfpos2d(:,2), unity_fac, &
291 commdata%Nnode_LCMeshFace, &
292 tmp_svec2d(:,1), tmp_svec2d(:,2), &
293 commdata%send_buf(:,varid), commdata%send_buf(:,varid+1) )
294 end do
295
296 deallocate( lcfpos2d, unity_fac, tmp_svec2d )
297 end if
298 end if
299
300 irs = ire + 1
301 end do
302 deallocate( fpos2d )
303 end do
304
305 !-----------------------
306
307 call meshfieldcommbase_exchange_core(this, this%commdata_list(:,:), do_wait )
308 if ( .not. this%call_wait_flag_sub_get ) &
309 call post_exchange_core( this )
310
311 return
312 end subroutine meshfieldcommcubedspheredom2d_exchange
313
314!- Private ---------------------------
315
316 subroutine post_exchange_core( this )
318 use scale_cubedsphere_coord_cnv, only: &
320 implicit none
321
322 class(meshfieldcommcubedspheredom2d), intent(inout), target :: this
323
324 integer :: n, f
325 integer :: varid
326
327 real(rp), allocatable :: fpos2d(:,:)
328 real(rp), allocatable :: lcfpos2d(:,:)
329 real(rp), allocatable :: unity_fac(:)
330
331 class(elementbase2d), pointer :: elem
332 type(localmesh2d), pointer :: lcmesh
333 type(localmeshcommdata), pointer :: commdata
334
335 integer :: irs, ire
336 !-----------------------------------------------------------------------------
337
338 do n=1, this%mesh%LOCAL_MESH_NUM
339 lcmesh => this%mesh2d%lcmesh_list(n)
340 elem => lcmesh%refElem2D
341
342 allocate( fpos2d(this%Nnode_LCMeshAllFace(n),2) )
343 call extract_boundary_data2d( lcmesh%pos_en(:,:,1), elem, lcmesh, fpos2d(:,1) )
344 call extract_boundary_data2d( lcmesh%pos_en(:,:,2), elem, lcmesh, fpos2d(:,2) )
345
346 irs = 1
347 do f=1, this%nfaces_comm
348 commdata => this%commdata_list(f,n)
349 ire = irs + commdata%Nnode_LCMeshFace - 1
350
351 if ( commdata%s_panelID /= lcmesh%panelID ) then
352 if ( this%hvfield_num > 0 ) then
353
354 allocate( lcfpos2d(commdata%Nnode_LCMeshFace,2), unity_fac(commdata%Nnode_LCMeshFace) )
355
356 call push_localsendbuf( lcfpos2d, &
357 fpos2d, f, this%is_f(f,n), commdata%Nnode_LCMeshFace, 2 )
358 unity_fac(:) = 1.0_rp
359
360 do varid=this%sfield_num+1, this%field_num_tot-1, 2
362 lcmesh%panelID, lcfpos2d(:,1), lcfpos2d(:,2), unity_fac(:), &
363 commdata%Nnode_LCMeshFace, &
364 commdata%recv_buf(:,varid), commdata%recv_buf(:,varid+1), &
365 this%recv_buf(irs:ire,varid,n), this%recv_buf(irs:ire,varid+1,n) )
366 end do
367 deallocate( lcfpos2d, unity_fac )
368 end if
369 end if
370
371 irs = ire + 1
372 end do
373
374 deallocate( fpos2d )
375 end do
376
377 return
378 end subroutine post_exchange_core
379
380 subroutine push_localsendbuf( lc_send_buf, send_buf, s_faceID, is, Nnode_LCMeshFace, var_num )
381 implicit none
382
383 integer, intent(in) :: var_num
384 integer, intent(in) :: nnode_lcmeshface
385 real(rp), intent(inout) :: lc_send_buf(nnode_lcmeshface,var_num)
386 real(rp), intent(in) :: send_buf(bufsize_per_field,var_num)
387 integer, intent(in) :: s_faceid, is
388
389 integer :: ie
390 !-----------------------------------------------------------------------------
391
392 ie = is + nnode_lcmeshface - 1
393 if ( s_faceid > 0 ) then
394 lc_send_buf(:,:) = send_buf(is:ie,:)
395 else
396 lc_send_buf(:,:) = send_buf(ie:is:-1,:)
397 end if
398
399 return
400 end subroutine push_localsendbuf
401
402 subroutine extract_boundary_data2d( var, refElem, mesh, buf )
403 implicit none
404
405 type(elementbase2d), intent(in) :: refelem
406 type(localmesh2d), intent(in) :: mesh
407 real(dp), intent(in) :: var(refelem%np * mesh%ne)
408 real(dp), intent(inout) :: buf(refelem%nfp * mesh%nex * 4)
409 !------------------------------------------------------------
410
411 buf(:) = var(mesh%VmapB(:))
412
413 return
414 end subroutine extract_boundary_data2d
415
416 subroutine set_boundary_data2d_u1u2( buf_U, buf_V, &
417 elem, mesh, G_ij, &
418 u1, u2)
419
420 implicit none
421
422 type(elementbase2d), intent(in) :: elem
423 type(localmesh2d), intent(in) :: mesh
424 real(dp), intent(in) :: buf_u(elem%nfp * mesh%nex * 4)
425 real(dp), intent(in) :: buf_v(elem%nfp * mesh%nex * 4)
426 real(dp), intent(in) :: g_ij(elem%np * mesh%ne,2,2)
427 real(dp), intent(inout) :: u1(elem%np * mesh%nea)
428 real(dp), intent(inout) :: u2(elem%np * mesh%nea)
429 !------------------------------------------------------------
430
431 u1(elem%Np*mesh%NeE+1:elem%Np*mesh%NeE+size(buf_u)) &
432 = g_ij(mesh%VmapB,1,1) * buf_u(:) + g_ij(mesh%VmapB,1,2) * buf_v(:)
433 u2(elem%Np*mesh%NeE+1:elem%Np*mesh%NeE+size(buf_u)) &
434 = g_ij(mesh%VmapB,2,1) * buf_u(:) + g_ij(mesh%VmapB,2,2) * buf_v(:)
435
436 return
437 end subroutine set_boundary_data2d_u1u2
438
Module common / Coordinate conversion with cubed-sphere projection.
subroutine, public cubedspherecoordcnv_cs2lonlatvec(panelid, alpha, beta, gam, np, vecalpha, vecbeta, veclon, veclat, lat)
Covert the components of a vector in local coordinates with an equiangular gnomonic cubed-sphere proj...
subroutine, public cubedspherecoordcnv_lonlat2csvec(panelid, alpha, beta, gam, np, veclon, veclat, vecalpha, vecbeta, lat)
Covert the components of a vector in longitude and latitude coordinates to those in local coordinates...
module FElib / Element / Base
module FElib / Mesh / Local 2D
module FElib / Mesh / Cubed-sphere 2D 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_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 in 2D cubed-sphere domain
Derived type representing a 2D reference element.
Derived type representing an arbitrary finite element.
Derived type representing a field with 2D 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 2D cubed-sphere domain.