libMesh
optimization_ex2.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 // <h1>Optimization Example 2 - Optimization with constraints</h1>
21 // \author David Knezevic
22 // \date 2015
23 //
24 // In this example we extend example 1 to demonstrate how to use
25 // OptimizationSystem's interface for imposing equality and inequality
26 // constraints.
27 
28 // C++ include files that we need
29 #include <iostream>
30 
31 // libMesh includes
32 #include "libmesh/libmesh.h"
33 #include "libmesh/mesh.h"
34 #include "libmesh/mesh_generation.h"
35 #include "libmesh/exodusII_io.h"
36 #include "libmesh/equation_systems.h"
37 #include "libmesh/fe.h"
38 #include "libmesh/quadrature_gauss.h"
39 #include "libmesh/dof_map.h"
40 #include "libmesh/sparse_matrix.h"
41 #include "libmesh/numeric_vector.h"
42 #include "libmesh/dense_matrix.h"
43 #include "libmesh/dense_vector.h"
44 #include "libmesh/elem.h"
45 #include "libmesh/boundary_info.h"
46 #include "libmesh/string_to_enum.h"
47 #include "libmesh/getpot.h"
48 #include "libmesh/zero_function.h"
49 #include "libmesh/dirichlet_boundaries.h"
50 #include "libmesh/optimization_system.h"
51 #include "libmesh/optimization_solver.h"
52 #include "libmesh/parallel.h"
53 
54 // Bring in everything from the libMesh namespace
55 using namespace libMesh;
56 
70 {
71 private:
72 
76  OptimizationSystem & _sys;
77 
78 public:
79 
84 
90  void assemble_A_and_F();
91 
95  virtual Number objective (const NumericVector<Number> & soln,
96  OptimizationSystem & /*sys*/);
97 
101  virtual void gradient (const NumericVector<Number> & soln,
102  NumericVector<Number> & grad_f,
103  OptimizationSystem & /*sys*/);
104 
108  virtual void hessian (const NumericVector<Number> & soln,
109  SparseMatrix<Number> & H_f,
110  OptimizationSystem & /*sys*/);
111 
115  virtual void equality_constraints (const NumericVector<Number> & X,
116  NumericVector<Number> & C_eq,
117  OptimizationSystem & /*sys*/);
118 
122  virtual void equality_constraints_jacobian (const NumericVector<Number> & X,
123  SparseMatrix<Number> & C_eq_jac,
124  OptimizationSystem & /*sys*/);
125 
129  virtual void inequality_constraints (const NumericVector<Number> & X,
130  NumericVector<Number> & C_ineq,
131  OptimizationSystem & /*sys*/);
132 
136  virtual void inequality_constraints_jacobian (const NumericVector<Number> & X,
137  SparseMatrix<Number> & C_ineq_jac,
138  OptimizationSystem & /*sys*/);
139 
143  virtual void lower_and_upper_bounds (OptimizationSystem & sys);
144 
150  SparseMatrix<Number> * A_matrix;
151 
156  NumericVector<Number> * F_vector;
157 };
158 
160  _sys(sys_in)
161 {}
162 
164 {
165  A_matrix->zero();
166  F_vector->zero();
167 
168  const MeshBase & mesh = _sys.get_mesh();
169 
170  const unsigned int dim = mesh.mesh_dimension();
171  const unsigned int u_var = _sys.variable_number ("u");
172 
173  const DofMap & dof_map = _sys.get_dof_map();
174  FEType fe_type = dof_map.variable_type(u_var);
175  std::unique_ptr<FEBase> fe (FEBase::build(dim, fe_type));
176  QGauss qrule (dim, fe_type.default_quadrature_order());
177  fe->attach_quadrature_rule (&qrule);
178 
179  const std::vector<Real> & JxW = fe->get_JxW();
180  const std::vector<std::vector<Real>> & phi = fe->get_phi();
181  const std::vector<std::vector<RealGradient>> & dphi = fe->get_dphi();
182 
183  std::vector<dof_id_type> dof_indices;
184 
187 
188  for (const auto & elem : mesh.active_local_element_ptr_range())
189  {
190  dof_map.dof_indices (elem, dof_indices);
191 
192  const unsigned int n_dofs = dof_indices.size();
193 
194  fe->reinit (elem);
195 
196  Ke.resize (n_dofs, n_dofs);
197  Fe.resize (n_dofs);
198 
199  for (unsigned int qp=0; qp<qrule.n_points(); qp++)
200  {
201  for (unsigned int dof_i=0; dof_i<n_dofs; dof_i++)
202  {
203  for (unsigned int dof_j=0; dof_j<n_dofs; dof_j++)
204  {
205  Ke(dof_i, dof_j) += JxW[qp] * (dphi[dof_j][qp]* dphi[dof_i][qp]);
206  }
207  Fe(dof_i) += JxW[qp] * phi[dof_i][qp];
208  }
209  }
210 
211  A_matrix->add_matrix (Ke, dof_indices);
212  F_vector->add_vector (Fe, dof_indices);
213  }
214 
215  A_matrix->close();
216  F_vector->close();
217 }
218 
220  OptimizationSystem & /*sys*/)
221 {
222  std::unique_ptr<NumericVector<Number>> AxU = soln.zero_clone();
223 
224  A_matrix->vector_mult(*AxU, soln);
225  Number UTxAxU = AxU->dot(soln);
226 
227  Number UTxF = F_vector->dot(soln);
228 
229  return 0.5 * UTxAxU - UTxF;
230 }
231 
233  NumericVector<Number> & grad_f,
234  OptimizationSystem & /*sys*/)
235 {
236  grad_f.zero();
237 
238  A_matrix->vector_mult(grad_f, soln);
239  grad_f.add(-1, *F_vector);
240 }
241 
242 
244  SparseMatrix<Number> & H_f,
245  OptimizationSystem & sys)
246 {
247  H_f.zero();
248  H_f.add(1., *A_matrix);
249 
250  // We also need to add the Hessian of the inequality and equality constraints,
251  // \sum_i^n_eq lambda_eq_i H_eq_i + \sum_i^n_ineq lambda_ineq_i H_ineq_i
252  // where lambda_eq and lambda_ineq denote Lagrange multipliers associated
253  // with the equality and inequality constraints, respectively.
254  // In this example, our equality constraints are linear, hence H_eq_i = 0.
255  // However, our inequality constraint is nonlinear, so it will contribute
256  // to the Hessian matrix.
257  sys.optimization_solver->get_dual_variables();
258  dof_id_type ineq_index = 0;
259  Number lambda_ineq_0 = 0.;
260  unsigned int lambda_rank = 0;
261  if ((sys.lambda_ineq->first_local_index() <= ineq_index) &&
262  (ineq_index < sys.lambda_ineq->last_local_index()))
263  {
264  lambda_ineq_0 = (*sys.lambda_ineq)(0);
265  lambda_rank = sys.comm().rank();
266  }
267 
268  // Sync lambda_rank across all processors.
269  sys.comm().sum(lambda_rank);
270  sys.comm().broadcast(lambda_rank, lambda_rank);
271 
272  if ((sys.get_dof_map().first_dof() <= 200) && (200 < sys.get_dof_map().end_dof()))
273  H_f.add(200, 200, 2. * lambda_ineq_0);
274 }
275 
277  NumericVector<Number> & C_eq,
278  OptimizationSystem & /*sys*/)
279 {
280  C_eq.zero();
281 
282  std::unique_ptr<NumericVector<Number>> X_localized =
284  X_localized->init(X.size(), false, SERIAL);
285  X.localize(*X_localized);
286 
287  std::vector<Number> constraint_values(3);
288  constraint_values[0] = (*X_localized)(17);
289  constraint_values[1] = (*X_localized)(23);
290  constraint_values[2] = (*X_localized)(98) + (*X_localized)(185);
291 
292  for (std::size_t i=0; i<constraint_values.size(); i++)
293  if ((C_eq.first_local_index() <= i) &&
294  (i < C_eq.last_local_index()))
295  C_eq.set(i, constraint_values[i]);
296 }
297 
299  SparseMatrix<Number> & C_eq_jac,
300  OptimizationSystem & sys)
301 {
302  C_eq_jac.zero();
303 
304  std::vector<std::vector<Number>> constraint_jac_values(3);
305  std::vector<std::vector<dof_id_type>> constraint_jac_indices(3);
306 
307  constraint_jac_values[0].resize(1);
308  constraint_jac_indices[0].resize(1);
309  constraint_jac_values[0][0] = 1.;
310  constraint_jac_indices[0][0] = 17;
311 
312  constraint_jac_values[1].resize(1);
313  constraint_jac_indices[1].resize(1);
314  constraint_jac_values[1][0] = 1.;
315  constraint_jac_indices[1][0] = 23;
316 
317  constraint_jac_values[2].resize(2);
318  constraint_jac_indices[2].resize(2);
319  constraint_jac_values[2][0] = 1.;
320  constraint_jac_values[2][1] = 1.;
321  constraint_jac_indices[2][0] = 98;
322  constraint_jac_indices[2][1] = 185;
323 
324  for (std::size_t i=0; i<constraint_jac_values.size(); i++)
325  for (std::size_t j=0; j<constraint_jac_values[i].size(); j++)
326  if ((sys.C_eq->first_local_index() <= i) &&
327  (i < sys.C_eq->last_local_index()))
328  {
329  dof_id_type col_index = constraint_jac_indices[i][j];
330  Number value = constraint_jac_values[i][j];
331  C_eq_jac.set(i, col_index, value);
332  }
333 }
334 
336  NumericVector<Number> & C_ineq,
337  OptimizationSystem & /*sys*/)
338 {
339  C_ineq.zero();
340 
341  std::unique_ptr<NumericVector<Number>> X_localized =
343  X_localized->init(X.size(), false, SERIAL);
344  X.localize(*X_localized);
345 
346  std::vector<Number> constraint_values(1);
347  constraint_values[0] = (*X_localized)(200)*(*X_localized)(200) + (*X_localized)(201) - 5.;
348 
349  for (std::size_t i=0; i<constraint_values.size(); i++)
350  if ((C_ineq.first_local_index() <= i) && (i < C_ineq.last_local_index()))
351  C_ineq.set(i, constraint_values[i]);
352 }
353 
355  SparseMatrix<Number> & C_ineq_jac,
356  OptimizationSystem & sys)
357 {
358  C_ineq_jac.zero();
359 
360  std::unique_ptr<NumericVector<Number>> X_localized =
362  X_localized->init(X.size(), false, SERIAL);
363  X.localize(*X_localized);
364 
365  std::vector<std::vector<Number>> constraint_jac_values(1);
366  std::vector<std::vector<dof_id_type>> constraint_jac_indices(1);
367 
368  constraint_jac_values[0].resize(2);
369  constraint_jac_indices[0].resize(2);
370  constraint_jac_values[0][0] = 2.* (*X_localized)(200);
371  constraint_jac_values[0][1] = 1.;
372  constraint_jac_indices[0][0] = 200;
373  constraint_jac_indices[0][1] = 201;
374 
375  for (std::size_t i=0; i<constraint_jac_values.size(); i++)
376  for (std::size_t j=0; j<constraint_jac_values[i].size(); j++)
377  if ((sys.C_ineq->first_local_index() <= i) &&
378  (i < sys.C_ineq->last_local_index()))
379  {
380  dof_id_type col_index = constraint_jac_indices[i][j];
381  Number value = constraint_jac_values[i][j];
382  C_ineq_jac.set(i, col_index, value);
383  }
384 
385 }
386 
388 {
389  for (unsigned int i=sys.get_dof_map().first_dof(); i<sys.get_dof_map().end_dof(); i++)
390  {
391  sys.get_vector("lower_bounds").set(i, -2.);
392  sys.get_vector("upper_bounds").set(i, 2.);
393  }
394 }
395 
396 int main (int argc, char ** argv)
397 {
398  LibMeshInit init (argc, argv);
399 
400 #ifndef LIBMESH_HAVE_PETSC_TAO
401 
402  libmesh_example_requires(false, "PETSc >= 3.5.0 with built-in TAO support");
403 
404 #elif LIBMESH_USE_COMPLEX_NUMBERS
405 
406  // According to
407  // http://www.mcs.anl.gov/research/projects/tao/documentation/installation.html
408  // TAO & PETSc-complex are currently mutually exclusive
409  libmesh_example_requires(false, "PETSc >= 3.5.0 with built-in TAO support & real-numbers only");
410 
411 #endif
412 
413  // We use a 2D domain.
414  libmesh_example_requires(LIBMESH_DIM > 1, "--disable-1D-only");
415 
416  // TAO is giving us problems in parallel?
417  if (init.comm().size() != 1)
418  {
419  libMesh::out << "This example can currently only be run in serial." << std::endl;
420  return 77;
421  }
422 
423  GetPot infile("optimization_ex2.in");
424  const std::string approx_order = infile("approx_order", "FIRST");
425  const std::string fe_family = infile("fe_family", "LAGRANGE");
426  const unsigned int n_elem = infile("n_elem", 10);
427 
428  Mesh mesh(init.comm());
430  n_elem,
431  n_elem,
432  -1., 1.,
433  -1., 1.,
434  QUAD9);
435 
436  mesh.print_info();
437 
438  EquationSystems equation_systems (mesh);
439 
440  OptimizationSystem & system =
441  equation_systems.add_system<OptimizationSystem> ("Optimization");
442 
443  // The default is to use PETSc/Tao solvers, but let the user change
444  // the optimization solver package on the fly.
445  {
446  const std::string optimization_solver_type = infile("optimization_solver_type",
447  "PETSC_SOLVERS");
448  SolverPackage sp = Utility::string_to_enum<SolverPackage>(optimization_solver_type);
449  std::unique_ptr<OptimizationSolver<Number>> new_solver =
451  system.optimization_solver.reset(new_solver.release());
452  }
453 
454  // Set tolerances and maximum iteration counts directly on the optimization solver.
455  system.optimization_solver->max_objective_function_evaluations = 128;
456  system.optimization_solver->objective_function_relative_tolerance = 1.e-4;
457  system.optimization_solver->verbose = true;
458 
459  AssembleOptimization assemble_opt(system);
460 
461  system.add_variable("u",
462  Utility::string_to_enum<Order> (approx_order),
463  Utility::string_to_enum<FEFamily>(fe_family));
464 
465  system.optimization_solver->objective_object = &assemble_opt;
466  system.optimization_solver->gradient_object = &assemble_opt;
467  system.optimization_solver->hessian_object = &assemble_opt;
468  system.optimization_solver->equality_constraints_object = &assemble_opt;
469  system.optimization_solver->equality_constraints_jacobian_object = &assemble_opt;
470  system.optimization_solver->inequality_constraints_object = &assemble_opt;
471  system.optimization_solver->inequality_constraints_jacobian_object = &assemble_opt;
472  system.optimization_solver->lower_and_upper_bounds_object = &assemble_opt;
473 
474  // system.matrix and system.rhs are used for the gradient and Hessian,
475  // so in this case we add an extra matrix and vector to store A and F.
476  // This makes it easy to write the code for evaluating the objective,
477  // gradient, and hessian.
478  system.add_matrix("A_matrix");
479  system.add_vector("F_vector");
480  assemble_opt.A_matrix = &system.get_matrix("A_matrix");
481  assemble_opt.F_vector = &system.get_vector("F_vector");
482 
483 
484  equation_systems.init();
485  equation_systems.print_info();
486 
487  assemble_opt.assemble_A_and_F();
488 
489  {
490  std::vector<std::set<numeric_index_type>> constraint_jac_sparsity;
491  std::set<numeric_index_type> sparsity_row;
492  sparsity_row.insert(17);
493  constraint_jac_sparsity.push_back(sparsity_row);
494  sparsity_row.clear();
495 
496  sparsity_row.insert(23);
497  constraint_jac_sparsity.push_back(sparsity_row);
498  sparsity_row.clear();
499 
500  sparsity_row.insert(98);
501  sparsity_row.insert(185);
502  constraint_jac_sparsity.push_back(sparsity_row);
503 
504  system.initialize_equality_constraints_storage(constraint_jac_sparsity);
505  }
506 
507  {
508  std::vector<std::set<numeric_index_type>> constraint_jac_sparsity;
509  std::set<numeric_index_type> sparsity_row;
510  sparsity_row.insert(200);
511  sparsity_row.insert(201);
512  constraint_jac_sparsity.push_back(sparsity_row);
513 
514  system.initialize_inequality_constraints_storage(constraint_jac_sparsity);
515  }
516 
517  // We need to close the matrix so that we can use it to store the
518  // Hessian during the solve.
519  system.matrix->close();
520  system.solve();
521 
522  // Print convergence information
523  system.optimization_solver->print_converged_reason();
524 
525  // Write the solution to file if the optimization solver converged
526  if (system.optimization_solver->get_converged_reason() > 0)
527  {
528  std::stringstream filename;
529  ExodusII_IO (mesh).write_equation_systems("optimization_soln.exo",
530  equation_systems);
531  }
532 
533  return 0;
534 }
libMesh::OptimizationSystem
This System subclass enables us to assemble an objective function, gradient, Hessian and bounds for o...
Definition: optimization_system.h:43
libMesh::SparseMatrix::close
virtual void close()=0
Calls the SparseMatrix's internal assembly routines, ensuring that the values are consistent across p...
libMesh::OptimizationSystem::C_ineq
std::unique_ptr< NumericVector< Number > > C_ineq
The vector that stores inequality constraints.
Definition: optimization_system.h:283
libMesh::NumericVector::zero
virtual void zero()=0
Set all entries to zero.
libMesh::SolverPackage
SolverPackage
Defines an enum for various linear solver packages.
Definition: enum_solver_package.h:34
libMesh::dof_id_type
uint8_t dof_id_type
Definition: id_types.h:67
libMesh::Number
Real Number
Definition: libmesh_common.h:195
libMesh::OptimizationSystem::initialize_inequality_constraints_storage
void initialize_inequality_constraints_storage(const std::vector< std::set< numeric_index_type >> &constraint_jac_sparsity)
Initialize storage for the inequality constraints, as per initialize_equality_constraints_storage.
Definition: optimization_system.C:127
libMesh::NumericVector::add
virtual void add(const numeric_index_type i, const T value)=0
Adds value to each entry of the vector.
libMesh::Mesh
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
Definition: mesh.h:50
AssembleOptimization::gradient
virtual void gradient(const NumericVector< Number > &soln, NumericVector< Number > &grad_f, OptimizationSystem &)
Evaluate the gradient.
Definition: optimization_ex1.C:216
libMesh::EquationSystems::add_system
virtual System & add_system(const std::string &system_type, const std::string &name)
Add the system of type system_type named name to the systems array.
Definition: equation_systems.C:345
libMesh::MeshTools::n_elem
dof_id_type n_elem(const MeshBase::const_element_iterator &begin, const MeshBase::const_element_iterator &end)
Count up the number of elements of a specific type (as defined by an iterator range).
Definition: mesh_tools.C:705
libMesh::QGauss
This class implements specific orders of Gauss quadrature.
Definition: quadrature_gauss.h:39
libMesh::DofMap::dof_indices
void dof_indices(const Elem *const elem, std::vector< dof_id_type > &di) const
Fills the vector di with the global degree of freedom indices for the element.
Definition: dof_map.C:1967
AssembleOptimization::AssembleOptimization
AssembleOptimization(OptimizationSystem &sys_in)
Constructor.
Definition: optimization_ex1.C:126
libMesh::OptimizationSystem::ComputeInequalityConstraintsJacobian
Abstract base class to be used to calculate the Jacobian of the inequality constraints.
Definition: optimization_system.h:187
AssembleOptimization::_sys
OptimizationSystem & _sys
Keep a reference to the OptimizationSystem.
Definition: optimization_ex1.C:76
libMesh::MeshBase::active_local_element_ptr_range
virtual SimpleRange< element_iterator > active_local_element_ptr_range()=0
libMesh::OptimizationSystem::solve
virtual void solve() override
Solves the optimization problem.
Definition: optimization_system.C:162
libMesh::SparseMatrix::vector_mult
void vector_mult(NumericVector< T > &dest, const NumericVector< T > &arg) const
Multiplies the matrix by the NumericVector arg and stores the result in NumericVector dest.
Definition: sparse_matrix.C:170
libMesh::SERIAL
Definition: enum_parallel_type.h:35
libMesh::NumericVector::last_local_index
virtual numeric_index_type last_local_index() const =0
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::OptimizationSystem::ComputeInequalityConstraints
Abstract base class to be used to calculate the inequality constraints.
Definition: optimization_system.h:169
AssembleOptimization::lower_and_upper_bounds
virtual void lower_and_upper_bounds(OptimizationSystem &sys)
Evaluate the lower and upper bounds vectors.
Definition: optimization_ex2.C:387
libMesh::NumericVector::close
virtual void close()=0
Calls the NumericVector's internal assembly routines, ensuring that the values are consistent across ...
libMesh::ImplicitSystem::add_matrix
SparseMatrix< Number > & add_matrix(const std::string &mat_name)
Adds the additional matrix mat_name to this system.
Definition: implicit_system.C:202
libMesh::NumericVector::init
virtual void init(const numeric_index_type n, const numeric_index_type n_local, const bool fast=false, const ParallelType ptype=AUTOMATIC)=0
Change the dimension of the vector to n.
AssembleOptimization
This class encapsulate all functionality required for assembling the objective function,...
Definition: optimization_ex1.C:66
AssembleOptimization::inequality_constraints
virtual void inequality_constraints(const NumericVector< Number > &X, NumericVector< Number > &C_ineq, OptimizationSystem &)
Evaluate the inequality constraints.
Definition: optimization_ex2.C:335
libMesh::ParallelObject::comm
const Parallel::Communicator & comm() const
Definition: parallel_object.h:94
libMesh::DenseMatrix< Number >
libMesh::OptimizationSystem::ComputeHessian
Abstract base class to be used to calculate the Hessian of an objective function.
Definition: optimization_system.h:116
mesh
MeshBase & mesh
Definition: mesh_communication.C:1257
libMesh::MeshBase::mesh_dimension
unsigned int mesh_dimension() const
Definition: mesh_base.C:135
libMesh::DofMap::first_dof
dof_id_type first_dof(const processor_id_type proc) const
Definition: dof_map.h:650
libMesh::ExodusII_IO
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
Definition: exodusII_io.h:51
libMesh::OptimizationSystem::lambda_ineq
std::unique_ptr< NumericVector< Number > > lambda_ineq
Definition: optimization_system.h:295
libMesh::SparseMatrix< Number >
dim
unsigned int dim
Definition: adaptivity_ex3.C:113
libMesh::DenseMatrix::resize
void resize(const unsigned int new_m, const unsigned int new_n)
Resize the matrix.
Definition: dense_matrix.h:822
libMesh::NumericVector::size
virtual numeric_index_type size() const =0
libMesh::NumericVector< Number >
libMesh::TriangleWrapper::init
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
AssembleOptimization::inequality_constraints_jacobian
virtual void inequality_constraints_jacobian(const NumericVector< Number > &X, SparseMatrix< Number > &C_ineq_jac, OptimizationSystem &)
Evaluate the inequality constraints Jacobian.
Definition: optimization_ex2.C:354
libMesh::OptimizationSystem::C_eq
std::unique_ptr< NumericVector< Number > > C_eq
The vector that stores equality constraints.
Definition: optimization_system.h:273
AssembleOptimization::F_vector
NumericVector< Number > * F_vector
Vector for storing F.
Definition: optimization_ex1.C:123
libMesh::MeshTools::Generation::build_square
void build_square(UnstructuredMesh &mesh, const unsigned int nx, const unsigned int ny, const Real xmin=0., const Real xmax=1., const Real ymin=0., const Real ymax=1., const ElemType type=INVALID_ELEM, const bool gauss_lobatto_grid=false)
A specialized build_cube() for 2D meshes.
Definition: mesh_generation.C:1501
libMesh::OptimizationSystem::ComputeEqualityConstraints
Abstract base class to be used to calculate the equality constraints.
Definition: optimization_system.h:135
libMesh::NumericVector::add_vector
virtual void add_vector(const T *v, const std::vector< numeric_index_type > &dof_indices)
Computes , where v is a pointer and each dof_indices[i] specifies where to add value v[i].
Definition: numeric_vector.C:363
libMesh::MeshBase
This is the MeshBase class.
Definition: mesh_base.h:78
libMesh::NumericVector::localize
virtual void localize(std::vector< T > &v_local) const =0
Creates a copy of the global vector in the local vector v_local.
AssembleOptimization::hessian
virtual void hessian(const NumericVector< Number > &soln, SparseMatrix< Number > &H_f, OptimizationSystem &)
Evaluate the Hessian.
Definition: optimization_ex1.C:230
libMesh::OptimizationSolver
This base class can be inherited from to provide interfaces to optimization solvers from different pa...
Definition: optimization_solver.h:60
libMesh::OptimizationSystem::ComputeLowerAndUpperBounds
Abstract base class to be used to calculate the lower and upper bounds for all dofs in the system.
Definition: optimization_system.h:204
AssembleOptimization::assemble_A_and_F
void assemble_A_and_F()
The optimization problem we consider here is: min_U 0.5*U^T A U - U^T F.
Definition: optimization_ex1.C:130
libMesh::System::get_mesh
const MeshBase & get_mesh() const
Definition: system.h:2083
libMesh::OptimizationSystem::ComputeEqualityConstraintsJacobian
Abstract base class to be used to calculate the Jacobian of the equality constraints.
Definition: optimization_system.h:153
libMesh::OptimizationSystem::optimization_solver
std::unique_ptr< OptimizationSolver< Number > > optimization_solver
The OptimizationSolver that is used for performing the optimization.
Definition: optimization_system.h:268
libMesh::System::add_variable
unsigned int add_variable(const std::string &var, const FEType &type, const std::set< subdomain_id_type > *const active_subdomains=nullptr)
Adds the variable var to the list of variables for this system.
Definition: system.C:1069
libMesh::EquationSystems::init
virtual void init()
Initialize all the systems.
Definition: equation_systems.C:96
libMesh::ImplicitSystem::matrix
SparseMatrix< Number > * matrix
The system matrix.
Definition: implicit_system.h:393
libMesh::DofMap::variable_type
const FEType & variable_type(const unsigned int c) const
Definition: dof_map.h:1924
libMesh::LibMeshInit
The LibMeshInit class, when constructed, initializes the dependent libraries (e.g.
Definition: libmesh.h:83
libMesh::OptimizationSystem::ComputeObjective
Abstract base class to be used to calculate the objective function for optimization.
Definition: optimization_system.h:74
libMesh::OptimizationSystem::ComputeGradient
Abstract base class to be used to calculate the gradient of an objective function.
Definition: optimization_system.h:95
libMesh::SparseMatrix::add
virtual void add(const numeric_index_type i, const numeric_index_type j, const T value)=0
Add value to the element (i,j).
libMesh::EquationSystems::print_info
void print_info(std::ostream &os=libMesh::out) const
Prints information about the equation systems, by default to libMesh::out.
Definition: equation_systems.C:1314
libMesh::EquationSystems
This is the EquationSystems class.
Definition: equation_systems.h:74
libMesh::MeshOutput::write_equation_systems
virtual void write_equation_systems(const std::string &, const EquationSystems &, const std::set< std::string > *system_names=nullptr)
This method implements writing a mesh with data to a specified file where the data is taken from the ...
Definition: mesh_output.C:31
libMesh::OptimizationSystem::initialize_equality_constraints_storage
void initialize_equality_constraints_storage(const std::vector< std::set< numeric_index_type >> &constraint_jac_sparsity)
Initialize storage for the equality constraints, and the corresponding Jacobian.
Definition: optimization_system.C:91
libMesh::DenseVector::resize
void resize(const unsigned int n)
Resize the vector.
Definition: dense_vector.h:355
libMesh::SparseMatrix::add_matrix
virtual void add_matrix(const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols)=0
Add the full matrix dm to the SparseMatrix.
libMesh::FEType
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:178
main
int main(int argc, char **argv)
Definition: optimization_ex2.C:396
AssembleOptimization::equality_constraints
virtual void equality_constraints(const NumericVector< Number > &X, NumericVector< Number > &C_eq, OptimizationSystem &)
Evaluate the equality constraints.
Definition: optimization_ex2.C:276
libMesh::NumericVector::set
virtual void set(const numeric_index_type i, const T value)=0
Sets v(i) = value.
value
static const bool value
Definition: xdr_io.C:56
libMesh::DofMap
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:176
libMesh::MeshBase::print_info
void print_info(std::ostream &os=libMesh::out) const
Prints relevant information about the mesh.
Definition: mesh_base.C:585
AssembleOptimization::objective
virtual Number objective(const NumericVector< Number > &soln, OptimizationSystem &)
Evaluate the objective function.
Definition: optimization_ex1.C:202
libMesh::System::variable_number
unsigned short int variable_number(const std::string &var) const
Definition: system.C:1232
AssembleOptimization::A_matrix
SparseMatrix< Number > * A_matrix
Sparse matrix for storing the matrix A.
Definition: optimization_ex1.C:117
libMesh::System::add_vector
NumericVector< Number > & add_vector(const std::string &vec_name, const bool projections=true, const ParallelType type=PARALLEL)
Adds the additional vector vec_name to this system.
Definition: system.C:661
libMesh::QUAD9
Definition: enum_elem_type.h:43
libMesh::System::get_dof_map
const DofMap & get_dof_map() const
Definition: system.h:2099
AssembleOptimization::equality_constraints_jacobian
virtual void equality_constraints_jacobian(const NumericVector< Number > &X, SparseMatrix< Number > &C_eq_jac, OptimizationSystem &)
Evaluate the equality constraints Jacobian.
Definition: optimization_ex2.C:298
libMesh::SparseMatrix::set
virtual void set(const numeric_index_type i, const numeric_index_type j, const T value)=0
Set the element (i,j) to value.
libMesh::FEType::default_quadrature_order
Order default_quadrature_order() const
Definition: fe_type.h:353
libMesh::out
OStreamProxy out
libMesh::ImplicitSystem::get_matrix
const SparseMatrix< Number > & get_matrix(const std::string &mat_name) const
Definition: implicit_system.C:262
libMesh::System::get_vector
const NumericVector< Number > & get_vector(const std::string &vec_name) const
Definition: system.C:774
libMesh::SparseMatrix::zero
virtual void zero()=0
Set all entries to 0.
libMesh::NumericVector::zero_clone
virtual std::unique_ptr< NumericVector< T > > zero_clone() const =0
libMesh::NumericVector::dot
virtual T dot(const NumericVector< T > &v) const =0
libMesh::DenseVector< Number >
libMesh::DofMap::end_dof
dof_id_type end_dof(const processor_id_type proc) const
Definition: dof_map.h:692
libMesh::NumericVector::first_local_index
virtual numeric_index_type first_local_index() const =0