OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
libeam2dnl.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 #include "../sm/Elements/Beams/libeam2dnl.h"
36 #include "../sm/Materials/structuralms.h"
37 #include "node.h"
38 #include "material.h"
39 #include "crosssection.h"
40 #include "gausspoint.h"
41 #include "gaussintegrationrule.h"
42 #include "floatmatrix.h"
43 #include "intarray.h"
44 #include "floatarray.h"
45 #include "mathfem.h"
46 #include "classfactory.h"
47 
48 #ifdef __OOFEG
49  #include "oofeggraphiccontext.h"
50 #endif
51 
52 namespace oofem {
53 REGISTER_Element(LIBeam2dNL);
54 
56 {
57  numberOfDofMans = 2;
58  length = 0.;
59  pitch = 10.; // a dummy value
60 }
61 
62 Interface *
64 {
65  if ( interface == LayeredCrossSectionInterfaceType ) {
66  return static_cast< LayeredCrossSectionInterface * >(this);
67  }
68 
69  return NULL;
70 }
71 
72 
73 void
75 // Returns the strain matrix of the receiver.
76 {
77  double l, ksi, n1, n2, n1x, n2x, n3x;
78 
79  l = this->computeLength();
80  ksi = gp->giveNaturalCoordinate(1);
81 
82  n1 = 0.5 * ( 1.0 - ksi );
83  n2 = 0.5 * ( 1.0 + ksi );
84  n1x = -1.0 / l;
85  n2x = 1.0 / l;
86  n3x = -4.0 * ksi / l;
87 
88  answer.resize(3, 6);
89  answer.zero();
90 
91  // eps_x
92  answer.at(1, 1) = n1x;
93  answer.at(1, 4) = n2x;
94  // kappa
95  answer.at(2, 3) = n1x;
96  answer.at(2, 6) = n2x;
97  // gamma_xz
98  answer.at(3, 2) = n1x;
99  answer.at(3, 3) = n1 - n3x * l / 8.0;
100  answer.at(3, 5) = n2x;
101  answer.at(3, 6) = n2 + n3x * l / 8.0;
102 }
103 
104 
105 void
107 {
108  //
109  // Returns nonlinear part of geometrical equations of the receiver at gp.
110  //
111  // Returns A matrix (see Bittnar & Sejnoha Num. Met. Mech. part II, chap 9)
112  double l, l8, ll88, ksi, n1x, n2x, n3x;
113 
114  l = this->computeLength();
115  ksi = gp->giveNaturalCoordinate(1);
116 
117  //n1 = 0.5*(1.0 - ksi);
118  //n2 = 0.5*(1.0 + ksi);
119  n1x = -1.0 / l;
120  n2x = 1.0 / l;
121  n3x = -4.0 * ksi / l;
122  l8 = l / 8.;
123  ll88 = l8 * l8;
124 
125  answer.resize(6, 6);
126  answer.zero();
127  /*
128  * if (i==1){
129  * answer.at(1,1)= n1x*n1x;
130  * answer.at(1,4)= n1x*n2x;
131  * answer.at(2,2)= n1x*n1x;
132  * answer.at(2,3)=-n1x*n3x*l8;
133  * answer.at(2,5)= n1x*n2x;
134  * answer.at(2,6)= n1x*n3x*l8;
135  * answer.at(3,2)=-n3x*l8*n1x;
136  * answer.at(3,3)= n3x*n3x*ll88;
137  * answer.at(3,5)=-n3x*l8*n2x;
138  * answer.at(3,6)=-n3x*n3x*ll88;
139  *
140  * answer.at(4,1)= n2x*n1x;
141  * answer.at(4,4)= n2x*n2x;
142  * answer.at(5,2)= n2x*n1x;
143  * answer.at(5,3)=-n2x*n3x*l8;
144  * answer.at(5,5)= n2x*n2x;
145  * answer.at(5,6)= n2x*n3x*l8;
146  * answer.at(6,2)= n3x*l8*n1x;
147  * answer.at(6,3)=-n3x*n3x*ll88;
148  * answer.at(6,5)= n3x*l8*n2x;
149  * answer.at(6,2)= n3x*n3x*ll88;
150  *
151  * }
152  *
153  * if (i==2){
154  * answer.at(1,2)=n1x*n1x;
155  * answer.at(1,5)=n1x*n2x;
156  * answer.at(2,1)=n1x*n1x;
157  * answer.at(2,4)=n1x*n2x;
158  * answer.at(4,2)=n2x*n1x;
159  * answer.at(4,5)=n2x*n2x;
160  * answer.at(5,1)=n2x*n1x;
161  * answer.at(5,4)=n2x*n2x;
162  *
163  * }
164  *
165  * if (i==3){
166  * answer.at(1,3)=n1x*n1;
167  * answer.at(1,6)=n1x*n2;
168  * answer.at(3,1)=n1*n1x;
169  * answer.at(3,4)=n1*n2x;
170  * answer.at(4,3)=n2x*n1;
171  * answer.at(4,6)=n2x*n2;
172  * answer.at(6,1)=n2*n1x;
173  * answer.at(6,4)=n2*n2x;
174  * }
175  */
176 
177  /* Working
178  *
179  * // large deflection
180  * // based on von Karman equation
181  *
182  * if (i==1) {
183  * answer.at(2,2)= n1x*n1x;
184  * answer.at(2,3)=-n1x*n3x*l8;
185  * answer.at(2,5)= n1x*n2x;
186  * answer.at(2,6)= n1x*n3x*l8;
187  * answer.at(3,2)=-n3x*l8*n1x;
188  * answer.at(3,3)= n3x*n3x*ll88;
189  * answer.at(3,5)=-n3x*l8*n2x;
190  * answer.at(3,6)=-n3x*n3x*ll88;
191  *
192  * answer.at(5,2)= n2x*n1x;
193  * answer.at(5,3)=-n2x*n3x*l8;
194  * answer.at(5,5)= n2x*n2x;
195  * answer.at(5,6)= n2x*n3x*l8;
196  * answer.at(6,2)= n3x*l8*n1x;
197  * answer.at(6,3)=-n3x*n3x*ll88;
198  * answer.at(6,5)= n3x*l8*n2x;
199  * answer.at(6,2)= n3x*n3x*ll88;
200  *
201  * }
202  */
203 
204  if ( i == 1 ) {
205  // answer.at(1,1)= n1x*n1x;
206  // answer.at(1,4)= n1x*n2x;
207  answer.at(2, 2) = n1x * n1x;
208  answer.at(2, 3) = -n1x * n3x * l8;
209  answer.at(2, 5) = n1x * n2x;
210  answer.at(2, 6) = n1x * n3x * l8;
211  answer.at(3, 2) = -n3x * l8 * n1x;
212  answer.at(3, 3) = n3x * n3x * ll88;
213  answer.at(3, 5) = -n3x * l8 * n2x;
214  answer.at(3, 6) = -n3x * n3x * ll88;
215 
216  // answer.at(4,1)= n2x*n1x;
217  // answer.at(4,4)= n2x*n2x;
218  answer.at(5, 2) = n2x * n1x;
219  answer.at(5, 3) = -n2x * n3x * l8;
220  answer.at(5, 5) = n2x * n2x;
221  answer.at(5, 6) = n2x * n3x * l8;
222  answer.at(6, 2) = n3x * l8 * n1x;
223  answer.at(6, 3) = -n3x * n3x * ll88;
224  answer.at(6, 5) = n3x * l8 * n2x;
225  answer.at(6, 2) = n3x * n3x * ll88;
226  }
227 
228  /*
229  * if (i==2){
230  * answer.at(1,2)=n1x*n1x;
231  * answer.at(1,5)=n1x*n2x;
232  * answer.at(2,1)=n1x*n1x;
233  * answer.at(2,4)=n1x*n2x;
234  * answer.at(4,2)=n2x*n1x;
235  * answer.at(4,5)=n2x*n2x;
236  * answer.at(5,1)=n2x*n1x;
237  * answer.at(5,4)=n2x*n2x;
238  *
239  * }
240  *
241  * if (i==3){
242  * answer.at(1,3)=n1x*n1;
243  * answer.at(1,6)=n1x*n2;
244  * answer.at(3,1)=n1*n1x;
245  * answer.at(3,4)=n1*n2x;
246  * answer.at(4,3)=n2x*n1;
247  * answer.at(4,6)=n2x*n2;
248  * answer.at(6,1)=n2*n1x;
249  * answer.at(6,4)=n2*n2x;
250  * }
251  */
252 }
253 
254 
255 void
257 {
258  double dV;
259  FloatArray stress;
260  FloatMatrix A;
261 
262  answer.resize(6, 6);
263  answer.zero();
264 
265  // assemble initial stress matrix
266  for ( GaussPoint *gp: *this->giveDefaultIntegrationRulePtr() ) {
267  dV = this->computeVolumeAround(gp);
268  stress = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() )->giveStressVector();
269  if ( stress.giveSize() ) {
270  for ( int j = 1; j <= stress.giveSize(); j++ ) {
271  // loop over each component of strain vector
272  this->computeNLBMatrixAt(A, gp, j);
273  if ( A.isNotEmpty() ) {
274  A.times(stress.at(j) * dV);
275  answer.add(A);
276  }
277  }
278  }
279  }
280 }
281 
282 
284 // Sets up the array of Gauss Points of the receiver.
285 {
286  if ( integrationRulesArray.size() == 0 ) {
287  integrationRulesArray.resize( 1 );
288  integrationRulesArray [ 0 ].reset( new GaussIntegrationRule(1, this, 1, 2) );
290  }
291 }
292 
293 
294 void
296 // Returns the lumped mass matrix of the receiver. This expression is
297 // valid in both local and global axes.
298 {
299  GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
300  double density = this->giveStructuralCrossSection()->give('d', gp);
301  double halfMass = density * this->giveCrossSection()->give(CS_Area, gp) * this->computeLength() / 2.;
302  answer.resize(6, 6);
303  answer.zero();
304  answer.at(1, 1) = halfMass;
305  answer.at(2, 2) = halfMass;
306  answer.at(4, 4) = halfMass;
307  answer.at(5, 5) = halfMass;
308 }
309 
310 
311 void
313 // Returns the displacement interpolation matrix {N} of the receiver, eva-
314 // luated at gp.
315 {
316  double ksi, n1, n2, n3;
317  double l = this->computeLength();
318 
319  ksi = iLocCoord.at(1);
320  n1 = ( 1. - ksi ) * 0.5;
321  n2 = ( 1. + ksi ) * 0.5;
322  n3 = ( 1. - ksi * ksi );
323 
324  answer.resize(3, 6);
325  answer.zero();
326 
327  // us
328  answer.at(1, 1) = n1;
329  answer.at(1, 4) = n2;
330  // w
331  answer.at(2, 2) = n1;
332  answer.at(2, 3) = -n3 * l / 8.;
333  answer.at(2, 5) = n2;
334  answer.at(2, 6) = n3 * l / 8.;
335  // fi_y
336  answer.at(3, 3) = n1;
337  answer.at(3, 6) = n2;
338 }
339 
340 
341 bool
343 {
344  double sine, cosine;
345 
346  sine = sin( this->givePitch() );
347  cosine = cos(pitch);
348 
349  answer.resize(6, 6);
350  answer.zero();
351 
352  answer.at(1, 1) = cosine;
353  answer.at(1, 2) = sine;
354  answer.at(2, 1) = -sine;
355  answer.at(2, 2) = cosine;
356  answer.at(3, 3) = 1.;
357  answer.at(4, 4) = cosine;
358  answer.at(4, 5) = sine;
359  answer.at(5, 4) = -sine;
360  answer.at(5, 5) = cosine;
361  answer.at(6, 6) = 1.;
362 
363  return true;
364 }
365 
366 
368 // Returns the length of the receiver. This method is valid only if 1
369 // Gauss point is used.
370 {
371  double weight = gp->giveWeight();
372  return weight * 0.5 * this->computeLength();
373 }
374 
375 
376 int
378 {
381  if ( type == IST_BeamForceMomentTensor ) {
382  answer = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() )->giveStressVector();
383  return 1;
384  } else if ( type == IST_BeamStrainCurvatureTensor ) {
385  answer = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() )->giveStrainVector();
386  return 1;
387  } else {
388  return StructuralElement :: giveIPValue(answer, gp, type, tStep);
389  }
390 }
391 
392 
393 void
395 {
396  this->giveStructuralCrossSection()->giveGeneralizedStress_Beam2d(answer, gp, strain, tStep);
397 }
398 
399 
400 void
402 {
403  this->giveStructuralCrossSection()->give2dBeamStiffMtrx(answer, rMode, gp, tStep);
404 }
405 
406 
407 void
408 LIBeam2dNL :: computeStrainVectorInLayer(FloatArray &answer, const FloatArray &masterGpStrain, GaussPoint *masterGp, GaussPoint *slaveGp, TimeStep *tStep)
409 //
410 // returns full 3d strain vector of given layer (whose z-coordinate from center-line is
411 // stored in slaveGp) for given tStep
412 //
413 {
414  double layerZeta, layerZCoord, top, bottom;
415 
416  top = this->giveCrossSection()->give(CS_TopZCoord, masterGp);
417  bottom = this->giveCrossSection()->give(CS_BottomZCoord, masterGp);
418  layerZeta = slaveGp->giveNaturalCoordinate(3);
419  layerZCoord = 0.5 * ( ( 1. - layerZeta ) * bottom + ( 1. + layerZeta ) * top );
420 
421  answer.resize(2); // {Exx,GMzx}
422 
423  answer.at(1) = masterGpStrain.at(1) + masterGpStrain.at(2) * layerZCoord;
424  answer.at(2) = masterGpStrain.at(3);
425 }
426 
427 
428 void
430 {
431  answer = {D_u, D_w, R_v};
432 }
433 
434 
435 int
437 {
438  double ksi, n1, n2;
439 
440  ksi = lcoords.at(1);
441  n1 = ( 1. - ksi ) * 0.5;
442  n2 = ( 1. + ksi ) * 0.5;
443 
444  answer.resize(3);
445  answer.at(1) = n1 * this->giveNode(1)->giveCoordinate(1) + n2 *this->giveNode(2)->giveCoordinate(1);
446  answer.at(3) = n1 * this->giveNode(1)->giveCoordinate(3) + n2 *this->giveNode(2)->giveCoordinate(3);
447 
448  return 1;
449 }
450 
451 
453 // Returns the length of the receiver.
454 {
455  double dx, dy;
456  Node *nodeA, *nodeB;
457 
458  if ( length == 0. ) {
459  nodeA = this->giveNode(1);
460  nodeB = this->giveNode(2);
461  dx = nodeB->giveCoordinate(1) - nodeA->giveCoordinate(1);
462  dy = nodeB->giveCoordinate(3) - nodeA->giveCoordinate(3);
463  length = sqrt(dx * dx + dy * dy);
464  }
465 
466  return length;
467 }
468 
469 
471 // Returns the pitch of the receiver.
472 {
473  double xA, xB, yA, yB;
474  Node *nodeA, *nodeB;
475 
476  if ( pitch == 10. ) { // 10. : dummy initialization value
477  nodeA = this->giveNode(1);
478  nodeB = this->giveNode(2);
479  xA = nodeA->giveCoordinate(1);
480  xB = nodeB->giveCoordinate(1);
481  yA = nodeA->giveCoordinate(3);
482  yB = nodeB->giveCoordinate(3);
483  pitch = atan2(yB - yA, xB - xA);
484  }
485 
486  return pitch;
487 }
488 
489 
492 {
494 }
495 
496 
497 void
499 {
500  /*
501  * provides dof mapping of local edge dofs (only nonzero are taken into account)
502  * to global element dofs
503  */
504 
505  if ( iEdge != 1 ) {
506  OOFEM_ERROR("wrong edge number");
507  }
508 
509 
510  answer.resize(6);
511  answer.at(1) = 1;
512  answer.at(2) = 2;
513  answer.at(3) = 3;
514  answer.at(4) = 4;
515  answer.at(5) = 5;
516  answer.at(6) = 6;
517 }
518 
519 double
521 {
522  if ( iEdge != 1 ) { // edge between nodes 1 2
523  OOFEM_ERROR("wrong egde number");
524  }
525 
526  double weight = gp->giveWeight();
527  return 0.5 * this->computeLength() * weight;
528 }
529 
530 
531 int
533 {
534  /*
535  * Returns transformation matrix from global coordinate system to local
536  * element coordinate system for element load vector components.
537  * If no transformation is necessary, answer is empty matrix (default);
538  */
539 
540  // f(elemLocal) = T * f (global)
541  double sine, cosine;
542 
543  answer.resize(3, 3);
544  answer.zero();
545 
546  sine = sin( this->givePitch() );
547  cosine = cos(pitch);
548 
549  answer.at(1, 1) = cosine;
550  answer.at(1, 2) = -sine;
551  answer.at(2, 1) = sine;
552  answer.at(2, 2) = cosine;
553  answer.at(3, 3) = 1.0;
554 
555  return 1;
556 }
557 
558 int
560 {
561  // returns transformation matrix from
562  // edge local coordinate system
563  // to element local c.s
564  // (same as global c.s in this case)
565  //
566  // i.e. f(element local) = T * f(edge local)
567  //
568  answer.clear();
569  return 0;
570 }
571 
572 
573 void
575 {
576  FloatArray lc(1);
577  StructuralElement :: computeBodyLoadVectorAt(answer, load, tStep, mode);
578  answer.times( this->giveCrossSection()->give(CS_Area, lc, this) );
579 }
580 
581 
582 #ifdef __OOFEG
584 {
585  GraphicObj *go;
586 
587  if ( !gc.testElementGraphicActivity(this) ) {
588  return;
589  }
590 
591  // if (!go) { // create new one
592  WCRec p [ 2 ]; /* poin */
593  EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
594  EASValsSetColor( gc.getElementColor() );
595  EASValsSetLayer(OOFEG_RAW_GEOMETRY_LAYER);
596  p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1);
597  p [ 0 ].y = 0.;
598  p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveCoordinate(3);
599  p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1);
600  p [ 1 ].y = 0.;
601  p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveCoordinate(3);
602  go = CreateLine3D(p);
603  EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
604  EGAttachObject(go, ( EObjectP ) this);
605  EMAddGraphicsToModel(ESIModel(), go);
606 }
607 
608 
610 {
611  GraphicObj *go;
612  if ( !gc.testElementGraphicActivity(this) ) {
613  return;
614  }
615 
616  double defScale = gc.getDefScale();
617  // if (!go) { // create new one
618  WCRec p [ 2 ]; /* poin */
619  EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
620  EASValsSetColor( gc.getDeformedElementColor() );
621  EASValsSetLayer(OOFEG_DEFORMED_GEOMETRY_LAYER);
622  p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(1, tStep, defScale);
623  p [ 0 ].y = 0.;
624  p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(3, tStep, defScale);
625 
626  p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(1, tStep, defScale);
627  p [ 1 ].y = 0.;
628  p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(3, tStep, defScale);
629  go = CreateLine3D(p);
630  EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
631  EMAddGraphicsToModel(ESIModel(), go);
632 }
633 #endif
634 } // end namespace oofem
CrossSection * giveCrossSection()
Definition: element.C:495
virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes lumped mass matrix of receiver.
Definition: libeam2dnl.C:295
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.
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
virtual bool computeGtoLRotationMatrix(FloatMatrix &answer)
Returns transformation matrix from global c.s.
Definition: libeam2dnl.C:342
virtual void giveGeneralizedStress_Beam2d(FloatArray &answer, GaussPoint *gp, const FloatArray &generalizedStrain, TimeStep *tStep)=0
Computes the generalized stress vector for given strain and integration point.
Class and object Domain.
Definition: domain.h:115
virtual IntegrationRule * giveDefaultIntegrationRulePtr()
Access method for default integration rule.
Definition: element.h:822
Bottom z coordinate.
Definition: crosssection.h:72
Abstract base class for "structural" finite elements with geometrical nonlinearities.
double & at(int i)
Coefficient access function.
Definition: floatarray.h:131
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
Definition: valuemodetype.h:78
#define OOFEG_RAW_GEOMETRY_LAYER
This class implements a structural material status information.
Definition: structuralms.h:65
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
virtual int computeLoadLEToLRotationMatrix(FloatMatrix &, int, GaussPoint *)
Returns transformation matrix from local edge c.s to element local coordinate system of load vector c...
Definition: libeam2dnl.C:559
Top z coordinate.
Definition: crosssection.h:71
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: libeam2dnl.C:491
virtual double giveCoordinate(int i)
Definition: node.C:82
Class implementing an array of integers.
Definition: intarray.h:61
int & at(int i)
Coefficient access function.
Definition: intarray.h:103
MatResponseMode
Describes the character of characteristic material matrix.
virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
Definition: libeam2dnl.C:583
virtual int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords)
Computes the global coordinates from given element&#39;s local coordinates.
Definition: libeam2dnl.C:436
#define OOFEG_DEFORMED_GEOMETRY_LAYER
double givePitch()
Definition: libeam2dnl.C:470
virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
Computes constitutive matrix of receiver.
Definition: libeam2dnl.C:401
virtual void giveEdgeDofMapping(IntArray &answer, int) const
Assembles edge dof mapping mask, which provides mapping between edge local DOFs and "global" element ...
Definition: libeam2dnl.C:498
bool isNotEmpty() const
Tests for empty matrix.
Definition: floatmatrix.h:162
double giveNaturalCoordinate(int i) const
Returns i-th natural element coordinate of receiver.
Definition: gausspoint.h:136
virtual void computeBodyLoadVectorAt(FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode)
Computes the load vector due to body load acting on receiver, at given time step. ...
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
Definition: libeam2dnl.C:283
virtual void computeBodyLoadVectorAt(FloatArray &answer, Load *load, TimeStep *tStep, ValueModeType mode)
Computes the load vector due to body load acting on receiver, at given time step. ...
Definition: libeam2dnl.C:574
#define OOFEM_ERROR(...)
Definition: error.h:61
virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
Computes interpolation matrix for element unknowns.
Definition: libeam2dnl.C:312
REGISTER_Element(LSpace)
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
Definition: libeam2dnl.C:377
#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
StructuralCrossSection * giveStructuralCrossSection()
Helper function which returns the structural cross-section for the element.
void times(double f)
Multiplies receiver by factor f.
Definition: floatmatrix.C:1594
virtual void computeStrainVectorInLayer(FloatArray &answer, const FloatArray &masterGpStrain, GaussPoint *masterGp, GaussPoint *slaveGp, TimeStep *tStep)
Computes full 3D strain vector in element layer.
Definition: libeam2dnl.C:408
virtual double giveUpdatedCoordinate(int ic, TimeStep *tStep, double scale=1.)
Returns updated ic-th coordinate of receiver.
Definition: node.C:245
virtual Interface * giveInterface(InterfaceType it)
Interface requesting service.
Definition: libeam2dnl.C:63
double at(int i, int j) const
Coefficient access function.
Definition: floatmatrix.h:176
void resize(int n)
Checks size of receiver towards requested bounds.
Definition: intarray.C:124
void computeNLBMatrixAt(FloatMatrix &answer, GaussPoint *gp, int)
Definition: libeam2dnl.C:106
virtual int setupIntegrationPoints(IntegrationRule &irule, int npoints, Element *element)
Sets up integration rule for the given element.
Definition: crosssection.C:54
LIBeam2dNL(int n, Domain *d)
Definition: libeam2dnl.C:55
Class representing vector of real numbers.
Definition: floatarray.h:82
Implementation of matrix containing floating point numbers.
Definition: floatmatrix.h:94
virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int=1, int=ALL_STRAINS)
Computes the geometrical matrix of receiver in given integration point.
Definition: libeam2dnl.C:74
IRResultType
Type defining the return values of InputRecord reading operations.
Definition: irresulttype.h:47
virtual double computeLength()
Computes the length (zero for all but 1D geometries)
Definition: libeam2dnl.C:452
virtual double give(CrossSectionProperty a, GaussPoint *gp)
Returns the value of cross section property at given point.
Definition: crosssection.C:151
IntegrationPointStatus * giveMaterialStatus()
Returns reference to associated material status (NULL if not defined).
Definition: gausspoint.h:205
virtual void computeInitialStressMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes the initial stiffness matrix of receiver.
Definition: libeam2dnl.C:256
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
void add(const FloatMatrix &a)
Adds matrix to the receiver.
Definition: floatmatrix.C:1023
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.
Definition: libeam2dnl.C:394
Class Interface.
Definition: interface.h:82
#define OOFEG_DEFORMED_GEOMETRY_WIDTH
void times(double s)
Multiplies receiver with scalar.
Definition: floatarray.C:818
virtual void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType)
Definition: libeam2dnl.C:609
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
virtual double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
Definition: libeam2dnl.C:367
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
void zero()
Zeroes all coefficient of receiver.
Definition: floatmatrix.C:1326
InterfaceType
Enumerative type, used to identify interface type.
Definition: interfacetype.h:43
virtual double computeEdgeVolumeAround(GaussPoint *, int)
Computes volume related to integration point on local edge.
Definition: libeam2dnl.C:520
Load is base abstract class for all loads.
Definition: load.h:61
The element interface required by LayeredCrossSection.
virtual void giveDofManDofIDMask(int inode, IntArray &) const
Returns dofmanager dof mask for node.
Definition: libeam2dnl.C:429
int giveSize() const
Returns the size of receiver.
Definition: floatarray.h:218
virtual int computeLoadGToLRotationMtrx(FloatMatrix &answer)
Returns transformation matrix from global coordinate system to local element coordinate system for el...
Definition: libeam2dnl.C:532
the oofem namespace is to define a context or scope in which all oofem names are defined.
void clear()
Sets size of receiver to be an empty matrix. It will have zero rows and zero columns size...
Definition: floatmatrix.h:516
Class implementing node in finite element mesh.
Definition: node.h:87
Node * giveNode(int i) const
Returns reference to the i-th node of element.
Definition: element.h:610
virtual void give2dBeamStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)=0
Computes the stiffness matrix for 2d beams.
Class representing integration point in finite element program.
Definition: gausspoint.h:93
Class representing solution step.
Definition: timestep.h:80
int numberOfDofMans
Number of dofmanagers.
Definition: element.h:149
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