[Nek5000-users] LES of a jet - Vreman model

nek5000-users at lists.mcs.anl.gov nek5000-users at lists.mcs.anl.gov
Wed Nov 25 11:05:06 CST 2015


Dear Neks,

I try to implement the Vreman LES model (see An eddy-viscosity 
subgrid-scale model for turbulent shear flow: Algebraic theory and 
applications, A. W. Vreman, Phys. Fl., 16, 10, 2004) to simulate a 
turbulent jet. To do so, I modified the turbChannel example provided in 
NEK. However, as soon as I impose in the uservp subroutine a 
space-varying viscosity, after some time steps, the simulation blows up 
even if this viscosity remains close from the standard 1/Re viscosity. 
Has anyone ever experienced this kind of issue? Does anyone have an idea 
on how to fix this problem? Below are the .usr file and the log of the 
simulation.

Thanks,
Best regards,

Léopold

-------------------------------------------
.usr file

c-----------------------------------------------------------------------
c
c     Based on the turbChannel example
c
c-----------------------------------------------------------------------
       subroutine uservp (ix,iy,iz,eg)
       include 'SIZE'
       include 'TOTAL'
       include 'NEKUSE'

       common /cdsmag/ ediff(lx1,ly1,lz1,lelv)
       integer ie,f,eg
       real aux
       ie = gllel(eg)

       udiff = ediff(ix,iy,iz,ie)
c      udiff = param(2)
       utrans=1.

       aux = ediff(ix,iy,iz,ie) - param(2)
       if (istep.eq.19) then
         write(6,*) aux
       endif
       return
       end
c-----------------------------------------------------------------------
       subroutine userf  (ix,iy,iz,eg)
       include 'SIZE'
       include 'TOTAL'
       include 'NEKUSE'

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

       ffx = 0.0
       ffy = 0.0
       ffz = 0.0

       return
       end
c-----------------------------------------------------------------------
       subroutine userq  (ix,iy,iz,ieg)
       include 'SIZE'
       include 'TOTAL'
       include 'NEKUSE'

       qvol = 0.

       return
       end
c-----------------------------------------------------------------------
       subroutine userchk
       include 'SIZE'
       include 'TOTAL'
       include 'RESTART'
       include 'ZPER'  ! for nelx,nely,nelz

       parameter(lt=lx1*ly1*lz1*lelv)
       real vort(lt,3),vort_pol(lt,3),w1(lt),w2(lt),rayon

       real x0(3)
       save x0

       integer icalld
       save    icalld
       data    icalld /0/

       character*1   snam1(80)
       character*1   f1nam1(80),f2nam1(80)
       character*80  f1name
       equivalence  (f1nam1,f1name)
       character*80  f2name
       equivalence  (f2nam1,f2name)

       real atime,timel
       save atime,timel

       integer icount
       save    icount

       common /cdsmag/ ediff(lx1,ly1,lz1,lelv)

       common /avg/    uravg(lx1,ly1,lz1,lelv)
      &             ,  utavg(lx1,ly1,lz1,lelv)
      &             ,  uzavg(lx1,ly1,lz1,lelv)

       common /lesleo/ sij (lx1*ly1*lz1,ldim,ldim)
      $              , dx2 (lx1,ly1,lz1,lelv,lp)
      $              , dy2 (lx1,ly1,lz1,lelv,lp)
      $              , dz2 (lx1,ly1,lz1,lelv,lp)
      $              , cdiff
       real sij,dx2,dy2,dz2,cdiff

       integer e
       logical ifverbose

       n=nx1*ny1*nz1*nelv

       pi    = 4.*atan(1.0)
       rho   = 1.
       dnu   = param(2)

       cdiff = 0.07

       ! add SGS term to RHS (explicit treatment in plan4)
       ! Explicit viscosity.
       if(ifsplit .and. param(30).gt.0) ifexplvis = .true.

c     Grid spacing.
c      if ((istep.lt.1).and.(nid.eq.0)) then
c      if (istep.lt.1) then
       if (nid.eq.(istep-1)) then
         call set_grid_spacing(dx2,dy2,dz2)
       endif

       ! Compute eddy viscosity using Vreman model
c     Waiting for several timesteps in order to be sure that gradients
c     are not zero
       if (istep.gt.17) then
           ifuservp = .true.
       else
           ifuservp = .false.
       endif
       if(istep.gt.16) then
         if(nid.eq.0) write(6,*) 'Calculating eddy visosity'
         do e=1,nelv
            call eddy_visc(ediff,e)
         enddo
c        call copy(t,ediff,n)
       endif

c
c     Below is just for postprocessing ...
c
c     Compute perturbation energy
       e2 = glsc3(vy,bm1,vy,n)+glsc3(vz,bm1,vz,n)
       e2 = e2/volvm1

       ifxyo = .true.                        ! Turn on xyz output
       if (istep.gt.iostep) ifxyo = .false.  ! Turn off xyz output after 
first dump

c     Average calculations.
       if(icalld.eq.0) then
         call rzero(uravg,n)
         call rzero(utavg,n)
         call rzero(uzavg,n)
         atime = 0.
         timel = time
         icalld = 1
       endif

       dtime = time - timel
       atime = atime + dtime

       if (atime.ne.0. .and. dtime.ne.0.) then
         beta      = dtime/atime
         alpha     = 1.-beta
         ifverbose = .false.

         call avg1(uravg,vx   ,alpha,beta,n,'uravg',ifverbose)
         call avg1(utavg,vy   ,alpha,beta,n,'utavg',ifverbose)
         call avg1(uzavg,vz   ,alpha,beta,n,'uzavg',ifverbose)
       endif

       timel = time



       ntot = nx1*ny1*nz1*nelv

c     Vorticity calculations.

       if (mod(istep,iostep).eq.0) then
           call comp_vort3(vort,w1,w2,vx,vy,vz)
