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