Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
RFLU_ModCellMapping.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 map cells.
26 !
27 ! Description: None.
28 !
29 ! Notes: None.
30 !
31 ! ******************************************************************************
32 !
33 ! $Id: RFLU_ModCellMapping.F90,v 1.16 2008/12/06 08:44:20 mtcampbe Exp $
34 !
35 ! Copyright: (c) 2002-2005 by the University of Illinois
36 !
37 ! ******************************************************************************
38 
40 
41  USE modglobal, ONLY: t_global
42  USE moddatatypes
43  USE modparameters
44  USE moderror
45  USE moddatastruct, ONLY: t_region
46  USE modgrid, ONLY: t_grid
47  USE modmpi
48 
49  IMPLICIT NONE
50 
51  PRIVATE
52  PUBLIC :: rflu_nullifycellmapping, &
59 
60 ! ******************************************************************************
61 ! Declarations and definitions
62 ! ******************************************************************************
63 
64  CHARACTER(CHRLEN) :: &
65  RCSIdentString = '$RCSfile: RFLU_ModCellMapping.F90,v $ $Revision: 1.16 $'
66 
67 ! ******************************************************************************
68 ! Procedures
69 ! ******************************************************************************
70 
71  CONTAINS
72 
73 
74 
75 
76 
77 
78 ! ******************************************************************************
79 !
80 ! Purpose: Nullify cell mapping.
81 !
82 ! Description: None.
83 !
84 ! Input:
85 ! pRegion Pointer to region
86 !
87 ! Output: None.
88 !
89 ! Notes: None.
90 !
91 ! ******************************************************************************
92 
93  SUBROUTINE rflu_nullifycellmapping(pRegion)
94 
95  IMPLICIT NONE
96 
97 ! ******************************************************************************
98 ! Declarations and definitions
99 ! ******************************************************************************
100 
101 ! ==============================================================================
102 ! Arguments
103 ! ==============================================================================
104 
105  TYPE(t_region), POINTER :: pregion
106 
107 ! ==============================================================================
108 ! Locals
109 ! ==============================================================================
110 
111  TYPE(t_grid), POINTER :: pgrid
112  TYPE(t_global), POINTER :: global
113 
114 ! ******************************************************************************
115 ! Start
116 ! ******************************************************************************
117 
118  global => pregion%global
119 
120  CALL registerfunction(global,'RFLU_NullifyCellMapping',&
121  'RFLU_ModCellMapping.F90')
122 
123  IF ( global%myProcid == masterproc .AND. &
124  global%verbLevel >= verbose_high ) THEN
125  WRITE(stdout,'(A,1X,A)') solver_name,'Nullifying cell mapping...'
126  END IF ! global%verbLevel
127 
128 ! ******************************************************************************
129 ! Set grid pointer
130 ! ******************************************************************************
131 
132  pgrid => pregion%grid
133 
134 ! ******************************************************************************
135 ! Nullify memory
136 ! ******************************************************************************
137 
138  nullify(pgrid%cellGlob2Loc)
139 
140  nullify(pgrid%tet2CellGlob)
141  nullify(pgrid%hex2CellGlob)
142  nullify(pgrid%pri2CellGlob)
143  nullify(pgrid%pyr2CellGlob)
144 
145 ! ******************************************************************************
146 ! End
147 ! ******************************************************************************
148 
149  IF ( global%myProcid == masterproc .AND. &
150  global%verbLevel >= verbose_high ) THEN
151  WRITE(stdout,'(A,1X,A)') solver_name,'Nullifying cell mapping done.'
152  END IF ! global%verbLevel
153 
154  CALL deregisterfunction(global)
155 
156  END SUBROUTINE rflu_nullifycellmapping
157 
158 
159 
160 
161 
162 
163 
164 
165 
166 ! ******************************************************************************
167 !
168 ! Purpose: Create cell mapping.
169 !
170 ! Description: None.
171 !
172 ! Input:
173 ! pRegion Pointer to region
174 !
175 ! Output: None.
176 !
177 ! Notes: None.
178 !
179 ! ******************************************************************************
180 
181  SUBROUTINE rflu_createcellmapping(pRegion)
182 
183  IMPLICIT NONE
184 
185 ! ******************************************************************************
186 ! Declarations and definitions
187 ! ******************************************************************************
188 
189 ! ==============================================================================
190 ! Arguments
191 ! ==============================================================================
192 
193  TYPE(t_region), POINTER :: pregion
194 
195 ! ==============================================================================
196 ! Locals
197 ! ==============================================================================
198 
199  INTEGER :: errorflag
200  TYPE(t_grid), POINTER :: pgrid
201  TYPE(t_global), POINTER :: global
202 
203 ! ******************************************************************************
204 ! Start
205 ! ******************************************************************************
206 
207  global => pregion%global
208 
209  CALL registerfunction(global,'RFLU_CreateCellMapping',&
210  'RFLU_ModCellMapping.F90')
211 
212  IF ( global%myProcid == masterproc .AND. &
213  global%verbLevel >= verbose_high ) THEN
214  WRITE(stdout,'(A,1X,A)') solver_name,'Creating cell mapping...'
215  END IF ! global%myProcid
216 
217 ! ******************************************************************************
218 ! Nullify
219 ! ******************************************************************************
220 
221  CALL rflu_nullifycellmapping(pregion)
222 
223 ! ******************************************************************************
224 ! Set grid pointer
225 ! ******************************************************************************
226 
227  pgrid => pregion%grid
228 
229 ! ******************************************************************************
230 ! Allocate memory
231 ! ******************************************************************************
232 
233  ALLOCATE(pgrid%cellGlob2Loc(2,pgrid%nCellsMax),stat=errorflag)
234  global%error = errorflag
235  IF ( global%error /= err_none ) THEN
236  CALL errorstop(global,err_allocate,__line__,'pGrid%cellGlob2Loc')
237  END IF ! global%error
238 
239  IF ( pgrid%nTetsMax > 0 ) THEN
240  ALLOCATE(pgrid%tet2CellGlob(pgrid%nTetsMax),stat=errorflag)
241  global%error = errorflag
242  IF ( global%error /= err_none ) THEN
243  CALL errorstop(global,err_allocate,__line__,'pGrid%tet2CellGlob')
244  END IF ! global%error
245  END IF ! pGrid%nTetsMax
246 
247  IF ( pgrid%nHexsMax > 0 ) THEN
248  ALLOCATE(pgrid%hex2CellGlob(pgrid%nHexsMax),stat=errorflag)
249  global%error = errorflag
250  IF ( global%error /= err_none ) THEN
251  CALL errorstop(global,err_allocate,__line__,'pGrid%hex2CellGlob')
252  END IF ! global%error
253  END IF ! pGrid%nHexsMax
254 
255  IF ( pgrid%nPrisMax > 0 ) THEN
256  ALLOCATE(pgrid%pri2CellGlob(pgrid%nPrisMax),stat=errorflag)
257  global%error = errorflag
258  IF ( global%error /= err_none ) THEN
259  CALL errorstop(global,err_allocate,__line__,'pGrid%pri2CellGlob')
260  END IF ! global%error
261  END IF ! pGrid%nPrisMax
262 
263  IF ( pgrid%nPyrsMax > 0 ) THEN
264  ALLOCATE(pgrid%pyr2CellGlob(pgrid%nPyrsMax),stat=errorflag)
265  global%error = errorflag
266  IF ( global%error /= err_none ) THEN
267  CALL errorstop(global,err_allocate,__line__,'pGrid%pyr2CellGlob')
268  END IF ! global%error
269  END IF ! pGrid%nPyrsMax
270 
271 ! ******************************************************************************
272 ! End
273 ! ******************************************************************************
274 
275  IF ( global%myProcid == masterproc .AND. &
276  global%verbLevel >= verbose_high ) THEN
277  WRITE(stdout,'(A,1X,A)') solver_name,'Creating cell mappings done.'
278  END IF ! global%myProcid
279 
280  CALL deregisterfunction(global)
281 
282  END SUBROUTINE rflu_createcellmapping
283 
284 
285 
286 
287 
288 
289 
290 
291 
292 ! ******************************************************************************
293 !
294 ! Purpose: Build global-to-local cell mapping.
295 !
296 ! Description: None.
297 !
298 ! Input:
299 ! pRegion Pointer to region
300 !
301 ! Output: None.
302 !
303 ! Notes: None.
304 !
305 ! ******************************************************************************
306 
307  SUBROUTINE rflu_buildglob2loccellmapping(pRegion)
308 
309  IMPLICIT NONE
310 
311 ! ******************************************************************************
312 ! Declarations and definitions
313 ! ******************************************************************************
314 
315 ! ==============================================================================
316 ! Arguments
317 ! ==============================================================================
318 
319  TYPE(t_region), POINTER :: pregion
320 
321 ! ==============================================================================
322 ! Locals
323 ! ==============================================================================
324 
325  CHARACTER(CHRLEN) :: errorstring
326  INTEGER :: icl,icg
327  TYPE(t_grid), POINTER :: pgrid
328  TYPE(t_global), POINTER :: global
329 
330 ! ******************************************************************************
331 ! Start
332 ! ******************************************************************************
333 
334  global => pregion%global
335 
336  CALL registerfunction(global,'RFLU_BuildGlob2LocCellMapping',&
337  'RFLU_ModCellMapping.F90')
338 
339  IF ( global%myProcid == masterproc .AND. &
340  global%verbLevel >= verbose_high ) THEN
341  WRITE(stdout,'(A,1X,A)') solver_name, &
342  'Building global-to-local cell mapping...'
343  END IF ! global%verbLevel
344 
345 ! ******************************************************************************
346 ! Set grid pointer
347 ! ******************************************************************************
348 
349  pgrid => pregion%grid
350 
351 ! ******************************************************************************
352 ! Build global-to-local cell mapping
353 ! ******************************************************************************
354 
355 ! ==============================================================================
356 ! Tetrahedra
357 ! ==============================================================================
358 
359  IF ( pgrid%nTetsTot > 0 ) THEN
360  IF ( global%myProcid == masterproc .AND. &
361  global%verbLevel >= verbose_high ) THEN
362  WRITE(stdout,'(A,3X,A)') solver_name,'Tetrahedra...'
363  END IF ! global%myProcid
364  END IF ! pGrid%nTetsTot
365 
366  DO icl = 1,pgrid%nTetsTot
367  icg = pgrid%tet2CellGlob(icl)
368 
369  pgrid%cellGlob2Loc(1,icg) = cell_type_tet
370  pgrid%cellGlob2Loc(2,icg) = icl
371  END DO ! icl
372 
373 ! ==============================================================================
374 ! Hexahedra
375 ! ==============================================================================
376 
377  IF ( pgrid%nHexsTot > 0 ) THEN
378  IF ( global%myProcid == masterproc .AND. &
379  global%verbLevel >= verbose_high ) THEN
380  WRITE(stdout,'(A,3X,A)') solver_name,'Hexahedra...'
381  END IF ! global%myProcid
382  END IF ! pGrid%nHexsTot
383 
384  DO icl = 1,pgrid%nHexsTot
385  icg = pgrid%hex2CellGlob(icl)
386 
387  pgrid%cellGlob2Loc(1,icg) = cell_type_hex
388  pgrid%cellGlob2Loc(2,icg) = icl
389  END DO ! icl
390 
391 ! ==============================================================================
392 ! Prisms
393 ! ==============================================================================
394 
395  IF ( pgrid%nPrisTot > 0 ) THEN
396  IF ( global%myProcid == masterproc .AND. &
397  global%verbLevel >= verbose_high ) THEN
398  WRITE(stdout,'(A,3X,A)') solver_name,'Prisms...'
399  END IF ! global%myProcid
400  END IF ! pGrid%nPrisTot
401 
402  DO icl = 1,pgrid%nPrisTot
403  icg = pgrid%pri2CellGlob(icl)
404 
405  pgrid%cellGlob2Loc(1,icg) = cell_type_pri
406  pgrid%cellGlob2Loc(2,icg) = icl
407  END DO ! icl
408 
409 ! ==============================================================================
410 ! Pyramids
411 ! ==============================================================================
412 
413  IF ( pgrid%nPyrsTot > 0 ) THEN
414  IF ( global%myProcid == masterproc .AND. &
415  global%verbLevel >= verbose_high ) THEN
416  WRITE(stdout,'(A,3X,A)') solver_name,'Pyramids...'
417  END IF ! global%myProcid
418  END IF ! pGrid%nPyrsTot
419 
420  DO icl = 1,pgrid%nPyrsTot
421  icg = pgrid%pyr2CellGlob(icl)
422 
423  pgrid%cellGlob2Loc(1,icg) = cell_type_pyr
424  pgrid%cellGlob2Loc(2,icg) = icl
425  END DO ! icl
426 
427 ! ******************************************************************************
428 ! End
429 ! ******************************************************************************
430 
431  IF ( global%myProcid == masterproc .AND. &
432  global%verbLevel >= verbose_high ) THEN
433  WRITE(stdout,'(A,1X,A)') solver_name, &
434  'Building global-to-local cell mapping done.'
435  END IF ! global%myProcid
436 
437  CALL deregisterfunction(global)
438 
439  END SUBROUTINE rflu_buildglob2loccellmapping
440 
441 
442 
443 
444 
445 
446 
447 
448 ! ******************************************************************************
449 !
450 ! Purpose: Build local-to-global cell mapping.
451 !
452 ! Description: None.
453 !
454 ! Input:
455 ! pRegion Pointer to region
456 !
457 ! Output: None.
458 !
459 ! Notes: None.
460 !
461 ! ******************************************************************************
462 
463  SUBROUTINE rflu_buildloc2globcellmapping(pRegion)
464 
465  IMPLICIT NONE
466 
467 ! ******************************************************************************
468 ! Declarations and definitions
469 ! ******************************************************************************
470 
471 ! ==============================================================================
472 ! Arguments
473 ! ==============================================================================
474 
475  TYPE(t_region), POINTER :: pregion
476 
477 ! ==============================================================================
478 ! Locals
479 ! ==============================================================================
480 
481  CHARACTER(CHRLEN) :: errorstring
482  INTEGER :: icl,iclsumactual,iclsumvirtual,icg
483  TYPE(t_grid), POINTER :: pgrid
484  TYPE(t_global), POINTER :: global
485 
486 ! ******************************************************************************
487 ! Start
488 ! ******************************************************************************
489 
490  global => pregion%global
491 
492  CALL registerfunction(global,'RFLU_BuildLoc2GlobCellMapping',&
493  'RFLU_ModCellMapping.F90')
494 
495  IF ( global%myProcid == masterproc .AND. &
496  global%verbLevel >= verbose_high ) THEN
497  WRITE(stdout,'(A,1X,A)') solver_name, &
498  'Building local-to-global cell mapping...'
499  END IF ! global%verbLevel
500 
501 ! ******************************************************************************
502 ! Set grid pointer
503 ! ******************************************************************************
504 
505  pgrid => pregion%grid
506 
507 ! ******************************************************************************
508 ! Build local-to-global cell mapping
509 ! ******************************************************************************
510 
511  iclsumactual = 0
512  iclsumvirtual = pgrid%nCells
513 
514 ! ==============================================================================
515 ! Tetrahedra
516 ! ==============================================================================
517 
518  IF ( pgrid%nTetsTot > 0 ) THEN
519  IF ( global%myProcid == masterproc .AND. &
520  global%verbLevel >= verbose_high ) THEN
521  WRITE(stdout,'(A,3X,A)') solver_name,'Tetrahedra...'
522  END IF ! global%myProcid
523  END IF ! pGrid%nTetsTot
524 
525  DO icl = 1,pgrid%nTets
526  iclsumactual = iclsumactual + 1
527  icg = iclsumactual
528 
529  pgrid%tet2CellGlob(icl) = icg
530  END DO ! icl
531 
532  DO icl = pgrid%nTets+1,pgrid%nTetsTot
533  iclsumvirtual = iclsumvirtual + 1
534  icg = iclsumvirtual
535 
536  pgrid%tet2CellGlob(icl) = icg
537  END DO ! icl
538 
539 ! ==============================================================================
540 ! Hexahedra
541 ! ==============================================================================
542 
543  IF ( pgrid%nHexsTot > 0 ) THEN
544  IF ( global%myProcid == masterproc .AND. &
545  global%verbLevel >= verbose_high ) THEN
546  WRITE(stdout,'(A,3X,A)') solver_name,'Hexahedra...'
547  END IF ! global%myProcid
548  END IF ! pGrid%nHexsTot
549 
550  DO icl = 1,pgrid%nHexs
551  iclsumactual = iclsumactual + 1
552  icg = iclsumactual
553 
554  pgrid%hex2CellGlob(icl) = icg
555  END DO ! icl
556 
557  DO icl = pgrid%nHexs+1,pgrid%nHexsTot
558  iclsumvirtual = iclsumvirtual + 1
559  icg = iclsumvirtual
560 
561  pgrid%hex2CellGlob(icl) = icg
562  END DO ! icl
563 
564 ! ==============================================================================
565 ! Prisms
566 ! ==============================================================================
567 
568  IF ( pgrid%nPrisTot > 0 ) THEN
569  IF ( global%myProcid == masterproc .AND. &
570  global%verbLevel >= verbose_high ) THEN
571  WRITE(stdout,'(A,3X,A)') solver_name,'Prisms...'
572  END IF ! global%myProcid
573  END IF ! pGrid%nPrisTot
574 
575  DO icl = 1,pgrid%nPris
576  iclsumactual = iclsumactual + 1
577  icg = iclsumactual
578 
579  pgrid%pri2CellGlob(icl) = icg
580  END DO ! icl
581 
582  DO icl = pgrid%nPris+1,pgrid%nPrisTot
583  iclsumvirtual = iclsumvirtual + 1
584  icg = iclsumvirtual
585 
586  pgrid%pri2CellGlob(icl) = icg
587  END DO ! icl
588 
589 ! ==============================================================================
590 ! Pyramids
591 ! ==============================================================================
592 
593  IF ( pgrid%nPyrsTot > 0 ) THEN
594  IF ( global%myProcid == masterproc .AND. &
595  global%verbLevel >= verbose_high ) THEN
596  WRITE(stdout,'(A,3X,A)') solver_name,'Pyramids...'
597  END IF ! global%myProcid
598  END IF ! pGrid%nPyrsTot
599 
600  DO icl = 1,pgrid%nPyrs
601  iclsumactual = iclsumactual + 1
602  icg = iclsumactual
603 
604  pgrid%pyr2CellGlob(icl) = icg
605  END DO ! icl
606 
607  DO icl = pgrid%nPyrs+1,pgrid%nPyrsTot
608  iclsumvirtual = iclsumvirtual + 1
609  icg = iclsumvirtual
610 
611  pgrid%pyr2CellGlob(icl) = icg
612  END DO ! icl
613 
614 ! ******************************************************************************
615 ! End
616 ! ******************************************************************************
617 
618  IF ( global%myProcid == masterproc .AND. &
619  global%verbLevel >= verbose_high ) THEN
620  WRITE(stdout,'(A,1X,A)') solver_name, &
621  'Building local-to-global cell mapping done.'
622  END IF ! global%myProcid
623 
624  CALL deregisterfunction(global)
625 
626  END SUBROUTINE rflu_buildloc2globcellmapping
627 
628 
629 
630 
631 
632 
633 ! ******************************************************************************
634 !
635 ! Purpose: Destroy cell mappings.
636 !
637 ! Description: None.
638 !
639 ! Input:
640 ! pRegion Pointer to region
641 !
642 ! Output: None.
643 !
644 ! Notes: None.
645 !
646 ! ******************************************************************************
647 
648  SUBROUTINE rflu_destroycellmapping(pRegion)
649 
650  IMPLICIT NONE
651 
652 ! ******************************************************************************
653 ! Declarations and definitions
654 ! ******************************************************************************
655 
656 ! ==============================================================================
657 ! Arguments
658 ! ==============================================================================
659 
660  TYPE(t_region), POINTER :: pregion
661 
662 ! ==============================================================================
663 ! Locals
664 ! ==============================================================================
665 
666  INTEGER :: errorflag
667  TYPE(t_grid), POINTER :: pgrid
668  TYPE(t_global), POINTER :: global
669 
670 ! ******************************************************************************
671 ! Start
672 ! ******************************************************************************
673 
674  global => pregion%global
675 
676  CALL registerfunction(global,'RFLU_DestroyCellMapping',&
677  'RFLU_ModCellMapping.F90')
678 
679  IF ( global%myProcid == masterproc .AND. &
680  global%verbLevel >= verbose_high ) THEN
681  WRITE(stdout,'(A,1X,A)') solver_name,'Destroying cell mapping...'
682  END IF ! global%myProcid
683 
684 ! ******************************************************************************
685 ! Set grid pointer
686 ! ******************************************************************************
687 
688  pgrid => pregion%grid
689 
690 ! ******************************************************************************
691 ! Deallocate memory
692 ! ******************************************************************************
693 
694  DEALLOCATE(pgrid%cellGlob2Loc,stat=errorflag)
695  global%error = errorflag
696  IF ( global%error /= err_none ) THEN
697  CALL errorstop(global,err_deallocate,__line__,'pGrid%cellGlob2Loc')
698  END IF ! global%error
699 
700  IF ( pgrid%nTetsMax > 0 ) THEN
701  DEALLOCATE(pgrid%tet2CellGlob,stat=errorflag)
702  IF ( global%error /= err_none ) THEN
703  CALL errorstop(global,err_deallocate,__line__,'pGrid%tet2CellGlob')
704  END IF ! global%error
705  END IF ! pGrid%nTetsMax
706 
707  IF ( pgrid%nHexsMax > 0 ) THEN
708  DEALLOCATE(pgrid%hex2CellGlob,stat=errorflag)
709  global%error = errorflag
710  IF ( global%error /= err_none ) THEN
711  CALL errorstop(global,err_deallocate,__line__,'pGrid%hex2CellGlob')
712  END IF ! global%error
713  END IF ! pGrid%nHexsMax
714 
715  IF ( pgrid%nPrisMax > 0 ) THEN
716  DEALLOCATE(pgrid%pri2CellGlob,stat=errorflag)
717  global%error = errorflag
718  IF ( global%error /= err_none ) THEN
719  CALL errorstop(global,err_deallocate,__line__,'pGrid%pri2CellGlob')
720  END IF ! global%error
721  END IF ! pGrid%nPrisMax
722 
723  IF ( pgrid%nPyrsMax > 0 ) THEN
724  DEALLOCATE(pgrid%pyr2CellGlob,stat=errorflag)
725  global%error = errorflag
726  IF ( global%error /= err_none ) THEN
727  CALL errorstop(global,err_deallocate,__line__,'pGrid%pyr2CellGlob')
728  END IF ! global%error
729  END IF ! pGrid%nPyrsMax
730 
731 ! ******************************************************************************
732 ! Nullify cell mapping
733 ! ******************************************************************************
734 
735  CALL rflu_nullifycellmapping(pregion)
736 
737 ! ******************************************************************************
738 ! End
739 ! ******************************************************************************
740 
741  IF ( global%myProcid == masterproc .AND. &
742  global%verbLevel >= verbose_high ) THEN
743  WRITE(stdout,'(A,1X,A)') solver_name,'Destroying cell mapping done.'
744  END IF ! global%myProcid
745 
746  CALL deregisterfunction(global)
747 
748  END SUBROUTINE rflu_destroycellmapping
749 
750 
751 
752 
753 
754 
755 
756 
757 
758 
759 ! ******************************************************************************
760 !
761 ! Purpose: Read local-to-global cell mapping.
762 !
763 ! Description: None.
764 !
765 ! Input:
766 ! pRegion Pointer to region
767 !
768 ! Output: None.
769 !
770 ! Notes: None.
771 !
772 ! ******************************************************************************
773 
774  SUBROUTINE rflu_readloc2globcellmapping(pRegion)
775 
777 
778  IMPLICIT NONE
779 
780 ! ******************************************************************************
781 ! Declarations and definitions
782 ! ******************************************************************************
783 
784 ! ==============================================================================
785 ! Arguments
786 ! ==============================================================================
787 
788  TYPE(t_region), POINTER :: pregion
789 
790 ! ==============================================================================
791 ! Local variables
792 ! ==============================================================================
793 
794  INTEGER :: errorflag,icg,ifile,loopcounter,nhexstot,npristot,npyrstot, &
795  ntetstot
796  CHARACTER(CHRLEN) :: ifilename,sectionstring
797  TYPE(t_grid), POINTER :: pgrid
798  TYPE(t_global), POINTER :: global
799 
800 ! ******************************************************************************
801 ! Start
802 ! ******************************************************************************
803 
804  global => pregion%global
805 
806  CALL registerfunction(global,'RFLU_ReadLoc2GlobCellMapping',&
807  'RFLU_ModCellMapping.F90')
808 
809  IF ( global%myProcid == masterproc .AND. &
810  global%verbLevel >= verbose_high ) THEN
811  WRITE(stdout,'(A,1X,A)') solver_name, &
812  'Reading local-to-global cell mapping...'
813  WRITE(stdout,'(A,3X,A,1X,I5.5)') solver_name,'Global region:', &
814  pregion%iRegionGlobal
815  END IF ! global%myProcid
816 
817  ifile = if_cell_maps
818 
819  CALL buildfilenamebasic(global,filedest_indir,'.cmp', &
820  pregion%iRegionGlobal,ifilename)
821 
822  OPEN(ifile,file=ifilename,form="FORMATTED",status="OLD",iostat=errorflag)
823  global%error = errorflag
824  IF ( global%error /= err_none ) THEN
825  CALL errorstop(global,err_file_open,__line__,ifilename)
826  END IF ! global%error
827 
828 ! ******************************************************************************
829 ! Header and general information
830 ! ******************************************************************************
831 
832  IF ( global%myProcid == masterproc .AND. &
833  global%verbLevel >= verbose_high) THEN
834  WRITE(stdout,'(A,3X,A)') solver_name,'Header information...'
835  END IF ! global%myProcid
836 
837  READ(ifile,'(A)') sectionstring
838  IF ( trim(sectionstring) /= '# ROCFLU cell mapping file' ) THEN
839  CALL errorstop(global,err_invalid_marker,__line__,sectionstring)
840  END IF ! TRIM
841 
842 ! ******************************************************************************
843 ! Dimensions
844 ! ******************************************************************************
845 
846  pgrid => pregion%grid
847 
848  READ(ifile,'(A)') sectionstring
849  IF ( trim(sectionstring) /= '# Dimensions' ) THEN
850  CALL errorstop(global,err_invalid_marker,__line__,sectionstring)
851  END IF ! TRIM
852 
853  READ(ifile,'(4(I8))') ntetstot,nhexstot,npristot,npyrstot
854 
855 ! ******************************************************************************
856 ! Check dimensions (against those read from dimensions file)
857 ! ******************************************************************************
858 
859  IF ( ntetstot /= pgrid%nTetsTot ) THEN
860  CALL errorstop(global,err_dimens_invalid,__line__)
861  END IF ! nTetsTot
862 
863  IF ( nhexstot /= pgrid%nHexsTot ) THEN
864  CALL errorstop(global,err_dimens_invalid,__line__)
865  END IF ! nHexsTot
866 
867  IF ( npristot /= pgrid%nPrisTot ) THEN
868  CALL errorstop(global,err_dimens_invalid,__line__)
869  END IF ! nPrisTot
870 
871  IF ( npyrstot /= pgrid%nPyrsTot ) THEN
872  CALL errorstop(global,err_dimens_invalid,__line__)
873  END IF ! nPyrsTot
874 
875 ! ******************************************************************************
876 ! Rest of file
877 ! ==============================================================================
878 
879  loopcounter = 0
880 
881  DO ! set up infinite loop
882  loopcounter = loopcounter + 1
883 
884  READ(ifile,'(A)') sectionstring
885 
886  SELECT CASE ( trim(sectionstring) )
887 
888 ! ==============================================================================
889 ! Tetrahedra
890 ! ==============================================================================
891 
892  CASE ( '# Tetrahedra' )
893  IF ( global%myProcid == masterproc .AND. &
894  global%verbLevel >= verbose_high) THEN
895  WRITE(stdout,'(A,3X,A)') solver_name,'Tetrahedra...'
896  END IF ! global%myProcid
897 
898  READ(ifile,'(10(I8))') (pgrid%tet2CellGlob(icg),icg=1,pgrid%nTetsTot)
899 
900 ! ==============================================================================
901 ! Hexahedra
902 ! ==============================================================================
903 
904  CASE ( '# Hexahedra' )
905  IF ( global%myProcid == masterproc .AND. &
906  global%verbLevel >= verbose_high) THEN
907  WRITE(stdout,'(A,3X,A)') solver_name,'Hexahedra...'
908  END IF ! global%myProcid
909 
910  READ(ifile,'(10(I8))') (pgrid%hex2CellGlob(icg),icg=1,pgrid%nHexsTot)
911 
912 ! ==============================================================================
913 ! Prisms
914 ! ==============================================================================
915 
916  CASE ( '# Prisms' )
917  IF ( global%myProcid == masterproc .AND. &
918  global%verbLevel >= verbose_high) THEN
919  WRITE(stdout,'(A,3X,A)') solver_name,'Prisms...'
920  END IF ! global%myProcid
921 
922  READ(ifile,'(10(I8))') (pgrid%pri2CellGlob(icg),icg=1,pgrid%nPrisTot)
923 
924 ! ==============================================================================
925 ! Pyramids
926 ! ==============================================================================
927 
928  CASE ( '# Pyramids' )
929  IF ( global%myProcid == masterproc .AND. &
930  global%verbLevel >= verbose_high) THEN
931  WRITE(stdout,'(A,3X,A)') solver_name,'Pyramids...'
932  END IF ! global%myProcid
933 
934  READ(ifile,'(10(I8))') (pgrid%pyr2CellGlob(icg),icg=1,pgrid%nPyrsTot)
935 
936 ! ==============================================================================
937 ! End marker
938 ! ==============================================================================
939 
940  CASE ( '# End' )
941  IF ( global%myProcid == masterproc .AND. &
942  global%verbLevel >= verbose_high) THEN
943  WRITE(stdout,'(A,3X,A)') solver_name,'End marker...'
944  END IF ! global%myProcid
945 
946  EXIT
947 
948 ! ==============================================================================
949 ! Invalid section string
950 ! ==============================================================================
951 
952  CASE default
953  IF ( global%myProcid == masterproc .AND. &
954  global%verbLevel >= verbose_high) THEN
955  WRITE(stdout,'(3X,A)') sectionstring
956  END IF ! global%myProcid
957 
958  CALL errorstop(global,err_invalid_marker,__line__,sectionstring)
959  END SELECT ! TRIM
960 
961 ! ==============================================================================
962 ! Guard against infinite loop - might be unnecessary because of read errors?
963 ! ==============================================================================
964 
965  IF ( loopcounter >= limit_infinite_loop ) THEN
966  CALL errorstop(global,err_infinite_loop,__line__)
967  END IF ! loopCounter
968  END DO ! <empty>
969 
970 ! ******************************************************************************
971 ! Close file
972 ! ******************************************************************************
973 
974  CLOSE(ifile,iostat=errorflag)
975  global%error = errorflag
976  IF ( global%error /= err_none ) THEN
977  CALL errorstop(global,err_file_close,__line__,ifilename)
978  END IF ! global%error
979 
980 ! ******************************************************************************
981 ! End
982 ! ******************************************************************************
983 
984  IF ( global%myProcid == masterproc .AND. &
985  global%verbLevel >= verbose_high ) THEN
986  WRITE(stdout,'(A,1X,A)') solver_name, &
987  'Reading local-to-global cell mapping done.'
988  END IF ! global%myProcid
989 
990  CALL deregisterfunction(global)
991 
992  END SUBROUTINE rflu_readloc2globcellmapping
993 
994 
995 
996 
997 
998 
999 
1000 ! ******************************************************************************
1001 !
1002 ! Purpose: Write local-to-global cell mapping.
1003 !
1004 ! Description: None.
1005 !
1006 ! Input:
1007 ! pRegion Pointer to region
1008 !
1009 ! Output: None.
1010 !
1011 ! Notes: None.
1012 !
1013 ! ******************************************************************************
1014 
1015  SUBROUTINE rflu_writeloc2globcellmapping(pRegion)
1016 
1018 
1019  IMPLICIT NONE
1020 
1021 ! ******************************************************************************
1022 ! Declarations and definitions
1023 ! ******************************************************************************
1024 
1025 ! ==============================================================================
1026 ! Arguments
1027 ! ==============================================================================
1028 
1029  TYPE(t_region), POINTER :: pregion
1030 
1031 ! ==============================================================================
1032 ! Local variables
1033 ! ==============================================================================
1034 
1035  INTEGER :: errorflag,icg,ifile
1036  CHARACTER(CHRLEN) :: ifilename,sectionstring
1037  TYPE(t_grid), POINTER :: pgrid
1038  TYPE(t_global), POINTER :: global
1039 
1040 ! ******************************************************************************
1041 ! Start
1042 ! ******************************************************************************
1043  global => pregion%global
1044 
1045  CALL registerfunction(global,'RFLU_WriteLoc2GlobCellMapping',&
1046  'RFLU_ModCellMapping.F90')
1047 
1048  IF ( global%myProcid == masterproc .AND. &
1049  global%verbLevel >= verbose_med ) THEN
1050  WRITE(stdout,'(A,1X,A)') solver_name, &
1051  'Writing local-to-global cell mapping...'
1052  WRITE(stdout,'(A,3X,A,1X,I5.5)') solver_name,'Global region:', &
1053  pregion%iRegionGlobal
1054  END IF ! global%myProcid
1055 
1056  ifile = if_cell_maps
1057 
1058  CALL buildfilenamebasic(global,filedest_indir,'.cmp', &
1059  pregion%iRegionGlobal,ifilename)
1060 
1061  OPEN(ifile,file=ifilename,form="FORMATTED",status="UNKNOWN", &
1062  iostat=errorflag)
1063  global%error = errorflag
1064  IF ( global%error /= err_none ) THEN
1065  CALL errorstop(global,err_file_open,__line__,ifilename)
1066  END IF ! global%error
1067 
1068 ! ******************************************************************************
1069 ! Header and general information
1070 ! ******************************************************************************
1071 
1072  IF ( global%myProcid == masterproc .AND. &
1073  global%verbLevel >= verbose_high) THEN
1074  WRITE(stdout,'(A,3X,A)') solver_name,'Header information...'
1075  END IF ! global%myProcid
1076 
1077  sectionstring = '# ROCFLU cell mapping file'
1078  WRITE(ifile,'(A)') trim(sectionstring)
1079 
1080 ! ******************************************************************************
1081 ! Dimensions
1082 ! ******************************************************************************
1083 
1084  pgrid => pregion%grid
1085 
1086  sectionstring = '# Dimensions'
1087  WRITE(ifile,'(A)') trim(sectionstring)
1088  WRITE(ifile,'(4(I8))') pgrid%nTetsTot,pgrid%nHexsTot,pgrid%nPrisTot, &
1089  pgrid%nPyrsTot
1090 
1091 ! ******************************************************************************
1092 ! Tetrahedra
1093 ! ******************************************************************************
1094 
1095  IF ( pgrid%nTetsTot > 0 ) THEN
1096  IF ( global%myProcid == masterproc .AND. &
1097  global%verbLevel >= verbose_high) THEN
1098  WRITE(stdout,'(A,3X,A)') solver_name,'Tetrahedra...'
1099  END IF ! global%myProcid
1100 
1101  sectionstring = '# Tetrahedra'
1102  WRITE(ifile,'(A)') trim(sectionstring)
1103  WRITE(ifile,'(10(I8))') (pgrid%tet2CellGlob(icg),icg=1,pgrid%nTetsTot)
1104  END IF ! pGrid%nTetsTot
1105 
1106 ! ******************************************************************************
1107 ! Hexahedra
1108 ! ******************************************************************************
1109 
1110  IF ( pgrid%nHexsTot > 0 ) THEN
1111  IF ( global%myProcid == masterproc .AND. &
1112  global%verbLevel >= verbose_high) THEN
1113  WRITE(stdout,'(A,3X,A)') solver_name,'Hexahedra...'
1114  END IF ! global%myProcid
1115 
1116  sectionstring = '# Hexahedra'
1117  WRITE(ifile,'(A)') trim(sectionstring)
1118  WRITE(ifile,'(10(I8))') (pgrid%hex2CellGlob(icg),icg=1,pgrid%nHexsTot)
1119  END IF ! pGrid%nHexsTot
1120 
1121 ! ******************************************************************************
1122 ! Prisms
1123 ! ******************************************************************************
1124 
1125  IF ( pgrid%nPrisTot > 0 ) THEN
1126  IF ( global%myProcid == masterproc .AND. &
1127  global%verbLevel >= verbose_high) THEN
1128  WRITE(stdout,'(A,3X,A)') solver_name,'Prisms...'
1129  END IF ! global%myProcid
1130 
1131  sectionstring = '# Prisms'
1132  WRITE(ifile,'(A)') trim(sectionstring)
1133  WRITE(ifile,'(10(I8))') (pgrid%pri2CellGlob(icg),icg=1,pgrid%nPrisTot)
1134  END IF ! pGrid%nPrisTot
1135 
1136 ! ******************************************************************************
1137 ! Pyramids
1138 ! ******************************************************************************
1139 
1140  IF ( pgrid%nPyrsTot > 0 ) THEN
1141  IF ( global%myProcid == masterproc .AND. &
1142  global%verbLevel >= verbose_high) THEN
1143  WRITE(stdout,'(A,3X,A)') solver_name,'Pyramids...'
1144  END IF ! global%myProcid
1145 
1146  sectionstring = '# Pyramids'
1147  WRITE(ifile,'(A)') trim(sectionstring)
1148  WRITE(ifile,'(10(I8))') (pgrid%pyr2CellGlob(icg),icg=1,pgrid%nPyrsTot)
1149  END IF ! pGrid%nPyrsTot
1150 
1151 ! ******************************************************************************
1152 ! End marker
1153 ! ******************************************************************************
1154 
1155  IF ( global%myProcid == masterproc .AND. &
1156  global%verbLevel >= verbose_high) THEN
1157  WRITE(stdout,'(A,3X,A)') solver_name,'End marker...'
1158  END IF ! global%myProcid
1159 
1160  sectionstring = '# End'
1161  WRITE(ifile,'(A)') trim(sectionstring)
1162 
1163 ! ******************************************************************************
1164 ! Close file
1165 ! ******************************************************************************
1166 
1167  CLOSE(ifile,iostat=errorflag)
1168  global%error = errorflag
1169  IF ( global%error /= err_none ) THEN
1170  CALL errorstop(global,err_file_close,__line__,ifilename)
1171  END IF ! global%error
1172 
1173 ! ******************************************************************************
1174 ! End
1175 ! ******************************************************************************
1176 
1177  IF ( global%myProcid == masterproc .AND. &
1178  global%verbLevel >= verbose_high ) THEN
1179  WRITE(stdout,'(A,1X,A)') solver_name, &
1180  'Writing local-to-global cell mapping done.'
1181  END IF ! global%myProcid
1182 
1183  CALL deregisterfunction(global)
1184 
1185  END SUBROUTINE rflu_writeloc2globcellmapping
1186 
1187 
1188 
1189 
1190 
1191 
1192 
1193 
1194 ! ******************************************************************************
1195 ! End
1196 ! ******************************************************************************
1197 
1198 END MODULE rflu_modcellmapping
1199 
1200 
1201 ! ******************************************************************************
1202 !
1203 ! RCS Revision history:
1204 !
1205 ! $Log: RFLU_ModCellMapping.F90,v $
1206 ! Revision 1.16 2008/12/06 08:44:20 mtcampbe
1207 ! Updated license.
1208 !
1209 ! Revision 1.15 2008/11/19 22:17:31 mtcampbe
1210 ! Added Illinois Open Source License/Copyright
1211 !
1212 ! Revision 1.14 2005/04/15 15:06:47 haselbac
1213 ! Small bug fix (use Max instead of Tot), cosmetics
1214 !
1215 ! Revision 1.13 2004/12/04 03:26:53 haselbac
1216 ! Substantial rewrite to maximize code reuse for partitioning
1217 !
1218 ! Revision 1.12 2004/11/03 17:01:36 haselbac
1219 ! Rewrite because of removal of vertex anc cell flags, cosmetics
1220 !
1221 ! Revision 1.11 2004/01/22 16:03:58 haselbac
1222 ! Made contents of modules PRIVATE, only procs PUBLIC, to avoid errors on ALC
1223 ! and titan
1224 !
1225 ! Revision 1.10 2003/11/25 21:03:21 haselbac
1226 ! Improved error reporting for invalid cell flags
1227 !
1228 ! Revision 1.9 2003/06/04 22:08:30 haselbac
1229 ! Added Nullify routines, some cosmetics
1230 !
1231 ! Revision 1.8 2003/03/15 18:03:21 haselbac
1232 ! Adaptation for || RFLU, removed Charm calls to separate routine
1233 !
1234 ! Revision 1.7 2003/01/28 16:26:03 haselbac
1235 ! Cosmetics only
1236 !
1237 ! Revision 1.6 2002/10/19 16:31:15 haselbac
1238 ! Bug fix in output for parallel runs - missing if statement
1239 !
1240 ! Revision 1.5 2002/10/08 15:49:21 haselbac
1241 ! {IO}STAT=global%error replaced by {IO}STAT=errorFlag - SGI problem
1242 !
1243 ! Revision 1.4 2002/09/10 20:30:27 haselbac
1244 ! Corrected bug in RFLU_BuildCellMapping
1245 !
1246 ! Revision 1.3 2002/09/09 15:02:31 haselbac
1247 ! global now under regions
1248 !
1249 ! Revision 1.2 2002/07/25 14:59:31 haselbac
1250 ! Only write out for MASTERPROC, commented out FEM calls, added RFLU_ModFEM
1251 !
1252 ! Revision 1.1 2002/06/27 15:48:16 haselbac
1253 ! Initial revision
1254 !
1255 ! ******************************************************************************
1256 
1257 
1258 
1259 
1260 
1261 
1262 
1263 
1264 
1265 
1266 
1267 
1268 
subroutine buildfilenamebasic(global, dest, ext, id, fileName)
subroutine registerfunction(global, funName, fileName)
Definition: ModError.F90:449
int status() const
Obtain the status of the attribute.
Definition: Attribute.h:240
subroutine, public rflu_writeloc2globcellmapping(pRegion)
subroutine, public rflu_nullifycellmapping(pRegion)
subroutine, public rflu_readloc2globcellmapping(pRegion)
subroutine, public rflu_buildglob2loccellmapping(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 form
subroutine, public rflu_destroycellmapping(pRegion)
subroutine errorstop(global, errorCode, errorLine, addMessage)
Definition: ModError.F90:483
subroutine, public rflu_createcellmapping(pRegion)
subroutine deregisterfunction(global)
Definition: ModError.F90:469
subroutine, public rflu_buildloc2globcellmapping(pRegion)