API
GalerkinToolkit.AbstractDomain — Type
abstract type AbstractDomain <: AbstractType endAbstract type representing a subset of $\mathbb{R}^d$, typically $d\in\{0,1,2,3\}$. Domains are defined using an underlying computational mesh. See also AbstractMesh.
Level
Beginner
Basic constructors
Basic queries
GalerkinToolkit.AbstractFaceDomain — Type
abstract type AbstractFaceDomain <: AbstractDomain endA domain defined on a single mesh face. Typically used as helper to identify cases that only make sense for a single mesh face.
Level
Advanced
Basic constructors
GalerkinToolkit.AbstractFaceSpace — Type
abstract type AbstractFaceSpace <: AbstractSpace end
Like AbstractSpace, but for a single mesh face. Typically used as helper to identify cases that only make sense for a single mesh face.
Level
Advanced
GalerkinToolkit.AbstractFaceTopology — Type
abstract type AbstractFaceTopology <: AbstractTopology end
Like AbstractTopology, but for a single mesh face. Typically used as helper to identify cases that only make sense for a single mesh face.
Level
Advanced
GalerkinToolkit.AbstractMesh — Type
abstract type AbstractMesh <: AbstractType endAbstract type representing a computational mesh.
Level
Beginner
Constructors
Iteration
Data retrival
GalerkinToolkit.AbstractMeshDomain — Type
abstract type AbstractMeshDomain <: AbstractDomain endGalerkinToolkit.AbstractQuadrature — Type
abstract type AbstractQuadratureBasic queries
Basic constructors
Supertype hierarchy
AbstractQuadrature <: GT.AbstractTypeGalerkinToolkit.AbstractSpace — Type
abstract type AbstractSpace <: AbstractType endAbstract type representing a finite element space.
Level
Basic
Basic constructors
lagrange_space raviart_thomas_space
Basic queries
domainnum_dofsface_dofsface_nodesface_reference_idreference_spacesgeometry_own_dofsgeometry_own_dofs_permutations
Additional queries
For spaces, used as reference spaces in AbstractMesh specializations.
GalerkinToolkit.AbstractTopology — Type
abstract type AbstractTopologyAbstract type representing the incidence relations in a face complex.
Level
Intermediate
Constructors
Queries
GalerkinToolkit.AbstractType — Type
abstract type AbstractType endParent of all types defined in GalerkinToolkit.
GalerkinToolkit.Mesh — Type
GalerkinToolkit.Options — Type
struct Options{...} <: AbstractTypeType of the objects returned by function options. All properties and type parameters are private.
Basic queries
GalerkinToolkit.UnitNCube — Type
struct UnitNCube{...} <: AbstractFaceDomainGalerkinToolkit.UnitSimplex — Type
struct UnitSimplex{...} <: AbstractFaceDomainGalerkinToolkit.ast_hoist! — Method
ast_hoist!(ast)Hoist loop invariant definitions in ast. The new AST is overwritten in ast.
Hypotheses:
- The root of the AST is a code block.
- The code block contains only for-loops, function calls, array indexing, definitions, and increments. No while loops, if statements, lambdas, etc.
- Functions that appear in a rhs of a definition are pure functions.
- RHSs contain a function call at most.
GalerkinToolkit.ast_lhs_indices — Method
ast_lhs_indices(ast)GalerkinToolkit.boundary — Function
GalerkinToolkit.bounding_box — Function
p0,p1 = bounding_box(x)Return a tuple of two vectors, where the vectors p0 and p1 define the span of the bounding box of x.
Level
Intermediate
GalerkinToolkit.cartesian_mesh — Method
cartesian_mesh(domain,cells_per_dir)Create a multi-dimensional Cartesian mesh. The dimension of the mesh is defined by the length of cells_per_dir. The number of cells in direction i is given by cells_per_dir[i]. The extends of the Cartesian mesh are given in domain. The range in direction i covered by the mesh is given by domain[2*i-1,2*i].
Keyword arguments
boundary=true[optional]: Include faces on the boundary and generate face groups identifying which faces are on which face of bounding box of the mesh. The groups are named$d-face-$ifor the faceiof dimensiondof the bounding box.complexify=true[optional]: Generate all low dimensional faces so that the mesh is a face complex.simplexify=false[optional]: Generate a mesh of simplex faces instead of hyper-cubes.
GalerkinToolkit.complexify — Function
complexify(x)Convert x into a mesh representing a face complex.
See also is_face_complex.
Level
Intermediate
GalerkinToolkit.coordinates — Function
GalerkinToolkit.create_chain — Method
create_chain(;kwargs...)Build an arbitrary mesh object, containing all faces of the same dimension. This function is similar to create_mesh but it only receives face arrays of one dimension.
See also create_mesh.
Level
Intermediate
Keyword arguments
node_coordinates: Like forcreate_mesh.face_nodes: A nested vector containing the node ids for each face in the mesh.node_coordinates[n]withn=face_nodes[i][k]is the global node coordinate vector for local node numberkin facei.reference_spaces: A tuple containing the reference spaces for faces.reference_spaces[i]is the reference space numberi.face_reference_id[optional]: A vector containing which reference space is assigned to each face.reference_sapces[r]withr=face_reference_id[i]is the reference space associated with face numberi. By default, all faces are assigned to the first reference space in its dimension.group_faces[optional]: A Dictionary containing labeled groups of faces.group_faces[group_name]is a vector of integers containing the ids of the faces in the group namedgroup_name. These groups might overlap. By default, no faces groups are created.normals=nothing[optinal]: Like forcreate_mesh.
GalerkinToolkit.create_mesh — Function
create_mesh(;kwargs...)Build an arbitrary mesh object.
See also cartesian_mesh, mesh_from_msh, and mesh_from_gmsh.
Level
Intermediate
Keyword arguments
node_coordinates: The vector containing the coordinates of all mesh nodes.node_coordinates[i]is the coordinate vector for global node numberi.face_nodes: A highly-nested vector containing the node ids for each face in the mesh.node_coordinates[n]withn=face_nodes[d+1][i][k]is the global node coordinate vector for local node numberkin faceiof dimensiond. The objectface_nodes[d+1]is a long vector of small vectors of integers. It is often represented using aJaggedArrayobject that uses continuous linear memory for performance.reference_spaces: A nested tuple containing the reference spaces for faces.reference_spaces[d+1][i]is the reference space numberiin dimensiond. Reference interpolation spaces are defined with functions likelagrange_space.face_reference_id[optional]: A nested vector containing which reference space is assigned to each face.reference_sapces[d+1][r]withr=face_reference_id[d+1][i]is the reference space associated with face numberiof dimensiond. By default, all faces are assigned to the first reference space in its dimension.group_faces[optional]: A vector of dictionaries containing labeled groups of faces.group_faces[d+1][group_name]is a vector of integers containing the ids of the faces of dimensiondin the group namedgroup_name. These groups might overlap. By default, no faces groups are created.is_face_complex=Val(false)[optional]:Val(true)if the input data represents a face complex,Val(false)otherwise.normals=nothing[optinal]: Vector containing the normal vectors for the faces of maximum dimension of the mesh. This is relevant for meshes of dimensiondembedded ind+1dimensions as there is no way to tell which should be the orientation of the normals from the other quantities in the mesh.normals[f]gives the normal vector of face numberfof dimensiond=length(face_nodes)-1.
GalerkinToolkit.diameter — Method
GalerkinToolkit.domain — Function
GalerkinToolkit.domain — Method
GalerkinToolkit.each_face — Function
GalerkinToolkit.face_around — Function
face_around(x)Return an integer that allows to break ties when faces in x need to point to faces around of one dimension higher. Return nothing if x does not break such ties.
Note: This function will eventually return a vector of integers.
Level
Intermediate
GalerkinToolkit.face_dofs — Function
GalerkinToolkit.face_incidence — Function
face_incidence(x,d)GalerkinToolkit.face_nodes — Function
GalerkinToolkit.face_permutation_ids — Function
GalerkinToolkit.face_reference_id — Function
GalerkinToolkit.faces — Function
faces(x)Return the subset of face ids in mesh(x) of dimension num_dims(x) defining the domain x. This is effectively the map from domain face id to mesh face id.
See also inverse_faces.
Level
Intermediate
GalerkinToolkit.geometries — Function
geometries(x,d)
geometries(x,Val(d))Return a vector of domains representing the geometrical entities of x of dimension d. The returned domains and x are defined on the same mesh. That is, faces(geometries(x,1)[2]) are the face ids in mesh(x) representing the second edge of x.
Notation
geometries(x,Val(0)) are referred to as the vertices of x. geometries(x,Val(1)) are referred to as the edges of x. geometries(x,Val(d)) are referred to as the d-faces of x.
Level
Advanced
GalerkinToolkit.geometry_nodes — Function
GalerkinToolkit.geometry_own_dofs — Function
GalerkinToolkit.global_int_type — Method
global_int_type(options::Options)Return the type of the integers used to enumerate global quantities.
GalerkinToolkit.group_boundary_faces! — Function
GalerkinToolkit.group_faces — Function
group_faces(mesh)
group_faces(mesh,d)Return the dictionary containing the faces in each group in dimension d. If d is omitted, it returns the dictionaries for all dimensions in a vector. I.e., calling group_faces(mesh,d) is equivalent to group_faces(mesh)[d+1].
The faces of dimension d in group group are group_faces(mesh,d)[group_name], where group_name is a string with the group name. One can create new groups by adding new keys to these dictionaries as long as the key is not already present. Calling group_faces(mesh,d)[new_group_name] = faces_in_newgroup will add a new group to dimension d with name equal to the string new_group_name with faces in vector faces_in_newgroup.
See also group_names.
Level
Beginner
GalerkinToolkit.group_faces_in_dim! — Function
GalerkinToolkit.group_interior_faces! — Function
GalerkinToolkit.group_names — Function
group_names(mesh)
group_names(mesh,d)Return the names of the groups in dimension d. Calling group_names(mesh,d) is equivalent to group_names(mesh)[d+1].
See also group_faces.
Level
Beginner
GalerkinToolkit.int_type — Method
int_type(options::Options)Return the default integer type used in the computation except for reference and global quantities.
GalerkinToolkit.interior — Function
GalerkinToolkit.interior_nodes — Function
GalerkinToolkit.inverse_faces — Function
inverse_faces(x)Return the inverse integer mas of faces(x). This is effectively the map from mesh face id to domain face id. Mesh faces not present in the domain, receive an invalid index id.
See also faces.
Level
Intermediate
GalerkinToolkit.is_axis_aligned — Function
is_axis_aligned(x)True if x is a unit simplex or a unit cube.
Level
Advanced
GalerkinToolkit.is_boundary — Function
is_boundary(x)True if x represent an (internal) boundary. Faces in an internal boundary "point" to only of the two faces around of a dimension higher.
See also face_around.
Level
Intermediate
GalerkinToolkit.is_face_complex — Function
GalerkinToolkit.is_n_cube — Function
is_n_cube(x)True if x is a n-cube (hypercube).
See also is_simplex, is_unit_simplex, is_unit_n_cube.
Level
Intermediate
GalerkinToolkit.is_simplex — Function
is_simplex(x)True if x is a simplex.
See also is_n_cube, is_unit_simplex, is_unit_n_cube.
Level
Intermediate
GalerkinToolkit.is_unit_n_cube — Function
is_unit_n_cube(x)True if x is a unit n-cube.
See also is_n_cube, is_unit_simplex, is_simplex.
Level
Intermediate
GalerkinToolkit.is_unit_simplex — Function
is_unit_simplex(x)True if x is a unit simplex.
See also is_n_cube, is_unit_n_cube, is_simplex.
Level
Intermediate
GalerkinToolkit.is_unitary — Function
is_unitary(x)True bounding_box(x) coincides with a unit n-cube.
Level
Advanced
GalerkinToolkit.lagrange_space — Function
GalerkinToolkit.mesh — Method
GalerkinToolkit.mesh_from_gmsh — Function
mesh_from_gmsh(gmsh::Module;complexify=true)Create a mesh objects from the current state of the gmsh module. If complexify==true, the mesh will be completed with all low dimensional faces into a face complex.
See also mesh_from_msh and with_gmsh.
Level
Beginner
GalerkinToolkit.mesh_from_msh — Function
mesh_from_msh(msh_file;kwargs...)Create a mesh object from a .msh file found in path msh_file.
See also mesh_from_gmsh and with_gmsh.
Keyword arguments
complexify=true[optional]: Ifcomplexify==true, the mesh will be completed with all low dimensional faces into a face complex.renumber=true[optional]: Ifrenumber==true, thengmsh.model.mesh.renumberNodes()andgmsh.model.mesh.renumberElements()will be called.- Any other keyword argument will be passed to function
with_gmsh.
Level
Beginner
GalerkinToolkit.mesh_from_space — Function
mesh_from_space(space)Return the mesh induced by space. For instance, a (high order) Lagrange space can be interpreted as a mesh using this function.
Level
Advanced
GalerkinToolkit.moebius_strip — Method
GalerkinToolkit.node_coordinates — Function
node_coordinates(x)Return the vector of node coordinates associated with `x.
GalerkinToolkit.node_quadrature — Function
GalerkinToolkit.nodes — Function
GalerkinToolkit.normals — Function
GalerkinToolkit.num_ambient_dims — Function
num_ambient_dims(x)Return the ambient dimension where object x lives.
See also num_codims, num_dims.
Level
Beginner
GalerkinToolkit.num_codims — Function
num_codims(x)Return num_ambient_dims(x)-num_dims(x).
See also num_ambient_dims, num_dims.
Level
Beginner
GalerkinToolkit.num_dims — Function
GalerkinToolkit.num_dofs — Function
GalerkinToolkit.num_faces — Function
num_faces(x)
num_faces(x,d)Return the number of faces of dimension d in mesh(x). If d is omitted, return a vector with the number of faces in each dimension, starting from dimension 0 up to num_dims(x).
Level
Beginner
GalerkinToolkit.num_nodes — Function
GalerkinToolkit.num_points — Function
GalerkinToolkit.options — Method
options(;kwargs...) -> OptionsCreate an object representing the default options for the current simulation. This object can be used as an optional argument in several object constructors in GalerkinToolkit, such as the mesh constructors cartesian_mesh and mesh_from_msh. In this case, the computations using the generated mesh, will use the given options by default.
GalerkinToolkit.periodic_nodes — Function
GalerkinToolkit.push — Function
push(a,ai)Like push!, but creates a new object to store the result. This function is used to push to immutable collections such as tuples.
GalerkinToolkit.quadrature — Method
GalerkinToolkit.raviart_thomas_space — Function
GalerkinToolkit.real_type — Method
real_type(options::Options)Return the default real type used in the computation.
GalerkinToolkit.reference_int_type — Method
reference_int_type(options::Options)Return the type of the integers used to enumerate reference quantities.
GalerkinToolkit.reference_quadratures — Function
GalerkinToolkit.reference_spaces — Function
GalerkinToolkit.reference_topologies — Function
reference_topologies(topo)
reference_topologies(topo,d)
reference_topologies(topo,Val(d))Return the list (a vector or a tuple) of reference topologies in topo of dimension d. If the second argument is omitted, then the function returns a collection such that reference_topologies(topo)[d+1] is equivalent to reference_topologies(topo,Val(d)).
The face reference topology of face f of dimension d, is accessed as reference_topologies(topo,d)[r] with r=face_reference_id(topo,d)[f].
See also face_reference_id.
Level
Intermediate
GalerkinToolkit.simplexify — Function
simplexify(x)Convert x into a mesh made of simplex cells.
Level
Intermediate
GalerkinToolkit.skeleton — Function
GalerkinToolkit.tabulator — Method
GalerkinToolkit.topology — Method
GalerkinToolkit.unit_n_cube — Method
unit_n_cube(d)
unit_n_cube(Val(d))Return an object representing a unit d-cube.
GalerkinToolkit.unit_simplex — Method
unit_simplex(d)
unit_simplex(Val(d))Return an object representing a unit simplex of dimension d.
GalerkinToolkit.val_parameter — Method
val_parameter(a)For a::Val{A} it returns A. Otherwise, it returns a.
GalerkinToolkit.vertex_permutations — Function
vertex_permutations(x)Return a list of permutations representing the admissible re-labelings of the vertices of x.
Level
Advanced
GalerkinToolkit.weights — Function
GalerkinToolkit.with_gmsh — Function
with_gmsh(f[;options])A safe way of initialize and finalize the gmsh module. The given function is called f(gmsh) on the gmsh module after is has been initialized. The module is finalized automatically when the function returns.
The optional keyword argument options is a vector for pairs k=>v containing gmesh options. Each of these options are set with gmsh.option.setNumber(k,v) just after gmsh has been initialized.
Level
Beginner
GalerkinToolkit.with_mesh_partitioner — Function
with_mesh_partitioner(mesher[,partitioner];[parts])Generate a mesh calling mesher() partition it, and distribute it over the part ids in parts.
Arguments
- Function
mesher()should have no arguments and returns a sequential mesh object. This function is called only on one process. partitioner[optional]: A function that takes a graph encoded as a sparse matrix, and returns a vector containing the part id of each node in the graph. Defaults toMetis.partition.
Keyword arguments
parts[optional]: A vector containing the part indices1:PwherePis the number of parts in the data distribution. By default,Pis the number of MPI ranks and1:Pis distributed one item per rank.