[Trilinos-Users] [EXTERNAL] Transforming Laplacians of Intrepid FEM basis functions

Roberts, Nathan V. nvroberts at alcf.anl.gov
Mon Mar 21 15:16:25 EDT 2016


Pavel and Kara,

Thanks for your replies.  As Pavel surmised, I am testing the Laplacian term with a $b \cdot grad v$ term, so that integration by parts will only shift the problem of computing second-order terms from the trial to the test space.  I'm attempting to implement "Method 5" from this paper (p. 18):

http://www.asc.tuwien.ac.at/~schoeberl/wiki/publications/MixedHybridDG.pdf

This does involve an element-wise definition of the stabilization parameter tau.  I do indeed intend to use higher-order elements in my implementation.

But I am open to rethinking the approach—this is part of a comparison of several methods, and for SUPG my hope would be to have a method that is reflective of what's usually done using SUPG.  Do you have a recommendation for an approach that avoids higher-order derivatives in the stabilization term, even in the context of higher-order basis functions?  (Kara, I would indeed be interested in seeing your SUPG stabilization examples using Intrepid.)

Thanks!
Nate

On Mar 21, 2016, at 12:54 PM, Bochev, Pavel B <pbboche at sandia.gov<mailto:pbboche at sandia.gov>> wrote:

Nate,

I would also advise against using strong form of the Laplacian. The term $-\Delta u, b\cdot\grad u$ could be potentially destabilizing if tau is nor carefully controlled on element by element basis.

Are you using higher order elements in your implementation? If you use P1/Q1 this term does not have to be included anyways.

If you must have it then, the proper way to transform it is to simply follow the change of variables formula for a single second order scalar derivative. Keep in mind, the HGRAD transformation serves a different purpose, as it is inteneded to transform vector field in a way that preserves the tangent component. You don’t need to treat D2 as a colective 2nd order tensor, just work out the change of variables formula for u_{xx} and u_{yy}, which I believe you have done in your writeup.

Pavel

From: Trilinos-Users <trilinos-users-bounces at trilinos.org<mailto:trilinos-users-bounces at trilinos.org>> on behalf of "Peterson, Kara J." <kjpeter at sandia.gov<mailto:kjpeter at sandia.gov>>
Date: Monday, March 21, 2016 at 10:54 AM
To: Nate Roberts <nvroberts at alcf.anl.gov<mailto:nvroberts at alcf.anl.gov>>, "<trilinos-users at trilinos.org<mailto:trilinos-users at trilinos.org>>" <trilinos-users at trilinos.org<mailto:trilinos-users at trilinos.org>>
Subject: Re: [Trilinos-Users] [EXTERNAL] Transforming Laplacians of Intrepid FEM basis functions

Hi Nate,

The transformation you included below does look correct, but there may be an easier way of implementing a stabilization term that involves a Laplacian than directly computing the second derivative of the basis functions.

Although I haven’t actually implemented this myself, I understand that hyperviscosity-like stabilizations that involve Laplacians are usually applied using a weak notion of the Laplacian in finite elements. So given a test function  v

\int (\nabla \cdot \nabla u) v d\Omega  = \int \nabla u \cdot \nabla v d\Omega + boundary terms.

In this way the Laplacian of a variable can be written as M^{-1}K u for a mass matrix M and stiffness matrix K.

Also, if it would be helpful to you I do have examples of SUPG stabilization using Intrepid, but the implementations did not involve any explicit Laplacian terms.

Best,
Kara

________________________________
Kara Peterson
Computational Mathematics
Sandia National Laboratories
505.844.9372


From: Trilinos-Users <trilinos-users-bounces at trilinos.org<mailto:trilinos-users-bounces at trilinos.org>> on behalf of "Roberts, Nathan V." <nvroberts at alcf.anl.gov<mailto:nvroberts at alcf.anl.gov>>
Date: Monday, March 21, 2016 at 8:15 AM
To: "<trilinos-users at trilinos.org<mailto:trilinos-users at trilinos.org>>" <trilinos-users at trilinos.org<mailto:trilinos-users at trilinos.org>>
Subject: [EXTERNAL] [Trilinos-Users] Transforming Laplacians of Intrepid FEM basis functions

Hello all,

I'm implementing an SUPG formulation for convection-diffusion, and the stabilization term for this involves a Laplacian.  So, I need to take second-order derivatives of the FEM basis functions in Intrepid.  I see how to do this in reference space; I can apply OPERATOR_D2 and use the Intrepid::getDkEnumeration() to navigate the returned values.

My question is this: how do I appropriately transform the values to physical space?  I.e., the second-order analogue to the FunctionSpaceTools::HGRADtransformGRAD()?

Below, I offer what I've worked out thus far—I'm hoping that someone can (a) confirm that what I have is reasonable, and (b) comment on the simplest/best way to approach this in terms of what's available in Intrepid.

Mathematically, what HGRADtransformGRAD() is doing goes something like this.  If the physical function is u and the reference-space function is \hat{u}, and there is a vector transformation function x = x(\xi) from reference to physical space, then the gradient of u can be computed as:


This can be computed by multiplying the inverse of the Jacobian of the reference-to-physical map by the reference space gradient values; this is what HGRADtransformGRAD does.

Proceeding to the second derivatives needed for the Laplacian, we have


I believe the first term on the right hand side is computable in terms of the reference-space OPERATOR_D2 values and the inverse of the Jacobian applied twice (I would be glad to see example code if anyone has that to offer).  The second term requires second derivatives of the reference-to-physical map (i.e. the Hessian), as well as first derivatives of the reference-space function.

Is that correct?  Can anyone comment on the best way to implement this?

Thanks!
Nate

P.S. In case my pasted PDF graphics do not come through, the LaTeX for the first equation is:
\frac{\partial u}{\partial x_i} = \frac{\partial \hat{u}}{\partial \xi_k} \frac{\partial \xi_k}{\partial x_i}

and the second is:
\frac{\partial^2 u}{\partial x_i^2} = \frac{\partial^2 \hat{u}}{\partial \xi_k \partial \xi_l} \frac{\partial \xi_k}{\partial x_i} \frac{\partial \xi_l}{\partial x_i} + \frac{\partial \hat{u}}{\partial \xi_k} \frac{\partial^2 \xi_k}{\partial x_i^2}

_______________________________________________
Trilinos-Users mailing list
Trilinos-Users at trilinos.org<mailto:Trilinos-Users at trilinos.org>
https://trilinos.org/mailman/listinfo/trilinos-users

-------------- next part --------------
An HTML attachment was scrubbed...
URL: <https://trilinos.org/pipermail/trilinos-users/attachments/20160321/13c1acca/attachment.html>


More information about the Trilinos-Users mailing list