[Nek5000-users] Parallel issue: Inflow conditions from external files

nek5000-users at lists.mcs.anl.gov nek5000-users at lists.mcs.anl.gov
Wed Sep 10 06:08:05 CDT 2014


Hi all,

I have a question on how to read in information from external files as an
inflow boundary condition when running in parallel.
Similar problem to this:
http://lists.mcs.anl.gov/pipermail/nek5000-users/2010-December/001151.html

My geometry is a simple 3D rectangle prism with the inlet on the Y-Z plane
@ x=0.

I have developed some code that reads in a 2D plane of velocity data from
another simulation (in this case fluent) as the boundary condition in
userbc(). I am happy to share the process and the code (Fortran
input/output from NEK, and MATLAB to do interpolation) however here I am
primarily focused on getting it working in parallel.

The process works fine in serial - it reads in 3 files (01.dat, 02.dat,
03.dat) that hold the x,y,z velocity (respectively) for each node on the
inlet.

However it does not work correctly when running in parallel - looking at
the results after a few iterations shows that the velocities are being
mapped to nodes that are not all on the inlet. It does compile and run
though.

I'd sincerely appreciate any the help in making it work.

Warm regards,

Tom

I have included my .usr file below.

C-----------------------------
------------------------------------------
C  nek5000 user-file template
C
C  user specified routines:
C     - userbc : boundary conditions
C     - useric : initial conditions
C     - uservp : variable properties
C     - userf  : local acceleration term for fluid
C     - userq  : local source term for scalars
C     - userchk: general purpose routine for checking errors etc.
C
C-----------------------------------------------------------------------
      subroutine uservp(ix,iy,iz,eg) ! set variable properties
      include 'SIZE'
      include 'TOTAL'
      include 'NEKUSE'

      integer e,f,eg
c     e = gllel(eg)

      udiff  = 0.0
      utrans = 0.0

      return
      end
c-----------------------------------------------------------------------
      subroutine userf(ix,iy,iz,eg) ! set acceleration term
c
c     Note: this is an acceleration term, NOT a force!
c     Thus, ffx will subsequently be multiplied by rho(x,t).
c
      include 'SIZE'
      include 'TOTAL'
      include 'NEKUSE'

      integer e,f,eg

      ffx=0
      ffy=0
      ffz=0


      return
      end
c-----------------------------------------------------------------------
      subroutine userq(ix,iy,iz,eg) ! set source term
      include 'SIZE'
      include 'TOTAL'
      include 'NEKUSE'

      integer e,f,eg
c     e = gllel(eg)

      qvol   = 0.0
      source = 0.0

      return
      end
c-----------------------------------------------------------------------
      subroutine userbc(ix,iy,iz,iside,eg) ! set up boundary conditions

c     NOTE: This routine may or may not be called by every processor

      include 'SIZE'
      include 'TOTAL'
      include 'NEKUSE'
      integer e,eg,ix,iy,iz,NELY,NELZ,egs

      common /myinlet/ Arrangedu(ly1,lz1,lelt),
     >                 Arrangedv(ly1,lz1,lelt),
     >                 Arrangedw(ly1,lz1,lelt)

      e  = gllel(eg)
      ux   = Arrangedu(iy,iz,eg)
      uy   = Arrangedv(iy,iz,eg)
      uz   = Arrangedw(iy,iz,eg)
      temp = 0.0

      return
      end
c-----------------------------------------------------------------------
      subroutine useric(ix,iy,iz,eg) ! set up initial conditions
      include 'SIZE'
      include 'TOTAL'
      include 'NEKUSE'
      integer e,eg

      ux   = 7.235
      uy   = 0.0
      uz   = 0.0
      temp = 0.0

      return
      end
c-----------------------------------------------------------------------
      subroutine userchk()
      include 'SIZE'
      include 'TOTAL'

      return
      end
c-----------------------------------------------------------------------
      subroutine usrdat()   ! This routine to modify element vertices
      include 'SIZE'
      include 'TOTAL'

      return
      end
