On Thu, Nov 6, 2008 at 11:03 AM, Roy Stogner <[EMAIL PROTECTED]> wrote:
> John Peterson wrote:
>
>> 1.) A possible refinement pattern for the pyramid involves 10
>> children: 4 sub-pyramids on the base, one in the apex, and one pyramid
>> upside-down below the apex pyramid.  Then, there are 4 tets filling in
>> the gaps.  (This is a little hard to visualize, but I can send you an
>> image if you want to see what I mean.)
>
> Sadly, this is probably the *easiest* alternative to visualize.
>
>> One bad thing about this
>> refinement pattern is that under repeated refinements, you end up with
>> a little "ribbon" of pyramids running through the original element.
>> So, we probably need to discuss other possibilities, but this one
>> would at least get us started.
>
> Replace those two central pyramids with four tets?  Should be just as easy
> and arguably just as good, especially if we do the same sort of
> geometrically adaptive axis selection that we do in tet refinement now.

Hey, that's a pretty good idea!  I will move that to the top of the
list of candidates.

>> Since embedding_matrix() is a
>> virtual function, it can be redefined in Pyramid sub-classes to do
>> whatever we want, e.g. we could place -1 entries in the actual tables
>> and just ensure, through embedding_matrix(), that those entries are
>> never accessed.  The real issue seems to be determining where in the
>> library, if anywhere, we've assumed that all the children have the
>> same geometric type.  I'd be glad to hear your comments there.
>
> There is an assumption in FEBase::coarsened_dof_values that elements and
> their children are oriented the same way - i.e. if child c has a side on
> side s of its parent, then the corresponding side number on c is also s.
>  That's actually an assumption that we don't *have* to break with either
> potential pyramid splitting, but in either case we'd have to be careful with
> the tet node ordering.

OK, interesting.  I wasn't aware of that particular function in FEBase
so I'm glad you mentioned it.  For non-homogeneous refinement
patterns, this assumption might not even make sense.  A pyramid has
five sides but its tet children will only have 4, for example.  It
probably wouldn't be too hard to remove that assumption from FEBase
but it would almost certainly be a performance hit to do so...


>> Anyone know anything about the
>> accuracy of quadrature for functions which are ratios of polynomials?
>
> We can derive custom quadrature rules which would integrate a mass matrix
> exactly... but would they then also integrate, say, a Laplacian matrix
> exactly?  The answer is an obvious "yes" for polynomial bases but I'd expect
> a "no" for pyramids.  That could be bad.
>
> What are we doing for them now?

The current quadrature rules have accuracies like you would expect for
1D elements, since they are conical products of Gauss-like rules.  So,
for example, a 2x2x2 rule will integrate exactly all monomials of the
form x^a y^b z^c, a+b+c <= 3.  I have no idea what will happen when we
try to integrate the rational basis functions...

-- 
John

-------------------------------------------------------------------------
This SF.Net email is sponsored by the Moblin Your Move Developer's challenge
Build the coolest Linux based applications with Moblin SDK & win great prizes
Grand prize is a trip for two to an Open Source event anywhere in the world
http://moblin-contest.org/redirect.php?banner_id=100&url=/
_______________________________________________
Libmesh-devel mailing list
[email protected]
https://lists.sourceforge.net/lists/listinfo/libmesh-devel

Reply via email to