ensightgold_geometry.F [SRC] [CPP] [JOB] [SCAN]
TOOLS / EXTERNAL



   1 | include(dom.inc)
   2 | !     ===============================================================
   3 | !     Copyright (c) CERFACS (all rights reserved)
   4 | !     ===============================================================
   5 | 
   6 | !     ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
   7 | !     
   8 | !     module: ensightgold_geometry.f
   9 | ! 
  10 | !     author: Andre Kaufmann, 10.08.2000
  11 | !     author: Charles Martin, 09.08.2004
  12 | !     
  13 | !     last change: J. Amaya, O. Cabrit, A. Roux - 11/2006
  14 | !     ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
  15 |       module ensightgold_geometry
  16 |       use avbp_coor
  17 |       use avbp_conn
  18 |       implicit none
  19 |  
  20 | !     variable declaration
  21 |  
  22 |       DOM_INT geometry_maxtype
  23 |       real*4, dimension(:,:), allocatable :: geometry_coor
  24 |       DOM_INT, dimension(:),allocatable   :: int_buf_array
  25 |       character*80, dimension(1:8)      :: geometry_eltype
  26 |       logical                           :: geometry_debug  = .false.
  27 |  
  28 | !     file format           1 binary
  29 | !                           2 ascii without line numbers
  30 | !                           3 ascii with line numbers
  31 |  
  32 | 
  33 |       DOM_INT  :: geometry_write_format = 1
  34 | 
  35 |       contains
  36 | !     ================================================================
  37 | !     geometry values
  38 | !     ================================================================
  39 |       subroutine intialisegeometry
  40 |          geometry_maxtype = 8
  41 |          geometry_eltype(1) = 'tria3'
  42 |          geometry_eltype(2) = 'quad4'
  43 |          geometry_eltype(3) = 'tetra4'
  44 |          geometry_eltype(4) = 'pyramid5'
  45 |          geometry_eltype(5) = 'penta6'
  46 |          geometry_eltype(6) = 'hexa8'
  47 |          geometry_eltype(7) = '     '
  48 |          geometry_eltype(8) = '     '
  49 |       end subroutine intialisegeometry
  50 | !     ================================================================
  51 | !     create binary geometry
  52 | !     ================================================================
  53 |       subroutine creategeometry(ensight_coorfile,                       &
  54 |      &                          ensight_connfile,                       &
  55 |      &                          ensight_asciiBoundfile,                 &
  56 |      &                          ensight_geometryfile)
  57 |       implicit none
  58 | !
  59 | !     variable declaration for subroutine
  60 | !
  61 |       character*80 ensight_coorfile,                                    &
  62 |      &             ensight_connfile,                                    &
  63 |      &             ensight_asciiBoundfile,                              &
  64 |      &             ensight_geometryfile
  65 |       character*80 buffer
  66 |       DOM_INT n,iel,nlen,nvert,nbeg,k
  67 |       DOM_INT nelem,buffer_int,lenofchar
  68 |       real z
  69 | !
  70 | !     debug information
  71 | !
  72 |       if(geometry_debug) then
  73 |          write(*,*) ' creating the binary coordinate file ensight.geo'
  74 |       endif
  75 | 
  76 | !
  77 | !     read the coor and conn file
  78 | !
  79 |       coorfile = ensight_coorfile
  80 |       connfile = ensight_connfile
  81 | !
  82 | !     read coorfile
  83 | !
  84 |       if(geometry_debug) write(*,*) ' reading *.coor file : ',coorfile
  85 |       call readcoor
  86 | 
  87 |       conn_ndim = coor_ndim
  88 | !
  89 | !     read connfile
  90 | !
  91 |       if(geometry_debug) write(*,*) ' reading *.conn file : ',connfile
  92 |       call readconn
  93 | 
  94 | !
  95 | !     allocate dynamik memory for simple precision 
  96 | !     coor data
  97 | !
  98 |       if(geometry_debug) write(*,*) ' copying grid from dbl to smpl pr'
  99 | 
 100 |       if(allocated(geometry_coor)) then
 101 |          deallocate(geometry_coor)
 102 |       endif
 103 |       allocate(geometry_coor(1:coor_nnode,1:3))
 104 | !
 105 | !     copy from double to simple precision
 106 | ! 
 107 |       do n=1,coor_nnode
 108 |         do k=1,coor_ndim
 109 |           geometry_coor(n,k) = coor(k,n)
 110 |         enddo
 111 |       enddo
 112 | !     2D treatment
 113 |       if (coor_ndim.eq.2) then 
 114 |         geometry_coor(:,3)=0
 115 |       endif
 116 | 
 117 | !
 118 | !     read the Boundary files
 119 | !
 120 | !
 121 | !     write ensight geometry (binary version)
 122 | !      
 123 |       n = len_trim(ensight_geometryfile)
 124 |       write(*,*) ' >>>> Writing geom to ', ensight_geometryfile(1:n)
 125 | 
 126 |         call ensightgold_write_geocoor_bin(ensight_geometryfile,n,      &
 127 |      &     coor_nnode,geometry_coor(:,1),                               &
 128 |      &                geometry_coor(:,2),                               &
 129 |      &                geometry_coor(:,3))
 130 |         nbeg    = 0
 131 |         
 132 |         if(geometry_debug) write(*,*) 'maxtype =',geometry_maxtype
 133 | 
 134 | !
 135 | !     loop over element types
 136 | !
 137 |         do iel=1,geometry_maxtype
 138 | 
 139 |           nlen = nclength(iel)
 140 | 
 141 |           if(geometry_debug) write(*,*) 'nlength = ',nlen
 142 | 
 143 |           if ( nlen.ne.0 ) then
 144 | 
 145 |             nvert = ieltype(2,iel)
 146 | 
 147 |             if(geometry_debug)  then
 148 |              write(*,*) 'geometry_eltype = ',geometry_eltype(iel)
 149 |              write(*,*) 'nlen            = ',nlen
 150 |             endif
 151 | 
 152 |             if(geometry_debug) then
 153 |                write(*,*) '     eltype',iel,nvert,nlen,                 &
 154 |      &              geometry_eltype(iel)
 155 |             endif
 156 | 
 157 |             if(allocated(int_buf_array)) then
 158 |                deallocate(int_buf_array)
 159 |             endif
 160 |             nelem=nlen*nvert
 161 |             allocate(int_buf_array(1:nelem))
 162 | 
 163 | !     are these statements about nvert = 6 or not really usefull?
 164 | !     ask andre kaufmann if not , int_buf_array can be supressed
 165 | 
 166 |             if ( nvert.ne.6 ) then  
 167 |                int_buf_array(1:nelem)=ielno(nbeg+1:nbeg+nelem)  
 168 |             else if ( nvert.eq.6 ) then
 169 |                do n=1,nlen
 170 |                   int_buf_array((n-1)*nvert+1)=ielno(nbeg+(n-1)*nvert+1)
 171 |                   int_buf_array((n-1)*nvert+2)=ielno(nbeg+(n-1)*nvert+4)
 172 |                   int_buf_array((n-1)*nvert+3)=ielno(nbeg+(n-1)*nvert+6)
 173 |                   int_buf_array((n-1)*nvert+4)=ielno(nbeg+(n-1)*nvert+2)
 174 |                   int_buf_array((n-1)*nvert+5)=ielno(nbeg+(n-1)*nvert+3)
 175 |                   int_buf_array((n-1)*nvert+6)=ielno(nbeg+(n-1)*nvert+5)
 176 |                enddo
 177 |             end if
 178 | 
 179 |             lenofchar = len_trim(geometry_eltype(iel))
 180 |             call ensightgold_write_geoelt_bin(ensight_geometryfile,n,   &
 181 |      &           geometry_eltype(iel),lenofchar,                        &
 182 |      &           nvert,nlen,int_buf_array)
 183 | 
 184 | 
 185 |             nbeg = nbeg+nlen*nvert
 186 |             
 187 |          end if
 188 | 
 189 |         end do
 190 | 
 191 |       close(1)
 192 | 
 193 |       if(geometry_debug) then
 194 |          write(*,*) ' >>>> Ensight geometry file written'
 195 |          write(*,*)
 196 |       endif
 197 | 
 198 |       if(allocated(geometry_coor)) then
 199 |          deallocate(geometry_coor)
 200 |       endif
 201 |       
 202 |       end subroutine creategeometry
 203 | !     ================================================================
 204 | !     create binary geometry_ascii
 205 | !     ================================================================
 206 |       subroutine creategeometry_ascii(ensight_coorfile,                 &
 207 |      &                                ensight_connfile,                 &
 208 |      &                                ensight_asciiBoundfile,           &
 209 |      &                                ensight_geometryfile)
 210 |       implicit none
 211 | !
 212 | !     variable declaration for subroutine
 213 | !
 214 |       character*80 ensight_coorfile,                                    &
 215 |      &             ensight_connfile,                                    &
 216 |      &             ensight_asciiBoundfile,                              &
 217 |      &             ensight_geometryfile
 218 |       character*80 buffer
 219 |       DOM_INT n,iel,nlen,nvert,nbeg,k
 220 |       real z
 221 | 
 222 |       if(geometry_debug) then
 223 |          write(*,*) ' --------- creategeometry_ascii ---------'
 224 |       endif
 225 | !
 226 | !     debug information
 227 | !
 228 |       if(geometry_debug) then
 229 |          write(*,*) ' creating the ascii coordinate file ensight.geo'
 230 |       endif
 231 | !
 232 | !     read the coor and conn file
 233 | !
 234 |       coorfile = ensight_coorfile
 235 |       connfile = ensight_connfile
 236 | !
 237 | !     read coorfile
 238 | !
 239 |       if(geometry_debug) write(*,*) ' reading *.coor file : ',coorfile
 240 |       call readcoor
 241 | 
 242 |       conn_ndim = coor_ndim
 243 | 
 244 | !
 245 | !     read connfile
 246 | !
 247 |       if(geometry_debug) write(*,*) ' reading *.conn file : ',connfile
 248 |       call readconn
 249 | !
 250 | !     allocate dynamik memory for simple precision 
 251 | !     coor data
 252 | !
 253 |       if(geometry_debug) write(*,*) ' copying grid dble to smple prec'
 254 | 
 255 |       if(allocated(geometry_coor)) then
 256 |          deallocate(geometry_coor)
 257 |       endif
 258 | 
 259 |       allocate(geometry_coor(1:coor_ndim,1:coor_nnode))
 260 | !
 261 | !     copy from double to simple precision
 262 | !
 263 |       do n=1,coor_nnode
 264 |         do k=1,coor_ndim
 265 |           geometry_coor(k,n) = coor(k,n)
 266 |         enddo
 267 |       enddo
 268 | !
 269 | !     read the Boundary files
 270 | !
 271 | 
 272 | !
 273 | !     write ensight geometry (ascii version)
 274 | !
 275 | !
 276 | !     see Ensight User Manual 5.5.2 p.3-49
 277 | !
 278 | !
 279 | !      
 280 |       n = len_trim(ensight_geometryfile)
 281 |       write(*,*) ' >>>> Writing ensight geometry to ',                  &
 282 |      &           ensight_geometryfile(1:n)
 283 |       open( unit=1, file=ensight_geometryfile, form='formatted' )
 284 | 
 285 |         buffer = 'Ensight Gold ASCII format'
 286 |         write(1,'(A)') buffer
 287 |         buffer = 'Generated by a2gold translator'
 288 |         write(1,'(A)') buffer
 289 | !
 290 | !     write header depending on file format
 291 | !        
 292 |         if(geometry_write_format.eq.1) STOP ' this is ascii write ...' 
 293 | 
 294 | !
 295 | !     file format 2
 296 | !
 297 |         if(geometry_write_format.eq.2) then
 298 | 
 299 |            buffer = 'node id assign'
 300 |            write(1,'(A)') buffer
 301 |            buffer = 'element id assign'
 302 |            write(1,'(A)') buffer
 303 |            buffer = 'part'
 304 |            write(1,'(A)') buffer
 305 |            
 306 |            write(1,'(i10)') 1
 307 |            buffer = 'First part'
 308 |            write(1,'(A)') buffer           
 309 |            buffer = 'coordinates'
 310 |            write(1,'(A)') buffer
 311 |            write(1,'(i10)') coor_nnode
 312 | !
 313 | !     file format 3
 314 | !
 315 |         else if(geometry_write_format.eq.3) then
 316 | 
 317 |            buffer = 'node id given'
 318 |            write(1,'(A)') buffer
 319 |            buffer = 'element id given'
 320 |            write(1,'(A)') buffer
 321 |            
 322 |            buffer = 'coordinates'
 323 |            write(1,'(A)') buffer
 324 |            write(1,'(i10)') coor_nnode
 325 | 
 326 |         endif
 327 | 
 328 |         if(geometry_write_format.eq.2) then
 329 | !
 330 | !     file format 2
 331 | !
 332 |            if (coor_ndim.eq.2 ) then
 333 |               z=0.0d0
 334 | !              do n=1,coor_nnode
 335 | !                 write(1,'(e12.5,e12.5,e12.5)') geometry_coor(1,n),
 336 | !     &                geometry_coor(2,n),
 337 | !     &                z 
 338 | !              enddo
 339 | 
 340 |                do n=1,coor_nnode
 341 |                   write(1,'(e12.5,e12.5,e12.5)') geometry_coor(1,n)
 342 |                enddo
 343 |                do n=1,coor_nnode
 344 |                   write(1,'(e12.5,e12.5,e12.5)') geometry_coor(2,n)
 345 |                enddo
 346 |                do n=1,coor_nnode
 347 |                   write(1,'(e12.5,e12.5,e12.5)') z
 348 |                enddo
 349 |            else
 350 |               
 351 |               do k=1,3
 352 |                  do n=1,coor_nnode
 353 |                     write(1,'(e12.5)') geometry_coor(k,n)    
 354 |                  enddo
 355 |               enddo
 356 |            endif
 357 | 
 358 | !
 359 | !     file format 3
 360 | !
 361 |         else if(geometry_write_format.eq.3) then
 362 | 
 363 |         if (coor_ndim.eq.2 ) then
 364 |          z=0.0d0
 365 |          do n=1,coor_nnode
 366 |            write(1,'(i10,e12.5,e12.5,e12.5)') n,                        &
 367 |      &                  geometry_coor(1,n),                             &
 368 |      &                  geometry_coor(2,n),                             &
 369 |      &                  z 
 370 |          enddo
 371 |         else
 372 |          do n=1,coor_nnode
 373 |            write(1,'(i10,e12.5,e12.5,e12.5)') n,                        &
 374 |      &                  geometry_coor(1,n),                             &
 375 |      &                  geometry_coor(2,n),                             &
 376 |      &                  geometry_coor(3,n)
 377 |          enddo
 378 |         endif
 379 |         
 380 |         endif
 381 | 
 382 | !
 383 | !     element part
 384 | !
 385 | 
 386 | !
 387 | !     file format 2
 388 | !
 389 |         if(geometry_write_format.eq.2) then
 390 | 
 391 |         nbeg    = 0
 392 |         
 393 |         if(geometry_debug) write(*,*) 'maxtype =',geometry_maxtype
 394 | 
 395 | !
 396 | !     loop over element types
 397 | !
 398 |         do iel=1,geometry_maxtype
 399 | 
 400 |           nlen = nclength(iel)
 401 | 
 402 |           if(geometry_debug) write(*,*) 'nlength = ',nlen
 403 | 
 404 |           if ( nlen.ne.0 ) then
 405 | 
 406 |             nvert = ieltype(2,iel)
 407 | 
 408 |             if(geometry_debug)  then
 409 |              write(*,*) 'geometry_eltype = ',geometry_eltype(iel)
 410 |              write(*,*) 'nlen            = ',nlen
 411 |             endif
 412 | !
 413 | !     geometry type information
 414 | !
 415 |             write(1,'(A)') geometry_eltype(iel)
 416 |             write(1,'(i10)') nlen
 417 | 
 418 |             if(geometry_debug) then
 419 |                write(*,*) '     eltype',iel,nvert,nlen,                 &
 420 |      &                    geometry_eltype(iel)
 421 |             endif
 422 | !
 423 | !     three vertices per element
 424 | !
 425 |             if ( nvert.eq.3 ) then
 426 |                
 427 |                do n=1,nlen
 428 |                 write(1,'(i10,i10,i10)') ielno(nbeg+(n-1)*nvert+1),     &
 429 |      &                       ielno(nbeg+(n-1)*nvert+2),                 &
 430 |      &                       ielno(nbeg+(n-1)*nvert+3)
 431 |                enddo
 432 | !
 433 | !     four vertices per element
 434 | !
 435 |             else if ( nvert.eq.4 ) then
 436 | 
 437 |                do n=1,nlen
 438 |                 write(1,'(i10,i10,i10,i10)') ielno(nbeg+(n-1)*nvert+1), &
 439 |      &                       ielno(nbeg+(n-1)*nvert+2),                 &
 440 |      &                       ielno(nbeg+(n-1)*nvert+3),                 &
 441 |      &                       ielno(nbeg+(n-1)*nvert+4)
 442 |                enddo
 443 | !
 444 | !     five vertices per element
 445 | !
 446 |             else if ( nvert.eq.5 ) then
 447 |                do n=1,nlen
 448 |                  write(1,'(i10,i10,i10,i10,i10)')                       &
 449 |      &                       ielno(nbeg+(n-1)*nvert+1),                 &
 450 |      &                       ielno(nbeg+(n-1)*nvert+2),                 &
 451 |      &                       ielno(nbeg+(n-1)*nvert+3),                 &
 452 |      &                       ielno(nbeg+(n-1)*nvert+4),                 &
 453 |      &                       ielno(nbeg+(n-1)*nvert+5)
 454 |                enddo
 455 | !
 456 | !     six vertices per element
 457 | !
 458 |             else if ( nvert.eq.6 ) then
 459 | 
 460 |                do n=1,nlen
 461 |                 write(1,'(i10,i10,i10,i10,i10,i10)')                    &
 462 |      &                       ielno(nbeg+(n-1)*nvert+1),                 &
 463 |      &                       ielno(nbeg+(n-1)*nvert+4),                 &
 464 |      &                       ielno(nbeg+(n-1)*nvert+6),                 &
 465 |      &                       ielno(nbeg+(n-1)*nvert+2),                 &
 466 |      &                       ielno(nbeg+(n-1)*nvert+3),                 &
 467 |      &                       ielno(nbeg+(n-1)*nvert+5)
 468 |                enddo
 469 | !
 470 | !     eight vertices per element
 471 | !
 472 |             else if ( nvert.eq.8 ) then
 473 | 
 474 |                do n=1,nlen
 475 |                 write(1,'(i10,i10,i10,i10,i10,i10,i10,i10)')            &
 476 |      &                       ielno(nbeg+(n-1)*nvert+1),                 &
 477 |      &                       ielno(nbeg+(n-1)*nvert+2),                 &
 478 |      &                       ielno(nbeg+(n-1)*nvert+3),                 &
 479 |      &                       ielno(nbeg+(n-1)*nvert+4),                 &
 480 |      &                       ielno(nbeg+(n-1)*nvert+5),                 &
 481 |      &                       ielno(nbeg+(n-1)*nvert+6),                 &
 482 |      &                       ielno(nbeg+(n-1)*nvert+7),                 &
 483 |      &                       ielno(nbeg+(n-1)*nvert+8)
 484 | 
 485 |                enddo
 486 |             end if
 487 | 
 488 |             nbeg = nbeg+nlen*nvert
 489 | 
 490 |           end if
 491 | 
 492 |         end do
 493 | 
 494 | !
 495 | !     file format 3
 496 | !
 497 |         else if(geometry_write_format.eq.3) then
 498 | 
 499 |         buffer = 'part 1'
 500 |         write(1,'(A)') buffer
 501 |         buffer = 'First part'
 502 |         write(1,'(A)') buffer
 503 | 
 504 |         nbeg    = 0
 505 |         
 506 |         if(geometry_debug) write(*,*) 'maxtype =',geometry_maxtype
 507 | !
 508 | !     loop over element types
 509 | !
 510 |         do iel=1,geometry_maxtype
 511 | 
 512 |           nlen = nclength(iel)
 513 | 
 514 |           if(geometry_debug) write(*,*) 'nlength = ',nlen
 515 | 
 516 |           if ( nlen.ne.0 ) then
 517 | 
 518 |             nvert = ieltype(2,iel)
 519 | 
 520 |             if(geometry_debug)  then
 521 |              write(*,*) 'geometry_eltype = ',geometry_eltype(iel)
 522 |              write(*,*) 'nlen            = ',nlen
 523 |             endif
 524 | !
 525 | !     geometry type information
 526 | !
 527 |             write(1,'(A)') geometry_eltype(iel)
 528 |             write(1,'(i10)') nlen
 529 | 
 530 |             if(geometry_debug) then
 531 |                write(*,*) '     eltype',iel,nvert,nlen,                 &
 532 |      &                    geometry_eltype(iel)
 533 |             endif
 534 | !
 535 | !     three vertices per element
 536 | !
 537 |             if ( nvert.eq.3 ) then
 538 |                
 539 |                do n=1,nlen
 540 |                 write(1,'(i10,i10,i10,i10)') n,                         &
 541 |      &                       ielno(nbeg+(n-1)*nvert+1),                 &
 542 |      &                       ielno(nbeg+(n-1)*nvert+2),                 &
 543 |      &                       ielno(nbeg+(n-1)*nvert+3)
 544 |                enddo
 545 | !
 546 | !     four vertices per element
 547 | !
 548 |             else if ( nvert.eq.4 ) then
 549 | 
 550 |                do n=1,nlen
 551 |                 write(1,'(i10,i10,i10,i10,i10)') n,                     &
 552 |      &                       ielno(nbeg+(n-1)*nvert+1),                 &
 553 |      &                       ielno(nbeg+(n-1)*nvert+2),                 &
 554 |      &                       ielno(nbeg+(n-1)*nvert+3),                 &
 555 |      &                       ielno(nbeg+(n-1)*nvert+4)
 556 |                enddo
 557 | !
 558 | !     five vertices per element
 559 | !
 560 |             else if ( nvert.eq.5 ) then
 561 |                do n=1,nlen
 562 |                 write(1,'(i10,i10,i10,i10,i10,i10)') n,                 &
 563 |      &                       ielno(nbeg+(n-1)*nvert+1),                 &
 564 |      &                       ielno(nbeg+(n-1)*nvert+2),                 &
 565 |      &                       ielno(nbeg+(n-1)*nvert+3),                 &
 566 |      &                       ielno(nbeg+(n-1)*nvert+4),                 &
 567 |      &                       ielno(nbeg+(n-1)*nvert+5)
 568 |                enddo
 569 | !
 570 | !     six vertices per element
 571 | !
 572 |             else if ( nvert.eq.6 ) then
 573 | 
 574 |                do n=1,nlen
 575 |                 write(1,'(i10,i10,i10,i10,i10,i10,i10)') n,             &
 576 |      &                       ielno(nbeg+(n-1)*nvert+1),                 &
 577 |      &                       ielno(nbeg+(n-1)*nvert+4),                 &
 578 |      &                       ielno(nbeg+(n-1)*nvert+6),                 &
 579 |      &                       ielno(nbeg+(n-1)*nvert+2),                 &
 580 |      &                       ielno(nbeg+(n-1)*nvert+3),                 &
 581 |      &                       ielno(nbeg+(n-1)*nvert+5)
 582 |                enddo
 583 | !
 584 | !     eight vertices per element
 585 | !
 586 |             else if ( nvert.eq.8 ) then
 587 | 
 588 |                do n=1,nlen
 589 |                 write(1,'(i10,i10,i10,i10,i10,i10,i10,i10,i10)') n,     &
 590 |      &                       ielno(nbeg+(n-1)*nvert+1),                 &
 591 |      &                       ielno(nbeg+(n-1)*nvert+2),                 &
 592 |      &                       ielno(nbeg+(n-1)*nvert+3),                 &
 593 |      &                       ielno(nbeg+(n-1)*nvert+4),                 &
 594 |      &                       ielno(nbeg+(n-1)*nvert+5),                 &
 595 |      &                       ielno(nbeg+(n-1)*nvert+6),                 &
 596 |      &                       ielno(nbeg+(n-1)*nvert+7),                 &
 597 |      &                       ielno(nbeg+(n-1)*nvert+8)
 598 | 
 599 |                enddo
 600 |             end if
 601 | 
 602 |             nbeg = nbeg+nlen*nvert
 603 | 
 604 |           end if
 605 | 
 606 |         end do
 607 | 
 608 |         endif
 609 | 
 610 |       close(1)
 611 | 
 612 |       if(geometry_debug) then
 613 |          write(*,*) ' >>>> Ensight geometry file written'
 614 |          write(*,*)
 615 |       endif
 616 | 
 617 |       if(allocated(geometry_coor)) then
 618 |          deallocate(geometry_coor)
 619 |       endif
 620 |       
 621 |       end subroutine creategeometry_ascii
 622 | !     ================================================================
 623 |       end module ensightgold_geometry


ensightgold_geometry.F could be called by:
creategeometry_gambit.F [TOOLS/VISU/SRC] - 13
Makefile [TOOLS/VISU] - 62
visual_ensight.F [TOOLS/VISU/SRC] - 16