55LIBeam3dBoundary :: LIBeam3dBoundary(
int n,
Domain *aDomain) :
LIBeam3d(n, aDomain)
65LIBeam3dBoundary :: initializeFrom(
InputRecord &ir,
int priority)
67 LIBeam3d :: initializeFrom(ir, priority);
75LIBeam3dBoundary :: postInitialize()
77 LIBeam3d :: postInitialize();
87LIBeam3dBoundary :: giveDofManDofIDMask(
int inode,
IntArray &answer)
const
90 answer = { E_xx, E_xy, E_xz, E_yx, E_yy, E_yz, E_zx, E_zy, E_zz };
92 answer = { D_u, D_v, D_w, R_u, R_v, R_w };
101 for (
int x = -1; x < 2; x++ ) {
102 for (
int y = -1; y < 2; y++ ) {
103 for (
int z = -1; z < 2; z++ ) {
104 if ( !( z == 0 && y == 0 && x == 0 ) ) {
120LIBeam3dBoundary :: recalculateCoordinates(
int nodeNumber,
FloatArray &coords)
124 unitCellSize.
at(1) = this->
giveNode(3)->giveCoordinate(1);
125 unitCellSize.
at(2) = this->
giveNode(3)->giveCoordinate(2);
126 unitCellSize.
at(3) = this->
giveNode(3)->giveCoordinate(3);
132 coords.
at(1) = this->
giveNode(nodeNumber)->giveCoordinate(1) + switches.
at(1) * unitCellSize.
at(1);
133 coords.
at(2) = this->
giveNode(nodeNumber)->giveCoordinate(2) + switches.
at(2) * unitCellSize.
at(2);
134 coords.
at(3) = this->
giveNode(nodeNumber)->giveCoordinate(3) + switches.
at(3) * unitCellSize.
at(3);
141LIBeam3dBoundary :: computeLength()
150 dx = nodeB.
at(1) - nodeA.
at(1);
151 dy = nodeB.
at(2) - nodeA.
at(2);
152 dz = nodeB.
at(3) - nodeA.
at(3);
153 length = sqrt(dx * dx + dy * dy + dz * dz);
167 n1 = ( 1. - ksi ) * 0.5;
168 n2 = ( 1. + ksi ) * 0.5;
173 answer.
at(1) = n1 * coordsNodeA.
at(1) + n2 * coordsNodeB.
at(1);
174 answer.
at(2) = n1 * coordsNodeA.
at(2) + n2 * coordsNodeB.
at(2);
175 answer.
at(3) = n1 * coordsNodeA.
at(3) + n2 * coordsNodeB.
at(3);
199 for (
int i = 1; i <= 3; i++ ) {
200 lx.at(i) = ( coordsNodeB.
at(i) - coordsNodeA.
at(i) ) /
length;
204 lz.beVectorProductOf(lx, help);
206 ly.beVectorProductOf(lz, lx);
209 for (
int i = 1; i <= 3; i++ ) {
210 answer.
at(1, i) = lx.at(i);
211 answer.
at(2, i) = ly.at(i);
212 answer.
at(3, i) = lz.at(i);
228 for (
int i = 1; i <= 3; i++ ) {
229 for (
int j = 1; j <= 3; j++ ) {
230 answer.
at(i, j) = lcs.
at(i, j);
231 answer.
at(i + 3, j + 3) = lcs.
at(i, j);
232 answer.
at(i + 6, j + 6) = lcs.
at(i, j);
233 answer.
at(i + 9, j + 9) = lcs.
at(i, j);
240 answer.
at(i, j) = 1.;
255 StructuralElement :: computeStiffnessMatrix(Korig, rMode, tStep);
279 unitCellSize.
at(1) = this->
giveNode(3)->giveCoordinate(1);
280 unitCellSize.
at(2) = this->
giveNode(3)->giveCoordinate(2);
281 unitCellSize.
at(3) = this->
giveNode(3)->giveCoordinate(3);
291 for (
int i = 1; i <= 3; i++ ) {
292 k1.
at(i, 3 * i - 2) = unitCellSize.
at(1) * switches1.
at(1);
293 k1.
at(i, 3 * i - 1) = unitCellSize.
at(2) * switches1.
at(2);
294 k1.
at(i, 3 * i) = unitCellSize.
at(3) * switches1.
at(3);
297 for (
int i = 1; i <= 3; i++ ) {
298 k2.
at(i, 3 * i - 2) = unitCellSize.
at(1) * switches2.
at(1);
299 k2.
at(i, 3 * i - 1) = unitCellSize.
at(2) * switches2.
at(2);
300 k2.
at(i, 3 * i) = unitCellSize.
at(3) * switches2.
at(3);
307 answer.
assemble(k1, { 1, 2, 3, 4, 5, 6 }, { 13, 14, 15, 16, 17, 18, 19, 20, 21 });
308 answer.
assemble(k2, { 7, 8, 9, 10, 11, 12 }, { 13, 14, 15, 16, 17, 18, 19, 20, 21 });
340LIBeam3dBoundary :: giveInternalForcesVector(
FloatArray &answer,
TimeStep *tStep,
int useUpdatedGpRecord)
343 FloatArray vStress, vStrain, fintsub, fintG, fintGT;
352 if ( useUpdatedGpRecord == 1 ) {
357 vStrain.
resize(StructuralMaterial :: giveSizeOfVoigtSymVector(gp->giveMaterialMode() ) );
369 StructuralMaterial :: giveReducedSymVectorForm(stressTemp, vStress, gp->giveMaterialMode() );
393 if ( type == IST_DisplacementVector ) {
402 return LIBeam3d :: giveIPValue(answer, gp, type, tStep);
#define REGISTER_Element(class)
double giveCoordinate(int i) const
Node * giveNode(int i) const
virtual bool giveRotationMatrix(FloatMatrix &answer)
virtual bool isActivated(TimeStep *tStep)
int numberOfDofMans
Number of dofmanagers.
void computeVectorOf(ValueModeType u, TimeStep *tStep, FloatArray &answer)
virtual IntegrationRule * giveDefaultIntegrationRulePtr()
Domain * giveDomain() const
int number
Component number.
void plusProduct(const FloatMatrix &b, const FloatArray &s, double dV)
Index giveSize() const
Returns the size of receiver.
void resizeWithValues(Index s, std::size_t allocChunk=0)
void zero()
Zeroes all coefficients of receiver.
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
void subtract(const FloatArray &src)
void rotatedWith(const FloatMatrix &r, char mode='n')
void resizeWithData(Index, Index)
void resize(Index rows, Index cols)
void beSubMatrixOf(const FloatMatrix &src, Index topRow, Index bottomRow, Index topCol, Index bottomCol)
void beProductOf(const FloatMatrix &a, const FloatMatrix &b)
void beTranspositionOf(const FloatMatrix &src)
void zero()
Zeroes all coefficient of receiver.
void assemble(const FloatMatrix &src, const IntArray &loc)
double at(std::size_t i, std::size_t j) const
void beUnitMatrix()
Sets receiver to unity matrix.
const FloatArray & giveSubPatchCoordinates() const
Returns local sub-patch coordinates of the receiver.
void recalculateCoordinates(int nodeNumber, FloatArray &coords) override
static ParamKey IPK_LIBeam3dBoundary_refnode
double computeLength() override
virtual void computeTransformationMatrix(FloatMatrix &answer, TimeStep *tStep)
void computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep) override
void giveSwitches(IntArray &answer, int location)
static ParamKey IPK_LIBeam3dBoundary_location
int giveLocalCoordinateSystem(FloatMatrix &answer) override
int computeNumberOfDofs() override
double computeVolumeAround(GaussPoint *gp) override
void computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) override
void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int=1, int=ALL_STRAINS) override
Node * giveReferenceNode(int refNode)
LIBeam3d(int n, Domain *d)
void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer) override
std::unique_ptr< FloatArray > initialDisplacements
Initial displacement vector, describes the initial nodal displacements when element has been casted.
const FloatArray & giveStressVector() const
Returns the const pointer to receiver's stress vector.
#define PM_ELEMENT_ERROR_IFNOTSET(_pm, _componentnum, _paramkey)
#define PM_UPDATE_PARAMETER(_val, _pm, _ir, _componentnum, _paramkey, _prio)