Nek5000
SEM for Incompressible NS
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
Data Types | Public Member Functions | Private Member Functions | List of all members
helmholtz Module Reference

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)
 

Detailed Description

Definition at line 15 of file navier4.F90.

Member Function/Subroutine Documentation

subroutine, private helmholtz::gensh ( real(dp), dimension (*), intent(inout)  v1,
real(dp), dimension (*), intent(in)  h1,
real(dp), dimension (*), intent(in)  h2,
real(dp), dimension(*), intent(in)  vml,
real(dp), dimension(*), intent(in)  vmk,
type(approx_space), intent(inout)  apx,
real(dp), dimension(:), intent(out ws 
)
private

Reconstruct the solution to the original problem by adding back the approximate solution, add solution to approximation space.

Parameters
[in,out]v1Full solution
[in]h1coefficient of A (stiffness)
[in]h2coefficient of M (mass)
[in]vmkmultiplicity array
[in]vmlmask array
[out]wssmall workspace to pass-through
[in,out]apxcurrent approximation space

Definition at line 215 of file navier4.F90.

References copy(), and hconj().

Referenced by hsolve().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

subroutine, private helmholtz::hconj ( type(approx_space), intent(inout)  apx,
integer, intent(in)  k,
real(dp), dimension(*), intent(in)  h1,
real(dp), dimension(*), intent(in)  h2,
real(dp), dimension(*), intent(in)  vml,
real(dp), dimension(*), intent(in)  vmk,
real(dp), dimension(*), intent(out ws 
)
private

Update the k-th row/column of H_red.

Parameters
[in,out]apxCurrent approximation space
[in]kIndex of new projector
[in]h1coefficient of A (stiffness)
[in]h2coefficient of M (mass)
[in]vmlmultiplicity array
[in]vmkmask array
[out]wssmall workspace

Definition at line 250 of file navier4.F90.

References axhelm(), ctimer::dnekclock(), dssum(), and gop().

Referenced by gensh(), and updrhsh().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

subroutine, private helmholtz::hmhzpf ( character(4)  name,
real(dp), dimension (lx1,ly1,lz1,1), intent(out u,
real(dp), dimension (lx1,ly1,lz1,1), intent(in)  r,
real(dp), dimension (lx1,ly1,lz1,1), intent(in)  h1,
real(dp), dimension (lx1,ly1,lz1,1), intent(in)  h2,
real(dp), dimension (lx1,ly1,lz1,1), intent(in)  mask,
real(dp), dimension (lx1,ly1,lz1,1), intent(in)  mult,
integer  imesh,
real(dp)  tli,
integer  maxit,
integer  isd,
real(dp), dimension (lx1,ly1,lz1,1), intent(in)  bi 
)
private
Parameters
[out]usolution vector
[in]rright hand side
[in]h1coefficient of A (stiffness)
[in]h2coefficient of M (mass)
[in]maskmask array
[in]multmultiplicity array
[in]biinverse of mass matrix

Definition at line 344 of file navier4.F90.

References cggo(), chktcg1(), and ctimer::dnekclock().

Referenced by hsolve().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

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.

Parameters
[in]namename of field we're solving for
[out]usolution vector
[in,out]rright hand side
[in]h1coefficient of A (stiffness)
[in]h2coefficient of M (mass)
[in]vmkmask array
[in]vmlmultiplicity array
[in]imshimesh?
[in]tolresidual tolerance
[in]maxitmaximum number of iterations
[in]isdsomething to do with axi-symmetric
[in,out]apxcurrent approximation space
[in]biinverse 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().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

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().

+ Here is the caller graph for this function:

subroutine, private helmholtz::projh ( real(dp), dimension(*), intent(inout)  r,
real(dp), dimension(*), intent(in)  h1,
real(dp), dimension(*), intent(in)  h2,
real(dp), dimension(*), intent(in)  bi,
real(dp), dimension(*), intent(in)  vml,
real(dp), dimension(*), intent(in)  vmk,
type(approx_space), intent(inout)  apx,
real(dp), dimension(*), intent(out wl,
real(dp), dimension(*), intent(out ws,
character(4), intent(in)  name4 
)
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)

Parameters
[in,out]rresidual
[in]h1coefficient of A (stiffness)
[in]h2coefficient of M (mass)
[in]vmlmultiplicity array
[in]vmkmask array
[in]biinverse mass matrix
[out]wllarge work array (size lx1*ly1*lz1*nelv)
[out]wssmall work array (size 2*max vecs)
[in,out]apxCurrent approx space
[in]name4Name of field for debug printing
Note
This dsyev call and the following dgemv are task-parallel!

Definition at line 62 of file navier4.F90.

References axhelm(), ctimer::dnekclock(), dssum(), glsc2(), gop(), and updrhsh().

Referenced by hsolve().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

subroutine, private helmholtz::updrhsh ( type(approx_space), intent(inout)  apx,
real(dp), dimension(*), intent(in)  h1,
real(dp), dimension(*), intent(in)  h2,
real(dp), dimension(*), intent(in)  vml,
real(dp), dimension(*), intent(in)  vmk,
real(dp), dimension(*), intent(out ws 
)
private

Recompute H_red if dt has changed.

Parameters
[in,out]apxcurrent approximation space
[in]h1coefficient of A (stiffness)
[in]h2coefficient of M (mass)
[in]vmlmultiplicity array
[in]vmkmask array
[out]wssmall workspace

Definition at line 301 of file navier4.F90.

References hconj().

Referenced by projh().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:


The documentation for this module was generated from the following file: