[Nek5000-users] dealias; mass matrix; interp_rstd
nek5000-users at lists.mcs.anl.gov
nek5000-users at lists.mcs.anl.gov
Wed May 23 07:37:14 CDT 2018
Dear Sandeep,
Look for routines intp_rst, grad_rstd, etc., and calls to them
in convect.f and induct.f
There are several routines that map onto the finer GL mesh
from the base GLL mesh.
The metrics, J*dr/dx, J*dr/dy, etc., are already on the fine
mesh and collocated with the quadrature weight on that fine
mesh when they are constructed (in set_dealias_rx in induct.f).
The mathematics of these procedure is described in Sec. 3 of
Stabilization of the spectral element method in convection dominated flows by recovery of skew-symmetry
* Johan Malm,
* Philipp Schlatter
* Paul F. Fischer
* Dan S. Henningson
NOTE: After you create your nonlinear source term, you
should divide the value each point x_{ijke} by bm1(i,j,k,e),
because Nek multiplies the entry in qvol by bm1(i,j,k,e).
In your case, however, you've already accounted for the mass
matrix by incorporated the mass matrix on the fine-GL points.
I would set up your u x nabla term via a call to userchk
and store the result (already divided by bm1), then just
dereference it inside qvol. via:
subroutine qvol(i,j,k,eg)
...
...
common /mynonlin/ s(lx1,ly1,lz1,lelt,ldim)
integer e,eg
e=gllel(eg)
qvol=s(i,j,k,e,ifield-2)
return
end
or something like that.
The time step will all be correct and you will get 2nd or
3rd-order accuracy in time via Nek's standard BDFk/EXTk
hth,
Paul
________________________________
From: Nek5000-users <nek5000-users-bounces at lists.mcs.anl.gov> on behalf of nek5000-users at lists.mcs.anl.gov <nek5000-users at lists.mcs.anl.gov>
Sent: Wednesday, May 23, 2018 4:50:19 AM
To: nek5000-users at lists.mcs.anl.gov
Subject: [Nek5000-users] dealias; mass matrix; interp_rstd
Dear Neks,
I solve
\partial A/\partial t = \nabla^2 A + u x Curl A,
where A = (s1,s2,s3)
I solve three scalar equations with respective components of u x Curl A as forcing terms in qvol.
I would like to compute the dealiased product of u x Curl A
I build my code based on kov_stokes_test.usr and subroutines in induct.f
I compute local curl, project curl and U onto fine mesh, perform product and project the product back to coarse mesh using interp_rstd
However, I am not sure about the multiplication of mass matrices for the product in these subroutines.
I attach my usr, rea and SIZE files.
Your suggestions are welcome.
Thanks,
Sandeep
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/nek5000-users/attachments/20180523/357c5dd3/attachment.html>
More information about the Nek5000-users
mailing list