Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
TURB_VFluxHybrid.F90
Go to the documentation of this file.
1 ! *********************************************************************
2 ! * Rocstar Simulation Suite *
3 ! * Copyright@2015, Illinois Rocstar LLC. All rights reserved. *
4 ! * *
5 ! * Illinois Rocstar LLC *
6 ! * Champaign, IL *
7 ! * www.illinoisrocstar.com *
8 ! * sales@illinoisrocstar.com *
9 ! * *
10 ! * License: See LICENSE file in top level of distribution package or *
11 ! * http://opensource.org/licenses/NCSA *
12 ! *********************************************************************
13 ! *********************************************************************
14 ! * THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, *
15 ! * EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES *
16 ! * OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND *
17 ! * NONINFRINGEMENT. IN NO EVENT SHALL THE CONTRIBUTORS OR *
18 ! * COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER *
19 ! * LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, *
20 ! * Arising FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE *
21 ! * USE OR OTHER DEALINGS WITH THE SOFTWARE. *
22 ! *********************************************************************
23 !******************************************************************************
24 !
25 ! Purpose: Compute total viscous fluxes of which turbulent stress contribution
26 ! is based on stress model as opposed to eddy viscosity model.
27 !
28 ! Description: The fluxes are computed in three consecutive sweeps (i, j and k
29 ! direction) at interior faces. The fluxes at the boundary patches
30 ! are then treated. The energy model contributions are also
31 ! treated on the fly.
32 !
33 ! Input: region = data of current region
34 !
35 ! Output: mixt%diss = viscous fluxes added to the dissipation.
36 !
37 ! Notes: This flux routine is used for the scale similarity and dyn.mixed model.
38 !
39 !******************************************************************************
40 !
41 ! $Id: TURB_VFluxHybrid.F90,v 1.10 2008/12/06 08:44:42 mtcampbe Exp $
42 !
43 ! Copyright: (c) 2001 by the University of Illinois
44 !
45 !******************************************************************************
46 
47 SUBROUTINE turb_vfluxhybrid( region )
48 
49 !$startcopy TURB_VisFluxEddy
50  USE moddatatypes
51  USE modbndpatch, ONLY : t_patch
52  USE moddatastruct, ONLY : t_region
53  USE modglobal, ONLY : t_global
54  USE modturbulence, ONLY : t_turb
55 !$endcopy
59 !$startcopy TURB_VisFluxEddy
60 #ifdef RFLO
61  USE modinterfaces, ONLY : rflo_getdimensphys, &
63 
64 #include "Indexing.h"
65 #endif
66  USE moderror
67  USE modparameters
69  IMPLICIT NONE
70 
71 ! ... parameters
72  TYPE(t_region), TARGET :: region
73 
74 ! ... loop variables
75  INTEGER :: i, j, k, ic, ipatch
76 
77 ! ... local variables
78  CHARACTER(CHRLEN) :: rcsidentstring
79  TYPE(t_patch), POINTER :: patch
80  TYPE(t_global), POINTER :: global
81  TYPE(t_turb), POINTER :: turb
82 
83  INTEGER :: ijkc0,ijkc1,ijkn, indcp, ibn,ien
84  REAL(RFREAL) :: one6th,beta,engmodel,rprt,cpprt,muel,tcol,tcot
85  REAL(RFREAL) :: muf,tx,ty,tz,tgradf
86  REAL(RFREAL) :: velf(3),fd(4),sface(3),sij(3,3)
87  REAL(RFREAL) :: velc0(3),velc1(3),tfd(4),fd4,ucitauij,taukk
88  REAL(RFREAL), POINTER :: diss(:,:),tv(:,:),gv(:,:),muet(:,:),srate(:,:)
89  REAL(RFREAL), POINTER :: trace(:),vol(:)
90  LOGICAL :: dowlm
91 
92 #ifdef RFLO
93  INTEGER :: ipcbeg, ipcend, jpcbeg, jpcend, kpcbeg, kpcend
94  INTEGER :: ilev,icoff,ijcoff,inoff,ijnoff
95  REAL(RFREAL), POINTER :: avgco(:,:), sf(:,:), dv(:,:), grad(:,:)
96 #endif
97 #ifdef RFLU
98  INTEGER, POINTER :: f2c(:,:)
99  REAL(RFREAL) :: rdens0, rdens1
100  REAL(RFREAL), POINTER :: fn(:,:), cv(:,:), grad(:,:,:)
101 #endif
102 !$endcopy
103 !******************************************************************************
104 
105  rcsidentstring = '$RCSfile: TURB_VFluxHybrid.F90,v $'
106 
107  global => region%global
108  CALL registerfunction( global,'TURB_VFluxHybrid',&
109  'TURB_VFluxHybrid.F90' )
110 
111 !$startcopy TURB_VisFluxEddy
112 ! get dimensions and pointers ------------------------------------------------
113 
114 #ifdef RFLO
115  ilev = region%currLevel
116  dv => region%levels(ilev)%mixt%dv
117  tv => region%levels(ilev)%mixt%tv
118  gv => region%levels(ilev)%mixt%gv
119  diss => region%levels(ilev)%mixt%diss
120  muet => region%levels(ilev)%turb%mueT
121  trace => region%levels(ilev)%turb%trace
122  vol => region%levels(ilev)%grid%vol
123  turb => region%levels(ilev)%turb
124 
125  CALL rflo_getdimensphys( region,ilev,ipcbeg,ipcend, &
126  jpcbeg,jpcend,kpcbeg,kpcend )
127  CALL rflo_getcelloffset( region,ilev,icoff,ijcoff )
128  CALL rflo_getnodeoffset( region,ilev,inoff,ijnoff )
129 #endif
130 #ifdef RFLU
131  cv => region%mixt%cv
132  tv => region%mixt%tv
133  gv => region%mixt%gv
134  diss => region%mixt%diss
135  muet => region%turb%mueT
136  trace => region%turb%trace
137  vol => region%grid%vol
138  turb => region%turb
139 #endif
140 
141 ! get coefficients and parameters
142 
143  one6th = 1._rfreal/6._rfreal
144 
145  IF (region%turbInput%engModel==off) THEN
146  engmodel = 0._rfreal
147  ELSE
148  engmodel = 1._rfreal
149  ENDIF
150  beta = region%mixtInput%betrk(region%irkStep)
151 
152 #ifdef RFLO
153  rprt = 1._rfreal/region%levels(ilev)%mixt%prTurb
154  indcp = region%levels(ilev)%mixt%indCp
155 #endif
156 #ifdef RFLU
157  rprt = 1._rfreal/region%mixtInput%prTurb
158  indcp = region%mixtInput%indCp
159 #endif
160 
161 ! compute interior fluxes through I, J and K faces
162 
163 #ifdef STATS
164 ! if desired, collect quantities of interest; for that, allocate work space
165  ibn = lbound( trace,1 )
166  ien = ubound( trace,1 )
167  ALLOCATE( turb%stwork(2,ibn:ien) )
168  turb%stwork = 0._rfreal
169 #endif
170 
171  CALL computefluxtot( diri )
172 #ifdef RFLO
173  CALL computefluxtot( dirj )
174  CALL computefluxtot( dirk )
175 #endif
176 
177 ! get fluxes through boundaries
178 #ifdef RFLO
179  DO ipatch=1,region%nPatches
180  patch => region%levels(ilev)%patches(ipatch)
181 #endif
182 #ifdef RFLU
183  DO ipatch=1,region%grid%nPatches
184  patch => region%patches(ipatch)
185 #endif
186 
187  dowlm = .false.
188 #ifdef RFLO
189  IF (patch%bcType>=bc_noslipwall .AND. &
190  patch%bcType<=bc_noslipwall+bc_range) THEN ! my boundary type
191  IF (patch%valBola%switches(wlm_input_model) /= wlm_model_nomodel) THEN
192  dowlm = .true.
193  ENDIF
194  ENDIF
195 #endif
196 !$endcopy
197 
198  IF (dowlm) THEN
199  CALL turb_wlmupdate( region,patch )
200  CALL turb_wlmtauwallmapping( region,patch )
201  CALL turb_wlmfluxpatch( region,patch )
202  ELSE
203  CALL turb_vfluxhybridpatch( region,patch )
204  ENDIF
205  ENDDO
206 
207 !$startcopy TURB_VisFluxEddy
208 #ifdef STATS
209 ! collect statistics and deallocate work space used
210 
211  IF ((region%turbInput%nSt > 0) .AND. (region%turbInput%engModel/=off)) THEN
212  turb%stwork(:,:) = one6th*turb%stwork(:,:)
213  CALL turb_statccollector( region,1,2,turb%stwork )
214  ENDIF
215  DEALLOCATE( turb%stwork )
216 #endif
217 
218 ! if active, get contribution of energy subgrid model alpha_4
219  IF (region%turbInput%engModel/=off) THEN
220  CALL turb_lesesgmodel4( region )
221  ENDIF
222 
223 ! finalize --------------------------------------------------------------------
224 
225  CALL deregisterfunction( global )
226 
227 ! =============================================================================
228 ! Flux computation subroutine
229 ! =============================================================================
230 
231 CONTAINS
232 
233  SUBROUTINE computefluxtot( ijk )
234 
235 ! ... parameters
236  INTEGER :: ijk
237 
238 ! ... local variables
239  INTEGER :: ibeg,iend,jbeg,jend,kbeg,kend, iadd,jadd,kadd
240  REAL(RFREAL) :: ac0, ac1
241 !$endcopy
242 !pert TURB
243  REAL(RFREAL) :: tij(3,3),sijsij,sijtij
244  REAL(RFREAL), POINTER :: tauij(:,:)
245 !endpert
246 
247 ! - Compute turbulent part tauij from scale-similarity; for Bardina model,
248 ! face values of cv have yet to be obtained.
249 
250 #ifdef RFLO
251  IF (ijk==diri) THEN
252  tauij => region%levels(ilev)%turb%fISij
253  ELSEIF (ijk==dirj) THEN
254  tauij => region%levels(ilev)%turb%fJSij
255  ELSEIF (ijk==dirk) THEN
256  tauij => region%levels(ilev)%turb%fKSij
257  ENDIF
258 #endif
259 #ifdef RFLU
260  tauij => region%turb%fISij
261 #endif
262 
263 !$startcopy TURB_VisFluxEddy
264 ! - Set limits and pointers ---------------------------------------------------
265 
266 #ifdef RFLO
267  IF (ijk==diri) THEN
268  ibeg = ipcbeg+1
269  iend = ipcend
270  jbeg = jpcbeg
271  jend = jpcend
272  kbeg = kpcbeg
273  kend = kpcend
274  iadd = -1
275  jadd = 0
276  kadd = 0
277  srate => region%levels(ilev)%turb%mISij
278  grad => region%levels(ilev)%mixt%gradi
279  sf => region%levels(ilev)%grid%si
280  avgco => region%levels(ilev)%grid%c2fCoI
281  ELSEIF (ijk==dirj) THEN
282  ibeg = ipcbeg
283  iend = ipcend
284  jbeg = jpcbeg+1
285  jend = jpcend
286  kbeg = kpcbeg
287  kend = kpcend
288  iadd = 0
289  jadd = -1
290  kadd = 0
291  srate => region%levels(ilev)%turb%mJSij
292  grad => region%levels(ilev)%mixt%gradj
293  sf => region%levels(ilev)%grid%sj
294  avgco => region%levels(ilev)%grid%c2fCoJ
295  ELSEIF (ijk==dirk) THEN
296  ibeg = ipcbeg
297  iend = ipcend
298  jbeg = jpcbeg
299  jend = jpcend
300  kbeg = kpcbeg+1
301  kend = kpcend
302  iadd = 0
303  jadd = 0
304  kadd = -1
305  srate => region%levels(ilev)%turb%mKSij
306  grad => region%levels(ilev)%mixt%gradk
307  sf => region%levels(ilev)%grid%sk
308  avgco => region%levels(ilev)%grid%c2fCoK
309  ENDIF
310 #endif
311 #ifdef RFLU
312  ibeg = 1
313  iend = region%grid%nFaces
314  f2c => region%grid%f2c
315  srate => region%turb%mISij
316  grad => region%mixt%gradFace
317  fn => region%grid%fn
318 #endif
319 
320 #ifdef RFLO
321  DO k=kbeg,kend
322  DO j=jbeg,jend
323  DO i=ibeg,iend
324 
325  ijkc0 = indijk(i ,j ,k ,icoff,ijcoff)
326  ijkc1 = indijk(i+iadd,j+jadd,k+kadd,icoff,ijcoff)
327  ijkn = indijk(i ,j ,k ,inoff,ijnoff)
328  sface(1)= sf(xcoord,ijkn)
329  sface(2)= sf(ycoord,ijkn)
330  sface(3)= sf(zcoord,ijkn)
331  ac0 = avgco(2,ijkn)
332  ac1 = avgco(1,ijkn)
333 #endif
334 #ifdef RFLU
335 ! - specific Rocflu, check the state of cv first
336  IF (region%mixt%cvState /= cv_mixt_state_cons) &
337  CALL errorstop(global,err_reached_default,__line__)
338  ac0 = 0.5_rfreal
339  ac1 = 0.5_rfreal
340  DO ic=ibeg,iend
341  ijkc0 = f2c(1,ic)
342  ijkc1 = f2c(2,ic)
343  ijkn = ic
344  sface(1)= fn(xcoord,ijkn)*fn(xyzmag,ijkn)
345  sface(2)= fn(ycoord,ijkn)*fn(xyzmag,ijkn)
346  sface(3)= fn(zcoord,ijkn)*fn(xyzmag,ijkn)
347  rdens0 = 1._rfreal/cv(cv_mixt_dens,ijkc0)
348  rdens1 = 1._rfreal/cv(cv_mixt_dens,ijkc1)
349 #endif
350  cpprt = (ac0*gv(gv_mixt_cp,ijkc0*indcp) + &
351  ac1*gv(gv_mixt_cp,ijkc1*indcp))*rprt
352  muel = ac0*tv(tv_mixt_muel,ijkc0)+ac1*tv(tv_mixt_muel,ijkc1)
353  tcol = ac0*tv(tv_mixt_tcol,ijkc0)+ac1*tv(tv_mixt_tcol,ijkc1)
354  tcot = cpprt*muet(ijk,ijkn)
355 #ifdef RFLO
356  velc0(1)= dv(dv_mixt_uvel,ijkc0)
357  velc0(2)= dv(dv_mixt_vvel,ijkc0)
358  velc0(3)= dv(dv_mixt_wvel,ijkc0)
359  velc1(1)= dv(dv_mixt_uvel,ijkc1)
360  velc1(2)= dv(dv_mixt_vvel,ijkc1)
361  velc1(3)= dv(dv_mixt_wvel,ijkc1)
362  velf(1) = ac0*dv(dv_mixt_uvel,ijkc0)+ac1*dv(dv_mixt_uvel,ijkc1)
363  velf(2) = ac0*dv(dv_mixt_vvel,ijkc0)+ac1*dv(dv_mixt_vvel,ijkc1)
364  velf(3) = ac0*dv(dv_mixt_wvel,ijkc0)+ac1*dv(dv_mixt_wvel,ijkc1)
365  tx = grad(gr_mixt_tx,ijkn)
366  ty = grad(gr_mixt_ty,ijkn)
367  tz = grad(gr_mixt_tz,ijkn)
368 #endif
369 #ifdef RFLU
370  velc0(1)= cv(cv_mixt_xmom,ijkc0)*rdens0
371  velc0(2)= cv(cv_mixt_ymom,ijkc0)*rdens0
372  velc0(3)= cv(cv_mixt_zmom,ijkc0)*rdens0
373  velc1(1)= cv(cv_mixt_xmom,ijkc1)*rdens1
374  velc1(2)= cv(cv_mixt_ymom,ijkc1)*rdens1
375  velc1(3)= cv(cv_mixt_zmom,ijkc1)*rdens1
376  velf(1) = ac0*velc0(1)+ac1*velc1(1)
377  velf(2) = ac0*velc0(2)+ac1*velc1(2)
378  velf(3) = ac0*velc0(3)+ac1*velc1(3)
379  tx = grad(xcoord,grf_mixt_temp,ijkn)
380  ty = grad(ycoord,grf_mixt_temp,ijkn)
381  tz = grad(zcoord,grf_mixt_temp,ijkn)
382 #endif
383 
384  tgradf= (tx*sface(1)+ty*sface(2)+tz*sface(3))
385 
386  sij(1,1) = srate(e11,ijkn)
387  sij(1,2) = srate(e12,ijkn)
388  sij(1,3) = srate(e13,ijkn)
389 
390  sij(2,1) = sij(1,2)
391  sij(2,2) = srate(e22,ijkn)
392  sij(2,3) = srate(e23,ijkn)
393 
394  sij(3,1) = sij(1,3)
395  sij(3,2) = sij(2,3)
396  sij(3,3) = srate(e33,ijkn)
397 
398  fd(1) = muel*(sij(1,1)*sface(1)+sij(1,2)*sface(2)+sij(1,3)*sface(3))
399  fd(2) = muel*(sij(2,1)*sface(1)+sij(2,2)*sface(2)+sij(2,3)*sface(3))
400  fd(3) = muel*(sij(3,1)*sface(1)+sij(3,2)*sface(2)+sij(3,3)*sface(3))
401  fd(4) = dot_product(fd(1:3),velf(1:3)) + tcol*tgradf
402 !$endcopy
403  tij(1,1) = tauij(e11,ijkn) - muet(ijk,ijkn)*sij(1,1)
404  tij(1,2) = tauij(e12,ijkn) - muet(ijk,ijkn)*sij(1,2)
405  tij(1,3) = tauij(e13,ijkn) - muet(ijk,ijkn)*sij(1,3)
406 
407  tij(2,1) = tauij(e12,ijkn) - muet(ijk,ijkn)*sij(2,1)
408  tij(2,2) = tauij(e22,ijkn) - muet(ijk,ijkn)*sij(2,2)
409  tij(2,3) = tauij(e23,ijkn) - muet(ijk,ijkn)*sij(2,3)
410 
411  tij(3,1) = tauij(e13,ijkn) - muet(ijk,ijkn)*sij(3,1)
412  tij(3,2) = tauij(e23,ijkn) - muet(ijk,ijkn)*sij(3,2)
413  tij(3,3) = tauij(e33,ijkn) - muet(ijk,ijkn)*sij(3,3)
414 
415  taukk = tij(1,1)+tij(2,2)+tij(3,3)
416  taukk = max( 0._rfreal,taukk )
417  trace(ijkc0) = trace(ijkc0) + taukk
418  trace(ijkc1) = trace(ijkc1) + taukk
419 
420  tfd(1) = tij(1,1)*sface(1)+tij(1,2)*sface(2)+tij(1,3)*sface(3)
421  tfd(2) = tij(2,1)*sface(1)+tij(2,2)*sface(2)+tij(2,3)*sface(3)
422  tfd(3) = tij(3,1)*sface(1)+tij(3,2)*sface(2)+tij(3,3)*sface(3)
423 
424  fd(1) = fd(1)-tfd(1)
425  fd(2) = fd(2)-tfd(2)
426  fd(3) = fd(3)-tfd(3)
427 
428 !$startcopy TURB_VisFluxEddy
429  fd4 = fd(4)
430  ucitauij = dot_product(tfd(1:3),velc0(1:3))
431  tfd(4) = ucitauij - tcot*tgradf
432  fd(4) = fd4 - engmodel*tfd(4)
433 
434 #ifdef STATS
435  turb%stwork(1,ijkc0)=turb%stwork(1,ijkc0)+ucitauij ! alp_1
436  turb%stwork(2,ijkc0)=turb%stwork(2,ijkc0)+tcot*tgradf ! alp_2+alp_3
437 #endif
438  diss(cv_mixt_xmom,ijkc0) = diss(cv_mixt_xmom,ijkc0) + fd(1)*beta
439  diss(cv_mixt_ymom,ijkc0) = diss(cv_mixt_ymom,ijkc0) + fd(2)*beta
440  diss(cv_mixt_zmom,ijkc0) = diss(cv_mixt_zmom,ijkc0) + fd(3)*beta
441  diss(cv_mixt_ener,ijkc0) = diss(cv_mixt_ener,ijkc0) + fd(4)*beta
442  global%esg1Sum = global%esg1Sum + vol(ijkc0)*ucitauij
443 
444  ucitauij = dot_product(tfd(1:3),velc1(1:3))
445  tfd(4) = ucitauij - tcot*tgradf
446  fd(4) = fd4 - engmodel*tfd(4)
447 
448 #ifdef STATS
449  turb%stwork(1,ijkc1)=turb%stwork(1,ijkc1)+ucitauij ! alp_1
450  turb%stwork(2,ijkc1)=turb%stwork(2,ijkc1)+tcot*tgradf ! alp_2+alp_3
451 #endif
452  diss(cv_mixt_xmom,ijkc1) = diss(cv_mixt_xmom,ijkc1) - fd(1)*beta
453  diss(cv_mixt_ymom,ijkc1) = diss(cv_mixt_ymom,ijkc1) - fd(2)*beta
454  diss(cv_mixt_zmom,ijkc1) = diss(cv_mixt_zmom,ijkc1) - fd(3)*beta
455  diss(cv_mixt_ener,ijkc1) = diss(cv_mixt_ener,ijkc1) - fd(4)*beta
456  global%esg1Sum = global%esg1Sum - vol(ijkc1)*ucitauij
457 
458 #ifdef RFLO
459  ENDDO ! i
460  ENDDO ! j
461  ENDDO ! k
462 #endif
463 #ifdef RFLU
464  ENDDO ! iC
465 #endif
466 
467 !$endcopy
468 
469 #ifdef RFLO
470 ! compute equivalent eddy viscosity from scale similarity
471 
472  IF (ijk==diri) THEN
473  ibeg = ipcbeg-1
474  iend = ipcend+2
475  jbeg = jpcbeg-1
476  jend = jpcend+1
477  kbeg = kpcbeg-1
478  kend = kpcend+1
479  ELSEIF (ijk==dirj) THEN
480  ibeg = ipcbeg-1
481  iend = ipcend+1
482  jbeg = jpcbeg-1
483  jend = jpcend+2
484  kbeg = kpcbeg-1
485  kend = kpcend+1
486  ELSEIF (ijk==dirk) THEN
487  ibeg = ipcbeg-1
488  iend = ipcend+1
489  jbeg = jpcbeg-1
490  jend = jpcend+1
491  kbeg = kpcbeg-1
492  kend = kpcend+2
493  ENDIF
494 
495  DO k=kbeg,kend
496  DO j=jbeg,jend
497  DO i=ibeg,iend
498  ijkn = indijk(i ,j ,k ,inoff,ijnoff)
499 #endif
500 #ifdef RFLU
501  ibeg = 1
502  iend = region%grid%nFaces
503 
504  DO ijkn=ibeg,iend
505 #endif
506  sij(1,1) = srate(e11,ijkn)
507  sij(1,2) = srate(e12,ijkn)
508  sij(1,3) = srate(e13,ijkn)
509 
510  sij(2,1) = sij(1,2)
511  sij(2,2) = srate(e22,ijkn)
512  sij(2,3) = srate(e23,ijkn)
513 
514  sij(3,1) = sij(1,3)
515  sij(3,2) = sij(2,3)
516  sij(3,3) = srate(e33,ijkn)
517 
518 ! ------- S_{ij}*S_{ij}:
519  sijsij = sij(1,1)*sij(1,1) &
520  +sij(2,2)*sij(2,2) &
521  +sij(3,3)*sij(3,3) &
522  +2.0d0*sij(1,2)*sij(1,2) &
523  +2.0d0*sij(1,3)*sij(1,3) &
524  +2.0d0*sij(2,3)*sij(2,3)
525 
526 ! ------- S_{ij}*Tau_{ij}:
527  sijtij = sij(1,1)*tauij(e11,ijkn) &
528  +sij(2,2)*tauij(e22,ijkn) &
529  +sij(3,3)*tauij(e33,ijkn) &
530  +2.0d0*sij(1,2)*tauij(e12,ijkn) &
531  +2.0d0*sij(1,3)*tauij(e13,ijkn) &
532  +2.0d0*sij(2,3)*tauij(e23,ijkn)
533 
534  muf = max( 0._rfreal, -(sijtij/max( sijsij,real_small )) )
535  muet(ijk,ijkn) = muet(ijk,ijkn) + muf
536 #ifdef RFLO
537  ENDDO ! i
538  ENDDO ! j
539  ENDDO ! k
540 #endif
541 #ifdef RFLU
542  ENDDO ! ijkN
543 #endif
544 
545  END SUBROUTINE computefluxtot
546 
547 END SUBROUTINE turb_vfluxhybrid
548 
549 !******************************************************************************
550 !
551 ! RCS Revision history:
552 !
553 ! $Log: TURB_VFluxHybrid.F90,v $
554 ! Revision 1.10 2008/12/06 08:44:42 mtcampbe
555 ! Updated license.
556 !
557 ! Revision 1.9 2008/11/19 22:17:54 mtcampbe
558 ! Added Illinois Open Source License/Copyright
559 !
560 ! Revision 1.8 2005/12/30 23:20:56 wasistho
561 ! exclude rocflu from WLM treatment
562 !
563 ! Revision 1.7 2004/10/22 23:20:33 wasistho
564 ! collect esg1 and esg2+3 into statistics
565 !
566 ! Revision 1.6 2004/08/02 19:33:58 wasistho
567 ! changed grid%avgCo to grid%c2fCo
568 !
569 ! Revision 1.5 2004/07/31 00:53:47 wasistho
570 ! replaced cell2face midpoint by linear averaging
571 !
572 ! Revision 1.4 2004/05/29 03:38:35 wasistho
573 ! extended computations of muet to first dummy layers
574 !
575 ! Revision 1.3 2004/03/27 02:16:42 wasistho
576 ! compiled with Rocflu
577 !
578 ! Revision 1.2 2004/03/24 03:37:02 wasistho
579 ! prepared for RFLU
580 !
581 ! Revision 1.1 2004/03/05 04:37:00 wasistho
582 ! changed nomenclature
583 !
584 ! Revision 1.7 2004/02/26 21:22:31 wasistho
585 ! changed turb%esg.. to global%esg..
586 !
587 ! Revision 1.6 2003/08/29 01:43:00 wasistho
588 ! Added TARGET attribute to region variable, since pointers are cached into it
589 !
590 ! Revision 1.5 2003/06/05 19:20:23 wasistho
591 ! implemented heat transfer model
592 !
593 ! Revision 1.4 2003/05/31 01:48:23 wasistho
594 ! installed turb. wall layer model
595 !
596 ! Revision 1.3 2003/05/15 02:57:06 jblazek
597 ! Inlined index function.
598 !
599 ! Revision 1.2 2003/02/04 04:00:11 wasistho
600 ! An Esgs-bug which caused fd(4) unsymmetry fixed
601 !
602 ! Revision 1.1 2002/10/14 23:55:30 wasistho
603 ! Install Rocturb
604 !
605 !******************************************************************************
606 
607 
608 
609 
610 
611 
612 
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE SUBROUTINE ibeg
j indices k indices k
Definition: Indexing.h:6
subroutine turb_vfluxhybrid(region)
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE SUBROUTINE kpcbeg
Vector_n max(const Array_n_const &v1, const Array_n_const &v2)
Definition: Vector_n.h:354
subroutine registerfunction(global, funName, fileName)
Definition: ModError.F90:449
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE SUBROUTINE ic
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE SUBROUTINE jpcbeg
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE SUBROUTINE ipcend
subroutine turb_wlmtauwallmapping(region, patch)
subroutine rflo_getnodeoffset(region, iLev, iNodeOffset, ijNodeOffset)
subroutine turb_wlmupdate(region, patch)
Definition: patch.h:74
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE SUBROUTINE knode iend
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE SUBROUTINE ipcbeg
blockLoc i
Definition: read.cpp:79
subroutine rflo_getcelloffset(region, iLev, iCellOffset, ijCellOffset)
subroutine turb_statccollector(region, iBegSt, iEndSt, colVar)
subroutine turb_lesesgmodel4(region)
Tfloat trace() const
Return the trace of the image, viewed as a matrix.
Definition: CImg.h:13217
j indices j
Definition: Indexing.h:6
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE SUBROUTINE jpcend
subroutine turb_vfluxhybridpatch(region, patch)
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE SUBROUTINE knode jend
subroutine computefluxtot(ijk)
subroutine errorstop(global, errorCode, errorLine, addMessage)
Definition: ModError.F90:483
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE SUBROUTINE knode jbeg
long double dot_product(pnt vec1, pnt vec2)
**********************************************************************Rocstar Simulation Suite Illinois Rocstar LLC All rights reserved ****Illinois Rocstar LLC IL **www illinoisrocstar com **sales illinoisrocstar com WITHOUT WARRANTY OF ANY **EXPRESS OR INCLUDING BUT NOT LIMITED TO THE WARRANTIES **OF FITNESS FOR A PARTICULAR PURPOSE AND **NONINFRINGEMENT IN NO EVENT SHALL THE CONTRIBUTORS OR **COPYRIGHT HOLDERS BE LIABLE FOR ANY DAMAGES OR OTHER WHETHER IN AN ACTION OF TORT OR **Arising OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE **USE OR OTHER DEALINGS WITH THE SOFTWARE **********************************************************************INTERFACE SUBROUTINE knode kbeg
subroutine deregisterfunction(global)
Definition: ModError.F90:469
subroutine turb_wlmfluxpatch(region, patch)
subroutine rflo_getdimensphys(region, iLev, ipcbeg, ipcend, jpcbeg, jpcend, kpcbeg, kpcend)