61LineDistributedSpring :: LineDistributedSpring(
int n,
Domain *aDomain) :
71LineDistributedSpring :: giveInterpolation(
DofIDItem id)
const
78LineDistributedSpring :: giveInterpolation()
const {
return &
interp_lin; }
82LineDistributedSpring :: computeGaussPoints()
94LineDistributedSpring :: computeBodyLoadVectorAt(
FloatArray &answer,
Load *forLoad,
TimeStep *tStep, ValueModeType mode)
96 OOFEM_ERROR(
"Body load not supported, use surface load instead");
106 int ndofs = this->
dofs.giveSize();
110 answer.
resize(ndofs, ndofs*2);
113 for (
int idof=1; idof<=ndofs; idof++) {
114 answer.
at(idof, idof) = n.
at(1);
115 answer.
at(idof, ndofs+idof) = n.
at(2);
123 int ndofs = this->
dofs.giveSize();
126 for (
int idof=1; idof<=ndofs; idof++) {
137 int ndofs = this->
dofs.giveSize();
138 answer.
resize(ndofs, ndofs);
145 TimeStep *tStep,
int useUpdatedGpRecord)
159LineDistributedSpring :: initializeFrom(
InputRecord &ir,
int priority)
168LineDistributedSpring :: postInitialize()
188void LineDistributedSpring :: printOutputAt(FILE *File,
TimeStep *tStep)
195 fprintf(File,
"\tGP 1.%d :\n", gp->giveNumber());
198 fprintf(File,
"\t\tStrain");
199 for (
int i=1; i<=strain.
giveSize(); i++) fprintf (File,
" %e", strain.
at(i));
200 fprintf(File,
"\n\t\tStress");
201 for (
int i=1; i<=stress.
giveSize(); i++) fprintf (File,
" %e", stress.
at(i));
210LineDistributedSpring :: giveDofManDofIDMask(
int inode,
IntArray &answer)
const
223 return detJ * weight;
238 return StructuralElement :: giveIPValue(answer, gp, type, tStep);
254LineDistributedSpring :: SPRNodalRecoveryMI_giveSPRAssemblyPoints(
IntArray &pap)
257 for (
int i = 1; i < 3; i++ ) {
263LineDistributedSpring :: SPRNodalRecoveryMI_giveDofMansDeterminedByPatch(
IntArray &answer,
int pap)
268 for (
int i = 1; i < 3; i++ ) {
#define REGISTER_Element(class)
Node * giveNode(int i) const
void initializeFrom(InputRecord &ir, int priority) override
int numberOfDofMans
Number of dofmanagers.
void computeVectorOf(ValueModeType u, TimeStep *tStep, FloatArray &answer)
void postInitialize() override
Performs post initialization steps.
std::vector< std ::unique_ptr< IntegrationRule > > integrationRulesArray
CrossSection * giveCrossSection()
virtual IntegrationRule * giveDefaultIntegrationRulePtr()
int checkConsistency() override
Domain * giveDomain() const
int number
Component number.
Index giveSize() const
Returns the size of receiver.
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
void resize(Index rows, Index cols)
void beDiagonal(const FloatArray &diag)
void zero()
Zeroes all coefficient of receiver.
double at(std::size_t i, std::size_t j) const
const FloatArray & giveNaturalCoordinates() const
Returns coordinate array of receiver.
double giveWeight()
Returns integration weight of receiver.
int checkConsistency() override
void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord) override
FloatArray springStiffnesses
static ParamKey IPK_LineDistributedSpring_Stiffnesses
void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) override
static ParamKey IPK_LineDistributedSpring_Dofs
static FEI3dLineLin interp_lin
SPRNodalRecoveryModelInterface()
Constructor.
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
StructuralElement(int n, Domain *d)
virtual void computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep)
ZZNodalRecoveryModelInterface(Element *element)
Constructor.
@ SPRNodalRecoveryModelInterfaceType
@ ZZNodalRecoveryModelInterfaceType
#define PM_ELEMENT_ERROR_IFNOTSET(_pm, _componentnum, _paramkey)
#define PM_UPDATE_PARAMETER(_val, _pm, _ir, _componentnum, _paramkey, _prio)