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