MFIX  2016-1
des_init_arrays.f
Go to the documentation of this file.
1 !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvv!
2 ! !
3 ! Subrourtine: DES_INIT_ARRAYS !
4 ! Author: Jay Boyalakuntla Date: 12-Jun-04 !
5 ! !
6 ! Purpose: Initialize arrays at the start of the simulation. Note that !
7 ! arrays based on the number of particles (MAX_PIP) should be added to !
8 ! the DES_INIT_PARTICLE_ARRAYS as they need to be reinitialized after !
9 ! particle arrays are grown (see PARTICLE_GROW). !
10 ! !
11 !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^!
12  SUBROUTINE des_init_arrays
13 
14  USE param
15  USE param1
16  USE discretelement
17  USE indices
18  USE geometry
19  USE compar
20  USE physprop
21  USE des_bc
22  USE run
23  use desgrid
24  use desmpi
25  USE des_thermo
26  USE des_rxns
27 
28  IMPLICIT NONE
29 
30 !-----------------------------------------------
31 
32  pinc(:) = zero
33 
34  p_force(:,:) = zero
35 
36  IF(allocated(drag_am)) drag_am = zero
37  IF(allocated(drag_bm)) drag_bm = zero
38 
39  f_gds = zero
40  vxf_gds = zero
41 
42  IF (des_continuum_hybrid) THEN
43  f_sds = zero
44  vxf_sds = zero
45  sdrag_am = zero
46  sdrag_bm = zero
47  ENDIF
48 
49  grav(:) = zero
50 
51  IF(allocated(avgdes_t_s)) avgdes_t_s = zero
52  IF(allocated(conv_sp)) conv_sp = zero
53  IF(allocated(conv_sc)) conv_sc = zero
54 
55  CALL des_init_particle_arrays(1,max_pip)
56 
57  RETURN
58  END SUBROUTINE des_init_arrays
59 
60 !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvv!
61 ! !
62 ! Subrourtine: DES_INIT_PARTICLE_ARRAYS !
63 ! Author: Jay Boyalakuntla Date: 12-Jun-04 !
64 ! !
65 ! Purpose: Initialize particle arrays. The upper and lower bounds are !
66 ! passed so that after resizing particle arrays (see GROW_PARTICLE) the !
67 ! new portions of the arrays can be initialized. !
68 ! !
69 !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^!
70  SUBROUTINE des_init_particle_arrays(LB,UB)
71 
72 !-----------------------------------------------
73 ! Modules
74 !-----------------------------------------------
75  use des_rxns
76  use des_thermo
77  use desgrid
78  use desmpi
79  use discretelement
80  use functions
81  use mfix_pic, only: avgsolvel_p, epg_p
82  use mfix_pic, only: mppic, des_stat_wt, ps_grad
83  use param1, only: zero
85  use particle_filter, only: filter_size
86  use run, only: any_species_eq
87  use run, only: energy_eq
88  use des_thermo_cond, only: des_qw_cond
89 
90  IMPLICIT NONE
91 
92  INTEGER, INTENT(IN) :: LB, UB
93  INTEGER :: II
94 
95  iglobal_id(lb:ub) = 0
96 
97 ! Physical properties:
98  des_radius(lb:ub) = zero
99  ro_sol(lb:ub) = zero
100  pvol(lb:ub) = huge(0.0)
101  pmass(lb:ub) = huge(0.0)
102  omoi(lb:ub) = zero
103 
104  ! Particle position, velocity, etc
105  des_pos_new(lb:ub,:) = zero
106  des_vel_new(lb:ub,:) = zero
107  omega_new(lb:ub,:) = zero
108  ppos(lb:ub,:) = zero
109  IF(particle_orientation) THEN
110  orientation(1,:) = init_orientation(1)
111  orientation(2,:) = init_orientation(2)
112  orientation(3,:) = init_orientation(3)
113  ENDIF
114 
115 ! Particle state flag
116  DO ii = lb, ub
117  call set_nonexistent(ii)
118  ENDDO
119  neighbor_index(:) = 0
120 
121 ! DES grid bin information
122  dg_pijk(lb:ub) = -1
123  dg_pijkprv(lb:ub) = -1
124  ighost_updated(lb:ub) = .false.
125 
126 ! Fluid cell bin information
127  pijk(lb:ub,:) = 0
128 
129 ! Translation and rotational forces
130  fc(lb:ub,:) = zero
131  tow(lb:ub,:) = zero
132 
133 ! Collision data
134  wall_collision_facet_id(:,lb:ub) = -1
135  wall_collision_pft(:,:,lb:ub) = zero
136 
137 ! Initializing user defined array
138  IF(des_usr_var_size > 0) &
139  des_usr_var(:,lb:ub) = zero
140 
141 ! Particle center drag coefficient and explicit drag force
142  f_gp(lb:ub) = zero
143  drag_fc(lb:ub,:) = zero
144 
145 ! Interpolation variables.
146  IF(filter_size > 0)THEN
147  filter_cell(:,lb:ub) = -1
148  filter_weight(:,lb:ub) = zero
149  ENDIF
150 
151 ! MPPIC variables
152  IF(mppic) THEN
153  des_stat_wt(lb:ub) = zero
154  ps_grad(:,lb:ub) = zero
155  avgsolvel_p(:,lb:ub) = zero
156  epg_p(lb:ub) = zero
157  ENDIF
158 
159 ! Higher order time integration variables.
160  IF (do_old) THEN
161  des_pos_old(lb:ub,:) = zero
162  des_vel_old(lb:ub,:) = zero
163  des_acc_old(lb:ub,:) = zero
164  omega_old(lb:ub,:) = zero
165  rot_acc_old(lb:ub,:) = zero
166  ENDIF
167 
168 ! Energy equation variables.
169  IF(energy_eq)THEN
170  des_t_s(lb:ub) = zero
171  des_c_ps(lb:ub) = zero
172  des_x_s(lb:ub,:) = zero
173 
174  IF(ALLOCATED(q_source)) q_source(lb:ub) = zero
175  IF(ALLOCATED(conv_qs)) conv_qs(lb:ub) = zero
176  IF(ALLOCATED(gammaxsa)) gammaxsa(lb:ub) = zero
177  IF (intg_adams_bashforth) &
178  q_source0(lb:ub) = zero
179  IF (ALLOCATED(des_qw_cond)) &
180  des_qw_cond(:,:) = zero
181  ENDIF
182 
183 ! Chemical reaction variables.
184  IF(any_species_eq)THEN
185  des_r_s(lb:ub,:) = zero
186  IF (intg_adams_bashforth) THEN
187  dmdt_old(lb:ub) = zero
188  dxdt_old(lb:ub,:) = zero
189  ENDIF
190  rxns_qs(lb:ub) = zero
191  ENDIF
192 
193 
194 ! Cohesion VDW forces
195  IF(use_cohesion) THEN
196  postcohesive(lb:ub) = zero
197  ENDIF
198 
199 
200  RETURN
201  END SUBROUTINE des_init_particle_arrays
202 
203 
204 
double precision, dimension(:), allocatable gammaxsa
subroutine des_init_arrays
double precision, dimension(:,:), allocatable des_r_s
Definition: des_rxns_mod.f:24
double precision, dimension(:), allocatable conv_sc
integer, dimension(:,:), allocatable filter_cell
double precision, dimension(:,:), allocatable avgsolvel_p
Definition: mfix_pic_mod.f:28
double precision, dimension(:,:), allocatable ps_grad
Definition: mfix_pic_mod.f:74
double precision, dimension(:), allocatable des_t_s
double precision, dimension(:), allocatable q_source0
double precision, dimension(:), allocatable avgdes_t_s
subroutine des_init_particle_arrays(LB, UB)
double precision, dimension(:,:), allocatable filter_weight
double precision, dimension(:,:), allocatable des_x_s
Definition: des_rxns_mod.f:21
logical any_species_eq
Definition: run_mod.f:118
double precision, dimension(:,:), allocatable dxdt_old
Definition: des_rxns_mod.f:33
Definition: run_mod.f:13
double precision, dimension(:), allocatable conv_sp
Definition: param_mod.f:2
double precision, dimension(:), allocatable conv_qs
logical energy_eq
Definition: run_mod.f:100
double precision, dimension(:), allocatable epg_p
Definition: mfix_pic_mod.f:31
logical mppic
Definition: mfix_pic_mod.f:14
double precision, dimension(:), allocatable q_source
double precision, dimension(:,:), allocatable des_qw_cond
double precision, dimension(:), allocatable rxns_qs
double precision, dimension(:), allocatable des_stat_wt
Definition: mfix_pic_mod.f:54
double precision, parameter zero
Definition: param1_mod.f:27
double precision, dimension(:), allocatable des_c_ps
double precision, dimension(:), allocatable dmdt_old
Definition: des_rxns_mod.f:31