implicit_system.C
Go to the documentation of this file.00001 // The libMesh Finite Element Library. 00002 // Copyright (C) 2002-2012 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 { 00212 libMesh::err << "ERROR: Too late. Cannot add matrices to the system after initialization" 00213 << std::endl 00214 << " any more. You should have done this earlier." 00215 << std::endl; 00216 libmesh_error(); 00217 } 00218 00219 // Return the matrix if it is already there. 00220 if (this->have_matrix(mat_name)) 00221 return *(_matrices[mat_name]); 00222 00223 // Otherwise build the matrix and return it. 00224 SparseMatrix<Number>* buf = SparseMatrix<Number>::build().release(); 00225 _matrices.insert (std::make_pair (mat_name, buf)); 00226 00227 return *buf; 00228 } 00229 00230 00231 00232 const SparseMatrix<Number> * ImplicitSystem::request_matrix (const std::string& mat_name) const 00233 { 00234 // Make sure the matrix exists 00235 const_matrices_iterator pos = _matrices.find (mat_name); 00236 00237 if (pos == _matrices.end()) 00238 return NULL; 00239 00240 return pos->second; 00241 } 00242 00243 00244 00245 SparseMatrix<Number> * ImplicitSystem::request_matrix (const std::string& mat_name) 00246 { 00247 // Make sure the matrix exists 00248 matrices_iterator pos = _matrices.find (mat_name); 00249 00250 if (pos == _matrices.end()) 00251 return NULL; 00252 00253 return pos->second; 00254 } 00255 00256 00257 00258 const SparseMatrix<Number> & ImplicitSystem::get_matrix (const std::string& mat_name) const 00259 { 00260 // Make sure the matrix exists 00261 const_matrices_iterator pos = _matrices.find (mat_name); 00262 00263 if (pos == _matrices.end()) 00264 { 00265 libMesh::err << "ERROR: matrix " 00266 << mat_name 00267 << " does not exist in this system!" 00268 << std::endl; 00269 libmesh_error(); 00270 } 00271 00272 return *(pos->second); 00273 } 00274 00275 00276 00277 SparseMatrix<Number> & ImplicitSystem::get_matrix (const std::string& mat_name) 00278 { 00279 // Make sure the matrix exists 00280 matrices_iterator pos = _matrices.find (mat_name); 00281 00282 if (pos == _matrices.end()) 00283 { 00284 libMesh::err << "ERROR: matrix " 00285 << mat_name 00286 << " does not exist in this system!" 00287 << std::endl; 00288 libmesh_error(); 00289 } 00290 00291 return *(pos->second); 00292 } 00293 00294 00295 00296 void ImplicitSystem::add_system_matrix () 00297 { 00298 // Possible that we cleared the _matrices but 00299 // forgot to NULL-out the matrix? 00300 if (_matrices.empty()) matrix = NULL; 00301 00302 00303 // Only need to add the matrix if it isn't there 00304 // already! 00305 if (matrix == NULL) 00306 matrix = &(this->add_matrix ("System Matrix")); 00307 00308 libmesh_assert(matrix); 00309 } 00310 00311 00312 00313 void ImplicitSystem::disable_cache () { 00314 this->assemble_before_solve = true; 00315 this->get_linear_solver()->reuse_preconditioner(false); 00316 } 00317 00318 00319 00320 std::pair<unsigned int, Real> 00321 ImplicitSystem::sensitivity_solve (const ParameterVector& parameters) 00322 { 00323 // Log how long the linear solve takes. 00324 START_LOG("sensitivity_solve()", "ImplicitSystem"); 00325 00326 // The forward system should now already be solved. 00327 // Now assemble the corresponding sensitivity system. 00328 00329 if (this->assemble_before_solve) 00330 { 00331 // Build the Jacobian 00332 this->assembly(false, true); 00333 this->matrix->close(); 00334 00335 // Reset and build the RHS from the residual derivatives 00336 this->assemble_residual_derivatives(parameters); 00337 } 00338 00339 // The sensitivity problem is linear 00340 LinearSolver<Number> *linear_solver = this->get_linear_solver(); 00341 00342 // Our iteration counts and residuals will be sums of the individual 00343 // results 00344 std::pair<unsigned int, Real> solver_params = 00345 this->get_linear_solve_parameters(); 00346 std::pair<unsigned int, Real> totalrval = std::make_pair(0,0.0); 00347 00348 // Solve the linear system. 00349 SparseMatrix<Number> *pc = this->request_matrix("Preconditioner"); 00350 for (unsigned int p=0; p != parameters.size(); ++p) 00351 { 00352 std::pair<unsigned int, Real> rval = 00353 linear_solver->solve (*matrix, pc, 00354 this->add_sensitivity_solution(p), 00355 this->get_sensitivity_rhs(p), 00356 solver_params.second, 00357 solver_params.first); 00358 00359 totalrval.first += rval.first; 00360 totalrval.second += rval.second; 00361 } 00362 00363 // The linear solver may not have fit our constraints exactly 00364 #ifdef LIBMESH_ENABLE_CONSTRAINTS 00365 for (unsigned int p=0; p != parameters.size(); ++p) 00366 this->get_dof_map().enforce_constraints_exactly 00367 (*this, &this->get_sensitivity_solution(p), 00368 /* homogeneous = */ true); 00369 #endif 00370 00371 this->release_linear_solver(linear_solver); 00372 00373 // Stop logging the nonlinear solve 00374 STOP_LOG("sensitivity_solve()", "ImplicitSystem"); 00375 00376 return totalrval; 00377 } 00378 00379 00380 00381 std::pair<unsigned int, Real> 00382 ImplicitSystem::adjoint_solve (const QoISet& qoi_indices) 00383 { 00384 // Log how long the linear solve takes. 00385 START_LOG("adjoint_solve()", "ImplicitSystem"); 00386 00387 if (this->assemble_before_solve) 00388 // Assemble the linear system 00389 this->assembly (/* get_residual = */ false, 00390 /* get_jacobian = */ true); 00391 00392 // The adjoint problem is linear 00393 LinearSolver<Number> *linear_solver = this->get_linear_solver(); 00394 00395 // Reset and build the RHS from the QOI derivative 00396 this->assemble_qoi_derivative(qoi_indices); 00397 00398 // Our iteration counts and residuals will be sums of the individual 00399 // results 00400 std::pair<unsigned int, Real> solver_params = 00401 this->get_linear_solve_parameters(); 00402 std::pair<unsigned int, Real> totalrval = std::make_pair(0,0.0); 00403 00404 for (unsigned int i=0; i != this->qoi.size(); ++i) 00405 if (qoi_indices.has_index(i)) 00406 { 00407 const std::pair<unsigned int, Real> rval = 00408 linear_solver->adjoint_solve (*matrix, this->add_adjoint_solution(i), 00409 this->get_adjoint_rhs(i), 00410 solver_params.second, 00411 solver_params.first); 00412 00413 totalrval.first += rval.first; 00414 totalrval.second += rval.second; 00415 } 00416 00417 this->release_linear_solver(linear_solver); 00418 00419 // The linear solver may not have fit our constraints exactly 00420 #ifdef LIBMESH_ENABLE_CONSTRAINTS 00421 for (unsigned int i=0; i != this->qoi.size(); ++i) 00422 if (qoi_indices.has_index(i)) 00423 this->get_dof_map().enforce_constraints_exactly 00424 (*this, &this->get_adjoint_solution(i), 00425 /* homogeneous = */ true); 00426 #endif 00427 00428 // Stop logging the nonlinear solve 00429 STOP_LOG("adjoint_solve()", "ImplicitSystem"); 00430 00431 return totalrval; 00432 } 00433 00434 00435 00436 std::pair<unsigned int, Real> 00437 ImplicitSystem::weighted_sensitivity_adjoint_solve (const ParameterVector& parameters, 00438 const ParameterVector& weights, 00439 const QoISet& qoi_indices) 00440 { 00441 // Log how long the linear solve takes. 00442 START_LOG("weighted_sensitivity_adjoint_solve()", "ImplicitSystem"); 00443 00444 // We currently get partial derivatives via central differencing 00445 const Real delta_p = TOLERANCE; 00446 00447 // The forward system should now already be solved. 00448 // The adjoint system should now already be solved. 00449 // Now we're assembling a weighted sum of adjoint-adjoint systems: 00450 // 00451 // dR/du (u, sum_l(w_l*z^l)) = sum_l(w_l*(Q''_ul - R''_ul (u, z))) 00452 00453 // We'll assemble the rhs first, because the R'' term will require 00454 // perturbing the jacobian 00455 00456 // We'll use temporary rhs vectors, because we haven't (yet) found 00457 // any good reasons why users might want to save these: 00458 00459 std::vector<NumericVector<Number> *> temprhs(this->qoi.size()); 00460 for (unsigned int i=0; i != this->qoi.size(); ++i) 00461 if (qoi_indices.has_index(i)) 00462 temprhs[i] = this->rhs->zero_clone().release(); 00463 00464 // We approximate the _l partial derivatives via a central 00465 // differencing perturbation in the w_l direction: 00466 // 00467 // sum_l(w_l*v_l) ~= (v(p + dp*w_l*e_l) - v(p - dp*w_l*e_l))/(2*dp) 00468 00469 // PETSc doesn't implement SGEMX, so neither does NumericVector, 00470 // so we want to avoid calculating f -= R'*z. We'll thus evaluate 00471 // the above equation by first adding -v(p+dp...), then multiplying 00472 // the intermediate result vectors by -1, then adding -v(p-dp...), 00473 // then finally dividing by 2*dp. 00474 00475 ParameterVector oldparameters, parameterperturbation; 00476 parameters.deep_copy(oldparameters); 00477 weights.deep_copy(parameterperturbation); 00478 parameterperturbation *= delta_p; 00479 parameters += parameterperturbation; 00480 00481 this->assembly(false, true); 00482 this->matrix->close(); 00483 00484 // Take the discrete adjoint, so that we can calculate R_u(u,z) with 00485 // a matrix-vector product of R_u and z. 00486 matrix->get_transpose(*matrix); 00487 00488 this->assemble_qoi_derivative(qoi_indices); 00489 for (unsigned int i=0; i != this->qoi.size(); ++i) 00490 if (qoi_indices.has_index(i)) 00491 { 00492 this->get_adjoint_rhs(i).close(); 00493 *(temprhs[i]) -= this->get_adjoint_rhs(i); 00494 this->matrix->vector_mult_add(*(temprhs[i]), this->get_adjoint_solution(i)); 00495 *(temprhs[i]) *= -1.0; 00496 } 00497 00498 oldparameters.value_copy(parameters); 00499 parameterperturbation *= -1.0; 00500 parameters += parameterperturbation; 00501 00502 this->assembly(false, true); 00503 this->matrix->close(); 00504 matrix->get_transpose(*matrix); 00505 00506 this->assemble_qoi_derivative(qoi_indices); 00507 for (unsigned int i=0; i != this->qoi.size(); ++i) 00508 if (qoi_indices.has_index(i)) 00509 { 00510 this->get_adjoint_rhs(i).close(); 00511 *(temprhs[i]) -= this->get_adjoint_rhs(i); 00512 this->matrix->vector_mult_add(*(temprhs[i]), this->get_adjoint_solution(i)); 00513 *(temprhs[i]) /= (2.0*delta_p); 00514 } 00515 00516 // Finally, assemble the jacobian at the non-perturbed parameter 00517 // values. Ignore assemble_before_solve; if we had a good 00518 // non-perturbed matrix before we've already overwritten it. 00519 oldparameters.value_copy(parameters); 00520 00521 // if (this->assemble_before_solve) 00522 { 00523 // Build the Jacobian 00524 this->assembly(false, true); 00525 this->matrix->close(); 00526 00527 // Take the discrete adjoint 00528 matrix->get_transpose(*matrix); 00529 } 00530 00531 // The weighted adjoint-adjoint problem is linear 00532 LinearSolver<Number> *linear_solver = this->get_linear_solver(); 00533 00534 // Our iteration counts and residuals will be sums of the individual 00535 // results 00536 std::pair<unsigned int, Real> solver_params = 00537 this->get_linear_solve_parameters(); 00538 std::pair<unsigned int, Real> totalrval = std::make_pair(0,0.0); 00539 00540 for (unsigned int i=0; i != this->qoi.size(); ++i) 00541 if (qoi_indices.has_index(i)) 00542 { 00543 const std::pair<unsigned int, Real> rval = 00544 linear_solver->solve (*matrix, this->add_weighted_sensitivity_adjoint_solution(i), 00545 *(temprhs[i]), 00546 solver_params.second, 00547 solver_params.first); 00548 00549 totalrval.first += rval.first; 00550 totalrval.second += rval.second; 00551 } 00552 00553 this->release_linear_solver(linear_solver); 00554 00555 for (unsigned int i=0; i != this->qoi.size(); ++i) 00556 if (qoi_indices.has_index(i)) 00557 delete temprhs[i]; 00558 00559 // The linear solver may not have fit our constraints exactly 00560 #ifdef LIBMESH_ENABLE_CONSTRAINTS 00561 for (unsigned int i=0; i != this->qoi.size(); ++i) 00562 if (qoi_indices.has_index(i)) 00563 this->get_dof_map().enforce_constraints_exactly 00564 (*this, &this->get_weighted_sensitivity_adjoint_solution(i), 00565 /* homogeneous = */ true); 00566 #endif 00567 00568 // Stop logging the nonlinear solve 00569 STOP_LOG("weighted_sensitivity_adjoint_solve()", "ImplicitSystem"); 00570 00571 return totalrval; 00572 } 00573 00574 00575 00576 std::pair<unsigned int, Real> 00577 ImplicitSystem::weighted_sensitivity_solve (const ParameterVector& parameters, 00578 const ParameterVector& weights) 00579 { 00580 // Log how long the linear solve takes. 00581 START_LOG("weighted_sensitivity_solve()", "ImplicitSystem"); 00582 00583 // We currently get partial derivatives via central differencing 00584 const Real delta_p = TOLERANCE; 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); 00609 this->rhs->close(); 00610 00611 AutoPtr<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); 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) 00661 { 00662 const unsigned int Np = libmesh_cast_int<unsigned int> 00663 (parameters.size()); 00664 Real deltap = TOLERANCE; 00665 00666 for (unsigned int p=0; p != Np; ++p) 00667 { 00668 NumericVector<Number> &sensitivity_rhs = this->add_sensitivity_rhs(p); 00669 00670 // Approximate -(partial R / partial p) by 00671 // (R(p-dp) - R(p+dp)) / (2*dp) 00672 00673 Number old_parameter = *parameters[p]; 00674 *parameters[p] -= deltap; 00675 00676 this->assembly(true, false); 00677 this->rhs->close(); 00678 sensitivity_rhs = *this->rhs; 00679 00680 *parameters[p] = old_parameter + deltap; 00681 00682 this->assembly(true, false); 00683 this->rhs->close(); 00684 00685 sensitivity_rhs -= *this->rhs; 00686 sensitivity_rhs /= (2*deltap); 00687 sensitivity_rhs.close(); 00688 00689 *parameters[p] = old_parameter; 00690 } 00691 } 00692 00693 00694 00695 void ImplicitSystem::adjoint_qoi_parameter_sensitivity 00696 (const QoISet& qoi_indices, 00697 const ParameterVector& parameters, 00698 SensitivityData& sensitivities) 00699 { 00700 const unsigned int Np = libmesh_cast_int<unsigned int> 00701 (parameters.size()); 00702 const unsigned int Nq = libmesh_cast_int<unsigned int> 00703 (qoi.size()); 00704 00705 // We currently get partial derivatives via central differencing 00706 const Real delta_p = TOLERANCE; 00707 00708 // An introduction to the problem: 00709 // 00710 // Residual R(u(p),p) = 0 00711 // partial R / partial u = J = system matrix 00712 // 00713 // This implies that: 00714 // d/dp(R) = 0 00715 // (partial R / partial p) + 00716 // (partial R / partial u) * (partial u / partial p) = 0 00717 00718 // We first do an adjoint solve: 00719 // J^T * z = (partial q / partial u) 00720 // if we havent already or dont have an initial condition for the adjoint 00721 if (!this->is_adjoint_already_solved()) 00722 { 00723 this->adjoint_solve(qoi_indices); 00724 } 00725 00726 // Get ready to fill in senstivities: 00727 sensitivities.allocate_data(qoi_indices, *this, parameters); 00728 00729 // We use the identities: 00730 // dq/dp = (partial q / partial p) + (partial q / partial u) * 00731 // (partial u / partial p) 00732 // dq/dp = (partial q / partial p) + (J^T * z) * 00733 // (partial u / partial p) 00734 // dq/dp = (partial q / partial p) + z * J * 00735 // (partial u / partial p) 00736 00737 // Leading to our final formula: 00738 // dq/dp = (partial q / partial p) - z * (partial R / partial p) 00739 00740 for (unsigned int j=0; j != Np; ++j) 00741 { 00742 // (partial q / partial p) ~= (q(p+dp)-q(p-dp))/(2*dp) 00743 // (partial R / partial p) ~= (rhs(p+dp) - rhs(p-dp))/(2*dp) 00744 00745 Number old_parameter = *parameters[j]; 00746 // Number old_qoi = this->qoi; 00747 00748 *parameters[j] = old_parameter - delta_p; 00749 this->assemble_qoi(qoi_indices); 00750 std::vector<Number> qoi_minus = this->qoi; 00751 00752 this->assembly(true, false); 00753 this->rhs->close(); 00754 00755 // FIXME - this can and should be optimized to avoid the clone() 00756 AutoPtr<NumericVector<Number> > partialR_partialp = this->rhs->clone(); 00757 *partialR_partialp *= -1; 00758 00759 *parameters[j] = old_parameter + delta_p; 00760 this->assemble_qoi(qoi_indices); 00761 std::vector<Number>& qoi_plus = this->qoi; 00762 00763 std::vector<Number> partialq_partialp(Nq, 0); 00764 for (unsigned int i=0; i != Nq; ++i) 00765 if (qoi_indices.has_index(i)) 00766 partialq_partialp[i] = (qoi_plus[i] - qoi_minus[i]) / (2.*delta_p); 00767 00768 this->assembly(true, false); 00769 this->rhs->close(); 00770 *partialR_partialp += *this->rhs; 00771 *partialR_partialp /= (2.*delta_p); 00772 00773 // Don't leave the parameter changed 00774 *parameters[j] = old_parameter; 00775 00776 for (unsigned int i=0; i != Nq; ++i) 00777 if (qoi_indices.has_index(i)) 00778 sensitivities[i][j] = partialq_partialp[i] - 00779 partialR_partialp->dot(this->get_adjoint_solution(i)); 00780 } 00781 00782 // All parameters have been reset. 00783 // We didn't cache the original rhs or matrix for memory reasons, 00784 // but we can restore them to a state consistent solution - 00785 // principle of least surprise. 00786 this->assembly(true, true); 00787 this->rhs->close(); 00788 this->matrix->close(); 00789 this->assemble_qoi(qoi_indices); 00790 } 00791 00792 00793 00794 void ImplicitSystem::forward_qoi_parameter_sensitivity 00795 (const QoISet& qoi_indices, 00796 const ParameterVector& parameters, 00797 SensitivityData& sensitivities) 00798 { 00799 const unsigned int Np = libmesh_cast_int<unsigned int> 00800 (parameters.size()); 00801 const unsigned int Nq = libmesh_cast_int<unsigned int> 00802 (qoi.size()); 00803 00804 // We currently get partial derivatives via central differencing 00805 const Real delta_p = TOLERANCE; 00806 00807 // An introduction to the problem: 00808 // 00809 // Residual R(u(p),p) = 0 00810 // partial R / partial u = J = system matrix 00811 // 00812 // This implies that: 00813 // d/dp(R) = 0 00814 // (partial R / partial p) + 00815 // (partial R / partial u) * (partial u / partial p) = 0 00816 00817 // We first solve for (partial u / partial p) for each parameter: 00818 // J * (partial u / partial p) = - (partial R / partial p) 00819 00820 this->sensitivity_solve(parameters); 00821 00822 // Get ready to fill in senstivities: 00823 sensitivities.allocate_data(qoi_indices, *this, parameters); 00824 00825 // We use the identity: 00826 // dq/dp = (partial q / partial p) + (partial q / partial u) * 00827 // (partial u / partial p) 00828 00829 // We get (partial q / partial u) from the user 00830 this->assemble_qoi_derivative(qoi_indices); 00831 00832 // We don't need these to be closed() in this function, but libMesh 00833 // standard practice is to have them closed() by the time the 00834 // function exits 00835 for (unsigned int i=0; i != this->qoi.size(); ++i) 00836 if (qoi_indices.has_index(i)) 00837 this->get_adjoint_rhs(i).close(); 00838 00839 for (unsigned int j=0; j != Np; ++j) 00840 { 00841 // (partial q / partial p) ~= (q(p+dp)-q(p-dp))/(2*dp) 00842 00843 Number old_parameter = *parameters[j]; 00844 00845 *parameters[j] = old_parameter - delta_p; 00846 this->assemble_qoi(); 00847 std::vector<Number> qoi_minus = this->qoi; 00848 00849 *parameters[j] = old_parameter + delta_p; 00850 this->assemble_qoi(); 00851 std::vector<Number>& qoi_plus = this->qoi; 00852 00853 std::vector<Number> partialq_partialp(Nq, 0); 00854 for (unsigned int i=0; i != Nq; ++i) 00855 if (qoi_indices.has_index(i)) 00856 partialq_partialp[i] = (qoi_plus[i] - qoi_minus[i]) / (2.*delta_p); 00857 00858 // Don't leave the parameter changed 00859 *parameters[j] = old_parameter; 00860 00861 for (unsigned int i=0; i != Nq; ++i) 00862 if (qoi_indices.has_index(i)) 00863 sensitivities[i][j] = partialq_partialp[i] + 00864 this->get_adjoint_rhs(i).dot(this->get_sensitivity_solution(j)); 00865 } 00866 00867 // All parameters have been reset. 00868 // We didn't cache the original rhs or matrix for memory reasons, 00869 // but we can restore them to a state consistent solution - 00870 // principle of least surprise. 00871 this->assembly(true, true); 00872 this->rhs->close(); 00873 this->matrix->close(); 00874 this->assemble_qoi(qoi_indices); 00875 } 00876 00877 00878 00879 void ImplicitSystem::qoi_parameter_hessian_vector_product 00880 (const QoISet& qoi_indices, 00881 const ParameterVector& parameters, 00882 const ParameterVector& vector, 00883 SensitivityData& sensitivities) 00884 { 00885 // We currently get partial derivatives via finite differencing 00886 const Real delta_p = TOLERANCE; 00887 00888 // We'll use a single temporary vector for matrix-vector-vector products 00889 AutoPtr<NumericVector<Number> > tempvec = this->solution->zero_clone(); 00890 00891 const unsigned int Np = libmesh_cast_int<unsigned int> 00892 (parameters.size()); 00893 const unsigned int Nq = libmesh_cast_int<unsigned int> 00894 (qoi.size()); 00895 00896 // For each quantity of interest q, the parameter sensitivity 00897 // Hessian is defined as q''_{kl} = {d^2 q}/{d p_k d p_l}. 00898 // Given a vector of parameter perturbation weights w_l, this 00899 // function evaluates the hessian-vector product sum_l(q''_{kl}*w_l) 00900 // 00901 // We calculate it from values and partial derivatives of the 00902 // quantity of interest function Q, solution u, adjoint solution z, 00903 // parameter sensitivity adjoint solutions z^l, and residual R, as: 00904 // 00905 // sum_l(q''_{kl}*w_l) = 00906 // sum_l(w_l * Q''_{kl}) + Q''_{uk}(u)*(sum_l(w_l u'_l)) - 00907 // R'_k(u, sum_l(w_l*z^l)) - R'_{uk}(u,z)*(sum_l(w_l u'_l) - 00908 // sum_l(w_l*R''_{kl}(u,z)) 00909 // 00910 // See the adjoints model document for more details. 00911 00912 // We first do an adjoint solve to get z for each quantity of 00913 // interest 00914 // if we havent already or dont have an initial condition for the adjoint 00915 if (!this->is_adjoint_already_solved()) 00916 { 00917 this->adjoint_solve(qoi_indices); 00918 } 00919 00920 // Get ready to fill in senstivities: 00921 sensitivities.allocate_data(qoi_indices, *this, parameters); 00922 00923 // We can't solve for all the solution sensitivities u'_l or for all 00924 // of the parameter sensitivity adjoint solutions z^l without 00925 // requiring O(Nq*Np) linear solves. So we'll solve directly for their 00926 // weighted sum - this is just O(Nq) solves. 00927 00928 // First solve for sum_l(w_l u'_l). 00929 this->weighted_sensitivity_solve(parameters, vector); 00930 00931 // Then solve for sum_l(w_l z^l). 00932 this->weighted_sensitivity_adjoint_solve(parameters, vector, qoi_indices); 00933 00934 for (unsigned int k=0; k != Np; ++k) 00935 { 00936 // We approximate sum_l(w_l * Q''_{kl}) with a central 00937 // differencing pertubation: 00938 // sum_l(w_l * Q''_{kl}) ~= 00939 // (Q(p + dp*w_l*e_l + dp*e_k) - Q(p - dp*w_l*e_l + dp*e_k) - 00940 // Q(p + dp*w_l*e_l - dp*e_k) + Q(p - dp*w_l*e_l - dp*e_k))/(4*dp^2) 00941 00942 // The sum(w_l*R''_kl) term requires the same sort of pertubation, 00943 // and so we subtract it in at the same time: 00944 // sum_l(w_l * R''_{kl}) ~= 00945 // (R(p + dp*w_l*e_l + dp*e_k) - R(p - dp*w_l*e_l + dp*e_k) - 00946 // R(p + dp*w_l*e_l - dp*e_k) + R(p - dp*w_l*e_l - dp*e_k))/(4*dp^2) 00947 00948 ParameterVector oldparameters, parameterperturbation; 00949 parameters.deep_copy(oldparameters); 00950 vector.deep_copy(parameterperturbation); 00951 parameterperturbation *= delta_p; 00952 parameters += parameterperturbation; 00953 00954 Number old_parameter = *parameters[k]; 00955 00956 *parameters[k] = old_parameter + delta_p; 00957 this->assemble_qoi(qoi_indices); 00958 this->assembly(true, false); 00959 this->rhs->close(); 00960 std::vector<Number> partial2q_term = this->qoi; 00961 std::vector<Number> partial2R_term(this->qoi.size()); 00962 for (unsigned int i=0; i != Nq; ++i) 00963 if (qoi_indices.has_index(i)) 00964 partial2R_term[i] = this->rhs->dot(this->get_adjoint_solution(i)); 00965 00966 *parameters[k] = old_parameter - delta_p; 00967 this->assemble_qoi(qoi_indices); 00968 this->assembly(true, false); 00969 this->rhs->close(); 00970 for (unsigned int i=0; i != Nq; ++i) 00971 if (qoi_indices.has_index(i)) 00972 { 00973 partial2q_term[i] -= this->qoi[i]; 00974 partial2R_term[i] -= this->rhs->dot(this->get_adjoint_solution(i)); 00975 } 00976 00977 oldparameters.value_copy(parameters); 00978 parameterperturbation *= -1.0; 00979 parameters += parameterperturbation; 00980 00981 // Re-center old_parameter, which may be affected by vector 00982 old_parameter = *parameters[k]; 00983 00984 *parameters[k] = old_parameter + delta_p; 00985 this->assemble_qoi(qoi_indices); 00986 this->assembly(true, false); 00987 this->rhs->close(); 00988 for (unsigned int i=0; i != Nq; ++i) 00989 if (qoi_indices.has_index(i)) 00990 { 00991 partial2q_term[i] -= this->qoi[i]; 00992 partial2R_term[i] -= this->rhs->dot(this->get_adjoint_solution(i)); 00993 } 00994 00995 *parameters[k] = old_parameter - delta_p; 00996 this->assemble_qoi(qoi_indices); 00997 this->assembly(true, false); 00998 this->rhs->close(); 00999 for (unsigned int i=0; i != Nq; ++i) 01000 if (qoi_indices.has_index(i)) 01001 { 01002 partial2q_term[i] += this->qoi[i]; 01003 partial2R_term[i] += this->rhs->dot(this->get_adjoint_solution(i)); 01004 } 01005 01006 for (unsigned int i=0; i != Nq; ++i) 01007 if (qoi_indices.has_index(i)) 01008 { 01009 partial2q_term[i] /= (4. * delta_p * delta_p); 01010 partial2R_term[i] /= (4. * delta_p * delta_p); 01011 } 01012 01013 for (unsigned int i=0; i != Nq; ++i) 01014 if (qoi_indices.has_index(i)) 01015 sensitivities[i][k] = partial2q_term[i] - partial2R_term[i]; 01016 01017 // We get (partial q / partial u), R, and 01018 // (partial R / partial u) from the user, but centrally 01019 // difference to get q_uk, R_k, and R_uk terms: 01020 // (partial R / partial k) 01021 // 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) 01022 // (partial^2 q / partial u partial k) 01023 // q_uk = (q_u(p+dp*e_k) - q_u(p-dp*e_k))/(2*dp) 01024 // (partial^2 R / partial u partial k) 01025 // 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) 01026 01027 // To avoid creating Nq temporary vectors for q_uk or R_uk, we add 01028 // subterms to the sensitivities output one by one. 01029 // 01030 // FIXME: this is probably a bad order of operations for 01031 // controlling floating point error. 01032 01033 *parameters[k] = old_parameter + delta_p; 01034 this->assembly(true, true); 01035 this->rhs->close(); 01036 this->matrix->close(); 01037 this->assemble_qoi_derivative(qoi_indices); 01038 01039 this->matrix->vector_mult(*tempvec, this->get_weighted_sensitivity_solution()); 01040 01041 for (unsigned int i=0; i != Nq; ++i) 01042 if (qoi_indices.has_index(i)) 01043 { 01044 this->get_adjoint_rhs(i).close(); 01045 sensitivities[i][k] += (this->get_adjoint_rhs(i).dot(this->get_weighted_sensitivity_solution()) - 01046 this->rhs->dot(this->get_weighted_sensitivity_adjoint_solution(i)) - 01047 this->get_adjoint_solution(i).dot(*tempvec)) / (2.*delta_p); 01048 } 01049 01050 *parameters[k] = old_parameter - delta_p; 01051 this->assembly(true, true); 01052 this->rhs->close(); 01053 this->matrix->close(); 01054 this->assemble_qoi_derivative(qoi_indices); 01055 01056 this->matrix->vector_mult(*tempvec, this->get_weighted_sensitivity_solution()); 01057 01058 for (unsigned int i=0; i != Nq; ++i) 01059 if (qoi_indices.has_index(i)) 01060 { 01061 this->get_adjoint_rhs(i).close(); 01062 sensitivities[i][k] += (-this->get_adjoint_rhs(i).dot(this->get_weighted_sensitivity_solution()) + 01063 this->rhs->dot(this->get_weighted_sensitivity_adjoint_solution(i)) + 01064 this->get_adjoint_solution(i).dot(*tempvec)) / (2.*delta_p); 01065 } 01066 } 01067 01068 // All parameters have been reset. 01069 // Don't leave the qoi or system changed - principle of least 01070 // surprise. 01071 this->assembly(true, true); 01072 this->rhs->close(); 01073 this->matrix->close(); 01074 this->assemble_qoi(qoi_indices); 01075 } 01076 01077 01078 01079 void ImplicitSystem::qoi_parameter_hessian 01080 (const QoISet& qoi_indices, 01081 const ParameterVector& parameters, 01082 SensitivityData& sensitivities) 01083 { 01084 // We currently get partial derivatives via finite differencing 01085 const Real delta_p = TOLERANCE; 01086 01087 // We'll use one temporary vector for matrix-vector-vector products 01088 AutoPtr<NumericVector<Number> > tempvec = this->solution->zero_clone(); 01089 01090 // And another temporary vector to hold a copy of the true solution 01091 // so we can safely perturb this->solution. 01092 AutoPtr<NumericVector<Number> > oldsolution = this->solution->clone(); 01093 01094 const unsigned int Np = libmesh_cast_int<unsigned int> 01095 (parameters.size()); 01096 const unsigned int Nq = libmesh_cast_int<unsigned int> 01097 (qoi.size()); 01098 01099 // For each quantity of interest q, the parameter sensitivity 01100 // Hessian is defined as q''_{kl} = {d^2 q}/{d p_k d p_l}. 01101 // 01102 // We calculate it from values and partial derivatives of the 01103 // quantity of interest function Q, solution u, adjoint solution z, 01104 // and residual R, as: 01105 // 01106 // q''_{kl} = 01107 // Q''_{kl} + Q''_{uk}(u)*u'_l + Q''_{ul}(u) * u'_k + 01108 // Q''_{uu}(u)*u'_k*u'_l - 01109 // R''_{kl}(u,z) - 01110 // R''_{uk}(u,z)*u'_l - R''_{ul}(u,z)*u'_k - 01111 // R''_{uu}(u,z)*u'_k*u'_l 01112 // 01113 // See the adjoints model document for more details. 01114 01115 // We first do an adjoint solve to get z for each quantity of 01116 // interest 01117 // if we havent already or dont have an initial condition for the adjoint 01118 if (!this->is_adjoint_already_solved()) 01119 { 01120 this->adjoint_solve(qoi_indices); 01121 } 01122 01123 // And a sensitivity solve to get u_k for each parameter 01124 this->sensitivity_solve(parameters); 01125 01126 // Get ready to fill in second derivatives: 01127 sensitivities.allocate_hessian_data(qoi_indices, *this, parameters); 01128 01129 for (unsigned int k=0; k != Np; ++k) 01130 { 01131 Number old_parameterk = *parameters[k]; 01132 01133 // The Hessian is symmetric, so we just calculate the lower 01134 // triangle and the diagonal, and we get the upper triangle from 01135 // the transpose of the lower 01136 01137 for (unsigned int l=0; l != k+1; ++l) 01138 { 01139 // The second partial derivatives with respect to parameters 01140 // are all calculated via a central finite difference 01141 // stencil: 01142 // F''_{kl} ~= (F(p+dp*e_k+dp*e_l) - F(p+dp*e_k-dp*e_l) - 01143 // F(p-dp*e_k+dp*e_l) + F(p-dp*e_k-dp*e_l))/(4*dp^2) 01144 // We will add Q''_{kl}(u) and subtract R''_{kl}(u,z) at the 01145 // same time. 01146 // 01147 // We have to be careful with the perturbations to handle 01148 // the k=l case 01149 01150 Number old_parameterl = *parameters[l]; 01151 01152 *parameters[k] += delta_p; 01153 *parameters[l] += delta_p; 01154 this->assemble_qoi(qoi_indices); 01155 this->assembly(true, false); 01156 this->rhs->close(); 01157 std::vector<Number> partial2q_term = this->qoi; 01158 std::vector<Number> partial2R_term(this->qoi.size()); 01159 for (unsigned int i=0; i != Nq; ++i) 01160 if (qoi_indices.has_index(i)) 01161 partial2R_term[i] = this->rhs->dot(this->get_adjoint_solution(i)); 01162 01163 *parameters[l] -= 2.*delta_p; 01164 this->assemble_qoi(qoi_indices); 01165 this->assembly(true, false); 01166 this->rhs->close(); 01167 for (unsigned int i=0; i != Nq; ++i) 01168 if (qoi_indices.has_index(i)) 01169 { 01170 partial2q_term[i] -= this->qoi[i]; 01171 partial2R_term[i] -= this->rhs->dot(this->get_adjoint_solution(i)); 01172 } 01173 01174 *parameters[k] -= 2.*delta_p; 01175 this->assemble_qoi(qoi_indices); 01176 this->assembly(true, false); 01177 this->rhs->close(); 01178 for (unsigned int i=0; i != Nq; ++i) 01179 if (qoi_indices.has_index(i)) 01180 { 01181 partial2q_term[i] += this->qoi[i]; 01182 partial2R_term[i] += this->rhs->dot(this->get_adjoint_solution(i)); 01183 } 01184 01185 *parameters[l] += 2.*delta_p; 01186 this->assemble_qoi(qoi_indices); 01187 this->assembly(true, false); 01188 this->rhs->close(); 01189 for (unsigned int i=0; i != Nq; ++i) 01190 if (qoi_indices.has_index(i)) 01191 { 01192 partial2q_term[i] -= this->qoi[i]; 01193 partial2R_term[i] -= this->rhs->dot(this->get_adjoint_solution(i)); 01194 partial2q_term[i] /= (4. * delta_p * delta_p); 01195 partial2R_term[i] /= (4. * delta_p * delta_p); 01196 } 01197 01198 for (unsigned int i=0; i != Nq; ++i) 01199 if (qoi_indices.has_index(i)) 01200 { 01201 Number current_terms = partial2q_term[i] - partial2R_term[i]; 01202 sensitivities.second_derivative(i,k,l) += current_terms; 01203 if (k != l) 01204 sensitivities.second_derivative(i,l,k) += current_terms; 01205 } 01206 01207 // Don't leave the parameters perturbed 01208 *parameters[l] = old_parameterl; 01209 *parameters[k] = old_parameterk; 01210 } 01211 01212 // We get (partial q / partial u) and 01213 // (partial R / partial u) from the user, but centrally 01214 // difference to get q_uk and R_uk terms: 01215 // (partial^2 q / partial u partial k) 01216 // q_uk*u'_l = (q_u(p+dp*e_k)*u'_l - q_u(p-dp*e_k)*u'_l)/(2*dp) 01217 // 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) 01218 // 01219 // To avoid creating Nq temporary vectors, we add these 01220 // subterms to the sensitivities output one by one. 01221 // 01222 // FIXME: this is probably a bad order of operations for 01223 // controlling floating point error. 01224 01225 *parameters[k] = old_parameterk + delta_p; 01226 this->assembly(false, true); 01227 this->matrix->close(); 01228 this->assemble_qoi_derivative(qoi_indices); 01229 01230 for (unsigned int l=0; l != Np; ++l) 01231 { 01232 this->matrix->vector_mult(*tempvec, this->get_sensitivity_solution(l)); 01233 for (unsigned int i=0; i != Nq; ++i) 01234 if (qoi_indices.has_index(i)) 01235 { 01236 this->get_adjoint_rhs(i).close(); 01237 Number current_terms = 01238 (this->get_adjoint_rhs(i).dot(this->get_sensitivity_solution(l)) - 01239 tempvec->dot(this->get_adjoint_solution(i))) / (2.*delta_p); 01240 sensitivities.second_derivative(i,k,l) += current_terms; 01241 01242 // We use the _uk terms twice; symmetry lets us reuse 01243 // these calculations for the _ul terms. 01244 01245 sensitivities.second_derivative(i,l,k) += current_terms; 01246 } 01247 } 01248 01249 *parameters[k] = old_parameterk - delta_p; 01250 this->assembly(false, true); 01251 this->matrix->close(); 01252 this->assemble_qoi_derivative(qoi_indices); 01253 01254 for (unsigned int l=0; l != Np; ++l) 01255 { 01256 this->matrix->vector_mult(*tempvec, this->get_sensitivity_solution(l)); 01257 for (unsigned int i=0; i != Nq; ++i) 01258 if (qoi_indices.has_index(i)) 01259 { 01260 this->get_adjoint_rhs(i).close(); 01261 Number current_terms = 01262 (-this->get_adjoint_rhs(i).dot(this->get_sensitivity_solution(l)) + 01263 tempvec->dot(this->get_adjoint_solution(i))) / (2.*delta_p); 01264 sensitivities.second_derivative(i,k,l) += current_terms; 01265 01266 // We use the _uk terms twice; symmetry lets us reuse 01267 // these calculations for the _ul terms. 01268 01269 sensitivities.second_derivative(i,l,k) += current_terms; 01270 } 01271 } 01272 01273 // Don't leave the parameter perturbed 01274 *parameters[k] = old_parameterk; 01275 01276 // Our last remaining terms are -R_uu(u,z)*u_k*u_l and 01277 // Q_uu(u)*u_k*u_l 01278 // 01279 // We take directional central finite differences of R_u and Q_u 01280 // to approximate these terms, e.g.: 01281 // 01282 // Q_uu(u)*u_k ~= (Q_u(u+dp*u_k) - Q_u(u-dp*u_k))/(2*dp) 01283 01284 *this->solution = this->get_sensitivity_solution(k); 01285 *this->solution *= delta_p; 01286 *this->solution += *oldsolution; 01287 this->assembly(false, true); 01288 this->matrix->close(); 01289 this->assemble_qoi_derivative(qoi_indices); 01290 01291 // The Hessian is symmetric, so we just calculate the lower 01292 // triangle and the diagonal, and we get the upper triangle from 01293 // the transpose of the lower 01294 // 01295 // Note that, because we took the directional finite difference 01296 // with respect to k and not l, we've added an O(delta_p^2) 01297 // error to any permutational symmetry in the Hessian... 01298 for (unsigned int l=0; l != k+1; ++l) 01299 { 01300 this->matrix->vector_mult(*tempvec, this->get_sensitivity_solution(l)); 01301 for (unsigned int i=0; i != Nq; ++i) 01302 if (qoi_indices.has_index(i)) 01303 { 01304 this->get_adjoint_rhs(i).close(); 01305 Number current_terms = 01306 (this->get_adjoint_rhs(i).dot(this->get_sensitivity_solution(l)) - 01307 tempvec->dot(this->get_adjoint_solution(i))) / (2.*delta_p); 01308 sensitivities.second_derivative(i,k,l) += current_terms; 01309 if (k != l) 01310 sensitivities.second_derivative(i,l,k) += current_terms; 01311 } 01312 } 01313 01314 *this->solution = this->get_sensitivity_solution(k); 01315 *this->solution *= -delta_p; 01316 *this->solution += *oldsolution; 01317 this->assembly(false, true); 01318 this->matrix->close(); 01319 this->assemble_qoi_derivative(qoi_indices); 01320 01321 for (unsigned int l=0; l != k+1; ++l) 01322 { 01323 this->matrix->vector_mult(*tempvec, this->get_sensitivity_solution(l)); 01324 for (unsigned int i=0; i != Nq; ++i) 01325 if (qoi_indices.has_index(i)) 01326 { 01327 this->get_adjoint_rhs(i).close(); 01328 Number current_terms = 01329 (-this->get_adjoint_rhs(i).dot(this->get_sensitivity_solution(l)) + 01330 tempvec->dot(this->get_adjoint_solution(i))) / (2.*delta_p); 01331 sensitivities.second_derivative(i,k,l) += current_terms; 01332 if (k != l) 01333 sensitivities.second_derivative(i,l,k) += current_terms; 01334 } 01335 } 01336 01337 // Don't leave the solution perturbed 01338 *this->solution = *oldsolution; 01339 } 01340 01341 // All parameters have been reset. 01342 // Don't leave the qoi or system changed - principle of least 01343 // surprise. 01344 this->assembly(true, true); 01345 this->rhs->close(); 01346 this->matrix->close(); 01347 this->assemble_qoi(qoi_indices); 01348 } 01349 01350 01351 01352 LinearSolver<Number>* ImplicitSystem::get_linear_solver() const 01353 { 01354 return LinearSolver<Number>::build().release(); 01355 } 01356 01357 01358 01359 std::pair<unsigned int, Real> ImplicitSystem::get_linear_solve_parameters() const 01360 { 01361 return std::make_pair(this->get_equation_systems().parameters.get<unsigned int>("linear solver maximum iterations"), 01362 this->get_equation_systems().parameters.get<Real>("linear solver tolerance")); 01363 } 01364 01365 01366 01367 void ImplicitSystem::release_linear_solver(LinearSolver<Number>* s) const 01368 { 01369 delete s; 01370 } 01371 01372 } // namespace libMesh 01373
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:47 UTC
Hosted By: