Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
RFLU_ModResidual.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 residual computation.
26 !
27 ! Description: None.
28 !
29 ! Notes: None.
30 !
31 ! ******************************************************************************
32 !
33 ! $Id: RFLU_ModResidual.F90,v 1.18 2008/12/06 08:44:24 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 modgrid, ONLY: t_grid
48  USE modbndpatch, ONLY: t_patch
49  USE modmixture, ONLY: t_mixt,t_mixt_input
50  USE modmpi
51 
53 
54  IMPLICIT NONE
55 
56  PRIVATE
57  PUBLIC :: rflu_res_computeresidual, &
63 
64 ! ******************************************************************************
65 ! Declarations and definitions
66 ! ******************************************************************************
67 
68  CHARACTER(CHRLEN) :: &
69  RCSIdentString = '$RCSfile: RFLU_ModResidual.F90,v $ $Revision: 1.18 $'
70 
71 ! ******************************************************************************
72 ! Routines
73 ! ******************************************************************************
74 
75  CONTAINS
76 
77 
78 
79 
80 
81 
82 
83 
84 ! ******************************************************************************
85 !
86 ! Purpose: Compute residual.
87 !
88 ! Description: None.
89 !
90 ! Input:
91 ! pRegion Pointer to region
92 !
93 ! Output: None.
94 !
95 ! Notes: None.
96 !
97 ! ******************************************************************************
98 
99 SUBROUTINE rflu_res_computeresidual(pRegion)
100 
105  USE rflu_modlimiters
107  USE rflu_modweno
108 
110 
111  IMPLICIT NONE
112 
113 ! ******************************************************************************
114 ! Declarations and definitions
115 ! ******************************************************************************
116 
117 ! ==============================================================================
118 ! Arguments
119 ! ==============================================================================
120 
121  TYPE(t_region), POINTER :: pregion
122 
123 ! ==============================================================================
124 ! Locals
125 ! ==============================================================================
126 
127  INTEGER :: ipatch
128  INTEGER :: varinfogradcells(cv_mixt_dens:cv_mixt_pres), &
129  varinfogradfaces(cv_mixt_xvel:cv_mixt_temp)
130  TYPE(t_global), POINTER :: global
131  TYPE(t_mixt), POINTER :: pmixt
132  TYPE(t_mixt_input), POINTER :: pmixtinput
133  TYPE(t_patch), POINTER :: ppatch
134 
135 ! ******************************************************************************
136 ! Start
137 ! ******************************************************************************
138 
139  global => pregion%global
140 
141  CALL registerfunction(global,'RFLU_RES_ComputeResidual',&
142  'RFLU_ModResidual.F90')
143 
144 ! ******************************************************************************
145 ! Set pointers
146 ! ******************************************************************************
147 
148  pmixt => pregion%mixt
149  pmixtinput => pregion%mixtInput
150 
151 ! ******************************************************************************
152 ! Set dissipation array
153 ! ******************************************************************************
154 
155  CALL rflu_res_setdiss(pregion)
156 
157 ! ******************************************************************************
158 ! Compute cell gradients. NOTE need to be done here so that back-up viscous
159 ! flux computation can make use of them.
160 ! ******************************************************************************
161 
162  IF ( pmixtinput%spaceOrder > discr_order_1 ) THEN
163  CALL rflu_convertcvcons2prim(pregion,cv_mixt_state_duvwp)
164 
165  varinfogradcells(cv_mixt_dens) = v_mixt_dens
166  varinfogradcells(cv_mixt_xvel) = v_mixt_xvel
167  varinfogradcells(cv_mixt_yvel) = v_mixt_yvel
168  varinfogradcells(cv_mixt_zvel) = v_mixt_zvel
169  varinfogradcells(cv_mixt_pres) = v_mixt_pres
170 
171  CALL rflu_computegradcellswrapper(pregion,cv_mixt_dens,cv_mixt_pres, &
172  grc_mixt_dens,grc_mixt_pres, &
173  varinfogradcells,pmixt%cv,pmixt%gradCell)
174 
175  SELECT CASE ( pregion%mixtInput%reconst )
176  CASE ( reconst_none )
177  CASE ( reconst_weno_simple )
178  CALL rflu_wenogradcellswrapper(pregion,grc_mixt_dens,grc_mixt_pres, &
179  pregion%mixt%gradCell)
180  CALL rflu_limitgradcellssimple(pregion,cv_mixt_dens,cv_mixt_pres, &
181  grc_mixt_dens,grc_mixt_pres, &
182  pregion%mixt%cv,pregion%mixt%cvInfo, &
183  pregion%mixt%gradCell)
184  CASE ( reconst_weno_xyz )
185  CALL rflu_wenogradcellsxyzwrapper(pregion,grc_mixt_dens,grc_mixt_pres, &
186  pregion%mixt%gradCell)
187  CALL rflu_limitgradcellssimple(pregion,cv_mixt_dens,cv_mixt_pres, &
188  grc_mixt_dens,grc_mixt_pres, &
189  pregion%mixt%cv,pregion%mixt%cvInfo, &
190  pregion%mixt%gradCell)
191  CASE ( reconst_lim_barthjesp )
192  CALL rflu_createlimiter(pregion,grc_mixt_dens,grc_mixt_pres, &
193  pregion%mixt%lim)
194  CALL rflu_computelimiterbarthjesp(pregion,cv_mixt_dens,cv_mixt_pres, &
195  grc_mixt_dens,grc_mixt_pres, &
196  pregion%mixt%cv, &
197  pregion%mixt%gradCell, &
198  pregion%mixt%lim)
199  CALL rflu_limitgradcells(pregion,grc_mixt_dens,grc_mixt_pres, &
200  pregion%mixt%gradCell,pregion%mixt%lim)
201  CALL rflu_destroylimiter(pregion,pregion%mixt%lim)
202  CASE ( reconst_lim_venkat )
203  CALL rflu_createlimiter(pregion,grc_mixt_dens,grc_mixt_pres, &
204  pregion%mixt%lim)
205  CALL rflu_computelimitervenkat(pregion,cv_mixt_dens,cv_mixt_pres, &
206  grc_mixt_dens,grc_mixt_pres, &
207  pregion%mixt%cv,pregion%mixt%gradCell, &
208  pregion%mixt%lim)
209  CALL rflu_limitgradcells(pregion,grc_mixt_dens,grc_mixt_pres, &
210  pregion%mixt%gradCell,pregion%mixt%lim)
211  CALL rflu_destroylimiter(pregion,pregion%mixt%lim)
212  CASE default
213  CALL errorstop(global,err_reached_default,__line__)
214  END SELECT ! pRegion%mixtInput%reconst
215 
216  CALL rflu_convertcvprim2cons(pregion,cv_mixt_state_cons)
217  END IF ! pMixtInput%spaceOrder
218 
219 ! ******************************************************************************
220 ! Viscous fluxes
221 ! ******************************************************************************
222 
223  IF ( pmixtinput%flowModel == flow_navst ) THEN
224  CALL rflu_convertcvcons2prim(pregion,cv_mixt_state_duvwt)
225 
226  varinfogradfaces(cv_mixt_xvel) = v_mixt_xvel
227  varinfogradfaces(cv_mixt_yvel) = v_mixt_yvel
228  varinfogradfaces(cv_mixt_zvel) = v_mixt_zvel
229  varinfogradfaces(cv_mixt_temp) = v_mixt_temp
230 
231  CALL rflu_computegradfaceswrapper(pregion,cv_mixt_xvel,cv_mixt_temp, &
232  grf_mixt_xvel,grf_mixt_temp, &
233  pregion%mixt%cv,pregion%mixt%gradFace)
234 
235  IF ( pregion%grid%nFacesConstr > 0 ) THEN
236  CALL rflu_computegradfacesconstr(pregion,cv_mixt_xvel,cv_mixt_temp, &
237  grf_mixt_xvel,grf_mixt_temp, &
238  varinfogradfaces,pregion%mixt%cv, &
239  pregion%mixt%gradFace)
240  END IF ! pRegion%grid%nFacesConstr
241 
242  DO ipatch = 1,pregion%grid%nPatches
243  ppatch => pregion%patches(ipatch)
244 
245  IF ( rflu_decideneedbgradface(pregion,ppatch) .EQV. .true. ) THEN
246  CALL rflu_computegradbfaceswrapper(pregion,ppatch,cv_mixt_xvel, &
247  cv_mixt_temp,grbf_mixt_xvel, &
248  grbf_mixt_temp,pregion%mixt%cv, &
249  ppatch%mixt%gradFace)
250  IF ( ppatch%cReconst /= constr_none ) THEN
251  CALL rflu_computebfgradconstrwrapper(pregion,ppatch,cv_mixt_xvel, &
252  cv_mixt_temp,grbf_mixt_xvel, &
253  grbf_mixt_temp,varinfogradfaces, &
254  pregion%mixt%cv, &
255  ppatch%mixt%gradFace)
256  END IF ! pPatch%cReconst
257  END IF ! RFLU_DecideNeedBGradFace
258  END DO ! iPatch
259 
260  IF ( pmixtinput%turbModel == turb_model_none ) THEN ! Only laminar
261  CALL rflu_enforceheatflux(pregion,pregion%mixt%tv,tv_mixt_tcol)
262  CALL rflu_viscousfluxes(pregion,pregion%mixt%tv,tv_mixt_muel,tv_mixt_tcol)
263  CALL rflu_viscousfluxespatches(pregion,pregion%mixt%tv,tv_mixt_muel, &
264  tv_mixt_tcol)
265  ELSE
266  CALL errorstop(global,err_reached_default,__line__)
267  END IF ! pMixtInput%turbModel
268 
269  CALL rflu_convertcvprim2cons(pregion,cv_mixt_state_cons)
270  END IF ! pMixtInput%flowModel
271 
272 ! ******************************************************************************
273 ! Set right-hand side array
274 ! ******************************************************************************
275 
276  CALL rflu_res_setrhs(pregion)
277 
278 ! ******************************************************************************
279 ! Inviscid fluxes
280 ! ******************************************************************************
281 
282  CALL rflu_computefluxinv(pregion,flux_part_both)
283 
284 ! ******************************************************************************
285 ! End
286 ! ******************************************************************************
287 
288  CALL deregisterfunction(global)
289 
290 END SUBROUTINE rflu_res_computeresidual
291 
292 
293 
294 
295 
296 
297 
298 
299 
300 
301 
302 ! ******************************************************************************
303 !
304 ! Purpose: Determine support of residual for first-order scheme
305 !
306 ! Description: None.
307 !
308 ! Input:
309 ! pRegion Pointer to region data
310 ! icg Cell index
311 ! rsSizeMax Maximum allowed size of residual support
312 !
313 ! Output:
314 ! rs Residual support
315 ! rsSize Actual size of residual support
316 !
317 ! Notes:
318 ! 1. At present, only take into account cell-to-cell stencil.
319 !
320 ! ******************************************************************************
321 
322 SUBROUTINE rflu_getresidualsupport1(pRegion,icg,rs,rsSizeMax,rsSize)
323 
324  USE moddatatypes
325  USE modparameters
326  USE moddatastruct, ONLY: t_region
327  USE modglobal, ONLY: t_global
328  USE modgrid, ONLY: t_grid
329  USE modbndpatch, ONLY: t_patch
330  USE moderror
331 
332  IMPLICIT NONE
333 
334 ! ******************************************************************************
335 ! Declarations and definitions
336 ! ******************************************************************************
337 
338 ! ==============================================================================
339 ! Arguments
340 ! ==============================================================================
341 
342  INTEGER, INTENT(IN) :: icg
343  INTEGER, INTENT(INOUT) :: rssizemax
344  INTEGER, INTENT(OUT) :: rssize
345  INTEGER, DIMENSION(:) :: rs
346  TYPE(t_region), POINTER :: pregion
347 
348 ! ==============================================================================
349 ! Locals
350 ! ==============================================================================
351 
352  INTEGER :: c1,c2,errorflag,icl,ict,ifg,ifl,ipatch,nfaces
353  INTEGER, DIMENSION(:,:), POINTER :: pc2f
354  TYPE(t_global), POINTER :: global
355  TYPE(t_grid), POINTER :: pgrid
356  TYPE(t_patch), POINTER :: ppatch
357 
358 ! ******************************************************************************
359 ! Start
360 ! ******************************************************************************
361 
362  global => pregion%global
363 
364  CALL registerfunction(global,'RFLU_GetResidualSupport1',&
365  'RFLU_ModResidual.F90')
366 
367  pgrid => pregion%grid
368 
369 ! ******************************************************************************
370 ! Set up cell-to-face pointer
371 ! ******************************************************************************
372 
373  ict = pgrid%cellGlob2Loc(1,icg)
374  icl = pgrid%cellGlob2Loc(2,icg)
375 
376  SELECT CASE ( ict )
377  CASE ( cell_type_tet )
378  pc2f => pgrid%tet2f(:,:,icl)
379  CASE ( cell_type_hex )
380  pc2f => pgrid%hex2f(:,:,icl)
381  CASE ( cell_type_pri )
382  pc2f => pgrid%pri2f(:,:,icl)
383  CASE ( cell_type_pyr )
384  pc2f => pgrid%pyr2f(:,:,icl)
385  CASE default
386  CALL errorstop(global,err_reached_default,__line__)
387  END SELECT ! ict
388 
389  nfaces = SIZE(pc2f,2)
390 
391 ! ******************************************************************************
392 ! Initialize
393 ! ******************************************************************************
394 
395  rssize = 1
396 
397  rs(rssize) = icg
398 
399 ! ******************************************************************************
400 ! Add face neighbors
401 ! ******************************************************************************
402 
403  DO ifl = 1,nfaces
404  ipatch = pc2f(1,ifl)
405  ifg = pc2f(2,ifl)
406 
407  IF ( ipatch == 0 ) THEN ! Interior face
408  c1 = pgrid%f2c(1,ifg)
409  c2 = pgrid%f2c(2,ifg)
410 
411  IF ( c1 == icg ) THEN
412  IF ( rssize < rssizemax ) THEN
413  rssize = rssize + 1
414  ELSE
415 ! TEMPORARY
416  WRITE(*,*) 'ERROR - About to exceed size of array rs!'
417  stop
418 ! END TEMPORARY
419  END IF ! rsSize
420 
421  rs(rssize) = c2
422  ELSE IF ( c2 == icg ) THEN
423  IF ( rssize < rssizemax ) THEN
424  rssize = rssize + 1
425  ELSE
426 ! TEMPORARY
427  WRITE(*,*) 'ERROR - About to exceed size of array rs!'
428  stop
429 ! END TEMPORARY
430  END IF ! rsSize
431 
432  rs(rssize) = c1
433  ELSE ! defensive programming
434  CALL errorstop(global,err_reached_default,__line__)
435  END IF ! c1
436  END IF ! iPatch
437  END DO ! ifl
438 
439 ! ******************************************************************************
440 ! End
441 ! ******************************************************************************
442 
443  CALL deregisterfunction(global)
444 
445 END SUBROUTINE rflu_getresidualsupport1
446 
447 
448 
449 
450 
451 
452 
453 
454 
455 ! ******************************************************************************
456 !
457 ! Purpose: Determine support of residual for second-order scheme
458 !
459 ! Description: None.
460 !
461 ! Input:
462 ! pRegion Pointer to region data
463 ! icg Cell index
464 ! rsSizeMax Maximum allowed size of residual support
465 !
466 ! Output:
467 ! rs Residual support
468 ! rsSize Actual size of residual support
469 !
470 ! Notes:
471 ! 1. At present, only take into account cell-to-cell stencil, i.e., only
472 ! inviscid flux contribution. This is because cell-to-cell stencil is
473 ! typically larger than face-to-cell stenci, at least for WENO scheme.
474 !
475 ! ******************************************************************************
476 
477 SUBROUTINE rflu_getresidualsupport2(pRegion,icg,rs,rsSizeMax,rsSize)
478 
479  USE moddatatypes
480  USE modparameters
481  USE moddatastruct, ONLY: t_region
482  USE modglobal, ONLY: t_global
483  USE modgrid, ONLY: t_grid
484  USE moderror
485 
486  USE modsortsearch
487 
488  IMPLICIT NONE
489 
490 ! ******************************************************************************
491 ! Declarations and definitions
492 ! ******************************************************************************
493 
494 ! ==============================================================================
495 ! Arguments
496 ! ==============================================================================
497 
498  INTEGER, INTENT(IN) :: icg
499  INTEGER, INTENT(INOUT) :: rssizemax
500  INTEGER, INTENT(OUT) :: rssize
501  INTEGER, DIMENSION(:) :: rs
502  TYPE(t_region), POINTER :: pregion
503 
504 ! ==============================================================================
505 ! Locals
506 ! ==============================================================================
507 
508  INTEGER :: errorflag,icg2,icg3,ilevel,iloc,irs,isl,nlevels,rssavesize, &
509  rssavesizemax,rstempsize,rstempsizemax
510  INTEGER, DIMENSION(:), ALLOCATABLE :: rssave,rstemp
511  TYPE(t_global), POINTER :: global
512  TYPE(t_grid), POINTER :: pgrid
513 
514 ! ******************************************************************************
515 ! Start
516 ! ******************************************************************************
517 
518  global => pregion%global
519 
520  CALL registerfunction(global,'RFLU_GetResidualSupport2',&
521  'RFLU_ModResidual.F90')
522 
523  pgrid => pregion%grid
524 
525 ! ******************************************************************************
526 ! Allocate temporary memory
527 ! ******************************************************************************
528 
529  rssavesizemax = 1000
530  rstempsizemax = 1000 ! Must be equal to or larger than <rsSaveSizeMax>
531 
532  ALLOCATE(rssave(rssavesizemax),stat=errorflag)
533  global%error = errorflag
534  IF ( global%error /= err_none ) THEN
535  CALL errorstop(global,err_allocate,__line__,'rsSave')
536  END IF ! global%errorFlag
537 
538  ALLOCATE(rstemp(rstempsizemax),stat=errorflag)
539  global%error = errorflag
540  IF ( global%error /= err_none ) THEN
541  CALL errorstop(global,err_allocate,__line__,'rsTemp')
542  END IF ! global%errorFlag
543 
544 ! ******************************************************************************
545 ! Initialize
546 ! ******************************************************************************
547 
548 ! TEMPORARY
549  nlevels = 3
550 ! END TEMPORARY
551 
552  rssize = 1
553 
554  rs(rssize) = icg
555 
556 ! ******************************************************************************
557 ! Cell-to-cell stencil for cell icg itself
558 ! ******************************************************************************
559 
560  DO isl = 1,pgrid%c2cs(icg)%nCellMembs
561  IF ( rssize < rssizemax ) THEN
562  rssize = rssize + 1
563 
564  rs(rssize) = pgrid%c2cs(icg)%cellMembs(isl)
565  ELSE
566 ! TEMPORARY
567  WRITE(*,*) 'ERROR - About to exceed size of array rs!'
568  stop
569 ! END TEMPORARY
570  END IF ! rsSize
571  END DO ! isl
572 
573 ! ==============================================================================
574 ! Save initial stencil and sort for searching
575 ! ==============================================================================
576 
577  rssavesize = rssize
578 
579  DO irs = 1,rssavesize
580  rssave(irs) = rs(irs)
581  END DO ! irs
582 
583  CALL quicksortinteger(rs(1:rssize),rssize)
584 
585 ! ******************************************************************************
586 ! Loop over levels
587 ! ******************************************************************************
588 
589  DO ilevel = 2,nlevels
590  rstempsize = 0
591 
592 ! ==============================================================================
593 ! Loop over members in saved stencil
594 ! ==============================================================================
595 
596  DO irs = 1,rssavesize
597  icg2 = rssave(irs)
598 
599 ! ------------------------------------------------------------------------------
600 ! Loop over stencil members
601 ! ------------------------------------------------------------------------------
602 
603  DO isl = 1,pgrid%c2cs(icg2)%nCellMembs
604  icg3 = pgrid%c2cs(icg2)%cellMembs(isl)
605 
606 ! ----- Add to temporary stencil if not already in original or saved stencil ---
607 
608  CALL binarysearchinteger(rs(1:rssize),rssize,icg3,iloc)
609 
610  IF ( iloc == element_not_found ) THEN
611  IF ( rstempsize > 1 ) THEN
612  CALL binarysearchinteger(rstemp(1:rstempsize),rstempsize,icg3,iloc)
613  ELSE
614  iloc = element_not_found
615  END IF ! rsTempSize
616 
617  IF ( iloc == element_not_found ) THEN
618  IF ( rstempsize < rstempsizemax ) THEN
619  rstempsize = rstempsize + 1
620 
621  rstemp(rstempsize) = icg3
622 
623  IF ( rstempsize > 1 ) THEN
624  CALL quicksortinteger(rstemp(1:rstempsize),rstempsize)
625  END IF ! rsTempSize
626  ELSE
627 ! TEMPORARY
628  WRITE(*,*) 'ERROR - About to exceed size of array rsTemp!'
629  stop
630 ! END TEMPORARY
631  END IF ! rsTempSize
632  END IF ! iLoc
633  END IF ! iLoc
634  END DO ! isl
635  END DO ! irs
636 
637 ! ==============================================================================
638 ! Loop over members in saved stencil, add to original stencil and sort
639 ! ==============================================================================
640 
641  DO irs = 1,rstempsize
642  rssave(irs) = rstemp(irs)
643 
644  IF ( rssize < rssizemax ) THEN
645  rssize = rssize + 1
646 
647  rs(rssize) = rstemp(irs)
648  ELSE
649 ! TEMPORARY
650  WRITE(*,*) 'ERROR - About to exceed size of array rs!'
651  stop
652 ! END TEMPORARY
653  END IF ! rsSize
654  END DO ! irs
655 
656  rssavesize = rstempsize
657 
658  CALL quicksortinteger(rs(1:rssize),rssize)
659  END DO ! iLevel
660 
661 ! DEBUG
662 ! WRITE(*,*) rs(1:rsSize)
663 ! END DEBUG
664 
665 ! ******************************************************************************
666 ! Deallocate temporary memory
667 ! ******************************************************************************
668 
669  DEALLOCATE(rssave,stat=errorflag)
670  global%error = errorflag
671  IF ( global%error /= err_none ) THEN
672  CALL errorstop(global,err_deallocate,__line__,'rsSave')
673  END IF ! global%errorFlag
674 
675  DEALLOCATE(rstemp,stat=errorflag)
676  global%error = errorflag
677  IF ( global%error /= err_none ) THEN
678  CALL errorstop(global,err_deallocate,__line__,'rsTemp')
679  END IF ! global%errorFlag
680 
681 ! ******************************************************************************
682 ! End
683 ! ******************************************************************************
684 
685  CALL deregisterfunction(global)
686 
687 END SUBROUTINE rflu_getresidualsupport2
688 
689 
690 
691 
692 
693 
694 
695 
696 
697 
698 
699 
700 
701 ! ******************************************************************************
702 !
703 ! Purpose: Set old values of conserved variables.
704 !
705 ! Description: None.
706 !
707 ! Input:
708 ! pRegion Pointer to region
709 !
710 ! Output: None.
711 !
712 ! Notes: None.
713 !
714 ! ******************************************************************************
715 
716 SUBROUTINE rflu_res_setcvold(pRegion)
717 
718  IMPLICIT NONE
719 
720 ! ******************************************************************************
721 ! Declarations
722 ! ******************************************************************************
723 
724 ! ==============================================================================
725 ! Arguments
726 ! ==============================================================================
727 
728  TYPE (t_region), POINTER :: pregion
729 
730 ! ==============================================================================
731 ! Locals
732 ! ==============================================================================
733 
734  INTEGER :: icg
735  TYPE(t_global), POINTER :: global
736  TYPE(t_grid), POINTER :: pgrid
737  TYPE(t_mixt), POINTER :: pmixt
738 
739 ! ******************************************************************************
740 ! Start
741 ! ******************************************************************************
742 
743  global => pregion%global
744 
745  CALL registerfunction(global,'RFLU_RES_SetCvOld',&
746  'RFLU_ModResidual.F90')
747 
748 ! ******************************************************************************
749 ! Set pointers
750 ! ******************************************************************************
751 
752  pgrid => pregion%grid
753  pmixt => pregion%mixt
754 
755 ! ******************************************************************************
756 ! Set old values
757 ! ******************************************************************************
758 
759  DO icg = 1,pgrid%nCellsTot
760  pmixt%cvOld(cv_mixt_dens,icg) = pmixt%cv(cv_mixt_dens,icg)
761  pmixt%cvOld(cv_mixt_xmom,icg) = pmixt%cv(cv_mixt_xmom,icg)
762  pmixt%cvOld(cv_mixt_ymom,icg) = pmixt%cv(cv_mixt_ymom,icg)
763  pmixt%cvOld(cv_mixt_zmom,icg) = pmixt%cv(cv_mixt_zmom,icg)
764  pmixt%cvOld(cv_mixt_ener,icg) = pmixt%cv(cv_mixt_ener,icg)
765  END DO ! icg
766 
767 ! ******************************************************************************
768 ! End
769 ! ******************************************************************************
770 
771  CALL deregisterfunction(global)
772 
773 END SUBROUTINE rflu_res_setcvold
774 
775 
776 
777 
778 
779 
780 
781 
782 
783 ! ******************************************************************************
784 !
785 ! Purpose: Set dissipation array.
786 !
787 ! Description: None.
788 !
789 ! Input:
790 ! pRegion Pointer to region
791 !
792 ! Output: None.
793 !
794 ! Notes: None.
795 !
796 ! ******************************************************************************
797 
798 SUBROUTINE rflu_res_setdiss(pRegion)
799 
800  IMPLICIT NONE
801 
802 ! ******************************************************************************
803 ! Declarations
804 ! ******************************************************************************
805 
806 ! ==============================================================================
807 ! Arguments
808 ! ==============================================================================
809 
810  TYPE (t_region), POINTER :: pregion
811 
812 ! ==============================================================================
813 ! Locals
814 ! ==============================================================================
815 
816  INTEGER :: icg
817  TYPE(t_global), POINTER :: global
818  TYPE(t_grid), POINTER :: pgrid
819  TYPE(t_mixt), POINTER :: pmixt
820 
821 ! ******************************************************************************
822 ! Start
823 ! ******************************************************************************
824 
825  global => pregion%global
826 
827  CALL registerfunction(global,'RFLU_RES_SetDiss',&
828  'RFLU_ModResidual.F90')
829 
830 ! ******************************************************************************
831 ! Set pointers
832 ! ******************************************************************************
833 
834  pgrid => pregion%grid
835  pmixt => pregion%mixt
836 
837 ! ******************************************************************************
838 ! Set old values
839 ! ******************************************************************************
840 
841  IF ( pregion%irkStep == 1 ) THEN
842  DO icg = 1,pgrid%nCellsTot
843  pmixt%diss(cv_mixt_dens,icg) = 0.0_rfreal
844  pmixt%diss(cv_mixt_xmom,icg) = 0.0_rfreal
845  pmixt%diss(cv_mixt_ymom,icg) = 0.0_rfreal
846  pmixt%diss(cv_mixt_zmom,icg) = 0.0_rfreal
847  pmixt%diss(cv_mixt_ener,icg) = 0.0_rfreal
848  END DO ! icg
849  END IF ! pRegion%irkStep
850 
851 ! ******************************************************************************
852 ! End
853 ! ******************************************************************************
854 
855  CALL deregisterfunction(global)
856 
857 END SUBROUTINE rflu_res_setdiss
858 
859 
860 
861 
862 
863 
864 
865 
866 
867 ! ******************************************************************************
868 !
869 ! Purpose: Set right-hand side.
870 !
871 ! Description: None.
872 !
873 ! Input:
874 ! pRegion Pointer to region
875 !
876 ! Output: None.
877 !
878 ! Notes: None.
879 !
880 ! ******************************************************************************
881 
882 SUBROUTINE rflu_res_setrhs(pRegion)
883 
884  IMPLICIT NONE
885 
886 ! ******************************************************************************
887 ! Declarations
888 ! ******************************************************************************
889 
890 ! ==============================================================================
891 ! Arguments
892 ! ==============================================================================
893 
894  TYPE (t_region), POINTER :: pregion
895 
896 ! ==============================================================================
897 ! Locals
898 ! ==============================================================================
899 
900  INTEGER :: icg
901  TYPE(t_global), POINTER :: global
902  TYPE(t_grid), POINTER :: pgrid
903  TYPE(t_mixt), POINTER :: pmixt
904 
905 ! ******************************************************************************
906 ! Start
907 ! ******************************************************************************
908 
909  global => pregion%global
910 
911  CALL registerfunction(global,'RFLU_RES_SetRhs',&
912  'RFLU_ModResidual.F90')
913 
914 ! ******************************************************************************
915 ! Set pointers
916 ! ******************************************************************************
917 
918  pgrid => pregion%grid
919  pmixt => pregion%mixt
920 
921 ! ******************************************************************************
922 ! Set right-hand side
923 ! ******************************************************************************
924 
925  DO icg = 1,pgrid%nCellsTot
926  pmixt%rhs(cv_mixt_dens,icg) = -pmixt%diss(cv_mixt_dens,icg)
927  pmixt%rhs(cv_mixt_xmom,icg) = -pmixt%diss(cv_mixt_xmom,icg)
928  pmixt%rhs(cv_mixt_ymom,icg) = -pmixt%diss(cv_mixt_ymom,icg)
929  pmixt%rhs(cv_mixt_zmom,icg) = -pmixt%diss(cv_mixt_zmom,icg)
930  pmixt%rhs(cv_mixt_ener,icg) = -pmixt%diss(cv_mixt_ener,icg)
931  END DO ! icg
932 
933 ! ******************************************************************************
934 ! End
935 ! ******************************************************************************
936 
937  CALL deregisterfunction(global)
938 
939 END SUBROUTINE rflu_res_setrhs
940 
941 
942 
943 
944 
945 
946 
947 ! ******************************************************************************
948 ! End
949 ! ******************************************************************************
950 
951 END MODULE rflu_modresidual
952 
953 
954 ! ******************************************************************************
955 !
956 ! RCS Revision history:
957 !
958 ! $Log: RFLU_ModResidual.F90,v $
959 ! Revision 1.18 2008/12/06 08:44:24 mtcampbe
960 ! Updated license.
961 !
962 ! Revision 1.17 2008/11/19 22:17:35 mtcampbe
963 ! Added Illinois Open Source License/Copyright
964 !
965 ! Revision 1.16 2006/08/19 15:39:16 mparmar
966 ! Moved region%mixt%bGradFace to patch%mixt%gradFace
967 !
968 ! Revision 1.15 2006/05/01 21:02:07 haselbac
969 ! Converted to new flux routine
970 !
971 ! Revision 1.14 2006/04/27 15:10:32 haselbac
972 ! Added VENKAT limiter option
973 !
974 ! Revision 1.13 2006/04/15 17:00:38 haselbac
975 ! Added RECONST_NONE and cReconst IF for bf grads
976 !
977 ! Revision 1.12 2006/04/07 15:19:20 haselbac
978 ! Removed tabs
979 !
980 ! Revision 1.11 2006/04/07 14:49:42 haselbac
981 ! Adapted to changes in f and bf diff modules and new WENO module
982 !
983 ! Revision 1.10 2006/01/06 22:13:11 haselbac
984 ! Added call to cell grad wrapper routine
985 !
986 ! Revision 1.9 2005/11/14 16:59:29 haselbac
987 ! Added support for pseudo-gas model
988 !
989 ! Revision 1.8 2005/11/10 02:28:34 haselbac
990 ! Adapted to changes in AUSM flux module
991 !
992 ! Revision 1.7 2005/10/27 19:18:46 haselbac
993 ! Adapted to changes in gradient routines
994 !
995 ! Revision 1.6 2005/10/14 14:08:09 haselbac
996 ! Added call to RFLU_EnforceHeatFlux
997 !
998 ! Revision 1.5 2005/10/05 14:05:51 haselbac
999 ! Adapted to changes in gradient computation
1000 !
1001 ! Revision 1.4 2005/08/17 20:25:22 hdewey2
1002 ! Moved RFLU_GetResidualSupport2 from RFLU_ModPETScNewtonKrylov, added RFLU_GetResidualSupport1
1003 !
1004 ! Revision 1.3 2005/07/14 21:44:28 haselbac
1005 ! Added limiters, new ENO scheme, and AUSM flux function
1006 !
1007 ! Revision 1.2 2005/07/11 19:32:40 mparmar
1008 ! Adapted call to RFLU_LimitGradCells to changes in modules
1009 !
1010 ! Revision 1.1 2005/05/19 18:17:52 haselbac
1011 ! Initial revision
1012 !
1013 ! ******************************************************************************
1014 
1015 
1016 
1017 
1018 
1019 
1020 
1021 
1022 
1023 
1024 
1025 
1026 
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)
subroutine rs(nm, n, a, w, matz, z, fv1, fv2, ierr)
subroutine, public rflu_res_setcvold(pRegion)
subroutine, public rflu_res_computeresidual(pRegion)
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_getresidualsupport1(pRegion, icg, rs, rsSizeMax, rsSize)
subroutine, public rflu_res_setrhs(pRegion)
subroutine quicksortinteger(a, n)
subroutine, public rflu_computegradfacesconstr(pRegion, iBegVar, iEndVar, iBegGrad, iEndGrad, varInfo, var, grad)
subroutine, public rflu_limitgradcells(pRegion, iBegGrad, iEndGrad, grad, lim)
subroutine binarysearchinteger(a, n, v, i, j)
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_convertcvcons2prim(pRegion, cvStateFuture)
subroutine, public rflu_res_setdiss(pRegion)
IndexType nfaces() const
Definition: Mesh.H:641
subroutine, public rflu_convertcvprim2cons(pRegion, cvStateFuture)
subroutine, public rflu_wenogradcellswrapper(pRegion, iBegGrad, iEndGrad, grad)
subroutine, public rflu_getresidualsupport2(pRegion, icg, rs, rsSizeMax, rsSize)
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)