In this section, it will be shown how to write the evolution of given fields to disk, e.g. the contact force as function of time. This is done with what is called “result curves”.
These curves are saved for each time step while archiving files (bfac
) are usually saved at a lower sampling frequency. This leads to a very good temporal resolution for these curves.
The respective files are saved in the folder of the test case (workspace/test_name
by default). The information of each curve is saved twice, once as a binary file used internally (values.v
), and once as a text file (valeurs.ascii
), easy to load in Matlab, Excel, …
The class “ValuesManager
” manages the interfac, allowing :
ValueExtractor
(extracting results), Depending on whether one or several values are extracted, the file on disk will be a vector or a matrix. The number of values saved for each time step depends on the extraction command.
The ValuesManager
manages the extraction of data at each time step. It is a container of ValueExtractors
defined in the analysis metafor
:
valuesmanager = metafor.getValuesManager() valuesmanager.add(nbr, extractor, name = "") # manages extraction, result storage... valuesmanager.add(nbr, extractor, v2sOp, name = "") # same + transformation from vector to scalar valuesmanager.getDataVector(nbr) # returns a reference of the vector storing values # (for example to create and view curves) valuesmanager.setOnFile(False) # Keep curves in memory instead of writing # .v and .ascii files
where
nbr | number of the curve (unique) |
extractor | reference of a ValueExtractor |
v2sOp | Operator such as VectorToScalarOperator |
name | name of the curve (of the file) |
can be generated automatically by the extractor (not always meaningful) |
The FacValuesManager
manages the extraction of data at each time a “bfac” file is written. It is exactly the same object as the ValuesManager
, but called at different times of the time integration.
To access the FacValuesManager
, just get it from metafor. The rest of the syntax is the same.
facValuesManager = metafor.getFacValuesManager()
The StageValuesManager
manages the extraction of data at the end of each stage. It is exactly the same object as the ValuesManager
, but called at different times of the time integration.
To access the StageValuesManager
, just get it from metafor. The rest of the syntax is the same.
stageValuesManager = metafor.getStageValuesManager()
ValueExtractor
defines which values are to be archived.
It manages:
Extract values stored in the data base, so fields which are indeed stored in nodes (positions, displacements, temperature, forces, velocities, …).
valueExtractor = DbNodalValueExtractor(gObject, dbField) valueExtractor = DbNodalValueExtractor(gObject, dbField, sOp = None, maxV = -1)
where
gObject | Geometric object containing the nodes at which the field is extracted |
dbField | Field1D describing the scalar to extract |
sOp | SortingOperator : sort nodes according to a geometric criterion def=None will give the order in which these nodes are generated (not reliable) |
maxV | Number of values saved after sorting def=-1 ⇒ keep all values |
Example:
valuesmanager.add(8, DbNodalValueExtractor( curveset(1), Field1D(TY,GF2)), SumOperator(), 'totalForce_X')
Extract a result curve numbered 8 defined by the sum of all internal forces along line 1. This curve is named totalForce_X
.
Extract Internal Fields Values at nodes. As these values only exists at integration point level of the elements, the extracted values are the result of the extrapolation of Integration points values to the node (a contribution of each element using the node is taken into account).
valueExtractor = IFNodalValueExtractor(gObject, InternalField, sOp=None, maxV=-1)
where
gObject | Geometric object containing the nodes at which the field is extracted |
InternalField | internal field of the element extrapolated in the node |
sOp | SortingOperator : sort nodes according to a geometric criterion def=None will give the order in which these nodes are generated (not reliable) |
maxV | Number of values saved after sorting def=-1 ⇒ keep all values |
A set of extractor can be defined on points that do not correspond to node position.
The value is computed from the interpolation of elementary values from the first element that belong the point. An outer tolerance is defined (OutsideTol) to allow detection of point just outside the mesh.
The point can be Lagrangian (fixed to material : default) or Eulerian (geometrically defined fixed or not) where element and reduced coordinate are computed at each extraction)
Extract values stored in the data base, On geometrically defined point (positions, displacements, temperature, forces, velocities, …).
valueExtractor = DbGeoPointValueExtractor (pointList, meshedObject, field) valueExtractor = DbGeoPointValueExtractor (pointList, meshedObjectList, field) valueExtractor.setEulerian(True) valueExtractor.setOutsideTol(tol)
where
listPt | List of points where the field is extracted |
meshObject | Meshed Object from where the point will try to find the element from which the value is interpolated |
meshObjectList | List of Mesh Meshed Object from where the point will try to find the element from which the value is interpolated |
dbField | Field1D describing the scalar to extract |
tol | Outside tolerance in terms of reduced coordinate (default = 0.2) |
Example:
pt102 = pointset.define(102,0.15*Lx,0.75*Ly,0) # geo point of extraction (must be in pointset) dyPt102 = DbGeoPointValueExtractor([pt102,], [volset(idx+1), volset(idx+2)], Field1D(TY, RE)) dyPt102.setEulerian(True) # Eulerian point (not moving with elements) dyPt102.setOutsideTol(0.5) # half an element outside detection tolerance
Extract Internal Fields Values on geometrically defined point(s). Value are interpolated from Integration Points Values.
valueExtractor = IFGeoPointValueExtractor (pointList, meshedObject, InternalField) valueExtractor = IFGeoPointValueExtractor (pointList, meshedObjectList, InternalField) valueExtractor.setEulerian(True) valueExtractor.setOutsideTol(tol)
where
listPt | List of points where the field is extracted |
meshObject | Meshed Object from where the point will try to find the element from which the value is interpolated |
meshObjectList | List of Mesh Meshed Object from where the point will try to find the element from which the value is interpolated |
InternalField | internal field of the element extrapolated in the node |
tol | Outside tolerance in terms of reduced coordinate (default = 0.2) |
Example:
pt102 = pointset.define(102,0.15*Lx,0.75*Ly,0) # geo point of extraction (must be in pointset) eplPt102 = IFGeoPointValueExtractor([pt102,], volset(101), IF_EPL) eplPt102.setEulerian(True) # Eulerian point (not moving with elements) eplPt102.setOutsideTol(0.5) # half an element outside detection tolerance
It is possible to modelize a strainGauge in Metafor using GeoPointValueExtractors. Selecting the 2 geometrical points (anywhere in the mesh) defining the Strain Gauge, the distance between the tho points can be computed during the simulation allowing to define a StrainGaufeValueExtractor. The strain measure define the kind of extractor :
valueExtractor = BiotStrainGaugeValueExtractor(listPt, meshedObject) valueExtractor = BiotStrainGaugeValueExtractor(listPt, meshedObjectList) valueExtractor = GLStrainGaugeValueExtractor(listPt, meshedObject) valueExtractor = GLStrainGaugeValueExtractor(listPt, meshedObjectList) valueExtractor = NatStrainGaugeValueExtractor(listPt, meshedObject) valueExtractor = NatStrainGaugeValueExtractor(listPt, meshedObjectList)
where
listPt | List of the 2 points defining the strain Gauge |
meshObject | Meshed Object from where the point will try to find the element from which the value is interpolated |
meshObjectList | List of Mesh Meshed Object from where the point will try to find the element from which the value is interpolated |
and if $\lambda = \frac{l}{l_0}$ the ratio of current to initial length of the Strain Gauge
BiotStrainGaugeValueExtractor | $\epsilon = \lambda-1 = \frac{l}{l_0}-1$ (Engineer Strain) |
GLStrainGaugeValueExtractor | $\epsilon = 0.5*(\lambda^2-1) $ |
NatStrainGaugeValueExtractor | $\epsilon = ln(\lambda)$ |
Extract thermodynamic fields :
valueExtractor = TdFieldValueExtractor(meta, gObject, tdFieldID) valueExtractor = TdFieldValueExtractor(meta, nodeSet, tdFieldID)
meta | Reference of the analysis Metafor |
gObject | Geometric object containing the nodes at which the thermodynamic field is extracted |
nodeset | Reference of Metafor NodeSet (calculation over the entire domain) |
tdFieldID | TdFieldID : thermodynamic field to extract |
THERMODYN_TRAV_FINT : Work of Internal Forces | |
THERMODYN_EN_CIN : Kinetic Energy | |
THERMODYN_EN_DIS : Dissipated Energy | |
THERMODYN_TRAV_FEXT : Work of External Forces | |
THERMODYN_POT_INT : Internal Potential | |
THERMODYN_LIN_MOM_X : X Linear Momentum | |
THERMODYN_LIN_MOM_Y : Y Linear Momentum | |
THERMODYN_LIN_MOM_Z : Z Linear Momentum | |
THERMODYN_ANG_MOM_X : X Angular Momentum | |
THERMODYN_ANG_MOM_Y : Y Angular Momentum | |
THERMODYN_ANG_MOM_Z : Z Angular Momentum |
momentVE = MomentValueExtractor(target, center, momentNature, targetVariant,sOp=None, maxV=-1)
Extract a component of the moment force vector (internal, external, inertial) with respect to a center of an object containing these nodes (meshed object, group, …).
target | Geometric object containing the nodes at which the field is extracted |
center | Point (GeoObject0D) around which the moment is computed |
momentNature | Component of the extracted moment vector (TX, TY or TZ) |
targetVariant | variant defining the type of force whose moment is computed (GF1, GF2 or GF3) |
sOp | SortingOperator : sort nodes according to a geometric criterion def=None will give the order in which these nodes are generated (not reliable) |
maxV | Number of values saved after sorting def=-1 ⇒ keep all values |
MomentValueExtractor
has all its meaning if used with VectorToScalarOperatormoment2AxeVE = Moment2AxeValueExtractor(target, axis, targetVariant,sOp=None, maxV=-1)
Extract the moment force (internal, external, inertial) with respect to an axis containing these nodes (meshed object, group, …). It computes the torque of forces around the axis.
target | Geometric object containing the nodes at which the field is extracted |
axis | Reference to a Axes object |
targetVariant | variant defining the type of force whose moment is computed (GF1, GF2 or GF3) |
sOp | SortingOperator : sort nodes according to a geometric criterion def=None will give the order in which these nodes are generated (not reliable) |
maxV | Number of values saved after sorting def=-1 ⇒ keep all values |
Moment2AxeValueExtractor
has all its meaning if used with VectorToScalarOperatorExtract values of internal fields at the integration points.
valueExtractor = IFGaussPointValueExtractor (topoCell, ifield)
where
topoCell | Topological object supporting the element |
2D : Topological face | |
3D : Topological volume | |
ifield | internal field of the element extrapolated to the node |
if ifield = TX, TY or TZ the extractor computes the current position of the Integration Point |
works, but it is difficult to select the topoCell
Extract values of internal fields at the integration points for all the elements of aninteraction.
valueExtractor = IFGPInteractionValueExtractor(inter, ifield)
where
inter | Reference to a volumic interaction (FieldApplicator) |
ifield | internal field of the element extrapolated to the node |
Extract a value of internal fields averaged over an element.
valueExtractor = IFElementValueExtractor (topoCell, ifield)
where
topoCell | Topological object supporting the element |
2D : Topological face | |
3D : Topological volume | |
ifield | internal field of the element averaged over the element |
works, but selection of the topoCell may be difficult
Extract values of internal fields averaged on each element of an interaction (1 mean value by element).
valueExtractor = IFElementsValueExtractor (inter, ifield)
where
inter | Reference to a volumic Interaction (FieldApplicator) |
ifield | internal field of the element averaged over the element |
Extract values associated to the computation (analysis values)
valueExtractor = MiscValueExtractor(meta, type)
meta | Reference of the analysis Metafor |
type | defines the value to extract: |
EXT_T : time |
|
EXT_DT : time step |
|
EXT_NT : time step number |
|
EXT_ITE : number of mechanical iterations over a time step |
|
EXT_ITE_TOT : total number of mechanical iterations (included failed steps !) |
|
EXT_ITE_TH : number of thermal iterations over a time step |
|
EXT_ITE_TH_TOT : total number of thermal iterations (included failed steps !) |
|
EXT_ITE_AR : number of mechanical iterations with update of tangent stiffness matrix over a time step |
|
EXT_ITE_SR : number of mechanical iterations without update of tangent stiffness matrix over a time step |
|
EXT_ITE_AR_TOT : total number of mechanical iterations with update of tangent stiffness matrix (included failed steps !) |
|
EXT_ITE_SR_TOT : total number of mechanical iterations without update of tangent stiffness matrix (included failed steps !) |
|
EXT_ITE_TH_AR : number of thermal iterations with update of tangent stiffness matrix over a time step |
|
EXT_ITE_TH_SR : number of thermal iterations without update of tangent stiffness matrix over a time step |
|
EXT_ITE_TH_AR_TOT : total number of thermal iterations with update of tangent stiffness matrix (included failed steps !) |
|
EXT_ITE_TH_SR_TOT : total number of thermal iterations without update of tangent stiffness matrix (included failed steps !) |
|
EXT_ITE_LS : number of mechanical line-search iterations over a time step |
|
EXT_ITE_LS_TOT : total number of mechanical line-search iterations (included failed steps !) |
|
EXT_ITE_TH_LS : number of thermal line-search iterations over a time step |
|
EXT_ITE_TH_LS_TOT : total number of thermal line-search iterations (included failed steps !) |
|
EXT_ITE_ALM : number of augmentations over a time step |
|
EXT_ITE_ALM_TOT : total number of augmentations (included failed steps !) |
|
EXT_RCOND : reciprocal condition number of the mechanical tangent stiffness matrix |
|
EXT_VM…. : Extractors of virtual memory (see table below) |
|
EXT_USER_CPU : User CPU |
|
EXT_REAL_CPU : Real CPU |
|
EXT_KERNEL_CPU : Kernel CPU |
Extractors of virtual memory in Linux and Windows
ExtractType | Linux description | Windows description |
EXT_VMPEAK | VmPeak : Peak Virtual Memory usage | PeakWorkingSetSize : The Peak Working Set size |
EXT_VMSIZE | VmSize : Current Virtual Memory usage | WorkingSetSizeCurrent : The current Working Set Size |
EXT_VMLCK | VMLCK : Current mlocked Memory | - |
EXT_VMHWM | VMHWM : Peak resident set size | - |
EXT_VMRSS | VmRSS : Virtual Memory ResidentSetSize | - |
EXT_VMDATA | VmData : size of “data” segment | - |
EXT_VMSTK | VmSTK : size of stack | - |
EXT_VMEXE | VmExe : size of “test” segment | - |
EXT_VMLIB | VmLib : Shared library usage | - |
EXT_VMPTE | VmPTE : pagetable entries size | PeakPagefileUsage : The peak value in bytes of the Commit Charge during te lifetime of this process |
EXT_VMSWAP | VmSwap : swap space used | PagefileUsage : The Commit Charge value in bytes for this process |
Example:
valuesmanager.add(1, MiscValueExtractor(metafor, EXT_T),'time')
Extract the time value for each time step, in the vector 1 named time
.
Extract a component of a vector computed by the elements of an interaction.
valueExtractor = InteractionValueExtractor (interaction, natureId, vectorId=GEN_EXT_FORC)
where
interaction | Reference of an interaction |
natureId | Fields at Nodes (TX , TY , TZ or TO ) (only scalar) |
vectorId | VectorID : ID defining the vector computed by the element |
GEN_INTER_FORC : Generalized internal forces |
|
GEN_EXT_FORC : Generalized external forces |
|
GEN_INERT_FORC : Generalized inertial forces (consistent mass matrix) |
|
GEN_DIAG_MASS_FORC : Generalized inertial forces (diagonalized mass matrix) |
|
GEN_DIS_FORC : Generalized dissipation forces |
Extract the time evolution of the position of the center of gravity and mass of an interaction. At each archiving step the extractor fills a four columns text file, where columns 1-3 contain x, y, and z coordinates of the center of gravity respectively, and column 4 contains the value of the mass.
valueExtractor = InteractionGravityCenterAndMassValueExtractor (interaction)
where
interaction
–> Reference of an interaction.
Series of extractors associated to contact interaction. They give values for contact elements, which is to say in each slave nodes.
The matricial output can be sorted and reduced using SortingOperators and/or can be made scalar (at the ValuesManager level) using a VectorToScalarOperator.
Extracted values have a filter allowing to archive the extractor name and a subset of slave nodes of the contact interaction.
valueExtractor = ContactStatusValueExtractor(contInt, sOp=None, maxV=-1) valueExtractor.setGeoFilter(gObject) gObject = meshed geometric object, included in the object on which the contact interaction is defined.
By default, all contact elements associated to the interaction are considered (for defo-defo contact in two stages, both master and slave nodes are considered).
Operations are done in the following order:
- Filter operation - Sorting operation - “resize” operation, so reducing the number of stored values after filtering and sorting.
Extract the status of normal contact in each slave node:
valueExtractor = ContactStatusValueExtractor(contInt, sOp=None, maxV=-1)
Inputs:
contInt | contact interaction |
sOp | SortingOperator |
defaut=None |
|
maxV | Number of values saved after sorting def=-1 ⇒ keep all values |
Outputs:
Status | Extracted value |
---|---|
Inactive contact | 0.0 |
Active contact | 1.0 |
The number of nodes in contact is easily deduced:
valuesmanager.add(no, ContactStatusValueExtractor(ci), SumOperator(), 'NodesInContact')
Extract the status of tangential contact in each slave node:
valueExtractor = SlidingStatusValueExtractor(contInt, sOp=None, maxV=-1)
Inputs:
contInt | contact interaction |
sOp | SortingOperator |
defaut=None |
|
maxV | Number of values saved after sorting def=-1 ⇒ keep all values |
Outputs :
Status | Extracted value |
---|---|
Inactive contact | 0.0 |
Sticking contact | 0.0 |
Sliding contact | 1.0 |
Extract the normal gap in each slave node:
valueExtractor = NormalGapValueExtractor(contInt, sOp=None, maxV=-1 , _nonContactGapsSetToZero=True)
Inputs :
contInt | contact interaction |
sOp | SortingOperator |
defaut=None |
|
maxV | Number of values saved after sorting def=-1 ⇒ keep all values |
_nonContactGapsSetToZero | Set the normal gap of inactive slave node to zero (def = True) |
Outputs :
Status | Extracted value |
---|---|
Inactive contact | 0.0 |
Active contact | normal gap (signed value, positive if overlap???) |
The maximal normal gap is easily deduced:
valuesmanager.add(no, NormalGapValueExtractor(ci), AbsMaxOperator(), 'GapMax')
Extract the norm of tangential gap in each slave node:
valueExtractor = TangentGapValueExtractor(contInt, sOp=None, maxV=-1)
Inputs:
contInt | contact interaction |
sOp | SortingOperator |
defaut=None |
|
maxV | Number of values saved after sorting def=-1 ⇒ keep all values |
The value corresponds to the total gap (sticking gap + sliding gap) with respect to the sticking point of the previous balanced configuration when the slave node is sliding/or sticking. But, when the slave node is sliding, this gap is considered as zero in this extractor and thus only the slave node in sticking contact receives a non-zero value.
It is possible to archive the sliding and the sticking gaps :
valueExtractor.setNonStickingGapsSetToZero(False)
Extract the component 1 of tangential gap in each slave node:
valueExtractor = TangentGap1ValueExtractor(contInt, sOp=None, maxV=-1)
Inputs:
contInt | contact interaction |
sOp | SortingOperator |
defaut=None |
|
maxV | Number of values saved after sorting def=-1 ⇒ keep all values |
By default, the components of tangential gap are related to the local tangent system defined by the sliding direction (The tangent 1 is aligned with this direction and the tangent 2 is obtained with the right hand rule.). It is possible to use another local tangent system related to the tool tangents (Pay attention when using with a defo-defo contact interaction !) :
valueExtractor.setUseToolLocalSystemAxes(True)
The value corresponds to the total gap (sticking gap + sliding gap) with respect to the sticking point of the previous balanced configuration when the slave node is sliding/or sticking. But, when the slave node is sliding, this gap is considered as zero in this extractor and thus only the slave node in sticking contact receives a non-zero value.
It is possible to archive the sliding and the sticking gaps :
valueExtractor.setNonStickingGapsSetToZero(False)
Extract the component 2 of tangential gap in each slave node:
valueExtractor = TangentGap2ValueExtractor(contInt, sOp=None, maxV=-1)
Inputs:
contInt | contact interaction |
sOp | SortingOperator |
defaut=None |
|
maxV | Number of values saved after sorting def=-1 ⇒ keep all values |
By default, the components of tangential gap are related to the local tangent system defined by the sliding direction (The tangent 1 is aligned with this direction and the tangent 2 is obtained with the right hand rule.). It is possible to use another local tangent system related to the tool tangents (Pay attention when using with a defo-defo contact interaction !) :
valueExtractor.setUseToolLocalSystemAxes(True)
The value corresponds to the total gap (sticking gap + sliding gap) with respect to the sticking point of the previous balanced configuration when the slave node is sliding/or sticking. But, when the slave node is sliding, this gap is considered as zero in this extractor and thus only the slave node in sticking contact receives a non-zero value.
It is possible to archive the sliding and the sticking gaps :
valueExtractor.setNonStickingGapsSetToZero(False)
Extract the nodal normal force in each slave node:
valueExtractor = NormalForceValueExtractor(contInt, sOp=None, maxV=-1)
Inputs:
contInt | contact interaction |
sOp | SortingOperator |
defaut=None |
|
maxV | Number of values saved after sorting def=-1 ⇒ keep all values |
Extract the contact pressure (negative or positive, if bodies are pushed against each other ???) in each slave node:
valueExtractor = PressureContactValueExtractor(contInt, sOp=None, maxV=-1)
Inputs:
contInt | contact interaction |
sOp | SortingOperator |
defaut=None |
|
maxV | Number of values saved after sorting def=-1 ⇒ keep all values |
Extract the normal augmented lagrangian in each slave node:
valueExtractor = NormalAugmentedLagrangianValueExtractor(contInt, sOp=None, maxV=-1)
Inputs:
contInt | contact interaction |
sOp | SortingOperator |
defaut=None |
|
maxV | Number of values saved after sorting def=-1 ⇒ keep all values |
Extract the norm of tangential force in each slave node:
valueExtractor = TangentForceValueExtractor(contInt, sOp=None, maxV=-1)
Inputs:
contInt | contact interaction |
sOp | SortingOperator |
defaut=None |
|
maxV | Number of values saved after sorting def=-1 ⇒ keep all values |
Extract the component 1 of tangential force in each slave node:
valueExtractor = Tangent1ForceValueExtractor(contInt, sOp=None, maxV=-1)
Inputs:
contInt | contact interaction |
sOp | SortingOperator |
defaut=None |
|
maxV | Number of values saved after sorting def=-1 ⇒ keep all values |
By default, the components of tangential force are related to the local tangent system defined by the sliding direction (The tangent 1 is aligned with this direction and the tangent 2 is obtained with the right hand rule.). It is possible to use another local tangent system related to the tool tangents (Pay attention when using with a defo-defo contact interaction !) :
valueExtractor.setUseToolLocalSystemAxes(True)
Extract the component 2 of tangential force in each slave node:
valueExtractor = Tangent2ForceValueExtractor(contInt, sOp=None, maxV=-1)
Inputs:
contInt | contact interaction |
sOp | SortingOperator |
defaut=None |
|
maxV | Number of values saved after sorting def=-1 ⇒ keep all values |
By default, the components of tangential force are related to the local tangent system defined by the sliding direction (The tangent 1 is aligned with this direction and the tangent 2 is obtained with the right hand rule.). It is possible to use another local tangent system related to the tool tangents (Pay attention when using with a defo-defo contact interaction !) :
valueExtractor.setUseToolLocalSystemAxes(True)
Extract the norm of shear stress in each slave node:
valueExtractor = ShearContactValueExtractor(contInt, sOp=None, maxV=-1)
Inputs:
contInt | contact interaction |
sOp | SortingOperator |
defaut=None |
|
maxV | Number of values saved after sorting def=-1 ⇒ keep all values |
Extract the component 1 of shear stress in each slave node:
valueExtractor = Shear1ContactValueExtractor(contInt, sOp=None, maxV=-1)
Inputs:
contInt | contact interaction |
sOp | SortingOperator |
defaut=None |
|
maxV | Number of values saved after sorting def=-1 ⇒ keep all values |
By default, the components of shear stress are related to the local tangent system defined by the sliding direction (The tangent 1 is aligned with this direction and the tangent 2 is obtained with the right hand rule.). It is possible to use another local tangent system related to the tool tangents (Pay attention when using with a defo-defo contact interaction !) :
valueExtractor.setUseToolLocalSystemAxes(True)
Extract the component 2 of shear stress in each slave node:
valueExtractor = Shear2ContactValueExtractor(contInt, sOp=None, maxV=-1)
Inputs:
contInt | contact interaction |
sOp | SortingOperator |
defaut=None |
|
maxV | Number of values saved after sorting def=-1 ⇒ keep all values |
By default, the components of shear stress are related to the local tangent system defined by the sliding direction (The tangent 1 is aligned with this direction and the tangent 2 is obtained with the right hand rule.). It is possible to use another local tangent system related to the tool tangents (Pay attention when using with a defo-defo contact interaction !) :
valueExtractor.setUseToolLocalSystemAxes(True)
Extract the norm of tangential augmented lagrangian in each slave node:
valueExtractor = TangentAugmentedLagrangianValueExtractor(contInt, sOp=None, maxV=-1)
Inputs:
contInt | contact interaction |
sOp | SortingOperator |
defaut=None |
|
maxV | Number of values saved after sorting def=-1 ⇒ keep all values |
Extract the component 1 of tangential augmented lagrangien in each slave node:
valueExtractor = TangentAugmentedLagrangian1ValueExtractor(contInt, sOp=None, maxV=-1)
Inputs:
contInt | contact interaction |
sOp | SortingOperator |
defaut=None |
|
maxV | Number of values saved after sorting def=-1 ⇒ keep all values |
By default, the components of tangential augmented lagrangian are related to the local tangent system defined by the sliding direction (The tangent 1 is aligned with this direction and the tangent 2 is obtained with the right hand rule.). It is possible to use another local tangent system related to the tool tangents (Pay attention when using with a defo-defo contact interaction !) :
valueExtractor.setUseToolLocalSystemAxes(True)
Extract the component 2 of tangential augmented lagrangien in each slave node:
valueExtractor = TangentAugmentedLagrangian2ValueExtractor(contInt, sOp=None, maxV=-1)
Inputs:
contInt | contact interaction |
sOp | SortingOperator |
defaut=None |
|
maxV | Number of values saved after sorting def=-1 ⇒ keep all values |
By default, the components of tangential augmented lagrangian are related to the local tangent system defined by the sliding direction (The tangent 1 is aligned with this direction and the tangent 2 is obtained with the right hand rule.). It is possible to use another local tangent system related to the tool tangents (Pay attention when using with a defo-defo contact interaction !) :
valueExtractor.setUseToolLocalSystemAxes(True)
Extract the slave node contact area in each slave node:
areaInContactValues=AreaInContactValueExtractor(contInt, sOp=None, maxV=-1)
Inputs:
contInt | contact interaction |
sOp | SortingOperator |
defaut=None |
|
maxV | Number of values saved after sorting def=-1 ⇒ keep all values |
Then, the area potentially in contact can be deduced:
valuesmanager.add(no, areaInContactValues, SumOperator(), 'AreaInContact')
Nodal contact area of slave nodes in contact (and master nodes in two stages defo-defo):
activeAreaInContactValues=AreaInContactValueExtractor(contInt, sOp=None, maxV=-1) activeAreaInContactValues.setOnlyInContactStatus()
Then, the active contact area (area which is truly in contact)
valuesmanager.add(no, activeAreaInContactValues, SumOperator(), 'ActiveAreaInContact')
Extract contact resultant forces along TX, TY or TZ of a contact interaction
contactForceValueExtractor=ContactForceValueExtractor(contInt, natureId, contactForceType=TOTAL_FORCE, sOp=None, maxV=-1)
Inputs:
contInt | contact interaction |
natureId | Fields at Nodes (TX , TY , TZ ) (only scalar) |
contactForceType | Type of contact force |
NORMAL_FORCE : Contact normal force |
|
TANGENTIAL_FORCE : Contact tangential force |
|
TOTAL_FORCE=default : Contact normal force + Contact tangential force |
|
sOp | SortingOperator |
default=None |
|
maxV | Number of values saved after sorting def=-1 ⇒ keep all values |
The contact force on a slave entity of a contact interaction is only taken into consideration in the extractor
Extract values associated to transfert during remeshing of an interaction (ALE or complete remeshing).
Extract the value of the field integral before or after transfer, the value of the difference between both or the value of the transfer error.
valueExtractor = TransferValueExtractor(transferRegion, field, type)
Inputs:
transferRegion | “TransferRegion” created on the studied interaction |
field | ScalarNatureID or Field1D |
type | BEFORE, AFTER, DIFF or TRANSFER_ERROR |
Value of the integral before transfer → type = BEFORE
Value of the integral after transfer → type = AFTER
Value of the difference between the integral before and after transfer → type = DIFF
value of the transfer error → type = TRANSFER_ERROR
The transfer error is defined as:
$$\frac{\int f^{old}-f^{new}}{\int f^{new} }$$
Number of active elements:
valueExtractor = NumberOfActiveElementsExtractor(interaction)
Number of inactive elements:
valueExtractor = NumberOfInactiveElementsExtractor(interaction)
where
interaction | reference to the interaction |
Example to generate a file brokenelements.ascii which contains the number of elements broken for interaction 1:
valuesmanager.add(1, NumberOfInactiveElementsExtractor(interactionset(1)),'brokenelements')
When a curve is defined by the ValuesManager
, a ''SortingOperator'' and the number of points to keep after sorting can be added. For example (see also apps/qs/cont2.py
):
valuesManager = metafor.getValuesManager() valuesManager .add(8, DbNodalValueExtractor(curveset(1), Field1D(TY,GF2), SortByDist0(0.5,0.0,0.0), 2), 'totalForce_X') valuesManager .add(9, IFNodalValueExtractor(curveset(1), IF_J2, SortByKsi0(curset(3), 8),'j2OnBaseLine')
The first curve displays the force field as a function of time, on the two nodes closest to (0.5, 0). This can be viewed in Matlab using :
load totalForce_X; mesh(totalForce_X);
The second curve displays the equivalent stress field as a function of time, on the 8 first nodes of line 1, sorted by their curvilinear abscissa when projected on curve 3.
Available sorting operators :
SortByX0 | sort along increasing X0 of the nodes |
SortByY0 | sort along increasing Y0 of the nodes |
SortByZ0 | sort along increasing Z0 of the nodes |
SortByDX0 | sort along decreasing X0 of the nodes |
SortByDY0 | sort along decreasing Y0 of the nodes |
SortByDZ0 | sort along decreasing Z0 of the nodes |
SortByKsi0 | sort along curvilinear abscissa of the projection on a given curve |
SortByDist0 | sort along the initial distance of the nodes from a given point |
SortByPos0 | sort along the position and a Lock - TX , TY or TZ |
SortByNo(num) | sort along the user number of the points (or the distance of the number of the points to a given number num ) |
In the previous example, note the parameter None
in the command defining values to extract. This parameter can be replaced by a mathematical operator which transforms the list of extracted value at a given time into a scalar.
Several operators exist :
None
: no operator (all component,s by default)MaxOperator
: maximal valueMinOperator
: minimal valueAbsMaxOperator
: compute the value whose absolute value is maximal (but keeps its sign)AbsMinOperator
: compute the value whose absolute value is minimal (but keeps its sign)MaxAbsOperator
: compute the maximal value of the vector containing absolute values (always positive)MinAbsOperator
: compute the minimal value of the vector containing absolute values (always positive)MaxOfNonZeroOperator
: maximal value excluding zero valuesMinOfNonZeroOperator
: minimal value excluding zero valuesAbsMaxNonZeroOperator
: compute the value whose absolute value is maximal (but keeps its sign), excluding zero valuesAbsMinNonZeroOperator
: compute the value whose absolute value is minimal (but keeps its sign), excluding zero valuesMaxAbsNonZeroOperator
: compute the maximal value of the vector containing absolute values (always positive), excluding zero valuesMinAbsNonZeroOperator
: compute the minimal value of the vector containing absolute values (always positive), excluding zero valuesMeanOperator
: compute the average valueSumOperator
: compute the sum (used toe extract the resultant force exercised on a given geometric entity)Use:
Extraction of maximal value of plastic deformation on a volume:
curves = metafor.getValuesManager() ext = IFNodalValueExtractor(curset(1), IF_EPL) curves.add(1, ext, MeanOperator(), 'meanEplOnCurve1')