Muscat.FE.FETools module¶
- CellDataToIntegrationPointsData(mesh, scalarFields, elementSet=None, relativeDimension=0, integrationRule=<Muscat.FE.IntegrationRules.MeshQuadrature object>)[source]¶
Change the representation of scalarFields from data constant by cell (elements) to data at integration points. (Lagrange isoparametric finite elements)
- Parameters:
mesh (Mesh) – mesh on which the function is applied
scalarFields (np.ndarray of size (nbe of fields, nbe of elements) or dict) – with “nbe of fields” as keys and np.ndarray of size (nbe of elements,) or floats as values fields whose representation is changed by the function
elementSet (elementSet defining the elements on which the function is) – computed. If None, takes all the elements of considered dimension
relativeDimension (int (0, -1 or -2)) – difference between the dimension of the elements on which the function is computed and the dimensionality of the mesh
- Returns:
of size (nbe of fields, numberOfIntegrationPoints)
- Return type:
np.ndarray
- ComputeGradPhiAtIntegPoint(mesh, elementSets=None, relativeDimension=0, integrationRule=<Muscat.FE.IntegrationRules.MeshQuadrature object>)[source]¶
Computes the components of the gradient of the shape functions on the integration points and the integration weights associated with the integration scheme
- Parameters:
mesh (Mesh) – mesh on which the function is applied
elementSets (list of strings) – element sets, support for the shape functions
relativeDimension (int) – 0, -1, or -2: the dimension of the element set relative to the dimension of the mesh
- Returns:
np.ndarray – integrationWeights
list (length dimensionality of the mesh) coo_matrix of size (NGauss, nbNodes) – gradPhiAtIntegPoint
- ComputeH10ScalarProductMatrix(mesh, numberOfComponents, integrationRule=<Muscat.FE.IntegrationRules.MeshQuadrature object>)[source]¶
Computes the H10 scalar product used to compute the correlations between the primal solution snapshots. The numberOfComponents depends on the solution type: 3 for solid mechanics in 3D, or 1 for thermal in any dimension (Lagrange isoparametric finite elements)
- Parameters:
mesh (Mesh) – mesh on which the function is applied
numberOfComponents (int) – the number of components of the primal variable snapshots
- Returns:
the sparse matrix of the H10 scalar product
- Return type:
scipy.sparse.csr_matrix
- ComputeIntegrationPointsTags(mesh, dimension=None, integrationRule=<Muscat.FE.IntegrationRules.MeshQuadrature object>)[source]¶
Returns a list of lists (of str) at each integration point (Lagrange isoparametric finite elements). Each list contains all the tags of the element of dimension “dimension” containg the considered integration points
- Parameters:
mesh (Mesh) – mesh on which the function is applied
elementSets (list if strings) – sets of elements on which the function is computed
- Returns:
of length numberOfIntegrationPoints
- Return type:
list of lists (of str)
- ComputeInterpolationMatrix_FE_GaussPoint(mesh, feSpace, integrationRule, feNumbering=None, ipNumbering=None, elementFilter=None)[source]¶
- ComputeJdetAtIntegPoint(mesh, elementSets=None, relativeDimension=0, integrationRule=<Muscat.FE.IntegrationRules.MeshQuadrature object>)[source]¶
Computes determinant of the Jacobian of the transformation of the transformation between the reference element and the mesh element, at the integration points. (Lagrange isoparametric finite elements)
- Parameters:
mesh (Mesh) – mesh on which the function is applied
elementSets (list if strings) – sets of elements on which the det-Jacobians are computed
relativeDimension (int (0, -1 or -2)) – difference between the dimension of the elements on which the function is computed and the dimensionality of the mesh
- Returns:
of size (NGauss,)
- Return type:
np.ndarray
- ComputeL2ScalarProductMatrix(mesh, numberOfComponents, elementFilter=None, integrationRule=<Muscat.FE.IntegrationRules.MeshQuadrature object>)[source]¶
Computes the L2 scalar product used to compute the correlations between the primal solution snapshots. The numberOfComponents depends on the solution type: 3 for solid mechanics in 3D, or 1 for thermal in any dimension. (Lagrange isoparametric finite elements)
- Parameters:
mesh (Mesh) – mesh on which the function is applied
numberOfComponents (int) – the number of components
elementFilter (ElementFilter, optional) – Filter on which the function is applied
- Returns:
the sparse matrix of the L2 scalar product
- Return type:
scipy.sparse.csr_matrix
- ComputeNormalsAtIntegPoint(mesh, elementSets, integrationRule=<Muscat.FE.IntegrationRules.MeshQuadrature object>)[source]¶
Computes the normals at the elements from the sets elementSets in the ambiant space (i.e. if mesh is of dimensionality 3, elementSets should contain 2D elements)
- Parameters:
mesh (Mesh) – mesh on which the function is applied
elementSets (list of strings) – sets of elements on which the function is computed
- Returns:
of size (dimensionality, NGauss)
- Return type:
np.ndarray
- ComputeNormalsAtPoints(mesh: Mesh, ef: ElementFilter | None = None) ndarray [source]¶
Compute the normal at the nodes selected by the element filter
- Parameters:
mesh (Mesh) – mesh to use, must have 2D element for a 3D mesh, or 1D elements for a 2D mesh
ef (ElementFilter, optional) – Selection of elements to compute the normal on , by default None
- Returns:
a (number of nodes, dimensionality) array with the normal at the nodes (zero filled)
- Return type:
np.ndarray
- ComputePhiAtIntegPoint(mesh, elementSets=None, relativeDimension=0, integrationRule=<Muscat.FE.IntegrationRules.MeshQuadrature object>)[source]¶
Computes the value of the finite element shape functions at the integration points. (Lagrange isoparametric finite elements)
- Parameters:
mesh (Mesh) – mesh on which the function is applied
elementSets (list if strings) – sets of elements on which the det-Jacobian are computed
relativeDimension (int (0, -1 or -2)) – difference between the dimension of the elements on which the function is computed and the dimensionality of the mesh
- Returns:
np.ndarray – of size (NGauss,)
coo_matrix – of size (NGauss, nbNodes)
- ComputePhiAtIntegPointFromElFilter(mesh, elFilter, space=<Muscat.FE.Spaces.FESpaces.FESpace object>, integrationRule=<Muscat.FE.IntegrationRules.MeshQuadrature object>)[source]¶
Computes the shape functions on the integration points and the integration weights associated with the integration scheme
- Parameters:
mesh (Mesh) –
elFilter (elementFilter) –
- Returns:
np.ndarray – integrationWeights
coo_matrix of size (NGauss, nbNodes) – phiAtIntegPointMatrix
- GetElementaryMatrixForFormulation(elemName, wform, unknownNames, space=Muscat.FE.Spaces.FESpaces.LagrangeSpaceP1)[source]¶
- IntegrationPointsToCellData(mesh, scalarFields, integrationRule=<Muscat.FE.IntegrationRules.MeshQuadrature object>)[source]¶
Change the representation of scalarFields from data at integration points to data constant by cell (elements) - taking the mean of values at the integration points in each elements. (Lagrange isoparametric finite elements)
- Parameters:
mesh (Mesh) – mesh on which the function is applied
scalarFields (np.ndarray of size (nbe of fields, nbe of elements) or dict) – with “nbe of fields” as keys and np.ndarray of size (nbe of elements,) as values fields whose representation in changed by the function
- Returns:
of size (nbe of elements,)
- Return type:
np.ndarray
- PrepareFEComputation(mesh, elementFilter=None, numberOfComponents=None, space=<Muscat.FE.Spaces.FESpaces.FESpace object>, integrationRule=<Muscat.FE.IntegrationRules.MeshQuadrature object>)[source]¶
Prepares a finite element computation with Lagrange isoparametric finite elements by computing FESpace, numberings, offset and NGauss
- Parameters:
mesh (Mesh) – mesh on which the function is applied
elementFilter (ElementFilter, optional) – Filter on which the function is applied
numberOfComponents (int) – the number of components
- Returns:
FESpace
numberings
offset (list of int)
NGauss (int)