That sounds good.  The whole process will be described in detail in my 
dissertation, which is due in a couple weeks.  If you send me directions 
on how to document it on the website, I'd be happy to add that in then.
-rd

On 10/22/2012 10:03 AM, Luca Antiga wrote:
> Hi Richard,
>   excellent, it's always nice to have good feedback about vmtk.
> Thanks a lot for offering to share. In fact, I think even documenting the 
> experience
> in a short page on the vmtk website with a few images would be great. Sharing
> code would be invaluable.
> Just think about it and let me know, in case you decide for it I'll give you 
> directions
> for contributing.
> Thanks a lot again
>
> Luca
>
>
> On Oct 22, 2012, at 4:54 PM, Richard Downe wrote:
>
>> Luca, thanks.
>> In the meantime, I noticed, digging around in the code, that I could pass 
>> source/target points into vmtkcenterlines, and get the same result; I might 
>> explore passing my own in at some point in the future, but the mesh 
>> generation is a tiny slice of the running time, so I may be satisfied where 
>> I am now.  (Moreover, my centerlines never actually meet, so using your code 
>> provides me with the necessary guarantee of intersections in appropriate 
>> places.)
>>
>> Adding the centerlines, however, has a) drastically reduced the odds of 
>> tetgen failures (previously on some of the more difficult meshes, tetgen 
>> would have to run 2 or 3 times before it got a convergence), and the nominal 
>> running time in fluent has dropped from 36 hours to 6.
>>
>> I did have to change edgelengthfactor to 0.2 (at 0.3, the mesh was a bit too 
>> coarse; when I tried 0.1 I got massive meshes that caused fluent to consume 
>> 24G ram), but at this point things are running very smoothly.
>>
>> Thanks again, if you think anything from my pipeline may be useful to others 
>> in the future, let me know and I'd be happy to share.
>> -rd
>>
>> On 10/22/2012 08:54 AM, Luca Antiga wrote:
>>> Hi Richard,
>>>   sorry for the long wait.
>>>
>>> Your pipeline looks good, but you'll need the radius at each point
>>> (you seem proficient in Python and VTK so I won't provide you with
>>> details on how to add a vtkDoubleArray to PointData, but let me know
>>> if you need directions).
>>>
>>> In any case, the mesh generator needs a scalar field over the input
>>> surface to be remeshed. That scalar field will be used (multiplied by
>>> a factor) as the nominal size of the triangle edges in the remeshed
>>> surface (as demonstrated in the mesh generation tutorial on the vmtk
>>> website).
>>>
>>> The way you generate that scalar field is really up to you, you could
>>> do it through curvatures or by evaluating the distance of the surface
>>> from centerlines, or the closest centerline radius. The latter is probably
>>> what you're aiming at now, so as soon as you have your dataset with
>>> centerlines and radii (in a point array named, say, Radius), you can do
>>>
>>> vmtkdistancetocenterlines -ifile surface.vtp -centerlinesfile 
>>> centerlines.vtp
>>> -radiusarray Radius -useradius 1 -centerlineradius 1 --pipe
>>> vmtkmeshgenerator -elementsizemode edgelengtharray -edgelengtharray
>>> Radius -edgelengthfactor 0.3 -ofile foo.vtu
>>>
>>> The first script will generate the centerline radius field on the surface
>>> and the second will remesh the surface taking that field into account for
>>> the local triangle size.
>>>
>>> Hope this helps
>>>
>>>
>>> Luca
>>>
>>>
>>> On Oct 12, 2012, at 6:43 PM, Richard Downe wrote:
>>>
>>>> Upon observing the behavior of my automated setup for a week, I do think
>>>> I'm going to have to include the centerlines.
>>>> I'm reasonably confident in the quality of my triangulated surfaces, but
>>>> I do think the fixed edge length passed to vmtkmeshgenerator is overly
>>>> limiting.
>>>>
>>>> My question is this: what arguments to I need to feed into
>>>> vmtkcenterlinemerge if I want to have an appropriate vtkPolyData for use
>>>> as such in subsequent processing?
>>>>
>>>> I have the code sample below for simply putting the centerlines into a
>>>> vtkPolyData; I have radius information available, but I'm curious what
>>>> it needs to look like internally for the rest of the vmtkpipeline to be
>>>> happy with it.
>>>> Thanks.
>>>> -rd
>>>>
>>>>      def ComputeCenterlinePolyData(self, fusmesher):
>>>>          vesselPoints = vtk.vtkPoints()
>>>>          vesselLines = vtk.vtkCellArray()
>>>>          vesselPoly = vtk.vtkPolyData()
>>>>
>>>>          vesselCenterline = fusmesher.GetVesselCenterline()
>>>>          vesselPoints.SetNumberOfPoints(len(vesselCenterline))
>>>>          vesselLines.InsertNextCell(len(vesselCenterline))
>>>>
>>>>          for i in range(len(vesselCenterline)):
>>>>              p = vesselCenterline[i]
>>>>              vesselPoints.SetPoint(i, p['x'], p['y'], p['z'])
>>>>              vesselLines.InsertCellPoint(i)
>>>>
>>>>          vesselPoly.SetPoints(vesselPoints)
>>>>          vesselPoly.SetLines(vesselLines)
>>>>
>>>>          branchPolys = []
>>>>
>>>>          for i in range(fusmesher.GetBranchCount()):
>>>>              branchPoints = vtk.vtkPoints()
>>>>              branchLines = vtk.vtkCellArray()
>>>>              branchPoly = vtk.vtkPolyData()
>>>>
>>>>              branchCenterline = fusmesher.GetBranchCenterline(i)
>>>>              branchPoints.SetNumberOfPoints(len(branchCenterline))
>>>>              branchLines.InsertNextCell(len(branchCenterline))
>>>>
>>>>              for j in range(len(branchCenterline)):
>>>>                  p = branchCenterline[j]
>>>>                  branchPoints.SetPoint(j, p['x'], p['y'], p['z'])
>>>>                  branchLines.InsertCellPoint(j)
>>>>
>>>>              branchPoly.SetPoints(branchPoints)
>>>>              branchPoly.SetLines(branchLines)
>>>>              branchPolys.append(branchPoly)
>>>>
>>>>          appender = vtk.vtkAppendPolyData()
>>>>          appender.AddInput(vesselPoly)
>>>>
>>>>          for i in range(fusmesher.GetBranchCount()):
>>>>              appender.AddInput(branchPolys[i])
>>>>
>>>>          return appender.GetOutput()
>>>>
>>>> On 10/09/2012 05:41 AM, Luca Antiga wrote:
>>>>> Hi Richard,
>>>>>
>>>>> On Oct 2, 2012, at 5:08 PM, Richard Downe wrote:
>>>>>
>>>>>> I have been attempting to use vmtk to generate meshes for use with 
>>>>>> fluent.
>>>>>> Because the nature of my workflow involves geometric transformations of
>>>>>> segmented borders, my starting point is a point set (well, more
>>>>>> accurately, several disjoint structured grids).
>>>>>> I have a program that generates a topologically sound triangular surface
>>>>>> with minimal triangle angle of 24 degrees, flow extensions added, and
>>>>>> caps removed at boundaries.
>>>>>>
>>>>>> My coordinates are also spaced in meters, so my edge lengths are
>>>>>> typically on the order of 0.0002.
>>>>>> After a few false starts, I have successfully been able to run
>>>>>> vmtkmeshgenerator on the STL file output by my initial triangulation,
>>>>>> and pipe it into vmtkmeshwriter to create a fluent file, and use a pipe
>>>>>> I wrote myself to identify the inlets based on proximity of the entities
>>>>>> to known locations of centerlines at boundaries.
>>>>> Sounds like great work.
>>>>>
>>>>>> The first surface I attempted ran in fluent fine.  2 subsequent, more
>>>>>> complex surfaces are causing fluent to fail in mysterious ways that are
>>>>>> almost certainly related to mesh quality.  I am invoking
>>>>>> vmtkmeshgenerator with -edgelength 0.0002 as the only argument.  Are
>>>>>> there known steps I can take to alleviate this/improve mesh quality?
>>>>>> There are clearly a large number of people using vmtk successfully, so I
>>>>>> can only assume I could be doing something better.  (When I examine the
>>>>>> mesh, fluent's quality rating is typically between 10 and 20, with a
>>>>>> maximum aspect ratio value of around 35).
>>>>> I'd really need to look into the mesh, would it be possible for you to 
>>>>> send
>>>>> me the stl file?
>>>>>
>>>>>> My current tack is to attempt to add the centerlines into the flow using
>>>>>> vmtkdistancetocenterlines, on the assumption that this will cause the
>>>>>> sizing function to generate more useful information, and thereby give me
>>>>>> better shaped tetrahedra -- but I could use advice, as I do somewhat
>>>>>> feel like I'm swinging in the dark.
>>>>> Depending on the nature of your domain, it is possible that 0.0002 as an
>>>>> element size ends up being inadequate for certain regions and overkill
>>>>> for others, so the use of centerline information almost certainly helps.
>>>>>
>>>>> However, I'd first take a look at the surface to see if it rings any bell.
>>>>>
>>>>> Best,
>>>>>
>>>>>
>>>>> Luca
>>>>>
>>>>>
>>>>>> Thanks.
>>>>>> --Richard
>>>>>>
>>>>>> ------------------------------------------------------------------------------
>>>>>> Don't let slow site performance ruin your business. Deploy New Relic APM
>>>>>> Deploy New Relic app performance management and know exactly
>>>>>> what is happening inside your Ruby, Python, PHP, Java, and .NET app
>>>>>> Try New Relic at no cost today and get our sweet Data Nerd shirt too!
>>>>>> http://p.sf.net/sfu/newrelic-dev2dev
>>>>>> _______________________________________________
>>>>>> vmtk-users mailing list
>>>>>> vmtk-users@lists.sourceforge.net
>>>>>> https://lists.sourceforge.net/lists/listinfo/vmtk-users
>>>> ------------------------------------------------------------------------------
>>>> Don't let slow site performance ruin your business. Deploy New Relic APM
>>>> Deploy New Relic app performance management and know exactly
>>>> what is happening inside your Ruby, Python, PHP, Java, and .NET app
>>>> Try New Relic at no cost today and get our sweet Data Nerd shirt too!
>>>> http://p.sf.net/sfu/newrelic-dev2dev
>>>> _______________________________________________
>>>> vmtk-users mailing list
>>>> vmtk-users@lists.sourceforge.net
>>>> https://lists.sourceforge.net/lists/listinfo/vmtk-users


------------------------------------------------------------------------------
Everyone hates slow websites. So do we.
Make your web apps faster with AppDynamics
Download AppDynamics Lite for free today:
http://p.sf.net/sfu/appdyn_sfd2d_oct
_______________________________________________
vmtk-users mailing list
vmtk-users@lists.sourceforge.net
https://lists.sourceforge.net/lists/listinfo/vmtk-users

Reply via email to