|
OOFEM 3.0
|
#include <tr1bubblestokes.h>
Public Member Functions | |
| Tr1BubbleStokes (int n, Domain *d) | |
| double | computeVolumeAround (GaussPoint *gp) override |
| void | computeGaussPoints () override |
| void | giveCharacteristicVector (FloatArray &answer, CharType type, ValueModeType mode, TimeStep *tStep) override |
| void | giveCharacteristicMatrix (FloatMatrix &answer, CharType type, TimeStep *tStep) override |
| void | computeInternalForcesVector (FloatArray &answer, TimeStep *tStep) |
| void | computeStiffnessMatrix (FloatMatrix &answer, MatResponseMode mode, TimeStep *tStep) |
| void | computeExternalForcesVector (FloatArray &answer, TimeStep *tStep) |
| void | computeLoadVector (FloatArray &answer, BodyLoad *load, CharType type, ValueModeType mode, TimeStep *tStep) override |
| void | computeBoundarySurfaceLoadVector (FloatArray &answer, BoundaryLoad *load, int boundary, CharType type, ValueModeType mode, TimeStep *tStep, bool global=true) override |
| const char * | giveClassName () const override |
| const char * | giveInputRecordName () const override |
| MaterialMode | giveMaterialMode () override |
| Element_Geometry_Type | giveGeometryType () const override |
| int | computeNumberOfDofs () override |
| int | giveNumberOfInternalDofManagers () const override |
| DofManager * | giveInternalDofManager (int i) const override |
| void | giveInternalDofManDofIDMask (int i, IntArray &answer) const override |
| FEInterpolation * | giveInterpolation () const override |
| FEInterpolation * | giveInterpolation (DofIDItem id) const override |
| void | giveDofManDofIDMask (int inode, IntArray &answer) const override |
| void | updateYourself (TimeStep *tStep) override |
| Interface * | giveInterface (InterfaceType it) override |
| void | computeField (ValueModeType u, TimeStep *tStep, const FloatArray &coords, FloatArray &answer) override |
| Public Member Functions inherited from oofem::FMElement | |
| FMElement (int n, Domain *aDomain) | |
| virtual void | updateStabilizationCoeffs (TimeStep *tStep) |
| void | computeVectorOfVelocities (ValueModeType mode, TimeStep *tStep, FloatArray &velocities) |
| void | computeVectorOfPressures (ValueModeType mode, TimeStep *tStep, FloatArray &pressures) |
| FloatArray | computeVectorOfVelocities (ValueModeType mode, TimeStep *tStep) |
| FloatArray | computeVectorOfPressures (ValueModeType mode, TimeStep *tStep) |
| Public Member Functions inherited from oofem::Element | |
| Element (int n, Domain *aDomain) | |
| Element (const Element &src)=delete | |
| Element & | operator= (const Element &src)=delete |
| virtual | ~Element () |
| Virtual destructor. | |
| void | giveLocationArray (IntArray &locationArray, const UnknownNumberingScheme &s, IntArray *dofIds=NULL) const |
| void | giveLocationArray (IntArray &locationArray, const IntArray &dofIDMask, const UnknownNumberingScheme &s, IntArray *dofIds=NULL) const |
| virtual void | giveBoundaryLocationArray (IntArray &locationArray, const IntArray &bNodes, const UnknownNumberingScheme &s, IntArray *dofIds=NULL) |
| virtual void | giveBoundaryLocationArray (IntArray &locationArray, const IntArray &bNodes, const IntArray &dofIDMask, const UnknownNumberingScheme &s, IntArray *dofIds=NULL) |
| virtual int | giveNumberOfDofs () |
| virtual void | setInternalDofManager (int num, std::unique_ptr< DofManager > dm) |
| virtual double | giveCharacteristicValue (CharType type, TimeStep *tStep) |
| virtual void | computeTangentFromSurfaceLoad (FloatMatrix &answer, BoundaryLoad *load, int boundary, MatResponseMode rmode, TimeStep *tStep) |
| virtual void | computeTangentFromEdgeLoad (FloatMatrix &answer, BoundaryLoad *load, int boundary, MatResponseMode rmode, TimeStep *tStep) |
| virtual void | computeBoundaryEdgeLoadVector (FloatArray &answer, BoundaryLoad *load, int edge, CharType type, ValueModeType mode, TimeStep *tStep, bool global=true) |
| const IntArray & | giveBodyLoadList () const |
| const IntArray & | giveBoundaryLoadList () const |
| void | computeVectorOf (ValueModeType u, TimeStep *tStep, FloatArray &answer) |
| void | computeVectorOf (const IntArray &dofIDMask, ValueModeType u, TimeStep *tStep, FloatArray &answer, bool padding=false) |
| void | computeBoundaryVectorOf (const IntArray &bNodes, const IntArray &dofIDMask, ValueModeType u, TimeStep *tStep, FloatArray &answer, bool padding=false) |
| void | computeVectorOf (PrimaryField &field, const IntArray &dofIDMask, ValueModeType u, TimeStep *tStep, FloatArray &answer, bool padding=false) |
| void | computeVectorOfPrescribed (ValueModeType u, TimeStep *tStep, FloatArray &answer) |
| void | computeVectorOfPrescribed (const IntArray &dofIDMask, ValueModeType type, TimeStep *tStep, FloatArray &answer) |
| virtual int | computeNumberOfGlobalDofs () |
| int | computeNumberOfPrimaryMasterDofs () |
| virtual bool | computeGtoLRotationMatrix (FloatMatrix &answer) |
| virtual bool | giveRotationMatrix (FloatMatrix &answer) |
| virtual bool | computeDofTransformationMatrix (FloatMatrix &answer, const IntArray &nodes, bool includeInternal) |
| virtual void | giveElementDofIDMask (IntArray &answer) const |
| virtual double | computeVolumeAreaOrLength () |
| Computes the volume, area or length of the element depending on its spatial dimension. | |
| double | computeMeanSize () |
| virtual double | computeVolume () |
| virtual double | computeArea () |
| virtual double | computeLength () |
| virtual IntArray | giveBoundaryEdgeNodes (int boundary, bool includeHierarchical=false) const |
| virtual IntArray | giveBoundarySurfaceNodes (int boundary, bool includeHierarchical=false) const |
| virtual IntArray | giveBoundaryNodes (int boundary) const |
| virtual std::unique_ptr< IntegrationRule > | giveBoundaryEdgeIntegrationRule (int order, int boundary) |
| virtual std::unique_ptr< IntegrationRule > | giveBoundarySurfaceIntegrationRule (int order, int boundary) |
| int | giveDofManagerNumber (int i) const |
| const IntArray & | giveDofManArray () const |
| void | addDofManager (DofManager *dMan) |
| DofManager * | giveDofManager (int i) const |
| Node * | giveNode (int i) const |
| virtual ElementSide * | giveSide (int i) const |
| virtual const FEInterpolation * | getGeometryInterpolation () const |
| virtual Material * | giveMaterial () |
| int | giveMaterialNumber () const |
| CrossSection * | giveCrossSection () |
| int | getActivityTimeFunctionNumber () |
| void | setActivityTimeFunctionNumber (int funcIndx) |
| void | setMaterial (int matIndx) |
| virtual void | setCrossSection (int csIndx) |
| virtual int | giveNumberOfDofManagers () const |
| void | setNumberOfDofManagers (int i) |
| Sets number of element dof managers. | |
| virtual int | giveNumberOfNodes () const |
| void | setDofManagers (const IntArray &dmans) |
| void | setDofManager (int id, int dm) |
| void | setBodyLoads (const IntArray &bodyLoads) |
| void | setIntegrationRules (std ::vector< std ::unique_ptr< IntegrationRule > > irlist) |
| virtual integrationDomain | giveIntegrationDomain () const |
| virtual int | giveIntegrationRuleLocalCodeNumbers (IntArray &answer, IntegrationRule &ie) |
| int | giveRegionNumber () |
| virtual void | updateInternalState (TimeStep *tStep) |
| virtual void | initializeYourself (TimeStep *timeStepWhenICApply) |
| int | checkConsistency () override |
| virtual bool | isActivated (TimeStep *tStep) |
| virtual bool | isCast (TimeStep *tStep) |
| virtual void | initForNewStep () |
| virtual Element_Geometry_Type | giveEdgeGeometryType (int id) const |
| Returns the receiver edge geometry type. | |
| virtual Element_Geometry_Type | giveSurfaceGeometryType (int id) const |
| Returns the receiver surface geometry type. | |
| virtual int | giveSpatialDimension () |
| virtual int | giveNumberOfBoundarySides () |
| Returns number of boundaries (entities of element_dimension-1: points, edges, surfaces). | |
| virtual int | giveNumberOfEdges () const |
| virtual int | giveNumberOfSurfaces () const |
| virtual int | giveDefaultIntegrationRule () const |
| virtual IntegrationRule * | giveDefaultIntegrationRulePtr () |
| int | giveNumberOfIntegrationRules () |
| virtual IntegrationRule * | giveIntegrationRule (int i) |
| std::vector< std ::unique_ptr< IntegrationRule > > & | giveIntegrationRulesArray () |
| virtual int | testElementExtension (ElementExtension ext) |
| virtual int | giveIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) |
| int | giveGlobalIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) |
| virtual double | giveLengthInDir (const FloatArray &normalToCrackPlane) |
| virtual double | giveCharacteristicLength (const FloatArray &normalToCrackPlane) |
| double | giveCharacteristicLengthForPlaneElements (const FloatArray &normalToCrackPlane) |
| double | giveCharacteristicLengthForAxisymmElements (const FloatArray &normalToCrackPlane) |
| virtual double | giveCharacteristicSize (GaussPoint *gp, FloatArray &normalToCrackPlane, ElementCharSizeMethod method) |
| virtual double | giveParentElSize () const |
| virtual void | updateBeforeNonlocalAverage (TimeStep *tStep) |
| virtual int | computeGlobalCoordinates (FloatArray &answer, const FloatArray &lcoords) |
| virtual bool | computeLocalCoordinates (FloatArray &answer, const FloatArray &gcoords) |
| virtual int | giveLocalCoordinateSystem (FloatMatrix &answer) |
| virtual void | giveLocalCoordinateSystemVector (InternalStateType isttype, FloatArray &answer) |
| virtual void | computeMidPlaneNormal (FloatArray &answer, const GaussPoint *gp) |
| virtual int | adaptiveMap (Domain *oldd, TimeStep *tStep) |
| virtual int | mapStateVariables (Domain &iOldDom, const TimeStep &iTStep) |
| virtual int | adaptiveUpdate (TimeStep *tStep) |
| virtual int | adaptiveFinish (TimeStep *tStep) |
| void | updateLocalNumbering (EntityRenumberingFunctor &f) override |
| template<class T> | |
| void | ipEvaluator (T *src, void(T ::*f)(GaussPoint *gp)) |
| Integration point evaluator, loops over receiver IP's and calls given function (passed as f parameter) on them. The IP is parameter to function f. | |
| template<class T, class S> | |
| void | ipEvaluator (T *src, void(T ::*f)(GaussPoint *, S &), S &_val) |
| Integration point evaluator, loops over receiver IP's and calls given function (passed as f parameter) on them. The IP is parameter to function f as well as additional array. | |
| virtual void | drawYourself (oofegGraphicContext &gc, TimeStep *tStep) |
| virtual void | drawAnnotation (oofegGraphicContext &gc, TimeStep *tStep) |
| virtual void | drawRawGeometry (oofegGraphicContext &gc, TimeStep *tStep) |
| virtual void | drawDeformedGeometry (oofegGraphicContext &gc, TimeStep *tStep, UnknownType) |
| virtual void | drawScalar (oofegGraphicContext &gc, TimeStep *tStep) |
| virtual void | drawSpecial (oofegGraphicContext &gc, TimeStep *tStep) |
| virtual void | giveLocalIntVarMaxMin (oofegGraphicContext &gc, TimeStep *tStep, double &emin, double &emax) |
| virtual int | giveInternalStateAtNode (FloatArray &answer, InternalStateType type, InternalStateMode mode, int node, TimeStep *tStep) |
| virtual int | giveInternalStateAtSide (FloatArray &answer, InternalStateType type, InternalStateMode mode, int side, TimeStep *tStep) |
| virtual void | showSparseMtrxStructure (CharType mtrx, oofegGraphicContext &gc, TimeStep *tStep) |
| Shows sparse structure. | |
| virtual void | showExtendedSparseMtrxStructure (CharType mtrx, oofegGraphicContext &gc, TimeStep *tStep) |
| Shows extended sparse structure (for example, due to nonlocal interactions for tangent stiffness). | |
| int | giveLabel () const |
| int | giveGlobalNumber () const |
| void | setGlobalNumber (int num) |
| elementParallelMode | giveParallelMode () const |
| void | setParallelMode (elementParallelMode _mode) |
| Sets parallel mode of element. | |
| virtual elementParallelMode | giveKnotSpanParallelMode (int) const |
| int | packUnknowns (DataStream &buff, TimeStep *tStep) |
| int | unpackAndUpdateUnknowns (DataStream &buff, TimeStep *tStep) |
| int | estimatePackSize (DataStream &buff) |
| const IntArray * | givePartitionList () const |
| void | setPartitionList (IntArray &pl) |
| virtual double | predictRelativeComputationalCost () |
| virtual double | giveRelativeSelfComputationalCost () |
| virtual double | predictRelativeRedistributionCost () |
| IntArray * | giveBodyLoadArray () |
| Returns array containing load numbers of loads acting on element. | |
| IntArray * | giveBoundaryLoadArray () |
| Returns array containing load numbers of boundary loads acting on element. | |
| void | initializeFrom (InputRecord &ir, int priority) override |
| void | initializeFinish () override |
| void | postInitialize () override |
| Performs post initialization steps. | |
| void | giveInputRecord (DynamicInputRecord &input) override |
| void | saveContext (DataStream &stream, ContextMode mode) override |
| void | restoreContext (DataStream &stream, ContextMode mode) override |
| void | printOutputAt (FILE *file, TimeStep *tStep) override |
| virtual const IntArray | giveLocation () |
| virtual void | recalculateCoordinates (int nodeNumber, FloatArray &coords) |
| void | setSharedEdgeID (int iedge, int globalID) |
| void | setSharedSurfaceID (int isurf, int globalID) |
| const IntArray * | giveSharedEdgeIDs () const |
| const IntArray * | giveSharedSurfaceIDs () const |
| Public Member Functions inherited from oofem::FEMComponent | |
| FEMComponent (int n, Domain *d) | |
| virtual | ~FEMComponent ()=default |
| Virtual destructor. | |
| Domain * | giveDomain () const |
| virtual void | setDomain (Domain *d) |
| int | giveNumber () const |
| void | setNumber (int num) |
| virtual void | initializeFrom (InputRecord &ir) |
| virtual void | printYourself () |
| Prints receiver state on stdout. Useful for debugging. | |
| std::string | errorInfo (const char *func) const |
| Returns string for prepending output (used by error reporting macros). | |
| Public Member Functions inherited from oofem::ZZNodalRecoveryModelInterface | |
| ZZNodalRecoveryModelInterface (Element *element) | |
| Constructor. | |
| virtual bool | ZZNodalRecoveryMI_computeNValProduct (FloatMatrix &answer, InternalStateType type, TimeStep *tStep) |
| virtual void | ZZNodalRecoveryMI_computeNNMatrix (FloatArray &answer, InternalStateType type) |
| Public Member Functions inherited from oofem::Interface | |
| Interface () | |
| Constructor. | |
| virtual | ~Interface () |
| Public Member Functions inherited from oofem::SpatialLocalizerInterface | |
| SpatialLocalizerInterface (Element *element) | |
| virtual int | SpatialLocalizerI_containsPoint (const FloatArray &coords) |
| int | SpatialLocalizerI_BBoxContainsPoint (const FloatArray &coords) |
| virtual void | SpatialLocalizerI_giveBBox (FloatArray &bb0, FloatArray &bb1) |
| virtual double | SpatialLocalizerI_giveClosestPoint (FloatArray &lcoords, FloatArray &closest, const FloatArray &gcoords) |
Protected Attributes | |
| std ::unique_ptr< ElementDofManager > | bubble |
| The extra dofs from the bubble. | |
| Protected Attributes inherited from oofem::Element | |
| int | numberOfDofMans |
| Number of dofmanagers. | |
| IntArray | dofManArray |
| Array containing dofmanager numbers. | |
| int | material |
| Number of associated material. | |
| int | crossSection |
| Number of associated cross section. | |
| IntArray | bodyLoadArray |
| IntArray | boundaryLoadArray |
| std::vector< std ::unique_ptr< IntegrationRule > > | integrationRulesArray |
| FloatMatrix | elemLocalCS |
| Transformation material matrix, used in orthotropic and anisotropic materials, global->local transformation. | |
| int | activityTimeFunction |
| Element activity time function. If defined, nonzero value indicates active receiver, zero value inactive element. | |
| int | globalNumber |
| int | numberOfGaussPoints |
| elementParallelMode | parallel_mode |
| Determines the parallel mode of the element. | |
| IntArray | partitions |
| IntArray | globalEdgeIDs |
| IntArray | globalSurfaceIDs |
| Protected Attributes inherited from oofem::FEMComponent | |
| int | number |
| Component number. | |
| Domain * | domain |
| Link to domain object, useful for communicating with other FEM components. | |
Static Protected Attributes | |
| static FEI2dTrLin | interp |
| Interpolation for pressure. | |
| static IntArray | momentum_ordering = {1, 2, 4, 5, 7, 8, 10, 11} |
| Ordering of dofs in element. Used to assemble the element stiffness. | |
| static IntArray | conservation_ordering = {3, 6, 9} |
| static IntArray | edge_ordering [3] |
| Ordering of dofs on edges. Used to assemble edge loads. | |
Additional Inherited Members | |
| Static Public Attributes inherited from oofem::Element | |
| static ParamKey | IPK_Element_mat |
| static ParamKey | IPK_Element_crosssect |
| static ParamKey | IPK_Element_nodes |
| static ParamKey | IPK_Element_bodyload |
| static ParamKey | IPK_Element_boundaryload |
| static ParamKey | IPK_Element_lcs |
| static ParamKey | IPK_Element_partitions |
| static ParamKey | IPK_Element_remote |
| static ParamKey | IPK_Element_activityTimeFunction |
| static ParamKey | IPK_Element_nip |
Triangular element for Stokes flow using Bubble basis function. Linear+Bubble interpolation of velocity, and linear interpolation of pressures. The element is exported as a linear triangle (bubble dofs are not exported). It can deal with nonlinear material models, but it is assumed that the fluid is without memory (which is usually the case).
Definition at line 57 of file tr1bubblestokes.h.
| oofem::Tr1BubbleStokes::Tr1BubbleStokes | ( | int | n, |
| Domain * | d ) |
Definition at line 66 of file tr1bubblestokes.C.
References bubble, oofem::FMElement::FMElement(), oofem::Element::numberOfDofMans, oofem::Element::numberOfGaussPoints, oofem::SpatialLocalizerInterface::SpatialLocalizerInterface(), and oofem::ZZNodalRecoveryModelInterface::ZZNodalRecoveryModelInterface().
|
overridevirtual |
Computes the contribution of the given load at the given boundary surface in global coordinate system. In general, the answer should include only relevant DOFs at the edge. The related is giveBoundaryLocationArray method, which should return corresponding code numbers.
| answer | Requested contribution of load. |
| load | Load to compute contribution from. |
| boundary | Boundary number. |
| type | Type of the contribution. |
| mode | Determines mode of answer. |
| tStep | Time step when answer is computed. |
| global | if true (default) then contribution is in global c.s., when false then contribution is in element local c.s. |
Reimplemented from oofem::Element.
Definition at line 256 of file tr1bubblestokes.C.
References oofem::FloatArray::assemble(), oofem::FloatArray::clear(), oofem::BoundaryLoad::computeValueAt(), edge_ordering, oofem::BoundaryLoad::giveApproxOrder(), oofem::Load::giveFormulationType(), oofem::BoundaryLoad::giveType(), interp, N, OOFEM_ERROR, oofem::FloatArray::resize(), oofem::GaussIntegrationRule::SetUpPointsOnLine(), oofem::TransmissionBC, and oofem::FloatArray::zero().
Referenced by computeExternalForcesVector().
| void oofem::Tr1BubbleStokes::computeExternalForcesVector | ( | FloatArray & | answer, |
| TimeStep * | tStep ) |
Definition at line 187 of file tr1bubblestokes.C.
References oofem::FloatArray::add(), oofem::Element::bodyLoadArray, oofem::BodyLoadBGT, oofem::Element::boundaryLoadArray, computeBoundarySurfaceLoadVector(), computeLoadVector(), oofem::FEMComponent::domain, oofem::EdgeLoadBGT, oofem::ForceLoadBVT, oofem::GeneralBoundaryCondition::giveBCGeoType(), oofem::GeneralBoundaryCondition::giveBCValType(), oofem::Element::giveBodyLoadArray(), oofem::FloatArray::resize(), and oofem::FloatArray::zero().
Referenced by giveCharacteristicVector().
|
overridevirtual |
Computes the unknown vector interpolated at the specified local coordinates. Used for exporting data and mapping fields.
| mode | Mode (total, increment, etc) of the output |
| tStep | Time step to evaluate at |
| lcoords | Local coordinates to evaluate at |
| answer | Results |
Reimplemented from oofem::Element.
Definition at line 389 of file tr1bubblestokes.C.
References oofem::FloatArray::at(), oofem::Element::computeVectorOf(), oofem::FloatArray::giveSize(), interp, oofem::FloatArray::resize(), and oofem::FloatArray::zero().
|
overridevirtual |
Initializes the array of integration rules member variable. Element can have multiple integration rules for different tasks. For example structural element family class uses this feature to implement transparent support for reduced and selective integration of some strain components. Must be defined by terminator classes.
Reimplemented from oofem::Element.
Definition at line 76 of file tr1bubblestokes.C.
References oofem::Element::giveCrossSection(), oofem::Element::integrationRulesArray, and oofem::Element::numberOfGaussPoints.
| void oofem::Tr1BubbleStokes::computeInternalForcesVector | ( | FloatArray & | answer, |
| TimeStep * | tStep ) |
Definition at line 140 of file tr1bubblestokes.C.
References oofem::FloatArray::assemble(), oofem::FluidDynamicMaterial::computeDeviatoricStress2D(), oofem::FMElement::computeVectorOfPressures(), oofem::FMElement::computeVectorOfVelocities(), conservation_ordering, oofem::dot(), oofem::Element::giveCrossSection(), oofem::Element::integrationRulesArray, interp, momentum_ordering, N, oofem::FloatArray::resize(), oofem::Tdot(), and oofem::FloatArray::zero().
Referenced by giveCharacteristicVector().
|
overridevirtual |
Computes the contribution of the given body load (volumetric).
| answer | Requested contribution of load. |
| load | Load to compute contribution from. |
| type | Type of the contribution. |
| mode | Determines mode of answer. |
| tStep | Time step when answer is computed. |
Reimplemented from oofem::Element.
Definition at line 221 of file tr1bubblestokes.C.
References oofem::FloatArray::assemble(), oofem::FloatArray::clear(), oofem::Load::computeComponentArrayAt(), oofem::Element::giveCrossSection(), oofem::FloatArray::giveSize(), oofem::Element::integrationRulesArray, interp, momentum_ordering, N, oofem::FloatArray::resize(), and oofem::FloatArray::zero().
Referenced by computeExternalForcesVector().
|
overridevirtual |
Computes or simply returns total number of element's local DOFs. Must be defined by particular element.
Reimplemented from oofem::Element.
Definition at line 85 of file tr1bubblestokes.C.
| void oofem::Tr1BubbleStokes::computeStiffnessMatrix | ( | FloatMatrix & | answer, |
| MatResponseMode | mode, | ||
| TimeStep * | tStep ) |
Definition at line 303 of file tr1bubblestokes.C.
References oofem::FloatMatrix::assemble(), oofem::FluidDynamicMaterial::computeTangents2D(), conservation_ordering, oofem::dot(), oofem::Element::giveCrossSection(), oofem::Element::integrationRulesArray, interp, momentum_ordering, N, oofem::FloatMatrix::resize(), oofem::Tdot(), oofem::transpose(), and oofem::FloatMatrix::zero().
Referenced by giveCharacteristicMatrix().
|
overridevirtual |
Returns volume related to given integration point. Used typically in subroutines, that perform integration over element volume. Should be implemented by particular elements.
| gp | Integration point for which volume is computed. |
Reimplemented from oofem::Element.
Definition at line 111 of file tr1bubblestokes.C.
References oofem::GaussPoint::giveNaturalCoordinates(), oofem::GaussPoint::giveWeight(), and interp.
|
overridevirtual |
Computes characteristic matrix of receiver of requested type in given time step.
| answer | Requested characteristic matrix (stiffness, tangent, ...). If element has no capability to compute requested type of characteristic matrix error function is invoked. |
| type | Id of characteristic component requested. |
| tStep | Time step when answer is computed. |
Reimplemented from oofem::Element.
Definition at line 130 of file tr1bubblestokes.C.
References computeStiffnessMatrix(), and OOFEM_ERROR.
|
overridevirtual |
Computes characteristic vector of receiver of requested type in given time step. If element has no capability to compute requested type of characteristic vector error function is invoked.
| answer | Requested characteristic vector. |
| type | Id of characteristic component requested. |
| mode | Determines mode of answer. |
| tStep | Time step when answer is computed. |
Reimplemented from oofem::Element.
Definition at line 117 of file tr1bubblestokes.C.
References computeExternalForcesVector(), computeInternalForcesVector(), and OOFEM_ERROR.
|
inlineoverridevirtual |
Reimplemented from oofem::Element.
Definition at line 91 of file tr1bubblestokes.h.
|
overridevirtual |
Returns dofmanager dof mask for node. This mask defines the dofs which are used by element in node. Mask influences the code number ordering for particular node. Code numbers are ordered according to node order and dofs belonging to particular node are ordered according to this mask. If element requests dofs using node mask which are not in node then error is generated. This masking allows node to be shared by different elements with different dofs in same node. Elements local code numbers are extracted from node using this mask. Must be defined by particular element.
| inode | Mask is computed for local dofmanager with inode number. |
| answer | Mask for node. |
Reimplemented from oofem::Element.
Definition at line 90 of file tr1bubblestokes.C.
|
inlineoverridevirtual |
Returns the element geometry type. This information is assumed to be of general interest, but it is required only for some specialized tasks.
Implements oofem::Element.
Definition at line 94 of file tr1bubblestokes.h.
|
inlineoverridevirtual |
Implements oofem::FEMComponent.
Definition at line 92 of file tr1bubblestokes.h.
References _IFT_Tr1BubbleStokes_Name.
|
overridevirtual |
Interface requesting service.
Reimplemented from oofem::FEMComponent.
Definition at line 375 of file tr1bubblestokes.C.
References oofem::SpatialLocalizerInterface::SpatialLocalizerInterface(), oofem::SpatialLocalizerInterfaceType, oofem::ZZNodalRecoveryModelInterface::ZZNodalRecoveryModelInterface(), and oofem::ZZNodalRecoveryModelInterfaceType.
|
overridevirtual |
Returns i-th internal element dof manager of the receiver
| i | Internal number of DOF. |
Reimplemented from oofem::Element.
Definition at line 106 of file tr1bubblestokes.C.
References bubble.
|
overridevirtual |
Returns internal dofmanager dof mask for node. This mask defines the dofs which are used by element in node. Mask influences the code number ordering for particular node. Code numbers are ordered according to node order and dofs belonging to particular node are ordered according to this mask. If element requests dofs using node mask which are not in node then error is generated. This masking allows node to be shared by different elements with different dofs in same node. Elements local code numbers are extracted from node using this mask. Must be defined by particular element.
| inode | Mask is computed for local dofmanager with inode number. |
| answer | mask for node. |
Reimplemented from oofem::Element.
Definition at line 95 of file tr1bubblestokes.C.
|
overridevirtual |
Reimplemented from oofem::Element.
Definition at line 358 of file tr1bubblestokes.C.
References interp.
|
overridevirtual |
Returns the interpolation for the specific dof id. Special elements which uses a mixed interpolation should reimplement this method.
| id | ID of the dof for the for the requested interpolation. |
Reimplemented from oofem::Element.
Definition at line 363 of file tr1bubblestokes.C.
References interp.
|
inlineoverridevirtual |
Returns material mode for receiver integration points. Should be specialized.
Reimplemented from oofem::Element.
Definition at line 93 of file tr1bubblestokes.h.
|
overridevirtual |
Reimplemented from oofem::Element.
Definition at line 100 of file tr1bubblestokes.C.
|
overridevirtual |
Updates element state after equilibrium in time step has been reached. Default implementation updates all integration rules defined by integrationRulesArray member variable. Doing this, all integration points and their material statuses are updated also. All temporary history variables, which now describe equilibrium state are copied into equilibrium ones. The existing internal state is used for update.
| tStep | Time step for newly reached state. |
Reimplemented from oofem::Element.
Definition at line 368 of file tr1bubblestokes.C.
|
protected |
The extra dofs from the bubble.
Definition at line 70 of file tr1bubblestokes.h.
Referenced by giveInternalDofManager(), and Tr1BubbleStokes().
|
staticprotected |
Definition at line 65 of file tr1bubblestokes.h.
Referenced by computeInternalForcesVector(), and computeStiffnessMatrix().
|
staticprotected |
Ordering of dofs on edges. Used to assemble edge loads.
Definition at line 67 of file tr1bubblestokes.h.
Referenced by computeBoundarySurfaceLoadVector().
|
staticprotected |
Interpolation for pressure.
Definition at line 63 of file tr1bubblestokes.h.
Referenced by computeBoundarySurfaceLoadVector(), computeField(), computeInternalForcesVector(), computeLoadVector(), computeStiffnessMatrix(), computeVolumeAround(), giveInterpolation(), and giveInterpolation().
|
staticprotected |
Ordering of dofs in element. Used to assemble the element stiffness.
Definition at line 65 of file tr1bubblestokes.h.
Referenced by computeInternalForcesVector(), computeLoadVector(), and computeStiffnessMatrix().