MAST
isotropic_element_property_card_3D.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 
25 
26 
27 namespace MAST {
28 
29  namespace IsotropicElementProperty3D {
30 
31 
33  public MAST::FieldFunction<RealMatrixX> {
34  public:
36 
37  virtual ~StiffnessMatrix() { }
38 
39  virtual void operator() (const libMesh::Point& p,
40  const Real t,
41  RealMatrixX& m) const;
42 
43 
44  virtual void derivative ( const MAST::FunctionBase& f,
45  const libMesh::Point& p,
46  const Real t,
47  RealMatrixX& m) const;
48 
49  protected:
50 
52  };
53 
54 
55 
56  class InertiaMatrix: public MAST::FieldFunction<RealMatrixX> {
57  public:
59 
60  virtual ~InertiaMatrix() { }
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 
78  class ThermalExpansionMatrix: public MAST::FieldFunction<RealMatrixX> {
79  public:
81  const MAST::FieldFunction<RealMatrixX>& mat_expansion);
82 
84 
85  virtual void operator() (const libMesh::Point& p,
86  const Real t,
87  RealMatrixX& m) const;
88 
89  virtual void derivative ( const MAST::FunctionBase& f,
90  const libMesh::Point& p,
91  const Real t,
92  RealMatrixX& m) const;
93 
94  protected:
95 
98  };
99 
100 
101 
102 
104  public MAST::FieldFunction<RealMatrixX> {
105  public:
107 
108  virtual ~PrestressAMatrix() { }
109 
110  virtual void operator() (const libMesh::Point& p,
111  const Real t,
112  RealMatrixX& m) const;
113 
114 
115  virtual void derivative ( const MAST::FunctionBase& f,
116  const libMesh::Point& p,
117  const Real t,
118  RealMatrixX& m) const;
119 
120  //virtual void convert_to_vector(const RealMatrixX& m, DenseRealVector& v) const;
121 
122  protected:
123 
125  };
126 
127 
128 
130  public MAST::FieldFunction<RealMatrixX> {
131 
132  public:
133 
135 
136  virtual ~ThermalConductanceMatrix();
137 
138  virtual void operator() (const libMesh::Point& p,
139  const Real t,
140  RealMatrixX& m) const;
141 
142 
143  virtual void derivative ( const MAST::FunctionBase& f,
144  const libMesh::Point& p,
145  const Real t,
146  RealMatrixX& m) const;
147 
148  protected:
149 
151  };
152 
153 
154 
155 
157  public MAST::FieldFunction<RealMatrixX> {
158 
159  public:
160 
162 
163  virtual ~ThermalCapacitanceMatrix();
164 
165  virtual void operator() (const libMesh::Point& p,
166  const Real t,
167  RealMatrixX& m) const;
168 
169 
170  virtual void derivative ( const MAST::FunctionBase& f,
171  const libMesh::Point& p,
172  const Real t,
173  RealMatrixX& m) const;
174 
175  protected:
176 
178  };
179 
180  }
181 
182 }
183 
184 
185 bool
187  return _material->depends_on(f) || // check if the material property depends on the function
188  MAST::ElementPropertyCardBase::depends_on(f); // check with this property card
189 }
190 
191 
192 
195 MAST::FieldFunction<RealMatrixX> ("StiffnessMatrix3D"),
196 _material_stiffness(mat) {
197 
198  _functions.insert(&mat);
199 }
200 
201 
202 void
204 StiffnessMatrix::operator() (const libMesh::Point& p,
205  const Real t,
206  RealMatrixX& m) const {
207  // this only returns the material stiffness
208  _material_stiffness(p, t, m);
209 }
210 
211 
212 
213 void
216  const libMesh::Point& p,
217  const Real t,
218  RealMatrixX& m) const {
219  // this only returns the material stiffness
220  _material_stiffness.derivative( f, p, t, m);
221 }
222 
223 
224 
225 
228 MAST::FieldFunction<RealMatrixX>("InertiaMatrix3D"),
229 _material_inertia(mat) {
230  _functions.insert(&mat);
231 }
232 
233 
234 
235 void
237 InertiaMatrix::operator() (const libMesh::Point& p,
238  const Real t,
239  RealMatrixX& m) const {
240 
241  _material_inertia(p, t, m);
242 }
243 
244 
245 
246 
247 void
250  const libMesh::Point& p,
251  const Real t,
252  RealMatrixX& m) const {
253 
254  _material_inertia.derivative( f, p, t, m);
255 }
256 
257 
258 
259 
262  const MAST::FieldFunction<RealMatrixX>& mat_expansion):
263 MAST::FieldFunction<RealMatrixX>("ThermalExpansionMatrix3D"),
264 _material_stiffness(mat_stiff),
265 _material_expansion(mat_expansion) {
266  _functions.insert(&mat_stiff);
267  _functions.insert(&mat_expansion);
268 }
269 
270 
271 
272 
273 void
275 ThermalExpansionMatrix::operator() (const libMesh::Point& p,
276  const Real t,
277  RealMatrixX& m) const {
278  RealMatrixX mat;
279  _material_stiffness(p, t, m);
280  _material_expansion(p, t, mat);
281  m *= mat;
282 }
283 
284 
285 
286 
287 
288 
289 void
292  const libMesh::Point& p,
293  const Real t,
294  RealMatrixX& m) const {
295  libmesh_error();
296 }
297 
298 
299 
300 
303 MAST::FieldFunction<RealMatrixX>("PrestressAMatrix3D"),
304 _prestress(prestress){
305  _functions.insert(&prestress);
306 }
307 
308 
309 
310 
311 void
313 PrestressAMatrix::operator() (const libMesh::Point& p,
314  const Real t,
315  RealMatrixX& m) const {
316  _prestress(p, t, m);
317 }
318 
319 
320 
321 
322 
323 
324 void
327  const libMesh::Point& p,
328  const Real t,
329  RealMatrixX& m) const {
330  _prestress.derivative( f, p, t, m);
331 }
332 
333 
334 
337 MAST::FieldFunction<RealMatrixX>("ThermalConductanceMatrix"),
338 _mat_cond(mat_cond) {
339 
340  _functions.insert(&mat_cond);
341 }
342 
343 
344 
345 
346 
349 
350 
351 void
353 operator() (const libMesh::Point& p,
354  const Real t,
355  RealMatrixX& m) const {
356 
357  _mat_cond(p, t, m);
358 }
359 
360 
361 
362 
363 void
365  const libMesh::Point& p,
366  const Real t,
367  RealMatrixX& m) const {
368  _mat_cond.derivative( f, p, t, m);
369 }
370 
371 
372 
373 
374 
375 
378 MAST::FieldFunction<RealMatrixX>("ThermalCapacitanceMatrix"),
379 _mat_cap(mat_cap) {
380 
381  _functions.insert(&mat_cap);
382 }
383 
384 
385 
386 
389 
390 
391 void
393 operator() (const libMesh::Point& p,
394  const Real t,
395  RealMatrixX& m) const {
396 
397  _mat_cap(p, t, m);
398 }
399 
400 
401 
402 void
404  const libMesh::Point& p,
405  const Real t,
406  RealMatrixX& m) const {
407  _mat_cap.derivative( f, p, t, m);
408 }
409 
410 
411 
412 //void
413 //MAST::IsotropicElementProperty3D::
414 //PrestressAMatrix::convert_to_vector(const RealMatrixX &m,
415 // RealVectorX &v) const {
416 // libmesh_error();
417 //}
418 
419 
420 
421 
422 std::unique_ptr<MAST::FieldFunction<RealMatrixX> >
425 
428  (_material->stiffness_matrix(3));
429 
430  return std::unique_ptr<MAST::FieldFunction<RealMatrixX> >(rval);
431 }
432 
433 
434 
435 std::unique_ptr<MAST::FieldFunction<RealMatrixX> >
438 
439  libmesh_assert(false);
440 
441  return std::unique_ptr<MAST::FieldFunction<RealMatrixX> >(nullptr);
442 }
443 
444 
445 
446 
447 std::unique_ptr<MAST::FieldFunction<RealMatrixX> >
450 
451  libmesh_assert(false);
452 
453  return std::unique_ptr<MAST::FieldFunction<RealMatrixX> >(nullptr);
454 }
455 
456 
457 
458 std::unique_ptr<MAST::FieldFunction<RealMatrixX> >
461 
462  libmesh_assert(false);
463 
464  return std::unique_ptr<MAST::FieldFunction<RealMatrixX> >(nullptr);
465 }
466 
467 
468 std::unique_ptr<MAST::FieldFunction<RealMatrixX> >
471 
474  (_material->inertia_matrix(3));
475 
476  return std::unique_ptr<MAST::FieldFunction<RealMatrixX> >(rval);
477 }
478 
479 
480 
481 std::unique_ptr<MAST::FieldFunction<RealMatrixX> >
484 
487  (_material->stiffness_matrix(3),
488  _material->thermal_expansion_matrix(3));
489 
490  return std::unique_ptr<MAST::FieldFunction<RealMatrixX> >(rval);
491 }
492 
493 
494 std::unique_ptr<MAST::FieldFunction<RealMatrixX> >
497 
498  // for 3D elements, there is no difference between the A and B matrices
499  return this->thermal_expansion_A_matrix(e);
500 }
501 
502 
503 
504 std::unique_ptr<MAST::FieldFunction<RealMatrixX> >
507 
508  libmesh_assert(false);
509 
510  return std::unique_ptr<MAST::FieldFunction<RealMatrixX> >(nullptr);
511 }
512 
513 
514 std::unique_ptr<MAST::FieldFunction<RealMatrixX> >
517 
518  libmesh_assert(false);
519 
520  return std::unique_ptr<MAST::FieldFunction<RealMatrixX> >(nullptr);
521 }
522 
523 
524 std::unique_ptr<MAST::FieldFunction<RealMatrixX> >
527 
528  libmesh_assert(false);
529 
530  return std::unique_ptr<MAST::FieldFunction<RealMatrixX> >(nullptr);
531 }
532 
533 
534 std::unique_ptr<MAST::FieldFunction<RealMatrixX> >
537 
540  (_material->conductance_matrix(3));
541 
542  return std::unique_ptr<MAST::FieldFunction<RealMatrixX> >(rval);
543 }
544 
545 
546 std::unique_ptr<MAST::FieldFunction<RealMatrixX> >
549 
552  (_material->capacitance_matrix(3));
553 
554  return std::unique_ptr<MAST::FieldFunction<RealMatrixX> >(rval);
555 }
556 
557 
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...
ThermalConductanceMatrix(const MAST::FieldFunction< RealMatrixX > &mat_cond)
InertiaMatrix(const MAST::FieldFunction< RealMatrixX > &mat)
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 > > stiffness_A_matrix(const MAST::ElementBase &e) const
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > thermal_conductance_matrix(const MAST::ElementBase &e) const
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > damping_matrix(const MAST::ElementBase &e) const
std::set< const MAST::FunctionBase * > _functions
set of functions that this function depends on
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > thermal_expansion_A_matrix(const MAST::ElementBase &e) const
libMesh::Real Real
StiffnessMatrix(const MAST::FieldFunction< RealMatrixX > &mat)
virtual void operator()(const libMesh::Point &p, const Real t, RealMatrixX &m) const
calculates the value of the function at the specified point, p, and time, t, and returns it in v...
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > stiffness_D_matrix(const MAST::ElementBase &e) const
virtual bool depends_on(const MAST::FunctionBase &f) const
returns true if the property card depends on the function f
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...
ThermalCapacitanceMatrix(const MAST::FieldFunction< RealMatrixX > &mat_cond)
ThermalExpansionMatrix(const MAST::FieldFunction< RealMatrixX > &mat_stiff, const MAST::FieldFunction< RealMatrixX > &mat_expansion)
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 > > thermal_expansion_B_matrix(const MAST::ElementBase &e) const
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > prestress_B_matrix(MAST::ElementBase &e) const
virtual void operator()(const libMesh::Point &p, const Real t, RealMatrixX &m) const
calculates the value of the function at the specified point, p, and time, t, and returns it in v...
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...
PrestressAMatrix(const MAST::FieldFunction< RealMatrixX > &prestress)
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...
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > prestress_A_matrix(MAST::ElementBase &e) const
virtual std::unique_ptr< MAST::FieldFunction< RealMatrixX > > stiffness_B_matrix(const MAST::ElementBase &e) const
const MAST::FieldFunction< RealMatrixX > & _material_inertia
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.
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...
const MAST::FieldFunction< RealMatrixX > & _material_stiffness
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_capacitance_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 std::unique_ptr< MAST::FieldFunction< RealMatrixX > > transverse_shear_stiffness_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 std::unique_ptr< MAST::FieldFunction< RealMatrixX > > inertia_matrix(const MAST::ElementBase &e) const
This is the base class for elements that implement calculation of finite element quantities over the ...
Definition: elem_base.h:72