Muscat.FE.Fields.FieldTools module

CheckIntegrity(GUI: bool = False)[source]
CheckIntegrity_ElementWiseFEToFETransferOp(GUI: bool = False)[source]
CheckIntegrity_FieldsAtIp(GUI: bool = False)[source]
CheckIntegrity_Legacy(GUI: bool = False)[source]
ComputeTransferOp(field1: FEField, field2: FEField, force: bool = False) Tuple[np.ndarray, np.ndarray][source]

Compute the transfer Operator from field1 to field2 Both fields must be compatibles fields: same mesh, same space but with different numberings

  • example of use:

    lIndex,rIndex = ComputeTransferOp(field1,field2):

    field2.data[lIndex] = field1.data[rIndex]

Parameters:
  • field1 (FEField) – field from where the information is extracted

  • field2 (FEField) – field to receive the information

  • force (bool, optional) – true to bypass the compatibility check, by default False

Returns:

index vector for field2, index vector for field1

Return type:

Tuple(np.ndarray,np.ndarray)

CreateFieldFromDescription(mesh: Mesh, fieldDefinition: Iterable[DescriptionUnit], fieldType: str = 'FE') FEField | IPField[source]

Create a FEField from a field description

Parameters:
  • mesh (Mesh) –

  • fieldDefinition (List[Tuple[Union[ElementFilter,NodeFilter],Union[float,Callable[[np.ndarray],float]]]]) – A field definition is a list of Tuples (with 2 element each ). The first element of the tuple is a ElementFilter or a NodeFilter The second element of the tuple is a float or float returning callable (argument is a point in the space)

  • fieldType (str, optional) – type of field created FEField or IPField (“FE”, “FE-P0”, “IP”) by default “FE” (isoparametric)

Returns:

created Field with values coming from the field description

Return type:

Union[FEField,IPField]

ElemFieldsToFEField(mesh: Mesh, elemFields: Dict[str, ArrayLike] | None = None, copy=False) Dict[str, FEField][source]

Create FEField(P0) from the elements field data. if elemFields is None the mesh.elemFields is used

Parameters:
  • mesh (Mesh) – The support for the FEFields

  • elemFields (Dict[str,ArrayLike], optional) – a dict containing the elements fields to be converted to FEFields if elemFields is None the mesh.elemFields is used, by default None

Returns:

A dictionary the keys are field names and the values are the FEFields

Return type:

Dict[str,FEField]

ElementWiseFEToFETransferOp(originSpace: FESpace, targetSpace: FESpace) Dict[str, np.ndarray][source]

Generate transfer operator (element wise) to pass information from a FE Field to a different FE field. This is done by solving a least-squares problem.

Parameters:
  • originSpace (FESpace) – The initial space

  • targetSpace (FESpace) – The target space

Returns:

the operator to pass the data from the initial space to the target FEField

Return type:

Dict[str,np.ndarray]

ElementWiseIpToFETransferOp(integrationRule: MeshQuadrature, space: FESpace) Dict[str, np.ndarray][source]

Generate transfer operator (element wise) to pass information from integration points to a FE field. This is done by solving a least-squares problem.

Parameters:
  • integrationRule (IntegrationRuleType) – The integration rule of the original data

  • space (FESpace) – The target space

Returns:

the operator to pass the data from the integration points to the FEField

Return type:

Dict[str,np.ndarray]

FEFieldsDataToVector(listOfFields: List[FEField], outVector: ArrayLike | None = None, copy: bool = False) np.ndarray[source]

Put FEFields data into a vector format

Parameters:
  • listOfFields (List[FEField]) – list of FEFields (the order is important)

  • outVector (Optional[ArrayLike], optional) – if provided the output vector, by default None

Returns:

_description_

Return type:

np.ndarray

FieldsAtIp(listOfFields: List[FEField | IPField | RestrictedIPField], rule, elementFilter: ElementFilter | None = None) List[IntegrationPointWrapper | IPField | RestrictedIPField][source]

Convert a list of Field (FEFields,IPField,RestrictedIPField) to a list of IPFields

Parameters:
  • listOfFields (List[Union[IPField,RestrictedIPField]]) – the list of fields to be treated

  • rule (_type_) – the integration rule to be used for the evaluation

  • elementFilter (Optional[ElementFilter], optional) – the filter to select the element to treat, by default None

Returns:

the list of IPFilters

Return type:

List[Union[IntegrationPointWrapper,IPField,RestrictedIPField]]

Raises:

Exception – in the case of a internal error

class FieldsEvaluator(fields=None)[source]

Bases: object

Helper to evaluate expression using FEField and IPFields. The user can add fields and constants. Then this class can be used to compute an expression (given using an callable object). The callable can capture some of the argument but the last must be **args (so extra argument are ignored).

AddConstant(name: str, val: float64)[source]

Add a constant value to the almanac

