MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
heat_conduction_transient_assembly.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 // MAST includes
25 #include "base/assembly_base.h"
26 #include "mesh/geom_elem.h"
27 
28 
32 
33 }
34 
35 
36 
37 
40 
41 }
42 
43 
44 
45 void
47 elem_calculations(bool if_jac,
48  RealVectorX& f_m,
49  RealVectorX& f_x,
50  RealMatrixX& f_m_jac_xdot,
51  RealMatrixX& f_m_jac,
52  RealMatrixX& f_x_jac) {
53 
54  libmesh_assert(_physics_elem);
55 
57  dynamic_cast<MAST::HeatConductionElementBase&>(*_physics_elem);
58 
59  f_m.setZero();
60  f_x.setZero();
61  f_m_jac_xdot.setZero();
62  f_m_jac.setZero();
63  f_x_jac.setZero();
64 
65  // assembly of the flux terms
66  e.internal_residual(if_jac, f_x, f_x_jac);
67  e.side_external_residual(if_jac, f_x, f_x_jac, _discipline->side_loads());
68  e.volume_external_residual(if_jac, f_x, f_x_jac, _discipline->volume_loads());
69 
70  //assembly of the capacitance term
71  e.velocity_residual(if_jac, f_m, f_m_jac_xdot, f_m_jac);
72 }
73 
74 
75 
76 void
79 
80  libmesh_error(); // to be implemented
81 }
82 
83 
84 
85 void
88  RealVectorX& f_m,
89  RealVectorX& f_x) {
90 
91  libmesh_assert(_physics_elem);
92 
94  dynamic_cast<MAST::HeatConductionElementBase&>(*_physics_elem);
95 
96  unsigned int
97  n = (unsigned int)f_m.size();
98 
100  dummy = RealMatrixX::Zero(n, n);
101 
102  f_m.setZero();
103  f_x.setZero();
104 
105  // assembly of the flux terms
109 
110  //assembly of the capacitance term
112 }
113 
114 
115 void
118 
119  libmesh_error(); // to be implemented
120 }
121 
122 
123 
124 void
126 set_elem_data(unsigned int dim,
127  const libMesh::Elem& ref_elem,
128  MAST::GeomElem& elem) const {
129 
130  libmesh_assert(!_physics_elem);
131 
132  if (dim == 1) {
133 
134  const MAST::ElementPropertyCard1D& p =
135  dynamic_cast<const MAST::ElementPropertyCard1D&>(_discipline->get_property_card(ref_elem));
136 
137  elem.set_local_y_vector(p.y_vector());
138  }
139 }
140 
141 
142 void
144 init(const MAST::GeomElem& elem) {
145 
146  libmesh_assert(!_physics_elem);
147  libmesh_assert(_system);
148  libmesh_assert(_assembly);
149 
151  dynamic_cast<const MAST::ElementPropertyCardBase&>(_discipline->get_property_card(elem));
152 
153  _physics_elem =
155 }
156 
virtual ~HeatConductionTransientAssemblyElemOperations()
destructor resets the association of this assembly object with the system
const MAST::ElementPropertyCardBase & get_property_card(const libMesh::Elem &elem) const
get property card for the specified element
virtual void internal_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac)
internal force contribution to system residual
RealVectorX & y_vector()
returns value of the property val.
void set_local_y_vector(const RealVectorX &y_vec)
for 1D elements the transformed coordinate system attached to the element defines the local x-axis al...
Definition: geom_elem.cpp:119
virtual void elem_calculations(bool if_jac, RealVectorX &f_m, RealVectorX &f_x, RealMatrixX &f_m_jac_x_dot, RealMatrixX &f_m_jac, RealMatrixX &f_x_jac)
performs the element calculations over elem, and returns the element vector and matrix quantities in ...
MAST::PhysicsDisciplineBase * _discipline
virtual void elem_sensitivity_calculations(const MAST::FunctionBase &f, RealVectorX &f_m, RealVectorX &f_x)
performs the element sensitivity calculations over elem, and returns the component of element residua...
virtual void set_elem_data(unsigned int dim, const libMesh::Elem &ref_elem, MAST::GeomElem &elem) const
sets the structural element y-vector if 1D element is used.
virtual void internal_residual_sensitivity(const MAST::FunctionBase &p, RealVectorX &f)
sensitivity of the internal force contribution to system residual
virtual void linearized_jacobian_solution_product(RealVectorX &f)
Calculates the product of Jacobian-solution, and Jacobian-velocity over the element for a system of t...
virtual void elem_second_derivative_dot_solution_assembly(RealMatrixX &mat)
calculates over elem, and returns the matrix in vec .
void side_external_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, std::multimap< libMesh::boundary_id_type, MAST::BoundaryConditionBase *> &bc)
side external force contribution to system residual
const MAST::SideBCMapType & side_loads() const
Matrix< Real, Dynamic, Dynamic > RealMatrixX
virtual void init(const MAST::GeomElem &elem)
initializes the object for the geometric element elem.
HeatConductionTransientAssemblyElemOperations()
constructor associates this assembly object with the system
Matrix< Real, Dynamic, 1 > RealVectorX
void side_external_residual_sensitivity(const MAST::FunctionBase &p, RealVectorX &f, std::multimap< libMesh::boundary_id_type, MAST::BoundaryConditionBase *> &bc)
sensitivity of the side external force contribution to system residual
This class acts as a wrapper around libMesh::Elem for the purpose of providing a uniform interface fo...
Definition: geom_elem.h:59
const MAST::VolumeBCMapType & volume_loads() const
void volume_external_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac, std::multimap< libMesh::subdomain_id_type, MAST::BoundaryConditionBase *> &bc)
volume external force contribution to system residual
void volume_external_residual_sensitivity(const MAST::FunctionBase &p, RealVectorX &f, std::multimap< libMesh::subdomain_id_type, MAST::BoundaryConditionBase *> &bc)
sensitivity of the volume external force contribution to system residual
This element implements the Galerkin discretization of the heat conduction problem with the flux pro...
virtual void velocity_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac_xdot, RealMatrixX &jac)
inertial force contribution to system residual
virtual void velocity_residual_sensitivity(const MAST::FunctionBase &p, RealVectorX &f)
sensitivity of the damping force contribution to system residual
MAST::SystemInitialization * _system