genericROM.FE.FETools

CellDataToIntegrationPointsData(mesh, scalarFields, elementSet=None, relativeDimension=0)[source]

Change the representation of scalarFields from data constant by cell (elements) to data at integration points. (Lagrange isoparametric finite elements)

Parameters
  • mesh (BasicToolsUnstructuredMesh) – high-dimensional mesh

  • 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

  • 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)[source]

Computes the components of the gradient of the shape functions at the integration points and the integration weights (Lagrange isoparametric finite elements)

Parameters
  • mesh (BasicToolsUnstructuredMesh) – high-dimensional mesh

  • elementSets (list of strings) – sets of elements on which the function is 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 (numberOfIntegrationPoints,)

  • list – of length dimensionality of the mesh, of scipy.sparse.coo_matrix of size (numberOfIntegrationPoints, numberOfModes)

ComputeH10ScalarProductMatrix(mesh, numberOfComponents)[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 (BasicToolsUnstructuredMesh) – high-dimensional mesh

  • 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

ComputeIndicesOfIntegPointsPerMaterial(listOfTags, keysConstitutiveLaws)[source]

From the list of lists containing all the tags of the element containing the integration points of the considered mesh, and the set of tags on which constitutive laws are defined, returns a dictionary containing the ranks of integration points for each tag (inverse the information contained in listOfTags).

Parameters
  • listOfTags (list of lists (of str)) – of length numberOfIntegrationPoints, containing all the tags of the element containing the integration points

  • keysConstitutiveLaws (set (of str)) – set of tags on which constitutive laws are defined

Returns

dict (str – dictionary with element tags (str) as keys and ranks of integration points (1D np.ndarray) as values

Return type

np.ndarray)

ComputeIntegrationPointsTags(mesh, dimension=None)[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” containing the considered integration points.

Parameters
  • mesh (BasicToolsUnstructuredMesh) – high-dimensional mesh

  • elementSets (list of strings) – sets of elements on which the function is computed

Returns

of length numberOfIntegrationPoints

Return type

list of lists (of str)

ComputeJdetAtIntegPoint(mesh, elementSets=None, relativeDimension=0)[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 (BasicToolsUnstructuredMesh) – high-dimensional mesh

  • elementSets (list of strings) – sets of elements on which the function is 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 (numberOfIntegrationPoints,)

Return type

np.ndarray

ComputeL2ScalarProducMatrix(mesh, numberOfComponents)[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 (BasicToolsUnstructuredMesh) – high-dimensional mesh

  • numberOfComponents (int) – the number of components of the primal variable snapshots

Returns

the sparse matrix of the L2 scalar product

Return type

scipy.sparse.csr

ComputeMaterialKeyPerIntegrationPoint(listOfTags, keysConstitutiveLaws)[source]

From the list of lists containing all the tags of the element containing the integration points of the considered mesh, and the set of tags on which constitutive laws are defined, returns the list containing the unique tag defining the constitutive laws key at each integration point.

Parameters
  • listOfTags (list of lists (of str)) – of length numberOfIntegrationPoints, containing all the tags of the element containing the integration points

  • keysConstitutiveLaws (set (of str)) – set of tags on which constitutive laws are defined

Returns

of length numberOfIntegrationPoints, containing the unique tag defining the constitutive laws key at each integration point

Return type

list (of str)

ComputeNormalsAtIntegPoint(mesh, elementSets=None)[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 (BasicToolsUnstructuredMesh) – high-dimensional mesh

  • elementSets (list of strings) – sets of elements on which the function is computed

Returns

of size (dimensionality, numberOfIntegrationPoints)

Return type

np.ndarray

ComputeNumberOfIntegrationPoints(mesh)[source]

Computes the number of integration points. (Lagrange isoparametric finite elements)

Parameters

mesh (BasicToolsUnstructuredMesh) – high-dimensional mesh

Returns

numberOfIntegrationPoints

Return type

int

ComputePhiAtIntegPoint(mesh, elementSets=None, relativeDimension=0)[source]

Computes the value of the finite element shape functions at the integration points and the integration weights (Lagrange isoparametric finite elements)

Parameters
  • mesh (BasicToolsUnstructuredMesh) – high-dimensional mesh

  • elementSets (list of strings) – sets of elements on which the function is 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 (numberOfIntegrationPoints,)

  • scipy.sparse.coo_matrix – of size (numberOfIntegrationPoints, numberOfModes)

ConvertMeshToUnstructuredMesh(mesh)[source]
IntegrationPointsToCellData(mesh, scalarFields)[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 (BasicToolsUnstructuredMesh) – high-dimensional mesh

  • 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