OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
oofem::CrossSection Class Referenceabstract

Base abstract class representing cross section in finite element mesh. More...

#include <crosssection.h>

+ Inheritance diagram for oofem::CrossSection:
+ Collaboration diagram for oofem::CrossSection:

Public Member Functions

 CrossSection (int n, Domain *d)
 Constructor. More...
 
virtual ~CrossSection ()
 Destructor. More...
 
int giveSetNumber () const
 
virtual bool hasProperty (CrossSectionProperty a)
 Returns true if the dictionary contains the requested property. More...
 
virtual double give (CrossSectionProperty a, GaussPoint *gp)
 Returns the value of cross section property at given point. More...
 
virtual double give (CrossSectionProperty a, const FloatArray &coords, Element *elem, bool local=true)
 Returns the value of cross section property at given point (belonging to given element). More...
 
virtual double give (int aProperty, GaussPoint *gp)
 Returns the value of cross section property. More...
 
virtual bool isCharacteristicMtrxSymmetric (MatResponseMode rMode)
 Check for symmetry of stiffness matrix. More...
 
virtual void printYourself ()
 Prints receiver state on stdout. Useful for debugging. More...
 
virtual int setupIntegrationPoints (IntegrationRule &irule, int npoints, Element *element)
 Sets up integration rule for the given element. More...
 
virtual int setupIntegrationPoints (IntegrationRule &irule, int npointsXY, int npointsZ, Element *element)
 Sets up integration rule for the given element. More...
 
virtual int testCrossSectionExtension (CrossSectExtension ext)
 Returns nonzero, if receiver implements required extension. More...
 
virtual int giveIPValue (FloatArray &answer, GaussPoint *ip, InternalStateType type, TimeStep *tStep)
 Returns the integration point corresponding value in Reduced form. More...
 
virtual int packUnknowns (DataStream &buff, TimeStep *tStep, GaussPoint *ip)=0
 Pack all necessary data of integration point (according to element parallel_mode) into given communication buffer. More...
 
virtual int unpackAndUpdateUnknowns (DataStream &buff, TimeStep *tStep, GaussPoint *ip)=0
 Unpack and updates all necessary data of given integration point (according to element parallel_mode) into given communication buffer. More...
 
virtual int estimatePackSize (DataStream &buff, GaussPoint *ip)=0
 Estimates the necessary pack size to hold all packed data of receiver. More...
 
virtual double predictRelativeComputationalCost (GaussPoint *ip)
 Returns the weight representing relative computational cost of receiver The reference cross section is integral model in plane stress. More...
 
virtual double giveRelativeSelfComputationalCost ()
 Returns the weight representing relative computational cost of receiver The reference element is integral model in plane stress. More...
 
virtual double predictRelativeRedistributionCost (GaussPoint *gp)
 
virtual IRResultType initializeFrom (InputRecord *ir)
 Initializes receiver according to object description stored in input record. More...
 
virtual void giveInputRecord (DynamicInputRecord &input)
 Setups the input record string of receiver. More...
 
virtual MaterialgiveMaterial (IntegrationPoint *ip)=0
 Returns the material associated with the GP. More...
 
virtual contextIOResultType saveIPContext (DataStream &stream, ContextMode mode, GaussPoint *gp)
 Stores integration point state to output stream. More...
 
virtual contextIOResultType restoreIPContext (DataStream &stream, ContextMode mode, GaussPoint *gp)
 Reads integration point state to output stream. 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 * giveClassName () const =0
 
virtual const char * giveInputRecordName () const =0
 
DomaingiveDomain () const
 
virtual void setDomain (Domain *d)
 Sets associated Domain. More...
 
int giveNumber () const
 
void setNumber (int num)
 Sets number of receiver. More...
 
virtual void updateLocalNumbering (EntityRenumberingFunctor &f)
 Local renumbering support. More...
 
virtual contextIOResultType saveContext (DataStream &stream, ContextMode mode, void *obj=NULL)
 Stores receiver state to output stream. More...
 
virtual contextIOResultType restoreContext (DataStream &stream, ContextMode mode, void *obj=NULL)
 Restores the receiver state previously written in stream. More...
 
