MATLAB/Octave API#
This page documents the API of the MOLE MATLAB/Octave module. Functions are organized by category.
Differential Operators#
Gradient Operators#
- grad(k, m, dx, dc, nc)#
PURPOSE#
Returns a one-dimensional mimetic gradient operator depending on whether
SYNOPSIS#
function grad
DESCRIPTION#
or not the operator will contain a periodic boundary condition type a0 U + b0 dU/dn = g, (optional) dc : a0 (2x1 vector for left and right vertices, resp.) (optional) nc : b0 (2x1 vector for left and right vertices, resp.) k : Order of accuracy m : Number of cells dx : Step size
CROSS-REFERENCE INFORMATION#
This function calls:
gradPeriodic()
Returns a m by m+2 one-dimensional mimetic gradient operatorgradNonPeriodic()
Returns a m+1 by m+2 one-dimensional mimetic gradient operatorThis function is called by:
grad3DCurv()
Returns a 3D curvilinear mimetic gradientrobinBC()
Returns a m+2 by m+2 one-dimensional mimetic boundary operator thatmixedBC()
Constructs a 1D mimetic mixed boundary conditions operatorlap()
Returns a one-dimensional mimetic Laplacian operator depending on whethergradNonUniform()
Returns a m+1 by m+2 one-dimensional non-uniform mimetic gradientmimeticB()
Returns a m+2 by m+1 one-dimensional mimetic boundary operatorweightsP()
Returns the m+1 weights of PaddScalarBC1Dlhs()
This functions uses geometry and boundary type conditions to create
- grad2D(k, m, dx, n, dy, dc, nc)#
PURPOSE#
Returns a two-dimensional mimetic gradient operator depending on whether
SYNOPSIS#
function grad2D
DESCRIPTION#
or not the operator will contain a periodic boundary condition type a0 U + b0 dU/dn = g, (optional) dc : a0 (4x1 vector for left, right, bottom, top boundaries, resp.) (optional) nc : b0 (4x1 vector for left, right, bottom, top boundaries, resp.) k : Order of accuracy m : Number of cells dx : Step size n : Number of cells along y-axis dy : Step size along y-axis
CROSS-REFERENCE INFORMATION#
This function calls:
gradNonPeriodic2D()
Returns a two-dimensional mimetic gradient operatorgradPeriodic()
Returns a m by m+2 one-dimensional mimetic gradient operatorgradNonPeriodic()
Returns a m+1 by m+2 one-dimensional mimetic gradient operatorThis function is called by:
grad2DCurv()
Returns a 2D curvilinear mimetic gradientlap2D()
Returns a two-dimensional mimetic Laplacian operator depending on whether
- grad2DCurv(k, X, Y)#
PURPOSE#
Returns a 2D curvilinear mimetic gradient
SYNOPSIS#
function grad2DCurv
DESCRIPTION#
CROSS-REFERENCE INFORMATION#
This function calls:
grad2D()
Returns a two-dimensional mimetic gradient operator depending on whetherjacobian2D()
Returns:
- grad2DNonUniform(k, xticks, yticks)#
PURPOSE#
Returns a two-dimensional non-uniform mimetic gradient operator
SYNOPSIS#
function grad2DNonUniform
DESCRIPTION#
(including the boundaries!) k : Order of accuracy xticks : Centers' ticks (x-axis) yticks : Centers' ticks (y-axis)
CROSS-REFERENCE INFORMATION#
This function calls:
gradNonUniform()
Returns a m+1 by m+2 one-dimensional non-uniform mimetic gradient
- grad3D(k, m, dx, n, dy, o, dz, dc, nc)#
PURPOSE#
Returns a three-dimensional mimetic gradient operator depending on whether
SYNOPSIS#
function grad3D
DESCRIPTION#
or not the operator will contain a periodic boundary condition type a0 U + b0 dU/dn = g, (optional) dc : a0 (6x1 vector for left, right, bottom, top, front, back boundary types, resp.) (optional) nc : b0 (6x1 vector for left, right, bottom, top, front, back boundary types, resp.) k : Order of accuracy m : Number of cells dx : Step size n : Number of cells along y-axis dy : Step size along y-axis o : Number of cells along z-axis dz : Step size along z-axis
CROSS-REFERENCE INFORMATION#
This function calls:
gradNonPeriodic3D()
Returns a three-dimensional mimetic gradient operatorgradPeriodic()
Returns a m by m+2 one-dimensional mimetic gradient operatorgradNonPeriodic()
Returns a m+1 by m+2 one-dimensional mimetic gradient operatorThis function is called by:
lap3D()
Returns a three-dimensional mimetic Laplacian operator depending on whethergrad3DCurv()
Returns a 3D curvilinear mimetic gradient
- grad3DCurv(k, X, Y, Z)#
PURPOSE#
Returns a 3D curvilinear mimetic gradient
SYNOPSIS#
function grad3DCurv
DESCRIPTION#
CROSS-REFERENCE INFORMATION#
This function calls:
grad()
Returns a one-dimensional mimetic gradient operator depending on whethergrad3D()
Returns a three-dimensional mimetic gradient operator depending on whetherjacobian3D()
Returns:
- grad3DNonUniform(k, xticks, yticks, zticks)#
PURPOSE#
Returns a three-dimensional non-uniform mimetic gradient operator
SYNOPSIS#
function grad3DNonUniform
DESCRIPTION#
(including the boundaries!) k : Order of accuracy xticks : Centers' ticks (x-axis) yticks : Centers' ticks (y-axis) zticks : Centers' ticks (z-axis)
CROSS-REFERENCE INFORMATION#
This function calls:
gradNonUniform()
Returns a m+1 by m+2 one-dimensional non-uniform mimetic gradient
- gradNonUniform(k, ticks)#
PURPOSE#
Returns a m+1 by m+2 one-dimensional non-uniform mimetic gradient
SYNOPSIS#
function gradNonUniform
DESCRIPTION#
operator (including the boundaries!) k : Order of accuracy ticks : Centers' ticks e.g. [0 0.5 1 3 5 7 8 9 9.5 10]
CROSS-REFERENCE INFORMATION#
This function calls:
grad()
Returns a one-dimensional mimetic gradient operator depending on whetherThis function is called by:
grad2DNonUniform()
Returns a two-dimensional non-uniform mimetic gradient operatorgrad3DNonUniform()
Returns a three-dimensional non-uniform mimetic gradient operator
Divergence Operators#
- div(k, m, dx, dc, nc)#
PURPOSE#
Returns a one-dimensional mimetic divergence operator depending on whether
SYNOPSIS#
function div
DESCRIPTION#
or not the operator will contain a periodic boundary condition type a0 U + b0 dU/dn = g, (optional) dc : a0 (2x1 vector for left and right vertices, resp.) (optional) nc : b0 (2x1 vector for left and right vertices, resp.) k : Order of accuracy m : Number of cells dx : Step size
CROSS-REFERENCE INFORMATION#
This function calls:
divNonPeriodic()
Returns a m+2 by m+1 one-dimensional mimetic divergence operatordivPeriodic()
Returns a m+2 by m one-dimensional mimetic divergence operatorThis function is called by:
divNonUniform()
Returns a m+2 by m+1 one-dimensional non-uniform mimetic divergencelap()
Returns a one-dimensional mimetic Laplacian operator depending on whethermimeticB()
Returns a m+2 by m+1 one-dimensional mimetic boundary operatorweightsQ()
Returns the m+2 weights of Q
- div2D(k, m, dx, n, dy, dc, nc)#
PURPOSE#
Returns a two-dimensional mimetic divergence operator depending on whether
SYNOPSIS#
function div2D
DESCRIPTION#
or not the operator will contain a periodic boundary condition type a0 U + b0 dU/dn = g, (optional) dc : a0 (4x1 vector for left, right, bottom, top boundaries, resp.) (optional) nc : b0 (4x1 vector for left, right, bottom, top boundaries, resp.) k : Order of accuracy m : Number of cells dx : Step size n : Number of cells along y-axis dy : Step size along y-axis
CROSS-REFERENCE INFORMATION#
This function calls:
divNonPeriodic()
Returns a m+2 by m+1 one-dimensional mimetic divergence operatordivPeriodic()
Returns a m+2 by m one-dimensional mimetic divergence operatordivNonPeriodic2D()
Returns a two-dimensional mimetic divergence operatorThis function is called by:
curl2D()
Returns a two-dimensional mimetic curl operatorlap2D()
Returns a two-dimensional mimetic Laplacian operator depending on whether
- div2DCurv(k, X, Y)#
- div2DNonUniform(k, xticks, yticks)#
PURPOSE#
Returns a two-dimensional non-uniform mimetic divergence operator
SYNOPSIS#
function div2DNonUniform
DESCRIPTION#
k : Order of accuracy xticks : Edges' ticks (x-axis) yticks : Edges' ticks (y-axis)
CROSS-REFERENCE INFORMATION#
This function calls:
divNonUniform()
Returns a m+2 by m+1 one-dimensional non-uniform mimetic divergence
- div3D(k, m, dx, n, dy, o, dz, dc, nc)#
PURPOSE#
Returns a three-dimensional mimetic divergence operator depending on whether
SYNOPSIS#
function div3D
DESCRIPTION#
or not the operator will contain a periodic boundary condition type a0 U + b0 dU/dn = g, (optional) dc : a0 (6x1 vector for left, right, bottom, top, front, back boundary types, resp.) (optional) nc : b0 (6x1 vector for left, right, bottom, top, front, back boundary types, resp.) k : Order of accuracy m : Number of cells dx : Step size n : Number of cells along y-axis dy : Step size along y-axis o : Number of cells along z-axis dz : Step size along z-axis
CROSS-REFERENCE INFORMATION#
This function calls:
divNonPeriodic()
Returns a m+2 by m+1 one-dimensional mimetic divergence operatordivNonPeriodic3D()
Returns a three-dimensional mimetic divergence operatordivPeriodic()
Returns a m+2 by m one-dimensional mimetic divergence operatorThis function is called by:
lap3D()
Returns a three-dimensional mimetic Laplacian operator depending on whether
- div3DCurv(k, X, Y, Z)#
- div3DNonUniform(k, xticks, yticks, zticks)#
PURPOSE#
Returns a three-dimensional non-uniform mimetic divergence operator
SYNOPSIS#
function div3DNonUniform
DESCRIPTION#
k : Order of accuracy xticks : Edges' ticks (x-axis) yticks : Edges' ticks (y-axis) zticks : Edges' ticks (z-axis)
CROSS-REFERENCE INFORMATION#
This function calls:
divNonUniform()
Returns a m+2 by m+1 one-dimensional non-uniform mimetic divergence
- divNonUniform(k, ticks)#
PURPOSE#
Returns a m+2 by m+1 one-dimensional non-uniform mimetic divergence
SYNOPSIS#
function divNonUniform
DESCRIPTION#
operator k : Order of accuracy ticks : Edges' ticks e.g. [0 0.1 0.15 0.2 0.3 0.4 0.45]
CROSS-REFERENCE INFORMATION#
This function calls:
div()
Returns a one-dimensional mimetic divergence operator depending on whetherThis function is called by:
div2DNonUniform()
Returns a two-dimensional non-uniform mimetic divergence operatordiv3DNonUniform()
Returns a three-dimensional non-uniform mimetic divergence operator
Curl Operators#
- curl2D(k, m, dx, n, dy, west, east, south, north, U, V)#
PURPOSE#
Returns a two-dimensional mimetic curl operator
SYNOPSIS#
function curl2D
DESCRIPTION#
west, east, south, north : west, east, south, north limits U(X,Y) must be defined as function handle V(X,Y) must be defined as function handle k : Order of accuracy m : Number of cells along x-axis dx : Step size along x-axis n : Number of cells along y-axis dy : Step size along y-axis U : Vector space function acting on x-direction V : Vector space function acting on y-direction
CROSS-REFERENCE INFORMATION#
This function calls:
div2D()
Returns a two-dimensional mimetic divergence operator depending on whether
Laplacian Operators#
- lap(k, m, dx, dc, nc)#
PURPOSE#
Returns a one-dimensional mimetic Laplacian operator depending on whether
SYNOPSIS#
function lap
DESCRIPTION#
or not the operator will contain a periodic boundary condition type a0 U + b0 dU/dn = g, (optional) dc : a0 (2x1 vector for left and right vertices, resp.) (optional) nc : b0 (2x1 vector for left and right vertices, resp.) k : Order of accuracy m : Number of cells dx : Step size
CROSS-REFERENCE INFORMATION#
This function calls:
div()
Returns a one-dimensional mimetic divergence operator depending on whethergrad()
Returns a one-dimensional mimetic gradient operator depending on whetherlapNonPeriodic()
Returns a m+2 by m+2 one-dimensional mimetic laplacian operator
- lap2D(k, m, dx, n, dy, dc, nc)#
PURPOSE#
Returns a two-dimensional mimetic Laplacian operator depending on whether
SYNOPSIS#
function lap2D
DESCRIPTION#
or not the operator will contain a periodic boundary condition type a0 U + b0 dU/dn = g, (optional) dc : a0 (4x1 vector for left, right, bottom, top boundaries, resp.) (optional) nc : b0 (4x1 vector for left, right, bottom, top boundaries, resp.) k : Order of accuracy m : Number of cells along x-axis dx : Step size along x-axis n : Number of cells along y-axis dy : Step size along y-axis
CROSS-REFERENCE INFORMATION#
This function calls:
div2D()
Returns a two-dimensional mimetic divergence operator depending on whethergrad2D()
Returns a two-dimensional mimetic gradient operator depending on whetherlapNonPeriodic2D()
Returns a two-dimensional mimetic laplacian operator
- lap3D(k, m, dx, n, dy, o, dz, dc, nc)#
PURPOSE#
Returns a three-dimensional mimetic Laplacian operator depending on whether
SYNOPSIS#
function lap3D
DESCRIPTION#
or not the operator will contain a periodic boundary condition type a0 U + b0 dU/dn = g, (optional) dc : a0 (6x1 vector for left, right, bottom, top, front, back boundary types, resp.) (optional) nc : b0 (6x1 vector for left, right, bottom, top, front, back boundary types, resp.) k : Order of accuracy m : Number of cells along x-axis dx : Step size along x-axis n : Number of cells along y-axis dy : Step size along y-axis o : Number of cells along z-axis dz : Step size along z-axis
CROSS-REFERENCE INFORMATION#
This function calls:
grad3D()
Returns a three-dimensional mimetic gradient operator depending on whetherlapNonPeriodic3D()
Returns a three-dimensional mimetic laplacian operatordiv3D()
Returns a three-dimensional mimetic divergence operator depending on whether
Interpolation Functions#
Node to Center Interpolation#
- interpolNodesToCenters1D(k, m)#
PURPOSE#
interpolation operator from nodal coordinates to staggered centers
SYNOPSIS#
function interpolNodesToCenters1D
DESCRIPTION#
m is the number of cells in the logical x-axis nodal logical coordinates are [1:1:m] centers logical coordinates [1,1.5:m-0.5,m]
CROSS-REFERENCE INFORMATION#
This function calls:
interpolFacesToCentersG1D()
1D interpolation from faces to centers
- interpolNodesToCenters2D(k, m, n)#
PURPOSE#
interpolation operator from nodal coordinates to staggered centers
SYNOPSIS#
function interpolNodesToCenters2D
DESCRIPTION#
m, n, are the number of cells in the logical x-, y- axes nodal logical coordinates are [1:1:m]x[1:1:n] centers logical coordinates [1,1.5:m-0.5,m]x[1,1.5:n-0.5,n]
CROSS-REFERENCE INFORMATION#
This function calls:
interpolFacesToCentersG1D()
1D interpolation from faces to centers
- interpolNodesToCenters3D(k, m, n, o)#
PURPOSE#
interpolation operator from nodal coordinates to staggered centers
SYNOPSIS#
function interpolNodesToCenters3D
DESCRIPTION#
m, n, o, are the number of cells in the logical x-, y-, z- axes nodal logical coordinates are [1:1:m]x[1:1:n]x[1:1:o] centers logical coordinates [1,1.5:m-0.5,m]x[1,1.5:n-0.5,n]x[1,1.5:o-0.5,o]
CROSS-REFERENCE INFORMATION#
This function calls:
interpolFacesToCentersG1D()
1D interpolation from faces to centers
Center to Node Interpolation#
- interpolCentersToNodes1D(k, m)#
PURPOSE#
interpolation operator from nodal coordinates to staggered centers
SYNOPSIS#
function interpolCentersToNodes1D
DESCRIPTION#
m is the number of cells in the logical x-axis nodal logical coordinates are [1:1:m] centers logical coordinates [1,1.5:m-0.5,m]
CROSS-REFERENCE INFORMATION#
This function calls:
interpolCentersToFacesD1D()
1D interpolation from centers to faces.
- interpolCentersToNodes2D(k, m, n)#
PURPOSE#
interpolation operator from staggered to nodes
SYNOPSIS#
function interpolCentersToNodes2D
DESCRIPTION#
m, n, are the number of cells in the logical x-, y- axes nodal logical coordinates are [1:1:m]x[1:1:n] centers logical coordinates [1,1.5:m-0.5,m]x[1,1.5:n-0.5,n]
CROSS-REFERENCE INFORMATION#
This function calls:
interpolCentersToFacesD1D()
1D interpolation from centers to faces.
- interpolCentersToNodes3D(k, m, n, o)#
PURPOSE#
interpolation operator from staggered to nodes
SYNOPSIS#
function interpolCentersToNodes3D
DESCRIPTION#
m, n, o, are the number of cells in the logical x-, y-, z- axes nodal logical coordinates are [1:1:m]x[1:1:n]x[1:1:o] centers logical coordinates [1,1.5:m-0.5,m]x[1,1.5:n-0.5,n]x[1,1.5:o-0.5,o]
CROSS-REFERENCE INFORMATION#
This function calls:
interpolCentersToFacesD1D()
1D interpolation from centers to faces.
Face Interpolation#
- interpolFacesToCentersG1D(k, m)#
PURPOSE#
1D interpolation from faces to centers
SYNOPSIS#
function interpolFacesToCentersG1D
DESCRIPTION#
centers logical coordinates [1,1.5:m-0.5,m] m is the number of cells in the logical x-axis
CROSS-REFERENCE INFORMATION#
This function is called by:
interpolNodesToCenters3D()
interpolation operator from nodal coordinates to staggered centersinterpolFacesToCentersG2D()
2D interpolation from faces to centersinterpolFacesToCentersG3D()
3D interpolation from faces to centersinterpolNodesToCenters1D()
interpolation operator from nodal coordinates to staggered centersinterpolNodesToCenters2D()
interpolation operator from nodal coordinates to staggered centers
- interpolFacesToCentersG2D(k, m, n)#
PURPOSE#
2D interpolation from faces to centers
SYNOPSIS#
function interpolFacesToCentersG2D
DESCRIPTION#
centers logical coordinates [1,1.5:m-0.5,m]x[1,1.5:n-0.5,n] m, n, are the number of cells in the logical x- and y- axes
CROSS-REFERENCE INFORMATION#
This function calls:
interpolFacesToCentersG1D()
1D interpolation from faces to centers
- interpolFacesToCentersG3D(k, m, n, o)#
PURPOSE#
3D interpolation from faces to centers
SYNOPSIS#
function interpolFacesToCentersG3D
DESCRIPTION#
centers logical coordinates [1,1.5:m-0.5,m]x[1,1.5:n-0.5,n]x[1,1.5:o-0.5,o] m, n, o, are the number of cells in the logical x-, y-, z- axes
CROSS-REFERENCE INFORMATION#
This function calls:
interpolFacesToCentersG1D()
1D interpolation from faces to centers
General Interpolation#
- interpol(m, c)#
PURPOSE#
Returns a m+1 by m+2 one-dimensional interpolator of 2nd-order
SYNOPSIS#
function interpol
DESCRIPTION#
m : Number of cells c : Left interpolation coeff.
CROSS-REFERENCE INFORMATION#
This function is called by:
interpol2D()
Returns a two-dimensional interpolator of 2nd-orderinterpol3D()
Returns a three-dimensional interpolator of 2nd-order
- interpol2D(m, n, c1, c2)#
PURPOSE#
Returns a two-dimensional interpolator of 2nd-order
SYNOPSIS#
function interpol2D
DESCRIPTION#
m : Number of cells along x-axis n : Number of cells along y-axis c1 : Left interpolation coeff. c2 : Bottom interpolation coeff.
CROSS-REFERENCE INFORMATION#
This function calls:
interpol()
Returns a m+1 by m+2 one-dimensional interpolator of 2nd-order
- interpol3D(m, n, o, c1, c2, c3)#
PURPOSE#
Returns a three-dimensional interpolator of 2nd-order
SYNOPSIS#
function interpol3D
DESCRIPTION#
m : Number of cells along x-axis n : Number of cells along y-axis o : Number of cells along z-axis c1 : Left interpolation coeff. c2 : Bottom interpolation coeff. c3 : Front interpolation coeff.
CROSS-REFERENCE INFORMATION#
This function calls:
interpol()
Returns a m+1 by m+2 one-dimensional interpolator of 2nd-order
- interpolD(m, c)#
PURPOSE#
Returns a m+2 by m+1 one-dimensional interpolator of 2nd-order
SYNOPSIS#
function interpolD
DESCRIPTION#
m : Number of cells c : Left interpolation coeff.
CROSS-REFERENCE INFORMATION#
This function is called by:
interpolD2D()
Returns a two-dimensional interpolator of 2nd-orderinterpolD3D()
Returns a three-dimensional interpolator of 2nd-order
- interpolD2D(m, n, c1, c2)#
PURPOSE#
Returns a two-dimensional interpolator of 2nd-order
SYNOPSIS#
function interpolD2D
DESCRIPTION#
m : Number of cells along x-axis n : Number of cells along y-axis c1 : Left interpolation coeff. c2 : Bottom interpolation coeff.
CROSS-REFERENCE INFORMATION#
This function calls:
interpolD()
Returns a m+2 by m+1 one-dimensional interpolator of 2nd-order
- interpolD3D(m, n, o, c1, c2, c3)#
PURPOSE#
Returns a three-dimensional interpolator of 2nd-order
SYNOPSIS#
function interpolD3D
DESCRIPTION#
m : Number of cells along x-axis n : Number of cells along y-axis o : Number of cells along z-axis c1 : Left interpolation coeff. c2 : Bottom interpolation coeff. c3 : Front interpolation coeff.
CROSS-REFERENCE INFORMATION#
This function calls:
interpolD()
Returns a m+2 by m+1 one-dimensional interpolator of 2nd-order
Boundary Conditions#
General Boundary Conditions#
- addScalarBC1D(A, b, k, m, dx, dc, nc, v)#
PURPOSE#
Separates cases non-periodic and periodic for dealing with boundary data
SYNOPSIS#
function addScalarBC1D
DESCRIPTION#
output input A0 : Linear operator with boundary conditions added b0 : Right hand side with boundary conditions added A : Linear operator without boundary conditions added b : Right hand side without boundary conditions added k : Order of accuracy m : Number of cells dx : Step size dc : a0 (2x1 vector for left and right vertices, resp.) nc : b0 (2x1 vector for left and right vertices, resp.) v : g (2x1 vector for left and right vertices, resp.)
CROSS-REFERENCE INFORMATION#
This function calls:
addScalarBC1Drhs()
This function uses the non-periodic boundary condition type of each vertexaddScalarBC1Dlhs()
This functions uses geometry and boundary type conditions to create
- addScalarBC1Dlhs(k, m, dx, dc, nc)#
PURPOSE#
This functions uses geometry and boundary type conditions to create
SYNOPSIS#
function addScalarBC1Dlhs
DESCRIPTION#
modifications of matrix A associated to each of the boundary faces. output input Al : modification of matrix A due to left boundary condition Ar : modification of matrix A due to right boundary condition k : Order of accuracy m : Number of cells dx : Step size dc : Dirichlet coefficient (2x1 vector for left and right vertices, resp.) nc : Neumann coefficient (2x1 vector for left and right vertices, resp.)
CROSS-REFERENCE INFORMATION#
This function calls:
grad()
Returns a one-dimensional mimetic gradient operator depending on whetherThis function is called by:
addScalarBC1D()
Separates cases non-periodic and periodic for dealing with boundary dataaddScalarBC2Dlhs()
This functions uses geometry and boundary type conditions to createaddScalarBC3Dlhs()
This functions uses geometry and boundary type conditions to create
- addScalarBC1Drhs(b, v, vec)#
PURPOSE#
This function uses the non-periodic boundary condition type of each vertex
SYNOPSIS#
function addScalarBC1Drhs
DESCRIPTION#
and the rhs b values associated to left, and right vertices to modify the rhs vector b. output input b : Right hand side with boundary conditions added b : Right hand side without boundary conditions added v : value (2x1 vector for left and right vertices, resp.) vec : vector with indices of rhs associated to bc
CROSS-REFERENCE INFORMATION#
This function is called by:
addScalarBC1D()
Separates cases non-periodic and periodic for dealing with boundary data
- addScalarBC2D(A, b, k, m, dx, n, dy, dc, nc, v)#
PURPOSE#
This function assumes that the unknown u, which represents the discrete
SYNOPSIS#
function addScalarBC2D
DESCRIPTION#
solution the continuous second-order 2D PDE operator L U = f, with continuous boundary condition a0 U + b0 dU/dn = g, are given at the 2D cell centers and boundary face centers. Furthermore, all discrete calculations are performed at the 2D cell centers and boundary face centers. The function receives as input quantities associated to the discrete analog of the continuous problem given by the squared linear system A u = b where A is the discrete analog of L and b is the discrete analog of g, both constructed by the user without boundary conditions. The function output is the modified square linear system A u = b where both A and b include boundary condition information. The boundary condition is always one of the following forms: For Dirichlet set: a0 not equal zero and b0 = 0. For Neumann set : a0 = 0 and b0 not equal zero. For Robin set : both a0 and b0 not equal zero. For Periodic set : both a0 = 0 and b0 = 0. For periodic bc, it is assumed that not only u but also du/dn are the same in both extremes of the domain since a second-order PDE is assumed. Periodic boundary conditions can be applied along some axes and non-periodic to some others. For consistence with the way boundary operators are calculated to avoid overwriting of the values v, the left and right boundary conditions are assumed to be column vectors of (m+2)*n components, and the bottom and top faces are assumed to be vectors of (m+2)*(n+2) components. The order of these components is as follows: For left and right edges, the ordering is the one given by columns vectors where x increases. For bottom and top faces, the ordering is the one given by columns vectors where y increases. The code assumes the following assertions: assert(k >= 2, 'k >= 2'); assert(mod(k, 2) == 0, 'k % 2 = 0'); assert(m >= 2*k+1, ['m >= ' num2str(2*k+1) ' for k = ' num2str(k)]); output input A : Linear operator with boundary conditions added b : Right hand side with boundary conditions added A : Linear operator without boundary conditions added b : Right hand side without boundary conditions added k : Order of accuracy m : Number of horizontal cells dx : Step size horizontal cells n : Number of vertical cells dy : Step size of vertical cells dc : a0 (4x1 vector for left, right, bottom, top boundaries, resp.) nc : b0 (4x1 vector for left, right, bottom, top boundaries, resp.) v : g (4x1 vector of arrays for left, right, bottom, top boundaries, resp.)
CROSS-REFERENCE INFORMATION#
This function calls:
addScalarBC2Dlhs()
This functions uses geometry and boundary type conditions to createaddScalarBC2Drhs()
function b = addBC2Drhs(b, m, n, dc, nc, v, vec)
- addScalarBC2Dlhs(k, m, dx, n, dy, dc, nc)#
PURPOSE#
This functions uses geometry and boundary type conditions to create
SYNOPSIS#
function addScalarBC2Dlhs
DESCRIPTION#
modifications of matrix A associated to each of the boundary edges. output input Abcl : Matrix coefficients associated to boundary conditions for left edge Abcr : Matrix coefficients associated to boundary conditions for right edge Abcb : Matrix coefficients associated to boundary conditions for bottom edge Abct : Matrix coefficients associated to boundary conditions for top edge k : Order of accuracy m : Number of the horizontal cells dx : Step size n : Number of the vertical cells dy : Horizontal cell size dc : a0 (4x1 vector for left, right, bottom, top boundaries, resp.) nc : b0 (4x1 vector for left, right, bottom, top boundaries resp.)
CROSS-REFERENCE INFORMATION#
This function calls:
addScalarBC1Dlhs()
This functions uses geometry and boundary type conditions to createThis function is called by:
addScalarBC2D()
This function assumes that the unknown u, which represents the discrete
- addScalarBC2Drhs(b, dc, nc, v, rl, rr, rb, rt)#
PURPOSE#
function b = addBC2Drhs(b, m, n, dc, nc, v, vec)
SYNOPSIS#
function b = addBC2Drhs(b, m, n, dc, nc, v, vec)
DESCRIPTION#
This function uses the boundary condition type of each face and the rhs b indices and values associated to left, right, bottom, top, front, back faces to modify the rhs vector b. output input b : Right hand side with boundary conditions added b : Right hand side without boundary conditions added dc : a0 (4x1 vector for left, right, bottom, top boundary types, resp.) nc : b0 (4x1 vector for left, right, bottom, top boundary types, resp.) v : g (4x1 vector of arrays for left, right, bottom, top boundaries, resp.) rl : indices of rhs left indices rr : indices of rhs right indices rb : indices of rhs bottom indices rt : indices of rhs top indices vec : vector with indices of rhs associated to bc
CROSS-REFERENCE INFORMATION#
This function is called by:
addScalarBC2D()
This function assumes that the unknown u, which represents the discrete
- addScalarBC3D(A, b, k, m, dx, n, dy, o, dz, dc, nc, v)#
PURPOSE#
This function assumes that the unknown u, which represents the discrete
SYNOPSIS#
function addScalarBC3D
DESCRIPTION#
solution the continuous second-order 3D PDE operator L U = f, with continuous boundary condition a0 U + b0 dU/dn = g, are given at the 3D cell centers and boundary face centers. Furthermore, all discrete calculations are performed at the 3D cell centers and boundary face centers. The function receives as input quantities associated to the discrete analog of the continuous problem given by the squared linear system A u = b where A is the discrete analog of L and b is the discrete analog of g, both constructed by the user without boundary conditions. The function output is the modified square linear system A u = b where both A and b include boundary condition information. The boundary condition is always one of the following forms: For Dirichlet set: a0 not equal zero and b0 = 0. For Neumann set : a0 = 0 and b0 not equal zero. For Robin set : both a0 and b0 not equal zero. For Periodic set : both a0 = 0 and b0 = 0. For periodic bc, it is assumed that not only u but also du/dn are the same in both extremes of the domain since a second-order PDE is assumed. Periodic boundary conditions can be applied along some axes and non-periodic to some others. For consistence with the way boundary operators are calculated to avoid overwriting of the values v, the left and right boundary conditions are assumed to be column vectors of n*o components, the bottom and top boundary conditions are assumed to be column vectors of (m+2)*o components, and the front and back faces are assumed to be vectors of (m+2)*(n+2) components. The order of these components is as follows: For left and right faces, the ordering is the one by columns of the matrix where y increase along rows, and z increase along columns. For bottom and top faces, the ordering is the one by columns of the matrix where x increase along rows, and z increase along columns. For front and back faces, the ordering is the one by columns of the matrix where x increase along rows, and y increase along columns. The code assumes the following assertions: assert(k >= 2, 'k >= 2'); assert(mod(k, 2) == 0, 'k % 2 = 0'); assert(m >= 2*k+1, ['m >= ' num2str(2*k+1) ' for k = ' num2str(k)]); output input A : Linear operator with boundary conditions added b : Right hand side with boundary conditions added A : Linear operator without boundary conditions added b : Right hand side without boundary conditions added k : Order of accuracy m : Number of horizontal cells dx : Step size of horizontal n : Number of vertical cells dy : Step size of vertical cells o : Number of depth cells dz : Step size of depth cells dc : a0 (6x1 vector for left, right, bottom, top, front, back boundary types, resp.) nc : b0 (6x1 vector for left, right, bottom, top, front, back boundary types, resp.) v : g (6x1 vector of arrays for left, right, bottom, top, front, back boundaries, resp.)
CROSS-REFERENCE INFORMATION#
This function calls:
addScalarBC3Drhs()
This function uses the boundary condition type of each face and the rhs baddScalarBC3Dlhs()
This functions uses geometry and boundary type conditions to create
- addScalarBC3Dlhs(k, m, dx, n, dy, o, dz, dc, nc)#
PURPOSE#
This functions uses geometry and boundary type conditions to create
SYNOPSIS#
function addScalarBC3Dlhs
DESCRIPTION#
modifications of matrix A associated to each of the boundary faces. output input Abcl : Matrix coefficients associated to boundary conditions for left face Abcr : Matrix coefficients associated to boundary conditions for right face Abcb : Matrix coefficients associated to boundary conditions for bottom face Abct : Matrix coefficients associated to boundary conditions for top face Abcf : Matrix coefficients associated to boundary conditions for front face Abcz : Matrix coefficients associated to boundary conditions for back face k : Order of accuracy m : Number of the horizontal cells dx : Step size n : Number of the vertical cells dy : Horizonttal cell size o : Number of the depth cells dz : Depth cell size dc : a0 (6x1 vector for left, right, bottom, top, front, back boundary types, resp.) nc : b0 (6x1 vector for left, right, bottom, top, front, back boundary types, resp.)
CROSS-REFERENCE INFORMATION#
This function calls:
addScalarBC1Dlhs()
This functions uses geometry and boundary type conditions to createThis function is called by:
addScalarBC3D()
This function assumes that the unknown u, which represents the discrete
- addScalarBC3Drhs(b, dc, nc, v, rl, rr, rb, rt, rf, rz)#
PURPOSE#
This function uses the boundary condition type of each face and the rhs b
SYNOPSIS#
function addScalarBC3Drhs
DESCRIPTION#
indices and values associated to left, right, bottom, top, front, back faces to modify the rhs vector b. output input b : Right hand side with boundary conditions added b : Right hand side without boundary conditions added dc : a0 (6x1 vector for left, right, bottom, top, front, back boundary types, resp.) nc : b0 (6x1 vector for left, right, bottom, top, front, back boundary types, resp.) v : g (6x1 vector of arrays for left, right, bottom, top, front, back boundaries, resp.) rl : indices of rhs left indices rr : indices of rhs right indices rb : indices of rhs bottom indices rt : indices of rhs top indices rf : indices of rhs front indices rz : indices of rhs back indices
CROSS-REFERENCE INFORMATION#
This function is called by:
addScalarBC3D()
This function assumes that the unknown u, which represents the discrete
Neumann Boundary Conditions#
Robin Boundary Conditions#
- robinBC(k, m, dx, a, b)#
PURPOSE#
Returns a m+2 by m+2 one-dimensional mimetic boundary operator that
SYNOPSIS#
function robinBC
DESCRIPTION#
imposes a boundary condition of Robin's type k : Order of accuracy m : Number of cells dx : Step size a : Dirichlet Coefficient b : Neumann Coefficient
CROSS-REFERENCE INFORMATION#
This function calls:
grad()
Returns a one-dimensional mimetic gradient operator depending on whetherThis function is called by:
robinBC2D()
Returns a two-dimensional mimetic boundary operator thatrobinBC3D()
Returns a three-dimensional mimetic boundary operator that
- robinBC2D(k, m, dx, n, dy, a, b)#
PURPOSE#
Returns a two-dimensional mimetic boundary operator that
SYNOPSIS#
function robinBC2D
DESCRIPTION#
imposes a boundary condition of Robin's type k : Order of accuracy m : Number of cells along x-axis dx : Step size along x-axis n : Number of cells along y-axis dy : Step size along y-axis a : Dirichlet Coefficient b : Neumann Coefficient
CROSS-REFERENCE INFORMATION#
This function calls:
robinBC()
Returns a m+2 by m+2 one-dimensional mimetic boundary operator that
- robinBC3D(k, m, dx, n, dy, o, dz, a, b)#
PURPOSE#
Returns a three-dimensional mimetic boundary operator that
SYNOPSIS#
function robinBC3D
DESCRIPTION#
imposes a boundary condition of Robin's type k : Order of accuracy m : Number of cells along x-axis dx : Step size along x-axis n : Number of cells along y-axis dy : Step size along y-axis o : Number of cells along z-axis dz : Step size along z-axis a : Dirichlet Coefficient b : Neumann Coefficient
CROSS-REFERENCE INFORMATION#
This function calls:
robinBC()
Returns a m+2 by m+2 one-dimensional mimetic boundary operator that
Mixed Boundary Conditions#
- mixedBC(k, m, dx, left, coeffs_left, right, coeffs_right)#
PURPOSE#
Constructs a 1D mimetic mixed boundary conditions operator
SYNOPSIS#
function mixedBC
DESCRIPTION#
k : Order of accuracy m : Number of cells dx : Step size left : Type of boundary condition at the left boundary ('Dirichlet', 'Neumann', 'Robin') coeffs_left : Coefficients for the left boundary condition (a, b for Robin, otherwise coeff. for Dirichlet or Neumann) right : Type of boundary condition at the right boundary ('Dirichlet', 'Neumann', 'Robin') coeffs_right : Coefficients for the right boundary condition (a, b for Robin, otherwise coeff. for Dirichlet or Neumann)
CROSS-REFERENCE INFORMATION#
This function calls:
grad()
Returns a one-dimensional mimetic gradient operator depending on whetherThis function is called by:
mixedBC3D()
Constructs a 3D mimetic mixed boundary conditions operatormixedBC2D()
Constructs a 2D mimetic mixed boundary conditions operator
- mixedBC2D(k, m, dx, n, dy, left, coeffs_left, right, coeffs_right, bottom, coeffs_bottom, top, coeffs_top)#
PURPOSE#
Constructs a 2D mimetic mixed boundary conditions operator
SYNOPSIS#
function mixedBC2D
DESCRIPTION#
k : Order of accuracy m : Number of cells in x-direction dx : Step size in x-direction n : Number of cells in y-direction dy : Step size in y-direction left : Type of boundary condition at the left boundary ('Dirichlet', 'Neumann', 'Robin') coeffs_left : Coefficients for the left boundary condition (a, b for Robin, otherwise coeff. for Dirichlet or Neumann) right : Type of boundary condition at the right boundary ('Dirichlet', 'Neumann', 'Robin') coeffs_right : Coefficients for the right boundary condition (a, b for Robin, otherwise coeff. for Dirichlet or Neumann) bottom : Type of boundary condition at the bottom boundary ('Dirichlet', 'Neumann', 'Robin') coeffs_bottom : Coefficients for the bottom boundary condition (a, b for Robin, otherwise coeff. for Dirichlet or Neumann) top : Type of boundary condition at the top boundary ('Dirichlet', 'Neumann', 'Robin') coeffs_top : Coefficients for the top boundary condition (a, b for Robin, otherwise coeff. for Dirichlet or Neumann)
CROSS-REFERENCE INFORMATION#
This function calls:
mixedBC()
Constructs a 1D mimetic mixed boundary conditions operator
- mixedBC3D(k, m, dx, n, dy, o, dz, left, coeffs_left, right, coeffs_right, bottom, coeffs_bottom, top, coeffs_top, front, coeffs_front, back, coeffs_back)#
PURPOSE#
Constructs a 3D mimetic mixed boundary conditions operator
SYNOPSIS#
function mixedBC3D
DESCRIPTION#
k : Order of accuracy m : Number of cells in x-direction dx : Step size in x-direction n : Number of cells in y-direction dy : Step size in y-direction o : Number of cells in z-direction dz : Step size in z-direction left : Type of boundary condition at the left boundary ('Dirichlet', 'Neumann', 'Robin') coeffs_left : Coefficients for the left boundary condition (a, b for Robin, otherwise coeff. for Dirichlet or Neumann) right : Type of boundary condition at the right boundary ('Dirichlet', 'Neumann', 'Robin') coeffs_right : Coefficients for the right boundary condition (a, b for Robin, otherwise coeff. for Dirichlet or Neumann) bottom : Type of boundary condition at the bottom boundary ('Dirichlet', 'Neumann', 'Robin') coeffs_bottom : Coefficients for the bottom boundary condition (a, b for Robin, otherwise coeff. for Dirichlet or Neumann) top : Type of boundary condition at the top boundary ('Dirichlet', 'Neumann', 'Robin') coeffs_top : Coefficients for the top boundary condition (a, b for Robin, otherwise coeff. for Dirichlet or Neumann) front : Type of boundary condition at the front boundary ('Dirichlet', 'Neumann', 'Robin') coeffs_front : Coefficients for the front boundary condition (a, b for Robin, otherwise coeff. for Dirichlet or Neumann) back : Type of boundary condition at the back boundary ('Dirichlet', 'Neumann', 'Robin') coeffs_back : Coefficients for the back boundary condition (a, b for Robin, otherwise coeff. for Dirichlet or Neumann)
CROSS-REFERENCE INFORMATION#
This function calls:
mixedBC()
Constructs a 1D mimetic mixed boundary conditions operator
Grid Generation and Transformation#
Grid Generation#
- gridGen(method, grid_name, m, n, plot_grid, varargin)#
PURPOSE#
Returns X and Y which are both m by n matrices that contains the physical
SYNOPSIS#
function gridGen
DESCRIPTION#
coordinates grid_name : String with the name of the grid folder m : Number of nodes along the horizontal axis n : Number of nodes along the vertical axis plot_grid : If true -> plot the grid varargin : Maximum number of iterations (Required for TTM)
CROSS-REFERENCE INFORMATION#
This function calls:
ttm()
https://www.sciencedirect.com/science/article/pii/0022247X78902172?via%3Dihubtfi()
https://en.wikipedia.org/wiki/Transfinite_interpolation
- tfi(grid_name, m, n, plot_grid)#
PURPOSE#
Returns X and Y which are both m by n matrices that contains the physical
SYNOPSIS#
function tfi
DESCRIPTION#
coordinates grid_name : String with the name of the grid folder m : Number of nodes along the horizontal axis n : Number of nodes along the vertical axis plot_grid : If defined -> grid will be plotted
CROSS-REFERENCE INFORMATION#
This function is called by:
gridGen()
Returns X and Y which are both m by n matrices that contains the physical
Jacobian Calculation#
- jacobian2D(k, X, Y)#
PURPOSE#
- returns:
Determinant of the Jacobian (XeYn - XnYe)
SYNOPSIS#
function jacobian2D
DESCRIPTION#
Xe : dx/de metric Xn : dx/dn metric Ye : dy/de metric Yn : dy/dn metric k : Order of accuracy X : x-coordinates (physical) of meshgrid Y : y-coordinates (physical) of meshgrid
CROSS-REFERENCE INFORMATION#
This function calls:
nodal2D()
Returns a two-dimensional operator that approximates the first-orderThis function is called by:
grad2DCurv()
Returns a 2D curvilinear mimetic gradient
- jacobian3D(k, X, Y, Z)#
PURPOSE#
- returns:
Determinant of the Jacobian
SYNOPSIS#
function jacobian3D
DESCRIPTION#
Xe : dx/de metric Xn : dx/dn metric Xc : dx/dc metric Ye : dy/de metric Yn : dy/dn metric Yc : dy/dc metric Ze : dz/de metric Zn : dz/dn metric Zc : dz/dc metric k : Order of accuracy X : x-coordinates (physical) of meshgrid Y : y-coordinates (physical) of meshgrid Z : z-coordinates (physical) of meshgrid
CROSS-REFERENCE INFORMATION#
This function calls:
nodal3D()
Returns a three-dimensional operator that approximates the first-orderThis function is called by:
grad3DCurv()
Returns a 3D curvilinear mimetic gradient
Nodal Operators#
- nodal(k, m, dx)#
PURPOSE#
Returns a m+1 by m+1 one-dimensional operator that approximates the
SYNOPSIS#
function nodal
DESCRIPTION#
first-order derivatives on a uniform nodal grid k : Order of accuracy m : Number of nodes dx : Step size
CROSS-REFERENCE INFORMATION#
This function is called by:
nodal2D()
Returns a two-dimensional operator that approximates the first-ordernodal3D()
Returns a three-dimensional operator that approximates the first-order
- nodal2D(k, m, dx, n, dy)#
PURPOSE#
Returns a two-dimensional operator that approximates the first-order
SYNOPSIS#
function nodal2D
DESCRIPTION#
derivatives on a uniform nodal grid k : Order of accuracy m : Number of nodes along x-axis dx : Step size along x-axis n : Number of nodes along y-axis dy : Step size along y-axis
CROSS-REFERENCE INFORMATION#
This function calls:
nodal()
Returns a m+1 by m+1 one-dimensional operator that approximates theThis function is called by:
jacobian2D()
Returns:
- nodal2DCurv(k, X, Y)#
- nodal3D(k, m, dx, n, dy, o, dz)#
PURPOSE#
Returns a three-dimensional operator that approximates the first-order
SYNOPSIS#
function nodal3D
DESCRIPTION#
derivatives on a uniform nodal grid k : Order of accuracy m : Number of nodes along x-axis dx : Step size along x-axis n : Number of nodes along y-axis dy : Step size along y-axis o : Number of nodes along z-axis dz : Step size along z-axis
CROSS-REFERENCE INFORMATION#
This function calls:
nodal()
Returns a m+1 by m+1 one-dimensional operator that approximates theThis function is called by:
jacobian3D()
Returns:
- nodal3DCurv(k, X, Y, Z)#
- sidedNodal(m, dx, type)#
PURPOSE#
Returns a m+1 by m+1 one-dimensional sided approximation for uniformly
SYNOPSIS#
function sidedNodal
DESCRIPTION#
spaced data points. This function is handy for advective terms. m : Number of cells dx : Step size type : 'backward', 'forward' or 'centered'
CROSS-REFERENCE INFORMATION#
No cross-reference information found. This typically means this function neither calls nor is called by other functions in the codebase.
Mimetic Weights#
- weightsP(k, m, dx)#
PURPOSE#
Returns the m+1 weights of P
SYNOPSIS#
function weightsP
DESCRIPTION#
k : Order of accuracy m : Number of cells dx : Step size
CROSS-REFERENCE INFORMATION#
This function calls:
grad()
Returns a one-dimensional mimetic gradient operator depending on whetherThis function is called by:
weightsP2D()
Returns the 2mn+m+n weights of P in 2-DmimeticB()
Returns a m+2 by m+1 one-dimensional mimetic boundary operator
- weightsP2D(k, m, dx, n, dy)#
PURPOSE#
Returns the 2mn+m+n weights of P in 2-D
SYNOPSIS#
function weightsP2D
DESCRIPTION#
k : Order of accuracy m : Number of cells along x-axis dx : Step size along x-axis n : Number of cells along y-axis dy : Step size along y-axis
CROSS-REFERENCE INFORMATION#
This function calls:
weightsP()
Returns the m+1 weights of P
- weightsQ(k, m, dx)#
PURPOSE#
Returns the m+2 weights of Q
SYNOPSIS#
function weightsQ
DESCRIPTION#
k : Order of accuracy m : Number of cells dx : Step size
CROSS-REFERENCE INFORMATION#
This function calls:
div()
Returns a one-dimensional mimetic divergence operator depending on whetherThis function is called by:
mimeticB()
Returns a m+2 by m+1 one-dimensional mimetic boundary operator
- weightsQ2D(m, n, d)#
PURPOSE#
Returns the (m+2)(n+2) weights of Q in 2-D
SYNOPSIS#
function weightsQ2D
DESCRIPTION#
Only works for 2nd-order 2-D Mimetic divergence operator m : Number of cells along x-axis n : Number of cells along y-axis d : Step size (assuming d = dx = dy)
CROSS-REFERENCE INFORMATION#
No cross-reference information found. This typically means this function neither calls nor is called by other functions in the codebase.
Utility Functions#
- amean(X)#
PURPOSE#
Returns the arithmetic mean for every two pairs in a column vector
SYNOPSIS#
function amean
DESCRIPTION#
And, Y(1) = X(1), Y(end) = X(end) X : Column vector
CROSS-REFERENCE INFORMATION#
No cross-reference information found. This typically means this function neither calls nor is called by other functions in the codebase.
- hmean(X)#
PURPOSE#
Returns the harmonic mean for every two pairs in a column vector
SYNOPSIS#
function hmean
DESCRIPTION#
And, Y(1) = X(1), Y(end) = X(end) X : Column vector
CROSS-REFERENCE INFORMATION#
No cross-reference information found. This typically means this function neither calls nor is called by other functions in the codebase.
- rk4(func, tspan, dt, y0)#
PURPOSE#
Explicit Runge-Kutta 4th-order method
SYNOPSIS#
function rk4
DESCRIPTION#
Returns : t (evaluation points) and y (solutions) of the specified ODE func : Function handler tspan : [t0 tf] dt : Step size y0 : Initial conditions
CROSS-REFERENCE INFORMATION#
No cross-reference information found. This typically means this function neither calls nor is called by other functions in the codebase.
- ttm(grid_name, m, n, iters, plot_grid)#
PURPOSE#
Returns X and Y which are both m by n matrices that contains the physical
SYNOPSIS#
function ttm
DESCRIPTION#
coordinates grid_name : String with the name of the grid folder m : Number of nodes along the horizontal axis n : Number of nodes along the vertical axis plot_grid : If defined -> grid will be plotted
CROSS-REFERENCE INFORMATION#
This function is called by:
gridGen()
Returns X and Y which are both m by n matrices that contains the physical
- boundaryIdx2D(m, n)#
PURPOSE#
Returns the indices of the nodes that lie on the boundary of a 2D nodal
SYNOPSIS#
function boundaryIdx2D
DESCRIPTION#
grid m : Number of nodes along x-axis n : Number of nodes along y-axis
CROSS-REFERENCE INFORMATION#
No cross-reference information found. This typically means this function neither calls nor is called by other functions in the codebase.
- DI2(m, n, type)#
- DI3(m, n, o, type)#
- GI1(M, m, n, type)#
- GI13(M, m, n, o, type)#
- GI2(M, m, n, type)#
- mimeticB(k, m)#
PURPOSE#
Returns a m+2 by m+1 one-dimensional mimetic boundary operator
SYNOPSIS#
function mimeticB
DESCRIPTION#
k : Order of accuracy m : Number of cells
CROSS-REFERENCE INFORMATION#
This function calls:
div()
Returns a one-dimensional mimetic divergence operator depending on whetherweightsQ()
Returns the m+2 weights of Qgrad()
Returns a one-dimensional mimetic gradient operator depending on whetherweightsP()
Returns the m+1 weights of P
- tensorGrad2D(K, G)#
PURPOSE#
Returns a two-dimensional flux operator
SYNOPSIS#
function tensorGrad2D
DESCRIPTION#
K : Tensor (e.g. diffusion tensor) G : 2D mimetic gradient operator
CROSS-REFERENCE INFORMATION#
No cross-reference information found. This typically means this function neither calls nor is called by other functions in the codebase.