OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
interfaceelem2dlin.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 - 2014 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 #include "interfaceelem2dlin.h"
36 #include "node.h"
37 #include "gausspoint.h"
38 #include "gaussintegrationrule.h"
39 #include "floatmatrix.h"
40 #include "floatarray.h"
41 #include "intarray.h"
42 #include "mathfem.h"
43 #include "fei2dlinelin.h"
44 #include "../sm/CrossSections/structuralinterfacecrosssection.h"
45 #include "classfactory.h"
46 
47 #ifdef __OOFEG
48  #include "oofeggraphiccontext.h"
49  #include <Emarkwd3d.h>
50 #endif
51 
52 namespace oofem {
53 REGISTER_Element(InterfaceElem2dLin);
54 
55 FEI2dLineLin InterfaceElem2dLin :: interp(1, 2);
56 
57 
59  StructuralElement(n, aDomain)
60 {
61  numberOfDofMans = 4;
62  axisymmode = false;
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(2, 8);
77  answer.zero();
78 
79  answer.at(1, 2) = answer.at(2, 1) = -n.at(1);
80  answer.at(1, 4) = answer.at(2, 3) = -n.at(2);
81 
82  answer.at(1, 6) = answer.at(2, 5) = n.at(1);
83  answer.at(1, 8) = answer.at(2, 7) = n.at(2);
84 }
85 
86 
87 void
89 // Sets up the array of Gauss Points of the receiver.
90 {
91  if ( integrationRulesArray.size() == 0 ) {
92  integrationRulesArray.resize( 1 );
93  integrationRulesArray [ 0 ].reset( new GaussIntegrationRule(1, this, 1, 2) );
94  integrationRulesArray [ 0 ]->SetUpPointsOnLine(2, _2dInterface);
95  }
96 }
97 
98 
99 double
101 // Returns the length of the receiver. This method is valid only if 1
102 // Gauss point is used.
103 {
104  double r = 1.0;
105  if (this->axisymmode) {
106  FloatArray n(2);
108  r = n.at(1)*this->giveNode(1)->giveCoordinate(1) + n.at(2)*this->giveNode(2)->giveCoordinate(1);
109  }
110 
111 
113  return result * gp->giveWeight() * r;
114 }
115 
116 
117 void
119 {
120  static_cast< StructuralInterfaceCrossSection* >(this->giveCrossSection())->giveEngTraction_2d(answer, gp, strain, tStep);
121 }
122 
123 
124 void
126 {
127  static_cast< StructuralInterfaceCrossSection* >(this->giveCrossSection())->give2dStiffnessMatrix_Eng(answer, rMode, gp, tStep);
128 }
129 
130 
133 {
136 }
137 
138 
139 void
141 {
142  answer = {D_u, D_v};
143 }
144 
145 
146 bool
148 {
149  FloatArray grad(2);
150 
151  // tangent
152  grad.at(1) = this->giveNode(2)->giveCoordinate(1) - this->giveNode(1)->giveCoordinate(1);
153  grad.at(2) = this->giveNode(2)->giveCoordinate(2) - this->giveNode(1)->giveCoordinate(2);
154  grad.normalize();
155 
156  answer.resize(8, 8);
157  for ( int i = 0; i < 4; i++ ) {
158  answer.at(i * 2 + 1, i * 2 + 1) = grad.at(1);
159  answer.at(i * 2 + 1, i * 2 + 2) = grad.at(2);
160  answer.at(i * 2 + 2, i * 2 + 1) = -grad.at(2);
161  answer.at(i * 2 + 2, i * 2 + 2) = grad.at(1);
162  }
163 
164  return 1;
165 }
166 
167 
170 {
171  return & interp;
172 }
173 
174 #ifdef __OOFEG
176 {
177  GraphicObj *go;
178  // if (!go) { // create new one
179  WCRec p [ 2 ]; /* poin */
180  if ( !gc.testElementGraphicActivity(this) ) {
181  return;
182  }
183 
184  EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
185  EASValsSetColor( gc.getElementColor() );
186  EASValsSetLayer(OOFEG_RAW_GEOMETRY_LAYER);
187  p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1);
188  p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveCoordinate(2);
189  p [ 0 ].z = 0.0;
190  p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1);
191  p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveCoordinate(2);
192  p [ 1 ].z = 0.0;
193  go = CreateLine3D(p);
194  EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
195  EGAttachObject(go, ( EObjectP ) this);
196  EMAddGraphicsToModel(ESIModel(), go);
197 }
198 
199 
201 {
202  GraphicObj *go;
203  // if (!go) { // create new one
204  WCRec p [ 2 ]; /* poin */
205  if ( !gc.testElementGraphicActivity(this) ) {
206  return;
207  }
208 
209  double defScale = gc.getDefScale();
210 
211  EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
212  EASValsSetColor( gc.getDeformedElementColor() );
213  EASValsSetLayer(OOFEG_DEFORMED_GEOMETRY_LAYER + 1);
214  p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(1, tStep, defScale);
215  p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(2, tStep, defScale);
216  p [ 0 ].z = 0.0;
217  p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(1, tStep, defScale);
218  p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(2, tStep, defScale);
219  p [ 1 ].z = 0.0;
220  go = CreateLine3D(p);
221  EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
222  EMAddGraphicsToModel(ESIModel(), go);
223 
224  p [ 0 ].x = ( FPNum ) this->giveNode(3)->giveUpdatedCoordinate(1, tStep, defScale);
225  p [ 0 ].y = ( FPNum ) this->giveNode(3)->giveUpdatedCoordinate(2, tStep, defScale);
226  p [ 0 ].z = 0.0;
227  p [ 1 ].x = ( FPNum ) this->giveNode(4)->giveUpdatedCoordinate(1, tStep, defScale);
228  p [ 1 ].y = ( FPNum ) this->giveNode(4)->giveUpdatedCoordinate(2, tStep, defScale);
229  p [ 1 ].z = 0.0;
230  go = CreateLine3D(p);
231  EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
232  EMAddGraphicsToModel(ESIModel(), go);
233 }
234 
235 
237 {
238  int indx, result = 0;
239  FloatArray gcoord(3), v1;
240  WCRec p [ 1 ];
241  GraphicObj *go;
242  double val [ 1 ];
243 
244  if ( !gc.testElementGraphicActivity(this) ) {
245  return;
246  }
247 
248  if ( gc.giveIntVarMode() == ISM_recovered ) {
249  return;
250  }
251 
252  for ( GaussPoint *gp: *this->giveDefaultIntegrationRulePtr() ) {
253  result = 0;
254  result += giveIPValue(v1, gp, gc.giveIntVarType(), tStep);
255  if ( result != 1 ) {
256  continue;
257  }
258 
259  indx = gc.giveIntVarIndx();
260 
261  result += this->computeGlobalCoordinates( gcoord, gp->giveNaturalCoordinates() );
262 
263  p [ 0 ].x = ( FPNum ) gcoord.at(1);
264  p [ 0 ].y = ( FPNum ) gcoord.at(2);
265  p [ 0 ].z = 0.;
266 
267  val [ 0 ] = v1.at(indx);
268  gc.updateFringeTableMinMax(val, 1);
269  //if (val[0] > 0.) {
270 
271  EASValsSetLayer(OOFEG_VARPLOT_PATTERN_LAYER);
272  EASValsSetMType(FILLED_CIRCLE_MARKER);
273  go = CreateMarkerWD3D(p, val [ 0 ]);
274  EGWithMaskChangeAttributes(LAYER_MASK | FILL_MASK | MTYPE_MASK, go);
275  EMAddGraphicsToModel(ESIModel(), go);
276  //}
277  }
278 }
279 
280 #endif
281 } // end namespace oofem
CrossSection * giveCrossSection()
Definition: element.C:495
virtual void giveDofManDofIDMask(int inode, IntArray &answer) const
Returns dofmanager dof mask for node.
int testElementGraphicActivity(Element *)
Test if particular element passed fulfills various filtering criteria for its graphics output...
virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
virtual double giveTransformationJacobian(const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the determinant of the transformation.
Definition: fei2dlinelin.C:120
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
Class and object Domain.
Definition: domain.h:115
virtual IntegrationRule * giveDefaultIntegrationRulePtr()
Access method for default integration rule.
Definition: element.h:822
double & at(int i)
Coefficient access function.
Definition: floatarray.h:131
#define OOFEG_RAW_GEOMETRY_LAYER
virtual void evalN(FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the array of interpolation functions (shape functions) at given point.
Definition: fei2dlinelin.C:42
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int=1, int=ALL_STRAINS)
Computes the geometrical matrix of receiver in given integration point.
virtual bool hasField(InputFieldType id)=0
Returns true if record contains field identified by idString keyword.
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.
#define OOFEG_DEFORMED_GEOMETRY_LAYER
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.
InterfaceElem2dLin(int n, Domain *d)
Class representing a general abstraction for finite element interpolation class.
Definition: feinterpol.h:132
InternalStateType giveIntVarType()
Abstract base class for all "structural" finite elements.
virtual FEInterpolation * giveInterpolation() const
static FEI2dLineLin interp
REGISTER_Element(LSpace)
#define OOFEG_RAW_GEOMETRY_WIDTH
virtual double giveWeight()
Returns integration weight of receiver.
Definition: gausspoint.h:181
UnknownType
Type representing particular unknown (its physical meaning).
Definition: unknowntype.h:55
Wrapper around element definition to provide FEICellGeometry interface.
Definition: feinterpol.h:95
virtual void drawScalar(oofegGraphicContext &gc, TimeStep *tStep)
virtual double giveUpdatedCoordinate(int ic, TimeStep *tStep, double scale=1.)
Returns updated ic-th coordinate of receiver.
Definition: node.C:245
Base class for all structural interface cross section models.
double at(int i, int j) const
Coefficient access function.
Definition: floatmatrix.h:176
bool axisymmode
Flag controlling axisymmetric mode (integration over unit circumferential angle)
InternalStateMode giveIntVarMode()
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
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 bool computeGtoLRotationMatrix(FloatMatrix &answer)
Returns transformation matrix from global c.s.
#define OOFEG_DEFORMED_GEOMETRY_WIDTH
virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
Computes constitutive matrix of receiver.
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
void updateFringeTableMinMax(double *s, int size)
#define _IFT_InterfaceElem2dLin_axisymmode
virtual void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType)
the oofem namespace is to define a context or scope in which all oofem names are defined.
double normalize()
Normalizes receiver.
Definition: floatarray.C:828
Node * giveNode(int i) const
Returns reference to the i-th node of element.
Definition: element.h:610
virtual double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
#define OOFEG_VARPLOT_PATTERN_LAYER
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 int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords)
Computes the global coordinates from given element&#39;s local coordinates.
Definition: element.C:1207
int numberOfDofMans
Number of dofmanagers.
Definition: element.h:149
const FloatArray & giveNaturalCoordinates()
Returns coordinate array of receiver.
Definition: gausspoint.h:138
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
Class representing Gaussian-quadrature integration rule.

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