virtual int checkConsistency ()
 Allows programmer to test some internal data, before computation begins. More...
 
virtual void printOutputAt (FILE *file, TimeStep *tStep)
 Prints output of receiver to stream, for given time step. More...
 
virtual InterfacegiveInterface (InterfaceType t)
 Interface requesting service. More...
 
std::string errorInfo (const char *func) const
 Returns string for prepending output (used by error reporting macros). More...
 

Protected Attributes

Dictionary propertyDictionary
 Dictionary for storing cross section parameters (like dimensions). More...
 
int setNumber
 
- Protected Attributes inherited from oofem::FEMComponent
int number
 Component number. More...
 
Domaindomain
 Link to domain object, useful for communicating with other FEM components. More...
 

Detailed Description

Base abstract class representing cross section in finite element mesh.

The main idea, why cross section has been introduced, is to hide all details of cross section description from particular element. Generally elements do not communicate directly with material but communicate through cross section interface, which therefore can perform necessary integration (for example over layers of fibers).

The cross section returns properties like thickness and area.

The derived classes are supposed to be base cross section classes for particular type of analysis. They should declare general interface methods necessary.

In particular cross section implementation, where is necessary to perform integration over cross section volume (over layers, fibers, ...) and therefore generally one must keep complete load history in these integration points, the concept of master-slave integration points should be used.

Integration point generally can contain list of slave integration points therefore is called as master point. Slaves are used for example to implement layered or fibered cross sections by cross section class. Then in one "macro" master Gauss point, cross section creates few slaves (one per layer) and puts them into master list. When cross sections completes requests for particular master integration point, it performs integration over layers. It therefore calls material class for each layer, sending corresponding slave as parameter and integrates results.

See also
GaussPoint class for more detail.

Definition at line 107 of file crosssection.h.

Constructor & Destructor Documentation

oofem::CrossSection::CrossSection ( int  n,
Domain d 
)

Constructor.

Creates cross section with number n belonging to domain d.

Parameters
nCross section number.
dDomain.

Definition at line 45 of file crosssection.C.

oofem::CrossSection::~CrossSection ( )
virtual

Destructor.

Definition at line 49 of file crosssection.C.

Member Function Documentation

virtual int oofem::CrossSection::estimatePackSize ( DataStream buff,
GaussPoint ip 
)
pure virtual

Estimates the necessary pack size to hold all packed data of receiver.

The corresponding material model service is invoked. The nature of packed data is typically material model dependent.

Parameters
buffCommunication buffer.
ipIntegration point.
Returns
Estimate of pack size.

Implemented in oofem::LayeredCrossSection, oofem::FiberedCrossSection, oofem::SimpleCrossSection, oofem::StructuralInterfaceCrossSection, oofem::SimpleTransportCrossSection, oofem::FluidCrossSection, and oofem::EmptyCS.

Referenced by oofem::Element::estimatePackSize().

double oofem::CrossSection::give ( CrossSectionProperty  a,
GaussPoint gp 
)
virtual

Returns the value of cross section property at given point.

The default implementation assumes constant properties stored in propertyDictionary.

Parameters
aId of requested property.
gpIntegration point
Returns
Property value.

Reimplemented in oofem::LayeredCrossSection, oofem::FiberedCrossSection, oofem::SimpleCrossSection, and oofem::VariableCrossSection.

Definition at line 151 of file crosssection.C.

References oofem::Dictionary::at(), oofem::Dictionary::includes(), OOFEM_ERROR, and propertyDictionary.

