Misc

Miscellaneous Function

This page documents a bunch of functions that are useful but don't fit into another catagory.

Actually, some of these functions do fit into another catagory and new files should be created to move them into.

A bunch of the things in euler_funcs.jl

EulerEquationMod.calcA0

This function calculates the A0 (ie. dq/dv, where q are the conservative and v are the entropy variables) for a node, and stores it A0

The formation of A0 is given in Hughes

Inputs: params : ParamType{Tdim, :entropy} q : vector of entropy variables, length Tdim + 2

Inputs/Outputs: A0 : (Tdim+2)x(Tdim+2) matrix populated with A0. Overwritten

Aliasing restrictions: none

source

EulerEquationMod.calcA0

This function calculates the A0 (ie. dq/dq, where q are the conservative variables at a node), and stores it A0. This function is provided for compatability purposes

Inputs: params : ParamType{2, :entropy} q : vector of conservative variables, length 4

Inputs/Outputs: A0 : 4x4 (or 5x5 in 3D) matrix populated with A0. Overwritten

Aliasing restrictions: none

source

EulerEquationMod.calcA0Inv

Calculates inv(A0), where A0 = dq/dv, where q are the conservative variables at a node and v are the entropy varaibles at a node, using the entropy variables.

Inputs: params : ParamType{Tdim, :entropy} q : vector (length 4 or 5) of entropy variables at a node

Inputs/Outputs: A0inv : matrix to be populated with inv(A0). Overwritten.

Aliasing restrictions: none

source

EulerEquationMod.calcA0Inv

Calculates inv(A0), where A0 = dq/dq, where q are the conservative variables at a node. This function is provided for compatability purposes

Inputs: params : ParamType{2, :entropy} q : vector (length Tdim + 2) of conservative variables at a node

Inputs/Outputs: A0inv : matrix to be populated with inv(A0). Overwritten.

Aliasing restrictions: none

source

EulerEquationMod.calcA1

This function calculates the A1 (ie. dF1/dq, where F1 is the first column of the Euler flux) for a node, aka the flux Jacobian of the Euler flux in the x direction. Methods are available for both conservative and entropy variables.

The formation of A1 is given in Hughes paper

Inputs: params : ParamType{2, :entropy} q : vector of entropy variables, length 4

Inputs/Outputs: A1 : 4x4 matrix to be populated. Overwritten

source

EulerEquationMod.calcA2

