$extrastylesheet
linear_implicit_system.C
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 // C++ includes
00021 
00022 // Local includes
00023 #include "libmesh/linear_implicit_system.h"
00024 #include "libmesh/linear_solver.h"
00025 #include "libmesh/equation_systems.h"
00026 #include "libmesh/numeric_vector.h" // for parameter sensitivity calcs
00027 //#include "libmesh/parameter_vector.h"
00028 #include "libmesh/sparse_matrix.h" // for get_transpose
00029 #include "libmesh/system_subset.h"
00030 
00031 namespace libMesh
00032 {
00033 
00034 
00035 // ------------------------------------------------------------
00036 // LinearImplicitSystem implementation
00037 LinearImplicitSystem::LinearImplicitSystem (EquationSystems& es,
00038                                             const std::string& name_in,
00039                                             const unsigned int number_in) :
00040 
00041   Parent                 (es, name_in, number_in),
00042   linear_solver          (LinearSolver<Number>::build(es.comm())),
00043   _n_linear_iterations   (0),
00044   _final_linear_residual (1.e20),
00045   _shell_matrix(NULL),
00046   _subset(NULL),
00047   _subset_solve_mode(SUBSET_ZERO)
00048 {
00049 }
00050 
00051 
00052 
00053 LinearImplicitSystem::~LinearImplicitSystem ()
00054 {
00055   // Clear data
00056   this->clear();
00057 }
00058 
00059 
00060 
00061 void LinearImplicitSystem::clear ()
00062 {
00063   // clear the linear solver
00064   linear_solver->clear();
00065 
00066   this->restrict_solve_to(NULL);
00067 
00068   // clear the parent data
00069   Parent::clear();
00070 }
00071 
00072 
00073 
00074 void LinearImplicitSystem::init_data ()
00075 {
00076   // initialize parent data
00077   Parent::init_data();
00078 
00079   // re-initialize the linear solver interface
00080   linear_solver->clear();
00081 }
00082 
00083 
00084 
00085 void LinearImplicitSystem::reinit ()
00086 {
00087   // re-initialize the linear solver interface
00088   linear_solver->clear();
00089 
00090   // initialize parent data
00091   Parent::reinit();
00092 }
00093 
00094 
00095 
00096 void LinearImplicitSystem::restrict_solve_to (const SystemSubset* subset,
00097                                               const SubsetSolveMode subset_solve_mode)
00098 {
00099   _subset = subset;
00100   _subset_solve_mode = subset_solve_mode;
00101   if(subset!=NULL)
00102     {
00103       libmesh_assert_equal_to (&subset->get_system(), this);
00104     }
00105 }
00106 
00107 
00108 
00109 void LinearImplicitSystem::solve ()
00110 {
00111   if (this->assemble_before_solve)
00112     // Assemble the linear system
00113     this->assemble ();
00114 
00115   // Log how long the linear solve takes.
00116   // This gets done by the LinearSolver classes now [RHS]
00117   // START_LOG("solve()", "System");
00118 
00119   // Get a reference to the EquationSystems
00120   const EquationSystems& es =
00121     this->get_equation_systems();
00122 
00123   // If the linear solver hasn't been initialized, we do so here.
00124   if (libMesh::on_command_line("--solver_system_names"))
00125     linear_solver->init((this->name()+"_").c_str());
00126   else
00127     linear_solver->init();
00128 
00129   // Get the user-specifiied linear solver tolerance
00130   const Real tol            =
00131     es.parameters.get<Real>("linear solver tolerance");
00132 
00133   // Get the user-specified maximum # of linear solver iterations
00134   const unsigned int maxits =
00135     es.parameters.get<unsigned int>("linear solver maximum iterations");
00136 
00137   if(_subset!=NULL)
00138     {
00139       linear_solver->restrict_solve_to(&_subset->dof_ids(),_subset_solve_mode);
00140     }
00141 
00142   // Solve the linear system.  Several cases:
00143   std::pair<unsigned int, Real> rval = std::make_pair(0,0.0);
00144   if(_shell_matrix)
00145     // 1.) Shell matrix with or without user-supplied preconditioner.
00146     rval = linear_solver->solve(*_shell_matrix, this->request_matrix("Preconditioner"), *solution, *rhs, tol, maxits);
00147   else
00148     // 2.) No shell matrix, with or without user-supplied preconditioner
00149     rval = linear_solver->solve (*matrix, this->request_matrix("Preconditioner"), *solution, *rhs, tol, maxits);
00150 
00151   if(_subset!=NULL)
00152     {
00153       linear_solver->restrict_solve_to(NULL);
00154     }
00155 
00156   // Store the number of linear iterations required to
00157   // solve and the final residual.
00158   _n_linear_iterations   = rval.first;
00159   _final_linear_residual = rval.second;
00160 
00161   // Stop logging the linear solve
00162   // This gets done by the LinearSolver classes now [RHS]
00163   // STOP_LOG("solve()", "System");
00164 
00165   // Update the system after the solve
00166   this->update();
00167 }
00168 
00169 
00170 
00171 void LinearImplicitSystem::attach_shell_matrix (ShellMatrix<Number>* shell_matrix)
00172 {
00173   _shell_matrix = shell_matrix;
00174 }
00175 
00176 
00177 /*
00178   void LinearImplicitSystem::sensitivity_solve (const ParameterVector& parameters)
00179   {
00180   if (this->assemble_before_solve)
00181   {
00182   // Assemble the linear system
00183   this->assemble ();
00184 
00185   // But now assemble right hand sides with the residual's
00186   // parameter derivatives
00187   this->assemble_residual_derivatives(parameters);
00188   }
00189 
00190   // Get a reference to the EquationSystems
00191   const EquationSystems& es =
00192   this->get_equation_systems();
00193 
00194   // Get the user-specifiied linear solver tolerance
00195   const Real tol            =
00196   es.parameters.get<Real>("sensitivity solver tolerance");
00197 
00198   // Get the user-specified maximum # of linear solver iterations
00199   const unsigned int maxits =
00200   es.parameters.get<unsigned int>("sensitivity solver maximum iterations");
00201 
00202   // Our iteration counts and residuals will be sums of the individual
00203   // results
00204   _n_linear_iterations = 0;
00205   _final_linear_residual = 0.0;
00206   std::pair<unsigned int, Real> rval = std::make_pair(0,0.0);
00207 
00208   // Solve the linear system.
00209   SparseMatrix<Number> *pc = this->request_matrix("Preconditioner");
00210   for (unsigned int p=0; p != parameters.size(); ++p)
00211   {
00212   rval = linear_solver->solve (*matrix, pc,
00213   this->get_sensitivity_solution(p),
00214   this->get_sensitivity_rhs(p), tol, maxits);
00215   _n_linear_iterations   += rval.first;
00216   _final_linear_residual += rval.second;
00217   }
00218 
00219   // Our matrix is the *negative* of the Jacobian for b-A*u, so our
00220   // solutions are all inverted
00221   for (unsigned int p=0; p != parameters.size(); ++p)
00222   {
00223   this->get_sensitivity_solution(p) *= -1.0;
00224   }
00225   }
00226 
00227 
00228 
00229   void LinearImplicitSystem::adjoint_solve (const QoISet &qoi_indices)
00230   {
00231   const unsigned int Nq = this->qoi.size();
00232 
00233   // We currently don't support adjoint solves of shell matrices
00234   // FIXME - we should let shell matrices support
00235   // vector_transpose_mult so that we can use them here.
00236   if(_shell_matrix!=NULL)
00237   libmesh_not_implemented();
00238 
00239   if (this->assemble_before_solve)
00240   {
00241   // Assemble the linear system
00242   this->assemble ();
00243 
00244   // And take the adjoint
00245   matrix->get_transpose(*matrix);
00246 
00247   // Including of any separate preconditioner
00248   SparseMatrix<Number> *pc = this->request_matrix("Preconditioner");
00249   if(pc)
00250   pc->get_transpose(*pc);
00251 
00252   // But now replace the right hand sides with the quantity of
00253   // interest functional's derivative
00254   this->assemble_qoi_derivative(qoi_indices);
00255   }
00256 
00257   // Get a reference to the EquationSystems
00258   const EquationSystems& es =
00259   this->get_equation_systems();
00260 
00261   // Get the user-specifiied linear solver tolerance
00262   const Real tol            =
00263   es.parameters.get<Real>("adjoint solver tolerance");
00264 
00265   // Get the user-specified maximum # of linear solver iterations
00266   const unsigned int maxits =
00267   es.parameters.get<unsigned int>("adjoint solver maximum iterations");
00268 
00269   // Our iteration counts and residuals will be sums of the individual
00270   // results
00271   _n_linear_iterations = 0;
00272   _final_linear_residual = 0.0;
00273   std::pair<unsigned int, Real> rval = std::make_pair(0,0.0);
00274 
00275   // Solve the linear system.
00276   SparseMatrix<Number> *pc = this->request_matrix("Preconditioner");
00277   for (unsigned int i=0; i != Nq; ++i)
00278   if (qoi_indices.has_index(i))
00279   {
00280   rval = linear_solver->solve (*matrix, pc,
00281   this->add_adjoint_solution(i),
00282   this->get_adjoint_rhs(i), tol, maxits);
00283   _n_linear_iterations   += rval.first;
00284   _final_linear_residual += rval.second;
00285   }
00286   }
00287 
00288 
00289 
00290   void LinearImplicitSystem::forward_qoi_parameter_sensitivity
00291   (const QoISet&          qoi_indices,
00292   const ParameterVector& parameters,
00293   SensitivityData&       sensitivities)
00294   {
00295   const unsigned int Np = parameters.size();
00296   const unsigned int Nq = this->qoi.size();
00297 
00298   // An introduction to the problem:
00299   //
00300   // A(p)*u(p) = b(p), where x is determined implicitly.
00301   // Residual R(u(p),p) := b(p) - A(p)*u(p)
00302   // partial R / partial u = -A
00303   //
00304   // This implies that:
00305   // d/dp(R) = 0
00306   // (partial b / partial p) -
00307   // (partial A / partial p) * u -
00308   // A * (partial u / partial p) = 0
00309   // A * (partial u / partial p) = (partial R / partial p)
00310   //   = (partial b / partial p) - (partial A / partial p) * u
00311 
00312   // We first solve for (partial u / partial p) for each parameter:
00313   // -A * (partial u / partial p) = - (partial R / partial p)
00314 
00315   this->sensitivity_solve(parameters);
00316 
00317   // Get ready to fill in senstivities:
00318   sensitivities.allocate_data(qoi_indices, *this, parameters);
00319 
00320   // We use the identity:
00321   // dq/dp = (partial q / partial p) + (partial q / partial u) *
00322   //         (partial u / partial p)
00323 
00324   // We get (partial q / partial u) from the user
00325   this->assemble_qoi_derivative(qoi_indices);
00326 
00327   for (unsigned int j=0; j != Np; ++j)
00328   {
00329   // We currently get partial derivatives via central differencing
00330   Number delta_p = 1e-6;
00331 
00332   // (partial q / partial p) ~= (q(p+dp)-q(p-dp))/(2*dp)
00333 
00334   Number old_parameter = *parameters[j];
00335 
00336   *parameters[j] = old_parameter - delta_p;
00337   this->assemble_qoi(qoi_indices);
00338   std::vector<Number> qoi_minus = this->qoi;
00339 
00340   *parameters[j] = old_parameter + delta_p;
00341   this->assemble_qoi(qoi_indices);
00342   std::vector<Number> &qoi_plus = this->qoi;
00343   std::vector<Number> partialq_partialp(Nq, 0);
00344   for (unsigned int i=0; i != Nq; ++i)
00345   if (qoi_indices.has_index(i))
00346   partialq_partialp[i] = (qoi_plus[i] - qoi_minus[i]) / (2.*delta_p);
00347 
00348   for (unsigned int i=0; i != Nq; ++i)
00349   if (qoi_indices.has_index(i))
00350   sensitivities[i][j] = partialq_partialp[i] +
00351   this->get_adjoint_rhs(i).dot(this->get_sensitivity_solution(i));
00352   }
00353 
00354   // All parameters have been reset.
00355   // Don't leave the qoi or system changed - principle of least
00356   // surprise.
00357   this->assemble();
00358   this->rhs->close();
00359   this->matrix->close();
00360   this->assemble_qoi(qoi_indices);
00361   }
00362 */
00363 
00364 
00365 
00366 LinearSolver<Number>* LinearImplicitSystem::get_linear_solver() const
00367 {
00368   return linear_solver.get();
00369 }
00370 
00371 
00372 
00373 void LinearImplicitSystem::release_linear_solver(LinearSolver<Number>*) const
00374 {
00375 }
00376 
00377 
00378 
00379 void LinearImplicitSystem::assembly(bool,
00380                                     bool,
00381                                     bool)
00382 {
00383   // Residual R(u(p),p) := A(p)*u(p) - b(p)
00384   // partial R / partial u = A
00385 
00386   this->assemble();
00387   this->rhs->close();
00388   this->matrix->close();
00389 
00390   *(this->rhs) *= -1.0;
00391   this->rhs->add_vector(*this->solution, *this->matrix);
00392 }
00393 
00394 } // namespace libMesh