MAST
solid_1d_section_element_property_card.cpp
Go to the documentation of this file.
1 /*
2  * MAST: Multidisciplinary-design Adaptation and Sensitivity Toolkit
3  * Copyright (C) 2013-2019 Manav Bhatia
4  *
5  * This library is free software; you can redistribute it and/or
6  * modify it under the terms of the GNU Lesser General Public
7  * License as published by the Free Software Foundation; either
8  * version 2.1 of the License, or (at your option) any later version.
9  *
10  * This library is distributed in the hope that it will be useful,
11  * but WITHOUT ANY WARRANTY; without even the implied warranty of
12  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
13  * Lesser General Public License for more details.
14  *
15  * You should have received a copy of the GNU Lesser General Public
16  * License along with this library; if not, write to the Free Software
17  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
18  */
19 
20 // MAST includes
24 #include "base/elem_base.h"
25 
26 
27 
28 namespace MAST {
29  namespace Solid1DSectionProperty {
30 
31  class Area: public MAST::FieldFunction<Real> {
32  public:
34  const MAST::FieldFunction<Real>& hz):
35  MAST::FieldFunction<Real>("Area"),
36  _hy(hy),
37  _hz(hz) {
38  _functions.insert(&hy);
39  _functions.insert(&hz);
40  }
41 
42  virtual ~Area() { }
43 
44  virtual void operator() (const libMesh::Point& p,
45  const Real t,
46  Real& m) const {
47  Real hy, hz;
48  _hy(p, t, hy);
49  _hz(p, t, hz);
50 
51  m = hy*hz;
52  }
53 
54  virtual void derivative ( const MAST::FunctionBase& f,
55  const libMesh::Point& p,
56  const Real t,
57  Real& m) const {
58  Real hy, hz, dhy, dhz;
59  _hy(p, t, hy); _hy.derivative( f, p, t, dhy);
60  _hz(p, t, hz); _hz.derivative( f, p, t, dhz);
61 
62  m = dhy*hz + hy*dhz;
63  }
64 
65  protected:
66 
68  };
69 
70 
71 
72  class TorsionalConstant: public MAST::FieldFunction<Real> {
73  public:
75  const MAST::FieldFunction<Real>& hz):
76  MAST::FieldFunction<Real>("TorsionalConstant"),
77  _hy(hy),
78  _hz(hz) {
79  _functions.insert(&hy);
80  _functions.insert(&hz);
81  }
82 
83 
84  virtual ~TorsionalConstant() { }
85 
86  virtual void operator() (const libMesh::Point& p,
87  const Real t,
88  Real& m) const {
89  Real hy, hz, a, b;
90  _hy(p, t, hy);
91  _hz(p, t, hz);
92 
93 
94  // shorter side is b, and longer side is a
95  if (hy > hz) {
96  a = hy;
97  b = hz;
98  }
99  else {
100  a = hz;
101  b = hy;
102  }
103 
104  m = a*pow(b,3)*(1./3.-.21*b/a*(1.-pow(b/a,4)/12.));
105  }
106 
107 
108  virtual void derivative ( const MAST::FunctionBase& f,
109  const libMesh::Point& p,
110  const Real t,
111  Real& m) const {
112  Real hy, hz, dhy, dhz, a, b, da, db;
113  _hy(p, t, hy); _hy.derivative( f, p, t, dhy);
114  _hz(p, t, hz); _hz.derivative( f, p, t, dhz);
115 
116  // shorter side is b, and longer side is a
117  if (hy > hz) {
118  a = hy; da = dhy;
119  b = hz; db = dhz;
120  }
121  else {
122  a = hz; da = dhz;
123  b = hy; db = dhy;
124  }
125 
126  m =
127  da*pow(b,3)*(1./3.-.21*b/a*(1.-pow(b/a,4)/12.)) +
128  a*3.*pow(b,2)*db*(1./3.-.21*b/a*(1.-pow(b/a,4)/12.)) +
129  a*pow(b,3)*(-.21*db/a*(1.-pow(b/a,4)/12.) +
130  (.21*b/pow(a,2)*da*(1.-pow(b/a,4)/12.)) +
131  (-.21*b/a*(-4.*pow(b,3)*db/pow(a,4)/12.+
132  4.*pow(b,4)/pow(a,5)*da/12.)));
133  }
134 
135  protected:
136 
138  };
139 
140 
141 
142  class PolarInertia: public MAST::FieldFunction<Real> {
143  public:
145  const MAST::FieldFunction<Real>& hz,
146  const MAST::FieldFunction<Real>& hy_offset,
147  const MAST::FieldFunction<Real>& hz_offset):
148  MAST::FieldFunction<Real>("PolarInertia"),
149  _hy(hy),
150  _hz(hz),
151  _hy_offset(hy_offset),
152  _hz_offset(hz_offset) {
153  _functions.insert(&hy);
154  _functions.insert(&hz);
155  _functions.insert(&hy_offset);
156  _functions.insert(&hz_offset);
157  }
158 
159  virtual ~PolarInertia() { }
160 
161  virtual void operator() (const libMesh::Point& p,
162  const Real t,
163  Real& m) const {
164  Real hy, hz, offy, offz;
165  _hy(p, t, hy);
166  _hz(p, t, hz);
167  _hy_offset(p, t, offy);
168  _hz_offset(p, t, offz);
169 
170  m = hy*hz*((pow(hy,2) + pow(hz,2))/12. +
171  pow(offy,2) + pow(offz,2));
172  }
173 
174  virtual void derivative ( const MAST::FunctionBase& f,
175  const libMesh::Point& p,
176  const Real t,
177  Real& m) const {
178  Real hy, hz, dhy, dhz, offy, offz, doffy, doffz;
179  _hy (p, t, hy); _hy.derivative( f, p, t, dhy);
180  _hz (p, t, hz); _hz.derivative( f, p, t, dhz);
181  _hy_offset (p, t, offy); _hy_offset.derivative( f, p, t, doffy);
182  _hz_offset (p, t, offz); _hz_offset.derivative( f, p, t, doffz);
183 
184 
185  m =
186  (dhy*hz + hy*dhz) * ((pow(hy,2) + pow(hz,2))/12. +
187  pow(offy,2) + pow(offz,2)) +
188  2.*hy*hz*((hy*dhy + hz*dhz)/12. +
189  offy*doffy + offz*doffz);
190  }
191 
192  protected:
193 
194  const MAST::FieldFunction<Real>& _hy, &_hz, &_hy_offset, &_hz_offset;
195  };
196 
197 
198 
203  class AreaYMoment: public MAST::FieldFunction<Real> {
204  public:
206  const MAST::FieldFunction<Real>& hz,
207  const MAST::FieldFunction<Real>& hz_offset):
208  MAST::FieldFunction<Real>("AreaYMoment"),
209  _hy(hy),
210  _hz(hz),
211  _hz_offset(hz_offset) {
212  _functions.insert(&hy);
213  _functions.insert(&hz);
214  _functions.insert(&hz_offset);
215  }
216 
217 
218  virtual ~AreaYMoment() { }
219 
220  virtual void operator() (const libMesh::Point& p,
221  const Real t,
222  Real& m) const {
223  Real hy, hz, off;
224  _hy(p, t, hy);
225  _hz(p, t, hz);
226  _hz_offset(p, t, off);
227 
228  m = hy*hz*off;
229  }
230 
231 
232  virtual void derivative ( const MAST::FunctionBase& f,
233  const libMesh::Point& p,
234  const Real t,
235  Real& m) const {
236  Real hy, hz, off, dhy, dhz, doff;
237  _hy (p, t, hy); _hy.derivative( f, p, t, dhy);
238  _hz (p, t, hz); _hz.derivative( f, p, t, dhz);
239  _hz_offset (p, t, off); _hz_offset.derivative( f, p, t, doff);
240 
241  m = dhy*hz*off + hy*dhz*off + hy*hz*doff;
242  }
243 
244  protected:
245 
247  };
248 
249 
250 
255  class AreaZMoment: public MAST::FieldFunction<Real> {
256  public:
258  const MAST::FieldFunction<Real>& hz,
259  const MAST::FieldFunction<Real>& hy_offset):
260  MAST::FieldFunction<Real>("AreaZMoment"),
261  _hy(hy),
262  _hz(hz),
263  _hy_offset(hy_offset) {
264  _functions.insert(&hy);
265  _functions.insert(&hz);
266  _functions.insert(&hy_offset);
267  }
268 
269  virtual ~AreaZMoment() { }
270 
271  virtual void operator() (const libMesh::Point& p,
272  const Real t,
273  Real& m) const {
274  Real hy, hz, off;
275  _hy(p, t, hy);
276  _hz(p, t, hz);
277  _hy_offset(p, t, off);
278 
279  m = hy*hz*off;
280  }
281 
282 
283  virtual void derivative ( const MAST::FunctionBase& f,
284  const libMesh::Point& p,
285  const Real t,
286  Real& m) const {
287  Real hy, hz, off, dhy, dhz, doff;
288  _hy(p, t, hy); _hy.derivative( f, p, t, dhy);
289  _hz(p, t, hz); _hz.derivative( f, p, t, dhz);
290  _hy_offset(p, t, off); _hy_offset.derivative( f, p, t, doff);
291 
292  m = dhy*hz*off + hy*dhz*off + hy*hz*doff;
293  }
294 
295  protected:
296 
297  const MAST::FieldFunction<Real>& _hy, &_hz, &_hy_offset;
298  };
299 
300 
301 
311  class AreaInertiaMatrix: public MAST::FieldFunction<RealMatrixX> {
312  public:
314  const MAST::FieldFunction<Real>& hz,
315  const MAST::FieldFunction<Real>& hy_offset,
316  const MAST::FieldFunction<Real>& hz_offset):
317  MAST::FieldFunction<RealMatrixX>("AreaInertiaMatrix"),
318  _hy(hy),
319  _hz(hz),
320  _hy_offset(hy_offset),
321  _hz_offset(hz_offset) {
322  _functions.insert(&hy);
323  _functions.insert(&hz);
324  _functions.insert(&hy_offset);
325  _functions.insert(&hz_offset);
326  }
327 
328  virtual ~AreaInertiaMatrix() { }
329 
330  virtual void operator() (const libMesh::Point& p,
331  const Real t,
332  RealMatrixX& m) const {
333  Real hy, hz, offy, offz;
334  m = RealMatrixX::Zero(2,2);
335  _hy(p, t, hy);
336  _hz(p, t, hz);
337  _hy_offset(p, t, offy);
338  _hz_offset(p, t, offz);
339 
340  m(0,0) = hz*pow(hy,3)/12. + hy*hz*pow(offy,2) ; // Izz for v-bending
341  m(0,1) = hy*hz*offy*offz;
342  m(1,0) = m(0,1);
343  m(1,1) = hy*pow(hz,3)/12. + hy*hz*pow(offz,2) ; // Iyy for w-bending
344  }
345 
346 
347  virtual void derivative ( const MAST::FunctionBase& f,
348  const libMesh::Point& p,
349  const Real t,
350  RealMatrixX& m) const {
351  Real hy, hz, offy, offz, dhy, dhz, doffy, doffz;
352  m = RealMatrixX::Zero(2,2);
353  _hy(p, t, hy); _hy.derivative( f, p, t, dhy);
354  _hz(p, t, hz); _hz.derivative( f, p, t, dhz);
355  _hy_offset(p, t, offy); _hy_offset.derivative( f, p, t, doffy);
356  _hz_offset(p, t, offz); _hz_offset.derivative( f, p, t, doffz);
357 
358 
359  m(0,0) = dhz*pow(hy,3)/12. + hz*pow(hy,2)/4.*dhy +
360  dhy*hz*pow(offy,2) + hy*dhz*pow(offy,2) + 2.*hy*hz*offy*doffy ;
361  m(0,1) = dhy*hz*offy*offz + hy*dhz*offy*offz +
362  hy*hz*doffy*offz + hy*hz*offy*doffz;
363  m(1,0) = m(0,1);
364  m(1,1) = dhy*pow(hz,3)/12. + hy*pow(hz,2)/4.*dhz +
365  dhy*hz*pow(offz,2) + hy*dhz*pow(offz,2) + 2.*hy*hz*offz*doffz ;
366  }
367 
368  protected:
369 
370  const MAST::FieldFunction<Real>& _hy, &_hz, &_hy_offset, &_hz_offset;
371  };
372 
373 
374  class ExtensionStiffnessMatrix: public MAST::FieldFunction<RealMatrixX> {
375  public:
377  const MAST::FieldFunction<Real>& A,
378  const MAST::FieldFunction<Real>& J);
379 
381 
382  virtual void operator() (const libMesh::Point& p,
383  const Real t,
384  RealMatrixX& m) const;
385 
386 
387 
388  virtual void derivative ( const MAST::FunctionBase& f,
389  const libMesh::Point& p,
390  const Real t,
391  RealMatrixX& m) const;
392 
393  protected:
394 
397  };
398 
399 
400 
402  public:
404  const MAST::FieldFunction<Real>& A_y_moment,
405  const MAST::FieldFunction<Real>& A_z_moment);
406 
408 
409  virtual void operator() (const libMesh::Point& p,
410  const Real t,
411  RealMatrixX& m) const;
412 
413 
414 
415  virtual void derivative ( const MAST::FunctionBase& f,
416  const libMesh::Point& p,
417  const Real t,
418  RealMatrixX& m) const;
419 
420  protected:
421 
424  };
425 
426 
427  class BendingStiffnessMatrix: public MAST::FieldFunction<RealMatrixX> {
428  public:
431 
433 
434  virtual void operator() (const libMesh::Point& p,
435  const Real t,
436  RealMatrixX& m) const;
437 
438 
439 
440  virtual void derivative ( const MAST::FunctionBase& f,
441  const libMesh::Point& p,
442  const Real t,
443  RealMatrixX& m) const;
444 
445  protected:
446 
449  };
450 
451 
452  class TransverseStiffnessMatrix: public MAST::FieldFunction<RealMatrixX> {
453  public:
455  const MAST::FieldFunction<Real>& A):
456  MAST::FieldFunction<RealMatrixX>("TransverseStiffnessMatrix1D"),
457  _material_stiffness(mat),
458  _A(A) {
459  _functions.insert(&mat);
460  _functions.insert(&A);
461  }
462 
464 
465  virtual void operator() (const libMesh::Point& p,
466  const Real t,
467  RealMatrixX& m) const {
468  Real A;
469  _A(p, t, A);
470  _material_stiffness(p, t, m);
471  m *= A;
472  }
473 
474  virtual void derivative ( const MAST::FunctionBase& f,
475  const libMesh::Point& p,
476  const Real t,
477  RealMatrixX& m) const {
478  RealMatrixX dm;
479  Real A, dA;
480  _A (p, t, A); _A.derivative( f, p, t, dA);
481  _material_stiffness (p, t, m); _material_stiffness.derivative( f, p, t, dm);
482 
483  m *= dA;
484  m += A*dm;
485  }
486 
487  protected:
488 
491  };
492 
493 
494  class InertiaMatrix: public MAST::FieldFunction<RealMatrixX> {
495  public:
497  const MAST::FieldFunction<Real>& A,
498  const MAST::FieldFunction<Real>& A_y_moment,
499  const MAST::FieldFunction<Real>& A_z_moment,
500  const MAST::FieldFunction<Real>& Ip,
502 
503  virtual ~InertiaMatrix() { }
504 
505  virtual void operator() (const libMesh::Point& p,
506  const Real t,
507  RealMatrixX& m) const;
508 
509 
510 
511  virtual void derivative ( const MAST::FunctionBase& f,
512  const libMesh::Point& p,
513  const Real t,
514  RealMatrixX& m) const;
515 
516  protected:
517 
518  const MAST::FieldFunction<Real>& _rho, &_A, &_A_y_moment, &_A_z_moment, &_Ip;
520  };
521 
522 
523 
524  class ThermalExpansionAMatrix: public MAST::FieldFunction<RealMatrixX> {
525  public:
527  const MAST::FieldFunction<RealMatrixX>& mat_expansion,
528  const MAST::FieldFunction<Real>& A);
529 
531 
532  virtual void operator() (const libMesh::Point& p,
533  const Real t,
534  RealMatrixX& m) const;
535 
536 
537 
538  virtual void derivative ( const MAST::FunctionBase& f,
539  const libMesh::Point& p,
540  const Real t,
541  RealMatrixX& m) const;
542 
543  protected:
544 
548  };
549 
550 
551 
552  class ThermalExpansionBMatrix: public MAST::FieldFunction<RealMatrixX> {
553  public:
555  const MAST::FieldFunction<RealMatrixX>& mat_expansion,
556  const MAST::FieldFunction<Real>& A_y_moment,
557  const MAST::FieldFunction<Real>& A_z_moment);
558 
560 
561  virtual void operator() (const libMesh::Point& p,
562  const Real t,
563  RealMatrixX& m) const;
564 
565 
566 
567  virtual void derivative ( const MAST::FunctionBase& f,
568  const libMesh::Point& p,
569  const Real t,
570  RealMatrixX& m) const;
571 
572  protected:
573 
577  };
578 
579 
580 
581 
583  public MAST::FieldFunction<RealMatrixX> {
584  public:
587  const MAST::FieldFunction<Real>& A);
588 
589  virtual ~PrestressAMatrix() { }
590 
591  virtual void operator() (const libMesh::Point& p,
592  const Real t,
593  RealMatrixX& m) const;
594 
595 
596 
597  virtual void derivative ( const MAST::FunctionBase& f,
598  const libMesh::Point& p,
599  const Real t,
600  RealMatrixX& m) const;
601 
602  //virtual void convert_to_vector(const RealMatrixX& m, DenseRealVector& v) const;
603 
604  protected:
605 
608  };
609 
610 
611 
613  public MAST::FieldFunction<RealMatrixX> {
614  public:
617  const MAST::FieldFunction<Real>& A_y_moment,
618  const MAST::FieldFunction<Real>& A_z_moment);
619 
620  virtual ~PrestressBMatrix() { }
621 
622  virtual void operator() (const libMesh::Point& p,
623  const Real t,
624  RealMatrixX& m) const;
625 
626 
627 
628  virtual void derivative ( const MAST::FunctionBase& f,
629  const libMesh::Point& p,
630  const Real t,
631  RealMatrixX& m) const;
632 
633  //virtual void convert_to_vector(const RealMatrixX& m, DenseRealVector& v) const;
634 
635  protected:
636 
639  };
640 
641 
642 
643 
645  public MAST::FieldFunction<RealMatrixX> {
646 
647  public:
648 
650  const MAST::FieldFunction<Real>& A);
651 
652  virtual ~ThermalConductanceMatrix();
653 
654  virtual void operator() (const libMesh::Point& p,
655  const Real t,
656  RealMatrixX& m) const;
657 
658 
659 
660  virtual void derivative ( const MAST::FunctionBase& f,
661  const libMesh::Point& p,
662  const Real t,
663  RealMatrixX& m) const;
664 
665  protected:
666 
668 
670  };
671 
672 
673 
674 
676  public MAST::FieldFunction<RealMatrixX> {
677 
678  public:
679 
681  const MAST::FieldFunction<Real>& h);
682 
683  virtual ~ThermalCapacitanceMatrix();
684 
685  virtual void operator() (const libMesh::Point& p,
686  const Real t,
687  RealMatrixX& m) const;
688 
689 
690 
691  virtual void derivative ( const MAST::FunctionBase& f,
692  const libMesh::Point& p,
693  const Real t,
694  RealMatrixX& m) const;
695 
696  protected:
697 
699 
701  };
702  }
703 
704 
705 }
706 
707 
711  const MAST::FieldFunction<Real>& A,
712  const MAST::FieldFunction<Real>& J):
713 MAST::FieldFunction<RealMatrixX> ("ExtensionStiffnessMatrix1D"),
714 _material_stiffness(mat),
715 _A(A),
716 _J(J) {
717  _functions.insert(&mat);
718  _functions.insert(&A);
719  _functions.insert(&J);
720 }
721 
722 
723 
724 
725 
726 void
729  const Real t,
730  RealMatrixX& m) const {
731  // [C]*h
732  Real A, J;
733  _A(p, t, A);
734  _J(p, t, J);
735  _material_stiffness(p, t, m);
736  m.row(0) *= A;
737  m.row(1) *= J;
738 }
739 
740 
741 
742 
743 void
746  const libMesh::Point& p,
747  const Real t,
748  RealMatrixX& m) const {
749  RealMatrixX dm;
750  m = RealMatrixX::Zero(2,2);
751  Real A, J, dA, dJ;
752  _A(p, t, A); _A.derivative( f, p, t, dA);
753  _J(p, t, J); _J.derivative( f, p, t, dJ);
754  _material_stiffness(p, t, m); _material_stiffness.derivative( f, p, t, dm);
755 
756  // [C]*dh
757  m.row(0) *= dA;
758  m.row(1) *= dJ;
759 
760  // += [dC]*h
761  dm.row(0) *= A;
762  dm.row(1) *= J;
763  m += dm;
764 }
765 
766 
767 
768 
769 
770 
773  const MAST::FieldFunction<Real>& A_y_moment,
774  const MAST::FieldFunction<Real>& A_z_moment):
775 MAST::FieldFunction<RealMatrixX> ("ExtensionBendingStiffnessMatrix1D"),
777 _A_y_moment(A_y_moment),
778 _A_z_moment(A_z_moment) {
779  _functions.insert(&mat);
780  _functions.insert(&A_y_moment);
781  _functions.insert(&A_z_moment);
782 }
783 
784 
785 
786 void
789  const Real t,
790  RealMatrixX& m) const {
791  Real Ay, Az;
792  _A_y_moment(p, t, Ay);
793  _A_z_moment(p, t, Az);
794  _material_stiffness(p, t, m);
795 
796  m(0,1) = m(0,0)*Ay; // coupling of u and w bending (== theta_y)
797  m(0,0) *= Az; // coupling of u and v bending (== theta_z)
798 
799  m.row(1) *= 0; // no coupling for torsion for symmetic sections
800 }
801 
802 
803 
804 
805 
806 void
809  const libMesh::Point& p,
810  const Real t,
811  RealMatrixX& m) const {
812  RealMatrixX dm;
813  Real Ay, Az, dAy, dAz;
814  _A_y_moment(p, t, Ay); _A_y_moment.derivative( f, p, t, dAy);
815  _A_z_moment(p, t, Az); _A_z_moment.derivative( f, p, t, dAz);
816  _material_stiffness(p, t, m); _material_stiffness.derivative( f, p, t, dm);
817 
818  m(0,1) = m(0,0)*dAy; // coupling of u and w bending (== theta_y)
819  m(0,0) *= dAz; // coupling of u and v bending (== theta_z)
820  m.row(1) *= 0; // no coupling for torsion for symmetic sections
821 
822  dm(0,1) = dm(0,0)*Ay;
823  dm(0,0) *= Az;
824  dm.row(1)*= 0.;
825  m += dm;
826 }
827 
828 
829 
830 
834 MAST::FieldFunction<RealMatrixX> ("BendingStiffnessMatrix1D"),
836 _I(I) {
837  _functions.insert(&mat);
838  _functions.insert(&I);
839 }
840 
841 
842 
843 void
845 BendingStiffnessMatrix::operator() (const libMesh::Point& p,
846  const Real t,
847  RealMatrixX& m) const {
848  RealMatrixX mat;
849  _I(p, t, m);
850  _material_stiffness(p, t, mat);
851 
852  // E*I
853  m *= mat(0,0); // scale the inertia matrix with modulus of elasticity
854 }
855 
856 
857 
858 
859 
860 
861 
862 
863 void
866  const libMesh::Point& p,
867  const Real t,
868  RealMatrixX& m) const {
869  RealMatrixX mat, dmat, dm;
870  _I(p, t, m); _I.derivative( f, p, t, dm);
871  _material_stiffness(p, t, mat); _material_stiffness.derivative( f, p, t, dmat);
872 
873  // dE*I
874  m *= dmat(0,0); // scale the inertia matrix with modulus of elasticity
875 
876  // E*dI
877  m += mat(0,0)*dm; // scale the inertia matrix with modulus of elasticity
878 }
879 
880 
881 
882 
883 
886  const MAST::FieldFunction<Real>& A,
887  const MAST::FieldFunction<Real>& A_y_moment,
888  const MAST::FieldFunction<Real>& A_z_moment,
889  const MAST::FieldFunction<Real>& Ip,
891 MAST::FieldFunction<RealMatrixX>("InertiaMatrix1D"),
892 _rho(rho),
893 _A(A),
894 _A_y_moment(A_y_moment),
895 _A_z_moment(A_z_moment),
896 _Ip(Ip),
897 _I(I) {
898  _functions.insert(&_rho);
899  _functions.insert(&_A);
900  _functions.insert(&_A_y_moment);
901  _functions.insert(&_A_z_moment);
902  _functions.insert(&_Ip);
903  _functions.insert(&_I);
904 }
905 
906 
907 
908 
909 void
911 InertiaMatrix::operator() (const libMesh::Point& p,
912  const Real t,
913  RealMatrixX& m) const {
914  m = RealMatrixX::Zero(6, 6);
915  RealMatrixX I;
916  Real rho, A, Ay, Az, Ip;
917  _rho(p, t, rho);
918  _A(p, t, A);
919  _A_y_moment(p, t, Ay);
920  _A_z_moment(p, t, Az);
921  _Ip(p, t, Ip);
922  _I(p, t, I);
923 
924  // translation velocities
925  m(0,0) = A; m(1,1) = A; m(2,2) = A;
926 
927  // torsion
928  m(3,3) = Ip;
929 
930  // rotational velocities
931  m(0,4) = Ay; m(4,0) = Ay; // w-displacement
932  m(0,5) = -Az; m(5,0) = -Az; // v-displacement
933 
934  // bending rotation inertia
935  for (unsigned int i=0; i<2; i++)
936  for (unsigned int j=0; j<2; j++)
937  m(4+i,4+j) = I(i,j)*1.e-6;
938 
939  m *= rho;
940 }
941 
942 
943 
944 
945 
946 
947 
948 
949 void
952  const libMesh::Point& p,
953  const Real t,
954  RealMatrixX& m) const {
955  RealMatrixX dm;
956  m = RealMatrixX::Zero(6, 6); dm = RealMatrixX::Zero(6, 6);
957  RealMatrixX I, dI;
958  Real rho, A, Ay, Az, Ip, drho, dA, dAy, dAz, dIp;
959  _rho(p, t, rho); _rho.derivative( f, p, t, drho);
960  _A(p, t, A); _A.derivative( f, p, t, dA);
961  _A_y_moment(p, t, Ay); _A_y_moment.derivative( f, p, t, dAy);
962  _A_z_moment(p, t, Az); _A_z_moment.derivative( f, p, t, dAz);
963  _Ip(p, t, Ip); _Ip.derivative( f, p, t, dIp);
964  _I(p, t, I); _I.derivative( f, p, t, dI);
965 
966  // translation velocities
967  m(0,0) = A; m(1,1) = A; m(2,2) = A;
968  dm(0,0) = dA; dm(1,1) = dA; dm(2,2) = dA;
969 
970  // torsion
971  m(3,3) = Ip;
972  dm(3,3) = dIp;
973 
974  // rotational velocities
975  m(0,4) = Ay; m(4,0) = Ay; // w-displacement
976  dm(0,4) = dAy; dm(4,0) = dAy; // w-displacement
977  m(0,5) = -Az; m(5,0) = -Az; // v-displacement
978  dm(0,5) = -dAz; m(5,0) = -dAz; // v-displacement
979 
980  // bending rotation inertia
981  for (unsigned int i=0; i<2; i++)
982  for (unsigned int j=0; j<2; j++) {
983  m(4+i,4+j) = I(i,j)*1.e-6;
984  dm(4+i,4+j) = dI(i,j)*1.e-6;
985  }
986 
987  m *= drho;
988  m += rho*dm;
989 }
990 
991 
992 
993 
996  const MAST::FieldFunction<RealMatrixX>& mat_expansion,
997  const MAST::FieldFunction<Real>& A):
998 MAST::FieldFunction<RealMatrixX>("ThermalExpansionAMatrix1D"),
999 _material_stiffness(mat_stiff),
1000 _material_expansion(mat_expansion),
1001 _A(A) {
1002  _functions.insert(&mat_stiff);
1003  _functions.insert(&mat_expansion);
1004  _functions.insert(&_A);
1005 }
1006 
1007 
1008 
1009 
1010 void
1013  const Real t,
1014  RealMatrixX& m) const {
1015  Real A;
1016  RealMatrixX at;
1017  _A(p, t, A);
1018  _material_stiffness(p, t, m);
1019  _material_expansion(p, t, at);
1020 
1021  m *= at;
1022  m *= A;
1023 }
1024 
1025 
1026 
1027 
1028 
1029 void
1032  const libMesh::Point& p,
1033  const Real t,
1034  RealMatrixX& m) const {
1035  Real A, dA;
1036  RealMatrixX m1, at, dat, dm;
1037  _A(p, t, A); _A.derivative( f, p, t, dA);
1038  _material_stiffness(p, t, m1); _material_stiffness.derivative( f, p, t, dm);
1039  _material_expansion(p, t, at); _material_expansion.derivative( f, p, t, dat);
1040 
1041  m=m1;
1042 
1043  m *= at;
1044  m *= dA;
1045 
1046  m1 *= dat;
1047  dm *= at;
1048  m1 += dm;
1049 
1050  m += A*m1;
1051 }
1052 
1053 
1054 
1055 
1058  const MAST::FieldFunction<RealMatrixX>& mat_expansion,
1059  const MAST::FieldFunction<Real>& A_y_moment,
1060  const MAST::FieldFunction<Real>& A_z_moment):
1061 MAST::FieldFunction<RealMatrixX>("ThermalExpansionBMatrix1D"),
1062 _material_stiffness(mat_stiff),
1063 _material_expansion(mat_expansion),
1064 _A_y_moment(A_y_moment),
1065 _A_z_moment(A_z_moment) {
1066  _functions.insert(&mat_stiff);
1067  _functions.insert(&mat_expansion);
1068  _functions.insert(&_A_y_moment);
1069  _functions.insert(&_A_z_moment);
1070 }
1071 
1072 
1073 
1074 
1075 void
1078  const Real t,
1079  RealMatrixX& m) const {
1080  Real Ay, Az;
1081  RealMatrixX at;
1082  _A_y_moment(p, t, Ay);
1083  _A_z_moment(p, t, Az);
1084  _material_stiffness(p, t, m);
1085  _material_expansion(p, t, at);
1086 
1087  m *= at;
1088  m(1,0) = Ay * m(0,0); // for w-displacement, area moment about Y-axis
1089  m(0,0) *= Az; // for v-displacement, area moment about Z-axis
1090 }
1091 
1092 
1093 
1094 
1095 
1096 
1097 void
1100  const libMesh::Point& p,
1101  const Real t,
1102  RealMatrixX& m) const {
1103  Real Ay, Az, dAy, dAz;
1104  RealMatrixX at, dat, m1, dm;
1105  _A_y_moment(p, t, Ay); _A_y_moment.derivative( f, p, t, dAy);
1106  _A_z_moment(p, t, Az); _A_z_moment.derivative( f, p, t, dAz);
1107  _material_stiffness(p, t, m1); _material_stiffness.derivative( f, p, t, dm);
1108  _material_expansion(p, t, at); _material_expansion.derivative( f, p, t, dat);
1109 
1110  m = m1;
1111  m *= at;
1112  m(1,0) = dAy * m(0,0);
1113  m(0,0) *= dAz;
1114 
1115  m1 *= dat;
1116  dm *= at;
1117  m1 += dm;
1118  m1(1,0) = Ay * m1(0,0);
1119  m1(0,0) *= Az;
1120 
1121  m += m1;
1122 }
1123 
1124 
1125 
1126 
1130  const MAST::FieldFunction<Real>& A):
1131 MAST::FieldFunction<RealMatrixX>("PrestressAMatrix1D"),
1132 _prestress(prestress),
1133 _T(T),
1134 _A(A) {
1135  _functions.insert(&prestress);
1136  _functions.insert(&T);
1137  _functions.insert(&A);
1138 }
1139 
1140 
1141 
1142 
1143 void
1145 PrestressAMatrix::operator() (const libMesh::Point& p,
1146  const Real t,
1147  RealMatrixX& m) const {
1148  RealMatrixX s, T;
1149  m = RealMatrixX::Zero(2, 2);
1150  Real A;
1151  _A(p, t, A);
1152  _prestress(p, t, s);
1153  _T(p, t, T);
1154 
1155  // convert the stress to the local coordinate
1156  s *= T;
1157  s = T.transpose() * s;
1158 
1159  m(0,0) = s(0,0)*A; // only sigma_xx is applied, and torsion is neglected
1160 }
1161 
1162 
1163 
1164 
1165 
1166 void
1169  const libMesh::Point& p,
1170  const Real t,
1171  RealMatrixX& m) const {
1172  RealMatrixX s, ds, T, dT;
1173  m = RealMatrixX::Zero(2, 2);
1174  Real A, dA;
1175  _A(p, t, A); _A.derivative( f, p, t, dA);
1176  _prestress(p, t, s); _prestress.derivative( f, p, t, ds);
1177  _T(p, t, T); _T.derivative( f, p, t, dT);
1178 
1179  // convert the stress to the local coordinate
1180  s *= T;
1181  s = T.transpose() * s;
1182 
1183  // ds = dT^T s T + T^T s dT + T^T ds T
1184  RealMatrixX tmp;
1185  ds * T;
1186  ds = T.transpose()*ds;
1187 
1188  tmp = T.transpose() * s * dT + dT.transpose() * s * T;
1189  ds += tmp;
1190 
1191  m(0,0) = s(0,0)*dA + ds(0,0)*A; // only sigma_xx is applied, and torsion is neglected
1192 }
1193 
1194 
1195 
1196 
1200  const MAST::FieldFunction<Real>& A_y_moment,
1201  const MAST::FieldFunction<Real>& A_z_moment):
1202 MAST::FieldFunction<RealMatrixX>("PrestressBMatrix1D"),
1203 _prestress(prestress),
1204 _T(T),
1205 _A_y_moment(A_y_moment),
1206 _A_z_moment(A_z_moment) {
1207  _functions.insert(&prestress);
1208  _functions.insert(&T);
1209  _functions.insert(&A_y_moment);
1210  _functions.insert(&A_z_moment);
1211 }
1212 
1213 
1214 
1215 
1216 void
1218 PrestressBMatrix::operator() (const libMesh::Point& p,
1219  const Real t,
1220  RealMatrixX& m) const {
1221  RealMatrixX s, T;
1222  m = RealMatrixX::Zero(2, 2);
1223  Real Ay, Az;
1224  _A_y_moment(p, t, Ay);
1225  _A_z_moment(p, t, Az);
1226  _prestress(p, t, s);
1227  _T(p, t, T);
1228 
1229  // convert the stress to the local coordinate
1230  s = T.transpose() * s * T;
1231 
1232  // only sigma_xx is applied, and torsion is neglected
1233  m(0,0) = s(0,0)*Az;
1234  m(0,1) = s(0,0)*Ay;
1235 }
1236 
1237 
1238 
1239 
1240 
1241 
1242 void
1245  const libMesh::Point& p,
1246  const Real t,
1247  RealMatrixX& m) const {
1248  RealMatrixX s, ds, T, dT;
1249  m = RealMatrixX::Zero(2, 2);
1250  Real Ay, Az, dAy, dAz;
1251  _A_y_moment(p, t, Ay); _A_y_moment.derivative( f, p, t, dAy);
1252  _A_z_moment(p, t, Az); _A_z_moment.derivative( f, p, t, dAz);
1253  _prestress(p, t, s); _prestress.derivative( f, p, t, ds);
1254  _T(p, t, T); _T.derivative( f, p, t, dT);
1255 
1256  // convert the stress to the local coordinate
1257  s = T.transpose() * s * T;
1258 
1259  // ds = dT^T s T + T^T s dT + T^T ds T
1260  RealMatrixX tmp;
1261  ds = (T.transpose() * ds * T +
1262  T.transpose() * s * dT +
1263  dT.transpose() * s * T);
1264 
1265  // only sigma_xx is applied, and torsion is neglected
1266  m(0,0) = (s(0,0)*dAz + ds(0,0)*Az);
1267  m(0,1) = s(0,0)*dAy + ds(0,0)*Ay;
1268 }
1269 
1270 
1271 
1272 
1273 
1274 
1277  const MAST::FieldFunction<Real>& A):
1278 MAST::FieldFunction<RealMatrixX>("ThermalConductanceMatrix"),
1279 _mat_cond(mat_cond),
1280 _A(A) {
1281  _functions.insert(&mat_cond);
1282  _functions.insert(&A);
1283 }
1284 
1285 
1288 
1289 
1290 void
1292 operator() (const libMesh::Point& p,
1293  const Real t,
1294  RealMatrixX& m) const {
1295 
1296  m = RealMatrixX::Zero(1, 1);
1297  Real A;
1298  _mat_cond(p, t, m);
1299  _A(p, t, A);
1300 
1301  m *= A;
1302 }
1303 
1304 
1305 
1306 void
1308  const libMesh::Point& p,
1309  const Real t,
1310  RealMatrixX& m) const {
1311  m = RealMatrixX::Zero(1, 1);
1312  RealMatrixX dm;
1313  Real A, dA;
1314  _mat_cond(p, t, m);
1315  _mat_cond.derivative( f, p, t, dm);
1316  _A(p, t, A);
1317  _A.derivative( f, p, t, dA);
1318 
1319  m *= dA;
1320  m += dm*A;
1321 }
1322 
1323 
1324 
1325 
1328  const MAST::FieldFunction<Real>& h):
1329 MAST::FieldFunction<RealMatrixX>("ThermalCapacitanceMatrix"),
1330 _mat_cap(mat_cap),
1331 _h(h) {
1332  _functions.insert(&mat_cap);
1333  _functions.insert(&h);
1334 }
1335 
1336 
1337 
1338 
1341 
1342 
1343 void
1345 operator() (const libMesh::Point& p,
1346  const Real t,
1347  RealMatrixX& m) const {
1348 
1349  m = RealMatrixX::Zero(1, 1);
1350  Real h;
1351  _mat_cap(p, t, m);
1352  _h(p, t, h);
1353 
1354  m *= h;
1355 }
1356 
1357 
1358 
1359 void
1361  const libMesh::Point& p,
1362  const Real t,
1363  RealMatrixX& m) const {
1364  m = RealMatrixX::Zero(1, 1);
1365  RealMatrixX dm;
1366  Real h, dh;
1367  _mat_cap(p, t, m);
1368  _mat_cap.derivative( f, p, t, dm);
1369  _h(p, t, h);
1370  _h.derivative( f, p, t, dh);
1371 
1372  m *= dh;
1373  m += dm*h;
1374 }
1375 
1376 
1377 
1378 
1379 bool
1382  return _material->depends_on(f) || // check if the material property depends on the function
1383  MAST::ElementPropertyCardBase::depends_on(f); // check with this property card
1384 }
1385 
1386 
1387 
1390 
1391  libmesh_assert(_initialized);
1392  return *_A;
1393 }
1394 
1395 
1396 
1399 
1400  libmesh_assert(_initialized);
1401  return *_J;
1402 }
1403 
1404 
1407 
1408  libmesh_assert(_initialized);
1409  return *_Ip;
1410 }
1411 
1412 
1415 
1416  libmesh_assert(_initialized);
1417  return *_Ay;
1418 }
1419 
1420 
1423 
1424  libmesh_assert(_initialized);
1425  return *_Az;
1426 }
1427 
1428 
1431 
1432  libmesh_assert(_initialized);
1433  return *_AI;
1434 }
1435 
1436 
1437 
1438 
1441 
1442  libmesh_assert(_initialized);
1443  return *_A;
1444 }
1445 
1446 
1447 
1450 
1451  libmesh_assert(_initialized);
1452  return *_J;
1453 }
1454 
1455 
1458 
1459  libmesh_assert(_initialized);
1460  return *_Ip;
1461 }
1462 
1463 
1466 
1467  libmesh_assert(_initialized);
1468  return *_Ay;
1469 }
1470 
1471 
1474 
1475  libmesh_assert(_initialized);
1476  return *_Az;
1477 }
1478 
1479 
1482 
1483  libmesh_assert(_initialized);
1484  return *_AI;
1485 }
1486 
1487 
1488 
1489 
1490 
1491 void
1493 
1494  libmesh_assert(!_initialized);
1495 
1496  _A.reset();
1497  _Ay.reset();
1498  _Az.reset();
1499  _J.reset();
1500  _Ip.reset();
1501  _AI.reset();
1502 
1503  _initialized = false;
1504 }
1505 
1506 
1507 
1508 void
1510 
1511  libmesh_assert(!_initialized);
1512 
1514  &hy = this->get<MAST::FieldFunction<Real> >("hy"),
1515  &hz = this->get<MAST::FieldFunction<Real> >("hz"),
1516  &hy_off = this->get<MAST::FieldFunction<Real> >("hy_off"),
1517  &hz_off = this->get<MAST::FieldFunction<Real> >("hz_off");
1518 
1519  _A.reset(new MAST::Solid1DSectionProperty::Area(hy,
1520  hz));
1522  hz,
1523  hz_off));
1525  hz,
1526  hy_off));
1528  hz));
1530  hz,
1531  hy_off,
1532  hz_off));
1534  hz,
1535  hy_off,
1536  hz_off));
1537 
1538  _initialized = true;
1539 }
1540 
1541 
1542 std::unique_ptr<MAST::FieldFunction<RealMatrixX> >
1545 
1546  // make sure that the init method has been called on the card
1547  libmesh_assert(_initialized);
1548 
1551  (_material->stiffness_matrix(1), *_A, *_J);
1552 
1553  return std::unique_ptr<MAST::FieldFunction<RealMatrixX> > (rval);
1554 }
1555 
1556 
1557 std::unique_ptr<MAST::FieldFunction<RealMatrixX> >
1560 
1561  // make sure that the init method has been called on the card
1562  libmesh_assert(_initialized);
1563 
1566  (_material->stiffness_matrix(1), *_Ay, *_Az);
1567 
1568  return std::unique_ptr<MAST::FieldFunction<RealMatrixX> > (rval);
1569 }
1570 
1571 
1572 std::unique_ptr<MAST::FieldFunction<RealMatrixX> >
1575 
1576  // make sure that the init method has been called on the card
1577  libmesh_assert(_initialized);
1578 
1581  (_material->stiffness_matrix(1), *_AI);
1582 
1583  return std::unique_ptr<MAST::FieldFunction<RealMatrixX> > (rval);
1584 }
1585 
1586 
1587 std::unique_ptr<MAST::FieldFunction<RealMatrixX> >
1590 
1591  libmesh_error();
1592 
1593  return std::unique_ptr<MAST::FieldFunction<RealMatrixX> > (nullptr);
1594 }
1595 
1596 
1597 
1598 std::unique_ptr<MAST::FieldFunction<RealMatrixX> >
1601 
1602  // make sure that the init method has been called on the card
1603  libmesh_assert(_initialized);
1604 
1607  (_material->get<FieldFunction<Real> >("rho"),
1608  *_A,
1609  *_Ay,
1610  *_Az,
1611  *_Ip,
1612  *_AI);
1613 
1614  return std::unique_ptr<MAST::FieldFunction<RealMatrixX> > (rval);
1615 }
1616 
1617 
1618 
1619 std::unique_ptr<MAST::FieldFunction<RealMatrixX> >
1622 
1623  // make sure that the init method has been called on the card
1624  libmesh_assert(_initialized);
1625 
1628  (_material->stiffness_matrix(1),
1629  _material->thermal_expansion_matrix(1),
1630  *_A);
1631 
1632  return std::unique_ptr<MAST::FieldFunction<RealMatrixX> > (rval);
1633 }
1634 
1635 
1636 
1637 std::unique_ptr<MAST::FieldFunction<RealMatrixX> >
1640 
1641 
1642  // make sure that the init method has been called on the card
1643  libmesh_assert(_initialized);
1644 
1647  (_material->stiffness_matrix(1),
1648  _material->thermal_expansion_matrix(1),
1649  *_Ay,
1650  *_Az);
1651 
1652  return std::unique_ptr<MAST::FieldFunction<RealMatrixX> > (rval);
1653 }
1654 
1655 
1656 std::unique_ptr<MAST::FieldFunction<RealMatrixX> >
1659 
1660 
1661  // make sure that the init method has been called on the card
1662  libmesh_assert(_initialized);
1663 
1666  (_material->transverse_shear_stiffness_matrix(),
1667  *_A);
1668 
1669  return std::unique_ptr<MAST::FieldFunction<RealMatrixX> > (rval);
1670 }
1671 
1672 
1673 std::unique_ptr<MAST::FieldFunction<RealMatrixX> >
1676 
1677  // make sure that the init method has been called on the card
1678  libmesh_assert(_initialized);
1679 
1681  // TODO: figure out the interface for prestress and T matrix
1682  libmesh_assert(false);
1683  // = new MAST::Solid1DSectionProperty::PrestressAMatrix
1684  //(this->get<MAST::FieldFunction<RealMatrixX> >("prestress"),
1685  // e.local_elem().T_matrix(),
1686  // *_A);
1687 
1688  return std::unique_ptr<MAST::FieldFunction<RealMatrixX> > (rval);
1689 }
1690 
1691 
1692 std::unique_ptr<MAST::FieldFunction<RealMatrixX> >
1695 
1696  // make sure that the init method has been called on the card
1697  libmesh_assert(_initialized);
1698 
1700  // TODO: figure out the interface for prestress and T matrix
1701  libmesh_assert(false);
1702  // = new MAST::Solid1DSectionProperty::PrestressBMatrix
1703  //(this->get<MAST::FieldFunction<RealMatrixX> >("prestress"),
1704  // e.local_elem().T_matrix(),
1705  // *_Ay,
1706  // *_Az);
1707 
1708  return std::unique_ptr<MAST::FieldFunction<RealMatrixX> > (rval);
1709 }
1710 
1711 
1712 
1713 std::unique_ptr<MAST::FieldFunction<RealMatrixX> >
1716 
1717  // make sure that the init method has been called on the card
1718  libmesh_assert(_initialized);
1719 
1722  (_material->conductance_matrix(1),
1723  *_A);
1724 
1725  return std::unique_ptr<MAST::FieldFunction<RealMatrixX> > (rval);
1726 }
1727 
1728 
1729 
1730 std::unique_ptr<MAST::FieldFunction<RealMatrixX> >
1733 
1734  // make sure that the init method has been called on the card
1735  libmesh_assert(_initialized);
1736 
1739  (_material->capacitance_matrix(1),
1740  *_A);
1741 
1742  return std::unique_ptr<MAST::FieldFunction<RealMatrixX> > (rval);
1743 }
1744 
1745 
1748 section(const MAST::ElementBase& e) const {
1749 
1750  return *_A;
1751 }
1752 
virtual void operator()(const libMesh::Point &p, const Real t, RealMatrixX &m) const
calculates the value of the function at the specified point, p, and time, t, and returns it in v...
calculates the area moment about the Y-axis due to an offset along the Z-axis
BendingStiffnessMatrix(const MAST::FieldFunction< RealMatrixX > &mat, const MAST::FieldFunction< RealMatrixX > &I)
ThermalCapacitanceMatrix(const MAST::FieldFunction< RealMatrixX > &mat_cond, const MAST::FieldFunction< Real > &h)
virtual const MAST::FieldFunction< Real > & section(const MAST::ElementBase &e) const
virtual void derivative(const MAST::FunctionBase &f, const libMesh::Point &p, const Real t, Real &m) const
calculates the value of the derivative of function with respect to the function f at the specified po...
virtual const MAST::FieldFunction< Real > & Ay() const
Area(const MAST::FieldFunction< Real > &hy, const MAST::FieldFunction< Real > &hz)
virtual void derivative(const MAST::FunctionBase &f, const libMesh::Point &p, const Real t, RealMatrixX &m) const
calculates the value of the derivative of function with respect to the function f at the specified po...
virtual void derivative(const MAST::FunctionBase &f, const libMesh::Point &p, const Real t, RealMatrixX &m) const
calculates the value of the derivative of function with respect to the function f at the specified po...
virtual void derivative(const MAST::FunctionBase &f, const libMesh::Point &p, const Real t, RealMatrixX &m) const
calculates the value of the derivative of function with respect to the function f at the specified po...
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > thermal_conductance_matrix(const MAST::ElementBase &e) const
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > thermal_expansion_A_matrix(const MAST::ElementBase &e) const
virtual void operator()(const libMesh::Point &p, const Real t, RealMatrixX &m) const
calculates the value of the function at the specified point, p, and time, t, and returns it in v...
virtual const MAST::FieldFunction< Real > & J() const
ExtensionBendingStiffnessMatrix(const MAST::FieldFunction< RealMatrixX > &mat, const MAST::FieldFunction< Real > &A_y_moment, const MAST::FieldFunction< Real > &A_z_moment)
ExtensionStiffnessMatrix(const MAST::FieldFunction< RealMatrixX > &mat, const MAST::FieldFunction< Real > &A, const MAST::FieldFunction< Real > &J)
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > stiffness_A_matrix(const MAST::ElementBase &e) const
AreaZMoment(const MAST::FieldFunction< Real > &hy, const MAST::FieldFunction< Real > &hz, const MAST::FieldFunction< Real > &hy_offset)
virtual void operator()(const libMesh::Point &p, const Real t, RealMatrixX &m) const
calculates the value of the function at the specified point, p, and time, t, and returns it in v...
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > thermal_capacitance_matrix(const MAST::ElementBase &e) const
virtual const MAST::FieldFunction< Real > & Az() const
PrestressBMatrix(const MAST::FieldFunction< RealMatrixX > &prestress, const MAST::FieldFunction< RealMatrixX > &T, const MAST::FieldFunction< Real > &A_y_moment, const MAST::FieldFunction< Real > &A_z_moment)
virtual void derivative(const MAST::FunctionBase &f, const libMesh::Point &p, const Real t, Real &m) const
calculates the value of the derivative of function with respect to the function f at the specified po...
std::set< const MAST::FunctionBase * > _functions
set of functions that this function depends on
virtual void derivative(const MAST::FunctionBase &f, const libMesh::Point &p, const Real t, RealMatrixX &m) const
calculates the value of the derivative of function with respect to the function f at the specified po...
InertiaMatrix(const MAST::FieldFunction< Real > &rho, const MAST::FieldFunction< Real > &A, const MAST::FieldFunction< Real > &A_y_moment, const MAST::FieldFunction< Real > &A_z_moment, const MAST::FieldFunction< Real > &Ip, const MAST::FieldFunction< RealMatrixX > &I)
libMesh::Real Real
virtual void operator()(const libMesh::Point &p, const Real t, Real &m) const
calculates the value of the function at the specified point, p, and time, t, and returns it in v...
PolarInertia(const MAST::FieldFunction< Real > &hy, const MAST::FieldFunction< Real > &hz, const MAST::FieldFunction< Real > &hy_offset, const MAST::FieldFunction< Real > &hz_offset)
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > transverse_shear_stiffness_matrix(const MAST::ElementBase &e) const
TransverseStiffnessMatrix(const MAST::FieldFunction< RealMatrixX > &mat, const MAST::FieldFunction< Real > &A)
virtual void derivative(const MAST::FunctionBase &f, const libMesh::Point &p, const Real t, RealMatrixX &m) const
calculates the value of the derivative of function with respect to the function f at the specified po...
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > thermal_expansion_B_matrix(const MAST::ElementBase &e) const
virtual bool depends_on(const MAST::FunctionBase &f) const
returns true if the property card depends on the function f
calculates the area moment about the Z-axis due to an offset along the Y-axis
virtual void derivative(const MAST::FunctionBase &f, ValType &v) const
calculates the value of the function derivative and returns it in v.
virtual const MAST::FieldFunction< Real > & Ip() const
virtual bool depends_on(const MAST::FunctionBase &f) const
returns true if the property card depends on the function f
ThermalExpansionAMatrix(const MAST::FieldFunction< RealMatrixX > &mat_stiff, const MAST::FieldFunction< RealMatrixX > &mat_expansion, const MAST::FieldFunction< Real > &A)
virtual void derivative(const MAST::FunctionBase &f, const libMesh::Point &p, const Real t, Real &m) const
calculates the value of the derivative of function with respect to the function f at the specified po...
TorsionalConstant(const MAST::FieldFunction< Real > &hy, const MAST::FieldFunction< Real > &hz)
virtual void derivative(const MAST::FunctionBase &f, const libMesh::Point &p, const Real t, RealMatrixX &m) const
calculates the value of the derivative of function with respect to the function f at the specified po...
virtual void operator()(const libMesh::Point &p, const Real t, RealMatrixX &m) const
calculates the value of the function at the specified point, p, and time, t, and returns it in v...
Matrix< Real, Dynamic, Dynamic > RealMatrixX
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > damping_matrix(const MAST::ElementBase &e) const
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > prestress_B_matrix(MAST::ElementBase &e) const
virtual void operator()(const libMesh::Point &p, const Real t, RealMatrixX &m) const
calculates the value of the function at the specified point, p, and time, t, and returns it in v...
AreaInertiaMatrix(const MAST::FieldFunction< Real > &hy, const MAST::FieldFunction< Real > &hz, const MAST::FieldFunction< Real > &hy_offset, const MAST::FieldFunction< Real > &hz_offset)
virtual void derivative(const MAST::FunctionBase &f, const libMesh::Point &p, const Real t, RealMatrixX &m) const
calculates the value of the derivative of function with respect to the function f at the specified po...
ThermalExpansionBMatrix(const MAST::FieldFunction< RealMatrixX > &mat_stiff, const MAST::FieldFunction< RealMatrixX > &mat_expansion, const MAST::FieldFunction< Real > &A_y_moment, const MAST::FieldFunction< Real > &A_z_moment)
This creates the base class for functions that have a saptial and temporal dependence, and provide sensitivity operations with respect to the functions and parameters.
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > stiffness_B_matrix(const MAST::ElementBase &e) const
virtual void derivative(const MAST::FunctionBase &f, const libMesh::Point &p, const Real t, Real &m) const
calculates the value of the derivative of function with respect to the function f at the specified po...
virtual void operator()(const libMesh::Point &p, const Real t, RealMatrixX &m) const
calculates the value of the function at the specified point, p, and time, t, and returns it in v...
calculates the 2x2 matrix of area inertia for the section with individual entries as ...
virtual void derivative(const MAST::FunctionBase &f, const libMesh::Point &p, const Real t, RealMatrixX &m) const
calculates the value of the derivative of function with respect to the function f at the specified po...
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > inertia_matrix(const MAST::ElementBase &e) const
virtual void derivative(const MAST::FunctionBase &f, const libMesh::Point &p, const Real t, RealMatrixX &m) const
calculates the value of the derivative of function with respect to the function f at the specified po...
ThermalConductanceMatrix(const MAST::FieldFunction< RealMatrixX > &mat_cond, const MAST::FieldFunction< Real > &A)
virtual const MAST::FieldFunction< Real > & A() const
virtual void derivative(const MAST::FunctionBase &f, const libMesh::Point &p, const Real t, RealMatrixX &m) const
calculates the value of the derivative of function with respect to the function f at the specified po...
virtual void operator()(const libMesh::Point &p, const Real t, RealMatrixX &m) const
calculates the value of the function at the specified point, p, and time, t, and returns it in v...
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > stiffness_D_matrix(const MAST::ElementBase &e) const
virtual void operator()(const libMesh::Point &p, const Real t, RealMatrixX &m) const
calculates the value of the function at the specified point, p, and time, t, and returns it in v...
virtual void derivative(const MAST::FunctionBase &f, const libMesh::Point &p, const Real t, RealMatrixX &m) const
calculates the value of the derivative of function with respect to the function f at the specified po...
virtual void derivative(const MAST::FunctionBase &f, const libMesh::Point &p, const Real t, RealMatrixX &m) const
calculates the value of the derivative of function with respect to the function f at the specified po...
virtual const MAST::FieldFunction< RealMatrixX > & I() const
PrestressAMatrix(const MAST::FieldFunction< RealMatrixX > &prestress, const MAST::FieldFunction< RealMatrixX > &T, const MAST::FieldFunction< Real > &A)
AreaYMoment(const MAST::FieldFunction< Real > &hy, const MAST::FieldFunction< Real > &hz, const MAST::FieldFunction< Real > &hz_offset)
virtual void operator()(const libMesh::Point &p, const Real t, RealMatrixX &m) const
calculates the value of the function at the specified point, p, and time, t, and returns it in v...
virtual void derivative(const MAST::FunctionBase &f, const libMesh::Point &p, const Real t, Real &m) const
calculates the value of the derivative of function with respect to the function f at the specified po...
virtual void operator()(const libMesh::Point &p, const Real t, RealMatrixX &m) const
calculates the value of the function at the specified point, p, and time, t, and returns it in v...
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > prestress_A_matrix(MAST::ElementBase &e) const
This is the base class for elements that implement calculation of finite element quantities over the ...
Definition: elem_base.h:72