This function calculates A2 (ie. dF2/dq, where F2 is the second column of the Euler flux, aka the flux jacobian in the y direction. Methods are available for both conservative and entropy variables.

The formation of A2 is given in Hughes

Inputs: params : ParamType{2, :entropy}, q : vector of entropy variables, length 4 Inputs/Outputs: A2 : 4x4 matrix to be populated. Overwritten

source

EulerEquationMod.calcEntropy

This function calculates the entropy at a node and returns it. Method are available for conservative and entropy variables, 2D or 3D

Inputs: params: ParamType{Tdim, var_type}, used to dispatch to the right method. q: vector of solution variables at a node.

Returns: entropy

This is a low level function

Aliasing Restrictions: none

source

This function calculates the entropy function U used to define the IR variablesat a node and returns it. It does not return the physical entropy s. This function is agnostic to the dimension of the equation.

Inputs: params: a ParamType q: a vector of conservative variables at a node

Outputs: U: the value of the entropy function

Aliasing restrictions: none.

source

EulerEquationMod.calcEulerFlux

This function calculates the Euler flux from the conservative variables at a single node in a particular direction. 2D only.

Inputs: params : ParamaterType{2, :conservative} q : vector of conservative variables aux_vars : vector of auxiliary variables dir : vector in direction to calculate the flux

Inputs/Outputs: F : vector to populate with the flux

The Tdim paramater of params determine whether this method or the 3D version is called.

This is a low level function

source

low level function

Calculates the Euler flux from entropy variables

Inputs:
params : ParameterType{Tdim, :entropy}
q : vector of entropy variables
aux_vars : vector of auxiliary variables
dir : vector specifying the direction to caculate the flux

Inputs/Outputs:
F  : vector to populate with the flux

This is a low level function.  The static parameters of
the ParameterType are used to dispatch to the right method for any
combination of variable type or equation dimension.
source

EulerEquationMod.calcEulerFlux

This is the 3D method. All arguments are same as the 2D version.

source

###EulerEquationMod.calcEulerFlux_revm

Compute the derivative of the euler flux in reverse mode w.r.t to unit vector flux direction.

source

###EulerEquationMod.calcEulerFlux_revq

Compute the derivative of the Euler flux in reverse mode w.r.t q

source

EulerEquationMod.calcMaxWaveSpeed

This function calculates the maximum wave speed (ie. acoustic wave speed) present in the domain and returns it. Methods are available for conservative and entropy variables. This function uses eqn.q_vec, not eqn.q to do the calculation.

This is a mid level function

source

calcMomentContribution_rev!

This is the reverse differentiated version of calcMomentContribution!. See docs of calcMomentContribution! for further details of the primal method. This function is differentiated with respect to the primal version's xsbp and dforce variables.

Inputs

  • sbpface: an SBP face operator type

  • xsbp: SBP-face nodes in physical space; [coord, sbp node, face]

  • dforce: scaled force at the sbpface nodes; [coord, sbp node, face]

  • xyz_about: point about which the moment is taken

  • moment_bar: left multiplies d(moment)/d(xsbp) and d(moment)/d(dforce)

InOuts

  • xsbp_bar: result of vector Jacobian product; [coord, sbp node, face]

  • dforce_bar: result of vector Jacobian product; [coord, sbp node, face]

source

EulerEquationMod.calcPressure

This function calculates the pressure from the conservative variables at a node in 2D. It returns a single value.

Inputs: params : ParamType{Tdim, var_type } q : vector of conservative variables

The parameter of params determines whether the 2D or 3D, conservative or entropy method is dispatched.

This is a low level function.

Aliasing restrictions: none

source

EulerEquationMod.calcPressure

This function calculates pressure using the entropy variables.

Inputs: params : ParamType{2, :entropy} q : vector of entropy varaibles

returns pressure

Aliasing restrictions: none

source

EulerEquationMod.calcPressure

3D method. See 2D method documentation

source

###EulerEquationMod.calcPressure_revq

Compute the gradient of pressure w.r.t q in the reverse mode

Arguments

  • params : Parameter object

  • q : Forward sweep solution variable

  • press_bar : Reverse pressure gradient

  • q_bar : Reverse mode solution gradient

`

source

EulerEquationMod.calcSpeedofSound

This function calculates the speed of sound at a node and returns it. Methods are available for both conservative and entropy variables.

Inputs: params: ParamType{Tdim, var_type} q vector of solution variables at a node

Returns: speed of sound

This is a low level function

Aliasing restrictions: none

source
source

Calculate (S .*F)1 where S is the skew symmetric part of sbp.Q and F is a symmetric numerical flux function. eqn.res is updated with the result. Methods are available for curvilinear and non-curvilinear meshes

Inputs: mesh sbp eqn opts functor: the numerical flux function F, of type FluxType

source

This function calculates I_S2F.'(S . F( uk(I_S2F wk), uk(I_S2F wk)))1. This is similar to the other method but for the staggered grid algorithm.

Inputs: mesh_s: the mesh object for the solution grid mesh_f: the mesh object for the flux grid sbp_s: SBP operator for solution grid sbp_f: SBP operator for flux grid functor: the numerical flux functor (FluxType) used to compute F

Inputs/Outputs: eqn: the equation object (implicitly on the solution grid). eqn.res is updated with the result

Aliasing Restrictions: none (the meshes and the sbp operators could alias, in which case this algorithm reduces to the non-staggered version

source

Calculate (S .* F)1, where S is the skew-symmetric part of sbp.Q and F is a symmetric numerical flux function. eqn.res is updated with the result. This function is used for curvilinear meshes.

source

Calculates the volume integrals for the weak form, computing the Euler flux as needed, rather than using eqn.flux_parametric

Inputs: mesh sbp eqn: eqn.res is updated with the result opts

source

This function calculates the vorticity at every node of an element (because vorticity requires derivatives of the velocity, it would be awkward to compute it at a node level).

3D, conservative variables only

Inputs: params: a ParamType q: numDofPerNode x numNodesPerElement array of conservative variables at each node dxidx: the 3 x 3 x numNodesPerElement scaled mapping jacobian at the node jac: numNodesPerElement vector of the mapping jacobian determinant at each node

Input/Outputs: vorticity: a 3 x numNodesPerElement array containing the vorticity in the x, y, and z directions at each node, overwritten

Aliasing restrictions: from params: dxidx_element, velocities, velocity_deriv, velocity_deriv_xy

source

EulerEquationMod.getAuxVars

This function calculates any extra variables that are stored across the mesh using the conservative variables eqn.q. Currently only calculates pressure.

Thi is a mid level function

source

EulerEquationMod.getEulerFlux

This function calculates the Euler flux across the entire mesh by passing pieces of the eqn.q, eqn.aux_vars, eqn.f_xi and eqn.params to a low level function. The flux is calculated in the xi and eta directions, scaled (mulitiplied) by the mapping jacobian (so that when performing the integral we don't have to explictly divide by the jacobian, it just cancels out with the jacobian factor introduced here.

Calls writeFlux to do any requested output.

This is a mid level function

source

EulerEquationMod.getEulerFlux2

This function calcules the euler flux over the entire mesh directly (ie. does not call a low level function. This function is deprecated, although useful for benchmarking purposes. 2D only.

This is a mid level function

source

This function does the interpolation from the solution grid to the flux grid for a single element

Inputs:

  • params: a ParamType. The qs_el1 and q_el1 fields are overwritten.

  • qs_el: the conservative variables for the element on the solution grid

Inputs/Outputs:

  • aux_vars: array to be populated with the aux vars on the staggered grid

  • qf_el: conservative variables on the flux grid

source

EulerEquationMod.matVecA0

This function is the same as matVecA0inv, except it multiplies by A0 not A0inv. See its documention.

source

EulerEquationMod.matVecA0inv

This function takes a 3D array and multiplies it in place by the inv(A0) matrix (calculated at each node), inplace, (where A0 =dq/dv, where q are the conservative variables and v are some other variables), inplace. Methods are available for conservative and entropy variables.

For conservative variables, A0 is the identity matrix and this is a no-op

Inputs: mesh sbp eqn opts

Inputs/Outputs: res_arr: the array to multiply

Aliasing restrictions: none

source

EulerEquationMod.writeFlux

This function writes the real part of Euler flux to a file named Fxi.dat, space delimited, controlled by the input options 'writeflux', of type Bool.

This is a high level function.

source

Calculates ( 1/(2V) )*integral(rho * omega dot omega dV), where V is the volume of the mesh omega is the vorticity vector, and rho is the density. 3D, conservative variables only. This should work for CG and DG, but has only been tested for the latter

Inputs: mesh: an AbstractMesh sbp: an SBP operator eqn: an EulerData object opts: options dictionary q_arr: a 3D array of conservative variables

Outputs: val: the value of the integral (over the entire domain, not just the part owned by this process)

Aliasing restrictions: see calcVorticity

source

This function calculates ( 1/(2V) )integral(rho * v dot v dV), where V is the volume of the mesh, v is the velocity vector, and rho is the density. This is the total kinetic energy normalized by the volume of the domain Conservative variables only.

This function contains an MPI blocking collective operation. It must be called by all processes at the same time.

This function relies on the sequential numbering of dofs on the same node

Inputs: mesh: a mesh sbp: an SBP Operator eqn: an EulerData object opts: options dictionary q_vec: the vector of conservative variables for the entire mesh

Outputs: val: the value of the integral (over the entire domain, not just teh part owned by this procss)

Aliasing restrictions: none

source

This function calclates the time derivative of calcKineticEnergy.

The idea is to expand the left hand side of d rho*u/dt = res using the product rule and solve for du/dt, then use it to compute the integral. Inputs: mesh: a mesh sbp: a SBP operator eqn: an EulerData opts: options dictionary q_vec: vector of conservative variables for the entire mesh res_vec: residual vector (dq/dt) of entire mesh

Aliasing restrictions: none

source

Like contractResEntropyVars, but uses a special summation technique

source