Hi.

I have done this in the past. I dug up the wrapper function I had for it, in 
case we decide that we should add it to the library:

https://github.com/dealii/dealii/pull/15838

The problem with this approach is that evaluating the reference space level set 
function gets very expensive. QuadratureGenerator uses many function calls and 
we have to use the Jacobian and its derivative in the transformation. This 
makes this approach significantly slower than using an interpolated level set 
function.

If you are working with a perfectly Cartesian background I would expect that it 
would actually be faster to generate the quadrature in real space and then 
transform the quadrature to reference space before passing it to FEValues. If you 
want to do this the functions cell->bounding_box() and 
BoundingBox::real_to_unit(..) are useful.

Best,
Simon

On 03/08/2023 23:11, Anton wrote:
Hello,

I would like to use the NonMatching::QuadratureGenerator<dim> class directly, 
as I need to generate multiple (as in many!, unfortunately) nonmatching quadratures 
for the same cell.  Therefore I would like to avoid the route of interpolating the 
level set function onto the whole mesh etc, as described in CutFEM tutorial.

Therefore I am thinking about calling the "generate" function, which only needs 
a level set function and a bounding box as inputs.  According to the documentation the 
level set function should provide the gradient and the Hessian.  I would like to do this 
in the reference cell coordinates, so that I can simply use this quadrature later via 
FEView as one normally does.
  * The bounding box I presume is the hypercube [0,1]^{dim} for the meshes I am 
interested in. (Can one request them for the reference cell and not the 
triangulation cell?)

  * The level set function should also be reasonably straightforward.  I know the 
function in the physical coordinates (say, f(x)) and all its derivatives.  I can also get 
the map from reference to physical coordinates via FEView, let us call this map x=M(y).  
The issue is how to extract the derivatives of this map so that the chain rule can be 
applied?  I am mostly working with Q1 maps, so of course I could figure out the 
derivatives "by hand".  I am just wondering if there is a more intelligent 
(code-wise) way of doing this?

--
The deal.II project is located at http://www.dealii.org/ 
<http://www.dealii.org/>
For mailing list/forum options, see https://groups.google.com/d/forum/dealii?hl=en 
<https://groups.google.com/d/forum/dealii?hl=en>
---
You received this message because you are subscribed to the Google Groups "deal.II 
User Group" group.
To unsubscribe from this group and stop receiving emails from it, send an email to 
dealii+unsubscr...@googlegroups.com 
<mailto:dealii+unsubscr...@googlegroups.com>.
To view this discussion on the web visit 
https://groups.google.com/d/msgid/dealii/4ac38563-e467-4f1d-8cb4-ee35379515a6n%40googlegroups.com
 
<https://groups.google.com/d/msgid/dealii/4ac38563-e467-4f1d-8cb4-ee35379515a6n%40googlegroups.com?utm_medium=email&utm_source=footer>.

--
The deal.II project is located at http://www.dealii.org/
For mailing list/forum options, see 
https://groups.google.com/d/forum/dealii?hl=en
--- You received this message because you are subscribed to the Google Groups "deal.II User Group" group.
To unsubscribe from this group and stop receiving emails from it, send an email 
to dealii+unsubscr...@googlegroups.com.
To view this discussion on the web visit 
https://groups.google.com/d/msgid/dealii/f1642630-6453-2028-1345-3db743d814f2%40gmail.com.

Reply via email to