Wolfgang,
> It doesn't have to be. For example, for FE_Q, we also build on
> TensorProductPolynomials but the ordering is not lexicographic. So it would
> still be of interest to document the order of shape functions if you end up
> finding out what it is!
Noted. So I have verified this with the following code.
#include <deal.II/base/quadrature_lib.h>
#include <deal.II/base/point.h>
#include <deal.II/fe/fe_dgq.h>
#include <iostream>
#include <cmath>
#include <vector>
using namespace dealii;
int main(int argc, char **argv)
{
const int degree = 4;
const int n_poly1 = degree+1;
const FE_DGQLegendre<1> fe1(degree);
const FE_DGQLegendre<3> fe3(degree);
const QGauss<3> quad(degree+1);
const int n_qp = quad.size();
const std::vector<Point<3>>& points = quad.get_points();
const std::vector<double>& weights = quad.get_weights();
for(int k=0; k<n_poly1; k++){
for(int j=0; j<n_poly1; j++){
for(int i=0; i<n_poly1; i++){
int index_3d = i + j*n_poly1 + k*n_poly1*n_poly1;
// compute error between 'index_3d'-th polynomial of 'fe3' and the product
of
// i-th, j-th and k-th polynomials of 'fe1'
double error = 0;
for(int q=0; q<n_qp; q++){
error += weights[q]*std::pow(
fe3.shape_value(index_3d, points[q]) - (
fe1.shape_value(i, Point<1>(points[q][0]))*
fe1.shape_value(j, Point<1>(points[q][1]))*
fe1.shape_value(k, Point<1>(points[q][2]))
),
2
);
}
std::cout << "Tensor indices: " << i << " " << j << " " << k << ", "
<< "3d index: " << index_3d << ", error: " << error << "\n";
}
}
}
}
And the output shows all errors to be 0! If this is ok, I will create a pr
with a patch to the documentation shortly.
I don't recall whether there is an easy way to achieve what you are looking
> for, but some elements definitely do something like what you are trying to
> achieve. For example, the difference between FE_Q and FE_DGQ is, in
> essence,
> just a change of basis where the basis functions are permuted. Similarly,
> the
> difference between FE_Q and FE_QHierarchical is similar to the nodal ->
> modal
> change you are interested in. Finally, there is also the case of the FE_Q
> constructor that receives a Quadrature object as argument and that then
> computes a basis change.
> You might want to look into how all of these are implemented. Most of this
> kind of functionality exists in some kind of helper function that might be
> useful to you.
Thanks a lot for providing these details. For now, I have done this by
myself, but I will keep this in mind.
Thanks very much!
On Tue, 25 May 2021 at 19:35, Wolfgang Bangerth <[email protected]>
wrote:
>
> Vachan
>
> > Thanks a lot for your reply. What I actually need is a change of basis
> from
> > Lagrange polynomials (nodal) to Legendre polynomials (modal). I then
> want to
> > know the coefficients of certain modes.
> >
> > So, if there is any straightforward way to do this in dealii, I would
> proceed
> > with that. I, could not find any such functions and hence planned on
> doing it
> > manually, using the shape functions. However, even using any inbuilt
> functions
> > would only address my issue partly if the ordering is not clear.
>
> I don't recall whether there is an easy way to achieve what you are
> looking
> for, but some elements definitely do something like what you are trying to
> achieve. For example, the difference between FE_Q and FE_DGQ is, in
> essence,
> just a change of basis where the basis functions are permuted. Similarly,
> the
> difference between FE_Q and FE_QHierarchical is similar to the nodal ->
> modal
> change you are interested in. Finally, there is also the case of the FE_Q
> constructor that receives a Quadrature object as argument and that then
> computes a basis change.
>
> You might want to look into how all of these are implemented. Most of this
> kind of functionality exists in some kind of helper function that might be
> useful to you.
>
> Best
> W.
> --
> ------------------------------------------------------------------------
> Wolfgang Bangerth email: [email protected]
> 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 [email protected].
> To view this discussion on the web visit
> https://groups.google.com/d/msgid/dealii/4d767eb0-cf66-b11a-8791-55deab3666de%40colostate.edu
> .
>
--
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 [email protected].
To view this discussion on the web visit
https://groups.google.com/d/msgid/dealii/CALAVa_yMw_i_hmj4GuB6iT0BOOcv2MC1DoUE1OEWfO9j%3D7s9Uw%40mail.gmail.com.