[Nek5000-users] Vorticity

nek5000-users at lists.mcs.anl.gov nek5000-users at lists.mcs.anl.gov
Thu Jul 14 12:43:09 CDT 2011



Hi Alireza,

The best option is to call the vorticity/s_ij routines, 
then interpolate the fields that they generate.

Attached is a routine that should do the interpolation ---
I just hacked into another routine that was available and
haven't tested, but it should work.

Note that if you want to change your point list subsequent
to the initial call, you should pass a handle so indicating.

Interpolation requires two steps (done by the attached code):


    1. Find the points (x,y,z)_i  in the (r,s,t,e) parameter space

           (done once, and expensive)

    2. Evaulate quantities of interest at (r,s,t,e)_i

           (potentially done often, and less expensive)


Please let me know if you have any questions on usage, etc.

Paul




On Thu, 14 Jul 2011, nek5000-users at lists.mcs.anl.gov wrote:

> Hi,
>
> Thanks for your quick reply. Just I wanna know if I can use the interpolated 
> values of velocity at some arbitrary points in the domain (not necessarily 
> the grid points) as the input arguments to theses subroutines.
>
> Regards,
> Alireza
>
>
>
>
>
> Quoting nek5000-users-request at lists.mcs.anl.gov:
>
>> Send Nek5000-users mailing list submissions to
>> 	nek5000-users at lists.mcs.anl.gov
>> 
>> To subscribe or unsubscribe via the World Wide Web, visit
>> 	https://lists.mcs.anl.gov/mailman/listinfo/nek5000-users
>> or, via email, send a message with subject or body 'help' to
>> 	nek5000-users-request at lists.mcs.anl.gov
>> 
>> You can reach the person managing the list at
>> 	nek5000-users-owner at lists.mcs.anl.gov
>> 
>> When replying, please edit your Subject line so it is more specific
>> than "Re: Contents of Nek5000-users digest..."
>> 
>> 
>> Today's Topics:
>>
>>   1. Vorticity (nek5000-users at lists.mcs.anl.gov)
>>   2. Re: Vorticity (nek5000-users at lists.mcs.anl.gov)
>>   3. Re: Vorticity (nek5000-users at lists.mcs.anl.gov)
>> 
>> 
>> ----------------------------------------------------------------------
>> 
>> Message: 1
>> Date: Wed, 13 Jul 2011 17:55:35 -0400
>> From: nek5000-users at lists.mcs.anl.gov
>> Subject: [Nek5000-users] Vorticity
>> To: nek5000-users at lists.mcs.anl.gov
>> Message-ID: <20110713175535.1042430oct81vz2v at webmail.vt.edu>
>> Content-Type: text/plain; charset=ISO-8859-1; DelSp="Yes";
>> 	format="flowed"
>> 
>> Hi,
>> 
>> I need to calculate the vorticity vector and the rate of strain tensor
>> at any arbitrary point in the domain. I wonder if these values are
>> computed somewhere in the code or if I should interpolate them in the
>> usr file. I appreciate any kind of help.
>> 
>> Thanks in advance,
>> 
>> 
>> Regards,
>> Alireza Karimi
>> 
>> 
>> 
>> ------------------------------
>> 
>> Message: 2
>> Date: Wed, 13 Jul 2011 17:09:48 -0500 (CDT)
>> From: nek5000-users at lists.mcs.anl.gov
>> Subject: Re: [Nek5000-users] Vorticity
>> To: nek5000-users at lists.mcs.anl.gov
>> Message-ID: <Pine.LNX.4.64.1107131704140.26989 at jacobi.mcs.anl.gov>
>> Content-Type: TEXT/PLAIN; charset=US-ASCII; format=flowed
>> 
>> Hi Alireza,
>> 
>> For vorticity computation, you can call
>>
>>       call comp_vort3(vort,work1,work2,vx,vy,vz)
>> 
>> from userchk with the appropriate definitions of vorticity and work arrays
>> like
>>
>>       common /scren/   vort (lx1,ly1,lz1,lelv,3) ! x y z components
>>      $               , work1(lx1,ly1,lz1,lelv)
>>      $               , work2(lx1,ly1,lz1,lelv)
>> 
>> Let me check what we have for the calculation of strain tensor.
>> Best,
>> Aleks
>> 
>> 
>> On Wed, 13 Jul 2011, nek5000-users at lists.mcs.anl.gov wrote:
>> 
>>> Hi,
>>> 
>>> I need to calculate the vorticity vector and the rate of strain tensor at 
>>> any
>>> arbitrary point in the domain. I wonder if these values are computed
>>> somewhere in the code or if I should interpolate them in the usr file. I
>>> appreciate any kind of help.
>>> 
>>> Thanks in advance,
>>> 
>>> 
>>> Regards,
>>> Alireza Karimi
>>> 
>>> _______________________________________________
>>> Nek5000-users mailing list
>>> Nek5000-users at lists.mcs.anl.gov
>>> https://lists.mcs.anl.gov/mailman/listinfo/nek5000-users
>> 
>> 
>> ------------------------------
>> 
>> Message: 3
>> Date: Wed, 13 Jul 2011 23:59:08 -0500 (CDT)
>> From: nek5000-users at lists.mcs.anl.gov
>> Subject: Re: [Nek5000-users] Vorticity
>> To: nek5000-users at lists.mcs.anl.gov
>> Message-ID:
>> 	<600285192.94128.1310619548370.JavaMail.root at zimbra.anl.gov>
>> Content-Type: text/plain; charset=utf-8
>> 
>> Hi Alireza,
>> 
>> 
>> For a strain stress computation you can call subroutine comp_sij
>>
>>       call comp_sij(sij,nij,vx,vy,vz
>>      $                     ,ur,us,ut,vr,vs,vt,wr,ws,wt)
>> 
>> 
>> (navier5.f:2225) with nij=6 for 3D somewhere from userchk once you declare 
>> the stress component array sij and arrays with local coordinate r-s-t 
>> derivatives of velocity components
>>
>>       common /csij/
>>      $     sij(lx1*ly1*lz1,6,lelv) ! 6 components
>>      $    ,ur (lx1*ly1*lz1)
>> ...
>>      $    ,wt (lx1*ly1*lz1)
>> 
>> 
>> Best,
>> Aleks
>> 
>> 
>> 
>> 
>> 
>> 
>> ----- Original Message -----
>> From: nek5000-users at lists.mcs.anl.gov
>> To: nek5000-users at lists.mcs.anl.gov
>> Sent: Wednesday, July 13, 2011 5:09:48 PM
>> Subject: Re: [Nek5000-users] Vorticity
>> 
>> Hi Alireza,
>> 
>> For vorticity computation, you can call
>>
>>       call comp_vort3(vort,work1,work2,vx,vy,vz)
>> 
>> from userchk with the appropriate definitions of vorticity and work arrays
>> like
>>
>>       common /scren/   vort (lx1,ly1,lz1,lelv,3) ! x y z components
>>      $               , work1(lx1,ly1,lz1,lelv)
>>      $               , work2(lx1,ly1,lz1,lelv)
>> 
>> Let me check what we have for the calculation of strain tensor.
>> Best,
>> Aleks
>> 
>> 
>> On Wed, 13 Jul 2011, nek5000-users at lists.mcs.anl.gov wrote:
>> 
>>> Hi,
>>> 
>>> I need to calculate the vorticity vector and the rate of strain tensor at 
>>> any
>>> arbitrary point in the domain. I wonder if these values are computed
>>> somewhere in the code or if I should interpolate them in the usr file. I
>>> appreciate any kind of help.
>>> 
>>> Thanks in advance,
>>> 
>>> 
>>> Regards,
>>> Alireza Karimi
>>> 
>>> _______________________________________________
>>> Nek5000-users mailing list
>>> Nek5000-users at lists.mcs.anl.gov
>>> https://lists.mcs.anl.gov/mailman/listinfo/nek5000-users
>> _______________________________________________
>> Nek5000-users mailing list
>> Nek5000-users at lists.mcs.anl.gov
>> https://lists.mcs.anl.gov/mailman/listinfo/nek5000-users
>> 
>> 
>> ------------------------------
>> 
>> _______________________________________________
>> Nek5000-users mailing list
>> Nek5000-users at lists.mcs.anl.gov
>> https://lists.mcs.anl.gov/mailman/listinfo/nek5000-users
>> 
>> 
>> End of Nek5000-users Digest, Vol 29, Issue 2
>> ********************************************
>> 
>
>
>
> -- 
> Alireza Karimi
> PhD Candidate
> Department of Engineering Science and Mechanics
> Virginia Polytechnic Institute and State University
> Blacksburg, VA 24061
> _______________________________________________
> Nek5000-users mailing list
> Nek5000-users at lists.mcs.anl.gov
> https://lists.mcs.anl.gov/mailman/listinfo/nek5000-users
-------------- next part --------------
C-----------------------------------------------------------------------
      subroutine interp_v3(uvw,xyz,v1,v2,v3,n)
