$extrastylesheet
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/dof_map.h"
00024 #include "libmesh/equation_systems.h"
00025 #include "libmesh/implicit_system.h"
00026 #include "libmesh/libmesh_logging.h"
00027 #include "libmesh/linear_solver.h"
00028 #include "libmesh/mesh_base.h"
00029 #include "libmesh/numeric_vector.h"
00030 #include "libmesh/parameters.h"
00031 #include "libmesh/parameter_vector.h"
00032 #include "libmesh/qoi_set.h"
00033 #include "libmesh/sensitivity_data.h"
00034 #include "libmesh/sparse_matrix.h"
00035 
00036 namespace libMesh
00037 {
00038 
00039 // ------------------------------------------------------------
00040 // ImplicitSystem implementation
00041 ImplicitSystem::ImplicitSystem (EquationSystems& es,
00042                                 const std::string& name_in,
00043                                 const unsigned int number_in) :
00044 
00045   Parent            (es, name_in, number_in),
00046   matrix            (NULL),
00047   _can_add_matrices (true)
00048 {
00049 }
00050 
00051 
00052 
00053 ImplicitSystem::~ImplicitSystem ()
00054 {
00055   // Clear data
00056   this->clear();
00057 }
00058 
00059 
00060 
00061 void ImplicitSystem::clear ()
00062 {
00063   // clear the parent data
00064   Parent::clear();
00065 
00066   // clear any user-added matrices
00067   {
00068     for (matrices_iterator pos = _matrices.begin();
00069          pos != _matrices.end(); ++pos)
00070       {
00071         pos->second->clear ();
00072         delete pos->second;
00073         pos->second = NULL;
00074       }
00075 
00076     _matrices.clear();
00077     _can_add_matrices = true;
00078   }
00079 
00080   // NULL the matrix.
00081   matrix = NULL;
00082 }
00083 
00084 
00085 
00086 void ImplicitSystem::init_data ()
00087 {
00088   // initialize parent data
00089   Parent::init_data();
00090 
00091   // Clear any existing matrices
00092   for (matrices_iterator pos = _matrices.begin();
00093        pos != _matrices.end(); ++pos)
00094     pos->second->clear();
00095 
00096   // Add the system matrix.
00097   this->add_system_matrix ();
00098 
00099   // Initialize the matrices for the system
00100   this->init_matrices ();
00101 }
00102 
00103 
00104 
00105 void ImplicitSystem::init_matrices ()
00106 {
00107   libmesh_assert(matrix);
00108 
00109   // Check for quick return in case the system matrix
00110   // (and by extension all the matrices) has already
00111   // been initialized
00112   if (matrix->initialized())
00113     return;
00114 
00115   // Get a reference to the DofMap
00116   DofMap& dof_map = this->get_dof_map();
00117 
00118   // no chance to add other matrices
00119   _can_add_matrices = false;
00120 
00121   // Tell the matrices about the dof map, and vice versa
00122   for (matrices_iterator pos = _matrices.begin();
00123        pos != _matrices.end(); ++pos)
00124     {
00125       SparseMatrix<Number> &m = *(pos->second);
00126       libmesh_assert (!m.initialized());
00127 
00128       // We want to allow repeated init() on systems, but we don't
00129       // want to attach the same matrix to the DofMap twice
00130       if (!dof_map.is_attached(m))
00131         dof_map.attach_matrix (m);
00132     }
00133 
00134   // Compute the sparsity pattern for the current
00135   // mesh and DOF distribution.  This also updates
00136   // additional matrices, \p DofMap now knows them
00137   dof_map.compute_sparsity (this->get_mesh());
00138 
00139   // Initialize matrices
00140   for (matrices_iterator pos = _matrices.begin();
00141        pos != _matrices.end(); ++pos)
00142     pos->second->init ();
00143 
00144   // Set the additional matrices to 0.
00145   for (matrices_iterator pos = _matrices.begin();
00146        pos != _matrices.end(); ++pos)
00147     pos->second->zero ();
00148 }
00149 
00150 
00151 
00152 void ImplicitSystem::reinit ()
00153 {
00154   // initialize parent data
00155   Parent::reinit();
00156 
00157   // Get a reference to the DofMap
00158   DofMap& dof_map = this->get_dof_map();
00159 
00160   // Clear the matrices
00161   for (matrices_iterator pos = _matrices.begin();
00162        pos != _matrices.end(); ++pos)
00163     {
00164       pos->second->clear();
00165       pos->second->attach_dof_map (dof_map);
00166     }
00167 
00168   // Clear the sparsity pattern
00169   this->get_dof_map().clear_sparsity();
00170 
00171   // Compute the sparsity pattern for the current
00172   // mesh and DOF distribution.  This also updates
00173   // additional matrices, \p DofMap now knows them
00174   dof_map.compute_sparsity (this->get_mesh());
00175 
00176   // Initialize matrices
00177   for (matrices_iterator pos = _matrices.begin();
00178        pos != _matrices.end(); ++pos)
00179     pos->second->init ();
00180 
00181   // Set the additional matrices to 0.
00182   for (matrices_iterator pos = _matrices.begin();
00183        pos != _matrices.end(); ++pos)
00184     pos->second->zero ();
00185 }
00186 
00187 
00188 
00189 void ImplicitSystem::assemble ()
00190 {
00191   libmesh_assert(matrix);
00192   libmesh_assert (matrix->initialized());
00193   libmesh_assert(rhs);
00194   libmesh_assert (rhs->initialized());
00195 
00196   // The user assembly gets to expect to accumulate on an initially
00197   // empty system
00198   matrix->zero ();
00199   rhs->zero ();
00200 
00201   // Call the base class assemble function
00202   Parent::assemble ();
00203 }
00204 
00205 
00206 
00207 SparseMatrix<Number> & ImplicitSystem::add_matrix (const std::string& mat_name)
00208 {
00209   // only add matrices before initializing...
00210   if (!_can_add_matrices)
00211     libmesh_error_msg("ERROR: Too late.  Cannot add matrices to the system after initialization"
00212                       << "\n any more.  You should have done this earlier.");
00213 
00214   // Return the matrix if it is already there.
00215   if (this->have_matrix(mat_name))
00216     return *(_matrices[mat_name]);
00217 
00218   // Otherwise build the matrix and return it.
00219   SparseMatrix<Number>* buf = SparseMatrix<Number>::build(this->comm()).release();
00220   _matrices.insert (std::make_pair (mat_name, buf));
00221 
00222   return *buf;
00223 }
00224 
00225 
00226 
00227 const SparseMatrix<Number> * ImplicitSystem::request_matrix (const std::string& mat_name) const
00228 {
00229   // Make sure the matrix exists
00230   const_matrices_iterator pos = _matrices.find (mat_name);
00231 
00232   if (pos == _matrices.end())
00233     return NULL;
00234 
00235   return pos->second;
00236 }
00237 
00238 
00239 
00240 SparseMatrix<Number> * ImplicitSystem::request_matrix (const std::string& mat_name)
00241 {
00242   // Make sure the matrix exists
00243   matrices_iterator pos = _matrices.find (mat_name);
00244 
00245   if (pos == _matrices.end())
00246     return NULL;
00247 
00248   return pos->second;
00249 }
00250 
00251 
00252 
00253 const SparseMatrix<Number> & ImplicitSystem::get_matrix (const std::string& mat_name) const
00254 {
00255   // Make sure the matrix exists
00256   const_matrices_iterator pos = _matrices.find (mat_name);
00257 
00258   if (pos == _matrices.end())
00259     libmesh_error_msg("ERROR: matrix " << mat_name << " does not exist in this system!");
00260 
00261   return *(pos->second);
00262 }
00263 
00264 
00265 
00266 SparseMatrix<Number> & ImplicitSystem::get_matrix (const std::string& mat_name)
00267 {
00268   // Make sure the matrix exists
00269   matrices_iterator pos = _matrices.find (mat_name);
00270 
00271   if (pos == _matrices.end())
00272     libmesh_error_msg("ERROR: matrix " << mat_name << " does not exist in this system!");
00273 
00274   return *(pos->second);
00275 }
00276 
00277 
00278 
00279 void ImplicitSystem::add_system_matrix ()
00280 {
00281   // Possible that we cleared the _matrices but
00282   // forgot to NULL-out the matrix?
00283   if (_matrices.empty()) matrix = NULL;
00284 
00285 
00286   // Only need to add the matrix if it isn't there
00287   // already!
00288   if (matrix == NULL)
00289     matrix = &(this->add_matrix ("System Matrix"));
00290 
00291   libmesh_assert(matrix);
00292 }
00293 
00294 
00295 
00296 void ImplicitSystem::disable_cache () {
00297   this->assemble_before_solve = true;
00298   this->get_linear_solver()->reuse_preconditioner(false);
00299 }
00300 
00301 
00302 
00303 std::pair<unsigned int, Real>
00304 ImplicitSystem::sensitivity_solve (const ParameterVector& parameters)
00305 {
00306   // Log how long the linear solve takes.
00307   START_LOG("sensitivity_solve()", "ImplicitSystem");
00308 
00309   // The forward system should now already be solved.
00310   // Now assemble the corresponding sensitivity system.
00311 
00312   if (this->assemble_before_solve)
00313     {
00314       // Build the Jacobian
00315       this->assembly(false, true);
00316       this->matrix->close();
00317 
00318       // Reset and build the RHS from the residual derivatives
00319       this->assemble_residual_derivatives(parameters);
00320     }
00321 
00322   // The sensitivity problem is linear
00323   LinearSolver<Number> *linear_solver = this->get_linear_solver();
00324 
00325   // Our iteration counts and residuals will be sums of the individual
00326   // results
00327   std::pair<unsigned int, Real> solver_params =
00328     this->get_linear_solve_parameters();
00329   std::pair<unsigned int, Real> totalrval = std::make_pair(0,0.0);
00330 
00331   // Solve the linear system.
00332   SparseMatrix<Number> *pc = this->request_matrix("Preconditioner");
00333   for (unsigned int p=0; p != parameters.size(); ++p)
00334     {
00335       std::pair<unsigned int, Real> rval =
00336         linear_solver->solve (*matrix, pc,
00337                               this->add_sensitivity_solution(p),
00338                               this->get_sensitivity_rhs(p),
00339                               solver_params.second,
00340                               solver_params.first);
00341 
00342       totalrval.first  += rval.first;
00343       totalrval.second += rval.second;
00344     }
00345 
00346   // The linear solver may not have fit our constraints exactly
00347 #ifdef LIBMESH_ENABLE_CONSTRAINTS
00348   for (unsigned int p=0; p != parameters.size(); ++p)
00349     this->get_dof_map().enforce_constraints_exactly
00350       (*this, &this->get_sensitivity_solution(p),
00351        /* homogeneous = */ true);
00352 #endif
00353 
00354   this->release_linear_solver(linear_solver);
00355 
00356   // Stop logging the nonlinear solve
00357   STOP_LOG("sensitivity_solve()", "ImplicitSystem");
00358 
00359   return totalrval;
00360 }
00361 
00362 
00363 
00364 std::pair<unsigned int, Real>
00365 ImplicitSystem::adjoint_solve (const QoISet& qoi_indices)
00366 {
00367   // Log how long the linear solve takes.
00368   START_LOG("adjoint_solve()", "ImplicitSystem");
00369 
00370   if (this->assemble_before_solve)
00371     // Assemble the linear system
00372     this->assembly (/* get_residual = */ false,
00373                     /* get_jacobian = */ true);
00374 
00375   // The adjoint problem is linear
00376   LinearSolver<Number> *linear_solver = this->get_linear_solver();
00377 
00378   // Reset and build the RHS from the QOI derivative
00379   this->assemble_qoi_derivative(qoi_indices,
00380                                 /* include_liftfunc = */ false,
00381                                 /* apply_constraints = */ true);
00382 
00383   // Our iteration counts and residuals will be sums of the individual
00384   // results
00385   std::pair<unsigned int, Real> solver_params =
00386     this->get_linear_solve_parameters();
00387   std::pair<unsigned int, Real> totalrval = std::make_pair(0,0.0);
00388 
00389   for (unsigned int i=0; i != this->qoi.size(); ++i)
00390     if (qoi_indices.has_index(i))
00391       {
00392         const std::pair<unsigned int, Real> rval =
00393           linear_solver->adjoint_solve (*matrix, this->add_adjoint_solution(i),
00394                                         this->get_adjoint_rhs(i),
00395                                         solver_params.second,
00396                                         solver_params.first);
00397 
00398         totalrval.first  += rval.first;
00399         totalrval.second += rval.second;
00400       }
00401 
00402   this->release_linear_solver(linear_solver);
00403 
00404   // The linear solver may not have fit our constraints exactly
00405 #ifdef LIBMESH_ENABLE_CONSTRAINTS
00406   for (unsigned int i=0; i != this->qoi.size(); ++i)
00407     if (qoi_indices.has_index(i))
00408       this->get_dof_map().enforce_adjoint_constraints_exactly
00409         (this->get_adjoint_solution(i), i);
00410 #endif
00411 
00412   // Stop logging the nonlinear solve
00413   STOP_LOG("adjoint_solve()", "ImplicitSystem");
00414 
00415   return totalrval;
00416 }
00417 
00418 
00419 
00420 std::pair<unsigned int, Real>
00421 ImplicitSystem::weighted_sensitivity_adjoint_solve (const ParameterVector& parameters_in,
00422                                                     const ParameterVector& weights,
00423                                                     const QoISet& qoi_indices)
00424 {
00425   // Log how long the linear solve takes.
00426   START_LOG("weighted_sensitivity_adjoint_solve()", "ImplicitSystem");
00427 
00428   // We currently get partial derivatives via central differencing
00429   const Real delta_p = TOLERANCE;
00430 
00431   ParameterVector& parameters =
00432     const_cast<ParameterVector&>(parameters_in);
00433 
00434   // The forward system should now already be solved.
00435   // The adjoint system should now already be solved.
00436   // Now we're assembling a weighted sum of adjoint-adjoint systems:
00437   //
00438   // dR/du (u, sum_l(w_l*z^l)) = sum_l(w_l*(Q''_ul - R''_ul (u, z)))
00439 
00440   // FIXME: The derivation here does not yet take adjoint boundary
00441   // conditions into account.
00442   for (unsigned int i=0; i != this->qoi.size(); ++i)
00443     if (qoi_indices.has_index(i))
00444       libmesh_assert(!this->get_dof_map().has_adjoint_dirichlet_boundaries(i));
00445 
00446   // We'll assemble the rhs first, because the R'' term will require
00447   // perturbing the jacobian
00448 
00449   // We'll use temporary rhs vectors, because we haven't (yet) found
00450   // any good reasons why users might want to save these:
00451 
00452   std::vector<NumericVector<Number> *> temprhs(this->qoi.size());
00453   for (unsigned int i=0; i != this->qoi.size(); ++i)
00454     if (qoi_indices.has_index(i))
00455       temprhs[i] = this->rhs->zero_clone().release();
00456 
00457   // We approximate the _l partial derivatives via a central
00458   // differencing perturbation in the w_l direction:
00459   //
00460   // sum_l(w_l*v_l) ~= (v(p + dp*w_l*e_l) - v(p - dp*w_l*e_l))/(2*dp)
00461 
00462   // PETSc doesn't implement SGEMX, so neither does NumericVector,
00463   // so we want to avoid calculating f -= R'*z.  We'll thus evaluate
00464   // the above equation by first adding -v(p+dp...), then multiplying
00465   // the intermediate result vectors by -1, then adding -v(p-dp...),
00466   // then finally dividing by 2*dp.
00467 
00468   ParameterVector oldparameters, parameterperturbation;
00469   parameters.deep_copy(oldparameters);
00470   weights.deep_copy(parameterperturbation);
00471   parameterperturbation *= delta_p;
00472   parameters += parameterperturbation;
00473 
00474   this->assembly(false, true);
00475   this->matrix->close();
00476 
00477   // Take the discrete adjoint, so that we can calculate R_u(u,z) with
00478   // a matrix-vector product of R_u and z.
00479   matrix->get_transpose(*matrix);
00480 
00481   this->assemble_qoi_derivative(qoi_indices,
00482                                 /* include_liftfunc = */ false,
00483                                 /* apply_constraints = */ true);
00484   for (unsigned int i=0; i != this->qoi.size(); ++i)
00485     if (qoi_indices.has_index(i))
00486       {
00487         this->get_adjoint_rhs(i).close();
00488         *(temprhs[i]) -= this->get_adjoint_rhs(i);
00489         this->matrix->vector_mult_add(*(temprhs[i]), this->get_adjoint_solution(i));
00490         *(temprhs[i]) *= -1.0;
00491       }
00492 
00493   oldparameters.value_copy(parameters);
00494   parameterperturbation *= -1.0;
00495   parameters += parameterperturbation;
00496 
00497   this->assembly(false, true);
00498   this->matrix->close();
00499   matrix->get_transpose(*matrix);
00500 
00501   this->assemble_qoi_derivative(qoi_indices,
00502                                 /* include_liftfunc = */ false,
00503                                 /* apply_constraints = */ true);
00504   for (unsigned int i=0; i != this->qoi.size(); ++i)
00505     if (qoi_indices.has_index(i))
00506       {
00507         this->get_adjoint_rhs(i).close();
00508         *(temprhs[i]) -= this->get_adjoint_rhs(i);
00509         this->matrix->vector_mult_add(*(temprhs[i]), this->get_adjoint_solution(i));
00510         *(temprhs[i]) /= (2.0*delta_p);
00511       }
00512 
00513   // Finally, assemble the jacobian at the non-perturbed parameter
00514   // values.  Ignore assemble_before_solve; if we had a good
00515   // non-perturbed matrix before we've already overwritten it.
00516   oldparameters.value_copy(parameters);
00517 
00518   // if (this->assemble_before_solve)
00519   {
00520     // Build the Jacobian
00521     this->assembly(false, true);
00522     this->matrix->close();
00523 
00524     // Take the discrete adjoint
00525     matrix->get_transpose(*matrix);
00526   }
00527 
00528   // The weighted adjoint-adjoint problem is linear
00529   LinearSolver<Number> *linear_solver = this->get_linear_solver();
00530 
00531   // Our iteration counts and residuals will be sums of the individual
00532   // results
00533   std::pair<unsigned int, Real> solver_params =
00534     this->get_linear_solve_parameters();
00535   std::pair<unsigned int, Real> totalrval = std::make_pair(0,0.0);
00536 
00537   for (unsigned int i=0; i != this->qoi.size(); ++i)
00538     if (qoi_indices.has_index(i))
00539       {
00540         const std::pair<unsigned int, Real> rval =
00541           linear_solver->solve (*matrix, this->add_weighted_sensitivity_adjoint_solution(i),
00542                                 *(temprhs[i]),
00543                                 solver_params.second,
00544                                 solver_params.first);
00545 
00546         totalrval.first  += rval.first;
00547         totalrval.second += rval.second;
00548       }
00549 
00550   this->release_linear_solver(linear_solver);
00551 
00552   for (unsigned int i=0; i != this->qoi.size(); ++i)
00553     if (qoi_indices.has_index(i))
00554       delete temprhs[i];
00555 
00556   // The linear solver may not have fit our constraints exactly
00557 #ifdef LIBMESH_ENABLE_CONSTRAINTS
00558   for (unsigned int i=0; i != this->qoi.size(); ++i)
00559     if (qoi_indices.has_index(i))
00560       this->get_dof_map().enforce_constraints_exactly
00561         (*this, &this->get_weighted_sensitivity_adjoint_solution(i),
00562          /* homogeneous = */ true);
00563 #endif
00564 
00565   // Stop logging the nonlinear solve
00566   STOP_LOG("weighted_sensitivity_adjoint_solve()", "ImplicitSystem");
00567 
00568   return totalrval;
00569 }
00570 
00571 
00572 
00573 std::pair<unsigned int, Real>
00574 ImplicitSystem::weighted_sensitivity_solve (const ParameterVector& parameters_in,
00575                                             const ParameterVector& weights)
00576 {
00577   // Log how long the linear solve takes.
00578   START_LOG("weighted_sensitivity_solve()", "ImplicitSystem");
00579 
00580   // We currently get partial derivatives via central differencing
00581   const Real delta_p = TOLERANCE;
00582 
00583   ParameterVector& parameters =
00584     const_cast<ParameterVector&>(parameters_in);
00585 
00586   // The forward system should now already be solved.
00587 
00588   // Now we're assembling a weighted sum of sensitivity systems:
00589   //
00590   // dR/du (u, v)(sum(w_l*u'_l)) = -sum_l(w_l*R'_l (u, v)) forall v
00591 
00592   // We'll assemble the rhs first, because the R' term will require
00593   // perturbing the system, and some applications may not be able to
00594   // assemble a perturbed residual without simultaneously constructing
00595   // a perturbed jacobian.
00596 
00597   // We approximate the _l partial derivatives via a central
00598   // differencing perturbation in the w_l direction:
00599   //
00600   // sum_l(w_l*v_l) ~= (v(p + dp*w_l*e_l) - v(p - dp*w_l*e_l))/(2*dp)
00601 
00602   ParameterVector oldparameters, parameterperturbation;
00603   parameters.deep_copy(oldparameters);
00604   weights.deep_copy(parameterperturbation);
00605   parameterperturbation *= delta_p;
00606   parameters += parameterperturbation;
00607 
00608   this->assembly(true, false, true);
00609   this->rhs->close();
00610 
00611   UniquePtr<NumericVector<Number> > temprhs = this->rhs->clone();
00612 
00613   oldparameters.value_copy(parameters);
00614   parameterperturbation *= -1.0;
00615   parameters += parameterperturbation;
00616 
00617   this->assembly(true, false, true);
00618   this->rhs->close();
00619 
00620   *temprhs -= *(this->rhs);
00621   *temprhs /= (2.0*delta_p);
00622 
00623   // Finally, assemble the jacobian at the non-perturbed parameter
00624   // values
00625   oldparameters.value_copy(parameters);
00626 
00627   // Build the Jacobian
00628   this->assembly(false, true);
00629   this->matrix->close();
00630 
00631   // The weighted sensitivity problem is linear
00632   LinearSolver<Number> *linear_solver = this->get_linear_solver();
00633 
00634   std::pair<unsigned int, Real> solver_params =
00635     this->get_linear_solve_parameters();
00636 
00637   const std::pair<unsigned int, Real> rval =
00638     linear_solver->solve (*matrix, this->add_weighted_sensitivity_solution(),
00639                           *temprhs,
00640                           solver_params.second,
00641                           solver_params.first);
00642 
00643   this->release_linear_solver(linear_solver);
00644 
00645   // The linear solver may not have fit our constraints exactly
00646 #ifdef LIBMESH_ENABLE_CONSTRAINTS
00647   this->get_dof_map().enforce_constraints_exactly
00648     (*this, &this->get_weighted_sensitivity_solution(),
00649      /* homogeneous = */ true);
00650 #endif
00651 
00652   // Stop logging the nonlinear solve
00653   STOP_LOG("weighted_sensitivity_solve()", "ImplicitSystem");
00654 
00655   return rval;
00656 }
00657 
00658 
00659 
00660 void ImplicitSystem::assemble_residual_derivatives(const ParameterVector& parameters_in)
00661 {
00662   Real deltap = TOLERANCE;
00663 
00664   ParameterVector& parameters =
00665     const_cast<ParameterVector&>(parameters_in);
00666 
00667   const unsigned int Np = cast_int<unsigned int>
00668     (parameters.size());
00669 
00670   for (unsigned int p=0; p != Np; ++p)
00671     {
00672       NumericVector<Number> &sensitivity_rhs = this->add_sensitivity_rhs(p);
00673 
00674       // Approximate -(partial R / partial p) by
00675       // (R(p-dp) - R(p+dp)) / (2*dp)
00676 
00677       Number old_parameter = *parameters[p];
00678       *parameters[p] -= deltap;
00679 
00680       this->assembly(true, false, true);
00681       this->rhs->close();
00682       sensitivity_rhs = *this->rhs;
00683 
00684       *parameters[p] = old_parameter + deltap;
00685 
00686       this->assembly(true, false, true);
00687       this->rhs->close();
00688 
00689       sensitivity_rhs -= *this->rhs;
00690       sensitivity_rhs /= (2*deltap);
00691       sensitivity_rhs.close();
00692 
00693       *parameters[p] = old_parameter;
00694     }
00695 }
00696 
00697 
00698 
00699 void ImplicitSystem::adjoint_qoi_parameter_sensitivity
00700 (const QoISet&          qoi_indices,
00701  const ParameterVector& parameters_in,
00702  SensitivityData&       sensitivities)
00703 {
00704   // We currently get partial derivatives via central differencing
00705   const Real delta_p = TOLERANCE;
00706 
00707   ParameterVector& parameters =
00708     const_cast<ParameterVector&>(parameters_in);
00709 
00710   const unsigned int Np = cast_int<unsigned int>
00711     (parameters.size());
00712   const unsigned int Nq = cast_int<unsigned int>
00713     (qoi.size());
00714 
00715   // An introduction to the problem:
00716   //
00717   // Residual R(u(p),p) = 0
00718   // partial R / partial u = J = system matrix
00719   //
00720   // This implies that:
00721   // d/dp(R) = 0
00722   // (partial R / partial p) +
00723   // (partial R / partial u) * (partial u / partial p) = 0
00724 
00725   // We first do an adjoint solve:
00726   // J^T * z = (partial q / partial u)
00727   // if we havent already or dont have an initial condition for the adjoint
00728   if (!this->is_adjoint_already_solved())
00729     {
00730       this->adjoint_solve(qoi_indices);
00731     }
00732 
00733   // Get ready to fill in senstivities:
00734   sensitivities.allocate_data(qoi_indices, *this, parameters);
00735 
00736   // We use the identities:
00737   // dq/dp = (partial q / partial p) + (partial q / partial u) *
00738   //         (partial u / partial p)
00739   // dq/dp = (partial q / partial p) + (J^T * z) *
00740   //         (partial u / partial p)
00741   // dq/dp = (partial q / partial p) + z * J *
00742   //         (partial u / partial p)
00743 
00744   // Leading to our final formula:
00745   // dq/dp = (partial q / partial p) - z * (partial R / partial p)
00746 
00747   // In the case of adjoints with heterogenous Dirichlet boundary
00748   // function phi, where
00749   // q := R(u,phi) + S(u)
00750   // the final formula works out to:
00751   // dq/dp = (partial S / partial p) - z * (partial R / partial p)
00752   // Because we currently have no direct access to
00753   // (partial S / partial p), we use the identity
00754   // (partial S / partial p) = (partial q / partial p) -
00755   //                           phi * (partial R / partial p)
00756   // to derive an equivalent equation:
00757   // dq/dp = (partial q / partial p) - (z+phi) * (partial R / partial p)
00758 
00759 
00760   // If we have non-zero adjoint dofs on Dirichlet constrained
00761   // boundary dofs, then we need the residual components
00762   // corresponding to those dofs when using r*z to compute R(u,z), so
00763   // we can't apply constraints.
00764   //
00765   // If we aren't in that situation we could apply constraints but
00766   // it will be faster not to.
00767 
00768   this->get_dof_map().stash_dof_constraints();
00769 
00770   for (unsigned int j=0; j != Np; ++j)
00771     {
00772       // (partial q / partial p) ~= (q(p+dp)-q(p-dp))/(2*dp)
00773       // (partial R / partial p) ~= (rhs(p+dp) - rhs(p-dp))/(2*dp)
00774 
00775       Number old_parameter = *parameters[j];
00776       // Number old_qoi = this->qoi;
00777 
00778       *parameters[j] = old_parameter - delta_p;
00779       this->assemble_qoi(qoi_indices);
00780       std::vector<Number> qoi_minus = this->qoi;
00781 
00782       this->assembly(true, false, true);
00783       this->rhs->close();
00784 
00785       // FIXME - this can and should be optimized to avoid the clone()
00786       UniquePtr<NumericVector<Number> > partialR_partialp = this->rhs->clone();
00787       *partialR_partialp *= -1;
00788 
00789       *parameters[j] = old_parameter + delta_p;
00790       this->assemble_qoi(qoi_indices);
00791       std::vector<Number>& qoi_plus = this->qoi;
00792 
00793       std::vector<Number> partialq_partialp(Nq, 0);
00794       for (unsigned int i=0; i != Nq; ++i)
00795         if (qoi_indices.has_index(i))
00796           partialq_partialp[i] = (qoi_plus[i] - qoi_minus[i]) / (2.*delta_p);
00797 
00798       this->assembly(true, false, true);
00799       this->rhs->close();
00800       *partialR_partialp += *this->rhs;
00801       *partialR_partialp /= (2.*delta_p);
00802 
00803       // Don't leave the parameter changed
00804       *parameters[j] = old_parameter;
00805 
00806       for (unsigned int i=0; i != Nq; ++i)
00807         if (qoi_indices.has_index(i))
00808           {
00809             sensitivities[i][j] = partialq_partialp[i] -
00810               partialR_partialp->dot(this->get_adjoint_solution(i));
00811 
00812             if (this->get_dof_map().has_adjoint_dirichlet_boundaries(i))
00813               {
00814                 UniquePtr<NumericVector<Number> > lift_func =
00815                   this->get_adjoint_solution(i).zero_clone();
00816                 this->get_dof_map().enforce_adjoint_constraints_exactly
00817                   (*lift_func.get(), i);
00818                 sensitivities[i][j] += partialR_partialp->dot(*lift_func);
00819               }
00820           }
00821     }
00822 
00823   // All parameters have been reset.
00824   // We didn't cache the original rhs or matrix for memory reasons,
00825   // but we can restore them to a state consistent solution -
00826   // principle of least surprise.
00827 
00828   this->get_dof_map().unstash_dof_constraints();
00829   this->assembly(true, true);
00830   this->rhs->close();
00831   this->matrix->close();
00832   this->assemble_qoi(qoi_indices);
00833 }
00834 
00835 
00836 
00837 void ImplicitSystem::forward_qoi_parameter_sensitivity
00838 (const QoISet&          qoi_indices,
00839  const ParameterVector& parameters_in,
00840  SensitivityData&       sensitivities)
00841 {
00842   // We currently get partial derivatives via central differencing
00843   const Real delta_p = TOLERANCE;
00844 
00845   ParameterVector& parameters =
00846     const_cast<ParameterVector&>(parameters_in);
00847 
00848   const unsigned int Np = cast_int<unsigned int>
00849     (parameters.size());
00850   const unsigned int Nq = cast_int<unsigned int>
00851     (qoi.size());
00852 
00853   // An introduction to the problem:
00854   //
00855   // Residual R(u(p),p) = 0
00856   // partial R / partial u = J = system matrix
00857   //
00858   // This implies that:
00859   // d/dp(R) = 0
00860   // (partial R / partial p) +
00861   // (partial R / partial u) * (partial u / partial p) = 0
00862 
00863   // We first solve for (partial u / partial p) for each parameter:
00864   // J * (partial u / partial p) = - (partial R / partial p)
00865 
00866   this->sensitivity_solve(parameters);
00867 
00868   // Get ready to fill in senstivities:
00869   sensitivities.allocate_data(qoi_indices, *this, parameters);
00870 
00871   // We use the identity:
00872   // dq/dp = (partial q / partial p) + (partial q / partial u) *
00873   //         (partial u / partial p)
00874 
00875   // We get (partial q / partial u) from the user
00876   this->assemble_qoi_derivative(qoi_indices,
00877                                 /* include_liftfunc = */ true,
00878                                 /* apply_constraints = */ false);
00879 
00880   // FIXME: what do we do with adjoint boundary conditions here?
00881 
00882   // We don't need these to be closed() in this function, but libMesh
00883   // standard practice is to have them closed() by the time the
00884   // function exits
00885   for (unsigned int i=0; i != this->qoi.size(); ++i)
00886     if (qoi_indices.has_index(i))
00887       this->get_adjoint_rhs(i).close();
00888 
00889   for (unsigned int j=0; j != Np; ++j)
00890     {
00891       // (partial q / partial p) ~= (q(p+dp)-q(p-dp))/(2*dp)
00892 
00893       Number old_parameter = *parameters[j];
00894 
00895       *parameters[j] = old_parameter - delta_p;
00896       this->assemble_qoi();
00897       std::vector<Number> qoi_minus = this->qoi;
00898 
00899       *parameters[j] = old_parameter + delta_p;
00900       this->assemble_qoi();
00901       std::vector<Number>& qoi_plus = this->qoi;
00902 
00903       std::vector<Number> partialq_partialp(Nq, 0);
00904       for (unsigned int i=0; i != Nq; ++i)
00905         if (qoi_indices.has_index(i))
00906           partialq_partialp[i] = (qoi_plus[i] - qoi_minus[i]) / (2.*delta_p);
00907 
00908       // Don't leave the parameter changed
00909       *parameters[j] = old_parameter;
00910 
00911       for (unsigned int i=0; i != Nq; ++i)
00912         if (qoi_indices.has_index(i))
00913           sensitivities[i][j] = partialq_partialp[i] +
00914             this->get_adjoint_rhs(i).dot(this->get_sensitivity_solution(j));
00915     }
00916 
00917   // All parameters have been reset.
00918   // We didn't cache the original rhs or matrix for memory reasons,
00919   // but we can restore them to a state consistent solution -
00920   // principle of least surprise.
00921   this->assembly(true, true);
00922   this->rhs->close();
00923   this->matrix->close();
00924   this->assemble_qoi(qoi_indices);
00925 }
00926 
00927 
00928 
00929 void ImplicitSystem::qoi_parameter_hessian_vector_product
00930 (const QoISet& qoi_indices,
00931  const ParameterVector& parameters_in,
00932  const ParameterVector& vector,
00933  SensitivityData& sensitivities)
00934 {
00935   // We currently get partial derivatives via finite differencing
00936   const Real delta_p = TOLERANCE;
00937 
00938   ParameterVector& parameters =
00939     const_cast<ParameterVector&>(parameters_in);
00940 
00941   // We'll use a single temporary vector for matrix-vector-vector products
00942   UniquePtr<NumericVector<Number> > tempvec = this->solution->zero_clone();
00943 
00944   const unsigned int Np = cast_int<unsigned int>
00945     (parameters.size());
00946   const unsigned int Nq = cast_int<unsigned int>
00947     (qoi.size());
00948 
00949   // For each quantity of interest q, the parameter sensitivity
00950   // Hessian is defined as q''_{kl} = {d^2 q}/{d p_k d p_l}.
00951   // Given a vector of parameter perturbation weights w_l, this
00952   // function evaluates the hessian-vector product sum_l(q''_{kl}*w_l)
00953   //
00954   // We calculate it from values and partial derivatives of the
00955   // quantity of interest function Q, solution u, adjoint solution z,
00956   // parameter sensitivity adjoint solutions z^l, and residual R, as:
00957   //
00958   // sum_l(q''_{kl}*w_l) =
00959   // sum_l(w_l * Q''_{kl}) + Q''_{uk}(u)*(sum_l(w_l u'_l)) -
00960   // R'_k(u, sum_l(w_l*z^l)) - R'_{uk}(u,z)*(sum_l(w_l u'_l) -
00961   // sum_l(w_l*R''_{kl}(u,z))
00962   //
00963   // See the adjoints model document for more details.
00964 
00965   // We first do an adjoint solve to get z for each quantity of
00966   // interest
00967   // if we havent already or dont have an initial condition for the adjoint
00968   if (!this->is_adjoint_already_solved())
00969     {
00970       this->adjoint_solve(qoi_indices);
00971     }
00972 
00973   // Get ready to fill in senstivities:
00974   sensitivities.allocate_data(qoi_indices, *this, parameters);
00975 
00976   // We can't solve for all the solution sensitivities u'_l or for all
00977   // of the parameter sensitivity adjoint solutions z^l without
00978   // requiring O(Nq*Np) linear solves.  So we'll solve directly for their
00979   // weighted sum - this is just O(Nq) solves.
00980 
00981   // First solve for sum_l(w_l u'_l).
00982   this->weighted_sensitivity_solve(parameters, vector);
00983 
00984   // Then solve for sum_l(w_l z^l).
00985   this->weighted_sensitivity_adjoint_solve(parameters, vector, qoi_indices);
00986 
00987   for (unsigned int k=0; k != Np; ++k)
00988     {
00989       // We approximate sum_l(w_l * Q''_{kl}) with a central
00990       // differencing pertubation:
00991       // sum_l(w_l * Q''_{kl}) ~=
00992       // (Q(p + dp*w_l*e_l + dp*e_k) - Q(p - dp*w_l*e_l + dp*e_k) -
00993       // Q(p + dp*w_l*e_l - dp*e_k) + Q(p - dp*w_l*e_l - dp*e_k))/(4*dp^2)
00994 
00995       // The sum(w_l*R''_kl) term requires the same sort of pertubation,
00996       // and so we subtract it in at the same time:
00997       // sum_l(w_l * R''_{kl}) ~=
00998       // (R(p + dp*w_l*e_l + dp*e_k) - R(p - dp*w_l*e_l + dp*e_k) -
00999       // R(p + dp*w_l*e_l - dp*e_k) + R(p - dp*w_l*e_l - dp*e_k))/(4*dp^2)
01000 
01001       ParameterVector oldparameters, parameterperturbation;
01002       parameters.deep_copy(oldparameters);
01003       vector.deep_copy(parameterperturbation);
01004       parameterperturbation *= delta_p;
01005       parameters += parameterperturbation;
01006 
01007       Number old_parameter = *parameters[k];
01008 
01009       *parameters[k] = old_parameter + delta_p;
01010       this->assemble_qoi(qoi_indices);
01011       this->assembly(true, false, true);
01012       this->rhs->close();
01013       std::vector<Number> partial2q_term = this->qoi;
01014       std::vector<Number> partial2R_term(this->qoi.size());
01015       for (unsigned int i=0; i != Nq; ++i)
01016         if (qoi_indices.has_index(i))
01017           partial2R_term[i] = this->rhs->dot(this->get_adjoint_solution(i));
01018 
01019       *parameters[k] = old_parameter - delta_p;
01020       this->assemble_qoi(qoi_indices);
01021       this->assembly(true, false, true);
01022       this->rhs->close();
01023       for (unsigned int i=0; i != Nq; ++i)
01024         if (qoi_indices.has_index(i))
01025           {
01026             partial2q_term[i] -= this->qoi[i];
01027             partial2R_term[i] -= this->rhs->dot(this->get_adjoint_solution(i));
01028           }
01029 
01030       oldparameters.value_copy(parameters);
01031       parameterperturbation *= -1.0;
01032       parameters += parameterperturbation;
01033 
01034       // Re-center old_parameter, which may be affected by vector
01035       old_parameter = *parameters[k];
01036 
01037       *parameters[k] = old_parameter + delta_p;
01038       this->assemble_qoi(qoi_indices);
01039       this->assembly(true, false, true);
01040       this->rhs->close();
01041       for (unsigned int i=0; i != Nq; ++i)
01042         if (qoi_indices.has_index(i))
01043           {
01044             partial2q_term[i] -= this->qoi[i];
01045             partial2R_term[i] -= this->rhs->dot(this->get_adjoint_solution(i));
01046           }
01047 
01048       *parameters[k] = old_parameter - delta_p;
01049       this->assemble_qoi(qoi_indices);
01050       this->assembly(true, false, true);
01051       this->rhs->close();
01052       for (unsigned int i=0; i != Nq; ++i)
01053         if (qoi_indices.has_index(i))
01054           {
01055             partial2q_term[i] += this->qoi[i];
01056             partial2R_term[i] += this->rhs->dot(this->get_adjoint_solution(i));
01057           }
01058 
01059       for (unsigned int i=0; i != Nq; ++i)
01060         if (qoi_indices.has_index(i))
01061           {
01062             partial2q_term[i] /= (4. * delta_p * delta_p);
01063             partial2R_term[i] /= (4. * delta_p * delta_p);
01064           }
01065 
01066       for (unsigned int i=0; i != Nq; ++i)
01067         if (qoi_indices.has_index(i))
01068           sensitivities[i][k] = partial2q_term[i] - partial2R_term[i];
01069 
01070       // We get (partial q / partial u), R, and
01071       // (partial R / partial u) from the user, but centrally
01072       // difference to get q_uk, R_k, and R_uk terms:
01073       // (partial R / partial k)
01074       // R_k*sum(w_l*z^l) = (R(p+dp*e_k)*sum(w_l*z^l) - R(p-dp*e_k)*sum(w_l*z^l))/(2*dp)
01075       // (partial^2 q / partial u partial k)
01076       // q_uk = (q_u(p+dp*e_k) - q_u(p-dp*e_k))/(2*dp)
01077       // (partial^2 R / partial u partial k)
01078       // R_uk*z*sum(w_l*u'_l) = (R_u(p+dp*e_k)*z*sum(w_l*u'_l) - R_u(p-dp*e_k)*z*sum(w_l*u'_l))/(2*dp)
01079 
01080       // To avoid creating Nq temporary vectors for q_uk or R_uk, we add
01081       // subterms to the sensitivities output one by one.
01082       //
01083       // FIXME: this is probably a bad order of operations for
01084       // controlling floating point error.
01085 
01086       *parameters[k] = old_parameter + delta_p;
01087       this->assembly(true, true);
01088       this->rhs->close();
01089       this->matrix->close();
01090       this->assemble_qoi_derivative(qoi_indices,
01091                                     /* include_liftfunc = */ true,
01092                                     /* apply_constraints = */ false);
01093 
01094       this->matrix->vector_mult(*tempvec, this->get_weighted_sensitivity_solution());
01095 
01096       for (unsigned int i=0; i != Nq; ++i)
01097         if (qoi_indices.has_index(i))
01098           {
01099             this->get_adjoint_rhs(i).close();
01100             sensitivities[i][k] += (this->get_adjoint_rhs(i).dot(this->get_weighted_sensitivity_solution()) -
01101                                     this->rhs->dot(this->get_weighted_sensitivity_adjoint_solution(i)) -
01102                                     this->get_adjoint_solution(i).dot(*tempvec)) / (2.*delta_p);
01103           }
01104 
01105       *parameters[k] = old_parameter - delta_p;
01106       this->assembly(true, true);
01107       this->rhs->close();
01108       this->matrix->close();
01109       this->assemble_qoi_derivative(qoi_indices,
01110                                     /* include_liftfunc = */ true,
01111                                     /* apply_constraints = */ false);
01112 
01113       this->matrix->vector_mult(*tempvec, this->get_weighted_sensitivity_solution());
01114 
01115       for (unsigned int i=0; i != Nq; ++i)
01116         if (qoi_indices.has_index(i))
01117           {
01118             this->get_adjoint_rhs(i).close();
01119             sensitivities[i][k] += (-this->get_adjoint_rhs(i).dot(this->get_weighted_sensitivity_solution()) +
01120                                     this->rhs->dot(this->get_weighted_sensitivity_adjoint_solution(i)) +
01121                                     this->get_adjoint_solution(i).dot(*tempvec)) / (2.*delta_p);
01122           }
01123     }
01124 
01125   // All parameters have been reset.
01126   // Don't leave the qoi or system changed - principle of least
01127   // surprise.
01128   this->assembly(true, true);
01129   this->rhs->close();
01130   this->matrix->close();
01131   this->assemble_qoi(qoi_indices);
01132 }
01133 
01134 
01135 
01136 void ImplicitSystem::qoi_parameter_hessian
01137 (const QoISet& qoi_indices,
01138  const ParameterVector& parameters_in,
01139  SensitivityData& sensitivities)
01140 {
01141   // We currently get partial derivatives via finite differencing
01142   const Real delta_p = TOLERANCE;
01143 
01144   ParameterVector& parameters =
01145     const_cast<ParameterVector&>(parameters_in);
01146 
01147   // We'll use one temporary vector for matrix-vector-vector products
01148   UniquePtr<NumericVector<Number> > tempvec = this->solution->zero_clone();
01149 
01150   // And another temporary vector to hold a copy of the true solution
01151   // so we can safely perturb this->solution.
01152   UniquePtr<NumericVector<Number> > oldsolution = this->solution->clone();
01153 
01154   const unsigned int Np = cast_int<unsigned int>
01155     (parameters.size());
01156   const unsigned int Nq = cast_int<unsigned int>
01157     (qoi.size());
01158 
01159   // For each quantity of interest q, the parameter sensitivity
01160   // Hessian is defined as q''_{kl} = {d^2 q}/{d p_k d p_l}.
01161   //
01162   // We calculate it from values and partial derivatives of the
01163   // quantity of interest function Q, solution u, adjoint solution z,
01164   // and residual R, as:
01165   //
01166   // q''_{kl} =
01167   // Q''_{kl} + Q''_{uk}(u)*u'_l + Q''_{ul}(u) * u'_k +
01168   // Q''_{uu}(u)*u'_k*u'_l -
01169   // R''_{kl}(u,z) -
01170   // R''_{uk}(u,z)*u'_l - R''_{ul}(u,z)*u'_k -
01171   // R''_{uu}(u,z)*u'_k*u'_l
01172   //
01173   // See the adjoints model document for more details.
01174 
01175   // We first do an adjoint solve to get z for each quantity of
01176   // interest
01177   // if we havent already or dont have an initial condition for the adjoint
01178   if (!this->is_adjoint_already_solved())
01179     {
01180       this->adjoint_solve(qoi_indices);
01181     }
01182 
01183   // And a sensitivity solve to get u_k for each parameter
01184   this->sensitivity_solve(parameters);
01185 
01186   // Get ready to fill in second derivatives:
01187   sensitivities.allocate_hessian_data(qoi_indices, *this, parameters);
01188 
01189   for (unsigned int k=0; k != Np; ++k)
01190     {
01191       Number old_parameterk = *parameters[k];
01192 
01193       // The Hessian is symmetric, so we just calculate the lower
01194       // triangle and the diagonal, and we get the upper triangle from
01195       // the transpose of the lower
01196 
01197       for (unsigned int l=0; l != k+1; ++l)
01198         {
01199           // The second partial derivatives with respect to parameters
01200           // are all calculated via a central finite difference
01201           // stencil:
01202           // F''_{kl} ~= (F(p+dp*e_k+dp*e_l) - F(p+dp*e_k-dp*e_l) -
01203           //              F(p-dp*e_k+dp*e_l) + F(p-dp*e_k-dp*e_l))/(4*dp^2)
01204           // We will add Q''_{kl}(u) and subtract R''_{kl}(u,z) at the
01205           // same time.
01206           //
01207           // We have to be careful with the perturbations to handle
01208           // the k=l case
01209 
01210           Number old_parameterl = *parameters[l];
01211 
01212           *parameters[k] += delta_p;
01213           *parameters[l] += delta_p;
01214           this->assemble_qoi(qoi_indices);
01215           this->assembly(true, false, true);
01216           this->rhs->close();
01217           std::vector<Number> partial2q_term = this->qoi;
01218           std::vector<Number> partial2R_term(this->qoi.size());
01219           for (unsigned int i=0; i != Nq; ++i)
01220             if (qoi_indices.has_index(i))
01221               partial2R_term[i] = this->rhs->dot(this->get_adjoint_solution(i));
01222 
01223           *parameters[l] -= 2.*delta_p;
01224           this->assemble_qoi(qoi_indices);
01225           this->assembly(true, false, true);
01226           this->rhs->close();
01227           for (unsigned int i=0; i != Nq; ++i)
01228             if (qoi_indices.has_index(i))
01229               {
01230                 partial2q_term[i] -= this->qoi[i];
01231                 partial2R_term[i] -= this->rhs->dot(this->get_adjoint_solution(i));
01232               }
01233 
01234           *parameters[k] -= 2.*delta_p;
01235           this->assemble_qoi(qoi_indices);
01236           this->assembly(true, false, true);
01237           this->rhs->close();
01238           for (unsigned int i=0; i != Nq; ++i)
01239             if (qoi_indices.has_index(i))
01240               {
01241                 partial2q_term[i] += this->qoi[i];
01242                 partial2R_term[i] += this->rhs->dot(this->get_adjoint_solution(i));
01243               }
01244 
01245           *parameters[l] += 2.*delta_p;
01246           this->assemble_qoi(qoi_indices);
01247           this->assembly(true, false, true);
01248           this->rhs->close();
01249           for (unsigned int i=0; i != Nq; ++i)
01250             if (qoi_indices.has_index(i))
01251               {
01252                 partial2q_term[i] -= this->qoi[i];
01253                 partial2R_term[i] -= this->rhs->dot(this->get_adjoint_solution(i));
01254                 partial2q_term[i] /= (4. * delta_p * delta_p);
01255                 partial2R_term[i] /= (4. * delta_p * delta_p);
01256               }
01257 
01258           for (unsigned int i=0; i != Nq; ++i)
01259             if (qoi_indices.has_index(i))
01260               {
01261                 Number current_terms = partial2q_term[i] - partial2R_term[i];
01262                 sensitivities.second_derivative(i,k,l) += current_terms;
01263                 if (k != l)
01264                   sensitivities.second_derivative(i,l,k) += current_terms;
01265               }
01266 
01267           // Don't leave the parameters perturbed
01268           *parameters[l] = old_parameterl;
01269           *parameters[k] = old_parameterk;
01270         }
01271 
01272       // We get (partial q / partial u) and
01273       // (partial R / partial u) from the user, but centrally
01274       // difference to get q_uk and R_uk terms:
01275       // (partial^2 q / partial u partial k)
01276       // q_uk*u'_l = (q_u(p+dp*e_k)*u'_l - q_u(p-dp*e_k)*u'_l)/(2*dp)
01277       // R_uk*z*u'_l = (R_u(p+dp*e_k)*z*u'_l - R_u(p-dp*e_k)*z*u'_l)/(2*dp)
01278       //
01279       // To avoid creating Nq temporary vectors, we add these
01280       // subterms to the sensitivities output one by one.
01281       //
01282       // FIXME: this is probably a bad order of operations for
01283       // controlling floating point error.
01284 
01285       *parameters[k] = old_parameterk + delta_p;
01286       this->assembly(false, true);
01287       this->matrix->close();
01288       this->assemble_qoi_derivative(qoi_indices,
01289                                     /* include_liftfunc = */ true,
01290                                     /* apply_constraints = */ false);
01291 
01292       for (unsigned int l=0; l != Np; ++l)
01293         {
01294           this->matrix->vector_mult(*tempvec, this->get_sensitivity_solution(l));
01295           for (unsigned int i=0; i != Nq; ++i)
01296             if (qoi_indices.has_index(i))
01297               {
01298                 this->get_adjoint_rhs(i).close();
01299                 Number current_terms =
01300                   (this->get_adjoint_rhs(i).dot(this->get_sensitivity_solution(l)) -
01301                    tempvec->dot(this->get_adjoint_solution(i))) / (2.*delta_p);
01302                 sensitivities.second_derivative(i,k,l) += current_terms;
01303 
01304                 // We use the _uk terms twice; symmetry lets us reuse
01305                 // these calculations for the _ul terms.
01306 
01307                 sensitivities.second_derivative(i,l,k) += current_terms;
01308               }
01309         }
01310 
01311       *parameters[k] = old_parameterk - delta_p;
01312       this->assembly(false, true);
01313       this->matrix->close();
01314       this->assemble_qoi_derivative(qoi_indices,
01315                                     /* include_liftfunc = */ true,
01316                                     /* apply_constraints = */ false);
01317 
01318       for (unsigned int l=0; l != Np; ++l)
01319         {
01320           this->matrix->vector_mult(*tempvec, this->get_sensitivity_solution(l));
01321           for (unsigned int i=0; i != Nq; ++i)
01322             if (qoi_indices.has_index(i))
01323               {
01324                 this->get_adjoint_rhs(i).close();
01325                 Number current_terms =
01326                   (-this->get_adjoint_rhs(i).dot(this->get_sensitivity_solution(l)) +
01327                    tempvec->dot(this->get_adjoint_solution(i))) / (2.*delta_p);
01328                 sensitivities.second_derivative(i,k,l) += current_terms;
01329 
01330                 // We use the _uk terms twice; symmetry lets us reuse
01331                 // these calculations for the _ul terms.
01332 
01333                 sensitivities.second_derivative(i,l,k) += current_terms;
01334               }
01335         }
01336 
01337       // Don't leave the parameter perturbed
01338       *parameters[k] = old_parameterk;
01339 
01340       // Our last remaining terms are -R_uu(u,z)*u_k*u_l and
01341       // Q_uu(u)*u_k*u_l
01342       //
01343       // We take directional central finite differences of R_u and Q_u
01344       // to approximate these terms, e.g.:
01345       //
01346       // Q_uu(u)*u_k ~= (Q_u(u+dp*u_k) - Q_u(u-dp*u_k))/(2*dp)
01347 
01348       *this->solution = this->get_sensitivity_solution(k);
01349       *this->solution *= delta_p;
01350       *this->solution += *oldsolution;
01351       this->assembly(false, true);
01352       this->matrix->close();
01353       this->assemble_qoi_derivative(qoi_indices,
01354                                     /* include_liftfunc = */ true,
01355                                     /* apply_constraints = */ false);
01356 
01357       // The Hessian is symmetric, so we just calculate the lower
01358       // triangle and the diagonal, and we get the upper triangle from
01359       // the transpose of the lower
01360       //
01361       // Note that, because we took the directional finite difference
01362       // with respect to k and not l, we've added an O(delta_p^2)
01363       // error to any permutational symmetry in the Hessian...
01364       for (unsigned int l=0; l != k+1; ++l)
01365         {
01366           this->matrix->vector_mult(*tempvec, this->get_sensitivity_solution(l));
01367           for (unsigned int i=0; i != Nq; ++i)
01368             if (qoi_indices.has_index(i))
01369               {
01370                 this->get_adjoint_rhs(i).close();
01371                 Number current_terms =
01372                   (this->get_adjoint_rhs(i).dot(this->get_sensitivity_solution(l)) -
01373                    tempvec->dot(this->get_adjoint_solution(i))) / (2.*delta_p);
01374                 sensitivities.second_derivative(i,k,l) += current_terms;
01375                 if (k != l)
01376                   sensitivities.second_derivative(i,l,k) += current_terms;
01377               }
01378         }
01379 
01380       *this->solution = this->get_sensitivity_solution(k);
01381       *this->solution *= -delta_p;
01382       *this->solution += *oldsolution;
01383       this->assembly(false, true);
01384       this->matrix->close();
01385       this->assemble_qoi_derivative(qoi_indices,
01386                                     /* include_liftfunc = */ true,
01387                                     /* apply_constraints = */ false);
01388 
01389       for (unsigned int l=0; l != k+1; ++l)
01390         {
01391           this->matrix->vector_mult(*tempvec, this->get_sensitivity_solution(l));
01392           for (unsigned int i=0; i != Nq; ++i)
01393             if (qoi_indices.has_index(i))
01394               {
01395                 this->get_adjoint_rhs(i).close();
01396                 Number current_terms =
01397                   (-this->get_adjoint_rhs(i).dot(this->get_sensitivity_solution(l)) +
01398                    tempvec->dot(this->get_adjoint_solution(i))) / (2.*delta_p);
01399                 sensitivities.second_derivative(i,k,l) += current_terms;
01400                 if (k != l)
01401                   sensitivities.second_derivative(i,l,k) += current_terms;
01402               }
01403         }
01404 
01405       // Don't leave the solution perturbed
01406       *this->solution = *oldsolution;
01407     }
01408 
01409   // All parameters have been reset.
01410   // Don't leave the qoi or system changed - principle of least
01411   // surprise.
01412   this->assembly(true, true);
01413   this->rhs->close();
01414   this->matrix->close();
01415   this->assemble_qoi(qoi_indices);
01416 }
01417 
01418 
01419 
01420 LinearSolver<Number>* ImplicitSystem::get_linear_solver() const
01421 {
01422   LinearSolver<Number>* new_solver =
01423     LinearSolver<Number>::build(this->comm()).release();
01424 
01425   if (libMesh::on_command_line("--solver_system_names"))
01426     new_solver->init((this->name()+"_").c_str());
01427   else
01428     new_solver->init();
01429 
01430   return new_solver;
01431 }
01432 
01433 
01434 
01435 std::pair<unsigned int, Real> ImplicitSystem::get_linear_solve_parameters() const
01436 {
01437   return std::make_pair(this->get_equation_systems().parameters.get<unsigned int>("linear solver maximum iterations"),
01438                         this->get_equation_systems().parameters.get<Real>("linear solver tolerance"));
01439 }
01440 
01441 
01442 
01443 void ImplicitSystem::release_linear_solver(LinearSolver<Number>* s) const
01444 {
01445   delete s;
01446 }
01447 
01448 } // namespace libMesh