FE-Project
Loading...
Searching...
No Matches
scale_meshfieldcomm_cubedspheredom3d.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: &
28
31
32
33 !-----------------------------------------------------------------------------
34 implicit none
35 private
36
37 !-----------------------------------------------------------------------------
38 !
39 !++ Public type & procedure
40 !
41
42 type :: veccovariantcomp
43 type(MeshField3D), pointer :: u1 => null()
44 type(MeshField3D), pointer :: u2 => null()
45 end type
46
48 class(meshcubedspheredom3d), pointer :: mesh3d
49 type(veccovariantcomp), allocatable :: vec_covariant_comp_ptrlist(:)
50 integer, allocatable :: nnode_lcmeshallface(:)
51 contains
52 procedure, public :: init => meshfieldcommcubedspheredom3d_init
53 procedure, public :: put => meshfieldcommcubedspheredom3d_put
54 procedure, public :: get => meshfieldcommcubedspheredom3d_get
55 procedure, public :: exchange => meshfieldcommcubedspheredom3d_exchange
56 procedure, public :: setcovariantvec => meshfieldcommcubedspheredom3d_set_covariantvec
57 procedure, public :: final => meshfieldcommcubedspheredom3d_final
59
60 !-----------------------------------------------------------------------------
61 !
62 !++ Public parameters & variables
63 !
64
65 !-----------------------------------------------------------------------------
66 !
67 !++ Private type & procedure
68 !
69 private :: post_exchange_core
70 private :: push_localsendbuf
71
72 !-----------------------------------------------------------------------------
73 !
74 !++ Private parameters & variables
75 !
77 integer, parameter :: comm_face_num = 6
78
79contains
80 subroutine meshfieldcommcubedspheredom3d_init( this, &
81 sfield_num, hvfield_num, htensorfield_num, mesh3d )
82
83 implicit none
84
85 class(meshfieldcommcubedspheredom3d), intent(inout), target :: this
86 integer, intent(in) :: sfield_num
87 integer, intent(in) :: hvfield_num
88 integer, intent(in) :: htensorfield_num
89 class(meshcubedspheredom3d), intent(in), target :: mesh3d
90
91 type(localmesh3d), pointer :: lcmesh
92 type(elementbase3d), pointer :: elem
93 integer :: n
94 integer :: nnode_lcmeshface(comm_face_num,mesh3d%local_mesh_num)
95 !-----------------------------------------------------------------------------
96
97 this%mesh3d => mesh3d
98 lcmesh => mesh3d%lcmesh_list(1)
99 elem => lcmesh%refElem3D
100
101 bufsize_per_field = 2*(lcmesh%NeX + lcmesh%NeY)*lcmesh%NeZ*elem%Nfp_h &
102 + 2*lcmesh%NeX*lcmesh%NeY*elem%Nfp_v
103
104 allocate( this%Nnode_LCMeshAllFace(mesh3d%LOCAL_MESH_NUM) )
105 do n=1, this%mesh3d%LOCAL_MESH_NUM
106 lcmesh => this%mesh3d%lcmesh_list(n)
107 nnode_lcmeshface(:,n) = &
108 (/ lcmesh%NeX, lcmesh%NeY, lcmesh%NeX, lcmesh%NeY, 0, 0 /) * lcmesh%NeZ * lcmesh%refElem3D%Nfp_h &
109 + (/ 0, 0, 0, 0, 1, 1 /) * lcmesh%NeX*lcmesh%NeY * lcmesh%refElem3D%Nfp_v
110 this%Nnode_LCMeshAllFace(n) = sum(nnode_lcmeshface(:,n))
111 end do
112
113 call meshfieldcommbase_init( this, sfield_num, hvfield_num, htensorfield_num, bufsize_per_field, comm_face_num, nnode_lcmeshface, mesh3d )
114
115 if (hvfield_num > 0) then
116 allocate( this%vec_covariant_comp_ptrlist(hvfield_num) )
117 end if
118
119 return
120 end subroutine meshfieldcommcubedspheredom3d_init
121
122 subroutine meshfieldcommcubedspheredom3d_final( this )
123
124 implicit none
125
126 class(meshfieldcommcubedspheredom3d), intent(inout) :: this
127 !-----------------------------------------------------------------------------
128
129 if ( this%hvfield_num > 0 ) then
130 deallocate( this%vec_covariant_comp_ptrlist )
131 end if
132
133 call meshfieldcommbase_final( this )
134
135 return
136 end subroutine meshfieldcommcubedspheredom3d_final
137
138 subroutine meshfieldcommcubedspheredom3d_set_covariantvec( &
139 this, hvfield_ID, u1, u2 )
140 implicit none
141 class(meshfieldcommcubedspheredom3d), intent(inout) :: this
142 integer, intent(in) :: hvfield_id
143 type(meshfield3d), intent(in), target :: u1
144 type(meshfield3d), intent(in), target :: u2
145 !--------------------------------------------------------------
146
147 this%vec_covariant_comp_ptrlist(hvfield_id)%u1 => u1
148 this%vec_covariant_comp_ptrlist(hvfield_id)%u2 => u2
149
150 return
151 end subroutine meshfieldcommcubedspheredom3d_set_covariantvec
152
153 subroutine meshfieldcommcubedspheredom3d_put(this, field_list, varid_s)
154 implicit none
155 class(meshfieldcommcubedspheredom3d), intent(inout) :: this
156 type(meshfieldcontainer), intent(in) :: field_list(:)
157 integer, intent(in) :: varid_s
158
159 integer :: i
160 integer :: n
161 type(localmesh3d), pointer :: lcmesh
162 !-----------------------------------------------------------------------------
163
164 do i=1, size(field_list)
165 do n=1, this%mesh%LOCAL_MESH_NUM
166 lcmesh => this%mesh3d%lcmesh_list(n)
167 call meshfieldcommbase_extract_bounddata( field_list(i)%field3d%local(n)%val, lcmesh%refElem, lcmesh, & ! (in)
168 this%send_buf(:,varid_s+i-1,n) ) ! (out)
169 end do
170 end do
171
172 return
173 end subroutine meshfieldcommcubedspheredom3d_put
174
175 subroutine meshfieldcommcubedspheredom3d_get(this, field_list, varid_s)
176 use scale_meshfieldcomm_base, only: &
178 implicit none
179
180 class(meshfieldcommcubedspheredom3d), intent(inout) :: this
181 type(meshfieldcontainer), intent(inout) :: field_list(:)
182 integer, intent(in) :: varid_s
183
184 integer :: i
185 integer :: n
186 integer :: ke
187 integer :: ke2d
188 type(localmesh3d), pointer :: lcmesh
189 type(elementbase3d), pointer :: elem
190
191 integer :: varnum
192 integer :: varid_e
193 integer :: varid_vec_s
194
195 real(rp), allocatable :: g_ij(:,:,:,:)
196 !-----------------------------------------------------------------------------
197
198 varnum = size(field_list)
199
200 !--
201 if ( this%call_wait_flag_sub_get ) then
202 call post_exchange_core( this )
203 call meshfieldcommbase_wait_core( this, this%commdata_list )
204 end if
205
206 !--
207 do i=1, varnum
208 do n=1, this%mesh3d%LOCAL_MESH_NUM
209 lcmesh => this%mesh3d%lcmesh_list(n)
210 call meshfieldcommbase_set_bounddata( this%recv_buf(:,varid_s+i-1,n), lcmesh%refElem, lcmesh, & !(in)
211 field_list(i)%field3d%local(n)%val ) !(out)
212 end do
213 end do
214
215 varid_e = varid_s + varnum - 1
216 if ( varid_e > this%sfield_num ) then
217 do i=1, this%hvfield_num
218
219 varid_vec_s = this%sfield_num + 2*i - 1
220 if ( varid_vec_s > varid_e ) exit
221
222 if ( associated(this%vec_covariant_comp_ptrlist(i)%u1 ) &
223 .and. associated(this%vec_covariant_comp_ptrlist(i)%u2 ) ) then
224
225 do n=1, this%mesh3d%LOCAL_MESH_NUM
226 lcmesh => this%mesh3d%lcmesh_list(n)
227 elem => lcmesh%refElem3D
228
229 allocate( g_ij(elem%Np,lcmesh%Ne,2,2) )
230 !$omp parallel do private(ke2D)
231 do ke=lcmesh%NeS, lcmesh%NeE
232 ke2d = lcmesh%EMap3Dto2D(ke)
233 g_ij(:,ke,1,1) = lcmesh%G_ij(elem%IndexH2Dto3D(:),ke2d,1,1)
234 g_ij(:,ke,2,1) = lcmesh%G_ij(elem%IndexH2Dto3D(:),ke2d,2,1)
235 g_ij(:,ke,1,2) = lcmesh%G_ij(elem%IndexH2Dto3D(:),ke2d,1,2)
236 g_ij(:,ke,2,2) = lcmesh%G_ij(elem%IndexH2Dto3D(:),ke2d,2,2)
237 end do
238
239 call set_boundary_data3d_u1u2( &
240 this%recv_buf(:,varid_vec_s,n), this%recv_buf(:,varid_vec_s+1,n), & ! (in)
241 lcmesh%refElem3D, lcmesh, g_ij(:,:,:,:), & ! (in)
242 this%vec_covariant_comp_ptrlist(i)%u1%local(n)%val, & ! (out)
243 this%vec_covariant_comp_ptrlist(i)%u2%local(n)%val ) ! (out)
244 end do
245 end if
246 end do
247 end if
248
249 return
250 end subroutine meshfieldcommcubedspheredom3d_get
251
252!OCL SERIAL
253 subroutine meshfieldcommcubedspheredom3d_exchange( this, do_wait )
254 use scale_meshfieldcomm_base, only: &
257
258 use scale_cubedsphere_coord_cnv, only: &
261
262 implicit none
263
264 class(meshfieldcommcubedspheredom3d), intent(inout), target :: this
265 logical, intent(in), optional :: do_wait
266
267 integer :: n, f
268 integer :: varid
269
270 real(rp), allocatable :: fpos3d(:,:)
271 real(rp), allocatable :: lcfpos3d(:,:)
272 real(rp), allocatable :: unity_fac(:)
273 real(rp), allocatable :: tmp_svec3d(:,:)
274 real(rp), allocatable :: tmp1_htensor3d(:,:,:)
275 real(rp), allocatable :: tmp2_htensor3d(:,:,:)
276
277 class(elementbase3d), pointer :: elem
278 type(localmesh3d), pointer :: lcmesh
279 type(localmeshcommdata), pointer :: commdata
280
281 integer :: irs, ire
282 !-----------------------------------------------------------------------------
283
284 do n=1, this%mesh%LOCAL_MESH_NUM
285 lcmesh => this%mesh3d%lcmesh_list(n)
286 elem => lcmesh%refElem3D
287
288 allocate( fpos3d(this%Nnode_LCMeshAllFace(n),2) )
289 call extract_boundary_data3d( lcmesh%pos_en(:,:,1), elem, lcmesh, fpos3d(:,1) )
290 call extract_boundary_data3d( lcmesh%pos_en(:,:,2), elem, lcmesh, fpos3d(:,2) )
291
292 do f=1, this%nfaces_comm
293 commdata => this%commdata_list(f,n)
294 call push_localsendbuf( commdata%send_buf(:,:), & ! (inout)
295 this%send_buf(:,:,n), commdata%s_faceID, this%is_f(f,n), & ! (in)
296 commdata%Nnode_LCMeshFace, this%field_num_tot, & ! (in)
297 lcmesh, elem ) ! (in)
298
299 if ( commdata%s_panelID /= lcmesh%panelID ) then
300 if ( this%hvfield_num > 0 ) then
301
302 allocate( lcfpos3d(commdata%Nnode_LCMeshFace,2), unity_fac(commdata%Nnode_LCMeshFace) )
303 allocate( tmp_svec3d(commdata%Nnode_LCMeshFace,2) )
304
305 call push_localsendbuf( lcfpos3d, &
306 fpos3d, commdata%s_faceID, this%is_f(f,n), commdata%Nnode_LCMeshFace, 2, &
307 lcmesh, elem )
308 unity_fac(:) = 1.0_rp
309
310 do varid=this%sfield_num+1, this%sfield_num+2*this%hvfield_num-1, 2
311 tmp_svec3d(:,1) = commdata%send_buf(:,varid )
312 tmp_svec3d(:,2) = commdata%send_buf(:,varid+1)
313
315 lcmesh%panelID, lcfpos3d(:,1), lcfpos3d(:,2), unity_fac(:), &
316 commdata%Nnode_LCMeshFace, &
317 tmp_svec3d(:,1), tmp_svec3d(:,2), &
318 commdata%send_buf(:,varid), commdata%send_buf(:,varid+1) )
319
320 end do
321 deallocate( lcfpos3d, unity_fac, tmp_svec3d )
322 end if
323
324 if ( this%htensorfield_num > 0 ) then
325 allocate( lcfpos3d(commdata%Nnode_LCMeshFace,2), unity_fac(commdata%Nnode_LCMeshFace) )
326 allocate( tmp1_htensor3d(commdata%Nnode_LCMeshFace,2,2) )
327 allocate( tmp2_htensor3d(commdata%Nnode_LCMeshFace,2,2) )
328
329 call push_localsendbuf( lcfpos3d, &
330 fpos3d, commdata%s_faceID, this%is_f(f,n), commdata%Nnode_LCMeshFace, 2, &
331 lcmesh, elem )
332 unity_fac(:) = 1.0_rp
333
334 do varid=this%sfield_num+2*this%hvfield_num+1, this%field_num_tot-3, 4
335 tmp1_htensor3d(:,1,1) = commdata%send_buf(:,varid )
336 tmp1_htensor3d(:,2,1) = commdata%send_buf(:,varid+1)
337 tmp1_htensor3d(:,1,2) = commdata%send_buf(:,varid+2)
338 tmp1_htensor3d(:,2,2) = commdata%send_buf(:,varid+3)
339
341 lcmesh%panelID, lcfpos3d(:,1), lcfpos3d(:,2), unity_fac(:), &
342 commdata%Nnode_LCMeshFace, &
343 tmp1_htensor3d(:,1,1), tmp1_htensor3d(:,2,1), &
344 tmp2_htensor3d(:,1,1), tmp2_htensor3d(:,2,1) )
346 lcmesh%panelID, lcfpos3d(:,1), lcfpos3d(:,2), unity_fac(:), &
347 commdata%Nnode_LCMeshFace, &
348 tmp1_htensor3d(:,1,2), tmp1_htensor3d(:,2,2), &
349 tmp2_htensor3d(:,1,2), tmp2_htensor3d(:,2,2) )
351 lcmesh%panelID, lcfpos3d(:,1), lcfpos3d(:,2), unity_fac(:), &
352 commdata%Nnode_LCMeshFace, &
353 tmp2_htensor3d(:,1,1), tmp2_htensor3d(:,1,2), &
354 commdata%send_buf(:,varid), commdata%send_buf(:,varid+2) )
356 lcmesh%panelID, lcfpos3d(:,1), lcfpos3d(:,2), unity_fac(:), &
357 commdata%Nnode_LCMeshFace, &
358 tmp2_htensor3d(:,2,1), tmp2_htensor3d(:,2,2), &
359 commdata%send_buf(:,varid+1), commdata%send_buf(:,varid+3) )
360 end do
361 deallocate( lcfpos3d, unity_fac, tmp1_htensor3d, tmp2_htensor3d )
362 end if
363
364 end if
365
366 end do
367 deallocate( fpos3d )
368 end do
369
370 !-----------------------
371
372 call meshfieldcommbase_exchange_core(this, this%commdata_list(:,:), do_wait )
373 if ( .not. this%call_wait_flag_sub_get ) &
374 call post_exchange_core( this )
375
376 return
377 end subroutine meshfieldcommcubedspheredom3d_exchange
378
379!----------------------------
380
381 subroutine post_exchange_core( this )
383 use scale_cubedsphere_coord_cnv, only: &
385 implicit none
386
387 class(meshfieldcommcubedspheredom3d), intent(inout), target :: this
388
389 integer :: n, f
390 integer :: varid
391
392 real(rp), allocatable :: fpos3d(:,:)
393 real(rp), allocatable :: lcfpos3d(:,:)
394 real(rp), allocatable :: unity_fac(:)
395 real(rp), allocatable :: tmp1_htensor3d(:,:,:)
396
397 class(elementbase3d), pointer :: elem
398 type(localmesh3d), pointer :: lcmesh
399 type(localmeshcommdata), pointer :: commdata
400
401 integer :: irs, ire
402 !-----------------------------------------------------------------------------
403
404 do n=1, this%mesh%LOCAL_MESH_NUM
405 lcmesh => this%mesh3d%lcmesh_list(n)
406 elem => lcmesh%refElem3D
407
408 allocate( fpos3d(this%Nnode_LCMeshAllFace(n),2) )
409 call extract_boundary_data3d( lcmesh%pos_en(:,:,1), elem, lcmesh, fpos3d(:,1) )
410 call extract_boundary_data3d( lcmesh%pos_en(:,:,2), elem, lcmesh, fpos3d(:,2) )
411
412 irs = 1
413 do f=1, this%nfaces_comm
414 commdata => this%commdata_list(f,n)
415 ire = irs + commdata%Nnode_LCMeshFace - 1
416
417 if ( commdata%s_panelID /= lcmesh%panelID ) then
418 if ( this%hvfield_num > 0 ) then
419
420 allocate( lcfpos3d(commdata%Nnode_LCMeshFace,2), unity_fac(commdata%Nnode_LCMeshFace) )
421 call push_localsendbuf( lcfpos3d, &
422 fpos3d, f, this%is_f(f,n), commdata%Nnode_LCMeshFace, 2, &
423 lcmesh, elem )
424 unity_fac(:) = 1.0_rp
425
426 do varid=this%sfield_num+1, this%sfield_num+2*this%hvfield_num-1, 2
428 lcmesh%panelID, lcfpos3d(:,1), lcfpos3d(:,2), unity_fac(:), &
429 commdata%Nnode_LCMeshFace, &
430 commdata%recv_buf(:,varid), commdata%recv_buf(:,varid+1), &
431 this%recv_buf(irs:ire,varid,n), this%recv_buf(irs:ire,varid+1,n) )
432 end do
433 deallocate( lcfpos3d, unity_fac )
434 end if
435
436 if ( this%htensorfield_num > 0 ) then
437 allocate( lcfpos3d(commdata%Nnode_LCMeshFace,2), unity_fac(commdata%Nnode_LCMeshFace) )
438 allocate( tmp1_htensor3d(commdata%Nnode_LCMeshFace,2,2) )
439 unity_fac(:) = 1.0_rp
440
441 call push_localsendbuf( lcfpos3d, &
442 fpos3d, f, this%is_f(f,n), commdata%Nnode_LCMeshFace, 2, &
443 lcmesh, elem )
444
445 do varid=this%sfield_num+2*this%hvfield_num+1, this%field_num_tot-3, 4
447 lcmesh%panelID, lcfpos3d(:,1), lcfpos3d(:,2), unity_fac(:), &
448 commdata%Nnode_LCMeshFace, &
449 commdata%recv_buf(:,varid), commdata%recv_buf(:,varid+1), &
450 tmp1_htensor3d(:,1,1), tmp1_htensor3d(:,2,1) )
452 lcmesh%panelID, lcfpos3d(:,1), lcfpos3d(:,2), unity_fac(:), &
453 commdata%Nnode_LCMeshFace, &
454 commdata%recv_buf(:,varid+2), commdata%recv_buf(:,varid+3), &
455 tmp1_htensor3d(:,1,2), tmp1_htensor3d(:,2,2) )
457 lcmesh%panelID, lcfpos3d(:,1), lcfpos3d(:,2), unity_fac(:), &
458 commdata%Nnode_LCMeshFace, &
459 tmp1_htensor3d(:,1,1), tmp1_htensor3d(:,1,2), &
460 this%recv_buf(irs:ire,varid,n), this%recv_buf(irs:ire,varid+2,n) )
462 lcmesh%panelID, lcfpos3d(:,1), lcfpos3d(:,2), unity_fac(:), &
463 commdata%Nnode_LCMeshFace, &
464 tmp1_htensor3d(:,2,1), tmp1_htensor3d(:,2,2), &
465 this%recv_buf(irs:ire,varid+1,n), this%recv_buf(irs:ire,varid+3,n) )
466 end do
467 deallocate( lcfpos3d, unity_fac, tmp1_htensor3d )
468 end if
469 end if
470
471 irs = ire + 1
472 end do
473
474 deallocate( fpos3d )
475 end do
476
477 return
478 end subroutine post_exchange_core
479
480!OCL SERIAL
481 subroutine push_localsendbuf( lc_send_buf, send_buf, s_faceID, is, Nnode_LCMeshFace, var_num, &
482 mesh3D, elem3D )
483 implicit none
484
485 integer, intent(in) :: var_num
486 integer, intent(in) :: nnode_lcmeshface
487 real(rp), intent(inout) :: lc_send_buf(nnode_lcmeshface,var_num)
488 real(rp), intent(in) :: send_buf(bufsize_per_field,var_num)
489 integer, intent(in) :: s_faceid, is
490 type(localmesh3d), pointer :: mesh3d
491 type(elementbase3d), intent(in) :: elem3d
492
493 integer :: ie
494 !-----------------------------------------------------------------------------
495
496 ie = is + nnode_lcmeshface - 1
497 if ( s_faceid > 0 ) then
498 lc_send_buf(:,:) = send_buf(is:ie,:)
499 else
500 call revert_hori( lc_send_buf, send_buf(is:ie,:), mesh3d, elem3d, mesh3d%lcmesh2D )
501 end if
502
503 return
504 contains
505 subroutine revert_hori( revert, ori, mesh, e3D, mesh2D )
506 implicit none
507 type(localmesh3d), intent(in) :: mesh
508 type(elementbase3d), intent(in) :: e3d
509 type(localmesh2d), intent(in) :: mesh2d
510 real(rp), intent(out) :: revert(e3d%nnode_h1d,e3d%nnode_v,mesh2d%nex,mesh%nez, var_num)
511 real(rp), intent(in) :: ori(e3d%nnode_h1d,e3d%nnode_v,mesh2d%nex,mesh%nez, var_num)
512
513 integer :: p1, p3, i, k, n
514 integer :: i_, p1_
515 !------------------------------------------------------------------------
516
517 do n=1, var_num
518 do k=1, mesh%NeZ
519 do i=1, mesh2d%NeX
520 i_ = mesh2d%NeX - i + 1
521 do p3=1, e3d%Nnode_v
522 do p1=1, e3d%Nnode_h1D
523 p1_ = e3d%Nnode_h1D - p1 + 1
524 revert(p1,p3,i,k,n) = ori(p1_,p3,i_,k,n)
525 end do
526 end do
527 end do
528 end do
529 end do
530
531 return
532 end subroutine revert_hori
533 end subroutine push_localsendbuf
534
535!OCL SERIAL
536 subroutine extract_boundary_data3d( var, elem, mesh, buf )
537 implicit none
538
539 type(elementbase3d), intent(in) :: elem
540 type(localmesh3d), intent(in) :: mesh
541 real(dp), intent(in) :: var(elem%np * mesh%ne)
542 real(dp), intent(inout) :: buf( 2*(mesh%nex + mesh%ney)*mesh%nez*elem%nfp_h &
543 + 2*mesh%nex*mesh%ney*elem%nfp_v )
544 !------------------------------------------------------------
545
546 buf(:) = var(mesh%VmapB(:))
547
548 return
549 end subroutine extract_boundary_data3d
550
551!OCL SERIAL
552 subroutine set_boundary_data3d_u1u2( buf_U, buf_V, &
553 elem, mesh, G_ij, &
554 u1, u2)
555
556 implicit none
557
558 type(elementbase3d), intent(in) :: elem
559 type(localmesh3d), intent(in) :: mesh
560 real(dp), intent(in) :: buf_u(2*(mesh%nex + mesh%ney)*mesh%nez*elem%nfp_h &
561 + 2*mesh%nex*mesh%ney*elem%nfp_v)
562 real(dp), intent(in) :: buf_v(2*(mesh%nex + mesh%ney)*mesh%nez*elem%nfp_h &
563 + 2*mesh%nex*mesh%ney*elem%nfp_v)
564 real(dp), intent(in) :: g_ij(elem%np * mesh%ne,2,2)
565 real(dp), intent(inout) :: u1(elem%np * mesh%nea)
566 real(dp), intent(inout) :: u2(elem%np * mesh%nea)
567 !------------------------------------------------------------
568
569 u1(elem%Np*mesh%NeE+1:elem%Np*mesh%NeE+size(buf_u)) &
570 = g_ij(mesh%VmapB,1,1) * buf_u(:) + g_ij(mesh%VmapB,1,2) * buf_v(:)
571 u2(elem%Np*mesh%NeE+1:elem%Np*mesh%NeE+size(buf_u)) &
572 = g_ij(mesh%VmapB,2,1) * buf_u(:) + g_ij(mesh%VmapB,2,2) * buf_v(:)
573
574 return
575 end subroutine set_boundary_data3d_u1u2
576
577end module scale_meshfieldcomm_cubedspheredom3d
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 / Local 3D
module FElib / Mesh / Cubed-sphere 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 in 3D cubed-sphere domain
Derived type to manage data communication between a face of local mesh.
Base derived type to manage data communication.