FE-Project
Loading...
Searching...
No Matches
scale_meshfieldcomm_rectdom2d.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
19 use scale_element_base, only: &
23 use scale_meshfieldcomm_base, only: &
30
31 !-----------------------------------------------------------------------------
32 implicit none
33 private
34
35 !-----------------------------------------------------------------------------
36 !
37 !++ Public type & procedure
38 !
39
40 type, public, extends(meshfieldcommbase) :: meshfieldcommrectdom2d
41 class(meshrectdom2d), pointer :: mesh2d
42 contains
43 procedure, public :: init => meshfieldcommrectdom2d_init
44 procedure, public :: put => meshfieldcommrectdom2d_put
45 procedure, public :: get => meshfieldcommrectdom2d_get
46 procedure, public :: exchange => meshfieldcommrectdom2d_exchange
47 procedure, public :: final => meshfieldcommrectdom2d_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 = 4
67
68contains
69 subroutine meshfieldcommrectdom2d_init( this, &
70 sfield_num, hvfield_num, htensorfield_num, mesh2d )
71 implicit none
72
73 class(meshfieldcommrectdom2d), intent(inout) :: this
74 integer, intent(in) :: sfield_num
75 integer, intent(in) :: hvfield_num
76 integer, intent(in) :: htensorfield_num
77 class(meshrectdom2d), intent(in), target :: mesh2d
78
79 type(localmesh2d), pointer :: lcmesh
80 integer :: n
81 integer :: nnode_lcmeshface(comm_face_num,mesh2d%local_mesh_num)
82 !-----------------------------------------------------------------------------
83
84 this%mesh2d => mesh2d
85 lcmesh => mesh2d%lcmesh_list(1)
86 bufsize_per_field = 2*(lcmesh%NeX + lcmesh%NeY)*lcmesh%refElem2D%Nfp
87
88 do n=1, this%mesh2d%LOCAL_MESH_NUM
89 lcmesh => this%mesh2d%lcmesh_list(n)
90 nnode_lcmeshface(:,n) = (/ lcmesh%NeX, lcmesh%NeY, lcmesh%NeX, lcmesh%NeY /) * lcmesh%refElem2D%Nfp
91 end do
92
93 call meshfieldcommbase_init( this, sfield_num, hvfield_num, htensorfield_num, bufsize_per_field, comm_face_num, nnode_lcmeshface, mesh2d)
94
95 return
96 end subroutine meshfieldcommrectdom2d_init
97
98 subroutine meshfieldcommrectdom2d_final( this )
99 implicit none
100
101 class(meshfieldcommrectdom2d), intent(inout) :: this
102 !-----------------------------------------------------------------------------
103
104 call meshfieldcommbase_final( this )
105
106 return
107 end subroutine meshfieldcommrectdom2d_final
108
109 subroutine meshfieldcommrectdom2d_put(this, field_list, varid_s)
110 implicit none
111 class(meshfieldcommrectdom2d), intent(inout) :: this
112 type(meshfieldcontainer), intent(in) :: field_list(:)
113 integer, intent(in) :: varid_s
114
115 integer :: i
116 integer :: n
117 type(localmesh2d), pointer :: lcmesh
118 !-----------------------------------------------------------------------------
119
120 do i=1, size(field_list)
121 do n=1, this%mesh%LOCAL_MESH_NUM
122 lcmesh => this%mesh2d%lcmesh_list(n)
123 call meshfieldcommbase_extract_bounddata( field_list(i)%field2d%local(n)%val, lcmesh%refElem, lcmesh, & ! (in)
124 this%send_buf(:,varid_s+i-1,n) ) ! (out)
125 end do
126 end do
127
128 return
129 end subroutine meshfieldcommrectdom2d_put
130
131 subroutine meshfieldcommrectdom2d_get(this, field_list, varid_s)
132 use scale_meshfieldcomm_base, only: &
134 implicit none
135
136 class(meshfieldcommrectdom2d), intent(inout) :: this
137 type(meshfieldcontainer), intent(inout) :: field_list(:)
138 integer, intent(in) :: varid_s
139
140 integer :: i
141 integer :: n
142 type(localmesh2d), pointer :: lcmesh
143 !-----------------------------------------------------------------------------
144
145 !--
146 if ( this%call_wait_flag_sub_get ) &
147 call meshfieldcommbase_wait_core( this, this%commdata_list )
148
149 do i=1, size(field_list)
150 do n=1, this%mesh2d%LOCAL_MESH_NUM
151 lcmesh => this%mesh2d%lcmesh_list(n)
152 call meshfieldcommbase_set_bounddata( this%recv_buf(:,varid_s+i-1,n), lcmesh%refElem, lcmesh, & !(in)
153 field_list(i)%field2d%local(n)%val ) !(out)
154 end do
155 end do
156
157 return
158 end subroutine meshfieldcommrectdom2d_get
159
160!OCL SERIAL
161 subroutine meshfieldcommrectdom2d_exchange( this, do_wait )
162 use scale_meshfieldcomm_base, only: &
165
166 implicit none
167
168 class(meshfieldcommrectdom2d), intent(inout), target :: this
169 logical, intent(in), optional :: do_wait
170
171 integer :: n, f
172 type(localmeshcommdata), pointer :: commdata
173 !-----------------------------------------------------------------------------
174
175 do n=1, this%mesh%LOCAL_MESH_NUM
176 do f=1, this%nfaces_comm
177 commdata => this%commdata_list(f,n)
178 call push_localsendbuf( commdata%send_buf(:,:), & ! (inout)
179 this%send_buf(:,:,n), commdata%s_faceID, this%is_f(f,n), & ! (in)
180 commdata%Nnode_LCMeshFace, this%field_num_tot) ! (in)
181 end do
182 end do
183
184 !-----------------------
185
186 call meshfieldcommbase_exchange_core( this, this%commdata_list(:,:), do_wait )
187
188 !---------------------
189
190 return
191 end subroutine meshfieldcommrectdom2d_exchange
192
193!----------------------------
194
195!OCL SERIAL
196 subroutine push_localsendbuf( lc_send_buf, send_buf, s_faceID, is, Nnode_LCMeshFace, var_num )
197 implicit none
198
199 integer, intent(in) :: var_num
200 integer, intent(in) :: nnode_lcmeshface
201 real(rp), intent(inout) :: lc_send_buf(nnode_lcmeshface,var_num)
202 real(rp), intent(in) :: send_buf(bufsize_per_field,var_num)
203 integer, intent(in) :: s_faceid, is
204
205 integer :: ie
206 !-----------------------------------------------------------------------------
207
208 ie = is + nnode_lcmeshface - 1
209 if ( s_faceid > 0 ) then
210 lc_send_buf(:,:) = send_buf(is:ie,:)
211 else
212 lc_send_buf(:,:) = send_buf(ie:is:-1,:)
213 end if
214
215 return
216 end subroutine push_localsendbuf
217
218end module scale_meshfieldcomm_rectdom2d
module FElib / Element / Base
module FElib / Mesh / Local 2D
module FElib / Mesh / Rectangle 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 2D rectangle domain
Derived type to manage data communication between a face of local mesh.
Base derived type to manage data communication.