Referenced by oofem::IntElLine1PF::computeAreaAround(), oofem::IntElLine1PhF::computeAreaAround(), oofem::IntElLine1::computeAreaAround(), oofem::TR1_2D_SUPG::computeBCRhsTerm_MB(), oofem::TrPlaneStrRot::computeBodyLoadVectorAt(), oofem::CCTPlate::computeBodyLoadVectorAt(), oofem::TrPlaneStrRot3d::computeBodyLoadVectorAt(), oofem::TrPlanestressRotAllman3d::computeBodyLoadVectorAt(), oofem::DKTPlate3d::computeBodyLoadVectorAt(), oofem::CCTPlate3d::computeBodyLoadVectorAt(), oofem::DKTPlate::computeBodyLoadVectorAt(), oofem::QDKTPlate::computeBodyLoadVectorAt(), oofem::Quad1Mindlin::computeBodyLoadVectorAt(), oofem::RerShell::computeBodyLoadVectorAt(), oofem::Quad1MindlinShell3D::computeBodyLoadVectorAt(), oofem::StructuralElement::computeBodyLoadVectorAt(), oofem::XfemStructuralElementInterface::computeCohesiveForces(), oofem::XfemStructuralElementInterface::computeCohesiveTangent(), oofem::Beam2d::computeConsistentMassMatrix(), oofem::Beam3d::computeConsistentMassMatrix(), oofem::StructuralElementEvaluator::computeConsistentMassMatrix(), oofem::StructuralElement::computeConsistentMassMatrix(), oofem::Shell7Base::computeConvectiveMassForce(), oofem::Tr1_ht::computeEdgeVolumeAround(), oofem::Quad1_ht::computeEdgeVolumeAround(), oofem::Beam3d::computeInternalForcesFromBodyLoadVectorAtPoint(), oofem::TR1_2D_SUPG::computeLoadVector(), oofem::SurfaceTensionBoundaryCondition::computeLoadVectorFromElement(), oofem::RerShell::computeLocalCoordinates(), oofem::CCTPlate3d::computeLocalCoordinates(), oofem::CCTPlate::computeLocalCoordinates(), oofem::DKTPlate3d::computeLocalCoordinates(), oofem::MITC4Shell::computeLocalCoordinates(), oofem::QDKTPlate::computeLocalCoordinates(), oofem::DKTPlate::computeLocalCoordinates(), oofem::LIBeam2dNL::computeLumpedMassMatrix(), oofem::LIBeam2d::computeLumpedMassMatrix(), oofem::LIBeam3d::computeLumpedMassMatrix(), oofem::Truss3d::computeLumpedMassMatrix(), oofem::Truss1d::computeLumpedMassMatrix(), oofem::LTRSpace::computeLumpedMassMatrix(), oofem::Truss2d::computeLumpedMassMatrix(), oofem::RerShell::computeLumpedMassMatrix(), oofem::LIBeam3dNL::computeLumpedMassMatrix(), oofem::LIBeam3dNL2::computeLumpedMassMatrix(), oofem::LIBeam3d2::computeLumpedMassMatrix(), oofem::Quad1Mindlin::computeLumpedMassMatrix(), oofem::Quad1MindlinShell3D::computeLumpedMassMatrix(), oofem::CCTPlate::computeLumpedMassMatrix(), oofem::DKTPlate::computeLumpedMassMatrix(), oofem::Quad1MindlinShell3D::computeStiffnessMatrix(), oofem::MITC4Shell::computeStiffnessMatrix(), oofem::TrPlanestressRotAllman::computeStiffnessMatrixZeroEnergyStabilization(), oofem::LIBeam2d::computeStrainVectorInLayer(), oofem::LIBeam2dNL::computeStrainVectorInLayer(), oofem::RerShell::computeStrainVectorInLayer(), oofem::TrPlanestressRotAllman::computeStrainVectorInLayer(), oofem::CCTPlate::computeStrainVectorInLayer(), oofem::Beam2d::computeStrainVectorInLayer(), oofem::QDKTPlate::computeStrainVectorInLayer(), oofem::DKTPlate::computeStrainVectorInLayer(), oofem::SurfaceTensionBoundaryCondition::computeTangentFromElement(), oofem::QClinearStatic::computeTotalVolumeOfInterpolationMesh(), oofem::PlaneStressStructuralElementEvaluator::computeVolumeAround(), oofem::Tr1_ht::computeVolumeAround(), oofem::Quad1_ht::computeVolumeAround(), oofem::InterfaceElement3dTrLin::computeVolumeAround(), oofem::InterfaceElem2dQuad::computeVolumeAround(), oofem::QTruss1d::computeVolumeAround(), oofem::Truss3d::computeVolumeAround(), oofem::TrPlaneStrRot3d::computeVolumeAround(), oofem::Truss1d::computeVolumeAround(), oofem::TrPlanestressRotAllman3d::computeVolumeAround(), oofem::Structural2DElement::computeVolumeAround(), oofem::Truss2d::computeVolumeAround(), oofem::SimpleCrossSection::give(), oofem::FiberedCrossSection::give(), oofem::LayeredCrossSection::give(), oofem::TrPlaneStrRot3d::giveCharacteristicTensor(), oofem::MITC4Shell::giveDirectorVectors(), oofem::Quad1MindlinShell3D::giveInternalForcesVector(), oofem::MITC4Shell::giveInternalForcesVector(), oofem::Shell7Base::giveMassFactorsAt(), oofem::MITC4Shell::giveMidplaneIPValue(), oofem::MITC4Shell::giveThickness(), oofem::Quad1_ht::giveThicknessAt(), oofem::Tr1_ht::giveThicknessAt(), oofem::Tr1Darcy::giveThicknessAt(), and oofem::XfemStructuralElementInterface::XfemElementInterface_computeConsistentMassMatrix().

