MFIX  2016-1
read_res0_des.f
Go to the documentation of this file.
1 !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvv!
2 ! !
3 ! Subroutine: DES_READ_RESTART !
4 ! Purpose : Reads either single restart file or multiple restart !
5 ! fles (based on bdist_io) flag. !
6 !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^!
7  SUBROUTINE read_res0_des
8 
9  use cdist
10  use compar
11  use des_allocate
12  use des_bc
13  use des_rxns
14  use des_thermo
15  use desmpi
16  use discretelement
17  use error_manager
18  use machine
19  use mfix_pic, only: mppic, des_stat_wt
20  use mpi_utility
21  use param, only: dimension_n_s
22  use read_res1_des
23  use run
24 
25  implicit none
26 
27  INTEGER :: LC1, LC2
28  INTEGER :: lDIMN, lNEXT_REC
29 
30  INTEGER :: lVAR_SIZE
31  DOUBLE PRECISION :: VERSION
32 
33  ldimn = merge(2,3,no_k)
34 
35  CALL init_read_res_des(trim(run_name), version, lnext_rec)
36 
37  CALL read_res_des(lnext_rec, vtp_findex)
38  CALL read_res_des(lnext_rec, tecplot_findex)
39  CALL read_res_des(lnext_rec, dtsolid)
40 
41 ! Position data is read and used to setup pARRAY reads.
42  CALL read_par_pos(lnext_rec)
43 
44  CALL read_res_parray(lnext_rec, iglobal_id)
45 
46  CALL read_res_parray(lnext_rec, particle_state)
47 
48  DO lc1 = 1, ldimn
49  CALL read_res_parray(lnext_rec, des_vel_new(:,lc1))
50  ENDDO
51 
52  DO lc1 = 1, merge(1,3,no_k)
53  CALL read_res_parray(lnext_rec, omega_new(:,lc1))
54  ENDDO
55 
56  CALL read_res_parray(lnext_rec, des_radius)
57  CALL read_res_parray(lnext_rec, ro_sol)
58 
59  IF(mppic) CALL read_res_parray(lnext_rec, des_stat_wt)
60  IF(energy_eq) CALL read_res_parray(lnext_rec, des_t_s)
61 
62  IF(any_species_eq) THEN
63  CALL read_res_parray(lnext_rec, pijk(:,5))
64  DO lc1=1, dimension_n_s
65  CALL read_res_parray(lnext_rec, des_x_s(:,lc1))
66  ENDDO
67  ENDIF
68 
69  IF(version >= 1.1) THEN
70  CALL read_res_des(lnext_rec, lvar_size)
71  DO lc1=1, lvar_size
72  if(lvar_size <= des_usr_var_size) &
73  CALL read_res_parray(lnext_rec, des_usr_var(lc1,:))
74  ENDDO
75  ENDIF
76 
77 ! RES2 does not need the collision of BC information.
78  IF(run_type == 'RESTART_2') RETURN
79 
80 ! Collision/neighbor data is read and used to setup cARRAY reads.
81  IF(.NOT.mppic) THEN
82  CALL read_par_col(lnext_rec)
83  DO lc1=1, ldimn
84  CALL read_res_carray(lnext_rec, pft_neighbor(lc1,:))
85  ENDDO
86  ENDIF
87 
88 ! Save the number of BCMI's read from input file, then read the
89 ! value from the restart file.
90  CALL read_res_des(lnext_rec, dem_bcmi)
91 
92 ! Allocation of MIs is done here to ignore changes to the mfix.dat
93 ! file during RES1.
94  IF(dem_bcmi > 0) CALL allocate_dem_mi
95 
96 ! Only save the number of mass inflows for RESTART_1. This allows
97 ! for mass inflows to be added/removed with RESTART_2.
98 ! Todo: Prune entering/exiting flagged particles for RESTART_2.
99  DO lc1=1, dem_bcmi
100  CALL read_res_des(lnext_rec, dem_mi_time(lc1))
101  CALL read_res_des(lnext_rec, dem_mi(lc1)%VACANCY)
102  CALL read_res_des(lnext_rec, dem_mi(lc1)%OCCUPANTS)
103  CALL read_res_des(lnext_rec, dem_mi(lc1)%WINDOW)
104  CALL read_res_des(lnext_rec, dem_mi(lc1)%OFFSET)
105  CALL read_res_des(lnext_rec, dem_mi(lc1)%L)
106 
107  lc2 = dem_mi(lc1)%OCCUPANTS
108 
109  allocate(dem_mi(lc1)%W(lc2))
110  CALL read_res_des(lnext_rec, dem_mi(lc1)%W(:))
111  allocate(dem_mi(lc1)%H(lc2))
112  CALL read_res_des(lnext_rec, dem_mi(lc1)%H(:))
113  allocate(dem_mi(lc1)%P(lc2))
114  CALL read_res_des(lnext_rec, dem_mi(lc1)%P(:))
115  allocate(dem_mi(lc1)%Q(lc2))
116  CALL read_res_des(lnext_rec, dem_mi(lc1)%Q(:))
117  ENDDO
118 
119  CALL finl_read_res_des
120 
121 
122  WRITE(err_msg,"('DES restart file read at Time = ',g12.5)") time
123  CALL flush_err_msg(header=.false., footer=.false.)
124 
125  RETURN
126  END SUBROUTINE read_res0_des
double precision, dimension(:), allocatable des_t_s
character(len=60) run_name
Definition: run_mod.f:24
subroutine, public allocate_dem_mi
subroutine, public finl_read_res_des
subroutine, public read_par_col(lNEXT_REC)
double precision, dimension(:), allocatable dem_mi_time
Definition: des_bc_mod.f:36
Definition: cdist_mod.f:2
character(len=16) run_type
Definition: run_mod.f:33
integer dem_bcmi
Definition: des_bc_mod.f:18
double precision, dimension(:,:), allocatable des_x_s
Definition: des_rxns_mod.f:21
logical any_species_eq
Definition: run_mod.f:118
subroutine, public read_par_pos(lNEXT_REC)
Definition: run_mod.f:13
Definition: param_mod.f:2
subroutine read_res0_des
Definition: read_res0_des.f:8
subroutine, public init_read_res_des(BASE, lVERSION, lNEXT_REC)
logical energy_eq
Definition: run_mod.f:100
character(len=line_length), dimension(line_count) err_msg
type(dem_mi_), dimension(:), allocatable, target dem_mi
Definition: des_bc_mod.f:90
logical mppic
Definition: mfix_pic_mod.f:14
integer dimension_n_s
Definition: param_mod.f:21
double precision time
Definition: run_mod.f:45
double precision, dimension(:), allocatable des_stat_wt
Definition: mfix_pic_mod.f:54
subroutine flush_err_msg(DEBUG, HEADER, FOOTER, ABORT, LOG, CALL_TREE)