I am trying to setup an adaptive quadrature rule that will use different
quadrature points/weights on different elements in the mesh. (In this
particular application, the number of quadrature points depends on the
size of the element.)
I am not sure how to get this to work, because QBase::init() seems to
default to using precomputed points/weights whenever you pass in the
same element type and p level, even for cases in which
QBase::shapes_need_reinit() returns true.
Modifying QBase::init() so that it does NOT skip initializing the
quadrature rule if shapes_need_reinit() seems to have the desired
effect: the quadrature rule appears to be reinitialized for each element.
I wonder if this (or something similar) should be the default behavior?
Is there already another / better way to accomplish this?
Thanks!
-- Boyce
Index: src/quadrature/quadrature.C
===================================================================
--- src/quadrature/quadrature.C (revision 4510)
+++ src/quadrature/quadrature.C (working copy)
@@ -31,7 +31,7 @@
{
// check to see if we have already
// done the work for this quadrature rule
- if (t == _type && p == _p_level)
+ if (!shapes_need_reinit() && t == _type && p == _p_level)
return;
else
{
------------------------------------------------------------------------------
All of the data generated in your IT infrastructure is seriously valuable.
Why? It contains a definitive record of application performance, security
threats, fraudulent activity, and more. Splunk takes this data and makes
sense of it. IT sense. And common sense.
http://p.sf.net/sfu/splunk-d2d-c2
_______________________________________________
Libmesh-devel mailing list
Libmesh-devel@lists.sourceforge.net
https://lists.sourceforge.net/lists/listinfo/libmesh-devel