44 procedure,
public :: init => meshfieldcommcubedom3d_init
45 procedure,
public :: put => meshfieldcommcubedom3d_put
46 procedure,
public :: get => meshfieldcommcubedom3d_get
47 procedure,
public :: exchange => meshfieldcommcubedom3d_exchange
48 procedure,
public :: final => meshfieldcommcubedom3d_final
60 private :: push_localsendbuf
67 integer,
parameter :: comm_face_num = 6
71 subroutine meshfieldcommcubedom3d_init( this, &
72 sfield_num, hvfield_num, htensorfield_num, mesh3d )
77 integer,
intent(in) :: sfield_num
78 integer,
intent(in) :: hvfield_num
79 integer,
intent(in) :: htensorfield_num
85 integer :: nnode_lcmeshface(comm_face_num,mesh3d%local_mesh_num)
89 lcmesh => mesh3d%lcmesh_list(1)
90 elem => lcmesh%refElem3D
92 + 2*lcmesh%NeX*lcmesh%NeY*elem%Nfp_v
94 do n=1, this%mesh3d%LOCAL_MESH_NUM
95 lcmesh => this%mesh3d%lcmesh_list(n)
96 nnode_lcmeshface(:,n) = &
97 (/ lcmesh%NeX, lcmesh%NeY, lcmesh%NeX, lcmesh%NeY, 0, 0 /) * lcmesh%NeZ * lcmesh%refElem3D%Nfp_h &
98 + (/ 0, 0, 0, 0, 1, 1 /) * lcmesh%NeX*lcmesh%NeY * lcmesh%refElem3D%Nfp_v
104 end subroutine meshfieldcommcubedom3d_init
107 subroutine meshfieldcommcubedom3d_final( this )
115 end subroutine meshfieldcommcubedom3d_final
118 subroutine meshfieldcommcubedom3d_put(this, field_list, varid_s)
122 integer,
intent(in) :: varid_s
143 end subroutine meshfieldcommcubedom3d_put
146 subroutine meshfieldcommcubedom3d_get(this, field_list, varid_s)
153 integer,
intent(in) :: varid_s
160 if ( this%call_wait_flag_sub_get )
then
163 field_list, 3, varid_s, this%mesh3d%lcmesh_list )
167 do i=1,
size(field_list)
168 do n=1, this%mesh3d%LOCAL_MESH_NUM
169 lcmesh => this%mesh3d%lcmesh_list(n)
171 field_list(i)%field3d%local(n)%val )
178 end subroutine meshfieldcommcubedom3d_get
184 subroutine meshfieldcommcubedom3d_exchange( this, do_wait )
192 logical,
intent(in),
optional :: do_wait
200 do n=1, this%mesh%LOCAL_MESH_NUM
201 lcmesh => this%mesh3d%lcmesh_list(n)
202 do f=1, this%nfaces_comm
203 commdata => this%commdata_list(f,n)
204 call push_localsendbuf( commdata%send_buf(:,:), &
205 this%send_buf(:,:,n), commdata%s_faceID, this%is_f(f,n), &
206 commdata%Nnode_LCMeshFace, this%field_num_tot, lcmesh )
216 end subroutine meshfieldcommcubedom3d_exchange
221 subroutine push_localsendbuf( lc_send_buf, send_buf, s_faceID, is, Nnode_LCMeshFace, var_num, lcmesh )
224 integer,
intent(in) :: nnode_lcmeshface
225 integer,
intent(in) :: var_num
227 real(rp),
intent(out) :: lc_send_buf(nnode_lcmeshface,var_num)
229 integer,
intent(in) :: s_faceid, is
236 if ( s_faceid > 0 )
then
239 do i=1, nnode_lcmeshface
240 lc_send_buf(i,vid) = send_buf(is+i-1,vid)
243 else if ( -5 < s_faceid .and. s_faceid < 0)
then
244 e3d => lcmesh%refElem3D
245 ie = is + nnode_lcmeshface - 1
246 call revert_hori( lc_send_buf(:,:), send_buf(is:ie,:), nnode_lcmeshface/lcmesh%NeZ, lcmesh%NeZ )
251 subroutine revert_hori(revert, ori, Ne_h1D, NeZ)
252 integer,
intent(in) :: ne_h1d, nez
253 real(rp),
intent(out) :: revert(e3d%nnode_h1d, e3d%nnode_v, ne_h1d, nez, var_num)
254 real(rp),
intent(in) :: ori(e3d%nnode_h1d, e3d%nnode_v, ne_h1d, nez, var_num)
256 integer :: p1, p3, i, k, n
265 do p1=1, e3d%Nnode_h1D
266 p1_ = e3d%Nnode_h1D - p1 + 1
267 revert(p1,p3,i,k,n) = ori(p1_,p3,i_,k,n)
273 end subroutine revert_hori
274 end subroutine push_localsendbuf
module FElib / Element / Base
module FElib / Mesh / Local 3D
module FElib / Mesh / Cubic 3D 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_extract_bounddata_2(field_list, dim, varid_s, lcmesh_list, buf)
Extract halo data from data array with MeshField object and set it to the recieving 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 3D cubic domain
integer bufsize_per_field
Buffer size per a field.
Derived type representing a 3D reference element.
Derived type representing an arbitrary finite element.
Derived type to manage a local 3D computational domain.
Derived type representing a field with 3D 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 3D cubic domain.