OOFEM
2.4
OOFEM.org - Object Oriented Finite Element Solver
|
This abstract class represent a general base element class for fluid dynamic problems solved using PFEM algorithm. More...
#include <pfemelement.h>
Public Member Functions | |
PFEMElement (int, Domain *) | |
Constructor. More... | |
~PFEMElement () | |
Destructor. More... | |
IRResultType | initializeFrom (InputRecord *ir) |
Initializes receiver acording to object description stored in input record. More... | |
void | giveCharacteristicMatrix (FloatMatrix &answer, CharType, TimeStep *) |
Computes characteristic matrix of receiver of requested type in given time step. More... | |
void | giveCharacteristicVector (FloatArray &answer, CharType, ValueModeType, TimeStep *) |
Computes characteristic vector of receiver of requested type in given time step. More... | |
virtual void | computeDiagonalMassMtrx (FloatArray &answer, TimeStep *)=0 |
Calculates diagonal mass matrix as vector. More... | |
virtual void | computeDiagonalMassMtrx (FloatMatrix &answer, TimeStep *)=0 |
Calculates diagonal mass matrix. More... | |
virtual double | computeCriticalTimeStep (TimeStep *tStep)=0 |
Calculates critical time step. More... | |
virtual void | updateInternalState (TimeStep *tStep) |
Updates element state after equilibrium in time step has been reached. More... | |
virtual void | printOutputAt (FILE *file, TimeStep *tStep) |
Prints output of receiver to stream, for given time step. More... | |
virtual int | checkConsistency () |
Performs consistency check. More... | |
virtual FEInterpolation * | giveVelocityInterpolation ()=0 |
Returns the interpolation for velocity. More... | |
virtual FEInterpolation * | givePressureInterpolation ()=0 |
Returns the interpolation for the pressure. More... | |
void | computeLoadVector (FloatArray &answer, BodyLoad *load, CharType type, ValueModeType mode, TimeStep *tStep) |
Computes the contribution of the given body load (volumetric). More... | |
virtual const char * | giveClassName () const |
virtual const IntArray & | giveVelocityDofMask () const =0 |
Returns mask of velocity Dofs. More... | |
virtual const IntArray & | givePressureDofMask () const =0 |
Returns mask of pressure Dofs. More... | |
virtual int | giveInternalStateAtNode (FloatArray &answer, InternalStateType type, InternalStateMode mode, int node, TimeStep *atTime) |
Returns internal state variable (like stress,strain) at node of element in Reduced form, the way how is obtained is dependent on InternalValueType. More... | |
virtual void | computeBodyLoadVectorAt (FloatArray &answer, BodyLoad *load, TimeStep *tStep, ValueModeType mode)=0 |
Computes the load vector due to body load acting on receiver, at given time step. More... | |
virtual void | computeDeviatoricStress (FloatArray &answer, GaussPoint *gp, TimeStep *)=0 |
Computes deviatoric stress vector in given integration point and solution step from given total strain vector. More... | |
virtual void | computeDeviatoricStressDivergence (FloatArray &answer, TimeStep *atTime)=0 |
Calculates the divergence of the deviatoric stress. More... | |
virtual void | computeBMatrix (FloatMatrix &answer, GaussPoint *gp)=0 |
Calculates the element shape function derivative matrix. More... | |
virtual void | computeStiffnessMatrix (FloatMatrix &answer, MatResponseMode mode, TimeStep *atTime)=0 |
Calculates the stiffness matrix. More... | |
virtual void | computePressureLaplacianMatrix (FloatMatrix &answer, TimeStep *atTime)=0 |
Calculates the pressure laplacian matrix. More... | |
virtual void | computeGradientMatrix (FloatMatrix &answer, TimeStep *atTime)=0 |
Calculates the pressure gradient matrix. More... | |
virtual void | computeDivergenceMatrix (FloatMatrix &answer, TimeStep *atTime)=0 |
Calculates the velocity divergence matrix. More... | |
virtual void | computePrescribedRhsVector (FloatArray &answer, TimeStep *tStep, ValueModeType mode)=0 |
Calculates the prescribed velocity vector for the right hand side of the pressure equation. More... | |
Public Member Functions inherited from oofem::FMElement | |
FMElement (int n, Domain *aDomain) | |
virtual | ~FMElement () |
virtual void | updateStabilizationCoeffs (TimeStep *tStep) |
Updates the stabilization coefficients used for CBS and SUPG algorithms. More... | |
virtual void | computeVectorOfVelocities (ValueModeType mode, TimeStep *tStep, FloatArray &velocities) |
virtual void | computeVectorOfPressures (ValueModeType mode, TimeStep *tStep, FloatArray &pressures) |
Public Member Functions inherited from oofem::Element | |
Element (int n, Domain *aDomain) | |
Constructor. More... | |
Element (const Element &src)=delete | |
Element & | operator= (const Element &src)=delete |
virtual | ~Element () |
Virtual destructor. More... | |
virtual void | drawYourself (oofegGraphicContext &gc, TimeStep *tStep) |
virtual void | drawAnnotation (oofegGraphicContext &gc, TimeStep *tStep) |
virtual void | drawRawGeometry (oofegGraphicContext &gc, TimeStep *tStep) |
virtual void | drawDeformedGeometry (oofegGraphicContext &gc, TimeStep *tStep, UnknownType) |
virtual void | drawScalar (oofegGraphicContext &gc, TimeStep *tStep) |
virtual void | drawSpecial (oofegGraphicContext &gc, TimeStep *tStep) |
virtual void | giveLocalIntVarMaxMin (oofegGraphicContext &gc, TimeStep *tStep, double &emin, double &emax) |
virtual int | giveInternalStateAtSide (FloatArray &answer, InternalStateType type, InternalStateMode mode, int side, TimeStep *tStep) |
Returns internal state variable (like stress,strain) at side of element in Reduced form If side is possessing DOFs, otherwise recover techniques will not work due to absence of side-shape functions. More... | |
virtual void | showSparseMtrxStructure (CharType mtrx, oofegGraphicContext &gc, TimeStep *tStep) |
Shows sparse structure. More... | |
virtual void | showExtendedSparseMtrxStructure (CharType mtrx, oofegGraphicContext &gc, TimeStep *tStep) |
Shows extended sparse structure (for example, due to nonlocal interactions for tangent stiffness) More... | |
int | giveLabel () const |
int | giveGlobalNumber () const |
void | setGlobalNumber (int num) |
Sets receiver globally unique number. More... | |
elementParallelMode | giveParallelMode () const |
Return elementParallelMode of receiver. More... | |
void | setParallelMode (elementParallelMode _mode) |
Sets parallel mode of element. More... | |
virtual elementParallelMode | giveKnotSpanParallelMode (int) const |
Returns the parallel mode for particular knot span of the receiver. More... | |
int | packUnknowns (DataStream &buff, TimeStep *tStep) |
Pack all necessary data of element (according to its parallel_mode) integration points into given communication buffer. More... | |
int | unpackAndUpdateUnknowns (DataStream &buff, TimeStep *tStep) |
Unpack and updates all necessary data of element (according to its parallel_mode) integration points into given communication buffer. More... | |
int | estimatePackSize (DataStream &buff) |
Estimates the necessary pack size to hold all packed data of receiver. More... | |
const IntArray * | givePartitionList () const |
Returns partition list of receiver. More... | |
void | setPartitionList (IntArray &pl) |
Sets partition list of receiver. More... | |
virtual double | predictRelativeComputationalCost () |
Returns the weight representing relative computational cost of receiver The reference element is triangular plane stress element with linear approximation, single integration point and linear isotropic material. More... | |
virtual double | giveRelativeSelfComputationalCost () |
Returns the weight representing relative computational cost of receiver The reference element is triangular plane stress element. More... | |
virtual double | predictRelativeRedistributionCost () |
Returns the relative redistribution cost of the receiver. More... | |
IntArray * | giveBodyLoadArray () |
Returns array containing load numbers of loads acting on element. More... | |
IntArray * | giveBoundaryLoadArray () |
Returns array containing load numbers of boundary loads acting on element. More... | |
virtual void | giveInputRecord (DynamicInputRecord &input) |
Setups the input record string of receiver. More... | |
virtual contextIOResultType | saveContext (DataStream &stream, ContextMode mode, void *obj=NULL) |
Stores receiver state to output stream. More... | |
virtual contextIOResultType | restoreContext (DataStream &stream, ContextMode mode, void *obj=NULL) |
Restores the receiver state previously written in stream. More... | |
void | giveLocationArray (IntArray &locationArray, const UnknownNumberingScheme &s, IntArray *dofIds=NULL) const |
Returns the location array (array of code numbers) of receiver for given numbering scheme. More... | |
void | giveLocationArray (IntArray &locationArray, const IntArray &dofIDMask, const UnknownNumberingScheme &s, IntArray *dofIds=NULL) const |
virtual void | giveBoundaryLocationArray (IntArray &locationArray, const IntArray &bNodes, const UnknownNumberingScheme &s, IntArray *dofIds=NULL) |
Returns the location array for the boundary of the element. More... | |
virtual void | giveBoundaryLocationArray (IntArray &locationArray, const IntArray &bNodes, const IntArray &dofIDMask, const UnknownNumberingScheme &s, IntArray *dofIds=NULL) |
virtual int | giveNumberOfDofs () |
virtual int | giveNumberOfInternalDofManagers () const |
virtual DofManager * | giveInternalDofManager (int i) const |
Returns i-th internal element dof manager of the receiver. More... | |
virtual double | giveCharacteristicValue (CharType type, TimeStep *tStep) |
Computes characteristic value of receiver of requested type in given time step. More... | |
virtual void | computeBoundarySurfaceLoadVector (FloatArray &answer, BoundaryLoad *load, int boundary, CharType type, ValueModeType mode, TimeStep *tStep, bool global=true) |
Computes the contribution of the given load at the given boundary surface in global coordinate system. More... | |
virtual void | computeTangentFromSurfaceLoad (FloatMatrix &answer, SurfaceLoad *load, int boundary, MatResponseMode rmode, TimeStep *tStep) |
Computes the tangent contribution of the given load at the given boundary. More... | |
virtual void | computeTangentFromEdgeLoad (FloatMatrix &answer, EdgeLoad *load, int boundary, MatResponseMode rmode, TimeStep *tStep) |
Computes the tangent contribution of the given load at the given boundary. More... | |
virtual void | computeBoundaryEdgeLoadVector (FloatArray &answer, BoundaryLoad *load, int edge, CharType type, ValueModeType mode, TimeStep *tStep, bool global=true) |
Computes the contribution of the given load at the given boundary edge. More... | |
const IntArray & | giveBodyLoadList () const |
Returns receiver list of bodyloads. More... | |
const IntArray & | giveBoundaryLoadList () const |
Returns receiver list of boundary loads. More... | |
void | computeVectorOf (ValueModeType u, TimeStep *tStep, FloatArray &answer) |
Returns local vector of unknowns. More... | |
void | computeVectorOf (const IntArray &dofIDMask, ValueModeType u, TimeStep *tStep, FloatArray &answer, bool padding=false) |
void | computeBoundaryVectorOf (const IntArray &bNodes, const IntArray &dofIDMask, ValueModeType u, TimeStep *tStep, FloatArray &answer, bool padding=false) |
Boundary version of computeVectorOf. More... | |
void | computeVectorOf (PrimaryField &field, const IntArray &dofIDMask, ValueModeType u, TimeStep *tStep, FloatArray &answer, bool padding=false) |
Returns local vector of unknowns. More... | |
void | computeVectorOfPrescribed (ValueModeType u, TimeStep *tStep, FloatArray &answer) |
Returns local vector of prescribed unknowns. More... | |
void | computeVectorOfPrescribed (const IntArray &dofIDMask, ValueModeType type, TimeStep *tStep, FloatArray &answer) |
Returns local vector of prescribed unknowns. More... | |
virtual int | computeNumberOfDofs () |
Computes or simply returns total number of element's local DOFs. More... | |
virtual int | computeNumberOfGlobalDofs () |
Computes the total number of element's global dofs. More... | |
int | computeNumberOfPrimaryMasterDofs () |
Computes the total number of element's primary master DOFs. More... | |
virtual bool | computeGtoLRotationMatrix (FloatMatrix &answer) |
Returns transformation matrix from global c.s. More... | |
virtual bool | giveRotationMatrix (FloatMatrix &answer) |
Transformation matrices updates rotation matrix between element-local and primary DOFs, taking into account nodal c.s. More... | |
virtual bool | computeDofTransformationMatrix (FloatMatrix &answer, const IntArray &nodes, bool includeInternal) |
Returns transformation matrix for DOFs from global coordinate system to local coordinate system in nodes. More... | |
virtual void | giveDofManDofIDMask (int inode, IntArray &answer) const |
Returns dofmanager dof mask for node. More... | |
virtual void | giveInternalDofManDofIDMask (int inode, IntArray &answer) const |
Returns internal dofmanager dof mask for node. More... | |
virtual void | giveElementDofIDMask (IntArray &answer) const |
Returns element dof mask for node. More... | |
virtual void | computeField (ValueModeType mode, TimeStep *tStep, const FloatArray &lcoords, FloatArray &answer) |
Computes the unknown vector interpolated at the specified local coordinates. More... | |
virtual double | computeVolumeAround (GaussPoint *gp) |
Returns volume related to given integration point. More... | |
virtual double | computeVolumeAreaOrLength () |
Computes the volume, area or length of the element depending on its spatial dimension. More... | |
double | computeMeanSize () |
Computes the size of the element defined as its length. More... | |
virtual double | computeVolume () |
Computes the volume. More... | |
virtual double | computeArea () |
Computes the area (zero for all but 2d geometries). More... | |
virtual double | computeLength () |
Computes the length (zero for all but 1D geometries) More... | |
virtual void | giveBoundaryEdgeNodes (IntArray &bNodes, int boundary) |
Returns list of receiver boundary nodes for given edge. More... | |
virtual void | giveBoundarySurfaceNodes (IntArray &bNodes, int boundary) |
Returns list of receiver boundary nodes for given surface. More... | |
virtual IntegrationRule * | giveBoundaryEdgeIntegrationRule (int order, int boundary) |
Returns boundary edge integration rule. More... | |
virtual IntegrationRule * | giveBoundarySurfaceIntegrationRule (int order, int boundary) |
Returns boundary surface integration rule. More... | |
int | giveDofManagerNumber (int i) const |
Translates local to global indices for dof managers. More... | |
const IntArray & | giveDofManArray () const |
void | addDofManager (DofManager *dMan) |
DofManager * | giveDofManager (int i) const |
Node * | giveNode (int i) const |
Returns reference to the i-th node of element. More... | |
virtual ElementSide * | giveSide (int i) const |
Returns reference to the i-th side of element. More... | |
virtual FEInterpolation * | giveInterpolation () const |
virtual FEInterpolation * | giveInterpolation (DofIDItem id) const |
Returns the interpolation for the specific dof id. More... | |
virtual Material * | giveMaterial () |
int | giveMaterialNumber () const |
CrossSection * | giveCrossSection () |
void | setMaterial (int matIndx) |
Sets the material of receiver. More... | |
virtual void | setCrossSection (int csIndx) |
Sets the cross section model of receiver. More... | |
virtual int | giveNumberOfDofManagers () const |
virtual int | giveNumberOfNodes () const |
Returns number of nodes of receiver. More... | |
void | setDofManagers (const IntArray &dmans) |
Sets receiver dofManagers. More... | |
void | setBodyLoads (const IntArray &bodyLoads) |
Sets receiver bodyLoadArray. More... | |
void | setIntegrationRules (std::vector< std::unique_ptr< IntegrationRule > > irlist) |
Sets integration rules. More... | |
virtual integrationDomain | giveIntegrationDomain () const |
Returns integration domain for receiver, used to initialize integration point over receiver volume. More... | |
virtual MaterialMode | giveMaterialMode () |
Returns material mode for receiver integration points. More... | |
virtual int | giveIntegrationRuleLocalCodeNumbers (IntArray &answer, IntegrationRule &ie) |
Assembles the code numbers of given integration element (sub-patch) This is done by obtaining list of nonzero shape functions and by collecting the code numbers of nodes corresponding to these shape functions. More... | |
int | giveRegionNumber () |
virtual void | postInitialize () |
Performs post initialization steps. More... | |
virtual void | updateYourself (TimeStep *tStep) |
Updates element state after equilibrium in time step has been reached. More... | |
virtual void | initializeYourself (TimeStep *timeStepWhenICApply) |
Initialization according to state given by initial conditions. More... | |
virtual bool | isActivated (TimeStep *tStep) |
virtual bool | isCast (TimeStep *tStep) |
virtual void | initForNewStep () |
Initializes receivers state to new time step. More... | |
virtual Element_Geometry_Type | giveGeometryType () const |
Returns the element geometry type. More... | |
virtual int | giveSpatialDimension () |
Returns the element spatial dimension (1, 2, or 3). More... | |
virtual int | giveNumberOfBoundarySides () |
virtual int | giveDefaultIntegrationRule () const |
Returns id of default integration rule. More... | |
virtual IntegrationRule * | giveDefaultIntegrationRulePtr () |
Access method for default integration rule. More... | |
int | giveNumberOfIntegrationRules () |
virtual IntegrationRule * | giveIntegrationRule (int i) |
std::vector< std::unique_ptr< IntegrationRule > > & | giveIntegrationRulesArray () |
virtual int | testElementExtension (ElementExtension ext) |
Tests if the element implements required extension. More... | |
virtual int | giveIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) |
Returns the integration point corresponding value in full form. More... | |
int | giveGlobalIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) |
virtual double | giveLengthInDir (const FloatArray &normalToCrackPlane) |
Default implementation returns length of element projection into specified direction. More... | |
virtual double | giveCharacteristicLength (const FloatArray &normalToCrackPlane) |
Returns the size of element in the given direction, in some cases adjusted (e.g. More... | |
double | giveCharacteristicLengthForPlaneElements (const FloatArray &normalToCrackPlane) |
Returns the size of element in the given direction if the direction is in the XY plane, otherwise gives the mean size defined as the square root of the element area. More... | |
double | giveCharacteristicLengthForAxisymmElements (const FloatArray &normalToCrackPlane) |
Returns the size of an axisymmetric element in the given direction if the direction is in the XY plane, otherwise gives the mean distance vrom the symmetry axis multiplied by pi. More... | |
virtual double | giveCharacteristicSize (GaussPoint *gp, FloatArray &normalToCrackPlane, ElementCharSizeMethod method) |
Returns characteristic element size for a given integration point and given direction. More... | |
virtual double | giveParentElSize () const |
Returns the size (length, area or volume depending on element type) of the parent element. More... | |
virtual void | updateBeforeNonlocalAverage (TimeStep *tStep) |
Updates internal element state (in all integration points of receiver) before nonlocal averaging takes place. More... | |
virtual int | computeGlobalCoordinates (FloatArray &answer, const FloatArray &lcoords) |
Computes the global coordinates from given element's local coordinates. More... | |
virtual bool | computeLocalCoordinates (FloatArray &answer, const FloatArray &gcoords) |
Computes the element local coordinates from given global coordinates. More... | |
virtual int | giveLocalCoordinateSystem (FloatMatrix &answer) |
Returns local coordinate system of receiver Required by material models with ortho- and anisotrophy. More... | |
virtual void | computeMidPlaneNormal (FloatArray &answer, const GaussPoint *gp) |
Computes mid-plane normal of receiver at integration point. More... | |
virtual int | adaptiveMap (Domain *oldd, TimeStep *tStep) |
Initializes the internal state variables stored in all IPs according to state in given domain. More... | |
virtual int | mapStateVariables (Domain &iOldDom, const TimeStep &iTStep) |
Maps the internal state variables stored in all IPs from the old domain to the new domain. More... | |
virtual int | adaptiveUpdate (TimeStep *tStep) |
Updates the internal state variables stored in all IPs according to already mapped state. More... | |
virtual int | adaptiveFinish (TimeStep *tStep) |
Finishes the mapping for given time step. More... | |
virtual void | updateLocalNumbering (EntityRenumberingFunctor &f) |
Local renumbering support. More... | |
template<class T > | |
void | ipEvaluator (T *src, void(T::*f)(GaussPoint *gp)) |
Integration point evaluator, loops over receiver IP's and calls given function (passed as f parameter) on them. The IP is parameter to function f. More... | |
template<class T , class S > | |
void | ipEvaluator (T *src, void(T::*f)(GaussPoint *, S &), S &_val) |
Integration point evaluator, loops over receiver IP's and calls given function (passed as f parameter) on them. The IP is parameter to function f as well as additional array. More... | |
Public Member Functions inherited from oofem::FEMComponent | |
FEMComponent (int n, Domain *d) | |
Regular constructor, creates component with given number and belonging to given domain. More... | |
virtual | ~FEMComponent () |
Virtual destructor. More... | |
virtual const char * | giveInputRecordName () const =0 |
Domain * | giveDomain () const |
virtual void | setDomain (Domain *d) |
Sets associated Domain. More... | |
int | giveNumber () const |
void | setNumber (int num) |
Sets number of receiver. More... | |
virtual void | printYourself () |
Prints receiver state on stdout. Useful for debugging. More... | |
virtual Interface * | giveInterface (InterfaceType t) |
Interface requesting service. More... | |
std::string | errorInfo (const char *func) const |
Returns string for prepending output (used by error reporting macros). More... | |
Additional Inherited Members | |
Protected Member Functions inherited from oofem::Element | |
virtual void | computeGaussPoints () |
Initializes the array of integration rules member variable. More... | |
Protected Attributes inherited from oofem::Element | |
int | numberOfDofMans |
Number of dofmanagers. More... | |
IntArray | dofManArray |
Array containing dofmanager numbers. More... | |
int | material |
Number of associated material. More... | |
int | crossSection |
Number of associated cross section. More... | |
IntArray | bodyLoadArray |
Array containing indexes of loads (body loads and boundary loads are kept separately), that apply on receiver. More... | |
IntArray | boundaryLoadArray |
std::vector< std::unique_ptr< IntegrationRule > > | integrationRulesArray |
List of integration rules of receiver (each integration rule contains associated integration points also). More... | |
FloatMatrix | elemLocalCS |
Transformation material matrix, used in orthotropic and anisotropic materials, global->local transformation. More... | |
int | activityTimeFunction |
Element activity time function. If defined, nonzero value indicates active receiver, zero value inactive element. More... | |
int | globalNumber |
In parallel mode, globalNumber contains globally unique DoFManager number. More... | |
int | numberOfGaussPoints |
Number of integration points as specified by nip. More... | |
elementParallelMode | parallel_mode |
Determines the parallel mode of the element. More... | |
IntArray | partitions |
List of partition sharing the shared element or remote partition containing remote element counterpart. More... | |
Protected Attributes inherited from oofem::FEMComponent | |
int | number |
Component number. More... | |
Domain * | domain |
Link to domain object, useful for communicating with other FEM components. More... | |
This abstract class represent a general base element class for fluid dynamic problems solved using PFEM algorithm.
Definition at line 66 of file pfemelement.h.
oofem::PFEMElement::PFEMElement | ( | int | n, |
Domain * | aDomain | ||
) |
Constructor.
Definition at line 58 of file pfemelement.C.
oofem::PFEMElement::~PFEMElement | ( | ) |
Destructor.
Definition at line 64 of file pfemelement.C.
|
virtual |
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::Element.
Reimplemented in oofem::TR1_2D_PFEM, and oofem::PFEMElement2d.
Definition at line 198 of file pfemelement.C.
Referenced by oofem::PFEMElement2d::checkConsistency().
|
pure virtual |
Calculates the element shape function derivative matrix.
Implemented in oofem::PFEMElement2d.
Referenced by giveClassName().
|
pure virtual |
Computes the load vector due to body load acting on receiver, at given time step.
Default implementation computes body load vector numerically as using default integration rule. Result is transformed to global c.s.
answer | Computed load vector due to body load |
load | Body load which contribution is computed. |
tStep | Time step. |
mode | determines the response mode |
Implemented in oofem::TR1_2D_PFEM.
Referenced by computeLoadVector(), and giveClassName().
|
pure virtual |
Calculates critical time step.
Implemented in oofem::TR1_2D_PFEM, and oofem::PFEMElement2d.
|
pure virtual |
Computes deviatoric stress vector in given integration point and solution step from given total strain vector.
Implemented in oofem::TR1_2D_PFEM, and oofem::PFEMElement2d.
Referenced by giveClassName(), and updateInternalState().
|
pure virtual |
Calculates the divergence of the deviatoric stress.
Implemented in oofem::TR1_2D_PFEM, and oofem::PFEMElement2d.
Referenced by giveClassName().
|
pure virtual |
Calculates diagonal mass matrix as vector.
Implemented in oofem::TR1_2D_PFEM.
Referenced by giveCharacteristicMatrix(), giveCharacteristicVector(), oofem::PFEMCorrectionRhsAssembler::vectorFromElement(), and oofem::PFEMMassVelocityAssembler::vectorFromElement().
|
pure virtual |
Calculates diagonal mass matrix.
Implemented in oofem::TR1_2D_PFEM.
|
pure virtual |
Calculates the velocity divergence matrix.
Implemented in oofem::PFEMElement2d.
Referenced by giveCharacteristicMatrix(), and giveClassName().
|
pure virtual |
Calculates the pressure gradient matrix.
Implemented in oofem::PFEMElement2d.
Referenced by giveCharacteristicMatrix(), giveClassName(), and oofem::PFEMCorrectionRhsAssembler::vectorFromElement().
|
virtual |
Computes the contribution of the given body load (volumetric).
answer | Requested contribution of load. |
load | Load to compute contribution from. |
type | Type of the contribution. |
mode | Determines mode of answer. |
tStep | Time step when answer is computed. |
Reimplemented from oofem::Element.
Definition at line 182 of file pfemelement.C.
References oofem::FloatArray::assemble(), oofem::FloatArray::clear(), computeBodyLoadVectorAt(), oofem::Element::computeNumberOfDofs(), giveVelocityDofMask(), oofem::FloatArray::resize(), and oofem::FloatArray::zero().
|
pure virtual |
Calculates the prescribed velocity vector for the right hand side of the pressure equation.
Implemented in oofem::PFEMElement2d.
Referenced by giveClassName(), and oofem::PFEMPressureRhsAssembler::vectorFromElement().
|
pure virtual |
Calculates the pressure laplacian matrix.
Implemented in oofem::PFEMElement2d.
Referenced by giveCharacteristicMatrix(), giveClassName(), and oofem::PFEMPressureLaplacianAssembler::matrixFromElement().
|
pure virtual |
Calculates the stiffness matrix.
Implemented in oofem::PFEMElement2d.
Referenced by giveCharacteristicMatrix(), giveClassName(), and oofem::PFEMLaplaceVelocityAssembler::vectorFromElement().
|
virtual |
Computes characteristic matrix of receiver of requested type in given time step.
answer | Requested characteristic matrix (stiffness, tangent, ...). If element has no capability to compute requested type of characteristic matrix error function is invoked. |
type | Id of characteristic component requested. |
tStep | Time step when answer is computed. |
Reimplemented from oofem::Element.
Definition at line 82 of file pfemelement.C.
References oofem::FloatMatrix::assemble(), computeDiagonalMassMtrx(), computeDivergenceMatrix(), computeGradientMatrix(), oofem::Element::computeNumberOfDofs(), computePressureLaplacianMatrix(), computeStiffnessMatrix(), givePressureDofMask(), giveVelocityDofMask(), OOFEM_ERROR, oofem::FloatMatrix::resize(), and oofem::FloatMatrix::zero().
Referenced by giveCharacteristicVector(), and oofem::PFEMPressureRhsAssembler::vectorFromElement().
|
virtual |
Computes characteristic vector of receiver of requested type in given time step.
If element has no capability to compute requested type of characteristic vector error function is invoked.
answer | Requested characteristic vector. |
type | Id of characteristic component requested. |
mode | Determines mode of answer. |
tStep | Time step when answer is computed. |
Reimplemented from oofem::Element.
Definition at line 118 of file pfemelement.C.
References oofem::FloatArray::assemble(), oofem::FloatArray::beProductOf(), computeDiagonalMassMtrx(), oofem::Element::computeNumberOfDofs(), oofem::Element::computeVectorOf(), oofem::FMElement::computeVectorOfPressures(), oofem::FMElement::computeVectorOfVelocities(), giveCharacteristicMatrix(), givePressureDofMask(), giveVelocityDofMask(), OOFEM_ERROR, oofem::FloatArray::resize(), and oofem::FloatArray::zero().
|
inlinevirtual |
Reimplemented from oofem::Element.
Reimplemented in oofem::TR1_2D_PFEM, and oofem::PFEMElement2d.
Definition at line 101 of file pfemelement.h.
References computeBMatrix(), computeBodyLoadVectorAt(), computeDeviatoricStress(), computeDeviatoricStressDivergence(), computeDivergenceMatrix(), computeGradientMatrix(), computePrescribedRhsVector(), computePressureLaplacianMatrix(), computeStiffnessMatrix(), giveInternalStateAtNode(), givePressureDofMask(), and giveVelocityDofMask().
|
virtual |
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::Element.
Reimplemented in oofem::TR1_2D_PFEM, and oofem::PFEMElement2d.
Definition at line 235 of file pfemelement.C.
References oofem::FloatArray::at(), oofem::DofManager::findDofWithDofId(), oofem::Element::giveInternalStateAtNode(), oofem::Element::giveNode(), oofem::Element::giveSpatialDimension(), and oofem::FloatArray::resize().
Referenced by giveClassName(), and oofem::TR1_2D_PFEM::giveInternalStateAtNode().
|
pure virtual |
Returns mask of pressure Dofs.
Implemented in oofem::TR1_2D_PFEM.
Referenced by giveCharacteristicMatrix(), giveCharacteristicVector(), giveClassName(), and oofem::PFEMPressureRhsAssembler::vectorFromElement().
|
pure virtual |
Returns the interpolation for the pressure.
Implemented in oofem::TR1_2D_PFEM, and oofem::PFEMElement2d.
|
pure virtual |
Returns mask of velocity Dofs.
Implemented in oofem::TR1_2D_PFEM.
Referenced by computeLoadVector(), giveCharacteristicMatrix(), giveCharacteristicVector(), and giveClassName().
|
pure virtual |
Returns the interpolation for velocity.
Implemented in oofem::TR1_2D_PFEM, and oofem::PFEMElement2d.
|
virtual |
Initializes receiver acording to object description stored in input record.
Reimplemented from oofem::Element.
Reimplemented in oofem::TR1_2D_PFEM, and oofem::PFEMElement2d.
Definition at line 69 of file pfemelement.C.
References oofem::Element::computeGaussPoints(), oofem::Element::initializeFrom(), and oofem::IRRT_OK.
Referenced by oofem::PFEMElement2d::initializeFrom(), and oofem::TR1_2D_PFEM::initializeFrom().
|
virtual |
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 219 of file pfemelement.C.
References oofem::Element::giveGlobalNumber(), oofem::FEMComponent::giveNumber(), oofem::Element::integrationRulesArray, and oofem::FEMComponent::number.
|
virtual |
Updates element state after equilibrium in time step has been reached.
Default implementation updates all integration rules defined by integrationRulesArray member variable. Doing this, all integration points and their material statuses are updated also. All temporary history variables, which now describe equilibrium state are copied into equilibrium ones. The existing internal state is used for update.
tStep | Time step for newly reached state. |
Reimplemented from oofem::Element.
Definition at line 207 of file pfemelement.C.
References computeDeviatoricStress(), and oofem::Element::integrationRulesArray.
Referenced by oofem::PFEM::applyIC().