Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
TURB_LesFluxScalSim.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: Obtain viscous fluxes based the Bardina scale similarity model.
26 !
27 ! Description: Get the scale-similarity tauij then add the viscous fluxes to
28 ! the dissipation residual, mixt%diss.
29 ! The ss tauij is obtained by calling LesLij then
30 ! the viscous fluxes by calling VFluxHybrid.
31 !
32 ! Input: region = data of current region
33 ! ibn,ien = begin and end node index
34 !
35 ! Output: tauij, eddy viscosity (mueT, derived from tauij) and viscous fluxes.
36 !
37 ! Notes: tauij in i, j and k faces are kept in fISij, fJSij and fKSij.
38 !
39 !******************************************************************************
40 !
41 ! $Id: TURB_LesFluxScalSim.F90,v 1.9 2008/12/06 08:44:41 mtcampbe Exp $
42 !
43 ! Copyright: (c) 2001 by the University of Illinois
44 !
45 !******************************************************************************
46 
47 SUBROUTINE turb_lesfluxscalsim( region,ibn,ien )
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
59 #endif
61  USE moderror
62  USE modparameters
64  IMPLICIT NONE
65 
66 ! ... parameters
67 #ifdef RFLO
68  TYPE(t_region), TARGET :: region
69 #endif
70 #ifdef RFLU
71  TYPE(t_region), POINTER :: region
72 #endif
73  INTEGER :: ibn, ien
74 
75 ! ... loop variables
76  INTEGER :: ipatch
77 
78 ! ... local variables
79  CHARACTER(CHRLEN) :: rcsidentstring
80  TYPE(t_global), POINTER :: global
81  TYPE(t_turb), POINTER :: turb
82 #ifdef RFLU
83  TYPE(t_patch), POINTER :: patch
84 #endif
85 
86  INTEGER :: errorflag
87 #ifdef RFLO
88  INTEGER :: ilev
89 #endif
90 #ifdef RFLU
91  INTEGER :: npatches, nbfaces
92 #endif
93 
94 !******************************************************************************
95 
96  rcsidentstring = '$RCSfile: TURB_LesFluxScalSim.F90,v $'
97 
98  global => region%global
99  CALL registerfunction( global,'TURB_LesFluxScalSim',&
100  'TURB_LesFluxScalSim.F90' )
101 
102 ! get pointers and parameters ------------------------------------------------
103 #ifdef RFLO
104  ilev = region%currLevel
105  turb => region%levels(ilev)%turb
106 #endif
107 #ifdef RFLU
108  npatches = region%grid%nPatches
109  nbfaces = 0
110 
111  DO ipatch = 1,npatches
112  patch => region%patches(ipatch)
113  nbfaces = nbfaces + patch%nBTris + patch%nBQuads
114  END DO ! iPatch
115  turb => region%turb
116 #endif
117 
118 ! get viscous fluxes from scale similarity contribution; first allocate arrays
119 
120  ALLOCATE( turb%fVar(cv_turb_nelm,ibn:ien),stat=errorflag )
121  global%error = errorflag
122  IF (global%error /= 0) CALL errorstop( global,err_allocate,__line__ )
123 
124  ALLOCATE( turb%ffVar(cv_turb_nelm,ibn:ien),stat=errorflag )
125  global%error = errorflag
126  IF (global%error /= 0) CALL errorstop( global,err_allocate,__line__ )
127 
128  ALLOCATE( turb%fISij(tensor_symm_nelm,ibn:ien),stat=errorflag )
129  global%error = errorflag
130  IF (global%error /= 0) CALL errorstop( global,err_allocate,__line__ )
131 
132 #ifdef RFLO
133  ALLOCATE( turb%fJSij(tensor_symm_nelm,ibn:ien),stat=errorflag )
134  global%error = errorflag
135  IF (global%error /= 0) CALL errorstop( global,err_allocate,__line__ )
136 
137  ALLOCATE( turb%fKSij(tensor_symm_nelm,ibn:ien),stat=errorflag )
138  global%error = errorflag
139  IF (global%error /= 0) CALL errorstop( global,err_allocate,__line__ )
140 #endif
141 #ifdef RFLU
142  ALLOCATE( turb%bfVar(cv_turb_nelm,nbfaces),stat=errorflag )
143  global%error = errorflag
144  IF (global%error /= 0) CALL errorstop( global,err_allocate,__line__ )
145 
146  ALLOCATE( turb%bffVar(cv_turb_nelm,nbfaces),stat=errorflag )
147  global%error = errorflag
148  IF (global%error /= 0) CALL errorstop( global,err_allocate,__line__ )
149 
150  ALLOCATE( turb%bfISij(tensor_symm_nelm,nbfaces),stat=errorflag )
151  global%error = errorflag
152  IF (global%error /= 0) CALL errorstop( global,err_allocate,__line__ )
153 #endif
154 
155 ! Leonard stress lij in I direction
156 #ifdef RFLO
157  CALL turb_flolesgenc2f( region,diri )
158 #endif
159 #ifdef RFLU
160  CALL turb_flulesc2f( region )
161 #endif
162 
163  CALL turb_leslij( region,diri,region%turbInput%filterWidth,turb%fISij )
164 #ifdef RFLU
165  CALL turb_flulesblij( region,region%turbInput%filterWidth,turb%bfISij )
166 #endif
167 
168 #ifdef RFLO
169 ! Leonard stress lij in J direction
170  CALL turb_flolesgenc2f( region,dirj )
171  CALL turb_leslij( region,dirj,region%turbInput%filterWidth,turb%fJSij )
172 
173 ! Leonard stress lij in K direction
174  CALL turb_flolesgenc2f( region,dirk )
175  CALL turb_leslij( region,dirk,region%turbInput%filterWidth,turb%fKSij )
176 #endif
177 
178 #ifdef RFLU
179 ! apply boundary conditions for LES variables
180 ! CALL TURB_FluLesBndConditions( targetFlag )
181 #endif
182 
183 ! get viscous fluxes
184  CALL turb_vfluxhybrid( region )
185 
186 ! deallocate retired arrays
187  DEALLOCATE( turb%fVar,turb%ffVar,turb%fISij )
188 #ifdef RFLO
189  DEALLOCATE( turb%fJSij,turb%fKSij )
190 #endif
191 #ifdef RFLU
192  DEALLOCATE( turb%bfVar,turb%bffVar,turb%bfISij )
193 #endif
194 
195 ! finalize --------------------------------------------------------------------
196 
197  CALL deregisterfunction( global )
198 
199 END SUBROUTINE turb_lesfluxscalsim
200 
201 !******************************************************************************
202 !
203 ! RCS Revision history:
204 !
205 ! $Log: TURB_LesFluxScalSim.F90,v $
206 ! Revision 1.9 2008/12/06 08:44:41 mtcampbe
207 ! Updated license.
208 !
209 ! Revision 1.8 2008/11/19 22:17:54 mtcampbe
210 ! Added Illinois Open Source License/Copyright
211 !
212 ! Revision 1.7 2004/07/30 22:34:45 wasistho
213 ! replaced floLesUniC2F by floLesGenC2F
214 !
215 ! Revision 1.6 2004/05/28 01:59:12 wasistho
216 ! update unstructured grid LES
217 !
218 ! Revision 1.5 2004/03/27 02:16:42 wasistho
219 ! compiled with Rocflu
220 !
221 ! Revision 1.4 2004/03/25 04:40:41 wasistho
222 ! prepared for RFLU
223 !
224 ! Revision 1.3 2004/03/19 02:50:09 wasistho
225 ! prepared for RFLU
226 !
227 ! Revision 1.2 2004/03/12 02:55:35 wasistho
228 ! changed rocturb routine names
229 !
230 ! Revision 1.1 2004/03/05 04:37:00 wasistho
231 ! changed nomenclature
232 !
233 ! Revision 1.3 2003/08/29 01:41:06 wasistho
234 ! Added TARGET attribute to region variable, since pointers are cached into it
235 !
236 ! Revision 1.2 2002/10/16 01:59:49 wasistho
237 ! Changed global%error flag
238 !
239 ! Revision 1.1 2002/10/14 23:55:29 wasistho
240 ! Install Rocturb
241 !
242 !
243 !******************************************************************************
244 
245 
246 
247 
248 
249 
250 
subroutine turb_flulesblij(region, nDel, lij)
subroutine turb_vfluxhybrid(region)
subroutine registerfunction(global, funName, fileName)
Definition: ModError.F90:449
subroutine turb_lesfluxscalsim(region, ibn, ien)
subroutine turb_flulesc2f(region)
Definition: patch.h:74
subroutine errorstop(global, errorCode, errorLine, addMessage)
Definition: ModError.F90:483
subroutine turb_leslij(region, ijk, nDel, lij)
Definition: TURB_LesLij.F90:57
subroutine deregisterfunction(global)
Definition: ModError.F90:469
subroutine turb_flolesgenc2f(region, ijk)