MAST
isotropic_material_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 
21 // MAST includes
24 
25 
26 namespace MAST {
27  namespace IsotropicMaterialProperty {
28 
29 
31  public MAST::FieldFunction<RealMatrixX> {
32  public:
33 
35  const MAST::FieldFunction<Real>& nu);
36 
37  virtual ~StiffnessMatrix1D() { }
38 
39  virtual void operator() (const libMesh::Point& p,
40  const Real t,
41  RealMatrixX& m) const;
42 
43  virtual void derivative( const MAST::FunctionBase& f,
44  const libMesh::Point& p,
45  const Real t,
46  RealMatrixX& m) const;
47 
48  protected:
49 
52  };
53 
54 
55 
57  public:
59  const MAST::FieldFunction<Real>& nu,
60  const MAST::FieldFunction<Real>& kappa);
61 
62 
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 
79  };
80 
81 
82  class StiffnessMatrix2D: public MAST::FieldFunction<RealMatrixX> {
83  public:
85  const MAST::FieldFunction<Real>& nu,
86  bool plane_stress);
87 
88 
89  virtual ~StiffnessMatrix2D() { }
90 
91 
92  virtual void operator() (const libMesh::Point& p,
93  const Real t,
94  RealMatrixX& m) const;
95 
96  virtual void derivative( const MAST::FunctionBase& f,
97  const libMesh::Point& p,
98  const Real t,
99  RealMatrixX& m) const;
100 
101  protected:
102 
106  };
107 
108 
109 
110  class StiffnessMatrix3D: public MAST::FieldFunction<RealMatrixX> {
111  public:
113  const MAST::FieldFunction<Real>& nu);
114 
115  virtual ~StiffnessMatrix3D() { }
116 
117  virtual void operator() (const libMesh::Point& p,
118  const Real t,
119  RealMatrixX& m) const;
120 
121  virtual void derivative( const MAST::FunctionBase& f,
122  const libMesh::Point& p,
123  const Real t,
124  RealMatrixX& m) const;
125 
126  protected:
127 
130  };
131 
132 
133 
134  class InertiaMatrix3D: public MAST::FieldFunction<RealMatrixX> {
135  public:
137 
138 
139  virtual ~InertiaMatrix3D() { }
140 
141  virtual void operator() (const libMesh::Point& p,
142  const Real t,
143  RealMatrixX& m) const;
144 
145  virtual void derivative( const MAST::FunctionBase& f,
146  const libMesh::Point& p,
147  const Real t,
148  RealMatrixX& m) const;
149 
150  protected:
151 
153 
154  };
155 
156 
157 
158  class ThermalExpansionMatrix: public MAST::FieldFunction<RealMatrixX> {
159  public:
160 
161  ThermalExpansionMatrix(unsigned int dim,
162  const MAST::FieldFunction<Real>& alpha):
163  MAST::FieldFunction<RealMatrixX>("ThermalExpansionMatrix"),
164  _dim(dim),
165  _alpha(alpha) {
166  _functions.insert(&_alpha);
167  }
168 
169 
171 
172  virtual void operator() (const libMesh::Point& p,
173  const Real t,
174  RealMatrixX& m) const;
175 
176  virtual void derivative( const MAST::FunctionBase& f,
177  const libMesh::Point& p,
178  const Real t,
179  RealMatrixX& m) const;
180 
181  protected:
182 
183  const unsigned int _dim;
184 
186  };
187 
188 
189 
190 
192  public MAST::FieldFunction<RealMatrixX> {
193  public:
194 
195  ThermalConductanceMatrix(unsigned int dim,
197  MAST::FieldFunction<RealMatrixX>("ThermalConductanceMatrix"),
198  _dim(dim),
199  _k(k) {
200  _functions.insert(&_k);
201  }
202 
203 
205 
206  virtual void operator() (const libMesh::Point& p,
207  const Real t,
208  RealMatrixX& m) const;
209 
210  virtual void derivative( const MAST::FunctionBase& f,
211  const libMesh::Point& p,
212  const Real t,
213  RealMatrixX& m) const;
214 
215  protected:
216 
217  const unsigned int _dim;
218 
220  };
221 
222 
223 
225  public MAST::FieldFunction<RealMatrixX> {
226  public:
227 
228  ThermalCapacitanceMatrix(unsigned int dim,
229  const MAST::FieldFunction<Real>& rho,
230  const MAST::FieldFunction<Real>& cp):
231  MAST::FieldFunction<RealMatrixX>("ThermalCapacitanceMatrix"),
232  _dim(dim),
233  _rho(rho),
234  _cp(cp) {
235 
236  _functions.insert(&_rho);
237  _functions.insert(&_cp);
238  }
239 
240 
242 
243  virtual void operator() (const libMesh::Point& p,
244  const Real t,
245  RealMatrixX& m) const;
246 
247  virtual void derivative( const MAST::FunctionBase& f,
248  const libMesh::Point& p,
249  const Real t,
250  RealMatrixX& m) const;
251 
252  protected:
253 
254  const unsigned int _dim;
255 
257 
259  };
260 
261 
262  }
263 }
264 
265 
266 
269  const MAST::FieldFunction<Real>& nu ):
270 MAST::FieldFunction<RealMatrixX>("StiffnessMatrix1D"),
271 _E(E),
272 _nu(nu)
273 {
274  _functions.insert(&E);
275  _functions.insert(&nu);
276 }
277 
278 
279 
280 
281 void
283 StiffnessMatrix1D::operator() (const libMesh::Point& p,
284  const Real t,
285  RealMatrixX& m) const {
286  m = RealMatrixX::Zero(2,2);
287  Real E, nu, G;
288  _E(p, t, E); _nu(p, t, nu);
289  G = E/2./(1.+nu);
290  m(0,0) = E;
291  m(1,1) = G;
292 }
293 
294 
295 void
298  const libMesh::Point &p,
299  const Real t,
300  RealMatrixX &m) const {
301 
302 
303  RealMatrixX dm;
304  m = RealMatrixX::Zero(2,2); dm = RealMatrixX::Zero(2,2);
305  Real E, nu, dEdf, dnudf;
306  _E(p, t, E); _E.derivative( f, p, t, dEdf);
307  _nu(p, t, nu); _nu.derivative( f, p, t, dnudf);
308 
309  // parM/parE * parE/parf
310  dm(0,0) = 1.;
311  dm(1,1) = 1./2./(1.+nu);
312  m += dEdf * dm;
313 
314 
315  // parM/parnu * parnu/parf
316  dm(0,0) = 0.;
317  dm(1,1) = -E/2./pow(1.+nu,2);
318  m+= dnudf*dm;
319 }
320 
321 
325  const MAST::FieldFunction<Real>& nu,
326  const MAST::FieldFunction<Real>& kappa):
327 MAST::FieldFunction<RealMatrixX>("TransverseShearStiffnessMatrix"),
328 _E(E),
329 _nu(nu),
330 _kappa(kappa)
331 {
332  _functions.insert(&E);
333  _functions.insert(&nu);
334  _functions.insert(&kappa);
335 }
336 
337 
338 
339 void
342  const Real t,
343  RealMatrixX& m) const {
344  m = RealMatrixX::Zero(2,2);
345  Real E, nu, kappa, G;
346  _E(p, t, E); _nu(p, t, nu); _kappa(p, t, kappa);
347  G = E/2./(1.+nu);
348  m(0,0) = G*kappa;
349  m(1,1) = m(0,0);
350 }
351 
352 
353 
354 void
357  const libMesh::Point& p,
358  const Real t,
359  RealMatrixX& m) const {
360  RealMatrixX dm;
361  m = RealMatrixX::Zero(2,2); dm = RealMatrixX::Zero(2, 2);
362  Real E, nu, kappa, dEdf, dnudf, dkappadf, G;
363  _E (p, t, E); _E.derivative( f, p, t, dEdf);
364  _nu (p, t, nu); _nu.derivative( f, p, t, dnudf);
365  _kappa(p, t, kappa); _kappa.derivative( f, p, t, dkappadf);
366  G = E/2./(1.+nu);
367 
368 
369  // parM/parE * parE/parf
370  dm(0,0) = 1./2./(1.+nu)*kappa;
371  dm(1,1) = dm(0,0);
372  m += dEdf * dm;
373 
374 
375  // parM/parnu * parnu/parf
376  dm(0,0) = -E/2./pow(1.+nu,2)*kappa;
377  dm(1,1) = dm(0,0);
378  m += dnudf*dm;
379 
380  // parM/parnu * parkappa/parf
381 
382  dm(0,0) = G; dm(1,1) = G;
383  m += dkappadf*dm;
384 }
385 
386 
387 
388 
391  const MAST::FieldFunction<Real>& nu ,
392  bool plane_stress ):
393 MAST::FieldFunction<RealMatrixX>("StiffnessMatrix2D"),
394 _E(E),
395 _nu(nu),
396 _plane_stress(plane_stress) {
397 
398  _functions.insert(&E);
399  _functions.insert(&nu);
400 }
401 
402 
403 
404 
405 void
407 StiffnessMatrix2D::operator() (const libMesh::Point& p,
408  const Real t,
409  RealMatrixX& m) const {
410  libmesh_assert(_plane_stress); // currently only implemented for plane stress
411  m = RealMatrixX::Zero(3,3);
412  Real E, nu;
413  _E(p, t, E); _nu(p, t, nu);
414  for (unsigned int i=0; i<2; i++) {
415  for (unsigned int j=0; j<2; j++)
416  if (i == j) // diagonal: direct stress
417  m(i,i) = E/(1.-nu*nu);
418  else // offdiagonal: direct stress
419  m(i,j) = E*nu/(1.-nu*nu);
420  }
421  m(2,2) = E/2./(1.+nu); // diagonal: shear stress
422 }
423 
424 
425 
426 
427 void
430  const libMesh::Point& p,
431  const Real t,
432  RealMatrixX& m) const {
433  libmesh_assert(_plane_stress); // currently only implemented for plane stress
434  RealMatrixX dm;
435  m = RealMatrixX::Zero(3,3); dm = RealMatrixX::Zero(3, 3);
436  Real E, nu, dEdf, dnudf;
437  _E (p, t, E); _E.derivative( f, p, t, dEdf);
438  _nu (p, t, nu); _nu.derivative( f, p, t, dnudf);
439 
440  // parM/parE * parE/parf
441  for (unsigned int i=0; i<2; i++) {
442  for (unsigned int j=0; j<2; j++)
443  if (i == j) // diagonal: direct stress
444  dm(i,i) = 1./(1.-nu*nu);
445  else // offdiagonal: direct stress
446  dm(i,j) = 1.*nu/(1.-nu*nu);
447  }
448  dm(2,2) = 1./2./(1.+nu); // diagonal: shear stress
449  m += dEdf * dm;
450 
451  // parM/parnu * parnu/parf
452  for (unsigned int i=0; i<2; i++) {
453  for (unsigned int j=0; j<2; j++)
454  if (i == j) // diagonal: direct stress
455  dm(i,i) = E/pow(1.-nu*nu, 2)*2.*nu;
456  else // offdiagonal: direct stress
457  dm(i,j) = E/(1.-nu*nu) + E*nu/pow(1.-nu*nu,2)*2.*nu;
458  }
459  dm(2,2) = -E/2./pow(1.+nu,2); // diagonal: shear stress
460  m+= dnudf*dm;
461 }
462 
463 
464 
467  const MAST::FieldFunction<Real>& nu):
468 MAST::FieldFunction<RealMatrixX>("StiffnessMatrix3D"),
469 _E(E),
470 _nu(nu) {
471 
472  _functions.insert(&E);
473  _functions.insert(&nu);
474 }
475 
476 
477 
478 
479 
480 void
482 StiffnessMatrix3D::operator() (const libMesh::Point& p,
483  const Real t,
484  RealMatrixX& m) const {
485  m = RealMatrixX::Zero(6,6);
486  Real E, nu;
487  _E(p, t, E); _nu(p, t, nu);
488  for (unsigned int i=0; i<3; i++) {
489  for (unsigned int j=0; j<3; j++)
490  if (i == j) // diagonal: direct stress
491  m(i,i) = E*(1.-nu)/(1.-nu-2.*nu*nu);
492  else // offdiagonal: direct stress
493  m(i,j) = E*nu/(1.-nu-2.*nu*nu);
494  m(i+3,i+3) = E/2./(1.+nu); // diagonal: shear stress
495  }
496 }
497 
498 
499 
500 
501 void
504  const libMesh::Point& p,
505  const Real t,
506  RealMatrixX& m) const {
507  RealMatrixX dm;
508  m = RealMatrixX::Zero(6,6); dm = RealMatrixX::Zero(6,6);
509  Real E, nu, dEdf, dnudf;
510  _E (p, t, E); _E.derivative( f, p, t, dEdf);
511  _nu (p, t, nu); _nu.derivative( f, p, t, dnudf);
512 
513  // parM/parE * parE/parf
514  for (unsigned int i=0; i<3; i++) {
515  for (unsigned int j=0; j<3; j++)
516  if (i == j) // diagonal: direct stress
517  dm(i,i) = (1.-nu)/(1.-nu-2.*nu*nu);
518  else // offdiagonal: direct stress
519  dm(i,j) = nu/(1.-nu-2.*nu*nu);
520  dm(i+3,i+3) = 1./2./(1.+nu); // diagonal: shear stress
521  }
522  m += dEdf * dm;
523 
524 
525  // parM/parnu * parnu/parf
526  for (unsigned int i=0; i<3; i++) {
527  for (unsigned int j=0; j<3; j++)
528  if (i == j) // diagonal: direct stress
529  dm(i,i) = -E/(1.-nu-2.*nu*nu) + E*(1.-nu)/pow(1.-nu-2.*nu*nu,2)*(1.+4.*nu);
530  else // offdiagonal: direct stress
531  dm(i,j) = E/(1.-nu-2.*nu*nu) + E*nu/pow(1.-nu-2.*nu*nu,2)*(1.+4.*nu);
532  dm(i+3,i+3) = -E/2./pow(1.+nu,2); // diagonal: shear stress
533  }
534  m+= dnudf*dm;
535 }
536 
537 
538 
541 MAST::FieldFunction<RealMatrixX>("InertiaMatrix3D"),
542 _rho(rho) {
543 
544  _functions.insert(&rho);
545 }
546 
547 
548 void
550 InertiaMatrix3D::operator() (const libMesh::Point& p,
551  const Real t,
552  RealMatrixX& m) const {
553  m = RealMatrixX::Zero(3,3);
554  Real rho;
555  _rho(p, t, rho);
556  m(0,0) = rho;
557  m(1,1) = rho;
558  m(2,2) = rho;
559 }
560 
561 
562 void
565  const MAST::FunctionBase &f,
566  const libMesh::Point &p,
567  const Real t,
568  RealMatrixX &m) const {
569 
570 
571  m = RealMatrixX::Zero(3,3);
572  Real rho;
573  _rho.derivative( f, p, t, rho);
574  m(0,0) = rho;
575  m(1,1) = rho;
576  m(2,2) = rho;
577 }
578 
579 
580 
581 void
583 operator() (const libMesh::Point& p,
584  const Real t,
585  RealMatrixX& m) const {
586 
587  Real alpha;
588  _alpha(p, t, alpha);
589  switch (_dim) {
590  case 1:
591  m = RealMatrixX::Zero(2,1);
592  break;
593 
594  case 2:
595  m = RealMatrixX::Zero(3,1);
596  break;
597 
598  case 3:
599  m = RealMatrixX::Zero(6,1);
600  break;
601  }
602 
603  for (unsigned int i=0; i<_dim; i++)
604  m(i,0) = alpha;
605 }
606 
607 
608 
609 
610 
611 void
614  const libMesh::Point& p,
615  const Real t,
616  RealMatrixX& m) const {
617 
618 
619  Real alpha;
620  _alpha.derivative( f, p, t, alpha);
621  switch (_dim) {
622  case 1:
623  m = RealMatrixX::Zero(2,1);
624  break;
625 
626  case 2:
627  m = RealMatrixX::Zero(3,1);
628  break;
629 
630  case 3:
631  m = RealMatrixX::Zero(6,1);
632  break;
633  }
634 
635  for (unsigned int i=0; i<_dim; i++)
636  m(i,0) = alpha;
637 
638 }
639 
640 
641 
642 
643 
644 void
646 operator() (const libMesh::Point& p,
647  const Real t,
648  RealMatrixX& m) const {
649 
650  Real cp, rho;
651  _cp (p, t, cp);
652  _rho (p, t, rho);
653 
654  m.setZero(1,1);
655 
656  m(0,0) = cp*rho;
657 }
658 
659 
660 
661 
662 
663 void
666  const libMesh::Point& p,
667  const Real t,
668  RealMatrixX& m) const {
669 
670 
671  Real cp, dcp, rho, drho;
672  _cp (p, t, cp); _cp.derivative( f, p, t, dcp);
673  _rho (p, t, rho); _rho.derivative( f, p, t, drho);
674 
675  m.setZero(1,1);
676 
677  m(0,0) = dcp*rho + cp*drho;
678 }
679 
680 
681 
682 
683 void
685 operator() (const libMesh::Point& p,
686  const Real t,
687  RealMatrixX& m) const {
688 
689  Real k;
690  _k (p, t, k);
691 
692  m.setIdentity(_dim, _dim);
693 
694  m *= k;
695 }
696 
697 
698 
699 
700 
701 void
704  const libMesh::Point& p,
705  const Real t,
706  RealMatrixX& m) const {
707 
708 
709  Real k;
710  _k.derivative( f, p, t, k);
711 
712  m.setIdentity(_dim, _dim);
713 
714  m *= k;
715 }
716 
717 
718 
719 
720 
723 _stiff_mat_1d (nullptr),
724 _stiff_mat_2d (nullptr),
725 _stiff_mat_3d (nullptr),
726 _damp_mat (nullptr),
727 _inertia_mat_3d (nullptr),
728 _thermal_exp_mat_1d (nullptr),
729 _thermal_exp_mat_2d (nullptr),
730 _thermal_exp_mat_3d (nullptr),
731 _transverse_shear_mat (nullptr),
732 _thermal_capacitance_mat_1d (nullptr),
733 _thermal_capacitance_mat_2d (nullptr),
734 _thermal_capacitance_mat_3d (nullptr),
735 _thermal_conductance_mat_1d (nullptr),
736 _thermal_conductance_mat_2d (nullptr),
737 _thermal_conductance_mat_3d (nullptr)
738 { }
739 
740 
741 
743 
744  if (_stiff_mat_1d) delete _stiff_mat_1d;
745  if (_stiff_mat_2d) delete _stiff_mat_2d;
746  if (_stiff_mat_3d) delete _stiff_mat_3d;
747 
748  if (_damp_mat) delete _damp_mat;
749  if (_inertia_mat_3d) delete _inertia_mat_3d;
754 
758 
762 }
763 
764 
767  const bool plane_stress) {
768 
769  MAST::FieldFunction<RealMatrixX> *rval = nullptr;
770 
771  switch (dim) {
772  case 1: {
773 
774  if (!_stiff_mat_1d)
776  (this->get<MAST::FieldFunction<Real> >("E"),
777  this->get<MAST::FieldFunction<Real> >("nu"));
778  return *_stiff_mat_1d;
779  }
780  break;
781 
782  case 2: {
783 
784  if (!_stiff_mat_2d)
786  (this->get<MAST::FieldFunction<Real> >("E"),
787  this->get<MAST::FieldFunction<Real> >("nu"),
788  plane_stress);
789 
790  return *_stiff_mat_2d;
791  }
792  break;
793 
794  case 3: {
795 
796  if (!_stiff_mat_3d)
798  (this->get<MAST::FieldFunction<Real> >("E"),
799  this->get<MAST::FieldFunction<Real> >("nu"));
800 
801  return *_stiff_mat_3d;
802  }
803  break;
804 
805  default:
806  // should not get here
807  libmesh_error();
808  }
809 }
810 
811 
812 
813 
816 
817  // make sure that this is not null
818  libmesh_assert(false);
819 
820  return *_damp_mat;
821 }
822 
823 
824 
825 
828 
829  switch (dim) {
830  case 3: {
831 
832  if (!_inertia_mat_3d)
834  (this->get<MAST::FieldFunction<Real> >("rho"));
835 
836  return *_inertia_mat_3d;
837  }
838  break;
839 
840  case 1:
841  case 2:
842  default:
843  // implemented only for 3D since the 2D and 1D elements
844  // calculate it themselves
845  libmesh_error();
846 
847  }
848 }
849 
850 
851 
852 
855 
856  switch (dim) {
857  case 3: {
858 
859  if (!_thermal_exp_mat_3d)
862  (dim,
863  this->get<MAST::FieldFunction<Real> >("alpha_expansion"));
864 
865  return *_thermal_exp_mat_3d;
866  }
867  break;
868 
869 
870  case 2: {
871 
872  if (!_thermal_exp_mat_2d)
875  (dim,
876  this->get<MAST::FieldFunction<Real> >("alpha_expansion"));
877 
878  return *_thermal_exp_mat_2d;
879  }
880  break;
881 
882 
883  case 1: {
884 
885  if (!_thermal_exp_mat_1d)
888  (dim,
889  this->get<MAST::FieldFunction<Real> >("alpha_expansion"));
890 
891  return *_thermal_exp_mat_1d;
892  }
893  break;
894 
895 
896  default:
897  libmesh_error();
898  }
899 
900 
901 }
902 
903 
904 
905 
908 
912  (this->get<MAST::FieldFunction<Real> >("E"),
913  this->get<MAST::FieldFunction<Real> >("nu"),
914  this->get<MAST::FieldFunction<Real> >("kappa"));
915 
916  return *_transverse_shear_mat;
917 }
918 
919 
920 
923 
924  switch (dim) {
925 
926  case 1: {
927 
931  (dim,
932  this->get<MAST::FieldFunction<Real> >("rho"),
933  this->get<MAST::FieldFunction<Real> >("cp"));
934 
936  }
937  break;
938 
939 
940  case 2: {
941 
945  (dim,
946  this->get<MAST::FieldFunction<Real> >("rho"),
947  this->get<MAST::FieldFunction<Real> >("cp"));
948 
950  }
951  break;
952 
953 
954  case 3: {
955 
959  (dim,
960  this->get<MAST::FieldFunction<Real> >("rho"),
961  this->get<MAST::FieldFunction<Real> >("cp"));
962 
964  }
965  break;
966 
967 
968  default:
969  // should not get here
970  libmesh_error();
971  }
972 }
973 
974 
975 
976 
977 
980 
981  switch (dim) {
982 
983  case 1: {
984 
988  (dim, this->get<MAST::FieldFunction<Real> >("k_th"));
989 
991  }
992  break;
993 
994 
995  case 2: {
996 
1000  (dim, this->get<MAST::FieldFunction<Real> >("k_th"));
1001 
1003  }
1004  break;
1005 
1006 
1007  case 3: {
1008 
1012  (dim, this->get<MAST::FieldFunction<Real> >("k_th"));
1013 
1015  }
1016  break;
1017 
1018 
1019  default:
1020  // should not get here
1021  libmesh_error();
1022  }
1023 }
1024 
1025 
1026 
1027 
MAST::FieldFunction< RealMatrixX > * _thermal_conductance_mat_3d
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 > & inertia_matrix(const unsigned int dim)
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 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...
StiffnessMatrix3D(const MAST::FieldFunction< Real > &E, const MAST::FieldFunction< Real > &nu)
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...
MAST::FieldFunction< RealMatrixX > * _thermal_exp_mat_3d
virtual const MAST::FieldFunction< RealMatrixX > & thermal_expansion_matrix(const unsigned int dim)
ThermalCapacitanceMatrix(unsigned int dim, const MAST::FieldFunction< Real > &rho, const MAST::FieldFunction< Real > &cp)
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...
std::set< const MAST::FunctionBase * > _functions
set of functions that this function depends on
StiffnessMatrix1D(const MAST::FieldFunction< Real > &E, const MAST::FieldFunction< Real > &nu)
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...
libMesh::Real Real
TransverseShearStiffnessMatrix(const MAST::FieldFunction< Real > &E, const MAST::FieldFunction< Real > &nu, const MAST::FieldFunction< Real > &kappa)
ThermalConductanceMatrix(unsigned int dim, MAST::FieldFunction< Real > &k)
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...
MAST::FieldFunction< RealMatrixX > * _inertia_mat_3d
virtual void derivative(const MAST::FunctionBase &f, ValType &v) const
calculates the value of the function derivative and returns it in v.
MAST::FieldFunction< RealMatrixX > * _damp_mat
Matrix< Real, Dynamic, Dynamic > RealMatrixX
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...
MAST::FieldFunction< RealMatrixX > * _transverse_shear_mat
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.
StiffnessMatrix2D(const MAST::FieldFunction< Real > &E, const MAST::FieldFunction< Real > &nu, bool plane_stress)
virtual const MAST::FieldFunction< RealMatrixX > & damping_matrix(const unsigned int dim)
virtual const MAST::FieldFunction< RealMatrixX > & conductance_matrix(const unsigned int dim)
virtual const MAST::FieldFunction< RealMatrixX > & capacitance_matrix(const unsigned int dim)
MAST::FieldFunction< RealMatrixX > * _thermal_exp_mat_1d
virtual const MAST::FieldFunction< RealMatrixX > & stiffness_matrix(const unsigned int dim, const bool plane_stress=true)
MAST::FieldFunction< RealMatrixX > * _stiff_mat_2d
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...
MAST::FieldFunction< RealMatrixX > * _thermal_capacitance_mat_3d
MAST::FieldFunction< RealMatrixX > * _stiff_mat_3d
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...
MAST::FieldFunction< RealMatrixX > * _thermal_conductance_mat_1d
ThermalExpansionMatrix(unsigned int dim, const MAST::FieldFunction< Real > &alpha)
virtual const MAST::FieldFunction< RealMatrixX > & transverse_shear_stiffness_matrix()
MAST::FieldFunction< RealMatrixX > * _thermal_exp_mat_2d
MAST::FieldFunction< RealMatrixX > * _thermal_capacitance_mat_1d
MAST::FieldFunction< RealMatrixX > * _thermal_capacitance_mat_2d
MAST::FieldFunction< RealMatrixX > * _thermal_conductance_mat_2d
MAST::FieldFunction< RealMatrixX > * _stiff_mat_1d
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...