MAST
assembly_base.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
21 #include "base/assembly_base.h"
24 #include "base/elem_base.h"
26 #include "base/nonlinear_system.h"
29 #include "mesh/fe_base.h"
30 #include "mesh/geom_elem.h"
31 #include "numerics/utility.h"
32 
33 
34 // libMesh includes
35 #include "libmesh/numeric_vector.h"
36 #include "libmesh/dof_map.h"
37 
38 
40 close_matrix (true),
41 _elem_ops (nullptr),
42 _discipline (nullptr),
43 _system (nullptr),
44 _sol_function (nullptr),
45 _solver_monitor (nullptr),
46 _param_dependence (nullptr) {
47 
48 }
49 
50 
51 
52 
54 
57 }
58 
59 
60 
63 
64  libmesh_assert_msg(_discipline,
65  "Error: Discipline not yet attached to Assembly.");
66  return *_discipline;
67 }
68 
69 
70 
73 
74  libmesh_assert_msg(_discipline,
75  "Error: Discipline not yet attached to Assembly.");
76  return *_discipline;
77 }
78 
79 
80 
81 
84 
85  libmesh_assert_msg(_discipline,
86  "Error: System not yet attached to Assembly.");
87  return _system->system();
88 }
89 
90 
93 
94  libmesh_assert_msg(_discipline,
95  "Error: System not yet attached to Assembly.");
96  return _system->system();
97 }
98 
99 
102 
103  libmesh_assert_msg(_elem_ops,
104  "Error: Not yet initialized.");
105  return *_elem_ops;
106 }
107 
108 
111 
112  libmesh_assert_msg(_system,
113  "Error: System not yet attached to Assembly.");
114  return *_system;
115 }
116 
117 
118 void
120 
121  libmesh_assert(!_solver_monitor);
122  _solver_monitor = &monitor;
123 }
124 
125 
126 
127 void
130 
131  libmesh_assert(!_param_dependence);
132 
133  _param_dependence = &dep;
134 }
135 
136 
137 
138 void
140 
141  _param_dependence = nullptr;
142 }
143 
144 
145 
148 
149  return _solver_monitor;
150 }
151 
152 
153 void
155 
156  _solver_monitor = nullptr;
157 }
158 
159 
160 void
164 
165  libmesh_assert_msg(!_discipline && !_system,
166  "Error: Assembly should be cleared before attaching System.");
167 
169  _system = &system;
170 }
171 
172 
173 
174 void
177 
178  close_matrix = true;
179  _discipline = nullptr;
180  _system = nullptr;
181  _param_dependence = nullptr;
182 }
183 
184 
185 
186 void
188 
189  libmesh_assert_msg(!_elem_ops,
190  "Error: Assembly should be cleared before attaching Elem.");
191 
192  _elem_ops = &elem_ops;
193  _elem_ops->set_assembly(*this);
194 }
195 
196 
197 
198 void
200 
201  if (_elem_ops) {
203  _elem_ops = nullptr;
204  }
205 }
206 
207 
208 
209 std::unique_ptr<libMesh::NumericVector<Real> >
211 build_localized_vector(const libMesh::System& sys,
212  const libMesh::NumericVector<Real>& global) const {
213 
214  libMesh::NumericVector<Real>* local =
215  libMesh::NumericVector<Real>::build(sys.comm()).release();
216 
217  const std::vector<libMesh::dof_id_type>& send_list =
218  sys.get_dof_map().get_send_list();
219 
220  local->init(sys.n_dofs(),
221  sys.n_local_dofs(),
222  send_list,
223  false,
224  libMesh::GHOSTED);
225  global.localize(*local, send_list);
226 
227  return std::unique_ptr<libMesh::NumericVector<Real> >(local);
228 }
229 
230 
231 
232 
233 void
235 
236  // make sure that no prior association is specified
237  libmesh_assert(!_sol_function);
238 
239  _sol_function = &f;
240 }
241 
242 
243 
244 
245 void
247  _sol_function = nullptr;
248 }
249 
250 
251 
252 
253 void
254 MAST::AssemblyBase::calculate_output(const libMesh::NumericVector<Real>& X,
256 
257  libmesh_assert(_discipline);
258  libmesh_assert(_system);
259 
260  MAST::NonlinearSystem& nonlin_sys = _system->system();
261  output.set_assembly(*this);
262 
263  output.zero_for_analysis();
264 
265  // iterate over each element, initialize it and get the relevant
266  // analysis quantities
267  RealVectorX sol;
268 
269  std::vector<libMesh::dof_id_type> dof_indices;
270  const libMesh::DofMap& dof_map = _system->system().get_dof_map();
271 
272  std::unique_ptr<libMesh::NumericVector<Real> > localized_solution;
273  localized_solution.reset(build_localized_vector(nonlin_sys,
274  X).release());
275 
276 
277  // if a solution function is attached, initialize it
278  //if (_sol_function)
279  // _sol_function->init( X);
280 
281  libMesh::MeshBase::const_element_iterator el =
282  nonlin_sys.get_mesh().active_local_elements_begin();
283  const libMesh::MeshBase::const_element_iterator end_el =
284  nonlin_sys.get_mesh().active_local_elements_end();
285 
286 
287  for ( ; el != end_el; ++el) {
288 
289  const libMesh::Elem* elem = *el;
290 
291  dof_map.dof_indices (elem, dof_indices);
292 
293  // get the solution
294  unsigned int ndofs = (unsigned int)dof_indices.size();
295  sol.setZero(ndofs);
296 
297  for (unsigned int i=0; i<dof_indices.size(); i++)
298  sol(i) = (*localized_solution)(dof_indices[i]);
299 
300 
301  //if (_sol_function)
302  // physics_elem->attach_active_solution_function(*_sol_function);
303 
304  MAST::GeomElem geom_elem;
305  output.set_elem_data(elem->dim(), *elem, geom_elem);
306  geom_elem.init(*elem, *_system);
307 
308  output.init(geom_elem);
309  output.set_elem_solution(sol);
310  output.evaluate();
311  output.clear_elem();
312 
313  //physics_elem->detach_active_solution_function();
314  }
315 
316  // if a solution function is attached, clear it
317  if (_sol_function)
318  _sol_function->clear();
319 
320  output.clear_assembly();
321 }
322 
323 
324 
325 
326 void
328 calculate_output_derivative(const libMesh::NumericVector<Real>& X,
330  libMesh::NumericVector<Real>& dq_dX) {
331 
332  libmesh_assert(_discipline);
333  libmesh_assert(_system);
334 
335  output.zero_for_sensitivity();
336 
337  MAST::NonlinearSystem& nonlin_sys = _system->system();
338  output.set_assembly(*this);
339 
340  dq_dX.zero();
341 
342  // iterate over each element, initialize it and get the relevant
343  // analysis quantities
344  RealVectorX vec, sol;
345 
346  std::vector<libMesh::dof_id_type> dof_indices;
347  const libMesh::DofMap& dof_map = _system->system().get_dof_map();
348 
349 
350  std::unique_ptr<libMesh::NumericVector<Real> > localized_solution;
351  localized_solution.reset(build_localized_vector(nonlin_sys,
352  X).release());
353 
354 
355  // if a solution function is attached, initialize it
356  if (_sol_function)
357  _sol_function->init( X);
358 
359 
360  libMesh::MeshBase::const_element_iterator el =
361  nonlin_sys.get_mesh().active_local_elements_begin();
362  const libMesh::MeshBase::const_element_iterator end_el =
363  nonlin_sys.get_mesh().active_local_elements_end();
364 
365 
366  for ( ; el != end_el; ++el) {
367 
368  const libMesh::Elem* elem = *el;
369 
370  dof_map.dof_indices (elem, dof_indices);
371 
372  // get the solution
373  unsigned int ndofs = (unsigned int)dof_indices.size();
374  sol.setZero(ndofs);
375  vec.setZero(ndofs);
376 
377  for (unsigned int i=0; i<dof_indices.size(); i++)
378  sol(i) = (*localized_solution)(dof_indices[i]);
379 
380 // if (_sol_function)
381 // physics_elem->attach_active_solution_function(*_sol_function);
382 
383  MAST::GeomElem geom_elem;
384  output.set_elem_data(elem->dim(), *elem, geom_elem);
385  geom_elem.init(*elem, *_system);
386 
387  output.init(geom_elem);
388  output.set_elem_solution(sol);
389  output.output_derivative_for_elem(vec);
390  output.clear_elem();
391 
392  DenseRealVector v;
393  MAST::copy(v, vec);
394  dof_map.constrain_element_vector(v, dof_indices);
395  dq_dX.add_vector(v, dof_indices);
396  dof_indices.clear();
397  }
398 
399  // if a solution function is attached, clear it
400  if (_sol_function)
401  _sol_function->clear();
402 
403  dq_dX.close();
404 
405  output.clear_assembly();
406 }
407 
408 
409 
410 void
412 calculate_output_direct_sensitivity(const libMesh::NumericVector<Real>& X,
413  const libMesh::NumericVector<Real>* dXdp,
414  const MAST::FunctionBase& p,
416 
417  libmesh_assert(_discipline);
418  libmesh_assert(_system);
419 
420  output.zero_for_sensitivity();
421 
422  MAST::NonlinearSystem& nonlin_sys = _system->system();
423  output.set_assembly(*this);
424 
425  // iterate over each element, initialize it and get the relevant
426  // analysis quantities
428  sol,
429  dsol;
430 
431  std::vector<libMesh::dof_id_type> dof_indices;
432  const libMesh::DofMap& dof_map = _system->system().get_dof_map();
433 
434 
435  std::unique_ptr<libMesh::NumericVector<Real> >
436  localized_solution,
437  localized_solution_sens;
438  localized_solution.reset(build_localized_vector(nonlin_sys,
439  X).release());
440  if (dXdp)
441  localized_solution_sens.reset(build_localized_vector(nonlin_sys,
442  *dXdp).release());
443 
444 
445  // if a solution function is attached, initialize it
446  if (_sol_function)
447  _sol_function->init( X);
448 
449 
450  libMesh::MeshBase::const_element_iterator el =
451  nonlin_sys.get_mesh().active_local_elements_begin();
452  const libMesh::MeshBase::const_element_iterator end_el =
453  nonlin_sys.get_mesh().active_local_elements_end();
454 
455 
456  for ( ; el != end_el; ++el) {
457 
458  const libMesh::Elem* elem = *el;
459 
460  // no sensitivity computation assembly is neeed in these cases
461  if (_param_dependence &&
462  // if object is specified and elem does not depend on it
464  // and if no sol_sens is given
465  (!dXdp ||
466  // or if it can be ignored for elem
467  (dXdp && _param_dependence->override_flag)))
468  continue;
469 
470  dof_map.dof_indices (elem, dof_indices);
471 
472  // get the solution
473  unsigned int ndofs = (unsigned int)dof_indices.size();
474  sol.setZero(ndofs);
475  dsol.setZero(ndofs);
476 
477  for (unsigned int i=0; i<dof_indices.size(); i++) {
478  sol(i) = (*localized_solution)(dof_indices[i]);
479  if (dXdp)
480  dsol(i) = (*localized_solution_sens)(dof_indices[i]);
481  }
482 
483  // if (_sol_function)
484  // physics_elem->attach_active_solution_function(*_sol_function);
485 
486  MAST::GeomElem geom_elem;
487  output.set_elem_data(elem->dim(), *elem, geom_elem);
488  geom_elem.init(*elem, *_system);
489 
490  output.init(geom_elem);
491  output.set_elem_solution(sol);
492  output.set_elem_solution_sensitivity(dsol);
493  output.evaluate_sensitivity(p);
494  output.clear_elem();
495 
496  // physics_elem->detach_active_solution_function();
497  }
498 
499  // if a solution function is attached, clear it
500  if (_sol_function)
501  _sol_function->clear();
502 
503  output.clear_assembly();
504 }
505 
506 
507 
508 
509 Real
511 calculate_output_adjoint_sensitivity(const libMesh::NumericVector<Real>& X,
512  const libMesh::NumericVector<Real>& dq_dX,
513  const MAST::FunctionBase& p,
516  const bool include_partial_sens) {
517 
518  libmesh_assert(_discipline);
519  libmesh_assert(_system);
520 
521  MAST::NonlinearSystem& nonlin_sys = _system->system();
522 
523  libMesh::NumericVector<Real>
524  &dres_dp = nonlin_sys.add_sensitivity_rhs();
525 
526  this->set_elem_operation_object(elem_ops);
527  this->sensitivity_assemble(p, dres_dp);
529 
530  Real
531  dq_dp = dq_dX.dot(dres_dp);
532 
533  if (include_partial_sens) {
534 
535  // calculate the partial sensitivity of the output, which is done
536  // with zero solution vector
537  this->calculate_output_direct_sensitivity(X, nullptr, p, output);
538 
539  dq_dp += output.output_sensitivity_total(p);
540  }
541 
542  return dq_dp;
543 }
544 
MAST::AssemblyElemOperations * _elem_ops
provides assembly elem operations for use by this class
MAST::NonlinearSystem & system()
virtual bool if_elem_depends_on_parameter(const libMesh::Elem &e, const MAST::FunctionBase &p) const =0
AssemblyBase()
constructor takes a reference to the discipline that provides the boundary conditions, volume loads, properties, etc.
void clear_solver_monitor()
clears the monitor object
This class implements a system for solution of nonlinear systems.
bool override_flag
if true, assume zero solution sensitivity when elem does not dependent on parameter.
Definition: assembly_base.h:97
virtual void clear_elem_operation_object()
clears the association of this object with the assembly element operation object. ...
virtual bool sensitivity_assemble(const MAST::FunctionBase &f, libMesh::NumericVector< Real > &sensitivity_rhs)
Assembly function.
This provides a wrapper FieldFunction compatible class that interpolates the solution using libMesh&#39;s...
void attach_solution_function(MAST::MeshFieldFunction &f)
tells the assembly object that this function is will need to be initialized before each residual eval...
virtual void calculate_output_derivative(const libMesh::NumericVector< Real > &X, MAST::OutputAssemblyElemOperations &output, libMesh::NumericVector< Real > &dq_dX)
calculates
void attach_elem_parameter_dependence_object(MAST::AssemblyBase::ElemParameterDependence &dep)
This object, if provided by user, will be used to reduce unnecessary computations in sensitivity anal...
virtual Real calculate_output_adjoint_sensitivity(const libMesh::NumericVector< Real > &X, const libMesh::NumericVector< Real > &dq_dX, const MAST::FunctionBase &p, MAST::AssemblyElemOperations &elem_ops, MAST::OutputAssemblyElemOperations &output, const bool include_partial_sens=true)
Evaluates the total sensitivity of output wrt p using the adjoint solution provided in dq_dX for a li...
virtual void set_elem_solution_sensitivity(const RealVectorX &sol)
sets the element solution sensitivity
void set_solver_monitor(MAST::AssemblyBase::SolverMonitor &monitor)
attaches the solver monitor, which is a user provided routine that is called each time ...
virtual void set_assembly(MAST::AssemblyBase &assembly)
sets the assembly object
This provides the base class for definitin of element level contribution of output quantity in an ana...
virtual void zero_for_sensitivity()=0
zeroes the output quantity values stored inside this object so that assembly process can begin...
bool close_matrix
flag to control the closing fo the Jacobian after assembly
virtual void set_elem_operation_object(MAST::AssemblyElemOperations &elem_ops)
attaches a element operation to this object, and associated this with the element operation object...
virtual void set_elem_solution(const RealVectorX &sol)
sets the element solution
libMesh::Real Real
virtual void evaluate()=0
this is the abstract interface to be implemented by derived classes.
MAST::SystemInitialization * _system
System for which this assembly is performed.
void detach_solution_function()
removes the attachment of the solution function
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...
virtual void clear_discipline_and_system()
clears association with a system to this discipline
virtual void calculate_output(const libMesh::NumericVector< Real > &X, MAST::OutputAssemblyElemOperations &output)
calculates the value of quantity .
libMesh::DenseVector< Real > DenseRealVector
virtual void clear_assembly()
clears the assembly object
void copy(DenseRealMatrix &m1, const RealMatrixX &m2)
Definition: utility.h:167
void clear()
clear the solution
virtual void calculate_output_direct_sensitivity(const libMesh::NumericVector< Real > &X, const libMesh::NumericVector< Real > *dXdp, const MAST::FunctionBase &p, MAST::OutputAssemblyElemOperations &output)
evaluates the sensitivity of the outputs in the attached discipline with respect to the parametrs in ...
Matrix< Real, Dynamic, 1 > RealVectorX
virtual Real output_sensitivity_total(const MAST::FunctionBase &p)=0
virtual void init(const MAST::GeomElem &elem)=0
initializes the object for calculation of element quantities for the specified elem.
virtual void output_derivative_for_elem(RealVectorX &dq_dX)=0
returns the output quantity derivative with respect to state vector in dq_dX.
void clear_elem_parameter_dependence_object()
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...
MAST::AssemblyBase::SolverMonitor * get_solver_monitor()
This class acts as a wrapper around libMesh::Elem for the purpose of providing a uniform interface fo...
Definition: geom_elem.h:59
Inherited objects from this class can be provided by the user provide assessment of whether or not an...
Definition: assembly_base.h:86
MAST::SystemInitialization & system_init()
virtual void zero_for_analysis()=0
zeroes the output quantity values stored inside this object so that assembly process can begin...
void init(const libMesh::NumericVector< Real > &sol, const libMesh::NumericVector< Real > *dsol=nullptr)
initializes the data structures to perform the interpolation function of sol.
MAST::AssemblyBase::ElemParameterDependence * _param_dependence
If provided by user, this object is used by sensitiivty analysis to check for whether or the current ...
virtual void set_discipline_and_system(MAST::PhysicsDisciplineBase &discipline, MAST::SystemInitialization &system)
attaches a system to this discipline
virtual void clear_elem()
clears the element initialization
virtual ~AssemblyBase()
virtual destructor
MAST::AssemblyElemOperations & get_elem_ops()
MAST::AssemblyBase::SolverMonitor * _solver_monitor
User provided solver monitor is attached to the linear nonlinear solvers, if provided.
const MAST::NonlinearSystem & system() const
virtual void evaluate_sensitivity(const MAST::FunctionBase &f)=0
this evaluates all relevant sensitivity components on the element.
MAST::MeshFieldFunction * _sol_function
system solution that will be initialized before each solution
const MAST::PhysicsDisciplineBase & discipline() const