|
OOFEM 3.0
|
#include <supgelement2.h>
Public Member Functions | |
| SUPGElement2 (int n, Domain *aDomain) | |
| void | initializeFrom (InputRecord &ir, int priority) override |
| void | giveInputRecord (DynamicInputRecord &input) override |
| void | giveCharacteristicMatrix (FloatMatrix &answer, CharType, TimeStep *tStep) override |
| void | giveCharacteristicVector (FloatArray &answer, CharType, ValueModeType, TimeStep *tStep) override |
| void | updateElementForNewInterfacePosition (TimeStep *tStep) override |
| void | computeAccelerationTerm_MB (FloatMatrix &answer, TimeStep *tStep) override |
| void | computeAdvectionTerm_MB (FloatArray &answer, TimeStep *tStep) override |
| void | computeAdvectionDerivativeTerm_MB (FloatMatrix &answer, TimeStep *tStep) override |
| void | computeDiffusionTerm_MB (FloatArray &answer, TimeStep *tStep) override |
| void | computeDiffusionDerivativeTerm_MB (FloatMatrix &answer, MatResponseMode mode, TimeStep *tStep) override |
| void | computePressureTerm_MB (FloatMatrix &answer, TimeStep *tStep) override |
| void | computeLSICStabilizationTerm_MB (FloatMatrix &answer, TimeStep *tStep) override |
| void | computeLinearAdvectionTerm_MC (FloatMatrix &answer, TimeStep *tStep) override |
| void | computeAdvectionTerm_MC (FloatArray &answer, TimeStep *tStep) override |
| void | computeAdvectionDerivativeTerm_MC (FloatMatrix &answer, TimeStep *tStep) override |
| void | computeDiffusionDerivativeTerm_MC (FloatMatrix &answer, TimeStep *tStep) override |
| void | computeDiffusionTerm_MC (FloatArray &answer, TimeStep *tStep) override |
| void | computeAccelerationTerm_MC (FloatMatrix &answer, TimeStep *tStep) override |
| void | computePressureTerm_MC (FloatMatrix &answer, TimeStep *tStep) override |
| void | computeBCRhsTerm_MB (FloatArray &answer, TimeStep *tStep) override |
| void | computeBCRhsTerm_MC (FloatArray &answer, TimeStep *tStep) override |
| void | computeLoadVector (FloatArray &answer, BodyLoad *load, CharType type, ValueModeType mode, TimeStep *tStep) override |
| double | computeCriticalTimeStep (TimeStep *tStep) override=0 |
| Computes the critical time increment. | |
| void | updateInternalState (TimeStep *tStep) override |
| int | checkConsistency () override |
| int | giveInternalStateAtNode (FloatArray &answer, InternalStateType type, InternalStateMode mode, int node, TimeStep *tStep) override |
| Public Member Functions inherited from oofem::SUPGElement | |
| SUPGElement (int n, Domain *aDomain) | |
| void | updateStabilizationCoeffs (TimeStep *tStep) override |
| virtual void | computeBCLhsPressureTerm_MC (FloatMatrix &answer, TimeStep *tStep) |
| virtual void | computeBCLhsTerm_MB (FloatMatrix &answer, TimeStep *tStep) |
| virtual void | computeBCLhsPressureTerm_MB (FloatMatrix &answer, TimeStep *tStep) |
| virtual void | computeSlipWithFrictionBCTerm_MB (FloatMatrix &answer, Load *load, int side, TimeStep *tStep) |
| virtual void | computePenetrationWithResistanceBCTerm_MB (FloatMatrix &answer, Load *load, int side, TimeStep *tStep) |
| virtual void | computeOutFlowBCTerm_MB (FloatMatrix &answer, int side, TimeStep *tStep) |
| virtual void | computeHomogenizedReinforceTerm_MB (FloatMatrix &answer, Load *load, TimeStep *tStep) |
| virtual void | computeHomogenizedReinforceTerm_MC (FloatMatrix &answer, Load *load, TimeStep *tStep) |
| virtual void | giveLocalVelocityDofMap (IntArray &map) |
| virtual void | giveLocalPressureDofMap (IntArray &map) |
| Public Member Functions inherited from oofem::FMElement | |
| FMElement (int n, Domain *aDomain) | |
| 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 int | giveNumberOfInternalDofManagers () const |
| virtual DofManager * | giveInternalDofManager (int i) const |
| virtual void | setInternalDofManager (int num, std::unique_ptr< DofManager > dm) |
| virtual double | giveCharacteristicValue (CharType type, TimeStep *tStep) |
| virtual void | computeBoundarySurfaceLoadVector (FloatArray &answer, BoundaryLoad *load, int boundary, CharType type, ValueModeType mode, TimeStep *tStep, bool global=true) |
| 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 | computeNumberOfDofs () |
| 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 | giveDofManDofIDMask (int inode, IntArray &answer) const |
| virtual void | giveInternalDofManDofIDMask (int inode, IntArray &answer) const |
| virtual void | giveElementDofIDMask (IntArray &answer) const |
| virtual void | computeField (ValueModeType mode, TimeStep *tStep, const FloatArray &lcoords, FloatArray &answer) |
| virtual double | computeVolumeAround (GaussPoint *gp) |
| 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 FEInterpolation * | giveInterpolation () const |
| virtual FEInterpolation * | giveInterpolation (DofIDItem id) 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 MaterialMode | giveMaterialMode () |
| virtual int | giveIntegrationRuleLocalCodeNumbers (IntArray &answer, IntegrationRule &ie) |
| int | giveRegionNumber () |
| virtual void | updateYourself (TimeStep *tStep) |
| virtual void | initializeYourself (TimeStep *timeStepWhenICApply) |
| virtual bool | isActivated (TimeStep *tStep) |
| virtual bool | isCast (TimeStep *tStep) |
| virtual void | initForNewStep () |
| virtual Element_Geometry_Type | giveGeometryType () const =0 |
| 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 | 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 | initializeFinish () override |
| void | postInitialize () override |
| Performs post initialization steps. | |
| void | saveContext (DataStream &stream, ContextMode mode) override |
| void | restoreContext (DataStream &stream, ContextMode mode) override |
| const char * | giveClassName () const 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. | |
| virtual const char * | giveInputRecordName () const =0 |
| 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. | |
| virtual Interface * | giveInterface (InterfaceType t) |
| std::string | errorInfo (const char *func) const |
| Returns string for prepending output (used by error reporting macros). | |
Protected Member Functions | |
| void | computeDeviatoricStrain (FloatArray &answer, GaussPoint *gp, TimeStep *tStep) override |
| virtual void | computeNuMatrix (FloatMatrix &answer, GaussPoint *gp)=0 |
| virtual void | computeUDotGradUMatrix (FloatMatrix &answer, GaussPoint *gp, TimeStep *tStep)=0 |
| virtual void | computeBMatrix (FloatMatrix &anwer, GaussPoint *gp)=0 |
| virtual void | computeDivUMatrix (FloatMatrix &answer, GaussPoint *gp)=0 |
| virtual void | computeNpMatrix (FloatMatrix &answer, GaussPoint *gp)=0 |
| virtual void | computeGradPMatrix (FloatMatrix &answer, GaussPoint *gp)=0 |
| virtual void | computeDivTauMatrix (FloatMatrix &answer, GaussPoint *gp, TimeStep *tStep)=0 |
| virtual void | computeGradUMatrix (FloatMatrix &answer, GaussPoint *gp, TimeStep *tStep)=0 |
| virtual int | giveNumberOfSpatialDimensions ()=0 |
| virtual void | computeEdgeLoadVector_MB (FloatArray &answer, Load *load, int id, TimeStep *tStep) |
| virtual void | computeSurfaceLoadVector_MB (FloatArray &answer, Load *load, int id, TimeStep *tStep) |
| virtual void | computeEdgeLoadVector_MC (FloatArray &answer, Load *load, int id, TimeStep *tStep) |
| virtual void | computeSurfaceLoadVector_MC (FloatArray &answer, Load *load, int id, TimeStep *tStep) |
| Protected Member Functions inherited from oofem::SUPGElement | |
| virtual void | computeDeviatoricStress (FloatArray &answer, const FloatArray &eps, GaussPoint *gp, TimeStep *tStep)=0 |
| virtual void | computeTangent (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)=0 |
| virtual void | computeGaussPoints () |
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 |
| Protected Attributes inherited from oofem::SUPGElement | |
| IntArray | boundarySides |
| Array of boundary sides. | |
| IntArray | boundaryCodes |
| Boundary sides codes. | |
| double | t_supg = 0. |
| double | t_pspg = 0. |
| double | t_lsic = 0. |
| 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 inherited from oofem::SUPGElement | |
| static ParamKey | IPK_SUPGElement_bsides |
| static ParamKey | IPK_SUPGElement_bcodes |
This abstract class represent a general base element class for fluid dynamic problems. It is derived from SUPGElement base, but it provides general algorithms how to obtain characteristic vectors & matrices by means of numerical integration. It is therefore suitable for high-order elements where exact integration would be difficult.
Definition at line 57 of file supgelement2.h.
| oofem::SUPGElement2::SUPGElement2 | ( | int | n, |
| Domain * | aDomain ) |
Definition at line 55 of file supgelement2.C.
References oofem::SUPGElement::SUPGElement().
Referenced by oofem::Quad10_2D_SUPG::Quad10_2D_SUPG(), oofem::Tet1_3D_SUPG::Tet1_3D_SUPG(), and oofem::TR21_2D_SUPG::TR21_2D_SUPG().
|
overridevirtual |
Performs consistency check. This method is called at startup for all elements in particular domain. This method is intended to check data compatibility. Particular element types should test if compatible material and crossSection both with required capabilities are specified. Derived classes should provide their own analysis specific tests. Some printed input if incompatibility is found should be provided (error or warning member functions). Method can be also used to initialize some variables, since this is invoked after all domain components are instanciated.
Reimplemented from oofem::SUPGElement.
Reimplemented in oofem::TR21_2D_SUPG.
Definition at line 189 of file supgelement2.C.
|
overridevirtual |
Computes acceleration terms (generalized mass matrix with stabilization terms) for momentum balance equations(s).
Implements oofem::SUPGElement.
Definition at line 253 of file supgelement2.C.
References oofem::FloatMatrix::clear(), computeNuMatrix(), computeUDotGradUMatrix(), oofem::Element::computeVolumeAround(), oofem::Element::giveCrossSection(), oofem::TimeStep::givePreviousStep(), oofem::Element::integrationRulesArray, oofem::FloatMatrix::plusProductUnsym(), and oofem::SUPGElement::t_supg.
|
overridevirtual |
Computes acceleration terms for mass conservation equation.
Implements oofem::SUPGElement.
Definition at line 552 of file supgelement2.C.
References oofem::FloatMatrix::clear(), computeGradPMatrix(), computeNuMatrix(), oofem::Element::computeVolumeAround(), oofem::Element::integrationRulesArray, oofem::FloatMatrix::plusProductUnsym(), and oofem::SUPGElement::t_pspg.
|
overridevirtual |
Computes the derivative of advection terms for momentum balance equations(s) with respect to nodal velocities.
Implements oofem::SUPGElement.
Definition at line 301 of file supgelement2.C.
References oofem::FloatMatrix::beProductOf(), oofem::FloatMatrix::clear(), computeGradUMatrix(), computeNuMatrix(), computeUDotGradUMatrix(), oofem::Element::computeVolumeAround(), oofem::Element::giveCrossSection(), oofem::TimeStep::givePreviousStep(), oofem::Element::integrationRulesArray, oofem::FloatMatrix::plusProductUnsym(), and oofem::SUPGElement::t_supg.
Referenced by giveCharacteristicMatrix().
|
overridevirtual |
Computes the derivative of advection terms for mass conservation equation with respect to nodal velocities.
Implements oofem::SUPGElement.
Definition at line 474 of file supgelement2.C.
References oofem::FloatMatrix::clear(), computeGradPMatrix(), computeUDotGradUMatrix(), oofem::Element::computeVolumeAround(), oofem::Element::integrationRulesArray, oofem::FloatMatrix::plusProductUnsym(), and oofem::SUPGElement::t_pspg.
Referenced by giveCharacteristicMatrix().
|
overridevirtual |
Computes nonlinear advection terms for momentum balance equations(s).
Implements oofem::SUPGElement.
Definition at line 274 of file supgelement2.C.
References oofem::FloatArray::beProductOf(), oofem::FloatArray::clear(), computeNuMatrix(), computeUDotGradUMatrix(), oofem::FMElement::computeVectorOfVelocities(), oofem::Element::computeVolumeAround(), oofem::Element::giveCrossSection(), oofem::TimeStep::givePreviousStep(), oofem::Element::integrationRulesArray, oofem::FloatArray::plusProduct(), and oofem::SUPGElement::t_supg.
Referenced by giveCharacteristicVector().
|
overridevirtual |
Computes advection terms for mass conservation equation.
Implements oofem::SUPGElement.
Definition at line 451 of file supgelement2.C.
References oofem::FloatArray::beProductOf(), oofem::FloatArray::clear(), computeGradPMatrix(), computeUDotGradUMatrix(), oofem::FMElement::computeVectorOfVelocities(), oofem::Element::computeVolumeAround(), oofem::Element::integrationRulesArray, oofem::FloatArray::plusProduct(), and oofem::SUPGElement::t_pspg.
Referenced by giveCharacteristicVector().
|
overridevirtual |
Computes Rhs terms due to boundary conditions.
Implements oofem::SUPGElement.
Definition at line 588 of file supgelement2.C.
References oofem::FloatArray::add(), oofem::Element::bodyLoadArray, oofem::BodyLoadBGT, oofem::Element::boundaryLoadArray, oofem::FloatArray::clear(), oofem::Load::computeComponentArrayAt(), computeEdgeLoadVector_MB(), computeNuMatrix(), computeSurfaceLoadVector_MB(), computeUDotGradUMatrix(), oofem::Element::computeVolumeAround(), oofem::FEMComponent::domain, oofem::EdgeLoadBGT, oofem::ForceLoadBVT, oofem::GeneralBoundaryCondition::giveBCGeoType(), oofem::GeneralBoundaryCondition::giveBCValType(), oofem::Element::giveBodyLoadArray(), oofem::Element::giveBoundaryLoadArray(), oofem::Element::giveCrossSection(), oofem::TimeStep::givePreviousStep(), oofem::FloatArray::giveSize(), oofem::Element::integrationRulesArray, OOFEM_ERROR, oofem::FloatArray::plusProduct(), oofem::SurfaceLoadBGT, and oofem::SUPGElement::t_supg.
Referenced by giveCharacteristicVector().
|
overridevirtual |
Computes Rhs terms due to boundary conditions.
Implements oofem::SUPGElement.
Definition at line 647 of file supgelement2.C.
References oofem::FloatArray::add(), oofem::Element::bodyLoadArray, oofem::BodyLoadBGT, oofem::Element::boundaryLoadArray, oofem::FloatArray::clear(), oofem::Load::computeComponentArrayAt(), computeEdgeLoadVector_MC(), computeGradPMatrix(), computeSurfaceLoadVector_MC(), oofem::Element::computeVolumeAround(), oofem::FEMComponent::domain, oofem::EdgeLoadBGT, oofem::ForceLoadBVT, oofem::GeneralBoundaryCondition::giveBCGeoType(), oofem::GeneralBoundaryCondition::giveBCValType(), oofem::Element::giveBodyLoadArray(), oofem::Element::giveBoundaryLoadArray(), oofem::FloatArray::giveSize(), oofem::Element::integrationRulesArray, OOFEM_ERROR, oofem::FloatArray::plusProduct(), oofem::SurfaceLoadBGT, and oofem::SUPGElement::t_pspg.
Referenced by giveCharacteristicVector().
|
protectedpure virtual |
Implemented in oofem::Quad10_2D_SUPG, oofem::Tet1_3D_SUPG, and oofem::TR21_2D_SUPG.
Referenced by computeDeviatoricStrain(), computeDiffusionDerivativeTerm_MB(), and computeDiffusionTerm_MB().
|
overridepure virtual |
Computes the critical time increment.
Implements oofem::SUPGElement.
Implemented in oofem::Quad10_2D_SUPG, oofem::Tet1_3D_SUPG, and oofem::TR21_2D_SUPG.
|
overrideprotectedvirtual |
Implements oofem::SUPGElement.
Definition at line 773 of file supgelement2.C.
References oofem::FloatArray::beProductOf(), computeBMatrix(), and oofem::FMElement::computeVectorOfVelocities().
Referenced by updateInternalState().
|
overridevirtual |
Computes the derivative of diffusion terms for momentum balance equations(s) with respect to nodal velocities.
Implements oofem::SUPGElement.
Definition at line 358 of file supgelement2.C.
References oofem::FloatMatrix::beProductOf(), oofem::FloatMatrix::clear(), computeBMatrix(), computeDivTauMatrix(), oofem::SUPGElement::computeTangent(), computeUDotGradUMatrix(), oofem::Element::computeVolumeAround(), oofem::FEMComponent::domain, oofem::TimeStep::givePreviousStep(), oofem::FluidModel::giveReynoldsNumber(), oofem::Element::integrationRulesArray, oofem::FloatMatrix::plusProductUnsym(), and oofem::SUPGElement::t_supg.
Referenced by giveCharacteristicMatrix().
|
overridevirtual |
Computes diffusion derivative terms for mass conservation equation.
Implements oofem::SUPGElement.
Definition at line 491 of file supgelement2.C.
References oofem::FloatMatrix::clear(), computeDivTauMatrix(), computeGradPMatrix(), oofem::Element::computeVolumeAround(), oofem::Element::giveCrossSection(), oofem::Element::integrationRulesArray, oofem::FloatMatrix::plusProductUnsym(), oofem::SUPGElement::t_pspg, and oofem::FloatMatrix::zero().
Referenced by giveCharacteristicMatrix().
|
overridevirtual |
Computes diffusion terms for momentum balance equations(s).
Implements oofem::SUPGElement.
Definition at line 330 of file supgelement2.C.
References oofem::FloatArray::beProductOf(), oofem::FloatArray::clear(), computeBMatrix(), oofem::SUPGElement::computeDeviatoricStress(), computeDivTauMatrix(), computeUDotGradUMatrix(), oofem::FMElement::computeVectorOfVelocities(), oofem::Element::computeVolumeAround(), oofem::FEMComponent::domain, oofem::TimeStep::givePreviousStep(), oofem::FluidModel::giveReynoldsNumber(), oofem::Element::integrationRulesArray, oofem::FloatArray::plusProduct(), and oofem::SUPGElement::t_supg.
Referenced by giveCharacteristicVector().
|
overridevirtual |
Computes diffusion terms for mass conservation equation.
Implements oofem::SUPGElement.
Definition at line 512 of file supgelement2.C.
References oofem::FloatArray::clear().
Referenced by giveCharacteristicVector().
|
protectedpure virtual |
Implemented in oofem::Quad10_2D_SUPG, oofem::Tet1_3D_SUPG, and oofem::TR21_2D_SUPG.
Referenced by computeDiffusionDerivativeTerm_MB(), computeDiffusionDerivativeTerm_MC(), and computeDiffusionTerm_MB().
|
protectedpure virtual |
Implemented in oofem::Quad10_2D_SUPG, oofem::Tet1_3D_SUPG, and oofem::TR21_2D_SUPG.
Referenced by computeLinearAdvectionTerm_MC(), computeLSICStabilizationTerm_MB(), and computePressureTerm_MB().
|
protectedvirtual |
Definition at line 748 of file supgelement2.C.
References OOFEM_ERROR.
Referenced by computeBCRhsTerm_MB().
|
protectedvirtual |
Definition at line 760 of file supgelement2.C.
References OOFEM_ERROR.
Referenced by computeBCRhsTerm_MC().
|
protectedpure virtual |
Implemented in oofem::Quad10_2D_SUPG, oofem::Tet1_3D_SUPG, and oofem::TR21_2D_SUPG.
Referenced by computeAccelerationTerm_MC(), computeAdvectionDerivativeTerm_MC(), computeAdvectionTerm_MC(), computeBCRhsTerm_MC(), computeDiffusionDerivativeTerm_MC(), computeLoadVector(), computePressureTerm_MB(), and computePressureTerm_MC().
|
protectedpure virtual |
Implemented in oofem::Quad10_2D_SUPG, oofem::Tet1_3D_SUPG, and oofem::TR21_2D_SUPG.
Referenced by computeAdvectionDerivativeTerm_MB().
|
overridevirtual |
Computes the linear advection term for mass conservation equation.
Implements oofem::SUPGElement.
Definition at line 433 of file supgelement2.C.
References oofem::FloatMatrix::clear(), computeDivUMatrix(), computeNpMatrix(), oofem::Element::computeVolumeAround(), oofem::Element::integrationRulesArray, and oofem::FloatMatrix::plusProductUnsym().
Referenced by giveCharacteristicMatrix(), and oofem::Quad10_2D_SUPG::updateStabilizationCoeffs().
|
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 702 of file supgelement2.C.
References oofem::FloatArray::assemble(), oofem::FloatArray::clear(), oofem::Load::computeComponentArrayAt(), computeGradPMatrix(), computeNuMatrix(), oofem::Element::computeNumberOfDofs(), computeUDotGradUMatrix(), oofem::FMElement::computeVectorOfVelocities(), oofem::Element::computeVolumeAround(), oofem::ForceLoadBVT, oofem::GeneralBoundaryCondition::giveBCValType(), oofem::Element::giveCrossSection(), oofem::SUPGElement::giveLocalPressureDofMap(), oofem::SUPGElement::giveLocalVelocityDofMap(), oofem::TimeStep::givePreviousStep(), oofem::Element::integrationRulesArray, oofem::FloatArray::plusProduct(), oofem::FloatArray::resize(), oofem::SUPGElement::t_pspg, oofem::SUPGElement::t_supg, and oofem::FloatArray::zero().
|
overridevirtual |
Computes SLIC stabilization term for momentum balance equation(s).
Implements oofem::SUPGElement.
Definition at line 414 of file supgelement2.C.
References oofem::FloatMatrix::clear(), computeDivUMatrix(), oofem::Element::computeVolumeAround(), oofem::Element::giveCrossSection(), oofem::Element::integrationRulesArray, oofem::FloatMatrix::plusProductSymmUpper(), oofem::FloatMatrix::symmetrized(), and oofem::SUPGElement::t_lsic.
Referenced by giveCharacteristicMatrix().
|
protectedpure virtual |
Implemented in oofem::Quad10_2D_SUPG, oofem::Tet1_3D_SUPG, and oofem::TR21_2D_SUPG.
Referenced by computeLinearAdvectionTerm_MC(), and computePressureTerm_MB().
|
protectedpure virtual |
Implemented in oofem::Quad10_2D_SUPG, oofem::Tet1_3D_SUPG, and oofem::TR21_2D_SUPG.
Referenced by computeAccelerationTerm_MB(), computeAccelerationTerm_MC(), computeAdvectionDerivativeTerm_MB(), computeAdvectionTerm_MB(), computeBCRhsTerm_MB(), and computeLoadVector().
|
overridevirtual |
Computes pressure terms for momentum balance equations(s).
Implements oofem::SUPGElement.
Definition at line 384 of file supgelement2.C.
References oofem::FloatMatrix::clear(), computeDivUMatrix(), computeGradPMatrix(), computeNpMatrix(), computeUDotGradUMatrix(), oofem::Element::computeVolumeAround(), oofem::TimeStep::givePreviousStep(), oofem::Element::integrationRulesArray, oofem::FloatMatrix::plusProductUnsym(), and oofem::SUPGElement::t_supg.
Referenced by giveCharacteristicMatrix().
|
overridevirtual |
Computes pressure terms for mass conservation equation.
Implements oofem::SUPGElement.
Definition at line 569 of file supgelement2.C.
References oofem::FloatMatrix::clear(), computeGradPMatrix(), oofem::Element::computeVolumeAround(), oofem::Element::giveCrossSection(), oofem::Element::integrationRulesArray, oofem::FloatMatrix::plusProductSymmUpper(), oofem::FloatMatrix::symmetrized(), and oofem::SUPGElement::t_pspg.
Referenced by giveCharacteristicMatrix().
|
protectedvirtual |
Definition at line 754 of file supgelement2.C.
References OOFEM_ERROR.
Referenced by computeBCRhsTerm_MB().
|
protectedvirtual |
Definition at line 766 of file supgelement2.C.
References OOFEM_ERROR.
Referenced by computeBCRhsTerm_MC().
|
protectedpure virtual |
Implemented in oofem::Quad10_2D_SUPG, oofem::Tet1_3D_SUPG, and oofem::TR21_2D_SUPG.
Referenced by computeAccelerationTerm_MB(), computeAdvectionDerivativeTerm_MB(), computeAdvectionDerivativeTerm_MC(), computeAdvectionTerm_MB(), computeAdvectionTerm_MC(), computeBCRhsTerm_MB(), computeDiffusionDerivativeTerm_MB(), computeDiffusionTerm_MB(), computeLoadVector(), and computePressureTerm_MB().
|
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::SUPGElement.
Definition at line 75 of file supgelement2.C.
References oofem::FloatMatrix::assemble(), computeAdvectionDerivativeTerm_MB(), computeAdvectionDerivativeTerm_MC(), computeDiffusionDerivativeTerm_MB(), computeDiffusionDerivativeTerm_MC(), computeLinearAdvectionTerm_MC(), computeLSICStabilizationTerm_MB(), oofem::Element::computeNumberOfDofs(), computePressureTerm_MB(), computePressureTerm_MC(), oofem::SUPGElement::giveLocalPressureDofMap(), oofem::SUPGElement::giveLocalVelocityDofMap(), OOFEM_ERROR, oofem::FloatMatrix::resize(), and oofem::FloatMatrix::zero().
Referenced by giveCharacteristicVector().
|
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::SUPGElement.
Definition at line 116 of file supgelement2.C.
References oofem::FloatArray::assemble(), oofem::FloatArray::beProductOf(), computeAdvectionTerm_MB(), computeAdvectionTerm_MC(), computeBCRhsTerm_MB(), computeBCRhsTerm_MC(), computeDiffusionTerm_MB(), computeDiffusionTerm_MC(), oofem::Element::computeNumberOfDofs(), oofem::FMElement::computeVectorOfPressures(), oofem::FMElement::computeVectorOfVelocities(), giveCharacteristicMatrix(), oofem::SUPGElement::giveLocalPressureDofMap(), oofem::SUPGElement::giveLocalVelocityDofMap(), OOFEM_ERROR, oofem::FloatArray::resize(), and oofem::FloatArray::zero().
|
overridevirtual |
Setups the input record string of receiver.
| input | Dynamic input record to be filled by receiver. |
Reimplemented from oofem::SUPGElement.
Definition at line 68 of file supgelement2.C.
|
overridevirtual |
Returns internal state variable (like stress,strain) at node of element in Reduced form, the way how is obtained is dependent on InternalValueType. The value may be local, or smoothed using some recovery technique. Returns zero if element is unable to respond to request.
| answer | Contains result, zero sized if not supported. |
| type | Determines the internal variable requested (physical meaning). |
| mode | Determines the mode of variable (recovered, local, ...). |
| node | Node number, for which variable is required. |
| tStep | Time step. |
Reimplemented from oofem::SUPGElement.
Reimplemented in oofem::TR21_2D_SUPG.
Definition at line 218 of file supgelement2.C.
References oofem::FloatArray::at(), oofem::DofManager::end(), oofem::DofManager::findDofWithDofId(), oofem::Element::giveNode(), oofem::Element::giveSpatialDimension(), and oofem::FloatArray::resize().
|
protectedpure virtual |
Implemented in oofem::Quad10_2D_SUPG, oofem::Tet1_3D_SUPG, and oofem::TR21_2D_SUPG.
|
overridevirtual |
Reimplemented from oofem::SUPGElement.
Definition at line 61 of file supgelement2.C.
|
inlineoverridevirtual |
Reimplemented from oofem::SUPGElement.
Definition at line 68 of file supgelement2.h.
|
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::SUPGElement.
Definition at line 202 of file supgelement2.C.
References computeDeviatoricStrain(), oofem::SUPGElement::computeDeviatoricStress(), and oofem::Element::integrationRulesArray.