MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
orthotropic_material_property_card.cpp
Go to the documentation of this file.
1 /*
2  * MAST: Multidisciplinary-design Adaptation and Sensitivity Toolkit
3  * Copyright (C) 2013-2020 Manav Bhatia and MAST authors
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 OrthotropicMaterialProperty {
28 
29 
30  class StiffnessMatrix1D: public MAST::FieldFunction<RealMatrixX> {
31  public:
32 
34  const MAST::FieldFunction<Real>& G12);
35 
36  virtual ~StiffnessMatrix1D() { }
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:
58 
59 
61 
62  virtual void operator() (const libMesh::Point& p,
63  const Real t,
64  RealMatrixX& m) const;
65 
66  virtual void derivative( const MAST::FunctionBase& f,
67  const libMesh::Point& p,
68  const Real t,
69  RealMatrixX& m) const;
70 
71  protected:
72 
74  };
75 
76 
77  class StiffnessMatrix2D: public MAST::FieldFunction<RealMatrixX> {
78  public:
80  const MAST::FieldFunction<Real>& E22,
81  const MAST::FieldFunction<Real>& E33,
82  const MAST::FieldFunction<Real>& nu12,
83  const MAST::FieldFunction<Real>& nu23,
84  const MAST::FieldFunction<Real>& nu31,
85  const MAST::FieldFunction<Real>& G12,
86  bool plane_stress);
87 
88  virtual ~StiffnessMatrix2D() { }
89 
90 
91  virtual void operator() (const libMesh::Point& p,
92  const Real t,
93  RealMatrixX& m) const;
94 
95  virtual void derivative( const MAST::FunctionBase& f,
96  const libMesh::Point& p,
97  const Real t,
98  RealMatrixX& m) const;
99 
100  protected:
101 
110  };
111 
112 
113 
114  class StiffnessMatrix3D: public MAST::FieldFunction<RealMatrixX> {
115  public:
117  const MAST::FieldFunction<Real>& E22,
118  const MAST::FieldFunction<Real>& E33,
119  const MAST::FieldFunction<Real>& nu12,
120  const MAST::FieldFunction<Real>& nu31,
121  const MAST::FieldFunction<Real>& nu23,
122  const MAST::FieldFunction<Real>& G12,
123  const MAST::FieldFunction<Real>& G13,
124  const MAST::FieldFunction<Real>& G23);
125 
126 
127  virtual ~StiffnessMatrix3D() { }
128 
129  virtual void operator() (const libMesh::Point& p,
130  const Real t,
131  RealMatrixX& m) const;
132 
133  virtual void derivative( const MAST::FunctionBase& f,
134  const libMesh::Point& p,
135  const Real t,
136  RealMatrixX& m) const;
137 
138  protected:
139 
141  const MAST::FieldFunction<Real> &_nu12, &_nu31, &_nu23;
143 
144  };
145 
146 
147 
148  class InertiaMatrix3D: public MAST::FieldFunction<RealMatrixX> {
149  public:
151 
152  virtual ~InertiaMatrix3D() { }
153 
154  virtual void operator() (const libMesh::Point& p,
155  const Real t,
156  RealMatrixX& m) const;
157 
158  virtual void derivative( const MAST::FunctionBase& f,
159  const libMesh::Point& p,
160  const Real t,
161  RealMatrixX& m) const;
162 
163  protected:
164 
166 
167  };
168 
169 
170 
171  class ThermalExpansionMatrix: public MAST::FieldFunction<RealMatrixX> {
172  public:
173 
175  MAST::FieldFunction<RealMatrixX>("ThermalExpansionMatrix"),
176  _dim(1),
177  _alpha(_dim) {
178 
179  _alpha[0] = &alpha11;
180  _functions.insert(_alpha[0]);
181  }
182 
183 
185  const MAST::FieldFunction<Real>& alpha22):
186  MAST::FieldFunction<RealMatrixX>("ThermalExpansionMatrix"),
187  _dim(2),
188  _alpha(_dim) {
189 
190  _alpha[0] = &alpha11;
191  _alpha[1] = &alpha22;
192  for (unsigned int i=0; i<_dim; i++)
193  _functions.insert(_alpha[i]);
194  }
195 
196 
198  const MAST::FieldFunction<Real>& alpha22,
199  const MAST::FieldFunction<Real>& alpha33):
200  MAST::FieldFunction<RealMatrixX>("ThermalExpansionMatrix"),
201  _dim(3),
202  _alpha(_dim) {
203 
204  _alpha[0] = &alpha11;
205  _alpha[1] = &alpha22;
206  _alpha[2] = &alpha33;
207 
208  for (unsigned int i=0; i<_dim; i++)
209  _functions.insert(_alpha[i]);
210  }
211 
212 
214 
215  virtual void operator() (const libMesh::Point& p,
216  const Real t,
217  RealMatrixX& m) const;
218 
219  virtual void derivative( const MAST::FunctionBase& f,
220  const libMesh::Point& p,
221  const Real t,
222  RealMatrixX& m) const;
223 
224  protected:
225 
226  const unsigned int _dim;
227 
228  std::vector<const MAST::FieldFunction<Real>*> _alpha;
229  };
230 
231 
232 
233 
235  public MAST::FieldFunction<RealMatrixX> {
236  public:
237 
239  MAST::FieldFunction<RealMatrixX>("ThermalConductanceMatrix"),
240  _dim(1),
241  _k(_dim) {
242 
243  _k[0] = &k11;
244 
245  _functions.insert(_k[0]);
246  }
247 
249  const MAST::FieldFunction<Real>& k22):
250  MAST::FieldFunction<RealMatrixX>("ThermalConductanceMatrix"),
251  _dim(2),
252  _k(_dim) {
253 
254  _k[0] = &k11;
255  _k[1] = &k22;
256 
257  for (unsigned int i=0; i<_dim; i++)
258  _functions.insert(_k[i]);
259  }
260 
262  const MAST::FieldFunction<Real>& k22,
263  const MAST::FieldFunction<Real>& k33):
264  MAST::FieldFunction<RealMatrixX>("ThermalConductanceMatrix"),
265  _dim(3),
266  _k(_dim) {
267 
268  _k[0] = &k11;
269  _k[1] = &k22;
270  _k[2] = &k33;
271 
272  for (unsigned int i=0; i<_dim; i++)
273  _functions.insert(_k[i]);
274  }
275 
276 
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  protected:
289 
290  const unsigned int _dim;
291 
292  std::vector<const MAST::FieldFunction<Real>*> _k;
293  };
294 
295 
296 
298  public MAST::FieldFunction<RealMatrixX> {
299  public:
300 
301  ThermalCapacitanceMatrix(unsigned int dim,
302  const MAST::FieldFunction<Real>& rho,
303  const MAST::FieldFunction<Real>& cp):
304  MAST::FieldFunction<RealMatrixX>("ThermalCapacitanceMatrix"),
305  _dim(dim),
306  _rho(rho),
307  _cp(cp) {
308 
309  _functions.insert(&_rho);
310  _functions.insert(&_cp);
311  }
312 
313 
315 
316  virtual void operator() (const libMesh::Point& p,
317  const Real t,
318  RealMatrixX& m) const;
319 
320  virtual void derivative( const MAST::FunctionBase& f,
321  const libMesh::Point& p,
322  const Real t,
323  RealMatrixX& m) const;
324 
325  protected:
326 
327  const unsigned int _dim;
328 
330 
332  };
333 
334 
335  }
336 }
337 
338 
339 
342  const MAST::FieldFunction<Real>& G12):
343 MAST::FieldFunction<RealMatrixX>("StiffnessMatrix1D"),
344 _E11(E11),
345 _G12(G12)
346 {
347  _functions.insert(&E11);
348  _functions.insert(&G12);
349 }
350 
351 
352 
353 
354 void
356 StiffnessMatrix1D::operator() (const libMesh::Point& p,
357  const Real t,
358  RealMatrixX& m) const {
359  m = RealMatrixX::Zero(2,2);
360  Real E11, G12;
361  _E11(p, t, E11);
362  _G12(p, t, G12);
363  m(0,0) = E11;
364  m(1,1) = G12;
365 }
366 
367 
368 void
371  const libMesh::Point &p,
372  const Real t,
373  RealMatrixX &m) const {
374 
375 
376  RealMatrixX dm;
377  m = RealMatrixX::Zero(2,2); dm = RealMatrixX::Zero(2,2);
378  Real E11, G12, dE11df,dG12df;
379  _E11(p, t, E11); _E11.derivative(f, p, t, dE11df);
380  _G12(p, t, G12); _G12.derivative(f, p, t, dG12df);
381 
382  m(0,0) = dE11df;
383  m(1,1) = dG12df;
384 }
385 
386 
389 MAST::FieldFunction<RealMatrixX>("TransverseShearStiffnessMatrix"),
390 _G12(G12)
391 {
392  _functions.insert(&G12);
393 }
394 
395 
396 
397 void
400  const Real t,
401  RealMatrixX& m) const {
402  // TODO: Allow for use of G1z and G2z for compatibility with Nastran input.
412  m = RealMatrixX::Zero(2,2);
413  Real G12;
414  _G12(p, t, G12);
415  m(0,0) = G12;
416  m(1,1) = G12;
417 }
418 
419 
420 
421 void
424  const libMesh::Point& p,
425  const Real t,
426  RealMatrixX& m) const {
427  RealMatrixX dm;
428  m = RealMatrixX::Zero(2,2); dm = RealMatrixX::Zero(2, 2);
429 
430  Real G12, dG12;
431  _G12(p, t, G12); _G12.derivative(f, p, t, dG12);
432  m(0,0) = dG12;
433  m(1,1) = dG12;
434 
435 
436 // // parM/parE * parE/parf
437 // dm(0,0) = 1./2./(1.+nu)*kappa;
438 // dm(1,1) = dm(0,0);
439 // m += dEdf * dm;
440 //
441 // // parM/parnu * parnu/parf
442 // dm(0,0) = -E/2./pow(1.+nu,2)*kappa;
443 // dm(1,1) = dm(0,0);
444 // m+= dnudf*dm;
445 //
446 // // parM/parnu * parkappa/parf
447 //
448 // dm(0,0) = G; dm(1,1) = G;
449 // dm += dkappadf*dm;
450 }
451 
452 
453 
454 
457  const MAST::FieldFunction<Real>& E22,
458  const MAST::FieldFunction<Real>& E33,
459  const MAST::FieldFunction<Real>& nu12,
460  const MAST::FieldFunction<Real>& nu23,
461  const MAST::FieldFunction<Real>& nu31,
462  const MAST::FieldFunction<Real>& G12,
463  bool plane_stress ):
464 MAST::FieldFunction<RealMatrixX>("StiffnessMatrix2D"),
465 _E11(E11),
466 _E22(E22),
467 _E33(E33),
468 _nu12(nu12),
469 _nu23(nu23),
470 _nu31(nu31),
471 _G12(G12),
472 _plane_stress(plane_stress)
473 {
474  _functions.insert(&E11);
475  _functions.insert(&E22);
476  _functions.insert(&nu12);
477  _functions.insert(&G12);
478 }
479 
480 
481 
482 
483 void
485 StiffnessMatrix2D::operator() (const libMesh::Point& p,
486  const Real t,
487  RealMatrixX& m) const {
488  // https://www.efunda.com/formulae/solid_mechanics/mat_mechanics/hooke_orthotropic.cfm
489  // Applying plane strain and plane stress assumptions to the 3D matrix at
490  // the website above, results in the same matrix for both. So, I think
491  // plane stress, and plane strain, may be equivalent for orthotropic
492  // materials.
493  // libmesh_assert(_plane_stress); // currently only implemented for plane stress
494  m = RealMatrixX::Zero(3,3);
495  Real E11, E22, E33, nu12, nu23, nu31, nu21, G12, D;
496  _E11 (p, t, E11);
497  _E22 (p, t, E22);
498  _E33 (p, t, E33);
499  _nu12 (p, t, nu12);
500  _nu23 (p, t, nu23);
501  _nu31 (p, t, nu31);
502  _G12 (p, t, G12);
503  nu21 = nu12 * E22/E11;
504 
505  if (_plane_stress) // plane stress
506  {
507  D = (1. - nu12*nu21);
508 
509  m(0,0) = E11;
510  m(0,1) = nu21*E11;
511 
512  m(1,0) = nu12*E22;
513  m(1,1) = E22;
514 
515  m.topLeftCorner(2, 2) /= D;
516 
517  m(2,2) = G12;
518  }
519  else // plane strain
520  {
521  m(0,0) = -((E11*E11)*E33*(E22-E33*(nu23*nu23)))/((E22*E22)*E33
522  *(nu12*nu12)+E11*(E33*E33)*(nu23*nu23)+(E11*E11)*E22
523  *(nu31*nu31)-E11*E22*E33+E11*E22*E33*nu12*nu23*nu31*2.0);
524 
525  m(1,1) = -(E11*(E22*E22)*(E33-E11*(nu31*nu31)))/((E22*E22)*E33
526  *(nu12*nu12)+E11*(E33*E33)*(nu23*nu23)+(E11*E11)*E22
527  *(nu31*nu31)-E11*E22*E33+E11*E22*E33*nu12*nu23*nu31*2.0);
528 
529  m(2,2) = G12;
530 
531  m(0,1) = m(1,0) = -(E11*E22*E33*(E22*nu12+E11*nu23*nu31))
532  /((E22*E22)*E33*(nu12*nu12)+E11*(E33*E33)*(nu23*nu23)
533  +(E11*E11)*E22*(nu31*nu31)-E11*E22*E33+E11*E22*E33
534  *nu12*nu23*nu31*2.0);
535  }
536 }
537 
538 
539 
540 
541 void
544  const libMesh::Point& p,
545  const Real t,
546  RealMatrixX& m) const {
547  // https://www.efunda.com/formulae/solid_mechanics/mat_mechanics/hooke_orthotropic.cfm
548  // Applying plane strain and plane stress assumptions to the 3D matrix at
549  // the website above, results in the same matrix for both. So, I think
550  // plane stress, and plane strain, may be equivalent for orthotropic
551  // materials.
552  // libmesh_assert(_plane_stress); // currently only implemented for plane stress
553  RealMatrixX dm;
554  m = RealMatrixX::Zero(3,3); dm = RealMatrixX::Zero(3, 3);
555  Real E11, E22, E33, nu12, nu23, nu31, nu21, D;
556  Real dE11df, dE22df, dE33df, dnu12df, dnu23df, dnu31df, dnu21df, dDdf, dG12df;
557  _E11 (p, t, E11); _E11.derivative( f, p, t, dE11df);
558  _E22 (p, t, E22); _E22.derivative( f, p, t, dE22df);
559  _E33 (p, t, E33); _E33.derivative( f, p, t, dE33df);
560  _nu12 (p, t, nu12); _nu12.derivative( f, p, t, dnu12df);
561  _nu23 (p, t, nu12); _nu23.derivative( f, p, t, dnu12df);
562  _nu31 (p, t, nu12); _nu31.derivative( f, p, t, dnu12df);
563  _G12.derivative( f, p, t, dG12df);
564 
565  if (_plane_stress)
566  {
567  nu21 = nu12 * E22/E11;
568  dnu21df = dnu12df * E22/E11 + nu12 * dE22df/E11 - nu12 * E22/E11/E11*dE11df;
569  D = (1. - nu12*nu21);
570  dDdf = (- dnu12df*nu21 - nu12*dnu21df);
571 
572  m(0,0) = E11;
573  m(0,1) = nu21*E11;
574 
575  m(1,0) = nu12*E22;
576  m(1,1) = E22;
577 
578  m.topLeftCorner(2, 2) *= -dDdf/D/D;
579  m(2,2) = dG12df;
580 
581  dm(0,0) = dE11df;
582  dm(0,1) = nu21*dE11df + dnu21df*E11;
583 
584  dm(1,0) = nu12*dE22df + dnu12df*E22;
585  dm(1,1) = dE22df;
586 
587  m.topLeftCorner(2, 2) += dm.topLeftCorner(3, 3)/D;
588  }
589  else
590  {
591  m(0,0) = 1.0/pow((E22*E22)*E33*(nu12*nu12)+E11*(E33*E33)*(nu23*nu23)
592  + (E11*E11)*E22*(nu31*nu31)-E11*E22*E33+E11*E22*E33*nu12*nu23*nu31*2.0,2.0)
593  *(E11*(E11*(E22*E22)*(E33*E33)*dE11df+E11*(E33*E33*E33*E33)*dE11df
594  *(nu23*nu23*nu23*nu23)-E11*E22*(E33*E33*E33)*dE11df*(nu23*nu23)*2.0)
595  +nu31*(E11*((E11*E11*E11)*(E22*E22)*E33*dnu31df*2.0+(E11*E11)*(E22*E22)
596  *(E33*E33)*dnu12df*nu23*2.0-(E11*E11)*E22*(E33*E33*E33)*dnu12df
597  *(nu23*nu23*nu23)*2.0-(E11*E11*E11)*E22*(E33*E33)*dnu31df
598  *(nu23*nu23)*2.0)+E11*nu12*((E11*E11)*(E22*E22)*(E33*E33)*dnu23df*2.0
599  -(E11*E11)*(E33*E33*E33)*dE22df*(nu23*nu23*nu23)*2.0+(E11*E11)*E22*(E33*E33)
600  *dE33df*(nu23*nu23*nu23)*2.0+(E11*E11)*E22*(E33*E33*E33)*dnu23df
601  *(nu23*nu23)*2.0-E11*(E22*E22)*(E33*E33)*dE11df*nu23*2.0+E11*E22
602  *(E33*E33*E33)*dE11df*(nu23*nu23*nu23)*2.0))+E11*nu12*(E11*(E22*E22*E22)
603  *(E33*E33)*dnu12df*2.0-E11*(E22*E22)*(E33*E33*E33)*dnu12df*(nu23*nu23)
604  *2.0+(E11*E11)*(E22*E22)*(E33*E33)*dnu31df*nu23*2.0-(E11*E11)*E22*(E33*E33*E33)
605  *dnu31df*(nu23*nu23*nu23)*2.0)-E11*(nu31*nu31)*((E11*E11*E11)*(E22*E22)
606  *dE33df+(E11*E11*E11)*(E33*E33)*dE22df*(nu23*nu23)-(E11*E11*E11)*E22*E33
607  *dE33df*(nu23*nu23)*2.0-(E11*E11*E11)*E22*(E33*E33)*dnu23df*nu23*2.0)+
608  E11*(nu12*nu12)*((E22*E22*E22)*(E33*E33)*dE11df*-2.0+E11*(E22*E22)*(E33*E33)
609  *dE22df+(E22*E22)*(E33*E33*E33)*dE11df*(nu23*nu23)*2.0+E11*(E22*E22)
610  *(E33*E33)*dE33df*(nu23*nu23)-E11*E22*(E33*E33*E33)*dE22df*(nu23*nu23)
611  *2.0+E11*(E22*E22)*(E33*E33*E33)*dnu23df*nu23*2.0));
612 
613  m(1,1) = 1.0/pow((E22*E22)*E33*(nu12*nu12)+E11*(E33*E33)*(nu23*nu23)+(E11*E11)
614  *E22*(nu31*nu31)-E11*E22*E33+E11*E22*E33*nu12*nu23*nu31*2.0,2.0)
615  *(E22*((E11*E11)*E22*(E33*E33)*dE22df+(E11*E11*E11*E11)*E22*dE22df
616  *(nu31*nu31*nu31*nu31)-(E11*E11*E11)*E22*E33*dE22df*(nu31*nu31)*2.0)
617  +nu23*(E22*((E11*E11)*E22*(E33*E33*E33)*dnu23df*2.0+(E11*E11)*(E22*E22)
618  *(E33*E33)*dnu12df*nu31*2.0-(E11*E11*E11)*(E22*E22)*E33*dnu12df
619  *(nu31*nu31*nu31)*2.0-(E11*E11*E11)*E22*(E33*E33)*dnu23df*(nu31*nu31)
620  *2.0)+E22*nu12*((E11*E11)*(E22*E22)*(E33*E33)*dnu31df*2.0-(E11*E11*E11)
621  *(E22*E22)*dE33df*(nu31*nu31*nu31)*2.0+(E11*E11)*(E22*E22)*E33*dE11df
622  *(nu31*nu31*nu31)*2.0+(E11*E11*E11)*(E22*E22)*E33*dnu31df*(nu31*nu31)*2.0
623  -(E11*E11)*E22*(E33*E33)*dE22df*nu31*2.0+(E11*E11*E11)*E22*E33*dE22df
624  *(nu31*nu31*nu31)*2.0))+E22*nu12*(E11*(E22*E22*E22)*(E33*E33)*dnu12df*2.0
625  -(E11*E11)*(E22*E22*E22)*E33*dnu12df*(nu31*nu31)*2.0+(E11*E11)*(E22*E22)
626  *(E33*E33)*dnu23df*nu31*2.0-(E11*E11*E11)*(E22*E22)*E33*dnu23df
627  *(nu31*nu31*nu31)*2.0)-E22*(nu12*nu12)*((E22*E22*E22)*(E33*E33)*dE11df
628  +(E11*E11)*(E22*E22*E22)*dE33df*(nu31*nu31)-E11*(E22*E22*E22)*E33*dE11df
629  *(nu31*nu31)*2.0-(E11*E11)*(E22*E22*E22)*E33*dnu31df*nu31*2.0)+E22
630  *(nu23*nu23)*((E11*E11)*(E33*E33*E33)*dE22df*-2.0+(E11*E11)*E22*(E33*E33)
631  *dE33df+(E11*E11*E11)*(E33*E33)*dE22df*(nu31*nu31)*2.0+(E11*E11)*E22
632  *(E33*E33)*dE11df*(nu31*nu31)-(E11*E11*E11)*E22*E33*dE33df*(nu31*nu31)
633  *2.0+(E11*E11*E11)*E22*(E33*E33)*dnu31df*nu31*2.0));
634 
635  m(2,2) = dG12df;
636 
637  m(0,1) = m(1,0) = -1.0/pow((E22*E22)*E33*(nu12*nu12)+E11*(E33*E33)
638  *(nu23*nu23)+(E11*E11)*E22*(nu31*nu31)-E11*E22*E33+E11*E22*E33*nu12*nu23
639  *nu31*2.0,2.0)*(-nu31*((E33*E33)*((E11*E11*E11)*(E22*E22)*dnu23df+(E11*E11)
640  *(E22*E22)*dE11df*nu23+(E11*E11*E11)*E22*dE33df*(nu23*nu23*nu23))
641  -(E33*E33*E33)*((E11*E11*E11)*dE22df*(nu23*nu23*nu23)+(E11*E11)*E22*dE11df
642  *(nu23*nu23*nu23)-(E11*E11*E11)*E22*dnu23df*(nu23*nu23)))
643  +(nu31*nu31*nu31)*((E11*E11*E11*E11)*(E22*E22)*E33*dnu23df+(E11*E11*E11*E11)
644  *(E22*E22)*dE33df*nu23)-(nu12*nu12)*((E33*E33)*(E11*(E22*E22*E22*E22)
645  *dnu12df+(E11*E11)*(E22*E22*E22)*dnu31df*nu23)-(E33*E33)*nu31*(-(E11*E11)
646  *(E22*E22*E22)*dnu23df+E11*(E22*E22*E22)*dE11df*nu23*2.0+(E11*E11)*(E22*E22)
647  *dE22df*nu23))+(nu31*nu31)*(E33*((E11*E11*E11)*(E22*E22*E22)*dnu12df
648  -(E11*E11*E11*E11)*(E22*E22)*dnu31df*nu23)-(E11*E11*E11)*(E22*E22)*(E33*E33)
649  *dnu12df*(nu23*nu23)*2.0)+nu12*((nu31*nu31)*(-E33*((E11*E11)
650  *(E22*E22*E22)*dE11df-(E11*E11*E11)*(E22*E22)*dE22df)+(E11*E11*E11)*(E22*E22*E22)
651  *dE33df+(E11*E11)*(E22*E22)*(E33*E33)*dE11df*(nu23*nu23)*2.0)-(E33*E33)
652  *((E11*E11)*(E22*E22)*dE22df+(E11*E11)*(E22*E22)*dE33df*(nu23*nu23))
653  +(E33*E33*E33)*((E11*E11)*E22*dE22df*(nu23*nu23)*2.0-(E11*E11)*(E22*E22)
654  *dnu23df*nu23*2.0)-nu31*((E11*E11*E11)*(E22*E22*E22)*E33*dnu31df*2.0
655  +(E11*E11)*(E22*E22*E22)*(E33*E33)*dnu12df*nu23*2.0))+(E33*E33*E33)
656  *((E11*E11*E11)*E22*dnu31df*(nu23*nu23*nu23)+(E11*E11)*(E22*E22)*dnu12df
657  *(nu23*nu23))-(E33*E33)*((E11*E11)*(E22*E22*E22)*dnu12df+(E11*E11*E11)
658  *(E22*E22)*dnu31df*nu23)+(E22*E22*E22*E22)*(E33*E33)*dE11df
659  *(nu12*nu12*nu12));
660  }
661 }
662 
663 
664 
667  const MAST::FieldFunction<Real>& E22,
668  const MAST::FieldFunction<Real>& E33,
669  const MAST::FieldFunction<Real>& nu12,
670  const MAST::FieldFunction<Real>& nu31,
671  const MAST::FieldFunction<Real>& nu23,
672  const MAST::FieldFunction<Real>& G12,
673  const MAST::FieldFunction<Real>& G13,
674  const MAST::FieldFunction<Real>& G23):
675 MAST::FieldFunction<RealMatrixX>("StiffnessMatrix3D"),
676 _E11(E11),
677 _E22(E22),
678 _E33(E33),
679 _nu12(nu12),
680 _nu31(nu31),
681 _nu23(nu23),
682 _G12(G12),
683 _G13(G13),
684 _G23(G23)
685 {
686  _functions.insert(&E11);
687  _functions.insert(&E22);
688  _functions.insert(&E33);
689  _functions.insert(&nu12);
690  _functions.insert(&nu31);
691  _functions.insert(&nu23);
692  _functions.insert(&G12);
693  _functions.insert(&G13);
694  _functions.insert(&G23);
695 }
696 
697 
698 
699 
700 
701 void
703 StiffnessMatrix3D::operator() (const libMesh::Point& p,
704  const Real t,
705  RealMatrixX& m) const {
706  m = RealMatrixX::Zero(6,6);
707  Real E11, E22, E33, G12, G13, G23, nu12, nu13, nu23, nu21, nu31, nu32, D;
708  _E11 (p, t, E11);
709  _E22 (p, t, E22);
710  _E33 (p, t, E33);
711  _nu12(p, t, nu12);
712  _nu31(p, t, nu31);
713  _nu23(p, t, nu23);
714  _G12 (p, t, G12);
715  _G13 (p, t, G13);
716  _G23 (p, t, G23);
717  nu21 = nu12 * E22/E11;
718  nu13 = nu31 * E11/E33;
719  nu32 = nu23 * E33/E22;
720  D = 1.- nu12*nu21 - nu13*nu31 - nu23*nu32 - nu12*nu23*nu31 - nu13*nu21*nu32;
721 
722 
723  m(0,0) = ( 1. - nu23*nu32)*E11;
724  m(0,1) = (nu21 + nu23*nu31)*E11;
725  m(0,2) = (nu31 + nu21*nu32)*E11;
726 
727  m(1,0) = (nu12 + nu13*nu32)*E22;
728  m(1,1) = ( 1. - nu13*nu31)*E22;
729  m(1,2) = (nu32 + nu12*nu31)*E22;
730 
731  m(2,0) = (nu13 + nu12*nu23)*E33;
732  m(2,1) = (nu23 + nu13*nu21)*E33;
733  m(2,2) = ( 1. - nu12*nu21)*E33;
734 
735  m.topLeftCorner(3, 3) /= D;
736  m(3,3) = G12;
737  m(4,4) = G23;
738  m(5,5) = G13;
739 }
740 
741 
742 
743 
744 void
747  const libMesh::Point& p,
748  const Real t,
749  RealMatrixX& m) const {
750  RealMatrixX dm;
751  m = RealMatrixX::Zero(6,6); dm = RealMatrixX::Zero(6,6);
752  Real
753  E11, dE11df,
754  E22, dE22df,
755  E33, dE33df,
756  dG12df,
757  dG13df,
758  dG23df,
759  nu12, dnu12df,
760  nu13, dnu13df,
761  nu23, dnu23df,
762  nu21, dnu21df,
763  nu31, dnu31df,
764  nu32, dnu32df,
765  D, dDdf;
766  _E11 (p, t, E11); _E11.derivative( f, p, t, dE11df);
767  _E22 (p, t, E22); _E22.derivative( f, p, t, dE22df);
768  _E33 (p, t, E33); _E33.derivative( f, p, t, dE33df);
769  _nu12 (p, t, nu12); _nu12.derivative( f, p, t, dnu12df);
770  _nu31 (p, t, nu31); _nu31.derivative( f, p, t, dnu31df);
771  _nu23 (p, t, nu23); _nu23.derivative( f, p, t, dnu23df);
772  _G12.derivative( f, p, t, dG12df);
773  _G13.derivative( f, p, t, dG13df);
774  _G23.derivative( f, p, t, dG23df);
775  nu21 = nu12 * E22/E11;
776  dnu21df = dnu12df * E22/E11 + nu12 * dE22df/E11 - nu12 * E22/E11/E11*dE11df;
777  nu13 = nu31 * E11/E33;
778  dnu13df = dnu31df * E11/E33 + nu31 * dE11df/E33 - nu31 * E11/E33/E33*dE33df;
779  nu32 = nu23 * E33/E22;
780  dnu32df = dnu23df * E33/E22 + nu23 * dE33df/E22 - nu23 * E33/E22/E22*dE22df;
781  D = 1.- nu12*nu21 - nu13*nu31 - nu23*nu32 - nu12*nu23*nu31 - nu13*nu21*nu32;
782  dDdf =
783  - dnu12df*nu21 - nu12*dnu21df
784  - dnu13df*nu31 - nu13*dnu31df
785  - dnu23df*nu32 - nu23*dnu32df
786  - dnu12df*nu23*nu31
787  - nu12*dnu23df*nu31
788  - nu12*nu23*dnu31df
789  - dnu13df*nu21*nu32
790  - nu13*dnu21df*nu32
791  - nu13*nu21*dnu32df;
792 
793  m(0,0) = ( 1. - nu23*nu32)*E11;
794  m(0,1) = (nu21 + nu23*nu31)*E11;
795  m(0,2) = (nu31 + nu21*nu32)*E11;
796 
797  m(1,0) = (nu12 + nu13*nu32)*E22;
798  m(1,1) = ( 1. - nu13*nu31)*E22;
799  m(1,2) = (nu32 + nu12*nu31)*E22;
800 
801  m(2,0) = (nu13 + nu12*nu23)*E33;
802  m(2,1) = (nu23 + nu13*nu21)*E33;
803  m(2,2) = ( 1. - nu12*nu21)*E33;
804  m *= -dDdf/D/D;
805 
806  m(3,3) = dG12df;
807  m(4,4) = dG23df;
808  m(5,5) = dG13df;
809 
810 
811  dm(0,0) = (- dnu23df*nu32 - nu23*dnu32df)*E11 + ( 1. - nu23*nu32)*dE11df;
812  dm(0,1) = (dnu21df + dnu23df*nu31 + nu23*dnu31df)*E11 + (nu21 + nu23*nu31)*dE11df;
813  dm(0,2) = (dnu31df + dnu21df*nu32 + nu21*dnu32df)*E11 + (nu31 + nu21*nu32)*dE11df;
814 
815  dm(1,0) = (dnu12df + dnu13df*nu32 + nu13*dnu32df)*E22 + (nu12 + nu13*nu32)*dE22df;
816  dm(1,1) = (- dnu13df*nu31 - nu13*dnu31df)*E22 + ( 1. - nu13*nu31)*dE22df;
817  dm(1,2) = (dnu32df + dnu12df*nu31 + nu12*dnu31df)*E22 + (nu32 + nu12*nu31)*dE22df;
818 
819  dm(2,0) = (dnu13df + dnu12df*nu23 + nu12*dnu23df)*E33 + (nu13 + nu12*nu23)*dE33df;
820  dm(2,1) = (dnu23df + dnu13df*nu21 + nu13*dnu21df)*E33 + (nu23 + nu13*nu21)*dE33df;
821  dm(2,2) = (- dnu12df*nu21 - nu12*dnu21df)*E33 + ( 1. - nu12*nu21)*dE33df;
822 
823  m += dm/D;
824 }
825 
826 
827 
830 MAST::FieldFunction<RealMatrixX>("InertiaMatrix3D"),
831 _rho(rho) {
832 
833  _functions.insert(&rho);
834 }
835 
836 
837 
838 void
840 InertiaMatrix3D::operator() (const libMesh::Point& p,
841  const Real t,
842  RealMatrixX& m) const {
843  m = RealMatrixX::Zero(3,3);
844  Real rho;
845  _rho(p, t, rho);
846  m(0,0) = rho;
847  m(1,1) = rho;
848  m(2,2) = rho;
849 }
850 
851 
852 void
855  const MAST::FunctionBase &f,
856  const libMesh::Point &p,
857  const Real t,
858  RealMatrixX &m) const {
859 
860 
861  m = RealMatrixX::Zero(3,3);
862  Real rho;
863  _rho.derivative( f, p, t, rho);
864  m(0,0) = rho;
865  m(1,1) = rho;
866  m(2,2) = rho;
867 }
868 
869 
870 void
872 operator() (const libMesh::Point& p,
873  const Real t,
874  RealMatrixX& m) const {
875 
876 
877  Real alpha;
878  switch (_dim) {
879  case 1:
880  m.setZero(2, 1);
881  break;
882 
883  case 2:
884  m.setZero(3, 1);
885  break;
886 
887  case 3:
888  m.setZero(6, 1);
889  break;
890  }
891 
892  for (unsigned int i=0; i<_dim; i++) {
893  (*_alpha[i]) (p, t, alpha);
894  m(i,0) = alpha;
895  }
896 }
897 
898 
899 
900 
901 
902 void
905  const libMesh::Point& p,
906  const Real t,
907  RealMatrixX& m) const {
908 
909 
910 
911  Real alpha;
912  switch (_dim) {
913  case 1:
914  m.setZero(2, 1);
915  break;
916 
917  case 2:
918  m.setZero(3, 1);
919  break;
920 
921  case 3:
922  m.setZero(6, 1);
923  break;
924  }
925 
926  for (unsigned int i=0; i<_dim; i++) {
927  _alpha[i]->derivative( f, p, t, alpha);
928  m(i,0) = alpha;
929  }
930 }
931 
932 
933 
934 
935 
936 void
938 operator() (const libMesh::Point& p,
939  const Real t,
940  RealMatrixX& m) const {
941 
942  Real cp, rho;
943  _cp (p, t, cp);
944  _rho (p, t, rho);
945 
946  m.setZero(1,1);
947 
948  m(0,0) = cp*rho;
949 }
950 
951 
952 
953 
954 
955 void
958  const libMesh::Point& p,
959  const Real t,
960  RealMatrixX& m) const {
961 
962 
963  Real cp, dcp, rho, drho;
964  _cp (p, t, cp); _cp.derivative( f, p, t, dcp);
965  _rho (p, t, rho); _rho.derivative( f, p, t, drho);
966 
967  m.setZero(1,1);
968 
969  m(0,0) = dcp*rho + cp*drho;
970 }
971 
972 
973 
974 
975 void
977 operator() (const libMesh::Point& p,
978  const Real t,
979  RealMatrixX& m) const {
980 
981  Real k;
982  m.setZero(_dim, _dim);
983  for (unsigned int i=0; i<_dim; i++) {
984  (*_k[i]) (p, t, k);
985  m(i,i) = k;
986  }
987 
988 }
989 
990 
991 
992 
993 
994 void
997  const libMesh::Point& p,
998  const Real t,
999  RealMatrixX& m) const {
1000 
1001  Real k;
1002  m.setZero(_dim, _dim);
1003  for (unsigned int i=0; i<_dim; i++) {
1004  _k[i]->derivative( f, p, t, k);
1005  m(i,i) = k;
1006  }
1007 }
1008 
1009 
1010 
1013 _stiff_mat_1d (nullptr),
1014 _stiff_mat_2d (nullptr),
1015 _stiff_mat_3d (nullptr),
1016 _damp_mat (nullptr),
1017 _inertia_mat_3d (nullptr),
1018 _thermal_exp_mat_1d (nullptr),
1019 _thermal_exp_mat_2d (nullptr),
1020 _thermal_exp_mat_3d (nullptr),
1021 _transverse_shear_mat (nullptr),
1022 _thermal_capacitance_mat_1d (nullptr),
1023 _thermal_capacitance_mat_2d (nullptr),
1024 _thermal_capacitance_mat_3d (nullptr),
1025 _thermal_conductance_mat_1d (nullptr),
1026 _thermal_conductance_mat_2d (nullptr),
1027 _thermal_conductance_mat_3d (nullptr)
1028 { }
1029 
1030 
1031 
1033 
1034  if (_stiff_mat_1d) delete _stiff_mat_1d;
1035  if (_stiff_mat_2d) delete _stiff_mat_2d;
1036  if (_stiff_mat_3d) delete _stiff_mat_3d;
1037 
1038  if (_damp_mat) delete _damp_mat;
1039  if (_inertia_mat_3d) delete _inertia_mat_3d;
1044 
1048 
1052 }
1053 
1054 
1055 
1058  const bool plane_stress) {
1059 
1060  switch (dim) {
1061 
1062  case 1: {
1063 
1064  if (!_stiff_mat_1d)
1066  (this->get<MAST::FieldFunction<Real> >("E11"),
1067  this->get<MAST::FieldFunction<Real> >("G12"));
1068 
1069  return *_stiff_mat_1d;
1070  }
1071  break;
1072 
1073  case 2: {
1074 
1075  if (!_stiff_mat_2d)
1077  (this->get<MAST::FieldFunction<Real> >("E11"),
1078  this->get<MAST::FieldFunction<Real> >("E22"),
1079  this->get<MAST::FieldFunction<Real> >("E33"),
1080  this->get<MAST::FieldFunction<Real> >("nu12"),
1081  this->get<MAST::FieldFunction<Real> >("nu23"),
1082  this->get<MAST::FieldFunction<Real> >("nu31"),
1083  this->get<MAST::FieldFunction<Real> >("G12"),
1084  plane_stress);
1085 
1086  return *_stiff_mat_2d;
1087  }
1088  break;
1089 
1090  case 3: {
1091 
1092  if (!_stiff_mat_3d)
1094  (this->get<MAST::FieldFunction<Real> >("E11"),
1095  this->get<MAST::FieldFunction<Real> >("E22"),
1096  this->get<MAST::FieldFunction<Real> >("E33"),
1097  this->get<MAST::FieldFunction<Real> >("nu12"),
1098  this->get<MAST::FieldFunction<Real> >("nu31"),
1099  this->get<MAST::FieldFunction<Real> >("nu23"),
1100  this->get<MAST::FieldFunction<Real> >("G12"),
1101  this->get<MAST::FieldFunction<Real> >("G13"),
1102  this->get<MAST::FieldFunction<Real> >("G23"));
1103 
1104  return *_stiff_mat_3d;
1105  }
1106  break;
1107 
1108  default:
1109  // should not get here
1110  libmesh_error();
1111  }
1112 }
1113 
1114 
1115 
1116 
1119 
1120  // make sure that this is not null
1121  libmesh_assert(false);
1122 
1123  return *_damp_mat;
1124 }
1125 
1126 
1127 
1128 
1131 
1132  switch (dim) {
1133  case 3: {
1134 
1135  if (!_inertia_mat_3d)
1137  (this->get<MAST::FieldFunction<Real> >("rho"));
1138 
1139  return *_inertia_mat_3d;
1140  }
1141  break;
1142 
1143  case 1:
1144  case 2:
1145  default:
1146  // implemented only for 3D since the 2D and 1D elements
1147  // calculate it themselves
1148  libmesh_error();
1149 
1150  }
1151 }
1152 
1153 
1154 
1155 
1158 
1159  switch (dim) {
1160  case 1: {
1161 
1162  if (!_thermal_exp_mat_1d)
1165  (this->get<MAST::FieldFunction<Real> >("alpha11_expansion"));
1166 
1167  return *_thermal_exp_mat_1d;
1168  }
1169  break;
1170 
1171 
1172  case 2: {
1173 
1174  if (!_thermal_exp_mat_2d)
1177  (this->get<MAST::FieldFunction<Real> >("alpha11_expansion"),
1178  this->get<MAST::FieldFunction<Real> >("alpha22_expansion"));
1179 
1180  return *_thermal_exp_mat_2d;
1181  }
1182  break;
1183 
1184 
1185  case 3: {
1186 
1187  if (!_thermal_exp_mat_3d)
1190  (this->get<MAST::FieldFunction<Real> >("alpha11_expansion"),
1191  this->get<MAST::FieldFunction<Real> >("alpha22_expansion"),
1192  this->get<MAST::FieldFunction<Real> >("alpha33_expansion"));
1193 
1194  return *_thermal_exp_mat_3d;
1195  }
1196  break;
1197 
1198 
1199  default:
1200  libmesh_error();
1201  }
1202 }
1203 
1204 
1205 
1206 
1209 
1210  if (!_transverse_shear_mat)
1213  (this->get<MAST::FieldFunction<Real> >("G12"));
1214 
1215  return *_transverse_shear_mat;
1216 
1217 }
1218 
1219 
1220 
1223 
1224  switch (dim) {
1225 
1226  case 1: {
1227 
1231  (dim,
1232  this->get<MAST::FieldFunction<Real> >("rho"),
1233  this->get<MAST::FieldFunction<Real> >("cp"));
1234 
1236  }
1237  break;
1238 
1239 
1240  case 2: {
1241 
1245  (dim,
1246  this->get<MAST::FieldFunction<Real> >("rho"),
1247  this->get<MAST::FieldFunction<Real> >("cp"));
1248 
1250  }
1251  break;
1252 
1253 
1254  case 3: {
1255 
1259  (dim,
1260  this->get<MAST::FieldFunction<Real> >("rho"),
1261  this->get<MAST::FieldFunction<Real> >("cp"));
1262 
1264  }
1265  break;
1266 
1267 
1268  default:
1269  // should not get here
1270  libmesh_error();
1271  }
1272 }
1273 
1274 
1275 
1276 
1277 
1280 
1281  switch (dim) {
1282 
1283  case 1: {
1284 
1288  (this->get<MAST::FieldFunction<Real> >("k11_th"));
1289 
1291  }
1292  break;
1293 
1294 
1295  case 2: {
1296 
1300  (this->get<MAST::FieldFunction<Real> >("k11_th"),
1301  this->get<MAST::FieldFunction<Real> >("k22_th"));
1302 
1304  }
1305  break;
1306 
1307 
1308  case 3: {
1309 
1313  (this->get<MAST::FieldFunction<Real> >("k11_th"),
1314  this->get<MAST::FieldFunction<Real> >("k22_th"),
1315  this->get<MAST::FieldFunction<Real> >("k33_th"));
1316 
1318  }
1319  break;
1320 
1321 
1322  default:
1323  // should not get here
1324  libmesh_error();
1325  }
1326 }
1327 
1328 
1329 
MAST::FieldFunction< RealMatrixX > * _thermal_exp_mat_1d
virtual const MAST::FieldFunction< RealMatrixX > & conductance_matrix(const unsigned int dim)
MAST::FieldFunction< RealMatrixX > * _thermal_exp_mat_3d
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 > * _thermal_capacitance_mat_3d
virtual const MAST::FieldFunction< RealMatrixX > & damping_matrix(const unsigned int dim)
ThermalExpansionMatrix(const MAST::FieldFunction< Real > &alpha11, const MAST::FieldFunction< Real > &alpha22, const MAST::FieldFunction< Real > &alpha33)
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...
StiffnessMatrix2D(const MAST::FieldFunction< Real > &E11, const MAST::FieldFunction< Real > &E22, const MAST::FieldFunction< Real > &E33, const MAST::FieldFunction< Real > &nu12, const MAST::FieldFunction< Real > &nu23, const MAST::FieldFunction< Real > &nu31, const MAST::FieldFunction< Real > &G12, bool plane_stress)
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...
ThermalExpansionMatrix(const MAST::FieldFunction< Real > &alpha11, const MAST::FieldFunction< Real > &alpha22)
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)
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...
virtual const MAST::FieldFunction< RealMatrixX > & inertia_matrix(const unsigned int dim)
MAST::FieldFunction< RealMatrixX > * _thermal_exp_mat_2d
ThermalConductanceMatrix(const MAST::FieldFunction< Real > &k11, const MAST::FieldFunction< Real > &k22)
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, ValType &v) const
calculates the value of the function derivative and returns it in v.
libMesh::Real Real
MAST::FieldFunction< RealMatrixX > * _thermal_capacitance_mat_1d
virtual const MAST::FieldFunction< RealMatrixX > & stiffness_matrix(const unsigned int dim, const bool plane_stress=true)
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...
MAST::FieldFunction< RealMatrixX > * _thermal_conductance_mat_3d
StiffnessMatrix1D(const MAST::FieldFunction< Real > &E11, const MAST::FieldFunction< Real > &G12)
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 > * _thermal_capacitance_mat_2d
MAST::FieldFunction< RealMatrixX > * _transverse_shear_mat
Matrix< Real, Dynamic, Dynamic > RealMatrixX
StiffnessMatrix3D(const MAST::FieldFunction< Real > &E11, const MAST::FieldFunction< Real > &E22, const MAST::FieldFunction< Real > &E33, const MAST::FieldFunction< Real > &nu12, const MAST::FieldFunction< Real > &nu31, const MAST::FieldFunction< Real > &nu23, const MAST::FieldFunction< Real > &G12, const MAST::FieldFunction< Real > &G13, const MAST::FieldFunction< Real > &G23)
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 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...
MAST::FieldFunction< RealMatrixX > * _thermal_conductance_mat_1d
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 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< RealMatrixX > & capacitance_matrix(const unsigned int dim)
ThermalConductanceMatrix(const MAST::FieldFunction< Real > &k11, const MAST::FieldFunction< Real > &k22, const MAST::FieldFunction< Real > &k33)
MAST::FieldFunction< RealMatrixX > * _inertia_mat_3d
MAST::FieldFunction< RealMatrixX > * _thermal_conductance_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...
virtual const MAST::FieldFunction< RealMatrixX > & transverse_shear_stiffness_matrix()
MAST::FieldFunction< RealMatrixX > * _stiff_mat_2d
MAST::FieldFunction< RealMatrixX > * _damp_mat
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...
MAST::FieldFunction< RealMatrixX > * _stiff_mat_3d