[vtk-developers] Pyramid and QuadraticPyramid formulation issues

Bill Lorensen bill.lorensen at gmail.com
Mon Jun 22 16:12:30 EDT 2015


Andy,

Thanks,

These are the sorts of things I am uncovering in the new unit test I'm
working  on.

Bill

On Mon, Jun 22, 2015 at 2:17 PM, Andy Bauer <andy.bauer at kitware.com> wrote:
> Hi Bill,
>
> It looks to me like the apex of vtkPyramid is at {.5, .5, 1}. I used this
> test code to play around with it:
> from vtk import *
> ugrid = vtkUnstructuredGrid()
> pts = vtkPoints()
> pts.InsertNextPoint(0, 0, 0)
> pts.InsertNextPoint(1, 0, 0)
> pts.InsertNextPoint(1, 1, 0)
> pts.InsertNextPoint(0, 1, 0)
> pts.InsertNextPoint(.5, .5, 1)
> ugrid.SetPoints(pts)
> ugrid.Allocate(10)
> p = [0, 1, 2, 3, 4]
> ugrid.InsertNextCell(VTK_PYRAMID, 5, p)
> c = ugrid.GetCell(0)
> pcoords = [0, 0, 0]
> x = [0, 0, 0]
> w = [0, 0, 0, 0, 0]
> subid = vtk.mutable(0)
> # EvaluateLocation() is from parametric to global coordinate
> c.EvaluateLocation(subid, pcoords, x, w)
> print x
> pcoords = [.5, .5, 1]
> c.EvaluateLocation(subid, pcoords, x, w)
>
> # EvaluatePosition() is from global to parametric coordinate
> x = [.5, .5, 1]
> closestpoint = [-1, -1, -1]
> dist2 = vtk.mutable(0)
> result = c.EvaluatePosition(x, closestpoint, subid, pcoords, dist2, w)
> print result, pcoords
>
>
>
>
>
> I didn't test vtkQuadraticPyramid though. Looking at the vtkPyramid.h though
> I see the following which looks wrong to me:
> //----------------------------------------------------------------------------
> inline int vtkPyramid::GetParametricCenter(double pcoords[3])
> {
>   pcoords[0] = pcoords[1] = 0.4;
>   pcoords[2] = 0.2;
>   return 0;
> }
>
>
> I would have thought that pcoords[0] and pcoords[1] should be 0.5 instead of
> 0.4. Since the volume of the above cell is 1/3 I would also think that
> pcoords[2] should be 1/3 instead of 0.2, but I may be wrong on that. It's
> also kind of weird that in vtkQuadraticPyramid.h that method is:
> //----------------------------------------------------------------------------
> // Return the center of the quadratic pyramid in parametric coordinates.
> //
> inline int vtkQuadraticPyramid::GetParametricCenter(double pcoords[3])
> {
>   pcoords[0] = pcoords[1] = 6./13;
>   pcoords[2] = 3./13;
>   return 0;
> }
>
> Just my 2 cents on the matter though since I haven't gone through the code
> that closely.
>
> Cheers,
> Andy
>
> On Mon, Jun 22, 2015 at 6:29 AM, Will Schroeder <will.schroeder at kitware.com>
> wrote:
>>
>> Bill it could be that the vtkPyramid functions were originally written (a
>> long time ago) as degenerate hex shape functions which might explain the
>> apex parametric coordinates. I believe Mathieu found and implemented an
>> alternative formulation and also added the quadratic version; I'm not sure
>> where to find that reference.
>> W
>>
>> On Sat, Jun 20, 2015 at 2:31 PM, Bill Lorensen <bill.lorensen at gmail.com>
>> wrote:
>>>
>>> Folks,
>>>
>>> As I continue to work through a comprehensive unit test for vtk cells,
>>> I have found an issue(s) with vtkPyramid and vtkQuadraticPyramid.
>>>
>>> For example, the parametric coordinates for the Pyramid (and Quadratic
>>> Pyramid) place the apex at 0, 0, 1. I believe it should be centered on
>>> the base (which would be .5, .5, 1).
>>>
>>> Changing that causes other issues with the shape functions.
>>>
>>> I am seeking a copy of the reference quoted in the QuadraticPyramid
>>> code so that I can verify the implementation.
>>> // The shape functions and derivatives could be implemented thanks to
>>> // the report Pyramid Solid Elements Linear and Quadratic Iso-P Models
>>> // From Center For Aerospace Structures
>>
>>
>>
>>
>> --
>> William J. Schroeder, PhD
>> Kitware, Inc.
>> 28 Corporate Drive
>> Clifton Park, NY 12065
>> will.schroeder at kitware.com
>> http://www.kitware.com
>> (518) 881-4902
>>
>> _______________________________________________
>> Powered by www.kitware.com
>>
>> Visit other Kitware open-source projects at
>> http://www.kitware.com/opensource/opensource.html
>>
>> Search the list archives at: http://markmail.org/search/?q=vtk-developers
>>
>> Follow this link to subscribe/unsubscribe:
>> http://public.kitware.com/mailman/listinfo/vtk-developers
>>
>>
>



-- 
Unpaid intern in BillsBasement at noware dot com


More information about the vtk-developers mailing list