API Reference

BEAST.DoubleLayerRotatedMW3DType
struct DoubleLayerRotatedMW3D{T,K} <: MaxwellOperator3D{T,K}

Bilinear form given by:

```math
    α ∬_{Γ^2} k(x) ⋅ [n̂(x) × (∇G_γ(x-y) × j(y))]
```

with ``G_γ = e^{-γ|x-y|} / 4π|x-y|``

Fields

  • alpha::T: Factor in front of bilinear form.
  • gamma::K: imaginary unit times the wavenumber.
source
BEAST.HH3DHyperSingularFDBIOType
∫_Γ dx ∫_Γ dy \left(α G g(x) n_x ⋅ n_y f(y) + β G \mbox{curl} g(x) ⋅ \mbox{curl} f(y) \right)

with $G(x,y) = \frac{e^{-γ |x-y|}}{4 π |x-y|}$

source
BEAST.HH3DLinearPotentialType
HH3DLinearPotential

A potential that linearly increases in direction with scaling coefficient amplitude. Its negative gradient will be a uniform vector field pointing in the opposite direction.

source
BEAST.MWFarField3DTDType

Operator to compute the far field of a current distribution in the time domain. In particular, given the current distribution $j$ this operator allows for the computation of

\[R = ̂x ⋅ y ffd = n × ∫_Γ j(r', t - R/c} dy\]

where $̂x$ is the unit vector in the direction of observation. Note that the assembly routing expects the observation directions to be normalised by the caller.

source
BEAST.MWSingleLayerField3DMethod
MWSingleLayerField3D(;gamma, wavenumber, alpha, beta)

Create the single layer near field operator, for use with potential.

source
BEAST.NCrossType
NCross <: LocalOperator

The identity operator where the trial function is rotated.

source
BEAST.NDRefSpaceType

Local shape function r is the one whose field lines straddle local edge r, which is the edge adjacent to vertex r.

This is not the edge starting at vertex r. The downside of this local numbering scheme is that it cannot be extended to cells that are not simplices because there is no well defined concept of adjacent-ness.

source
BEAST.RungeKuttaConvolutionQuadratureType
RungeKuttaConvolutionQuadrature{T,N,NN}

T: the value type of the basis function. N: the number of stages. NN: N*N.

Performs a convolution quadrature on a laplaceKernel to represent an operator in time domain using an implicit Runge-Kutta method.

laplaceKernel: function of the Laplace variable s that returns an IntegralOperator. A, b: Coefficient matrix and vectors from the Butcher tableau. Δt: time step. zTransformedTermCount: Number of terms in the inverse Z-transform. contourRadius: radius of circle used as integration contour for the inverse Z-transform.

source
BEAST.SingleLayerTraceType

Describe a single layer operator from the surface to a line.

\[<v, Su> = ∫_γ dx v(x) ∫_Γ dy rac{e^{-ikR}}{4πR} u(y)\]

source
BEAST.StagedTimeStepType
StagedTimeStep{T,N,NN}

T: the value type of the basis function. N: the number of stages. NN: the number of stages squared Each time step has intermediary stages given by the vertor c in a Butcher tableau (A,b,c)

source
BEAST.TimeBasisFunctionType
TimeBasisFunction{N,D}

T: the value type of the time basis function N: the number of intervals in the support (this included the semi infinite interval stretching to +∞) D1: the degree of the TBF restricted to each of the intervals plus one

source
BEAST.VIEFarField3DType

Operator to compute the far field of a current distribution. In particular, given the current distribution $j$ this operator allows for the computation of

\[A j = n × ∫_Ω j e^{γ ̂x ⋅ y} dy\]

where $̂x$ is the unit vector in the direction of observation. Note that the assembly routing expects the observation directions to be normalised by the caller.

source
BEAST.assembleMethod
assemble(operator, test_functions, trial_functions;
    storage_policy = Val{:bandedstorage},
    threading = Threading{:multi},
    quadstrat=defaultquadstrat(operator, test_functions, trial_functions))

Assemble the system matrix corresponding to the operator operator tested with the test functions test_functions and the trial functions trial_functions.

source
BEAST.assembleMethod
assemble(fn, tfs)

Assemble the vector of test coefficients corresponding to functional fn and test functions tfs.

source
BEAST.assemblechunk!Method
assemblechunk!(biop::IntegralOperator, tfs, bfs, store)

Computes the matrix of operator biop wrt the finite element spaces tfs and bfs

source
BEAST.assemblydataMethod
charts, admap, act_to_global = assemblydata(basis; onlyactives=true)

Given a basis this function returns a data structure containing the information required for matrix assemble, that is, the vector charts containing Simplex elements, a variable admap of type AssemblyData, and a mapping from indices of actively used simplices to global simplices.

When onlyactives is true, another layer of indices is introduced to filter out all cells of the mesh that are not in the union of the support of the basis functions (i.e., when the basis functions are defined only on a part of the mesh).

admap is, in essence, a three-dimensional array of named tuples, which, by wrapping it in the struct AssemblyData, allows the definition of iterators. The tuple consists of the two entries

admap[i,r,c].globalindex
admap[i,r,c].coefficient

Here, c and r are indices in the iterable set of (active) simplices and the set of shape functions on each cell/simplex: r ranges from 1 to the number of shape functions on a cell/simplex, c ranges from 1 to the number of active simplices, and i ranges from 1 to the number of maximal number of basis functions, where any of the shape functions contributes to.

For example, for continuous piecewise linear lagrange functions (c0d1), each of the three shape functions on a triangle are associated with exactly one Lagrange function, and therefore i is limited to 1.

Note: When onlyactives=false, the indices c correspond to the position of the corresponding cell/simplex whilst iterating over geometry(basis). When onlyactives=true, then act_to_global(c) correspond to the position of the corresponding cell/simplex whilst iterating over geometry(basis).

For a triplet (i,r,c), globalindex is the index in the basis of the ith basis function that has a contribution from shape function r on (active) cell/simplex c. coefficient is the coefficient of that contribution in the linear combination defining that basis function in terms of shape function.

source
BEAST.basisfunctionMethod
basisfunction(basis, i)

Returns a vector of the shape functions defining the ith function of the basis.

source
BEAST.blockassemblerFunction
blockassembler(operator, test_space, trial_space) -> assembler

Return a callable object for the creation of blocks within a BEM matrix.

This function performs all tasks common to the assembly of several blocks within a single boundary element matrix. The return value can be used to generate blocks by calling it as follows:

assembler(I,J,storefn)

where I and J are arrays of indices in test_space and trial_space, respectively, corresponding to the rows and columns of the desired block.

Note that the block will be constructed in compressed form, i.e. the rows and columns of the store that are written into are the positions within I and J (as opposed to the positions within 1:numfunctions(test_space) and 1:numfunctions(trial_space)). In particular the size of the constructed block will be (length(I), length(J)).

This last property allows the assembly of permutations of the BEM matrix by supplying for I and J permutations of 1:numfunctions(test_space) and 1:numfunctions(trial_space).

source
BEAST.buffachristiansenFunction
buffachristiansen(Γ, γ)

Construct the set of Buffa-Christiansen functions subject to mesh Γ and only enforcing zero normal components on ∂Γ ∖ γ.

source
BEAST.curlMethod
curl(X)

Compute the curl of a finite element basis. The resulting set of functions might be linearly dependent because of the kernel of the curl operator.

source
BEAST.dipolemw3dMethod
dipolemw3d(;location, orientation, wavenumber)

Create an electric dipole solution to Maxwell's equations representing the electric field part. Implementation is based on (9.18) of Jackson's “Classical electrodynamics”, with the notable difference that the $xp(ikr)$ is used.

source
BEAST.duallagrangec0d1Method
duallagrangec0d1(originalmesh, refinedmesh)

It is the user responsibility to provide two meshes representing the same object. The second mesh needs to be obtained using "barycentric_refinement(originalmesh)". This basis function creats the dual Lagrange basis function and return an object that contains array of shapes [fns] It also return a gemoetry containing the refined mesh.

source
BEAST.duallagrangecxd0Function
duallagrangecxd0(mesh, jct) -> basis

Build dual Lagrange piecewise constant elements. Boundary nodes are only considered if they are in the interior of jct.

The default dual function (interpolatory=false) is similar to the one depicted in Figure 3 of Buffa et al (doi: 10.1090/S0025-5718-07-01965-5), with the difference that each individual shape function is normalized with respect to the area so that overall the integral over the dual function is one.

When interpolatory=true is used, the function value is one on the support, and thus, it gives rise to a partition of unity.

source
BEAST.elementsMethod

elements(geo)

Create an iterable collection of the elements stored in geo. The order in which this collection produces the elements determines the index used for lookup in the data structures returned by assemblydata and quaddata.

source
BEAST.facecurrentsMethod
fcr, geo = facecurrents(coeffs, basis)

Compute the value of the function with the given collection of coeffient in the provided basis in all the centroids of the mesh underlying the basis. The mesh is returned together with the currents.

source
BEAST.gamma_wavenumber_handlerMethod
gamma_wavenumber_handler(gamma, wavenumber)

This function handles the input of gamma and wavenumber. It throws an error if both gamma and wavenumber are provided. If neither is provided, it assumes a static problem and returns Val(0) for gamma and wavenumber.

Arguments

  • gamma: im * wavenumber or nothing.
  • wavenumber: wavenumber or nothing.

Returns

  • gamma and wavenumber: Appropriate pair gamma and wavenumber.
source
BEAST.geometryMethod
geometry(basis)

Returns an iterable collection of geometric elements on which the functions in basis are defined. The order the elements are encountered needs correspond to the element indices used in the data structure returned by assemblydata.

source
BEAST.getcommonedgeMethod
getcommonedge(cell1, cell2) -> e1, e2, edge

Returns in edge the common vertices of cell1 and cell2. e1 contains the index of the vertex of cell1 opposite to this common edge, and with a plus or minus sign depending on whether the orientation of the common edge is along or against the internal orientation of cell1. Similar for e2.

source
BEAST.getindex_rtgMethod
getindex_rtg(RT::RTBasis)

Returns the indices of the global half RWGs present in RT. RT is typically gotten from rt_ports

source
BEAST.index_activesMethod
active, index_among_actives, num_active_cells, act_to_global =
    index_actives(num_cells, celltonum)

Given the number of cells of the mesh num_cells and an array indicating, which shape functions of each cell are used in the definition of the basis, celltonum, index_actives(num_cells, celltonum) computes a Boolean vector active that indicates for each cell index if the cell is in the support of any basis function. For this reduced set of active cells, a new index set is introduced, where index_among_actives is a vector providing the mapping from the cells of the mesh to the indices of the active cells, num_active_cells provides the number of all active cells, and act_to_global is a vector providing the mapping from the indices of the active cells to the indices of cells of the underlying mesh of the basis function spce.

source
BEAST.indices_splitfemglobalMethod
indices_splitfemglobal(lincombgfs)

Given a linear combination of FEMFunctions and/or GlobalFunctions, return the indices of the fem and the global functions

source
BEAST.interpolateMethod
interpolate(interpolant::RefSpace, chart1, interpolee::RefSpace, chart2)

Computes by interpolation approximations of the local shape functions for interpolee on chart2 in terms of the local shape functions for interpolant on chart1. The returned value is a matrix Q such that

\[\phi_i \approx \sum_j Q_{ij} \psi_j\]

with $\phi_i$ the i-th local shape function for interpolee and $\psi_j$ the j-th local shape function for interpolant.

source
BEAST.inverse_z_transformMethod
inverse_z_transform(k, rho, N, X)

Returns the k-th term of the inverse z-transform. X is an array of the z-transform evaluated in the points z=rhoexp(2impin/N) for n in 0:(N-1).

source
BEAST.isstaticMethod
isstatic(gamma)

This function checks if the provided gamma value represents a static problem. It returns true if gamma is of type Val{0} indicating a static problem.

Arguments

  • gamma: gamma value.

Returns

  • A boolean indicating whether the problem is static or not.
source
BEAST.lagdimensionFunction

The dimension of the space of Lagrange shape functions of degree d over a simplex of dimension n is binom(n+d,d) == binom(n+d,n)

source
BEAST.lagrangec0d1Method
lagrangec0d1(mesh; dirichlet=[true|false]) -> basis

Build lagrangec0d1 elements, including (dirichlet=false) or excluding (dirichlet=true) those attached to boundary vertices.

source
BEAST.lagrangec0d1_dirichletMethod
lagrangec0d1(mesh[, bnd])

Construct the basis of continuous, piecewise linear basis functions subordinate to mesh mesh. Basis functions are constructed at vertices in the interionr of the mesh and on the closure of 'bnd'. In particular, leaving out the second argument creates a finite element space subject to homogeneous Dirichlet boundary conditions.

source
BEAST.laplace_to_zMethod
laplace_to_z(rho, n, N, dt, A, b)

Returns the complex matrix valued Laplace variable s that correspond to the variable z = rhoexp(2impin/N) for a given Butcher tableau (A,b,c) and a time step dt.

source
BEAST.linearpotentialvieMethod

linearpotentialvie(; direction = error("missing argument direction"), amplitude = 1, )

Linear potential for volume integral equations.

source
BEAST.localindicesMethod
localindices(dof, chart, i)

Returns a vector of indices into the vector of local shape functions that correspond to global degrees of freedom supported on sub-entity i, where the type of entity (nodes, edge, face) is encoded in the type of 'dof'.

source
BEAST.make_celltonumMethod
celltonum = make_celltonum(num_cells, num_refs, num_bfs, basis)

Computes the array celltonum[c,r] that, given the index c of a cell of the mesh associated with the basis and the index r of the shape function in the refspace associated with basis, provides the number of basis functions where the rth shape function on cell c is used in their definition.

source
BEAST.marchonintimeMethod
marchonintime(W0,Z,B,I; convhist=false)

Solve by marching-on-in-time the causal convolution problem defined by (W0,Z,B) up to timestep I. Here, Z is an array of order 3 that contains a discretisation of a time translation invariant retarded potential operator. W0 is the inverse of the slice Z[:,:,1].

Keyword arguments: - 'convhist': when true, return in addition to the space-time data for the solution also the vector of convergence histories as returned each time step by the supplied solver W0.

source
BEAST.momintegrals!Method

momintegrals!(biop, tshs, bshs, tcell, bcell, interactions, strat)

Function for the computation of moment integrals using simple double quadrature.

source
BEAST.nedelec2Function
nedelec2(mesh, edges)

Constructs the 2nd degree Nedelec basis of the first kind i.e. H(curl).

Returns an object of type 'ND2Basis'.

source
BEAST.ntraceMethod
ntrace(refspace, element, localindex, face)

Compute the normal trace of all local shape functions on elements belonging to refspace on face. This function returns a matrix expressing the traces of local shape functions in refspace as linear combinations of functions in the local trace space. Cf. restrict. localindex is the index of face in the enumeration of faces of elements. In many special cases knowing this index allows for highly optimised implementations.

source
BEAST.ntraceMethod
ntrace(X::Space, γ::Mesh)

Compute the normal trace of basis X on mesh γ. γ is assumed to be part of the boundary of geometry(X).

source
BEAST.planewaveMethod
planewave(polarisation,direction,amplitude,speedoflight)

Time-domain plane wave.

source
BEAST.planewavemw3dMethod
planewavemw3d(;direction, polarization, wavenumber, gamma[, amplitude=1])

Create a plane wave solution to Maxwell's equations.

source
BEAST.planewavevieMethod

planewavevie(; direction = error("missing arguement direction"), polarization = error("missing arguement polarization"), wavenumber = error("missing arguement wavenumber"), amplitude = 1, )

For volume integral equations

source
BEAST.portcellsMethod
portcells(Γ::Mesh, γ::Mesh)

returns an array containing cell pairs of mesh Γ around a boundary edge that overlaps with mesh γ

source
BEAST.potentialMethod
potential(op, points, coeffs, basis)

Evaluate operator for a given bases and expansion coefficients at the given points.

source
BEAST.qhMethod
Q = qd(T,dh,::Val{N})

Q[k] is the factor in front resulting from differentiating t^(k-1) dh times.

source
BEAST.quaddataFunction
quaddata(operator, test_refspace, trial_refspace, test_elements, trial_elements)

Returns an object cashing data required for the computation of boundary element interactions. It is up to the client programmer to decide what (if any) data is cached. For double numberical quadrature, storing the integration points for example can significantly speed up matrix assembly.

  • operator is an integration kernel.
  • test_refspace and trial_refspace are reference space objects. quadata

is typically overloaded on the type of these local spaces of shape functions. (See the implementation in maxwell.jl for an example).

  • test_elements and trial_elements are iterable collections of the geometric

elements on which the finite element space are defined. These are provided to allow computation of the actual integrations points - as opposed to only their coordinates.

source
BEAST.quadruleFunction
quadrule(operator, test_refspace, trial_refspace, test_index, test_chart, trial_index, trial_chart, quad_data)

Based on the operator kernel and the test and trial elements, this function builds an object whose type and data fields specify the quadrature rule that needs to be used to accurately compute the interaction integrals. The quad_data object created by quaddata is passed to allow reuse of any precomputed data such as quadrature points and weights, geometric quantities, etc.

The type of the returned quadrature rule will help in deciding which method of momintegrals to dispatch to.

source
BEAST.raviartthomasMethod
raviartthomas(mesh)

Conducts pre-processing on the input mesh by extracting the cell edges, cell pairs and indices required to construct the RT basis on the mesh.

Calls raviartthomas(mesh::Mesh, cellpairs::Array{Int,2}), which constructs the RT basis on the mesh, using the cell pairs identified.

Returns the RT basis object.

source
BEAST.raviartthomasMethod
raviartthomas(mesh, cellpairs::Array{Int,2})

Constructs the RT basis on the input mesh. The i-th RT basis function will represent a current distribution flowing from cell cellpairs[1,i] to cellpairs[2,i] on the mesh.

Returns an object of type RTBasis, which comprises both the mesh and pairs of Shape objects which corresponds to the cell pairs, containing the necsessary coefficients and indices to compute the exact basis functions when required by the solver.

source
BEAST.raviartthomas2Method
raviartthomas2(mesh)

Conducts pre-processing on the input mesh by extracting the cell edges, cell pairs and indices required to construct the RT2 basis on the mesh.

Calls raviartthomas2(mesh::Mesh, cellpairs::Array{Int,2}), which constructs the RT2 basis on the mesh, using the cell pairs identified.

Returns the RT2 basis object.

source
BEAST.raviartthomas2Method
raviartthomas2(mesh, cellpairs::Array{Int,2})

Constructs the RT2 basis on the input mesh. The i-th RT2 basis function will represent a current distribution flowing from cell cellpairs[1,i] to cellpairs[2,i] on the mesh.

Returns an object of type RT2Basis, which comprises both the mesh and pairs of Shape objects which corresponds to the cell pairs, containing the necsessary coefficients and indices to compute the exact basis functions when required by the solver.

source
BEAST.real_inverse_z_transformMethod
real_inverse_z_transform(k, rho, N, X)

Returns the k-th term of the inverse z-transform. It is assumed that X[n+1] = conj(X[N-n]) for each n in 1:(N-1) so that Nmax = N/2+1 or (N+1)/2 (resp. if N%2==0 or N%2==1) terms are used in X X is an array of the z-transform evaluated in the points z=rhoexp(2impin/N) for n in 0:(Nmax-1).

source
BEAST.refspaceFunction
refspace(basis)

Returns the ReferenceSpace of local shape functions on which the basis is built.

source
BEAST.restrictFunction
restrict(refspace, element1, element2)

Computes the restriction of a set of local shape functions on element1 as linear combinations of the set of local shape functions on element2. More precisely restrict returns an NxM matrix P such that the i-th local shape $g_i$ function on element2 can be written as:

$g_i = sum_{j=1}^{M} P_{ij} f_j$

source
BEAST.rt_cedgeMethod
rt_cedge(cps::Array{Int,2}, weight)

Computes single basis function with equally distributed constant current leaving or entering port defined by cellpairs cps. weight defines the total current over the port and its direction (+ve = out, -ve = in)

source
BEAST.rt_portsMethod
rt_ports(Γ::Mesh, γ::Mesh ...)

Constructs the RT space on Γ, relative to boundary pairs in γ. γ expects any number of pairs-of-ports as arguments and accepts tuples, arrays, vectors etc. e.g rt_ports(Γ, a, b ...); where a = [γ₁ γ₂], b = (γ₃,γ₄) etc. rt_ports with no pair of ports supplied i.e rt_ports(Γ) reduces to the raviartthomas(Γ) function. The RT space ensures current continuity in each pair of ports. i.e. current leaving mesh Γ through γ₁ is accounted for in γ₂.

Returns the RT basis object.

source
BEAST.rt_vedgeMethod
rt_vedge(cps::Array{Int,2}, weight)

Computes n-1 basis function with oscillating current in(leaving and entering) pairs of half triangles defined over port specified by cellpairs cps. weight defines the magnitude of individual current in and out the half triangles, and it's polarity simply defines whether to start with in or out

source
BEAST.scalartypeFunction
scalartype(x)

The scalar field over which the values of a global or local basis function, or an operator are defined. This should always be a scalar type, even if the basis or operator takes on values in a vector or tensor space. This data type is used to determine the eltype of assembled discrete operators.

source
BEAST.scalartypeMethod
scalartype(s)

The scalar field over which the argument to a basis function or operator's integration kernel are defined. This is always a salar data type, even if the function or kernel is defined over a multi-dimensional space.

source
BEAST.singleduallagd0Method
singleduallagd0(fine, F, v; interpolatory=false)

Build a single dual constant Lagrange element a mesh fine. F contains the indices to cells in the support and v is the index in the vertex list of the defining vertex.

The default dual function (interpolatory=false) is similar to the one depicted in Figure 3 of Buffa et al (doi: 10.1090/S0025-5718-07-01965-5), with the difference that each individual shape function is normalized with respect to the area so that overall the integral over the dual function is one.

When interpolatory=true is used, the function value is one on the support, and thus, it gives rise to a partition of unity.

source
BEAST.solveMethod

Solves a variational equation by simply creating the full system matrix and calling a traditional lu decomposition.

source
BEAST.timebasisc0d1Function
timebasisc0d1(type, timestep, numfunctions)

Build the space of continuous, piecewise linear time basis functions. The DoFs are the time steps. numfunctions basis functions will be built in total.

source
BEAST.timebasiscxd0Function
timebasiscxd0(timestep, numfunctions, T::Type=Float64)

Create a temporal basis based on shifted copies of the nodal continuous, piecewise linear interpolant.

source
BEAST.timebasisspline2Function
timebasisspline2(timestep, numfunctions, T::Type=Float64)

Create a temporal basis based on shifted copies of the quadratic spline. The spline is the convolution of a cxd0 and a c0d1 basis function.

source
BEAST.ttraceMethod

Does not give the correct result for an imput basis with non-vanishing input basis on the boundary

source
BEAST.unitfunctionc0d1Method
unitfunctionc0d1(mesh)

Constructs a constant function with value 1 on mesh consisting of linear shapes. For dirichlet=true goes to zero on the boundary.

source
BEAST.@discretiseMacro
discr(eq, pairs...)

This macro provides syntactical sugar for the definition of a discretisation of a varational formulation. Given a variational equation EQ: Find j ∈ X such that for all k ∈ Y a(k,j) = f(k) can be discretised by stating:

eq = @discretise EQ j∈X k∈Y
source
BEAST.VIE.hhboundaryMethod

hhboundary(;gamma, alpha, tau) hhboundary(;wavenumber, alpha, tau)

Bilinear form given by:

\[ α ∬_{∂Ω×Ω} n̂(x) ⋅ j(x) G_{γ}(x,y) τ(y) (grad k(y))\]

with $G_{γ} = e^{-γ|x-y|} / 4π|x-y|$ and $τ(y)$ contrast function and $n̂(x)$ normal vector

source
BEAST.VIE.hhvolumeMethod
hhvolume(;gamma, alpha, tau)
hhvolume(;wavenumber, alpha, tau)

Bilinear form given by:

\[ α ∬_{Ω×Ω} (grad j(x)) ⋅ G_{γ}(x,y)τ(y) (grad k(y))\]

with $G_{γ} = e^{-γ|x-y|} / 4π|x-y|$ and $τ(y)$ contrast function

source
BEAST.VIE.hhvolumegradGMethod
hhvolumegradG(;gamma, alpha, tau)
hhvolumegradG(;wavenumber, alpha, tau)

Bilinear form given by:

\[ α ∬_{Ω×Ω} j(x) grad_y(G_{γ}(x,y)) τ(y) ⋅ (grad k(y))\]

with $G_{γ} = e^{-γ|x-y|} / 4π|x-y|$ and $τ(y)$ constant function

source
BEAST.VIE.hhvolumek0Method

hhvolumek0(;gamma, alpha, tau) hhvolumek0(;wavenumber, alpha, tau)

Bilinear form given by:

\[ α ∬_{Ω×Ω} j(x) G_{γ}(x,y)τ(y) k(y)\]

with $G_{γ} = e^{-γ|x-y|} / 4π|x-y|$ and $τ(y)$ contrast function

source
BEAST.Maxwell3D.doublelayerMethod
doublelayer(;gamma)
doublelayer(;wavenumber)

Maxwell double layer operator.

Either gamma or the wavenumber must be provided. Optionally, also alpha can be provided.

If alpha is not provided explitly, it is set to $α = 1$.

source
BEAST.Maxwell3D.hypersingularMethod
hypersingular(;wavenumber)

Hyper singular part of the Maxwell 3D single layer operator.

$β = -1/\mathrm{j} k$ is set with the wavenumber $k$.

source
BEAST.Maxwell3D.planewaveMethod
planewave(;
        direction    = error("missing arguement `direction`"),
        polarization = error("missing arguement `polarization`"),
        wavenumber   = error("missing arguement `wavenumber`"),
        amplitude    = one(real(typeof(wavenumber))))

Time-harmonic plane wave.

source
BEAST.Maxwell3D.singlelayerMethod
singlelayer(;gamma, alpha, beta)
singlelayer(;wavenumber, alpha, beta)

Maxwell 3D single layer operator.

Either gamma, or the wavenumber, or α and β must be provided.

If α and β are not provided explitly, they are set to $α = -γ$ and $β = -1/γ$ with $γ=\mathrm{j} k$.

source
BEAST.Maxwell3D.weaklysingularMethod
weaklysingular(;wavenumber)

Weakly singular part of the Maxwell 3D single layer operator.

$α = -\mathrm{j} k$ is set with the wavenumber $k$.

source
BEAST.Variational.depthfirstMethod
depthfirst(xp)

Returns an iterator that visits all nodes in an Expr depth first. head and args are visited before the Expr they define.

source
Base.:==Method
==(lhs::BilForm, rhs::LinForm)

Build an equation from a left hand and right hand side

source
Base.getindexMethod
getindex(A, v::HilbertVector, u::HilbertVector)

Create a BilForm corresponding to A[v,u]

source
Base.getindexMethod
getindex(f, v::HilbertVector)
f[v]

Return a LinForm corresponding to f[v]

source
BEAST.Variational.@varformMacro
@varform <form-definition>

The Julia form compiler uses the Julia parser and meta-programming based traversal of the AST to create a structure containing all information required for the description of a variational formulation from an Expr that follows closely widespread mathematical convention.

E.g:

EFIE = @varform T[k,j] = e[k]
MFIE = @varform 0.5*I[k,j] + K[k,j] = h[k]
PMCH = @varform M[k,j] - η*T[k,m] + 1/η*T[l,j] + M[l,m] = e[k] + h[l]
source