OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
beam2d.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/beam2d.h"
36 #include "../sm/Materials/structuralms.h"
37 #include "fei2dlinelin.h"
38 #include "fei2dlinehermite.h"
39 #include "node.h"
40 #include "material.h"
41 #include "crosssection.h"
42 #include "gausspoint.h"
43 #include "gaussintegrationrule.h"
44 #include "floatmatrix.h"
45 #include "intarray.h"
46 #include "floatarray.h"
47 #include "engngm.h"
48 #include "bodyload.h"
49 #include "boundaryload.h"
50 #include "mathfem.h"
51 #include "classfactory.h"
52 #include "elementinternaldofman.h"
53 #include "masterdof.h"
54 
55 #ifdef __OOFEG
56  #include "oofeggraphiccontext.h"
57 #endif
58 
59 namespace oofem {
60 REGISTER_Element(Beam2d);
61 
62 // Set up interpolation coordinates
63 FEI2dLineLin Beam2d :: interp_geom(1, 3);
64 FEI2dLineHermite Beam2d :: interp_beam(1, 3);
65 
67 {
68  numberOfDofMans = 2;
69 
70  kappa = -1; // set kappa to undef value (should be always > 0.)
71  length = 0.;
72  pitch = 10.; // a dummy value
74 
75  ghostNodes [ 0 ] = ghostNodes [ 1 ] = NULL;
77 }
78 
79 
81 {
82  delete ghostNodes [ 0 ];
83  delete ghostNodes [ 1 ];
84 }
85 
86 
88 
89 
90 Interface *
92 {
93  if ( interface == LayeredCrossSectionInterfaceType ) {
94  return static_cast< LayeredCrossSectionInterface * >(this);
95  }
96 
97  return NULL;
98 }
99 
100 
101 void
102 Beam2d :: computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int li, int ui)
103 {
104  double l, ksi, kappa, c1;
106 
107  l = this->computeLength();
108  ksi = 0.5 + 0.5 * gp->giveNaturalCoordinate(1);
109  kappa = this->giveKappaCoeff(tStep);
110  c1 = 1. + 2. * kappa;
111 
112  answer.resize(3, 6);
113  answer.zero();
114 
115  answer.at(1, 1) = -1. / l;
116  answer.at(1, 4) = 1. / l;
117  answer.at(2, 2) = ( 6. - 12. * ksi ) / ( l * l * c1 );
118  answer.at(2, 3) = ( -2. * ( 2. + kappa ) + 6. * ksi ) / ( l * c1 );
119  answer.at(2, 5) = ( -6. + 12. * ksi ) / ( l * l * c1 );
120  answer.at(2, 6) = ( -2. * ( 1. - kappa ) + 6. * ksi ) / ( l * c1 );
121  answer.at(3, 2) = ( -2. * kappa ) / ( l * c1 );
122  answer.at(3, 3) = kappa / ( c1 );
123  answer.at(3, 5) = 2. * kappa / ( l * c1 );
124  answer.at(3, 6) = kappa / ( c1 );
125 }
126 
127 
128 void
130 {
131  if ( integrationRulesArray.size() == 0 ) {
132  // the gauss point is used only when methods from crosssection and/or material
133  // classes are requested
134  integrationRulesArray.resize(1);
135  integrationRulesArray [ 0 ].reset( new GaussIntegrationRule(1, this, 1, 3) );
137  }
138 }
139 
140 
141 void
143 // Returns the displacement interpolation matrix {N} of the receiver, eva-
144 // luated at gp. Used for numerical calculation of consistent mass
145 // matrix. Must contain only interpolation for displacement terms,
146 // not for any rotations. (Inertia forces do not work on rotations).
147 // r = {u1,w1,fi_y1,u2,w2,fi_y2}^T
148 {
149  double l, ksi, ksi2, ksi3, kappa, c1;
150  TimeStep *tStep = this->domain->giveEngngModel()->giveCurrentStep();
151 
152  l = this->computeLength();
153  ksi = 0.5 + 0.5 * iLocCoord.at(1);
154  kappa = this->giveKappaCoeff(tStep);
155  c1 = 1. + 2. * kappa;
156  ksi2 = ksi * ksi;
157  ksi3 = ksi2 * ksi;
158 
159  answer.resize(3, 6);
160  answer.zero();
161 
162  answer.at(1, 1) = 1. - ksi;
163  answer.at(1, 4) = ksi;
164  answer.at(2, 2) = ( c1 - 2. * kappa * ksi - 3. * ksi2 + 2. * ksi3 ) / c1;
165  answer.at(2, 3) = l * ( -( 1. + kappa ) * ksi + ( 2. + kappa ) * ksi2 - ksi3 ) / c1;
166  answer.at(2, 5) = ( 2. * kappa * ksi + 3. * ksi2 - 2. * ksi3 ) / c1;
167  answer.at(2, 6) = l * ( kappa * ksi + ( 1. - kappa ) * ksi2 - ksi3 ) / c1;
168  answer.at(3, 2) = ( 6. * ksi - 6. * ksi2 ) / ( l * c1 );
169  answer.at(3, 3) = ( c1 - 2. * ( 2. + kappa ) * ksi + 3. * ksi2 ) / c1;
170  answer.at(3, 5) = ( -6. * ksi + 6. * ksi2 ) / ( l * c1 );
171  answer.at(3, 6) = ( -2. * ( 1. - kappa ) * ksi + 3. * ksi2 ) / c1;
172 }
173 
174 
175 void
177 {
178  double l = this->computeLength();
179  FloatMatrix B, d, DB;
180  answer.clear();
181  for ( auto &gp : *this->giveDefaultIntegrationRulePtr() ) {
182  this->computeBmatrixAt(gp, B);
183  this->computeConstitutiveMatrixAt(d, rMode, gp, tStep);
184  double dV = gp->giveWeight() * 0.5 * l;
185  DB.beProductOf(d, B);
186  answer.plusProductSymmUpper(B, DB, dV);
187  }
188  answer.symmetrized();
189 }
190 
191 
192 void
194 {
195  this->giveStructuralCrossSection()->give2dBeamStiffMtrx(answer, rMode, gp, tStep);
196 }
197 
198 
199 void
201 {
202  this->giveStructuralCrossSection()->giveGeneralizedStress_Beam2d(answer, gp, strain, tStep);
203 }
204 
205 
206 bool
208 {
209  double sine, cosine;
210 
211  int ndofs = computeNumberOfGlobalDofs();
212  answer.resize(ndofs, ndofs);
213  answer.zero();
214 
215  sine = sin( this->givePitch() );
216  cosine = cos(pitch);
217  answer.at(1, 1) = cosine;
218  answer.at(1, 2) = sine;
219  answer.at(2, 1) = -sine;
220  answer.at(2, 2) = cosine;
221  answer.at(3, 3) = 1.;
222  answer.at(4, 4) = cosine;
223  answer.at(4, 5) = sine;
224  answer.at(5, 4) = -sine;
225  answer.at(5, 5) = cosine;
226  answer.at(6, 6) = 1.;
227 
228  for ( int i = 7; i <= ndofs; i++ ) {
229  answer.at(i, i) = 1.0;
230  }
231 
232  if ( this->hasDofs2Condense() ) {
233  int condensedDofCounter = 0;
234  DofIDItem dofids[] = {
235  D_u, D_w, R_v
236  };
237  FloatMatrix l2p(6, ndofs); // local -> primary
238  l2p.zero();
239  // loop over nodes
240  for ( int inode = 0; inode < 2; inode++ ) {
241  // loop over DOFs
242  for ( int idof = 0; idof < 3; idof++ ) {
243  int eq = inode * 3 + idof + 1;
244  if ( ghostNodes [ inode ] ) {
245  if ( ghostNodes [ inode ]->hasDofID(dofids [ idof ]) ) {
246  condensedDofCounter++;
247  l2p.at(eq, 6 + condensedDofCounter) = 1.0;
248  continue;
249  }
250  }
251  l2p.at(eq, eq) = 1.0;
252  }
253  }
254 
255  FloatMatrix g2l(answer);
256  answer.beProductOf(l2p, g2l);
257  }
258 
259  return true;
260 }
261 
262 
263 double
265 {
266  return 0.5 * this->computeLength() * gp->giveWeight();
267 }
268 
269 
270 void
271 Beam2d :: computeStrainVectorInLayer(FloatArray &answer, const FloatArray &masterGpStrain, GaussPoint *masterGp, GaussPoint *slaveGp, TimeStep *tStep)
272 //
273 // returns full 3d strain vector of given layer (whose z-coordinate from center-line is
274 // stored in slaveGp) for given tStep
275 //
276 {
277  double layerZeta, layerZCoord, top, bottom;
278 
279  top = this->giveCrossSection()->give(CS_TopZCoord, masterGp);
280  bottom = this->giveCrossSection()->give(CS_BottomZCoord, masterGp);
281  layerZeta = slaveGp->giveNaturalCoordinate(3);
282  layerZCoord = 0.5 * ( ( 1. - layerZeta ) * bottom + ( 1. + layerZeta ) * top );
283 
284  answer.resize(2); // {Exx,GMzx}
285 
286  answer.at(1) = masterGpStrain.at(1) + masterGpStrain.at(2) * layerZCoord;
287  answer.at(2) = masterGpStrain.at(3);
288 }
289 
290 
291 void
292 Beam2d :: giveDofManDofIDMask(int inode, IntArray &answer) const
293 {
294  answer = {
295  D_u, D_w, R_v
296  };
297 }
298 
299 
300 double
302 {
303  double dx, dy;
304  Node *nodeA, *nodeB;
305 
306  if ( length == 0. ) {
307  nodeA = this->giveNode(1);
308  nodeB = this->giveNode(2);
309  dx = nodeB->giveCoordinate(1) - nodeA->giveCoordinate(1);
310  dy = nodeB->giveCoordinate(3) - nodeA->giveCoordinate(3);
311  length = sqrt(dx * dx + dy * dy);
312  }
313 
314  return length;
315 }
316 
317 
318 double
320 {
321  double xA, xB, yA, yB;
322  Node *nodeA, *nodeB;
323 
324  if ( pitch == 10. ) { // 10. : dummy initialization value
325  nodeA = this->giveNode(1);
326  nodeB = this->giveNode(2);
327  xA = nodeA->giveCoordinate(1);
328  xB = nodeB->giveCoordinate(1);
329  yA = nodeA->giveCoordinate(3);
330  yB = nodeB->giveCoordinate(3);
331  pitch = atan2(yB - yA, xB - xA);
332  }
333 
334  return pitch;
335 }
336 
337 
338 double
340 {
341  // kappa = (6*E*I)/(k*G*A*l^2)
342 
343  if ( kappa < 0. ) {
344  FloatMatrix d;
345  double l = this->computeLength();
346 
347  this->computeConstitutiveMatrixAt(d, ElasticStiffness, integrationRulesArray [ 0 ]->getIntegrationPoint(0), tStep);
348  kappa = 6. * d.at(2, 2) / ( d.at(3, 3) * l * l );
349  }
350 
351  return kappa;
352 }
353 
354 
355 int
357 //
358 // returns a unit vectors of local coordinate system at element
359 // stored rowwise (mainly used by some materials with ortho and anisotrophy)
360 //
361 {
362  double sine, cosine;
363 
364  answer.resize(3, 3);
365  answer.zero();
366 
367  sine = sin( this->givePitch() );
368  cosine = cos(pitch);
369 
370  answer.at(1, 1) = cosine;
371  answer.at(1, 2) = sine;
372  answer.at(2, 1) = -sine;
373  answer.at(2, 2) = cosine;
374  answer.at(3, 3) = 1.0;
375 
376  return 1;
377 }
378 
379 
382 {
383  IRResultType result; // Required by IR_GIVE_FIELD macro
384  // first call parent
386 
388  IntArray val;
390  if ( val.giveSize() >= 6 ) {
391  OOFEM_WARNING("wrong input data for condensed dofs");
392  return IRRT_BAD_FORMAT;
393  }
394 
395  DofIDItem mask[] = {
396  D_u, D_w, R_v
397  };
398  this->numberOfCondensedDofs = val.giveSize();
399  for ( int i = 1; i <= val.giveSize(); i++ ) {
400  if ( val.at(i) <= 3 ) {
401  if ( ghostNodes [ 0 ] == NULL ) {
402  ghostNodes [ 0 ] = new ElementDofManager(1, giveDomain(), this);
403  }
404  ghostNodes [ 0 ]->appendDof( new MasterDof(ghostNodes [ 0 ], mask [ val.at(i) - 1 ]) );
405  } else {
406  if ( ghostNodes [ 1 ] == NULL ) {
407  ghostNodes [ 1 ] = new ElementDofManager(1, giveDomain(), this);
408  }
409  ghostNodes [ 1 ]->appendDof( new MasterDof(ghostNodes [ 1 ], mask [ val.at(i) - 4 ]) );
410  }
411  }
412 
413  }
414  return IRRT_OK;
415 }
416 
417 
418 
419 void
420 Beam2d :: giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord)
421 {
422  BeamBaseElement :: giveInternalForcesVector(answer, tStep, useUpdatedGpRecord);
423 }
424 
425 
426 void
428 {
429  // stress equivalent vector in nodes (vector of internal forces)
430  FloatArray load;
431 
432  this->giveInternalForcesVector(answer, tStep, false);
433 
434  // subtract exact end forces due to nonnodal loading
435  this->computeLocalForceLoadVector(load, tStep, VM_Total);
436  if ( load.isNotEmpty() ) {
437  answer.subtract(load);
438  }
439 
440 }
441 
442 
443 void
445 {
446  answer.clear();
447 
448  if ( edge != 1 ) {
449  OOFEM_ERROR("Beam2D only has 1 edge (the midline) that supports loads. Attempted to apply load to edge %d", edge);
450  }
451 
452  if ( type != ExternalForcesVector ) {
453  return;
454  }
455 
456  double l = this->computeLength();
457  FloatArray coords, t;
458  FloatMatrix N, T;
459 
460  answer.clear();
461  for ( auto &gp : *this->giveDefaultIntegrationRulePtr() ) {
462  const FloatArray &lcoords = gp->giveNaturalCoordinates();
463  this->computeNmatrixAt(lcoords, N);
464  if ( load ) {
465  this->computeGlobalCoordinates(coords, lcoords);
466  load->computeValues(t, tStep, coords, { D_u, D_w, R_v }, mode);
467  } else {
468  load->computeValues(t, tStep, lcoords, { D_u, D_w, R_v }, mode);
469  }
470 
471  if ( load->giveCoordSystMode() == Load :: CST_Global ) {
472  if ( this->computeLoadGToLRotationMtrx(T) ) {
473  t.rotatedWith(T, 'n');
474  }
475  }
476 
477  double dl = gp->giveWeight() * 0.5 * l;
478  answer.plusProduct(N, t, dl);
479  }
480 
481  if (global) {
482  // Loads from sets expects global c.s.
483  this->computeGtoLRotationMatrix(T);
484  answer.rotatedWith(T, 't');
485  }
486 }
487 
488 
489 void
491 {
492  FloatArray lc(1);
493  BeamBaseElement :: computeBodyLoadVectorAt(answer, load, tStep, mode);
494  answer.times( this->giveCrossSection()->give(CS_Area, lc, this) );
495 }
496 
497 
498 int
500 {
501  if ( type == IST_BeamForceMomentTensor ) {
502  answer = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() )->giveStressVector();
503  return 1;
504  } else if ( type == IST_BeamStrainCurvatureTensor ) {
505  answer = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() )->giveStrainVector();
506  return 1;
507  } else if ( type == IST_ShellForceTensor || type == IST_ShellStrainTensor ) {
508  // Order in generalized strain is:
509  // {\eps_x, \gamma_xz, \kappa_y}
510  const FloatArray &help = type == IST_ShellForceTensor ?
511  static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() )->giveStressVector() :
512  static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() )->giveStrainVector();
513 
514  answer.resize(6);
515  answer.at(1) = help.at(1); // nx
516  answer.at(2) = 0.0; // ny
517  answer.at(3) = 0.0; // nz
518  answer.at(4) = 0.0; // vyz
519  answer.at(5) = help.at(2); // vxz
520  answer.at(6) = 0.0; // vxy
521  return 1;
522  } else if ( type == IST_ShellMomentTensor || type == IST_CurvatureTensor ) {
523  const FloatArray &help = type == IST_ShellMomentTensor ?
524  static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() )->giveStressVector() :
525  static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() )->giveStrainVector();
526  answer.resize(6);
527  answer.at(1) = 0.0; // mx
528  answer.at(2) = 0.0; // my
529  answer.at(3) = 0.0; // mz
530  answer.at(4) = 0.0; // myz
531  answer.at(5) = 0.0; // mxz
532  answer.at(6) = help.at(3); // mxy
533  return 1;
534  } else {
535  return BeamBaseElement :: giveIPValue(answer, gp, type, tStep);
536  }
537 }
538 
539 
540 void
542 {
543  FloatArray rl, Fl;
544 
545  fprintf(File, "beam element %d (%8d) :\n", this->giveLabel(), this->giveNumber() );
546 
547  // ask for global element displacement vector
548  this->computeVectorOf(VM_Total, tStep, rl);
549 
550  // ask for global element end forces vector
551  this->giveEndForcesVector(Fl, tStep);
552 
553  fprintf(File, " local displacements ");
554  for ( auto &val : rl ) {
555  fprintf(File, " %.4e", val);
556  }
557 
558  fprintf(File, "\n local end forces ");
559  for ( auto &val : Fl ) {
560  fprintf(File, " %.4e", val);
561  }
562 
563  fprintf(File, "\n");
564 
565  for ( auto &iRule : integrationRulesArray ) {
566  iRule->printOutputAt(File, tStep);
567  }
568 }
569 
570 
571 void
572 Beam2d :: computeConsistentMassMatrix(FloatMatrix &answer, TimeStep *tStep, double &mass, const double *ipDensity)
573 {
574  // computes mass matrix of the receiver
575 
576  FloatMatrix stiff;
577  GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
578 
579  /*
580  * StructuralElement::computeMassMatrix(answer, tStep);
581  * answer.times(this->giveCrossSection()->give('A'));
582  */
583  double l = this->computeLength();
584  double kappa = this->giveKappaCoeff(tStep);
585  double kappa2 = kappa * kappa;
586 
587  double density = this->giveStructuralCrossSection()->give('d', gp); // constant density assumed
588  if ( ipDensity != NULL ) {
589  // Override density if desired
590  density = * ipDensity;
591  }
592 
593  double area = this->giveCrossSection()->give(CS_Area, gp); // constant area assumed
594  double c2 = ( area * density ) / ( ( 1. + 2. * kappa ) * ( 1. + 2. * kappa ) );
595  double c1 = ( area * density );
596 
597  answer.resize(6, 6);
598  answer.zero();
599 
600  answer.at(1, 1) = c1 * l / 3.0;
601  answer.at(1, 4) = c1 * l / 6.0;
602  answer.at(2, 2) = c2 * l * ( 13. / 35. + 7. * kappa / 5. + 4. * kappa2 / 3. );
603  answer.at(2, 3) = -c2 * l * l * ( 11. / 210. + kappa * 11. / 60. + kappa2 / 6. );
604  answer.at(2, 5) = c2 * l * ( 9. / 70. + kappa * 3. / 5. + kappa2 * 2. / 3. );
605  answer.at(2, 6) = c2 * l * l * ( 13. / 420. + kappa * 3. / 20. + kappa2 / 6. );
606  answer.at(3, 3) = c2 * l * l * l * ( 1. / 105. + kappa / 30. + kappa2 / 30. );
607  answer.at(3, 5) = -c2 * l * l * ( 13. / 420. + kappa * 3. / 20. + kappa2 / 6. );
608  answer.at(3, 6) = -c2 * l * l * l * ( 1. / 140. + kappa / 30. + kappa2 / 30. );
609 
610  answer.at(4, 4) = c1 * l / 3.0;
611  answer.at(5, 5) = c2 * l * ( 13. / 35. + kappa * 7. / 5. + kappa2 * 4. / 3. );
612  answer.at(5, 6) = c2 * l * l * ( 11. / 210. + kappa * 11. / 60. + kappa2 / 6. );
613  answer.at(6, 6) = c2 * l * l * l * ( 1. / 105. + kappa / 30. + kappa2 / 30. );
614 
615  answer.symmetrized();
616 
617  mass = l * area * density;
618 }
619 
620 
621 void
623 {
624  // computes initial stress matrix of receiver (or geometric stiffness matrix)
625 
626  FloatMatrix stiff;
627  FloatArray endForces;
628 
629  double l = this->computeLength();
630  double kappa = this->giveKappaCoeff(tStep);
631  double kappa2 = kappa * kappa;
632  double N;
633 
634  answer.resize(6, 6);
635  answer.zero();
636 
637  answer.at(2, 2) = 4. * kappa2 + 4. * kappa + 6. / 5.;
638  answer.at(2, 3) = -l / 10.;
639  answer.at(2, 5) = -4. * kappa2 - 4. * kappa - 6. / 5.;
640  answer.at(2, 6) = -l / 10.;
641  answer.at(3, 3) = l * l * ( kappa2 / 3. + kappa / 3. + 2. / 15. );
642  answer.at(3, 5) = l / 10.;
643  answer.at(3, 6) = -l * l * ( kappa2 / 3. + kappa / 3. + 1. / 30. );
644 
645  answer.at(5, 5) = 4. * kappa2 + 4. * kappa + 6. / 5.;
646  answer.at(5, 6) = l / 10.;
647  answer.at(6, 6) = l * l * ( kappa2 / 3. + kappa / 3. + 2. / 15. );
648 
649  answer.at(1, 1) = min( fabs( answer.at(2, 2) ), fabs( answer.at(3, 3) ) ) / 1000.;
650  answer.at(1, 4) = -answer.at(1, 1);
651  answer.at(4, 4) = answer.at(1, 1);
652 
653  answer.symmetrized();
654  // ask end forces in g.c.s
655  this->giveEndForcesVector(endForces, tStep);
656 
657  N = ( -endForces.at(1) + endForces.at(4) ) / 2.;
658  answer.times( N / ( l * ( 1. + 2. * kappa ) * ( 1. + 2. * kappa ) ) );
659 
660  //answer.beLumpedOf (mass);
661 }
662 
663 
664 #ifdef __OOFEG
666 {
667  GraphicObj *go;
668 
669  if ( !gc.testElementGraphicActivity(this) ) {
670  return;
671  }
672 
673  // if (!go) { // create new one
674  WCRec p [ 2 ]; /* poin */
675  EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
676  EASValsSetColor( gc.getElementColor() );
677  EASValsSetLayer(OOFEG_RAW_GEOMETRY_LAYER);
678  p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1);
679  p [ 0 ].y = 0.;
680  p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveCoordinate(3);
681  p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1);
682  p [ 1 ].y = 0.;
683  p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveCoordinate(3);
684  go = CreateLine3D(p);
685  EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
686  EGAttachObject(go, ( EObjectP ) this);
687  EMAddGraphicsToModel(ESIModel(), go);
688 }
689 
690 
692 {
693  GraphicObj *go;
694 
695  if ( !gc.testElementGraphicActivity(this) ) {
696  return;
697  }
698 
699  double defScale = gc.getDefScale();
700  // if (!go) { // create new one
701  WCRec p [ 2 ]; /* poin */
702  EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
703  EASValsSetColor( gc.getDeformedElementColor() );
704  EASValsSetLayer(OOFEG_DEFORMED_GEOMETRY_LAYER);
705  p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(1, tStep, defScale);
706  p [ 0 ].y = 0.;
707  p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(3, tStep, defScale);
708 
709  p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(1, tStep, defScale);
710  p [ 1 ].y = 0.;
711  p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(3, tStep, defScale);
712  go = CreateLine3D(p);
713  EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
714  EMAddGraphicsToModel(ESIModel(), go);
715 }
716 #endif
717 } // end namespace oofem
CrossSection * giveCrossSection()
Definition: element.C:495
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
Definition: beam2d.C:129
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.
void subtract(const FloatArray &src)
Subtracts array src to receiver.
Definition: floatarray.C:258
bool hasDofs2Condense()
Definition: beam2d.h:169
virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &)
Computes interpolation matrix for element unknowns.
Definition: beam2d.C:142
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
virtual Interface * giveInterface(InterfaceType)
Interface requesting service.
Definition: beam2d.C:91
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
void computeVectorOf(ValueModeType u, TimeStep *tStep, FloatArray &answer)
Returns local vector of unknowns.
Definition: element.C:86
virtual IntegrationRule * giveDefaultIntegrationRulePtr()
Access method for default integration rule.
Definition: element.h:822
Bottom z coordinate.
Definition: crosssection.h:72
double kappa
Definition: beam2d.h:65
Domain * domain
Link to domain object, useful for communicating with other FEM components.
Definition: femcmpnn.h:82
Class implementing internal element dof manager having some DOFs.
virtual bool computeGtoLRotationMatrix(FloatMatrix &answer)
Returns transformation matrix from global c.s.
Definition: beam2d.C:207
virtual void computeBoundaryEdgeLoadVector(FloatArray &answer, BoundaryLoad *load, int edge, CharType type, ValueModeType mode, TimeStep *tStep, bool global=true)
Computes the contribution of the given load at the given boundary edge.
Definition: beam2d.C:444
Load is specified in global c.s.
Definition: load.h:70
virtual void computeStrainVectorInLayer(FloatArray &answer, const FloatArray &masterGpStrain, GaussPoint *masterGp, GaussPoint *slaveGp, TimeStep *tStep)
Computes full 3D strain vector in element layer.
Definition: beam2d.C:271
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
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: beam2d.C:200
#define OOFEG_RAW_GEOMETRY_LAYER
This class implements a structural material status information.
Definition: structuralms.h:65
void clear()
Clears receiver (zero size).
Definition: floatarray.h:206
virtual void giveDofManDofIDMask(int inode, IntArray &) const
Returns dofmanager dof mask for node.
Definition: beam2d.C:292
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]
virtual int computeNumberOfGlobalDofs()
Computes the total number of element&#39;s global dofs.
Definition: beam2d.h:98
EngngModel * giveEngngModel()
Returns engineering model to which receiver is associated.
Definition: domain.C:433
virtual CoordSystType giveCoordSystMode()
Returns receiver&#39;s coordinate system.
Definition: boundaryload.h:151
virtual bool hasField(InputFieldType id)=0
Returns true if record contains field identified by idString keyword.
double pitch
Definition: beam2d.h:65
Top z coordinate.
Definition: crosssection.h:71
virtual double computeLength()
Computes the length (zero for all but 1D geometries)
Definition: beam2d.C:301
int numberOfCondensedDofs
number of condensed DOFs
Definition: beam2d.h:74
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 computeLocalForceLoadVector(FloatArray &answer, TimeStep *tStep, ValueModeType mode)
Computes element end force vector from applied loading in local coordinate system.
#define OOFEG_DEFORMED_GEOMETRY_LAYER
Beam2d(int n, Domain *aDomain)
Definition: beam2d.C:66
double givePitch()
Definition: beam2d.C:319
void rotatedWith(FloatMatrix &r, char mode)
Returns the receiver a rotated according the change-of-base matrix r.
Definition: floatarray.C:799
virtual FEInterpolation * giveInterpolation() const
Definition: beam2d.C:87
Class representing a general abstraction for finite element interpolation class.
Definition: feinterpol.h:132
Class representing "master" degree of freedom.
Definition: masterdof.h:92
virtual void computeValues(FloatArray &answer, TimeStep *tStep, const FloatArray &coords, const IntArray &dofids, ValueModeType mode)
Computes components values for specified dof ids.
Definition: load.C:57
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. ...
Abstract base class representing a boundary load (force, momentum, ...) that acts directly on a bound...
Definition: boundaryload.h:110
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
Computes numerically stiffness matrix of receiver.
Definition: beam2d.C:176
virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
Computes constitutive matrix of receiver.
Definition: beam2d.C:193
DofManager * ghostNodes[2]
Ghost nodes are used to introduce additional DOFs at element.
Definition: beam2d.h:72
#define OOFEM_ERROR(...)
Definition: error.h:61
REGISTER_Element(LSpace)
DofIDItem
Type representing particular dof type.
Definition: dofiditem.h:86
virtual int computeLoadGToLRotationMtrx(FloatMatrix &answer)
Returns transformation matrix from global coordinate system to local element coordinate system for el...
#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 double giveUpdatedCoordinate(int ic, TimeStep *tStep, double scale=1.)
Returns updated ic-th coordinate of receiver.
Definition: node.C:245
#define N(p, q)
Definition: mdm.C:367
void plusProductSymmUpper(const FloatMatrix &a, const FloatMatrix &b, double dV)
Adds to the receiver the product .
Definition: floatmatrix.C:698
double at(int i, int j) const
Coefficient access function.
Definition: floatmatrix.h:176
void appendDof(Dof *dof)
Adds the given Dof into the receiver.
Definition: dofmanager.C:134
virtual int setupIntegrationPoints(IntegrationRule &irule, int npoints, Element *element)
Sets up integration rule for the given element.
Definition: crosssection.C:54
int numberOfGaussPoints
Number of integration points as specified by nip.
Definition: element.h:188
virtual void printOutputAt(FILE *file, TimeStep *tStep)
Prints output of receiver to stream, for given time step.
Definition: beam2d.C:541
Class representing vector of real numbers.
Definition: floatarray.h:82
Implementation of matrix containing floating point numbers.
Definition: floatmatrix.h:94
virtual void giveEndForcesVector(FloatArray &answer, TimeStep *tStep)
Definition: beam2d.C:427
IRResultType
Type defining the return values of InputRecord reading operations.
Definition: irresulttype.h:47
virtual void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord=0)
Returns equivalent nodal forces vectors.
Definition: beam2d.C:420
virtual double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
Definition: beam2d.C:264
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
CharType
Definition: chartype.h:87
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 giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord=0)
Returns equivalent nodal forces vectors.
Class Interface.
Definition: interface.h:82
virtual ~Beam2d()
Definition: beam2d.C:80
#define OOFEG_DEFORMED_GEOMETRY_WIDTH
void times(double s)
Multiplies receiver with scalar.
Definition: floatarray.C:818
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
double length
Definition: beam2d.h:65
static FEI2dLineLin interp_geom
Definition: beam2d.h:76
void zero()
Zeroes all coefficient of receiver.
Definition: floatmatrix.C:1326
Domain * giveDomain() const
Definition: femcmpnn.h:100
virtual void computeConsistentMassMatrix(FloatMatrix &answer, TimeStep *tStep, double &mass, const double *ipDensity=NULL)
Computes consistent mass matrix of receiver using numerical integration over element volume...
Definition: beam2d.C:572
InterfaceType
Enumerative type, used to identify interface type.
Definition: interfacetype.h:43
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: beam2d.C:381
virtual void computeInitialStressMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes initial stress matrix for linear stability problem.
Definition: beam2d.C:622
int min(int i, int j)
Returns smaller value from two given decimals.
Definition: mathfem.h:59
Load is base abstract class for all loads.
Definition: load.h:61
The element interface required by LayeredCrossSection.
void beProductOf(const FloatMatrix &a, const FloatMatrix &b)
Assigns to the receiver product of .
Definition: floatmatrix.C:337
virtual TimeStep * giveCurrentStep(bool force=false)
Returns current time step.
Definition: engngm.h:683
int giveSize() const
Definition: intarray.h:203
int giveLabel() const
Definition: element.h:1055
the oofem namespace is to define a context or scope in which all oofem names are defined.
bool isNotEmpty() const
Returns true if receiver is not empty.
Definition: floatarray.h:220
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
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in full form.
Definition: beam2d.C:499
#define IR_GIVE_FIELD(__ir, __value, __id)
Macro facilitating the use of input record reading methods.
Definition: inputrecord.h:69
int giveNumber() const
Definition: femcmpnn.h:107
Node * giveNode(int i) const
Returns reference to the i-th node of element.
Definition: element.h:610
static FEI2dLineHermite interp_beam
Definition: beam2d.h:77
virtual void computeBmatrixAt(GaussPoint *, FloatMatrix &, int=1, int=ALL_STRAINS)
Computes the geometrical matrix of receiver in given integration point.
Definition: beam2d.C:102
void symmetrized()
Initializes the lower half of the receiver according to the upper half.
Definition: floatmatrix.C:1576
#define _IFT_Beam2d_dofstocondense
Definition: beam2d.h:45
virtual int giveLocalCoordinateSystem(FloatMatrix &answer)
Returns local coordinate system of receiver Required by material models with ortho- and anisotrophy...
Definition: beam2d.C:356
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
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
#define OOFEM_WARNING(...)
Definition: error.h:62
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
virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
Definition: beam2d.C:665
This class implements a base beam intented to be a base class for beams based on lagrangian interpola...
virtual void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType)
Definition: beam2d.C:691
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: beam2d.C:490
Class representing Gaussian-quadrature integration rule.
void resize(int s)
Resizes receiver towards requested size.
Definition: floatarray.C:631
double giveKappaCoeff(TimeStep *tStep)
Definition: beam2d.C:339
void plusProduct(const FloatMatrix &b, const FloatArray &s, double dV)
Adds the product .
Definition: floatarray.C:226

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:27 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011