double oofem::CrossSection::give ( CrossSectionProperty  a,
const FloatArray coords,
Element elem,
bool  local = true 
)
virtual

Returns the value of cross section property at given point (belonging to given element).

the point coordinates can be specified using its local element coordinates or global coordinates (one of these two can be set to NULL) The default implementation assumes constant properties stored in propertyDictionary.

Parameters
aId of requested property.
coordslocal or global coordinates (determined by local parameter) of point of interest
elemreference to underlying element containing given point
gpIntegration point
Returns
Property value.

Reimplemented in oofem::LayeredCrossSection, oofem::SimpleCrossSection, and oofem::VariableCrossSection.

Definition at line 164 of file crosssection.C.

References oofem::Dictionary::at(), oofem::Dictionary::includes(), OOFEM_ERROR, and propertyDictionary.

virtual double oofem::CrossSection::give ( int  aProperty,
GaussPoint gp 
)
inlinevirtual

Returns the value of cross section property.

Parameters
aPropertyId of requested property.
gpIntegration point.
Returns
Property value.

Reimplemented in oofem::LayeredCrossSection, oofem::SimpleCrossSection, and oofem::FiberedCrossSection.

Definition at line 163 of file crosssection.h.

int oofem::CrossSection::giveIPValue ( FloatArray answer,
GaussPoint ip,
InternalStateType  type,
TimeStep tStep 
)
virtual

Returns the integration point corresponding value in Reduced form.

Parameters
answercontain corresponding ip value, zero sized if not available
ipIntegration point.
typeDetermines the type of internal variable.
tStepTime step.
Returns
Nonzero if o.k, zero otherwise.

Reimplemented in oofem::LayeredCrossSection, oofem::SimpleCrossSection, oofem::FiberedCrossSection, oofem::StructuralInterfaceCrossSection, oofem::SimpleTransportCrossSection, and oofem::FluidCrossSection.

Definition at line 89 of file crosssection.C.

References oofem::FloatArray::at(), oofem::Material::giveIPValue(), oofem::GaussPoint::giveMaterial(), oofem::FEMComponent::giveNumber(), and oofem::FloatArray::resize().

Referenced by oofem::FiberedCrossSection::giveIPValue(), oofem::Element::giveIPValue(), and oofem::Tr1Darcy::NodalAveragingRecoveryMI_computeNodalValue().

virtual double oofem::CrossSection::giveRelativeSelfComputationalCost ( )
inlinevirtual

Returns the weight representing relative computational cost of receiver The reference element is integral model in plane stress.

Its weight is equal to 1.0. The other cross section models should compare to this reference.

Returns
Relative computational cost of self.

Definition at line 260 of file crosssection.h.

Referenced by predictRelativeComputationalCost().

int oofem::CrossSection::giveSetNumber ( ) const
inline
bool oofem::CrossSection::hasProperty ( CrossSectionProperty  a)
virtual

