MAST
heat_conduction_nonlinear_assembly.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
23 #include "base/assembly_base.h"
28 
29 
33 
34 }
35 
36 
37 
40 
41 }
42 
43 
44 void
46 set_elem_data(unsigned int dim,
47  const libMesh::Elem& ref_elem,
48  MAST::GeomElem& elem) const {
49 
50  libmesh_assert(!_physics_elem);
51 
52  if (dim == 1) {
53 
55  dynamic_cast<const MAST::ElementPropertyCard1D&>(_discipline->get_property_card(ref_elem));
56 
57  elem.set_local_y_vector(p.y_vector());
58  }
59 }
60 
61 
62 void
64 init(const MAST::GeomElem& elem) {
65 
66  libmesh_assert(!_physics_elem);
67  libmesh_assert(_system);
68  libmesh_assert(_assembly);
69 
71  dynamic_cast<const MAST::ElementPropertyCardBase&>
73 
76 }
77 
78 
79 
80 
81 void
83 elem_calculations(bool if_jac,
84  RealVectorX& vec,
85  RealMatrixX& mat) {
86 
87  libmesh_assert(_physics_elem);
88 
90  dynamic_cast<MAST::HeatConductionElementBase&>(*_physics_elem);
91 
92  vec.setZero();
93  mat.setZero();
94 
95  e.internal_residual(if_jac, vec, mat);
96  e.side_external_residual(if_jac, vec, mat, _discipline->side_loads());
97  e.volume_external_residual(if_jac, vec, mat, _discipline->volume_loads());
98 }
99 
100 
101 
102 
103 void
106  RealVectorX& vec) {
107 
108  libmesh_assert(_physics_elem);
109 
111  dynamic_cast<MAST::HeatConductionElementBase&>(*_physics_elem);
112 
113  vec.setZero();
114 
118 }
119 
120 
121 void
125  RealVectorX& vec) {
126 
127  libmesh_assert(_physics_elem);
128  libmesh_assert(f.is_topology_parameter());
129 
131  &elem = dynamic_cast<const MAST::LevelSetIntersectedElem&>(_physics_elem->elem());
132 
133  // sensitivity only exists at the boundary. So, we proceed with calculation
134  // only if this element has an intersection in the interior, or with a side.
135  if (elem.if_elem_has_level_set_boundary() &&
136  elem.if_subelem_has_side_on_level_set_boundary()) {
137 
139  dynamic_cast<MAST::HeatConductionElementBase&>(*_physics_elem);
140 
141  vec.setZero();
143  dummy = RealMatrixX::Zero(vec.size(), vec.size());
144 
146  elem.get_subelem_side_on_level_set_boundary(),
147  vel);
149  elem.get_subelem_side_on_level_set_boundary(),
150  vel,
152  /*e.side_external_residual_sensitivity(f, false,
153  vec,
154  dummy,
155  dummy,
156  _discipline->side_loads());*/
157  }
158 }
159 
160 
161 void
164 
165  libmesh_error(); // to be implemented
166 }
167 
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::VolumeBCMapType & volume_loads() const
virtual void elem_second_derivative_dot_solution_assembly(RealMatrixX &mat)
calculates over elem, and returns the matrix in vec .
virtual void elem_sensitivity_calculations(const MAST::FunctionBase &f, RealVectorX &vec)
performs the element sensitivity calculations over elem, and returns the element residual sensitivity...
virtual void internal_residual(bool request_jacobian, RealVectorX &f, RealMatrixX &jac)
internal force contribution to system residual
virtual void elem_topology_sensitivity_calculations(const MAST::FunctionBase &f, const MAST::FieldFunction< RealVectorX > &vel, RealVectorX &vec)
performs the element topology sensitivity calculations over elem, and returns the element residual se...
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
RealVectorX & y_vector()
returns value of the property val.
const MAST::GeomElem & elem() const
Definition: elem_base.h:117
virtual void init(const MAST::GeomElem &elem)
initializes the object for the geometric element elem.
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
void volume_external_residual_boundary_velocity(const MAST::FunctionBase &p, RealVectorX &f, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f, std::multimap< libMesh::subdomain_id_type, MAST::BoundaryConditionBase * > &bc)
boundary velocity contribution of volume external force.
MAST::PhysicsDisciplineBase * _discipline
HeatConductionNonlinearAssemblyElemOperations()
constructor associates this assembly object with the system
virtual ~HeatConductionNonlinearAssemblyElemOperations()
destructor resets the association of this assembly object with the system
virtual void elem_calculations(bool if_jac, RealVectorX &vec, RealMatrixX &mat)
performs the element calculations over elem, and returns the element vector and matrix quantities in ...
virtual void internal_residual_sensitivity(const MAST::FunctionBase &p, RealVectorX &f)
sensitivity of the internal force contribution to system residual
Matrix< Real, Dynamic, Dynamic > RealMatrixX
This class inherits from MAST::GeomElem and provides an interface to initialize FE objects on sub-ele...
Matrix< Real, Dynamic, 1 > RealVectorX
const MAST::ElementPropertyCardBase & get_property_card(const libMesh::Elem &elem) const
get property card for the specified element
This class acts as a wrapper around libMesh::Elem for the purpose of providing a uniform interface fo...
Definition: geom_elem.h:59
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.
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
This element implements the Galerkin discretization of the heat conduction problem with the flux pro...
const MAST::SideBCMapType & side_loads() const
virtual void internal_residual_boundary_velocity(const MAST::FunctionBase &p, RealVectorX &f, const unsigned int s, const MAST::FieldFunction< RealVectorX > &vel_f)
sensitivity of the internal force contribution to system residual
virtual bool is_topology_parameter() const
Definition: function_base.h:97
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
MAST::SystemInitialization * _system