Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
RFLU_ModStencilsVert.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 construct vertex stencils.
26 !
27 ! Description: None.
28 !
29 ! Notes: None.
30 !
31 ! ******************************************************************************
32 !
33 ! $Id: RFLU_ModStencilsVert.F90,v 1.5 2008/12/06 08:44:24 mtcampbe Exp $
34 !
35 ! Copyright: (c) 2005 by the University of Illinois
36 !
37 ! ******************************************************************************
38 
40 
41  USE modglobal, ONLY: t_global
42  USE modparameters
43  USE moddatatypes
44  USE moderror
45  USE modbndpatch, ONLY: t_patch
46  USE moddatastruct, ONLY: t_region
47  USE modgrid, ONLY: t_grid
48  USE modstencil, ONLY: t_stencil
49  USE modsortsearch
50  USE modmpi
51 
54 
55  IMPLICIT NONE
56 
57  PRIVATE
58  PUBLIC :: rflu_buildstencilvert2cell, &
63 
64 
65 ! ******************************************************************************
66 ! Declarations and definitions
67 ! ******************************************************************************
68 
69  CHARACTER(CHRLEN) :: RCSIdentString = &
70  '$RCSfile: RFLU_ModStencilsVert.F90,v $ $Revision: 1.5 $'
71 
72 ! ******************************************************************************
73 ! Routines
74 ! ******************************************************************************
75 
76  CONTAINS
77 
78 
79 
80 
81 
82 ! *******************************************************************************
83 !
84 ! Purpose: Build vertex-to-cell stencil.
85 !
86 ! Description: None.
87 !
88 ! Input:
89 ! pRegion Pointer to region
90 !
91 ! Output: None.
92 !
93 ! Notes: None.
94 !
95 ! ******************************************************************************
96 
97  SUBROUTINE rflu_buildstencilvert2cell(pRegion)
98 
99  IMPLICIT NONE
100 
101 ! ******************************************************************************
102 ! Declarations and definitions
103 ! ******************************************************************************
104 
105 ! ==============================================================================
106 ! Arguments
107 ! ==============================================================================
108 
109  TYPE(t_region), POINTER :: pregion
110 
111 ! ==============================================================================
112 ! Locals
113 ! ==============================================================================
114 
115  INTEGER :: degr,errorflag,icg,icl,ict,ilayer,iloc,isl,ivg,ivl,iv2c, &
116  nbfacemembsmax,nbfacemembsmaxtemp,ncellmembsinfomax, &
117  ncellmembsinfomaxloc,ncellmembsinfomin,ncellmembsinfominloc, &
118  nlayersinfomax,nlayersinfomaxloc,nlayersinfomin, &
119  nlayersinfominloc,nlayersmax,nrows,order,ordernominal,scount, &
120  stencilsizemax,stencilsizemin,v2csbeg,v2csend
121  INTEGER, DIMENSION(:), ALLOCATABLE :: v2cs
122  INTEGER, DIMENSION(:,:), ALLOCATABLE :: layerinfo
123  REAL(RFREAL), DIMENSION(:,:), ALLOCATABLE :: dr,wts
124  TYPE(t_grid), POINTER :: pgrid
125  TYPE(t_global), POINTER :: global
126 
127 ! ******************************************************************************
128 ! Start
129 ! ******************************************************************************
130 
131  global => pregion%global
132 
133  CALL registerfunction(global,'RFLU_BuildStencilVert2Cell',&
134  'RFLU_ModStencilsVert.F90')
135 
136  IF ( global%myProcid == masterproc .AND. &
137  global%verbLevel > verbose_none ) THEN
138  WRITE(stdout,'(A,1X,A)') solver_name, &
139  'Building vertex-to-cell stencil...'
140  END IF ! global%verbLevel
141 
142 ! ******************************************************************************
143 ! Set grid pointer
144 ! ******************************************************************************
145 
146  pgrid => pregion%grid
147 
148 ! ******************************************************************************
149 ! Set variables
150 ! ******************************************************************************
151 
152  ordernominal = pgrid%v2csInfo%orderNominal
153  nlayersmax = pgrid%v2csInfo%nLayersMax
154  nbfacemembsmax = pgrid%v2csInfo%nBFaceMembsMax
155  stencilsizemax = pgrid%v2csInfo%nCellMembsMax
156  stencilsizemin = pgrid%v2csInfo%nCellMembsMin
157 
158  ncellmembsinfomax = 0
159  ncellmembsinfomin = huge(1)
160 
161  nlayersinfomax = 0
162  nlayersinfomin = huge(1)
163 
164  nbfacemembsmaxtemp = 2*nbfacemembsmax
165 
166 ! ******************************************************************************
167 ! Allocate temporary memory
168 ! ******************************************************************************
169 
170  ALLOCATE(v2cs(stencilsizemax),stat=errorflag)
171  global%error = errorflag
172  IF ( global%error /= err_none ) THEN
173  CALL errorstop(global,err_allocate,__line__,'v2cs')
174  END IF ! global%error
175 
176  ALLOCATE(layerinfo(x2cs_layer_beg:x2cs_layer_end,nlayersmax), &
177  stat=errorflag)
178  global%error = errorflag
179  IF ( global%error /= err_none ) THEN
180  CALL errorstop(global,err_allocate,__line__,'layerInfo')
181  END IF ! global%error
182 
183 ! ******************************************************************************
184 ! Loop over vertices
185 ! ******************************************************************************
186 
187  DO ivg = 1,pgrid%nVertTot
188 
189 ! ==============================================================================
190 ! Initialize
191 ! ==============================================================================
192 
193  degr = 0
194 
195  DO isl = 1,stencilsizemax
196  v2cs(isl) = 0
197  END DO ! isl
198 
199  DO ilayer = 1,nlayersmax
200  layerinfo(x2cs_layer_beg,ilayer) = 0
201  layerinfo(x2cs_layer_end,ilayer) = 0
202  END DO ! iLayer
203 
204  pgrid%v2cs(ivg)%nLayers = 1
205 
206 ! ==============================================================================
207 ! Build basic stencil
208 ! ==============================================================================
209 
210  DO iv2c = pgrid%v2cInfo(v2c_beg,ivg),pgrid%v2cInfo(v2c_end,ivg)
211  degr = degr + 1
212  v2cs(degr) = pgrid%v2c(iv2c)
213  END DO ! iv2c
214 
215  layerinfo(x2cs_layer_beg,1) = 1
216  layerinfo(x2cs_layer_end,1) = degr
217 
218 ! ==============================================================================
219 ! Extend basic stencil
220 ! ==============================================================================
221 
222  DO ilayer = 2,nlayersmax
223  order = ordernominal
224  scount = 0
225 
226 ! ------------------------------------------------------------------------------
227 ! Check whether stencil weights are singular
228 ! ------------------------------------------------------------------------------
229 
230  IF ( degr >= 4 ) THEN
231  nrows = degr
232 
233  ALLOCATE(dr(xcoord:zcoord,nrows),stat=errorflag)
234  global%error = errorflag
235  IF ( global%error /= err_none ) THEN
236  CALL errorstop(global,err_allocate,__line__,'dr')
237  END IF ! global%error
238 
239  ALLOCATE(wts(1,nrows),stat=errorflag)
240  global%error = errorflag
241  IF ( global%error /= err_none ) THEN
242  CALL errorstop(global,err_allocate,__line__,'wts')
243  END IF ! global%error
244 
245  DO isl = 1,degr
246  icg = v2cs(isl)
247 
248  dr(xcoord,isl) = pgrid%cofg(xcoord,icg) - pgrid%xyz(xcoord,ivg)
249  dr(ycoord,isl) = pgrid%cofg(ycoord,icg) - pgrid%xyz(ycoord,ivg)
250  dr(zcoord,isl) = pgrid%cofg(zcoord,icg) - pgrid%xyz(zcoord,ivg)
251  END DO ! isl
252 
253  CALL rflu_computestencilweights(global,pregion%mixtInput%dimens, &
254  compwts_mode_fixed, &
255  compwts_scal_invdist,deriv_degree_0, &
256  order,nrows,dr,wts,scount)
257 
258  DEALLOCATE(dr,stat=errorflag)
259  global%error = errorflag
260  IF ( global%error /= err_none ) THEN
261  CALL errorstop(global,err_deallocate,__line__,'dr')
262  END IF ! global%error
263 
264  DEALLOCATE(wts,stat=errorflag)
265  global%error = errorflag
266  IF ( global%error /= err_none ) THEN
267  CALL errorstop(global,err_deallocate,__line__,'wts')
268  END IF ! global%error
269  END IF ! degr
270 
271 ! ------------------------------------------------------------------------------
272 ! Check whether to reject or accept stencil. If singular or too small,
273 ! add layer of cells. Pass 0 instead of ivg because want to prevent
274 ! present cell from being added, not present vertex. If pass present
275 ! vertex, could actually prevent a proper cell from being added in
276 ! rare cases...
277 ! ------------------------------------------------------------------------------
278 
279  IF ( scount /= 0 .OR. degr <= stencilsizemin ) THEN
280  v2csbeg = layerinfo(x2cs_layer_beg,ilayer-1)
281  v2csend = layerinfo(x2cs_layer_end,ilayer-1)
282 
283  CALL rflu_addcelllayer(global,pgrid,stencilsizemax,0,degr, &
284  v2csbeg,v2csend,v2cs)
285 
286  pgrid%v2cs(ivg)%nLayers = pgrid%v2cs(ivg)%nLayers + 1
287 
288  layerinfo(x2cs_layer_beg,ilayer) = &
289  layerinfo(x2cs_layer_end,ilayer-1) + 1
290  layerinfo(x2cs_layer_end,ilayer) = degr
291  ELSE
292  EXIT
293  END IF ! sCount
294  END DO ! iLayer
295 
296 ! ==============================================================================
297 ! Store stencil
298 ! ==============================================================================
299 
300  pgrid%v2cs(ivg)%nCellMembs = degr
301 
302  ALLOCATE(pgrid%v2cs(ivg)%cellMembs(pgrid%v2cs(ivg)%nCellMembs), &
303  stat=errorflag)
304  global%error = errorflag
305  IF ( global%error /= err_none ) THEN
306  CALL errorstop(global,err_allocate,__line__,'pGrid%v2cs%cellMembs')
307  END IF ! global%error
308 
309  DO isl = 1,pgrid%v2cs(ivg)%nCellMembs
310  pgrid%v2cs(ivg)%cellMembs(isl) = v2cs(isl)
311  END DO ! isl
312 
313  ALLOCATE(pgrid%v2cs(ivg)%layerInfo(x2cs_layer_beg:x2cs_layer_end, &
314  pgrid%v2cs(ivg)%nLayers),stat=errorflag)
315  global%error = errorflag
316  IF ( global%error /= err_none ) THEN
317  CALL errorstop(global,err_allocate,__line__,'pGrid%v2cs%layerInfo')
318  END IF ! global%error
319 
320  DO ilayer = 1,pgrid%v2cs(ivg)%nLayers
321  pgrid%v2cs(ivg)%layerInfo(x2cs_layer_beg,ilayer) = &
322  layerinfo(x2cs_layer_beg,ilayer)
323  pgrid%v2cs(ivg)%layerInfo(x2cs_layer_end,ilayer) = &
324  layerinfo(x2cs_layer_end,ilayer)
325  END DO ! iLayer
326 
327 ! ==============================================================================
328 ! Extract information for later printing
329 ! ==============================================================================
330 
331  IF ( pgrid%v2cs(ivg)%nLayers < nlayersinfomin ) THEN
332  nlayersinfomin = pgrid%v2cs(ivg)%nLayers
333  nlayersinfominloc = ivg
334  END IF ! pGrid%v2cs(ivg)%nLayers
335 
336  IF ( pgrid%v2cs(ivg)%nLayers > nlayersinfomax ) THEN
337  nlayersinfomax = pgrid%v2cs(ivg)%nLayers
338  nlayersinfomaxloc = ivg
339  END IF ! pGrid%v2cs(ivg)%nLayers
340 
341  IF ( pgrid%v2cs(ivg)%nCellMembs < ncellmembsinfomin ) THEN
342  ncellmembsinfomin = pgrid%v2cs(ivg)%nCellMembs
343  ncellmembsinfominloc = ivg
344  END IF ! pGrid%v2cs(ivg)%nCellMembs
345 
346  IF ( pgrid%v2cs(ivg)%nCellMembs > ncellmembsinfomax ) THEN
347  ncellmembsinfomax = pgrid%v2cs(ivg)%nCellMembs
348  ncellmembsinfomaxloc = ivg
349  END IF ! pGrid%v2cs(ivg)%nCellMembs
350  END DO ! ivg
351 
352 ! ******************************************************************************
353 ! Write out information on stencils
354 ! ******************************************************************************
355 
356  IF ( global%myProcid == masterproc .AND. &
357  global%verbLevel > verbose_low ) THEN
358  WRITE(stdout,'(A,3X,A)') solver_name,'Statistics:'
359  WRITE(stdout,'(A,5X,A,2(1X,I3),2(1X,I9))') solver_name, &
360  'Minimum/maximum number of cell layers: ',nlayersinfomin, &
361  nlayersinfomax,nlayersinfominloc,nlayersinfomaxloc
362  WRITE(stdout,'(A,5X,A,2(1X,I3),2(1X,I9))') solver_name, &
363  'Minimum/maximum number of cell members:',ncellmembsinfomin, &
364  ncellmembsinfomax,ncellmembsinfominloc,ncellmembsinfomaxloc
365  END IF ! global%myProcid
366 
367 ! ******************************************************************************
368 ! Deallocate temporary memory
369 ! ******************************************************************************
370 
371  DEALLOCATE(v2cs,stat=errorflag)
372  global%error = errorflag
373  IF ( global%error /= err_none ) THEN
374  CALL errorstop(global,err_deallocate,__line__,'v2cs')
375  END IF ! global%error
376 
377  DEALLOCATE(layerinfo,stat=errorflag)
378  global%error = errorflag
379  IF ( global%error /= err_none ) THEN
380  CALL errorstop(global,err_deallocate,__line__,'layerInfo')
381  END IF ! global%error
382 
383 #ifdef CHECK_DATASTRUCT
384 ! ******************************************************************************
385 ! Data structure output for checking
386 ! ******************************************************************************
387 
388  WRITE(stdout,'(A)') solver_name
389  WRITE(stdout,'(A,1X,A)') solver_name,'### START CHECK OUTPUT ###'
390  WRITE(stdout,'(A,1X,A)') solver_name,'Vertex-to-cell stencils'
391  WRITE(stdout,'(A,1X,A,1X,I6)') solver_name,'Maximum number of layers:', &
392  pgrid%v2csInfo%nLayersMax
393  WRITE(stdout,'(A,1X,A,1X,I6)') solver_name,'Minimum stencil size:', &
394  pgrid%v2csInfo%nCellMembsMin
395 
396  DO ivg = 1,pgrid%nVertTot
397  WRITE(stdout,'(A,1X,I6,2(1X,I3),3X,20(1X,I6))') solver_name,ivg, &
398  pgrid%v2cs(ivg)%nLayers,pgrid%v2cs(ivg)%nCellMembs, &
399  pgrid%v2cs(ivg)%cellMembs(1:pgrid%v2cs(ivg)%nCellMembs)
400  END DO ! ivg
401 
402  WRITE(stdout,'(A,1X,A)') solver_name,'### END CHECK OUTPUT ###'
403  WRITE(stdout,'(A)') solver_name
404 #endif
405 
406 ! ******************************************************************************
407 ! End
408 ! ******************************************************************************
409 
410  IF ( global%myProcid == masterproc .AND. &
411  global%verbLevel > verbose_none ) THEN
412  WRITE(stdout,'(A,1X,A)') solver_name, &
413  'Building vertex-to-cell stencil done.'
414  END IF ! global%verbLevel
415 
416  CALL deregisterfunction(global)
417 
418  END SUBROUTINE rflu_buildstencilvert2cell
419 
420 
421 
422 
423 
424 
425 
426 ! *******************************************************************************
427 !
428 ! Purpose: Create vertex-to-cell stencil.
429 !
430 ! Description: None.
431 !
432 ! Input:
433 ! pRegion Pointer to region
434 !
435 ! Output: None.
436 !
437 ! Notes: None.
438 !
439 ! ******************************************************************************
440 
441  SUBROUTINE rflu_createstencilvert2cell(pRegion)
442 
443  IMPLICIT NONE
444 
445 ! ******************************************************************************
446 ! Declarations and definitions
447 ! ******************************************************************************
448 
449 ! ==============================================================================
450 ! Arguments
451 ! ==============================================================================
452 
453  TYPE(t_region), POINTER :: pregion
454 
455 ! ==============================================================================
456 ! Locals
457 ! ==============================================================================
458 
459  INTEGER :: errorflag,ivg
460  TYPE(t_grid), POINTER :: pgrid
461  TYPE(t_global), POINTER :: global
462 
463 ! ******************************************************************************
464 ! Start
465 ! ******************************************************************************
466 
467  global => pregion%global
468 
469  CALL registerfunction(global,'RFLU_BuildStencilVert2Cell',&
470  'RFLU_ModStencilsVert.F90')
471 
472  IF ( global%myProcid == masterproc .AND. &
473  global%verbLevel > verbose_none ) THEN
474  WRITE(stdout,'(A,1X,A)') solver_name, &
475  'Creating vertex-to-cell stencil...'
476  END IF ! global%verbLevel
477 
478 ! ******************************************************************************
479 ! Nullify memory
480 ! ******************************************************************************
481 
482  CALL rflu_nullifystencilvert2cell(pregion)
483 
484 ! ******************************************************************************
485 ! Set grid pointer
486 ! ******************************************************************************
487 
488  pgrid => pregion%grid
489 
490 ! ******************************************************************************
491 ! Allocate memory and initialize
492 ! ******************************************************************************
493 
494  ALLOCATE(pgrid%v2cs(pgrid%nVertTot),stat=errorflag)
495  global%error = errorflag
496  IF ( global%error /= err_none ) THEN
497  CALL errorstop(global,err_allocate,__line__,'pGrid%v2cs')
498  END IF ! global%error
499 
500  DO ivg = 1,pgrid%nVertTot
501  pgrid%v2cs(ivg)%nCellMembs = 0
502  pgrid%v2cs(ivg)%nBFaceMembs = 0
503  END DO ! ivg
504 
505 ! ******************************************************************************
506 ! End
507 ! ******************************************************************************
508 
509  IF ( global%myProcid == masterproc .AND. &
510  global%verbLevel > verbose_none ) THEN
511  WRITE(stdout,'(A,1X,A)') solver_name, &
512  'Creating vertex-to-cell stencil done.'
513  END IF ! global%verbLevel
514 
515  CALL deregisterfunction(global)
516 
517  END SUBROUTINE rflu_createstencilvert2cell
518 
519 
520 
521 
522 
523 
524 
525 
526 ! *******************************************************************************
527 !
528 ! Purpose: Destroy vertex-to-cell stencil.
529 !
530 ! Description: None.
531 !
532 ! Input:
533 ! pRegion Pointer to region
534 !
535 ! Output: None.
536 !
537 ! Notes: None.
538 !
539 ! ******************************************************************************
540 
541  SUBROUTINE rflu_destroystencilvert2cell(pRegion)
542 
543  IMPLICIT NONE
544 
545 ! ******************************************************************************
546 ! Declarations and definitions
547 ! ******************************************************************************
548 
549 ! ==============================================================================
550 ! Arguments
551 ! ==============================================================================
552 
553  TYPE(t_region), POINTER :: pregion
554 
555 ! ==============================================================================
556 ! Locals
557 ! ==============================================================================
558 
559  INTEGER :: errorflag,ivg
560  TYPE(t_grid), POINTER :: pgrid
561  TYPE(t_global), POINTER :: global
562 
563 ! ******************************************************************************
564 ! Start
565 ! ******************************************************************************
566 
567  global => pregion%global
568 
569  CALL registerfunction(global,'RFLU_DestroyStencilVert2Cell',&
570  'RFLU_ModStencilsVert.F90')
571 
572  IF ( global%myProcid == masterproc .AND. &
573  global%verbLevel > verbose_none ) THEN
574  WRITE(stdout,'(A,1X,A)') solver_name, &
575  'Destroying vertex-to-cell stencil...'
576  END IF ! global%verbLevel
577 
578 ! ******************************************************************************
579 ! Set grid pointer
580 ! ******************************************************************************
581 
582  pgrid => pregion%grid
583 
584 ! ******************************************************************************
585 ! Deallocate memory
586 ! ******************************************************************************
587 
588  DO ivg = 1,pgrid%nVertTot
589  DEALLOCATE(pgrid%v2cs(ivg)%cellMembs,stat=errorflag)
590  global%error = errorflag
591  IF ( global%error /= err_none ) THEN
592  CALL errorstop(global,err_deallocate,__line__,'pGrid%v2cs%cellMembs')
593  END IF ! global%error
594  END DO ! ivg
595 
596  DEALLOCATE(pgrid%v2cs,stat=errorflag)
597  global%error = errorflag
598  IF ( global%error /= err_none ) THEN
599  CALL errorstop(global,err_deallocate,__line__,'pGrid%v2cs')
600  END IF ! global%error
601 
602 ! ******************************************************************************
603 ! Nullify memory
604 ! ******************************************************************************
605 
606  CALL rflu_nullifystencilvert2cell(pregion)
607 
608 ! ******************************************************************************
609 ! End
610 ! ******************************************************************************
611 
612  IF ( global%myProcid == masterproc .AND. &
613  global%verbLevel > verbose_none ) THEN
614  WRITE(stdout,'(A,1X,A)') solver_name, &
615  'Destroying vertex-to-cell stencil done.'
616  END IF ! global%verbLevel
617 
618  CALL deregisterfunction(global)
619 
620  END SUBROUTINE rflu_destroystencilvert2cell
621 
622 
623 
624 
625 
626 
627 ! *******************************************************************************
628 !
629 ! Purpose: Nullify vertex-to-cell stencil.
630 !
631 ! Description: None.
632 !
633 ! Input:
634 ! pRegion Pointer to region
635 !
636 ! Output: None.
637 !
638 ! Notes: None.
639 !
640 ! ******************************************************************************
641 
642  SUBROUTINE rflu_nullifystencilvert2cell(pRegion)
643 
644  IMPLICIT NONE
645 
646 ! ******************************************************************************
647 ! Declarations and definitions
648 ! ******************************************************************************
649 
650 ! ==============================================================================
651 ! Arguments
652 ! ==============================================================================
653 
654  TYPE(t_region), POINTER :: pregion
655 
656 ! ==============================================================================
657 ! Locals
658 ! ==============================================================================
659 
660  TYPE(t_grid), POINTER :: pgrid
661  TYPE(t_global), POINTER :: global
662 
663 ! ******************************************************************************
664 ! Start
665 ! ******************************************************************************
666 
667  global => pregion%global
668 
669  CALL registerfunction(global,'RFLU_NullifyStencilVert2Cell',&
670  'RFLU_ModStencilsVert.F90')
671 
672  IF ( global%myProcid == masterproc .AND. &
673  global%verbLevel > verbose_none ) THEN
674  WRITE(stdout,'(A,1X,A)') solver_name, &
675  'Nullifying vertex-to-cell stencil...'
676  END IF ! global%verbLevel
677 
678 ! ******************************************************************************
679 ! Set grid pointer
680 ! ******************************************************************************
681 
682  pgrid => pregion%grid
683 
684 ! ******************************************************************************
685 ! Nullify memory
686 ! ******************************************************************************
687 
688  nullify(pgrid%v2cs)
689 
690 ! ******************************************************************************
691 ! End
692 ! ******************************************************************************
693 
694  IF ( global%myProcid == masterproc .AND. &
695  global%verbLevel > verbose_none ) THEN
696  WRITE(stdout,'(A,1X,A)') solver_name, &
697  'Nullifying vertex-to-cell stencil done.'
698  END IF ! global%verbLevel
699 
700  CALL deregisterfunction(global)
701 
702  END SUBROUTINE rflu_nullifystencilvert2cell
703 
704 
705 
706 
707 
708 ! *******************************************************************************
709 !
710 ! Purpose: Set vertex-to-cell stencil information.
711 !
712 ! Description: None.
713 !
714 ! Input:
715 ! pRegion Pointer to region
716 ! orderNominal Nominal order of accuracy
717 !
718 ! Output: None.
719 !
720 ! Notes: None.
721 !
722 ! ******************************************************************************
723 
724  SUBROUTINE rflu_setinfostencilvert2cell(pRegion,orderNominal)
725 
726  IMPLICIT NONE
727 
728 ! ******************************************************************************
729 ! Declarations and definitions
730 ! ******************************************************************************
731 
732 ! ==============================================================================
733 ! Arguments
734 ! ==============================================================================
735 
736  INTEGER, INTENT(IN) :: ordernominal
737  TYPE(t_region), POINTER :: pregion
738 
739 ! ==============================================================================
740 ! Locals
741 ! ==============================================================================
742 
743  INTEGER :: nbfacemembsmax,nlayersmax,stencilsizemax,stencilsizemin
744  TYPE(t_grid), POINTER :: pgrid
745  TYPE(t_global), POINTER :: global
746 
747 ! ******************************************************************************
748 ! Start
749 ! ******************************************************************************
750 
751  global => pregion%global
752 
753  CALL registerfunction(global,'RFLU_SetInfoStencilVert2Cell',&
754  'RFLU_ModStencilsVert.F90')
755 
756  IF ( global%myProcid == masterproc .AND. &
757  global%verbLevel > verbose_none ) THEN
758  WRITE(stdout,'(A,1X,A)') solver_name, &
759  'Setting vert-to-cell stencil information...'
760  END IF ! global%verbLevel
761 
762 ! ******************************************************************************
763 ! Set grid pointer
764 ! ******************************************************************************
765 
766  pgrid => pregion%grid
767 
768 ! ******************************************************************************
769 ! Set stencil information. NOTE nBFaceMembsMax must be one less than the
770 ! number of unknowns (columns), otherwise LAPACK routine used for constrained
771 ! least-squares problem always gives trivial solution.
772 ! ******************************************************************************
773 
774  nlayersmax = 6
775  nbfacemembsmax = rflu_computestencilsize(global,pregion%mixtInput%dimens, &
776  1,ordernominal) - 1
777  stencilsizemin = rflu_computestencilsize(global,pregion%mixtInput%dimens, &
778  2,ordernominal)
779  stencilsizemax = 10*stencilsizemin
780 
781  pgrid%v2csInfo%orderNominal = ordernominal
782  pgrid%v2csInfo%nLayersMax = nlayersmax
783  pgrid%v2csInfo%nBFaceMembsMax = nbfacemembsmax
784  pgrid%v2csInfo%nCellMembsMax = stencilsizemax
785  pgrid%v2csInfo%nCellMembsMin = stencilsizemin
786 
787 ! ******************************************************************************
788 ! Print stencil information
789 ! ******************************************************************************
790 
791  IF ( global%myProcid == masterproc .AND. &
792  global%verbLevel > verbose_low ) THEN
793  WRITE(stdout,'(A,3X,A,1X,I3)') solver_name, &
794  'Maximum allowed number of cell layers: ',nlayersmax
795  WRITE(stdout,'(A,3X,A,1X,I3)') solver_name, &
796  'Minimum required number of cell members:',stencilsizemin
797  WRITE(stdout,'(A,3X,A,1X,I3)') solver_name, &
798  'Maximum allowed number of cell members: ',stencilsizemax
799  END IF ! global%myProcid
800 
801 ! ******************************************************************************
802 ! End
803 ! ******************************************************************************
804 
805  IF ( global%myProcid == masterproc .AND. &
806  global%verbLevel > verbose_none ) THEN
807  WRITE(stdout,'(A,1X,A)') solver_name, &
808  'Setting vertex-to-cell stencil information done.'
809  END IF ! global%verbLevel
810 
811  CALL deregisterfunction(global)
812 
813  END SUBROUTINE rflu_setinfostencilvert2cell
814 
815 
816 
817 
818 
819 ! ******************************************************************************
820 ! End
821 ! ******************************************************************************
822 
823 END MODULE rflu_modstencilsvert
824 
825 
826 ! ******************************************************************************
827 !
828 ! RCS Revision history:
829 !
830 ! $Log: RFLU_ModStencilsVert.F90,v $
831 ! Revision 1.5 2008/12/06 08:44:24 mtcampbe
832 ! Updated license.
833 !
834 ! Revision 1.4 2008/11/19 22:17:35 mtcampbe
835 ! Added Illinois Open Source License/Copyright
836 !
837 ! Revision 1.3 2006/12/15 13:26:36 haselbac
838 ! Fixed bug in format statement, found by ifort
839 !
840 ! Revision 1.2 2006/04/07 15:19:20 haselbac
841 ! Removed tabs
842 !
843 ! Revision 1.1 2005/10/05 14:33:44 haselbac
844 ! Initial revision
845 !
846 ! ******************************************************************************
847 
848 
849 
850 
851 
852 
853 
854 
855 
856 
857 
subroutine, public rflu_computestencilweights(global, dimens, wtsMode, scalMode, derivDegree, orderNominal, nRows, dr, wts, sCount)
subroutine, public rflu_addcelllayer(global, pGrid, stencilSizeMax, ixg, degr, x2csBeg, x2csEnd, x2cs)
Size order() const
Degree of the element. 1 for linear and 2 for quadratic.
subroutine, public rflu_setinfostencilvert2cell(pRegion, orderNominal)
subroutine registerfunction(global, funName, fileName)
Definition: ModError.F90:449
subroutine, public rflu_createstencilvert2cell(pRegion)
subroutine, public rflu_buildstencilvert2cell(pRegion)
subroutine, public rflu_destroystencilvert2cell(pRegion)
INTEGER function rflu_computestencilsize(global, factor, order)
subroutine, public rflu_nullifystencilvert2cell(pRegion)
subroutine errorstop(global, errorCode, errorLine, addMessage)
Definition: ModError.F90:483
subroutine deregisterfunction(global)
Definition: ModError.F90:469