Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
TURB_ReadTurbSection.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: Read in user input within TURB section (done on all processors).
26 !
27 ! Description: none.
28 !
29 ! Input: regions = user input file of all regions.
30 !
31 ! Output: regions = turbulence input parameters.
32 !
33 ! Notes: Mother routine = ReadInputFile.
34 !
35 !******************************************************************************
36 !
37 ! $Id: TURB_ReadTurbSection.F90,v 1.7 2008/12/06 08:44:42 mtcampbe Exp $
38 !
39 ! Copyright: (c) 2001 by the University of Illinois
40 !
41 !******************************************************************************
42 
43 SUBROUTINE turb_readturbsection( regions )
44 
45  USE moddatatypes
46  USE moddatastruct, ONLY : t_region
47  USE modglobal, ONLY : t_global
48 #ifdef RFLO
50 #endif
51 #ifdef RFLU
52  USE modinterfaces, ONLY : readsection
53 #endif
54  USE modturbulence
55  USE moderror
56  USE modparameters
58 
59  IMPLICIT NONE
60 
61 ! ... parameters
62  TYPE(t_region), POINTER :: regions(:)
63 
64 ! ... local variables
65  CHARACTER(CHRLEN) :: rcsidentstring
66  TYPE(t_global), POINTER :: global
67 
68 
69  INTEGER :: nvals
70  INTEGER :: brbeg, brend
71  INTEGER, PARAMETER :: nvals_max = 25
72 
73  LOGICAL :: defined(nvals_max)
74  REAL(RFREAL) :: vals(nvals_max)
75  CHARACTER(20) :: keys(nvals_max)
76 
77 !******************************************************************************
78 
79  rcsidentstring = '$RCSfile: TURB_ReadTurbSection.F90,v $'
80 
81  global => regions(1)%global
82  CALL registerfunction( global,'TURB_ReadTurbSection',&
83  'TURB_ReadTurbSection.F90' )
84 
85 ! specify keywords and search for them
86 
87  nvals = nvals_max
88 
89  defined(:) = .false.
90 
91  keys( 1) = 'TURBMODEL'
92  keys( 2) = 'OUTPUTNUMBER'
93  keys( 3) = 'CSMAGORINSKY'
94  keys( 4) = 'XTRIPLOC'
95  keys( 5) = 'YTRIPLOC'
96  keys( 6) = 'ZTRIPLOC'
97  keys( 7) = 'FILTERTYPE'
98  keys( 8) = 'DELTATYPE'
99  keys( 9) = 'IFILTERWIDTH'
100  keys(10) = 'JFILTERWIDTH'
101  keys(11) = 'KFILTERWIDTH'
102  keys(12) = 'IHOMOGENDIR'
103  keys(13) = 'JHOMOGENDIR'
104  keys(14) = 'KHOMOGENDIR'
105  keys(15) = 'ENERGYMODEL'
106  keys(16) = 'CALCVORTIC'
107  keys(17) = 'WALLDISTMETHOD'
108  keys(18) = 'WALLDISTFREQ'
109  keys(19) = 'VISCFUNCTION'
110  keys(20) = 'CDES'
111  keys(21) = 'SMOOCF'
112  keys(22) = 'DISCR'
113  keys(23) = 'K2'
114  keys(24) = '1/K4'
115  keys(25) = 'ORDER'
116 
117 #ifdef RFLO
118  CALL readregionsection( global,if_input,nvals,keys(1:nvals),vals(1:nvals), &
119  brbeg,brend,defined(1:nvals) )
120 #endif
121 #ifdef RFLU
122  CALL readsection( global,if_input,nvals,keys(1:nvals),vals(1:nvals), &
123  defined(1:nvals) )
124 
125  brbeg = lbound(regions,1) ! temporary for now before zonal modeling apply
126  brend = ubound(regions,1)
127 #endif
128  IF (defined(1) .eqv. .true.) THEN
129  regions(brbeg:brend)%mixtInput%turbModel = int(vals(1)+0.5_rfreal)
130  ENDIF
131 
132  IF (defined(2) .eqv. .true.) THEN
133  regions(brbeg:brend)%turbInput%nOutField = int(vals(2)+0.5_rfreal)
134  ENDIF
135 
136  IF (defined(3) .eqv. .true.) regions(brbeg:brend)%turbInput%cSmag = abs(vals(3))
137 
138  IF (defined(4) .eqv. .true.) THEN
139  regions(brbeg:brend)%turbInput%xyzSmag(xcoord) = abs(vals(4))
140  ENDIF
141 
142  IF (defined(5) .eqv. .true.) THEN
143  regions(brbeg:brend)%turbInput%xyzSmag(ycoord) = abs(vals(5))
144  ENDIF
145 
146  IF (defined(6) .eqv. .true.) THEN
147  regions(brbeg:brend)%turbInput%xyzSmag(zcoord) = abs(vals(6))
148  ENDIF
149 
150 #ifdef RFLO
151  IF (defined(7) .eqv. .true.) THEN
152  regions(brbeg:brend)%turbInput%filterType = int(vals(7)+0.5_rfreal)
153  ENDIF
154  IF (defined(8) .eqv. .true.) THEN
155  regions(brbeg:brend)%turbInput%deltaType = int(vals(8)+0.5_rfreal)
156  ENDIF
157 #endif
158  IF (defined(9) .eqv. .true.) THEN
159  regions(brbeg:brend)%turbInput%filterWidth(diri) = int(vals(9)+0.5_rfreal)
160  ENDIF
161 #ifdef RFLO
162  IF (defined(10) .eqv. .true.) THEN
163  regions(brbeg:brend)%turbInput%filterWidth(dirj) = int(vals(10)+0.5_rfreal)
164  ENDIF
165  IF (defined(11) .eqv. .true.) THEN
166  regions(brbeg:brend)%turbInput%filterWidth(dirk) = int(vals(11)+0.5_rfreal)
167  ENDIF
168  IF (defined(12) .eqv. .true.) THEN
169  regions(brbeg:brend)%turbInput%homDir(diri) = int(vals(12)+0.5_rfreal)
170  ENDIF
171  IF (defined(13) .eqv. .true.) THEN
172  regions(brbeg:brend)%turbInput%homDir(dirj) = int(vals(13)+0.5_rfreal)
173  ENDIF
174  IF (defined(14) .eqv. .true.) THEN
175  regions(brbeg:brend)%turbInput%homDir(dirk) = int(vals(14)+0.5_rfreal)
176  ENDIF
177 #endif
178  IF (defined(15) .eqv. .true.) THEN
179  regions(brbeg:brend)%turbInput%engModel = int(vals(15)+0.5_rfreal)
180  ENDIF
181  IF (defined(16) .eqv. .true.) THEN
182  regions(brbeg:brend)%turbInput%calcVort = int(vals(16)+0.5_rfreal)
183  ENDIF
184  IF (defined(17) .eqv. .true.) THEN
185  regions(brbeg:brend)%turbInput%wDistMethod = int(vals(17)+0.5_rfreal)
186  ENDIF
187  IF (defined(18) .eqv. .true.) THEN
188  global%turbCalcWDistFreq = max( global%turbCalcWDistFreq, &
189  int(vals(18)+0.5_rfreal) )
190  ENDIF
191  IF (defined(19) .eqv. .true.) THEN
192  regions(brbeg:brend)%turbInput%functV1 = int(vals(19)+0.5_rfreal)
193  ENDIF
194 
195  IF (defined(20) .eqv. .true.) regions(brbeg:brend)%turbInput%cDes = abs(vals(20))
196  IF (defined(21) .eqv. .true.) regions(brbeg:brend)%turbInput%smoocf = vals(21)
197 
198  IF (defined(22) .eqv. .true.) THEN
199  regions(brbeg:brend)%turbInput%spaceDiscr = int(vals(22)+0.5_rfreal)
200  ENDIF
201 
202  IF (defined(23) .eqv. .true.) regions(brbeg:brend)%turbInput%vis2 = abs(vals(23))
203 
204  IF (defined(24) .eqv. .true.) THEN
205  IF (vals(24) > 1.e-10_rfreal) THEN
206  regions(brbeg:brend)%turbInput%vis4 = 1._rfreal/vals(24)
207  ELSEIF (vals(24) > 0._rfreal .AND. vals(24) <= 1.e-10_rfreal) THEN
208  regions(brbeg:brend)%turbInput%vis4 = 1.e+10_rfreal
209  ELSEIF (vals(24) <= 0._rfreal ) THEN
210  regions(brbeg:brend)%turbInput%vis4 = 0.0_rfreal
211  ENDIF
212  ENDIF
213  IF (defined(25) .eqv. .true.) THEN
214  regions(brbeg:brend)%turbInput%spaceOrder = int(vals(25)+0.5_rfreal)
215  ENDIF
216 
217 ! finalize --------------------------------------------------------------------
218 
219  CALL deregisterfunction( global )
220 
221 END SUBROUTINE turb_readturbsection
222 
223 !******************************************************************************
224 !
225 ! RCS Revision history:
226 !
227 ! $Log: TURB_ReadTurbSection.F90,v $
228 ! Revision 1.7 2008/12/06 08:44:42 mtcampbe
229 ! Updated license.
230 !
231 ! Revision 1.6 2008/11/19 22:17:54 mtcampbe
232 ! Added Illinois Open Source License/Copyright
233 !
234 ! Revision 1.5 2008/10/23 18:20:57 mtcampbe
235 ! Crazy number of changes to track and fix initialization and
236 ! restart bugs. Many improperly formed logical expressions
237 ! were fixed, and bug in allocation for data associated with
238 ! the BC_INFLOWVELTEMP boundary condition squashed in
239 ! RFLO_ReadBcInflowVelSection.F90.
240 !
241 ! Revision 1.4 2006/01/12 09:48:15 wasistho
242 ! enabled tripping fixed Smagorinsky
243 !
244 ! Revision 1.3 2004/04/20 20:49:57 wasistho
245 ! added user option for frequency in computing wall distance
246 !
247 ! Revision 1.2 2004/03/19 02:48:11 wasistho
248 ! prepared for RFLU
249 !
250 ! Revision 1.1 2004/03/05 04:37:00 wasistho
251 ! changed nomenclature
252 !
253 ! Revision 1.11 2004/02/19 04:03:34 wasistho
254 ! added new rans/SA parameter VISCFUNCTION
255 !
256 ! Revision 1.10 2004/02/11 03:23:50 wasistho
257 ! added feature: variable number of turbulence output fields
258 !
259 ! Revision 1.9 2003/10/27 23:12:12 wasistho
260 ! bug fixed in reading vis2
261 !
262 ! Revision 1.8 2003/10/26 00:08:26 wasistho
263 ! added multiple discr.types and order
264 !
265 ! Revision 1.7 2003/10/15 03:40:59 wasistho
266 ! added 2nd order dissipation coeff. k2
267 !
268 ! Revision 1.6 2003/10/09 20:49:38 wasistho
269 ! added DES lengthscale coefficient CDES
270 !
271 ! Revision 1.5 2003/10/07 02:07:06 wasistho
272 ! initial installation of RaNS-SA and DES
273 !
274 ! Revision 1.4 2003/08/06 15:55:50 wasistho
275 ! added vorticities computation
276 !
277 ! Revision 1.3 2003/08/01 22:17:52 wasistho
278 ! prepared rocturb for Genx
279 !
280 ! Revision 1.2 2003/07/22 02:59:55 wasistho
281 ! prepare more accurate rocturb restart
282 !
283 ! Revision 1.1 2002/10/14 23:55:30 wasistho
284 ! Install Rocturb
285 !
286 !
287 !******************************************************************************
288 
289 
290 
291 
292 
293 
294 
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 brbeg
subroutine readsection(global, fileID, nvals, keys, vals, defined)
subroutine turb_readturbsection(regions)
subroutine readregionsection(global, fileID, nvals, keys, vals, brbeg, brend, defined)
**********************************************************************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 USE ModDataTypes USE nvals
subroutine deregisterfunction(global)
Definition: ModError.F90:469