[Nek5000-users] Temperature field values in mesh 1 only

nek5000-users at lists.mcs.anl.gov nek5000-users at lists.mcs.anl.gov
Thu May 19 06:09:28 CDT 2016


Hi Padyumna,

I just checked that gradm1() operates on all elements (1:nelt)

If you want a field that is zero in the solid elements, here is one approach:


      real tm1(lx1,ly1,lz1,lelt,ldimt)   [ note, use lelt here ]

      nv = nx1*ny1*nz1*nelv
      nt = nx1*ny1*nz1*nelt

      do i=1,ldimt
          call rzero(tm1(1,1,1,1,i),nt)
          call copy (tm1(1,1,1,1,i),t(1,1,1,1,i),nv)
      enddo

Note that if you're trying to compute (say)   grad V_x * T, you can do so via:

      call col2(vxt,vx,t,nt)
      call gradm1(vxtx,vxty,vxtz,vxt)

Note that vxt will be zero in the solid elements because vx is zero there.

Best, Paul

________________________________
From: nek5000-users-bounces at lists.mcs.anl.gov [nek5000-users-bounces at lists.mcs.anl.gov] on behalf of nek5000-users at lists.mcs.anl.gov [nek5000-users at lists.mcs.anl.gov]
Sent: Thursday, May 19, 2016 5:27 AM
To: nek5000-users at lists.mcs.anl.gov
Subject: [Nek5000-users] Temperature field values in mesh 1 only

Hi Neks!

I am running a Conjugate Heat Transfer simulation with a fluid box (velocity mesh) and two solid boxes (all boxes- temperature mesh). I am interested in computing the heat flux budget for this simulation. While doing this, on multiple occasions, I need to multiply the temperature or passive scalars throughout the fluid domain field with the velocity component values at the respective locations. I also need to calculate the gradients of these variables throughout the fluid domain. For this I use the sub-routine “gradm1” present in navier5.f . But I believe gradm1 calculates the gradients for field variables in the velocity mesh only. For all these functions, I think it’s best if I could isolate an array with the values of the temperature and multiple passive scalars in the fluid domain (or velocity mesh) only, which should be a subset of the whole temperature array.

Thus my question is, is there any way to get this? I am using multiple processors which complicates things. I want to create an array that contains the values of the temperature and the passive scalars at all the velocity mesh elements and this I have to do across all processors. In other words, from the standard NEK array t(lx2,ly2,lz2,lelt,ldimt) which contains the temperature and the passive scalars, I have to create an array like, say, tm1(lx1,ly1,lz1,lelv,ldimt) which would have these values. I am not sure if I have been able to explain my query clearly enough, but I would greatly appreciate if anybody could help me with this or if you could think of an easier way to calculate the gradient of the temperature and passive scalars throughout the velocity field. Let me know if you want more clarity with the question, I will try to explain better.

Thanks in advance and cheers!

Pradyumna


-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/nek5000-users/attachments/20160519/19a7e3a7/attachment.html>


More information about the Nek5000-users mailing list