Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
RADI_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 set based on user input parameters.
28 ! Derived variables are components of radiInput data type
29 ! (region%radiInput).
30 !
31 ! Input: regions = input parameters for all regions.
32 !
33 ! Output: regions = derived variables stored as part of radiInput data.
34 !
35 ! Notes: Unlike mixture, derived parameters/variables are not necessarily dv.
36 !
37 !******************************************************************************
38 !
39 ! $Id: RADI_DerivedInputValues.F90,v 1.4 2008/12/06 08:44:37 mtcampbe Exp $
40 !
41 ! Copyright: (c) 2001 by the University of Illinois
42 !
43 !******************************************************************************
44 
45 #ifdef RFLO
46 SUBROUTINE radi_derivedinputvalues( regions )
47 #endif
48 #ifdef RFLU
49 SUBROUTINE radi_derivedinputvalues
50 #endif
51 
52  USE moddatatypes
53 #ifdef RFLO
54  USE moddatastruct, ONLY : t_region
55 #endif
56  USE modglobal, ONLY : t_global
57  USE modradiation, ONLY : t_radi_input
58  USE moderror
59  USE modparameters
61  IMPLICIT NONE
62 
63 #ifdef RFLO
64 ! ... parameters
65  TYPE(t_region), POINTER :: regions(:)
66 
67 ! ... loop variables
68  INTEGER :: ireg, l, m, n
69 #endif
70 
71 ! ... local variables
72  TYPE(t_global), POINTER :: global
73  TYPE(t_radi_input), POINTER :: input
74 
75  INTEGER :: errorflag, angles(3,2)
76  REAL(RFREAL) :: pi, twopi
77 
78 !******************************************************************************
79 
80  global => regions(1)%global
81  CALL registerfunction( global,'RADI_DerivedInputValues',&
82  'RADI_DerivedInputValues.F90' )
83 
84 ! set local constants ---------------------------------------------------------
85 
86  pi = global%pi
87  twopi = 2._rfreal*pi
88 
89 #ifdef RFLO
90 ! global values ---------------------------------------------------------------
91 
92 ! region related data (all levels) --------------------------------------------
93 
94  DO ireg=1,global%nRegions
95 
96  IF (regions(ireg)%procid==global%myProcid .AND. & ! region active and
97  regions(ireg)%active==active) THEN ! on my processor
98 
99  input => regions(ireg)%radiInput
100 
101 ! --- Stefan-Boltzman constant
102  input%stBoltz = 5.67e-8_rfreal
103 
104 ! --- set logical value of radiUsed
105 
106  IF (input%radiModel /= radi_model_none) THEN
107  regions(ireg)%mixtInput%radiUsed = .true.
108  ENDIF
109 
110 ! --- set number of cv, dv and grad components
111 
112  input%nCv = 0
113  input%nGrad = 0
114  IF (input%radiModel == radi_model_fldtran) THEN
115  input%nCv = 1
116  input%nGrad = 3
117  ENDIF
118 
119  input%nDv = 0
120  IF ((input%radiModel /= radi_model_ross) .AND. &
121  (input%radiModel /= radi_model_fldsrc) .AND. &
122  (input%radiModel /= radi_model_fldtran)) THEN
123  input%nDv = 1
124  ENDIF
125 
126 ! --- discrete ordinates and/or intensity angles of RTE models
127 
128  IF ((input%radiModel == radi_model_rtegray) .OR. &
129  (input%radiModel == radi_model_rteband)) THEN
130 
131  IF (input%solMethod == radi_num_dom4) THEN
132  input%nOrdin = 4
133  input%nAng = input%nOrdin
134  ELSEIF (input%solMethod == radi_num_dom8) THEN
135  input%nOrdin = 8
136  input%nAng = input%nOrdin
137  ELSEIF (input%solMethod == radi_num_dom16) THEN
138  input%nOrdin = 16
139  input%nAng = input%nOrdin
140  ELSEIF (input%solMethod == radi_num_fvm) THEN
141  input%nAng = (input%nPol+1)*(input%nAzi+1)
142  ENDIF ! solMethod
143  ENDIF ! radiModel
144 
145 ! --- assign angles read from input file to data structure
146 
147  ALLOCATE( input%angles(input%nAng,radi_angle_ncomp),stat=errorflag )
148  global%error = errorflag
149  IF (global%error /= 0) CALL errorstop( global,err_allocate,__line__ )
150 
151  IF ((input%radiModel == radi_model_ross) .OR. &
152  (input%radiModel == radi_model_fldsrc) .OR. &
153  (input%radiModel == radi_model_fldtran)) THEN
154  READ(input%line(1),*,err=10,end=20) (input%angles(l,radi_angle_polar), &
155  l=1,input%nAng )
156  READ(input%line(2),*,err=10,end=20) (input%angles(l,radi_angle_azimu), &
157  l=1,input%nAng )
158  input%angles = input%angles*global%rad
159 
160  ELSEIF ((input%radiModel == radi_model_rtegray) .OR. &
161  (input%radiModel == radi_model_rteband)) THEN
162  IF (input%solMethod == radi_num_fvm) THEN
163  DO m = 1,input%nPol+1
164  DO l = 1,input%nAzi+1
165  n = (m-1)*(input%nAzi+1) + l
166  input%angles(n,radi_angle_azimu) = pi*dble(l-1)/dble(input%nAzi)
167  input%angles(n,radi_angle_polar) =twopi*dble(m-1)/dble(input%nPol)
168  ENDDO
169  ENDDO
170  ELSE ! DOM
171 
172  ENDIF ! solMethod
173  ENDIF ! radiModel
174 
175  ENDIF ! region active and my processor
176  ENDDO ! iReg
177 
178 #endif
179 
180 ! finalize --------------------------------------------------------------------
181 
182  CALL deregisterfunction( global )
183  goto 999
184 
185 10 CONTINUE
186  CALL errorstop( global,err_file_read,__line__, &
187  'Error in reading real numbers from string' )
188 20 CONTINUE
189  CALL errorstop( global,err_file_read,__line__, &
190  'Number of intensity angles is inconsistent' )
191 
192 999 CONTINUE
193 
194 END SUBROUTINE radi_derivedinputvalues
195 
196 !******************************************************************************
197 !
198 ! RCS Revision history:
199 !
200 ! $Log: RADI_DerivedInputValues.F90,v $
201 ! Revision 1.4 2008/12/06 08:44:37 mtcampbe
202 ! Updated license.
203 !
204 ! Revision 1.3 2008/11/19 22:17:49 mtcampbe
205 ! Added Illinois Open Source License/Copyright
206 !
207 ! Revision 1.2 2004/09/30 17:10:04 wasistho
208 ! prepared for full FLD radiation model
209 !
210 ! Revision 1.1 2004/09/22 02:35:49 wasistho
211 ! changed file nomenclature from lower to upper case
212 !
213 ! Revision 1.6 2004/09/22 01:31:23 wasistho
214 ! switch LFD to FLD for flux limited diffusion
215 !
216 ! Revision 1.5 2004/09/18 17:41:18 wasistho
217 ! install Limited Flux Diffusion radiation
218 !
219 ! Revision 1.4 2003/08/01 22:16:10 wasistho
220 ! prepared rocrad for Genx
221 !
222 ! Revision 1.3 2003/07/23 03:13:49 wasistho
223 ! cured baby illness
224 !
225 ! Revision 1.2 2003/07/22 03:05:41 wasistho
226 ! include logical write-parameter
227 !
228 ! Revision 1.1 2003/07/17 01:16:59 wasistho
229 ! initial activation rocrad
230 !
231 !
232 !
233 !******************************************************************************
234 
235 
236 
237 
238 
239 
240 
FT m(int i, int j) const
subroutine registerfunction(global, funName, fileName)
Definition: ModError.F90:449
static const double pi
Definition: smooth_medial.C:43
subroutine input(X, NNODE, NDC, NCELL, NFCE, NBPTS, NBFACE, ITYP, NPROP, XBNDY, XFAR, YFAR, ZFAR)
const NT & n
subroutine errorstop(global, errorCode, errorLine, addMessage)
Definition: ModError.F90:483
subroutine deregisterfunction(global)
Definition: ModError.F90:469