11#include "scaleFElib.h"
20 use scale_file_history,
only: file_history_reg
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
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
86 nullify( this%ptr_comm )
92 subroutine modelvarmanager_final( this )
97 call this%list%Traverse( field_release )
99 call this%list%Final()
101 nullify( this%ptr_comm )
102 if (
allocated(this%comm_list) )
deallocate(this%comm_list)
105 end subroutine modelvarmanager_final
108 subroutine field_release( key, pField, done)
110 class(*),
intent(in) :: key
111 class(*),
pointer :: pField
112 logical,
intent(out) :: done
115 select type( pfield )
125 end subroutine field_release
128 subroutine modelvarmanager_regist1d( this, &
129 varinfo, mesh, field, reg_file_history_flag, &
130 monitor_flag, fill_zero )
139 logical,
intent(in) :: reg_file_history_flag
140 logical,
intent(in),
optional :: monitor_flag
141 logical,
intent(in),
optional :: fill_zero
143 class(*),
pointer :: ptr_field
146 logical :: fill_zero_
149 if (
present(fill_zero))
then
150 fill_zero_ = fill_zero
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 )
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', &
163 field%monitor_id, dim_type=
'ATM1D', is_tendency=.false. )
168 call this%list%AddByPointer( varinfo%keyID, ptr_field )
170 if ( fill_zero_ )
then
171 do domid=1, mesh%LOCAL_MESH_NUM
173 do ke=mesh%lcmesh_list(domid)%NeS, mesh%lcmesh_list(domid)%NeE
174 field%local(domid)%val(:,ke) = 0.0_rp
180 end subroutine modelvarmanager_regist1d
183 subroutine modelvarmanager_regist2d( this, &
184 varinfo, mesh, field, reg_file_history_flag, &
185 monitor_flag, fill_zero )
193 logical,
intent(in) :: reg_file_history_flag
194 logical,
intent(in),
optional :: monitor_flag
195 logical,
intent(in),
optional :: fill_zero
197 class(*),
pointer :: ptr_field
200 logical :: fill_zero_
203 if (
present(fill_zero))
then
204 fill_zero_ = fill_zero
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)
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', &
217 field%monitor_id, dim_type=
'ATM2D', is_tendency=.false. )
222 call this%list%AddByPointer( varinfo%keyID, ptr_field )
224 if ( fill_zero_ )
then
225 do domid=1, mesh%LOCAL_MESH_NUM
227 do ke=mesh%lcmesh_list(domid)%NeS, mesh%lcmesh_list(domid)%NeE
228 field%local(domid)%val(:,ke) = 0.0_rp
234 end subroutine modelvarmanager_regist2d
237 subroutine modelvarmanager_regist3d( this, &
238 varinfo, mesh, field, reg_file_history_flag, &
239 monitor_flag, fill_zero )
247 logical,
intent(in) :: reg_file_history_flag
248 logical,
intent(in),
optional :: monitor_flag
249 logical,
intent(in),
optional :: fill_zero
251 class(*),
pointer :: ptr_field
254 logical :: fill_zero_
257 if (
present(fill_zero))
then
258 fill_zero_ = fill_zero
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 )
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', &
272 field%monitor_id, dim_type=
'ATM3D', is_tendency=.false. )
277 call this%list%AddByPointer( varinfo%keyID, ptr_field )
279 if ( fill_zero_ )
then
280 do domid=1, mesh%LOCAL_MESH_NUM
282 do ke=mesh%lcmesh_list(domid)%NeS, mesh%lcmesh_list(domid)%NeE
283 field%local(domid)%val(:,ke) = 0.0_rp
289 end subroutine modelvarmanager_regist3d
292 subroutine modelvarmanager_get( this, keyID, pField )
297 integer,
intent(in) :: keyID
300 class(*),
pointer :: ptr_field
303 call this%list%Get(keyid, ptr_field)
306 select type( ptr_field )
316 end subroutine modelvarmanager_get
319 subroutine modelvarmanager_get2d( this, keyID, pField )
324 integer,
intent(in) :: keyID
327 class(*),
pointer :: ptr_field
330 call this%list%Get(keyid, ptr_field)
333 select type( ptr_field )
339 end subroutine modelvarmanager_get2d
342 subroutine modelvarmanager_get3d( this, keyID, pField )
347 integer,
intent(in) :: keyID
350 class(*),
pointer :: ptr_field
353 call this%list%Get(keyid, ptr_field)
356 select type( ptr_field )
362 end subroutine modelvarmanager_get3d
365 subroutine modelvarmanager_getlocalmeshfield( this, keyID, domID, pField_lc )
371 integer,
intent(in) :: keyID
372 integer,
intent(in) :: domID
378 call this%Get(keyid, pfield)
379 call pfield%GetLocalMeshField(domid, pfield_lc)
382 end subroutine modelvarmanager_getlocalmeshfield
385 subroutine modelvarmanager_getlocalmeshfieldlist( this, keyID_list, domID, lcfield_list )
391 integer,
intent(in) :: keyID_list(:)
392 integer,
intent(in) :: domID
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 )
408 end subroutine modelvarmanager_getlocalmeshfieldlist
411 subroutine modelvarmanager_meshfiled_comm_prepare( this, &
416 class(meshfieldbase),
target,
intent(in) :: fields(:)
420 class(meshfieldbase),
pointer :: ptr_field
423 this%ptr_comm => comm
425 nfields =
size(fields)
426 allocate( this%comm_list(nfields) )
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
441 end subroutine modelvarmanager_meshfiled_comm_prepare
444 subroutine modelvarmanager_meshfiled_comm_exchange( this, do_wait )
447 logical,
intent(in),
optional :: do_wait
452 call this%ptr_comm%Put( this%comm_list, 1 )
453 call this%ptr_comm%Exchange( do_wait )
455 if (
present(do_wait) )
then
461 call this%ptr_comm%Get( this%comm_list, 1 )
464 end subroutine modelvarmanager_meshfiled_comm_exchange
467 subroutine modelvarmanager_meshfiled_comm_get( this )
471 call this%ptr_comm%Get( this%comm_list, 1 )
473 end subroutine modelvarmanager_meshfiled_comm_get
474end module scale_model_var_manager
module FElib / File / Monitor
module common / data collection
module FElib / Data / base
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.
subroutine modelvarmanager_init(this)
module FElib / Data / base
Base derived type to manage data communication.