Returns true if the dictionary contains the requested property.

Parameters
aId of requested property.

Definition at line 145 of file crosssection.C.

References oofem::Dictionary::includes(), and propertyDictionary.

IRResultType oofem::CrossSection::initializeFrom ( InputRecord ir)
virtual

Initializes receiver according to object description stored in input record.

This function is called immediately after creating object using constructor. Input record can be imagined as data record in component database belonging to receiver. Receiver may use value-name extracting functions to extract particular field from record.

See also
IR_GIVE_FIELD
IR_GIVE_OPTIONAL_FIELD
Parameters
irInput record to initialize from.
Returns
IRResultType

Reimplemented from oofem::FEMComponent.

Reimplemented in oofem::FiberedCrossSection, oofem::SimpleCrossSection, oofem::VariableCrossSection, oofem::LayeredCrossSection, oofem::StructuralInterfaceCrossSection, oofem::SimpleTransportCrossSection, oofem::FluidCrossSection, and oofem::WarpingCrossSection.

Definition at line 67 of file crosssection.C.

References _IFT_CrossSection_SetNumber, IR_GIVE_OPTIONAL_FIELD, oofem::IRRT_OK, and setNumber.

Referenced by oofem::Quasicontinuum::addCrosssectionToInterpolationElements(), oofem::Subdivision::createMesh(), oofem::FluidCrossSection::initializeFrom(), oofem::SimpleTransportCrossSection::initializeFrom(), oofem::StructuralInterfaceCrossSection::initializeFrom(), oofem::LayeredCrossSection::initializeFrom(), oofem::VariableCrossSection::initializeFrom(), oofem::SimpleCrossSection::initializeFrom(), and oofem::T3DInterface::t3d_2_OOFEM().

virtual int oofem::CrossSection::packUnknowns ( DataStream buff,
TimeStep tStep,
GaussPoint ip 
)
pure virtual

Pack all necessary data of integration point (according to element parallel_mode) into given communication buffer.

The corresponding material model service for particular integration point is invoked. The nature of packed data is material model dependent. Typically, for material of "local" response (response depends only on integration point local state) no data are exchanged. For "nonlocal" constitutive models the send/receive of local values which undergo averaging is performed between local and corresponding remote elements.

Parameters
buffCommunication buffer.
tStepSolution step.
ipIntegration point.
Returns
Nonzero if successful.

Implemented in oofem::LayeredCrossSection, oofem::FiberedCrossSection, oofem::SimpleCrossSection, oofem::StructuralInterfaceCrossSection, oofem::SimpleTransportCrossSection, oofem::FluidCrossSection, and oofem::EmptyCS.

Referenced by oofem::Element::packUnknowns().

double oofem::CrossSection::predictRelativeComputationalCost ( GaussPoint ip)
virtual

Returns the weight representing relative computational cost of receiver The reference cross section is integral model in plane stress.

Its weight is equal to 1.0 Default implementation computes average computational cost of material model and multiplies it by cross section type weight (obtained by giveRelativeSelfComputationalCost()) The other cross section models should compare to this reference.

Parameters
ipIntegration point.
Returns
Prediction of the computational cost.

Definition at line 178 of file crosssection.C.

References giveMaterial(), giveRelativeSelfComputationalCost(), and oofem::Material::predictRelativeComputationalCost().

Referenced by oofem::Element::predictRelativeComputationalCost().

virtual double oofem::CrossSection::predictRelativeRedistributionCost ( GaussPoint gp)
inlinevirtual
Returns
Relative redistribution cost of the receiver.

Definition at line 264 of file crosssection.h.

void oofem::CrossSection::printYourself ( )
virtual

Prints receiver state on stdout. Useful for debugging.

Reimplemented from oofem::FEMComponent.

Reimplemented in oofem::LayeredCrossSection, and oofem::FiberedCrossSection.

Definition at line 101 of file crosssection.C.

References oofem::Dictionary::printYourself(), and propertyDictionary.

contextIOResultType oofem::CrossSection::restoreIPContext ( DataStream stream,
ContextMode  mode,
GaussPoint gp 
)
virtual

