Elements
Relevant module: edelweissfe.config.elementlibrary
EdelweissFE currently supports finite element implementations provided by the Marmot library. In future, elements by other providers or elements directly implemented in EdelweissFE may be added here.
*element, type=C3D8, provider=marmot
** el_label, node1, node2, node3, node4, ...
1000, 1, 2, 3, 4, ...
- edelweissfe.config.elementlibrary.getElementClass(provider=None)[source]
Get the class type of the requested element provider.
- Parameters
provider (str) – The name of the element provider ot load.
- Returns
The element provider class type.
- Return type
type
Provider marmot
Relevant module: edelweissfe.elements.marmotelement.element
- class edelweissfe.elements.marmotelement.element.MarmotElementWrapper
This element serves as a wrapper for MarmotElements.
For the documentation of MarmotElements, please refer to Marmot.
- Parameters
elementType – The Marmot element which should be represented, e.g., CPE4.
elNumber – The number of the element.
- Attributes
- dofIndicesPermutation
- elNumber
- elType
- ensightType
- fields
- hasMaterial
- nDof
- nNodes
- nSpatialDimensions
- nodeCoordinates
- nodes
- visualizationNodes
Methods
Accept the computed state (in nonlinear iteration schemes).
computeBodyForce
(P, K, load, U, time, dTime)Evaluate residual and stiffness for given time, field, and field increment due to a volume load.
computeDistributedLoad
(loadType, P, K, ...)Evaluate residual and stiffness for given time, field, and field increment due to a surface load.
computeYourself
(Ke, Pe, U, dU, time, dTime)Evaluate residual and stiffness for given time, field, and field increment.
Compute the underlying MarmotElement centroid coordinates.
Compute the underlying MarmotElement qp coordinates.
Compute the underlying MarmotElement qp coordinates.
getResultArray
(result, quadraturePoint[, ...])Get the array of a result, possibly as a persistent view which is continiously updated by the underlying MarmotElement.
Let the underlying MarmotElement initialize itself
Reset to the last valid state.
setInitialCondition
(stateType, values)Assign initial conditions to the underlying Marmot element
setMaterial
(materialName, materialProperties)Assign a material and material properties to the underlying MarmotElement.
setNodes
(nodes)Assign the nodes coordinates to the underyling MarmotElement
setProperties
(elementProperties)Assign a set of properties to the underyling MarmotElement
- acceptLastState()
Accept the computed state (in nonlinear iteration schemes).
- computeBodyForce(P, K, load, U, time, dTime)
Evaluate residual and stiffness for given time, field, and field increment due to a volume load.
- computeDistributedLoad(loadType, P, K, faceID, load, U, time, dTime)
Evaluate residual and stiffness for given time, field, and field increment due to a surface load.
- computeYourself(Ke, Pe, U, dU, time, dTime)
Evaluate residual and stiffness for given time, field, and field increment.
- getCoordinatesAtCenter()
Compute the underlying MarmotElement centroid coordinates.
- getCoordinatesAtQuadraturePoints()
Compute the underlying MarmotElement qp coordinates.
- getNumberOfQuadraturePoints()
Compute the underlying MarmotElement qp coordinates.
- getResultArray(result, quadraturePoint, getPersistentView=True)
Get the array of a result, possibly as a persistent view which is continiously updated by the underlying MarmotElement.
- initializeElement()
Let the underlying MarmotElement initialize itself
- resetToLastValidState()
Reset to the last valid state.
- setInitialCondition(stateType, values)
Assign initial conditions to the underlying Marmot element
- setMaterial(materialName, materialProperties)
Assign a material and material properties to the underlying MarmotElement. Furthermore, create two sets of state vars:
the actual set,
and a temporary set for backup in nonlinear iteration schemes.
- setNodes(nodes)
Assign the nodes coordinates to the underyling MarmotElement
- setProperties(elementProperties)
Assign a set of properties to the underyling MarmotElement
Provider marmotsingleqpelement
Relevant module: edelweissfe.elements.marmotsingleqpelement.element
- class edelweissfe.elements.marmotsingleqpelement.element.MarmotMaterialWrappingElement(materialType, elNumber)[source]
This element serves as a wrapper for MarmotMaterials, cf. Marmot.
It has a single quadrature point, and one (dummy) node. For interfacing with specific Marmot materials, specialized material wrappers are used.
The element allows to run quadrature point simulations investigating materials for development purposes.
- Parameters
materialType (str) – The Marmot material class which should be represented, e.g., MarmotMaterialHypoElastic.
elNumber (int) – The number of the element.
- Attributes
dofIndicesPermutation
If the provides computes residual vectors and stiffness matrices not nodewise, but e.g., fieldwise, this permutation pattern is used to aggregate all entries in order to resemble the defined fields nodewise.
elNumber
The unique number of this element
ensightType
The shape of the element in Ensight Gold notation.
fields
The list of fields per node which this entity couples.
hasMaterial
Flag to check if a material was assigned to this element.
nDof
The total number of degrees of freedom this entity has.
nNodes
The number of nodes this entity couples.
- nSpatialDimensions
nodes
The list of nodes this currently entity holds.
visualizationNodes
The nodes for visualization.
Methods
Accept the computed state.
computeBodyForce
(P, K, load, U, time, dTime)Not implemented for this wrapper.
computeDistributedLoad
(loadType, P, K, ...)Not implemented for this wrapper.
computeYourself
(Ke, Pe, U, dU, time, dTime)Evaluate the residual and stiffness for given time, field, and field increment due to a surface load.
Return the only node's coordinates.
Return the only qp's coordinates.
Return the only qp's coordinates.
getResultArray
(result, quadraturePoint[, ...])Get the array of a result, possibly as a persistent view which is continiously updated by the element.
Not used by this wrapper
Rest to the last valid state.
setInitialCondition
(stateType, values)Assign initial conditions.
setMaterial
(materialName, materialProperties)Assign a material and material properties to the underlying Wrapper.
setNodes
(nodes)Assign the nodes.
setProperties
(elementProperties)Not used by this wrapper
- computeDistributedLoad(loadType, P, K, faceID, load, U, time, dTime)[source]
Not implemented for this wrapper.
- computeYourself(Ke, Pe, U, dU, time, dTime)[source]
Evaluate the residual and stiffness for given time, field, and field increment due to a surface load.
- Parameters
P – The external load vector to be defined.
K – The stiffness matrix to be defined.
U – The current solution vector.
dU – The current solution vector increment.
time – Array of step time and total time.
dTime – The time increment.
- getCoordinatesAtCenter()[source]
Return the only node’s coordinates.
- Returns
The node’s coordinates.
- Return type
np.ndarray
- getCoordinatesAtQuadraturePoints()[source]
Return the only qp’s coordinates.
- Returns
The qp’s coordinates.
- Return type
np.ndarray
- getNumberOfQuadraturePoints()[source]
Return the only qp’s coordinates.
- Returns
The qp’s coordinates.
- Return type
np.ndarray
- getResultArray(result, quadraturePoint, getPersistentView=True)[source]
Get the array of a result, possibly as a persistent view which is continiously updated by the element.
- Parameters
result (str) – The name of the result.
quadraturePoint (int) – The number of the quadrature point.
getPersistentView (bool) – If true, the returned array should be continiously updated by the element.
- Returns
The result.
- Return type
np.ndarray
- setInitialCondition(stateType, values)[source]
Assign initial conditions.
- Parameters
stateType – The type of initial state.
values – The numpy array describing the initial state.
- setMaterial(materialName, materialProperties)[source]
Assign a material and material properties to the underlying Wrapper. Furthermore, create two sets of state vars:
the actual set,
and a temporary set for backup in nonlinear iteration schemes
- Parameters
materialName (str) – The name of the requested material.
materialProperties (ndarray) – The properties for he requested material.
- setNodes(nodes)[source]
Assign the nodes.
Only the first node is considered.
- Parameters
nodes (list) – The list of node instances.
- property dofIndicesPermutation
If the provides computes residual vectors and stiffness matrices not nodewise, but e.g., fieldwise, this permutation pattern is used to aggregate all entries in order to resemble the defined fields nodewise.
For performance reasons, this is computed once and stored for later use.
- property elNumber
The unique number of this element
- property ensightType
The shape of the element in Ensight Gold notation. Valid types are: point bar2 bar3 tria3 tria6 quad4 quad8 tetra4 tetra10 pyramid5 pyramid13 penta6 penta15 hexa8 hexa20 nsided nfaced
- property fields
The list of fields per node which this entity couples.
- property hasMaterial
Flag to check if a material was assigned to this element.
- property nDof
The total number of degrees of freedom this entity has. Not that this information is redundent, as it results from the list of fields on nodes, but we keep this redundant information for the sake of performance.
- property nNodes
The number of nodes this entity couples.
- property nodes
The list of nodes this currently entity holds.
- property visualizationNodes: list[edelweissfe.points.node.Node]
The nodes for visualization. Commonly, these are the same as the nodes of the entity. However, in some cases, the visualization nodes are different from the nodes of the entity, e.g., in case of mixed formulations.
testfiles/QPMarmotMaterialHypoElastic/test.inp
*material, name=LinearElastic, id=myMaterial
**Isotropic
**E | nu |
1.8e4, 0.22
*section, name=section1, material=myMaterial, type=solid
all
*node
** this is our single dummy node, which carries the strain 'field'
1, 0, 0, 0
*element, provider=MarmotSingleQpElement, type=MarmotMaterialHypoElastic
** This is our mateiral wrapping element. It needs a dummy node to be attached
** to a spatial position.
** The type of the element actually specifies the type of requested MarmotMaterial
1, 1
*job, name=mySingleElementJob, domain=3d
*solver, solver=NIST, name=theSolver
*step, solver=theSolver, maxInc=1e0, minInc=1e0, maxNumInc=100, maxIter=25, stepLength=1
setinitialconditions, property=characteristic element length, values='100,'
dirichlet, name=prescribed_strain, nSet=all, field=strain symmetric, components='[ 0, x, x, 0, 0, 0 ]',
nodeforces, name=prescribed_stress, nSet=all, field=strain symmetric, components='[ x, 0, 0, x, x, x ]'
*step, solver=theSolver, maxInc=1e-1, minInc=1e0, maxNumInc=10, maxIter=25, stepLength=1
dirichlet, name=prescribed_strain, components='[ -0.2, x, x, 0, 0, 0 ]'
nodeforces, name=prescribed_stress, components='[ x, .1, .2, x, x, x ]'
*fieldOutput
create=perElement, name=stress, elSet=all, result=stress, quadraturePoint=0, saveHistory=True
create=perElement, name=strain, elSet=all, result=strain, quadraturePoint=0, saveHistory=True
create=perElement, name=dStress_dStrain, elSet=all, result=dStress_dStrain, quadraturePoint=0, saveHistory=True, f(x)="np.reshape(x, (6,-1), order='F')"
*output, type=monitor
fieldOutput=stress, f(x)="'\n' + np.array2string(x, formatter={'float_kind':'{:>8.2f}'.format})"
fieldOutput=strain, f(x)="'\n' + np.array2string(x, formatter={'float_kind':'{:>8.2f}'.format})"
fieldOutput=dStress_dStrain, f(x)="'\n' + np.array2string(x, formatter={'float_kind':'{:>8.2f}'.format})"
Provider displacementelement
Relevant module: edelweissfe.elements.displacementelement.element
- class edelweissfe.elements.displacementelement.element.DisplacementElement(elementType, elNumber)[source]
Finite elements in EdelweissFE should be derived from this base class in order to follow the general interface.
EdelweissFE expects the layout of the internal and external load vectors, P, PExt, (and the stiffness) to be of the form
[ node 1 - dofs field 1, node 1 - dofs field 2, node 1 - ... , node 1 - dofs field n, node 2 - dofs field 1, ... , node N - dofs field n].
- Parameters
elementType (str) – A string identifying the requested element formulation.
elNumber (int) – A unique integer label used for all kinds of purposes.
- Attributes
dofIndicesPermutation
The permutation pattern for the residual vector and the stiffness matrix to aggregate all entries in order to resemble the defined fields nodewise.
elNumber
The unique number of this element
ensightType
The shape of the element in Ensight Gold notation.
fields
The list of fields per nodes.
hasMaterial
Flag to check if a material was assigned to this element.
nDof
The total number of degrees of freedom this element has
nNodes
The number of nodes this element requires
nodes
The list of nodes this element holds
visualizationNodes
The nodes for visualization.
Methods
Accept the computed state (in nonlinear iteration schemes).
computeBodyForce
(P, K, load, U, time, dTime)Evaluate residual and stiffness for given time, field, and field increment due to a body force load.
computeDistributedLoad
(loadType, P, K, ...)Evaluate residual and stiffness for given time, field, and field increment due to a surface load.
computeYourself
(K_, P, U, dU, time, dTime)Evaluate the residual and stiffness matrix for given time, field, and field increment due to a displacement or load.
Compute the underlying MarmotElement centroid coordinates.
Compute the underlying MarmotElement qp coordinates.
Get the number of Quadrature points the element has.
getResultArray
(result, quadraturePoint[, ...])Get the array of a result, possibly as a persistent view which is continiously updated by the element.
Initalize the element to be ready for computing.
Reset to the last valid state.
setInitialCondition
(stateType, values)Assign initial conditions.
setMaterial
(material)Assign a material.
setNodes
(nodes)Assign the nodes to the element.
setProperties
(elementProperties)Assign a set of properties to the element.
- computeBodyForce(P, K, load, U, time, dTime)[source]
Evaluate residual and stiffness for given time, field, and field increment due to a body force load.
- Parameters
P (ndarray) – The external load vector to be defined.
K (ndarray) – The stiffness matrix to be defined.
load (ndarray) – The magnitude (or vector) describing the load.
U (ndarray) – The current solution vector.
time (ndarray) – Array of step time and total time.
dTime (float) – The time increment.
- computeDistributedLoad(loadType, P, K, faceID, load, U, time, dTime)[source]
Evaluate residual and stiffness for given time, field, and field increment due to a surface load.
- Parameters
loadType (str) – The type of load.
P (ndarray) – The external load vector to be defined.
K (ndarray) – The stiffness matrix to be defined.
faceID (int) – The number of the elements face this load acts on.
load (ndarray) – The magnitude (or vector) describing the load.
U (ndarray) – The current solution vector.
time (ndarray) – Array of step time and total time.
dTime (float) – The time increment.
- computeYourself(K_, P, U, dU, time, dTime)[source]
Evaluate the residual and stiffness matrix for given time, field, and field increment due to a displacement or load.
- Parameters
P (ndarray) – The external load vector gets calculated.
K – The stiffness matrix gets calculated.
U (ndarray) – The current solution vector.
dU (ndarray) – The current solution vector increment.
time (ndarray) – Array of step time and total time.
dTime (float) – The time increment.
K_ (ndarray) –
- getCoordinatesAtCenter()[source]
Compute the underlying MarmotElement centroid coordinates.
- Returns
The element’s central coordinates.
- Return type
np.ndarray
- getCoordinatesAtQuadraturePoints()[source]
Compute the underlying MarmotElement qp coordinates.
- Returns
The element’s qp coordinates.
- Return type
np.ndarray
- getNumberOfQuadraturePoints()[source]
Get the number of Quadrature points the element has.
- Returns
The number of Quadrature points.
- Return type
nInt
- getResultArray(result, quadraturePoint, getPersistentView=True)[source]
Get the array of a result, possibly as a persistent view which is continiously updated by the element.
- Parameters
result (str) – The name of the result.
quadraturePoint (int) – The number of the quadrature point.
getPersistentView (bool) – If true, the returned array should be continiously updated by the element.
- Returns
The result.
- Return type
np.ndarray
- setInitialCondition(stateType, values)[source]
Assign initial conditions.
- Parameters
stateType (str) – The type of initial state.
values (ndarray) – The numpy array describing the initial state.
- setMaterial(material)[source]
Assign a material.
- Parameters
material (type) – An initalized instance of a material.
- setNodes(nodes)[source]
Assign the nodes to the element.
- Parameters
nodes (list[edelweissfe.points.node.Node]) – A list of nodes.
- setProperties(elementProperties)[source]
Assign a set of properties to the element.
- Parameters
elementProperties (ndarray) – A numpy array containing the element properties.
- thickness
Thickness of 2D elements.
- property dofIndicesPermutation: ndarray
The permutation pattern for the residual vector and the stiffness matrix to aggregate all entries in order to resemble the defined fields nodewise. In this case it stays the same because we use the nodes exactly like they are.
- property elNumber: int
The unique number of this element
- property ensightType: str
The shape of the element in Ensight Gold notation.
- property fields: list[list[str]]
The list of fields per nodes.
- property hasMaterial: str
Flag to check if a material was assigned to this element.
- property nDof: int
The total number of degrees of freedom this element has
- property nNodes: int
The number of nodes this element requires
- property nodes: int
The list of nodes this element holds
- property visualizationNodes: str
The nodes for visualization.
The elementType definition has to be of the following form:
Elements
One of the following element types needs to be included in the definition (elementType[0:5]
).
Quad4 - quadrilateral 2D element with 4 nodes.
Quad8 - quadrilateral 2D element with 8 nodes.
Hexa8 - hexahedron 3D element with 8 nodes.
additional Parameters
The following optional Parameters are also included in the element type definition.
R - reduced integration for element, in
elementType[5]
, currently possible for Q8 and Q4.E - extended integration for element, in
elementType[5]
, currently possible for Q4.PE - use plane strain for 2D elements, in
elementType[6:8]
orelementType[5:7]
.PS - use plane stress for 2D elements, in
elementType[6:8]
orelementType[5:7]
.
Note
If R or E is not given by the user, we assume regular integration.
If PE or PS is not given by the user, we assume PE.
Example of the usage of a quadrilateral element with 4 nodes, regular integration and plane strain.
examples/CantileverBeam/test4.inp
** cantilever beam, steel, l=10000mm, h=100mm, b=40mm, calculation in [N/mm²]
*material, name=linearelastic, id=linearelastic,
** elasticity module and poisson's ratio
210000.0, 0.15
*section, name=section1, thickness=40.0, material=linearelastic, type=plane
** assign section with 40mm thickness and linear elastic to all elements
all
** make job
*job, name=cps4job, domain=2d
** choose solver
*solver, solver=NIST, name=theSolver
*fieldOutput
create=perNode, nSet=gen_bottom, result=P, field=displacement, name=RF
** displacement output for all elements
create=perNode, elSet=all, field=displacement, result=U, name=displacement
*modelGenerator, generator=planeRectQuad, name=gen
** create model: length 10000mm, height 100mm with 300 elements in l direction and 4 in h direction
x0=0, l=10000
y0=0, h=100
elType=Quad4NPE
** quadrilateral element with 4 nodes, normal integration and plane strain
** set provider
elProvider=displacementelement
nX=300
nY=4
*output, type=monitor
fieldOutput=RF, f(x)='mean(x[:,1])'
** create output data to open with paraview
*output, type=ensight, name=esExport
configuration, overwrite=yes
create=perNode, fieldOutput=displacement
** set steps with size between 1 and 1e-8, 1000 iterations and max 25 iterations
*step, solver=theSolver, maxInc=1e-0, minInc=1e-8, maxNumInc=1000, maxIter=25, stepLength=100
** fix one side (fixed support)
dirichlet, name=left, nSet=gen_left, field=displacement, 2=0, 1=0
** set load -1000N in one point down
nodeforces, name=cloadTop, nSet=gen_rightTop, components='0,-1000', field=displacement
Example of the usage of a quadrilateral element with 8 nodes, regular integration and plane strain for the same example.
examples/CantileverBeam/test8.inp
** cantilever beam, steel, l=10000mm, h=100mm, b=40mm, calculation in [N/mm²]
*material, name=linearelastic, id=linearelastic,
** elasticity module and poisson's ratio
210000.0, 0.15
*section, name=section1, thickness=40.0, material=linearelastic, type=plane
** assign section with 40mm thickness and linear elastic to all elements
all
** make job
*job, name=cps4job, domain=2d
** choose solver
*solver, solver=NIST, name=theSolver
*fieldOutput
create=perNode, nSet=gen_bottom, result=P, field=displacement, name=RF
** displacement output for all elements
create=perNode, elSet=all, field=displacement, result=U, name=displacement
*modelGenerator, generator=planeRectQuad, name=gen
** create model: length 10000mm, height 100mm with 300 elements in l direction and 4 in h direction
x0=0, l=10000
y0=0, h=100
elType=Quad8NPE
** quadrilateral element with 4 nodes, normal integration and plane strain
** set provider
elProvider=displacementelement
nX=300
nY=4
*output, type=monitor
fieldOutput=RF, f(x)='mean(x[:,1])'
** create output data to open with paraview
*output, type=ensight, name=esExport
configuration, overwrite=yes
create=perNode, fieldOutput=displacement
** set steps with size between 1 and 1e-8, 1000 iterations and max 25 iterations with step length 100
*step, solver=theSolver, maxInc=1e-0, minInc=1e-8, maxNumInc=1000, maxIter=25, stepLength=100
** fix one side (fixed support)
dirichlet, name=left, nSet=gen_left, field=displacement, 2=0, 1=0
** set load -1000N in one point down
nodeforces, name=cloadTop, nSet=gen_rightTop, components='0,-1000', field=displacement
Example of the usage of a 3D hexahedron element with 8 nodes.
examples/CantileverBeam/test3D.inp
*material, name=LinearElastic, id=linearelastic
**Isotropic material
**E nu -> elasticity module and poisson's ratio
1.8e4, 0.22
** make job
*job, name=c3d8job, domain=3d
** choose solver
*solver, solver=NIST, name=theSolver
*modelGenerator, generator=boxGen, name=gen
** set number of elements for directions: 400, 4, 4
nX =400
nY =4
nZ =4
** create model: length 500mm, with section of 5mm x 5mm
lX =500
lY =5
lZ =5
** use provider: displacementelement
elProvider =displacementelement
** use hexahedron element with 8 nodes and normal integration
elType =Hexa8
*section, name=section1, material=linearelastic, type=solid
** assign linear elastic section to all elements
all
** displacement output for all elements
*fieldOutput
create=perNode, elSet=all, field=displacement, result=U, name=displacement
** create output data to open with paraview
*output, type=ensight, name=test
create=perNode, fieldOutput=displacement
configuration, overwrite=yes
*output, type=monitor
** set step size between 1 and 1e-9, 1000 iterations, 25 max iterations and step length 1
*step, solver=theSolver, maxInc=1e0, minInc=1e-9, maxNumInc=1000, maxIter=25, stepLength=1
options, category=NISTSolver, extrapolation=off
** fix one side (fixed support)
dirichlet, name=Left, nSet=gen_left, field=displacement, 1=0.0, 2=0.0, 3=.0
** set load in z to 0.012N for all the points on right back
nodeforces, name=cloadTop, nSet=gen_rightBack, components='0, 0, .012', field=displacement
Implementing your own elements
Relevant module: edelweissfe.elements.base.baseelement
Implementing your own finite elements can be done easily by subclassing from
the abstract base class BaseElement
.
- class edelweissfe.elements.base.baseelement.BaseElement(elementType, elNumber)[source]
Finite elements in EdelweissFE should be derived from this base class in order to follow the general interface.
EdelweissFE expects the layout of the internal and external load vectors, P, PExt, (and the stiffness) to be of the form
[ node 1 - dofs field 1, node 1 - dofs field 2, node 1 - ... , node 1 - dofs field n, node 2 - dofs field 1, ... , node N - dofs field n].
- Parameters
elementType (str) – A string identifying the requested element formulation.
elNumber (int) – A unique integer label used for all kinds of purposes.
- Attributes
dofIndicesPermutation
If the provides computes residual vectors and stiffness matrices not nodewise, but e.g., fieldwise, this permutation pattern is used to aggregate all entries in order to resemble the defined fields nodewise.
elNumber
The unique number of this element
ensightType
The shape of the element in Ensight Gold notation.
fields
The list of fields per node which this entity couples.
hasMaterial
Flag to check if a material was assigned to this element.
nDof
The total number of degrees of freedom this entity has.
nNodes
The number of nodes this entity couples.
nodes
The list of nodes this currently entity holds.
visualizationNodes
The nodes for visualization.
Methods
Accept the computed state (in nonlinear iteration schemes).
computeBodyForce
(P, K, load, U, time, dTime)Evaluate residual and stiffness for given time, field, and field increment due to a body force load.
computeDistributedLoad
(loadType, P, K, ...)Evaluate residual and stiffness for given time, field, and field increment due to a surface load.
computeYourself
(P, K, U, dU, time, dT)Evaluate the residual and stiffness for given time, field, and field increment due to a surface load.
Compute the underlying MarmotElement centroid coordinates.
Compute the underlying MarmotElement qp coordinates.
Compute the underlying MarmotElement qp coordinates.
getResultArray
(result, quadraturePoint[, ...])Get the array of a result, possibly as a persistent view which is continiously updated by the element.
Initalize the element to be ready for computing.
Rest to the last valid state.
setInitialCondition
(stateType, values)Assign initial conditions.
setMaterial
(materialName, materialProperties)Assign a material and material properties.
setNodes
(nodes)Assign the nodes to the element.
setProperties
(elementProperties)Assign a set of properties to the element.
- abstract computeBodyForce(P, K, load, U, time, dTime)[source]
Evaluate residual and stiffness for given time, field, and field increment due to a body force load.
- Parameters
P (ndarray) – The external load vector to be defined.
K (ndarray) – The stiffness matrix to be defined.
load (ndarray) – The magnitude (or vector) describing the load.
U (ndarray) – The current solution vector.
time (ndarray) – Array of step time and total time.
dTime (float) – The time increment.
- abstract computeDistributedLoad(loadType, P, K, faceID, load, U, time, dT)[source]
Evaluate residual and stiffness for given time, field, and field increment due to a surface load.
- Parameters
loadType (str) – The type of load.
P (ndarray) – The external load vector to be defined.
K (ndarray) – The stiffness matrix to be defined.
faceID (int) – The number of the elements face this load acts on.
load (ndarray) – The magnitude (or vector) describing the load.
U (ndarray) – The current solution vector.
time (ndarray) – Array of step time and total time.
dTime – The time increment.
dT (float) –
- abstract computeYourself(P, K, U, dU, time, dT)[source]
Evaluate the residual and stiffness for given time, field, and field increment due to a surface load.
- Parameters
P (ndarray) – The external load vector to be defined.
K (ndarray) – The stiffness matrix to be defined.
U (ndarray) – The current solution vector.
dU (ndarray) – The current solution vector increment.
time (ndarray) – Array of step time and total time.
dTime – The time increment.
dT (float) –
- abstract getCoordinatesAtCenter()[source]
Compute the underlying MarmotElement centroid coordinates.
- Returns
The element’s central coordinates.
- Return type
np.ndarray
- abstract getCoordinatesAtQuadraturePoints()[source]
Compute the underlying MarmotElement qp coordinates.
- Returns
The element’s qp coordinates.
- Return type
np.ndarray
- abstract getNumberOfQuadraturePoints()[source]
Compute the underlying MarmotElement qp coordinates.
- Returns
The element’s qp coordinates.
- Return type
np.ndarray
- abstract getResultArray(result, quadraturePoint, getPersistentView=True)[source]
Get the array of a result, possibly as a persistent view which is continiously updated by the element.
- Parameters
result (str) – The name of the result.
quadraturePoint (int) – The number of the quadrature point.
getPersistentView (bool) – If true, the returned array should be continiously updated by the element.
- Returns
The result.
- Return type
np.ndarray
- abstract setInitialCondition(stateType, values)[source]
Assign initial conditions.
- Parameters
stateType (str) – The type of initial state.
values (ndarray) – The numpy array describing the initial state.
- abstract setMaterial(materialName, materialProperties)[source]
Assign a material and material properties.
- Parameters
materialName (str) – The name of the requested material.
materialProperties (ndarray) – The numpy array containing the material properties.
- abstract setNodes(nodes)[source]
Assign the nodes to the element.
- Parameters
nodes (list[edelweissfe.points.node.Node]) – A list of nodes.
- abstract setProperties(elementProperties)[source]
Assign a set of properties to the element.
- Parameters
elementProperties (ndarray) – A numpy array containing the element properties.
- abstract property elNumber: int
The unique number of this element
- abstract property hasMaterial: str
Flag to check if a material was assigned to this element.