c         Calculations done in polar coordinates.
           do i=1,ntot
               rayon = SQRT(xm1(i,1,1,1)**2 + ym1(i,1,1,1)**2)
               if (rayon>0.0) then
                   vort_pol(i,1)=(xm1(i,1,1,1)*vort(i,1) +
      $                ym1(i,1,1,1)*vort(i,2))/rayon
                   vort_pol(i,2)=(-ym1(i,1,1,1)*vort(i,1) +
      $                xm1(i,1,1,1)*vort(i,2))/rayon
                   vort_pol(i,3)=vort(i,3)
               else
                   vort_pol(i,1)=0.0
                   vort_pol(i,2)=0.0
                   vort_pol(i,3)=vort(i,3)
               endif
           enddo
c          call outpost(vort(1,1),vort(1,2),vort(1,3),
c     $         pr,t(1,1,1,1,1),'vrt')
           call outpost(vort_pol(1,1),vort_pol(1,2),vort_pol(1,3),
      $         pr,t(1,1,1,1,1),'pol')
       endif

       if (istep.eq.1) then

       endif

c     We export the grid spacing and the average values of velocities.
       if(istep.eq.18) then
           ifxyo = .true.
           call outpost(dx2(1,1,1,1,1),dy2(1,1,1,1,1),dz2(1,1,1,1,1),
      $         ediff(1,1,1,1),t(1,1,1,1,1),'les')
           call outpost(uravg(1,1,1,1),utavg(1,1,1,1),uzavg(1,1,1,1),
      $         pr,t(1,1,1,1,1),'avg')
           if (nid.eq.0) then
           open(223,file='diffusion.dat',status='old',action='write',
      $      form='formatted',position="append")
           do i=1,ntot
c              write(223,*) xm1(i,1,1,1),ym1(i,1,1,1),zm1(i,1,1,1),nid,
c     $ ediff(i,1,1,1),dx2(i,1,1,1,nid+1),dy2(i,1,1,1,nid+1),
c     $            dz2(i,1,1,1,nid+1)
               write(223,*) nid
           enddo
           close(223)
           endif
           if (nid.eq.1) then
open(225,file='diffusion1.dat',status='old',action='write',
      $      form='formatted',position="append")
           do i=1,ntot
c              write(223,*) xm1(i,1,1,1),ym1(i,1,1,1),zm1(i,1,1,1),nid,
c     $ ediff(i,1,1,1),dx2(i,1,1,1,nid+1),dy2(i,1,1,1,nid+1),
c     $            dz2(i,1,1,1,nid+1)
               write(225,*) nid
           enddo
           close(225)
           endif

       endif
       return
       end
c-----------------------------------------------------------------------
       subroutine userbc (ix,iy,iz,iside,ieg)
c     NOTE ::: This subroutine MAY NOT be called by every process
       include 'SIZE'
       include 'TOTAL'
       include 'NEKUSE'
c
       ux=0.0
       uy=0.0
       uz=0.0
       temp=0
c
       z0 = 0.05
       rayon = SQRT(x**2 + y**2) ! Radius in polar coordinates
       R = 1.0                   ! Radius of the jet
       theta0 = 0.05*R           ! Thickness of the mixing layer
       axialamp = 0.05        ! Turbulence amplitude (axial excitation)
       helicamp = 0.05        ! Turbulence amplitude (helical direction)
       faxial = 0.55             ! Axial frequency = St*U0/D et D=1.0
       fhelic = faxial/2.0       ! Helical frequency = faxial/2.0
       phihelic = 0.0            ! Phase lag between axial and helical
c
       if (z.le.z0) then
           ux=0.0
           uy=0.0
           if (rayon.gt.(0.0)) then
               uz=0.5*(1-tanh(0.25*(rayon/R - R/rayon)*R/theta0))
           else
               uz=1.0
           endif
       endif
c
       return
       end
c-----------------------------------------------------------------------
       subroutine useric (ix,iy,iz,ieg)
       include 'SIZE'
       include 'TOTAL'
       include 'NEKUSE'

       rayon = SQRT(x**2 + y**2)  ! Radius in polar coordinates
       R = 0.5                    ! Radius of the jet
       theta0 = 0.05*R            ! Thickness of the mixing layer

       ux=0.0
       uy=0.0
       if (rayon.gt.(0.0)) then
           uz=0.5*(1-tanh(0.25*(rayon/R - R/rayon)*R/theta0))
       else
           uz=1.0
       endif
       temp=0
       return
       end
c-----------------------------------------------------------------------
       subroutine usrdat
       include 'SIZE'
       include 'TOTAL'
       common /cdsmag/ ediff(lx1,ly1,lz1,lelv)

       n = nx1*ny1*nz1*nelt
       call cfill(ediff,param(2),n)  ! initialize viscosity

       ! enable stress formulation if we use PN-PN-2
       if(lx1.ne.lx2 .and. param(30).gt.0) IFSTRS = .true.

       return
       end
c-----------------------------------------------------------------------
       subroutine usrdat2
       include 'SIZE'
       include 'TOTAL'


       return
       end
c-----------------------------------------------------------------------
       subroutine usrdat3
       include 'SIZE'
       include 'TOTAL'
c
       return
       end
c-----------------------------------------------------------------------
       subroutine set_obj  ! define objects for surface integrals
c
       include 'SIZE'
       include 'TOTAL'

       integer e,f,eg

       nobj = 1
       iobj = 0
       do ii=nhis+1,nhis+nobj
          iobj = iobj+1
          hcode(10,ii) = 'I'
          hcode( 1,ii) = 'F'
          hcode( 2,ii) = 'F'
          hcode( 3,ii) = 'F'
          lochis(1,ii) = iobj
       enddo
       nhis = nhis + nobj

       if (maxobj.lt.nobj) call exitti('increase maxobj in SIZE$',nobj)
       nxyz  = nx1*ny1*nz1
       nface = 2*ndim

       do e=1,nelv
       do f=1,nface
          if (cbc(f,e,1).eq.'W  ') then
             iobj  = 1
             if (iobj.gt.0) then
                nmember(iobj) = nmember(iobj) + 1
                mem = nmember(iobj)
                eg  = lglel(e)
                object(iobj,mem,1) = eg
                object(iobj,mem,2) = f
c              write(6,1) iobj,mem,f,eg,e,nid,' OBJ'
c   1          format(6i9,a4)

             endif
          endif
       enddo
       enddo
