MFIX  2016-1
qmomk_bc_mod.f
Go to the documentation of this file.
1 !vvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvvC
2 ! C
3 ! Module name: QMOMK_BC C
4 ! Author: Alberto Passalacqua (A.P.) Date: C
5 ! Reviewer: Date: C
6 ! C
7 !^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^C
8 
9 MODULE qmomk_bc
10 
11  USE param
12  USE param1
13  USE constant
14  USE physprop
15  USE fldvar
16  USE geometry
17  USE compar
18  USE indices
19  USE bc
22  USE functions
23 
24  IMPLICIT NONE
25 
26  PRIVATE
27 
28  PUBLIC :: qmomk_reflective_wall_bc
29  PUBLIC :: qmomk_outlet_bc
30  PUBLIC :: qmomk_inlet_bc
31  PUBLIC :: qmomk_cyclic_bc
32 
33 CONTAINS
34 
35  ! A.P. Purely reflective boundary condition with restitution coefficien e_w
36  SUBROUTINE qmomk_reflective_wall_bc(L, I1, I2, J1, J2, K1, K2, INIT)
37  IMPLICIT NONE
38 
39  INTEGER, INTENT(IN) :: L, I1, I2, J1, J2, K1, K2
40  LOGICAL, INTENT(IN) :: INIT
41 
42  INTEGER :: I, J, K, M, IJK, LFLUID
43 
44  IF (init) THEN
45  DO k = k1, k2
46  DO j = j1, j2
47  DO i = i1, i2
48  !//SP Check if current i,j,k resides on this PE
49  IF (.NOT.is_on_mype_plus2layers(i,j,k)) cycle
50  ijk = funijk(i,j,k)
51  !
52  ! Fluid cell at West
53  !
54  IF (fluid_at(im_of(ijk))) THEN
55  lfluid = im_of(ijk)
56  DO m = 1, mmax
57  qmomk_n0(:,ijk,m) = qmomk_n0(:,lfluid,m)/e_w
58  qmomk_u0(:,ijk,m) = -e_w*qmomk_u0(:,lfluid,m)
59  qmomk_v0(:,ijk,m) = qmomk_v0(:,lfluid,m)
60  qmomk_w0(:,ijk,m) = qmomk_w0(:,lfluid,m)
61  END DO
62  ENDIF
63 
64  !
65  ! Fluid cell at East
66  !
67  IF (fluid_at(ip_of(ijk))) THEN
68  lfluid = ip_of(ijk)
69  DO m = 1, mmax
70  qmomk_n0(:,ijk,m) = qmomk_n0(:,lfluid,m)/e_w
71  qmomk_u0(:,ijk,m) = -e_w*qmomk_u0(:,lfluid,m)
72  qmomk_v0(:,ijk,m) = qmomk_v0(:,lfluid,m)
73  qmomk_w0(:,ijk,m) = qmomk_w0(:,lfluid,m)
74  END DO
75  END IF
76 
77  !
78  ! Fluid cell at South
79  !
80  IF (fluid_at(jm_of(ijk))) THEN
81  lfluid = jm_of(ijk)
82  DO m = 1, mmax
83  qmomk_n0(:,ijk,m) = qmomk_n0(:,lfluid,m)/e_w
84  qmomk_u0(:,ijk,m) = qmomk_u0(:,lfluid,m)
85  qmomk_v0(:,ijk,m) = -e_w*qmomk_v0(:,lfluid,m)
86  qmomk_w0(:,ijk,m) = qmomk_w0(:,lfluid,m)
87  END DO
88  ENDIF
89  !
90  ! Fluid cell at North
91  !
92  IF (fluid_at(jp_of(ijk))) THEN
93  lfluid = jp_of(ijk)
94  DO m = 1, mmax
95  qmomk_n0(:,ijk,m) = qmomk_n0(:,lfluid,m)/e_w
96  qmomk_u0(:,ijk,m) = qmomk_u0(:,lfluid,m)
97  qmomk_v0(:,ijk,m) = -e_w*qmomk_v0(:,lfluid,m)
98  qmomk_w0(:,ijk,m) = qmomk_w0(:,lfluid,m)
99  END DO
100  ENDIF
101  !
102  ! Fluid cell at Bottom
103  !
104  IF (fluid_at(km_of(ijk))) THEN
105  lfluid = km_of(ijk)
106  DO m = 1, mmax
107  qmomk_n0(:,ijk,m) = qmomk_n0(:,lfluid,m)/e_w
108  qmomk_u0(:,ijk,m) = qmomk_u0(:,lfluid,m)
109  qmomk_v0(:,ijk,m) = qmomk_v0(:,lfluid,m)
110  qmomk_w0(:,ijk,m) = -e_w*qmomk_w0(:,lfluid,m)
111  END DO
112  ENDIF
113  !
114  ! Fluid cell at Top
115  !
116  IF (fluid_at(kp_of(ijk))) THEN
117  lfluid = kp_of(ijk)
118  DO m = 1, mmax
119  qmomk_n0(:,ijk,m) = qmomk_n0(:,lfluid,m)/e_w
120  qmomk_u0(:,ijk,m) = qmomk_u0(:,lfluid,m)
121  qmomk_v0(:,ijk,m) = qmomk_v0(:,lfluid,m)
122  qmomk_w0(:,ijk,m) = -e_w*qmomk_w0(:,lfluid,m)
123  END DO
124  ENDIF
125  END DO
126  END DO
127  END DO
128 
129  ! Running...
130  ELSE
131  DO k = k1, k2
132  DO j = j1, j2
133  DO i = i1, i2
134  !//SP Check if current i,j,k resides on this PE
135  IF (.NOT.is_on_mype_plus2layers(i,j,k)) cycle
136  ijk = funijk(i,j,k)
137  !
138  ! Fluid cell at West
139  !
140  IF (fluid_at(im_of(ijk))) THEN
141  lfluid = im_of(ijk)
142  DO m = 1, mmax
143  qmomk_n1(:,ijk,m) = qmomk_n1(:,lfluid,m)/e_w
144  qmomk_u1(:,ijk,m) = -e_w*qmomk_u1(:,lfluid,m)
145  qmomk_v1(:,ijk,m) = qmomk_v1(:,lfluid,m)
146  qmomk_w1(:,ijk,m) = qmomk_w1(:,lfluid,m)
147  END DO
148  ENDIF
149 
150  !
151  ! Fluid cell at East
152  !
153  IF (fluid_at(ip_of(ijk))) THEN
154  lfluid = ip_of(ijk)
155  DO m = 1, mmax
156  qmomk_n1(:,ijk,m) = qmomk_n1(:,lfluid,m)/e_w
157  qmomk_u1(:,ijk,m) = -e_w*qmomk_u1(:,lfluid,m)
158  qmomk_v1(:,ijk,m) = qmomk_v1(:,lfluid,m)
159  qmomk_w1(:,ijk,m) = qmomk_w1(:,lfluid,m)
160  END DO
161  END IF
162 
163  !
164  ! Fluid cell at South
165  !
166  IF (fluid_at(jm_of(ijk))) THEN
167  lfluid = jm_of(ijk)
168  DO m = 1, mmax
169  qmomk_n1(:,ijk,m) = qmomk_n1(:,lfluid,m)/e_w
170  qmomk_u1(:,ijk,m) = qmomk_u1(:,lfluid,m)
171  qmomk_v1(:,ijk,m) = -e_w*qmomk_v1(:,lfluid,m)
172  qmomk_w1(:,ijk,m) = qmomk_w1(:,lfluid,m)
173  END DO
174  ENDIF
175  !
176  ! Fluid cell at North
177  !
178  IF (fluid_at(jp_of(ijk))) THEN
179  lfluid = jp_of(ijk)
180  DO m = 1, mmax
181  qmomk_n1(:,ijk,m) = qmomk_n1(:,lfluid,m)/e_w
182  qmomk_u1(:,ijk,m) = qmomk_u1(:,lfluid,m)
183  qmomk_v1(:,ijk,m) = -e_w*qmomk_v1(:,lfluid,m)
184  qmomk_w1(:,ijk,m) = qmomk_w1(:,lfluid,m)
185  END DO
186  ENDIF
187  !
188  ! Fluid cell at Bottom
189  !
190  IF (fluid_at(km_of(ijk))) THEN
191  lfluid = km_of(ijk)
192  DO m = 1, mmax
193  qmomk_n1(:,ijk,m) = qmomk_n1(:,lfluid,m)/e_w
194  qmomk_u1(:,ijk,m) = qmomk_u1(:,lfluid,m)
195  qmomk_v1(:,ijk,m) = qmomk_v1(:,lfluid,m)
196  qmomk_w1(:,ijk,m) = -e_w*qmomk_w1(:,lfluid,m)
197  END DO
198  ENDIF
199  !
200  ! Fluid cell at Top
201  !
202  IF (fluid_at(kp_of(ijk))) THEN
203  lfluid = kp_of(ijk)
204  DO m = 1, mmax
205  qmomk_n1(:,ijk,m) = qmomk_n1(:,lfluid,m)/e_w
206  qmomk_u1(:,ijk,m) = qmomk_u1(:,lfluid,m)
207  qmomk_v1(:,ijk,m) = qmomk_v1(:,lfluid,m)
208  qmomk_w1(:,ijk,m) = -e_w*qmomk_w1(:,lfluid,m)
209  END DO
210  ENDIF
211  END DO
212  END DO
213  END DO
214  END IF
215  RETURN
216 
217  END SUBROUTINE qmomk_reflective_wall_bc
218 
219  ! A.P. Zero-gradient outlet
220  SUBROUTINE qmomk_outlet_bc(L, I1, I2, J1, J2, K1, K2, INIT)
221  IMPLICIT NONE
222 
223  INTEGER, INTENT(IN) :: I1, I2, J1, J2, K1, K2, L
224  LOGICAL, INTENT(IN) :: INIT
225 
226  INTEGER :: I, J, K, M, IJK, LFLUID
227 
228  IF (init) THEN
229  DO k = k1, k2
230  DO j = j1, j2
231  DO i = i1, i2
232 !//SP Check if current i,j,k resides on this PE
233  IF (.NOT.is_on_mype_plus2layers(i,j,k)) cycle
234  ijk = funijk(i,j,k)
235 !
236 ! Fluid cell at West
237 !
238  IF (fluid_at(im_of(ijk))) THEN
239  lfluid = im_of(ijk)
240  DO m = 1, mmax
241  qmomk_n0(:,ijk,m) = qmomk_n0(:,lfluid,m)
242  qmomk_u0(:,ijk,m) = qmomk_u0(:,lfluid,m)
243  qmomk_v0(:,ijk,m) = qmomk_v0(:,lfluid,m)
244  qmomk_w0(:,ijk,m) = qmomk_w0(:,lfluid,m)
245  END DO
246  ENDIF
247 
248 !
249 ! Fluid cell at East
250 !
251  IF (fluid_at(ip_of(ijk))) THEN
252  lfluid = ip_of(ijk)
253  DO m = 1, mmax
254  qmomk_n0(:,ijk,m) = qmomk_n0(:,lfluid,m)
255  qmomk_u0(:,ijk,m) = qmomk_u0(:,lfluid,m)
256  qmomk_v0(:,ijk,m) = qmomk_v0(:,lfluid,m)
257  qmomk_w0(:,ijk,m) = qmomk_w0(:,lfluid,m)
258  END DO
259  END IF
260 
261 !
262 ! Fluid cell at South
263 !
264  IF (fluid_at(jm_of(ijk))) THEN
265  lfluid = jm_of(ijk)
266  DO m = 1, mmax
267  qmomk_n0(:,ijk,m) = qmomk_n0(:,lfluid,m)
268  qmomk_u0(:,ijk,m) = qmomk_u0(:,lfluid,m)
269  qmomk_v0(:,ijk,m) = qmomk_v0(:,lfluid,m)
270  qmomk_w0(:,ijk,m) = qmomk_w0(:,lfluid,m)
271  END DO
272  ENDIF
273 !
274 ! Fluid cell at North
275 !
276  IF (fluid_at(jp_of(ijk))) THEN
277  lfluid = jp_of(ijk)
278  DO m = 1, mmax
279  qmomk_n0(:,ijk,m) = qmomk_n0(:,lfluid,m)
280  qmomk_u0(:,ijk,m) = qmomk_u0(:,lfluid,m)
281  qmomk_v0(:,ijk,m) = qmomk_v0(:,lfluid,m)
282  qmomk_w0(:,ijk,m) = qmomk_w0(:,lfluid,m)
283  END DO
284  ENDIF
285 !
286 ! Fluid cell at Bottom
287 !
288  IF (fluid_at(km_of(ijk))) THEN
289  lfluid = km_of(ijk)
290  DO m = 1, mmax
291  qmomk_n0(:,ijk,m) = qmomk_n0(:,lfluid,m)
292  qmomk_u0(:,ijk,m) = qmomk_u0(:,lfluid,m)
293  qmomk_v0(:,ijk,m) = qmomk_v0(:,lfluid,m)
294  qmomk_w0(:,ijk,m) = qmomk_w0(:,lfluid,m)
295  END DO
296  ENDIF
297 !
298 ! Fluid cell at Top
299 !
300  IF (fluid_at(kp_of(ijk))) THEN
301  lfluid = kp_of(ijk)
302  DO m = 1, mmax
303  qmomk_n0(:,ijk,m) = qmomk_n0(:,lfluid,m)
304  qmomk_u0(:,ijk,m) = qmomk_u0(:,lfluid,m)
305  qmomk_v0(:,ijk,m) = qmomk_v0(:,lfluid,m)
306  qmomk_w0(:,ijk,m) = qmomk_w0(:,lfluid,m)
307  END DO
308  ENDIF
309  END DO
310  END DO
311  END DO
312 
313  ELSE
314  DO k = k1, k2
315  DO j = j1, j2
316  DO i = i1, i2
317 !//SP Check if current i,j,k resides on this PE
318  IF (.NOT.is_on_mype_plus2layers(i,j,k)) cycle
319  ijk = funijk(i,j,k)
320 !
321 ! Fluid cell at West
322 !
323  IF (fluid_at(im_of(ijk))) THEN
324  lfluid = im_of(ijk)
325  DO m = 1, mmax
326  qmomk_n1(:,ijk,m) = qmomk_n1(:,lfluid,m)
327  qmomk_u1(:,ijk,m) = qmomk_u1(:,lfluid,m)
328  qmomk_v1(:,ijk,m) = qmomk_v1(:,lfluid,m)
329  qmomk_w1(:,ijk,m) = qmomk_w1(:,lfluid,m)
330  END DO
331  ENDIF
332 
333 !
334 ! Fluid cell at East
335 !
336  IF (fluid_at(ip_of(ijk))) THEN
337  lfluid = ip_of(ijk)
338  DO m = 1, mmax
339  qmomk_n1(:,ijk,m) = qmomk_n1(:,lfluid,m)
340  qmomk_u1(:,ijk,m) = qmomk_u1(:,lfluid,m)
341  qmomk_v1(:,ijk,m) = qmomk_v1(:,lfluid,m)
342  qmomk_w1(:,ijk,m) = qmomk_w1(:,lfluid,m)
343  END DO
344  END IF
345 
346 !
347 ! Fluid cell at South
348 !
349  IF (fluid_at(jm_of(ijk))) THEN
350  lfluid = jm_of(ijk)
351  DO m = 1, mmax
352  qmomk_n1(:,ijk,m) = qmomk_n1(:,lfluid,m)
353  qmomk_u1(:,ijk,m) = qmomk_u1(:,lfluid,m)
354  qmomk_v1(:,ijk,m) = qmomk_v1(:,lfluid,m)
355  qmomk_w1(:,ijk,m) = qmomk_w1(:,lfluid,m)
356  END DO
357  ENDIF
358 !
359 ! Fluid cell at North
360 !
361  IF (fluid_at(jp_of(ijk))) THEN
362  lfluid = jp_of(ijk)
363  DO m = 1, mmax
364  qmomk_n1(:,ijk,m) = qmomk_n1(:,lfluid,m)
365  qmomk_u1(:,ijk,m) = qmomk_u1(:,lfluid,m)
366  qmomk_v1(:,ijk,m) = qmomk_v1(:,lfluid,m)
367  qmomk_w1(:,ijk,m) = qmomk_w1(:,lfluid,m)
368  END DO
369  ENDIF
370 !
371 ! Fluid cell at Bottom
372 !
373  IF (fluid_at(km_of(ijk))) THEN
374  lfluid = km_of(ijk)
375  DO m = 1, mmax
376  qmomk_n1(:,ijk,m) = qmomk_n1(:,lfluid,m)
377  qmomk_u1(:,ijk,m) = qmomk_u1(:,lfluid,m)
378  qmomk_v1(:,ijk,m) = qmomk_v1(:,lfluid,m)
379  qmomk_w1(:,ijk,m) = qmomk_w1(:,lfluid,m)
380  END DO
381  ENDIF
382 !
383 ! Fluid cell at Top
384 !
385  IF (fluid_at(kp_of(ijk))) THEN
386  lfluid = kp_of(ijk)
387  DO m = 1, mmax
388  qmomk_n1(:,ijk,m) = qmomk_n1(:,lfluid,m)
389  qmomk_u1(:,ijk,m) = qmomk_u1(:,lfluid,m)
390  qmomk_v1(:,ijk,m) = qmomk_v1(:,lfluid,m)
391  qmomk_w1(:,ijk,m) = qmomk_w1(:,lfluid,m)
392  END DO
393  ENDIF
394  END DO
395  END DO
396  END DO
397  END IF
398 
399  RETURN
400 
401  END SUBROUTINE qmomk_outlet_bc
402 
403 
404  ! A.P. Velocity inlet
405  SUBROUTINE qmomk_inlet_bc(L, INIT)
406  IMPLICIT NONE
407 
408  INTEGER, INTENT(IN) :: L
409  LOGICAL, INTENT(IN) :: INIT
410 
411  INTEGER :: I, J, K, M, IJK, IJK2
412  DOUBLE PRECISION :: InitVal
413 
414  IF (init) THEN
415  DO k = bc_k_b(l), bc_k_t(l)
416  DO j = bc_j_s(l), bc_j_n(l)
417  DO i = bc_i_w(l), bc_i_e(l)
418  IF (.NOT.is_on_mype_plus2layers(i,j,k)) cycle
419  ijk = funijk(i,j,k)
420  SELECT CASE (trim(bc_plane(l)))
421  CASE ('W')
422  ijk2 = im_of(ijk)
423  qmomk_n0(:, ijk2, m) = bc_rop_s(l, m)/(qmomk_nn * ro_s(ijk2,m))
424 
425  initval = max(sqrt(bc_theta_m(l,m)), minimum_theta)
426 
427  qmomk_u0(1, ijk2, m) = -initval + bc_u_s(l, m)
428  qmomk_u0(2, ijk2, m) = +initval + bc_u_s(l, m)
429  qmomk_u0(3, ijk2, m) = -initval + bc_u_s(l, m)
430  qmomk_u0(4, ijk2, m) = +initval + bc_u_s(l, m)
431  qmomk_u0(5, ijk2, m) = -initval + bc_u_s(l, m)
432  qmomk_u0(6, ijk2, m) = +initval + bc_u_s(l, m)
433  qmomk_u0(7, ijk2, m) = -initval + bc_u_s(l, m)
434  qmomk_u0(8, ijk2, m) = +initval + bc_u_s(l, m)
435 
436  qmomk_v0(1, ijk2, m) = -initval
437  qmomk_v0(2, ijk2, m) = -initval
438  qmomk_v0(3, ijk2, m) = +initval
439  qmomk_v0(4, ijk2, m) = +initval
440  qmomk_v0(5, ijk2, m) = -initval
441  qmomk_v0(6, ijk2, m) = -initval
442  qmomk_v0(7, ijk2, m) = +initval
443  qmomk_v0(8, ijk2, m) = +initval
444 
445  qmomk_w0(1, ijk2, m) = -initval
446  qmomk_w0(2, ijk2, m) = -initval
447  qmomk_w0(3, ijk2, m) = -initval
448  qmomk_w0(4, ijk2, m) = -initval
449  qmomk_w0(5, ijk2, m) = +initval
450  qmomk_w0(6, ijk2, m) = +initval
451  qmomk_w0(7, ijk2, m) = +initval
452  qmomk_w0(8, ijk2, m) = +initval
453 
454  CASE ('E')
455  qmomk_n0(:, ijk, m) = bc_rop_s(l, m)/(qmomk_nn * ro_s(ijk,m))
456 
457  initval = max(sqrt(bc_theta_m(l,m)), minimum_theta)
458 
459  qmomk_u0(1, ijk, m) = -initval + bc_u_s(ijk, m)
460  qmomk_u0(2, ijk, m) = +initval + bc_u_s(ijk, m)
461  qmomk_u0(3, ijk, m) = -initval + bc_u_s(ijk, m)
462  qmomk_u0(4, ijk, m) = +initval + bc_u_s(ijk, m)
463  qmomk_u0(5, ijk, m) = -initval + bc_u_s(ijk, m)
464  qmomk_u0(6, ijk, m) = +initval + bc_u_s(ijk, m)
465  qmomk_u0(7, ijk, m) = -initval + bc_u_s(ijk, m)
466  qmomk_u0(8, ijk, m) = +initval + bc_u_s(ijk, m)
467 
468  qmomk_v0(1, ijk, m) = -initval
469  qmomk_v0(2, ijk, m) = -initval
470  qmomk_v0(3, ijk, m) = +initval
471  qmomk_v0(4, ijk, m) = +initval
472  qmomk_v0(5, ijk, m) = -initval
473  qmomk_v0(6, ijk, m) = -initval
474  qmomk_v0(7, ijk, m) = +initval
475  qmomk_v0(8, ijk, m) = +initval
476 
477  qmomk_w0(1, ijk, m) = -initval
478  qmomk_w0(2, ijk, m) = -initval
479  qmomk_w0(3, ijk, m) = -initval
480  qmomk_w0(4, ijk, m) = -initval
481  qmomk_w0(5, ijk, m) = +initval
482  qmomk_w0(6, ijk, m) = +initval
483  qmomk_w0(7, ijk, m) = +initval
484  qmomk_w0(8, ijk, m) = +initval
485  CASE ('S')
486  ijk2 = jm_of(ijk)
487  qmomk_n0(:, ijk2, m) = bc_rop_s(l, m)/(qmomk_nn * ro_s(ijk,m))
488 
489  initval = max(sqrt(bc_theta_m(l,m)), minimum_theta)
490 
491  qmomk_u0(1, ijk2, m) = -initval
492  qmomk_u0(2, ijk2, m) = +initval
493  qmomk_u0(3, ijk2, m) = -initval
494  qmomk_u0(4, ijk2, m) = +initval
495  qmomk_u0(5, ijk2, m) = -initval
496  qmomk_u0(6, ijk2, m) = +initval
497  qmomk_u0(7, ijk2, m) = -initval
498  qmomk_u0(8, ijk2, m) = +initval
499 
500  qmomk_v0(1, ijk2, m) = -initval + bc_v_s(l, m)
501  qmomk_v0(2, ijk2, m) = -initval + bc_v_s(l, m)
502  qmomk_v0(3, ijk2, m) = +initval + bc_v_s(l, m)
503  qmomk_v0(4, ijk2, m) = +initval + bc_v_s(l, m)
504  qmomk_v0(5, ijk2, m) = -initval + bc_v_s(l, m)
505  qmomk_v0(6, ijk2, m) = -initval + bc_v_s(l, m)
506  qmomk_v0(7, ijk2, m) = +initval + bc_v_s(l, m)
507  qmomk_v0(8, ijk2, m) = +initval + bc_v_s(l, m)
508 
509  qmomk_w0(1, ijk2, m) = -initval
510  qmomk_w0(2, ijk2, m) = -initval
511  qmomk_w0(3, ijk2, m) = -initval
512  qmomk_w0(4, ijk2, m) = -initval
513  qmomk_w0(5, ijk2, m) = +initval
514  qmomk_w0(6, ijk2, m) = +initval
515  qmomk_w0(7, ijk2, m) = +initval
516  qmomk_w0(8, ijk2, m) = +initval
517 
518  CASE ('N')
519  qmomk_n0(:, ijk, m) = bc_rop_s(l, m)/(qmomk_nn * ro_s(ijk,m))
520 
521  initval = max(sqrt(bc_theta_m(l,m)), minimum_theta)
522 
523  qmomk_u0(1, ijk, m) = -initval
524  qmomk_u0(2, ijk, m) = +initval
525  qmomk_u0(3, ijk, m) = -initval
526  qmomk_u0(4, ijk, m) = +initval
527  qmomk_u0(5, ijk, m) = -initval
528  qmomk_u0(6, ijk, m) = +initval
529  qmomk_u0(7, ijk, m) = -initval
530  qmomk_u0(8, ijk, m) = +initval
531 
532  qmomk_v0(1, ijk, m) = -initval + bc_v_s(ijk, m)
533  qmomk_v0(2, ijk, m) = -initval + bc_v_s(ijk, m)
534  qmomk_v0(3, ijk, m) = +initval + bc_v_s(ijk, m)
535  qmomk_v0(4, ijk, m) = +initval + bc_v_s(ijk, m)
536  qmomk_v0(5, ijk, m) = -initval + bc_v_s(ijk, m)
537  qmomk_v0(6, ijk, m) = -initval + bc_v_s(ijk, m)
538  qmomk_v0(7, ijk, m) = +initval + bc_v_s(ijk, m)
539  qmomk_v0(8, ijk, m) = +initval + bc_v_s(ijk, m)
540 
541  qmomk_w0(1, ijk, m) = -initval
542  qmomk_w0(2, ijk, m) = -initval
543  qmomk_w0(3, ijk, m) = -initval
544  qmomk_w0(4, ijk, m) = -initval
545  qmomk_w0(5, ijk, m) = +initval
546  qmomk_w0(6, ijk, m) = +initval
547  qmomk_w0(7, ijk, m) = +initval
548  qmomk_w0(8, ijk, m) = +initval
549 
550  CASE ('B')
551  ijk2 = km_of(ijk)
552  qmomk_n0(:, ijk2, m) = bc_rop_s(l, m)/(qmomk_nn * ro_s(ijk2,m))
553 
554  initval = max(sqrt(bc_theta_m(l,m)), minimum_theta)
555 
556  qmomk_u0(1, ijk2, m) = -initval
557  qmomk_u0(2, ijk2, m) = +initval
558  qmomk_u0(3, ijk2, m) = -initval
559  qmomk_u0(4, ijk2, m) = +initval
560  qmomk_u0(5, ijk2, m) = -initval
561  qmomk_u0(6, ijk2, m) = +initval
562  qmomk_u0(7, ijk2, m) = -initval
563  qmomk_u0(8, ijk2, m) = +initval
564 
565  qmomk_v0(1, ijk2, m) = -initval
566  qmomk_v0(2, ijk2, m) = -initval
567  qmomk_v0(3, ijk2, m) = +initval
568  qmomk_v0(4, ijk2, m) = +initval
569  qmomk_v0(5, ijk2, m) = -initval
570  qmomk_v0(6, ijk2, m) = -initval
571  qmomk_v0(7, ijk2, m) = +initval
572  qmomk_v0(8, ijk2, m) = +initval
573 
574  qmomk_w0(1, ijk2, m) = -initval + bc_w_s(l, m)
575  qmomk_w0(2, ijk2, m) = -initval + bc_w_s(l, m)
576  qmomk_w0(3, ijk2, m) = -initval + bc_w_s(l, m)
577  qmomk_w0(4, ijk2, m) = -initval + bc_w_s(l, m)
578  qmomk_w0(5, ijk2, m) = +initval + bc_w_s(l, m)
579  qmomk_w0(6, ijk2, m) = +initval + bc_w_s(l, m)
580  qmomk_w0(7, ijk2, m) = +initval + bc_w_s(l, m)
581  qmomk_w0(8, ijk2, m) = +initval + bc_w_s(l, m)
582 
583  CASE ('T')
584  qmomk_n0(:, ijk, m) = bc_rop_s(l, m)/(qmomk_nn * ro_s(ijk,m))
585 
586  initval = max(sqrt(bc_theta_m(l,m)), minimum_theta)
587 
588  qmomk_u0(1, ijk, m) = -initval
589  qmomk_u0(2, ijk, m) = +initval
590  qmomk_u0(3, ijk, m) = -initval
591  qmomk_u0(4, ijk, m) = +initval
592  qmomk_u0(5, ijk, m) = -initval
593  qmomk_u0(6, ijk, m) = +initval
594  qmomk_u0(7, ijk, m) = -initval
595  qmomk_u0(8, ijk, m) = +initval
596 
597  qmomk_v0(1, ijk, m) = -initval
598  qmomk_v0(2, ijk, m) = -initval
599  qmomk_v0(3, ijk, m) = +initval
600  qmomk_v0(4, ijk, m) = +initval
601  qmomk_v0(5, ijk, m) = -initval
602  qmomk_v0(6, ijk, m) = -initval
603  qmomk_v0(7, ijk, m) = +initval
604  qmomk_v0(8, ijk, m) = +initval
605 
606  qmomk_w0(1, ijk, m) = -initval + bc_w_s(l, m)
607  qmomk_w0(2, ijk, m) = -initval + bc_w_s(l, m)
608  qmomk_w0(3, ijk, m) = -initval + bc_w_s(l, m)
609  qmomk_w0(4, ijk, m) = -initval + bc_w_s(l, m)
610  qmomk_w0(5, ijk, m) = +initval + bc_w_s(l, m)
611  qmomk_w0(6, ijk, m) = +initval + bc_w_s(l, m)
612  qmomk_w0(7, ijk, m) = +initval + bc_w_s(l, m)
613  qmomk_w0(8, ijk, m) = +initval + bc_w_s(l, m)
614  END SELECT
615  END DO
616  END DO
617  END DO
618 
619  ELSE
620  DO k = bc_k_b(l), bc_k_t(l)
621  DO j = bc_j_s(l), bc_j_n(l)
622  DO i = bc_i_w(l), bc_i_e(l)
623  IF (.NOT.is_on_mype_plus2layers(i,j,k)) cycle
624  ijk = funijk(i,j,k)
625  SELECT CASE (trim(bc_plane(l)))
626 
627  CASE ('W')
628  ijk2 = im_of(ijk)
629  qmomk_n1(:, ijk2, m) = bc_rop_s(l, m)/(qmomk_nn * ro_s(ijk2,m))
630 
631  initval = max(sqrt(bc_theta_m(l,m)), minimum_theta)
632 
633  qmomk_u1(1, ijk2, m) = -initval + bc_u_s(l, m)
634  qmomk_u1(2, ijk2, m) = +initval + bc_u_s(l, m)
635  qmomk_u1(3, ijk2, m) = -initval + bc_u_s(l, m)
636  qmomk_u1(4, ijk2, m) = +initval + bc_u_s(l, m)
637  qmomk_u1(5, ijk2, m) = -initval + bc_u_s(l, m)
638  qmomk_u1(6, ijk2, m) = +initval + bc_u_s(l, m)
639  qmomk_u1(7, ijk2, m) = -initval + bc_u_s(l, m)
640  qmomk_u1(8, ijk2, m) = +initval + bc_u_s(l, m)
641 
642  qmomk_v1(1, ijk2, m) = -initval
643  qmomk_v1(2, ijk2, m) = -initval
644  qmomk_v1(3, ijk2, m) = +initval
645  qmomk_v1(4, ijk2, m) = +initval
646  qmomk_v1(5, ijk2, m) = -initval
647  qmomk_v1(6, ijk2, m) = -initval
648  qmomk_v1(7, ijk2, m) = +initval
649  qmomk_v1(8, ijk2, m) = +initval
650 
651  qmomk_w1(1, ijk2, m) = -initval
652  qmomk_w1(2, ijk2, m) = -initval
653  qmomk_w1(3, ijk2, m) = -initval
654  qmomk_w1(4, ijk2, m) = -initval
655  qmomk_w1(5, ijk2, m) = +initval
656  qmomk_w1(6, ijk2, m) = +initval
657  qmomk_w1(7, ijk2, m) = +initval
658  qmomk_w1(8, ijk2, m) = +initval
659 
660  CASE ('E')
661  qmomk_n1(:, ijk, m) = bc_rop_s(l, m)/(qmomk_nn * ro_s(ijk,m))
662 
663  initval = max(sqrt(bc_theta_m(l,m)), minimum_theta)
664 
665  qmomk_u1(1, ijk, m) = -initval + bc_u_s(ijk, m)
666  qmomk_u1(2, ijk, m) = +initval + bc_u_s(ijk, m)
667  qmomk_u1(3, ijk, m) = -initval + bc_u_s(ijk, m)
668  qmomk_u1(4, ijk, m) = +initval + bc_u_s(ijk, m)
669  qmomk_u1(5, ijk, m) = -initval + bc_u_s(ijk, m)
670  qmomk_u1(6, ijk, m) = +initval + bc_u_s(ijk, m)
671  qmomk_u1(7, ijk, m) = -initval + bc_u_s(ijk, m)
672  qmomk_u1(8, ijk, m) = +initval + bc_u_s(ijk, m)
673 
674  qmomk_v1(1, ijk, m) = -initval
675  qmomk_v1(2, ijk, m) = -initval
676  qmomk_v1(3, ijk, m) = +initval
677  qmomk_v1(4, ijk, m) = +initval
678  qmomk_v1(5, ijk, m) = -initval
679  qmomk_v1(6, ijk, m) = -initval
680  qmomk_v1(7, ijk, m) = +initval
681  qmomk_v1(8, ijk, m) = +initval
682 
683  qmomk_w1(1, ijk, m) = -initval
684  qmomk_w1(2, ijk, m) = -initval
685  qmomk_w1(3, ijk, m) = -initval
686  qmomk_w1(4, ijk, m) = -initval
687  qmomk_w1(5, ijk, m) = +initval
688  qmomk_w1(6, ijk, m) = +initval
689  qmomk_w1(7, ijk, m) = +initval
690  qmomk_w1(8, ijk, m) = +initval
691 
692  CASE ('S')
693  ijk2 = jm_of(ijk)
694  qmomk_n1(:, ijk2, m) = bc_rop_s(l, m)/(qmomk_nn * ro_s(ijk2,m))
695 
696  initval = max(sqrt(bc_theta_m(l,m)), minimum_theta)
697 
698  qmomk_u1(1, ijk2, m) = -initval
699  qmomk_u1(2, ijk2, m) = +initval
700  qmomk_u1(3, ijk2, m) = -initval
701  qmomk_u1(4, ijk2, m) = +initval
702  qmomk_u1(5, ijk2, m) = -initval
703  qmomk_u1(6, ijk2, m) = +initval
704  qmomk_u1(7, ijk2, m) = -initval
705  qmomk_u1(8, ijk2, m) = +initval
706 
707  qmomk_v1(1, ijk2, m) = -initval + bc_v_s(l, m)
708  qmomk_v1(2, ijk2, m) = -initval + bc_v_s(l, m)
709  qmomk_v1(3, ijk2, m) = +initval + bc_v_s(l, m)
710  qmomk_v1(4, ijk2, m) = +initval + bc_v_s(l, m)
711  qmomk_v1(5, ijk2, m) = -initval + bc_v_s(l, m)
712  qmomk_v1(6, ijk2, m) = -initval + bc_v_s(l, m)
713  qmomk_v1(7, ijk2, m) = +initval + bc_v_s(l, m)
714  qmomk_v1(8, ijk2, m) = +initval + bc_v_s(l, m)
715 
716  qmomk_w1(1, ijk2, m) = -initval
717  qmomk_w1(2, ijk2, m) = -initval
718  qmomk_w1(3, ijk2, m) = -initval
719  qmomk_w1(4, ijk2, m) = -initval
720  qmomk_w1(5, ijk2, m) = +initval
721  qmomk_w1(6, ijk2, m) = +initval
722  qmomk_w1(7, ijk2, m) = +initval
723  qmomk_w1(8, ijk2, m) = +initval
724 
725  CASE ('N')
726  qmomk_n1(:, ijk, m) = bc_rop_s(l, m)/(qmomk_nn * ro_s(ijk,m))
727 
728  initval = max(sqrt(bc_theta_m(l,m)), minimum_theta)
729 
730  qmomk_u1(1, ijk, m) = -initval
731  qmomk_u1(2, ijk, m) = +initval
732  qmomk_u1(3, ijk, m) = -initval
733  qmomk_u1(4, ijk, m) = +initval
734  qmomk_u1(5, ijk, m) = -initval
735  qmomk_u1(6, ijk, m) = +initval
736  qmomk_u1(7, ijk, m) = -initval
737  qmomk_u1(8, ijk, m) = +initval
738 
739  qmomk_v1(1, ijk, m) = -initval + bc_v_s(ijk, m)
740  qmomk_v1(2, ijk, m) = -initval + bc_v_s(ijk, m)
741  qmomk_v1(3, ijk, m) = +initval + bc_v_s(ijk, m)
742  qmomk_v1(4, ijk, m) = +initval + bc_v_s(ijk, m)
743  qmomk_v1(5, ijk, m) = -initval + bc_v_s(ijk, m)
744  qmomk_v1(6, ijk, m) = -initval + bc_v_s(ijk, m)
745  qmomk_v1(7, ijk, m) = +initval + bc_v_s(ijk, m)
746  qmomk_v1(8, ijk, m) = +initval + bc_v_s(ijk, m)
747 
748  qmomk_w1(1, ijk, m) = -initval
749  qmomk_w1(2, ijk, m) = -initval
750  qmomk_w1(3, ijk, m) = -initval
751  qmomk_w1(4, ijk, m) = -initval
752  qmomk_w1(5, ijk, m) = +initval
753  qmomk_w1(6, ijk, m) = +initval
754  qmomk_w1(7, ijk, m) = +initval
755  qmomk_w1(8, ijk, m) = +initval
756 
757  CASE ('B')
758  ijk2 = km_of(ijk)
759  qmomk_n1(:, ijk2, m) = bc_rop_s(l, m)/(qmomk_nn * ro_s(ijk2,m))
760 
761  initval = max(sqrt(bc_theta_m(l,m)), minimum_theta)
762 
763  qmomk_u1(1, ijk2, m) = -initval
764  qmomk_u1(2, ijk2, m) = +initval
765  qmomk_u1(3, ijk2, m) = -initval
766  qmomk_u1(4, ijk2, m) = +initval
767  qmomk_u1(5, ijk2, m) = -initval
768  qmomk_u1(6, ijk2, m) = +initval
769  qmomk_u1(7, ijk2, m) = -initval
770  qmomk_u1(8, ijk2, m) = +initval
771 
772  qmomk_v1(1, ijk2, m) = -initval
773  qmomk_v1(2, ijk2, m) = -initval
774  qmomk_v1(3, ijk2, m) = +initval
775  qmomk_v1(4, ijk2, m) = +initval
776  qmomk_v1(5, ijk2, m) = -initval
777  qmomk_v1(6, ijk2, m) = -initval
778  qmomk_v1(7, ijk2, m) = +initval
779  qmomk_v1(8, ijk2, m) = +initval
780 
781  qmomk_w1(1, ijk2, m) = -initval + bc_w_s(l, m)
782  qmomk_w1(2, ijk2, m) = -initval + bc_w_s(l, m)
783  qmomk_w1(3, ijk2, m) = -initval + bc_w_s(l, m)
784  qmomk_w1(4, ijk2, m) = -initval + bc_w_s(l, m)
785  qmomk_w1(5, ijk2, m) = +initval + bc_w_s(l, m)
786  qmomk_w1(6, ijk2, m) = +initval + bc_w_s(l, m)
787  qmomk_w1(7, ijk2, m) = +initval + bc_w_s(l, m)
788  qmomk_w1(8, ijk2, m) = +initval + bc_w_s(l, m)
789 
790  CASE ('T')
791  qmomk_n1(:, ijk, m) = bc_rop_s(l, m)/(qmomk_nn * ro_s(ijk,m))
792 
793  initval = max(sqrt(bc_theta_m(l,m)), minimum_theta)
794 
795  qmomk_u1(1, ijk, m) = -initval
796  qmomk_u1(2, ijk, m) = +initval
797  qmomk_u1(3, ijk, m) = -initval
798  qmomk_u1(4, ijk, m) = +initval
799  qmomk_u1(5, ijk, m) = -initval
800  qmomk_u1(6, ijk, m) = +initval
801  qmomk_u1(7, ijk, m) = -initval
802  qmomk_u1(8, ijk, m) = +initval
803 
804  qmomk_v1(1, ijk, m) = -initval
805  qmomk_v1(2, ijk, m) = -initval
806  qmomk_v1(3, ijk, m) = +initval
807  qmomk_v1(4, ijk, m) = +initval
808  qmomk_v1(5, ijk, m) = -initval
809  qmomk_v1(6, ijk, m) = -initval
810  qmomk_v1(7, ijk, m) = +initval
811  qmomk_v1(8, ijk, m) = +initval
812 
813  qmomk_w1(1, ijk, m) = -initval + bc_w_s(l, m)
814  qmomk_w1(2, ijk, m) = -initval + bc_w_s(l, m)
815  qmomk_w1(3, ijk, m) = -initval + bc_w_s(l, m)
816  qmomk_w1(4, ijk, m) = -initval + bc_w_s(l, m)
817  qmomk_w1(5, ijk, m) = +initval + bc_w_s(l, m)
818  qmomk_w1(6, ijk, m) = +initval + bc_w_s(l, m)
819  qmomk_w1(7, ijk, m) = +initval + bc_w_s(l, m)
820  qmomk_w1(8, ijk, m) = +initval + bc_w_s(l, m)
821  END SELECT
822  END DO
823  END DO
824  END DO
825  END IF
826  RETURN
827  END SUBROUTINE qmomk_inlet_bc
828 
829  ! A.P. Cyclic boundary conditions
830  SUBROUTINE qmomk_cyclic_bc(INIT)
831  IMPLICIT NONE
832 
833  LOGICAL, INTENT(IN) :: INIT
834  INTEGER :: IJK, IJK_CYCLIC, I, J, K, IJKN, IJKS, IJKE, IJKW, IJKT, IJKB, M
835  DOUBLE PRECISION, DIMENSION(QMOMK_NN) :: QMOMK_N_TMP, QMOMK_U_TMP
836  DOUBLE PRECISION, DIMENSION(QMOMK_NN) :: QMOMK_V_TMP, QMOMK_W_TMP
837 
838  IF (init) THEN
839  DO m = 1, mmax
840  DO ijk = ijkstart3, ijkend3
841  IF (fluid_at(ijk)) THEN
842  ijkn = north_of(ijk)
843  ijks = south_of(ijk)
844  ijke = east_of(ijk)
845  ijkw = west_of(ijk)
846  ijkt = top_of(ijk)
847  ijkb = bottom_of(ijk)
848 
849  ! x direction cyclic
850  IF (cyclic_x .OR. cyclic_x_pd) THEN
851  IF (cyclic_at_e(ijk)) THEN
852  i = i_of(ijke)
853  j = j_of(ijke)
854  k = k_of(ijke)
855 
856  ijk_cyclic = funijk(ip1(i), j, k)
857 
858  qmomk_n_tmp(:) = qmomk_n0(:,ijke,m)
859  qmomk_u_tmp(:) = qmomk_u0(:,ijke,m)
860  qmomk_v_tmp(:) = qmomk_v0(:,ijke,m)
861  qmomk_w_tmp(:) = qmomk_w0(:,ijke,m)
862 
863  qmomk_n0(:,ijke,m) = qmomk_n0(:,ijk_cyclic,m)
864  qmomk_u0(:,ijke,m) = qmomk_u0(:,ijk_cyclic,m)
865  qmomk_v0(:,ijke,m) = qmomk_v0(:,ijk_cyclic,m)
866  qmomk_w0(:,ijke,m) = qmomk_w0(:,ijk_cyclic,m)
867 
868  qmomk_n0(:,ijk_cyclic,m) = qmomk_n_tmp(:)
869  qmomk_u0(:,ijk_cyclic,m) = qmomk_u_tmp(:)
870  qmomk_v0(:,ijk_cyclic,m) = qmomk_v_tmp(:)
871  qmomk_w0(:,ijk_cyclic,m) = qmomk_w_tmp(:)
872  END IF
873  END IF
874  ! y direction cyclic
875  IF (cyclic_y .OR. cyclic_y_pd) THEN
876  IF (cyclic_at_n(ijk)) THEN
877  i = i_of(ijkn)
878  j = j_of(ijkn)
879  k = k_of(ijkn)
880 
881  ijk_cyclic = funijk(i, jp1(j), k)
882 
883  qmomk_n_tmp(:) = qmomk_n0(:,ijkn,m)
884  qmomk_u_tmp(:) = qmomk_u0(:,ijkn,m)
885  qmomk_v_tmp(:) = qmomk_v0(:,ijkn,m)
886  qmomk_w_tmp(:) = qmomk_w0(:,ijkn,m)
887 
888  qmomk_n0(:,ijkn,m) = qmomk_n0(:,ijk_cyclic,m)
889  qmomk_u0(:,ijkn,m) = qmomk_u0(:,ijk_cyclic,m)
890  qmomk_v0(:,ijkn,m) = qmomk_v0(:,ijk_cyclic,m)
891  qmomk_w0(:,ijkn,m) = qmomk_w0(:,ijk_cyclic,m)
892 
893  qmomk_n0(:,ijk_cyclic,m) = qmomk_n_tmp(:)
894  qmomk_u0(:,ijk_cyclic,m) = qmomk_u_tmp(:)
895  qmomk_v0(:,ijk_cyclic,m) = qmomk_v_tmp(:)
896  qmomk_w0(:,ijk_cyclic,m) = qmomk_w_tmp(:)
897  END IF
898  END IF
899  ! z direction cyclic
900  IF (cyclic_z .OR. cyclic_z_pd) THEN
901  IF (cyclic_at_t(ijk)) THEN
902  i = i_of(ijkt)
903  j = j_of(ijkt)
904  k = k_of(ijkt)
905 
906  ijk_cyclic = funijk(i, j, kp1(k))
907 
908  qmomk_n_tmp(:) = qmomk_n0(:,ijkt,m)
909  qmomk_u_tmp(:) = qmomk_u0(:,ijkt,m)
910  qmomk_v_tmp(:) = qmomk_v0(:,ijkt,m)
911  qmomk_w_tmp(:) = qmomk_w0(:,ijkt,m)
912 
913  qmomk_n0(:,ijkt,m) = qmomk_n0(:,ijk_cyclic,m)
914  qmomk_u0(:,ijkt,m) = qmomk_u0(:,ijk_cyclic,m)
915  qmomk_v0(:,ijkt,m) = qmomk_v0(:,ijk_cyclic,m)
916  qmomk_w0(:,ijkt,m) = qmomk_w0(:,ijk_cyclic,m)
917 
918  qmomk_n0(:,ijk_cyclic,m) = qmomk_n_tmp(:)
919  qmomk_u0(:,ijk_cyclic,m) = qmomk_u_tmp(:)
920  qmomk_v0(:,ijk_cyclic,m) = qmomk_v_tmp(:)
921  qmomk_w0(:,ijk_cyclic,m) = qmomk_w_tmp(:)
922  END IF
923  END IF
924  END IF
925  END DO
926  END DO
927  ELSE
928  DO m = 1, mmax
929  DO ijk = ijkstart3, ijkend3
930  IF (fluid_at(ijk)) THEN
931  ijkn = north_of(ijk)
932  ijks = south_of(ijk)
933  ijke = east_of(ijk)
934  ijkw = west_of(ijk)
935  ijkt = top_of(ijk)
936  ijkb = bottom_of(ijk)
937 
938  ! x direction cyclic
939  IF (cyclic_x .OR. cyclic_x_pd) THEN
940  IF (cyclic_at_e(ijk)) THEN
941  i = i_of(ijke)
942  j = j_of(ijke)
943  k = k_of(ijke)
944 
945  ijk_cyclic = funijk(ip1(i), j, k)
946 
947  qmomk_n_tmp(:) = qmomk_n1(:,ijke,m)
948  qmomk_u_tmp(:) = qmomk_u1(:,ijke,m)
949  qmomk_v_tmp(:) = qmomk_v1(:,ijke,m)
950  qmomk_w_tmp(:) = qmomk_w1(:,ijke,m)
951 
952  qmomk_n1(:,ijke,m) = qmomk_n1(:,ijk_cyclic,m)
953  qmomk_u1(:,ijke,m) = qmomk_u1(:,ijk_cyclic,m)
954  qmomk_v1(:,ijke,m) = qmomk_v1(:,ijk_cyclic,m)
955  qmomk_w1(:,ijke,m) = qmomk_w1(:,ijk_cyclic,m)
956 
957  qmomk_n1(:,ijk_cyclic,m) = qmomk_n_tmp(:)
958  qmomk_u1(:,ijk_cyclic,m) = qmomk_u_tmp(:)
959  qmomk_v1(:,ijk_cyclic,m) = qmomk_v_tmp(:)
960  qmomk_w1(:,ijk_cyclic,m) = qmomk_w_tmp(:)
961  END IF
962  END IF
963  ! y direction cyclic
964  IF (cyclic_y .OR. cyclic_y_pd) THEN
965  IF (cyclic_at_n(ijk)) THEN
966  i = i_of(ijkn)
967  j = j_of(ijkn)
968  k = k_of(ijkn)
969 
970  ijk_cyclic = funijk(i, jp1(j), k)
971 
972  qmomk_n_tmp(:) = qmomk_n1(:,ijkn,m)
973  qmomk_u_tmp(:) = qmomk_u1(:,ijkn,m)
974  qmomk_v_tmp(:) = qmomk_v1(:,ijkn,m)
975  qmomk_w_tmp(:) = qmomk_w1(:,ijkn,m)
976 
977  qmomk_n1(:,ijkn,m) = qmomk_n1(:,ijk_cyclic,m)
978  qmomk_u1(:,ijkn,m) = qmomk_u1(:,ijk_cyclic,m)
979  qmomk_v1(:,ijkn,m) = qmomk_v1(:,ijk_cyclic,m)
980  qmomk_w1(:,ijkn,m) = qmomk_w1(:,ijk_cyclic,m)
981 
982  qmomk_n1(:,ijk_cyclic,m) = qmomk_n_tmp(:)
983  qmomk_u1(:,ijk_cyclic,m) = qmomk_u_tmp(:)
984  qmomk_v1(:,ijk_cyclic,m) = qmomk_v_tmp(:)
985  qmomk_w1(:,ijk_cyclic,m) = qmomk_w_tmp(:)
986  END IF
987  END IF
988  ! z direction cyclic
989  IF (cyclic_z .OR. cyclic_z_pd) THEN
990  IF (cyclic_at_t(ijk)) THEN
991  i = i_of(ijkt)
992  j = j_of(ijkt)
993  k = k_of(ijkt)
994 
995  ijk_cyclic = funijk(i, j, kp1(k))
996 
997  qmomk_n_tmp(:) = qmomk_n1(:,ijkt,m)
998  qmomk_u_tmp(:) = qmomk_u1(:,ijkt,m)
999  qmomk_v_tmp(:) = qmomk_v1(:,ijkt,m)
1000  qmomk_w_tmp(:) = qmomk_w1(:,ijkt,m)
1001 
1002  qmomk_n1(:,ijkt,m) = qmomk_n1(:,ijk_cyclic,m)
1003  qmomk_u1(:,ijkt,m) = qmomk_u1(:,ijk_cyclic,m)
1004  qmomk_v1(:,ijkt,m) = qmomk_v1(:,ijk_cyclic,m)
1005  qmomk_w1(:,ijkt,m) = qmomk_w1(:,ijk_cyclic,m)
1006 
1007  qmomk_n1(:,ijk_cyclic,m) = qmomk_n_tmp(:)
1008  qmomk_u1(:,ijk_cyclic,m) = qmomk_u_tmp(:)
1009  qmomk_v1(:,ijk_cyclic,m) = qmomk_v_tmp(:)
1010  qmomk_w1(:,ijk_cyclic,m) = qmomk_w_tmp(:)
1011  END IF
1012  END IF
1013  END IF
1014  END DO
1015  END DO
1016  END IF
1017  END SUBROUTINE qmomk_cyclic_bc
1018 
1019 END MODULE qmomk_bc
integer, dimension(:), allocatable ip1
Definition: indices_mod.f:50
integer, dimension(dimension_bc) bc_k_b
Definition: bc_mod.f:70
double precision e_w
Definition: constant_mod.f:85
integer, dimension(:), allocatable i_of
Definition: indices_mod.f:45
integer ijkend3
Definition: compar_mod.f:80
integer, dimension(dimension_bc) bc_i_w
Definition: bc_mod.f:54
integer, dimension(dimension_bc) bc_j_n
Definition: bc_mod.f:66
double precision, dimension(:,:,:), allocatable qmomk_w1
double precision, dimension(dimension_bc, dim_m) bc_w_s
Definition: bc_mod.f:129
double precision, dimension(:,:,:), allocatable qmomk_v0
double precision, dimension(:,:,:), allocatable qmomk_w0
logical cyclic_z
Definition: geometry_mod.f:153
logical cyclic_z_pd
Definition: geometry_mod.f:159
character, dimension(dimension_bc) bc_plane
Definition: bc_mod.f:217
integer, dimension(:), allocatable k_of
Definition: indices_mod.f:47
subroutine, public qmomk_outlet_bc(L, I1, I2, J1, J2, K1, K2, INIT)
Definition: qmomk_bc_mod.f:221
integer, dimension(dimension_bc) bc_k_t
Definition: bc_mod.f:74
integer mmax
Definition: physprop_mod.f:19
double precision, dimension(:,:,:), allocatable qmomk_n1
subroutine, public qmomk_cyclic_bc(INIT)
Definition: qmomk_bc_mod.f:831
integer, dimension(:), allocatable j_of
Definition: indices_mod.f:46
logical cyclic_y_pd
Definition: geometry_mod.f:157
integer, dimension(dimension_bc) bc_j_s
Definition: bc_mod.f:62
double precision, dimension(:,:,:), allocatable qmomk_n0
logical cyclic_y
Definition: geometry_mod.f:151
integer, dimension(:), allocatable jp1
Definition: indices_mod.f:51
double precision, dimension(:,:,:), allocatable qmomk_u1
integer, dimension(:), allocatable kp1
Definition: indices_mod.f:52
logical cyclic_x
Definition: geometry_mod.f:149
Definition: param_mod.f:2
double precision, dimension(:,:), allocatable ro_s
Definition: fldvar_mod.f:45
double precision, dimension(dimension_bc, dim_m) bc_v_s
Definition: bc_mod.f:121
subroutine, public qmomk_reflective_wall_bc(L, I1, I2, J1, J2, K1, K2, INIT)
Definition: qmomk_bc_mod.f:37
logical cyclic_x_pd
Definition: geometry_mod.f:155
integer ijkstart3
Definition: compar_mod.f:80
double precision, dimension(dimension_bc, dim_m) bc_u_s
Definition: bc_mod.f:113
double precision, dimension(dimension_bc, dim_m) bc_theta_m
Definition: bc_mod.f:105
subroutine, public qmomk_inlet_bc(L, INIT)
Definition: qmomk_bc_mod.f:406
double precision, dimension(:,:,:), allocatable qmomk_v1
integer, dimension(dimension_bc) bc_i_e
Definition: bc_mod.f:58
double precision, dimension(dimension_bc, dim_m) bc_rop_s
Definition: bc_mod.f:92
Definition: bc_mod.f:23
double precision, dimension(:,:,:), allocatable qmomk_u0