|
OOFEM 3.0
|
#include <tr21_2d_supg.h>
Public Member Functions | |
| TR21_2D_SUPG (int n, Domain *aDomain) | |
| FEInterpolation * | giveInterpolation () const override |
| FEInterpolation * | giveInterpolation (DofIDItem id) const override |
| const char * | giveClassName () const override |
| const char * | giveInputRecordName () const override |
| MaterialMode | giveMaterialMode () override |
| Element_Geometry_Type | giveGeometryType () const override |
| void | giveDofManDofIDMask (int inode, IntArray &answer) const override |
| int | computeNumberOfDofs () override |
| void | updateYourself (TimeStep *tStep) override |
| int | checkConsistency () override |
| Used to check consistency and initialize some element geometry data (area,b,c). | |
| void | saveContext (DataStream &stream, ContextMode mode) override |
| void | restoreContext (DataStream &stream, ContextMode mode) override |
| double | LS_PCS_computeF (LevelSetPCS *ls, TimeStep *tStep) override |
| void | LS_PCS_computedN (FloatMatrix &answer) override |
| double | LS_PCS_computeVolume () override |
| Returns receiver's volume. | |
| void | LS_PCS_computeVolume (double &answer, const FloatArray **coordinates) |
| double | LS_PCS_computeS (LevelSetPCS *ls, TimeStep *tStep) override |
| void | LS_PCS_computeVOFFractions (FloatArray &answer, FloatArray &fi) override |
| void | NodalAveragingRecoveryMI_computeNodalValue (FloatArray &answer, int node, InternalStateType type, TimeStep *tStep) override |
| Public Member Functions inherited from oofem::SUPGElement2 | |
| 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 |
| void | updateInternalState (TimeStep *tStep) override |
| Public Member Functions inherited from oofem::SUPGElement | |
| SUPGElement (int n, Domain *aDomain) | |
| 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) |
| 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 | 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 | 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 | 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 | initializeYourself (TimeStep *timeStepWhenICApply) |
| 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) |
| 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 | drawDeformedGeometry (oofegGraphicContext &gc, TimeStep *tStep, UnknownType) |
| 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 | 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::LevelSetPCSElementInterface | |
| LevelSetPCSElementInterface () | |
| Public Member Functions inherited from oofem::Interface | |
| Interface () | |
| Constructor. | |
| virtual | ~Interface () |
| 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::NodalAveragingRecoveryModelInterface | |
| NodalAveragingRecoveryModelInterface () | |
| Constructor. | |
Static Protected Attributes | |
| static FEI2dTrQuad | velocityInterpolation |
| static FEI2dTrLin | pressureInterpolation |
| Static Protected Attributes inherited from oofem::SUPGElement | |
| static ParamKey | IPK_SUPGElement_bsides |
| static ParamKey | IPK_SUPGElement_bcodes |
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 Member Functions inherited from oofem::SUPGElement2 | |
| void | computeDeviatoricStrain (FloatArray &answer, GaussPoint *gp, TimeStep *tStep) override |
| 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 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. | |
Class representing 2d triangular element with quadratic velocity and linear pressure approximation for solving incompressible fluid problems with SUPG solver.
Definition at line 57 of file tr21_2d_supg.h.
| oofem::TR21_2D_SUPG::TR21_2D_SUPG | ( | int | n, |
| Domain * | aDomain ) |
Definition at line 66 of file tr21_2d_supg.C.
References oofem::Element::numberOfDofMans, oofem::SUPGElement2::SUPGElement2(), and oofem::ZZNodalRecoveryModelInterface::ZZNodalRecoveryModelInterface().
|
overridevirtual |
Used to check consistency and initialize some element geometry data (area,b,c).
Reimplemented from oofem::SUPGElement2.
Definition at line 1332 of file tr21_2d_supg.C.
| void oofem::TR21_2D_SUPG::computeAdvectionDeltaTerm | ( | FloatMatrix & | answer, |
| TimeStep * | tStep ) |
Definition at line 331 of file tr21_2d_supg.C.
References oofem::FloatMatrix::clear(), computeNuMatrix(), computeUDotGradUMatrix(), computeVolumeAround(), oofem::Element::giveCrossSection(), oofem::Element::integrationRulesArray, and oofem::FloatMatrix::plusProductUnsym().
Referenced by updateStabilizationCoeffs().
| void oofem::TR21_2D_SUPG::computeAdvectionTerm | ( | FloatMatrix & | answer, |
| TimeStep * | tStep ) |
Definition at line 313 of file tr21_2d_supg.C.
References oofem::FloatMatrix::clear(), computeNuMatrix(), computeUDotGradUMatrix(), computeVolumeAround(), oofem::Element::giveCrossSection(), oofem::Element::integrationRulesArray, and oofem::FloatMatrix::plusProductUnsym().
Referenced by updateStabilizationCoeffs().
|
overrideprotectedvirtual |
Implements oofem::SUPGElement2.
Definition at line 172 of file tr21_2d_supg.C.
References oofem::FloatMatrix::at(), oofem::GaussPoint::giveNaturalCoordinates(), oofem::FloatMatrix::resize(), velocityInterpolation, and oofem::FloatMatrix::zero().
| void oofem::TR21_2D_SUPG::computeCenterOf | ( | FloatArray & | C, |
| FloatArray | c, | ||
| int | dim ) |
Definition at line 1157 of file tr21_2d_supg.C.
References oofem::FloatArray::at().
Referenced by computeMiddlePointOnParabolicArc(), and LS_PCS_computeVOFFractions().
| void oofem::TR21_2D_SUPG::computeCoordsOfEdge | ( | FloatArray & | answer, |
| int | iedge ) |
Definition at line 1193 of file tr21_2d_supg.C.
References oofem::FloatArray::at(), oofem::Element::giveNode(), and velocityInterpolation.
Referenced by computeQuadraticFunct().
|
overridevirtual |
Computes the critical time increment.
Implements oofem::SUPGElement2.
Definition at line 395 of file tr21_2d_supg.C.
|
overrideprotectedvirtual |
Implements oofem::SUPGElement.
Definition at line 129 of file tr21_2d_supg.C.
References oofem::Element::giveCrossSection().
|
overrideprotectedvirtual |
Implements oofem::SUPGElement2.
Definition at line 246 of file tr21_2d_supg.C.
References oofem::FloatMatrix::at(), oofem::Element::giveCrossSection(), oofem::GaussPoint::giveNaturalCoordinates(), oofem::Element::integrationRulesArray, oofem::FloatMatrix::resize(), velocityInterpolation, and oofem::FloatMatrix::zero().
|
overrideprotectedvirtual |
Implements oofem::SUPGElement2.
Definition at line 189 of file tr21_2d_supg.C.
References oofem::FloatMatrix::at(), oofem::GaussPoint::giveNaturalCoordinates(), oofem::FloatMatrix::resize(), velocityInterpolation, and oofem::FloatMatrix::zero().
Referenced by computeLSICTerm().
|
overrideprotectedvirtual |
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 105 of file tr21_2d_supg.C.
References oofem::Element::giveCrossSection(), and oofem::Element::integrationRulesArray.
|
overrideprotectedvirtual |
Implements oofem::SUPGElement2.
Definition at line 235 of file tr21_2d_supg.C.
References oofem::FloatMatrix::beTranspositionOf(), oofem::GaussPoint::giveNaturalCoordinates(), and pressureInterpolation.
|
overrideprotectedvirtual |
Implements oofem::SUPGElement2.
Definition at line 219 of file tr21_2d_supg.C.
References oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::FloatMatrix::beProductOf(), oofem::FMElement::computeVectorOfVelocities(), oofem::GaussPoint::giveNaturalCoordinates(), and velocityInterpolation.
| void oofem::TR21_2D_SUPG::computeIntersection | ( | int | iedge, |
| FloatArray & | intcoords, | ||
| FloatArray & | fi ) |
Definition at line 1062 of file tr21_2d_supg.C.
References oofem::FloatArray::at(), computeQuadraticRoots(), oofem::FloatArray::resize(), velocityInterpolation, and oofem::FloatArray::zero().
Referenced by LS_PCS_computeVOFFractions().
| void oofem::TR21_2D_SUPG::computeLSICTerm | ( | FloatMatrix & | answer, |
| TimeStep * | tStep ) |
Definition at line 369 of file tr21_2d_supg.C.
References oofem::FloatMatrix::clear(), computeDivUMatrix(), computeVolumeAround(), oofem::Element::giveCrossSection(), oofem::Element::integrationRulesArray, oofem::FloatMatrix::plusProductSymmUpper(), and oofem::FloatMatrix::symmetrized().
Referenced by updateStabilizationCoeffs().
| void oofem::TR21_2D_SUPG::computeMassDeltaTerm | ( | FloatMatrix & | answer, |
| TimeStep * | tStep ) |
Definition at line 351 of file tr21_2d_supg.C.
References oofem::FloatMatrix::clear(), computeNuMatrix(), computeUDotGradUMatrix(), computeVolumeAround(), oofem::Element::giveCrossSection(), oofem::Element::integrationRulesArray, and oofem::FloatMatrix::plusProductUnsym().
Referenced by updateStabilizationCoeffs().
| void oofem::TR21_2D_SUPG::computeMiddlePointOnParabolicArc | ( | FloatArray & | answer, |
| int | iedge, | ||
| FloatArray | borderpoints ) |
Definition at line 1132 of file tr21_2d_supg.C.
References oofem::FloatArray::at(), computeCenterOf(), and computeQuadraticFunct().
Referenced by LS_PCS_computeVOFFractions().
|
overrideprotectedvirtual |
Implements oofem::SUPGElement2.
Definition at line 204 of file tr21_2d_supg.C.
References oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::GaussPoint::giveNaturalCoordinates(), pressureInterpolation, oofem::FloatMatrix::resize(), and oofem::FloatMatrix::zero().
|
overrideprotectedvirtual |
Implements oofem::SUPGElement2.
Definition at line 143 of file tr21_2d_supg.C.
References oofem::FloatMatrix::beNMatrixOf(), oofem::GaussPoint::giveNaturalCoordinates(), and velocityInterpolation.
Referenced by computeAdvectionDeltaTerm(), computeAdvectionTerm(), computeMassDeltaTerm(), computeUDotGradUMatrix(), and LS_PCS_computeF().
|
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 89 of file tr21_2d_supg.C.
| void oofem::TR21_2D_SUPG::computeQuadraticFunct | ( | FloatArray & | answer, |
| FloatArray | line ) |
Definition at line 1267 of file tr21_2d_supg.C.
References oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::FloatMatrix::giveDeterminant(), and oofem::FloatArray::resize().
| void oofem::TR21_2D_SUPG::computeQuadraticFunct | ( | FloatArray & | answer, |
| int | iedge ) |
Definition at line 1213 of file tr21_2d_supg.C.
References oofem::FloatArray::at(), oofem::FloatMatrix::at(), computeCoordsOfEdge(), oofem::FloatMatrix::giveDeterminant(), and oofem::FloatArray::resize().
Referenced by computeMiddlePointOnParabolicArc(), and LS_PCS_computeVOFFractions().
| void oofem::TR21_2D_SUPG::computeQuadraticRoots | ( | FloatArray | Coeff, |
| double & | r1, | ||
| double & | r2 ) |
Definition at line 1178 of file tr21_2d_supg.C.
References oofem::FloatArray::at().
Referenced by computeIntersection(), and LS_PCS_computeVOFFractions().
|
overrideprotectedvirtual |
Implements oofem::SUPGElement.
Definition at line 136 of file tr21_2d_supg.C.
References oofem::Element::giveCrossSection().
|
overrideprotectedvirtual |
Implements oofem::SUPGElement2.
Definition at line 153 of file tr21_2d_supg.C.
References oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::FloatArray::beProductOf(), computeNuMatrix(), oofem::FMElement::computeVectorOfVelocities(), oofem::GaussPoint::giveNaturalCoordinates(), oofem::FloatMatrix::resize(), velocityInterpolation, and oofem::FloatMatrix::zero().
Referenced by computeAdvectionDeltaTerm(), computeAdvectionTerm(), and computeMassDeltaTerm().
|
overrideprotectedvirtual |
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 1362 of file tr21_2d_supg.C.
References oofem::GaussPoint::giveNaturalCoordinates(), oofem::GaussPoint::giveWeight(), and velocityInterpolation.
Referenced by computeAdvectionDeltaTerm(), computeAdvectionTerm(), computeLSICTerm(), computeMassDeltaTerm(), LS_PCS_computeF(), LS_PCS_computeS(), and LS_PCS_computeVolume().
|
overridevirtual |
Reimplemented from oofem::Element.
Definition at line 1419 of file tr21_2d_supg.C.
References gc, oofem::Element::giveNode(), OOFEG_RAW_GEOMETRY_LAYER, and OOFEG_RAW_GEOMETRY_WIDTH.
|
overridevirtual |
Reimplemented from oofem::Element.
Definition at line 1449 of file tr21_2d_supg.C.
References oofem::FloatArray::at(), oofem::Polygon::draw(), gc, giveInternalStateAtNode(), giveIPValue(), oofem::Element::giveNode(), oofem::Element::integrationRulesArray, oofem::ISM_local, oofem::ISM_recovered, OOFEG_DEFORMED_GEOMETRY_WIDTH, OOFEG_VARPLOT_PATTERN_LAYER, oofem::SA_COLORZPROFILE, oofem::SA_ISO_SURF, and oofem::SA_ZPROFILE.
|
inlineoverridevirtual |
Reimplemented from oofem::Element.
Definition at line 70 of file tr21_2d_supg.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 95 of file tr21_2d_supg.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 73 of file tr21_2d_supg.h.
|
inlineoverridevirtual |
Implements oofem::FEMComponent.
Definition at line 71 of file tr21_2d_supg.h.
References _IFT_TR21_2D_SUPG_Name.
|
overridevirtual |
Interface requesting service.
Reimplemented from oofem::FEMComponent.
Definition at line 1378 of file tr21_2d_supg.C.
References oofem::LevelSetPCSElementInterface::LevelSetPCSElementInterface(), oofem::LevelSetPCSElementInterfaceType, oofem::NodalAveragingRecoveryModelInterface::NodalAveragingRecoveryModelInterface(), oofem::NodalAveragingRecoveryModelInterfaceType, oofem::ZZNodalRecoveryModelInterface::ZZNodalRecoveryModelInterface(), and oofem::ZZNodalRecoveryModelInterfaceType.
|
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::SUPGElement2.
Definition at line 1410 of file tr21_2d_supg.C.
Referenced by drawScalar().
|
overridevirtual |
Reimplemented from oofem::Element.
Definition at line 73 of file tr21_2d_supg.C.
References velocityInterpolation.
|
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 79 of file tr21_2d_supg.C.
References pressureInterpolation, and velocityInterpolation.
|
overridevirtual |
Returns the integration point corresponding value in full form.
| answer | Contain corresponding integration point value, zero sized if not available. |
| gp | Integration point to check. |
| type | Determines the type of internal variable. |
| tStep | Time step. |
Reimplemented from oofem::Element.
Definition at line 1345 of file tr21_2d_supg.C.
Referenced by drawScalar(), and NodalAveragingRecoveryMI_computeNodalValue().
|
overrideprotectedvirtual |
Reimplemented from oofem::SUPGElement.
Definition at line 1399 of file tr21_2d_supg.C.
|
overrideprotectedvirtual |
Reimplemented from oofem::SUPGElement.
Definition at line 1393 of file tr21_2d_supg.C.
|
inlineoverridevirtual |
Returns material mode for receiver integration points. Should be specialized.
Reimplemented from oofem::Element.
Definition at line 72 of file tr21_2d_supg.h.
|
overrideprotectedvirtual |
Implements oofem::SUPGElement2.
Definition at line 388 of file tr21_2d_supg.C.
|
protected |
Definition at line 1327 of file tr21_2d_supg.C.
|
overridevirtual |
Returns gradient of shape functions.
Implements oofem::LevelSetPCSElementInterface.
Definition at line 429 of file tr21_2d_supg.C.
References oofem::FloatMatrix::add(), oofem::FloatMatrix::clear(), oofem::Element::integrationRulesArray, and velocityInterpolation.
|
overridevirtual |
Evaluates F in level set equation of the form
\[ \phi_t + F(\nabla\phi, x) |\nabla\phi| = 0 \]
where for interface position driven by flow with speed u:
\[ F = u\cdot \frac{\nabla\phi}{|\nabla\phi|} \]
Implements oofem::LevelSetPCSElementInterface.
Definition at line 402 of file tr21_2d_supg.C.
References oofem::FloatArray::beProductOf(), oofem::FloatArray::beTProductOf(), oofem::FloatArray::computeNorm(), computeNuMatrix(), oofem::FMElement::computeVectorOfVelocities(), computeVolumeAround(), oofem::Element::dofManArray, oofem::FloatArray::dotProduct(), oofem::LevelSetPCS::giveLevelSetDofManValue(), oofem::Element::integrationRulesArray, and velocityInterpolation.
|
overridevirtual |
Evaluates S in level set equation of the form
\[ \phi_t = S(\phi) (1-|\nabla\phi|) = 0 \]
where
\[ S=\frac{\phi}{\sqrt{\phi^2+\epsilon^2}} \]
Implements oofem::LevelSetPCSElementInterface.
Definition at line 477 of file tr21_2d_supg.C.
References computeVolumeAround(), oofem::Element::dofManArray, oofem::FloatArray::dotProduct(), oofem::LevelSetPCS::giveLevelSetDofManValue(), oofem::Element::integrationRulesArray, S, and velocityInterpolation.
|
overridevirtual |
Returns VOF fractions for each material on element according to nodal values of level set function (passed as parameter)
Implements oofem::LevelSetPCSElementInterface.
Definition at line 501 of file tr21_2d_supg.C.
References oofem::FloatArray::at(), computeCenterOf(), computeIntersection(), computeMiddlePointOnParabolicArc(), computeQuadraticFunct(), computeQuadraticRoots(), oofem::FloatArray::dotProduct(), oofem::Element::giveNode(), oofem::FEMComponent::giveNumber(), LS_PCS_computeVolume(), OOFEM_LOG_INFO, oofem::FloatArray::resize(), velocityInterpolation, oofem::FloatArray::zero(), and oofem::zero().
|
overridevirtual |
Returns receiver's volume.
Implements oofem::LevelSetPCSElementInterface.
Definition at line 465 of file tr21_2d_supg.C.
References computeVolumeAround(), and oofem::Element::integrationRulesArray.
Referenced by LS_PCS_computeVOFFractions().
| void oofem::TR21_2D_SUPG::LS_PCS_computeVolume | ( | double & | answer, |
| const FloatArray ** | coordinates ) |
Definition at line 444 of file tr21_2d_supg.C.
References oofem::Element::integrationRulesArray, and velocityInterpolation.
|
overridevirtual |
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 1318 of file tr21_2d_supg.C.
References giveIPValue(), and oofem::Element::integrationRulesArray.
|
overridevirtual |
Restores the receiver state previously written in stream.
| stream | Input stream. |
| mode | Determines amount of info available in stream (state, definition, ...). |
| throws | an ContextIOERR exception if error encountered. |
Reimplemented from oofem::Element.
Definition at line 1355 of file tr21_2d_supg.C.
|
overridevirtual |
Stores receiver state to output stream.
| stream | Output stream. |
| mode | Determines amount of info required in stream (state, definition, ...). |
| throws | an ContextIOERR exception if error encountered. |
Reimplemented from oofem::Element.
Definition at line 1350 of file tr21_2d_supg.C.
|
overrideprotectedvirtual |
Updates the stabilization coefficients used for CBS and SUPG algorithms.
| tStep | Active time step. |
Reimplemented from oofem::SUPGElement.
Definition at line 266 of file tr21_2d_supg.C.
References computeAdvectionDeltaTerm(), computeAdvectionTerm(), oofem::FloatMatrix::computeFrobeniusNorm(), computeLSICTerm(), computeMassDeltaTerm(), oofem::FMElement::computeVectorOfVelocities(), oofem::Element::giveCrossSection(), oofem::TimeStep::givePreviousStep(), oofem::Element::integrationRulesArray, N, oofem::SUPGElement::t_lsic, oofem::SUPGElement::t_pspg, and oofem::SUPGElement::t_supg.
|
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 1339 of file tr21_2d_supg.C.
|
staticprotected |
Definition at line 61 of file tr21_2d_supg.h.
Referenced by computeGradPMatrix(), computeNpMatrix(), and giveInterpolation().
|
staticprotected |
Definition at line 60 of file tr21_2d_supg.h.
Referenced by computeBMatrix(), computeCoordsOfEdge(), computeDivTauMatrix(), computeDivUMatrix(), computeGradUMatrix(), computeIntersection(), computeNuMatrix(), computeUDotGradUMatrix(), computeVolumeAround(), giveInterpolation(), giveInterpolation(), LS_PCS_computedN(), LS_PCS_computeF(), LS_PCS_computeS(), LS_PCS_computeVOFFractions(), and LS_PCS_computeVolume().