OOFEM
2.4
OOFEM.org - Object Oriented Finite Element Solver
|
Tetrahedral Taylor-Hood element for Stokes flow. More...
#include <tet21stokes.h>
Public Member Functions | |
Tet21Stokes (int n, Domain *d) | |
virtual | ~Tet21Stokes () |
virtual void | computeGaussPoints () |
Initializes the array of integration rules member variable. More... | |
virtual void | giveCharacteristicVector (FloatArray &answer, CharType type, ValueModeType mode, TimeStep *tStep) |
Computes characteristic vector of receiver of requested type in given time step. More... | |
virtual void | giveCharacteristicMatrix (FloatMatrix &answer, CharType type, TimeStep *tStep) |
Computes characteristic matrix of receiver of requested type in given time step. More... | |
void | computeInternalForcesVector (FloatArray &answer, TimeStep *tStep) |
void | computeStiffnessMatrix (FloatMatrix &answer, MatResponseMode mode, TimeStep *tStep) |
void | computeExternalForcesVector (FloatArray &answer, TimeStep *tStep) |
virtual void | computeLoadVector (FloatArray &answer, BodyLoad *load, CharType type, ValueModeType mode, TimeStep *tStep) |
Computes the contribution of the given body load (volumetric). More... | |
virtual void | computeBoundarySurfaceLoadVector (FloatArray &answer, BoundaryLoad *load, int boundary, CharType type, ValueModeType mode, TimeStep *tStep, bool global=true) |
Computes the contribution of the given load at the given boundary surface in global coordinate system. More... | |
virtual const char * | giveClassName () const |
virtual const char * | giveInputRecordName () const |
virtual MaterialMode | giveMaterialMode () |
Returns material mode for receiver integration points. More... | |
virtual int | computeNumberOfDofs () |
Computes or simply returns total number of element's local DOFs. More... | |
virtual FEInterpolation * | giveInterpolation () const |
virtual FEInterpolation * | giveInterpolation (DofIDItem id) const |
Returns the interpolation for the specific dof id. More... | |
virtual void | giveDofManDofIDMask (int inode, IntArray &answer) const |
Gives the dof ID mask for the element. More... | |
virtual void | updateYourself (TimeStep *tStep) |
Updates element state after equilibrium in time step has been reached. More... | |
virtual Interface * | giveInterface (InterfaceType it) |
Interface requesting service. More... | |
virtual void | EIPrimaryUnknownMI_computePrimaryUnknownVectorAtLocal (ValueModeType u, TimeStep *tStep, const FloatArray &coords, FloatArray &answer) |
Computes the element vector of primary unknowns at given point in the local coordinate system. More... | |
virtual void | NodalAveragingRecoveryMI_computeNodalValue (FloatArray &answer, int node, InternalStateType type, TimeStep *tStep) |
Computes the element value in given node. More... | |
Public Member Functions inherited from oofem::FMElement | |
FMElement (int n, Domain *aDomain) | |
virtual | ~FMElement () |
virtual void | updateStabilizationCoeffs (TimeStep *tStep) |
Updates the stabilization coefficients used for CBS and SUPG algorithms. More... | |
virtual void | computeVectorOfVelocities (ValueModeType mode, TimeStep *tStep, FloatArray &velocities) |
virtual void | computeVectorOfPressures (ValueModeType mode, TimeStep *tStep, FloatArray &pressures) |
Public Member Functions inherited from oofem::Element | |
Element (int n, Domain *aDomain) | |
Constructor. More... | |
Element (const Element &src)=delete | |
Element & | operator= (const Element &src)=delete |
virtual | ~Element () |
Virtual destructor. More... | |
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) |
Returns internal state variable (like stress,strain) at node of element in Reduced form, the way how is obtained is dependent on InternalValueType. More... | |
virtual int | giveInternalStateAtSide (FloatArray &answer, InternalStateType type, InternalStateMode mode, int side, TimeStep *tStep) |
Returns internal state variable (like stress,strain) at side of element in Reduced form If side is possessing DOFs, otherwise recover techniques will not work due to absence of side-shape functions. More... | |
virtual void | showSparseMtrxStructure (CharType mtrx, oofegGraphicContext &gc, TimeStep *tStep) |
Shows sparse structure. More... | |
virtual void | showExtendedSparseMtrxStructure (CharType mtrx, oofegGraphicContext &gc, TimeStep *tStep) |
Shows extended sparse structure (for example, due to nonlocal interactions for tangent stiffness) More... | |
int | giveLabel () const |
int | giveGlobalNumber () const |
void | setGlobalNumber (int num) |
Sets receiver globally unique number. More... | |
elementParallelMode | giveParallelMode () const |
Return elementParallelMode of receiver. More... | |
void | setParallelMode (elementParallelMode _mode) |
Sets parallel mode of element. More... | |
virtual elementParallelMode | giveKnotSpanParallelMode (int) const |
Returns the parallel mode for particular knot span of the receiver. More... | |
int | packUnknowns (DataStream &buff, TimeStep *tStep) |
Pack all necessary data of element (according to its parallel_mode) integration points into given communication buffer. More... | |
int | unpackAndUpdateUnknowns (DataStream &buff, TimeStep *tStep) |
Unpack and updates all necessary data of element (according to its parallel_mode) integration points into given communication buffer. More... | |
int | estimatePackSize (DataStream &buff) |
Estimates the necessary pack size to hold all packed data of receiver. More... | |
const IntArray * | givePartitionList () const |
Returns partition list of receiver. More... | |
void | setPartitionList (IntArray &pl) |
Sets partition list of receiver. More... | |
virtual double | predictRelativeComputationalCost () |
Returns the weight representing relative computational cost of receiver The reference element is triangular plane stress element with linear approximation, single integration point and linear isotropic material. More... | |
virtual double | giveRelativeSelfComputationalCost () |
Returns the weight representing relative computational cost of receiver The reference element is triangular plane stress element. More... | |
virtual double | predictRelativeRedistributionCost () |
Returns the relative redistribution cost of the receiver. More... | |
IntArray * | giveBodyLoadArray () |
Returns array containing load numbers of loads acting on element. More... | |
IntArray * | giveBoundaryLoadArray () |
Returns array containing load numbers of boundary loads acting on element. More... | |
virtual IRResultType | initializeFrom (InputRecord *ir) |
Initializes receiver according to object description stored in input record. More... | |
virtual void | giveInputRecord (DynamicInputRecord &input) |
Setups the input record string of receiver. More... | |
virtual contextIOResultType | saveContext (DataStream &stream, ContextMode mode, void *obj=NULL) |
Stores receiver state to output stream. More... | |
virtual contextIOResultType | restoreContext (DataStream &stream, ContextMode mode, void *obj=NULL) |
Restores the receiver state previously written in stream. More... | |
virtual void | printOutputAt (FILE *file, TimeStep *tStep) |
Prints output of receiver to stream, for given time step. More... | |
void | giveLocationArray (IntArray &locationArray, const UnknownNumberingScheme &s, IntArray *dofIds=NULL) const |
Returns the location array (array of code numbers) of receiver for given numbering scheme. More... | |
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) |
Returns the location array for the boundary of the element. More... | |
virtual void | giveBoundaryLocationArray (IntArray &locationArray, const IntArray &bNodes, const IntArray &dofIDMask, const UnknownNumberingScheme &s, IntArray *dofIds=NULL) |
virtual int | giveNumberOfDofs () |
virtual int | giveNumberOfInternalDofManagers () const |
virtual DofManager * | giveInternalDofManager (int i) const |
Returns i-th internal element dof manager of the receiver. More... | |
virtual double | giveCharacteristicValue (CharType type, TimeStep *tStep) |
Computes characteristic value of receiver of requested type in given time step. More... | |
virtual void | computeTangentFromSurfaceLoad (FloatMatrix &answer, SurfaceLoad *load, int boundary, MatResponseMode rmode, TimeStep *tStep) |
Computes the tangent contribution of the given load at the given boundary. More... | |
virtual void | computeTangentFromEdgeLoad (FloatMatrix &answer, EdgeLoad *load, int boundary, MatResponseMode rmode, TimeStep *tStep) |
Computes the tangent contribution of the given load at the given boundary. More... | |
virtual void | computeBoundaryEdgeLoadVector (FloatArray &answer, BoundaryLoad *load, int edge, CharType type, ValueModeType mode, TimeStep *tStep, bool global=true) |
Computes the contribution of the given load at the given boundary edge. More... | |
const IntArray & | giveBodyLoadList () const |
Returns receiver list of bodyloads. More... | |
const IntArray & | giveBoundaryLoadList () const |
Returns receiver list of boundary loads. More... | |
void | computeVectorOf (ValueModeType u, TimeStep *tStep, FloatArray &answer) |
Returns local vector of unknowns. More... | |
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) |
Boundary version of computeVectorOf. More... | |
void | computeVectorOf (PrimaryField &field, const IntArray &dofIDMask, ValueModeType u, TimeStep *tStep, FloatArray &answer, bool padding=false) |
Returns local vector of unknowns. More... | |
void | computeVectorOfPrescribed (ValueModeType u, TimeStep *tStep, FloatArray &answer) |
Returns local vector of prescribed unknowns. More... | |
void | computeVectorOfPrescribed (const IntArray &dofIDMask, ValueModeType type, TimeStep *tStep, FloatArray &answer) |
Returns local vector of prescribed unknowns. More... | |
virtual int | computeNumberOfGlobalDofs () |
Computes the total number of element's global dofs. More... | |
int | computeNumberOfPrimaryMasterDofs () |
Computes the total number of element's primary master DOFs. More... | |
virtual bool | computeGtoLRotationMatrix (FloatMatrix &answer) |
Returns transformation matrix from global c.s. More... | |
virtual bool | giveRotationMatrix (FloatMatrix &answer) |
Transformation matrices updates rotation matrix between element-local and primary DOFs, taking into account nodal c.s. More... | |
virtual bool | computeDofTransformationMatrix (FloatMatrix &answer, const IntArray &nodes, bool includeInternal) |
Returns transformation matrix for DOFs from global coordinate system to local coordinate system in nodes. More... | |
virtual void | giveInternalDofManDofIDMask (int inode, IntArray &answer) const |
Returns internal dofmanager dof mask for node. More... | |
virtual void | giveElementDofIDMask (IntArray &answer) const |
Returns element dof mask for node. More... | |
virtual void | computeField (ValueModeType mode, TimeStep *tStep, const FloatArray &lcoords, FloatArray &answer) |
Computes the unknown vector interpolated at the specified local coordinates. More... | |
virtual double | computeVolumeAround (GaussPoint *gp) |
Returns volume related to given integration point. More... | |
virtual double | computeVolumeAreaOrLength () |
Computes the volume, area or length of the element depending on its spatial dimension. More... | |
double | computeMeanSize () |
Computes the size of the element defined as its length. More... | |
virtual double | computeVolume () |
Computes the volume. More... | |
virtual double | computeArea () |
Computes the area (zero for all but 2d geometries). More... | |
virtual double | computeLength () |
Computes the length (zero for all but 1D geometries) More... | |
virtual void | giveBoundaryEdgeNodes (IntArray &bNodes, int boundary) |
Returns list of receiver boundary nodes for given edge. More... | |
virtual void | giveBoundarySurfaceNodes (IntArray &bNodes, int boundary) |
Returns list of receiver boundary nodes for given surface. More... | |
virtual IntegrationRule * | giveBoundaryEdgeIntegrationRule (int order, int boundary) |
Returns boundary edge integration rule. More... | |
virtual IntegrationRule * | giveBoundarySurfaceIntegrationRule (int order, int boundary) |
Returns boundary surface integration rule. More... | |
int | giveDofManagerNumber (int i) const |
Translates local to global indices for dof managers. More... | |
const IntArray & | giveDofManArray () const |
void | addDofManager (DofManager *dMan) |
DofManager * | giveDofManager (int i) const |
Node * | giveNode (int i) const |
Returns reference to the i-th node of element. More... | |
virtual ElementSide * | giveSide (int i) const |
Returns reference to the i-th side of element. More... | |
virtual Material * | giveMaterial () |
int | giveMaterialNumber () const |
CrossSection * | giveCrossSection () |
void | setMaterial (int matIndx) |
Sets the material of receiver. More... | |
virtual void | setCrossSection (int csIndx) |
Sets the cross section model of receiver. More... | |
virtual int | giveNumberOfDofManagers () const |
virtual int | giveNumberOfNodes () const |
Returns number of nodes of receiver. More... | |
void | setDofManagers (const IntArray &dmans) |
Sets receiver dofManagers. More... | |
void | setBodyLoads (const IntArray &bodyLoads) |
Sets receiver bodyLoadArray. More... | |
void | setIntegrationRules (std::vector< std::unique_ptr< IntegrationRule > > irlist) |
Sets integration rules. More... | |
virtual integrationDomain | giveIntegrationDomain () const |
Returns integration domain for receiver, used to initialize integration point over receiver volume. More... | |
virtual int | giveIntegrationRuleLocalCodeNumbers (IntArray &answer, IntegrationRule &ie) |
Assembles the code numbers of given integration element (sub-patch) This is done by obtaining list of nonzero shape functions and by collecting the code numbers of nodes corresponding to these shape functions. More... | |
int | giveRegionNumber () |
virtual void | postInitialize () |
Performs post initialization steps. More... | |
virtual void | updateInternalState (TimeStep *tStep) |
Updates element state after equilibrium in time step has been reached. More... | |
virtual void | initializeYourself (TimeStep *timeStepWhenICApply) |
Initialization according to state given by initial conditions. More... | |
virtual int | checkConsistency () |
Performs consistency check. More... | |
virtual bool | isActivated (TimeStep *tStep) |
virtual bool | isCast (TimeStep *tStep) |
virtual void | initForNewStep () |
Initializes receivers state to new time step. More... | |
virtual Element_Geometry_Type | giveGeometryType () const |
Returns the element geometry type. More... | |
virtual int | giveSpatialDimension () |
Returns the element spatial dimension (1, 2, or 3). More... | |
virtual int | giveNumberOfBoundarySides () |
virtual int | giveDefaultIntegrationRule () const |
Returns id of default integration rule. More... | |
virtual IntegrationRule * | giveDefaultIntegrationRulePtr () |
Access method for default integration rule. More... | |
int | giveNumberOfIntegrationRules () |
virtual IntegrationRule * | giveIntegrationRule (int i) |
std::vector< std::unique_ptr< IntegrationRule > > & | giveIntegrationRulesArray () |
virtual int | testElementExtension (ElementExtension ext) |
Tests if the element implements required extension. More... | |
virtual int | giveIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) |
Returns the integration point corresponding value in full form. More... | |
int | giveGlobalIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) |
virtual double | giveLengthInDir (const FloatArray &normalToCrackPlane) |
Default implementation returns length of element projection into specified direction. More... | |
virtual double | giveCharacteristicLength (const FloatArray &normalToCrackPlane) |
Returns the size of element in the given direction, in some cases adjusted (e.g. More... | |
double | giveCharacteristicLengthForPlaneElements (const FloatArray &normalToCrackPlane) |
Returns the size of element in the given direction if the direction is in the XY plane, otherwise gives the mean size defined as the square root of the element area. More... | |
double | giveCharacteristicLengthForAxisymmElements (const FloatArray &normalToCrackPlane) |
Returns the size of an axisymmetric element in the given direction if the direction is in the XY plane, otherwise gives the mean distance vrom the symmetry axis multiplied by pi. More... | |
virtual double | giveCharacteristicSize (GaussPoint *gp, FloatArray &normalToCrackPlane, ElementCharSizeMethod method) |
Returns characteristic element size for a given integration point and given direction. More... | |
virtual double | giveParentElSize () const |
Returns the size (length, area or volume depending on element type) of the parent element. More... | |
virtual void | updateBeforeNonlocalAverage (TimeStep *tStep) |
Updates internal element state (in all integration points of receiver) before nonlocal averaging takes place. More... | |
virtual int | computeGlobalCoordinates (FloatArray &answer, const FloatArray &lcoords) |
Computes the global coordinates from given element's local coordinates. More... | |
virtual bool | computeLocalCoordinates (FloatArray &answer, const FloatArray &gcoords) |
Computes the element local coordinates from given global coordinates. More... | |
virtual int | giveLocalCoordinateSystem (FloatMatrix &answer) |
Returns local coordinate system of receiver Required by material models with ortho- and anisotrophy. More... | |
virtual void | computeMidPlaneNormal (FloatArray &answer, const GaussPoint *gp) |
Computes mid-plane normal of receiver at integration point. More... | |
virtual int | adaptiveMap (Domain *oldd, TimeStep *tStep) |
Initializes the internal state variables stored in all IPs according to state in given domain. More... | |
virtual int | mapStateVariables (Domain &iOldDom, const TimeStep &iTStep) |
Maps the internal state variables stored in all IPs from the old domain to the new domain. More... | |
virtual int | adaptiveUpdate (TimeStep *tStep) |
Updates the internal state variables stored in all IPs according to already mapped state. More... | |
virtual int | adaptiveFinish (TimeStep *tStep) |
Finishes the mapping for given time step. More... | |
virtual void | updateLocalNumbering (EntityRenumberingFunctor &f) |
Local renumbering support. More... | |
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. More... | |
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. More... | |
Public Member Functions inherited from oofem::FEMComponent | |
FEMComponent (int n, Domain *d) | |
Regular constructor, creates component with given number and belonging to given domain. More... | |
virtual | ~FEMComponent () |
Virtual destructor. More... | |
Domain * | giveDomain () const |
virtual void | setDomain (Domain *d) |
Sets associated Domain. More... | |
int | giveNumber () const |
void | setNumber (int num) |
Sets number of receiver. More... | |
virtual void | printYourself () |
Prints receiver state on stdout. Useful for debugging. More... | |
std::string | errorInfo (const char *func) const |
Returns string for prepending output (used by error reporting macros). More... | |
Public Member Functions inherited from oofem::NodalAveragingRecoveryModelInterface | |
NodalAveragingRecoveryModelInterface () | |
Constructor. More... | |
Public Member Functions inherited from oofem::Interface | |
Interface () | |
Constructor. More... | |
virtual | ~Interface () |
Public Member Functions inherited from oofem::SpatialLocalizerInterface | |
SpatialLocalizerInterface (Element *element) | |
virtual int | SpatialLocalizerI_containsPoint (const FloatArray &coords) |
Checks if element contains specified coordinate. More... | |
int | SpatialLocalizerI_BBoxContainsPoint (const FloatArray &coords) |
Creates a bounding box of the nodes and checks if it includes the given coordinate. More... | |
virtual void | SpatialLocalizerI_giveBBox (FloatArray &bb0, FloatArray &bb1) |
Creates a bounding box of the nodes and checks if it includes the given coordinate. More... | |
virtual double | SpatialLocalizerI_giveClosestPoint (FloatArray &lcoords, FloatArray &closest, const FloatArray &gcoords) |
Gives the closest point on the element. More... | |
Public Member Functions inherited from oofem::EIPrimaryUnknownMapperInterface | |
EIPrimaryUnknownMapperInterface () | |
Static Protected Attributes | |
static FEI3dTetLin | interpolation_lin |
Interpolation for pressure. More... | |
static FEI3dTetQuad | interpolation_quad |
Interpolation for geometry and velocity. More... | |
static IntArray | momentum_ordering |
Ordering of momentum balance dofs in element. Used to assemble the element stiffness. More... | |
static IntArray | conservation_ordering = {4, 8, 12, 16} |
Ordering of conservation dofs in element. Used to assemble the element stiffness. More... | |
static IntArray | surf_ordering [4] |
Ordering of dofs on surfaces. Used to assemble edge loads (only momentum balance) More... | |
Additional Inherited Members | |
Protected Attributes inherited from oofem::Element | |
int | numberOfDofMans |
Number of dofmanagers. More... | |
IntArray | dofManArray |
Array containing dofmanager numbers. More... | |
int | material |
Number of associated material. More... | |
int | crossSection |
Number of associated cross section. More... | |
IntArray | bodyLoadArray |
Array containing indexes of loads (body loads and boundary loads are kept separately), that apply on receiver. More... | |
IntArray | boundaryLoadArray |
std::vector< std::unique_ptr< IntegrationRule > > | integrationRulesArray |
List of integration rules of receiver (each integration rule contains associated integration points also). More... | |
FloatMatrix | elemLocalCS |
Transformation material matrix, used in orthotropic and anisotropic materials, global->local transformation. More... | |
int | activityTimeFunction |
Element activity time function. If defined, nonzero value indicates active receiver, zero value inactive element. More... | |
int | globalNumber |
In parallel mode, globalNumber contains globally unique DoFManager number. More... | |
int | numberOfGaussPoints |
Number of integration points as specified by nip. More... | |
elementParallelMode | parallel_mode |
Determines the parallel mode of the element. More... | |
IntArray | partitions |
List of partition sharing the shared element or remote partition containing remote element counterpart. More... | |
Protected Attributes inherited from oofem::FEMComponent | |
int | number |
Component number. More... | |
Domain * | domain |
Link to domain object, useful for communicating with other FEM components. More... | |
Tetrahedral Taylor-Hood element for Stokes flow.
Quadratic interpolation of geometry and velocity, and linear interpolation of pressures.
Definition at line 55 of file tet21stokes.h.
oofem::Tet21Stokes::Tet21Stokes | ( | int | n, |
Domain * | d | ||
) |
Definition at line 71 of file tet21stokes.C.
References oofem::Element::numberOfDofMans, and oofem::Element::numberOfGaussPoints.
|
virtual |
Definition at line 77 of file tet21stokes.C.
|
virtual |
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 235 of file tet21stokes.C.
References oofem::FloatArray::assemble(), oofem::FloatArray::clear(), oofem::BoundaryLoad::computeValueAt(), oofem::Load::FT_Entity, oofem::BoundaryLoad::giveApproxOrder(), oofem::Load::giveFormulationType(), oofem::FloatArray::giveSize(), oofem::BoundaryLoad::giveType(), interpolation_quad, N, OOFEM_ERROR, oofem::FloatArray::resize(), oofem::GaussIntegrationRule::SetUpPointsOnTriangle(), surf_ordering, oofem::FEI3dTetQuad::surfaceEvalN(), oofem::FEI3dTetQuad::surfaceGiveTransformationJacobian(), oofem::FEI3dTetQuad::surfaceLocal2global(), oofem::TransmissionBC, and oofem::FloatArray::zero().
Referenced by computeExternalForcesVector().
void oofem::Tet21Stokes::computeExternalForcesVector | ( | FloatArray & | answer, |
TimeStep * | tStep | ||
) |
Definition at line 166 of file tet21stokes.C.
References oofem::FloatArray::add(), oofem::IntArray::at(), oofem::Element::bodyLoadArray, oofem::BodyLoadBGT, oofem::Element::boundaryLoadArray, oofem::FloatArray::clear(), computeBoundarySurfaceLoadVector(), computeLoadVector(), oofem::FEMComponent::domain, oofem::ForceLoadBVT, oofem::GeneralBoundaryCondition::giveBCGeoType(), oofem::GeneralBoundaryCondition::giveBCValType(), oofem::Element::giveBodyLoadArray(), oofem::Domain::giveLoad(), oofem::IntArray::giveSize(), and oofem::SurfaceLoadBGT.
Referenced by giveCharacteristicVector().
|
virtual |
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 80 of file tet21stokes.C.
References oofem::Element::giveCrossSection(), oofem::Element::integrationRulesArray, oofem::Element::numberOfGaussPoints, and oofem::CrossSection::setupIntegrationPoints().
void oofem::Tet21Stokes::computeInternalForcesVector | ( | FloatArray & | answer, |
TimeStep * | tStep | ||
) |
Definition at line 127 of file tet21stokes.C.
References oofem::FloatArray::add(), oofem::FloatArray::assemble(), oofem::FloatArray::beProductOf(), oofem::FluidDynamicMaterial::computeDeviatoricStressVector(), oofem::FMElement::computeVectorOfPressures(), oofem::FMElement::computeVectorOfVelocities(), conservation_ordering, oofem::FloatArray::dotProduct(), oofem::FEI3dTetQuad::evaldNdx(), oofem::FEI3dTetLin::evalN(), oofem::Element::giveCrossSection(), oofem::FloatMatrix::giveNumberOfRows(), oofem::Element::integrationRulesArray, interpolation_lin, interpolation_quad, momentum_ordering, oofem::FloatArray::plusProduct(), oofem::FloatArray::resize(), oofem::FloatArray::zero(), and oofem::FloatMatrix::zero().
Referenced by giveCharacteristicVector().
|
virtual |
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 202 of file tet21stokes.C.
References oofem::FloatArray::assemble(), oofem::FloatArray::clear(), oofem::Load::computeComponentArrayAt(), oofem::FEI3dTetQuad::evalN(), oofem::Element::giveCrossSection(), oofem::FloatArray::giveSize(), oofem::FEInterpolation::giveTransformationJacobian(), oofem::Element::integrationRulesArray, interpolation_quad, momentum_ordering, N, oofem::FloatArray::resize(), and oofem::FloatArray::zero().
Referenced by computeExternalForcesVector().
|
virtual |
Computes or simply returns total number of element's local DOFs.
Must be defined by particular element.
Reimplemented from oofem::Element.
Definition at line 89 of file tet21stokes.C.
Referenced by giveMaterialMode().
void oofem::Tet21Stokes::computeStiffnessMatrix | ( | FloatMatrix & | answer, |
MatResponseMode | mode, | ||
TimeStep * | tStep | ||
) |
Definition at line 283 of file tet21stokes.C.
References oofem::FloatMatrix::add(), oofem::FloatMatrix::assemble(), oofem::FloatMatrix::beProductOf(), oofem::FloatArray::beTProductOf(), oofem::FloatMatrix::beTranspositionOf(), conservation_ordering, oofem::FEI3dTetQuad::evaldNdx(), oofem::FEI3dTetLin::evalN(), oofem::Element::giveCrossSection(), oofem::FloatMatrix::giveNumberOfRows(), oofem::FluidDynamicMaterial::giveStiffnessMatrices(), oofem::Element::integrationRulesArray, interpolation_lin, interpolation_quad, momentum_ordering, oofem::FloatMatrix::plusDyadSymmUpper(), oofem::FloatMatrix::plusDyadUnsym(), oofem::FloatMatrix::plusProductSymmUpper(), oofem::FloatMatrix::resize(), oofem::FloatMatrix::symmetrized(), oofem::FloatArray::zero(), and oofem::FloatMatrix::zero().
Referenced by giveCharacteristicMatrix().
|
virtual |
Computes the element vector of primary unknowns at given point in the local coordinate system.
mode | Identifies mode of unknown (eg. total value or velocity of unknown). |
tStep | Time step, when vector of unknowns is requested. |
lcoords | Local coordinates of point of interest. |
answer | Vector of unknowns. |
Reimplemented from oofem::EIPrimaryUnknownMapperInterface.
Definition at line 380 of file tet21stokes.C.
References oofem::FloatArray::at(), oofem::FEI3dTetLin::evalN(), oofem::FEI3dTetQuad::evalN(), oofem::DofManager::giveDofWithID(), oofem::Element::giveNode(), oofem::FloatArray::giveSize(), oofem::Dof::giveUnknown(), interpolation_lin, interpolation_quad, oofem::FloatArray::resize(), and oofem::FloatArray::zero().
Referenced by giveMaterialMode().
|
virtual |
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 116 of file tet21stokes.C.
References computeStiffnessMatrix(), and OOFEM_ERROR.
|
virtual |
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 103 of file tet21stokes.C.
References computeExternalForcesVector(), computeInternalForcesVector(), and OOFEM_ERROR.
|
inlinevirtual |
Reimplemented from oofem::Element.
Definition at line 87 of file tet21stokes.h.
|
virtual |
Gives the dof ID mask for the element.
This element (Taylor-Hood) has V_u, V_v, P_f in node corner and V_u, V_v in edge nodes.
inode | Node to check. |
answer | List of dof IDs. |
Reimplemented from oofem::Element.
Definition at line 94 of file tet21stokes.C.
Referenced by giveMaterialMode().
|
inlinevirtual |
Implements oofem::FEMComponent.
Definition at line 88 of file tet21stokes.h.
References _IFT_Tet21Stokes_Name.
|
virtual |
Interface requesting service.
Reimplemented from oofem::FEMComponent.
Definition at line 363 of file tet21stokes.C.
References oofem::EIPrimaryUnknownMapperInterfaceType, oofem::FEMComponent::giveInterface(), oofem::NodalAveragingRecoveryModelInterfaceType, and oofem::SpatialLocalizerInterfaceType.
Referenced by giveMaterialMode().
|
virtual |
Reimplemented from oofem::Element.
Definition at line 342 of file tet21stokes.C.
References interpolation_quad.
Referenced by giveMaterialMode().
|
virtual |
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 347 of file tet21stokes.C.
References interpolation_lin, and interpolation_quad.
|
inlinevirtual |
Returns material mode for receiver integration points.
Should be specialized.
Reimplemented from oofem::Element.
Definition at line 89 of file tet21stokes.h.
References computeNumberOfDofs(), EIPrimaryUnknownMI_computePrimaryUnknownVectorAtLocal(), giveDofManDofIDMask(), giveInterface(), giveInterpolation(), NodalAveragingRecoveryMI_computeNodalValue(), and updateYourself().
|
virtual |
Computes the element value in given node.
answer | Contains the result. |
node | Element node number. |
type | Determines the type of internal variable to be recovered. |
tStep | Time step. |
Implements oofem::NodalAveragingRecoveryModelInterface.
Definition at line 399 of file tet21stokes.C.
References oofem::IntArray::at(), oofem::FloatArray::at(), oofem::FloatArray::clear(), oofem::FEI3dTetQuad::computeLocalEdgeMapping(), oofem::DofManager::giveDofWithID(), oofem::Element::giveNode(), oofem::Dof::giveUnknown(), interpolation_quad, and oofem::FloatArray::resize().
Referenced by giveMaterialMode().
|
virtual |
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 356 of file tet21stokes.C.
References oofem::Element::updateYourself().
Referenced by giveMaterialMode().
|
staticprotected |
Ordering of conservation dofs in element. Used to assemble the element stiffness.
Definition at line 68 of file tet21stokes.h.
Referenced by computeInternalForcesVector(), and computeStiffnessMatrix().
|
staticprotected |
Interpolation for pressure.
Definition at line 62 of file tet21stokes.h.
Referenced by computeInternalForcesVector(), computeStiffnessMatrix(), EIPrimaryUnknownMI_computePrimaryUnknownVectorAtLocal(), and giveInterpolation().
|
staticprotected |
Interpolation for geometry and velocity.
Definition at line 64 of file tet21stokes.h.
Referenced by computeBoundarySurfaceLoadVector(), computeInternalForcesVector(), computeLoadVector(), computeStiffnessMatrix(), EIPrimaryUnknownMI_computePrimaryUnknownVectorAtLocal(), giveInterpolation(), and NodalAveragingRecoveryMI_computeNodalValue().
|
staticprotected |
Ordering of momentum balance dofs in element. Used to assemble the element stiffness.
Definition at line 66 of file tet21stokes.h.
Referenced by computeInternalForcesVector(), computeLoadVector(), and computeStiffnessMatrix().
|
staticprotected |
Ordering of dofs on surfaces. Used to assemble edge loads (only momentum balance)
Definition at line 70 of file tet21stokes.h.
Referenced by computeBoundarySurfaceLoadVector().