FE-Project
Loading...
Searching...
No Matches
scale_meshfield_base.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
26
27 use scale_localmeshfield_base, only: &
30
31 use scale_element_base, only: &
33
34
35 !-----------------------------------------------------------------------------
36 implicit none
37 private
38
39 !-----------------------------------------------------------------------------
40 !
41 !++ Public type & procedure
42 !
43
44 type, abstract, public :: meshfieldbase
45 character(len=H_SHORT) :: varname
46 character(len=H_SHORT) :: unit
47 integer :: hist_id
48 integer :: monitor_id
49 contains
50 procedure(meshfieldbase_get_localmeshfield), deferred, public :: getlocalmeshfield
51 end type meshfieldbase
52
53 interface
54 subroutine meshfieldbase_get_localmeshfield( this, domID, ptr_lcmeshField )
55 import meshfieldbase
57 class(meshfieldbase), target, intent(in) :: this
58 integer, intent(in) :: domID
59 class(localmeshfieldbase), pointer, intent(out) :: ptr_lcmeshfield
61 end interface
62
63 !------
64
65 type, extends(meshfieldbase), public :: meshfield1d
66 type(localmeshfield1d), allocatable :: local(:)
67 class(meshbase1d), pointer :: mesh
68 contains
69 procedure :: init => meshfield1d_init
70 procedure :: final => meshfield1d_final
71 procedure :: getlocalmeshfield => meshfield1d_get_localmeshfield
72 end type meshfield1d
73
74 type, extends(meshfieldbase), public :: meshfield2d
75 type(localmeshfield2d), allocatable :: local(:)
76 class(meshbase2d), pointer :: mesh
77 contains
78 procedure :: init => meshfield2d_init
79 procedure :: final => meshfield2d_final
80 procedure :: getlocalmeshfield => meshfield2d_get_localmeshfield
81 end type meshfield2d
82
83 type, extends(meshfieldbase), public :: meshfield3d
84 type(localmeshfield3d), allocatable :: local(:)
85 class(meshbase3d), pointer :: mesh
86 contains
87 procedure :: init => meshfield3d_init
88 procedure :: final => meshfield3d_final
89 procedure :: getlocalmeshfield => meshfield3d_get_localmeshfield
90 end type meshfield3d
91
92 type, public :: meshfield3dlist
93 class(meshfield3d), pointer :: ptr
94 end type meshfield3dlist
95
96 !-----------------------------------------------------------------------------
97 !
98 !++ Public parameters & variables
99 !
100
101 !-----------------------------------------------------------------------------
102 !
103 !++ Private procedure
104 !
105
106 !-----------------------------------------------------------------------------
107 !
108 !++ Private parameters & variables
109 !
110
111contains
112 !**** 1D **********************************************************************
113
114!OCL SERIAL
115 subroutine meshfield1d_init( this, varname, units, mesh, data_type )
116 implicit none
117 class(meshfield1d), intent(inout) :: this
118 character(len=*), intent(in) :: varname
119 character(len=*), intent(in) :: units
120 class(meshbase1d), target, intent(in) :: mesh
121 integer, intent(in), optional :: data_type
122
123 integer :: n
124 !-----------------------------------------------------------------------------
125
126 this%varname = varname
127 this%unit = units
128 this%mesh => mesh
129 this%hist_id = -1
130 this%monitor_id = -1
131
132 allocate( this%local(mesh%LOCAL_MESH_NUM) )
133 do n=1, mesh%LOCAL_MESH_NUM
134 call this%local(n)%Init( mesh%lcmesh_list(n), data_type )
135 end do
136
137 return
138 end subroutine meshfield1d_init
139
140!OCL SERIAL
141 subroutine meshfield1d_final(this)
142 implicit none
143 class(meshfield1d), intent(inout) :: this
144
145 integer :: n
146 !-----------------------------------------------------------------------------
147
148 do n=1, size(this%local)
149 call this%local(n)%Final()
150 end do
151 deallocate( this%local )
152
153 return
154 end subroutine meshfield1d_final
155
156!OCL SERIAL
157 subroutine meshfield1d_get_localmeshfield( this, domID, ptr_lcmeshField )
158 implicit none
159
160 class(meshfield1d), target, intent(in) :: this
161 integer, intent(in) :: domID
162 class(localmeshfieldbase), pointer, intent(out) :: ptr_lcmeshfield
163 !-----------------------------------------------------------------------------
164
165 ptr_lcmeshfield => this%local(domid)
166 return
167 end subroutine meshfield1d_get_localmeshfield
168
169 !**** 2D **********************************************************************
170
171!OCL SERIAL
172 subroutine meshfield2d_init( this, varname, units, mesh, data_type )
173 implicit none
174
175 class(meshfield2d), intent(inout) :: this
176 character(len=*), intent(in) :: varname
177 character(len=*), intent(in) :: units
178 class(meshbase2d), target, intent(in) :: mesh
179 integer, intent(in), optional :: data_type
180
181 integer :: n
182 !-----------------------------------------------------------------------------
183
184 this%varname = varname
185 this%unit = units
186 this%mesh => mesh
187 this%hist_id = -1
188
189 allocate( this%local(mesh%LOCAL_MESH_NUM) )
190 do n=1, mesh%LOCAL_MESH_NUM
191 call this%local(n)%Init( mesh%lcmesh_list(n), data_type )
192 end do
193
194 return
195 end subroutine meshfield2d_init
196
197!OCL SERIAL
198 subroutine meshfield2d_final(this)
199 implicit none
200
201 class(meshfield2d), intent(inout) :: this
202
203 integer :: n
204 !-----------------------------------------------------------------------------
205
206 do n=1, size(this%local)
207 call this%local(n)%Final()
208 end do
209 deallocate( this%local )
210
211 return
212 end subroutine meshfield2d_final
213
214!OCL SERIAL
215 subroutine meshfield2d_get_localmeshfield( this, domID, ptr_lcmeshField )
216 implicit none
217
218 class(meshfield2d), target, intent(in) :: this
219 integer, intent(in) :: domID
220 class(localmeshfieldbase), pointer, intent(out) :: ptr_lcmeshfield
221 !-----------------------------------------------------------------------------
222
223 ptr_lcmeshfield => this%local(domid)
224 return
225 end subroutine meshfield2d_get_localmeshfield
226
227 !**** 3D **********************************************************************
228
229!OCL SERIAL
230 subroutine meshfield3d_init( this, varname, units, mesh, data_type )
231 implicit none
232
233 class(meshfield3d), intent(inout) :: this
234 character(len=*), intent(in) :: varname
235 character(len=*), intent(in) :: units
236 class(meshbase3d), target, intent(in) :: mesh
237 integer, intent(in), optional :: data_type
238
239 integer :: n
240 !-----------------------------------------------------------------------------
241
242 this%varname = varname
243 this%unit = units
244 this%mesh => mesh
245 this%hist_id = -1
246
247 allocate( this%local(mesh%LOCAL_MESH_NUM) )
248 do n=1, mesh%LOCAL_MESH_NUM
249 call this%local(n)%Init( mesh%lcmesh_list(n), data_type )
250 end do
251
252 return
253 end subroutine meshfield3d_init
254
255!OCL SERIAL
256 subroutine meshfield3d_final(this)
257 implicit none
258
259 class(meshfield3d), intent(inout) :: this
260
261 integer :: n
262 !-----------------------------------------------------------------------------
263
264 do n=1, size(this%local)
265 call this%local(n)%Final()
266 end do
267 deallocate( this%local )
268
269 return
270 end subroutine meshfield3d_final
271
272 subroutine meshfield3d_get_localmeshfield( this, domID, ptr_lcmeshField )
273 implicit none
274
275 class(meshfield3d), target, intent(in) :: this
276 integer, intent(in) :: domID
277 class(localmeshfieldbase), pointer, intent(out) :: ptr_lcmeshfield
278 !-----------------------------------------------------------------------------
279
280 ptr_lcmeshfield => this%local(domid)
281 return
282 end subroutine meshfield3d_get_localmeshfield
283
284end module scale_meshfield_base
module FElib / Element / Base
module FElib / Mesh / Local 1D
module FElib / Mesh / Local 2D
module FElib / Mesh / Local 3D
module FElib / Mesh / Base 1D
module FElib / Mesh / Base 2D
module FElib / Mesh / Base 3D
module FElib / Data / base
subroutine meshfield1d_init(this, varname, units, mesh, data_type)