MFIX  2016-1
write_res0_des.f
Go to the documentation of this file.
1 !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvc
2 !
3 ! module name: des_write_restart
4 ! purpose: writing des data for restart
5 !
6 ! Author : Pradeep G
7 ! Purpose : Reads either single restart file or multiple restart files
8 ! (based on bdist_io) flag
9 !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^c
10  subroutine write_res0_des
11 
12  use compar
13  use des_bc
14  use des_rxns
15  use des_thermo
16  use discretelement
17  use error_manager
18  use mfix_pic, only: des_stat_wt
19  use mfix_pic, only: mppic
20  use mpi_utility
21  use param, only: dimension_n_s
22  use param1
23  use run
24  use write_res1_des
25 
26  implicit none
27 !-----------------------------------------------
28 ! local variables
29 !-----------------------------------------------
30  INTEGER :: LC1
31  INTEGER :: lNEXT_REC
32  INTEGER :: lDIMN
33 
34  DOUBLE PRECISION :: VERSION
35 
36 ! Set the version of the DES RES file.
37  version = 1.1
38 
39 ! Set the output dimension.
40  ldimn = merge(2,3,no_k)
41 
42  CALL init_write_res_des(trim(run_name), version, lnext_rec)
43 
44  CALL write_res_des(lnext_rec, vtp_findex)
45  CALL write_res_des(lnext_rec, tecplot_findex)
46  CALL write_res_des(lnext_rec, dtsolid)
47 
48  DO lc1 = 1, ldimn
49  CALL write_res_parray(lnext_rec, des_pos_new(:,lc1))
50  ENDDO
51 
52  CALL write_res_parray(lnext_rec, iglobal_id)
53 
54  CALL write_res_parray(lnext_rec, particle_state)
55 
56  DO lc1 = 1, ldimn
57  CALL write_res_parray(lnext_rec, des_vel_new(:,lc1))
58  ENDDO
59 
60  DO lc1 = 1, merge(1,3,no_k)
61  CALL write_res_parray(lnext_rec, omega_new(:,lc1))
62  ENDDO
63 
64  CALL write_res_parray(lnext_rec, des_radius)
65  CALL write_res_parray(lnext_rec, ro_sol)
66 
67  IF(mppic) &
68  CALL write_res_parray(lnext_rec, des_stat_wt)
69 
70  IF(energy_eq) &
71  CALL write_res_parray(lnext_rec, des_t_s)
72 
73  IF(any_species_eq) THEN
74  CALL write_res_parray(lnext_rec, pijk(:,5))
75  DO lc1=1, dimension_n_s
76  CALL write_res_parray(lnext_rec, des_x_s(:,lc1))
77  ENDDO
78  ENDIF
79 
80 ! DES User defined variable :: added for VERSION>= 1.1
81  CALL write_res_des(lnext_rec, des_usr_var_size)
82  DO lc1=1,des_usr_var_size
83  CALL write_res_parray(lnext_rec, des_usr_var(lc1,:))
84  ENDDO
85 
86 
87  IF(.NOT.mppic) THEN
88  CALL write_res_carray(lnext_rec, neighbors(:), ploc2glb=.true.)
89  CALL write_res_carray(lnext_rec, neighbor_index(:), ploc2glb=.true.)
90  ENDIF
91 
92  DO lc1=1, ldimn
93  CALL write_res_carray(lnext_rec,pft_neighbor(lc1,:))
94  ENDDO
95 
96  CALL write_res_des(lnext_rec, dem_bcmi)
97  DO lc1=1, dem_bcmi
98  CALL write_res_des(lnext_rec, dem_mi_time(lc1))
99  CALL write_res_des(lnext_rec, dem_mi(lc1)%VACANCY)
100  CALL write_res_des(lnext_rec, dem_mi(lc1)%OCCUPANTS)
101  CALL write_res_des(lnext_rec, dem_mi(lc1)%WINDOW)
102  CALL write_res_des(lnext_rec, dem_mi(lc1)%OFFSET)
103  CALL write_res_des(lnext_rec, dem_mi(lc1)%L)
104  CALL write_res_des(lnext_rec, dem_mi(lc1)%W(:))
105  CALL write_res_des(lnext_rec, dem_mi(lc1)%H(:))
106  CALL write_res_des(lnext_rec, dem_mi(lc1)%P(:))
107  CALL write_res_des(lnext_rec, dem_mi(lc1)%Q(:))
108  ENDDO
109 
110  CALL finl_write_res_des
111 
112  RETURN
113  END SUBROUTINE write_res0_des
double precision, dimension(:), allocatable des_t_s
character(len=60) run_name
Definition: run_mod.f:24
subroutine, public init_write_res_des(BASE, lVERSION, lNEXT_REC)
double precision, dimension(:), allocatable dem_mi_time
Definition: des_bc_mod.f:36
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 write_res0_des
Definition: run_mod.f:13
Definition: param_mod.f:2
logical energy_eq
Definition: run_mod.f:100
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, dimension(:), allocatable des_stat_wt
Definition: mfix_pic_mod.f:54
subroutine, public finl_write_res_des