[Nek5000-users] Volume of an element
nek5000-users at lists.mcs.anl.gov
nek5000-users at lists.mcs.anl.gov
Tue Feb 27 06:53:52 CST 2018
Sorry for an empty email before, it was a mistake. Here is my problem I was trying to make the blending coefficient(h) depending upon the volume of elements to have a mesh deoformation rigid near the walls (beasev = 1) and 0 at boundaries (beasev = 0). I had also seen a way to make blending coefficient from the papper that you had sent to me and it was there in an example called ocyl2.usr . But there is variable called delta which is an unknown scale to me and it depends up on the mesh.
How ever, I have also seen another equation which deals with volumes of elements in mesh (you had mentioned in you papper also). I was successfully implemented that function on 2D cases by the code shown at the end. But when it becomes to 3D case, I realise that the volume of elements that i tried to capture is in a wrong way. So could you please help me to find the volume of element corresponding to all the grid points (I would like to have the element volume not the spectral elements volume for each grid points) for example if the spectral element order is 4. then it would have 25 small spectral elements in 1 big element right? What I would like to have what ever the grid point in that element I would like to have the element volume.
that is if i tried like below,,
n = nx1*ny1*nz1*nelv
do e=1,n
volume = bm1(e,1,1,1)
enddo
it will throw spectral elements volume right?
and if i tried this,
do e=1,nelt
vol(e) = volel(e)
enddo
I will get element volume,
but i would like to create smoothing varibale for all grid points. That is i need h(n) not h(nelt) (n and nelt are the array sizes)
Code which works for 2D
v = volume(1,1,1,1)
a=glmax(v,1)
b=glmin(v,1)
do i=1,nelv
denom = 1/a
h1(i) = ((1-(b/a))/(denom*v))
enddo
Regards
Sijo
De: "nek5000-users" <nek5000-users at lists.mcs.anl.gov>
À: "nek5000-users" <nek5000-users at lists.mcs.anl.gov>
Envoyé: Mardi 27 Février 2018 13:28:54
Objet: Re: [Nek5000-users] Volume of an element
Dear Sijo,
The volume per grid point is given by the entries in bm1(i,j,k,e)
For each element e, there are points [i,j,k] with index ranges i=1,...,nx1
j=1,...,ny1, k=1,...,nz1
You must reference elements by the 4th index (...,e).
Fortunately, all the data is contiguous, so, if accessing _all_ points,
you can use a loop like
do i=1,n
volume_of_point_i = bm1(i,1,1,1)
enddo
You certainly do not need to do any summing if you only want the
point wise volume.
Perhaps it'd be best if you indicate what exactly you are planning to
do with the volume.
Finally, I see your use of glmax and glmin. How many processors
are you using? If more than one, are you intending that they are
all writing to the file hcoeff.dat?
hth,
Paul
From: Nek5000-users <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: Tuesday, February 27, 2018 4:45:30 AM
To: nek5000-users
Subject: Re: [Nek5000-users] Volume of an element
Thanks Paul for your Kind replay. You were right. I want to get the volume associated with each grid point. And I tried the code below, but at each time when I run I get different answers for volume so I think my ocde is not rubust. Could you please verify my code below.
nxyz = nx1*ny1*nz1
n = nx1*ny1*nz1*nelv
vol = vlsum(bm1(1,1,1,1),nxyz)
a=glmax(vol,1)
b=glmin(vol,1)
do e=1,n
volume = vlsum(bm1(e,1,1,1),nxyz)
open (unit=10,file="hcoeff.dat")
write (10,7020) a,b,volume
7020 format(7f20.14)
enddo
Regards
Sijo GEORGE
De: "nek5000-users" <nek5000-users at lists.mcs.anl.gov>
À: "nek5000-users" <nek5000-users at lists.mcs.anl.gov>
Envoyé: Lundi 26 Février 2018 21:46:02
Objet: Re: [Nek5000-users] Volume of an element
Dear Sijo,
We typically refer to elements as the structures that have (N+1)^3
points. (Here, N+1 = lx1.)
They are enumerated in nek5000 by the element number, e:
integer e
nxyz = nx1*ny1*nz1
do e=1,nelt
volume = volel(e)
volume = vlsum(bm1(1,1,1,e),nxyz)
enddo
would both yield the same value for the volume of element e.
It appears to me that you are after something different, however.
It seems like you want to associate a volume with each grid point.
If that is the case, then the entries of the diagonal mass matrix, bm1,
are indeed the quantities that you want.
Each entry of bm1(i,1,1,1) for i=1 to n:=nxyz*nelt, corresponds to the
integral of the underlying basis function.
That is, if we think about a function phi_i(x) such that
u(x) = sum_i=1^n u(i,1,1,1)*phi_i(x)
then
bm1(i,1,1,1) = \int phi_i(x) dx
where the integral (\int) is taken over the entire computational domain.
[ Note that we really don't have such a global phi_i(x) ... I use it simply
for illustration, here.]
hth
Paul
From: Nek5000-users <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: Monday, February 26, 2018 11:59:19 AM
To: nek5000-users
Subject: [Nek5000-users] Volume of an element
Hi Neks,
How to calculate volume of an element in a 3D mesh? If i write
n = nx1*ny1*nz1*nelv
do i=1,n
v(i) = volume(1,1,1,i)
enddo
Can i store the volume of each element into a varibale called v?
What happens if i use bm1(1,1,1,i)? Instead of volume(1,1,1,i)
Thanks
Sijo
_______________________________________________
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
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/nek5000-users/attachments/20180227/7b625b25/attachment.html>
More information about the Nek5000-users
mailing list