New issue
Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.
By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.
Already on GitHub? Sign in to your account
Request for input: Laplace-based linear and bilinear integrators #4238
Comments
Implementation effort can be reduced by using the |
Is there a need for a |
This sounds great to me!
I'm not sure we really need to target exact integration of Laplacian terms, even if that is possible in some cases. I think using the lowest quadrature order (e.g., ignoring the transformation order) that works robustly is best. It is hard to test all possible use cases for "works robustly", so go with what you find in your tests or what is recommended in the literature, if there are such recommendations. How do you plan to compute the Laplacian on a general element? Do you plan to use exact computation at quadrature points using the Hessian and the Jacobian of the transformation? An approximate alternative that only requires the Jacobian is to project the gradient (e.g. using nodal interpolation) back into the same FE space and then take the divergence. I think @bslazarov and/or @vladotomov wanted to try something like that but I'm not 100% certain. |
The Laplacian/Hessian is already implemented - in principle. The Laplacian is compute as the trace of the Hessian which sounds expensive but is the easiest way to have the correct effect of the mapping . There is a check of to see if the Transformation is affine, |
@michi002 My apologies, I pushed it to an other local repo. I will push it tomorrow, unfortunately I can not ssh into my machine due to overly zealous security... |
@IdoAkkerman thank you. I had a look on your commit. |
@michi002 Convection-diffusion only requires the scalar implementation. However, the SUPG/PSPG terms for Stokes/Navier-Stokes involve terms where the velocity and pressure spaces are combined, or is there a way to bypass that? Regarding, the term you are interested in, that can be obtained by by using the 4th term, viz. |
@IdoAkkerman Concerning the coefficient: |
Implementing Do you have a clear application in mind? |
@IdoAkkerman I do not know if it is an important case for the general userbase, but I think a mixed version ( |
So you are converting Why not test with |
@IdoAkkerman If you are really interested why we do that, you can look at our preprint: https://arxiv.org/abs/2404.10350 I hope that answers your question. |
I am interested in implementing
LinearFormIntegrator
andBilinearFormIntegrator
involving the Laplacian.Any heads up on potential conflicting or duplicate effort is highly appreciated.
Any tips or warning for pitfalls are also welcome.
Last but not least, if I am skipping some useful functionality for other use cases please let me know.
My personal interest for these terms are stabilized methods. So first would be the scalar convection-diffusion case which would require the following:
Coefficient f
Coefficient q
VectorCoefficient a
Coefficient q
VectorCoefficient a
Coefficient q
Next would be vector based case, such as navier-stokes. Which would require the same as above but all the scalar shape functions become vector shape function. The coefficients will in principle be of Matrix form, probably more efficient routines are available for scalar or vector (= diagonal of matrix). I will consider implementing all three.
To enable the pressure and PSPG terms two mixed terms are required.
MatrixCoefficient Q
MatrixCoefficient Q
Hessian based integrators will not be considered in this PR.
The text was updated successfully, but these errors were encountered: