Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
TURB_DerivedInputValues.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: Set values derived from user input.
26 !
27 ! Description: Derived variables are mostly components of turbInput data type.
28 !
29 ! Input: regions = Input parameters for all regions.
30 !
31 ! Output: regions = Derived variables stored as turbInput data.
32 !
33 ! Notes: none.
34 !
35 !******************************************************************************
36 !
37 ! $Id: TURB_DerivedInputValues.F90,v 1.11 2008/12/06 08:44:41 mtcampbe Exp $
38 !
39 ! Copyright: (c) 2001 by the University of Illinois
40 !
41 !******************************************************************************
42 
43 SUBROUTINE turb_derivedinputvalues( regions )
44 
45  USE moddatatypes
46  USE moddatastruct, ONLY : t_region
47  USE modglobal, ONLY : t_global
48  USE modturbulence, ONLY : t_turb_input
49  USE moderror
50  USE modmpi
51  USE modparameters
53  IMPLICIT NONE
54 
55 ! ... parameters
56  TYPE(t_region), POINTER :: regions(:)
57 
58 ! ... loop variables
59  INTEGER :: ireg, m
60 
61 ! ... local variables
62  CHARACTER(CHRLEN) :: rcsidentstring
63  TYPE(t_global), POINTER :: global
64  TYPE(t_turb_input), POINTER :: input
65 
66  INTEGER :: flowmodel, turbmodel, ilev, nsv, nst, errorflag
67 #ifdef RFLO
68  REAL(RFREAL) :: one3rd, filtlscale, rndel(3)
69 #endif
70 
71 !******************************************************************************
72 
73  rcsidentstring = '$RCSfile: TURB_DerivedInputValues.F90,v $'
74 
75  global => regions(1)%global
76  CALL registerfunction( global,'TURB_DerivedInputValues',&
77  'TURB_DerivedInputValues.F90' )
78 
79  IF ( global%myProcid == masterproc .AND. &
80  global%verbLevel > verbose_none ) THEN
81  WRITE(stdout,'(A,1X,A)') solver_name,'Entering TURB_DerivedInputValues...'
82  END IF ! global%verbLevel
83 
84 ! set local constants ---------------------------------------------------------
85 
86 #ifdef RFLO
87  one3rd = 1._rfreal/3._rfreal
88 #endif
89 
90 ! region related data (all levels) --------------------------------------------
91 
92 #ifdef RFLO
93  DO ireg = 1,global%nRegions
94 #endif
95 #ifdef RFLU
96  DO ireg = lbound(regions,1),ubound(regions,1)
97 #endif
98 
99  flowmodel = regions(ireg)%mixtInput%flowModel
100  turbmodel = regions(ireg)%mixtInput%turbModel
101  input => regions(ireg)%turbInput
102 
103 ! - number of N-S transport variables
104 
105 #ifdef RFLO
106  DO ilev=1,regions(ireg)%nGridLevels
107  regions(ireg)%levels(ilev)%mixt%nTv = 4
108  ENDDO
109 #endif
110 #ifdef RFLU
111  regions(ireg)%mixtInput%nTv = 4
112 #endif
113 
114 ! - number of dummy cells
115 
116 #ifdef RFLO
117  IF (regions(ireg)%nDumCells < 3) THEN
118  IF ((turbmodel==turb_model_scalsim).OR. &
119  (turbmodel==turb_model_dynsmag)) THEN ! LES with filtering
120  IF ((input%filterWidth(diri) > 1) .OR. &
121  (input%filterWidth(dirj) > 1) .OR. &
122  (input%filterWidth(dirk) > 1)) THEN
123  regions(ireg)%nDumCells = 3
124  ENDIF
125  ELSEIF (turbmodel==turb_model_dynmixd .OR. &
126  turbmodel==turb_model_hdessa) THEN
127  regions(ireg)%nDumCells = 3
128  ELSE
129  ENDIF
130  ENDIF
131 #endif
132 
133 ! - GENERAL TURB --------------------------------------------------------------
134 ! - class of turbulence model
135 
136  IF ((turbmodel==turb_model_fixsmag).OR. &
137  (turbmodel==turb_model_scalsim).OR. &
138  (turbmodel==turb_model_dynsmag).OR. &
139  (turbmodel==turb_model_dynmixd)) THEN
140  input%modelClass = model_les
141  ELSEIF ((turbmodel==turb_model_sa).OR. &
142  (turbmodel==turb_model_dessa).OR. &
143  (turbmodel==turb_model_hdessa)) THEN
144  input%modelClass = model_rans
145  ENDIF
146 
147 ! - number of conservative and derived variables
148 
149  input%nCv = 0
150  input%nDv = 0
151  IF (input%modelClass == model_les) THEN
152  input%nDv = 1
153  ELSEIF (input%modelClass == model_rans) THEN
154  IF ((turbmodel==turb_model_sa).OR. &
155  (turbmodel==turb_model_dessa).OR. &
156  (turbmodel==turb_model_hdessa)) THEN
157  input%nCv = 1
158  ENDIF
159  ENDIF
160 
161 ! - number of zero-one switch fields
162 
163  IF (turbmodel==turb_model_fixsmag .OR. &
164  turbmodel==turb_model_dynsmag .OR. &
165  turbmodel==turb_model_dynmixd) THEN
166  input%nZof = input%nZof + 1
167  global%calcFaceCtr = .true.
168  ENDIF
169 
170 ! - number of permanent and collected variables to be time averaged
171 
172  input%nFixSt = 3
173  nst = 0
174 #ifdef STATS
175  IF ((global%flowType == flow_unsteady) .AND. &
176  (global%doStat == active) .AND. &
177  (global%turbNStat > 0)) THEN
178  nst = 3
179  ENDIF
180 #endif
181  input%nSt = 0 ! = nSt to activate
182 
183 ! - number of stress components
184 
185  nsv = 6
186  input%nSv = 0 ! = nSv to activate
187 
188 ! - number of gradient variables
189 
190  IF ((turbmodel == turb_model_dynsmag).OR. &
191  (turbmodel == turb_model_dynmixd)) THEN
192 #ifdef RFLO
193  input%nGrad = 9
194 #endif
195 #ifdef RFLU
196  input%nGrad = 3
197 #endif
198  ELSEIF ((turbmodel == turb_model_sa).OR. &
199  (turbmodel == turb_model_dessa).OR. &
200  (turbmodel == turb_model_hdessa)) THEN
201 #ifdef RFLO
202  input%nGrad = 3
203 #endif
204 #ifdef RFLU
205  input%nGrad = 1
206 #endif
207  ELSE
208  input%nGrad = 0
209  ENDIF
210 
211 ! - LES SPECIFIC --------------------------------------------------------------
212 ! - filter width scaling factor
213 
214 #ifdef RFLO
215  IF (input%deltaType == deltype_cbrt) THEN
216  filtlscale = 1._rfreal
217  ELSEIF (input%deltaType == deltype_sqrt) THEN
218  filtlscale = 3._rfreal
219  ENDIF
220 
221  rndel(1) = REAL(input%filterWidth(DIRI))
222  rndel(2) = REAL(input%filterWidth(DIRJ))
223  rndel(3) = REAL(input%filterWidth(DIRK))
224  input%delFac2 = (rndel(1)*rndel(2)*rndel(3))**one3rd
225  IF (input%delFac2 < real_small) THEN
226  input%delFac2 = sqrt(rndel(1)*rndel(2))+ &
227  sqrt(rndel(1)*rndel(3))+ &
228  sqrt(rndel(2)*rndel(3))
229  ENDIF
230  IF ((input%filterWidth(diri)/=0.AND.input%filterWidth(dirj)==0 &
231  .AND.input%filterWidth(dirk)==0).OR. &
232  (input%filterWidth(diri)==0.AND.input%filterWidth(dirj)/=0 &
233  .AND.input%filterWidth(dirk)==0).OR. &
234  (input%filterWidth(diri)==0.AND.input%filterWidth(dirj)==0 &
235  .AND.input%filterWidth(dirk)/=0)) THEN
236  input%delFac2 = rndel(1)+rndel(2)+rndel(3)
237  ENDIF
238  input%delFac2 = filtlscale*input%delFac2**2
239 #endif
240 #ifdef RFLU
241  input%delFac2 = REAL(input%filterwidth(diri))
242  input%delFac2 = input%delFac2**2
243 #endif
244 
245  IF (turbmodel==turb_model_fixsmag) THEN
246  input%delFac2 = 1._rfreal
247  ENDIF
248 
249 ! - RANS SPECIFIC -------------------------------------------------------------
250 ! - model constants
251 
252 #ifdef RFLO
253  IF (global%flowType == flow_unsteady) THEN
254  input%smoocf = -1._rfreal
255  ENDIF
256 #endif
257 
258  IF ((turbmodel == turb_model_sa).OR. &
259  (turbmodel == turb_model_dessa).OR. &
260  (turbmodel == turb_model_hdessa)) THEN
261  ALLOCATE( input%const(mc_sa_nelm),stat=errorflag )
262  global%error = errorflag
263  IF (global%error /= 0) CALL errorstop( global,err_allocate,__line__ )
264 
265  input%const(mc_sa_cb1) = 0.1355_rfreal
266  input%const(mc_sa_cb2) = 0.622_rfreal
267  input%const(mc_sa_cw2) = 0.3_rfreal
268  input%const(mc_sa_cw3) = 2.0_rfreal
269  input%const(mc_sa_cv1) = 7.1_rfreal
270  input%const(mc_sa_rsig) = 3.0_rfreal/2.0_rfreal
271  input%const(mc_sa_rkap) = 1.0_rfreal/0.41_rfreal
272  input%const(mc_sa_cw1) = &
273  input%const(mc_sa_cb1)*input%const(mc_sa_rkap)**2 + &
274  (1._rfreal+input%const(mc_sa_cb2))*input%const(mc_sa_rsig)
275  ELSE
276  nullify( input%const )
277  ENDIF
278 
279  ENDDO ! iReg
280 
281 ! finalize --------------------------------------------------------------------
282 
283  IF ( global%myProcid == masterproc .AND. &
284  global%verbLevel > verbose_none ) THEN
285  WRITE(stdout,'(A,1X,A)') solver_name,'Leaving TURB_DerivedInputValues.'
286  END IF ! global%verbLevel
287 
288  CALL deregisterfunction( global )
289 
290 END SUBROUTINE turb_derivedinputvalues
291 
292 !******************************************************************************
293 !
294 ! RCS Revision history:
295 !
296 ! $Log: TURB_DerivedInputValues.F90,v $
297 ! Revision 1.11 2008/12/06 08:44:41 mtcampbe
298 ! Updated license.
299 !
300 ! Revision 1.10 2008/11/19 22:17:53 mtcampbe
301 ! Added Illinois Open Source License/Copyright
302 !
303 ! Revision 1.9 2006/02/04 05:00:05 wasistho
304 ! added enter and leave statements
305 !
306 ! Revision 1.8 2006/01/17 17:25:22 wasistho
307 ! applied tripping to all eddy viscosity models
308 !
309 ! Revision 1.7 2006/01/12 09:49:56 wasistho
310 ! enabled tripping fixed Smagorinsky
311 !
312 ! Revision 1.6 2005/04/18 21:40:02 wasistho
313 ! use nDumCell=3 for hybrid DES
314 !
315 ! Revision 1.5 2005/03/09 06:34:51 wasistho
316 ! incorporated HDESSA
317 !
318 ! Revision 1.4 2004/10/22 23:16:35 wasistho
319 ! set max. number of collected extra statistics, nSt to 3
320 !
321 ! Revision 1.3 2004/08/09 18:38:40 wasistho
322 ! fixed input%delFac2 for RFLU
323 !
324 ! Revision 1.2 2004/03/19 02:45:21 wasistho
325 ! prepared for RFLU
326 !
327 ! Revision 1.1 2004/03/05 04:37:00 wasistho
328 ! changed nomenclature
329 !
330 ! Revision 1.10 2003/10/26 00:24:52 wasistho
331 ! another slightly more accurate model constant
332 !
333 ! Revision 1.9 2003/10/24 20:56:49 wasistho
334 ! corrected SA model constant cb2
335 !
336 ! Revision 1.8 2003/10/07 23:56:05 wasistho
337 ! set RaNS smoocf to -1.0 for unsteady flow
338 !
339 ! Revision 1.6 2003/08/01 22:17:25 wasistho
340 ! prepared rocturb for Genx
341 !
342 ! Revision 1.5 2003/07/22 02:59:09 wasistho
343 ! prepare more accurate rocturb restart
344 !
345 ! Revision 1.4 2003/05/24 02:41:56 wasistho
346 ! sv collector turned off
347 !
348 ! Revision 1.3 2003/05/24 02:10:51 wasistho
349 ! turbulence statistics expanded
350 !
351 ! Revision 1.2 2002/10/15 00:00:00 wasistho
352 ! Removed temporary safety
353 !
354 ! Revision 1.1 2002/10/14 23:55:29 wasistho
355 ! Install Rocturb
356 !
357 !
358 !******************************************************************************
359 
360 
361 
362 
363 
364 
365 
subroutine turb_derivedinputvalues(regions)
FT m(int i, int j) const
subroutine registerfunction(global, funName, fileName)
Definition: ModError.F90:449
double sqrt(double d)
Definition: double.h:73
subroutine input(X, NNODE, NDC, NCELL, NFCE, NBPTS, NBFACE, ITYP, NPROP, XBNDY, XFAR, YFAR, ZFAR)
subroutine errorstop(global, errorCode, errorLine, addMessage)
Definition: ModError.F90:483
subroutine deregisterfunction(global)
Definition: ModError.F90:469