c     write(6,*) 'number',(nmember(k),k=1,4)
c
       return
       end
c-----------------------------------------------------------------------
       subroutine eddy_visc(ediff,e)
c
c     Compute eddy viscosity using Vreman's model.
c
       include 'SIZE'
       include 'TOTAL'
       include 'ZPER'

       real ediff(nx1*ny1*nz1,nelv)
       real betaij(lx1*ly1*lz1,ldim,ldim)
       real Bbeta(lx1*ly1*lz1)
       real Sommeij2(lx1*ly1*lz1)
       integer e,i,ldim1,ldim2

       common /lesleo/ sij (lx1*ly1*lz1,ldim,ldim)
      $              , dx2 (lx1*ly1*lz1,lelv,lp)
      $              , dy2 (lx1*ly1*lz1,lelv,lp)
      $              , dz2 (lx1*ly1*lz1,lelv,lp)
      $              , cdiff
       real sij,dx2,dy2,dz2,cdiff
       integer ntot

       ntot = nx1*ny1*nz1

c     Velocity gradient.
       call comp_gije(sij,vx(1,1,1,e),vy(1,1,1,e),vz(1,1,1,e),e)

c     Calculation of Beta.
       do ldim2=1,ldim
         do ldim1=1,ldim
           do i=1,ntot
             betaij(i,ldim1,ldim2) =
      $        dx2(i,e,nid+1)*sij(i,1,ldim1)*sij(i,1,ldim2)
      $        + dy2(i,e,nid+1)*sij(i,2,ldim1)*sij(i,2,ldim2)
      $        + dz2(i,e,nid+1)*sij(i,3,ldim1)*sij(i,3,ldim2)
           enddo
         enddo
       enddo


c     Calculation of Bbeta.
       do i=1,ntot
         Bbeta(i) = betaij(i,1,1)*betaij(i,2,2) - (betaij(i,1,2)**2)
      $           + betaij(i,2,2)*betaij(i,3,3) - (betaij(i,2,3)**2)
      $           + betaij(i,3,3)*betaij(i,1,1) - (betaij(i,3,1)**2)
       enddo

c     Calculons Sommeij2 = Sum_{i,j} (alpha_{ij}^2) pour le modèle de
c     Vreman. On met d'abord Sommeij2 à zero.
       call rzero(Sommeij2,ntot)
       do ldim2=1,ldim
         do ldim1=1,ldim
           do i=1,ntot
             Sommeij2(i) = Sommeij2(i) + (sij(i,ldim1,ldim2)**2)
           enddo
         enddo
       enddo

c     ediff calculation.
       do i=1,ntot
         if (Sommeij2(i).gt.0.0) then
           ediff(i,e) = param(2) + cdiff*sqrt(Bbeta(i)/Sommeij2(i))
c          ediff(i,e) = param(2)
         else
           ediff(i,e) = param(2)
         endif
       enddo

       return
       end
c-----------------------------------------------------------------------
       subroutine set_grid_spacing(dx2,dy2,dz2)
c
c     Compute D^2, the grid spacing used in the DS sgs model.
c
       include 'SIZE'
       include 'TOTAL'


       real dx2(nx1,ny1,nz1,nelv,np)
       real dy2(nx1,ny1,nz1,nelv,np)
       real dz2(nx1,ny1,nz1,nelv,np)
       real longueur,largeur,dhlong,dhlarg,deltai,deltaj

       integer e,i,j,k,im,jm,km,ip,jp,kp

c      return               ! Comment this line for a non-trivial Delta defn

       write(6,*) 'LadHyX'

       open (224,file='dx2dy2dz2.dat',status='old',action='write',
      $      form='formatted',position="append")
       do e=1,nelv
          do k=1,nz1
            km = max(1  ,k-1)
            kp = min(nz1,k+1)
            do j=1,ny1
              jm = max(1  ,j-1)
              jp = min(ny1,j+1)
              do i=1,nx1
                im = max(1  ,i-1)
                ip = min(nx1,i+1)
                longueur = 0.5*(abs(xm1(ip,jp,k,e) - xm1(im,jp,k,e))
      $           + abs(xm1(ip,jm,k,e) - xm1(im,jm,k,e)))/(ip - im)
                largeur = 0.5*(abs(ym1(ip,jp,k,e) - ym1(ip,jm,k,e))
      $           + abs(ym1(im,jp,k,e) - ym1(im,jm,k,e)))/(jp - jm)
                dhlong = 0.5*(abs(ym1(ip,jp,k,e) - ym1(im,jp,k,e))
      $           + abs(ym1(ip,jm,k,e) - ym1(im,jm,k,e)))/(ip - im)
                dhlarg = 0.5*(abs(xm1(ip,jp,k,e) - xm1(ip,jm,k,e))
      $           + abs(xm1(im,jp,k,e) - xm1(im,jm,k,e)))/(jp - jm)
                deltai = sqrt(longueur**2 + dhlong**2)
                deltaj = sqrt(largeur**2 + dhlarg**2)

                dx2(i,j,k,e,nid+1) = ((longueur**2)/deltai + (dhlarg**2)/
      $                              deltaj)**2
                dy2(i,j,k,e,nid+1) = ((dhlong**2)/deltai + (largeur**2)/
      $                              deltaj)**2
                dz2(i,j,k,e,nid+1) = ((zm1(i,j,kp,e)-zm1(i,j,km,e))/
      $                              (kp-km))**2
                write(224,*) i,j,k,e,nid,dx2(i,j,k,e,nid+1),
      $                      dy2(i,j,k,e,nid+1),dz2(i,j,k,e,nid+1)
              enddo
            enddo
          enddo
       enddo
       close(224)

       return
       end
c-----------------------------------------------------------------------
c
c automatically added by makenek
       subroutine usrsetvert(glo_num,nel,nx,ny,nz) ! to modify glo_num
       integer*8 glo_num(1)
       return
       end



-------------------------------------------
log file

