MAST
interface_dof_handler.h
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 #ifndef __mast_level_set_interface_dof_handler_h__
21 #define __mast_level_set_interface_dof_handler_h__
22 
23 
24 // MAST includes
25 #include "base/mast_data_types.h"
26 
27 
28 // libMesh includes
29 #include "libmesh/elem.h"
30 #include "libmesh/numeric_vector.h"
31 
32 
33 namespace MAST {
34 
35  // Forward declerations
36  class SystemInitialization;
37  template <typename ValType> class FieldFunction;
38  class LevelSetIntersection;
39  class NonlinearImplicitAssemblyElemOperations;
40 
42 
43  public:
44 
46 
48 
49  void init(const MAST::SystemInitialization& sys_init,
50  MAST::LevelSetIntersection& intersection,
52 
55  libmesh_assert(_phi);
56  return *_phi;
57  }
58 
64  bool if_factor_element(const libMesh::Elem& elem) const;
65 
74  void partition_local_elem_rows(const libMesh::Elem& elem,
75  std::vector<libMesh::dof_id_type>& material_dofs,
76  std::vector<libMesh::dof_id_type>& void_dofs);
77 
78 
83  void partition_global_elem_rows(const libMesh::Elem& elem,
84  std::vector<libMesh::dof_id_type>& material_dofs,
85  std::vector<libMesh::dof_id_type>& void_dofs);
86 
87 
92  void solution_of_factored_element(const libMesh::Elem& elem,
93  RealVectorX& elem_sol);
94 
95 
104  void element_factored_jacobian(const libMesh::Elem& elem,
105  const RealMatrixX& jac,
106  std::vector<libMesh::dof_id_type>& material_dof_ids,
107  RealMatrixX& jac_factored_uu);
108 
118  void
119  element_factored_residual_and_jacobian(const libMesh::Elem& elem,
120  const RealMatrixX& jac,
121  const RealVectorX& res,
122  std::vector<libMesh::dof_id_type>& material_dof_ids,
123  RealMatrixX& jac_factored_uu,
124  RealVectorX& res_factored_u);
125 
133  void element_factored_jacobian(const libMesh::Elem& elem,
134  const RealMatrixX& jac,
135  std::vector<libMesh::dof_id_type>& material_dof_ids,
136  std::vector<libMesh::dof_id_type>& void_dof_ids,
137  RealMatrixX& jac_uu,
138  RealMatrixX& jac_uf,
139  RealMatrixX& jac_fu,
140  RealMatrixX& jac_ff,
141  RealMatrixX& jac_factored_uu);
142 
143  void update_factored_element_solution(const libMesh::Elem& elem,
144  const RealMatrixX& res,
145  const RealMatrixX& jac,
146  const RealMatrixX& sol,
147  const RealMatrixX& dsol,
148  RealVectorX& updated_sol);
149 
150 
151  protected:
152 
153  void _compute_matrix_inverse(const RealMatrixX& mat,
154  RealMatrixX& mat_inv);
155 
158 
159  std::map<const libMesh::Elem*, RealVectorX> _elem_sol;
160 
165  std::map<const libMesh::Elem*, std::set<const libMesh::Node*>> _elem_void_nodes;
166 
170  std::map<const libMesh::Node*, std::set<const libMesh::Elem*>> _void_node_elems;
171 
175  std::map<const libMesh::Elem*, std::map<libMesh::dof_id_type, libMesh::dof_id_type>> _dof_ids;
176  };
177 }
178 
179 
180 #endif // __mast_level_set_interface_dof_handler_h__
void element_factored_jacobian(const libMesh::Elem &elem, const RealMatrixX &jac, std::vector< libMesh::dof_id_type > &material_dof_ids, RealMatrixX &jac_factored_uu)
a wrapper around the second element_factored_jacobian.
void update_factored_element_solution(const libMesh::Elem &elem, const RealMatrixX &res, const RealMatrixX &jac, const RealMatrixX &sol, const RealMatrixX &dsol, RealVectorX &updated_sol)
std::map< const libMesh::Elem *, RealVectorX > _elem_sol
void solution_of_factored_element(const libMesh::Elem &elem, RealVectorX &elem_sol)
updates the components of the solution vector in elem_sol for the void domain using the stored soluti...
std::map< const libMesh::Node *, std::set< const libMesh::Elem * > > _void_node_elems
map of elements that each node will provide a new dof for
void _compute_matrix_inverse(const RealMatrixX &mat, RealMatrixX &mat_inv)
MAST::FieldFunction< Real > * _phi
Matrix< Real, Dynamic, Dynamic > RealMatrixX
std::map< const libMesh::Elem *, std::set< const libMesh::Node * > > _elem_void_nodes
map of nodes on each element that will be added as independent dofs.
void element_factored_residual_and_jacobian(const libMesh::Elem &elem, const RealMatrixX &jac, const RealVectorX &res, std::vector< libMesh::dof_id_type > &material_dof_ids, RealMatrixX &jac_factored_uu, RealVectorX &res_factored_u)
factorizes the residual and jacobian into the components for the dofs on material nodes...
Matrix< Real, Dynamic, 1 > RealVectorX
void partition_global_elem_rows(const libMesh::Elem &elem, std::vector< libMesh::dof_id_type > &material_dofs, std::vector< libMesh::dof_id_type > &void_dofs)
fills the material_dofs and void_dofs with the dofs_ids in the global system corresponding to these d...
void init(const MAST::SystemInitialization &sys_init, MAST::LevelSetIntersection &intersection, MAST::FieldFunction< Real > &phi)
std::map< const libMesh::Elem *, std::map< libMesh::dof_id_type, libMesh::dof_id_type > > _dof_ids
new dof ids for each elem/old_dof pair.
const MAST::SystemInitialization * _sys_init
bool if_factor_element(const libMesh::Elem &elem) const
void partition_local_elem_rows(const libMesh::Elem &elem, std::vector< libMesh::dof_id_type > &material_dofs, std::vector< libMesh::dof_id_type > &void_dofs)
identifies which rows in the element residual vector and rows/columns in the jacobian matrix correspo...
MAST::FieldFunction< Real > & get_level_set_function()