File: N:\mfix\model\des\read_res0_des.f
1
2
3
4
5
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
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
78 IF(RUN_TYPE == 'RESTART_2') RETURN
79
80
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
89
90 CALL READ_RES_DES(lNEXT_REC, DEM_BCMI)
91
92
93
94 IF(DEM_BCMI > 0) CALL ALLOCATE_DEM_MI
95
96
97
98
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
127