/----------------------------------------------------------\\
|      _   __ ______ __ __  ______  ____    ____    ____   |
|     / | / // ____// //_/ / ____/ / __ \\ / __ \\ / __ \\ |
|    /  |/ // __/  / ,<   /___ \\  / / / // / / // / / /   |
|   / /|  // /___ / /| | ____/ / / /_/ // /_/ // /_/ /     |
|  /_/ |_//_____//_/ |_|/_____/  \\____/ \\____/ \\____/   |
|                                                          |
|----------------------------------------------------------|
|                                                          |
| NEK5000:  Open Source Spectral Element Solver            |
| COPYRIGHT (c) 2008-2010 UCHICAGO ARGONNE, LLC            |
| Version:  1.0rc1 / SVN  r1088                            |
| Web:      http://nek5000.mcs.anl.gov                     |
|                                                          |
\\----------------------------------------------------------/


  Number of processors:          16
  REAL    wdsize      :           8
  INTEGER wdsize      :           4


   Beginning session:
/home/leopold/NEK5000/turbessai_2511_1714/jet3d.rea


  timer accuracy:   0.0000000E+00 sec

  read .rea file
  nelgt/nelgv/lelt:        1200        1200          75
  lx1  /lx2  /lx3 :           8           6           6

  mapping elements to processors
            1          75          75        1200        1200  NELV
            2          75          75        1200        1200  NELV
            3          75          75        1200        1200  NELV
            4          75          75        1200        1200  NELV
            5          75          75        1200        1200  NELV
            6          75          75        1200        1200  NELV
            7          75          75        1200        1200  NELV
            8          75          75        1200        1200  NELV
            9          75          75        1200        1200  NELV
           10          75          75        1200        1200  NELV
           11          75          75        1200        1200  NELV
           12          75          75        1200        1200  NELV
           13          75          75        1200        1200  NELV
           14          75          75        1200        1200  NELV
           15          75          75        1200        1200  NELV
            0          75          75        1200        1200  NELV
  RANK     0 IEG     851     861     871     951     952     953 954     961
                     962     963     964     965     971     972 973     974
                     975     981     982     983     984     985 991     992
                     993     994     995    1051    1052    1053 1054    
1061
                    1062    1063    1064    1065    1071    1072 1073    
1074
                    1075    1081    1082    1083    1084    1085 1091    
1092
                    1093    1094    1095    1151    1152    1153 1154    
1161
                    1162    1163    1164    1165    1171    1172 1173    
1174
                    1175    1181    1182    1183    1184    1185 1191    
1192
                    1193    1194    1195
  element load imbalance:            0          75          75
  done :: mapping elements to processors

            0  objects found
  done :: read .rea file   0.32928E-01 sec

  setup mesh topology
    Right-handed check complete for    1200 elements. OK.
    setvert3d:   8      169285      428485      169285      169285
  call usrsetvert
  done :: usrsetvert

gs_setup: 30185 unique labels shared
    pairwise times (avg, min, max): 8.32036e-05 8.0204e-05 8.54015e-05
    crystal router                : 0.00022305 0.000220394 0.000225186
    all reduce                    : 0.000695395 0.000693202 0.000698304
    used all_to_all method: pairwise
    handle bytes (avg, min, max): 297659 277636 321556
    buffer bytes (avg, min, max): 65870 50368 84512
    setupds time 1.1183E-01 seconds   0  8      169285        1200
            8  max multiplicity
  done :: setup mesh topology

  call usrdat
  done :: usrdat

  generate geometry data
  vol_t,vol_v:   76.800000000007586        76.800000000007586
  done :: generate geometry data

  call usrdat2
  done :: usrdat2

  regenerate geometry data           1
  vol_t,vol_v:   76.800000000007373        76.800000000007373
  done :: regenerate geometry data           1

  verify mesh topology
   -2.0000000000000018        2.0000000000000018       Xrange
   -2.0000000000000018        2.0000000000000018       Yrange
    0.0000000000000000        4.8000000000000069       Zrange
  done :: verify mesh topology

  103   Parameters from 
file:/home/leopold/NEK5000/turbessai_2511_1714/jet3d.rea
    1          1.000000         p01 DENSITY
    2         -600.000         p02 VISCOS
    7          1.000000            p07 RHOCP
    8         -1492.000         p08 CONDUCT
   11         300.00           p11 NSTEPS
   12         4.0000000E-02     p12 DT
   15         25.000         p15 IOSTEP
   21         1.0000000E-07     p21 DIVERGENCE
   22         8.0000000E-09     p22 HELMHOLTZ
   24         5.0000000E-03     p24 TOLREL
   25         1.0000000E-02     p25 TOLABS
   26          2.000000         p26 COURANT/NTAU
   27          2.000000         p27 TORDER
   28         0.0000000E+00     p28 TORDER: mesh velocity (0: p28=p27)
   30         1.0000000E+00     p30 > 0 ==> properties set in uservp()
   65          1.000000         p65 #iofile(eg,0 or 64); <0 --> sep. dirs
   66          6.000000         p66 write fmt:ONLY postx uses rea value
   67          6.000000         p67 read fmt: same modes as p66
   68          1000.000         p68 iastep: freq for avg_all
   93          30.00000         p93 Numbr of prev pressure solns saved
   99          3.000000         p99    dealiasing:if <0 disable
  103         1.0000000E-03     p103   weight of stabilizing filter (.01)

  IFTRAN   = T
  IFFLOW   = T
  IFHEAT   = F
  IFSPLIT  = F
  IFLOMACH = F
  IFUSERVP = T
  IFUSERMV = F
  IFSTRS   = T
  IFCHAR   = T
  IFCYCLIC = F
  IFAXIS   = F
  IFMVBD   = F
  IFMELT   = F
  IFMODEL  = F
  IFKEPS   = F
  IFMOAB   = F
  IFNEKNEK = F
  IFSYNC   = T

  IFVCOR   = F
  IFINTQ   = F
  IFCWUZ   = F
  IFSWALL  = F
  IFGEOM   = F
  IFSURT   = F
  IFWCNO   = F
  IFCMT    = F
  IFVISC   = F
  IFFLTR   = F

  IFTMSH for field           1    =  F
  IFADVC for field           1    =  T
  IFNONL for field           1    =  F

  Dealiasing enabled, lxd=          12

  Estimated eigenvalues
  EIGAA =   0.55402293223707511
  EIGGA =    43231.464629758724
  EIGAE =   0.42836824657505768
  EIGAS =    6.4432989690721601E-002
  EIGGE =    43231.464629758724
  EIGGS =    2.0000000000000000

  verify mesh topology
   -2.0000000000000018        2.0000000000000018       Xrange
   -2.0000000000000018        2.0000000000000018       Yrange
    0.0000000000000000        4.8000000000000069       Zrange
  done :: verify mesh topology

   E-solver strategy:  1 itr
  mg_nx:           1           3           7
  mg_ny:           1           3           7
  mg_nz:           1           3           7
  call usrsetvert
  done :: usrsetvert

