MFIX  2016-1
make_arrays_des.f
Go to the documentation of this file.
1 !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvv!
2 ! Module name: MAKE_ARRAYS_DES !
3 ! Author: Jay Boyalakuntla Date: 12-Jun-04 !
4 ! !
5 ! Purpose: DES - allocating DES arrays
6 ! !
7 !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^!
8  SUBROUTINE make_arrays_des
9 
11  USE compar
12  use constant, only: pi
13  USE cutcell
14  USE des_rxns
15  USE des_thermo
16  USE desgrid
17  USE discretelement
18  USE error_manager
19  USE functions
20  USE funits
21  use fldvar, only: ro_s
23  USE geometry
24  use mfix_pic, only: mppic
26  USE mpi_utility
27  USE param1
28  use physprop, only: mmax, ro_s0
29  USE run
30  USE stl
32  use stl_preproc_des, only: add_facet
33 
34  IMPLICIT NONE
35 !-----------------------------------------------
36 ! Local variables
37 !-----------------------------------------------
38  INTEGER :: I, J, K, L, IJK
39  INTEGER :: I1, I2, J1, J2, K1, K2, II, JJ, KK, IJK2
40  INTEGER :: lcurpar, lpip_all(0:numpes-1), lglobal_id
41  INTEGER :: CELL_ID, I_CELL, J_CELL, K_CELL, COUNT, NF
42  INTEGER :: IMINUS1, IPLUS1, JMINUS1, JPLUS1, KMINUS1, KPLUS1
43 
44 ! MPPIC related quantities
45  CALL init_err_msg("MAKE_ARRAYS_DES")
46 
47 ! Check interpolation input.
48  CALL set_filter_des
49 
50 ! cfassign and des_init_bc called before reading the particle info
51  CALL cfassign
52 
53  vol_surr(:) = zero
54 
55  ! initialize VOL_SURR array
56  DO k = kstart2, kend1
57  DO j = jstart2, jend1
58  DO i = istart2, iend1
59  IF (dead_cell_at(i,j,k)) cycle ! skip dead cells
60  ijk = funijk(i,j,k)
61  i1 = i
62  i2 = i+1
63  j1 = j
64  j2 = j+1
65  k1 = k
66  k2 = merge(k, k+1, no_k)
67 
68 ! looping over stencil points (node values)
69  DO kk = k1, k2
70  DO jj = j1, j2
71  DO ii = i1, i2
72  IF (dead_cell_at(ii,jj,kk)) cycle ! skip dead cells
73  ijk2 = funijk_map_c(ii, jj, kk)
74  IF(fluid_at(ijk2)) vol_surr(ijk) = &
75  vol_surr(ijk)+vol(ijk2)
76  ENDDO
77  ENDDO
78  ENDDO
79  ENDDO
80  ENDDO
81  ENDDO
82 
83 
84 
85 ! Set the initial particle data.
86  IF(run_type == 'NEW') THEN
87  IF(particles /= 0) THEN
88  IF(gener_part_config) THEN
90  ELSE
91  CALL read_par_input
92  ENDIF
93  ENDIF
94 
95 ! Set the global ID for the particles and set the ghost cnt
96  ighost_cnt = 0
97  lpip_all = 0
98  lpip_all(mype) = pip
99  call global_all_sum(lpip_all)
100  lglobal_id = sum(lpip_all(0:mype-1))
101  imax_global_id = 0
102  do lcurpar = 1,pip
103  lglobal_id = lglobal_id + 1
104  iglobal_id(lcurpar) = lglobal_id
105  imax_global_id = iglobal_id(pip)
106  end do
107  call global_all_max(imax_global_id)
108 
109 ! Initialize old values
110  omega_new(:,:) = zero
111 
112 ! Particle orientation
113  IF(particle_orientation) THEN
114  orientation(1,:) = init_orientation(1)
115  orientation(2,:) = init_orientation(2)
116  orientation(3,:) = init_orientation(3)
117  ENDIF
118 
119 
120  IF (do_old) THEN
121  omega_old(:,:) = zero
122  des_pos_old(:,:) = des_pos_new(:,:)
123  des_vel_old(:,:) = des_vel_new(:,:)
124  ENDIF
125 
126 ! Read the restart file.
127  ELSEIF(run_type == 'RESTART_1' .OR. run_type == 'RESTART_2') THEN
128 
129  CALL read_res0_des
130  imax_global_id = maxval(iglobal_id(1:pip))
131  call global_all_max(imax_global_id)
132 
133 ! Initizlie the old values.
134  IF (do_old) THEN
135  omega_old(:,:) = omega_new(:,:)
136  des_pos_old(:,:) = des_pos_new(:,:)
137  des_vel_old(:,:) = des_vel_new(:,:)
138  ENDIF
139 
140  ELSE
141 
142  WRITE(err_msg, 1100)
143  CALL flush_err_msg(abort=.true.)
144  1100 FORMAT('Error 1100: Unsupported RUN_TYPE for DES.')
145 
146  ENDIF
147 
148  IF(run_type == 'RESTART_2') vtp_findex=0
149 
150 ! setting additional particle properties now that the particles
151 ! have been identified
152  DO l = 1, max_pip
153 ! Skip 'empty' locations when populating the particle property arrays.
154  IF(is_nonexistent(l)) cycle
155  IF(is_ghost(l) .OR. is_entering_ghost(l) .OR. is_exiting_ghost(l)) cycle
156  pvol(l) = (4.0d0/3.0d0)*pi*des_radius(l)**3
157  pmass(l) = pvol(l)*ro_sol(l)
158  omoi(l) = 2.5d0/(pmass(l)*des_radius(l)**2) !ONE OVER MOI
159  ENDDO
160 
161  CALL set_phase_index
163 
164 ! do_nsearch should be set before calling particle in cell
165  do_nsearch =.true.
166 ! Bin the particles to the DES grid.
167  CALL desgrid_pic(plocate=.true.)
168  CALL des_par_exchange
169  CALL particles_in_cell
170 
171  IF(dem_solids) THEN
172  CALL neighbour
173  CALL init_settling_dem
174  ENDIF
175 
176  IF(run_type == 'NEW') CALL set_ic_dem
177 
178 
179 ! Calculate interpolation weights
181 ! Calculate mean fields using either interpolation or cell averaging.
182  CALL comp_mean_fields
183 
184  IF(mppic) CALL calc_dtpic
185 
186 
187  IF(run_type /= 'RESTART_1' .AND. print_des_data) THEN
188  s_time = time
189  CALL write_des_data
190  ENDIF
191 
192  CALL finl_err_msg
193 
194  RETURN
195  END SUBROUTINE make_arrays_des
subroutine read_par_input
subroutine comp_mean_fields
logical dem_solids
Definition: run_mod.f:257
double precision, dimension(:), allocatable vol_surr
Definition: geometry_mod.f:215
integer kend1
Definition: compar_mod.f:80
subroutine write_des_data
subroutine finl_err_msg
integer iend1
Definition: compar_mod.f:80
subroutine init_settling_dem
integer istart2
Definition: compar_mod.f:80
subroutine desgrid_pic(plocate)
Definition: desgrid_mod.f:711
subroutine cfassign
Definition: cfassign.f:21
subroutine set_phase_index
subroutine des_par_exchange()
subroutine set_filter_des
integer kstart2
Definition: compar_mod.f:80
subroutine add_facet(IJK, FACET_ID)
integer numpes
Definition: compar_mod.f:24
subroutine calc_interp_weights
subroutine init_err_msg(CALLER)
integer mmax
Definition: physprop_mod.f:19
logical, dimension(:,:,:), allocatable dead_cell_at
Definition: compar_mod.f:127
Definition: stl_mod.f:1
integer jstart2
Definition: compar_mod.f:80
character(len=16) run_type
Definition: run_mod.f:33
subroutine set_ic_dem
Definition: set_ic_dem.f:14
Definition: run_mod.f:13
subroutine particles_in_cell
subroutine calc_dtpic
Definition: calc_dtpic.f:9
double precision, dimension(:,:), allocatable ro_s
Definition: fldvar_mod.f:45
logical no_k
Definition: geometry_mod.f:28
subroutine read_res0_des
Definition: read_res0_des.f:8
subroutine neighbour
Definition: neighbour.f:16
integer mype
Definition: compar_mod.f:24
subroutine init_particles_in_cell
character(len=line_length), dimension(line_count) err_msg
logical mppic
Definition: mfix_pic_mod.f:14
double precision, dimension(dim_m) ro_s0
Definition: physprop_mod.f:28
double precision, parameter pi
Definition: constant_mod.f:158
double precision time
Definition: run_mod.f:45
double precision, dimension(:), allocatable vol
Definition: geometry_mod.f:212
double precision, parameter zero
Definition: param1_mod.f:27
subroutine make_arrays_des
integer jend1
Definition: compar_mod.f:80
subroutine flush_err_msg(DEBUG, HEADER, FOOTER, ABORT, LOG, CALL_TREE)