c-----------------------------------------------------------------------
      subroutine usrdat2()  ! This routine to modify mesh coordinates
      include 'SIZE'
      include 'TOTAL'

      common /myinlet/ Arrangedu(ly1,lz1,lelt),
     >                 Arrangedv(ly1,lz1,lelt),
     >                 Arrangedw(ly1,lz1,lelt)

      call read_inlet(Arrangedu,Arrangedv,Arrangedw) ! Read in inlet data
for BCs

      return
      end
c-----------------------------------------------------------------------
      subroutine usrdat3()
      include 'SIZE'
      include 'TOTAL'

      return
      end
C-------------------------------------------------------------------------------
      subroutine read_inlet(Arrangedu,Arrangedv,Arrangedw) ! Read in the
inlet data

      include 'SIZE'
      include 'TOTAL'

      Integer    Max_Nodes,          Lun_Dat
      Parameter( Max_Nodes = 200000, Lun_Dat = 3)

      real Var_1( Max_Nodes,lx1), Var_2( Max_Nodes), Var_3( Max_Nodes)
      real Var_4( Max_Nodes), Var_5( Max_Nodes), Var_6( Max_Nodes)
      real Arrangedu(ly1,lz1,lelt)
      real Arrangedv(ly1,lz1,lelt)
      real Arrangedw(ly1,lz1,lelt)

      integer i_row( Max_Nodes), i_eg( Max_Nodes), i_y( Max_Nodes)

      real junk

      integer i, j, k, l, m, n, imin, imax, N_Nodes,ze, iyz, egs,egz
      integer fcount

      Character Line* 72

      Character  Directory* 57
      Parameter( Directory=
     >   '/home/tom/NEK/meshControl/MultiBox/MultiBox_10/NewInlet4/')
c         1234567890123456789012345678901234567890123456789012345678

      Character  Filename* 4
      Parameter( Filename= '.dat')
c                           12345
      character(2) fc

C     +-----------+
C     | Execution |
C     +-----------+

1       format( 1x, a)

      do fcount=1,3 !for each file 01=u;02=v;03=w.

37    format (I2.2)
      write (fc,37) fcount ! converting integer to string
      print *,'Readining file #: ',fc

      open( unit = Lun_Dat, status = 'old', form='formatted', err = 2,
     >      file = Directory//trim(fc)//Filename)
      go to 3
2     write( 6, 1) 'Error Opening file - EXIT'


3     read( Lun_Dat, 1)  Line

      N_Nodes = 0
      do i = 1, Max_Nodes
         read( Lun_Dat, 11, end=5) i_row(i),  i_eg(i),  i_y(i),
     >                      Var_1(i,1), Var_1(i,2), Var_1(i,3),
     >                      Var_1(i,4), Var_1(i,5), Var_1(i,6),
     >                      Var_1(i,7), Var_1(i,8), Var_1(i,9)
11       format( i6, 1x, i5, 1x, i4, 9( 1x, 1P E16.9))
         N_Nodes = N_Nodes + 1

       if( fcount .EQ. 1) THEN
         iyz=i_y(i)
         egz=i_eg(i)
         do ze=1,lx1
         Arrangedu(iyz,ze,egz)= Var_1(i,ze)
         end do
       else if( fcount .EQ. 2) THEN
         iyz=i_y(i)
         egz=i_eg(i)
         do ze=1,lx1
         Arrangedv(iyz,ze,egz)= Var_1(i,ze)
         end do
      else if( fcount .EQ. 3) THEN
         iyz=i_y(i)
         egz=i_eg(i)
         do ze=1,lx1
         Arrangedw(iyz,ze,egz)= Var_1(i,ze)
         end do
      end if

      end do
      write( 6, 1) ' '
      write( 6, 1) ' WARNING - Max nodes exceeeded'
      write( 6, 1) ' '
5     write( 6, 6) N_Nodes
6     format( ' Number of nodes read = ', I8)

      end do

      end
C-------------------------------------------------------------------------------
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/nek5000-users/attachments/20140910/4e069ca4/attachment-0001.html>


More information about the Nek5000-users mailing list