gs_setup: 641 unique labels shared
    pairwise times (avg, min, max): 5.16847e-05 5.04971e-05 5.30005e-05
    crystal router                : 2.23234e-05 2.19107e-05 2.26021e-05
    all reduce                    : 4.48287e-05 4.42982e-05 4.52042e-05
    used all_to_all method: crystal router
    handle bytes (avg, min, max): 22825 20364 25796
    buffer bytes (avg, min, max): 3360 3040 3744
    setupds time 3.4559E-03 seconds   1  2        1573        1200
    setvert3d:   4       25957       35557       25957       25957
  call usrsetvert
  done :: usrsetvert

gs_setup: 5593 unique labels shared
    pairwise times (avg, min, max): 4.74274e-05 4.41074e-05 4.97103e-05
    crystal router                : 5.89028e-05 5.86033e-05 5.95093e-05
    all reduce                    : 0.000171191 0.000170612 0.000171995
    used all_to_all method: pairwise
    handle bytes (avg, min, max): 58603 54564 63444
    buffer bytes (avg, min, max): 13614 10304 17696
    setupds time 9.9120E-03 seconds   2  4       25957        1200
    setvert3d:   4       25957       35557       25957       25957
  call usrsetvert
  done :: usrsetvert

gs_setup: 5593 unique labels shared
    pairwise times (avg, min, max): 4.033e-05 3.91006e-05 4.20094e-05
    crystal router                : 5.31495e-05 5.26905e-05 5.39064e-05
    all reduce                    : 0.000171348 0.000170612 0.000172019
    used all_to_all method: pairwise
    handle bytes (avg, min, max): 58603 54564 63444
    buffer bytes (avg, min, max): 13614 10304 17696
    setupds time 9.2380E-03 seconds   3  4       25957        1200
    setvert3d:   6       81861      158661       81861       81861
  call usrsetvert
  done :: usrsetvert

gs_setup: 15441 unique labels shared
    pairwise times (avg, min, max): 6.1132e-05 6.02007e-05 6.27041e-05
    crystal router                : 0.00013117 0.00013051 0.000131893
    all reduce                    : 0.000425158 0.000423312 0.000426292
    used all_to_all method: pairwise
    handle bytes (avg, min, max): 154899 144404 167444
    buffer bytes (avg, min, max): 34846 26560 44896
    setupds time 2.1295E-02 seconds   4  6       81861        1200
  setup h1 coarse grid, nx_crs=           2
  call usrsetvert
  done :: usrsetvert

gs_setup: 641 unique labels shared
    pairwise times (avg, min, max): 5.51209e-05 5.09977e-05 5.77927e-05
    crystal router                : 1.62676e-05 1.59979e-05 1.65939e-05
    all reduce                    : 5.0886e-05 4.99964e-05 5.14984e-05
    used all_to_all method: crystal router
    handle bytes (avg, min, max): 22825 20364 25796
    buffer bytes (avg, min, max): 3360 3040 3744
  done :: setup h1 coarse grid    2.9180049896240234E-002  sec

  call usrdat3
  done :: usrdat3

  set initial conditions
  nekuic (1) for ifld            1
  call nekuic for vel
  xyz min    -2.0000      -2.0000       0.0000
  uvwpt min   0.0000       0.0000       0.0000       0.0000 0.0000
  PS min      0.0000       0.0000      0.99000E+22
  xyz max     2.0000       2.0000       4.8000
  uvwpt max   0.0000       0.0000       1.0000       0.0000 0.0000
  PS max      0.0000       0.0000     -0.99000E+22
  done :: set initial conditions

  call userchk
  schfile:/home/leopold/NEK5000/turbessai_2511_1714/jet3d.sch
  call outfld: ifpsco: F

         0  0.0000E+00 Write checkpoint:
        0       0 OPEN: poljet3d0.f00001

         0  0.0000E+00 done :: Write checkpoint
                               file size =      17.    MB
                               avg data-throughput =   520.9MB/s
                               io-nodes =     1

  done :: userchk

gridpoints unique/tot:        428485       614400
   dofs:               423444               259200

  Initialization successfully completed   0.99081     sec

Starting time loop ...

      DT/DTCFL/DTFS/DTINIT   0.400E-01   0.513E-01   0.000E+00 0.400E-01
Step      1, t= 4.0000000E-02, DT= 4.0000000E-02, C=  1.559 0.0000E+00 
0.0000E+00
  Solving for fluid F T T
   0  8.0000E-09  1.0338E-06  8.0000E-09 tol,matmod
            1 Helmholtz3/fluid:     11   0.5267E-08   0.8000E-08 0.1797E-01
     1 1.00000E-07 4.40911E-01 9.30604E-01 4.73790E-01       1 Divergence
     2 1.00000E-07 1.76488E-01 9.30604E-01 1.89649E-01       1 Divergence
     3 1.00000E-07 5.55944E-02 9.30604E-01 5.97401E-02       1 Divergence
     4 1.00000E-07 1.75339E-02 9.30604E-01 1.88415E-02       1 Divergence
     5 1.00000E-07 6.73386E-03 9.30604E-01 7.23600E-03       1 Divergence
     6 1.00000E-07 2.27044E-03 9.30604E-01 2.43975E-03       1 Divergence
     7 1.00000E-07 7.55678E-04 9.30604E-01 8.12029E-04       1 Divergence
     8 1.00000E-07 2.71533E-04 9.30604E-01 2.91781E-04       1 Divergence
     9 1.00000E-07 9.26116E-05 9.30604E-01 9.95177E-05       1 Divergence
    10 1.00000E-07 2.89531E-05 9.30604E-01 3.11121E-05       1 Divergence
    11 1.00000E-07 1.00088E-05 9.30604E-01 1.07552E-05       1 Divergence
    12 1.00000E-07 4.02671E-06 9.30604E-01 4.32698E-06       1 Divergence
    13 1.00000E-07 1.70162E-06 9.30604E-01 1.82852E-06       1 Divergence
    14 1.00000E-07 6.61159E-07 9.30604E-01 7.10462E-07       1 Divergence
    15 1.00000E-07 2.33259E-07 9.30604E-01 2.50654E-07       1 Divergence
    16 1.00000E-07 7.88484E-08 9.30604E-01 8.47282E-08       1 Divergence
          1 U-PRES gmres:     16  7.8848E-08  1.0000E-07  9.3060E-01 
