$extrastylesheet
implicit_system.h
Go to the documentation of this file.
00001 // The libMesh Finite Element Library.
00002 // Copyright (C) 2002-2014 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
00003 
00004 // This library is free software; you can redistribute it and/or
00005 // modify it under the terms of the GNU Lesser General Public
00006 // License as published by the Free Software Foundation; either
00007 // version 2.1 of the License, or (at your option) any later version.
00008 
00009 // This library is distributed in the hope that it will be useful,
00010 // but WITHOUT ANY WARRANTY; without even the implied warranty of
00011 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00012 // Lesser General Public License for more details.
00013 
00014 // You should have received a copy of the GNU Lesser General Public
00015 // License along with this library; if not, write to the Free Software
00016 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00017 
00018 
00019 
00020 #ifndef LIBMESH_IMPLICIT_SYSTEM_H
00021 #define LIBMESH_IMPLICIT_SYSTEM_H
00022 
00023 // Local Includes
00024 #include "libmesh/explicit_system.h"
00025 
00026 // C++ includes
00027 #include <cstddef>
00028 
00029 namespace libMesh
00030 {
00031 
00032 // Forward declarations
00033 template <typename T> class LinearSolver;
00034 template <typename T> class SparseMatrix;
00035 
00036 
00037 
00046 // ------------------------------------------------------------
00047 // ImplicitSystem class definition
00048 
00049 class ImplicitSystem : public ExplicitSystem
00050 {
00051 public:
00052 
00057   ImplicitSystem (EquationSystems& es,
00058                   const std::string& name,
00059                   const unsigned int number);
00060 
00064   virtual ~ImplicitSystem ();
00065 
00069   typedef ImplicitSystem sys_type;
00070 
00074   sys_type & system () { return *this; }
00075 
00079   typedef ExplicitSystem Parent;
00080 
00085   virtual void clear ();
00086 
00091   virtual void reinit ();
00092 
00098   virtual void assemble ();
00099 
00100   //   /**
00101   //    * Assembles & solves the linear system Ax=b.
00102   //    */
00103   //   virtual void solve ();
00104 
00105 
00110   virtual void disable_cache ();
00111 
00116   virtual std::string system_type () const { return "Implicit"; }
00117 
00126   virtual LinearSolver<Number> *get_linear_solver() const;
00127 
00133   virtual std::pair<unsigned int, Real>
00134   get_linear_solve_parameters() const;
00135 
00140   virtual void release_linear_solver(LinearSolver<Number> *) const;
00141 
00149   virtual void assembly
00150     (bool /* get_residual */,
00151      bool /* get_jacobian */,
00152      bool /* apply_heterogeneous_constraints */ = false)
00153   { libmesh_not_implemented(); }
00154 
00166   virtual void assemble_residual_derivatives (const ParameterVector& parameters);
00167 
00175   virtual std::pair<unsigned int, Real>
00176   sensitivity_solve (const ParameterVector& parameters);
00177 
00186   virtual std::pair<unsigned int, Real>
00187   weighted_sensitivity_solve (const ParameterVector& parameters,
00188                               const ParameterVector& weights);
00189 
00199   virtual std::pair<unsigned int, Real>
00200   adjoint_solve (const QoISet& qoi_indices = QoISet());
00201 
00214   virtual std::pair<unsigned int, Real>
00215   weighted_sensitivity_adjoint_solve (const ParameterVector& parameters,
00216                                       const ParameterVector& weights,
00217                                       const QoISet& qoi_indices = QoISet());
00218 
00230   virtual void adjoint_qoi_parameter_sensitivity (const QoISet& qoi_indices,
00231                                                   const ParameterVector& parameters,
00232                                                   SensitivityData& sensitivities);
00233 
00245   virtual void forward_qoi_parameter_sensitivity (const QoISet& qoi_indices,
00246                                                   const ParameterVector& parameters,
00247                                                   SensitivityData& sensitivities);
00248 
00257   virtual void qoi_parameter_hessian(const QoISet& qoi_indices,
00258                                      const ParameterVector& parameters,
00259                                      SensitivityData& hessian);
00260 
00271   virtual void qoi_parameter_hessian_vector_product(const QoISet& qoi_indices,
00272                                                     const ParameterVector& parameters,
00273                                                     const ParameterVector& vector,
00274                                                     SensitivityData& product);
00275 
00279   typedef std::map<std::string, SparseMatrix<Number>* >::iterator        matrices_iterator;
00280   typedef std::map<std::string, SparseMatrix<Number>* >::const_iterator  const_matrices_iterator;
00281 
00290   SparseMatrix<Number> & add_matrix (const std::string& mat_name);
00291 
00296   bool have_matrix (const std::string& mat_name) const;
00297 
00303   const SparseMatrix<Number> * request_matrix (const std::string& mat_name) const;
00304 
00310   SparseMatrix<Number> * request_matrix (const std::string& mat_name);
00311 
00318   const SparseMatrix<Number> & get_matrix (const std::string& mat_name) const;
00319 
00326   SparseMatrix<Number> & get_matrix (const std::string& mat_name);
00327 
00331   virtual unsigned int n_matrices () const;
00332 
00338   SparseMatrix<Number> * matrix;
00339 
00340 
00341 
00342 protected:
00343 
00348   virtual void init_data ();
00349 
00353   virtual void init_matrices ();
00354 
00355 
00356 
00357 private:
00358 
00363   void add_system_matrix ();
00364 
00368   std::map<std::string, SparseMatrix<Number>* > _matrices;
00369 
00373   bool _can_add_matrices;
00374 };
00375 
00376 
00377 
00378 // ------------------------------------------------------------
00379 // ImplicitSystem inline methods
00380 inline
00381 bool ImplicitSystem::have_matrix (const std::string& mat_name) const
00382 {
00383   return (_matrices.count(mat_name));
00384 }
00385 
00386 
00387 inline
00388 unsigned int ImplicitSystem::n_matrices () const
00389 {
00390   return cast_int<unsigned int>(_matrices.size());
00391 }
00392 
00393 } // namespace libMesh
00394 
00395 #endif // LIBMESH_IMPLICIT_SYSTEM_H