Reads integration point state to output stream.

Parameters
streamOutput stream.
modeDetermines amount of info required in stream (state, definition, ...).
gpintegration point.
Returns
contextIOResultType.
Exceptions
throwsan ContextIOERR exception if error encountered.

Reimplemented in oofem::LayeredCrossSection, and oofem::FiberedCrossSection.

Definition at line 128 of file crosssection.C.

References oofem::CIO_OK, giveMaterial(), oofem::Material::restoreIPContext(), and THROW_CIOERR.

Referenced by oofem::IntegrationRule::restoreContext(), oofem::FiberedCrossSection::restoreIPContext(), and oofem::LayeredCrossSection::restoreIPContext().

contextIOResultType oofem::CrossSection::saveIPContext ( DataStream stream,
ContextMode  mode,
GaussPoint gp 
)
virtual

Stores integration point state to output stream.

Parameters
streamOutput stream.
modeDetermines amount of info required in stream (state, definition, ...).
gpintegration point.
Returns
contextIOResultType.
Exceptions
throwsan ContextIOERR exception if error encountered.

Reimplemented in oofem::LayeredCrossSection, and oofem::FiberedCrossSection.

Definition at line 110 of file crosssection.C.

References oofem::CIO_OK, giveMaterial(), oofem::Material::saveIPContext(), and THROW_CIOERR.

Referenced by oofem::IntegrationRule::saveContext(), oofem::FiberedCrossSection::saveIPContext(), and oofem::LayeredCrossSection::saveIPContext().

int oofem::CrossSection::setupIntegrationPoints ( IntegrationRule irule,
int  npoints,
Element element 
)
virtual

Sets up integration rule for the given element.

Default behavior is just to call the Gauss integration rule, but for example the layered and fibered crosssections need to do their own thing.

Parameters
iruleIntegration rule to set up.
npointsNumber of integration points.
elementElement which the integration rule belongs to.
Returns
Number of integration points.

Reimplemented in oofem::LayeredCrossSection.

Definition at line 54 of file crosssection.C.

References oofem::Element::giveIntegrationDomain(), oofem::Element::giveMaterialMode(), and oofem::IntegrationRule::setUpIntegrationPoints().

Referenced by oofem::QTrPlaneStrainGrad::computeGaussPoints(), oofem::TrPlaneStrRot::computeGaussPoints(), oofem::QPlaneStressGrad::computeGaussPoints(), oofem::QWedgeGrad::computeGaussPoints(), oofem::QPlaneStrainGrad::computeGaussPoints(), oofem::QTrPlaneStressGrad::computeGaussPoints(), oofem::QTRSpaceGrad::computeGaussPoints(), oofem::SolidShell::computeGaussPoints(), oofem::QSpaceGrad::computeGaussPoints(), oofem::Tr1Darcy::computeGaussPoints(), oofem::Hexa21Stokes::computeGaussPoints(), oofem::QTruss1dGrad::computeGaussPoints(), oofem::Tet21Stokes::computeGaussPoints(), oofem::Quad2PlateSubSoil::computeGaussPoints(), oofem::Tr21Stokes::computeGaussPoints(), oofem::Tr1BubbleStokes::computeGaussPoints(), oofem::CCTPlate::computeGaussPoints(), oofem::QBrick1_ht::computeGaussPoints(), oofem::QTruss1d::computeGaussPoints(), oofem::QDKTPlate::computeGaussPoints(), oofem::Tet1_3D_SUPG::computeGaussPoints(), oofem::Tet1BubbleStokes::computeGaussPoints(), oofem::Tr1_ht::computeGaussPoints(), oofem::Quad1_ht::computeGaussPoints(), oofem::Tetrah1_ht::computeGaussPoints(), oofem::Structural3DElement::computeGaussPoints(), oofem::Wedge_ht::computeGaussPoints(), oofem::QWedge_ht::computeGaussPoints(), oofem::Brick1_ht::computeGaussPoints(), oofem::Structural2DElement::computeGaussPoints(), oofem::TrPlanestressRotAllman::computeGaussPoints(), oofem::Tria1PlateSubSoil::computeGaussPoints(), oofem::CohesiveSurface3d::computeGaussPoints(), oofem::Quad1PlateSubSoil::computeGaussPoints(), oofem::LineDistributedSpring::computeGaussPoints(), oofem::TR1_2D_SUPG_AXI::computeGaussPoints(), oofem::LIBeam2d::computeGaussPoints(), oofem::Tr_Warp::computeGaussPoints(), oofem::Quad1Mindlin::computeGaussPoints(), oofem::LIBeam2dNL::computeGaussPoints(), oofem::LIBeam3d::computeGaussPoints(), oofem::Truss3d::computeGaussPoints(), oofem::tet21ghostsolid::computeGaussPoints(), oofem::MITC4Shell::computeGaussPoints(), oofem::Truss2d::computeGaussPoints(), oofem::Truss1d::computeGaussPoints(), oofem::LIBeam3dNL2::computeGaussPoints(), oofem::LIBeam3dNL::computeGaussPoints(), oofem::RerShell::computeGaussPoints(), oofem::TR21_2D_SUPG::computeGaussPoints(), oofem::Quad1MindlinShell3D::computeGaussPoints(), oofem::LIBeam3d2::computeGaussPoints(), oofem::DKTPlate::computeGaussPoints(), oofem::Quad10_2D_SUPG::computeGaussPoints(), oofem::AxisymElement::computeGaussPoints(), oofem::TR1_2D_CBS::computeGaussPoints(), oofem::Beam2d::computeGaussPoints(), oofem::TR1_2D_SUPG2::computeGaussPoints(), oofem::TR1_2D_SUPG::computeGaussPoints(), oofem::Beam3d::computeGaussPoints(), and oofem::TR1_2D_SUPG2_AXI::updateIntegrationRules().