7.7317E-02  3.0742E-01
           1   4.0000E-02  5.1049E-01 Fluid done
filt amp 0.0000 0.0000 0.0000 0.0000 0.0000 0.0000 0.0000 0.0010
filt trn 1.0000 1.0000 1.0000 1.0000 1.0000 1.0000 1.0000 0.9990
  LadHyX
      DT/DTCFL/DTFS/DTINIT   0.329E-01   0.329E-01   0.000E+00 0.400E-01
Step      2, t= 7.2887311E-02, DT= 3.2887311E-02, C=  1.676 1.1179E+00 
1.1179E+00
  Solving for fluid F T T
            2 Helmholtz3/fluid:     11   0.2403E-08   0.8000E-08 0.3008E+00
          2 U-PRES gmres:     16  6.7074E-08  1.0000E-07  6.6007E-01 
7.7229E-02  2.9921E-01
           2   7.2887E-02  5.0106E-01 Fluid done
  LadHyX
Step      3, t= 1.0577462E-01, DT= 3.2887311E-02, C=  1.686 2.5436E+00 
1.4257E+00
  Solving for fluid F T T
            3 Helmholtz3/fluid:     10   0.3690E-08   0.8000E-08 0.9894E-01
          3 U-PRES gmres:     15  6.4098E-08  1.0000E-07  2.0530E-01 
6.6209E-02  2.5669E-01
           3   1.0577E-01  4.3805E-01 Fluid done
  LadHyX
Step      4, t= 1.3866193E-01, DT= 3.2887311E-02, C=  1.693 3.8424E+00 
1.2988E+00
  Solving for fluid F T T
            4 Helmholtz3/fluid:      9   0.7802E-08   0.8000E-08 0.3080E-01
          4 U-PRES gmres:     12  8.7915E-08  1.0000E-07  4.0701E-02 
5.3055E-02  2.0200E-01
           4   1.3866E-01  3.6896E-01 Fluid done
  LadHyX
Step      5, t= 1.7154924E-01, DT= 3.2887311E-02, C=  1.699 5.0340E+00 
1.1916E+00
  Solving for fluid F T T
            5 Helmholtz3/fluid:      9   0.7287E-08   0.8000E-08 0.2859E-01
          5 U-PRES gmres:     12  6.0219E-08  1.0000E-07  2.5639E-02 
5.3179E-02  2.0134E-01
           5   1.7155E-01  3.6834E-01 Fluid done
  LadHyX
Step      6, t= 2.0443656E-01, DT= 3.2887311E-02, C=  1.695 6.2218E+00 
1.1877E+00
  Solving for fluid F T T
            6 Helmholtz3/fluid:      9   0.6680E-08   0.8000E-08 0.2745E-01
          6 U-PRES gmres:     12  6.6203E-08  1.0000E-07  3.2161E-02 
5.2715E-02  2.0061E-01
           6   2.0444E-01  3.6690E-01 Fluid done
  LadHyX
Step      7, t= 2.3732387E-01, DT= 3.2887311E-02, C=  1.697 7.4158E+00 
1.1941E+00
  Solving for fluid F T T
            7 Helmholtz3/fluid:      9   0.6353E-08   0.8000E-08 0.2680E-01
          7 U-PRES gmres:     12  5.1712E-08  1.0000E-07  2.6367E-02 
5.2906E-02  2.0242E-01
           7   2.3732E-01  3.6932E-01 Fluid done
  LadHyX
Step      8, t= 2.7021118E-01, DT= 3.2887311E-02, C=  1.698 8.6288E+00 
1.2129E+00
  Solving for fluid F T T
            8 Helmholtz3/fluid:      9   0.6673E-08   0.8000E-08 0.2685E-01
          8 U-PRES gmres:     12  4.4087E-08  1.0000E-07  2.3825E-02 
5.2858E-02  2.0136E-01
           8   2.7021E-01  3.6842E-01 Fluid done
  LadHyX
Step      9, t= 3.0309849E-01, DT= 3.2887311E-02, C=  1.684 9.8241E+00 
1.1954E+00
  Solving for fluid F T T
            9 Helmholtz3/fluid:      9   0.7127E-08   0.8000E-08 0.2677E-01
          9 U-PRES gmres:     12  4.5466E-08  1.0000E-07  2.2156E-02 
5.2914E-02  2.0164E-01
           9   3.0310E-01  3.6852E-01 Fluid done
  LadHyX
Step     10, t= 3.3598580E-01, DT= 3.2887311E-02, C=  1.645 1.1023E+01 
1.1992E+00
  Solving for fluid F T T
           10 Helmholtz3/fluid:      9   0.6816E-08   0.8000E-08 0.2582E-01
         10 U-PRES gmres:     12  4.4157E-08  1.0000E-07  1.8823E-02 
5.2885E-02  2.0112E-01
          10   3.3599E-01  3.6791E-01 Fluid done
  LadHyX
Step     11, t= 3.6887311E-01, DT= 3.2887311E-02, C=  1.632 1.2209E+01 
1.1855E+00
  Solving for fluid F T T
           11 Helmholtz3/fluid:      9   0.6294E-08   0.8000E-08 0.2484E-01
         11 U-PRES gmres:     12  4.1769E-08  1.0000E-07  1.6670E-02 
