Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
RFLU_ModAllocateMemory.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: Suite of routines to allocate memory.
26 !
27 ! Description: None.
28 !
29 ! Notes: None.
30 !
31 ! ******************************************************************************
32 !
33 ! $Id: RFLU_ModAllocateMemory.F90,v 1.16 2008/12/06 08:44:20 mtcampbe Exp $
34 !
35 ! Copyright: (c) 2004-2005 by the University of Illinois
36 !
37 ! ******************************************************************************
38 
40 
41  USE moddatatypes
42  USE moderror
43  USE modglobal, ONLY: t_global
44  USE modparameters
45  USE modgrid, ONLY: t_grid
46  USE modbndpatch, ONLY: t_patch
47  USE moddatastruct, ONLY: t_region
48  USE modmixture, ONLY: t_mixt_input
49  USE modmpi
50 
52 
53  IMPLICIT NONE
54 
55  PRIVATE
56  PUBLIC :: rflu_allocatememorygspeeds, &
65 
66 ! ******************************************************************************
67 ! Declarations and definitions
68 ! ******************************************************************************
69 
70  CHARACTER(CHRLEN) :: RCSIdentString = &
71  '$RCSfile: RFLU_ModAllocateMemory.F90,v $ $Revision: 1.16 $'
72 
73 ! ******************************************************************************
74 ! Routines
75 ! ******************************************************************************
76 
77  CONTAINS
78 
79 
80 
81 
82 
83 ! ******************************************************************************
84 !
85 ! Purpose: Allocate memory for grid speeds.
86 !
87 ! Description: None.
88 !
89 ! Input:
90 ! pRegion Region pointer
91 !
92 ! Output: None.
93 !
94 ! Notes: None.
95 !
96 ! ******************************************************************************
97 
98 SUBROUTINE rflu_allocatememorygspeeds(pRegion)
99 
101 
102  IMPLICIT NONE
103 
104 ! ******************************************************************************
105 ! Definitions and declarations
106 ! ******************************************************************************
107 
108 ! ==============================================================================
109 ! Arguments
110 ! ==============================================================================
111 
112  TYPE(t_region), POINTER :: pregion
113 
114 ! ==============================================================================
115 ! Locals
116 ! ==============================================================================
117 
118  INTEGER :: errorflag,ifc,ipatch
119  TYPE(t_grid), POINTER :: pgrid
120  TYPE(t_mixt_input), POINTER :: pmixtinput
121  TYPE(t_patch), POINTER :: ppatch
122  TYPE(t_global), POINTER :: global
123 
124 ! ******************************************************************************
125 ! Start
126 ! ******************************************************************************
127 
128  global => pregion%global
129 
130  CALL registerfunction(global,'RFLU_AllocateMemoryGSpeeds',&
131  'RFLU_ModAllocateMemory.F90')
132 
133 ! ******************************************************************************
134 ! Set pointers
135 ! ******************************************************************************
136 
137  pgrid => pregion%grid
138  pmixtinput => pregion%mixtInput
139 
140 ! ******************************************************************************
141 ! Allocate memory
142 ! ******************************************************************************
143 
144 ! ==============================================================================
145 ! Require grid speeds
146 ! ==============================================================================
147 
148  IF ( rflu_decideneedgridspeeds(pregion) .EQV. .true. ) THEN
149 
150 ! ------------------------------------------------------------------------------
151 ! Interior faces
152 ! ------------------------------------------------------------------------------
153 
154 #ifndef GENX
155  ALLOCATE(pgrid%gs(pgrid%nFaces),stat=errorflag)
156  global%error = errorflag
157  IF ( global%error /= err_none ) THEN
158  CALL errorstop(global,err_allocate,__line__,'pGrid%gs')
159  END IF ! global%error
160 
161  DO ifc = 1,pgrid%nFaces ! Explicit loop to avoid Frost problem
162  pgrid%gs(ifc) = 0.0_rfreal
163  END DO ! ifc
164 #else
165  ALLOCATE(pgrid%gs(pgrid%nFacesTot),stat=errorflag)
166  global%error = errorflag
167  IF ( global%error /= err_none ) THEN
168  CALL errorstop(global,err_allocate,__line__,'pGrid%gs')
169  END IF ! global%error
170 
171  DO ifc = 1,pgrid%nFacesTot ! Explicit loop to avoid Frost problem
172  pgrid%gs(ifc) = 0.0_rfreal
173  END DO ! ifc
174 #endif
175 
176 ! ------------------------------------------------------------------------------
177 ! Patch faces
178 ! ------------------------------------------------------------------------------
179 
180  IF ( pgrid%nPatches > 0 ) THEN
181  DO ipatch = 1,pgrid%nPatches
182  ppatch => pregion%patches(ipatch)
183 
184 #ifndef GENX
185  ALLOCATE(ppatch%gs(ppatch%nBFaces),stat=errorflag)
186  global%error = errorflag
187  IF ( global%error /= err_none ) THEN
188  CALL errorstop(global,err_allocate,__line__,'pPatch%gs')
189  END IF ! global%error
190 
191  DO ifc = 1,ppatch%nBFaces ! Explicit loop to avoid Frost problem
192  ppatch%gs(ifc) = 0.0_rfreal
193  END DO ! ifc
194 #else
195  ALLOCATE(ppatch%gs(ppatch%nBFacesTot),stat=errorflag)
196  global%error = errorflag
197  IF ( global%error /= err_none ) THEN
198  CALL errorstop(global,err_allocate,__line__,'pPatch%gs')
199  END IF ! global%error
200 
201  DO ifc = 1,ppatch%nBFaces ! Explicit loop to avoid Frost problem
202  ppatch%gs(ifc) = 0.0_rfreal
203  END DO ! ifc
204 
205  DO ifc = ppatch%nBFaces+1,ppatch%nBFacesTot
206  ppatch%gs(ifc) = REAL(crazy_value_int,kind=rfreal)
207  END DO ! ifc
208 #endif
209  END DO ! iPatch
210  END IF ! pGrid%nPatches
211 
212 ! ==============================================================================
213 ! Do not require grid speeds
214 ! ==============================================================================
215 
216  ELSE
217 
218 ! ------------------------------------------------------------------------------
219 ! Interior faces
220 ! ------------------------------------------------------------------------------
221 
222  ALLOCATE(pgrid%gs(0:1),stat=errorflag)
223  global%error = errorflag
224  IF ( global%error /= err_none ) THEN
225  CALL errorstop(global,err_allocate,__line__,'pGrid%gs')
226  END IF ! global%error
227 
228  DO ifc = 0,1 ! Explicit loop to avoid Frost problem
229  pgrid%gs(ifc) = 0.0_rfreal
230  END DO ! ifc
231 
232 ! ------------------------------------------------------------------------------
233 ! Patch faces
234 ! ------------------------------------------------------------------------------
235 
236  IF ( pgrid%nPatches > 0 ) THEN
237  DO ipatch = 1,pgrid%nPatches
238  ppatch => pregion%patches(ipatch)
239 
240  ALLOCATE(ppatch%gs(0:1),stat=errorflag)
241  global%error = errorflag
242  IF ( global%error /= err_none ) THEN
243  CALL errorstop(global,err_allocate,__line__,'pPatch%gs')
244  END IF ! global%error
245 
246  DO ifc = 0,1 ! Explicit loop to avoid Frost problem
247  ppatch%gs(ifc) = 0.0_rfreal
248  END DO ! ifc
249  END DO ! iPatch
250  END IF ! pGrid%nPatches
251  END IF ! pMixtInput%moveGrid
252 
253 ! ******************************************************************************
254 ! End
255 ! ******************************************************************************
256 
257  CALL deregisterfunction(global)
258 
259 END SUBROUTINE rflu_allocatememorygspeeds
260 
261 
262 
263 
264 
265 
266 ! ******************************************************************************
267 !
268 ! Purpose: Allocate memory for mixture solution.
269 !
270 ! Description: None.
271 !
272 ! Input:
273 ! pRegion Region pointer
274 !
275 ! Output: None.
276 !
277 ! Notes: None.
278 !
279 ! ******************************************************************************
280 
281 SUBROUTINE rflu_allocatememorysol(pRegion)
282 
283  IMPLICIT NONE
284 
285 ! ******************************************************************************
286 ! Definitions and declarations
287 ! ******************************************************************************
288 
289 ! ==============================================================================
290 ! Arguments
291 ! ==============================================================================
292 
293  TYPE(t_region), POINTER :: pregion
294 
295 ! ==============================================================================
296 ! Locals
297 ! ==============================================================================
298 
299  TYPE(t_global), POINTER :: global
300 
301 ! ******************************************************************************
302 ! Start
303 ! ******************************************************************************
304 
305  global => pregion%global
306 
307  CALL registerfunction(global,'RFLU_AllocateMemorySol',&
308  'RFLU_ModAllocateMemory.F90')
309 
310 ! ******************************************************************************
311 ! Allocate memory
312 ! ******************************************************************************
313 
314  CALL rflu_allocatememorysolcv(pregion)
315  CALL rflu_allocatememorysoldv(pregion)
316  CALL rflu_allocatememorysolgv(pregion)
317  CALL rflu_allocatememorysoltv(pregion)
318 
319 ! ******************************************************************************
320 ! End
321 ! ******************************************************************************
322 
323  CALL deregisterfunction(global)
324 
325 END SUBROUTINE rflu_allocatememorysol
326 
327 
328 
329 
330 
331 
332 ! ******************************************************************************
333 !
334 ! Purpose: Allocate memory for conserved variables of mixture.
335 !
336 ! Description: None.
337 !
338 ! Input:
339 ! pRegion Region pointer
340 !
341 ! Output: None.
342 !
343 ! Notes: None.
344 !
345 ! ******************************************************************************
346 
347 SUBROUTINE rflu_allocatememorysolcv(pRegion)
348 
349  IMPLICIT NONE
350 
351 ! ******************************************************************************
352 ! Definitions and declarations
353 ! ******************************************************************************
354 
355 ! ==============================================================================
356 ! Arguments
357 ! ==============================================================================
358 
359  TYPE(t_region), POINTER :: pregion
360 
361 ! ==============================================================================
362 ! Locals
363 ! ==============================================================================
364 
365  INTEGER :: errorflag
366  TYPE(t_global), POINTER :: global
367  TYPE(t_grid), POINTER :: pgrid
368  TYPE(t_mixt_input), POINTER :: pmixtinput
369 
370 ! ******************************************************************************
371 ! Start
372 ! ******************************************************************************
373 
374  global => pregion%global
375 
376  CALL registerfunction(global,'RFLU_AllocateMemorySolCv',&
377  'RFLU_ModAllocateMemory.F90')
378 
379 ! ******************************************************************************
380 ! Set pointers
381 ! ******************************************************************************
382 
383  pgrid => pregion%grid
384  pmixtinput => pregion%mixtInput
385 
386 ! ******************************************************************************
387 ! Allocate memory
388 ! ******************************************************************************
389 
390  ALLOCATE(pregion%mixt%cv(pmixtinput%nCv,pgrid%nCellsTot),stat=errorflag)
391  global%error = errorflag
392  IF (global%error /= err_none) THEN
393  CALL errorstop(global,err_allocate,__line__,'pRegion%mixt%cv')
394  END IF ! global%error
395 
396  ALLOCATE(pregion%mixt%cvInfo(pmixtinput%nCv),stat=errorflag)
397  global%error = errorflag
398  IF (global%error /= err_none) THEN
399  CALL errorstop(global,err_allocate,__line__,'pRegion%mixt%cvInfo')
400  END IF ! global%error
401 
402 ! ******************************************************************************
403 ! End
404 ! ******************************************************************************
405 
406  CALL deregisterfunction(global)
407 
408 END SUBROUTINE rflu_allocatememorysolcv
409 
410 
411 
412 
413 
414 ! ******************************************************************************
415 !
416 ! Purpose: Allocate memory for dependent variables of mixture.
417 !
418 ! Description: None.
419 !
420 ! Input:
421 ! pRegion Region pointer
422 !
423 ! Output: None.
424 !
425 ! Notes: None.
426 !
427 ! ******************************************************************************
428 
429 SUBROUTINE rflu_allocatememorysoldv(pRegion)
430 
431  IMPLICIT NONE
432 
433 ! ******************************************************************************
434 ! Definitions and declarations
435 ! ******************************************************************************
436 
437 ! ==============================================================================
438 ! Arguments
439 ! ==============================================================================
440 
441  TYPE(t_region), POINTER :: pregion
442 
443 ! ==============================================================================
444 ! Locals
445 ! ==============================================================================
446 
447  INTEGER :: errorflag
448  TYPE(t_global), POINTER :: global
449  TYPE(t_grid), POINTER :: pgrid
450  TYPE(t_mixt_input), POINTER :: pmixtinput
451 
452 ! ******************************************************************************
453 ! Start
454 ! ******************************************************************************
455 
456  global => pregion%global
457 
458  CALL registerfunction(global,'RFLU_AllocateMemorySolDv',&
459  'RFLU_ModAllocateMemory.F90')
460 
461 ! ******************************************************************************
462 ! Set pointers
463 ! ******************************************************************************
464 
465  pgrid => pregion%grid
466  pmixtinput => pregion%mixtInput
467 
468 ! ******************************************************************************
469 ! Allocate memory
470 ! ******************************************************************************
471 
472  ALLOCATE(pregion%mixt%dv(pmixtinput%nDv,pgrid%nCellsTot),stat=errorflag)
473  global%error = errorflag
474  IF (global%error /= err_none) THEN
475  CALL errorstop(global,err_allocate,__line__,'pRegion%mixt%dv')
476  END IF ! global%error
477 
478 ! ******************************************************************************
479 ! End
480 ! ******************************************************************************
481 
482  CALL deregisterfunction(global)
483 
484 END SUBROUTINE rflu_allocatememorysoldv
485 
486 
487 
488 
489 
490 ! ******************************************************************************
491 !
492 ! Purpose: Allocate memory for gas variables of mixture.
493 !
494 ! Description: None.
495 !
496 ! Input:
497 ! pRegion Region pointer
498 !
499 ! Output: None.
500 !
501 ! Notes: None.
502 !
503 ! ******************************************************************************
504 
505 SUBROUTINE rflu_allocatememorysolgv(pRegion)
506 
507  IMPLICIT NONE
508 
509 ! ******************************************************************************
510 ! Definitions and declarations
511 ! ******************************************************************************
512 
513 ! ==============================================================================
514 ! Arguments
515 ! ==============================================================================
516 
517  TYPE(t_region), POINTER :: pregion
518 
519 ! ==============================================================================
520 ! Locals
521 ! ==============================================================================
522 
523  INTEGER :: errorflag
524  TYPE(t_global), POINTER :: global
525  TYPE(t_grid), POINTER :: pgrid
526  TYPE(t_mixt_input), POINTER :: pmixtinput
527 
528 ! ******************************************************************************
529 ! Start
530 ! ******************************************************************************
531 
532  global => pregion%global
533 
534  CALL registerfunction(global,'RFLU_AllocateMemorySolGv',&
535  'RFLU_ModAllocateMemory.F90')
536 
537 ! ******************************************************************************
538 ! Set pointers
539 ! ******************************************************************************
540 
541  pgrid => pregion%grid
542  pmixtinput => pregion%mixtInput
543 
544 ! ******************************************************************************
545 ! Allocate memory
546 ! ******************************************************************************
547 
548  IF ( pmixtinput%gasModel == gas_model_tcperf ) THEN
549  ALLOCATE(pregion%mixt%gv(pmixtinput%nGv,0:1),stat=errorflag)
550  ELSE
551  ALLOCATE(pregion%mixt%gv(pmixtinput%nGv,pgrid%nCellsTot),stat=errorflag)
552  END IF ! pRegion
553  global%error = errorflag
554  IF (global%error /= err_none) THEN
555  CALL errorstop(global,err_allocate,__line__,'pRegion%mixt%gv')
556  END IF ! global%error
557 
558 ! ******************************************************************************
559 ! End
560 ! ******************************************************************************
561 
562  CALL deregisterfunction(global)
563 
564 END SUBROUTINE rflu_allocatememorysolgv
565 
566 
567 
568 
569 
570 ! ******************************************************************************
571 !
572 ! Purpose: Allocate memory for transport variables of mixture.
573 !
574 ! Description: None.
575 !
576 ! Input:
577 ! pRegion Region pointer
578 !
579 ! Output: None.
580 !
581 ! Notes: None.
582 !
583 ! ******************************************************************************
584 
585 SUBROUTINE rflu_allocatememorysoltv(pRegion)
586 
587  IMPLICIT NONE
588 
589 ! ******************************************************************************
590 ! Definitions and declarations
591 ! ******************************************************************************
592 
593 ! ==============================================================================
594 ! Arguments
595 ! ==============================================================================
596 
597  TYPE(t_region), POINTER :: pregion
598 
599 ! ==============================================================================
600 ! Locals
601 ! ==============================================================================
602 
603  INTEGER :: errorflag
604  TYPE(t_global), POINTER :: global
605  TYPE(t_grid), POINTER :: pgrid
606  TYPE(t_mixt_input), POINTER :: pmixtinput
607 
608 ! ******************************************************************************
609 ! Start
610 ! ******************************************************************************
611 
612  global => pregion%global
613 
614  CALL registerfunction(global,'RFLU_AllocateMemorySolTv',&
615  'RFLU_ModAllocateMemory.F90')
616 
617 ! ******************************************************************************
618 ! Set pointers
619 ! ******************************************************************************
620 
621  pgrid => pregion%grid
622  pmixtinput => pregion%mixtInput
623 
624 ! ******************************************************************************
625 ! Allocate memory
626 ! ******************************************************************************
627 
628  IF ( pmixtinput%computeTv .EQV. .true. ) THEN
629  ALLOCATE(pregion%mixt%tv(pmixtinput%nTv,pgrid%nCellsTot),stat=errorflag)
630  global%error = errorflag
631  IF (global%error /= err_none) THEN
632  CALL errorstop(global,err_allocate,__line__,'pRegion%mixt%tv')
633  END IF ! global%error
634  ELSE
635  nullify(pregion%mixt%tv)
636  END IF ! pMixtInput%computeTv
637 
638 ! ******************************************************************************
639 ! End
640 ! ******************************************************************************
641 
642  CALL deregisterfunction(global)
643 
644 END SUBROUTINE rflu_allocatememorysoltv
645 
646 
647 
648 
649 
650 ! ******************************************************************************
651 !
652 ! Purpose: Allocate memory for mixture time-stepping.
653 !
654 ! Description: None.
655 !
656 ! Input:
657 ! pRegion Region pointer
658 !
659 ! Output: None.
660 !
661 ! Notes: None.
662 !
663 ! ******************************************************************************
664 
665 SUBROUTINE rflu_allocatememorytstep(pRegion)
666 
667  USE rflu_modoles
668 
669  IMPLICIT NONE
670 
671 ! ******************************************************************************
672 ! Definitions and declarations
673 ! ******************************************************************************
674 
675 ! ==============================================================================
676 ! Arguments
677 ! ==============================================================================
678 
679  TYPE(t_region), POINTER :: pregion
680 
681 ! ==============================================================================
682 ! Locals
683 ! ==============================================================================
684 
685  INTEGER :: arraylimlow,arraylimupp,errorflag,ic,icmp,ifc,iv,ipatch,nbfaces, &
686  nbfacestot
687  TYPE(t_grid), POINTER :: pgrid,pgridold,pgridold2
688  TYPE(t_patch), POINTER :: ppatch
689  TYPE(t_global), POINTER :: global
690  TYPE(t_mixt_input), POINTER :: pmixtinput
691 
692 ! ******************************************************************************
693 ! Start
694 ! ******************************************************************************
695 
696  global => pregion%global
697 
698  CALL registerfunction(global,'RFLU_AllocateMemoryTStep',&
699  'RFLU_ModAllocateMemory.F90')
700 
701 ! ******************************************************************************
702 ! Set pointers
703 ! ******************************************************************************
704 
705  pgrid => pregion%grid
706  pgridold => pregion%gridOld
707  pgridold2 => pregion%gridOld2
708  pmixtinput => pregion%mixtInput
709 
710 ! ******************************************************************************
711 ! Allocate memory
712 ! ******************************************************************************
713 
714 ! ==============================================================================
715 ! Compute number of boundary faces for later use
716 ! ==============================================================================
717 
718  nbfaces = 0
719  nbfacestot = 0
720 
721  DO ipatch = 1,pgrid%nPatches
722  ppatch => pregion%patches(ipatch)
723 
724  nbfaces = nbfaces + ppatch%nBTris + ppatch%nBQuads
725  nbfacestot = nbfacestot + ppatch%nBTrisTot + ppatch%nBQuadsTot
726  END DO ! iPatch
727 
728 ! ==============================================================================
729 ! Old solutions
730 ! ==============================================================================
731 
732  ALLOCATE(pregion%mixt%cvOld(pmixtinput%nCv,pgrid%nCellsTot),stat=errorflag)
733  global%error = errorflag
734  IF (global%error /= err_none) THEN
735  CALL errorstop(global,err_allocate,__line__,'pRegion%mixt%cvOld')
736  END IF ! global%error
737 
738  IF ( global%solverType == solv_implicit_nk ) THEN
739  ALLOCATE(pregion%mixt%cvOld1(pmixtinput%nCv,pgrid%nCellsTot),stat=errorflag)
740  global%error = errorflag
741  IF (global%error /= err_none) THEN
742  CALL errorstop(global,err_allocate,__line__,'pRegion%mixt%cvOld1')
743  END IF ! global%error
744 
745  ALLOCATE(pregion%mixt%cvOld2(pmixtinput%nCv,pgrid%nCellsTot),stat=errorflag)
746  global%error = errorflag
747  IF (global%error /= err_none) THEN
748  CALL errorstop(global,err_allocate,__line__,'pRegion%mixt%cvOld2')
749  END IF ! global%error
750  END IF ! global%solverType
751 
752 ! ==============================================================================
753 ! Time step
754 ! ==============================================================================
755 
756  ALLOCATE(pregion%dt(pgrid%nCellsTot),stat=errorflag)
757  global%error = errorflag
758  IF ( global%error /= err_none ) THEN
759  CALL errorstop(global,err_allocate,__line__,'pRegion%dt')
760  END IF ! global%error
761 
762 ! ==============================================================================
763 ! Residuals
764 ! ==============================================================================
765 
766  ALLOCATE(pregion%mixt%rhs(pmixtinput%nCv,pgrid%nCellsTot),stat=errorflag)
767  global%error = errorflag
768  IF ( global%error /= err_none ) THEN
769  CALL errorstop(global,err_allocate,__line__,'pRegion%mixt%rhs')
770  END IF ! global%error
771 
772  ALLOCATE(pregion%mixt%diss(pmixtinput%nCv,pgrid%nCellsTot),stat=errorflag)
773  global%error = errorflag
774  IF ( global%error /= err_none ) THEN
775  CALL errorstop(global,err_allocate,__line__,'pRegion%mixt%diss')
776  END IF ! global%error
777 
778  IF ( global%flowType == flow_unsteady ) THEN
779  ALLOCATE(pregion%mixt%rhsSum(pmixtinput%nCv,pgrid%nCellsTot),stat=errorflag)
780  global%error = errorflag
781  IF ( global%error /= err_none ) THEN
782  CALL errorstop(global,err_allocate,__line__,'pRegion%mixt%rhsSum')
783  END IF ! global%error
784  ELSE
785  nullify(pregion%mixt%rhsSum)
786  END IF ! global%flowType
787 
788 ! ==============================================================================
789 ! Gradients
790 ! ==============================================================================
791 
792 ! ------------------------------------------------------------------------------
793 ! Cell gradients
794 ! ------------------------------------------------------------------------------
795 
796  IF ( (pmixtinput%spaceDiscr == discr_upw_roe ) .OR. &
797  (pmixtinput%spaceDiscr == discr_upw_hllc ) .OR. &
798  (pmixtinput%spaceDiscr == discr_upw_ausmplus) ) THEN
799  IF ( pmixtinput%spaceOrder > 1 ) THEN
800  ALLOCATE(pregion%mixt%gradCell(xcoord:zcoord, &
801  grc_mixt_dens:grc_mixt_pres, &
802  pgrid%nCellsTot),stat=errorflag)
803  global%error = errorflag
804  IF ( global%error /= err_none ) THEN
805  CALL errorstop(global,err_allocate,__line__,'pRegion%mixt%gradCell')
806  END IF ! global%error
807  ELSE
808  nullify(pregion%mixt%gradCell)
809  END IF ! pMixtInput%spaceOrder
810  ELSE IF ( pmixtinput%spaceDiscr == discr_opt_les ) THEN
811  ALLOCATE(pregion%mixt%gradCell(xcoord:zcoord, &
812  grc_mixt_xvel:grc_mixt_zvel, &
813  pgrid%nCellsTot),stat=errorflag)
814  global%error = errorflag
815  IF ( global%error /= err_none ) THEN
816  CALL errorstop(global,err_allocate,__line__,'pRegion%mixt%gradCell')
817  END IF ! global%error
818  ELSE
819  nullify(pregion%mixt%gradCell)
820  END IF ! pMixtInput%spaceDiscr
821 
822 ! ------------------------------------------------------------------------------
823 ! Face gradients
824 ! ------------------------------------------------------------------------------
825 
826  IF ( pmixtinput%flowModel == flow_navst ) THEN
827  ALLOCATE(pregion%mixt%gradFace(xcoord:zcoord, &
828  grf_mixt_xvel:grf_mixt_temp, &
829  pgrid%nFaces),stat=errorflag)
830  global%error = errorflag
831  IF ( global%error /= err_none ) THEN
832  CALL errorstop(global,err_allocate,__line__,'pRegion%mixt%gradFace')
833  END IF ! global%error
834  ELSE
835  nullify(pregion%mixt%gradFace)
836  END IF ! pMixtInput%flowModel
837 
838  IF ( pgrid%nPatches > 0 ) THEN
839  DO ipatch = 1,pregion%grid%nPatches
840  ppatch => pregion%patches(ipatch)
841 
842  IF ( rflu_decideneedbgradface(pregion,ppatch) .EQV. .true. ) THEN
843  ALLOCATE(ppatch%mixt%gradFace(xcoord:zcoord, &
844  grbf_mixt_dens:grbf_mixt_pres, &
845  ppatch%nBFaces),stat=errorflag)
846  global%error = errorflag
847  IF ( global%error /= err_none ) THEN
848  CALL errorstop(global,err_allocate,__line__,'pPatch%mixt%gradFace')
849  END IF ! global%error
850  ELSE
851  nullify(ppatch%mixt%gradFace)
852  END IF ! RFLU_DecideNeedBGradFace
853  END DO ! iPatch
854  END IF ! pGrid%nPatches
855 
856 ! ==============================================================================
857 ! Grid motion. NOTE grid speeds are allocated separately because they are
858 ! written into grid file, and hence they need to be allocated in pre- and
859 ! postprocessors also.
860 ! ==============================================================================
861 
862  IF ( pmixtinput%moveGrid .EQV. .true. ) THEN
863 
864 ! ------------------------------------------------------------------------------
865 ! Residual
866 ! ------------------------------------------------------------------------------
867 
868  ALLOCATE(pgrid%rhs(xcoord:zcoord,pgrid%nVertTot),stat=errorflag)
869  global%error = errorflag
870  IF ( global%error /= err_none ) THEN
871  CALL errorstop(global,err_allocate,__line__,'pGrid%rhs')
872  END IF ! global%error
873 
874  DO icmp = xcoord,zcoord ! Explicit loop to avoid Frost problem
875  DO iv = 1,pgrid%nVertTot
876  pgrid%rhs(icmp,iv) = 0.0_rfreal
877  END DO ! iv
878  END DO ! icmp
879 
880 ! ------------------------------------------------------------------------------
881 ! Displacement
882 ! ------------------------------------------------------------------------------
883 
884  IF ( pmixtinput%moveGridType /= movegrid_type_xyz ) THEN
885  ALLOCATE(pgrid%disp(xcoord:zcoord,pgrid%nVertTot),stat=errorflag)
886  global%error = errorflag
887  IF ( global%error /= err_none ) THEN
888  CALL errorstop(global,err_allocate,__line__,'pGrid%disp')
889  END IF ! global%error
890 
891  DO icmp = xcoord,zcoord ! Explicit loop to avoid Frost problem
892  DO iv = 1,pgrid%nVertTot
893  pgrid%disp(icmp,iv) = 0.0_rfreal
894  END DO ! iv
895  END DO ! icmp
896  END IF ! pMixtInput%moveGridType
897 
898 ! ------------------------------------------------------------------------------
899 ! Old coordinates
900 ! ------------------------------------------------------------------------------
901 
902  ALLOCATE(pgridold%xyz(xcoord:zcoord,pgrid%nVertTot),stat=errorflag)
903  global%error = errorflag
904  IF ( global%error /= err_none ) THEN
905  CALL errorstop(global,err_allocate,__line__,'pRegion%gridOld%xyz')
906  END IF ! global%error
907 
908  DO icmp = xcoord,zcoord ! Explicit loop to avoid Frost problem
909  DO iv = 1,pgrid%nVertTot
910  pgridold%xyz(icmp,iv) = 0.0_rfreal
911  END DO ! iv
912  END DO ! icmp
913 
914  IF ( global%solverType == solv_implicit_nk ) THEN
915  ALLOCATE(pgridold2%xyz(xcoord:zcoord,pgrid%nVertTot),stat=errorflag)
916  global%error = errorflag
917  IF ( global%error /= err_none ) THEN
918  CALL errorstop(global,err_allocate,__line__,'pRegion%gridOld2%xyz')
919  END IF ! global%error
920 
921  DO icmp = xcoord,zcoord ! Explicit loop to avoid Frost problem
922  DO iv = 1,pgrid%nVertTot
923  pgridold2%xyz(icmp,iv) = 0.0_rfreal
924  END DO ! iv
925  END DO ! icmp
926  END IF ! global%solverType
927 
928 ! ------------------------------------------------------------------------------
929 ! Old volume
930 ! ------------------------------------------------------------------------------
931 
932  ALLOCATE(pgridold%vol(pgrid%nCellsTot),stat=errorflag)
933  global%error = errorflag
934  IF ( global%error /= err_none ) THEN
935  CALL errorstop(global,err_allocate,__line__,'pRegion%gridOld%vol')
936  END IF ! global%error
937 
938  DO ic = 1,pgrid%nCellsTot ! Explicit loop to avoid Frost problem
939  pgridold%vol(ic) = 0.0_rfreal
940  END DO ! ic
941 
942  IF ( global%solverType == solv_implicit_nk ) THEN
943  ALLOCATE(pgridold2%vol(pgrid%nCellsTot),stat=errorflag)
944  global%error = errorflag
945  IF ( global%error /= err_none ) THEN
946  CALL errorstop(global,err_allocate,__line__,'pRegion%gridOld2%vol')
947  END IF ! global%error
948 
949  DO ic = 1,pgrid%nCellsTot ! Explicit loop to avoid Frost problem
950  pgridold2%vol(ic) = 0.0_rfreal
951  END DO ! ic
952  END IF ! global%solverType
953 
954 #ifndef GENX
955 ! ------------------------------------------------------------------------------
956 ! Patch displacements. NOTE allocate here only if not running inside GENX,
957 ! because when running inside GENX allocate displacements also for virtual
958 ! vertices.
959 ! ------------------------------------------------------------------------------
960 
961  IF ( pgrid%nPatches > 0 ) THEN
962  DO ipatch = 1,pgrid%nPatches
963  ppatch => pregion%patches(ipatch)
964 
965  ALLOCATE(ppatch%dXyz(xcoord:zcoord,ppatch%nBVert),stat=errorflag)
966  global%error = errorflag
967  IF ( global%error /= err_none ) THEN
968  CALL errorstop(global,err_allocate,__line__,'pPatch%dXyz')
969  END IF ! global%error
970 
971  DO icmp = xcoord,zcoord ! Explicit loop to avoid Frost problem
972  DO iv = 1,ppatch%nBVert
973  ppatch%dXyz(icmp,iv) = 0.0_rfreal
974  END DO ! iv
975  END DO ! icmp
976  END DO ! iPatch
977  END IF ! pGrid%nPatches
978 #endif
979  END IF ! pMixtInput%moveGrid
980 
981 #ifdef STATS
982 ! ==============================================================================
983 ! Time averaged statistics
984 ! ==============================================================================
985 
986  IF ( (global%flowType == flow_unsteady) .AND. &
987  (global%doStat == active) ) THEN
988  ALLOCATE(pregion%mixt%tav(global%mixtNStat,pgrid%nCellsTot),stat=errorflag)
989  global%error = errorflag
990  IF ( global%error /= err_none ) THEN
991  CALL errorstop(global,err_allocate,__line__,'pRegion%mixt%tav')
992  END IF ! global%error
993  ELSE
994  nullify(pregion%mixt%tav)
995  END IF ! global%flowType
996 #endif
997 
998 ! ==============================================================================
999 ! Optimal LES
1000 ! ==============================================================================
1001 
1002  IF ( pmixtinput%spaceDiscr == discr_opt_les ) THEN
1003  CALL rflu_createstencilsweightsoles(pregion)
1004  CALL rflu_createintegralsoles(pregion)
1005  END IF ! pMixtInput%spaceDiscr
1006 
1007 ! ==============================================================================
1008 ! Substantial derivative. NOTE only needed for Equilibrium Eulerian method,
1009 ! but allocate anyway, because would need IF statement inside loops otherwise.
1010 ! ==============================================================================
1011 
1012  SELECT CASE ( pmixtinput%indSd )
1013  CASE ( 0 )
1014  arraylimlow = 0
1015  arraylimupp = 1
1016  CASE ( 1 )
1017  arraylimlow = 1
1018  arraylimupp = pgrid%nCellsTot
1019  CASE default
1020  CALL errorstop(global,err_reached_default,__line__)
1021  END SELECT ! pMixtInput%indSd
1022 
1023  ALLOCATE(pregion%mixt%sd(sd_xmom:sd_zmom,arraylimlow:arraylimupp), &
1024  stat=errorflag)
1025  global%error = errorflag
1026  IF ( global%error /= err_none ) THEN
1027  CALL errorstop(global,err_allocate,__line__,'pRegion%mixt%sd')
1028  END IF ! global%error
1029 
1030 ! ******************************************************************************
1031 ! End
1032 ! ******************************************************************************
1033 
1034  CALL deregisterfunction(global)
1035 
1036 END SUBROUTINE rflu_allocatememorytstep
1037 
1038 
1039 
1040 
1041 
1042 
1043 ! ******************************************************************************
1044 !
1045 ! Purpose: Allocate memory for mixture time-stepping for compressible fluid.
1046 !
1047 ! Description: None.
1048 !
1049 ! Input:
1050 ! pRegion Region pointer
1051 !
1052 ! Output: None.
1053 !
1054 ! Notes: None.
1055 !
1056 ! ******************************************************************************
1057 
1058 SUBROUTINE rflu_allocatememorytstep_c(pRegion)
1059 
1060  USE rflu_modoles
1061 
1062  IMPLICIT NONE
1063 
1064 ! ******************************************************************************
1065 ! Definitions and declarations
1066 ! ******************************************************************************
1067 
1068 ! ==============================================================================
1069 ! Arguments
1070 ! ==============================================================================
1071 
1072  TYPE(t_region), POINTER :: pregion
1073 
1074 ! ==============================================================================
1075 ! Locals
1076 ! ==============================================================================
1077 
1078  INTEGER :: arraylimlow,arraylimupp,errorflag,ipatch
1079  TYPE(t_grid), POINTER :: pgrid
1080  TYPE(t_patch), POINTER :: ppatch
1081  TYPE(t_global), POINTER :: global
1082  TYPE(t_mixt_input), POINTER :: pmixtinput
1083 
1084 ! ******************************************************************************
1085 ! Start
1086 ! ******************************************************************************
1087 
1088  global => pregion%global
1089 
1090  CALL registerfunction(global,'RFLU_AllocateMemoryTStep_C',&
1091  'RFLU_ModAllocateMemory.F90')
1092 
1093 ! ******************************************************************************
1094 ! Set pointers
1095 ! ******************************************************************************
1096 
1097  pgrid => pregion%grid
1098  pmixtinput => pregion%mixtInput
1099 
1100 ! ******************************************************************************
1101 ! Allocate memory
1102 ! ******************************************************************************
1103 
1104 ! ==============================================================================
1105 ! Mass fluxes. NOTE only needed when solving additional scalar transport
1106 ! equations, but allocate anyway, because would need IF statement inside flux
1107 ! loops otherwise.
1108 ! ==============================================================================
1109 
1110  IF ( pmixtinput%indMfMixt == 1 ) THEN
1111  arraylimlow = 1
1112  arraylimupp = pgrid%nFaces
1113  ELSE
1114  arraylimlow = 0
1115  arraylimupp = 1
1116  END IF ! pMixtInput%indMfMixt
1117 
1118  ALLOCATE(pregion%mixt%mfMixt(arraylimlow:arraylimupp),stat=errorflag)
1119  global%error = errorflag
1120  IF ( global%error /= err_none ) THEN
1121  CALL errorstop(global,err_allocate,__line__,'pRegion%mixt%mfMixt')
1122  END IF ! global%error
1123 
1124  IF ( pregion%grid%nPatches > 0 ) THEN
1125  DO ipatch = 1,pregion%grid%nPatches
1126  ppatch => pregion%patches(ipatch)
1127 
1128  IF ( pmixtinput%indMfMixt == 1 ) THEN
1129  arraylimlow = 1
1130  arraylimupp = ppatch%nBFaces
1131  ELSE
1132  arraylimlow = 0
1133  arraylimupp = 1
1134  END IF ! pMixtInput%indMfMixt
1135 
1136  ALLOCATE(ppatch%mfMixt(arraylimlow:arraylimupp),stat=errorflag)
1137  global%error = errorflag
1138  IF ( global%error /= err_none ) THEN
1139  CALL errorstop(global,err_allocate,__line__,'pPatch%mfMixt')
1140  END IF ! global%error
1141  END DO ! iPatch
1142  END IF ! pRegion%grid%nPatches
1143 
1144 ! ******************************************************************************
1145 ! End
1146 ! ******************************************************************************
1147 
1148  CALL deregisterfunction(global)
1149 
1150 END SUBROUTINE rflu_allocatememorytstep_c
1151 
1152 
1153 
1154 
1155 
1156 ! ******************************************************************************
1157 !
1158 ! Purpose: Allocate memory for mixture time-stepping for incompressible fluid.
1159 !
1160 ! Description: None.
1161 !
1162 ! Input:
1163 ! pRegion Region pointer
1164 !
1165 ! Output: None.
1166 !
1167 ! Notes: None.
1168 !
1169 ! ******************************************************************************
1170 
1171 SUBROUTINE rflu_allocatememorytstep_i(pRegion)
1172 
1173  IMPLICIT NONE
1174 
1175 ! ******************************************************************************
1176 ! Definitions and declarations
1177 ! ******************************************************************************
1178 
1179 ! ==============================================================================
1180 ! Arguments
1181 ! ==============================================================================
1182 
1183  TYPE(t_region), POINTER :: pregion
1184 
1185 ! ==============================================================================
1186 ! Locals
1187 ! ==============================================================================
1188 
1189  INTEGER :: errorflag,ipatch
1190  TYPE(t_grid), POINTER :: pgrid
1191  TYPE(t_patch), POINTER :: ppatch
1192  TYPE(t_global), POINTER :: global
1193 
1194 ! ******************************************************************************
1195 ! Start
1196 ! ******************************************************************************
1197 
1198  global => pregion%global
1199 
1200  CALL registerfunction(global,'RFLU_AllocateMemoryTStep_I',&
1201  'RFLU_ModAllocateMemory.F90')
1202 
1203 ! ******************************************************************************
1204 ! Set pointers
1205 ! ******************************************************************************
1206 
1207  pgrid => pregion%grid
1208 
1209 ! ******************************************************************************
1210 ! Allocate memory
1211 ! ******************************************************************************
1212 
1213 ! ==============================================================================
1214 ! Face velocities
1215 ! ==============================================================================
1216 
1217  ALLOCATE(pregion%mixt%vfMixt(pgrid%nFaces),stat=errorflag)
1218  global%error = errorflag
1219  IF ( global%error /= err_none ) THEN
1220  CALL errorstop(global,err_allocate,__line__,'pRegion%mixt%vfMixt')
1221  END IF ! global%error
1222 
1223  IF ( pregion%grid%nPatches > 0 ) THEN
1224  DO ipatch = 1,pregion%grid%nPatches
1225  ppatch => pregion%patches(ipatch)
1226 
1227  ALLOCATE(ppatch%vfMixt(ppatch%nBFaces),stat=errorflag)
1228  global%error = errorflag
1229  IF ( global%error /= err_none ) THEN
1230  CALL errorstop(global,err_allocate,__line__,'pPatch%vfMixt')
1231  END IF ! global%error
1232  END DO ! iPatch
1233  END IF ! pRegion%grid%nPatches
1234 
1235 ! ******************************************************************************
1236 ! End
1237 ! ******************************************************************************
1238 
1239  CALL deregisterfunction(global)
1240 
1241 END SUBROUTINE rflu_allocatememorytstep_i
1242 
1243 
1244 
1245 
1246 
1247 ! ******************************************************************************
1248 ! End
1249 ! ******************************************************************************
1250 
1251 END MODULE rflu_modallocatememory
1252 
1253 
1254 ! ******************************************************************************
1255 !
1256 ! RCS Revision history:
1257 !
1258 ! $Log: RFLU_ModAllocateMemory.F90,v $
1259 ! Revision 1.16 2008/12/06 08:44:20 mtcampbe
1260 ! Updated license.
1261 !
1262 ! Revision 1.15 2008/11/19 22:17:31 mtcampbe
1263 ! Added Illinois Open Source License/Copyright
1264 !
1265 ! Revision 1.14 2006/11/01 15:50:00 haselbac
1266 ! Changed so implicit-solver arrays dealt with properly
1267 !
1268 ! Revision 1.13 2006/10/20 21:17:34 mparmar
1269 ! Fixed a bug in allocation of pPatch%mixt%gradFace
1270 !
1271 ! Revision 1.12 2006/08/19 15:38:59 mparmar
1272 ! Moved region%mixt%bGradFace to patch%mixt%gradFace
1273 !
1274 ! Revision 1.11 2006/02/08 21:03:16 hdewey2
1275 ! Added old2 quantities
1276 !
1277 ! Revision 1.10 2006/01/07 10:16:33 wasistho
1278 ! STATS tav allocation to nCellsTot
1279 !
1280 ! Revision 1.9 2005/10/31 21:09:36 haselbac
1281 ! Changed specModel and SPEC_MODEL_NONE
1282 !
1283 ! Revision 1.8 2005/09/22 17:09:21 hdewey2
1284 ! Added allocation of cvOld1 and cvOld2 for transient implicit solver.
1285 !
1286 ! Revision 1.7 2005/07/14 21:42:17 haselbac
1287 ! Added AUSM flux function to memory allocation IF statement
1288 !
1289 ! Revision 1.6 2005/03/31 17:00:19 haselbac
1290 ! Changed allocation of sd
1291 !
1292 ! Revision 1.5 2004/12/19 15:46:56 haselbac
1293 ! Added memory allocation for incompressible solver
1294 !
1295 ! Revision 1.4 2004/11/02 02:30:22 haselbac
1296 ! Removed use of CV_MIXT_NEQS and init of cvInfo
1297 !
1298 ! Revision 1.3 2004/10/19 19:27:44 haselbac
1299 ! Modified allocation criteria for grid speeds
1300 !
1301 ! Revision 1.2 2004/07/30 22:47:36 jferry
1302 ! Implemented Equilibrium Eulerian method for Rocflu
1303 !
1304 ! Revision 1.1 2004/03/19 21:15:19 haselbac
1305 ! Initial revision
1306 !
1307 ! ******************************************************************************
1308 
1309 
1310 
1311 
1312 
1313 
1314 
1315 
1316 
1317 
1318 
1319 
1320 
1321 
1322 
LOGICAL function rflu_decideneedbgradface(pRegion, pPatch)
subroutine registerfunction(global, funName, fileName)
Definition: ModError.F90:449
subroutine, public rflu_allocatememorysolcv(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 ic
LOGICAL function, public rflu_decideneedgridspeeds(pRegion)
subroutine, public rflu_allocatememorysoldv(pRegion)
subroutine, public rflu_allocatememorysoltv(pRegion)
subroutine, public rflu_createintegralsoles(pRegion)
subroutine, public rflu_allocatememorysol(pRegion)
subroutine, public rflu_allocatememorytstep_i(pRegion)
subroutine errorstop(global, errorCode, errorLine, addMessage)
Definition: ModError.F90:483
subroutine, public rflu_allocatememorygspeeds(pRegion)
subroutine deregisterfunction(global)
Definition: ModError.F90:469
subroutine, public rflu_allocatememorytstep_c(pRegion)
subroutine, public rflu_createstencilsweightsoles(pRegion)
subroutine, public rflu_allocatememorysolgv(pRegion)
subroutine, public rflu_allocatememorytstep(pRegion)