[Nek5000-users] Spectral coefficients for velocity
nek5000-users at lists.mcs.anl.gov
nek5000-users at lists.mcs.anl.gov
Wed May 25 10:06:16 CDT 2016
Hi Andrey,
I don't know if you got an answer yet but here are my suggestions.
If you are asking about the theory behind spectral transforms, I
suggest that you look at appendix B in the book "High-Order Methods
for Incompressible Fluid Flow" by Fischer, Deville and Mund. In
particular, Section B.3.1 gives the formula for the spectral
coefficients uh_k
uh_k = 1/gamma_k \int_-1^1 w(x) u(x) p_k(x) dx
where gamma_k=||p_k||^2 is the discrete norm and w(x) is the weight
function, both of them associated to polynomial p_k.
If you want to compute these coefficients in Nek5000, then you can
look at the routine "local_err_est" in the file navier5.f. There, one
goes to the Legendre space by using
call tensr3(uh,nx,u,nx,Lj,Ljt,Ljt,w)
where uh is the array of spectral elements and the operators Lj and
Ljt are computed in the routine "build_legend_transform" (also in
navier5.f) from the GLL points.
I hope that helps.
Best regards,
Nicolas Offermans
Quoting nek5000-users at lists.mcs.anl.gov:
> Hi Neks,
>
> Priviously, the parameters p101 and p103 have been discussed:
>
> "In each direction, on each element, the solution is represented as a
> polynomial of degree N in the reference element (i.e., on the interval
> [-1,1]).
>
> For any Nth-order polynomial there is a modal expansion of the form:
>
> u(x) = \sum_k=0^N uh_k phi_k(x)
>
> where phi_0 = 1, phi_1 = x, and phi_k(x)=L_k(x) - L_{k-1} (x) for k > 1,
> where L_k(x)=kth Legendre polynomial."
>
> How can I obtain these coefficients "uh_k" on the several elements?
>
> Thanks in advance for any suggestion,
>
> Andrey
More information about the Nek5000-users
mailing list