Nek5000
SEM for Incompressible NS
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
Functions/Subroutines
plan4.F90 File Reference

Routine plan4 and supporting routines. More...

Go to the source code of this file.

Functions/Subroutines

subroutine plan4 ()
 Splitting scheme [2]. More...
 
subroutine crespsp (respr, vext)
 Compute start residual/right-hand-side in the pressure. More...
 
subroutine cresvsp (resv1, resv2, resv3, h1, h2)
 Compute the residual for the velocity. More...
 
subroutine op_curl (w1, w2, w3, u1, u2, u3, ifavg, work1, work2)
 Compute curl of U. More...
 
subroutine v_extrap (vext)
 Extrapolate the velocity forward in time with AB(k) More...
 

Detailed Description

Routine plan4 and supporting routines.

This file contains high-level routines for computing the right hand sides for the Helmholtz and Poisson solves in the P(N)-P(N) formalism (colocation).

Definition in file plan4.F90.

Function/Subroutine Documentation

subroutine crespsp ( real(dp), dimension (lx2,ly2,lz2,lelv), intent(out respr,
real(dp), dimension (lx1*ly1*lz1*lelv,3), intent(in)  vext 
)

Compute start residual/right-hand-side in the pressure.

\( R = \vec{F_b} - \nabla \times \nabla \times \vec{v} - \nabla \left(P + \frac{1}{3} H_1 \nabla \cdot Q \right) \)

where \( F_b \) is the boundary force, \( A \) is the stiffness matrix, \( B \) is the mass matrix, \( \vec{v} \) is the velocity, \( P \) is the pressure, and \( \nabla \cdot Q \) is the "thermal divergence"

Definition at line 219 of file plan4.F90.

References axhelm(), cdtp(), ctimer::dnekclock(), faccl2(), faccl3(), mxm(), op_curl(), opdssum(), and ortho().

Referenced by plan4().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

subroutine cresvsp ( real(dp), dimension(lx1,ly1,lz1,lelv), intent(out resv1,
real(dp), dimension(lx1,ly1,lz1,lelv), intent(out resv2,
real(dp), dimension(lx1,ly1,lz1,lelv), intent(out resv3,
real(dp), dimension (lx1,ly1,lz1,lelv), intent(out h1,
real(dp), dimension (lx1,ly1,lz1,lelv), intent(out h2 
)

Compute the residual for the velocity.

Definition at line 418 of file plan4.F90.

References ctimer::dnekclock(), opgrad(), ophx(), and sethlm().

Referenced by plan4().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

subroutine op_curl ( real(dp), dimension(lx1,ly1,lz1,lelv), intent(out w1,
real(dp), dimension(lx1,ly1,lz1,lelv), intent(out w2,
real(dp), dimension(lx1,ly1,lz1,lelv), intent(out w3,
real(dp), dimension(lx1,ly1,lz1,lelv), intent(in)  u1,
real(dp), dimension(lx1,ly1,lz1,lelv), intent(in)  u2,
real(dp), dimension(lx1,ly1,lz1,lelv), intent(in)  u3,
logical, intent(in)  ifavg,
real(dp), dimension(lx1,ly1,lz1,lelv), intent(out work1,
real(dp), dimension(lx1,ly1,lz1,lelv), intent(out work2 
)

Compute curl of U.

\( (w_1, w_2, w_3) = \nabla \times (u_1, u_2, u_3) \)

Parameters
[out]w11st component of curl U
[out]w22nd component of curl U
[out]w33rd component of curl U
[in]u11st component of U
[in]u22nd component of U
[in]u33rd component of U
[out]work1work array
[out]work2work array
[in]ifavgAverage at boundary?

Definition at line 480 of file plan4.F90.

References copy(), dudxyz(), mxm(), and opdssum().

Referenced by crespsp().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

subroutine plan4 ( )

Splitting scheme [2].

NOTE: QTL denotes the so called thermal divergence and has to be provided by an external subroutine e.g qthermal

Note
These three calls are task-parallel

Definition at line 14 of file plan4.F90.

References bcdirvc(), crespsp(), cresvsp(), ctolspl(), ctimer::dnekclock(), dssum(), helmholtz::hsolve(), helmholtz::init_approx_space(), lagvel(), makef(), opdiv(), ortho(), and v_extrap().

Referenced by fluid().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

subroutine v_extrap ( real(dp), dimension(lx1,ly1,lz1,lelv,3), intent(out vext)

Extrapolate the velocity forward in time with AB(k)

This is the first half of (6.5.8) in HOMfIFF [1] : \( \hat{v} = \sum_{j=1}^k \beta_{k-j} v^{n+1-j} + ... \)

Parameters
[out]vextvelocity (3 components) extrapolated forward in time

Definition at line 604 of file plan4.F90.

Referenced by plan4().

+ Here is the caller graph for this function: