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

nek5000-users at lists.mcs.anl.gov nek5000-users at lists.mcs.anl.gov
Sun Jul 11 17:00:26 CDT 2010


Hi Markus,

It's possible that it's returning the local (to the processor)
element number.  You would get the global (to the computational
domain, defined by the .rea file) element number via:

        integer e,eg

        eg = lglel(e)


Paul


On Sun, 11 Jul 2010, nek5000-users at lists.mcs.anl.gov wrote:

> Hi again,
>
> Regarding finding element numbers (such as what gllel would return) for 
> arbitrary history points, I was experimenting with the history point routine 
> hpts a little and it seems that this function call:
>
>        call findpts(inth_hpts,rcode,1,
>     &                 proc,1,
>     &                 elid,1,
>     &                 rst,ndim,
>     &                 dist,1,
>     &                 pts(1,1),ndim,
>     &                 pts(2,1),ndim,
>     &                 pts(3,1),ndim,npoints)
>
> does not return the same element number (elid) as it can be found in the .rea 
> file, for instance.
> For a handful of manually set points, I checked that the returned element 
> number does not provide the right coordinates to enclose the desired point.
>
> Is there a way to transform the function return or am I missing something 
> else?
>
> Thanks,
> Markus
>
> _______________________________________________
> 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