int oofem::CrossSection::setupIntegrationPoints ( IntegrationRule irule,
int  npointsXY,
int  npointsZ,
Element element 
)
virtual

Sets up integration rule for the given element.

Default behavior is just to call the Gauss integration rule, but for example the layered and fibered crosssections need to do their own thing.

Parameters
iruleIntegration rule to set up.
npointsNumber of integration points.
elementElement which the integration rule belongs to.
Returns
Number of integration points.

Reimplemented in oofem::LayeredCrossSection.

Definition at line 61 of file crosssection.C.

References oofem::Element::giveIntegrationDomain(), oofem::Element::giveMaterialMode(), and oofem::IntegrationRule::setUpIntegrationPoints().

virtual int oofem::CrossSection::testCrossSectionExtension ( CrossSectExtension  ext)
inlinevirtual

Returns nonzero, if receiver implements required extension.

Parameters
extRequired extension.
Returns
Nonzero, if supported, zero otherwise.

Reimplemented in oofem::StructuralCrossSection, and oofem::StructuralInterfaceCrossSection.

Definition at line 198 of file crosssection.h.

virtual int oofem::CrossSection::unpackAndUpdateUnknowns ( DataStream buff,
TimeStep tStep,
GaussPoint ip 
)
pure virtual

Unpack and updates all necessary data of given integration point (according to element parallel_mode) into given communication buffer.

See also
packUnknowns service.
Parameters
buffCommunication buffer.
tStepSolution step.
ipIntegration point.
Returns
Nonzero if successful.

Implemented in oofem::LayeredCrossSection, oofem::FiberedCrossSection, oofem::SimpleCrossSection, oofem::StructuralInterfaceCrossSection, oofem::SimpleTransportCrossSection, oofem::FluidCrossSection, and oofem::EmptyCS.

Referenced by oofem::Element::unpackAndUpdateUnknowns().

Member Data Documentation

Dictionary oofem::CrossSection::propertyDictionary
protected
int oofem::CrossSection::setNumber
protected

Definition at line 117 of file crosssection.h.

Referenced by giveInputRecord(), and initializeFrom().


The documentation for this class was generated from the following files:

This page is part of the OOFEM documentation. Copyright (c) 2011 Borek Patzak
Project e-mail: info@oofem.org
Generated at Tue Jan 2 2018 20:07:34 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011