5.2812E-02  2.0193E-01
          11   3.6887E-01  3.6934E-01 Fluid done
  LadHyX
Step     12, t= 4.0176042E-01, DT= 3.2887311E-02, C=  1.641 1.3405E+01 
1.1956E+00
  Solving for fluid F T T
           12 Helmholtz3/fluid:      9   0.6230E-08   0.8000E-08 0.2475E-01
         12 U-PRES gmres:     12  4.2914E-08  1.0000E-07  1.6074E-02 
5.2975E-02  2.0160E-01
          12   4.0176E-01  3.6909E-01 Fluid done
  LadHyX
Step     13, t= 4.3464773E-01, DT= 3.2887311E-02, C=  1.651 1.4600E+01 
1.1958E+00
  Solving for fluid F T T
           13 Helmholtz3/fluid:      9   0.6215E-08   0.8000E-08 0.2488E-01
         13 U-PRES gmres:     12  3.9520E-08  1.0000E-07  1.4803E-02 
5.2862E-02  2.0123E-01
          13   4.3465E-01  3.6785E-01 Fluid done
  LadHyX
Step     14, t= 4.6753504E-01, DT= 3.2887311E-02, C=  1.661 1.5797E+01 
1.1967E+00
  Solving for fluid F T T
           14 Helmholtz3/fluid:      9   0.6096E-08   0.8000E-08 0.2437E-01
         14 U-PRES gmres:     12  3.6261E-08  1.0000E-07  1.3651E-02 
5.3142E-02  2.0131E-01
          14   4.6754E-01  3.6796E-01 Fluid done
  LadHyX
Step     15, t= 5.0042236E-01, DT= 3.2887311E-02, C=  1.671 1.6990E+01 
1.1930E+00
  Solving for fluid F T T
           15 Helmholtz3/fluid:      9   0.5814E-08   0.8000E-08 0.2304E-01
         15 U-PRES gmres:     12  3.8622E-08  1.0000E-07  1.2689E-02 
5.2549E-02  2.0080E-01
          15   5.0042E-01  3.6748E-01 Fluid done
  LadHyX
Step     16, t= 5.3330967E-01, DT= 3.2887311E-02, C=  1.686 1.8185E+01 
1.1949E+00
  Solving for fluid F T T
           16 Helmholtz3/fluid:      9   0.5392E-08   0.8000E-08 0.2159E-01
         16 U-PRES gmres:     12  4.0755E-08  1.0000E-07  1.1259E-02 
5.2697E-02  2.0102E-01
          16   5.3331E-01  3.6759E-01 Fluid done
  LadHyX
Step     17, t= 5.6619698E-01, DT= 3.2887311E-02, C=  1.706 1.9387E+01 
1.2019E+00
  Solving for fluid F T T
           17 Helmholtz3/fluid:      9   0.5216E-08   0.8000E-08 0.2097E-01
         17 U-PRES gmres:     12  3.9763E-08  1.0000E-07  1.0706E-02 
5.2939E-02  2.0170E-01
          17   5.6620E-01  3.6837E-01 Fluid done
  Calculating eddy visosity
Step     18, t= 5.9908429E-01, DT= 3.2887311E-02, C=  1.725 2.0349E+01 
9.6170E-01
  Solving for fluid F T T
           18 Helmholtz3/fluid:      9   0.5243E-08   0.8000E-08 0.2106E-01
         18 U-PRES gmres:     12  3.4908E-08  1.0000E-07  1.0308E-02 
5.3016E-02  2.0171E-01
          18   5.9908E-01  3.6860E-01 Fluid done
  Calculating eddy visosity
  call outfld: ifpsco: F

        18  5.9908E-01 Write checkpoint:
        0      18 OPEN: lesjet3d0.f00001

        18  5.9908E-01 done :: Write checkpoint
                               file size =      17.    MB
                               avg data-throughput =   617.1MB/s
                               io-nodes =     1

  call outfld: ifpsco: F

        18  5.9908E-01 Write checkpoint:
        0      18 OPEN: avgjet3d0.f00001

        18  5.9908E-01 done :: Write checkpoint
                               file size =      17.    MB
                               avg data-throughput =   647.7MB/s
                               io-nodes =     1

Step     19, t= 6.3197160E-01, DT= 3.2887311E-02, C=  1.739 2.1382E+01 
1.0337E+00
  Solving for fluid F T T
           19 Helmholtz3/fluid:     10   0.1805E-08   0.8000E-08 0.2278E-01
         19 U-PRES gmres:     11  9.1061E-08  1.0000E-07  9.2230E-03 
5.2888E-02  2.0412E-01
          19   6.3197E-01  4.0092E-01 Fluid done
  Calculating eddy visosity
Step     20, t= 6.6485891E-01, DT= 3.2887311E-02, C=  1.732 2.5370E+01 
3.9873E+00
  Solving for fluid F T T
           20 Helmholtz3/fluid:     10   0.1673E-08   0.8000E-08 0.2231E-01
         20 U-PRES gmres:     11  8.8236E-08  1.0000E-07  8.5941E-03 
5.2522E-02  2.0177E-01
          20   6.6486E-01  3.8955E-01 Fluid done
  Calculating eddy visosity
Step     21, t= 6.9774622E-01, DT= 3.2887311E-02, C=  1.700 2.6416E+01 
1.0463E+00
  Solving for fluid F T T
           21 Helmholtz3/fluid:     10   0.1618E-08   0.8000E-08 0.2169E-01
         21 U-PRES gmres:     11  8.8444E-08  1.0000E-07  8.0928E-03 
5.2019E-02  1.9286E-01
          21   6.9775E-01  3.8227E-01 Fluid done
  Calculating eddy visosity

......................................................

Step    245, t= 8.6498981E+00, DT= 3.9464773E-02, C=  1.878 2.4108E+02 
9.7130E-01
  Solving for fluid F T T
          245 Helmholtz3/fluid:     11   0.1876E-08   0.8000E-08 0.3329E-01
        245 U-PRES gmres:     12  4.0605E-08  1.0000E-07  1.2480E-02 
5.3032E-02  2.0306E-01
         245   8.6499E+00  3.9528E-01 Fluid done
  Calculating eddy visosity
