OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
libeam3d2.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/libeam3d2.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 "timestep.h"
46 #include "contextioerr.h"
47 #include "mathfem.h"
48 #include "classfactory.h"
49 #include "engngm.h"
50 
51 #ifdef __OOFEG
52  #include "oofeggraphiccontext.h"
53  #include "oofegutils.h"
54  #include "connectivitytable.h"
55 #endif
56 
57 namespace oofem {
58 REGISTER_Element(LIBeam3d2);
59 
60 LIBeam3d2 :: LIBeam3d2(int n, Domain *aDomain) : NLStructuralElement(n, aDomain), tc(), tempTc()
61 {
62  numberOfDofMans = 2;
63  referenceNode = 0;
64  length = 0.;
65 
66  tempTcCounter = 0;
67 }
68 
69 
70 void
72 // Returns the strain matrix of the receiver.
73 // eeps = {\eps_x, \gamma_xz, \gamma_xy, \der{phi_x}{x}, \kappa_y, \kappa_z}^T
74 {
75  double l, ksi, n1, n2, n1x, n2x;
76 
77  if ( this->nlGeometry ) {
79  } else {
80  l = this->computeLength();
81  }
82 
83  ksi = gp->giveNaturalCoordinate(1);
84 
85  answer.resize(6, 12);
86  answer.zero();
87 
88  n1 = 0.5 * ( 1 - ksi );
89  n2 = 0.5 * ( 1. + ksi );
90  n1x = -1.0 / l;
91  n2x = 1.0 / l;
92 
93  answer.at(1, 1) = -1. / l;
94  answer.at(1, 7) = 1. / l;
95 
96  answer.at(2, 3) = n1x;
97  answer.at(2, 5) = n1;
98  answer.at(2, 9) = n2x;
99  answer.at(2, 11) = n2;
100 
101  answer.at(3, 2) = n1x;
102  answer.at(3, 6) = -n1;
103  answer.at(3, 8) = n2x;
104  answer.at(3, 12) = -n2;
105 
106  answer.at(4, 4) = -1. / l;
107  answer.at(4, 10) = 1. / l;
108 
109  answer.at(5, 5) = n1x;
110  answer.at(5, 11) = n2x;
111 
112  answer.at(6, 6) = n1x;
113  answer.at(6, 12) = n2x;
114 }
115 
116 void
118 // Sets up the array of Gauss Points of the receiver.
119 {
120  if ( integrationRulesArray.size() == 0 ) {
121  integrationRulesArray.resize( 1 );
122  integrationRulesArray [ 0 ].reset( new GaussIntegrationRule(1, this, 1, 2) );
124  }
125 }
126 
127 
128 
129 void
131 // Returns the lumped mass matrix of the receiver. This expression is
132 // valid in both local and global axes.
133 {
134  GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
135  double density = this->giveStructuralCrossSection()->give('d', gp);
136  double halfMass = density * this->giveCrossSection()->give(CS_Area, gp) * this->computeLength() / 2.;
137  answer.resize(12, 12);
138  answer.zero();
139  answer.at(1, 1) = answer.at(2, 2) = answer.at(3, 3) = halfMass;
140  answer.at(7, 7) = answer.at(8, 8) = answer.at(9, 9) = halfMass;
141 }
142 
143 
144 void
146 // Returns the displacement interpolation matrix {N} of the receiver, eva-
147 // luated at gp.
148 {
149  double ksi, n1, n2;
150 
151  ksi = iLocCoord.at(1);
152  n1 = ( 1. - ksi ) * 0.5;
153  n2 = ( 1. + ksi ) * 0.5;
154 
155  answer.resize(6, 12);
156  answer.zero();
157 
158  answer.at(1, 1) = n1;
159  answer.at(1, 7) = n2;
160  answer.at(2, 2) = n1;
161  answer.at(2, 8) = n2;
162  answer.at(3, 3) = n1;
163  answer.at(3, 9) = n2;
164 
165  answer.at(4, 4) = n1;
166  answer.at(4, 10) = n2;
167  answer.at(5, 5) = n1;
168  answer.at(5, 11) = n2;
169  answer.at(6, 6) = n1;
170  answer.at(6, 12) = n2;
171 }
172 
173 
174 void
176 // Returns the stiffness matrix of the receiver, expressed in the global
177 // axes.
178 {
179  StructuralElement :: computeStiffnessMatrix(answer, rMode, tStep);
180 }
181 
182 
183 bool
185 {
186  answer.resize(12, 12);
187  answer.zero();
188 
189  if ( this->nlGeometry == 0 ) {
190  FloatMatrix lcs;
191 
192  this->giveLocalCoordinateSystem(lcs);
193  for ( int i = 1; i <= 3; i++ ) {
194  for ( int j = 1; j <= 3; j++ ) {
195  answer.at(i, j) = lcs.at(i, j);
196  answer.at(i + 3, j + 3) = lcs.at(i, j);
197  answer.at(i + 6, j + 6) = lcs.at(i, j);
198  answer.at(i + 9, j + 9) = lcs.at(i, j);
199  }
200  }
201  } else {
203 
204  for ( int i = 1; i <= 3; i++ ) {
205  for ( int j = 1; j <= 3; j++ ) {
206  answer.at(i, j) = tempTc.at(j, i);
207  answer.at(i + 3, j + 3) = tempTc.at(j, i);
208  answer.at(i + 6, j + 6) = tempTc.at(j, i);
209  answer.at(i + 9, j + 9) = tempTc.at(j, i);
210  }
211  }
212  }
213 
214  return 1;
215 }
216 
217 
218 double
220 // Returns the length of the receiver. This method is valid only if 1
221 // Gauss point is used.
222 {
223  double weight = gp->giveWeight();
224  return weight * 0.5 * this->computeLength();
225 }
226 
227 
228 void
230 {
231  answer = {D_u, D_v, D_w, R_u, R_v, R_w};
232 }
233 
234 
235 int
237 {
238  double ksi, n1, n2;
239 
240  ksi = lcoords.at(1);
241  n1 = ( 1. - ksi ) * 0.5;
242  n2 = ( 1. + ksi ) * 0.5;
243 
244  answer.resize(3);
245  answer.at(1) = n1 * this->giveNode(1)->giveCoordinate(1) + n2 *this->giveNode(2)->giveCoordinate(1);
246  answer.at(2) = n1 * this->giveNode(1)->giveCoordinate(2) + n2 *this->giveNode(2)->giveCoordinate(2);
247  answer.at(3) = n1 * this->giveNode(1)->giveCoordinate(3) + n2 *this->giveNode(2)->giveCoordinate(3);
248 
249  return 1;
250 }
251 
252 
253 void
255 {
256  this->giveStructuralCrossSection()->give3dBeamStiffMtrx(answer, rMode, gp, tStep);
257 }
258 
259 
260 void
262 {
263  this->giveStructuralCrossSection()->giveGeneralizedStress_Beam3d(answer, gp, strain, tStep);
264 }
265 
266 
267 int
269 {
270  return ( ( ext == Element_EdgeLoadSupport ) ? 1 : 0 );
271 }
272 
273 
274 double
276 // Returns the length of the receiver.
277 {
278  double dx, dy, dz;
279  Node *nodeA, *nodeB;
280 
281  if ( length == 0. ) {
282  nodeA = this->giveNode(1);
283  nodeB = this->giveNode(2);
284  dx = nodeB->giveCoordinate(1) - nodeA->giveCoordinate(1);
285  dy = nodeB->giveCoordinate(2) - nodeA->giveCoordinate(2);
286  dz = nodeB->giveCoordinate(3) - nodeA->giveCoordinate(3);
287  length = sqrt(dx * dx + dy * dy + dz * dz);
288  }
289 
290  return length;
291 }
292 
293 
296 {
297  IRResultType result; // Required by IR_GIVE_FIELD macro
298 
299  // first call parent
301  if ( result != IRRT_OK ) {
302  return result;
303  }
304 
306  if ( referenceNode == 0 ) {
307  OOFEM_ERROR("wrong reference node specified");
308  }
309 
310  // if (this->hasString (initString, "dofstocondense")) {
311  // dofsToCondense = this->ReadIntArray (initString, "dofstocondense");
312  // if (dofsToCondense->giveSize() >= 12)
313  // OOFEM_ERROR("wrong input data for condensed dofs");
314  // } else {
315  // dofsToCondense = NULL;
316  // }
317  // compute initial triad at centre - requires nodal coordinates
318  FloatMatrix lcs;
319  this->giveLocalCoordinateSystem(lcs);
320  this->tc.beTranspositionOf(lcs);
321 
322  return IRRT_OK;
323 }
324 
325 
326 void
327 LIBeam3d2 :: giveEdgeDofMapping(IntArray &answer, int iEdge) const
328 {
329  /*
330  * provides dof mapping of local edge dofs (only nonzero are taken into account)
331  * to global element dofs
332  */
333  if ( iEdge != 1 ) {
334  OOFEM_ERROR("wrong edge number");
335  }
336 
337 
338  answer.resize(12);
339  for ( int i = 1; i <= 12; i++ ) {
340  answer.at(i) = i;
341  }
342 }
343 
344 
345 double
347 {
348  if ( iEdge != 1 ) { // edge between nodes 1 2
349  OOFEM_ERROR("wrong egde number");
350  }
351 
352  double weight = gp->giveWeight();
353  return 0.5 * this->computeLength() * weight;
354 }
355 
356 
357 int
359 {
360  /*
361  * Returns transformation matrix from global coordinate system to local
362  * element coordinate system for element load vector components.
363  * If no transformation is necessary, answer is empty matrix (default);
364  */
365 
366  // f(elemLocal) = T * f (global)
367 
368  FloatMatrix lcs;
369 
370  answer.resize(6, 6);
371  answer.zero();
372 
373  this->giveLocalCoordinateSystem(lcs);
374  for ( int i = 1; i <= 3; i++ ) {
375  for ( int j = 1; j <= 3; j++ ) {
376  answer.at(i, j) = lcs.at(i, j);
377  answer.at(3 + i, 3 + j) = lcs.at(i, j);
378  }
379  }
380 
381  return 1;
382 }
383 
384 
385 void
387 {
388  FloatArray lc(1);
389  NLStructuralElement :: computeBodyLoadVectorAt(answer, load, tStep, mode);
390  answer.times( this->giveCrossSection()->give(CS_Area, lc, this) );
391 }
392 
393 
394 int
396 {
397  // returns transformation matrix from
398  // edge local coordinate system
399  // to element local c.s
400  // (same as global c.s in this case)
401  //
402  // i.e. f(element local) = T * f(edge local)
403  //
404  answer.clear();
405  return 0;
406 }
407 
408 
409 int
411 //
412 // returns a unit vectors of local coordinate system at element
413 // stored rowwise (mainly used by some materials with ortho and anisotrophy)
414 //
415 {
416  FloatArray lx(3), ly(3), lz(3), help(3);
417  double length = this->computeLength();
418  Node *nodeA, *nodeB, *refNode;
419 
420  answer.resize(3, 3);
421  answer.zero();
422  nodeA = this->giveNode(1);
423  nodeB = this->giveNode(2);
424  refNode = this->giveDomain()->giveNode(this->referenceNode);
425 
426  for ( int i = 1; i <= 3; i++ ) {
427  lx.at(i) = ( nodeB->giveCoordinate(i) - nodeA->giveCoordinate(i) ) / length;
428  help.at(i) = ( refNode->giveCoordinate(i) - nodeA->giveCoordinate(i) );
429  }
430 
431  lz.beVectorProductOf(lx, help);
432  lz.normalize();
433  ly.beVectorProductOf(lz, lx);
434  ly.normalize();
435 
436  for ( int i = 1; i <= 3; i++ ) {
437  answer.at(1, i) = lx.at(i);
438  answer.at(2, i) = ly.at(i);
439  answer.at(3, i) = lz.at(i);
440  }
441 
442  return 1;
443 }
444 
445 
446 void
448 {
449  // test if not previously done
450  if ( tStep->giveSolutionStateCounter() == tempTcCounter ) {
451  return;
452  }
453 
454  FloatArray u, centreSpin(3);
455  FloatMatrix dR(3, 3);
456 
457  // ask element's displacement increments
458  this->computeVectorOf(VM_Incremental, tStep, u);
459 
460  // interpolate spin at the centre
461  centreSpin.at(1) = 0.5 * ( u.at(4) + u.at(10) );
462  centreSpin.at(2) = 0.5 * ( u.at(5) + u.at(11) );
463  centreSpin.at(3) = 0.5 * ( u.at(6) + u.at(12) );
464 
465  // compute rotation matrix from centreSpin pseudovector
466  this->computeRotMtrx(dR, centreSpin);
467 
468  // update triad
469  tempTc.beProductOf(dR, tc);
470  // remember timestamp
472 }
473 
474 
475 void
477 {
478  FloatMatrix S(3, 3), SS(3, 3);
479  double psiSize;
480 
481  if ( psi.giveSize() != 3 ) {
482  OOFEM_ERROR("psi param size mismatch");
483  }
484 
485  answer.resize(3, 3);
486  answer.zero();
487 
488  psiSize = psi.computeNorm();
489  answer.at(1, 1) = answer.at(2, 2) = answer.at(3, 3) = 1.;
490 
491  if ( psiSize <= 1.e-40 ) {
492  return;
493  }
494 
495  this->computeSMtrx(S, psi);
496  SS.beProductOf(S, S);
497  S.times(sin(psiSize) / psiSize);
498  SS.times( ( 1. - cos(psiSize) ) / ( psiSize * psiSize ) );
499 
500  answer.add(S);
501  answer.add(SS);
502 }
503 
504 
505 void
507 {
508  if ( vec.giveSize() != 3 ) {
509  OOFEM_ERROR("vec param size mismatch");
510  }
511 
512  answer.resize(3, 3);
513 
514  answer.at(1, 1) = answer.at(2, 2) = answer.at(3, 3) = 0.;
515  answer.at(1, 2) = -vec.at(3);
516  answer.at(1, 3) = vec.at(2);
517  answer.at(2, 1) = vec.at(3);
518  answer.at(2, 3) = -vec.at(1);
519  answer.at(3, 1) = -vec.at(2);
520  answer.at(3, 2) = vec.at(1);
521 }
522 
523 
524 void
526 {
527  FloatArray ui, PrevEpsilon;
528  FloatMatrix b;
529 
530  this->computeVectorOf(VM_Incremental, tStep, ui);
531 
532  this->computeBmatrixAt(gp, b);
533  // increment of strains
534  answer.beProductOf(b, ui);
535 
536  if ( this->nlGeometry ) {
537  // add increment to previous total value
538  PrevEpsilon = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() )->giveStrainVector();
539  if ( PrevEpsilon.giveSize() ) {
540  answer.add(PrevEpsilon);
541  }
542  }
543 }
544 
545 
546 double
548 // Returns the length of the receiver.
549 {
550  double dx, dy, dz;
551  FloatArray u;
552 
553  this->computeVectorOf(VM_Total, tStep, u);
554 
555  dx = ( this->giveNode(2)->giveCoordinate(1) + u.at(7) ) -
556  ( this->giveNode(1)->giveCoordinate(1) + u.at(1) );
557  dy = ( this->giveNode(2)->giveCoordinate(2) + u.at(8) ) -
558  ( this->giveNode(1)->giveCoordinate(2) + u.at(2) );
559  dz = ( this->giveNode(2)->giveCoordinate(3) + u.at(9) ) -
560  ( this->giveNode(1)->giveCoordinate(3) + u.at(3) );
561 
562  return sqrt(dx * dx + dy * dy + dz * dz);
563 }
564 
565 
567 // Updates the receiver at end of step.
568 {
570 
571  // update triad
572  this->updateTempTriad(tStep);
573  tc = tempTc;
574 
575  // update curvature
576  //FloatArray curv;
577  //this->computeTempCurv (curv, tStep);
578  //kappa = curv;
579 }
580 
581 
582 void
584 // initializes receiver to new time step or can be used
585 // if current time step must be restarted
586 {
588  tempTc = tc;
589 }
590 
591 
592 
594 //
595 // saves full element context (saves state variables, that completely describe
596 // current state)
597 //
598 {
599  contextIOResultType iores;
600  if ( ( iores = NLStructuralElement :: saveContext(stream, mode, obj) ) != CIO_OK ) {
601  THROW_CIOERR(iores);
602  }
603 
604  if ( ( iores = tc.storeYourself(stream) ) != CIO_OK ) {
605  THROW_CIOERR(iores);
606  }
607 
608  return CIO_OK;
609 }
610 
611 
613 //
614 // restores full element context (saves state variables, that completely describe
615 // current state)
616 //
617 {
618  contextIOResultType iores;
619  if ( ( iores = NLStructuralElement :: restoreContext(stream, mode, obj) ) != CIO_OK ) {
620  THROW_CIOERR(iores);
621  }
622 
623  if ( ( iores = tc.restoreYourself(stream) ) != CIO_OK ) {
624  THROW_CIOERR(iores);
625  }
626 
627  return CIO_OK;
628 }
629 
630 
631 void
633  GaussPoint *slaveGp, TimeStep *tStep)
634 {
635  double layerYCoord, layerZCoord;
636 
637  layerZCoord = slaveGp->giveNaturalCoordinate(2);
638  layerYCoord = slaveGp->giveNaturalCoordinate(1);
639 
640  answer.resize(3); // {Exx,GMzx,GMxy}
641 
642  answer.at(1) = masterGpStrain.at(1) + masterGpStrain.at(5) * layerZCoord - masterGpStrain.at(6) * layerYCoord;
643  answer.at(2) = masterGpStrain.at(2);
644  answer.at(3) = masterGpStrain.at(3);
645 }
646 
647 
648 Interface *
650 {
651  if ( interface == FiberedCrossSectionInterfaceType ) {
652  return static_cast< FiberedCrossSectionInterface * >(this);
653  }
654 
655  return NULL;
656 }
657 
658 
659 #ifdef __OOFEG
660 void
662 {
663  GraphicObj *go;
664 
665  if ( !gc.testElementGraphicActivity(this) ) {
666  return;
667  }
668 
669  // if (!go) { // create new one
670  WCRec p [ 2 ]; /* poin */
671  EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
672  EASValsSetColor( gc.getElementColor() );
673  EASValsSetLayer(OOFEG_RAW_GEOMETRY_LAYER);
674  p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1);
675  p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveCoordinate(2);
676  p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveCoordinate(3);
677  p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1);
678  p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveCoordinate(2);
679  p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveCoordinate(3);
680  go = CreateLine3D(p);
681  EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
682  EGAttachObject(go, ( EObjectP ) this);
683  EMAddGraphicsToModel(ESIModel(), go);
684 }
685 
686 
687 void
689 {
690  GraphicObj *go;
691 
692  if ( !gc.testElementGraphicActivity(this) ) {
693  return;
694  }
695 
696  double defScale = gc.getDefScale();
697  // if (!go) { // create new one
698  WCRec p [ 2 ]; /* poin */
699  char const *colors[] = {
700  "red", "green", "blue"
701  };
702 
703  EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
704  EASValsSetColor( gc.getDeformedElementColor() );
705  EASValsSetLayer(OOFEG_DEFORMED_GEOMETRY_LAYER);
706  p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(1, tStep, defScale);
707  p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(2, tStep, defScale);
708  p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(3, tStep, defScale);
709 
710  p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(1, tStep, defScale);
711  p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(2, tStep, defScale);
712  p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(3, tStep, defScale);
713  go = CreateLine3D(p);
714  EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
715  EMAddGraphicsToModel(ESIModel(), go);
716 
717  // draw centre triad
718  int i, succ;
719  double coeff = this->computeLength() / 3.;
720  //this->updateTempTriad (tStep);
721 
722  p [ 0 ].x = 0.5 * ( ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(1, tStep, defScale) + ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(1, tStep, defScale) );
723  p [ 0 ].y = 0.5 * ( ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(2, tStep, defScale) + ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(2, tStep, defScale) );
724  p [ 0 ].z = 0.5 * ( ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(3, tStep, defScale) + ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(3, tStep, defScale) );
725 
726  // draw t1
727  for ( i = 1; i <= 3; i++ ) {
728  p [ 1 ].x = p [ 0 ].x + coeff *tc.at(1, i);
729  p [ 1 ].y = p [ 0 ].y + coeff *tc.at(2, i);
730  p [ 1 ].z = p [ 0 ].z + coeff *tc.at(3, i);
731 
732  EASValsSetColor( ColorGetPixelFromString(const_cast< char * >(colors [ i - 1 ]), & succ) );
733 
734  go = CreateLine3D(p);
735  EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
736  EMAddGraphicsToModel(ESIModel(), go);
737  }
738 }
739 
740 
741 void
743 {
744  WCRec p [ 2 ];
745  GraphicObj *go;
746  FloatArray v;
747  double defScale;
748 
749  if ( !gc.testElementGraphicActivity(this) ) {
750  return;
751  }
752 
753  EASValsSetLayer(OOFEG_VARPLOT_PATTERN_LAYER);
754  if ( gc.getInternalVarsDefGeoFlag() ) {
755  // use deformed geometry
756  defScale = gc.getDefScale();
757  p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(1, tStep, defScale);
758  p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(2, tStep, defScale);
759  p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(3, tStep, defScale);
760  p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(1, tStep, defScale);
761  p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(2, tStep, defScale);
762  p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(3, tStep, defScale);
763  } else {
764  p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1);
765  p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveCoordinate(2);
766  p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveCoordinate(3);
767  p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1);
768  p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveCoordinate(2);
769  p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveCoordinate(3);
770  }
771 
772  GaussPoint *gp;
773  double s;
774  gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
775  if ( giveIPValue(v, gp, gc.giveIntVarType(), tStep) == 0 ) {
776  return;
777  }
778 
779  s = v.at(1);
780  gc.updateFringeTableMinMax(& s, 1);
781 
782  go = CreateLine3D(p);
783  EASValsSetColor( gc.getElementColor() );
784  EASValsSetLayer(OOFEG_VARPLOT_PATTERN_LAYER);
785  EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
786  EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
787  EMAddGraphicsToModel(ESIModel(), go);
788 
789  EASValsSetMType(FILLED_CIRCLE_MARKER);
790  go = CreateMarkerWD3D( p, v.at(1) );
791  EGWithMaskChangeAttributes(LAYER_MASK | FILL_MASK | MTYPE_MASK, go);
792  EMAddGraphicsToModel(ESIModel(), go);
793 }
794 
795 #endif
796 } // 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 giveEdgeDofMapping(IntArray &answer, int iEdge) const
Assembles edge dof mapping mask, which provides mapping between edge local DOFs and "global" element ...
Definition: libeam3d2.C:327
void computeSMtrx(FloatMatrix &answer, FloatArray &vec)
Definition: libeam3d2.C:506
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: libeam3d2.C:295
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
virtual void updateYourself(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
Definition: libeam3d2.C:566
int nlGeometry
Flag indicating if geometrical nonlinearities apply.
virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
Computes constitutive matrix of receiver.
Definition: libeam3d2.C:254
Class and object Domain.
Definition: domain.h:115
void computeVectorOf(ValueModeType u, TimeStep *tStep, FloatArray &answer)
Returns local vector of unknowns.
Definition: element.C:86
#define _IFT_LIBeam3d2_refnode
Definition: libeam3d2.h:44
void initForNewStep()
Initializes receivers state to new time step.
Definition: libeam3d2.C:583
Domain * domain
Link to domain object, useful for communicating with other FEM components.
Definition: femcmpnn.h:82
Abstract base class for "structural" finite elements with geometrical nonlinearities.
The purpose of DataStream abstract class is to allow to store/restore context to different streams...
Definition: datastream.h:54
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
Definition: libeam3d2.C:117
virtual void drawScalar(oofegGraphicContext &gc, TimeStep *tStep)
Definition: libeam3d2.C:742
virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
Definition: libeam3d2.C:661
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]
EngngModel * giveEngngModel()
Returns engineering model to which receiver is associated.
Definition: domain.C:433
virtual void initForNewStep()
Initializes receivers state to new time step.
Definition: element.C:846
virtual void giveGeneralizedStress_Beam3d(FloatArray &answer, GaussPoint *gp, const FloatArray &generalizedStrain, TimeStep *tStep)=0
#define S(p)
Definition: mdm.C:481
virtual int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords)
Computes the global coordinates from given element&#39;s local coordinates.
Definition: libeam3d2.C:236
virtual double giveCoordinate(int i)
Definition: node.C:82
virtual double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
Definition: libeam3d2.C:219
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: libeam3d2.C:386
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 computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes lumped mass matrix of receiver.
Definition: libeam3d2.C:130
virtual double computeEdgeVolumeAround(GaussPoint *gp, int iEdge)
Computes volume related to integration point on local edge.
Definition: libeam3d2.C:346
#define THROW_CIOERR(e)
Definition: contextioerr.h:61
#define OOFEG_DEFORMED_GEOMETRY_LAYER
contextIOResultType storeYourself(DataStream &stream) const
Definition: floatmatrix.C:1852
virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Restores the receiver state previously written in stream.
Definition: element.C:970
StateCounterType tempTcCounter
Time stamp of temporary centre triad.
Definition: libeam3d2.h:71
InternalStateType giveIntVarType()
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 int computeLoadLEToLRotationMatrix(FloatMatrix &answer, int iEdge, GaussPoint *gp)
Returns transformation matrix from local edge c.s to element local coordinate system of load vector c...
Definition: libeam3d2.C:395
StateCounterType giveSolutionStateCounter()
Returns current solution state counter.
Definition: timestep.h:188
virtual double computeLength()
Computes the length (zero for all but 1D geometries)
Definition: libeam3d2.C:275
#define OOFEM_ERROR(...)
Definition: error.h:61
REGISTER_Element(LSpace)
ElementExtension
Type representing element extension.
#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
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Receiver becomes the result of the product of aMatrix and anArray.
Definition: floatarray.C:676
virtual bool computeGtoLRotationMatrix(FloatMatrix &answer)
Returns transformation matrix from global c.s.
Definition: libeam3d2.C:184
virtual double giveUpdatedCoordinate(int ic, TimeStep *tStep, double scale=1.)
Returns updated ic-th coordinate of receiver.
Definition: node.C:245
double at(int i, int j) const
Coefficient access function.
Definition: floatmatrix.h:176
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
Computes numerically stiffness matrix of receiver.
void resize(int n)
Checks size of receiver towards requested bounds.
Definition: intarray.C:124
virtual void give3dBeamStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)=0
Computes the stiffness matrix for 2d beams.
virtual int computeLoadGToLRotationMtrx(FloatMatrix &answer)
Returns transformation matrix from global coordinate system to local element coordinate system for el...
Definition: libeam3d2.C:358
virtual int setupIntegrationPoints(IntegrationRule &irule, int npoints, Element *element)
Sets up integration rule for the given element.
Definition: crosssection.C:54
FloatMatrix tc
Last equilibrium triad at the centre.
Definition: libeam3d2.h:67
double giveCurrentLength(TimeStep *tStep)
Definition: libeam3d2.C:547
contextIOResultType restoreYourself(DataStream &stream)
Definition: floatmatrix.C:1877
Class representing vector of real numbers.
Definition: floatarray.h:82
virtual Interface * giveInterface(InterfaceType it)
Interface requesting service.
Definition: libeam3d2.C:649
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 void updateYourself(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
void computeRotMtrx(FloatMatrix &answer, FloatArray &psi)
Definition: libeam3d2.C:476
IntegrationPointStatus * giveMaterialStatus()
Returns reference to associated material status (NULL if not defined).
Definition: gausspoint.h:205
double computeNorm() const
Computes the norm (or length) of the vector.
Definition: floatarray.C:840
void resize(int rows, int cols)
Checks size of receiver towards requested bounds.
Definition: floatmatrix.C:1358
virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
Computes interpolation matrix for element unknowns.
Definition: libeam3d2.C:145
Class representing the general Input Record.
Definition: inputrecord.h:101
virtual void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType)
Definition: libeam3d2.C:688
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: libeam3d2.C:261
void add(const FloatMatrix &a)
Adds matrix to the receiver.
Definition: floatmatrix.C:1023
virtual void computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep)
Compute strain vector of receiver evaluated at given integration point at given time step from elemen...
Definition: libeam3d2.C:525
virtual int giveLocalCoordinateSystem(FloatMatrix &answer)
Returns local coordinate system of receiver Required by material models with ortho- and anisotrophy...
Definition: libeam3d2.C:410
Class Interface.
Definition: interface.h:82
virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Stores receiver state to output stream.
Definition: element.C:885
#define OOFEG_DEFORMED_GEOMETRY_WIDTH
void times(double s)
Multiplies receiver with scalar.
Definition: floatarray.C:818
void beTranspositionOf(const FloatMatrix &src)
Assigns to the receiver the transposition of parameter.
Definition: floatmatrix.C:323
FloatMatrix tempTc
Temporary triad at the centre.
Definition: libeam3d2.h:69
long ContextMode
Context mode (mask), defining the type of information written/read to/from context.
Definition: contextmode.h:43
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 IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
The element interface required by FiberedCrossSection.
Definition: fiberedcs.h:220
void zero()
Zeroes all coefficient of receiver.
Definition: floatmatrix.C:1326
Domain * giveDomain() const
Definition: femcmpnn.h:100
InterfaceType
Enumerative type, used to identify interface type.
Definition: interfacetype.h:43
void updateFringeTableMinMax(double *s, int size)
Load is base abstract class for all loads.
Definition: load.h:61
virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj)
Stores receiver state to output stream.
Definition: libeam3d2.C:593
void beProductOf(const FloatMatrix &a, const FloatMatrix &b)
Assigns to the receiver product of .
Definition: floatmatrix.C:337
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
Computes the stiffness matrix of receiver.
Definition: libeam3d2.C:175
virtual TimeStep * giveCurrentStep(bool force=false)
Returns current time step.
Definition: engngm.h:683
int giveSize() const
Returns the size of receiver.
Definition: floatarray.h:218
Node * giveNode(int n)
Service for accessing particular domain node.
Definition: domain.h:371
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
#define IR_GIVE_FIELD(__ir, __value, __id)
Macro facilitating the use of input record reading methods.
Definition: inputrecord.h:69
virtual void computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int=1, int=ALL_STRAINS)
Computes the geometrical matrix of receiver in given integration point.
Definition: libeam3d2.C:71
Node * giveNode(int i) const
Returns reference to the i-th node of element.
Definition: element.h:610
virtual void giveDofManDofIDMask(int inode, IntArray &) const
Returns dofmanager dof mask for node.
Definition: libeam3d2.C:229
#define OOFEG_VARPLOT_PATTERN_LAYER
LIBeam3d2(int n, Domain *d)
Definition: libeam3d2.C:60
virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj)
Restores the receiver state previously written in stream.
Definition: libeam3d2.C:612
Class representing integration point in finite element program.
Definition: gausspoint.h:93
void updateTempTriad(TimeStep *tStep)
Definition: libeam3d2.C:447
void FiberedCrossSectionInterface_computeStrainVectorInFiber(FloatArray &answer, const FloatArray &masterGpStrain, GaussPoint *slaveGp, TimeStep *tStep)
Computes full 3d strain vector in element fiber.
Definition: libeam3d2.C:632
Class representing solution step.
Definition: timestep.h:80
virtual int testElementExtension(ElementExtension ext)
Tests if the element implements required extension.
Definition: libeam3d2.C:268
int numberOfDofMans
Number of dofmanagers.
Definition: element.h:149
void add(const FloatArray &src)
Adds array src to receiver.
Definition: floatarray.C:156
Element extension for edge loads.
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