MFIX  2016-1
mpi_node_des_mod.f
Go to the documentation of this file.
1 !------------------------------------------------------------------------
2 ! Module : desmpi
3 ! Purpose : Contains wrapper class for mpi communications- send,recv
4 !
5 ! Author : Pradeep.G
6 !
7 ! Purpose : Module contains subroutines and variables related to
8 ! des mpi communication.
9 !
10 ! Comments : do_nsearch flag should be set to true before calling
11 ! des_par_exchange; when do_nsearch is true ghost particles of the
12 ! system will be updated, which will be later used to generate
13 ! neighbour list.
14 !
15 ! Contains following subroutines:
16 ! des_addnodevalues, des_addnodevalues2, des_addnodevalues_mean_fields
17 !------------------------------------------------------------------------
18  module mpi_node_des
19 
20 !-----------------------------------------------
21 ! Modules
22 !-----------------------------------------------
23  use parallel_mpi
24  use mpi_utility
25  use discretelement
26  use desgrid
27  use compar
28  use physprop
29  use sendrecv
30  use des_bc
31  use desmpi_wrapper
32  use sendrecvnode
33  use mfix_pic
34  use des_thermo
35  use run, only: energy_eq,any_species_eq
36  use param, only: dimension_n_s
37  use des_rxns
38  use desmpi
39 
42 
43  contains
44 
45 
46 !------------------------------------------------------------------------
47 ! Subroutine : des_addnodevalues_mean_fields
48 ! Purpose : This routine is specially used for computing mean
49 ! fields by backward interpolation.
50 !
51 ! Parameters : None
52 !------------------------------------------------------------------------
54 
55  use functions
56 !-----------------------------------------------
57  implicit none
58 !-----------------------------------------------
59 ! local variables
60 !-----------------------------------------------
61  integer :: li, lj, lk
62  integer :: lm,lijkmin,lijkmax
63 !-----------------------------------------------
64 
65 ! fill the temporary buffer
66  DO lm = 1,des_mmax+mmax
67  CALL des_exchangenode(des_rops_node(:,lm),padd=.true.)
68  DO li =1,dimn
69  CALL des_exchangenode(des_vel_node(:,li,lm),padd=.true.)
70  END DO
71  END DO
72 
73 ! adjust for periodic boundaries with no domain decomposition
74  if (des_periodic_walls_x .and. nodesi.eq.1) then
75  do lk = kstart2,kend2
76  do lj = jstart2,jend2
77  lijkmin = funijk(1,lj,lk)
78  lijkmax = funijk(imax1,lj,lk)
79  des_rops_node(lijkmin,:) = des_rops_node(lijkmin,:)+des_rops_node(lijkmax,:)
80  des_vel_node(lijkmin,:,:) = des_vel_node(lijkmin,:,:)+des_vel_node(lijkmax,:,:)
81  des_rops_node(lijkmax,:) = des_rops_node(lijkmin,:)
82  des_vel_node(lijkmax,:,:) = des_vel_node(lijkmin,:,:)
83  end do
84  end do
85  end if
86  if (des_periodic_walls_y .and. nodesj.eq.1) then
87  do lk = kstart2,kend2
88  do li = istart2,iend2
89  lijkmin = funijk(li,1,lk)
90  lijkmax = funijk(li,jmax1,lk)
91  des_rops_node(lijkmin,:) = des_rops_node(lijkmin,:)+des_rops_node(lijkmax,:)
92  des_vel_node(lijkmin,:,:) = des_vel_node(lijkmin,:,:)+des_vel_node(lijkmax,:,:)
93  des_rops_node(lijkmax,:) = des_rops_node(lijkmin,:)
94  des_vel_node(lijkmax,:,:) = des_vel_node(lijkmin,:,:)
95  end do
96  end do
97  end if
98  if (des_periodic_walls_z .and. nodesk.eq.1 .and. do_k) then
99  do li = istart2,iend2
100  do lj = jstart2,jend2
101  lijkmin = funijk(li,lj,1)
102  lijkmax = funijk(li,lj,kmax1)
103  des_rops_node(lijkmin,:) = des_rops_node(lijkmin,:)+des_rops_node(lijkmax,:)
104  des_vel_node(lijkmin,:,:) = des_vel_node(lijkmin,:,:)+des_vel_node(lijkmax,:,:)
105  des_rops_node(lijkmax,:) = des_rops_node(lijkmin,:)
106  des_vel_node(lijkmax,:,:) = des_vel_node(lijkmin,:,:)
107  end do
108  end do
109  end if
110 
111  return
112 
113  end subroutine des_addnodevalues_mean_fields
114 
115 
116 
117 !------------------------------------------------------------------------
118 ! Subroutine : des_addnodevalues
119 ! Purpose : This routine is specially used for des_drag_gs
120 ! The backward interpolation in des_drag_gs computes
121 ! the grid node values of drag_am and drag_bm
122 ! node values are from istart2 to iend1;
123 ! hence a separate module is created to exchange
124 ! node values
125 !
126 ! Parameters : None
127 !------------------------------------------------------------------------
128  subroutine des_addnodevalues()
130  use functions
131 !-----------------------------------------------
132  implicit none
133 !-----------------------------------------------
134 ! local variables
135 !-----------------------------------------------
136  integer :: li, lj, lk
137  integer :: lijkmin,lijkmax
138 !-----------------------------------------------
139 
140 ! fill the temporary buffer
141  call des_exchangenode(drag_am, padd=.true.)
142  do li =1,dimn
143  call des_exchangenode(drag_bm(:,li), padd=.true.)
144  end do
145 
146 ! adjust for periodic boundaries with no domain decomposition
147  if (des_periodic_walls_x .and. nodesi.eq.1) then
148  do lk = kstart2,kend2
149  do lj = jstart2,jend2
150  lijkmin = funijk(1,lj,lk)
151  lijkmax = funijk(imax1,lj,lk)
152  drag_am(lijkmin) = drag_am(lijkmin)+drag_am(lijkmax)
153  drag_bm(lijkmin,:) = drag_bm(lijkmin,:)+drag_bm(lijkmax,:)
154  drag_am(lijkmax) = drag_am(lijkmin)
155  drag_bm(lijkmax,:) = drag_bm(lijkmin,:)
156  end do
157  end do
158  end if
159  if (des_periodic_walls_y .and. nodesj.eq.1) then
160  do lk = kstart2,kend2
161  do li = istart2,iend2
162  lijkmin = funijk(li,1,lk)
163  lijkmax = funijk(li,jmax1,lk)
164  drag_am(lijkmin) = drag_am(lijkmin)+drag_am(lijkmax)
165  drag_bm(lijkmin,:) = drag_bm(lijkmin,:)+drag_bm(lijkmax,:)
166  drag_am(lijkmax) = drag_am(lijkmin)
167  drag_bm(lijkmax,:) = drag_bm(lijkmin,:)
168  end do
169  end do
170  end if
171  if (des_periodic_walls_z .and. nodesk.eq.1 .and. do_k) then
172  do li = istart2,iend2
173  do lj = jstart2,jend2
174  lijkmin = funijk(li,lj,1)
175  lijkmax = funijk(li,lj,kmax1)
176  drag_am(lijkmin) = drag_am(lijkmin)+drag_am(lijkmax)
177  drag_bm(lijkmin,:) = drag_bm(lijkmin,:)+drag_bm(lijkmax,:)
178  drag_am(lijkmax) = drag_am(lijkmin)
179  drag_bm(lijkmax,:) = drag_bm(lijkmin,:)
180  end do
181  end do
182  end if
183 
184  return
185 
186  end subroutine des_addnodevalues
187 
188 
189 !------------------------------------------------------------------------
190 ! Subroutine : des_addnodevalues2
191 ! Purpose : This routine is specially used for calc_des_rop_s
192 ! The backward interpolation in calc_des_rop_s computes
193 ! the grid node values of des_rops_node
194 ! node values are from istart2 to iend1;
195 ! hence a separate module is created to exchange
196 ! node values
197 !
198 ! Parameters : None
199 !------------------------------------------------------------------------
200  subroutine des_addnodevalues2()
202  use functions
203 !-----------------------------------------------
204  implicit none
205 !-----------------------------------------------
206 ! local variables
207 !-----------------------------------------------
208  integer :: li, lj, lk
209  integer :: lm,lijkmin,lijkmax
210 !-----------------------------------------------
211 
212 ! fill the temporary buffer
213  do lm = 1,des_mmax+mmax
214  call des_exchangenode(des_rops_node(:,lm),padd=.true.)
215  end do
216 
217 ! adjust for periodic boundaries with no domain decomposition
218  if (des_periodic_walls_x .and. nodesi.eq.1) then
219  do lk = kstart2,kend2
220  do lj = jstart2,jend2
221  lijkmin = funijk(1,lj,lk)
222  lijkmax = funijk(imax1,lj,lk)
223  des_rops_node(lijkmin,:) = des_rops_node(lijkmin,:)+des_rops_node(lijkmax,:)
224  des_rops_node(lijkmax,:) = des_rops_node(lijkmin,:)
225  end do
226  end do
227  end if
228  if (des_periodic_walls_y .and. nodesj.eq.1) then
229  do lk = kstart2,kend2
230  do li = istart2,iend2
231  lijkmin = funijk(li,1,lk)
232  lijkmax = funijk(li,jmax1,lk)
233  des_rops_node(lijkmin,:) = des_rops_node(lijkmin,:)+des_rops_node(lijkmax,:)
234  des_rops_node(lijkmax,:) = des_rops_node(lijkmin,:)
235  end do
236  end do
237  end if
238  if (des_periodic_walls_z .and. nodesk.eq.1 .and. do_k) then
239  do li = istart2,iend2
240  do lj = jstart2,jend2
241  lijkmin = funijk(li,lj,1)
242  lijkmax = funijk(li,lj,kmax1)
243  des_rops_node(lijkmin,:) = des_rops_node(lijkmin,:)+des_rops_node(lijkmax,:)
244  des_rops_node(lijkmax,:) = des_rops_node(lijkmin,:)
245  end do
246  end do
247  end if
248 
249  return
250 
251  end subroutine des_addnodevalues2
252 
253 
254 
255  end module mpi_node_des
integer jend2
Definition: compar_mod.f:80
subroutine des_addnodevalues_mean_fields()
subroutine desmpi_sendrecv_init(pface, pdebug)
integer istart2
Definition: compar_mod.f:80
subroutine desmpi_sendrecv_wait(pface, pdebug)
subroutine des_addnodevalues()
integer iend2
Definition: compar_mod.f:80
subroutine, public des_exchangenode(pvar, padd)
integer kend2
Definition: compar_mod.f:80
integer kstart2
Definition: compar_mod.f:80
subroutine des_addnodevalues2()
integer mmax
Definition: physprop_mod.f:19
integer jstart2
Definition: compar_mod.f:80
logical any_species_eq
Definition: run_mod.f:118
Definition: run_mod.f:13
Definition: param_mod.f:2
integer nodesj
Definition: compar_mod.f:37
logical energy_eq
Definition: run_mod.f:100
integer nodesk
Definition: compar_mod.f:37
integer nodesi
Definition: compar_mod.f:37
integer dimension_n_s
Definition: param_mod.f:21