Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
RFLU_ModExplicitMultiStage.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: Collection of routines related to explicit multistage schemes.
26 !
27 ! Description: None.
28 !
29 ! Notes: None.
30 !
31 ! ******************************************************************************
32 !
33 ! $Id: RFLU_ModExplicitMultiStage.F90,v 1.19 2008/12/06 08:44:21 mtcampbe Exp $
34 !
35 ! Copyright: (c) 2005-2006 by the University of Illinois
36 !
37 ! ******************************************************************************
38 
40 
41  USE moddatatypes
42  USE modglobal, ONLY: t_global
43  USE modparameters
44  USE moderror
45  USE modbndpatch, ONLY: t_patch
46  USE moddatastruct, ONLY: t_region
47  USE modbndpatch, ONLY: t_patch
48  USE modgrid, ONLY: t_grid
49  USE modmixture, ONLY: t_mixt,t_mixt_input
50  USE modmpi
51 
53 
54  IMPLICIT NONE
55 
56  PRIVATE
57  PUBLIC :: rflu_ems_computeresidual, &
60  rflu_ems_setrhs, &
62 
63 ! ******************************************************************************
64 ! Declarations and definitions
65 ! ******************************************************************************
66 
67  CHARACTER(CHRLEN) :: &
68  RCSIdentString = '$RCSfile: RFLU_ModExplicitMultiStage.F90,v $ $Revision: 1.19 $'
69 
70 ! ******************************************************************************
71 ! Routines
72 ! ******************************************************************************
73 
74  CONTAINS
75 
76 
77 
78 
79 
80 
81 
82 
83 ! ******************************************************************************
84 !
85 ! Purpose: Compute residual.
86 !
87 ! Description: None.
88 !
89 ! Input:
90 ! pRegion Pointer to region
91 !
92 ! Output: None.
93 !
94 ! Notes: None.
95 !
96 ! ******************************************************************************
97 
98 SUBROUTINE rflu_ems_computeresidual(pRegion)
99 
104  USE rflu_modlimiters
106  USE rflu_modweno
107 
109 
110  IMPLICIT NONE
111 
112 ! ******************************************************************************
113 ! Declarations and definitions
114 ! ******************************************************************************
115 
116 ! ==============================================================================
117 ! Arguments
118 ! ==============================================================================
119 
120  TYPE(t_region), POINTER :: pregion
121 
122 ! ==============================================================================
123 ! Locals
124 ! ==============================================================================
125 
126  INTEGER :: ipatch
127  INTEGER :: varinfogradcells(cv_mixt_dens:cv_mixt_pres), &
128  varinfogradfaces(cv_mixt_xvel:cv_mixt_temp)
129  TYPE(t_global), POINTER :: global
130  TYPE(t_mixt), POINTER :: pmixt
131  TYPE(t_mixt_input), POINTER :: pmixtinput
132  TYPE(t_patch), POINTER :: ppatch
133 
134 ! ******************************************************************************
135 ! Start
136 ! ******************************************************************************
137 
138  global => pregion%global
139 
140  CALL registerfunction(global,'RFLU_EMS_ComputeResidual',&
141  'RFLU_ModExplicitMultiStage.F90')
142 
143 ! ******************************************************************************
144 ! Set pointers
145 ! ******************************************************************************
146 
147  pmixt => pregion%mixt
148  pmixtinput => pregion%mixtInput
149 
150 ! ******************************************************************************
151 ! Set dissipation array
152 ! ******************************************************************************
153 
154  CALL rflu_ems_setdiss(pregion)
155 
156 ! ******************************************************************************
157 ! Compute cell gradients
158 ! ******************************************************************************
159 
160  IF ( pmixtinput%spaceOrder > discr_order_1 ) THEN
161  CALL rflu_convertcvcons2prim(pregion,cv_mixt_state_duvwp)
162 
163  varinfogradcells(cv_mixt_dens) = v_mixt_dens
164  varinfogradcells(cv_mixt_xvel) = v_mixt_xvel
165  varinfogradcells(cv_mixt_yvel) = v_mixt_yvel
166  varinfogradcells(cv_mixt_zvel) = v_mixt_zvel
167  varinfogradcells(cv_mixt_pres) = v_mixt_pres
168 
169  CALL rflu_computegradcellswrapper(pregion,cv_mixt_dens,cv_mixt_pres, &
170  grc_mixt_dens,grc_mixt_pres, &
171  varinfogradcells,pmixt%cv,pmixt%gradCell)
172 
173  SELECT CASE ( pregion%mixtInput%reconst )
174  CASE ( reconst_none )
175  CASE ( reconst_weno_simple )
176  CALL rflu_wenogradcellswrapper(pregion,grc_mixt_dens,grc_mixt_pres, &
177  pregion%mixt%gradCell)
178  CALL rflu_limitgradcellssimple(pregion,cv_mixt_dens,cv_mixt_pres, &
179  grc_mixt_dens,grc_mixt_pres, &
180  pregion%mixt%cv,pregion%mixt%cvInfo, &
181  pregion%mixt%gradCell)
182  CASE ( reconst_weno_xyz )
183  CALL rflu_wenogradcellsxyzwrapper(pregion,grc_mixt_dens,grc_mixt_pres, &
184  pregion%mixt%gradCell)
185  CALL rflu_limitgradcellssimple(pregion,cv_mixt_dens,cv_mixt_pres, &
186  grc_mixt_dens,grc_mixt_pres, &
187  pregion%mixt%cv,pregion%mixt%cvInfo, &
188  pregion%mixt%gradCell)
189  CASE ( reconst_lim_barthjesp )
190  CALL rflu_createlimiter(pregion,grc_mixt_dens,grc_mixt_pres, &
191  pregion%mixt%lim)
192  CALL rflu_computelimiterbarthjesp(pregion,cv_mixt_dens,cv_mixt_pres, &
193  grc_mixt_dens,grc_mixt_pres, &
194  pregion%mixt%cv, &
195  pregion%mixt%gradCell, &
196  pregion%mixt%lim)
197  CALL rflu_limitgradcells(pregion,grc_mixt_dens,grc_mixt_pres, &
198  pregion%mixt%gradCell,pregion%mixt%lim)
199  CALL rflu_destroylimiter(pregion,pregion%mixt%lim)
200  CASE ( reconst_lim_venkat )
201  CALL rflu_createlimiter(pregion,grc_mixt_dens,grc_mixt_pres, &
202  pregion%mixt%lim)
203  CALL rflu_computelimitervenkat(pregion,cv_mixt_dens,cv_mixt_pres, &
204  grc_mixt_dens,grc_mixt_pres, &
205  pregion%mixt%cv, &
206  pregion%mixt%gradCell, &
207  pregion%mixt%lim)
208  CALL rflu_limitgradcells(pregion,grc_mixt_dens,grc_mixt_pres, &
209  pregion%mixt%gradCell,pregion%mixt%lim)
210  CALL rflu_destroylimiter(pregion,pregion%mixt%lim)
211  CASE default
212  CALL errorstop(global,err_reached_default,__line__)
213  END SELECT ! pRegion%mixtInput%reconst
214 
215  CALL rflu_convertcvprim2cons(pregion,cv_mixt_state_cons)
216  END IF ! pMixtInput%spaceOrder
217 
218 ! ******************************************************************************
219 ! Convective fluxes (dissipative part)
220 ! ******************************************************************************
221 
222  IF ( pmixtinput%ldiss(pregion%irkStep) /=0 ) THEN
223  CALL rflu_computefluxinv(pregion,flux_part_dissip)
224  END IF ! pMixtInput%ldiss
225 
226 ! ******************************************************************************
227 ! Viscous fluxes
228 ! ******************************************************************************
229 
230  IF ( (pmixtinput%flowModel == flow_navst) .AND. &
231  (pmixtinput%ldiss(pregion%irkStep) /=0) ) THEN
232  CALL rflu_convertcvcons2prim(pregion,cv_mixt_state_duvwt)
233 
234  varinfogradfaces(cv_mixt_xvel) = v_mixt_xvel
235  varinfogradfaces(cv_mixt_yvel) = v_mixt_yvel
236  varinfogradfaces(cv_mixt_zvel) = v_mixt_zvel
237  varinfogradfaces(cv_mixt_temp) = v_mixt_temp
238 
239  CALL rflu_computegradfaceswrapper(pregion,cv_mixt_xvel,cv_mixt_temp, &
240  grf_mixt_xvel,grf_mixt_temp, &
241  pregion%mixt%cv, &
242  pregion%mixt%gradFace)
243 
244  IF ( pregion%grid%nFacesConstr > 0 ) THEN
245  CALL rflu_computegradfacesconstr(pregion,cv_mixt_xvel,cv_mixt_temp, &
246  grf_mixt_xvel,grf_mixt_temp, &
247  varinfogradfaces,pregion%mixt%cv, &
248  pregion%mixt%gradFace)
249  END IF ! pRegion%grid%nFacesConstr
250 
251  DO ipatch = 1,pregion%grid%nPatches
252  ppatch => pregion%patches(ipatch)
253 
254  IF ( rflu_decideneedbgradface(pregion,ppatch) .EQV. .true. ) THEN
255  CALL rflu_computegradbfaceswrapper(pregion,ppatch,cv_mixt_xvel, &
256  cv_mixt_temp,grbf_mixt_xvel, &
257  grbf_mixt_temp,pregion%mixt%cv, &
258  ppatch%mixt%gradFace)
259  IF ( ppatch%cReconst /= constr_none ) THEN
260  CALL rflu_computebfgradconstrwrapper(pregion,ppatch,cv_mixt_xvel, &
261  cv_mixt_temp,grbf_mixt_xvel, &
262  grbf_mixt_temp,varinfogradfaces, &
263  pregion%mixt%cv, &
264  ppatch%mixt%gradFace)
265  END IF ! pPatch%cReconst
266  END IF ! RFLU_DecideNeedBGradFace
267  END DO ! iPatch
268 
269  IF ( pmixtinput%turbModel == turb_model_none ) THEN ! Only laminar
270  CALL rflu_enforceheatflux(pregion,pregion%mixt%tv,tv_mixt_tcol)
271  CALL rflu_viscousfluxes(pregion,pregion%mixt%tv,tv_mixt_muel,tv_mixt_tcol)
272  CALL rflu_viscousfluxespatches(pregion,pregion%mixt%tv,tv_mixt_muel, &
273  tv_mixt_tcol)
274  ELSE
275  CALL errorstop(global,err_reached_default,__line__)
276  END IF ! pMixtInput%turbModel
277 
278  CALL rflu_convertcvprim2cons(pregion,cv_mixt_state_cons)
279  END IF ! pMixtInput%flowModel
280 
281 ! ******************************************************************************
282 ! Set right-hand side array
283 ! ******************************************************************************
284 
285  CALL rflu_ems_setrhs(pregion)
286 
287 ! ******************************************************************************
288 ! Convective fluxes (central part)
289 ! ******************************************************************************
290 
291  CALL rflu_computefluxinv(pregion,flux_part_central)
292 
293 ! ******************************************************************************
294 ! End
295 ! ******************************************************************************
296 
297  CALL deregisterfunction(global)
298 
299 END SUBROUTINE rflu_ems_computeresidual
300 
301 
302 
303 
304 
305 
306 
307 
308 
309 
310 
311 ! ******************************************************************************
312 !
313 ! Purpose: Set old values of conserved variables.
314 !
315 ! Description: None.
316 !
317 ! Input:
318 ! pRegion Pointer to region
319 !
320 ! Output: None.
321 !
322 ! Notes: None.
323 !
324 ! ******************************************************************************
325 
326 SUBROUTINE rflu_ems_setcvold(pRegion)
327 
328  IMPLICIT NONE
329 
330 ! ******************************************************************************
331 ! Declarations
332 ! ******************************************************************************
333 
334 ! ==============================================================================
335 ! Arguments
336 ! ==============================================================================
337 
338  TYPE (t_region), POINTER :: pregion
339 
340 ! ==============================================================================
341 ! Locals
342 ! ==============================================================================
343 
344  INTEGER :: icg
345  TYPE(t_global), POINTER :: global
346  TYPE(t_grid), POINTER :: pgrid
347  TYPE(t_mixt), POINTER :: pmixt
348 
349 ! ******************************************************************************
350 ! Start
351 ! ******************************************************************************
352 
353  global => pregion%global
354 
355  CALL registerfunction(global,'RFLU_EMS_SetCvOld',&
356  'RFLU_ModExplicitMultiStage.F90')
357 
358 ! ******************************************************************************
359 ! Set pointers
360 ! ******************************************************************************
361 
362  pgrid => pregion%grid
363  pmixt => pregion%mixt
364 
365 ! ******************************************************************************
366 ! Set old values
367 ! ******************************************************************************
368 
369  DO icg = 1,pgrid%nCellsTot
370  pmixt%cvOld(cv_mixt_dens,icg) = pmixt%cv(cv_mixt_dens,icg)
371  pmixt%cvOld(cv_mixt_xmom,icg) = pmixt%cv(cv_mixt_xmom,icg)
372  pmixt%cvOld(cv_mixt_ymom,icg) = pmixt%cv(cv_mixt_ymom,icg)
373  pmixt%cvOld(cv_mixt_zmom,icg) = pmixt%cv(cv_mixt_zmom,icg)
374  pmixt%cvOld(cv_mixt_ener,icg) = pmixt%cv(cv_mixt_ener,icg)
375  END DO ! icg
376 
377 ! ******************************************************************************
378 ! End
379 ! ******************************************************************************
380 
381  CALL deregisterfunction(global)
382 
383 END SUBROUTINE rflu_ems_setcvold
384 
385 
386 
387 
388 
389 
390 
391 
392 
393 ! ******************************************************************************
394 !
395 ! Purpose: Set dissipation array.
396 !
397 ! Description: None.
398 !
399 ! Input:
400 ! pRegion Pointer to region
401 !
402 ! Output: None.
403 !
404 ! Notes: None.
405 !
406 ! ******************************************************************************
407 
408 SUBROUTINE rflu_ems_setdiss(pRegion)
409 
410  IMPLICIT NONE
411 
412 ! ******************************************************************************
413 ! Declarations
414 ! ******************************************************************************
415 
416 ! ==============================================================================
417 ! Arguments
418 ! ==============================================================================
419 
420  TYPE (t_region), POINTER :: pregion
421 
422 ! ==============================================================================
423 ! Locals
424 ! ==============================================================================
425 
426  INTEGER :: icg
427  REAL(RFREAL) :: term
428  TYPE(t_global), POINTER :: global
429  TYPE(t_grid), POINTER :: pgrid
430  TYPE(t_mixt), POINTER :: pmixt
431  TYPE(t_mixt_input), POINTER :: pmixtinput
432 
433 ! ******************************************************************************
434 ! Start
435 ! ******************************************************************************
436 
437  global => pregion%global
438 
439  CALL registerfunction(global,'RFLU_EMS_SetDiss',&
440  'RFLU_ModExplicitMultiStage.F90')
441 
442 ! ******************************************************************************
443 ! Set pointers
444 ! ******************************************************************************
445 
446  pgrid => pregion%grid
447  pmixt => pregion%mixt
448  pmixtinput => pregion%mixtInput
449 
450 ! ******************************************************************************
451 ! Set old values
452 ! ******************************************************************************
453 
454  IF ( pregion%irkStep == 1 ) THEN
455  DO icg = 1,pgrid%nCellsTot
456  pmixt%diss(cv_mixt_dens,icg) = 0.0_rfreal
457  pmixt%diss(cv_mixt_xmom,icg) = 0.0_rfreal
458  pmixt%diss(cv_mixt_ymom,icg) = 0.0_rfreal
459  pmixt%diss(cv_mixt_zmom,icg) = 0.0_rfreal
460  pmixt%diss(cv_mixt_ener,icg) = 0.0_rfreal
461  END DO ! icg
462  ELSE
463  IF ( pmixtinput%ldiss(pregion%irkStep) /=0 ) THEN
464  term = 1.0_rfreal - pmixtinput%betrk(pregion%irkStep)
465 
466  DO icg = 1,pgrid%nCellsTot
467  pmixt%diss(cv_mixt_dens,icg) = term*pmixt%diss(cv_mixt_dens,icg)
468  pmixt%diss(cv_mixt_xmom,icg) = term*pmixt%diss(cv_mixt_xmom,icg)
469  pmixt%diss(cv_mixt_ymom,icg) = term*pmixt%diss(cv_mixt_ymom,icg)
470  pmixt%diss(cv_mixt_zmom,icg) = term*pmixt%diss(cv_mixt_zmom,icg)
471  pmixt%diss(cv_mixt_ener,icg) = term*pmixt%diss(cv_mixt_ener,icg)
472  END DO ! icg
473  END IF ! pMixtInput%ldiss
474  END IF ! pRegion%irkStep
475 
476 ! ******************************************************************************
477 ! End
478 ! ******************************************************************************
479 
480  CALL deregisterfunction(global)
481 
482 END SUBROUTINE rflu_ems_setdiss
483 
484 
485 
486 
487 
488 
489 
490 
491 
492 ! ******************************************************************************
493 !
494 ! Purpose: Set right-hand side.
495 !
496 ! Description: None.
497 !
498 ! Input:
499 ! pRegion Pointer to region
500 !
501 ! Output: None.
502 !
503 ! Notes: None.
504 !
505 ! ******************************************************************************
506 
507 SUBROUTINE rflu_ems_setrhs(pRegion)
508 
509  IMPLICIT NONE
510 
511 ! ******************************************************************************
512 ! Declarations
513 ! ******************************************************************************
514 
515 ! ==============================================================================
516 ! Arguments
517 ! ==============================================================================
518 
519  TYPE (t_region), POINTER :: pregion
520 
521 ! ==============================================================================
522 ! Locals
523 ! ==============================================================================
524 
525  INTEGER :: icg
526  TYPE(t_global), POINTER :: global
527  TYPE(t_grid), POINTER :: pgrid
528  TYPE(t_mixt), POINTER :: pmixt
529 
530 ! ******************************************************************************
531 ! Start
532 ! ******************************************************************************
533 
534  global => pregion%global
535 
536  CALL registerfunction(global,'RFLU_EMS_SetRhs',&
537  'RFLU_ModExplicitMultiStage.F90')
538 
539 ! ******************************************************************************
540 ! Set pointers
541 ! ******************************************************************************
542 
543  pgrid => pregion%grid
544  pmixt => pregion%mixt
545 
546 ! ******************************************************************************
547 ! Set right-hand side
548 ! ******************************************************************************
549 
550  DO icg = 1,pgrid%nCellsTot
551  pmixt%rhs(cv_mixt_dens,icg) = -pmixt%diss(cv_mixt_dens,icg)
552  pmixt%rhs(cv_mixt_xmom,icg) = -pmixt%diss(cv_mixt_xmom,icg)
553  pmixt%rhs(cv_mixt_ymom,icg) = -pmixt%diss(cv_mixt_ymom,icg)
554  pmixt%rhs(cv_mixt_zmom,icg) = -pmixt%diss(cv_mixt_zmom,icg)
555  pmixt%rhs(cv_mixt_ener,icg) = -pmixt%diss(cv_mixt_ener,icg)
556  END DO ! icg
557 
558 ! ******************************************************************************
559 ! End
560 ! ******************************************************************************
561 
562  CALL deregisterfunction(global)
563 
564 END SUBROUTINE rflu_ems_setrhs
565 
566 
567 
568 
569 
570 
571 
572 
573 ! ******************************************************************************
574 !
575 ! Purpose: Update conserved variables.
576 !
577 ! Description: None.
578 !
579 ! Input:
580 ! pRegion Pointer to region
581 !
582 ! Output: None.
583 !
584 ! Notes: None.
585 !
586 ! ******************************************************************************
587 
588 SUBROUTINE rflu_ems_updateconservedvars(pRegion)
589 
590  IMPLICIT NONE
591 
592 ! ******************************************************************************
593 ! Declarations
594 ! ******************************************************************************
595 
596 ! ==============================================================================
597 ! Arguments
598 ! ==============================================================================
599 
600  TYPE (t_region), POINTER :: pregion
601 
602 ! ==============================================================================
603 ! Locals
604 ! ==============================================================================
605 
606  INTEGER :: icg
607  REAL(RFREAL) :: term
608  TYPE(t_global), POINTER :: global
609  TYPE(t_grid), POINTER :: pgrid
610  TYPE(t_mixt), POINTER :: pmixt
611  TYPE(t_mixt_input), POINTER :: pmixtinput
612 
613 ! ******************************************************************************
614 ! Start
615 ! ******************************************************************************
616 
617  global => pregion%global
618 
619  CALL registerfunction(global,'RFLU_EMS_UpdateConservedVars',&
620  'RFLU_ModExplicitMultiStage.F90')
621 
622 ! ******************************************************************************
623 ! Set pointers and values
624 ! ******************************************************************************
625 
626  pgrid => pregion%grid
627  pmixt => pregion%mixt
628  pmixtinput => pregion%mixtInput
629 
630 ! ******************************************************************************
631 ! Set old values
632 ! ******************************************************************************
633 
634  DO icg = 1,pgrid%nCells
635  term = pmixtinput%ark(pregion%irkStep)*pmixtinput%cfl*pregion%dt(icg)/pgrid%vol(icg)
636 
637  pmixt%cv(cv_mixt_dens,icg) = pmixt%cvOld(cv_mixt_dens,icg) - term*pmixt%rhs(cv_mixt_dens,icg)
638  pmixt%cv(cv_mixt_xmom,icg) = pmixt%cvOld(cv_mixt_xmom,icg) - term*pmixt%rhs(cv_mixt_xmom,icg)
639  pmixt%cv(cv_mixt_ymom,icg) = pmixt%cvOld(cv_mixt_ymom,icg) - term*pmixt%rhs(cv_mixt_ymom,icg)
640  pmixt%cv(cv_mixt_zmom,icg) = pmixt%cvOld(cv_mixt_zmom,icg) - term*pmixt%rhs(cv_mixt_zmom,icg)
641  pmixt%cv(cv_mixt_ener,icg) = pmixt%cvOld(cv_mixt_ener,icg) - term*pmixt%rhs(cv_mixt_ener,icg)
642  END DO ! icg
643 
644 ! ******************************************************************************
645 ! End
646 ! ******************************************************************************
647 
648  CALL deregisterfunction(global)
649 
650 END SUBROUTINE rflu_ems_updateconservedvars
651 
652 
653 
654 
655 
656 
657 
658 
659 ! ******************************************************************************
660 ! End
661 ! ******************************************************************************
662 
664 
665 
666 ! ******************************************************************************
667 !
668 ! RCS Revision history:
669 !
670 ! $Log: RFLU_ModExplicitMultiStage.F90,v $
671 ! Revision 1.19 2008/12/06 08:44:21 mtcampbe
672 ! Updated license.
673 !
674 ! Revision 1.18 2008/11/19 22:17:32 mtcampbe
675 ! Added Illinois Open Source License/Copyright
676 !
677 ! Revision 1.17 2006/08/19 15:39:07 mparmar
678 ! Moved region%mixt%bGradFace to patch%mixt%gradFace
679 !
680 ! Revision 1.16 2006/05/01 21:01:19 haselbac
681 ! Converted to new flux routine
682 !
683 ! Revision 1.15 2006/04/27 15:10:32 haselbac
684 ! Added VENKAT limiter option
685 !
686 ! Revision 1.14 2006/04/15 16:59:24 haselbac
687 ! Added RECONST_NONE and cReconst IF for bf grads
688 !
689 ! Revision 1.13 2006/04/07 15:19:19 haselbac
690 ! Removed tabs
691 !
692 ! Revision 1.12 2006/04/07 14:49:41 haselbac
693 ! Adapted to changes in f and bf diff modules and new WENO module
694 !
695 ! Revision 1.11 2006/01/06 22:12:44 haselbac
696 ! Added call to cell grad wrapper routine
697 !
698 ! Revision 1.10 2005/11/14 16:59:29 haselbac
699 ! Added support for pseudo-gas model
700 !
701 ! Revision 1.9 2005/11/10 02:28:49 haselbac
702 ! Adapted to changes in AUSM flux module
703 !
704 ! Revision 1.8 2005/10/27 19:18:35 haselbac
705 ! Adapted to changes in gradient routines
706 !
707 ! Revision 1.7 2005/10/14 14:06:07 haselbac
708 ! Added call to RFLU_EnforceHeatFlux
709 !
710 ! Revision 1.6 2005/10/05 13:57:06 haselbac
711 ! Adapted to changes in gradient computation
712 !
713 ! Revision 1.5 2005/07/14 21:43:57 haselbac
714 ! Added limiters, new ENO scheme, and AUSM flux function
715 !
716 ! Revision 1.4 2005/07/11 19:32:26 mparmar
717 ! Adapted call to RFLU_LimitGradCells to changes in modules
718 !
719 ! Revision 1.3 2005/07/07 22:45:01 haselbac
720 ! Added profiling calls
721 !
722 ! Revision 1.2 2005/05/19 18:18:44 haselbac
723 ! Cosmetics only
724 !
725 ! Revision 1.1 2005/05/16 20:36:29 haselbac
726 ! Initial revision
727 !
728 ! ******************************************************************************
729 
730 
731 
732 
733 
734 
735 
736 
737 
738 
739 
740 
subroutine, public rflu_destroylimiter(pRegion, lim)
subroutine, public rflu_computegradfaceswrapper(pRegion, iBegVar, iEndVar, iBegGrad, iEndGrad, var, grad)
subroutine, public rflu_computelimitervenkat(pRegion, iBegVar, iEndVar, iBegGrad, iEndGrad, var, grad, lim)
subroutine, public rflu_computelimiterbarthjesp(pRegion, iBegVar, iEndVar, iBegGrad, iEndGrad, var, grad, lim)
subroutine, public rflu_createlimiter(pRegion, iBegGrad, iEndGrad, lim)
subroutine, public rflu_viscousfluxespatches(pRegion, tv, tvIndxVisc, tvIndxCond)
LOGICAL function rflu_decideneedbgradface(pRegion, pPatch)
subroutine, public rflu_computebfgradconstrwrapper(pRegion, pPatch, iBegVar, iEndVar, iBegGrad, iEndGrad, varInfo, var, grad)
subroutine rflu_computefluxinv(pRegion, fluxPart)
subroutine registerfunction(global, funName, fileName)
Definition: ModError.F90:449
subroutine, public rflu_enforceheatflux(pRegion, tv, tvIndxCond)
subroutine, public rflu_ems_computeresidual(pRegion)
subroutine, public rflu_computegradfacesconstr(pRegion, iBegVar, iEndVar, iBegGrad, iEndGrad, varInfo, var, grad)
subroutine, public rflu_limitgradcells(pRegion, iBegGrad, iEndGrad, grad, lim)
subroutine, public rflu_limitgradcellssimple(pRegion, iBegVar, iEndVar, iBegGrad, iEndGrad, var, varInfo, grad)
subroutine, public rflu_computegradbfaceswrapper(pRegion, pPatch, iBegVar, iEndVar, iBegGrad, iEndGrad, var, grad)
subroutine, public rflu_ems_setcvold(pRegion)
subroutine, public rflu_convertcvcons2prim(pRegion, cvStateFuture)
subroutine, public rflu_ems_setdiss(pRegion)
subroutine, public rflu_convertcvprim2cons(pRegion, cvStateFuture)
subroutine, public rflu_wenogradcellswrapper(pRegion, iBegGrad, iEndGrad, grad)
subroutine, public rflu_ems_setrhs(pRegion)
subroutine, public rflu_ems_updateconservedvars(pRegion)
subroutine, public rflu_wenogradcellsxyzwrapper(pRegion, iBegGrad, iEndGrad, grad)
subroutine, public rflu_computegradcellswrapper(pRegion, iBegVar, iEndVar, iBegGrad, iEndGrad, varInfo, var, grad)
subroutine errorstop(global, errorCode, errorLine, addMessage)
Definition: ModError.F90:483
subroutine deregisterfunction(global)
Definition: ModError.F90:469
subroutine, public rflu_viscousfluxes(pRegion, tv, tvIndxVisc, tvIndxCond)