FE-Project
Loading...
Searching...
No Matches
scale_model_var_manager.F90
Go to the documentation of this file.
1!-------------------------------------------------------------------------------
10!-------------------------------------------------------------------------------
11#include "scaleFElib.h"
13 !-----------------------------------------------------------------------------
14 !
15 !++ Used modules
16 !
17 use scale_precision
18 use scale_io
19
20 use scale_file_history, only: file_history_reg
21 use scale_file_monitor_meshfield, only: file_monitor_meshfield_reg
22
23
24 use scale_meshfield_base, only: &
26 use scale_linkedlist, only: &
28 use scale_meshfieldcomm_base, only: &
30
31 use scale_variableinfo, only: &
33
34 !-----------------------------------------------------------------------------
35 implicit none
36 private
37 !-----------------------------------------------------------------------------
38 !
39 !++ Public type & procedures
40 !
41 type, public :: modelvarmanager
42 type(linkedlist) :: list
43 class(meshfieldcommbase), pointer :: ptr_comm
44 type(meshfieldcontainer), allocatable :: comm_list(:)
45 contains
46 procedure, public :: init => modelvarmanager_init
47 procedure, public :: final => modelvarmanager_final
48 procedure, private :: regist1d => modelvarmanager_regist1d
49 procedure, private :: regist2d => modelvarmanager_regist2d
50 procedure, private :: regist3d => modelvarmanager_regist3d
51 generic, public :: regist => regist1d, regist2d, regist3d
52 procedure, public :: get => modelvarmanager_get
53 procedure, public :: get2d => modelvarmanager_get2d
54 procedure, public :: get3d => modelvarmanager_get3d
55 procedure, public :: getlocalmeshfield => modelvarmanager_getlocalmeshfield
56 procedure, public :: getlocalmeshfieldlist => modelvarmanager_getlocalmeshfieldlist
57
58 !--
59 procedure, public :: meshfieldcomm_prepair => modelvarmanager_meshfiled_comm_prepare
60 procedure, public :: meshfieldcomm_exchange => modelvarmanager_meshfiled_comm_exchange
61 procedure, public :: meshfieldcomm_get => modelvarmanager_meshfiled_comm_get
62 end type modelvarmanager
63
64 public :: variableinfo ! Cascade
65
66 !-----------------------------------------------------------------------------
67 !
68 !++ Public parameters & variables
69 !
70
71 !-----------------------------------------------------------------------------
72 !
73 !++ Private procedures & variables
74 !
75 !------------------
76
77contains
78
79!OCL SERIAL
80 subroutine modelvarmanager_init( this )
81 implicit none
82 class(modelvarmanager), intent(inout) :: this
83 !------------------------------------------------
84
85 call this%list%Init()
86 nullify( this%ptr_comm )
87
88 return
89 end subroutine modelvarmanager_init
90
91!OCL SERIAL
92 subroutine modelvarmanager_final( this )
93 implicit none
94 class(modelvarmanager), intent(inout) :: this
95 !------------------------------------------------
96
97 call this%list%Traverse( field_release )
98
99 call this%list%Final()
100
101 nullify( this%ptr_comm )
102 if ( allocated(this%comm_list) ) deallocate(this%comm_list)
103
104 return
105 end subroutine modelvarmanager_final
106
107!OCL SERIAL
108 subroutine field_release( key, pField, done)
109 implicit none
110 class(*), intent(in) :: key
111 class(*), pointer :: pField
112 logical, intent(out) :: done
113 !------------------------------------------------
114
115 select type( pfield )
116 type is ( meshfield1d )
117 call pfield%Final()
118 type is ( meshfield2d )
119 call pfield%Final()
120 type is ( meshfield3d )
121 call pfield%Final()
122 end select
123
124 return
125 end subroutine field_release
126
127!OCL SERIAL
128 subroutine modelvarmanager_regist1d( this, &
129 varinfo, mesh, field, reg_file_history_flag, &
130 monitor_flag, fill_zero )
131
133 implicit none
134
135 class(modelvarmanager), intent(inout) :: this
136 type(variableinfo), intent(in) :: varinfo
137 class(meshbase1d), intent(in) :: mesh
138 class(meshfield1d), intent(inout), target :: field
139 logical, intent(in) :: reg_file_history_flag
140 logical, intent(in), optional :: monitor_flag
141 logical, intent(in), optional :: fill_zero
142
143 class(*), pointer :: ptr_field
144
145 integer :: domID, ke
146 logical :: fill_zero_
147 !------------------------------------------------
148
149 if (present(fill_zero)) then
150 fill_zero_ = fill_zero
151 else
152 fill_zero_ = .true.
153 end if
154
155 call field%Init( varinfo%NAME, varinfo%UNIT, mesh )
156 if (reg_file_history_flag) then
157 call file_history_reg( varinfo%NAME, varinfo%DESC, varinfo%UNIT, field%hist_id, dim_type=varinfo%dim_type )
158 end if
159 if ( present(monitor_flag) ) then
160 if (monitor_flag) then
161 call file_monitor_meshfield_reg( &
162 varinfo%NAME, varinfo%DESC, trim(varinfo%UNIT)//'*m', & ! (in)
163 field%monitor_id, dim_type='ATM1D', is_tendency=.false. ) ! (in)
164 end if
165 end if
166
167 ptr_field => field
168 call this%list%AddByPointer( varinfo%keyID, ptr_field )
169
170 if ( fill_zero_ ) then
171 do domid=1, mesh%LOCAL_MESH_NUM
172 !$omp parallel do
173 do ke=mesh%lcmesh_list(domid)%NeS, mesh%lcmesh_list(domid)%NeE
174 field%local(domid)%val(:,ke) = 0.0_rp
175 end do
176 end do
177 end if
178
179 return
180 end subroutine modelvarmanager_regist1d
181
182!OCL SERIAL
183 subroutine modelvarmanager_regist2d( this, &
184 varinfo, mesh, field, reg_file_history_flag, &
185 monitor_flag, fill_zero )
187 implicit none
188
189 class(modelvarmanager), intent(inout) :: this
190 type(variableinfo), intent(in) :: varinfo
191 class(meshbase2d), intent(in) :: mesh
192 class(meshfield2d), intent(inout), target :: field
193 logical, intent(in) :: reg_file_history_flag
194 logical, intent(in), optional :: monitor_flag
195 logical, intent(in), optional :: fill_zero
196
197 class(*), pointer :: ptr_field
198
199 integer :: domID, ke
200 logical :: fill_zero_
201 !------------------------------------------------
202
203 if (present(fill_zero)) then
204 fill_zero_ = fill_zero
205 else
206 fill_zero_ = .false.
207 end if
208
209 call field%Init( varinfo%NAME, varinfo%UNIT, mesh )
210 if (reg_file_history_flag) then
211 call file_history_reg( varinfo%NAME, varinfo%DESC, varinfo%UNIT, field%hist_id, dim_type=varinfo%dim_type)
212 end if
213 if ( present(monitor_flag) ) then
214 if (monitor_flag) then
215 call file_monitor_meshfield_reg( &
216 varinfo%NAME, varinfo%DESC, trim(varinfo%UNIT)//'*m2', & ! (in)
217 field%monitor_id, dim_type='ATM2D', is_tendency=.false. ) ! (in)
218 end if
219 end if
220
221 ptr_field => field
222 call this%list%AddByPointer( varinfo%keyID, ptr_field )
223
224 if ( fill_zero_ ) then
225 do domid=1, mesh%LOCAL_MESH_NUM
226 !$omp parallel do
227 do ke=mesh%lcmesh_list(domid)%NeS, mesh%lcmesh_list(domid)%NeE
228 field%local(domid)%val(:,ke) = 0.0_rp
229 end do
230 end do
231 end if
232
233 return
234 end subroutine modelvarmanager_regist2d
235
236!OCL SERIAL
237 subroutine modelvarmanager_regist3d( this, &
238 varinfo, mesh, field, reg_file_history_flag, &
239 monitor_flag, fill_zero )
241 implicit none
242
243 class(modelvarmanager), intent(inout) :: this
244 type(variableinfo), intent(in) :: varinfo
245 class(meshbase3d), intent(in) :: mesh
246 class(meshfield3d), intent(inout), target :: field
247 logical, intent(in) :: reg_file_history_flag
248 logical, intent(in), optional :: monitor_flag
249 logical, intent(in), optional :: fill_zero
250
251 class(*), pointer :: ptr_field
252
253 integer :: domID, ke
254 logical :: fill_zero_
255 !------------------------------------------------
256
257 if (present(fill_zero)) then
258 fill_zero_ = fill_zero
259 else
260 fill_zero_ = .false.
261 end if
262
263 call field%Init( varinfo%NAME, varinfo%UNIT, mesh )
264 if (reg_file_history_flag) then
265 call file_history_reg( field%varname, varinfo%DESC, field%unit, field%hist_id, &
266 dim_type=varinfo%dim_type )
267 end if
268 if ( present(monitor_flag) ) then
269 if (monitor_flag) then
270 call file_monitor_meshfield_reg( &
271 varinfo%NAME, varinfo%DESC, trim(varinfo%UNIT)//'*m3', & ! (in)
272 field%monitor_id, dim_type='ATM3D', is_tendency=.false. ) ! (in)
273 end if
274 end if
275
276 ptr_field => field
277 call this%list%AddByPointer( varinfo%keyID, ptr_field )
278
279 if ( fill_zero_ ) then
280 do domid=1, mesh%LOCAL_MESH_NUM
281 !$omp parallel do
282 do ke=mesh%lcmesh_list(domid)%NeS, mesh%lcmesh_list(domid)%NeE
283 field%local(domid)%val(:,ke) = 0.0_rp
284 end do
285 end do
286 end if
287
288 return
289 end subroutine modelvarmanager_regist3d
290
291!OCL SERIAL
292 subroutine modelvarmanager_get( this, keyID, pField )
294 implicit none
295
296 class(modelvarmanager), intent(inout) :: this
297 integer, intent(in) :: keyID
298 class(meshfieldbase), pointer, intent(out) :: pField
299
300 class(*), pointer :: ptr_field
301 !------------------------------------------------
302
303 call this%list%Get(keyid, ptr_field)
304
305 nullify( pfield )
306 select type( ptr_field )
307 type is (meshfield1d)
308 pfield => ptr_field
309 type is (meshfield2d)
310 pfield => ptr_field
311 type is (meshfield3d)
312 pfield => ptr_field
313 end select
314
315 return
316 end subroutine modelvarmanager_get
317
318!OCL SERIAL
319 subroutine modelvarmanager_get2d( this, keyID, pField )
321 implicit none
322
323 class(modelvarmanager), intent(inout) :: this
324 integer, intent(in) :: keyID
325 class(meshfield2d), pointer, intent(out) :: pField
326
327 class(*), pointer :: ptr_field
328 !------------------------------------------------
329
330 call this%list%Get(keyid, ptr_field)
331
332 nullify( pfield )
333 select type( ptr_field )
334 type is (meshfield2d)
335 pfield => ptr_field
336 end select
337
338 return
339 end subroutine modelvarmanager_get2d
340
341!OCL SERIAL
342 subroutine modelvarmanager_get3d( this, keyID, pField )
344 implicit none
345
346 class(modelvarmanager), intent(inout) :: this
347 integer, intent(in) :: keyID
348 class(meshfield3d), pointer, intent(out) :: pField
349
350 class(*), pointer :: ptr_field
351 !------------------------------------------------
352
353 call this%list%Get(keyid, ptr_field)
354
355 nullify( pfield )
356 select type( ptr_field )
357 type is (meshfield3d)
358 pfield => ptr_field
359 end select
360
361 return
362 end subroutine modelvarmanager_get3d
363
364!OCL SERIAL
365 subroutine modelvarmanager_getlocalmeshfield( this, keyID, domID, pField_lc )
368 implicit none
369
370 class(modelvarmanager), intent(inout) :: this
371 integer, intent(in) :: keyID
372 integer, intent(in) :: domID
373 class(localmeshfieldbase), pointer, intent(out) :: pField_lc
374
375 class(meshfieldbase), pointer :: pField
376 !------------------------------------------------
377
378 call this%Get(keyid, pfield)
379 call pfield%GetLocalMeshField(domid, pfield_lc)
380
381 return
382 end subroutine modelvarmanager_getlocalmeshfield
383
384!OCL SERIAL
385 subroutine modelvarmanager_getlocalmeshfieldlist( this, keyID_list, domID, lcfield_list )
388 implicit none
389
390 class(modelvarmanager), intent(inout) :: this
391 integer, intent(in) :: keyID_list(:)
392 integer, intent(in) :: domID
393 type(localmeshfieldbaselist), intent(out) :: lcfield_list(size(keyID_list))
394
395 integer :: i
396 integer :: keyID
397
398 class(meshfieldbase), pointer :: pField
399 !------------------------------------------------
400
401 do i=1, size(keyid_list)
402 keyid = keyid_list(i)
403 call this%Get(keyid, pfield)
404 call pfield%GetLocalMeshField( domid, lcfield_list(i)%ptr )
405 end do
406
407 return
408 end subroutine modelvarmanager_getlocalmeshfieldlist
409
410!OCL SERIAL
411 subroutine modelvarmanager_meshfiled_comm_prepare( this, &
412 comm, fields )
413 implicit none
414 class(modelvarmanager), intent(inout) :: this
415 class(meshfieldcommbase), target, intent(in) :: comm
416 class(meshfieldbase), target, intent(in) :: fields(:)
417
418 integer :: v
419 integer :: nFields
420 class(meshfieldbase), pointer :: ptr_field
421 !------------------------------------------------
422
423 this%ptr_comm => comm
424
425 nfields = size(fields)
426 allocate( this%comm_list(nfields) )
427
428 do v = 1, size(fields)
429 ptr_field => fields(v)
430 select type( ptr_field )
431 type is (meshfield1d)
432 this%comm_list(v)%field1d => ptr_field
433 type is (meshfield2d)
434 this%comm_list(v)%field2d => ptr_field
435 type is (meshfield3d)
436 this%comm_list(v)%field3d => ptr_field
437 end select
438 end do
439
440 return
441 end subroutine modelvarmanager_meshfiled_comm_prepare
442
443!OCL SERIAL
444 subroutine modelvarmanager_meshfiled_comm_exchange( this, do_wait )
445 implicit none
446 class(modelvarmanager), intent(inout) :: this
447 logical, intent(in), optional :: do_wait
448
449 logical :: do_wait_
450 !------------------------------------------------
451
452 call this%ptr_comm%Put( this%comm_list, 1 )
453 call this%ptr_comm%Exchange( do_wait )
454
455 if ( present(do_wait) ) then
456 do_wait_ = do_wait
457 else
458 do_wait_ = .true.
459 end if
460 if ( do_wait_ ) &
461 call this%ptr_comm%Get( this%comm_list, 1 )
462
463 return
464 end subroutine modelvarmanager_meshfiled_comm_exchange
465
466!OCL SERIAL
467 subroutine modelvarmanager_meshfiled_comm_get( this )
468 implicit none
469 class(modelvarmanager), intent(inout) :: this
470 !------------------------------------------------
471 call this%ptr_comm%Get( this%comm_list, 1 )
472 return
473 end subroutine modelvarmanager_meshfiled_comm_get
474end module scale_model_var_manager
module common / data collection
module FElib / Mesh / Base 1D
module FElib / Mesh / Base 2D
module FElib / Mesh / Base 3D
module FElib / Data / base
module FElib / Data / Communication base
FElib / model framework / variable manager.
module FElib / Data / base
Base derived type to manage data communication.