OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
interfaceelem3dtrlin.C
Go to the documentation of this file.
1 /*
2  *
3  * ##### ##### ###### ###### ### ###
4  * ## ## ## ## ## ## ## ### ##
5  * ## ## ## ## #### #### ## # ##
6  * ## ## ## ## ## ## ## ##
7  * ## ## ## ## ## ## ## ##
8  * ##### ##### ## ###### ## ##
9  *
10  *
11  * OOFEM : Object Oriented Finite Element Code
12  *
13  * Copyright (C) 1993 - 2013 Borek Patzak
14  *
15  *
16  *
17  * Czech Technical University, Faculty of Civil Engineering,
18  * Department of Structural Mechanics, 166 29 Prague, Czech Republic
19  *
20  * This library is free software; you can redistribute it and/or
21  * modify it under the terms of the GNU Lesser General Public
22  * License as published by the Free Software Foundation; either
23  * version 2.1 of the License, or (at your option) any later version.
24  *
25  * This program is distributed in the hope that it will be useful,
26  * but WITHOUT ANY WARRANTY; without even the implied warranty of
27  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
28  * Lesser General Public License for more details.
29  *
30  * You should have received a copy of the GNU Lesser General Public
31  * License along with this library; if not, write to the Free Software
32  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
33  */
34 
35 
36 #include "interfaceelem3dtrlin.h"
37 #include "fei2dtrlin.h"
38 #include "node.h"
39 #include "gausspoint.h"
40 #include "gaussintegrationrule.h"
41 #include "floatmatrix.h"
42 #include "floatarray.h"
43 #include "intarray.h"
44 #include "mathfem.h"
45 #include "../sm/CrossSections/structuralinterfacecrosssection.h"
46 #include "classfactory.h"
47 
48 #ifdef __OOFEG
49  #include "oofeggraphiccontext.h"
50 
51  #include <Emarkwd3d.h>
52 #endif
53 
54 namespace oofem {
55 REGISTER_Element(InterfaceElement3dTrLin);
56 
58 
60  StructuralElement(n, aDomain)
61 {
62  numberOfDofMans = 6;
63 }
64 
65 
66 void
68 //
69 // Returns linear part of geometrical equations of the receiver at gp.
70 // Returns the linear part of the B matrix
71 //
72 {
73  FloatArray n;
75 
76  answer.resize(3, 18);
77  answer.zero();
78 
79  answer.at(1, 10) = answer.at(2, 11) = answer.at(3, 12) = n.at(1);
80  answer.at(1, 1) = answer.at(2, 2) = answer.at(3, 3) = -n.at(1);
81 
82  answer.at(1, 13) = answer.at(2, 14) = answer.at(3, 15) = n.at(2);
83  answer.at(1, 4) = answer.at(2, 5) = answer.at(3, 6) = -n.at(2);
84 
85  answer.at(1, 16) = answer.at(2, 17) = answer.at(3, 18) = n.at(3);
86  answer.at(1, 7) = answer.at(2, 8) = answer.at(3, 9) = -n.at(3);
87 }
88 
89 
90 void
92 // Sets up the array of Gauss Points of the receiver.
93 {
94  if ( integrationRulesArray.size() == 0 ) {
95  integrationRulesArray.resize( 1 );
96  //integrationRulesArray[0].reset( new LobattoIntegrationRule (1,domain, 1, 2) );
97  integrationRulesArray [ 0 ].reset( new GaussIntegrationRule(1, this, 1, 3) );
98  integrationRulesArray [ 0 ]->SetUpPointsOnTriangle(4, _3dInterface);
99  }
100 }
101 
102 
103 int
105 {
106  FloatArray n;
107 
108  this->interpolation.evalN( n, lcoords, FEIElementGeometryWrapper(this) );
109 
110  answer.resize(3);
111  answer.zero();
112  for ( int i = 1; i <= 3; i++ ) {
113  answer.at(1) += n.at(i) * this->giveNode(i)->giveCoordinate(1);
114  answer.at(2) += n.at(i) * this->giveNode(i)->giveCoordinate(2);
115  answer.at(3) += n.at(i) * this->giveNode(i)->giveCoordinate(3);
116  }
117 
118  return 1;
119 }
120 
121 
122 bool
124 {
125  OOFEM_ERROR("Not implemented");
126  return false;
127 }
128 
129 
130 double
132 // Returns the length of the receiver. This method is valid only if 1
133 // Gauss point is used.
134 {
135  double determinant, weight, thickness, volume;
136  // first compute local nodal coordinates in element plane
137  std::vector< FloatArray > lncp(3);
138  FloatMatrix lcs(3, 3);
139  this->computeLCS(lcs);
140  for ( int i = 1; i <= 3; i++ ) {
141  lncp[ i - 1 ].beProductOf(lcs, *this->giveNode(i)->giveCoordinates());
142  }
143 
145  weight = gp->giveWeight();
146  thickness = this->giveCrossSection()->give(CS_Thickness, gp);
147  volume = determinant * weight * thickness;
148 
149  return volume;
150 }
151 
152 
153 void
155 {
156  static_cast< StructuralInterfaceCrossSection* >(this->giveCrossSection())->giveEngTraction_3d(answer, gp, strain, tStep);
157 }
158 
159 
160 void
162 {
163  static_cast< StructuralInterfaceCrossSection* >(this->giveCrossSection())->give3dStiffnessMatrix_Eng(answer, rMode, gp, tStep);
164 }
165 
166 
169 {
171 }
172 
173 
174 void
176 {
177  answer = {D_u, D_v, D_w};
178 }
179 
180 
181 void
183 {
184  // computes local coordinate system unit vectors (expressed in global cs)
185  // unit vectors are stored rowwise
186  FloatArray xl(3), yl(3), zl(3), t2(3);
187 
188  // compute local x-axis xl (node(2)-node(1))
189  xl.at(1) = this->giveNode(2)->giveCoordinate(1) - this->giveNode(1)->giveCoordinate(1);
190  xl.at(2) = this->giveNode(2)->giveCoordinate(2) - this->giveNode(1)->giveCoordinate(2);
191  xl.at(3) = this->giveNode(2)->giveCoordinate(3) - this->giveNode(1)->giveCoordinate(3);
192 
193  xl.normalize();
194 
195  // compute another in-plane tangent vector t2 (node(3)-node(1))
196  t2.at(1) = this->giveNode(3)->giveCoordinate(1) - this->giveNode(1)->giveCoordinate(1);
197  t2.at(2) = this->giveNode(3)->giveCoordinate(2) - this->giveNode(1)->giveCoordinate(2);
198  t2.at(3) = this->giveNode(3)->giveCoordinate(3) - this->giveNode(1)->giveCoordinate(3);
199 
200  // compute local z axis as product of xl and t2
201  zl.beVectorProductOf(xl, t2);
202  zl.normalize();
203 
204  // compute local y axis as product of zl x xl
205  yl.beVectorProductOf(zl, xl);
206 
207  answer.resize(3, 3);
208  for ( int i = 1; i <= 3; i++ ) {
209  answer.at(1, i) = xl.at(i);
210  answer.at(2, i) = yl.at(i);
211  answer.at(3, i) = zl.at(i);
212  }
213 }
214 
215 
216 bool
218 {
219  // planar geometry is assumed
220  FloatMatrix lcs(3, 3);
221  this->computeLCS(lcs);
222 
223  answer.resize(18, 18);
224  for ( int i = 0; i < 6; i++ ) {
225  for ( int j = 1; j <= 3; j++ ) {
226  answer.at(i * 3 + 1, i * 3 + j) = lcs.at(3, j);
227  answer.at(i * 3 + 2, i * 3 + j) = lcs.at(1, j);
228  answer.at(i * 3 + 3, i * 3 + j) = lcs.at(2, j);
229  }
230  }
231 
232  return 1;
233 }
234 
235 
236 #ifdef __OOFEG
238 {
239  GraphicObj *go;
240  // if (!go) { // create new one
241  WCRec p [ 3 ]; /* triangle */
242  if ( !gc.testElementGraphicActivity(this) ) {
243  return;
244  }
245 
246  EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
247  EASValsSetColor( gc.getElementColor() );
248  EASValsSetEdgeColor( gc.getElementEdgeColor() );
249  EASValsSetEdgeFlag(true);
250  EASValsSetLayer(OOFEG_RAW_GEOMETRY_LAYER);
251  p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1);
252  p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveCoordinate(2);
253  p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveCoordinate(3);
254  p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1);
255  p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveCoordinate(2);
256  p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveCoordinate(3);
257  p [ 2 ].x = ( FPNum ) this->giveNode(3)->giveCoordinate(1);
258  p [ 2 ].y = ( FPNum ) this->giveNode(3)->giveCoordinate(2);
259  p [ 2 ].z = ( FPNum ) this->giveNode(3)->giveCoordinate(3);
260 
261  go = CreateTriangle3D(p);
262  EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | EDGE_COLOR_MASK | EDGE_FLAG_MASK | LAYER_MASK, go);
263  EGAttachObject(go, ( EObjectP ) this);
264  EMAddGraphicsToModel(ESIModel(), go);
265 }
266 
267 
269 { }
270 
271 
273 { }
274 
275 #endif
276 } // end namespace oofem
CrossSection * giveCrossSection()
Definition: element.C:495
int testElementGraphicActivity(Element *)
Test if particular element passed fulfills various filtering criteria for its graphics output...
virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
Class and object Domain.
Definition: domain.h:115
virtual void evalN(FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the array of interpolation functions (shape functions) at given point.
Definition: fei2dtrlin.C:43
virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
Computes constitutive matrix of receiver.
Wrapper around cell with vertex coordinates stored in FloatArray**.
Definition: feinterpol.h:115
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.
double & at(int i)
Coefficient access function.
Definition: floatarray.h:131
#define OOFEG_RAW_GEOMETRY_LAYER
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
virtual double giveCoordinate(int i)
Definition: node.C:82
Class implementing an array of integers.
Definition: intarray.h:61
MatResponseMode
Describes the character of characteristic material matrix.
virtual void drawScalar(oofegGraphicContext &gc, TimeStep *tStep)
Abstract base class for all "structural" finite elements.
virtual void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType)
virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int=1, int=ALL_STRAINS)
Computes the geometrical matrix of receiver in given integration point.
virtual void giveDofManDofIDMask(int inode, IntArray &answer) const
Returns dofmanager dof mask for node.
#define OOFEM_ERROR(...)
Definition: error.h:61
REGISTER_Element(LSpace)
#define OOFEG_RAW_GEOMETRY_WIDTH
virtual double giveWeight()
Returns integration weight of receiver.
Definition: gausspoint.h:181
virtual double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
UnknownType
Type representing particular unknown (its physical meaning).
Definition: unknowntype.h:55
Wrapper around element definition to provide FEICellGeometry interface.
Definition: feinterpol.h:95
Base class for all structural interface cross section models.
void computeLCS(FloatMatrix &answer)
double at(int i, int j) const
Coefficient access function.
Definition: floatmatrix.h:176
Class representing vector of real numbers.
Definition: floatarray.h:82
Implementation of matrix containing floating point numbers.
Definition: floatmatrix.h:94
IRResultType
Type defining the return values of InputRecord reading operations.
Definition: irresulttype.h:47
virtual double give(CrossSectionProperty a, GaussPoint *gp)
Returns the value of cross section property at given point.
Definition: crosssection.C:151
virtual double giveTransformationJacobian(const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the determinant of the transformation.
Definition: fei2dtrlin.C:147
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
void resize(int rows, int cols)
Checks size of receiver towards requested bounds.
Definition: floatmatrix.C:1358
Class representing the general Input Record.
Definition: inputrecord.h:101
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
void zero()
Zeroes all coefficients of receiver.
Definition: floatarray.C:658
virtual bool computeLocalCoordinates(FloatArray &answer, const FloatArray &gcoords)
Computes the element local coordinates from given global coordinates.
std::vector< std::unique_ptr< IntegrationRule > > integrationRulesArray
List of integration rules of receiver (each integration rule contains associated integration points a...
Definition: element.h:170
void zero()
Zeroes all coefficient of receiver.
Definition: floatmatrix.C:1326
the oofem namespace is to define a context or scope in which all oofem names are defined.
Node * giveNode(int i) const
Returns reference to the i-th node of element.
Definition: element.h:610
virtual int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords)
Computes the global coordinates from given element&#39;s local coordinates.
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Class representing integration point in finite element program.
Definition: gausspoint.h:93
Class representing solution step.
Definition: timestep.h:80
virtual bool computeGtoLRotationMatrix(FloatMatrix &answer)
Returns transformation matrix from global c.s.
int numberOfDofMans
Number of dofmanagers.
Definition: element.h:149
const FloatArray & giveNaturalCoordinates()
Returns coordinate array of receiver.
Definition: gausspoint.h:138
Class representing Gaussian-quadrature integration rule.
void resize(int s)
Resizes receiver towards requested size.
Definition: floatarray.C:631

This page is part of the OOFEM documentation. Copyright (c) 2011 Borek Patzak
Project e-mail: info@oofem.org
Generated at Tue Jan 2 2018 20:07:29 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011