35 #include "../sm/Elements/3D/ltrspace.h" 36 #include "../sm/CrossSections/structuralcrosssection.h" 115 mss1 = dV * density / 4.;
117 for (
int i = 1; i <= 12; i++ ) {
118 answer.
at(i, i) = mss1;
156 for (
int i = 1; i <= 4; i++ ) {
186 int &localNodeId,
int &localElemId,
int &localBcId,
190 FloatArray *corner [ 4 ], midSide [ 6 ], midFace [ 4 ], midNode;
191 double x = 0.0, y = 0.0, z = 0.0;
192 int inode, nodes = 4, iside, sides = 6, iface, faces = 4, nd, nd1, nd2;
194 static int sideNode [ 6 ] [ 2 ] = { { 1, 2 }, { 2, 3 }, { 3, 1 }, { 1, 4 }, { 2, 4 }, { 3, 4 } };
195 static int faceNode [ 4 ] [ 3 ] = { { 1, 2, 3 }, { 1, 2, 4 }, { 2, 3, 4 }, { 3, 1, 4 } };
201 int hexaSideNode [ 4 ] [ 3 ] = { { 1, 3, 4 }, { 2, 1, 5 }, { 3, 2, 6 }, { 4, 6, 5 } };
202 int hexaFaceNode [ 4 ] [ 3 ] = { { 1, 2, 4 }, { 1, 3, 2 }, { 1, 4, 3 }, { 4, 2, 3 } };
206 for ( inode = 0; inode < nodes; inode++ ) {
209 x += corner [ inode ]->
at(1);
210 y += corner [ inode ]->
at(2);
211 z += corner [ inode ]->
at(3);
214 for ( iside = 0; iside < sides; iside++ ) {
215 midSide [ iside ].
resize(3);
217 nd1 = sideNode [ iside ] [ 0 ] - 1;
218 nd2 = sideNode [ iside ] [ 1 ] - 1;
220 midSide [ iside ].
at(1) = ( corner [ nd1 ]->
at(1) + corner [ nd2 ]->
at(1) ) / 2.0;
221 midSide [ iside ].
at(2) = ( corner [ nd1 ]->
at(2) + corner [ nd2 ]->
at(2) ) / 2.0;
222 midSide [ iside ].
at(3) = ( corner [ nd1 ]->
at(3) + corner [ nd2 ]->
at(3) ) / 2.0;
227 midNode.
at(1) = x / nodes;
228 midNode.
at(2) = y / nodes;
229 midNode.
at(3) = z / nodes;
231 for ( iface = 0; iface < faces; iface++ ) {
233 for ( inode = 0; inode < 3; inode++ ) {
234 nd = faceNode [ iface ] [ inode ] - 1;
235 x += corner [ nd ]->
at(1);
236 y += corner [ nd ]->
at(2);
237 z += corner [ nd ]->
at(3);
240 midFace [ iface ].
resize(3);
242 midFace [ iface ].
at(1) = x / 3;
243 midFace [ iface ].
at(2) = y / 3;
244 midFace [ iface ].
at(3) = z / 3;
249 sMode, tStep, nodes, corner, midSide, midFace, midNode,
250 localNodeId, localElemId, localBcId, hexaSideNode, hexaFaceNode,
251 controlNode, controlDof, aMode,
"LSpace");
260 #define TR_LENGHT_REDUCT 0.3333 274 EASValsSetEdgeFlag(
true);
276 EASValsSetFillStyle(FILL_SOLID);
277 #ifdef __PARALLEL_MODE 297 EGWithMaskChangeAttributes(WIDTH_MASK | FILL_MASK | COLOR_MASK | EDGE_COLOR_MASK | EDGE_FLAG_MASK | LAYER_MASK, go);
298 EGAttachObject(go, ( EObjectP )
this);
299 EMAddGraphicsToModel(ESIModel(), go);
316 EASValsSetEdgeFlag(
true);
318 EASValsSetFillStyle(FILL_SOLID);
333 EGWithMaskChangeAttributes(WIDTH_MASK | FILL_MASK | COLOR_MASK | EDGE_COLOR_MASK | EDGE_FLAG_MASK | LAYER_MASK, go);
334 EMAddGraphicsToModel(ESIModel(), go);
340 int i, indx, result = 0;
344 double s [ 4 ], defScale = 0.0;
351 for ( i = 1; i <= 4; i++ ) {
371 for ( i = 1; i <= 4; i++ ) {
372 s [ i - 1 ] = v [ i - 1 ].
at(indx);
376 EASValsSetEdgeFlag(
true);
379 for ( i = 0; i < 4; i++ ) {
393 tr = CreateTetraWD(p, s);
394 EGWithMaskChangeAttributes(LAYER_MASK | EDGE_COLOR_MASK | EDGE_FLAG_MASK, tr);
395 EMAddGraphicsToModel(ESIModel(), tr);
414 double xc, yc, zc, length;
424 if ( this->
giveIPValue(cf, gp, IST_CrackedFlag, tStep) == 0 ) {
428 if ( (
int ) cf.
at(1) == 0 ) {
436 for ( i = 0; i < 4; i++ ) {
453 if ( this->
giveIPValue(crackDir, gp, IST_CrackDirs, tStep) ) {
454 this->
giveIPValue(crackStatuses, gp, IST_CrackStatuses, tStep);
457 for ( i = 1; i <= 3; i++ ) {
458 crackStatus = ( int ) crackStatuses.
at(i);
466 }
else if ( i == 2 ) {
474 q [ 0 ].x = ( FPNum ) xc + 0.5 * crackDir.
at(0 + j) * length + 0.5 * crackDir.
at(0 + k) * length;
475 q [ 0 ].y = ( FPNum ) yc + 0.5 * crackDir.
at(3 + j) * length + 0.5 * crackDir.
at(3 + k) * length;
476 q [ 0 ].z = ( FPNum ) zc + 0.5 * crackDir.
at(6 + j) * length + 0.5 * crackDir.
at(6 + k) * length;
477 q [ 1 ].x = ( FPNum ) xc + 0.5 * crackDir.
at(0 + j) * length - 0.5 * crackDir.
at(0 + k) * length;
478 q [ 1 ].y = ( FPNum ) yc + 0.5 * crackDir.
at(3 + j) * length - 0.5 * crackDir.
at(3 + k) * length;
479 q [ 1 ].z = ( FPNum ) zc + 0.5 * crackDir.
at(6 + j) * length - 0.5 * crackDir.
at(6 + k) * length;
480 q [ 2 ].x = ( FPNum ) xc - 0.5 * crackDir.
at(0 + j) * length - 0.5 * crackDir.
at(0 + k) * length;
481 q [ 2 ].y = ( FPNum ) yc - 0.5 * crackDir.
at(3 + j) * length - 0.5 * crackDir.
at(3 + k) * length;
482 q [ 2 ].z = ( FPNum ) zc - 0.5 * crackDir.
at(6 + j) * length - 0.5 * crackDir.
at(6 + k) * length;
483 q [ 3 ].x = ( FPNum ) xc - 0.5 * crackDir.
at(0 + j) * length + 0.5 * crackDir.
at(0 + k) * length;
484 q [ 3 ].y = ( FPNum ) yc - 0.5 * crackDir.
at(3 + j) * length + 0.5 * crackDir.
at(3 + k) * length;
485 q [ 3 ].z = ( FPNum ) zc - 0.5 * crackDir.
at(6 + j) * length + 0.5 * crackDir.
at(6 + k) * length;
496 tr = CreateQuad3D(q);
497 EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, tr);
498 EMAddGraphicsToModel(ESIModel(), tr);
virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes lumped mass matrix of receiver.
virtual IntegrationRule * GetSurfaceIntegrationRule(int)
int testElementGraphicActivity(Element *)
Test if particular element passed fulfills various filtering criteria for its graphics output...
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
The element interface required by NodalAvergagingRecoveryModel.
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 IntegrationRule * giveDefaultIntegrationRulePtr()
Access method for default integration rule.
virtual void drawSpecial(oofegGraphicContext &gc, TimeStep *tStep)
ScalarAlgorithmType getScalarAlgo()
The element interface required by ZZNodalRecoveryModel.
const FloatArray & giveSubPatchCoordinates()
Returns local sub-patch coordinates of the receiver.
virtual void NodalAveragingRecoveryMI_computeNodalValue(FloatArray &answer, int node, InternalStateType type, TimeStep *tStep)
Computes the element value in given node.
AnalysisMode
Mode of analysis.
double & at(int i)
Coefficient access function.
virtual int giveInternalStateAtNode(FloatArray &answer, InternalStateType type, InternalStateMode mode, int node, TimeStep *tStep)
Returns internal state variable (like stress,strain) at node of element in Reduced form...
#define OOFEG_RAW_GEOMETRY_LAYER
EPixel getRemoteElementColor()
void clear()
Clears receiver (zero size).
EPixel getElementEdgeColor()
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
void setupRefinedElementProblem3D(Element *element, RefinedElement *refinedElement, int level, int nodeId, IntArray &localNodeIdArray, IntArray &globalNodeIdArray, HuertaErrorEstimatorInterface::SetupMode mode, TimeStep *tStep, int nodes, FloatArray **corner, FloatArray *midSide, FloatArray *midFace, FloatArray &midNode, int &localNodeId, int &localElemId, int &localBcId, int hexaSideNode[1][3], int hexaFaceNode[1][3], IntArray &controlNode, IntArray &controlDof, HuertaErrorEstimator::AnalysisMode aMode, const char *hexatype)
EPixel getRemoteElementEdgeColor()
virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
virtual double giveCoordinate(int i)
virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
Computes interpolation matrix for element unknowns.
Class implementing an array of integers.
int & at(int i)
Coefficient access function.
virtual void SPRNodalRecoveryMI_giveDofMansDeterminedByPatch(IntArray &answer, int pap)
virtual void HuertaErrorEstimatorI_setupRefinedElementProblem(RefinedElement *refinedElement, int level, int nodeId, IntArray &localNodeIdArray, IntArray &globalNodeIdArray, HuertaErrorEstimatorInterface::SetupMode sMode, TimeStep *tStep, int &localNodeId, int &localElemId, int &localBcId, IntArray &controlNode, IntArray &controlDof, HuertaErrorEstimator::AnalysisMode aMode)
#define OOFEG_DEFORMED_GEOMETRY_LAYER
LTRSpace(int n, Domain *d)
Abstract base class representing integration rule.
int getInternalVarsDefGeoFlag()
EPixel getDeformedElementColor()
Class representing a general abstraction for finite element interpolation class.
InternalStateType giveIntVarType()
#define OOFEG_CRACK_PATTERN_LAYER
EPixel getActiveCrackColor()
The element interface corresponding to ZZErrorEstimator.
The element interface corresponding to HuertaErrorEstimator.
#define OOFEG_RAW_GEOMETRY_WIDTH
virtual int SPRNodalRecoveryMI_giveNumberOfIP()
UnknownType
Type representing particular unknown (its physical meaning).
StructuralCrossSection * giveStructuralCrossSection()
Helper function which returns the structural cross-section for the element.
#define OOFEG_CRACK_PATTERN_WIDTH
SetupMode
Mode for problem setup.
virtual double giveUpdatedCoordinate(int ic, TimeStep *tStep, double scale=1.)
Returns updated ic-th coordinate of receiver.
EPixel getCrackPatternColor()
double at(int i, int j) const
Coefficient access function.
void resize(int n)
Checks size of receiver towards requested bounds.
virtual Interface * giveInterface(InterfaceType it)
Interface requesting service.
virtual SPRPatchType SPRNodalRecoveryMI_givePatchType()
int numberOfGaussPoints
Number of integration points as specified by nip.
InternalStateMode giveIntVarMode()
Class representing vector of real numbers.
elementParallelMode giveParallelMode() const
Return elementParallelMode of receiver.
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.
void resize(int rows, int cols)
Checks size of receiver towards requested bounds.
virtual int getRequiredNumberOfIntegrationPoints(integrationDomain dType, int approxOrder)
Abstract service.
int giveNumberOfIntegrationPoints() const
Returns number of integration points of receiver.
virtual FEInterpolation * giveInterpolation() const
#define OOFEG_DEFORMED_GEOMETRY_WIDTH
virtual void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType)
The spatial localizer element interface associated to spatial localizer.
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.
Element in active domain is only mirror of some remote element.
virtual FloatArray * giveCoordinates()
void zero()
Zeroes all coefficient of receiver.
InterfaceType
Enumerative type, used to identify interface type.
virtual double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
void updateFringeTableMinMax(double *s, int size)
the oofem namespace is to define a context or scope in which all oofem names are defined.
virtual void SPRNodalRecoveryMI_giveSPRAssemblyPoints(IntArray &pap)
virtual void drawScalar(oofegGraphicContext &gc, TimeStep *tStep)
Node * giveNode(int i) const
Returns reference to the i-th node of element.
static FEI3dTetLin interpolation
#define OOFEG_VARPLOT_PATTERN_LAYER
Class representing integration point in finite element program.
Class representing solution step.
int numberOfDofMans
Number of dofmanagers.
virtual void HuertaErrorEstimatorI_computeNmatrixAt(GaussPoint *gp, FloatMatrix &answer)
Class representing Gaussian-quadrature integration rule.
void resize(int s)
Resizes receiver towards requested size.