MAST
solid_2d_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 namespace MAST {
28  namespace Solid2DSectionProperty {
29 
31  public MAST::FieldFunction<RealMatrixX> {
32  public:
34  const MAST::FieldFunction<Real>& h);
35 
37 
38  virtual void operator() (const libMesh::Point& p,
39  const Real t,
40  RealMatrixX& m) const;
41 
42  virtual void derivative ( const MAST::FunctionBase& f,
43  const libMesh::Point& p,
44  const Real t,
45  RealMatrixX& m) const;
46 
47  protected:
48 
51  };
52 
53 
54 
56  public MAST::FieldFunction<RealMatrixX> {
57  public:
60  const MAST::FieldFunction<Real>& off);
61 
63 
64 
65  virtual void operator() (const libMesh::Point& p,
66  const Real t,
67  RealMatrixX& m) const;
68 
69  virtual void derivative ( const MAST::FunctionBase& f,
70  const libMesh::Point& p,
71  const Real t,
72  RealMatrixX& m) const;
73 
74  protected:
75 
78  };
79 
80 
81  class BendingStiffnessMatrix: public MAST::FieldFunction<RealMatrixX> {
82  public:
85  const MAST::FieldFunction<Real>& off);
86 
87 
89 
90  virtual void operator() (const libMesh::Point& p,
91  const Real t,
92  RealMatrixX& m) const;
93 
94  virtual void derivative ( const MAST::FunctionBase& f,
95  const libMesh::Point& p,
96  const Real t,
97  RealMatrixX& m) const;
98 
99  protected:
100 
103  };
104 
105 
106 
107 
108  class TransverseStiffnessMatrix: public MAST::FieldFunction<RealMatrixX> {
109  public:
111  const MAST::FieldFunction<Real>& h):
112  MAST::FieldFunction<RealMatrixX>("TransverseStiffnessMatrix2D"),
113  _material_stiffness(mat),
114  _h(h) {
115  _functions.insert(&mat);
116  _functions.insert(&h);
117  }
118 
119 
121 
122  virtual void operator() (const libMesh::Point& p,
123  const Real t,
124  RealMatrixX& m) const {
125  Real h;
126  _h(p, t, h);
127  _material_stiffness(p, t, m);
128  m *= h;
129  }
130 
131  virtual void derivative ( const MAST::FunctionBase& f,
132  const libMesh::Point& p,
133  const Real t,
134  RealMatrixX& m) const {
135  RealMatrixX dm;
136  Real h, dh;
137  _h(p, t, h); _h.derivative( f, p, t, dh);
138  _material_stiffness(p, t, m); _material_stiffness.derivative( f, p, t, dm);
139 
140  m *= dh;
141  m += h*dm;
142  }
143 
144  protected:
145 
148  };
149 
150 
151  class InertiaMatrix: public MAST::FieldFunction<RealMatrixX> {
152  public:
154  const MAST::FieldFunction<Real>& h,
155  const MAST::FieldFunction<Real>& off);
156 
157 
158  virtual ~InertiaMatrix() { }
159 
160 
161  virtual void operator() (const libMesh::Point& p,
162  const Real t,
163  RealMatrixX& m) const;
164 
165  virtual void derivative ( const MAST::FunctionBase& f,
166  const libMesh::Point& p,
167  const Real t,
168  RealMatrixX& m) const;
169 
170  protected:
171 
173  };
174 
175 
176 
177  class ThermalExpansionAMatrix: public MAST::FieldFunction<RealMatrixX> {
178  public:
180  const MAST::FieldFunction<RealMatrixX>& mat_expansion,
181  const MAST::FieldFunction<Real>& h);
182 
183 
185 
186 
187 
188  virtual void operator() (const libMesh::Point& p,
189  const Real t,
190  RealMatrixX& m) const;
191 
192 
193  virtual void derivative ( const MAST::FunctionBase& f,
194  const libMesh::Point& p,
195  const Real t,
196  RealMatrixX& m) const;
197 
198  protected:
199 
203  };
204 
205 
206 
207  class ThermalExpansionBMatrix: public MAST::FieldFunction<RealMatrixX> {
208  public:
210  const MAST::FieldFunction<RealMatrixX>& mat_expansion,
211  const MAST::FieldFunction<Real>& h,
212  const MAST::FieldFunction<Real>& off);
213 
214 
216 
217 
218  virtual void operator() (const libMesh::Point& p,
219  const Real t,
220  RealMatrixX& m) const;
221 
222 
223  virtual void derivative ( const MAST::FunctionBase& f,
224  const libMesh::Point& p,
225  const Real t,
226  RealMatrixX& m) const;
227 
228 
229  protected:
230 
234  };
235 
236 
237 
238 
240  public MAST::FieldFunction<RealMatrixX> {
241  public:
244  const MAST::FieldFunction<Real>& h);
245 
246 
247  virtual ~PrestressAMatrix() { }
248 
249  virtual void operator() (const libMesh::Point& p,
250  const Real t,
251  RealMatrixX& m) const;
252 
253  virtual void derivative ( const MAST::FunctionBase& f,
254  const libMesh::Point& p,
255  const Real t,
256  RealMatrixX& m) const;
257 
258  //virtual void convert_to_vector(const RealMatrixX& m, DenseRealVector& v) const;
259 
260  protected:
261 
264  };
265 
266 
267 
269  public MAST::FieldFunction<RealMatrixX> {
270  public:
273  const MAST::FieldFunction<Real>& h,
274  const MAST::FieldFunction<Real>& off);
275 
276 
277  virtual ~PrestressBMatrix() { }
278 
279  virtual void operator() (const libMesh::Point& p,
280  const Real t,
281  RealMatrixX& m) const;
282 
283  virtual void derivative ( const MAST::FunctionBase& f,
284  const libMesh::Point& p,
285  const Real t,
286  RealMatrixX& m) const;
287 
288  //virtual void convert_to_vector(const RealMatrixX& m, DenseRealVector& v) const;
289 
290  protected:
291 
294  };
295 
296 
298  public MAST::FieldFunction<RealMatrixX> {
299 
300  public:
301 
303  const MAST::FieldFunction<Real>& h);
304 
305 
306  virtual ~ThermalConductanceMatrix();
307 
308  virtual void operator() (const libMesh::Point& p,
309  const Real t,
310  RealMatrixX& m) const;
311 
312 
313 
314  virtual void derivative ( const MAST::FunctionBase& f,
315  const libMesh::Point& p,
316  const Real t,
317  RealMatrixX& m) const;
318 
319  protected:
320 
322 
324  };
325 
326 
327 
328 
330  public MAST::FieldFunction<RealMatrixX> {
331 
332  public:
333 
335  const MAST::FieldFunction<Real>& h);
336 
337 
338  virtual ~ThermalCapacitanceMatrix();
339 
340  virtual void operator() (const libMesh::Point& p,
341  const Real t,
342  RealMatrixX& m) const;
343 
344 
345 
346  virtual void derivative ( const MAST::FunctionBase& f,
347  const libMesh::Point& p,
348  const Real t,
349  RealMatrixX& m) const;
350 
351  protected:
352 
354 
356  };
357  }
358 }
359 
360 
361 
362 
363 bool
365 
366  return _material->depends_on(f) || // check if the material property depends on the function
367  MAST::ElementPropertyCardBase::depends_on(f); // check with this property card
368 }
369 
370 
371 
372 
375  const MAST::FieldFunction<Real>& h):
376 MAST::FieldFunction<RealMatrixX> ("ExtensionStiffnessMatrix2D"),
378 _h(h) {
379  _functions.insert(&mat);
380  _functions.insert(&h);
381 }
382 
383 
384 
385 void
388  const Real t,
389  RealMatrixX& m) const {
390  // [C]*h
391  Real h;
392  _h(p, t, h);
393  _material_stiffness(p, t, m);
394  m *= h;
395 }
396 
397 
398 
399 
400 void
403  const libMesh::Point& p,
404  const Real t,
405  RealMatrixX& m) const {
406  RealMatrixX dm;
407  Real h, dhdf;
408  _h(p, t, h); _h.derivative( f, p, t, dhdf);
409  _material_stiffness(p, t, m); _material_stiffness.derivative( f, p, t, dm);
410 
411  // [C]*dh
412  m *= dhdf;
413 
414  // += [dC]*h
415  m += h*dm;
416 }
417 
418 
419 
420 
421 
422 
425  const MAST::FieldFunction<Real>& h,
426  const MAST::FieldFunction<Real>& off):
427 MAST::FieldFunction<RealMatrixX> ("ExtensionBendingStiffnessMatrix2D"),
429 _h(h),
430 _off(off) {
431  _functions.insert(&mat);
432  _functions.insert(&h);
433  _functions.insert(&off);
434 }
435 
436 
437 
438 void
441  const Real t,
442  RealMatrixX& m) const {
443  // [C]*h
444  Real h, off;
445  _h(p, t, h);
446  _off(p, t, off);
447  _material_stiffness(p, t, m);
448  m *= h*off;
449 }
450 
451 
452 
453 
454 void
457  const libMesh::Point& p,
458  const Real t,
459  RealMatrixX& m) const {
460  RealMatrixX dm;
461  m = RealMatrixX::Zero(3,3); dm = RealMatrixX::Zero(3, 3);
462  Real h, off, dh, doff;
463 
464  _h(p, t, h); _h.derivative( f, p, t, dh);
465  _off(p, t, off); _off.derivative( f, p, t, doff);
466  _material_stiffness(p, t, m); _material_stiffness.derivative( f, p, t, dm);
467  m *= dh*off + h*doff;
468  m += h*off*dm;
469 }
470 
471 
472 
473 
476  const MAST::FieldFunction<Real>& h,
477  const MAST::FieldFunction<Real>& off):
478 MAST::FieldFunction<RealMatrixX> ("BendingStiffnessMatrix2D"),
480 _h(h),
481 _off(off) {
482  _functions.insert(&mat);
483  _functions.insert(&h);
484  _functions.insert(&off);
485 }
486 
487 
488 
489 void
491 BendingStiffnessMatrix::operator() (const libMesh::Point& p,
492  const Real t,
493  RealMatrixX& m) const {
494  // [C]*h
495  Real h, off;
496  _h(p, t, h);
497  _off(p, t, off);
498  _material_stiffness(p, t, m);
499  m *= (pow(h,3)/12. + h*pow(off,2));
500 }
501 
502 
503 
504 
505 void
508  const libMesh::Point& p,
509  const Real t,
510  RealMatrixX& m) const {
511  RealMatrixX dm;
512  m = RealMatrixX::Zero(3,3); dm = RealMatrixX::Zero(3, 3);
513  Real h, dhdf, off, doff;
514  _h(p, t, h); _h.derivative( f, p, t, dhdf);
515  _off(p, t, off); _off.derivative( f, p, t, doff);
516  _material_stiffness(p, t, m); _material_stiffness.derivative( f, p, t, dm);
517 
518  // [C]*dh
519  m *= (pow(h,2)/4.*dhdf + dhdf*pow(off,2) + h*2.*off*doff);
520 
521  // += [dC]*h
522  m += (pow(h,3)/12. + h*pow(off, 2))* dm;
523 }
524 
525 
526 
527 
528 
531  const MAST::FieldFunction<Real>& h,
532  const MAST::FieldFunction<Real>& off):
533 MAST::FieldFunction<RealMatrixX>("InertiaMatrix2D"),
534 _rho(rho),
535 _h(h),
536 _off(off) {
537  _functions.insert(&rho);
538  _functions.insert(&h);
539  _functions.insert(&off);
540 }
541 
542 
543 
544 
545 void
547 InertiaMatrix::operator() (const libMesh::Point& p,
548  const Real t,
549  RealMatrixX& m) const {
550  m = RealMatrixX::Zero(6, 6);
551  Real h, rho, off;
552  _h(p, t, h);
553  _off(p, t, off);
554  _rho(p, t, rho);
555 
556  for (unsigned int i=0; i<3; i++)
557  m(i,i) = h;
558 
559  m(0,4) = off*h; m(4,0) = m(0,4); // extension-bending coupling
560  m(1,3) = -off*h; m(3,1) = m(1,3); // extension-bending coupling
561  m(3,3) = pow(h,3)/12. + h*pow(off,2); // rotary inertia
562  m(4,4) = pow(h,3)/12. + h*pow(off,2); // rotary inertia
563  m(5,5) = pow(h,3)/12.*1.0e-6; // neglect the rotary inertia wrt theta_z
564 
565  m *= rho;
566 }
567 
568 
569 
570 
571 
572 void
575  const libMesh::Point& p,
576  const Real t,
577  RealMatrixX& m) const {
578  m = RealMatrixX::Zero(6,6);
579  Real h, dhdf, rho, drhodf, off, doff;
580  _h(p, t, h); _h.derivative( f, p, t, dhdf);
581  _off(p, t, off); _off.derivative( f, p, t, doff);
582  _rho(p, t, rho); _rho.derivative( f, p, t, drhodf);
583 
584  for (unsigned int i=0; i<3; i++)
585  m(i,i) = drhodf*h + rho*dhdf;
586 
587  m(0,4) = doff*h+off*dhdf; m(4,0) = m(0,4); // extension-bending coupling
588  m(1,3) = -doff*h-off*dhdf; m(3,1) = m(1,3); // extension-bending coupling
589  m(3,3) = drhodf*pow(h,3)/12.+rho*pow(h,2)/4.*dhdf; // rotary inertia
590  m(4,4) = drhodf*pow(h,3)/12.+rho*pow(h,2)/4.*dhdf; // rotary inertia
591  m(5,5) = (drhodf*pow(h,3)/12.+rho*pow(h,2)/4.*dhdf)*1.0e-6; // neglect the rotary inertia wrt theta_z
592 }
593 
594 
595 
596 
599  const MAST::FieldFunction<RealMatrixX>& mat_expansion,
600  const MAST::FieldFunction<Real>& h):
601 MAST::FieldFunction<RealMatrixX>("ThermalExpansionAMatrix2D"),
602 _material_stiffness(mat_stiff),
603 _material_expansion(mat_expansion),
604 _h(h) {
605  _functions.insert(&mat_stiff);
606  _functions.insert(&mat_expansion);
607  _functions.insert(&h);
608 }
609 
610 
611 
612 
613 void
615 ThermalExpansionAMatrix::operator() (const libMesh::Point& p,
616  const Real t,
617  RealMatrixX& m) const {
618  RealMatrixX at;
619  Real h;
620  _h(p, t, h);
621  _material_stiffness(p, t, m);
622  _material_expansion(p, t, at);
623 
624  m *= at;
625  m *= h;
626 }
627 
628 
629 
630 
631 
632 
633 void
636  const libMesh::Point& p,
637  const Real t,
638  RealMatrixX& m) const {
639  RealMatrixX m1, at, dm, dat;
640  Real h, dh;
641  _h(p, t, h); _h.derivative( f, p, t, dh);
642  _material_stiffness(p, t, m1); _material_stiffness.derivative( f, p, t, dm);
643  _material_expansion(p, t, at); _material_expansion.derivative( f, p, t, dat);
644 
645  m = m1;
646 
647  m *= at;
648  m *= dh;
649 
650  m1 *= dat;
651  dm *= at;
652  m1 += dm;
653 
654  m += h*m1;
655 }
656 
657 
658 
659 
662  const MAST::FieldFunction<RealMatrixX>& mat_expansion,
663  const MAST::FieldFunction<Real>& h,
664  const MAST::FieldFunction<Real>& off):
665 MAST::FieldFunction<RealMatrixX>("ThermalExpansionBMatrix2D"),
666 _material_stiffness(mat_stiff),
667 _material_expansion(mat_expansion),
668 _h(h),
669 _off(off) {
670  _functions.insert(&mat_stiff);
671  _functions.insert(&mat_expansion);
672  _functions.insert(&h);
673  _functions.insert(&off);
674 }
675 
676 
677 
678 
679 void
681 ThermalExpansionBMatrix::operator() (const libMesh::Point& p,
682  const Real t,
683  RealMatrixX& m) const {
684  RealMatrixX at;
685  Real h, off;
686  _h(p, t, h);
687  _off(p, t, off);
688  _material_stiffness(p, t, m);
689  _material_expansion(p, t, at);
690 
691  m *= at;
692  m *= h*off;
693 }
694 
695 
696 
697 
698 
699 void
702  const libMesh::Point& p,
703  const Real t,
704  RealMatrixX& m) const {
705  RealMatrixX m1, at, dm, dat;
706  Real h, dh, off, doff;
707  _h(p, t, h); _h.derivative( f, p, t, dh);
708  _off(p, t, off); _off.derivative( f, p, t, doff);
709  _material_stiffness(p, t, m1); _material_stiffness.derivative( f, p, t, dm);
710  _material_expansion(p, t, at); _material_expansion.derivative( f, p, t, dat);
711 
712  m = m1;
713 
714  m *= at;
715  m *= (dh*off+h*doff);
716 
717  m1 *= dat;
718  dm *= at;
719  m1 += dm;
720 
721  m += h*off*m1;
722 }
723 
724 
725 
726 
730  const MAST::FieldFunction<Real>& h):
731 MAST::FieldFunction<RealMatrixX>("PrestressAMatrix2D"),
732 _prestress(prestress),
733 _T(T),
734 _h(h) {
735  _functions.insert(&prestress);
736  _functions.insert(&T);
737  _functions.insert(&h);
738 }
739 
740 
741 
742 
743 void
745 PrestressAMatrix::operator() (const libMesh::Point& p,
746  const Real t,
747  RealMatrixX& m) const {
748  RealMatrixX s, T;
749  m = RealMatrixX::Zero(2, 2);
750  Real h;
751  _h(p, t, h);
752  _prestress(p, t, s);
753  _T(p, t, T);
754 
755  // convert the stress to the local coordinate
756  s *= T;
757  s = T.transpose() * s;
758 
759  for (unsigned int i=0; i<2; i++)
760  for (unsigned int j=0; j<2; j++)
761  m(i,j) = s(i,j)*h;
762 }
763 
764 
765 
766 
767 
768 
769 void
772  const libMesh::Point& p,
773  const Real t,
774  RealMatrixX& m) const {
775  RealMatrixX s, ds, T, dT;
776  m = RealMatrixX::Zero(2, 2);
777  Real h, dh;
778  _h(p, t, h); _h.derivative( f, p, t, dh);
779  _prestress(p, t, s); _prestress.derivative( f, p, t, ds);
780  _T(p, t, T); _T.derivative( f, p, t, dT);
781 
782  // convert the stress to the local coordinate
783  s *= T;
784  s = T.transpose() * s;
785 
786  // ds = dT^T s T + T^T s dT + T^T ds T
787  RealMatrixX tmp;
788  ds *= T;
789  ds = T.transpose()*ds;
790 
791  tmp = s;
792  tmp *= dT;
793  tmp = T.transpose() * tmp;
794  ds += tmp;
795 
796  tmp = s;
797  tmp *= T;
798  tmp = dT.transpose() * tmp;
799  ds += tmp;
800 
801 
802 
803  for (unsigned int i=0; i<2; i++)
804  for (unsigned int j=0; j<2; j++)
805  m(i,j) = ds(i,j)*h + s(i,j)*dh;
806 }
807 
808 
809 
810 
811 //void
812 //MAST::Solid2DSectionProperty::
813 //PrestressAMatrix::convert_to_vector(const RealMatrixX &m,
814 // RealVectorX &v) const {
815 // libmesh_assert_equal_to(m.rows(), 2);
816 // libmesh_assert_equal_to(m.cols(), 2);
817 // v.resize(3);
818 // v(0) = m(0,0); // sigma x
819 // v(2) = m(1,1); // sigma y
820 // v(2) = m(0,1); // tau xy
821 //}
822 
823 
824 
825 
829  const MAST::FieldFunction<Real>& h,
830  const MAST::FieldFunction<Real>& off):
831 MAST::FieldFunction<RealMatrixX>("PrestressBMatrix2D"),
832 _prestress(prestress),
833 _T(T),
834 _h(h),
835 _off(off) {
836  _functions.insert(&prestress);
837  _functions.insert(&T);
838  _functions.insert(&h);
839  _functions.insert(&off);
840 }
841 
842 
843 
844 
845 void
847 PrestressBMatrix::operator() (const libMesh::Point& p,
848  const Real t,
849  RealMatrixX& m) const {
850  RealMatrixX s, T;
851  m = RealMatrixX::Zero(2, 2);
852  Real h, off;
853  _h(p, t, h);
854  _off(p, t, off);
855  _prestress(p, t, s);
856  _T(p, t, T);
857 
858  // convert the stress to the local coordinate
859  s *= T;
860  s = T.transpose() * s;
861 
862  for (unsigned int i=0; i<2; i++)
863  for (unsigned int j=0; j<2; j++)
864  m(i,j) = s(i,j)*(h*off);
865 }
866 
867 
868 
869 
870 
871 void
874  const libMesh::Point& p,
875  const Real t,
876  RealMatrixX& m) const {
877  RealMatrixX s, ds, T, dT;
878  m = RealMatrixX::Zero(2, 2);
879  Real h, dh, off, doff;
880  _h(p, t, h); _h.derivative( f, p, t, dh);
881  _off(p, t, off); _off.derivative( f, p, t, doff);
882  _prestress(p, t, s); _prestress.derivative( f, p, t, ds);
883  _T(p, t, T); _T.derivative( f, p, t, dT);
884 
885  // convert the stress to the local coordinate
886  s *= T;
887  s = T.transpose() * s;
888 
889  // ds = dT^T s T + T^T s dT + T^T ds T
890  RealMatrixX tmp;
891  ds *= T;
892  ds = T.transpose() * ds;
893 
894  tmp = s;
895  tmp *= dT;
896  tmp = dT.transpose() * tmp;
897  ds += tmp;
898 
899  tmp = s;
900  tmp *= T;
901  tmp = T.transpose()*tmp;
902  ds += tmp;
903 
904 
905 
906  for (unsigned int i=0; i<2; i++)
907  for (unsigned int j=0; j<2; j++)
908  m(i,j) = ds(i,j)*(h*off) + s(i,j)*(dh*off+h*doff);
909 }
910 
911 
912 
913 
916  const MAST::FieldFunction<Real>& h):
917 MAST::FieldFunction<RealMatrixX>("ThermalConductanceMatrix"),
918 _mat_cond(mat_cond),
919 _h(h) {
920  _functions.insert(&mat_cond);
921  _functions.insert(&h);
922 }
923 
924 
925 
928 
929 
930 void
932 operator() (const libMesh::Point& p,
933  const Real t,
934  RealMatrixX& m) const {
935 
936  m = RealMatrixX::Zero(2, 2);
937  Real h;
938  _mat_cond(p, t, m);
939  _h(p, t, h);
940 
941  m *= h;
942 }
943 
944 
945 
946 
947 void
949  const libMesh::Point& p,
950  const Real t,
951  RealMatrixX& m) const {
952  m = RealMatrixX::Zero(2, 2);
953  RealMatrixX dm;
954  Real h, dh;
955  _mat_cond(p, t, m);
956  _mat_cond.derivative( f, p, t, dm);
957  _h(p, t, h);
958  _h.derivative( f, p, t, dh);
959 
960  m *= dh;
961  m += dm*h;
962 }
963 
964 
965 
966 
967 
970  const MAST::FieldFunction<Real>& h):
971 MAST::FieldFunction<RealMatrixX>("ThermalCapacitanceMatrix"),
972 _mat_cap(mat_cap),
973 _h(h) {
974  _functions.insert(&mat_cap);
975  _functions.insert(&h);
976 }
977 
978 
979 
982 
983 
984 void
986 operator() (const libMesh::Point& p,
987  const Real t,
988  RealMatrixX& m) const {
989 
990  m = RealMatrixX::Zero(1, 1);
991  Real h;
992  _mat_cap(p, t, m);
993  _h(p, t, h);
994 
995  m *= h;
996 }
997 
998 
999 
1000 
1001 void
1004  const libMesh::Point& p,
1005  const Real t,
1006  RealMatrixX& m) const {
1007 
1008  m = RealMatrixX::Zero(1, 1);
1009  RealMatrixX dm;
1010  Real h, dh;
1011  _mat_cap(p, t, m);
1012  _mat_cap.derivative( f, p, t, dm);
1013  _h(p, t, h);
1014  _h.derivative( f, p, t, dh);
1015 
1016  m *= dh;
1017  m += dm*h;
1018 }
1019 
1020 
1021 
1022 //void
1023 //MAST::Solid2DSectionProperty::
1024 //PrestressBMatrix::convert_to_vector(const RealMatrixX &m,
1025 // RealVectorX &v) const {
1026 // // nothing to be done for a symmetric section
1027 // v.resize(3);
1028 //}
1029 
1030 
1031 
1032 
1033 std::unique_ptr<MAST::FieldFunction<RealMatrixX> >
1036 
1039  (_material->stiffness_matrix(2),
1040  this->get<const FieldFunction<Real> >("h"));
1041 
1042  return std::unique_ptr<MAST::FieldFunction<RealMatrixX> > (rval);
1043 }
1044 
1045 
1046 std::unique_ptr<MAST::FieldFunction<RealMatrixX> >
1049 
1052  (_material->stiffness_matrix(2),
1053  this->get<FieldFunction<Real> >("h"),
1054  this->get<FieldFunction<Real> >("off"));
1055 
1056  return std::unique_ptr<MAST::FieldFunction<RealMatrixX> > (rval);
1057 }
1058 
1059 
1060 std::unique_ptr<MAST::FieldFunction<RealMatrixX> >
1063 
1064 
1067  (_material->stiffness_matrix(2),
1068  this->get<FieldFunction<Real> >("h"),
1069  this->get<FieldFunction<Real> >("off"));
1070 
1071  return std::unique_ptr<MAST::FieldFunction<RealMatrixX> > (rval);
1072 }
1073 
1074 
1075 std::unique_ptr<MAST::FieldFunction<RealMatrixX> >
1078 
1079  libmesh_error();
1080 
1081  return std::unique_ptr<MAST::FieldFunction<RealMatrixX> > (nullptr);
1082 }
1083 
1084 
1085 
1086 std::unique_ptr<MAST::FieldFunction<RealMatrixX> >
1089 
1090 
1093  (_material->get<FieldFunction<Real> >("rho"),
1094  this->get<FieldFunction<Real> >("h"),
1095  this->get<FieldFunction<Real> >("off"));
1096 
1097  return std::unique_ptr<MAST::FieldFunction<RealMatrixX> > (rval);
1098 }
1099 
1100 
1101 
1102 std::unique_ptr<MAST::FieldFunction<RealMatrixX> >
1105 
1108  (_material->stiffness_matrix(2),
1109  _material->thermal_expansion_matrix(2),
1110  this->get<FieldFunction<Real> >("h"));
1111 
1112  return std::unique_ptr<MAST::FieldFunction<RealMatrixX> > (rval);
1113 }
1114 
1115 
1116 
1117 std::unique_ptr<MAST::FieldFunction<RealMatrixX> >
1120 
1121 
1124  (_material->stiffness_matrix(2),
1125  _material->thermal_expansion_matrix(2),
1126  this->get<FieldFunction<Real> >("h"),
1127  this->get<FieldFunction<Real> >("off"));
1128 
1129  return std::unique_ptr<MAST::FieldFunction<RealMatrixX> > (rval);
1130 }
1131 
1132 
1133 std::unique_ptr<MAST::FieldFunction<RealMatrixX> >
1136 
1137 
1140  (_material->transverse_shear_stiffness_matrix(),
1141  this->get<FieldFunction<Real> >("h"));
1142 
1143  return std::unique_ptr<MAST::FieldFunction<RealMatrixX> > (rval);
1144 }
1145 
1146 
1147 std::unique_ptr<MAST::FieldFunction<RealMatrixX> >
1150 
1152  // TODO: figure out the interface for prestress and T matrix
1153  libmesh_assert(false);
1154  // = new MAST::Solid2DSectionProperty::PrestressAMatrix
1155  // (this->get<MAST::FieldFunction<RealMatrixX> >("prestress"),
1156  // e.local_elem().T_matrix(),
1157  // this->get<FieldFunction<Real> >("h"));
1158 
1159  return std::unique_ptr<MAST::FieldFunction<RealMatrixX> > (rval);
1160 }
1161 
1162 
1163 std::unique_ptr<MAST::FieldFunction<RealMatrixX> >
1166 
1168  // TODO: figure out the interface for prestress and T matrix
1169  libmesh_assert(false);
1170  // = new MAST::Solid2DSectionProperty::PrestressBMatrix
1171  // (this->get<MAST::FieldFunction<RealMatrixX> >("prestress"),
1172  // e.local_elem().T_matrix(),
1173  // this->get<FieldFunction<Real> >("h"),
1174  // this->get<FieldFunction<Real> >("off"));
1175 
1176  return std::unique_ptr<MAST::FieldFunction<RealMatrixX> > (rval);
1177 }
1178 
1179 
1180 
1181 std::unique_ptr<MAST::FieldFunction<RealMatrixX> >
1184 
1185 
1188  (_material->conductance_matrix(2),
1189  this->get<FieldFunction<Real> >("h"));
1190 
1191  return std::unique_ptr<MAST::FieldFunction<RealMatrixX> > (rval);
1192 }
1193 
1194 
1195 
1196 std::unique_ptr<MAST::FieldFunction<RealMatrixX> >
1199 
1200 
1203  (_material->capacitance_matrix(2),
1204  this->get<FieldFunction<Real> >("h"));
1205 
1206  return std::unique_ptr<MAST::FieldFunction<RealMatrixX> > (rval);
1207 }
1208 
1209 
1212 section(const MAST::ElementBase& e) const {
1213 
1214  return this->get<FieldFunction<Real>>("h");
1215 }
1216 
ExtensionBendingStiffnessMatrix(const MAST::FieldFunction< RealMatrixX > &mat, const MAST::FieldFunction< Real > &h, const MAST::FieldFunction< Real > &off)
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > prestress_A_matrix(MAST::ElementBase &e) const
BendingStiffnessMatrix(const MAST::FieldFunction< RealMatrixX > &mat, const MAST::FieldFunction< Real > &h, const MAST::FieldFunction< Real > &off)
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 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...
ExtensionStiffnessMatrix(const MAST::FieldFunction< RealMatrixX > &mat, const MAST::FieldFunction< Real > &h)
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...
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< Real > & section(const MAST::ElementBase &e) const
TransverseStiffnessMatrix(const MAST::FieldFunction< RealMatrixX > &mat, const MAST::FieldFunction< Real > &h)
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 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 > > inertia_matrix(const MAST::ElementBase &e) const
std::set< const MAST::FunctionBase * > _functions
set of functions that this function depends on
libMesh::Real Real
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...
PrestressBMatrix(const MAST::FieldFunction< RealMatrixX > &prestress, const MAST::FieldFunction< RealMatrixX > &T, const MAST::FieldFunction< Real > &h, const MAST::FieldFunction< Real > &off)
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > stiffness_D_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
virtual void derivative(const MAST::FunctionBase &f, ValType &v) const
calculates the value of the function derivative and returns it in v.
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > transverse_shear_stiffness_matrix(const MAST::ElementBase &e) const
Matrix< Real, Dynamic, Dynamic > RealMatrixX
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > thermal_conductance_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...
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.
PrestressAMatrix(const MAST::FieldFunction< RealMatrixX > &prestress, const MAST::FieldFunction< RealMatrixX > &T, const MAST::FieldFunction< Real > &h)
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > thermal_capacitance_matrix(const MAST::ElementBase &e) const
ThermalConductanceMatrix(const MAST::FieldFunction< RealMatrixX > &mat_cond, const MAST::FieldFunction< Real > &h)
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_expansion_B_matrix(const MAST::ElementBase &e) const
ThermalExpansionAMatrix(const MAST::FieldFunction< RealMatrixX > &mat_stiff, const MAST::FieldFunction< RealMatrixX > &mat_expansion, const MAST::FieldFunction< Real > &h)
ThermalCapacitanceMatrix(const MAST::FieldFunction< RealMatrixX > &mat_cond, const MAST::FieldFunction< Real > &h)
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_B_matrix(const MAST::ElementBase &e) const
InertiaMatrix(const MAST::FieldFunction< Real > &rho, const MAST::FieldFunction< Real > &h, const MAST::FieldFunction< Real > &off)
ThermalExpansionBMatrix(const MAST::FieldFunction< RealMatrixX > &mat_stiff, const MAST::FieldFunction< RealMatrixX > &mat_expansion, const MAST::FieldFunction< Real > &h, const MAST::FieldFunction< Real > &off)
virtual bool depends_on(const MAST::FunctionBase &f) const
returns true if the property card depends on the function f
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > damping_matrix(const MAST::ElementBase &e) const
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > stiffness_A_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...
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 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 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...
This is the base class for elements that implement calculation of finite element quantities over the ...
Definition: elem_base.h:72