OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
libeam3dnl2.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/libeam3dnl2.h"
36 #include "../sm/CrossSections/structuralcrosssection.h"
37 #include "../sm/Materials/structuralms.h"
38 #include "node.h"
39 #include "material.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 "timestep.h"
47 #include "contextioerr.h"
48 #include "classfactory.h"
49 
50 #ifdef __OOFEG
51  #include "oofeggraphiccontext.h"
52 #endif
53 
54 namespace oofem {
55 REGISTER_Element(LIBeam3dNL2);
56 
57 LIBeam3dNL2 :: LIBeam3dNL2(int n, Domain *aDomain) : NLStructuralElement(n, aDomain), q(4), tempQ(4) //, kappa (3)
58 {
59  numberOfDofMans = 2;
60  l0 = 0.;
61  tempQCounter = 0;
62  referenceNode = 0;
63  // init kappa vector at centre
64  // kappa.zero();
65 }
66 
67 
68 void
70 {
71  if ( vec.giveSize() != 3 ) {
72  OOFEM_ERROR("vec param size mismatch");
73  }
74 
75  answer.resize(3, 3);
76 
77  answer.at(1, 1) = answer.at(2, 2) = answer.at(3, 3) = 0.;
78  answer.at(1, 2) = -vec.at(3);
79  answer.at(1, 3) = vec.at(2);
80  answer.at(2, 1) = vec.at(3);
81  answer.at(2, 3) = -vec.at(1);
82  answer.at(3, 1) = -vec.at(2);
83  answer.at(3, 2) = vec.at(1);
84 }
85 
86 
87 void
89 {
90  FloatMatrix S(3, 3), SS(3, 3);
91  double psiSize;
92 
93  if ( psi.giveSize() != 3 ) {
94  OOFEM_ERROR("psi param size mismatch");
95  }
96 
97  answer.resize(3, 3);
98  answer.zero();
99 
100  psiSize = psi.computeNorm();
101  answer.at(1, 1) = answer.at(2, 2) = answer.at(3, 3) = 1.;
102 
103  if ( psiSize <= 1.e-40 ) {
104  return;
105  }
106 
107  this->computeSMtrx(S, psi);
108  SS.beProductOf(S, S);
109  S.times(sin(psiSize) / psiSize);
110  SS.times( ( 1. - cos(psiSize) ) / ( psiSize * psiSize ) );
111 
112  answer.add(S);
113  answer.add(SS);
114 }
115 
116 
117 void
119 {
120  // test if not previously updated temporary quaternion
121  if ( tStep->giveSolutionStateCounter() != tempQCounter ) {
122  // update temporary quaternion
123  FloatArray u, centreSpin(3), q2(4);
124  double centreSpinSize;
125 
126  // ask element's displacement increments
127  this->computeVectorOf(VM_Incremental, tStep, u);
128 
129  // interpolate spin at the centre
130  centreSpin.at(1) = 0.5 * ( u.at(4) + u.at(10) );
131  centreSpin.at(2) = 0.5 * ( u.at(5) + u.at(11) );
132  centreSpin.at(3) = 0.5 * ( u.at(6) + u.at(12) );
133 
134  centreSpinSize = centreSpin.computeNorm();
135  if ( centreSpinSize > 1.e-30 ) {
136  centreSpin.normalize();
137  q2.at(1) = sin(centreSpinSize / 2.) * centreSpin.at(1);
138  q2.at(2) = sin(centreSpinSize / 2.) * centreSpin.at(2);
139  q2.at(3) = sin(centreSpinSize / 2.) * centreSpin.at(3);
140  q2.at(4) = cos(centreSpinSize / 2.);
141 
142  // update temporary quaternion at the center
143  tempQ.at(1) = q.at(4) * q2.at(1) + q2.at(4) * q.at(1) - ( q.at(2) * q2.at(3) - q.at(3) * q2.at(2) );
144  tempQ.at(2) = q.at(4) * q2.at(2) + q2.at(4) * q.at(2) - ( q.at(3) * q2.at(1) - q.at(1) * q2.at(3) );
145  tempQ.at(3) = q.at(4) * q2.at(3) + q2.at(4) * q.at(3) - ( q.at(1) * q2.at(2) - q.at(2) * q2.at(1) );
146  tempQ.at(4) = q.at(4) * q2.at(4) - ( q.at(1) * q2.at(1) + q.at(2) * q2.at(2) + q.at(3) * q2.at(3) );
147  } else {
148  tempQ = q;
149  }
150 
151  // remember timestamp
153  }
154 }
155 
156 
157 void
159 {
160  answer.resize(3, 3);
161 
162  answer.at(1, 1) = q.at(4) * q.at(4) + q.at(1) * q.at(1) - 0.5;
163  answer.at(1, 2) = q.at(1) * q.at(2) - q.at(3) * q.at(4);
164  answer.at(1, 3) = q.at(1) * q.at(3) + q.at(2) * q.at(4);
165 
166  answer.at(2, 1) = q.at(2) * q.at(1) + q.at(3) * q.at(4);
167  answer.at(2, 2) = q.at(4) * q.at(4) + q.at(2) * q.at(2) - 0.5;
168  answer.at(2, 3) = q.at(2) * q.at(3) - q.at(1) * q.at(4);
169 
170  answer.at(3, 1) = q.at(3) * q.at(1) - q.at(2) * q.at(4);
171  answer.at(3, 2) = q.at(3) * q.at(2) + q.at(1) * q.at(4);
172  answer.at(3, 3) = q.at(4) * q.at(4) + q.at(3) * q.at(3) - 0.5;
173 
174  answer.times(2.);
175 }
176 
177 
178 void
180 {
181  // Spurrier's algorithm
182 
183  int i, ii;
184  double a, trR;
185 
186  answer.resize(4);
187 
188  trR = R.at(1, 1) + R.at(2, 2) + R.at(3, 3);
189  a = trR;
190  ii = 0;
191  for ( i = 1; i <= 3; i++ ) {
192  if ( R.at(i, i) > a ) {
193  a = R.at(i, i);
194  ii = i;
195  }
196  }
197 
198  if ( a == trR ) {
199  //printf (".");
200  answer.at(4) = 0.5 * sqrt(1. + a);
201  answer.at(1) = ( R.at(3, 2) - R.at(2, 3) ) / ( 4. * answer.at(4) );
202  answer.at(2) = ( R.at(1, 3) - R.at(3, 1) ) / ( 4. * answer.at(4) );
203  answer.at(3) = ( R.at(2, 1) - R.at(1, 2) ) / ( 4. * answer.at(4) );
204  } else {
205  //printf (":");
206  int jj, kk;
207  if ( ii == 1 ) {
208  jj = 2;
209  kk = 3;
210  } else if ( ii == 2 ) {
211  jj = 3;
212  kk = 1;
213  } else {
214  jj = 1;
215  kk = 2;
216  }
217 
218  answer.at(ii) = sqrt( 0.5 * a + 0.25 * ( 1. - trR ) );
219  answer.at(4) = 0.25 * ( R.at(kk, jj) - R.at(jj, kk) ) / answer.at(ii);
220  answer.at(jj) = 0.25 * ( R.at(jj, ii) + R.at(ii, jj) ) / answer.at(ii);
221  answer.at(kk) = 0.25 * ( R.at(kk, ii) + R.at(ii, kk) ) / answer.at(ii);
222  }
223 }
224 
225 
226 void
228 {
229  FloatArray xd(3), eps(3), curv(3);
230  FloatMatrix tempTc;
231 
232  // update temp triad
233  this->updateTempQuaternion(tStep);
234  this->computeRotMtrxFromQuaternion(tempTc, this->tempQ);
235 
236  this->computeXdVector(xd, tStep);
237 
238  // compute eps_xl, gamma_l2, gamma_l3
239  eps.beTProductOf(tempTc, xd);
240  eps.times(1. / this->l0);
241  eps.at(1) -= 1.0;
242 
243  // update curvature at midpoint
244  this->computeTempCurv(curv, tStep);
245 
246  answer.resize(6);
247  answer.at(1) = eps.at(1); // eps_xl
248  answer.at(2) = eps.at(2); // gamma_l2
249  answer.at(3) = eps.at(3); // gamma_l3
250  answer.at(4) = curv.at(1); // kappa_1
251  answer.at(5) = curv.at(2); // kappa_2
252  answer.at(6) = curv.at(3); // kappa_3
253 }
254 
255 
256 void
258 {
259  FloatArray xd(3);
260  FloatMatrix s(3, 3);
261 
262  this->computeXdVector(xd, tStep);
263  this->computeSMtrx(s, xd);
264 
265  answer.resize(12, 6);
266  answer.zero();
267 
268  for ( int i = 1; i < 4; i++ ) {
269  answer.at(i, i) = -1.0;
270  answer.at(i + 6, i) = 1.0;
271  answer.at(i + 3, i + 3) = -1.0;
272  answer.at(i + 9, i + 3) = 1.0;
273 
274  for ( int j = 1; j < 4; j++ ) {
275  answer.at(i + 3, j) = answer.at(i + 9, j) = 0.5 * s.at(j, i);
276  }
277  }
278 }
279 
280 
281 void
282 LIBeam3dNL2 :: giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord)
283 {
285  FloatArray nm(6), stress, strain;
286  FloatMatrix x, tempTc;
287  double s1, s2;
288 
289  // update temp triad
290  this->updateTempQuaternion(tStep);
291  this->computeRotMtrxFromQuaternion(tempTc, this->tempQ);
292 
293  if ( useUpdatedGpRecord == 1 ) {
294  stress = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() )->giveStressVector();
295  } else {
296  this->computeStrainVector(strain, gp, tStep);
297  this->computeStressVector(stress, strain, gp, tStep);
298  }
299 
300  for ( int i = 1; i <= 3; i++ ) {
301  s1 = s2 = 0.0;
302  for ( int j = 1; j <= 3; j++ ) {
303  s1 += tempTc.at(i, j) * stress.at(j);
304  s2 += tempTc.at(i, j) * stress.at(j + 3);
305  }
306 
307  nm.at(i) = s1;
308  nm.at(i + 3) = s2;
309  }
310 
311  this->computeXMtrx(x, tStep);
312  answer.beProductOf(x, nm);
313 }
314 
315 
316 void
318 {
319  FloatArray u(3);
320 
321  answer.resize(3);
322  // ask element's displacements
323  this->computeVectorOf(VM_Total, tStep, u);
324 
325  answer.at(1) = ( this->giveNode(2)->giveCoordinate(1) + u.at(7) ) -
326  ( this->giveNode(1)->giveCoordinate(1) + u.at(1) );
327  answer.at(2) = ( this->giveNode(2)->giveCoordinate(2) + u.at(8) ) -
328  ( this->giveNode(1)->giveCoordinate(2) + u.at(2) );
329  answer.at(3) = ( this->giveNode(2)->giveCoordinate(3) + u.at(9) ) -
330  ( this->giveNode(1)->giveCoordinate(3) + u.at(3) );
331 }
332 
333 
334 void
336 {
337  double s1, s2;
338  FloatMatrix d, x, xt(12, 6), dxt, sn, sm, sxd, y, tempTc;
339  FloatArray n(3), m(3), xd(3), stress, strain;
341 
342  answer.clear();
343 
344  // linear part
345 
346  this->updateTempQuaternion(tStep);
347  this->computeRotMtrxFromQuaternion(tempTc, this->tempQ);
348  this->computeXMtrx(x, tStep);
349  xt.zero();
350  for ( int i = 1; i <= 12; i++ ) {
351  for ( int j = 1; j <= 3; j++ ) {
352  for ( int k = 1; k <= 3; k++ ) {
353  // compute x*Tbar, taking into account sparsity of Tbar
354  xt.at(i, j) += x.at(i, k) * tempTc.at(k, j);
355  xt.at(i, j + 3) += x.at(i, k + 3) * tempTc.at(k, j);
356  }
357  }
358  }
359 
360  this->computeConstitutiveMatrixAt(d, rMode, gp, tStep);
361  dxt.beProductTOf(d, xt);
362  answer.beProductOf(xt, dxt);
363  answer.times(1. / this->l0);
364 
365  // geometric stiffness ks = ks1+ks2
366  // ks1
367  this->computeStrainVector(strain, gp, tStep);
368  this->computeStressVector(stress, strain, gp, tStep);
369 
370  for ( int i = 1; i <= 3; i++ ) {
371  s1 = s2 = 0.0;
372  for ( int j = 1; j <= 3; j++ ) {
373  s1 += tempTc.at(i, j) * stress.at(j);
374  s2 += tempTc.at(i, j) * stress.at(j + 3);
375  }
376 
377  n.at(i) = s1;
378  m.at(i) = s2;
379  }
380 
381  this->computeSMtrx(sn, n);
382  this->computeSMtrx(sm, m);
383 
384  for ( int i = 1; i <= 3; i++ ) {
385  for ( int j = 1; j <= 3; j++ ) {
386  answer.at(i, j + 3) += sn.at(i, j);
387  answer.at(i, j + 9) += sn.at(i, j);
388  answer.at(i + 3, j + 3) += sm.at(i, j);
389  answer.at(i + 3, j + 9) += sm.at(i, j);
390 
391  answer.at(i + 6, j + 3) -= sn.at(i, j);
392  answer.at(i + 6, j + 9) -= sn.at(i, j);
393  answer.at(i + 9, j + 3) -= sm.at(i, j);
394  answer.at(i + 9, j + 9) -= sm.at(i, j);
395  }
396  }
397 
398  // ks2
399  this->computeXdVector(xd, tStep);
400  this->computeSMtrx(sxd, xd);
401 
402  y.beProductOf(sxd, sn);
403  y.times(0.5);
404 
405  for ( int i = 1; i <= 3; i++ ) {
406  for ( int j = 1; j <= 3; j++ ) {
407  answer.at(i + 3, j) -= sn.at(i, j);
408  answer.at(i + 3, j + 3) += y.at(i, j);
409  answer.at(i + 3, j + 6) += sn.at(i, j);
410  answer.at(i + 3, j + 9) += y.at(i, j);
411 
412  answer.at(i + 9, j) -= sn.at(i, j);
413  answer.at(i + 9, j + 3) += y.at(i, j);
414  answer.at(i + 9, j + 6) += sn.at(i, j);
415  answer.at(i + 9, j + 9) += y.at(i, j);
416  }
417  }
418 }
419 
420 
421 void
423 // Sets up the array of Gauss Points of the receiver.
424 {
425  if ( integrationRulesArray.size() == 0 ) {
426  integrationRulesArray.resize( 1 );
427  integrationRulesArray [ 0 ].reset( new GaussIntegrationRule(1, this, 1, 2) );
429  }
430 }
431 
432 
433 void
435 {
436  this->giveStructuralCrossSection()->give3dBeamStiffMtrx(answer, rMode, gp, tStep);
437 }
438 
439 
440 void
442 {
443  this->giveStructuralCrossSection()->giveGeneralizedStress_Beam3d(answer, gp, strain, tStep);
444 }
445 
446 
449 {
450  IRResultType result; // Required by IR_GIVE_FIELD macro
451 
452  // first call parent
454  if ( result != IRRT_OK ) {
455  return result;
456  }
457 
459  if ( referenceNode == 0 ) {
460  OOFEM_ERROR("wrong reference node specified");
461  }
462 
463  /*
464  * if (this->hasString (initString, "dofstocondense")) {
465  * dofsToCondense = this->ReadIntArray (initString, "dofstocondense");
466  * if (dofsToCondense->giveSize() >= 12)
467  * OOFEM_ERROR("wrong input data for condensed dofs");
468  * } else {
469  * dofsToCondense = NULL;
470  * }
471  */
472 
474  // compute initial triad at centre - requires nodal coordinates
475  FloatMatrix lcs, tc;
476  this->giveLocalCoordinateSystem(lcs);
477  tc.beTranspositionOf(lcs);
478 
479  this->computeQuaternionFromRotMtrx(q, tc);
480  return IRRT_OK;
481 }
482 
483 
484 double
486 // Returns the original length (l0) of the receiver.
487 {
488  double dx, dy, dz;
489  Node *nodeA, *nodeB;
490 
491  if ( l0 == 0. ) {
492  nodeA = this->giveNode(1);
493  nodeB = this->giveNode(2);
494  dx = nodeB->giveCoordinate(1) - nodeA->giveCoordinate(1);
495  dy = nodeB->giveCoordinate(2) - nodeA->giveCoordinate(2);
496  dz = nodeB->giveCoordinate(3) - nodeA->giveCoordinate(3);
497  l0 = sqrt(dx * dx + dy * dy + dz * dz);
498  }
499 
500  return l0;
501 }
502 
503 
504 void
506 // Returns the lumped mass matrix of the receiver. This expression is
507 // valid in both local and global axes.
508 {
509  GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
510  double density = this->giveStructuralCrossSection()->give('d', gp);
511  double halfMass = density * this->giveCrossSection()->give(CS_Area, gp) * this->computeLength() / 2.;
512  answer.resize(12, 12);
513  answer.zero();
514  answer.at(1, 1) = answer.at(2, 2) = answer.at(3, 3) = halfMass;
515  answer.at(7, 7) = answer.at(8, 8) = answer.at(9, 9) = halfMass;
516 }
517 
518 
519 void
521 // Returns the displacement interpolation matrix {N} of the receiver, eva-
522 // luated at gp.
523 {
524  double ksi, n1, n2;
525 
526  ksi = iLocCoord.at(1);
527  n1 = ( 1. - ksi ) * 0.5;
528  n2 = ( 1. + ksi ) * 0.5;
529 
530  answer.resize(6, 12);
531  answer.zero();
532 
533  // u
534  answer.at(1, 1) = n1;
535  answer.at(1, 7) = n2;
536  // v
537  answer.at(2, 2) = n1;
538  answer.at(2, 8) = n2;
539  // w
540  answer.at(3, 3) = n1;
541  answer.at(3, 9) = n2;
542  // fi_x
543  answer.at(4, 4) = n1;
544  answer.at(4, 10) = n2;
545  // fi_y
546  answer.at(5, 5) = n1;
547  answer.at(5, 11) = n2;
548  // fi_z
549  answer.at(6, 6) = n1;
550  answer.at(6, 12) = n2;
551 }
552 
553 
554 double
556 // Returns the length of the receiver. This method is valid only if 1
557 // Gauss point is used.
558 {
559  double weight = gp->giveWeight();
560  return weight * 0.5 * this->computeLength();
561 }
562 
563 
564 void
566 {
567  answer = {D_u, D_v, D_w, R_u, R_v, R_w};
568 }
569 
570 int
572 {
573  double ksi, n1, n2;
574 
575  ksi = lcoords.at(1);
576  n1 = ( 1. - ksi ) * 0.5;
577  n2 = ( 1. + ksi ) * 0.5;
578 
579  answer.resize(3);
580  answer.at(1) = n1 * this->giveNode(1)->giveCoordinate(1) + n2 *this->giveNode(2)->giveCoordinate(1);
581  answer.at(2) = n1 * this->giveNode(1)->giveCoordinate(2) + n2 *this->giveNode(2)->giveCoordinate(2);
582  answer.at(3) = n1 * this->giveNode(1)->giveCoordinate(3) + n2 *this->giveNode(2)->giveCoordinate(3);
583 
584  return 1;
585 }
586 
587 
588 void
590 {
591  /*
592  * provides dof mapping of local edge dofs (only nonzero are taken into account)
593  * to global element dofs
594  */
595  if ( iEdge != 1 ) {
596  OOFEM_ERROR("wrong edge number");
597  }
598 
599  answer.resize(12);
600  for ( int i = 1; i <= 12; i++ ) {
601  answer.at(i) = i;
602  }
603 }
604 
605 
606 double
608 {
609  if ( iEdge != 1 ) { // edge between nodes 1 2
610  OOFEM_ERROR("wrong egde number");
611  }
612 
613  double weight = gp->giveWeight();
614  return 0.5 * this->computeLength() * weight;
615 }
616 
617 
618 int
620 //
621 // returns a unit vectors of local coordinate system at element
622 // stored rowwise (mainly used by some materials with ortho and anisotrophy)
623 //
624 {
625  FloatArray lx(3), ly(3), lz(3), help(3);
626  double length = this->computeLength();
627  Node *nodeA, *nodeB, *refNode;
628 
629  answer.resize(3, 3);
630  answer.zero();
631  nodeA = this->giveNode(1);
632  nodeB = this->giveNode(2);
633  refNode = this->giveDomain()->giveNode(this->referenceNode);
634 
635  for ( int i = 1; i <= 3; i++ ) {
636  lx.at(i) = ( nodeB->giveCoordinate(i) - nodeA->giveCoordinate(i) ) / length;
637  help.at(i) = ( refNode->giveCoordinate(i) - nodeA->giveCoordinate(i) );
638  }
639 
640  lz.beVectorProductOf(lx, help);
641  lz.normalize();
642  ly.beVectorProductOf(lz, lx);
643  ly.normalize();
644 
645  for ( int i = 1; i <= 3; i++ ) {
646  answer.at(1, i) = lx.at(i);
647  answer.at(2, i) = ly.at(i);
648  answer.at(3, i) = lz.at(i);
649  }
650 
651  return 1;
652 }
653 
654 
655 int
657 {
658  /*
659  * Returns transformation matrix from global coordinate system to local
660  * element coordinate system for element load vector components.
661  * If no transformation is necessary, answer is empty matrix (default);
662  *
663  * Does not support follower load
664  */
665 
666  FloatMatrix lcs;
667 
668  answer.resize(6, 6);
669  answer.zero();
670 
671  this->giveLocalCoordinateSystem(lcs);
672 
673  for ( int i = 1; i <= 3; i++ ) {
674  for ( int j = 1; j <= 3; j++ ) {
675  answer.at(i, j) = lcs.at(i, j);
676  answer.at(3 + i, 3 + j) = lcs.at(i, j);
677  }
678  }
679 
680  return 1;
681 }
682 
683 
684 int
686 {
687  // returns transformation matrix from
688  // edge local coordinate system
689  // to element local c.s
690  // (same as global c.s in this case)
691  //
692  // i.e. f(element local) = T * f(edge local)
693  //
694  answer.clear();
695  return 0;
696 }
697 
698 
699 void
701 {
702  FloatArray lc(1);
703  NLStructuralElement :: computeBodyLoadVectorAt(answer, load, tStep, mode);
704  answer.times( this->giveCrossSection()->give(CS_Area, lc, this) );
705 }
706 
707 
708 void
710 // Updates the receiver at end of step.
711 {
713 
714  // update quaternion
715  this->updateTempQuaternion(tStep);
716  q = tempQ;
717 
718  // update curvature
719  //FloatArray curv;
720  //this->computeTempCurv (curv, tStep);
721  //kappa = curv;
722 }
723 
724 void
726 // initializes receiver to new time step or can be used
727 // if current time step must be restarted
728 {
730  tempQ = q;
731 }
732 
733 
734 void
736 {
738  FloatArray ui(3), ac(3), PrevEpsilon;
739 
740  answer.resize(3);
741 
742  /*
743  * // update curvature at midpoint
744  * // first, compute Tmid
745  * // ask increments
746  * this -> computeVectorOf(UnknownMode_Incremental,tStep, ui) ;
747  * ac.at(1) = 0.5*(ui.at(10) - ui.at(4));
748  * ac.at(2) = 0.5*(ui.at(11) - ui.at(5));
749  * ac.at(3) = 0.5*(ui.at(12) - ui.at(6));
750  * this->computeSMtrx (sc, ac);
751  * sc.times (1./2.);
752  * // compute I+sc
753  * sc.at(1,1) += 1.0;
754  * sc.at(2,2) += 1.0;
755  * sc.at(3,3) += 1.0;
756  *
757  * this->computeRotMtrxFromQuaternion(tc, this->q);
758  * tmid.beProductOf (sc, tc);
759  *
760  * // update curvature at centre
761  * ac.at(1) = (ui.at(10) - ui.at(4));
762  * ac.at(2) = (ui.at(11) - ui.at(5));
763  * ac.at(3) = (ui.at(12) - ui.at(6));
764  *
765  * answer.beTProductOf (tmid, ac);
766  * answer.times (1/this->l0);
767  * // ask for previous kappa
768  * PrevEpsilon = ((StructuralMaterialStatus*) mat->giveStatus(gp)) -> giveStrainVector ();
769  * if (PrevEpsilon.giveSize()) {
770  * answer.at(1) += PrevEpsilon.at(4);
771  * answer.at(2) += PrevEpsilon.at(5);
772  * answer.at(3) += PrevEpsilon.at(6);
773  * }
774  */
775  // update curvature
776  // exact procedure due to Simo & Vu Quoc
777  FloatMatrix dR, Rn, Ro;
778  FloatArray om, omp, acp(3), kapgn1(3);
779  double acSize, coeff;
780 
781  this->computeVectorOf(VM_Incremental, tStep, ui);
782 
783  ac.at(1) = 0.5 * ( ui.at(10) - ui.at(4) );
784  ac.at(2) = 0.5 * ( ui.at(11) - ui.at(5) );
785  ac.at(3) = 0.5 * ( ui.at(12) - ui.at(6) );
786 
787  this->computeRotMtrx(dR, ac);
788  this->computeRotMtrxFromQuaternion(Ro, this->q);
789  Rn.beProductOf(dR, Ro);
790 
791  acSize = ac.computeSquaredNorm();
792 
793  if ( acSize > 1.e-30 ) {
794  FloatMatrix h(3, 3);
795  ac.normalize();
796  om = ac;
797  om.times( 2. * tan(acSize / 2.) );
798 
799  coeff = ( 1. - ( acSize / sin(acSize) ) );
800  for ( int i = 1; i <= 3; i++ ) {
801  for ( int j = 1; j <= 3; j++ ) {
802  h.at(i, j) = coeff * ac.at(i) * ac.at(j);
803  }
804  }
805 
806  for ( int i = 1; i <= 3; i++ ) {
807  h.at(i, i) = 1. - h.at(i, i);
808  }
809 
810  // compute dPsi/ds
811  acp.at(1) = ( ui.at(10) - ui.at(4) ) / this->l0;
812  acp.at(2) = ( ui.at(11) - ui.at(5) ) / this->l0;
813  acp.at(3) = ( ui.at(12) - ui.at(6) ) / this->l0;
814 
815  omp.beProductOf(h, acp);
816  omp.times(2. * tan(acSize / 2.) / acSize);
817 
818  kapgn1.beVectorProductOf(om, omp);
819  kapgn1.times(1. / 2.);
820  kapgn1.add(omp);
821 
822  coeff = 1. / ( 1. + 0.25 * om.computeSquaredNorm() );
823  kapgn1.times(coeff);
824 
825  answer.beTProductOf(Rn, kapgn1);
826  }
827 
828  // ask for previous kappa
829  StructuralMaterialStatus *matStat = static_cast< StructuralMaterialStatus * >( gp->giveMaterialStatus() );
830  if ( matStat ) {
831  PrevEpsilon = matStat->giveStrainVector();
832  }
833  if ( PrevEpsilon.giveSize() ) {
834  answer.at(1) += PrevEpsilon.at(4);
835  answer.at(2) += PrevEpsilon.at(5);
836  answer.at(3) += PrevEpsilon.at(6);
837  }
838 }
839 
840 
843 //
844 // saves full element context (saves state variables, that completely describe
845 // current state)
846 //
847 {
848  contextIOResultType iores;
849  if ( ( iores = NLStructuralElement :: saveContext(stream, mode, obj) ) != CIO_OK ) {
850  THROW_CIOERR(iores);
851  }
852 
853  if ( ( iores = q.storeYourself(stream) ) != CIO_OK ) {
854  THROW_CIOERR(iores);
855  }
856 
857  return CIO_OK;
858 }
859 
860 
863 //
864 // restores full element context (saves state variables, that completely describe
865 // current state)
866 //
867 {
868  contextIOResultType iores;
869  if ( ( iores = NLStructuralElement :: restoreContext(stream, mode, obj) ) != CIO_OK ) {
870  THROW_CIOERR(iores);
871  }
872 
873  if ( ( iores = q.restoreYourself(stream) ) != CIO_OK ) {
874  THROW_CIOERR(iores);
875  }
876 
877  return CIO_OK;
878 }
879 
880 
881 #ifdef __OOFEG
883 {
884  GraphicObj *go;
885 
886  if ( !gc.testElementGraphicActivity(this) ) {
887  return;
888  }
889 
890  // if (!go) { // create new one
891  WCRec p [ 2 ]; /* poin */
892  EASValsSetLineWidth(OOFEG_RAW_GEOMETRY_WIDTH);
893  EASValsSetColor( gc.getElementColor() );
894  EASValsSetLayer(OOFEG_RAW_GEOMETRY_LAYER);
895  p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveCoordinate(1);
896  p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveCoordinate(2);
897  p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveCoordinate(3);
898  p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveCoordinate(1);
899  p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveCoordinate(2);
900  p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveCoordinate(3);
901  go = CreateLine3D(p);
902  EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
903  EGAttachObject(go, ( EObjectP ) this);
904  EMAddGraphicsToModel(ESIModel(), go);
905 }
906 
907 
909 {
910  GraphicObj *go;
911 
912  if ( !gc.testElementGraphicActivity(this) ) {
913  return;
914  }
915 
916  double defScale = gc.getDefScale();
917  // if (!go) { // create new one
918  WCRec p [ 2 ]; /* poin */
919  const char *colors[] = {
920  "red", "green", "blue"
921  };
922 
923  EASValsSetLineWidth(OOFEG_DEFORMED_GEOMETRY_WIDTH);
924  EASValsSetColor( gc.getDeformedElementColor() );
925  EASValsSetLayer(OOFEG_DEFORMED_GEOMETRY_LAYER);
926  p [ 0 ].x = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(1, tStep, defScale);
927  p [ 0 ].y = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(2, tStep, defScale);
928  p [ 0 ].z = ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(3, tStep, defScale);
929 
930  p [ 1 ].x = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(1, tStep, defScale);
931  p [ 1 ].y = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(2, tStep, defScale);
932  p [ 1 ].z = ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(3, tStep, defScale);
933  go = CreateLine3D(p);
934  EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
935  EMAddGraphicsToModel(ESIModel(), go);
936 
937  // draw centre triad
938  FloatMatrix tc;
939  int i, succ;
940  double coeff = this->l0 / 3.;
941  this->computeRotMtrxFromQuaternion(tc, this->q);
942 
943  p [ 0 ].x = 0.5 * ( ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(1, tStep, defScale) + ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(1, tStep, defScale) );
944  p [ 0 ].y = 0.5 * ( ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(2, tStep, defScale) + ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(2, tStep, defScale) );
945  p [ 0 ].z = 0.5 * ( ( FPNum ) this->giveNode(1)->giveUpdatedCoordinate(3, tStep, defScale) + ( FPNum ) this->giveNode(2)->giveUpdatedCoordinate(3, tStep, defScale) );
946 
947  // draw t1
948  for ( i = 1; i <= 3; i++ ) {
949  p [ 1 ].x = p [ 0 ].x + coeff *tc.at(1, i);
950  p [ 1 ].y = p [ 0 ].y + coeff *tc.at(2, i);
951  p [ 1 ].z = p [ 0 ].z + coeff *tc.at(3, i);
952 
953  EASValsSetColor( ColorGetPixelFromString(const_cast< char * >(colors [ i - 1 ]), & succ) );
954 
955  go = CreateLine3D(p);
956  EGWithMaskChangeAttributes(WIDTH_MASK | COLOR_MASK | LAYER_MASK, go);
957  EMAddGraphicsToModel(ESIModel(), go);
958  }
959 }
960 #endif
961 } // 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 computeStrainVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep)
Compute strain vector of receiver evaluated at given integration point at given time step from elemen...
Definition: libeam3dnl2.C:227
void beVectorProductOf(const FloatArray &v1, const FloatArray &v2)
Computes vector product (or cross product) of vectors given as parameters, , and stores the result in...
Definition: floatarray.C:415
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
virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
Computes interpolation matrix for element unknowns.
Definition: libeam3dnl2.C:520
int referenceNode
Reference node.
Definition: libeam3dnl2.h:69
void computeQuaternionFromRotMtrx(FloatArray &answer, FloatMatrix &R)
Computes the normalized quaternion from the given rotation matrix.
Definition: libeam3dnl2.C:179
Abstract base class for "structural" finite elements with geometrical nonlinearities.
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: libeam3dnl2.C:441
virtual int computeLoadGToLRotationMtrx(FloatMatrix &answer)
Returns transformation matrix from global coordinate system to local element coordinate system for el...
Definition: libeam3dnl2.C:656
The purpose of DataStream abstract class is to allow to store/restore context to different streams...
Definition: datastream.h:54
contextIOResultType storeYourself(DataStream &stream) const
Definition: floatarray.C:872
double & at(int i)
Coefficient access function.
Definition: floatarray.h:131
StateCounterType tempQCounter
Time stamp of temporary centre quaternion.
Definition: libeam3dnl2.h:67
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 void drawDeformedGeometry(oofegGraphicContext &gc, TimeStep *tStep, UnknownType)
Definition: libeam3dnl2.C:908
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
virtual int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords)
Computes the global coordinates from given element&#39;s local coordinates.
Definition: libeam3dnl2.C:571
#define S(p)
Definition: mdm.C:481
FloatArray tempQ
Temporary quaternion at the center.
Definition: libeam3dnl2.h:63
virtual double giveCoordinate(int i)
Definition: node.C:82
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: libeam3dnl2.C:700
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.
#define THROW_CIOERR(e)
Definition: contextioerr.h:61
virtual double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
Definition: libeam3dnl2.C:555
#define OOFEG_DEFORMED_GEOMETRY_LAYER
virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Restores the receiver state previously written in stream.
Definition: element.C:970
virtual void giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord=0)
Evaluates nodal representation of real internal forces.
Definition: libeam3dnl2.C:282
virtual void computeLumpedMassMatrix(FloatMatrix &answer, TimeStep *tStep)
Computes lumped mass matrix of receiver.
Definition: libeam3dnl2.C:505
void computeSMtrx(FloatMatrix &answer, FloatArray &vec)
Evaluates the S matrix from given vector vec.
Definition: libeam3dnl2.C:69
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
Computes the stiffness matrix of receiver.
Definition: libeam3dnl2.C:335
virtual void giveDofManDofIDMask(int inode, IntArray &answer) const
Returns dofmanager dof mask for node.
Definition: libeam3dnl2.C:565
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. ...
void updateTempQuaternion(TimeStep *tStep)
Updates the temporary triad at the centre to the state identified by given solution step...
Definition: libeam3dnl2.C:118
virtual double computeLength()
Computes the length (zero for all but 1D geometries)
Definition: libeam3dnl2.C:485
StateCounterType giveSolutionStateCounter()
Returns current solution state counter.
Definition: timestep.h:188
#define OOFEM_ERROR(...)
Definition: error.h:61
REGISTER_Element(LSpace)
double computeSquaredNorm() const
Computes the square of the norm.
Definition: floatarray.C:846
void computeRotMtrx(FloatMatrix &answer, FloatArray &psi)
Evaluates the rotation matrix for large rotations according to Rodrigues formula for given pseudovect...
Definition: libeam3dnl2.C:88
#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
#define _IFT_LIBeam3dNL2_refnode
Definition: libeam3dnl2.h:43
contextIOResultType restoreYourself(DataStream &stream)
Definition: floatarray.C:895
FloatArray q
Quaternion at the center (last equilibrated).
Definition: libeam3dnl2.h:61
virtual contextIOResultType restoreContext(DataStream &stream, ContextMode mode, void *obj)
Restores the receiver state previously written in stream.
Definition: libeam3dnl2.C:862
void computeTempCurv(FloatArray &answer, TimeStep *tStep)
Compute the temporary curvature at the centre to the state identified by given solution step...
Definition: libeam3dnl2.C:735
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Receiver becomes the result of the product of aMatrix and anArray.
Definition: floatarray.C:676
virtual double giveUpdatedCoordinate(int ic, TimeStep *tStep, double scale=1.)
Returns updated ic-th coordinate of receiver.
Definition: node.C:245
virtual void computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
Computes constitutive matrix of receiver.
Definition: libeam3dnl2.C:434
void beTProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Receiver becomes the result of the product of aMatrix^T and anArray.
Definition: floatarray.C:708
virtual void updateYourself(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
Definition: libeam3dnl2.C:709
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
virtual double computeEdgeVolumeAround(GaussPoint *gp, int iEdge)
Computes volume related to integration point on local edge.
Definition: libeam3dnl2.C:607
virtual void give3dBeamStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)=0
Computes the stiffness matrix for 2d beams.
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: libeam3dnl2.C:685
virtual int setupIntegrationPoints(IntegrationRule &irule, int npoints, Element *element)
Sets up integration rule for the given element.
Definition: crosssection.C:54
Class representing vector of real numbers.
Definition: floatarray.h:82
Implementation of matrix containing floating point numbers.
Definition: floatmatrix.h:94
IRResultType
Type defining the return values of InputRecord reading operations.
Definition: irresulttype.h:47
virtual double give(CrossSectionProperty a, GaussPoint *gp)
Returns the value of cross section property at given point.
Definition: crosssection.C:151
GaussPoint * getIntegrationPoint(int n)
Access particular integration point of receiver.
virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj)
Stores receiver state to output stream.
Definition: libeam3dnl2.C:842
virtual void updateYourself(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
virtual int giveLocalCoordinateSystem(FloatMatrix &answer)
Returns local coordinate system of receiver Required by material models with ortho- and anisotrophy...
Definition: libeam3dnl2.C:619
virtual void drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
Definition: libeam3dnl2.C:882
virtual void giveEdgeDofMapping(IntArray &answer, int iEdge) const
Assembles edge dof mapping mask, which provides mapping between edge local DOFs and "global" element ...
Definition: libeam3dnl2.C:589
IntegrationPointStatus * giveMaterialStatus()
Returns reference to associated material status (NULL if not defined).
Definition: gausspoint.h:205
void computeRotMtrxFromQuaternion(FloatMatrix &answer, FloatArray &q)
Computes rotation matrix from given quaternion.
Definition: libeam3dnl2.C:158
LIBeam3dNL2(int n, Domain *d)
Definition: libeam3dnl2.C:57
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
Class representing the general Input Record.
Definition: inputrecord.h:101
double l0
Initial length.
Definition: libeam3dnl2.h:59
void computeXdVector(FloatArray &answer, TimeStep *tStep)
Computes x_21&#39; vector for given solution state.
Definition: libeam3dnl2.C:317
void add(const FloatMatrix &a)
Adds matrix to the receiver.
Definition: floatmatrix.C:1023
virtual contextIOResultType saveContext(DataStream &stream, ContextMode mode, void *obj=NULL)
Stores receiver state to output stream.
Definition: element.C:885
virtual void initForNewStep()
Initializes receivers state to new time step.
Definition: libeam3dnl2.C:725
#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
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.
void beProductTOf(const FloatMatrix &a, const FloatMatrix &b)
Assigns to the receiver product of .
Definition: floatmatrix.C:397
void zero()
Zeroes all coefficient of receiver.
Definition: floatmatrix.C:1326
Domain * giveDomain() const
Definition: femcmpnn.h:100
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: libeam3dnl2.C:448
void computeXMtrx(FloatMatrix &answer, TimeStep *tStep)
Computes X matrix at given solution state.
Definition: libeam3dnl2.C:257
Load is base abstract class for all loads.
Definition: load.h:61
void beProductOf(const FloatMatrix &a, const FloatMatrix &b)
Assigns to the receiver product of .
Definition: floatmatrix.C:337
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
Definition: libeam3dnl2.C:422
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
Node * giveNode(int i) const
Returns reference to the i-th node of element.
Definition: element.h:610
const FloatArray & giveStrainVector() const
Returns the const pointer to receiver&#39;s strain vector.
Definition: structuralms.h:105
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
void add(const FloatArray &src)
Adds array src to receiver.
Definition: floatarray.C:156
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