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