Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
ConvectiveFluxes.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 convective fluxes, add them to dissipation in order
26 ! to obtain the residual.
27 !
28 ! Description: none.
29 !
30 ! Input: region = data of current region.
31 !
32 ! Output: region%levels%mixt%rhs = convective fluxes + dissipation.
33 !
34 ! Notes: none.
35 !
36 ! ******************************************************************************
37 !
38 ! $Id: ConvectiveFluxes.F90,v 1.10 2008/12/06 08:44:08 mtcampbe Exp $
39 !
40 ! Copyright: (c) 2001-2006 by the University of Illinois
41 !
42 ! ******************************************************************************
43 
44 SUBROUTINE convectivefluxes( region )
45 
46  USE moddatatypes
47  USE modglobal, ONLY: t_global
48  USE moddatastruct, ONLY : t_region
49  USE moderror
50  USE modparameters
51 
52 #ifdef RFLU
53  USE rflu_modnscbc
54 
56 #endif
57 
58 
59 #ifdef RFLO
62 #include "Indexing.h"
63 #endif
64 
65  IMPLICIT NONE
66 
67 ! ******************************************************************************
68 ! Definitions and declarations
69 ! ******************************************************************************
70 
71 ! ==============================================================================
72 ! Arguments
73 ! ==============================================================================
74 
75  TYPE(t_region), TARGET :: region
76 
77 ! ==============================================================================
78 ! Locals
79 ! ==============================================================================
80 
81 #ifdef RFLO
82  INTEGER :: ibc, ic, idcbeg, idcend, iec, jdcbeg, jdcend, kdcbeg, kdcend
83  INTEGER :: ilev, icoff, ijcoff
84 #endif
85 #ifdef RFLU
86  INTEGER :: icg
87 #endif
88  INTEGER :: spacediscr, spaceorder
89  REAL(RFREAL), POINTER :: diss(:,:), rhs(:,:)
90  TYPE(t_global), POINTER :: global
91 #ifdef RFLU
92  TYPE(t_region), POINTER :: pregion
93 #endif
94 
95 ! ******************************************************************************
96 ! Start
97 ! ******************************************************************************
98 
99  global => region%global
100 
101  CALL registerfunction( global,'ConvectiveFluxes',&
102  'ConvectiveFluxes.F90' )
103 
104 ! ******************************************************************************
105 ! Get dimensions and pointers
106 ! ******************************************************************************
107 
108 #ifdef RFLO
109  ilev = region%currLevel
110 
111  CALL rflo_getdimensdummy( region,ilev,idcbeg,idcend,jdcbeg,jdcend, &
112  kdcbeg,kdcend )
113  CALL rflo_getcelloffset( region,ilev,icoff,ijcoff )
114  ibc = indijk(idcbeg,jdcbeg,kdcbeg,icoff,ijcoff)
115  iec = indijk(idcend,jdcend,kdcend,icoff,ijcoff)
116 
117  diss => region%levels(ilev)%mixt%diss
118  rhs => region%levels(ilev)%mixt%rhs
119 #endif
120 #ifdef RFLU
121  pregion => region
122 
123  diss => region%mixt%diss
124  rhs => region%mixt%rhs
125 #endif
126 
127  spacediscr = region%mixtInput%spaceDiscr
128  spaceorder = region%mixtInput%spaceOrder
129 
130 #ifdef RFLO
131 ! ******************************************************************************
132 ! Initialize residual (set = dissipation)
133 ! ******************************************************************************
134 
135  DO ic=ibc,iec
136  rhs(cv_mixt_dens,ic) = -diss(cv_mixt_dens,ic)
137  rhs(cv_mixt_xmom,ic) = -diss(cv_mixt_xmom,ic)
138  rhs(cv_mixt_ymom,ic) = -diss(cv_mixt_ymom,ic)
139  rhs(cv_mixt_zmom,ic) = -diss(cv_mixt_zmom,ic)
140  rhs(cv_mixt_ener,ic) = -diss(cv_mixt_ener,ic)
141  ENDDO
142 
143 ! ******************************************************************************
144 ! Compute fluxes
145 ! ******************************************************************************
146 
147 ! ==============================================================================
148 ! 2nd-order central scheme
149 ! ==============================================================================
150 
151  IF (spacediscr==discr_cen_scal .AND. &
152  (spaceorder==discr_order_1 .OR. spaceorder==discr_order_2)) THEN
153  CALL rflo_centralflux( region )
154  ENDIF
155 
156 ! ==============================================================================
157 ! Roe upwind scheme
158 ! ==============================================================================
159 
160  IF (spacediscr == discr_upw_roe) THEN
161  IF (spaceorder == discr_order_1) THEN
162  CALL rflo_roefluxfirst( region )
163  ELSE IF (spaceorder == discr_order_2) THEN
164  CALL rflo_roefluxsecond( region )
165  ENDIF
166  ENDIF
167 #endif
168 
169 #ifdef RFLU
170 ! ******************************************************************************
171 ! Initialize residual (set = dissipation). NOTE only need to do this for
172 ! unsteady flows because for steady flows this is now done by a dedicated
173 ! routine (RFLU_EMS_SetRhs).
174 ! ******************************************************************************
175 
176  IF ( global%flowType == flow_unsteady ) THEN ! Unsteady
177  DO icg = 1,region%grid%nCellsTot
178  rhs(cv_mixt_dens,icg) = -diss(cv_mixt_dens,icg)
179  rhs(cv_mixt_xmom,icg) = -diss(cv_mixt_xmom,icg)
180  rhs(cv_mixt_ymom,icg) = -diss(cv_mixt_ymom,icg)
181  rhs(cv_mixt_zmom,icg) = -diss(cv_mixt_zmom,icg)
182  rhs(cv_mixt_ener,icg) = -diss(cv_mixt_ener,icg)
183  END DO ! icg
184  END IF ! global%flowType
185 
186 ! ******************************************************************************
187 ! Compute fluxes
188 ! ******************************************************************************
189 
190 ! ==============================================================================
191 ! Check state of conserved variable vector
192 ! ==============================================================================
193 
194  IF ( pregion%mixt%cvState /= cv_mixt_state_cons ) THEN
195  CALL errorstop(pregion%global,err_cv_state_invalid,__line__)
196  END IF ! pRegion%mixt%cvState
197 
198 ! ==============================================================================
199 ! Call flux functions
200 ! ==============================================================================
201 
202  CALL rflu_computefluxinv(pregion,flux_part_both)
203 
204 ! ******************************************************************************
205 ! Compute Boundary Patch Rhs
206 ! ******************************************************************************
207 
208  IF ( rflu_nscbc_decidehavenscbc(pregion) .EQV. .true. ) THEN
209  CALL rflu_nscbc_comprhs(pregion)
210  END IF ! RFLU_NSCBC_DecideHaveNSCBC(pRegion)
211 #endif
212 
213 ! ******************************************************************************
214 ! End
215 ! ******************************************************************************
216 
217  CALL deregisterfunction( region%global )
218 
219 END SUBROUTINE convectivefluxes
220 
221 ! ******************************************************************************
222 !
223 ! RCS Revision history:
224 !
225 ! $Log: ConvectiveFluxes.F90,v $
226 ! Revision 1.10 2008/12/06 08:44:08 mtcampbe
227 ! Updated license.
228 !
229 ! Revision 1.9 2008/11/19 22:17:22 mtcampbe
230 ! Added Illinois Open Source License/Copyright
231 !
232 ! Revision 1.8 2006/08/19 15:38:26 mparmar
233 ! Added computing of fluxes for boundary variables
234 !
235 ! Revision 1.7 2006/05/01 20:57:43 haselbac
236 ! Converted to new flux routine
237 !
238 ! Revision 1.6 2006/03/26 20:21:12 haselbac
239 ! Rewrite bcos of of GL model
240 !
241 ! Revision 1.5 2005/11/14 16:53:37 haselbac
242 ! Added AUSM flux computation for pseudo-gas
243 !
244 ! Revision 1.4 2005/11/10 01:55:16 haselbac
245 ! Added new logic for AUSM scheme
246 !
247 ! Revision 1.3 2005/07/14 21:39:28 haselbac
248 ! Added AUSM flux function
249 !
250 ! Revision 1.2 2005/05/16 20:37:56 haselbac
251 ! Changed args to pRegion, rmvd OLES call, changed init of rhs, cosmetics
252 !
253 ! Revision 1.1 2004/12/01 16:48:25 haselbac
254 ! Initial revision after changing case
255 !
256 ! Revision 1.17 2004/02/13 02:56:23 haselbac
257 ! Added calls (commented out for now) to fast Roe routines
258 !
259 ! Revision 1.16 2004/01/29 22:52:38 haselbac
260 ! Rewrote for positive species update
261 !
262 ! Revision 1.15 2003/12/04 03:22:57 haselbac
263 ! Added second-order schemes
264 !
265 ! Revision 1.14 2003/11/20 16:40:35 mdbrandy
266 ! Backing out RocfluidMP changes from 11-17-03
267 !
268 ! Revision 1.11 2003/05/29 17:28:42 jblazek
269 ! Implemented Roe scheme.
270 !
271 ! Revision 1.10 2003/05/16 22:05:40 haselbac
272 ! Added HLLC flux
273 !
274 ! Revision 1.9 2003/05/15 02:57:02 jblazek
275 ! Inlined index function.
276 !
277 ! Revision 1.8 2002/09/09 13:56:51 haselbac
278 ! mixtInput under region, bug fix for iec, added state check
279 !
280 ! Revision 1.7 2002/09/05 17:40:20 jblazek
281 ! Variable global moved into regions().
282 !
283 ! Revision 1.6 2002/07/29 17:13:08 jblazek
284 ! Clean up after RFLU and TURB.
285 !
286 ! Revision 1.5 2002/07/25 14:46:49 haselbac
287 ! Added optimal LES flux routine
288 !
289 ! Revision 1.4 2002/05/04 16:20:24 haselbac
290 ! Added RFLU statements
291 !
292 ! Revision 1.3 2002/02/21 23:25:05 jblazek
293 ! Blocks renamed as regions.
294 !
295 ! Revision 1.2 2002/01/28 23:55:22 jblazek
296 ! Added flux computation (central scheme).
297 !
298 ! Revision 1.1 2002/01/23 03:51:24 jblazek
299 ! Added low-level time-stepping routines.
300 !
301 ! ******************************************************************************
302 
303 
304 
305 
306 
307 
308 
subroutine rflo_roefluxfirst(region)
**********************************************************************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 idcend
NT rhs
subroutine rflu_computefluxinv(pRegion, fluxPart)
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 ic
subroutine, public rflu_nscbc_comprhs(pRegion)
subroutine rflo_getdimensdummy(region, iLev, idcbeg, idcend, jdcbeg, jdcend, kdcbeg, kdcend)
LOGICAL function, public rflu_nscbc_decidehavenscbc(pRegion)
**********************************************************************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 kdcbeg
subroutine rflo_getcelloffset(region, iLev, iCellOffset, ijCellOffset)
**********************************************************************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 idcbeg
subroutine rflo_centralflux(region)
**********************************************************************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 jdcend
subroutine convectivefluxes(region)
**********************************************************************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 jdcbeg
subroutine errorstop(global, errorCode, errorLine, addMessage)
Definition: ModError.F90:483
subroutine deregisterfunction(global)
Definition: ModError.F90:469
subroutine rflo_roefluxsecond(region)