Parameters:
  • name (str) – name of the value

  • val (MuscatFloat) – Value

AddField(field: FieldBase)[source]

Add a field to the internal almanac

Parameters:

field (FieldBase) – A field to be added to the almanac

Compute(func: Callable, on: str, useSympy: bool = False)[source]

Compute the function func using the user fields on the support “on”

Parameters:
  • func (Callable) – the function to evaluate

  • on (str) – [“IPField”, “Centroids”, “FEField”, “Nodes”]

  • useSympy (bool, optional) – if true the function is compiled using sympy , by default False

Returns:

The callable evaluated

Return type:

Any

GetFieldsAt(on: str) Dict[source]

Get the internal representations of the user fields ant different support

Parameters:

on (str) – [“IPField”, “Centroids”, “FEField”, “Nodes”]

Returns:

A dictionary containing the different representation of the user fields

Return type:

Dict

Raises:

ValueError – In the case the “on” parameter is not valid

GetOptimizedFunction(func)[source]
Update(what: str = 'all')[source]

Function to update the field calculated at different support (a FEField can be internally evaluated at the integration points) The objective of this function is to update this evaluation

Parameters:

what (str, optional) – [“all”, “IPField”, “Centroids”], by default “all”

class FieldsMeshTransportation(oldMesh, newMesh)[source]

Bases: object

Class to help the transfer of field (FEField and IPField) between old and new meshes

GetNumbering(mesh, space, fromConnectivity=False, discontinuous=False)[source]
ResetCacheData()[source]
TransportFieldToNewMesh(inputField: FEField | IPField) FEField | IPField[source]

Function to create a Field on the new mesh, the new mesh must be a transformation of the mesh in the inputField. This means the new mesh originalIds (for nodes and elements) must be with respect to the mesh of the inputField

Parameters:

inputField (FEField or a IPField) – the FEField to be transported to the new mesh

Returns:

a FEField or a IPField (depending on the inputField Type) on the new mesh filled with the values of the input field

Return type:

FEField or IPField

TransportFieldToOldMesh(inputField: FEField | IPField, fillValue: float64 = 0.0) FEField | IPField[source]

Function to transport a FEField on the old mesh, the inputField mesh must be a transformation of the old mesh. This means the inputField mesh originalIds (for nodes and elements) must be with respect to the old mesh.

Parameters:
  • inputField (FEField or a IPField) – the field to be transported

  • fillValue (MuscatFloat, optional) – Value to fill the values of the dofs not present in the inFEField , by default 0.

Returns:

a FEField or a IPField (depending on the inputField Type) on the old mesh filled with the values of the input field

Return type:

FEField or IPField

FillFEField(field: FEField, fieldDefinition) None[source]

Fill a FEField using a definition

Parameters:
  • field (FEField) – FEField to treat

  • fieldDefinition (List[Tuple[Union[ElementFilter, NodeFilter], Union[float, Callable[[np.ndarray], float]]]]) – A field definition is a list of Tuples (with 2 element each ). The first element of the tuple is a ElementFilter or a NodeFilter The second element of the tuple is a float or float returning callable (argument is a point in the space)

FillIPField(field: IPField, fieldDefinition) None[source]

Fill a IPField using a definition

Parameters:
  • field (IPField) – IPField to treat

  • fieldDefinition (List[Tuple[ElementFilter,Union[float,Callable[[np.ndarray],float]]]]) – A field definition is a list of Tuples (with 2 element each ). The first element of the tuple is a ElementFilter or a NodeFilter The second element of the tuple is a float or float returning callable (argument is a point in the space)

GetCellRepresentation(listOfFEFields: List[FEField], fillValue: MuscatFloat = 0.0, method='mean') np.ndarray[source]

Get a np.ndarray compatible with the points of the mesh for a list of FEFields. Each field in the list is treated as a scalar component for the output field

Parameters:
  • listOfFEFields (list) – list of FEFields to extract the information

  • fillValue (float) – value to use in the case some field are not defined over all the nodes, by default 0.

  • method (str | mean) – Read FEField.GetCellRepresentation for more information, by default “mean”.

Returns:

Array of size (number of elements, number of fields)

Return type:

np.array

GetPointRepresentation(listOfFEFields: List[FEField], fillValue: MuscatFloat = 0.0) np.ndarray[source]

Get a np.ndarray compatible with the points of the mesh for a list of FEFields. Each field in the list is treated as a scalar component for the output field

Parameters:
  • listOfFEFields (List[FEField]) – list of FEFields to extract the information

  • fillValue (MuscatFloat, optional) – value to use in the case some field are not defined over all the nodes, by default 0.

Returns:

the output and nd array of size (nb points, nb of FEField)

Return type:

np.ndarray

