Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
TURB_LesCalcEddyVis.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 turbulent viscosity mu_t.
26 !
27 ! Description: The eddy viscosity is computed at faces using the model
28 ! coefficient depending on the model selected by user.
29 !
30 ! Input: region = data of current region
31 ! ibn,ien = begin and end nodes index
32 ! ijk = integer identifying wich face is being treated
33 !
34 ! Output: mueT(ijk,:) = eddy viscosity at ijk face
35 !
36 ! Notes: The face values contribution are collected here for model coeffs.
37 ! at cell centers.
38 !
39 !******************************************************************************
40 !
41 ! $Id: TURB_LesCalcEddyVis.F90,v 1.17 2009/08/12 04:15:59 mtcampbe Exp $
42 !
43 ! Copyright: (c) 2001 by the University of Illinois
44 !
45 !******************************************************************************
46 
47 SUBROUTINE turb_lescalceddyvis( region,ibn,ien,ijk )
48 
49  USE moddatatypes
50  USE moddatastruct, ONLY : t_region
51  USE modturbulence, ONLY : t_turb
52  USE modglobal, ONLY : t_global
53 #ifdef RFLU
54  USE modbndpatch, ONLY : t_patch
56 #endif
57 #ifdef RFLO
58  USE modinterfaces, ONLY : rflo_getdimensphys, &
61 
62 #include "Indexing.h"
63 #endif
65  USE moderror
66  USE modparameters
68  IMPLICIT NONE
69 
70 ! ... parameters
71 #ifdef RFLO
72  TYPE(t_region), TARGET :: region
73 #endif
74 #ifdef RFLU
75  TYPE(t_region), POINTER :: region
76 #endif
77  INTEGER :: ibn,ien,ijk
78 
79 ! ... loop variables
80  INTEGER :: i, j, k, m, ijkn, ipatch, ifl
81 
82 ! ... local variables
83  CHARACTER(CHRLEN) :: rcsidentstring
84  TYPE(t_global), POINTER :: global
85  TYPE(t_turb), POINTER :: turb
86 #ifdef RFLU
87  TYPE(t_patch), POINTER :: patch
88 #endif
89 
90  INTEGER :: turbmodel, zofid, xco, yco, zco, errorflag
91  INTEGER :: ibeg, iend, jbeg, jend, kbeg, kend, ijkc0, ijkc1
92  REAL(RFREAL) :: two3rd,csmag,delfac2,delta2, zoffac,bzoffac
93  REAL(RFREAL) :: vormag,volface,rhoface,muturb,rhomin,mutmin
94  REAL(RFREAL) :: sxx,sxy,sxz,syy,syz,szz
95  REAL(RFREAL), POINTER :: fvol(:),fvar(:,:),cmodel(:,:),muet(:,:)
96  REAL(RFREAL), POINTER :: tdv(:,:),sij(:,:),zof(:,:,:)
97 #ifdef RFLU
98  REAL(RFREAL), POINTER :: bfvol(:),bfvar(:,:),bcmodel(:,:),bmuet(:,:),bsij(:,:)
99  REAL(RFREAL), POINTER :: bzof(:,:,:)
100 #endif
101 
102 #ifdef RFLO
103  INTEGER :: ipcbeg, ipcend, jpcbeg, jpcend, kpcbeg, kpcend, iadd, jadd, kadd
104  INTEGER :: ilev, icoff, ijcoff, inoff, ijnoff, filtertype
105 #endif
106 #ifdef RFLU
107  INTEGER :: bctype
108  INTEGER :: ifg,ifgt,ifgbeg,ifgtbeg
109  INTEGER :: npatches,nfaces,nfacestot,nbfaces,nbfacestot
110  REAL(RFREAL) :: bvalfactor
111 #endif
112 
113 !******************************************************************************
114 
115  rcsidentstring = '$RCSfile: TURB_LesCalcEddyVis.F90,v $'
116 
117  global => region%global
118  CALL registerfunction( global,'Turb_LesCalcEddyVis',&
119  'TURB_LesCalcEddyVis.F90' )
120 
121 ! get parameters and coefficients -------------------------------------------
122 
123  xco = xcoord
124  yco = ycoord
125  zco = zcoord
126  zofid = zof_les_eddyvis
127  turbmodel = region%mixtInput%turbModel
128  csmag = region%turbInput%cSmag
129  delfac2 = region%turbInput%delFac2
130 #ifdef RFLO
131  filtertype = region%turbInput%filterType
132  ilev = region%currLevel
133 #endif
134 
135  rhomin = 10000._rfreal
136  mutmin = 10000._rfreal
137  two3rd = 2._rfreal/3._rfreal
138 
139 ! allocate required arrays within this scope
140 
141 #ifdef RFLO
142  turb => region%levels(ilev)%turb
143 
144  ALLOCATE( turb%coef(1,ibn:ien),stat=errorflag )
145  ALLOCATE( turb%fVar(cv_turb_nelm,ibn:ien),stat=errorflag )
146  global%error = errorflag
147  IF (global%error /= 0) CALL errorstop( global,err_allocate,__line__ )
148 #endif
149 #ifdef RFLU
150  npatches = region%grid%nPatches
151  nfaces = region%grid%nFaces
152  nfacestot = region%grid%nFacesTot
153  nbfaces = 0
154  nbfacestot = 0
155 
156  DO ipatch = 1,npatches
157  patch => region%patches(ipatch)
158 
159  nbfaces = nbfaces + patch%nBTris + patch%nBQuads
160  nbfacestot = nbfacestot + patch%nBTrisTot + patch%nBQuadsTot
161  END DO ! iPatch
162 
163  IF (ijk /= diri) THEN
164  CALL errorstop( global,err_turb_fixparam,__line__,'RFLU but ijk/=DIRI' )
165  ENDIF
166  turb => region%turb
167 
168  ALLOCATE( turb%coef( 1,ibn:ien),stat=errorflag )
169  ALLOCATE( turb%bCoef(1,nbfaces),stat=errorflag )
170  global%error = errorflag
171  IF (global%error /= 0) CALL errorstop( global,err_allocate,__line__ )
172 
173  ALLOCATE( turb%fVar( cv_turb_nelm, nfaces),stat=errorflag )
174  ALLOCATE( turb%bfVar(cv_turb_nelm,nbfaces),stat=errorflag )
175  global%error = errorflag
176  IF (global%error /= 0) CALL errorstop( global,err_allocate,__line__ )
177 #endif
178 
179 ! get model coefficients -----------------------------------------------
180 
181  IF (turbmodel==turb_model_fixsmag) THEN
182 
183 #ifdef RFLO
184  CALL turb_flolesgenc2f( region,ijk )
185  DO ijkn = ibn,ien
186  turb%fVar(cv_turb_dens,ijkn) = 1._rfreal/turb%fVar(cv_turb_dens,ijkn)
187  turb%coef(1,ijkn) = csmag*csmag
188  ENDDO
189 #endif
190 #ifdef RFLU
191  CALL turb_flulesc2f( region )
192  DO ijkn = ibn,ien
193  turb%coef(1,ijkn) = csmag*csmag
194  turb%fVar(cv_turb_dens,ijkn) = 1._rfreal/turb%fVar(cv_turb_dens,ijkn)
195  ENDDO
196  DO ijkn = 1,nbfaces
197  turb%bCoef(1,ijkn) = csmag*csmag
198  turb%bfVar(cv_turb_dens,ijkn) = 1._rfreal/turb%bfVar(cv_turb_dens,ijkn)
199  ENDDO
200 #endif
201  ELSEIF (turbmodel==turb_model_dynsmag) THEN
202  CALL turb_lescoefdynsmag( region,ibn,ien,ijk )
203  ELSEIF (turbmodel==turb_model_dynmixd) THEN
204  CALL turb_lescoefdynmixd( region,ibn,ien,ijk )
205  ENDIF
206 
207 ! apply spatial tripping to eddy viscosity
208 
209 #ifdef RFLO
210  IF (ijk==diri) THEN
211  zof => turb%zofi
212  ELSEIF (ijk==dirj) THEN
213  zof => turb%zofj
214  ELSEIF (ijk==dirk) THEN
215  zof => turb%zofk
216  ENDIF
217  DO ijkn = ibn,ien
218  zoffac = zof(xco,zofid,ijkn)*zof(yco,zofid,ijkn)*zof(zco,zofid,ijkn)
219  turb%coef(1,ijkn) = zoffac*turb%coef(1,ijkn)
220  ENDDO
221 #endif
222 #ifdef RFLU
223  zof => turb%zofi
224  bzof => turb%bZofi
225  DO ijkn = ibn,ien
226  zoffac = zof(xco,zofid,ijkn)*zof(yco,zofid,ijkn)*zof(zco,zofid,ijkn)
227  turb%coef(1,ijkn) = zoffac*turb%coef(1,ijkn)
228  ENDDO
229  DO ijkn = 1,nbfaces
230  bzoffac = bzof(xco,zofid,ijkn)*bzof(yco,zofid,ijkn)*bzof(zco,zofid,ijkn)
231  turb%bCoef(1,ijkn) = bzoffac*turb%bCoef(1,ijkn)
232  ENDDO
233 #endif
234 
235 ! get dimensions and pointers
236 
237  tdv => turb%dv
238  fvar => turb%fVar
239  cmodel => turb%coef
240  muet => turb%mueT
241 
242 #ifdef RFLO
243 
244  CALL rflo_getdimensphys( region,ilev,ipcbeg,ipcend, &
245  jpcbeg,jpcend,kpcbeg,kpcend )
246  CALL rflo_getcelloffset( region,ilev,icoff,ijcoff )
247  CALL rflo_getnodeoffset( region,ilev,inoff,ijnoff )
248 
249  IF (ijk==diri) THEN
250 ! - pointers and indices for i-direction
251  ibeg = ipcbeg-1
252  iend = ipcend+2
253  jbeg = jpcbeg-1
254  jend = jpcend+1
255  kbeg = kpcbeg-1
256  kend = kpcend+1
257  iadd = -1
258  jadd = 0
259  kadd = 0
260  fvol => turb%fvolI
261  sij => turb%mISij
262  ELSEIF (ijk==dirj) THEN
263 ! - pointers and indices for j-direction
264  ibeg = ipcbeg-1
265  iend = ipcend+1
266  jbeg = jpcbeg-1
267  jend = jpcend+2
268  kbeg = kpcbeg-1
269  kend = kpcend+1
270  iadd = 0
271  jadd = -1
272  kadd = 0
273  fvol => turb%fvolJ
274  sij => turb%mJSij
275  ELSEIF (ijk==dirk) THEN
276 ! - pointers and indices for k-direction
277  ibeg = ipcbeg-1
278  iend = ipcend+1
279  jbeg = jpcbeg-1
280  jend = jpcend+1
281  kbeg = kpcbeg-1
282  kend = kpcend+2
283  iadd = 0
284  jadd = 0
285  kadd = -1
286  fvol => turb%fvolK
287  sij => turb%mKSij
288  END IF
289 
290 ! compute eddy viscosity; note that face density is stored in fVar as 1/rho
291 ! which is leftover from lesLij; so we devide by fVar(1,:) instead of multiply
292 
293  DO k=kbeg,kend
294  DO j=jbeg,jend
295  DO i=ibeg,iend
296 
297  ijkc0 = indijk(i ,j ,k ,icoff,ijcoff)
298  ijkc1 = ijkc0 + iadd + jadd*icoff + kadd*ijcoff
299  ijkn = indijk(i ,j ,k ,inoff,ijnoff)
300 #endif
301 #ifdef RFLU
302 
303  fvol => turb%fvolI
304  sij => turb%mISij
305 
306  DO ijkn =ibn,ien
307 #endif
308 
309  sxx = sij(e11,ijkn)
310  sxy = sij(e12,ijkn)
311  sxz = sij(e13,ijkn)
312 
313  syy = sij(e22,ijkn)
314  syz = sij(e23,ijkn)
315 
316  szz = sij(e33,ijkn)
317 
318  vormag = 0.5_rfreal*(sxx*sxx+syy*syy+szz*szz)+sxy*sxy+sxz*sxz+syz*syz
319  vormag = sqrt(vormag)
320 
321  volface= fvol(ijkn)
322  delta2 = volface**two3rd ! another option: delta2=3.*volFace**two3rd
323 
324 ! rhoFace= 1._RFREAL/fVar(CV_TURB_DENS,ijkN)
325  muturb = cmodel(1,ijkn)/fvar(cv_turb_dens,ijkn)*delfac2*delta2*vormag
326 
327 ! Temporary clipping fix
328  muturb = max(muturb,real_small)
329 !
330 
331 ! rhoMin = MIN(rhoFace,rhoMin)
332 ! muTMin = MIN(muTurb ,muTMin)
333 
334  muet(ijk,ijkn) = muturb
335 
336 #ifdef RFLO
337  tdv(dv_turb_cdyn,ijkc0) = tdv(dv_turb_cdyn,ijkc0)+cmodel(1,ijkn)
338  tdv(dv_turb_cdyn,ijkc1) = tdv(dv_turb_cdyn,ijkc1)+cmodel(1,ijkn)
339  ENDDO ! i
340  ENDDO ! j
341  ENDDO ! k
342 #endif
343 #ifdef RFLU
344  ENDDO ! ijkN
345 
346 ! sum model coefficient from faces; the division by nfaces is in GetTvCell
347  DO ijkn = 1,region%grid%nFaces
348  ijkc0 = region%grid%f2c(1,ijkn)
349  ijkc1 = region%grid%f2c(2,ijkn)
350 
351  tdv(dv_turb_cdyn,ijkc0) = tdv(dv_turb_cdyn,ijkc0)+cmodel(1,ijkn)
352  tdv(dv_turb_cdyn,ijkc1) = tdv(dv_turb_cdyn,ijkc1)+cmodel(1,ijkn)
353  ENDDO ! ijkN
354 #endif
355 
356 ! RFLU boundary treatment =====================================================
357 
358 #ifdef RFLU
359  bfvol => turb%bfVolI
360  bfvar => turb%bfVar
361  bcmodel => turb%bCoef
362  bmuet => turb%bMueT
363  bsij => turb%bmISij
364 
365  DO ipatch = 1,region%grid%nPatches
366  patch => region%patches(ipatch)
367 ! TEMPORARY : removing usage of bf2bg from everywhere
368 ! ifgBeg = patch%bf2bg(BF2BG_BEG)
369 ! ifgtBeg= patch%bf2bgTot(BF2BG_BEG)
370 
371  bctype = patch%bcType
372  IF (bctype>=bc_noslipwall .AND. bctype<=bc_noslipwall+bc_range .OR. &
373  bctype>=bc_injection .AND. bctype<=bc_injection+bc_range) THEN
374  IF (turbmodel==turb_model_fixsmag) THEN
375  bvalfactor = 1._rfreal
376  ELSE
377  bvalfactor = 0._rfreal
378  ENDIF ! turbModel
379  ELSE
380  bvalfactor = 1._rfreal
381  ENDIF
382 
383 ! - sum model coefficient from bFaces; the division by nfaces is in GetTvCell
384 ! Note, bValFactor is needed for Rocflu to set eddy viscosity explicitly
385 ! to zero at injection and noslip-walls if desired.
386 
387  DO ifl = 1,patch%nBFaces
388  ijkc0 = patch%bf2c(ifl)
389  ifg = ifl + ifgbeg-1
390 
391  sxx = bsij(e11,ifg)
392  sxy = bsij(e12,ifg)
393  sxz = bsij(e13,ifg)
394 
395  syy = bsij(e22,ifg)
396  syz = bsij(e23,ifg)
397 
398  szz = bsij(e33,ifg)
399 
400  vormag = 0.5_rfreal*(sxx*sxx+syy*syy+szz*szz)+sxy*sxy+sxz*sxz+syz*syz
401  vormag = sqrt(vormag)
402 
403  volface= bfvol(ifg)
404  delta2 = volface**two3rd ! another option: delta2=3.*volFace**two3rd
405 
406 ! rhoFace= 1._RFREAL/bfVar(CV_TURB_DENS,ifg)
407  muturb = cmodel(1,ifg)/bfvar(cv_turb_dens,ifg)*delfac2*delta2*vormag
408 ! Temporary clipping fix
409  muturb = max(muturb,real_small)
410 !
411 
412 ! rhoMin = MIN(rhoFace,rhoMin)
413 ! muTMin = MIN(muTurb ,muTMin)
414 
415  bmuet(ijk,ifg) = bvalfactor*muturb
416  tdv(dv_turb_cdyn,ijkc0) = tdv(dv_turb_cdyn,ijkc0)+bcmodel(1,ifg)
417  ENDDO ! ifl
418  ENDDO ! iPatch
419 #endif
420 
421 ! check minimum density
422 
423  IF (rhomin < 0._rfreal) THEN
424  write(*,*) 'TURB_LesCalcEddyVis: rhoMin < 0'
425  ENDIF
426  IF (mutmin < 0._rfreal) THEN
427  write(*,*) 'TURB_LesCalcEddyVis: muTMin < 0'
428  ENDIF
429 
430 ! deallocate temporary arrays
431 
432  DEALLOCATE( turb%coef, turb%fVar )
433 #ifdef RFLU
434  DEALLOCATE( turb%bCoef, turb%bfVar )
435 #endif
436 
437 ! finalize --------------------------------------------------------------------
438 
439  CALL deregisterfunction( global )
440 
441 END SUBROUTINE turb_lescalceddyvis
442 
443 !******************************************************************************
444 !
445 ! RCS Revision history:
446 !
447 ! $Log: TURB_LesCalcEddyVis.F90,v $
448 ! Revision 1.17 2009/08/12 04:15:59 mtcampbe
449 ! Major update, bugfix from Abe development, more propagation compatibility,
450 ! some Rocstar IO changes, Ju's temporary clipping fix for turbulence. A bug
451 ! fix for initialization IO.
452 !
453 ! Revision 1.16 2008/12/06 08:44:41 mtcampbe
454 ! Updated license.
455 !
456 ! Revision 1.15 2008/11/19 22:17:53 mtcampbe
457 ! Added Illinois Open Source License/Copyright
458 !
459 ! Revision 1.14 2006/08/19 15:40:59 mparmar
460 ! Removed bf2bg, bf2bgTot
461 !
462 ! Revision 1.13 2006/01/17 17:51:55 wasistho
463 ! applied tripping to all eddy viscosity models
464 !
465 ! Revision 1.12 2006/01/13 07:13:49 wasistho
466 ! set rflu bValFactor=1 for FixSM and 0 else
467 !
468 ! Revision 1.11 2006/01/13 06:53:42 wasistho
469 ! set rflu bValFactor=1 for more stable sim
470 !
471 ! Revision 1.10 2006/01/12 09:48:31 wasistho
472 ! enabled tripping fixed Smagorinsky
473 !
474 ! Revision 1.9 2006/01/01 00:11:03 wasistho
475 ! multiplied muturb by bValFactor for rocflu
476 !
477 ! Revision 1.8 2005/12/29 19:51:00 wasistho
478 ! modified face indexing for sij
479 !
480 ! Revision 1.7 2004/07/30 22:33:54 wasistho
481 ! replaced floLesUniC2F by floLesGenC2F
482 !
483 ! Revision 1.6 2004/05/28 02:00:17 wasistho
484 ! update unstructured grid LES
485 !
486 ! Revision 1.5 2004/05/17 20:47:39 wasistho
487 ! compute first layer dummy mu_t instead of extrapolated
488 !
489 ! Revision 1.4 2004/03/27 02:16:42 wasistho
490 ! compiled with Rocflu
491 !
492 ! Revision 1.3 2004/03/19 02:48:36 wasistho
493 ! prepared for RFLU
494 !
495 ! Revision 1.2 2004/03/12 02:55:35 wasistho
496 ! changed rocturb routine names
497 !
498 ! Revision 1.1 2004/03/05 04:37:00 wasistho
499 ! changed nomenclature
500 !
501 ! Revision 1.1 2003/10/09 23:10:45 wasistho
502 ! renamed CalcEddyVis to LesCalcEddyVis
503 !
504 ! Revision 1.5 2003/08/29 01:42:51 wasistho
505 ! Added TARGET attribute to region variable, since pointers are cached into it
506 !
507 ! Revision 1.4 2003/05/15 02:57:06 jblazek
508 ! Inlined index function.
509 !
510 ! Revision 1.3 2002/10/16 07:47:41 wasistho
511 ! Enable Fix Smagorinsky
512 !
513 ! Revision 1.2 2002/10/16 01:59:03 wasistho
514 ! Changed global%error flag
515 !
516 ! Revision 1.1 2002/10/14 23:55:29 wasistho
517 ! Install Rocturb
518 !
519 !
520 !******************************************************************************
521 
522 
523 
524 
525 
526 
527 
**********************************************************************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
FT m(int i, int j) const
j indices k indices k
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 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
subroutine turb_flulesc2f(region)
double sqrt(double d)
Definition: double.h:73
**********************************************************************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
subroutine turb_lescoefdynsmag(region, ibn, ien, ijk)
**********************************************************************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 rflo_getnodeoffset(region, iLev, iNodeOffset, ijNodeOffset)
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
IndexType nfaces() const
Definition: Mesh.H:641
blockLoc i
Definition: read.cpp:79
subroutine rflo_getcelloffset(region, iLev, iCellOffset, ijCellOffset)
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
**********************************************************************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 turb_lescalceddyvis(region, ibn, ien, 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
subroutine turb_lescoefdynmixd(region, ibn, ien, ijk)
**********************************************************************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 rflo_getdimensphys(region, iLev, ipcbeg, ipcend, jpcbeg, jpcend, kpcbeg, kpcend)
subroutine turb_flolesgenc2f(region, ijk)