avbp2dom.F [SRC] [CPP] [JOB] [SCAN]
TOOLS / PREDATAS / INOUT



   1 | include(dom.inc)
   2 | 
   3 |       SUBROUTINE avbp2dom(file1, file2, file3 )
   4 | 
   5 | !       ===============================================================!
   6 | !                                                                      !
   7 | !       avbp2dom.F : Reads AVBP mesh files and starts data structure   !
   8 | !                                                                      !
   9 | !       in         : The .coor file name 'file1'                       !
  10 | !                    The .conn file name 'file2'                       !
  11 | !                    The .exBound filename 'file3'                     !
  12 | !                                                                      !
  13 | !       author     : J. AMAYA (avril 2007)                             !
  14 | !                                                                      !
  15 | !       Nota       : i_nnodes and i_nfacesmax must be initialized in   !
  16 | !                    this subroutine                                   !
  17 | !                                                                      !
  18 | !       ===============================================================!
  19 |  
  20 |         USE avbp_coor
  21 |         USE avbp_conn
  22 |         USE avbp_exBound
  23 |         USE datas
  24 | 
  25 |         IMPLICIT NONE
  26 | 
  27 |         INCLUDE 'dom_constants.h'
  28 | 
  29 | !       IN
  30 |         CHARACTER*80 file1
  31 |         CHARACTER*80 file2
  32 |         CHARACTER*80 file3
  33 | 
  34 | !       LOCAL
  35 |         DOM_INT    :: i, j, iel, k, beg, iend, idx
  36 |         DOM_INT    :: nbeg, nlen, nend, pnode
  37 |         DOM_INT    :: nfacevert, nfaces, nvert
  38 |         DOM_INT    :: cellid
  39 |         DOM_INT    :: facelist6(6), facelist4(4), facelist3(3)
  40 |         DOM_INT    :: nodelist4(4), nodelist3(3), nodelist2(2)
  41 |         DOM_INT    :: nbpatches, facepatch
  42 | 
  43 |         character*80                        :: patchname, patchkey
  44 | 
  45 |         type(cell), pointer                 :: new_cell
  46 | 
  47 | !-- CL are treated in initbc
  48 | !       DOM_REAL, allocatable, dimension(:) :: bndy_epsil
  49 | !       DOM_REAL, allocatable, dimension(:) :: bndy_Tw
  50 | !--
  51 | 
  52 | !       --------------------!
  53 | !       Read AVBP coor file !
  54 | !       --------------------!
  55 |         coorfile = file1
  56 |         CALL readcoor
  57 | 
  58 |         i_nnodes = coor_nnode
  59 |         conn_ndim = coor_ndim
  60 | 
  61 | !       ---------------------!
  62 | !       Initialise the nodes !
  63 | !       ---------------------!
  64 | 
  65 |         ALLOCATE(node_list(3,i_nnodes))
  66 |         ALLOCATE(facesatnode(i_nnodes))
  67 | 
  68 |         DO i=1, i_nnodes
  69 |           NULLIFY(facesatnode(i)%fatnode_ptr)
  70 |         ENDDO
  71 | 
  72 |         IF (coor_ndim.eq.3) THEN
  73 |           node_list(1,:) = coor(1,:)
  74 |           node_list(2,:) = coor(2,:)
  75 |           node_list(3,:) = coor(3,:)
  76 |         ELSE
  77 |           node_list(1,:) = coor(1,:)
  78 |           node_list(2,:) = coor(2,:)
  79 |           node_list(3,:) = 0.0d0
  80 |         ENDIF
  81 | 
  82 |         CALL destroycoormemory
  83 | 
  84 | !       ---------------------!
  85 | !       Read connection file !
  86 | !       ---------------------!
  87 | 
  88 |         connfile = file2
  89 |         CALL readconn
  90 | 
  91 | !       -------------------------------------------!
  92 | !       Read .exBound file and boundary properties !
  93 | !       -------------------------------------------!
  94 | 
  95 |         exBound_ndim = coor_ndim
  96 |         exBoundfile = file3
  97 |         CALL readexBound
  98 |         CALL getExBoundinfo 
  99 | 
 100 |         ALLOCATE(nodepatch(0:exBound_npbound, coor_nnode))
 101 |         nodepatch = 0
 102 | 
 103 |         do j=1, exBound_npbound 
 104 |           WRITE(*,*) "Reading nodes at patch: ", j
 105 |           nbeg = exBound_ibound1(j)
 106 |           nlen = exBound_ibound2(j)
 107 |           nend = nbeg + nlen - 1
 108 | 
 109 |           do i = nbeg, nend
 110 |             pnode              = exBound_ibound(i)
 111 |             nodepatch(0,pnode) = nodepatch(0,pnode) + 1
 112 |             idx                = nodepatch(0,pnode)
 113 |             nodepatch(idx,pnode) = j
 114 |           enddo
 115 | 
 116 |         enddo
 117 | 
 118 | !-- CL are donne in initbc
 119 | !       -------------------------!
 120 | !       Read boundary properties !
 121 | !       -------------------------!
 122 | !
 123 | !       IF (ALLOCATED(bndy_epsil)) DEALLOCATE(bndy_epsil)
 124 | !       IF (ALLOCATED(bndy_Tw))    DEALLOCATE(bndy_Tw)
 125 | !
 126 | !       ALLOCATE(bndy_epsil(0:exBound_mnbpatch))
 127 | !       ALLOCATE(   bndy_Tw(0:exBound_mnbpatch))
 128 | !
 129 | !       ------------------------------------------------!
 130 | !       Setup default non-boundary faces epsilon and Tw !
 131 | !       ------------------------------------------------!
 132 | !
 133 | !       bndy_epsil(0) = -1.
 134 | !       bndy_Tw   (0) = 0.
 135 | !
 136 | !       --------------------------------------!
 137 | !       Read BC properties from boundary file !
 138 | !       --------------------------------------!
 139 | !
 140 | !       OPEN(1, FILE=bndyfile, FORM='formatted')
 141 | !
 142 | !       READ(1,*)
 143 | !       READ(1,*) nbpatches
 144 | !
 145 | !       IF (nbpatches.ne.exBound_mnbpatch) THEN
 146 | !         CLOSE(1)
 147 | !         WRITE(*,*) " Fatal Error: The number of patches in input "
 148 | !         WRITE(*,*) "              file input_boundary.dat does not"
 149 | !         WRITE(*,*) "              match the number of patches of "
 150 | !         WRITE(*,*) "              the exBound mesh file."
 151 | !         STOP
 152 | !       ENDIF
 153 | !
 154 | !       DO i=1, nbpatches
 155 | !         READ(1,*)
 156 | !         READ(1,*)
 157 | !         READ(1,*) patchname
 158 | !         READ(1,*) patchkey
 159 | !          print*,   trim(patchkey)
 160 | !
 161 | !         -----------------------------------!
 162 | !         Detect patch type and set variales !
 163 | !         -----------------------------------!
 164 | !
 165 | !         IF (trim(patchkey).eq.'EMISSIVE_WALL') THEN
 166 | !
 167 | !           READ(1,*) bndy_epsil(i)
 168 | !           READ(1,*) bndy_Tw(i)
 169 | !
 170 | !         ELSEIF(trim(patchkey).eq.'BLACK_WALL') THEN
 171 | !
 172 | !           bndy_epsil(i) = 1.
 173 | !           READ(1,*) bndy_Tw(i)
 174 | !
 175 | !         ELSE
 176 | !
 177 | !           WRITE(*,*) " Fatal Error: Unknown patch type"
 178 | !           WRITE(*,*) "              ", trim(patchkey)
 179 | !           STOP
 180 | !
 181 | !         ENDIF
 182 | !
 183 | !       ENDDO
 184 | !
 185 | !        print*, " => Boundary epsil:",bndy_epsil
 186 | !        print*, " => Boundary Tw:",bndy_Tw
 187 | !
 188 | !       CLOSE(1)
 189 | !--
 190 | 
 191 | !       --------------------------------------!
 192 | !       Allocating memory for face indexation !
 193 | !       --------------------------------------!
 194 | 
 195 |         ALLOCATE(faceidx(MAX_NFACES_CELL*conn_ntcell))
 196 | 
 197 | !       -----------------------!
 198 | !       Create faces and cells !
 199 | !       -----------------------!
 200 | 
 201 |         WRITE(*,*) " Building the data structure..."
 202 |         cellid = 1
 203 | 
 204 |         DO iel=1, conn_maxtype
 205 |           DO i=1, nclength(iel)
 206 | 
 207 | !           -----------------------------------------------------------!
 208 | !           Allocate memory space for the new cell (needed in addface) !
 209 | !           -----------------------------------------------------------!
 210 |             ALLOCATE(new_cell)
 211 |             NULLIFY(new_cell%next_cell)
 212 |             current_cell => new_cell
 213 |             current_cell%cell_id = cellid
 214 | 
 215 |             SELECT CASE (iel)
 216 | 
 217 |               CASE (EL_TRI)
 218 | 
 219 |                 nfacevert = 2
 220 |                 nvert     = 3
 221 |                 nfaces    = 3
 222 |                 beg       = ncbeg(iel) +((i-1)*nvert)
 223 |                 iend      = beg + nvert
 224 |                 
 225 |                 IF (nfaces.gt.i_nfacesmax) i_nfacesmax = nfaces
 226 | 
 227 | !               -----------------!
 228 | !               Create Tri faces !
 229 | !               -----------------!
 230 | 
 231 |                 nodelist2=(/ ielno(beg+2), ielno(beg+3) /)
 232 |                 CALL detectpatch(nfacevert, nodelist2, facepatch,       &
 233 |      &                           nodepatch,coor_nnode)
 234 |                 CALL addface(cellid, nfacevert, nodelist2, facelist3(1),&
 235 |      &                       facepatch)
 236 | 
 237 |                 nodelist2=(/ ielno(beg+3), ielno(beg+1) /)
 238 |                 CALL detectpatch(nfacevert, nodelist2, facepatch,       &
 239 |      &                           nodepatch,coor_nnode)
 240 |                 CALL addface(cellid, nfacevert, nodelist2, facelist3(2),&
 241 |      &                       facepatch)
 242 | 
 243 |                 nodelist2=(/ ielno(beg+1), ielno(beg+2) /)
 244 |                 CALL detectpatch(nfacevert, nodelist2, facepatch,       &
 245 |      &                           nodepatch,coor_nnode)
 246 |                 CALL addface(cellid, nfacevert, nodelist2, facelist3(3),&
 247 |      &                       facepatch)
 248 | 
 249 | !               ------------------------!
 250 | !               Create the new TRI cell !
 251 | !               ------------------------!
 252 | 
 253 |                 CALL addcell(cellid, EL_TRI, 3, facelist3)
 254 | !               print*, "  >> Cell ", cellid, " added"
 255 |                 cellid = cellid + 1
 256 | 
 257 |               CASE (EL_QUAD)
 258 | 
 259 |                 nfacevert = 2
 260 |                 nvert     = 4
 261 |                 nfaces    = 4
 262 |                 beg       = ncbeg(iel) +((i-1)*nvert)
 263 |                 iend      = beg + nvert
 264 | 
 265 |                 IF (nfaces.gt.i_nfacesmax) i_nfacesmax = nfaces
 266 | 
 267 | !               ------------------!
 268 | !               Create Quad faces !
 269 | !               ------------------!
 270 | 
 271 |                 nodelist2=(/ ielno(beg+1), ielno(beg+2) /)
 272 |                 CALL detectpatch(nfacevert, nodelist2, facepatch,       &
 273 |      &                           nodepatch,coor_nnode)
 274 |                 CALL addface(cellid, nfacevert, nodelist2, facelist4(1),&
 275 |      &                       facepatch)
 276 | 
 277 |                 nodelist2=(/ ielno(beg+2), ielno(beg+3) /)
 278 |                 CALL detectpatch(nfacevert, nodelist2, facepatch,       &
 279 |      &                           nodepatch,coor_nnode)
 280 |                 CALL addface(cellid, nfacevert, nodelist2, facelist4(2),&
 281 |      &                       facepatch)
 282 | 
 283 |                 nodelist2=(/ ielno(beg+3), ielno(beg+4) /)
 284 |                 CALL detectpatch(nfacevert, nodelist2, facepatch,       &
 285 |      &                           nodepatch,coor_nnode)
 286 |                 CALL addface(cellid, nfacevert, nodelist2, facelist4(3),&
 287 |      &                       facepatch)
 288 | 
 289 |                 nodelist2=(/ ielno(beg+4), ielno(beg+1) /)
 290 |                 CALL detectpatch(nfacevert, nodelist2, facepatch,       &
 291 |      &                           nodepatch,coor_nnode)
 292 |                 CALL addface(cellid, nfacevert, nodelist2, facelist4(4),&
 293 |      &                       facepatch)
 294 | 
 295 | !               -------------------------!
 296 | !               Create the new QUAD cell !
 297 | !               -------------------------!
 298 | 
 299 |                 CALL addcell(cellid, EL_QUAD, 4, facelist4)
 300 | !               print*, "  >> Cell ", cellid, " added"
 301 |                 cellid = cellid + 1
 302 | 
 303 |               CASE (EL_TETRA)
 304 | 
 305 |                 nfacevert = 3
 306 |                 nvert     = 4
 307 |                 nfaces    = 4
 308 |                 beg       = ncbeg(iel) +((i-1)*nvert)
 309 |                 iend      = beg + nvert
 310 | 
 311 |                 IF (nfaces.gt.i_nfacesmax) i_nfacesmax = nfaces
 312 | 
 313 | !               -------------------!
 314 | !               Create Tetra faces !
 315 | !               -------------------!
 316 | 
 317 |                 nodelist3=(/ ielno(beg+2), ielno(beg+4), ielno(beg+3) /)
 318 |                 CALL detectpatch(nfacevert, nodelist3, facepatch,       &
 319 |      &                           nodepatch,coor_nnode)
 320 |                 CALL addface(cellid, nfacevert, nodelist3, facelist4(1),&
 321 |      &                       facepatch)
 322 | 
 323 |                 nodelist3=(/ ielno(beg+1), ielno(beg+3), ielno(beg+4) /)
 324 |                 CALL detectpatch(nfacevert, nodelist3, facepatch,       &
 325 |      &                           nodepatch,coor_nnode)
 326 |                 CALL addface(cellid, nfacevert, nodelist3, facelist4(2),&
 327 |      &                       facepatch)
 328 | 
 329 |                 nodelist3=(/ ielno(beg+1), ielno(beg+4), ielno(beg+2) /)
 330 |                 CALL detectpatch(nfacevert, nodelist3, facepatch,       &
 331 |      &                           nodepatch,coor_nnode)
 332 |                 CALL addface(cellid, nfacevert, nodelist3, facelist4(3),&
 333 |      &                       facepatch)
 334 | 
 335 |                 nodelist3=(/ ielno(beg+1), ielno(beg+2), ielno(beg+3) /)
 336 |                 CALL detectpatch(nfacevert, nodelist3, facepatch,       &
 337 |      &                           nodepatch,coor_nnode)
 338 |                 CALL addface(cellid, nfacevert, nodelist3, facelist4(4),&
 339 |      &                       facepatch)
 340 | 
 341 | !               --------------------------!
 342 | !               Create the new TETRA cell !
 343 | !               --------------------------!
 344 | 
 345 |                 CALL addcell(cellid, EL_TETRA, 4, facelist4)
 346 | !               print*, "  >> Cell ", cellid, " added"
 347 |                 cellid = cellid + 1
 348 | 
 349 |               CASE (EL_HEXA)
 350 | 
 351 |                 nfacevert = 4
 352 |                 nvert     = 8
 353 |                 nfaces    = 6
 354 |                 beg       = ncbeg(iel) +((i-1)*nvert) 
 355 |                 iend      = beg + nvert
 356 | 
 357 |                 IF (nfaces.gt.i_nfacesmax) i_nfacesmax = nfaces
 358 | 
 359 | !               ------------------!
 360 | !               Create Hexa faces !
 361 | !               ------------------!
 362 |                 nodelist4 = (/ ielno(beg+1), ielno(beg+2),              &
 363 |      &                        ielno(beg+6), ielno(beg+5) /)
 364 |                 CALL detectpatch(nfacevert, nodelist4, facepatch,       &
 365 |      &                           nodepatch,coor_nnode)
 366 |                 CALL addface(cellid, nfacevert, nodelist4, facelist6(1),&
 367 |      &                       facepatch)
 368 | 
 369 |                 nodelist4 = (/ ielno(beg+6), ielno(beg+2),              &
 370 |      &                        ielno(beg+3), ielno(beg+7) /)
 371 |                 CALL detectpatch(nfacevert, nodelist4, facepatch,       &
 372 |      &                           nodepatch,coor_nnode)
 373 |                 CALL addface(cellid, nfacevert, nodelist4, facelist6(2),&
 374 |      &                       facepatch)
 375 | 
 376 |                 nodelist4 = (/ ielno(beg+8), ielno(beg+7),              &
 377 |      &                        ielno(beg+3), ielno(beg+4) /)
 378 |                 CALL detectpatch(nfacevert, nodelist4, facepatch,       &
 379 |      &                           nodepatch,coor_nnode)
 380 |                 CALL addface(cellid, nfacevert, nodelist4, facelist6(3),&
 381 |      &                       facepatch)
 382 | 
 383 |                 nodelist4 = (/ ielno(beg+1), ielno(beg+5),              &
 384 |      &                        ielno(beg+8), ielno(beg+4) /)
 385 |                 CALL detectpatch(nfacevert, nodelist4, facepatch,       &
 386 |      &                           nodepatch,coor_nnode)
 387 |                 CALL addface(cellid, nfacevert, nodelist4, facelist6(4),&
 388 |      &                       facepatch)
 389 | 
 390 |                 nodelist4 = (/ ielno(beg+1), ielno(beg+4),              &
 391 |      &                        ielno(beg+3), ielno(beg+2) /)
 392 |                 CALL detectpatch(nfacevert, nodelist4, facepatch,       &
 393 |      &                           nodepatch,coor_nnode)
 394 |                 CALL addface(cellid, nfacevert, nodelist4, facelist6(5),&
 395 |      &                       facepatch)
 396 | 
 397 |                 nodelist4 = (/ ielno(beg+5), ielno(beg+6),              &
 398 |      &                        ielno(beg+7), ielno(beg+8) /)
 399 |                 CALL detectpatch(nfacevert, nodelist4, facepatch,       &
 400 |      &                           nodepatch,coor_nnode)
 401 |                 CALL addface(cellid, nfacevert, nodelist4, facelist6(6),&
 402 |      &                       facepatch)
 403 | 
 404 | !               -------------------------!
 405 | !               Create the new HEXA cell !
 406 | !               -------------------------!
 407 | 
 408 |                 CALL addcell(cellid, EL_HEXA, 6, facelist6)
 409 | !               print*, "  >> Cell ", cellid, " added"
 410 |                 cellid = cellid + 1
 411 | 
 412 |             END SELECT
 413 | 
 414 |           ENDDO
 415 |         ENDDO
 416 | 
 417 |         WRITE(*,*) " Total cells added: ", i_ncells
 418 |         WRITE(*,*) " Total faces added: ", i_nfaces
 419 |         CALL destroyconnmemory
 420 | 
 421 |       END SUBROUTINE avbp2dom


avbp2dom.F could be called by:
avbptest.F [TOOLS/PREDATAS/EXAMPLES] - 37
Makefile [TOOLS/PREDATAS] - 101
predatas.F [TOOLS/PREDATAS/SRC] - 84