GetTransferOpToIPField(inputField: FEField, ruleName: str | None = None, rule: MeshQuadrature | None = None, der: int = -1, elementFilter: ElementFilter | None = None) FeToIPOp[source]
Compute the transfer operator (.dot()) to construct a IPField (of RestrictedIPField)

IPField = FeToIPOp.dot(FEField)

Parameters:
  • inputField (FEField) – the input FEField

  • ruleName (Optional[str], optional) – the ruleName of the integration rule to use, by default None (see Muscat.FE.IntegrationsRules:GetRule for more info)

  • rule (Optional[Tuple[np.ndarray,np.ndarray]], optional) – the integration rule to use, by default None (see Muscat.FE.IntegrationsRules:GetRule for more info)

  • der (int, optional) – the coordinate to be derived -1 no derivations only evaluation at IP [0,1,2] the coordinate number to compute derivative of the FEField by default -1

  • elementFilter (Optional[ElementFilter], optional) – An element filter to restrict the operator. In this case a RestrictedIPField is generated _description_, by default None

Returns:

An instance of an object with the .dot operator

Return type:

FeToIPOp

class IntegrationPointWrapper(field: FEField, rule, elementFilter: ElementFilter | None = None)[source]

Bases: FieldBase

Class to generate a FEField at the integration points

Two important function are available : diff and GetIpField.

Parameters:
  • field (FEField) – the field to evaluate

  • rule (_type_) – the integration rule

  • elementFilter (ElementFilter, optional) – in this case a RestrictedIPField is generated, by default None

GetIpField() IPField | RestrictedIPField[source]

Get the integration point representation of the field

Returns:

Depending on the presence of an elementFilter

Return type:

Union[IPField, RestrictedIPField]

ResetCache()[source]
binaryOp(other, op)[source]
property data
diff(compName: str) IPField | RestrictedIPField[source]

Get the IPField representation of this field with a derivation in the direction of compName

Parameters:

compName (str) – compName are available in : from Muscat.FE.SymWeakForm import space

Returns:

Depending on the presence of an elementFilter

Return type:

(IPField | RestrictedIPField)

property name
unaryOp(op)[source]
Maximum(A, B)[source]
Minimum(A, B)[source]
NodeFieldToFEField(mesh: Mesh, nodeFields: Dict[str, ArrayLike] | None = None, copy: bool = False) Dict[str, FEField][source]

Create FEField(P isoparametric) from the node field data. if nodesFields is None the mesh.nodeFields is used

Parameters:
  • mesh (Mesh) – The support for the FEFields

  • nodeFields (Dict[str,ArrayLike], optional) – the dictionary containing the nodes fields to be converted to FEFields, if None the mesh.nodeFields is used, by default None

Returns:

A dictionary the keys are field names and the values are the FEFields

Return type:

Dict[str,FEField]

TransferFEFieldToIPField(inFEField: FEField, ruleName: str | None = None, rule=None, der: int = -1, elementFilter: ElementFilter | None = None, op=None) IPField | RestrictedIPField[source]

Transfer FEField to a IPField

Parameters:
  • inFEField (FEField) – the input FEField

  • ruleName (Optional[str], optional) – the ruleName of the integration rule to use, by default None (see Muscat.FE.IntegrationsRules:GetRule for more info)

  • rule (Optional[Tuple[np.ndarray,np.ndarray]], optional) – the integration rule to use, by default None (see Muscat.FE.IntegrationsRules:GetRule for more info)

  • der (int, optional) – the coordinate to be derived -1 no derivations only evaluation at IP [0,1,2] the coordinate number to compute derivative of the FEField by default -1

  • elementFilter (Optional[ElementFilter], optional) – An element filter to restrict the operator. In this case a RestrictedIPField is generated _description_, by default None

  • op (_type_, optional) – an object returned by the GetTransferOpToIPField function if not given then this op is calculated internally by default None

Returns:

The field evaluated at the integration points

Return type:

Union[IPField,RestrictedIPField]

TransferPosToIPField(mesh: Mesh, ruleName: str | None = None, rule: Tuple[np.ndarray, np.ndarray] | None = None, elementFilter: ElementFilter | None = None) List[IPField][source]

Create (2 or 3) IPFields with the values of the space coordinates

Parameters:
  • mesh (Mesh) –

  • ruleName (str, optional) – The Integration rule name, by default None (“LagrangeIsoParamQuadrature”)

  • rule (Tuple[np.ndarray,np.ndarray], optional) – The Integration Rule, by default None (“LagrangeIsoParamQuadrature”)

  • elementFilter (ElementFilter, optional) – the zone where the transfer must be applied, by default None (all the elements)

Returns:

The IPFields containing the coordinates

Return type:

List[IPField]

VectorToFEFieldsData(inVector: ArrayLike, listOfFields: List[FEField])[source]

Put vector data in FEFields

Parameters:
  • inVector (ArrayLike) – vector with the data

  • listOfFields (List[FEField]) – list of FEFields to store the data