OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
simplecrosssection.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/CrossSections/simplecrosssection.h"
36 #include "../sm/Materials/structuralmaterial.h"
37 #include "../sm/Materials/structuralms.h"
38 #include "../sm/Elements/structuralelement.h"
39 #include "gausspoint.h"
40 #include "floatarray.h"
41 #include "classfactory.h"
42 #include "dynamicinputrecord.h"
43 #include "engngm.h"
44 
45 namespace oofem {
46 REGISTER_CrossSection(SimpleCrossSection);
47 
48 
49 void
51 {
52  StructuralMaterial *mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
53  mat->giveRealStressVector_3d(answer, gp, strain, tStep);
54 }
55 
56 void
58 {
59  StructuralMaterial *mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
60  IntArray strainControl = {
61  1, 2, 4, 5, 6
62  };
63  mat->giveRealStressVector_ShellStressControl(answer, gp, strain, strainControl, tStep);
64 }
65 
66 
67 void
69 {
70  StructuralMaterial *mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
71  mat->giveRealStressVector_PlaneStrain(answer, gp, strain, tStep);
72 }
73 
74 
75 void
77 {
78  StructuralMaterial *mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
79  mat->giveRealStressVector_PlaneStress(answer, gp, strain, tStep);
80 }
81 
82 
83 void
85 {
86  StructuralMaterial *mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
87  mat->giveRealStressVector_1d(answer, gp, strain, tStep);
88 }
89 
90 
91 void
93 {
94  StructuralMaterial *mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
95  mat->giveRealStressVector_Warping(answer, gp, strain, tStep);
96 }
97 
98 
99 void
101 {
102  StructuralMaterial *mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
103  mat->give3dMaterialStiffnessMatrix(answer, rMode, gp, tStep);
104 }
105 
106 
107 void
109 {
110  StructuralMaterial *mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
111  mat->givePlaneStressStiffMtrx(answer, rMode, gp, tStep);
112 }
113 
114 
115 void
117 {
118  StructuralMaterial *mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
119  mat->givePlaneStrainStiffMtrx(answer, rMode, gp, tStep);
120 }
121 
122 
123 void
125 {
126  StructuralMaterial *mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
127  mat->give1dStressStiffMtrx(answer, rMode, gp, tStep);
128 }
129 
130 
131 void
133 {
141  StructuralMaterial *mat = static_cast< StructuralMaterial * >( this->giveMaterial(gp) );
142  FloatArray elasticStrain, et, e0;
143  FloatMatrix tangent;
144  elasticStrain = strain;
145  this->giveTemperatureVector(et, gp, tStep);
146  if ( et.giveSize() > 0 ) {
147  mat->giveThermalDilatationVector(e0, gp, tStep);
148  double thick = this->give(CS_Thickness, gp);
149  elasticStrain.at(1) -= e0.at(1) * ( et.at(1) - mat->giveReferenceTemperature() );
150  if ( et.giveSize() > 1 ) {
151  elasticStrain.at(2) -= e0.at(1) * et.at(2) / thick; // kappa_x
152  }
153  }
154  this->give2dBeamStiffMtrx(tangent, ElasticStiffness, gp, tStep);
155  answer.beProductOf(tangent, elasticStrain);
156  StructuralMaterialStatus *status = static_cast< StructuralMaterialStatus * >( mat->giveStatus(gp) );
157  status->letTempStrainVectorBe(strain);
158  status->letTempStressVectorBe(answer);
159 }
160 
161 
162 void
164 {
170  StructuralMaterial *mat = static_cast< StructuralMaterial * >( this->giveMaterial(gp) );
171  FloatArray elasticStrain, et, e0;
172  FloatMatrix tangent;
173  elasticStrain = strain;
174  this->giveTemperatureVector(et, gp, tStep);
175  if ( et.giveSize() > 0 ) {
176  double thick = this->give(CS_Thickness, gp);
177  double width = this->give(CS_Width, gp);
178  mat->giveThermalDilatationVector(e0, gp, tStep);
179  elasticStrain.at(1) -= e0.at(1) * ( et.at(1) - mat->giveReferenceTemperature() );
180  if ( et.giveSize() > 1 ) {
181  elasticStrain.at(5) -= e0.at(1) * et.at(2) / thick; // kappa_y
182  if ( et.giveSize() > 2 ) {
183  elasticStrain.at(6) -= e0.at(1) * et.at(3) / width; // kappa_z
184  }
185  }
186  }
187  this->give3dBeamStiffMtrx(tangent, ElasticStiffness, gp, tStep);
188  answer.beProductOf(tangent, elasticStrain);
189  StructuralMaterialStatus *status = static_cast< StructuralMaterialStatus * >( mat->giveStatus(gp) );
190  status->letTempStrainVectorBe(strain);
191  status->letTempStressVectorBe(answer);
192 }
193 
194 
195 void
197 {
204  StructuralMaterial *mat = static_cast< StructuralMaterial * >( this->giveMaterial(gp) );
205  FloatArray elasticStrain, et, e0;
206  FloatMatrix tangent;
207  elasticStrain = strain;
208  this->giveTemperatureVector(et, gp, tStep);
209  if ( et.giveSize() > 0 ) {
210  double thick = this->give(CS_Thickness, gp);
211  mat->giveThermalDilatationVector(e0, gp, tStep);
212  if ( et.giveSize() > 1 ) {
213  elasticStrain.at(1) -= e0.at(1) * et.at(2) / thick; // kappa_x
214  elasticStrain.at(2) -= e0.at(2) * et.at(2) / thick; // kappa_y
215  }
216  }
217  this->give2dPlateStiffMtrx(tangent, ElasticStiffness, gp, tStep);
218  answer.beProductOf(tangent, elasticStrain);
219  StructuralMaterialStatus *status = static_cast< StructuralMaterialStatus * >( mat->giveStatus(gp) );
220  status->letTempStrainVectorBe(strain);
221  status->letTempStressVectorBe(answer);
222 }
223 
224 
225 void
227 {
234  StructuralMaterial *mat = static_cast< StructuralMaterial * >( this->giveMaterial(gp) );
235  FloatArray elasticStrain, et, e0;
236  FloatMatrix tangent;
237  elasticStrain = strain;
238  this->giveTemperatureVector(et, gp, tStep);
239  if ( et.giveSize() ) {
240  double thick = this->give(CS_Thickness, gp);
241  mat->giveThermalDilatationVector(e0, gp, tStep);
242  elasticStrain.at(1) -= e0.at(1) * ( et.at(1) - mat->giveReferenceTemperature() );
243  elasticStrain.at(2) -= e0.at(2) * ( et.at(1) - mat->giveReferenceTemperature() );
244  if ( et.giveSize() > 1 ) {
245  elasticStrain.at(4) -= e0.at(1) * et.at(2) / thick; // kappa_x
246  elasticStrain.at(5) -= e0.at(2) * et.at(2) / thick; // kappa_y
247  }
248  }
249  this->give3dShellStiffMtrx(tangent, ElasticStiffness, gp, tStep);
250  answer.beProductOf(tangent, elasticStrain);
251  StructuralMaterialStatus *status = static_cast< StructuralMaterialStatus * >( mat->giveStatus(gp) );
252  status->letTempStrainVectorBe(strain);
253  status->letTempStressVectorBe(answer);
254 }
255 
256 
257 void
259 {
260  FloatMatrix tangent;
261  this->giveMembraneRotStiffMtrx(tangent, ElasticStiffness, gp, tStep);
262  answer.beProductOf(tangent, strain);
263 
264  StructuralMaterialStatus *status = static_cast< StructuralMaterialStatus * >( this->giveMaterial(gp)->giveStatus(gp) );
265  status->letTempStrainVectorBe(strain);
266  status->letTempStressVectorBe(answer);
267 
268 
271 }
272 
273 void
275 {
276  StructuralMaterial *mat = static_cast< StructuralMaterial * >( this->giveMaterial(gp) );
277  return mat->giveRealStressVector_2dPlateSubSoil(answer, gp, generalizedStrain, tStep);
278 }
279 
280 
281 void
283 {
284  MaterialMode mode = gp->giveMaterialMode();
285  if ( mode == _2dBeam ) {
286  this->give2dBeamStiffMtrx(answer, rMode, gp, tStep);
287  } else if ( mode == _3dBeam ) {
288  this->give3dBeamStiffMtrx(answer, rMode, gp, tStep);
289  } else if ( mode == _2dPlate ) {
290  this->give2dPlateStiffMtrx(answer, rMode, gp, tStep);
291  } else if ( mode == _3dShell ) {
292  this->give3dShellStiffMtrx(answer, rMode, gp, tStep);
293  } else if ( mode == _3dDegeneratedShell ) {
294  this->give3dDegeneratedShellStiffMtrx(answer, rMode, gp, tStep);
295  } else {
296  StructuralMaterial *mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
297 
298  if ( mode == _3dMat ) {
299  mat->give3dMaterialStiffnessMatrix(answer, rMode, gp, tStep);
300  } else if ( mode == _PlaneStress ) {
301  mat->givePlaneStressStiffMtrx(answer, rMode, gp, tStep);
302  } else if ( mode == _PlaneStrain ) {
303  mat->givePlaneStrainStiffMtrx(answer, rMode, gp, tStep);
304  } else if ( mode == _1dMat ) {
305  mat->give1dStressStiffMtrx(answer, rMode, gp, tStep);
306  } else {
307  mat->giveStiffnessMatrix(answer, rMode, gp, tStep);
308  }
309  }
310 }
311 
312 
313 void
315 {
316  StructuralMaterial *mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
317 
318  FloatMatrix mat3d;
319  double area, Iy, shearAreaz;
320 
321  mat->give1dStressStiffMtrx(mat3d, rMode, gp, tStep);
322  area = this->give(CS_Area, gp);
323  Iy = this->give(CS_InertiaMomentY, gp);
324  shearAreaz = this->give(CS_ShearAreaZ, gp);
325 
326  answer.resize(3, 3);
327  answer.zero();
328 
329  answer.at(1, 1) = mat3d.at(1, 1) * area;
330  answer.at(2, 2) = mat3d.at(1, 1) * Iy;
331  answer.at(3, 3) = shearAreaz * mat->give('G', gp);
332 }
333 
334 
335 void
337 {
338  StructuralMaterial *mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
339 
340  FloatMatrix mat3d;
341  double area, E, G, Iy, Iz, Ik;
342  double shearAreay, shearAreaz;
343 
344  mat->give1dStressStiffMtrx(mat3d, rMode, gp, tStep);
345  E = mat3d.at(1, 1);
346  G = mat->give('G', gp);
347  area = this->give(CS_Area, gp);
348  Iy = this->give(CS_InertiaMomentY, gp);
349  Iz = this->give(CS_InertiaMomentZ, gp);
350  Ik = this->give(CS_TorsionMomentX, gp);
351 
352  //shearCoeff = this->give(CS_BeamShearCoeff);
353  shearAreay = this->give(CS_ShearAreaY, gp);
354  shearAreaz = this->give(CS_ShearAreaZ, gp);
355 
356  answer.resize(6, 6);
357  answer.zero();
358 
359  answer.at(1, 1) = E * area;
361  answer.at(2, 2) = shearAreay * G;
362  answer.at(3, 3) = shearAreaz * G;
363  //answer.at(2, 2) = shearCoeff * G * area;
364  //answer.at(3, 3) = shearCoeff * G * area;
365  answer.at(4, 4) = G * Ik;
366  answer.at(5, 5) = E * Iy;
367  answer.at(6, 6) = E * Iz;
368 }
369 
370 
371 void
373 {
374  StructuralMaterial *mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
375 
376  FloatMatrix mat3d;
377  double thickness3, thickness;
378 
379  mat->givePlaneStressStiffMtrx(mat3d, rMode, gp, tStep);
380  thickness = this->give(CS_Thickness, gp);
381  thickness3 = thickness * thickness * thickness;
382 
383  answer.resize(5, 5);
384  answer.zero();
385 
386  for ( int i = 1; i <= 2; i++ ) {
387  for ( int j = 1; j <= 2; j++ ) {
388  answer.at(i, j) = mat3d.at(i, j) * thickness3 / 12.;
389  }
390  }
391 
392  answer.at(3, 3) = mat3d.at(3, 3) * thickness3 / 12.;
393  answer.at(4, 4) = mat3d.at(3, 3) * thickness * ( 5. / 6. );
394  answer.at(5, 5) = answer.at(4, 4);
395 }
396 
397 
398 void
400 {
401  StructuralMaterial *mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
402 
403  FloatMatrix mat3d;
404 
405  double thickness = this->give(CS_Thickness, gp);
406  double thickness3 = thickness * thickness * thickness;
407 
408  mat->givePlaneStressStiffMtrx(mat3d, rMode, gp, tStep);
409 
410  answer.resize(8, 8);
411  answer.zero();
412 
413  for ( int i = 1; i <= 3; i++ ) {
414  for ( int j = 1; j <= 3; j++ ) {
415  answer.at(i, j) = mat3d.at(i, j) * thickness;
416  }
417  }
418  for ( int i = 1; i <= 3; i++ ) {
419  for ( int j = 1; j <= 3; j++ ) {
420  answer.at(i + 3, j + 3) = mat3d.at(i, j) * thickness3 / 12.0;
421  }
422  }
423 
424  answer.at(8, 8) = answer.at(7, 7) = mat3d.at(3, 3) * thickness * ( 5. / 6. );
425 }
426 
427 
428 void
430 {
431  StructuralMaterial *mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
432  mat->give3dMaterialStiffnessMatrix(answer, rMode, gp, tStep);
433 
434  answer.at(1, 1) -= answer.at(1, 3) * answer.at(3, 1) / answer.at(3, 3);
435  answer.at(2, 1) -= answer.at(2, 3) * answer.at(3, 1) / answer.at(3, 3);
436  answer.at(1, 2) -= answer.at(1, 3) * answer.at(3, 2) / answer.at(3, 3);
437  answer.at(2, 2) -= answer.at(2, 3) * answer.at(3, 2) / answer.at(3, 3);
438 
439  answer.at(3, 1) = 0.0;
440  answer.at(3, 2) = 0.0;
441  answer.at(3, 3) = 0.0;
442  answer.at(2, 3) = 0.0;
443  answer.at(1, 3) = 0.0;
444 }
445 
446 void
448 {
449  StructuralMaterial *mat;
450  mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
451  mat->givePlaneStressStiffMtrx(answer, ElasticStiffness, gp, tStep);
452  answer.resizeWithData(4, 4);
453  answer.at(4, 4) = this->give(CS_DrillingStiffness, gp);
454 }
455 
456 void
458 {
459  StructuralMaterial *mat;
460  mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
461  mat->give2dPlateSubSoilStiffMtrx(answer, ElasticStiffness, gp, tStep);
462 }
463 
464 
467 //
468 // instanciates receiver from input record
469 //
470 {
471  IRResultType result; // Required by IR_GIVE_FIELD macro
472  double value;
473 
474  double thick = 0.0;
478  }
479 
480  double width = 0.0;
484  }
485 
486  double area = 0.0;
489  } else {
490  area = thick * width;
491  }
493 
494  value = 0.0;
497 
498  value = 0.0;
501 
502  value = 0.0;
505 
506  double beamshearcoeff = 0.0;
508  propertyDictionary.add(CS_BeamShearCoeff, beamshearcoeff);
509 
510  value = 0.0;
512  if ( value == 0.0 ) {
513  value = beamshearcoeff * area;
514  }
516 
517  value = 0.0;
519  if ( value == 0.0 ) {
520  value = beamshearcoeff * area;
521  }
523 
524  value = 0.0;
527 
528  value = 0.0;
531 
532  value = 0.0;
535 
536  this->materialNumber = 0;
538 
540  value = 0.0;
543 
544  value = 0.0;
547 
548  value = 1.0;
551  }
552 
554 }
555 
556 
557 void
559 {
561 
564  }
565 
566  if ( this->propertyDictionary.includes(CS_Width) ) {
568  }
569 
570  if ( this->propertyDictionary.includes(CS_Area) ) {
572  }
573 
576  }
577 
580  }
581 
584  }
585 
588  }
589 
592  }
593 
596  }
597 
599 
602  }
603 
606  }
607 
610  }
611 }
612 
613 void
615 {
617  MaterialStatus *matStat = mat->CreateStatus(& iGP);
618  iGP.setMaterialStatus(matStat);
619 }
620 
621 
622 bool
624 {
625  if ( this->giveMaterialNumber() ) {
626  return this->domain->giveMaterial( this->giveMaterialNumber() )->isCharacteristicMtrxSymmetric(rMode);
627  } else {
628  return false; // Bet false...
629  }
630 }
631 
632 
633 Material
635 {
636  if ( this->giveMaterialNumber() ) {
637  return this->giveDomain()->giveMaterial( this->giveMaterialNumber() );
638  } else {
639  return ip->giveElement()->giveMaterial();
640  }
641 }
642 
643 
644 double
646 {
647  return this->giveMaterial(gp)->give(aProperty, gp);
648 }
649 
650 
651 int
653 {
654  if ( type == IST_CrossSectionNumber ) {
655  answer.resize(1);
656  answer.at(1) = this->giveNumber();
657  return 1;
658  }
659  return this->giveMaterial(ip)->giveIPValue(answer, ip, type, tStep);
660 }
661 
662 
663 
664 int
666 //
667 // check internal consistency
668 // mainly tests, whether material and crossSection data
669 // are safe for conversion to "Structural" versions
670 //
671 {
672  int result = 1;
673  Material *mat = this->giveDomain()->giveMaterial(this->materialNumber);
674  if ( !dynamic_cast< StructuralMaterial * >(mat) ) {
675  OOFEM_WARNING( "material %s without structural support", mat->giveClassName() );
676  result = 0;
677  }
678 
679  return result;
680 }
681 
682 
683 Interface
685 {
686  return this->giveMaterial(ip)->giveInterface(t);
687 }
688 
689 
690 
691 
692 void
694 {
695  // This function returns the first Piola-Kirchoff stress in vector format
696  // corresponding to a given deformation gradient according to the stress-deformation
697  // mode stored in the each gp.
698 
699  MaterialMode mode = gp->giveMaterialMode();
700  StructuralMaterial *mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
701 
702  if ( mode == _3dMat ) {
703  mat->giveFirstPKStressVector_3d(answer, gp, reducedvF, tStep);
704  } else if ( mode == _PlaneStrain ) {
705  mat->giveFirstPKStressVector_PlaneStrain(answer, gp, reducedvF, tStep);
706  } else if ( mode == _PlaneStress ) {
707  mat->giveFirstPKStressVector_PlaneStress(answer, gp, reducedvF, tStep);
708  } else if ( mode == _1dMat ) {
709  mat->giveFirstPKStressVector_1d(answer, gp, reducedvF, tStep);
710  } else {
711  OOFEM_ERROR( "unknown mode (%s)", __MaterialModeToString(mode) );
712  }
713 }
714 
715 
716 void
718 {
719  // This function returns the Cauchy stress in vector format
720  // corresponding to a given deformation gradient according to the stress-deformation
721  // mode stored in the each gp.
722 
723  MaterialMode mode = gp->giveMaterialMode();
724  StructuralMaterial *mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
725 
726  if ( mode == _3dMat ) {
727  mat->giveCauchyStressVector_3d(answer, gp, reducedvF, tStep);
728  } else if ( mode == _PlaneStrain ) {
729  mat->giveCauchyStressVector_PlaneStrain(answer, gp, reducedvF, tStep);
730  } else if ( mode == _PlaneStress ) {
731  mat->giveCauchyStressVector_PlaneStress(answer, gp, reducedvF, tStep);
732  } else if ( mode == _1dMat ) {
733  mat->giveCauchyStressVector_1d(answer, gp, reducedvF, tStep);
734  }
735 }
736 
737 void
739 {
740  MaterialMode mode = gp->giveMaterialMode();
741  StructuralMaterial *mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
742 
743  if ( mode == _3dMat ) {
744  // mat->giveCauchyStressVector_3d(answer, gp, reducedvF, tStep);
745  } else if ( mode == _PlaneStrain ) {
746  mat->giveEshelbyStressVector_PlaneStrain(answer, gp, reducedvF, tStep);
747  } else if ( mode == _PlaneStress ) {
748  // mat->giveCauchyStressVector_PlaneStress(answer, gp, reducedvF, tStep);
749  } else if ( mode == _1dMat ) {
750  // mat->giveCauchyStressVector_1d(answer, gp, reducedvF, tStep);
751  }
752 }
753 
754 
755 void
757  MatResponseMode rMode, GaussPoint *gp,
758  TimeStep *tStep)
759 {
760  StructuralMaterial *mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
761 
762  MaterialMode mode = gp->giveMaterialMode();
763  if ( mode == _3dMat ) {
764  mat->give3dMaterialStiffnessMatrix_dPdF(answer, rMode, gp, tStep);
765  } else if ( mode == _PlaneStress ) {
766  mat->givePlaneStressStiffMtrx_dPdF(answer, rMode, gp, tStep);
767  } else if ( mode == _PlaneStrain ) {
768  mat->givePlaneStrainStiffMtrx_dPdF(answer, rMode, gp, tStep);
769  } else if ( mode == _1dMat ) {
770  mat->give1dStressStiffMtrx_dPdF(answer, rMode, gp, tStep);
771  } else {
772  OOFEM_ERROR( "unknown mode (%s)", __MaterialModeToString(mode) );
773  }
774 }
775 
776 
777 void
779  MatResponseMode rMode, GaussPoint *gp,
780  TimeStep *tStep)
781 {
782  StructuralMaterial *mat = dynamic_cast< StructuralMaterial * >( this->giveMaterial(gp) );
783 
784  MaterialMode mode = gp->giveMaterialMode();
785  if ( mode == _3dMat ) {
786  mat->give3dMaterialStiffnessMatrix_dCde(answer, rMode, gp, tStep);
787  } else if ( mode == _PlaneStress ) {
788  mat->givePlaneStressStiffMtrx_dCde(answer, rMode, gp, tStep);
789  } else if ( mode == _PlaneStrain ) {
790  mat->givePlaneStrainStiffMtrx_dCde(answer, rMode, gp, tStep);
791  } else if ( mode == _1dMat ) {
792  mat->give1dStressStiffMtrx_dCde(answer, rMode, gp, tStep);
793  } else {
794  OOFEM_ERROR( "unknown mode (%s)", __MaterialModeToString(mode) );
795  }
796 }
797 
798 
799 void
801 {
802  Element *elem = gp->giveElement();
803  answer.clear();
804  //sum up all prescribed temperatures over an element
805  StructuralElement *selem = dynamic_cast< StructuralElement * >(elem);
806  selem->computeResultingIPTemperatureAt(answer, tStep, gp, VM_Total);
807 
808  /* add external source, if provided */
810  FieldPtr tf;
811 
812  if ( ( tf = fm->giveField(FT_Temperature) ) ) {
813  // temperature field registered
814  FloatArray gcoords, et2;
815  int err;
816  elem->computeGlobalCoordinates( gcoords, gp->giveNaturalCoordinates() );
817  if ( ( err = tf->evaluateAt(et2, gcoords, VM_Total, tStep) ) ) {
818  OOFEM_ERROR("tf->evaluateAt failed, element %d, error code %d", elem->giveNumber(), err);
819  }
820  if ( et2.isNotEmpty() ) {
821  if ( answer.isEmpty() ) {
822  answer = et2;
823  } else {
824  answer.at(1) += et2.at(1);
825  }
826  }
827  }
828 }
829 
830 int
832 {
833  return this->giveMaterial(gp)->packUnknowns(buff, tStep, gp);
834 }
835 
836 int
838 {
839  return this->giveMaterial(gp)->unpackAndUpdateUnknowns(buff, tStep, gp);
840 }
841 
842 int
844 {
845  return this->giveMaterial(gp)->estimatePackSize(buff, gp);
846 }
847 } // end namespace oofem
virtual void giveFirstPKStressVector_1d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
Default implementation relies on giveFirstPKStressVector_3d.
virtual int unpackAndUpdateUnknowns(DataStream &buff, TimeStep *tStep, GaussPoint *gp)
Unpack and updates all necessary data of given integration point (according to element parallel_mode)...
Relative penalty stiffness for drilling DOFs.
Definition: crosssection.h:69
Shear area in z direction.
Definition: crosssection.h:67
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
virtual void giveTemperatureVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep)
void setField(int item, InputFieldType id)
virtual void giveEshelbyStressVector_PlaneStrain(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
Prototype for computation of Eshelby stress.
MaterialMode giveMaterialMode()
Returns corresponding material mode of receiver.
Definition: gausspoint.h:191
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: crosssection.C:67
virtual void computeResultingIPTemperatureAt(FloatArray &answer, TimeStep *tStep, GaussPoint *gp, ValueModeType mode)
Computes at given time (tStep) the the resulting temperature component array.
void letTempStrainVectorBe(const FloatArray &v)
Assigns tempStrainVector to given vector v.
Definition: structuralms.h:137
virtual void giveStiffnessMatrix_3d(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Method for computing the stiffness matrix.
#define _IFT_SimpleCrossSection_ik
Torsion moment x.
virtual void give1dStressStiffMtrx(FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
Method for computing 1d stiffness matrix of receiver.
virtual MaterialStatus * giveStatus(GaussPoint *gp) const
Returns material status of receiver in given integration point.
Definition: material.C:244
std::shared_ptr< Field > FieldPtr
Definition: field.h:72
virtual void giveFirstPKStressVector_3d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress.
virtual void giveRealStress_PlaneStress(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep)
virtual void giveRealStressVector_PlaneStrain(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
Default implementation relies on giveRealStressVector_3d.
virtual void giveGeneralizedStress_Beam2d(FloatArray &answer, GaussPoint *gp, const FloatArray &generalizedStrain, TimeStep *tStep)
Computes the generalized stress vector for given strain and integration point.
virtual int packUnknowns(DataStream &buff, TimeStep *tStep, GaussPoint *gp)
Pack all necessary data of integration point (according to element parallel_mode) into given communic...
IntegrationPointStatus * setMaterialStatus(IntegrationPointStatus *ptr, int n)
Sets Material status managed by receiver.
Definition: gausspoint.h:213
virtual void give1dStressStiffMtrx_dPdF(FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
Domain * domain
Link to domain object, useful for communicating with other FEM components.
Definition: femcmpnn.h:82
Director vector component in y-axis.
Definition: crosssection.h:75
virtual void giveFirstPKStressVector_PlaneStrain(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
Default implementation relies on giveFirstPKStressVector_3d.
#define _IFT_SimpleCrossSection_iz
Inertia moment z.
virtual Interface * giveMaterialInterface(InterfaceType t, IntegrationPoint *ip)
virtual void giveCauchyStressVector_1d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
virtual int packUnknowns(DataStream &buff, TimeStep *tStep, GaussPoint *ip)
Pack all necessary data of integration point (according to element parallel_mode) into given communic...
Definition: material.h:285
double & at(int aKey)
Returns the value of the pair which key is aKey.
Definition: dictionary.C:106
virtual void giveCauchyStressVector_PlaneStrain(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
The purpose of DataStream abstract class is to allow to store/restore context to different streams...
Definition: datastream.h:54
#define _IFT_SimpleCrossSection_shearareaz
Shear area z direction.
virtual void giveGeneralizedStress_MembraneRot(FloatArray &answer, GaussPoint *gp, const FloatArray &generalizedStrain, TimeStep *tStep)
double & at(int i)
Coefficient access function.
Definition: floatarray.h:131
virtual void givePlaneStrainStiffMtrx(FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
Method for computing plane strain stiffness matrix of receiver.
#define _IFT_SimpleCrossSection_area
#define _IFT_SimpleCrossSection_drillType
Type of artificially added stiffnes for drilling DOFs.
This class implements a structural material status information.
Definition: structuralms.h:65
virtual void giveEshelbyStresses(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedvF, TimeStep *tStep)
Computes the Eshelby stress vector.
void clear()
Clears receiver (zero size).
Definition: floatarray.h:206
bool includes(int aKey)
Checks if dictionary includes given key.
Definition: dictionary.C:126
virtual void giveRealStress_3dDegeneratedShell(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep)
#define _IFT_SimpleCrossSection_iy
Inertia moment y.
double giveReferenceTemperature()
Returns the reference temperature of receiver.
virtual void give3dMaterialStiffnessMatrix_dCde(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
#define _IFT_SimpleCrossSection_MaterialNumber
Material number for the bulk material.
EngngModel * giveEngngModel()
Returns engineering model to which receiver is associated.
Definition: domain.C:433
virtual bool hasField(InputFieldType id)=0
Returns true if record contains field identified by idString keyword.
virtual double give(int aProperty, GaussPoint *gp)
Returns the value of cross section property.
REGISTER_CrossSection(EmptyCS)
Penalty stiffness for drilling DOFs.
Definition: crosssection.h:68
Abstract base class for all finite elements.
Definition: element.h:145
Moment of inertia around z-axis.
Definition: crosssection.h:64
Element * giveElement()
Returns corresponding element to receiver.
Definition: gausspoint.h:188
#define _IFT_SimpleCrossSection_width
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver acording to object description stored in input record.
virtual void give2dBeamStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Computes the stiffness matrix for 2d beams.
#define _IFT_SimpleCrossSection_directorz
virtual void giveInputRecord(DynamicInputRecord &input)
Setups the input record string of receiver.
Definition: crosssection.C:82
MaterialMode
Type representing material mode of integration point.
Definition: materialmode.h:89
FieldManager * giveFieldManager()
Definition: engngm.h:131
Class implementing an array of integers.
Definition: intarray.h:61
MatResponseMode
Describes the character of characteristic material matrix.
virtual void giveInputRecord(DynamicInputRecord &input)
Setups the input record string of receiver.
#define _IFT_SimpleCrossSection_directorx
#define _IFT_SimpleCrossSection_shearcoeff
virtual void giveRealStress_PlaneStrain(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep)
virtual void givePlaneStressStiffMtrx_dPdF(FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
virtual void givePlaneStrainStiffMtrx_dCde(FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
const char * __MaterialModeToString(MaterialMode _value)
Definition: cltypes.C:314
virtual void giveRealStressVector_2dPlateSubSoil(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
Default implementation is not provided.
virtual double give(int aProperty, GaussPoint *gp)
Returns the value of material property &#39;aProperty&#39;.
Definition: material.C:52
virtual void giveRealStress_Warping(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep)
virtual void giveCauchyStressVector_PlaneStress(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
Abstract base class for all "structural" finite elements.
virtual void giveStiffnessMatrix_PlaneStrain(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
virtual void giveRealStress_1d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep)
virtual void giveStiffnessMatrix(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Computes the stiffness matrix for giveRealStressVector of receiver in given integration point...
Shear coefficient of beam.
Definition: crosssection.h:61
#define _IFT_SimpleCrossSection_directory
Dictionary propertyDictionary
Dictionary for storing cross section parameters (like dimensions).
Definition: crosssection.h:115
Director vector component in z-axis.
Definition: crosssection.h:76
#define E(p)
Definition: mdm.C:368
Material * giveMaterial(int n)
Service for accessing particular domain material model.
Definition: domain.C:281
#define OOFEM_ERROR(...)
Definition: error.h:61
virtual void givePlaneStressStiffMtrx_dCde(FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
EngngModelContext * giveContext()
Context requesting service.
Definition: engngm.h:1078
virtual void giveGeneralizedStress_Plate(FloatArray &answer, GaussPoint *gp, const FloatArray &generalizedStrain, TimeStep *tStep)
virtual void giveFirstPKStresses(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedFIncrement, TimeStep *tStep)
Computes the First Piola-Kirchoff stress vector for a given deformation gradient and integration poin...
virtual bool isCharacteristicMtrxSymmetric(MatResponseMode rMode)
Returns true if stiffness matrix of receiver is symmetric Default implementation returns true...
Definition: material.h:129
virtual void give3dMaterialStiffnessMatrix(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Computes full 3d material stiffness matrix at given integration point, time, respecting load history ...
void resizeWithData(int, int)
Checks size of receiver towards requested bounds.
Definition: floatmatrix.C:1369
Moment of inertia around x-axis.
Definition: crosssection.h:65
virtual void giveStiffnessMatrix_dPdF(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Computes the material stiffness matrix dPdF of receiver in a given integration point, respecting its history.
virtual void give2dPlateStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Method for computing 2d plate stiffness matrix.
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 isCharacteristicMtrxSymmetric(MatResponseMode mode)
Check for symmetry of stiffness matrix.
virtual int giveIPValue(FloatArray &answer, GaussPoint *ip, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in Reduced form.
bool isEmpty() const
Returns true if receiver is empty.
Definition: floatarray.h:222
Abstract base class for all material models.
Definition: material.h:95
double at(int i, int j) const
Coefficient access function.
Definition: floatmatrix.h:176
virtual void giveRealStressVector_ShellStressControl(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, const IntArray &strainControl, TimeStep *tStep)
virtual void givePlaneStrainStiffMtrx_dPdF(FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
FieldPtr giveField(FieldType key)
Returns the previously registered field under given key; NULL otherwise.
Definition: fieldmanager.C:63
virtual void giveRealStressVector_3d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress.
Abstract base class representing a material status information.
Definition: matstatus.h:84
Pair * add(int aKey, double value)
Adds a new Pair with given keyword and value into receiver.
Definition: dictionary.C:81
virtual void giveFirstPKStressVector_PlaneStress(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
Default implementation relies on giveFirstPKStressVector_3d.
virtual void giveCauchyStresses(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedFIncrement, TimeStep *tStep)
Computes the Cauchy stress vector for a given increment of deformation gradient and given integration...
Class representing vector of real numbers.
Definition: floatarray.h:82
virtual void give2dPlateSubSoilStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Method for computing subsoil stiffness matrix for plates.
virtual void giveRealStressVector_Warping(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
Default implementation relies on giveRealStressVector_StressControl.
Implementation of matrix containing floating point numbers.
Definition: floatmatrix.h:94
Shear area in y direction.
Definition: crosssection.h:66
IRResultType
Type defining the return values of InputRecord reading operations.
Definition: irresulttype.h:47
void letTempStressVectorBe(const FloatArray &v)
Assigns tempStressVector to given vector v.
Definition: structuralms.h:135
virtual void give3dDegeneratedShellStiffMtrx(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
Method for computing 3d shell stiffness matrix on degenerated shell elements.
virtual const char * giveClassName() const =0
#define _IFT_SimpleCrossSection_relDrillStiffness
Relative penalty term for drilling stiffness.
virtual void give3dShellStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Method for computing 3d shell stiffness matrix.
void resize(int rows, int cols)
Checks size of receiver towards requested bounds.
Definition: floatmatrix.C:1358
virtual void giveThermalDilatationVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep)
Returns a vector of coefficients of thermal dilatation in direction of each material principal (local...
Class representing the general Input Record.
Definition: inputrecord.h:101
virtual void giveGeneralizedStress_Shell(FloatArray &answer, GaussPoint *gp, const FloatArray &generalizedStrain, TimeStep *tStep)
virtual Material * giveMaterial(IntegrationPoint *ip)
Returns the material associated with the GP.
Class Interface.
Definition: interface.h:82
virtual void giveGeneralizedStress_Beam3d(FloatArray &answer, GaussPoint *gp, const FloatArray &generalizedStrain, TimeStep *tStep)
virtual Interface * giveInterface(InterfaceType t)
Interface requesting service.
Definition: femcmpnn.h:179
Class representing the a dynamic Input Record.
virtual void give2dPlateSubSoilStiffMtrx(FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
Method for computing stiffness matrix of plate subsoil model.
virtual void giveStiffnessMatrix_1d(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
virtual void giveRealStressVector_1d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
Default implementation relies on giveRealStressVector_StressControl.
Abstract base class for all "structural" constitutive models.
virtual void giveStiffnessMatrix_PlaneStress(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
virtual int estimatePackSize(DataStream &buff, GaussPoint *ip)
Estimates the necessary pack size to hold all packed data of receiver.
Definition: material.h:298
void zero()
Zeroes all coefficient of receiver.
Definition: floatmatrix.C:1326
virtual void giveRealStressVector_PlaneStress(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
Default implementation relies on giveRealStressVector_StressControl.
Domain * giveDomain() const
Definition: femcmpnn.h:100
Type of artificially added drilling stiffness for drilling DOFs.
Definition: crosssection.h:70
InterfaceType
Enumerative type, used to identify interface type.
Definition: interfacetype.h:43
virtual void giveRealStress_3d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep)
virtual void giveGeneralizedStress_PlateSubSoil(FloatArray &answer, GaussPoint *gp, const FloatArray &generalizedStrain, TimeStep *tStep)
#define _IFT_SimpleCrossSection_shearareay
Shear area y direction.
#define _IFT_SimpleCrossSection_thick
virtual void giveCharMaterialStiffnessMatrix(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Computes the stiffness matrix of receiver in given integration point, respecting its history...
#define IR_GIVE_OPTIONAL_FIELD(__ir, __value, __id)
Macro facilitating the use of input record reading methods.
Definition: inputrecord.h:78
virtual void give3dMaterialStiffnessMatrix_dPdF(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
int giveSize() const
Returns the size of receiver.
Definition: floatarray.h:218
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
virtual int unpackAndUpdateUnknowns(DataStream &buff, TimeStep *tStep, GaussPoint *ip)
Unpack and updates all necessary data of given integration point (according to element parallel_mode)...
Definition: material.h:294
virtual void give3dBeamStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Computes the stiffness matrix for 2d beams.
#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
virtual void givePlaneStressStiffMtrx(FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
Method for computing plane stress stiffness matrix of receiver.
Director vector component in x-axis.
Definition: crosssection.h:74
Moment of inertia around y-axis.
Definition: crosssection.h:63
Class representing integration point in finite element program.
Definition: gausspoint.h:93
#define OOFEM_WARNING(...)
Definition: error.h:62
virtual void giveMembraneRotStiffMtrx(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Method for computing membrane stiffness matrix with added drilling stiffness.
Class representing solution step.
Definition: timestep.h:80
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in Reduced form.
Definition: material.C:142
virtual int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords)
Computes the global coordinates from given element&#39;s local coordinates.
Definition: element.C:1207
virtual int checkConsistency()
Allows programmer to test some internal data, before computation begins.
virtual void createMaterialStatus(GaussPoint &iGP)
const FloatArray & giveNaturalCoordinates()
Returns coordinate array of receiver.
Definition: gausspoint.h:138
virtual int estimatePackSize(DataStream &buff, GaussPoint *gp)
Estimates the necessary pack size to hold all packed data of receiver.
virtual Material * giveMaterial()
Definition: element.C:484
virtual void give1dStressStiffMtrx_dCde(FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
#define _IFT_SimpleCrossSection_drillStiffness
Penalty term for drilling stiffness.
void resize(int s)
Resizes receiver towards requested size.
Definition: floatarray.C:631
virtual MaterialStatus * CreateStatus(GaussPoint *gp) const
Creates new copy of associated status and inserts it into given integration point.
Definition: material.h:316
virtual void giveCauchyStressVector_3d(FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
virtual void giveStiffnessMatrix_dCde(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Computes the material stiffness matrix dCde of receiver in a given integration point, respecting its history.

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