MAST
eigenproblem_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 // MAST includes
22 #include "base/nonlinear_system.h"
24 #include "base/elem_base.h"
27 #include "mesh/geom_elem.h"
28 #include "numerics/utility.h"
29 
30 
31 // libMesh includes
32 #include "libmesh/numeric_vector.h"
33 #include "libmesh/sparse_matrix.h"
34 #include "libmesh/dof_map.h"
35 
36 
37 
40 _base_sol (nullptr),
41 _base_sol_sensitivity (nullptr) {
42 
43 }
44 
45 
46 
48 
49 }
50 
51 
52 
53 void
54 MAST::EigenproblemAssembly::set_base_solution(const libMesh::NumericVector<Real>& sol,
55  bool if_sens) {
56 
57  if (!if_sens) {
58 
59  // make sure the pointer has been cleared
60  libmesh_assert(!_base_sol);
61 
62  _base_sol = &sol;
63  }
64  else {
65 
66  // make sure the pointer has been cleared
67  libmesh_assert(!_base_sol_sensitivity);
68 
69  _base_sol_sensitivity = &sol;
70  }
71 }
72 
73 
74 
75 
76 void
78 
79  if (!if_sens)
80  _base_sol = nullptr;
81  else
82  _base_sol_sensitivity = nullptr;
83 }
84 
85 
86 
87 
88 void
90 eigenproblem_assemble(libMesh::SparseMatrix<Real>* A,
91  libMesh::SparseMatrix<Real>* B) {
92 
93  libmesh_assert(_system);
94  libmesh_assert(_discipline);
95  libmesh_assert(_elem_ops);
96 
97  MAST::NonlinearSystem& eigen_sys =
98  dynamic_cast<MAST::NonlinearSystem&>(_system->system());
99 
100  libMesh::SparseMatrix<Real>
101  &matrix_A = *A,
102  &matrix_B = *B;
103 
104  matrix_A.zero();
105  matrix_B.zero();
106 
107 
108  // build localized solutions if needed
109  std::unique_ptr<libMesh::NumericVector<Real> >
110  localized_solution;
111 
112  if (_base_sol)
113  localized_solution.reset(build_localized_vector(eigen_sys,
114  *_base_sol).release());
115 
116 
117  // iterate over each element, initialize it and get the relevant
118  // analysis quantities
119  RealVectorX sol;
120  RealMatrixX mat_A, mat_B;
121  std::vector<libMesh::dof_id_type> dof_indices;
122  const libMesh::DofMap& dof_map = eigen_sys.get_dof_map();
123 
124 
125  libMesh::MeshBase::const_element_iterator el =
126  eigen_sys.get_mesh().active_local_elements_begin();
127  const libMesh::MeshBase::const_element_iterator end_el =
128  eigen_sys.get_mesh().active_local_elements_end();
129 
131  &ops = dynamic_cast<MAST::EigenproblemAssemblyElemOperations&>(*_elem_ops);
132 
133  for ( ; el != end_el; ++el) {
134 
135  const libMesh::Elem* elem = *el;
136 
137  dof_map.dof_indices (elem, dof_indices);
138 
139  MAST::GeomElem geom_elem;
140  ops.set_elem_data(elem->dim(), *elem, geom_elem);
141  geom_elem.init(*elem, *_system);
142 
143  ops.init(geom_elem);
144 
145  // get the solution
146  unsigned int ndofs = (unsigned int)dof_indices.size();
147  sol.setZero(ndofs);
148  mat_A.setZero(ndofs, ndofs);
149  mat_B.setZero(ndofs, ndofs);
150 
151  // if the base solution is provided, then tell the element about it
152  if (_base_sol) {
153  for (unsigned int i=0; i<dof_indices.size(); i++)
154  sol(i) = (*localized_solution)(dof_indices[i]);
155  }
156 
157  ops.set_elem_solution(sol);
158  ops.elem_calculations(mat_A, mat_B);
159  ops.clear_elem();
160 
161  // copy to the libMesh matrix for further processing
162  DenseRealMatrix A, B;
163  MAST::copy(A, mat_A);
164  MAST::copy(B, mat_B);
165 
166  // constrain the element matrices.
167  dof_map.constrain_element_matrix(A, dof_indices);
168  dof_map.constrain_element_matrix(B, dof_indices);
169 
170  matrix_A.add_matrix (A, dof_indices); // load independent
171  matrix_B.add_matrix (B, dof_indices); // load dependent
172  }
173 
174  // finalize the data structures
175  A->close();
176  B->close();
177 }
178 
179 
180 
181 
182 
183 bool
186  libMesh::SparseMatrix<Real>* sensitivity_A,
187  libMesh::SparseMatrix<Real>* sensitivity_B) {
188 
189  libmesh_assert(_system);
190  libmesh_assert(_discipline);
191  libmesh_assert(_elem_ops);
192 
193  MAST::NonlinearSystem& eigen_sys =
194  dynamic_cast<MAST::NonlinearSystem&>(_system->system());
195 
196  libMesh::SparseMatrix<Real>& matrix_A = *sensitivity_A;
197  libMesh::SparseMatrix<Real>& matrix_B = *sensitivity_B;
198 
199  matrix_A.zero();
200  matrix_B.zero();
201 
202  // build localized solutions if needed
203  std::unique_ptr<libMesh::NumericVector<Real> >
204  localized_solution,
205  localized_solution_sens;
206 
207  if (_base_sol) {
208 
209  localized_solution.reset(build_localized_vector(eigen_sys,
210  *_base_sol).release());
211 
212  // make sure that the sensitivity was also provided
213  libmesh_assert(_base_sol_sensitivity);
214  localized_solution_sens.reset(build_localized_vector(eigen_sys,
215  *_base_sol_sensitivity).release());
216  }
217 
218 
219  // iterate over each element, initialize it and get the relevant
220  // analysis quantities
221  RealVectorX sol;
222  RealMatrixX mat_A, mat_B;
223  std::vector<libMesh::dof_id_type> dof_indices;
224  const libMesh::DofMap& dof_map = eigen_sys.get_dof_map();
225 
226 
227  libMesh::MeshBase::const_element_iterator el =
228  eigen_sys.get_mesh().active_local_elements_begin();
229  const libMesh::MeshBase::const_element_iterator end_el =
230  eigen_sys.get_mesh().active_local_elements_end();
231 
233  &ops = dynamic_cast<MAST::EigenproblemAssemblyElemOperations&>(*_elem_ops);
234 
235  for ( ; el != end_el; ++el) {
236 
237  const libMesh::Elem* elem = *el;
238 
239  dof_map.dof_indices (elem, dof_indices);
240 
241  MAST::GeomElem geom_elem;
242  ops.set_elem_data(elem->dim(), *elem, geom_elem);
243  geom_elem.init(*elem, *_system);
244 
245  ops.init(geom_elem);
246 
247  // get the solution
248  unsigned int ndofs = (unsigned int)dof_indices.size();
249  sol.setZero(ndofs);
250  mat_A.setZero(ndofs, ndofs);
251  mat_B.setZero(ndofs, ndofs);
252 
253  // if the base solution is provided, tell the element about it
254  if (_base_sol) {
255 
256  // set the element's base solution
257  for (unsigned int i=0; i<dof_indices.size(); i++)
258  sol(i) = (*localized_solution)(dof_indices[i]);
259  }
260 
261  ops.set_elem_solution(sol);
262 
263  // set the element's base solution sensitivity
264  if (_base_sol) {
265 
266  for (unsigned int i=0; i<dof_indices.size(); i++)
267  sol(i) = (*localized_solution_sens)(dof_indices[i]);
268  }
269 
272  _base_sol!=nullptr,
273  mat_A,
274  mat_B);
275  ops.clear_elem();
276 
277  // copy to the libMesh matrix for further processing
278  DenseRealMatrix A, B;
279  MAST::copy(A, mat_A);
280  MAST::copy(B, mat_B);
281 
282  // constrain the element matrices.
283  dof_map.constrain_element_matrix(A, dof_indices);
284  dof_map.constrain_element_matrix(B, dof_indices);
285 
286  matrix_A.add_matrix (A, dof_indices);
287  matrix_B.add_matrix (B, dof_indices);
288  }
289 
290  // finalize the data structures
291  sensitivity_A->close();
292  sensitivity_B->close();
293 
294  return true;
295 }
296 
297 
298 
MAST::AssemblyElemOperations * _elem_ops
provides assembly elem operations for use by this class
MAST::NonlinearSystem & system()
virtual void eigenproblem_assemble(libMesh::SparseMatrix< Real > *A, libMesh::SparseMatrix< Real > *B)
assembles the matrices for eigenproblem depending on the analysis type
libMesh::DenseMatrix< Real > DenseRealMatrix
This class implements a system for solution of nonlinear systems.
virtual ~EigenproblemAssembly()
destructor resets the association with the eigen system from this assembly object ...
void set_base_solution(const libMesh::NumericVector< Real > &sol, bool if_sens=false)
if the eigenproblem is defined about a non-zero base solution, then this method provides the object w...
virtual void set_elem_solution_sensitivity(const RealVectorX &sol)
sets the element solution sensitivity
virtual void elem_calculations(RealMatrixX &mat_A, RealMatrixX &mat_B)=0
performs the element calculations over elem, and returns the element matrices for the eigenproblem ...
virtual void set_elem_solution(const RealVectorX &sol)
sets the element solution
MAST::SystemInitialization * _system
System for which this assembly is performed.
MAST::PhysicsDisciplineBase * _discipline
PhysicsDisciplineBase object for which this class is assembling.
virtual void set_elem_data(unsigned int dim, const libMesh::Elem &ref_elem, MAST::GeomElem &elem) const =0
some analyses may want to set additional element data before initialization of the GeomElem...
void clear_base_solution(bool if_sens=false)
Clears the pointer to the solution.
virtual bool eigenproblem_sensitivity_assemble(const MAST::FunctionBase &f, libMesh::SparseMatrix< Real > *sensitivity_A, libMesh::SparseMatrix< Real > *sensitivity_B)
Assembly function.
Matrix< Real, Dynamic, Dynamic > RealMatrixX
const libMesh::NumericVector< Real > * _base_sol
base solution about which this eigenproblem is defined.
void copy(DenseRealMatrix &m1, const RealMatrixX &m2)
Definition: utility.h:167
Matrix< Real, Dynamic, 1 > RealVectorX
virtual void init(const MAST::GeomElem &elem)=0
initializes the object for calculation of element quantities for the specified elem.
std::unique_ptr< libMesh::NumericVector< Real > > build_localized_vector(const libMesh::System &sys, const libMesh::NumericVector< Real > &global) const
localizes the parallel vector so that the local copy stores all values necessary for calculation of t...
This class acts as a wrapper around libMesh::Elem for the purpose of providing a uniform interface fo...
Definition: geom_elem.h:59
const libMesh::NumericVector< Real > * _base_sol_sensitivity
sensitivity of base solution may be needed for sensitivity analysis.
virtual void elem_sensitivity_calculations(const MAST::FunctionBase &f, bool base_sol, RealMatrixX &mat_A, RealMatrixX &mat_B)=0
performs the element sensitivity calculations over elem, and returns the element matrices for the eig...
virtual void clear_elem()
clears the element initialization
virtual void init(const libMesh::Elem &elem, const MAST::SystemInitialization &sys_init)
initialize the object for the specified reference elem.
Definition: geom_elem.cpp:128
EigenproblemAssembly()
constructor associates the eigen system with this assembly object