Step    246, t= 8.6893628E+00, DT= 3.9464773E-02, C=  1.889 2.4207E+02 
9.8505E-01
  Solving for fluid F T T
          246 Helmholtz3/fluid:     11   0.1944E-08   0.8000E-08 0.3391E-01
        246 U-PRES gmres:     12  4.5293E-08  1.0000E-07  1.4596E-02 
5.3483E-02  2.0330E-01
         246   8.6894E+00  3.9605E-01 Fluid done
  Calculating eddy visosity
Step    247, t= 8.7288276E+00, DT= 3.9464773E-02, C=  1.898 2.4305E+02 
9.8502E-01
  Solving for fluid F T T
          247 Helmholtz3/fluid:     11   0.2644E-08   0.8000E-08 0.3535E-01
        247 U-PRES gmres:     12  8.9201E-08  1.0000E-07  2.1450E-02 
5.7194E-02  2.1392E-01
         247   8.7288E+00  4.1121E-01 Fluid done
  Calculating eddy visosity
Step    248, t= 8.7682924E+00, DT= 3.9464773E-02, C=  1.906 2.4411E+02 
1.0537E+00
  Solving for fluid F T T
          248 Helmholtz3/fluid:     11   0.4365E-08   0.8000E-08 0.4022E-01
        248 U-PRES gmres:     12  8.4866E-08  1.0000E-07  2.8833E-02 
5.6473E-02  2.1303E-01
         248   8.7683E+00  4.1053E-01 Fluid done
  Calculating eddy visosity
Step    249, t= 8.8077571E+00, DT= 3.9464773E-02, C=  2.072 2.4516E+02 
1.0495E+00
  Solving for fluid F T T
          249 Helmholtz3/fluid:     12   0.1785E-08   0.8000E-08 0.6760E-01
        249 U-PRES gmres:     13  8.3312E-08  1.0000E-07  1.0799E-01 
5.7963E-02  2.2072E-01
         249   8.8078E+00  4.2749E-01 Fluid done
  Calculating eddy visosity
Step    250, t= 8.8393290E+00, DT= 3.1571819E-02, C=  3.334 2.4622E+02 
1.0646E+00
  Solving for fluid F T T
          250 Helmholtz3/fluid:     13   0.3976E-08   0.8000E-08 0.3083E+00
        250 U-PRES gmres:     15  8.8405E-08  1.0000E-07  6.2205E-01 
6.7180E-02  2.5722E-01
         250   8.8393E+00  4.7414E-01 Fluid done
  Calculating eddy visosity
  call outfld: ifpsco: F

       250  8.8393E+00 Write checkpoint:
        0     250 OPEN: poljet3d0.f00011

       250  8.8393E+00 done :: Write checkpoint
                               file size =      9.4    MB
                               avg data-throughput =   774.5MB/s
                               io-nodes =     1

  call outfld: ifpsco: F

       250  8.8393E+00 Write checkpoint:
        0     250 OPEN: jet3d0.f00010

       250  8.8393E+00 done :: Write checkpoint
                               file size =      9.4    MB
                               avg data-throughput =   847.0MB/s
                               io-nodes =     1

Step    251, t= 8.8645864E+00, DT= 2.5257455E-02, C=  6.125 2.4732E+02 
1.1000E+00
  Solving for fluid F T T
          251 Helmholtz3/fluid:     20   0.3689E-08   0.8000E-08 0.4576E+02
        251 U-PRES gmres:     19  7.6441E-08  1.0000E-07  4.7293E+01 
8.9832E-02  3.4811E-01
         251   8.8646E+00  6.6288E-01 Fluid done
  Calculating eddy visosity
Step    252, t= 8.8847924E+00, DT= 2.0205964E-02, C=304.531 2.4858E+02 
1.2546E+00
  Solving for fluid F T T
    252   201 Unconverged Helmholtz3/Fluid: rbnorm = 0.140100E+18 
0.800000E-08
        252 U-PRES gmres:     76  9.1868E-08  1.0000E-07  2.8678E+21 
3.4409E-01  1.4581E+00
         252   8.8848E+00  4.0444E+00 Fluid done
  Calculating eddy visosity
Step    253, t= 8.9009572E+00, DT= 1.6164771E-02, C=******* 2.5326E+02 
4.6784E+00
  Solving for fluid F T T
    253   201 Unconverged Helmholtz3/Fluid: rbnorm =          NaN 
0.800000E-08
        253 U-PRES gmres:    120         NaN  1.0000E-07         NaN 
5.4340E-01  2.3310E+00
         253   8.9010E+00  4.8603E+00 Fluid done
  Calculating eddy visosity
Step    254, t= 8.9203549E+00, DT= 1.9397725E-02, C=******* 2.5870E+02 
5.4484E+00
  Solving for fluid F T T
    254   201 Unconverged Helmholtz3/Fluid: rbnorm =          NaN 
0.800000E-08
        254 U-PRES gmres:    120         NaN  1.0000E-07         NaN 
5.4476E-01  2.3302E+00
         254   8.9204E+00  4.8624E+00 Fluid done
  Calculating eddy visosity
Step    255, t= 8.9436322E+00, DT= 2.3277270E-02, C=******* 2.6416E+02 
5.4521E+00
  Solving for fluid F T T
    255   201 Unconverged Helmholtz3/Fluid: rbnorm =          NaN 
0.800000E-08
        255 U-PRES gmres:    120         NaN  1.0000E-07         NaN 
5.4309E-01  2.3275E+00
         255   8.9436E+00  4.8631E+00 Fluid done
  Calculating eddy visosity
Step    256, t= 8.9715649E+00, DT= 2.7932725E-02, C=******* 2.6961E+02 
5.4533E+00
  Solving for fluid F T T
    256   201 Unconverged Helmholtz3/Fluid: rbnorm =          NaN 
0.800000E-08
        256 U-PRES gmres:    120         NaN  1.0000E-07         NaN 
5.8407E-01  2.4309E+00
         256   8.9716E+00  4.9608E+00 Fluid done
  Calculating eddy visosity





More information about the Nek5000-users mailing list