[Nek5000-users] r,s,t index for given point
nek5000-users at lists.mcs.anl.gov
nek5000-users at lists.mcs.anl.gov
Wed Jul 14 14:55:06 CDT 2010
Howdy,
I am not sure if it is of use for anybody else, but the element number
as returned from findpts (vector elid), as for example used in hpts(),
can be converted to a global element number (as it would appear in the
rea file) by:
"
lglel(elid(i)+1)
"
I guess the jl/jl2 routines start counting the (local) element index
from 0 instead of 1 - sounds simple here, but I probably pulled out a
third of my scalp hair to figure it out :-)
Also, the variables returned from findpts, even though this function
needs to be called on every node in parallel, are only known to node0,
which makes it necessary to do something like this if one really wants
to know the global number:
"
call bcast(proc,lhis*isize)
call bcast(elid,lhis*isize)
call bcast(npoints,isize) !Don't do this before jl2 functions are done
!Zero out globelnum
do i=1,npoints
globelnum(i)=0
enddo
do i=1,npoints
if (nid.eq.proc(i)) then
globelnum(i)=lglel(elid(i)+1)
endif
enddo
call igop(globelnum,w1,'+ ',lhis)
"
Markus
nek5000-users at lists.mcs.anl.gov wrote:
>
> Hi Markus,
>
> lglel reads "local to global element number" so input is
> local element number and output is global element number.
>
> I'll be seeing Stefan this week so I'm sure we can straighten
> out the interpolation routines.
>
> You're right about the other points - the ntotv wasn't flagged
> by my compiler and the .his file was taken care of by my own
> scripts, so I didn't catch either of them. Sorry about that.
> I've updated the source so that those issues should be
> resolved.
>
> Paul
>
>
> On Mon, 12 Jul 2010, nek5000-users at lists.mcs.anl.gov wrote:
>
>> Hi,
>>
>> that was not it, I'm afraid (though I thought gllel transforms from
>> local to global).
>>
>> But I also realized that the hpts routine does not deliver the correct
>> values for me. To test this, I set up a simple 3d box with inflow,
>> outflow and symmetry all around, Re very low, initial
>> condition=boundary condition and 4 points.
>> Please find all files attached, the output for the first point is OK,
>> but then the values seem to shift positions and don't make sense.
>>
>> Also, when running it for 1 step or more, I get the error message:
>> "
>> dump history points
>> done :: dump history points
>> hisfile:
>> /Users/m0s1978/nektonruns/rectest/rectest.his
>> ERROR: .sch file already exists. ierr= 1
>> "
>>
>> And lastly, the newest repo won't compile, here is the error message:
>> "
>> Error: /Users/m0s1978/nek5trial/trunk/nek/drive2.f, line 1551: This
>> name has already been assigned a data type. [NVTOT]
>> integer*8 ntot,nvtot,nptot
>> ---------------------^
>> "
>>
>> Thanks for any help,
>> Markus
>>
>>
>> nek5000-users at lists.mcs.anl.gov wrote:
>>>
>>> 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
>>>>
>>> _______________________________________________
>>> 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