| 1 | c-----------------------------------------------------------------------
|
|---|
| 2 | subroutine read_re2_data(ifbswap, ifxyz, ifcur, ifbc) ! .re2 reader
|
|---|
| 3 |
|
|---|
| 4 | include 'SIZE'
|
|---|
| 5 | include 'TOTAL'
|
|---|
| 6 | include 'RESTART'
|
|---|
| 7 | include 'CTIMER'
|
|---|
| 8 |
|
|---|
| 9 | logical ifbswap
|
|---|
| 10 | logical ifxyz, ifcur, ifbc
|
|---|
| 11 | integer idummy(100)
|
|---|
| 12 |
|
|---|
| 13 | common /nekmpi/ nidd,npp,nekcomm,nekgroup,nekreal
|
|---|
| 14 |
|
|---|
| 15 |
|
|---|
| 16 | etime0 = dnekclock_sync()
|
|---|
| 17 |
|
|---|
| 18 | ibc = 2
|
|---|
| 19 | if (ifflow) ibc = 1
|
|---|
| 20 |
|
|---|
| 21 | nfldt = 1
|
|---|
| 22 | if (ifheat) nfldt = 2+npscal
|
|---|
| 23 | if (ifmhd ) nfldt = 2+npscal+1
|
|---|
| 24 |
|
|---|
| 25 | ! first field to read
|
|---|
| 26 | if (param(33).gt.0) ibc = int(param(33))
|
|---|
| 27 |
|
|---|
| 28 | ! number of fields to read
|
|---|
| 29 | if (param(32).gt.0) nfldt = ibc + int(param(32)) - 1
|
|---|
| 30 |
|
|---|
| 31 | call blank(cbc,3*size(cbc))
|
|---|
| 32 | call rzero(bc ,size(bc))
|
|---|
| 33 |
|
|---|
| 34 | #ifndef NOMPIIO
|
|---|
| 35 | call fgslib_crystal_setup(cr_re2,nekcomm,np)
|
|---|
| 36 |
|
|---|
| 37 | call byte_open_mpi(re2fle,fh_re2,.true.,ierr)
|
|---|
| 38 | call err_chk(ierr,' Cannot open .re2 file!$')
|
|---|
| 39 |
|
|---|
| 40 | call readp_re2_mesh (ifbswap,ifxyz)
|
|---|
| 41 | call readp_re2_curve(ifbswap,ifcur)
|
|---|
| 42 | do ifield = ibc,nfldt
|
|---|
| 43 | call readp_re2_bc(cbc(1,1,ifield),bc(1,1,1,ifield),
|
|---|
| 44 | & ifbswap,ifbc)
|
|---|
| 45 | enddo
|
|---|
| 46 |
|
|---|
| 47 | call fgslib_crystal_free(cr_re2)
|
|---|
| 48 | call byte_close_mpi(fh_re2,ierr)
|
|---|
| 49 | #else
|
|---|
| 50 | call byte_open(re2fle,ierr)
|
|---|
| 51 | call byte_read(idummy,21,ierr) ! skip hdr+endian code
|
|---|
| 52 |
|
|---|
| 53 | call bin_rd1_mesh (ifbswap)
|
|---|
| 54 | call bin_rd1_curve(ifbswap)
|
|---|
| 55 | do ifield = ibc,nfldt
|
|---|
| 56 | call bin_rd1_bc (cbc(1,1,ifield),bc(1,1,1,ifield),ifbswap)
|
|---|
| 57 | enddo
|
|---|
| 58 |
|
|---|
| 59 | call byte_close(ierr)
|
|---|
| 60 | #endif
|
|---|
| 61 |
|
|---|
| 62 | etime_t = dnekclock_sync() - etime0
|
|---|
| 63 | if(nio.eq.0) write(6,'(A,1(1g9.2),A,/)')
|
|---|
| 64 | & ' done :: read .re2 file ',
|
|---|
| 65 | & etime_t, ' sec'
|
|---|
| 66 |
|
|---|
| 67 | return
|
|---|
| 68 | end
|
|---|
| 69 | c-----------------------------------------------------------------------
|
|---|
| 70 | subroutine readp_re2_mesh(ifbswap,ifread) ! version 2 of .re2 reader
|
|---|
| 71 |
|
|---|
| 72 | include 'SIZE'
|
|---|
| 73 | include 'TOTAL'
|
|---|
| 74 |
|
|---|
| 75 | logical ifbswap
|
|---|
| 76 | logical ifread
|
|---|
| 77 |
|
|---|
| 78 | parameter(nrmax = lelt) ! maximum number of records
|
|---|
| 79 | parameter(lrs = 1+ldim*(2**ldim)) ! record size: group x(:,c) ...
|
|---|
| 80 | parameter(li = 2*lrs+2)
|
|---|
| 81 |
|
|---|
| 82 | integer bufr(li-2,nrmax)
|
|---|
| 83 | common /scrns/ bufr
|
|---|
| 84 |
|
|---|
| 85 | integer vi (li ,nrmax)
|
|---|
| 86 | common /ctmp1/ vi
|
|---|
| 87 |
|
|---|
| 88 | integer*8 lre2off_b,dtmp8
|
|---|
| 89 | integer*8 nrg
|
|---|
| 90 |
|
|---|
| 91 | nrg = nelgt
|
|---|
| 92 | nr = nelt
|
|---|
| 93 | irankoff = igl_running_sum(nr) - nr
|
|---|
| 94 | dtmp8 = irankoff
|
|---|
| 95 | re2off_b = 84 ! set initial offset (hdr + endian)
|
|---|
| 96 | lre2off_b = re2off_b + dtmp8*lrs*wdsizi
|
|---|
| 97 | lrs4 = lrs*wdsizi/4
|
|---|
| 98 |
|
|---|
| 99 | ! read coordinates from file
|
|---|
| 100 | nwds4r = nr*lrs4
|
|---|
| 101 | call byte_set_view(lre2off_b,fh_re2)
|
|---|
| 102 | call byte_read_mpi(bufr,nwds4r,-1,fh_re2,ierr)
|
|---|
| 103 | re2off_b = re2off_b + nrg*4*lrs4
|
|---|
| 104 | if (ierr.gt.0) goto 100
|
|---|
| 105 |
|
|---|
| 106 | if (.not.ifread) return
|
|---|
| 107 |
|
|---|
| 108 | if (nio.eq.0) write(6,*) 'reading mesh '
|
|---|
| 109 |
|
|---|
| 110 | ! pack buffer
|
|---|
| 111 | do i = 1,nr
|
|---|
| 112 | jj = (i-1)*lrs4 + 1
|
|---|
| 113 | ielg = irankoff + i ! elements are stored in global order
|
|---|
| 114 | vi(1,i) = gllnid(ielg)
|
|---|
| 115 | vi(2,i) = ielg
|
|---|
| 116 | call icopy(vi(3,i),bufr(jj,1),lrs4)
|
|---|
| 117 | enddo
|
|---|
| 118 |
|
|---|
| 119 | ! crystal route nr real items of size lrs to rank vi(key,1:nr)
|
|---|
| 120 | n = nr
|
|---|
| 121 | key = 1
|
|---|
| 122 | call fgslib_crystal_tuple_transfer(cr_re2,n,nrmax,vi,li,
|
|---|
| 123 | & vl,0,vr,0,key)
|
|---|
| 124 |
|
|---|
| 125 | ! unpack buffer
|
|---|
| 126 | ierr = 0
|
|---|
| 127 | if (n.gt.nrmax) then
|
|---|
| 128 | ierr = 1
|
|---|
| 129 | goto 100
|
|---|
| 130 | endif
|
|---|
| 131 |
|
|---|
| 132 | do i = 1,n
|
|---|
| 133 | iel = gllel(vi(2,i))
|
|---|
| 134 | call icopy (bufr,vi(3,i),lrs4)
|
|---|
| 135 | call buf_to_xyz(bufr,iel,ifbswap,ierr)
|
|---|
| 136 | enddo
|
|---|
| 137 |
|
|---|
| 138 | 100 call err_chk(ierr,'Error reading .re2 mesh$')
|
|---|
| 139 |
|
|---|
| 140 | return
|
|---|
| 141 | end
|
|---|
| 142 | c-----------------------------------------------------------------------
|
|---|
| 143 | subroutine readp_re2_curve(ifbswap,ifread)
|
|---|
| 144 |
|
|---|
| 145 | include 'SIZE'
|
|---|
| 146 | include 'TOTAL'
|
|---|
| 147 |
|
|---|
| 148 | logical ifbswap
|
|---|
| 149 | logical ifread
|
|---|
| 150 |
|
|---|
| 151 | common /nekmpi/ nidd,npp,nekcomm,nekgroup,nekreal
|
|---|
| 152 |
|
|---|
| 153 | parameter(nrmax = 12*lelt) ! maximum number of records
|
|---|
| 154 | parameter(lrs = 2+1+5) ! record size: eg iside curve(5) ccurve
|
|---|
| 155 | parameter(li = 2*lrs+1)
|
|---|
| 156 |
|
|---|
| 157 | integer bufr(li-1,nrmax)
|
|---|
| 158 | common /scrns/ bufr
|
|---|
| 159 |
|
|---|
| 160 | integer vi (li ,nrmax)
|
|---|
| 161 | common /ctmp1/ vi
|
|---|
| 162 |
|
|---|
| 163 | integer*8 lre2off_b,dtmp8
|
|---|
| 164 | integer*8 nrg
|
|---|
| 165 | integer*4 nrg4(2)
|
|---|
| 166 |
|
|---|
| 167 |
|
|---|
| 168 | ! read total number of records
|
|---|
| 169 | nwds4r = 1*wdsizi/4
|
|---|
| 170 | lre2off_b = re2off_b
|
|---|
| 171 | call byte_set_view(lre2off_b,fh_re2)
|
|---|
| 172 | call byte_read_mpi(nrg4,nwds4r,-1,fh_re2,ierr)
|
|---|
| 173 | if(ierr.gt.0) goto 100
|
|---|
| 174 |
|
|---|
| 175 | if(wdsizi.eq.8) then
|
|---|
| 176 | if(ifbswap) call byte_reverse8(nrg4,nwds4r,ierr)
|
|---|
| 177 | call copy(dnrg,nrg4,1)
|
|---|
| 178 | nrg = dnrg
|
|---|
| 179 | else
|
|---|
| 180 | if(ifbswap) call byte_reverse (nrg4,nwds4r,ierr)
|
|---|
| 181 | nrg = nrg4(1)
|
|---|
| 182 | endif
|
|---|
| 183 | re2off_b = re2off_b + 4*nwds4r
|
|---|
| 184 |
|
|---|
| 185 | if(nrg.eq.0) return
|
|---|
| 186 |
|
|---|
| 187 | ! read data from file
|
|---|
| 188 | dtmp8 = np
|
|---|
| 189 | nr = nrg/dtmp8
|
|---|
| 190 | do i = 0,mod(nrg,dtmp8)-1
|
|---|
| 191 | if(i.eq.nid) nr = nr + 1
|
|---|
| 192 | enddo
|
|---|
| 193 | irankoff = igl_running_sum(nr) - nr
|
|---|
| 194 | dtmp8 = irankoff
|
|---|
| 195 | lre2off_b = re2off_b + dtmp8*lrs*wdsizi
|
|---|
| 196 | lrs4 = lrs*wdsizi/4
|
|---|
| 197 |
|
|---|
| 198 | re2off_b = re2off_b + nrg*4*lrs4
|
|---|
| 199 |
|
|---|
| 200 | if (.not.ifread) return
|
|---|
| 201 | if(nio.eq.0) write(6,*) 'reading curved sides '
|
|---|
| 202 |
|
|---|
| 203 | nwds4r = nr*lrs4
|
|---|
| 204 | call byte_set_view(lre2off_b,fh_re2)
|
|---|
| 205 | call byte_read_mpi(bufr,nwds4r,-1,fh_re2,ierr)
|
|---|
| 206 | if(ierr.gt.0) goto 100
|
|---|
| 207 |
|
|---|
| 208 | ! pack buffer
|
|---|
| 209 | do i = 1,nr
|
|---|
| 210 | jj = (i-1)*lrs4 + 1
|
|---|
| 211 |
|
|---|
| 212 | if(ifbswap) then
|
|---|
| 213 | lrs4s = lrs4 - wdsizi/4 ! words to swap (last is char)
|
|---|
| 214 | if(wdsizi.eq.8) call byte_reverse8(bufr(jj,1),lrs4s,ierr)
|
|---|
| 215 | if(wdsizi.eq.4) call byte_reverse (bufr(jj,1),lrs4s,ierr)
|
|---|
| 216 | endif
|
|---|
| 217 |
|
|---|
| 218 | ielg = bufr(jj,1)
|
|---|
| 219 | if(wdsizi.eq.8) call copyi4(ielg,bufr(jj,1),1)
|
|---|
| 220 |
|
|---|
| 221 | if(ielg.le.0 .or. ielg.gt.nelgt) goto 100
|
|---|
| 222 | vi(1,i) = gllnid(ielg)
|
|---|
| 223 |
|
|---|
| 224 | call icopy (vi(2,i),bufr(jj,1),lrs4)
|
|---|
| 225 | enddo
|
|---|
| 226 |
|
|---|
| 227 | ! crystal route nr real items of size lrs to rank vi(key,1:nr)
|
|---|
| 228 | n = nr
|
|---|
| 229 | key = 1
|
|---|
| 230 | call fgslib_crystal_tuple_transfer(cr_re2,n,nrmax,vi,li,vl,0,vr,0,
|
|---|
| 231 | & key)
|
|---|
| 232 |
|
|---|
| 233 | ! unpack buffer
|
|---|
| 234 | if(n.gt.nrmax) goto 100
|
|---|
| 235 | do i = 1,n
|
|---|
| 236 | call icopy (bufr,vi(2,i),lrs4)
|
|---|
| 237 | call buf_to_curve(bufr)
|
|---|
| 238 | enddo
|
|---|
| 239 |
|
|---|
| 240 | return
|
|---|
| 241 |
|
|---|
| 242 | 100 ierr = 1
|
|---|
| 243 | call err_chk(ierr,'Error reading .re2 curved data$')
|
|---|
| 244 |
|
|---|
| 245 | end
|
|---|
| 246 | c-----------------------------------------------------------------------
|
|---|
| 247 | subroutine readp_re2_bc(cbl,bl,ifbswap,ifread)
|
|---|
| 248 |
|
|---|
| 249 | include 'SIZE'
|
|---|
| 250 | include 'TOTAL'
|
|---|
| 251 |
|
|---|
| 252 | character*3 cbl( 6,lelt)
|
|---|
| 253 | real bl (5,6,lelt)
|
|---|
| 254 | logical ifbswap
|
|---|
| 255 | logical ifread
|
|---|
| 256 |
|
|---|
| 257 | parameter(nrmax = 6*lelt) ! maximum number of records
|
|---|
| 258 | parameter(lrs = 2+1+5) ! record size: eg iside bl(5) cbl
|
|---|
| 259 | parameter(li = 2*lrs+1)
|
|---|
| 260 |
|
|---|
| 261 | integer bufr(li-1,nrmax)
|
|---|
| 262 | common /scrns/ bufr
|
|---|
| 263 |
|
|---|
| 264 | integer vi (li ,nrmax)
|
|---|
| 265 | common /ctmp1/ vi
|
|---|
| 266 |
|
|---|
| 267 | integer*8 lre2off_b,dtmp8
|
|---|
| 268 | integer*8 nrg
|
|---|
| 269 | integer*4 nrg4(2)
|
|---|
| 270 |
|
|---|
| 271 | ! read total number of records
|
|---|
| 272 | nwds4r = 1*wdsizi/4
|
|---|
| 273 | lre2off_b = re2off_b
|
|---|
| 274 | call byte_set_view(lre2off_b,fh_re2)
|
|---|
| 275 | call byte_read_mpi(nrg4,nwds4r,-1,fh_re2,ierr)
|
|---|
| 276 | if(ierr.gt.0) goto 100
|
|---|
| 277 |
|
|---|
| 278 | if(wdsizi.eq.8) then
|
|---|
| 279 | if(ifbswap) call byte_reverse8(nrg4,nwds4r,ierr)
|
|---|
| 280 | call copy(dnrg,nrg4,1)
|
|---|
| 281 | nrg = dnrg
|
|---|
| 282 | else
|
|---|
| 283 | if(ifbswap) call byte_reverse (nrg4,nwds4r,ierr)
|
|---|
| 284 | nrg = nrg4(1)
|
|---|
| 285 | endif
|
|---|
| 286 | re2off_b = re2off_b + 4*nwds4r
|
|---|
| 287 |
|
|---|
| 288 | if(nrg.eq.0) return
|
|---|
| 289 |
|
|---|
| 290 | ! read data from file
|
|---|
| 291 | dtmp8 = np
|
|---|
| 292 | nr = nrg/dtmp8
|
|---|
| 293 | do i = 0,mod(nrg,dtmp8)-1
|
|---|
| 294 | if(i.eq.nid) nr = nr + 1
|
|---|
| 295 | enddo
|
|---|
| 296 | irankoff = igl_running_sum(nr) - nr
|
|---|
| 297 | dtmp8 = irankoff
|
|---|
| 298 | lre2off_b = re2off_b + dtmp8*lrs*wdsizi
|
|---|
| 299 | lrs4 = lrs*wdsizi/4
|
|---|
| 300 |
|
|---|
| 301 | re2off_b = re2off_b + nrg*4*lrs4
|
|---|
| 302 |
|
|---|
| 303 | if (.not.ifread) return
|
|---|
| 304 | if(nio.eq.0) write(6,*) 'reading bc for ifld',ifield
|
|---|
| 305 |
|
|---|
| 306 | nwds4r = nr*lrs4
|
|---|
| 307 | call byte_set_view(lre2off_b,fh_re2)
|
|---|
| 308 | call byte_read_mpi(bufr,nwds4r,-1,fh_re2,ierr)
|
|---|
| 309 | if(ierr.gt.0) goto 100
|
|---|
| 310 |
|
|---|
| 311 | ! pack buffer
|
|---|
| 312 | do i = 1,nr
|
|---|
| 313 | jj = (i-1)*lrs4 + 1
|
|---|
| 314 |
|
|---|
| 315 | if(ifbswap) then
|
|---|
| 316 | lrs4s = lrs4 - wdsizi/4 ! words to swap (last is char)
|
|---|
| 317 | if(wdsizi.eq.8) call byte_reverse8(bufr(jj,1),lrs4s,ierr)
|
|---|
| 318 | if(wdsizi.eq.4) call byte_reverse (bufr(jj,1),lrs4s,ierr)
|
|---|
| 319 | endif
|
|---|
| 320 |
|
|---|
| 321 | ielg = bufr(jj,1)
|
|---|
| 322 | if(wdsizi.eq.8) call copyi4(ielg,bufr(jj,1),1)
|
|---|
| 323 |
|
|---|
| 324 | if(ielg.le.0 .or. ielg.gt.nelgt) goto 100
|
|---|
| 325 | vi(1,i) = gllnid(ielg)
|
|---|
| 326 |
|
|---|
| 327 | call icopy (vi(2,i),bufr(jj,1),lrs4)
|
|---|
| 328 | enddo
|
|---|
| 329 |
|
|---|
| 330 | ! crystal route nr real items of size lrs to rank vi(key,1:nr)
|
|---|
| 331 | n = nr
|
|---|
| 332 | key = 1
|
|---|
| 333 |
|
|---|
| 334 | call fgslib_crystal_tuple_transfer(cr_re2,n,nrmax,vi,li,vl,0,vr,0,
|
|---|
| 335 | & key)
|
|---|
| 336 |
|
|---|
| 337 | ! fill up with default
|
|---|
| 338 | do iel=1,nelt
|
|---|
| 339 | do k=1,6
|
|---|
| 340 | cbl(k,iel) = 'E '
|
|---|
| 341 | enddo
|
|---|
| 342 | enddo
|
|---|
| 343 |
|
|---|
| 344 | ! unpack buffer
|
|---|
| 345 | if(n.gt.nrmax) goto 100
|
|---|
| 346 | do i = 1,n
|
|---|
| 347 | call icopy (bufr,vi(2,i),lrs4)
|
|---|
| 348 | call buf_to_bc(cbl,bl,bufr)
|
|---|
| 349 | enddo
|
|---|
| 350 |
|
|---|
| 351 | return
|
|---|
| 352 |
|
|---|
| 353 | 100 ierr = 1
|
|---|
| 354 | call err_chk(ierr,'Error reading .re2 boundary data$')
|
|---|
| 355 |
|
|---|
| 356 | end
|
|---|
| 357 | c-----------------------------------------------------------------------
|
|---|
| 358 | subroutine buf_to_xyz(buf,e,ifbswap,ierr)! version 1 of binary reader
|
|---|
| 359 |
|
|---|
| 360 | include 'SIZE'
|
|---|
| 361 | include 'TOTAL'
|
|---|
| 362 | logical ifbswap
|
|---|
| 363 |
|
|---|
| 364 | c integer e,eg,buf(0:49)
|
|---|
| 365 | integer e,eg,buf(0:49)
|
|---|
| 366 |
|
|---|
| 367 | nwds = (1 + ldim*(2**ldim))*(wdsizi/4) ! group + 2x4 for 2d, 3x8 for 3d
|
|---|
| 368 |
|
|---|
| 369 | if (ifbswap.and.ierr.eq.0.and.wdsizi.eq.8) then
|
|---|
| 370 | call byte_reverse8(buf,nwds,ierr)
|
|---|
| 371 | elseif (ifbswap.and.ierr.eq.0.and.wdsizi.eq.4) then
|
|---|
| 372 | call byte_reverse (buf,nwds,ierr)
|
|---|
| 373 | endif
|
|---|
| 374 | if(ierr.ne.0) return
|
|---|
| 375 |
|
|---|
| 376 | if(wdsizi.eq.8) then
|
|---|
| 377 | call copyi4(igroup(e),buf(0),1) !0-1
|
|---|
| 378 | if (ldim.eq.3) then
|
|---|
| 379 | call copy (xc(1,e),buf( 2),8) !2 --17
|
|---|
| 380 | call copy (yc(1,e),buf(18),8) !18--33
|
|---|
| 381 | call copy (zc(1,e),buf(34),8) !34--49
|
|---|
| 382 | else
|
|---|
| 383 | call copy (xc(1,e),buf( 2),4) !2 --9
|
|---|
| 384 | call copy (yc(1,e),buf(10),4) !10--17
|
|---|
| 385 | endif
|
|---|
| 386 | else
|
|---|
| 387 | igroup(e) = buf(0)
|
|---|
| 388 | if (if3d) then
|
|---|
| 389 | call copy4r(xc(1,e),buf( 1),8)
|
|---|
| 390 | call copy4r(yc(1,e),buf( 9),8)
|
|---|
| 391 | call copy4r(zc(1,e),buf(17),8)
|
|---|
| 392 | else
|
|---|
| 393 | call copy4r(xc(1,e),buf( 1),4)
|
|---|
| 394 | call copy4r(yc(1,e),buf( 5),4)
|
|---|
| 395 | endif
|
|---|
| 396 | endif
|
|---|
| 397 |
|
|---|
| 398 | return
|
|---|
| 399 | end
|
|---|
| 400 | c-----------------------------------------------------------------------
|
|---|
| 401 | subroutine buf_to_curve(buf) ! version 1 of binary reader
|
|---|
| 402 |
|
|---|
| 403 | include 'SIZE'
|
|---|
| 404 | include 'TOTAL'
|
|---|
| 405 |
|
|---|
| 406 | integer e,eg,f,buf(30)
|
|---|
| 407 |
|
|---|
| 408 | if(wdsizi.eq.8) then
|
|---|
| 409 | call copyi4(eg,buf(1),1) !1-2
|
|---|
| 410 | e = gllel(eg)
|
|---|
| 411 |
|
|---|
| 412 | call copyi4(f,buf(3),1) !3-4
|
|---|
| 413 |
|
|---|
| 414 | call copy ( curve(1,f,e),buf(5) ,5) !5--14
|
|---|
| 415 | call chcopy(ccurve( f,e),buf(15),1)!15
|
|---|
| 416 | else
|
|---|
| 417 | eg = buf(1)
|
|---|
| 418 | e = gllel(eg)
|
|---|
| 419 | f = buf(2)
|
|---|
| 420 |
|
|---|
| 421 | call copy4r( curve(1,f,e),buf(3),5)
|
|---|
| 422 | call chcopy(ccurve(f,e) ,buf(8),1)
|
|---|
| 423 | endif
|
|---|
| 424 |
|
|---|
| 425 | c write(6,1) eg,e,f,(curve(k,f,e),k=1,5),ccurve(f,e)
|
|---|
| 426 | c 1 format(2i7,i3,5f10.3,1x,a1,'ccurve')
|
|---|
| 427 |
|
|---|
| 428 | return
|
|---|
| 429 | end
|
|---|
| 430 | c-----------------------------------------------------------------------
|
|---|
| 431 | subroutine buf_to_bc(cbl,bl,buf) ! version 1 of binary reader
|
|---|
| 432 |
|
|---|
| 433 | include 'SIZE'
|
|---|
| 434 | include 'TOTAL'
|
|---|
| 435 |
|
|---|
| 436 | character*3 cbl(6,lelt)
|
|---|
| 437 | real bl(5,6,lelt)
|
|---|
| 438 |
|
|---|
| 439 | integer e,eg,f,buf(30)
|
|---|
| 440 |
|
|---|
| 441 | if(wdsizi.eq.8) then
|
|---|
| 442 | call copyi4(eg,buf(1),1) !1-2
|
|---|
| 443 | e = gllel(eg)
|
|---|
| 444 |
|
|---|
| 445 | call copyi4(f,buf(3),1) !3-4
|
|---|
| 446 |
|
|---|
| 447 | call copy (bl(1,f,e),buf(5),5) !5--14
|
|---|
| 448 | call chcopy(cbl( f,e),buf(15),3)!15-16
|
|---|
| 449 |
|
|---|
| 450 | if(nelt.ge.1000000.and.cbl(f,e).eq.'P ')
|
|---|
| 451 | $ call copyi4(bl(1,f,e),buf(5),1) !Integer assign connecting P element
|
|---|
| 452 |
|
|---|
| 453 | else
|
|---|
| 454 | eg = buf(1)
|
|---|
| 455 | e = gllel(eg)
|
|---|
| 456 | f = buf(2)
|
|---|
| 457 |
|
|---|
| 458 | call copy4r ( bl(1,f,e),buf(3),5)
|
|---|
| 459 | call chcopy (cbl( f,e),buf(8),3)
|
|---|
| 460 |
|
|---|
| 461 | if (nelgt.ge.1000000.and.cbl(f,e).eq.'P ')
|
|---|
| 462 | $ bl(1,f,e) = buf(3) ! Integer assign of connecting periodic element
|
|---|
| 463 | endif
|
|---|
| 464 |
|
|---|
| 465 | c write(6,1) eg,e,f,cbl(f,e),' CBC',nid
|
|---|
| 466 | c 1 format(2i8,i4,2x,a3,a4,i8)
|
|---|
| 467 |
|
|---|
| 468 | return
|
|---|
| 469 | end
|
|---|
| 470 | c-----------------------------------------------------------------------
|
|---|
| 471 | subroutine bin_rd1_mesh(ifbswap) ! version 1 of binary reader
|
|---|
| 472 |
|
|---|
| 473 | include 'SIZE'
|
|---|
| 474 | include 'TOTAL'
|
|---|
| 475 | logical ifbswap
|
|---|
| 476 |
|
|---|
| 477 | integer e,eg,buf(55)
|
|---|
| 478 |
|
|---|
| 479 | if (nio.eq.0) write(6,*) ' reading mesh '
|
|---|
| 480 |
|
|---|
| 481 | nwds = (1 + ldim*(2**ldim))*(wdsizi/4) ! group + 2x4 for 2d, 3x8 for 3d
|
|---|
| 482 | len = 4*nwds ! 4 bytes / wd
|
|---|
| 483 |
|
|---|
| 484 | if (nwds.gt.55.or.isize.gt.4) then
|
|---|
| 485 | write(6,*) nid,' Error in bin_rd1_mesh: buf size',nwds,isize
|
|---|
| 486 | call exitt
|
|---|
| 487 | endif
|
|---|
| 488 |
|
|---|
| 489 | call nekgsync()
|
|---|
| 490 |
|
|---|
| 491 | niop = 10
|
|---|
| 492 | do k=1,8
|
|---|
| 493 | if (nelgt/niop .lt. 100) goto 10
|
|---|
| 494 | niop = niop*10
|
|---|
| 495 | enddo
|
|---|
| 496 | 10 continue
|
|---|
| 497 |
|
|---|
| 498 | ierr = 0
|
|---|
| 499 | ierr2 = 0
|
|---|
| 500 | len1 = 4
|
|---|
| 501 | do eg=1,nelgt ! sync NOT needed here
|
|---|
| 502 |
|
|---|
| 503 | mid = gllnid(eg)
|
|---|
| 504 | e = gllel (eg)
|
|---|
| 505 | #ifdef DEBUG
|
|---|
| 506 | if (nio.eq.0.and.mod(eg,niop).eq.0) write(6,*) eg,' mesh read'
|
|---|
| 507 | #endif
|
|---|
| 508 | if (mid.ne.nid.and.nid.eq.0) then ! read & send
|
|---|
| 509 |
|
|---|
| 510 | if(ierr.eq.0) then
|
|---|
| 511 | call byte_read (buf,nwds,ierr)
|
|---|
| 512 | call csend(e,ierr,len1,mid,0)
|
|---|
| 513 | if(ierr.eq.0) call csend(e,buf,len,mid,0)
|
|---|
| 514 | else
|
|---|
| 515 | call csend(e,ierr,len1,mid,0)
|
|---|
| 516 | endif
|
|---|
| 517 |
|
|---|
| 518 | elseif (mid.eq.nid.and.nid.ne.0) then ! recv & process
|
|---|
| 519 |
|
|---|
| 520 | call crecv (e,ierr,len1)
|
|---|
| 521 | if(ierr.eq.0) then
|
|---|
| 522 | call crecv (e,buf,len)
|
|---|
| 523 | call buf_to_xyz (buf,e,ifbswap,ierr2)
|
|---|
| 524 | endif
|
|---|
| 525 |
|
|---|
| 526 | elseif (mid.eq.nid.and.nid.eq.0) then ! read & process
|
|---|
| 527 |
|
|---|
| 528 | if(ierr.eq.0) then
|
|---|
| 529 | call byte_read (buf,nwds,ierr)
|
|---|
| 530 | call buf_to_xyz (buf,e,ifbswap,ierr2)
|
|---|
| 531 | endif
|
|---|
| 532 | endif
|
|---|
| 533 |
|
|---|
| 534 | enddo
|
|---|
| 535 | ierr = ierr + ierr2
|
|---|
| 536 | call err_chk(ierr,'Error reading .re2 mesh. Abort. $')
|
|---|
| 537 |
|
|---|
| 538 | return
|
|---|
| 539 | end
|
|---|
| 540 | c-----------------------------------------------------------------------
|
|---|
| 541 | subroutine bin_rd1_curve (ifbswap) ! v. 1 of curve side reader
|
|---|
| 542 |
|
|---|
| 543 | include 'SIZE'
|
|---|
| 544 | include 'TOTAL'
|
|---|
| 545 | logical ifbswap
|
|---|
| 546 |
|
|---|
| 547 | integer e,eg,buf(55)
|
|---|
| 548 | real rcurve
|
|---|
| 549 |
|
|---|
| 550 | nwds = (2 + 1 + 5)*(wdsizi/4) !eg+iside+ccurve+curve(6,:,:) !only 5 in rea
|
|---|
| 551 | len = 4*nwds ! 4 bytes / wd
|
|---|
| 552 |
|
|---|
| 553 | if (nwds.gt.55.or.isize.gt.4) then
|
|---|
| 554 | write(6,*)nid,' Error in bin_rd1_curve: buf size',nwds,isize
|
|---|
| 555 | call exitt
|
|---|
| 556 | endif
|
|---|
| 557 |
|
|---|
| 558 | call nekgsync()
|
|---|
| 559 |
|
|---|
| 560 | ierr = 0
|
|---|
| 561 | len1 = 4
|
|---|
| 562 | if (nid.eq.0) then ! read & send/process
|
|---|
| 563 |
|
|---|
| 564 | if(wdsizi.eq.8) then
|
|---|
| 565 | call byte_read(rcurve,2,ierr)
|
|---|
| 566 | if (ifbswap) call byte_reverse8(rcurve,2,ierr)
|
|---|
| 567 | ncurve = rcurve
|
|---|
| 568 | else
|
|---|
| 569 | call byte_read(ncurve,1,ierr)
|
|---|
| 570 | if (ifbswap) call byte_reverse(ncurve,1,ierr)
|
|---|
| 571 | endif
|
|---|
| 572 |
|
|---|
| 573 | if(ncurve.ne.0) write(6,*) ' reading curved sides '
|
|---|
| 574 | do k=1,ncurve
|
|---|
| 575 | if(ierr.eq.0) then
|
|---|
| 576 | call byte_read(buf,nwds,ierr)
|
|---|
| 577 | if(wdsizi.eq.8) then
|
|---|
| 578 | if(ifbswap) call byte_reverse8(buf,nwds-2,ierr)
|
|---|
| 579 | call copyi4(eg,buf(1),1) !1,2
|
|---|
| 580 | else
|
|---|
| 581 | if (ifbswap) call byte_reverse(buf,nwds-1,ierr) ! last is char
|
|---|
| 582 | eg = buf(1)
|
|---|
| 583 | endif
|
|---|
| 584 |
|
|---|
| 585 | mid = gllnid(eg)
|
|---|
| 586 | if (mid.eq.0.and.ierr.eq.0) then
|
|---|
| 587 | call buf_to_curve(buf)
|
|---|
| 588 | else
|
|---|
| 589 | if(ierr.eq.0) then
|
|---|
| 590 | call csend(mid,buf,len,mid,0)
|
|---|
| 591 | else
|
|---|
| 592 | goto 98
|
|---|
| 593 | endif
|
|---|
| 594 | endif
|
|---|
| 595 | else
|
|---|
| 596 | goto 98
|
|---|
| 597 | endif
|
|---|
| 598 | enddo
|
|---|
| 599 | 98 call buf_close_out ! notify all procs: no more data
|
|---|
| 600 |
|
|---|
| 601 | else ! wait for data from node 0
|
|---|
| 602 |
|
|---|
| 603 | ncurve_mx = 12*nelt
|
|---|
| 604 | do k=1,ncurve_mx+1 ! +1 to make certain we receive the close-out
|
|---|
| 605 |
|
|---|
| 606 | call crecv(nid,buf,len)
|
|---|
| 607 | if(wdsizi.eq.8) then
|
|---|
| 608 | call copyi4(ichk,buf(1),1)
|
|---|
| 609 | if(ichk.eq.0) goto 99
|
|---|
| 610 | call buf_to_curve(buf)
|
|---|
| 611 | elseif (buf(1).eq.0) then
|
|---|
| 612 | goto 99
|
|---|
| 613 | else
|
|---|
| 614 | call buf_to_curve(buf)
|
|---|
| 615 | endif
|
|---|
| 616 |
|
|---|
| 617 | enddo
|
|---|
| 618 | 99 call buf_close_out
|
|---|
| 619 |
|
|---|
| 620 | endif
|
|---|
| 621 | call err_chk(ierr,'Error reading .re2 curved data. Abort.$')
|
|---|
| 622 |
|
|---|
| 623 |
|
|---|
| 624 | return
|
|---|
| 625 | end
|
|---|
| 626 | c-----------------------------------------------------------------------
|
|---|
| 627 | subroutine bin_rd1_bc (cbl,bl,ifbswap) ! v. 1 of bc reader
|
|---|
| 628 |
|
|---|
| 629 | include 'SIZE'
|
|---|
| 630 | include 'TOTAL'
|
|---|
| 631 | logical ifbswap
|
|---|
| 632 |
|
|---|
| 633 | character*3 cbl(6,lelt)
|
|---|
| 634 | real bl(5,6,lelt)
|
|---|
| 635 |
|
|---|
| 636 | integer e,eg,buf(55)
|
|---|
| 637 | real rbc_max
|
|---|
| 638 |
|
|---|
| 639 | nwds = (2 + 1 + 5)*(wdsizi/4) ! eg + iside + cbc + bc(5,:,:)
|
|---|
| 640 | len = 4*nwds ! 4 bytes / wd
|
|---|
| 641 |
|
|---|
| 642 | if (nwds.gt.55.or.isize.gt.4) then
|
|---|
| 643 | write(6,*) nid,' Error in bin_rd1_bc: buf size',nwds,isize
|
|---|
| 644 | call exitt
|
|---|
| 645 | endif
|
|---|
| 646 |
|
|---|
| 647 | do e=1,nelt ! fill up cbc w/ default
|
|---|
| 648 | do k=1,6
|
|---|
| 649 | cbl(k,e) = 'E '
|
|---|
| 650 | enddo
|
|---|
| 651 | enddo
|
|---|
| 652 |
|
|---|
| 653 | call nekgsync()
|
|---|
| 654 | ierr=0
|
|---|
| 655 | len1=4
|
|---|
| 656 | if (nid.eq.0) then ! read & send/process
|
|---|
| 657 |
|
|---|
| 658 | if(wdsizi.eq.8) then
|
|---|
| 659 | call byte_read(rbc_max,2,ierr)
|
|---|
| 660 | if (ifbswap) call byte_reverse8(rbc_max,2,ierr) ! last is char
|
|---|
| 661 | nbc_max = rbc_max
|
|---|
| 662 | else
|
|---|
| 663 | call byte_read(nbc_max,1,ierr)
|
|---|
| 664 | if (ifbswap) call byte_reverse(nbc_max,1,ierr) ! last is char
|
|---|
| 665 | endif
|
|---|
| 666 |
|
|---|
| 667 | if(nbc_max.ne.0) write(6,*) ' reading bc for ifld',ifield
|
|---|
| 668 | do k=1,nbc_max
|
|---|
| 669 | c write(6,*) k,' dobc1 ',nbc_max
|
|---|
| 670 | if(ierr.eq.0) then
|
|---|
| 671 | call byte_read(buf,nwds,ierr)
|
|---|
| 672 | if(wdsizi.eq.8) then
|
|---|
| 673 | if (ifbswap) call byte_reverse8(buf,nwds-2,ierr)
|
|---|
| 674 | call copyi4(eg,buf(1),1) !1&2 of buf
|
|---|
| 675 | else
|
|---|
| 676 | if (ifbswap) call byte_reverse(buf,nwds-1,ierr) ! last is char
|
|---|
| 677 | eg = buf(1)
|
|---|
| 678 | endif
|
|---|
| 679 | mid = gllnid(eg)
|
|---|
| 680 | c write(6,*) k,' dobc3 ',eg,mid
|
|---|
| 681 |
|
|---|
| 682 | if (mid.eq.0.and.ierr.eq.0) then
|
|---|
| 683 | call buf_to_bc(cbl,bl,buf)
|
|---|
| 684 | else
|
|---|
| 685 | c write(6,*) mid,' sendbc1 ',eg
|
|---|
| 686 | if(ierr.eq.0) then
|
|---|
| 687 | call csend(mid,buf,len,mid,0)
|
|---|
| 688 | else
|
|---|
| 689 | goto 98
|
|---|
| 690 | endif
|
|---|
| 691 | c write(6,*) mid,' sendbc2 ',eg
|
|---|
| 692 | endif
|
|---|
| 693 | c write(6,*) k,' dobc2 ',nbc_max,eg
|
|---|
| 694 | else
|
|---|
| 695 | goto 98
|
|---|
| 696 | endif
|
|---|
| 697 | enddo
|
|---|
| 698 | c write(6,*) mid,' bclose ',eg,nbc_max
|
|---|
| 699 | 98 call buf_close_outv ! notify all procs: no more data
|
|---|
| 700 |
|
|---|
| 701 | else ! wait for data from node 0
|
|---|
| 702 |
|
|---|
| 703 | nbc_max = 2*ldim*nelt
|
|---|
| 704 | do k=1,nbc_max+1 ! Need one extra !
|
|---|
| 705 |
|
|---|
| 706 | c write(6,*) nid,' recvbc1',k
|
|---|
| 707 | call crecv(nid,buf,len)
|
|---|
| 708 | c write(6,*) nid,' recvbc2',k,buf(1)
|
|---|
| 709 |
|
|---|
| 710 | if(wdsizi.eq.8) then
|
|---|
| 711 | call copyi4(ichk,buf(1),1)
|
|---|
| 712 | if(ichk.eq.0) goto 99
|
|---|
| 713 | call buf_to_bc(cbl,bl,buf)
|
|---|
| 714 | elseif (buf(1).eq.0) then
|
|---|
| 715 | goto 99
|
|---|
| 716 | else
|
|---|
| 717 | call buf_to_bc(cbl,bl,buf)
|
|---|
| 718 | endif
|
|---|
| 719 |
|
|---|
| 720 | enddo
|
|---|
| 721 | 99 call buf_close_outv
|
|---|
| 722 |
|
|---|
| 723 | endif
|
|---|
| 724 |
|
|---|
| 725 | call err_chk(ierr,'Error reading boundary data for re2. Abort.$')
|
|---|
| 726 |
|
|---|
| 727 | return
|
|---|
| 728 | end
|
|---|
| 729 | c-----------------------------------------------------------------------
|
|---|
| 730 | subroutine buf_close_outv ! this is the stupid O(P) formulation
|
|---|
| 731 |
|
|---|
| 732 | include 'SIZE'
|
|---|
| 733 | include 'PARALLEL'
|
|---|
| 734 | integer*4 zero
|
|---|
| 735 | real rzero
|
|---|
| 736 |
|
|---|
| 737 | len = wdsizi
|
|---|
| 738 | rzero = 0
|
|---|
| 739 | zero = 0
|
|---|
| 740 | c write(6,*) nid,' bufclose'
|
|---|
| 741 | if (nid.eq.0) then
|
|---|
| 742 | do mid=1,np-1
|
|---|
| 743 | if(wdsizi.eq.8)call csend(mid,rzero,len,mid,0)
|
|---|
| 744 | if(wdsizi.eq.4)call csend(mid, zero,len,mid,0)
|
|---|
| 745 | c write(6,*) mid,' sendclose'
|
|---|
| 746 | enddo
|
|---|
| 747 | endif
|
|---|
| 748 |
|
|---|
| 749 | return
|
|---|
| 750 | end
|
|---|
| 751 | c-----------------------------------------------------------------------
|
|---|
| 752 | subroutine buf_close_out ! this is the stupid O(P) formulation
|
|---|
| 753 |
|
|---|
| 754 | include 'SIZE'
|
|---|
| 755 | include 'PARALLEL'
|
|---|
| 756 | integer*4 zero
|
|---|
| 757 | real rzero
|
|---|
| 758 |
|
|---|
| 759 | c len = 4
|
|---|
| 760 | len = wdsizi
|
|---|
| 761 | zero = 0
|
|---|
| 762 | rzero = 0
|
|---|
| 763 | if (nid.eq.0) then
|
|---|
| 764 | do mid=1,np-1
|
|---|
| 765 | if(wdsizi.eq.8)call csend(mid,rzero,len,mid,0)
|
|---|
| 766 | if(wdsizi.eq.4)call csend(mid, zero,len,mid,0)
|
|---|
| 767 | enddo
|
|---|
| 768 | endif
|
|---|
| 769 |
|
|---|
| 770 | return
|
|---|
| 771 | end
|
|---|
| 772 | c-----------------------------------------------------------------------
|
|---|
| 773 | subroutine read_re2_hdr(ifbswap, ifverbose) ! open file & chk for byteswap
|
|---|
| 774 |
|
|---|
| 775 | include 'SIZE'
|
|---|
| 776 | include 'TOTAL'
|
|---|
| 777 |
|
|---|
| 778 | logical ifbswap, ifverbose
|
|---|
| 779 | logical if_byte_swap_test
|
|---|
| 780 |
|
|---|
| 781 | integer fnami (33)
|
|---|
| 782 | character*132 fname
|
|---|
| 783 | equivalence (fname,fnami)
|
|---|
| 784 |
|
|---|
| 785 | character*132 hdr
|
|---|
| 786 | character*5 version
|
|---|
| 787 | real*4 test
|
|---|
| 788 |
|
|---|
| 789 | logical iffound
|
|---|
| 790 |
|
|---|
| 791 | ierr=0
|
|---|
| 792 |
|
|---|
| 793 | if (nid.eq.0) then
|
|---|
| 794 | if (ifverbose) write(6,'(A,A)') ' Reading ', re2fle
|
|---|
| 795 | call izero(fnami,33)
|
|---|
| 796 | m = indx2(re2fle,132,' ',1)-1
|
|---|
| 797 | call chcopy(fname,re2fle,m)
|
|---|
| 798 |
|
|---|
| 799 | inquire(file=fname, exist=iffound)
|
|---|
| 800 | if(.not.iffound) ierr = 1
|
|---|
| 801 | endif
|
|---|
| 802 | call err_chk(ierr,' Cannot find re2 file!$')
|
|---|
| 803 |
|
|---|
| 804 | if (nid.eq.0) then
|
|---|
| 805 | call byte_open(fname,ierr)
|
|---|
| 806 | if(ierr.ne.0) goto 100
|
|---|
| 807 | call byte_read(hdr,20,ierr)
|
|---|
| 808 | if(ierr.ne.0) goto 100
|
|---|
| 809 |
|
|---|
| 810 | read (hdr,1) version,nelgt,ldimr,nelgv
|
|---|
| 811 | 1 format(a5,i9,i3,i9)
|
|---|
| 812 |
|
|---|
| 813 | wdsizi = 4
|
|---|
| 814 | if(version.eq.'#v002') wdsizi = 8
|
|---|
| 815 | if(version.eq.'#v003') then
|
|---|
| 816 | wdsizi = 8
|
|---|
| 817 | param(32) = 1
|
|---|
| 818 | endif
|
|---|
| 819 |
|
|---|
| 820 | call byte_read(test,1,ierr)
|
|---|
| 821 | if(ierr.ne.0) goto 100
|
|---|
| 822 | ifbswap = if_byte_swap_test(test,ierr)
|
|---|
| 823 | if(ierr.ne.0) goto 100
|
|---|
| 824 |
|
|---|
| 825 | call byte_close(ierr)
|
|---|
| 826 | endif
|
|---|
| 827 |
|
|---|
| 828 | 100 call err_chk(ierr,'Error reading re2 header$')
|
|---|
| 829 |
|
|---|
| 830 | call bcast(wdsizi, ISIZE)
|
|---|
| 831 | call bcast(ifbswap,LSIZE)
|
|---|
| 832 | call bcast(nelgv ,ISIZE)
|
|---|
| 833 | call bcast(nelgt ,ISIZE)
|
|---|
| 834 | call bcast(ldimr ,ISIZE)
|
|---|
| 835 | call bcast(param(32),WDSIZE)
|
|---|
| 836 |
|
|---|
| 837 | if(wdsize.eq.4.and.wdsizi.eq.8)
|
|---|
| 838 | $ call exitti('wdsize=4 & wdsizi(re2)=8 not compatible$',wdsizi)
|
|---|
| 839 |
|
|---|
| 840 | return
|
|---|
| 841 | end
|
|---|