Hi Martin, Thanks a lot for the elucidation! I found out after your answer that there is already an AnisotropicPolynomials class. It was remiss of me not even doing a quick search...
I noticed that FE_DGQ class is derived from FE_Poly with a polynomial space argument: *TensorProductPolynomials <https://dealii.org/current/doxygen/deal.II/classTensorProductPolynomials.html><dim>(Polynomials <https://dealii.org/current/doxygen/deal.II/namespacePolynomials.html>::generate_complete_Lagrange_basis(internal <https://dealii.org/current/doxygen/deal.II/namespaceinternal.html>::FE_DGQ <https://dealii.org/current/doxygen/deal.II/classFE__DGQ.html>::get_QGaussLobatto_points(degree))* It seems that if we write a FE_DGQAnistropic class, we'd probably want to replace this argument with a counterpart of the AnisotropicPolynomials class. Do you think this replacement alone would cut it? We are not looking for any more features out of the fe space than what FE_DGQ is already doing. FE_DGQ involves implementing quite a bunch of virtual functions <https://dealii.org/current/doxygen/deal.II/fe__dgq_8h_source.html#l00375>. Apart from, as far as I can parse, FE_DGQ<dim, spacedim>::rotate_indices <https://dealii.org/current/doxygen/deal.II/classFE__DGQ.html#ab7491f8ec8f7081f1f989f1e3e5da5f7>, they seem to be coordinate direction agnostic. Would really appreciate if you can shed more light on this! Best, Greg On Sunday, March 12, 2023 at 8:33:02 AM UTC Martin Kronbichler wrote: > Dear Greg, > > Just to extend on what Peter said: It should be possible to define > anisotropic degrees in a completely discontinuous finite element, say > FE_DGP or FE_DGQ. In that case, you would derive from `FE_Poly` and insert > an anisotropic polynomial space, doing whatever you think is appropriate. > For elements imposing some continuity, include FE_FaceP/FE_FaceQ or the > regular H1 functions (FE_Q), you can't do that as Peter said. However, one > could work around this by constructing constraints that constrain the > polynomial space of uniform high order in one of the directions to the > lower polynomial degree. It would be a bit of work to realize this (say > 200-400 lines of code to set up the constraints), but the easiest solution > I could come up. > > Best, > Martin > > > On 11.03.23 17:26, Greg Wang wrote: > > Hi Peter, > > Thanks a lot for the info! > > We have recently implemented in deal.II a space-time HDG method for the > advection-diffusion problem on deforming domains. Our analysis allows the > spatial polynomial order p_s to differ from the temporal one p_t. When p_t > = p_s, we just call FE_DGP, FE_FaceP constructors with nothing extra to do > and the rates come out nicely matching the a priori rates of convergence. > And now we are trying to have a p_t = p_s +1 test case to further validate > our analysis, which led to this post. > > Best regards, > Greg > > On Saturday, March 11, 2023 at 2:20:37 PM UTC peterr...@gmail.com wrote: > >> Hi Greg, >> >> unfortunately we do not support anisotropic polynomials: we make the >> assumption that all faces (of the same type) and lines have the same number >> of DoFs. For what do you need such polynomial spaces? >> >> Best, >> Peter >> >> On Saturday, 11 March 2023 at 15:05:25 UTC+1 ygre...@gmail.com wrote: >> >>> Hi there, >>> >>> I was wondering if finite elements with anisotropic polynomial degrees >>> are possible in deal.II. As an example, for a 3D element, can we construct >>> a tensor product polynomial space of {1,x,y,z,xy,xz,yz, z^2,xz^2,yz^2}, >>> i.e., order 1 in x- and y-directions and 2 in z-direction? >>> >>> I was looking at, for example, constructors of FE_DGQ() class [1]. The >>> second constructor [2] takes an arbitrary vector of polynomials to build >>> the tensor product polynomial space. This is close to what I want but it >>> seems the argument can only be one-dimensional polynomials, which means >>> equal order on all dimensions. >>> >>> Would really appreciate any insights and/or tips! >>> >>> Best, >>> Greg >>> >>> ------------------------------------------ >>> [1] https://dealii.org/current/doxygen/deal.II/classFE__DGQ.html >>> [1] >>> https://dealii.org/current/doxygen/deal.II/fe__dgq_8cc_source.html#l00100 >>> >> -- > 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+un...@googlegroups.com. > To view this discussion on the web visit > https://groups.google.com/d/msgid/dealii/3dcc5829-7218-4c51-9c86-ae8fccb4ba82n%40googlegroups.com > > <https://groups.google.com/d/msgid/dealii/3dcc5829-7218-4c51-9c86-ae8fccb4ba82n%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/7951a3b7-6b00-40fd-95d8-94cba5f4968an%40googlegroups.com.