[Nek5000-users] r,s,t index for given point

nek5000-users at lists.mcs.anl.gov nek5000-users at lists.mcs.anl.gov
Sat Jul 10 15:32:04 CDT 2010


Hi,

yes, that helped, thank you.
That also explains why I saw the original point coordinates when I 
replaced vx, vy&vz with xm1,ym1,zm1 in the interpolation routine - the 
hpts function does not write out values at the nearest grid point (as 
the old one did), but instead interpolates all field values at the 
desired, arbitrary points.

I think what I'll do for the sake of simplicity is to find the element 
number from findpts and then use the routine from postx (loops over 
lx1*ly1*lz1 for the specified element) to find the gll indices. These 
values are needed only once and thus I don't have to pay too close 
attention to speed.

Thanks again,
Markus


nek5000-users at lists.mcs.anl.gov wrote:
> 
> Hi Markus,
> 
> Functions in the spectral element method are expressed as:
> 
>        e                  e
> [1]   u (r,s,t) = sum    u     h (r) h (s) h (t)
>                    ijk    ijk   i     j     k
> 
> 
> Where u^e_ijk are the basis coefficients for element e, with
> i,j,k ranging from 0 to N  (1 to lx1 in the code), and the
> h_i's are polynomials of degree N on the interval from [-1,1].
> Given r,s,t and the set of basis coefficients, one can compute
> u^e(r,s,t) using [1].
> 
> Note that the h_i's are Lagrangian interpolants on the 
> Gauss-Lobatto-Legendre points, {xi_i}, i=0,...,N, such that
> 
>       h_i (xi_j)  = delta_ij,   ( = 1 for i=j, 0 otherwise ).
> 
> Thus for any point (r,s,t) = (xi_i,xi_j,xi_k), we can associate basis 
> coefficient u^e_ijk with u^e at that point.  More specifically,
> we can associated the index set (i,j,k) with that point.
> 
> However, if the continuous variables (r,s,t) do not coincide with
> the quadrature points then there is no associated index --- that is,
> the value of u^e(r,s,t) would then depend on _all_ the basis
> coefficients, using Eq. [1] above.
> 
> The interpolation routine returns the continuous variables - as
> are needed for the general interpolation formula [1].
> 
> To get the indices, you could check the rst values against the
> GLL coordinate points that are stored in WZ in the array zgm1(),
> for instance.
> 
> Does this help ?
> 
> Paul
> 
> 
> 
> 
> 
> 
> 
> 
> On Fri, 9 Jul 2010, nek5000-users at lists.mcs.anl.gov wrote:
> 
>> Hi,
>>
>> what is the difference between "index" and "value"?
>> Essentially, I am looking for the information postx provides when 
>> figuring out history points. I browsed the postx source code, located 
>> the routine and implemented it in .usr, but I expect trouble when 
>> moving to bigger domains. I also figured there must be something in nek.
>>
>> So I took the new hpts routine (see attached my modified .usr file), 
>> what I am getting for r,s,t for a specific point is:
>> 1.00000000000000      -5.551115123125783E-017  -1.00000000000000
>>
>> The real deal for the same point (from postx) is 4    1    1.
>>
>> What am I missing? Is there a way to go from the output I have to the 
>> r,s,t indices?
>>
>> Thanks for any help,
>> Markus
>>
>>
>> nek5000-users at lists.mcs.anl.gov wrote:
>>>
>>> Markus,
>>>
>>> those r,s,t values (not indices) are valid  - do they not correspond
>>> to your expectations?  Off-hand one might expect that the associated
>>> indices would be 1,NY1/2,NZ1
>>>
>>> Paul
>>>
>>>
>>> On Fri, 9 Jul 2010, nek5000-users at lists.mcs.anl.gov wrote:
>>>
>>>> Howdy Neks,
>>>>
>>>> is there a function in Nek that for a given point in x,y,z would 
>>>> return processor number, element number and r,s,t index of the 
>>>> nearest GLL grid point? I was desperately trying to figure out how 
>>>> the interpolation routines (f.e. for history points) would help me, 
>>>> but I only got garbage (-1, 0 and 1) back for the variables I would 
>>>> expect contain r,s,t. It worked with the element number, however.
>>>>
>>>> Thanks,
>>>> Markus
>>>>
>>>> _______________________________________________
>>>> 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
> 



More information about the Nek5000-users mailing list