Nek5000
SEM for Incompressible NS
|
Data Types | |
type | approx_space |
Type to hold the approximation space. Should not be modified outside this module, so more of a handle. More... | |
Public Member Functions | |
subroutine, public | init_approx_space (apx, n_max, ntot) |
Initialize approximation space object. More... | |
subroutine, public | hsolve (name, u, r, h1, h2, vmk, vml, imsh, tol, maxit, isd, apx, bi) |
Either std. Helmholtz solve, or a projection + Helmholtz solve. More... | |
Private Member Functions | |
subroutine, private | projh (r, h1, h2, bi, vml, vmk, apx, wl, ws, name4) |
Project out the part of the residual in the approx space. More... | |
subroutine, private | gensh (v1, h1, h2, vml, vmk, apx, ws) |
Reconstruct the solution to the original problem by adding back the approximate solution, add solution to approximation space. More... | |
subroutine, private | hconj (apx, k, h1, h2, vml, vmk, ws) |
Update the k-th row/column of H_red. More... | |
subroutine, private | updrhsh (apx, h1, h2, vml, vmk, ws) |
Recompute H_red if dt has changed. More... | |
subroutine, private | hmhzpf (name, u, r, h1, h2, mask, mult, imesh, tli, maxit, isd, bi) |
Definition at line 15 of file navier4.F90.
|
private |
Reconstruct the solution to the original problem by adding back the approximate solution, add solution to approximation space.
[in,out] | v1 | Full solution |
[in] | h1 | coefficient of A (stiffness) |
[in] | h2 | coefficient of M (mass) |
[in] | vmk | multiplicity array |
[in] | vml | mask array |
[out] | ws | small workspace to pass-through |
[in,out] | apx | current approximation space |
Definition at line 215 of file navier4.F90.
References copy(), and hconj().
Referenced by hsolve().
|
private |
Update the k-th row/column of H_red.
[in,out] | apx | Current approximation space |
[in] | k | Index of new projector |
[in] | h1 | coefficient of A (stiffness) |
[in] | h2 | coefficient of M (mass) |
[in] | vml | multiplicity array |
[in] | vmk | mask array |
[out] | ws | small workspace |
Definition at line 250 of file navier4.F90.
References axhelm(), ctimer::dnekclock(), dssum(), and gop().
Referenced by gensh(), and updrhsh().
|
private |
[out] | u | solution vector |
[in] | r | right hand side |
[in] | h1 | coefficient of A (stiffness) |
[in] | h2 | coefficient of M (mass) |
[in] | mask | mask array |
[in] | mult | multiplicity array |
[in] | bi | inverse of mass matrix |
Definition at line 344 of file navier4.F90.
References cggo(), chktcg1(), and ctimer::dnekclock().
Referenced by hsolve().
subroutine, public helmholtz::hsolve | ( | character(4), intent(in) | name, |
real(dp), dimension (lx1,ly1,lz1,lelv), intent(out) | u, | ||
real(dp), dimension (lx1,ly1,lz1,lelv), intent(inout) | r, | ||
real(dp), dimension (lx1,ly1,lz1,lelv), intent(in) | h1, | ||
real(dp), dimension (lx1,ly1,lz1,lelv), intent(in) | h2, | ||
real(dp), dimension (lx1,ly1,lz1,lelv), intent(in) | vmk, | ||
real(dp), dimension (lx1,ly1,lz1,lelv), intent(in) | vml, | ||
integer, intent(in) | imsh, | ||
real(dp), intent(in) | tol, | ||
integer, intent(in) | maxit, | ||
integer, intent(in) | isd, | ||
type(approx_space), intent(inout) | apx, | ||
real(dp), dimension (lx1,ly1,lz1,*), intent(in) | bi | ||
) |
Either std. Helmholtz solve, or a projection + Helmholtz solve.
[in] | name | name of field we're solving for |
[out] | u | solution vector |
[in,out] | r | right hand side |
[in] | h1 | coefficient of A (stiffness) |
[in] | h2 | coefficient of M (mass) |
[in] | vmk | mask array |
[in] | vml | multiplicity array |
[in] | imsh | imesh? |
[in] | tol | residual tolerance |
[in] | maxit | maximum number of iterations |
[in] | isd | something to do with axi-symmetric |
[in,out] | apx | current approximation space |
[in] | bi | inverse of mass matrix |
Definition at line 396 of file navier4.F90.
References string::capit(), chcopy(), dssum(), gensh(), hmholtz(), hmhzpf(), and projh().
Referenced by cdscal(), and plan4().
subroutine, public helmholtz::init_approx_space | ( | type(approx_space), intent(out) | apx, |
integer, intent(in) | n_max, | ||
integer, intent(in) | ntot | ||
) |
Initialize approximation space object.
Simple assigns and allocations
Definition at line 41 of file navier4.F90.
Referenced by plan4().
|
private |
Project out the part of the residual in the approx space.
Starts by finding the EVD H_red to get orthogonal projectors. Then, computes the overlap of the residual with each projector and mixes them based on the EVD to get an orthogonal projection. Next, multiplies by 1/ to take the inverse and expands back into the full space to populate the approximate solution: projectors(:,0)
[in,out] | r | residual |
[in] | h1 | coefficient of A (stiffness) |
[in] | h2 | coefficient of M (mass) |
[in] | vml | multiplicity array |
[in] | vmk | mask array |
[in] | bi | inverse mass matrix |
[out] | wl | large work array (size lx1*ly1*lz1*nelv) |
[out] | ws | small work array (size 2*max vecs) |
[in,out] | apx | Current approx space |
[in] | name4 | Name of field for debug printing |
Definition at line 62 of file navier4.F90.
References axhelm(), ctimer::dnekclock(), dssum(), glsc2(), gop(), and updrhsh().
Referenced by hsolve().
|
private |
Recompute H_red if dt has changed.
[in,out] | apx | current approximation space |
[in] | h1 | coefficient of A (stiffness) |
[in] | h2 | coefficient of M (mass) |
[in] | vml | multiplicity array |
[in] | vmk | mask array |
[out] | ws | small workspace |
Definition at line 301 of file navier4.F90.
References hconj().
Referenced by projh().