|
OOFEM 3.0
|
#include <tr1_2d_cbs.h>
Public Member Functions | |
| TR1_2D_CBS (int n, Domain *aDomain) | |
| FEInterpolation * | giveInterpolation () const override |
| void | computeConsistentMassMtrx (FloatMatrix &answer, TimeStep *tStep) override |
| void | computeDiagonalMassMtrx (FloatArray &answer, TimeStep *tStep) override |
| void | computeConvectionTermsI (FloatArray &answer, TimeStep *tStep) override |
| void | computeDiffusionTermsI (FloatArray &answer, TimeStep *tStep) override |
| void | computeDensityRhsVelocityTerms (FloatArray &answer, TimeStep *tStep) override |
| Computes velocity terms on RHS for density equation. | |
| void | computeDensityRhsPressureTerms (FloatArray &answer, TimeStep *tStep) override |
| Computes pressure terms on RHS for density equation. | |
| void | computePrescribedTractionPressure (FloatArray &answer, TimeStep *tStep) override |
| Computes prescribed pressure due to applied tractions. | |
| void | computeNumberOfNodalPrescribedTractionPressureContributions (FloatArray &answer, TimeStep *tStep) override |
| Computes number of edges/sides with prescribed traction contributing to node with prescribed pressure. | |
| void | computePressureLhs (FloatMatrix &answer, TimeStep *tStep) override |
| Calculates the pressure LHS. | |
| void | computeCorrectionRhs (FloatArray &answer, TimeStep *tStep) override |
| Calculates the RHS of velocity correction step. | |
| double | computeCriticalTimeStep (TimeStep *tStep) override |
| Calculates critical time step. | |
| 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 | initializeFrom (InputRecord &ir, int priority) override |
| void | giveInputRecord (DynamicInputRecord &input) 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 |
| Interface * | giveInterface (InterfaceType) override |
| int | EIPrimaryFieldI_evaluateFieldVectorAt (FloatArray &answer, PrimaryField &pf, const FloatArray &coords, IntArray &dofId, ValueModeType mode, TimeStep *tStep) override |
| double | computeLEPLICVolumeFraction (const FloatArray &n, const double p, LEPlic *matInterface, bool updFlag) override |
| Computes corresponding volume fraction to given interface position. | |
| void | formMaterialVolumePoly (Polygon &matvolpoly, LEPlic *matInterface, const FloatArray &normal, const double p, bool updFlag) override |
| Assembles the true element material polygon (takes receiver vof into accout). | |
| void | formVolumeInterfacePoly (Polygon &matvolpoly, LEPlic *matInterface, const FloatArray &normal, const double p, bool updFlag) override |
| Assembles receiver material polygon based solely on given interface line. | |
| double | truncateMatVolume (const Polygon &matvolpoly, double &volume) override |
| Truncates given material polygon to receiver. | |
| void | giveElementCenter (LEPlic *mat_interface, FloatArray ¢er, bool upd) override |
| Computes the receiver center (in updated Lagrangian configuration). | |
| void | formMyVolumePoly (Polygon &myPoly, LEPlic *mat_interface, bool updFlag) override |
| Assembles receiver volume. | |
| Element * | giveElement () override |
| Return number of receiver's element. | |
| double | computeMyVolume (LEPlic *matInterface, bool updFlag) override |
| Computes the volume of receiver. | |
| double | computeCriticalLEPlicTimeStep (TimeStep *tStep) override |
| Computes critical time step. | |
| void | NodalAveragingRecoveryMI_computeNodalValue (FloatArray &answer, int node, InternalStateType type, TimeStep *tStep) override |
| void | SPRNodalRecoveryMI_giveSPRAssemblyPoints (IntArray &pap) override |
| void | SPRNodalRecoveryMI_giveDofMansDeterminedByPatch (IntArray &answer, int pap) override |
| int | SPRNodalRecoveryMI_giveNumberOfIP () override |
| SPRPatchType | SPRNodalRecoveryMI_givePatchType () override |
| int | giveIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override |
| void | printOutputAt (FILE *file, TimeStep *tStep) override |
| int | giveInternalStateAtNode (FloatArray &answer, InternalStateType type, InternalStateMode mode, int node, TimeStep *tStep) override |
| void | drawRawGeometry (oofegGraphicContext &gc, TimeStep *tStep) override |
| void | drawScalar (oofegGraphicContext &gc, TimeStep *tStep) override |
| Public Member Functions inherited from oofem::CBSElement | |
| CBSElement (int n, Domain *aDomain) | |
| void | initializeFinish () override |
| void | giveCharacteristicMatrix (FloatMatrix &answer, CharType type, TimeStep *tStep) override |
| void | giveCharacteristicVector (FloatArray &answer, CharType type, ValueModeType mode, TimeStep *tStep) override |
| virtual void | computePrescribedTermsI (FloatArray &answer, TimeStep *tStep) |
| void | updateInternalState (TimeStep *tStep) override |
| Public Member Functions inherited from oofem::FMElement | |
| FMElement (int n, Domain *aDomain) | |
| virtual void | updateStabilizationCoeffs (TimeStep *tStep) |
| void | computeVectorOfVelocities (ValueModeType mode, TimeStep *tStep, FloatArray &velocities) |
| void | computeVectorOfPressures (ValueModeType mode, TimeStep *tStep, FloatArray &pressures) |
| FloatArray | computeVectorOfVelocities (ValueModeType mode, TimeStep *tStep) |
| FloatArray | computeVectorOfPressures (ValueModeType mode, TimeStep *tStep) |
| Public Member Functions inherited from oofem::Element | |
| Element (int n, Domain *aDomain) | |
| Element (const Element &src)=delete | |
| Element & | operator= (const Element &src)=delete |
| virtual | ~Element () |
| Virtual destructor. | |
| void | giveLocationArray (IntArray &locationArray, const UnknownNumberingScheme &s, IntArray *dofIds=NULL) const |
| void | giveLocationArray (IntArray &locationArray, const IntArray &dofIDMask, const UnknownNumberingScheme &s, IntArray *dofIds=NULL) const |
| virtual void | giveBoundaryLocationArray (IntArray &locationArray, const IntArray &bNodes, const UnknownNumberingScheme &s, IntArray *dofIds=NULL) |
| virtual void | giveBoundaryLocationArray (IntArray &locationArray, const IntArray &bNodes, const IntArray &dofIDMask, const UnknownNumberingScheme &s, IntArray *dofIds=NULL) |
| virtual int | giveNumberOfDofs () |
| virtual 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 | computeLoadVector (FloatArray &answer, BodyLoad *load, CharType type, ValueModeType mode, 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 | 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 (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 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. | |
| 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::SpatialLocalizerInterface | |
| SpatialLocalizerInterface (Element *element) | |
| virtual int | SpatialLocalizerI_containsPoint (const FloatArray &coords) |
| int | SpatialLocalizerI_BBoxContainsPoint (const FloatArray &coords) |
| virtual void | SpatialLocalizerI_giveBBox (FloatArray &bb0, FloatArray &bb1) |
| virtual double | SpatialLocalizerI_giveClosestPoint (FloatArray &lcoords, FloatArray &closest, const FloatArray &gcoords) |
| Public Member Functions inherited from oofem::Interface | |
| Interface () | |
| Constructor. | |
| virtual | ~Interface () |
| Public Member Functions inherited from oofem::EIPrimaryFieldInterface | |
| EIPrimaryFieldInterface () | |
| 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. | |
| Public Member Functions inherited from oofem::SPRNodalRecoveryModelInterface | |
| SPRNodalRecoveryModelInterface () | |
| Constructor. | |
| Public Member Functions inherited from oofem::LEPlicElementInterface | |
| LEPlicElementInterface () | |
| bool | isBoundary () |
| Returns true if cell is boundary. | |
| void | setTempLineConstant (double tp) |
| void | setTempInterfaceNormal (FloatArray tg) |
| void | setTempVolumeFraction (double v) |
| void | setPermanentVolumeFraction (double v) |
| void | addTempVolumeFraction (double v) |
| double | giveVolumeFraction () |
| double | giveTempVolumeFraction () |
| void | giveTempInterfaceNormal (FloatArray &n) |
| double | giveTempLineConstant () |
| void | updateYourself (TimeStep *tStep) |
| void | saveContext (DataStream &stream, ContextMode mode) |
| void | restoreContext (DataStream &stream, ContextMode mode) |
Protected Member Functions | |
| void | computeGaussPoints () override |
| void | computeDeviatoricStress (FloatArray &answer, GaussPoint *gp, TimeStep *tStep) override |
| Computes deviatoric stress vector in given integration point and solution step from given total strain vector. | |
Protected Attributes | |
| double | b [3] |
| double | c [3] |
| double | area = 0. |
| Protected Attributes inherited from oofem::CBSElement | |
| IntArray | boundarySides |
| Array of boundary sides. | |
| IntArray | boundaryCodes |
| Boundary sides codes. | |
| 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. | |
| Protected Attributes inherited from oofem::LEPlicElementInterface | |
| bool | permanentVofFlag |
| double | vof |
| Volume fraction of reference fluid in element. | |
| double | temp_vof |
| double | p |
| Line constant of line segment representing interface. | |
| double | temp_p |
| FloatArray | normal |
| Interface segment normal. | |
| FloatArray | temp_normal |
Static Protected Attributes | |
| static FEI2dTrLin | interp |
| static ParamKey | IPK_TR1_2D_CBS_vof |
| static ParamKey | IPK_TR1_2D_CBS_pvof |
| Static Protected Attributes inherited from oofem::CBSElement | |
| static ParamKey | IPK_CBSElement_bsides |
| static ParamKey | IPK_CBSElement_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 |
This class is the implementation of triangular CFD element with linear (and equal order) interpolation of velocity and pressure fields. Should be used with CBS solution algorithm.
Definition at line 67 of file tr1_2d_cbs.h.
| oofem::TR1_2D_CBS::TR1_2D_CBS | ( | int | n, |
| Domain * | aDomain ) |
Definition at line 75 of file tr1_2d_cbs.C.
References oofem::CBSElement::CBSElement(), oofem::LEPlicElementInterface::LEPlicElementInterface(), oofem::Element::numberOfDofMans, oofem::SpatialLocalizerInterface::SpatialLocalizerInterface(), and oofem::ZZNodalRecoveryModelInterface::ZZNodalRecoveryModelInterface().
|
overridevirtual |
Used to check consistency and initialize some element geometry data (area,b,c).
Reimplemented from oofem::CBSElement.
Definition at line 676 of file tr1_2d_cbs.C.
References area, b, c, and oofem::Element::giveNode().
|
overridevirtual |
Calculates consistent mass matrix.
Implements oofem::CBSElement.
Definition at line 148 of file tr1_2d_cbs.C.
References area, oofem::FloatMatrix::at(), oofem::Element::giveCrossSection(), oofem::Element::integrationRulesArray, oofem::FloatMatrix::resize(), and oofem::FloatMatrix::zero().
|
overridevirtual |
Calculates convection component for (*) velocities.
Implements oofem::CBSElement.
Definition at line 185 of file tr1_2d_cbs.C.
References area, oofem::FloatArray::at(), b, oofem::Element::bodyLoadArray, oofem::BodyLoadBGT, c, oofem::FMElement::computeVectorOfVelocities(), oofem::FEMComponent::domain, oofem::ForceLoadBVT, oofem::Element::giveBodyLoadArray(), oofem::Element::giveCrossSection(), oofem::TimeStep::givePreviousStep(), oofem::FloatArray::giveSize(), oofem::TimeStep::giveTimeIncrement(), oofem::Element::integrationRulesArray, oofem::FloatArray::resize(), oofem::FloatArray::times(), and oofem::FloatArray::zero().
|
overridevirtual |
Calculates the RHS of velocity correction step.
Implements oofem::CBSElement.
Definition at line 589 of file tr1_2d_cbs.C.
References area, oofem::FloatArray::at(), b, c, oofem::FMElement::computeVectorOfPressures(), oofem::FMElement::computeVectorOfVelocities(), oofem::TimeStep::givePreviousStep(), oofem::TimeStep::giveTimeIncrement(), oofem::FloatArray::resize(), and oofem::FloatArray::zero().
|
inlineoverridevirtual |
Computes critical time step.
Implements oofem::LEPlicElementInterface.
Definition at line 135 of file tr1_2d_cbs.h.
|
overridevirtual |
Calculates critical time step.
Implements oofem::CBSElement.
Definition at line 706 of file tr1_2d_cbs.C.
References oofem::FloatArray::at(), b, c, oofem::FMElement::computeVectorOfVelocities(), oofem::FEMComponent::domain, oofem::FluidModel::giveReynoldsNumber(), oofem::max(), and oofem::min().
|
overridevirtual |
Computes pressure terms on RHS for density equation.
Implements oofem::CBSElement.
Definition at line 549 of file tr1_2d_cbs.C.
References area, oofem::FloatArray::at(), b, c, oofem::FMElement::computeVectorOfPressures(), oofem::FEMComponent::domain, oofem::TimeStep::givePreviousStep(), oofem::CBS::giveTheta1(), oofem::TimeStep::giveTimeIncrement(), oofem::FloatArray::resize(), and oofem::FloatArray::zero().
|
overridevirtual |
Computes velocity terms on RHS for density equation.
Implements oofem::CBSElement.
Definition at line 373 of file tr1_2d_cbs.C.
References oofem::FloatArray::add(), area, oofem::FloatArray::at(), b, oofem::CBSElement::boundaryCodes, oofem::CBSElement::boundarySides, c, oofem::Element::computeVectorOfPrescribed(), oofem::FMElement::computeVectorOfVelocities(), oofem::FEMComponent::domain, FMElement_PrescribedPressureBC, FMElement_PrescribedTractionBC, FMElement_PrescribedUnBC, oofem::Element::giveCrossSection(), oofem::Element::giveNode(), oofem::CBS::giveTheta1(), oofem::Element::integrationRulesArray, oofem::FloatArray::resize(), oofem::FloatArray::times(), and oofem::FloatArray::zero().
|
overrideprotectedvirtual |
Computes deviatoric stress vector in given integration point and solution step from given total strain vector.
Implements oofem::CBSElement.
Definition at line 661 of file tr1_2d_cbs.C.
References oofem::FloatArray::at(), b, c, oofem::FMElement::computeVectorOfVelocities(), and oofem::Element::giveCrossSection().
|
overridevirtual |
Calculates diagonal mass matrix.
Implements oofem::CBSElement.
Definition at line 172 of file tr1_2d_cbs.C.
References area, oofem::FloatArray::at(), oofem::Element::giveCrossSection(), oofem::Element::integrationRulesArray, oofem::FloatArray::resize(), and oofem::FloatArray::zero().
|
overridevirtual |
Calculates contribution from diffusion terms for (*) velocities.
Implements oofem::CBSElement.
Definition at line 270 of file tr1_2d_cbs.C.
References area, oofem::FloatArray::at(), b, oofem::Element::bodyLoadArray, oofem::BodyLoadBGT, oofem::CBSElement::boundaryCodes, oofem::Element::boundaryLoadArray, oofem::CBSElement::boundarySides, c, oofem::FEMComponent::domain, FMElement_PrescribedTractionBC, FMElement_PrescribedUnBC, FMElement_PrescribedUsBC, oofem::ForceLoadBVT, oofem::Element::giveBodyLoadArray(), oofem::Element::giveBoundaryLoadArray(), oofem::Element::giveCrossSection(), oofem::FluidDynamicMaterialStatus::giveDeviatoricStressVector(), oofem::GaussPoint::giveMaterialStatus(), oofem::Element::giveNode(), oofem::FluidModel::giveReynoldsNumber(), oofem::FloatArray::giveSize(), oofem::Element::integrationRulesArray, oofem::FloatArray::resize(), oofem::FloatArray::times(), and oofem::FloatArray::zero().
|
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 137 of file tr1_2d_cbs.C.
References oofem::Element::giveCrossSection(), and oofem::Element::integrationRulesArray.
|
overridevirtual |
Computes corresponding volume fraction to given interface position.
Implements oofem::LEPlicElementInterface.
Definition at line 739 of file tr1_2d_cbs.C.
References computeMyVolume(), oofem::Polygon::computeVolume(), formVolumeInterfacePoly(), and oofem::LEPlicElementInterface::p.
|
overridevirtual |
Computes the volume of receiver.
Implements oofem::LEPlicElementInterface.
Definition at line 934 of file tr1_2d_cbs.C.
References area, oofem::Element::giveNode(), oofem::FEMComponent::giveNumber(), oofem::LEPlic::giveUpdatedXCoordinate(), and oofem::LEPlic::giveUpdatedYCoordinate().
Referenced by computeLEPLICVolumeFraction().
|
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 91 of file tr1_2d_cbs.C.
|
overridevirtual |
Computes number of edges/sides with prescribed traction contributing to node with prescribed pressure.
Implements oofem::CBSElement.
Definition at line 523 of file tr1_2d_cbs.C.
References oofem::FloatArray::at(), oofem::CBSElement::boundaryCodes, oofem::CBSElement::boundarySides, FMElement_PrescribedTractionBC, oofem::FloatArray::resize(), and oofem::FloatArray::zero().
|
overridevirtual |
Computes prescribed pressure due to applied tractions.
Implements oofem::CBSElement.
Definition at line 453 of file tr1_2d_cbs.C.
References oofem::FloatArray::at(), oofem::CBSElement::boundaryCodes, oofem::Element::boundaryLoadArray, oofem::CBSElement::boundarySides, oofem::BoundaryLoad::computeValueAt(), oofem::FEMComponent::domain, FMElement_PrescribedTractionBC, oofem::Element::giveBoundaryLoadArray(), oofem::FluidDynamicMaterialStatus::giveDeviatoricStressVector(), oofem::GaussPoint::giveMaterialStatus(), oofem::Element::giveNode(), oofem::FluidModel::giveReynoldsNumber(), oofem::Element::integrationRulesArray, oofem::FloatArray::resize(), oofem::FloatArray::times(), and oofem::FloatArray::zero().
|
overridevirtual |
Calculates the pressure LHS.
Implements oofem::CBSElement.
Definition at line 572 of file tr1_2d_cbs.C.
References area, oofem::FloatMatrix::at(), b, c, and oofem::FloatMatrix::resize().
|
overridevirtual |
Reimplemented from oofem::Element.
Definition at line 1148 of file tr1_2d_cbs.C.
References gc, oofem::Element::giveNode(), OOFEG_RAW_GEOMETRY_LAYER, OOFEG_RAW_GEOMETRY_WIDTH, and oofem::LEPlicElementInterface::p.
|
overridevirtual |
Reimplemented from oofem::Element.
Definition at line 1178 of file tr1_2d_cbs.C.
References oofem::FloatArray::at(), gc, giveInternalStateAtNode(), giveIPValue(), oofem::Element::giveNode(), oofem::Element::integrationRulesArray, oofem::ISM_local, oofem::ISM_recovered, OOFEG_DEFORMED_GEOMETRY_WIDTH, OOFEG_VARPLOT_PATTERN_LAYER, oofem::LEPlicElementInterface::p, oofem::SA_COLORZPROFILE, oofem::SA_ISO_SURF, and oofem::SA_ZPROFILE.
|
overridevirtual |
Evaluates the value of field at given point of interest (should be located inside receiver's volume) using element interpolation.
This should use local coordinates instead of having all elements search for it manually.
Shouldn't this just be part of Element? It's very much the core of functionality for elements.
| answer | Field evaluated at coordinate. |
| pf | Field to use for evaluation. |
| coords | Coordinate. |
| dofId | IDs of DOFs to evaluate. |
| mode | Mode of field. |
| tStep | Time step to evaluate at. |
Implements oofem::EIPrimaryFieldInterface.
Definition at line 977 of file tr1_2d_cbs.C.
References oofem::FloatArray::at(), oofem::IntArray::at(), oofem::Element::computeLocalCoordinates(), oofem::Element::computeVectorOf(), oofem::IntArray::findFirstIndexOf(), oofem::Element::giveElementDofIDMask(), oofem::IntArray::giveSize(), OOFEM_ERROR, oofem::FloatArray::resize(), oofem::sum(), and oofem::FloatArray::zero().
|
overridevirtual |
Assembles the true element material polygon (takes receiver vof into accout).
Implements oofem::LEPlicElementInterface.
Definition at line 753 of file tr1_2d_cbs.C.
References oofem::Polygon::addVertex(), oofem::Polygon::clear(), formVolumeInterfacePoly(), oofem::Element::giveNode(), oofem::FEMComponent::giveNumber(), oofem::LEPlic::giveUpdatedXCoordinate(), oofem::LEPlic::giveUpdatedYCoordinate(), oofem::LEPlicElementInterface::normal, oofem::LEPlicElementInterface::p, oofem::Vertex::setCoords(), TRSUPG_ZERO_VOF, and oofem::LEPlicElementInterface::vof.
|
overridevirtual |
Assembles receiver volume.
Implements oofem::LEPlicElementInterface.
Definition at line 911 of file tr1_2d_cbs.C.
References oofem::Polygon::addVertex(), oofem::Polygon::clear(), oofem::Element::giveNode(), oofem::FEMComponent::giveNumber(), oofem::LEPlic::giveUpdatedXCoordinate(), oofem::LEPlic::giveUpdatedYCoordinate(), and oofem::Vertex::setCoords().
Referenced by truncateMatVolume().
|
overridevirtual |
Assembles receiver material polygon based solely on given interface line.
Implements oofem::LEPlicElementInterface.
Definition at line 785 of file tr1_2d_cbs.C.
References oofem::Polygon::addVertex(), oofem::Polygon::clear(), oofem::Element::giveNode(), oofem::FEMComponent::giveNumber(), oofem::LEPlic::giveUpdatedXCoordinate(), oofem::LEPlic::giveUpdatedYCoordinate(), oofem::LEPlicElementInterface::normal, oofem::LEPlicElementInterface::p, and oofem::Vertex::setCoords().
Referenced by computeLEPLICVolumeFraction(), and formMaterialVolumePoly().
|
inlineoverridevirtual |
Reimplemented from oofem::Element.
Definition at line 101 of file tr1_2d_cbs.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 97 of file tr1_2d_cbs.C.
|
inlineoverridevirtual |
Return number of receiver's element.
Implements oofem::LEPlicElementInterface.
Definition at line 133 of file tr1_2d_cbs.h.
|
overridevirtual |
Computes the receiver center (in updated Lagrangian configuration).
Implements oofem::LEPlicElementInterface.
Definition at line 952 of file tr1_2d_cbs.C.
References oofem::FloatArray::add(), oofem::FloatArray::at(), oofem::Element::giveNode(), oofem::FEMComponent::giveNumber(), oofem::LEPlic::giveUpdatedCoordinate(), oofem::FloatArray::resize(), oofem::FloatArray::times(), and oofem::FloatArray::zero().
|
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 104 of file tr1_2d_cbs.h.
|
overridevirtual |
Setups the input record string of receiver.
| input | Dynamic input record to be filled by receiver. |
Reimplemented from oofem::CBSElement.
Definition at line 125 of file tr1_2d_cbs.C.
References IPK_TR1_2D_CBS_pvof, IPK_TR1_2D_CBS_vof, oofem::LEPlicElementInterface::permanentVofFlag, oofem::DynamicInputRecord::setField(), and oofem::LEPlicElementInterface::vof.
|
inlineoverridevirtual |
Implements oofem::FEMComponent.
Definition at line 102 of file tr1_2d_cbs.h.
References _IFT_TR1_2D_CBS_Name.
|
overridevirtual |
Interface requesting service.
Reimplemented from oofem::FEMComponent.
Definition at line 637 of file tr1_2d_cbs.C.
References oofem::EIPrimaryFieldInterface::EIPrimaryFieldInterface(), oofem::EIPrimaryFieldInterfaceType, oofem::LEPlicElementInterface::LEPlicElementInterface(), oofem::LEPlicElementInterfaceType, oofem::NodalAveragingRecoveryModelInterface::NodalAveragingRecoveryModelInterface(), oofem::NodalAveragingRecoveryModelInterfaceType, oofem::SpatialLocalizerInterface::SpatialLocalizerInterface(), oofem::SpatialLocalizerInterfaceType, oofem::SPRNodalRecoveryModelInterface::SPRNodalRecoveryModelInterface(), oofem::SPRNodalRecoveryModelInterfaceType, 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::CBSElement.
Definition at line 1128 of file tr1_2d_cbs.C.
References oofem::FloatArray::at(), oofem::Element::giveCrossSection(), oofem::LEPlicElementInterface::giveTempVolumeFraction(), oofem::Element::integrationRulesArray, and oofem::FloatArray::resize().
Referenced by drawScalar().
|
overridevirtual |
Reimplemented from oofem::Element.
Definition at line 88 of file tr1_2d_cbs.C.
References interp.
|
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 1041 of file tr1_2d_cbs.C.
References oofem::FloatArray::at(), oofem::LEPlicElementInterface::giveTempVolumeFraction(), and oofem::FloatArray::resize().
Referenced by drawScalar(), and NodalAveragingRecoveryMI_computeNodalValue().
|
inlineoverridevirtual |
Returns material mode for receiver integration points. Should be specialized.
Reimplemented from oofem::Element.
Definition at line 103 of file tr1_2d_cbs.h.
|
overridevirtual |
Reimplemented from oofem::CBSElement.
Definition at line 104 of file tr1_2d_cbs.C.
References oofem::FEMComponent::giveDomain(), IPK_TR1_2D_CBS_pvof, IPK_TR1_2D_CBS_vof, oofem::FEMComponent::number, PM_UPDATE_PARAMETER_AND_REPORT, oofem::LEPlicElementInterface::setPermanentVolumeFraction(), oofem::LEPlicElementInterface::temp_vof, and oofem::LEPlicElementInterface::vof.
|
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 1053 of file tr1_2d_cbs.C.
References giveIPValue(), and oofem::Element::integrationRulesArray.
|
overridevirtual |
Prints output of receiver to stream, for given time step. This is used for output into the standard output file.
| file | File pointer to print to. |
| tStep | Time step to write for. |
Reimplemented from oofem::Element.
Definition at line 1096 of file tr1_2d_cbs.C.
References oofem::Element::giveCrossSection(), oofem::LEPlicElementInterface::giveVolumeFraction(), 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 1117 of file tr1_2d_cbs.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 1108 of file tr1_2d_cbs.C.
|
overridevirtual |
Implements oofem::SPRNodalRecoveryModelInterface.
Definition at line 1070 of file tr1_2d_cbs.C.
References oofem::IntArray::at(), oofem::Element::giveNode(), oofem::FEMComponent::giveNumber(), OOFEM_ERROR, and oofem::IntArray::resize().
|
overridevirtual |
Implements oofem::SPRNodalRecoveryModelInterface.
Definition at line 1083 of file tr1_2d_cbs.C.
|
overridevirtual |
Implements oofem::SPRNodalRecoveryModelInterface.
Definition at line 1088 of file tr1_2d_cbs.C.
References oofem::SPRPatchType_2dxy.
|
overridevirtual |
Implements oofem::SPRNodalRecoveryModelInterface.
Definition at line 1061 of file tr1_2d_cbs.C.
References oofem::IntArray::at(), oofem::Element::giveNode(), and oofem::IntArray::resize().
|
overridevirtual |
Truncates given material polygon to receiver.
Implements oofem::LEPlicElementInterface.
Definition at line 894 of file tr1_2d_cbs.C.
References area, oofem::Graph::clip(), oofem::Polygon::computeVolume(), oofem::Polygon::draw(), formMyVolumePoly(), gc, and OOFEG_DEBUG_LAYER.
|
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 1032 of file tr1_2d_cbs.C.
|
protected |
Definition at line 78 of file tr1_2d_cbs.h.
Referenced by checkConsistency(), computeConsistentMassMtrx(), computeConvectionTermsI(), computeCorrectionRhs(), computeDensityRhsPressureTerms(), computeDensityRhsVelocityTerms(), computeDiagonalMassMtrx(), computeDiffusionTermsI(), computeMyVolume(), computePressureLhs(), and truncateMatVolume().
|
protected |
Definition at line 76 of file tr1_2d_cbs.h.
Referenced by checkConsistency(), computeConvectionTermsI(), computeCorrectionRhs(), computeCriticalTimeStep(), computeDensityRhsPressureTerms(), computeDensityRhsVelocityTerms(), computeDeviatoricStress(), computeDiffusionTermsI(), and computePressureLhs().
|
protected |
Definition at line 77 of file tr1_2d_cbs.h.
Referenced by checkConsistency(), computeConvectionTermsI(), computeCorrectionRhs(), computeCriticalTimeStep(), computeDensityRhsPressureTerms(), computeDensityRhsVelocityTerms(), computeDeviatoricStress(), computeDiffusionTermsI(), and computePressureLhs().
|
staticprotected |
Definition at line 74 of file tr1_2d_cbs.h.
Referenced by giveInterpolation().
|
staticprotected |
Definition at line 81 of file tr1_2d_cbs.h.
Referenced by giveInputRecord(), and initializeFrom().
|
staticprotected |
Definition at line 80 of file tr1_2d_cbs.h.
Referenced by giveInputRecord(), and initializeFrom().