Re: [deal.II] Non-tensor polynomials with bilinear elements
On 5/3/19 4:14 PM, Doug wrote: > Ok, so my initial understanding of deal.II was accurate. It all makes sense > now. Although FE_DGP is spanning {1,x,y} in reference space, it is not > spanning those bases in physical space. Instead, it is spanning some > combination of {1,x,y,xy} with 3 degrees of freedom, that only spans > {1,x,y} in physical space if the mapping is affine. As a result, loss of > convergence order entails since it cannot properly represent linear > functions in physical space. > > On the other hand, a bilinear basis on the reference space can indeed > represent physical linear functions if a bilinear mapping is used! Hence > why FE_DGQ work. Yes. You could also use the FE_DGPNonparametric class in that case. Best W. -- Wolfgang Bangerth email: bange...@colostate.edu www: http://www.math.colostate.edu/~bangerth/ -- 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. For more options, visit https://groups.google.com/d/optout.
Re: [deal.II] Non-tensor polynomials with bilinear elements
Prof. Bangerth, Ok, so my initial understanding of deal.II was accurate. It all makes sense now. Although FE_DGP is spanning {1,x,y} in reference space, it is not spanning those bases in physical space. Instead, it is spanning some combination of {1,x,y,xy} with 3 degrees of freedom, that only spans {1,x,y} in physical space if the mapping is affine. As a result, loss of convergence order entails since it cannot properly represent linear functions in physical space. On the other hand, a bilinear basis on the reference space can indeed represent physical linear functions if a bilinear mapping is used! Hence why FE_DGQ work. If we were working with triangles, the Legendre basis {1,x,y} would be able to represent linear functions in physical space since the mapping is automatically affine for straight sided elements, hence my wrong assumption that Legendre polynomials should produce optimal orders of convergence. Thank you for the response, it was immensely helpful. Doug -- 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. For more options, visit https://groups.google.com/d/optout.
Re: [deal.II] Non-tensor polynomials with bilinear elements
On 5/2/19 7:26 PM, Doug wrote: > > Based on this, I assume that FEValues, which takes a FiniteElement and a > Mapping, uses them to evaluate the Jacobian field instead of using a bilinear > mapping as the documentation of Mapping suggests. So FE_DGP would only be > able > to parametrize the cell geometry with the polynomial space {1, x, y} instead > of {1, x, y, xy}. No, that's not correct. The *element* may only have 3 basis functions, but the mapping is still bilinear (with four basis functions). These two aspects are separate -- in fact, the mapping never knows anything about the element anyway. What this means is that on cells that are not parallelograms, the *mapped* basis functions may not be linear. You can probably try this out easily if you just do a single cell that's not a parallelogram, output a solution vector with one 1 and the rest all zeros, and use DataOut::build_patches(N) with a sufficiently large N. > The documentation of FE_DGP does mention similar to this, but I thought this > only applied to the solution, not the geometric mapping. I would expect the > geometric mapping to always use a tensor product basis in order to be able to > represent any straight-sided quads, especially since deal.II only handles > quads. Yes, and it does. > Since the FEValues seems to evaluate Jacobian based on the provided > FiniteElement, FE_DGP is basically useless since it can only handle truly > linear elements. I don't know if that's true (that it's useless), but it's worth keeping the Jacobian of the mapping and the Jacobians of the shape functions separate -- they're not the same. Best Wolfgang -- Wolfgang Bangerth email: bange...@colostate.edu www: http://www.math.colostate.edu/~bangerth/ -- 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. For more options, visit https://groups.google.com/d/optout.
[deal.II] Non-tensor polynomials with bilinear elements
Hello, I am not sure if FE_DGP and FE_DGPMonomial is behaving as expected. I had recently commented on a thread about distorted meshes ( https://groups.google.com/forum/#!topic/dealii/ZSHIDbp7yfE). I am using DG to solve a linear advection problem with a manufactured solution. A 2D square with perfectly straight cartesian elements, MappingQ, FE_DGP (Legendre polynomials), would give me optimal (p+1) order of convergence with h-refinement. However, if the h-refinement is followed by a small random distortion through (distort_random), the optimal convergence is lost. If a distorted mesh is h-refined subsequently (such that the distortion is only on the coarse grid), then optimal convergence orders are achieved. As Prof. Bangerth mentionned, this refinement results in elements that tend to parallelograms, which have a linear Jacobian. Now, if I use FE_DGQ (tensor Lagrange) or FE_DGQLegendre (tensor-product Legendre), I do obtain optimal orders of convergence, no matter the distortion. Based on this, I assume that FEValues, which takes a FiniteElement and a Mapping, uses them to evaluate the Jacobian field instead of using a bilinear mapping as the documentation of Mapping suggests. So FE_DGP would only be able to parametrize the cell geometry with the polynomial space {1, x, y} instead of {1, x, y, xy}. The documentation of FE_DGP does mention similar to this, but I thought this only applied to the solution, not the geometric mapping. I would expect the geometric mapping to always use a tensor product basis in order to be able to represent any straight-sided quads, especially since deal.II only handles quads. Since the FEValues seems to evaluate Jacobian based on the provided FiniteElement, FE_DGP is basically useless since it can only handle truly linear elements. It's not really an issue for me since I will just switch to a tensor-product basis. However, I am curious whether or not my understanding of what is happening is correct, and if this behaviour is expected. Best Regards, Doug -- 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. For more options, visit https://groups.google.com/d/optout.