35 #include "../sm/Elements/linedistributedspring.h" 36 #include "../sm/Materials/structuralms.h" 37 #include "../sm/CrossSections/structuralcrosssection.h" 91 OOFEM_ERROR(
"Body load not supported, use surface load instead");
105 answer.
resize(ndofs, ndofs*2);
108 for (
int idof=1; idof<=ndofs; idof++) {
109 answer.
at(idof, idof) = n.
at(1);
110 answer.
at(idof, ndofs+idof) = n.
at(2);
121 for (
int idof=1; idof<=ndofs; idof++) {
133 answer.
resize(ndofs, ndofs);
140 TimeStep *tStep,
int useUpdatedGpRecord)
182 fprintf(File,
"\tGP 1.%d :\n", gp->giveNumber());
185 fprintf(File,
"\t\tStrain");
186 for (
int i=1; i<=strain.
giveSize(); i++) fprintf (File,
" %e", strain.
at(i));
187 fprintf(File,
"\n\t\tStress");
188 for (
int i=1; i<=stress.
giveSize(); i++) fprintf (File,
" %e", stress.
at(i));
210 return detJ * weight;
244 for (
int i = 1; i < 3; i++ ) {
255 for (
int i = 1; i < 3; i++ ) {
CrossSection * giveCrossSection()
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
The element interface required by ZZNodalRecoveryModel.
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
LineDistributedSpring(int n, Domain *d)
void computeVectorOf(ValueModeType u, TimeStep *tStep, FloatArray &answer)
Returns local vector of unknowns.
virtual IntegrationRule * giveDefaultIntegrationRulePtr()
Access method for default integration rule.
static FEI3dLineLin interp_lin
int checkConsistency()
Performs consistency check.
The element interface required by ZZNodalRecoveryModel.
virtual void giveDofManDofIDMask(int inode, IntArray &) const
Returns dofmanager dof mask for node.
double & at(int i)
Coefficient access function.
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
Computes constitutive matrix of receiver.
#define _IFT_LineDistributedSpring_Stifnesses
void beDiagonal(const FloatArray &diag)
Modifies receiver to be a diagonal matrix with the components specified in diag.
virtual double giveTransformationJacobian(const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the determinant of the transformation.
FloatArray springStiffnesses
Class implementing an array of integers.
int & at(int i)
Coefficient access function.
MatResponseMode
Describes the character of characteristic material matrix.
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
virtual void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)
Computes the stress vector of receiver at given integration point, at time step tStep.
Class representing a general abstraction for finite element interpolation class.
Abstract base class for all "structural" finite elements.
virtual void SPRNodalRecoveryMI_giveSPRAssemblyPoints(IntArray &pap)
DofIDItem
Type representing particular dof type.
virtual void evalN(FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the array of interpolation functions (shape functions) at given point.
virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int=1, int=ALL_STRAINS)
Computes the geometrical matrix of receiver in given integration point.
virtual double giveWeight()
Returns integration weight of receiver.
virtual Interface * giveInterface(InterfaceType it)
Interface requesting service.
Wrapper around element definition to provide FEICellGeometry interface.
virtual int checkConsistency()
Performs consistency check.
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Receiver becomes the result of the product of aMatrix and anArray.
void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord)
Returns equivalent nodal forces vectors.
virtual void computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep)
Compute strain vector of receiver evaluated at given integration point at given time step from elemen...
double at(int i, int j) const
Coefficient access function.
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
Computes numerically stiffness matrix of receiver.
void resize(int n)
Checks size of receiver towards requested bounds.
virtual int setupIntegrationPoints(IntegrationRule &irule, int npoints, Element *element)
Sets up integration rule for the given element.
int numberOfGaussPoints
Number of integration points as specified by nip.
Class representing vector of real numbers.
Implementation of matrix containing floating point numbers.
IRResultType
Type defining the return values of InputRecord reading operations.
virtual FEInterpolation * giveInterpolation() const
void resize(int rows, int cols)
Checks size of receiver towards requested bounds.
std::vector< std::unique_ptr< IntegrationRule > > integrationRulesArray
List of integration rules of receiver (each integration rule contains associated integration points a...
void zero()
Zeroes all coefficient of receiver.
virtual void computeBodyLoadVectorAt(FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode)
Computes the load vector due to body load acting on receiver, at given time step. ...
void printOutputAt(FILE *File, TimeStep *tStep)
Prints output of receiver to stream, for given time step.
InterfaceType
Enumerative type, used to identify interface type.
#define _IFT_LineDistributedSpring_Dofs
virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes lumped mass matrix of receiver.
Load is base abstract class for all loads.
int giveSize() const
Returns the size of receiver.
the oofem namespace is to define a context or scope in which all oofem names are defined.
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
Node * giveNode(int i) const
Returns reference to the i-th node of element.
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Class representing integration point in finite element program.
virtual double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
Class representing solution step.
virtual void SPRNodalRecoveryMI_giveDofMansDeterminedByPatch(IntArray &answer, int pap)
int numberOfDofMans
Number of dofmanagers.
const FloatArray & giveNaturalCoordinates()
Returns coordinate array of receiver.
Class representing Gaussian-quadrature integration rule.
void resize(int s)
Resizes receiver towards requested size.