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