c
c     evaluate velocity for list of points xyz
c
c     Note:  -- modify
c     intpts to get rid off " WARNING: point on boundary or ..."

      include 'SIZE'
      include 'TOTAL'

      real uvw(ldim,n),xyz(ldim,n),v1(1),v2(1),v3(1)
      logical ifjac,ifpts

      parameter(nmax=lpart,nfldmax=ldim) 
      common /rv_intp/ pts(ldim*nmax)
      common /iv_intp/ ihandle
      common /outtmp/ wrk(lx1*ly1*lz1*lelt,nfldmax)

      integer icalld,e
      save    icalld
      data    icalld /0/

      nxyz  = nx1*ny1*nz1
      ntot  = nxyz*nelt

      if (n.gt.nmax) call exitti ('ABORT: interp_v() n > nmax!$',n)
      
      if (nelgt.ne.nelgv) call exitti
     $   ('ABORT: interp_v() nelgt.ne.nelgv not yet supported!$',nelgv)

      do i=1,n                          ! ? not moving -> save?
         pts(i)     = xyz(1,i)
         pts(i + n) = xyz(2,i)
         if (if3d)  pts(i + n*2) = xyz(3,i)
      enddo

      if (icalld.eq.0) then             ! interpolation setu !? intpts_done(ih_intp_v)?
        icalld = 1
        tolin  = 1.e-8
        call intpts_setup(tolin,ihandle)
      endif

      nflds  = ndim ! number of fields to interpolate

      ! pack working array
      call opcopy(wrk(1,1),wrk(1,2),wrk(1,3),v1,v2,v3)

      ! interpolate
      ifjac  = .true.           ! output transpose (of Jacobian)
      ifpts  = .true.            ! find points
      call intpts(wrk,nflds,pts,n,uvw,ifjac,ifpts,ihandle) ! copy array instead?

      return
      end
c-----------------------------------------------------------------------


More information about the Nek5000-users mailing list