35 #include "../sm/Elements/Shells/cct3d.h" 36 #include "../sm/Materials/structuralms.h" 64 OOFEM_ERROR(
"cannot transform coordinates - size mismatch");
79 double &y1,
double &y2,
double &y3,
80 double &z1,
double &z2,
double &z3)
106 answer = {D_u, D_v, D_w, R_u, R_v, R_w};
118 std::vector< FloatArray > lc(3);
121 for (
int _i = 0; _i < 3; _i++ ) {
127 answer.
at(1) = inputCoords_ElCS.
at(1);
128 answer.
at(2) = inputCoords_ElCS.
at(2);
129 GaussPoint _gp(NULL, 1, answer, 2.0, _2dPlate);
133 return inplane && outofplane;
140 double l1 = lcoords.
at(1);
141 double l2 = lcoords.
at(2);
142 double l3 = 1. - l2 - l1;
145 for (
int _i = 1; _i <= 3; _i++ ) {
185 for (
int i = 1; i <= 3; i++ ) {
207 for (
int i = 1; i <= 3; i++ ) {
231 answer.
at(1, 3) = charVect.
at(4);
232 answer.
at(3, 1) = charVect.
at(4);
233 answer.
at(2, 3) = charVect.
at(5);
234 answer.
at(3, 2) = charVect.
at(5);
239 answer.
at(1, 1) = charVect.
at(1);
240 answer.
at(2, 2) = charVect.
at(2);
241 answer.
at(1, 2) = charVect.
at(3);
242 answer.
at(2, 1) = charVect.
at(3);
247 answer.
at(1, 3) = charVect.
at(4) / 2.;
248 answer.
at(3, 1) = charVect.
at(4) / 2.;
249 answer.
at(2, 3) = charVect.
at(5) / 2.;
250 answer.
at(3, 2) = charVect.
at(5) / 2.;
255 answer.
at(1, 1) = charVect.
at(1);
256 answer.
at(2, 2) = charVect.
at(2);
257 answer.
at(1, 2) = charVect.
at(3) / 2.;
258 answer.
at(2, 1) = charVect.
at(3) / 2.;
280 if ( type == IST_CurvatureTensor || type == IST_ShellStrainTensor ) {
281 if ( type == IST_CurvatureTensor ) {
289 answer.
at(1) = globTensor.
at(1, 1);
290 answer.
at(2) = globTensor.
at(2, 2);
291 answer.
at(3) = globTensor.
at(3, 3);
292 answer.
at(4) = 2 * globTensor.
at(2, 3);
293 answer.
at(5) = 2 * globTensor.
at(1, 3);
294 answer.
at(6) = 2 * globTensor.
at(1, 2);
297 }
else if ( type == IST_ShellMomentTensor || type == IST_ShellForceTensor ) {
298 if ( type == IST_ShellMomentTensor ) {
306 answer.
at(1) = globTensor.
at(1, 1);
307 answer.
at(2) = globTensor.
at(2, 2);
308 answer.
at(3) = globTensor.
at(3, 3);
309 answer.
at(4) = globTensor.
at(2, 3);
310 answer.
at(5) = globTensor.
at(1, 3);
311 answer.
at(6) = globTensor.
at(1, 2);
329 for (
int i = 1; i <= 3; i++ ) {
346 double dens, dV, load;
372 load = force.
at(3) * dens * dV / 3.0;
377 load = force.
at(4) * dens * dV / 3.0;
382 load = force.
at(5) * dens * dV / 3.0;
408 2, 3, 4, 8, 9, 10, 14, 15, 16
411 for (
int i = 0; i < 3; i++ ) {
412 for (
int j = 0; j < 9; j++ ) {
413 answer(ri [ i ], ci [ j ]) = ne(i, j);
473 fprintf( file,
" GP %2d.%-2d :", i + 1, gp->giveNumber() );
475 this->
giveIPValue(v, gp, IST_ShellStrainTensor, tStep);
476 fprintf(file,
" strains ");
477 for (
auto &val : v ) fprintf(file,
" %.4e", val);
479 this->
giveIPValue(v, gp, IST_CurvatureTensor, tStep);
480 fprintf(file,
"\n curvatures ");
481 for (
auto &val : v ) fprintf(file,
" %.4e", val);
484 this->
giveIPValue(v, gp, IST_ShellForceTensor, tStep);
485 fprintf(file,
"\n stresses ");
486 for (
auto &val : v ) fprintf(file,
" %.4e", val);
488 this->
giveIPValue(v, gp, IST_ShellMomentTensor, tStep);
489 fprintf(file,
"\n moments ");
490 for (
auto &val : v ) fprintf(file,
" %.4e", val);
CrossSection * giveCrossSection()
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
void subtract(const FloatArray &src)
Subtracts array src to receiver.
void giveLocalCoordinates(FloatArray &answer, FloatArray &global)
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
void beVectorProductOf(const FloatArray &v1, const FloatArray &v2)
Computes vector product (or cross product) of vectors given as parameters, , and stores the result in...
virtual void giveDofManDofIDMask(int inode, IntArray &) const
Returns dofmanager dof mask for node.
Wrapper around cell with vertex coordinates stored in FloatArray**.
void zero()
Sets all component to zero.
double & at(int i)
Coefficient access function.
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
virtual double computeSurfaceVolumeAround(GaussPoint *gp, int iSurf)
Computes volume related to integration point on local surface.
This class implements a structural material status information.
void clear()
Clears receiver (zero size).
virtual int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords)
Computes the global coordinates from given element's local coordinates.
virtual int computeLoadGToLRotationMtrx(FloatMatrix &answer)
Returns transformation matrix from global coordinate system to local element coordinate system for el...
virtual void computeComponentArrayAt(FloatArray &answer, TimeStep *tStep, ValueModeType mode)
Computes boundary condition value - its components values at given time.
void giveCharacteristicTensor(FloatMatrix &answer, CharTensor type, GaussPoint *gp, TimeStep *tStep)
virtual double giveCoordinate(int i)
Class implementing an array of integers.
int & at(int i)
Coefficient access function.
void rotatedWith(FloatMatrix &r, char mode)
Returns the receiver a rotated according the change-of-base matrix r.
Abstract base class representing integration rule.
virtual int computeLoadLSToLRotationMatrix(FloatMatrix &answer, int iSurf, GaussPoint *gp)
Returns transformation matrix from local surface c.s to element local coordinate system of load vecto...
void beDifferenceOf(const FloatArray &a, const FloatArray &b)
Sets receiver to be a - b.
CCTPlate3d(int n, Domain *d)
This class implements an triangular three-node plate CCT finite element.
bool isNotEmpty() const
Tests for empty matrix.
virtual void printOutputAt(FILE *file, TimeStep *tStep)
Prints output of receiver to stream, for given time step.
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
Class representing a 2d triangular linear interpolation based on area coordinates.
virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
Computes interpolation matrix for element unknowns.
FloatMatrix GtoLRotationMatrix
Transformation Matrix form GtoL(3,3) is stored at the element level for computation efficiency...
StructuralCrossSection * giveStructuralCrossSection()
Helper function which returns the structural cross-section for the element.
void rotatedWith(const FloatMatrix &r, char mode= 'n')
Returns the receiver 'a' transformed using give transformation matrix r.
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Receiver becomes the result of the product of aMatrix and anArray.
virtual int global2local(FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Default implementation using Newton's method to find the local coordinates.
double at(int i, int j) const
Coefficient access function.
void resize(int n)
Checks size of receiver towards requested bounds.
Class representing vector of real numbers.
virtual void computeBodyLoadVectorAt(FloatArray &answer, Load *forLoad, TimeStep *tStep, ValueModeType mode)
Computes the load vector due to body load acting on receiver, at given time step. ...
Implementation of matrix containing floating point numbers.
virtual double give(CrossSectionProperty a, GaussPoint *gp)
Returns the value of cross section property at given point.
GaussPoint * getIntegrationPoint(int n)
Access particular integration point of receiver.
const FloatArray & giveStressVector() const
Returns the const pointer to receiver's stress vector.
IntegrationPointStatus * giveMaterialStatus()
Returns reference to associated material status (NULL if not defined).
void resize(int rows, int cols)
Checks size of receiver towards requested bounds.
virtual int getRequiredNumberOfIntegrationPoints(integrationDomain dType, int approxOrder)
Abstract service.
void zero()
Zeroes all coefficients of receiver.
virtual bcGeomType giveBCGeoType() const
Returns geometry character of boundary condition.
virtual void giveNodeCoordinates(double &x1, double &x2, double &x3, double &y1, double &y2, double &y3, double &z1, double &z2, double &z3)
virtual const FloatMatrix * computeGtoLRotationMatrix()
std::vector< std::unique_ptr< IntegrationRule > > integrationRulesArray
List of integration rules of receiver (each integration rule contains associated integration points a...
virtual int SetUpPointsOnTriangle(int, MaterialMode mode)
Sets up receiver's integration points on triangular (area coords) integration domain.
void zero()
Zeroes all coefficient of receiver.
virtual void giveSurfaceDofMapping(IntArray &answer, int iSurf) const
Assembles surface dof mapping mask, which provides mapping between surface local DOFs and "global" el...
Load is base abstract class for all loads.
virtual void computeSurfaceNMatrixAt(FloatMatrix &answer, int iSurf, GaussPoint *gp)
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 double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
virtual bcValType giveBCValType() const
Returns receiver load type.
double normalize()
Normalizes receiver.
Node * giveNode(int i) const
Returns reference to the i-th node of element.
const FloatArray & giveStrainVector() const
Returns the const pointer to receiver's strain vector.
Class representing integration point in finite element program.
Class representing solution step.
virtual bool computeLocalCoordinates(FloatArray &answer, const FloatArray &gcoords)
Computes the element local coordinates from given global coordinates.
virtual IntegrationRule * GetSurfaceIntegrationRule(int iSurf)
const FloatArray & giveNaturalCoordinates()
Returns coordinate array of receiver.
virtual int SetUpPointsOnTriangle(int nPoints, MaterialMode mode)
Sets up receiver's integration points on triangular (area coords) integration domain.
Class representing Gaussian-quadrature integration rule.
void resize(int s)
Resizes receiver towards requested size.