$extrastylesheet
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 // LibMesh includes 00019 #include "libmesh/equation_systems.h" 00020 #include "libmesh/continuation_system.h" 00021 #include "libmesh/linear_solver.h" 00022 #include "libmesh/time_solver.h" 00023 #include "libmesh/newton_solver.h" 00024 #include "libmesh/sparse_matrix.h" 00025 00026 namespace libMesh 00027 { 00028 00029 ContinuationSystem::ContinuationSystem (EquationSystems& es, 00030 const std::string& name_in, 00031 const unsigned int number_in) 00032 : Parent(es, name_in, number_in), 00033 continuation_parameter(NULL), 00034 quiet(true), 00035 continuation_parameter_tolerance(1.e-6), 00036 solution_tolerance(1.e-6), 00037 initial_newton_tolerance(0.01), 00038 old_continuation_parameter(0.), 00039 min_continuation_parameter(0.), 00040 max_continuation_parameter(0.), 00041 Theta(1.), 00042 Theta_LOCA(1.), 00043 //tau(1.), 00044 n_backtrack_steps(5), 00045 n_arclength_reductions(5), 00046 ds_min(1.e-8), 00047 predictor(Euler), 00048 newton_stepgrowth_aggressiveness(1.), 00049 newton_progress_check(true), 00050 rhs_mode(Residual), 00051 linear_solver(LinearSolver<Number>::build(es.comm())), 00052 tangent_initialized(false), 00053 newton_solver(NULL), 00054 dlambda_ds(0.707), 00055 ds(0.1), 00056 ds_current(0.1), 00057 previous_dlambda_ds(0.), 00058 previous_ds(0.), 00059 newton_step(0) 00060 { 00061 // Warn about using untested code 00062 libmesh_experimental(); 00063 00064 if (libMesh::on_command_line("--solver_system_names")) 00065 linear_solver->init((this->name()+"_").c_str()); 00066 else 00067 linear_solver->init(); 00068 } 00069 00070 00071 00072 00073 ContinuationSystem::~ContinuationSystem () 00074 { 00075 this->clear(); 00076 } 00077 00078 00079 00080 00081 void ContinuationSystem::clear() 00082 { 00083 // FIXME: Do anything here, e.g. zero vectors, etc? 00084 00085 // Call the Parent's clear function 00086 Parent::clear(); 00087 } 00088 00089 00090 00091 void ContinuationSystem::init_data () 00092 { 00093 // Add a vector which stores the tangent "du/ds" to the system and save its pointer. 00094 du_ds = &(add_vector("du_ds")); 00095 00096 // Add a vector which stores the tangent "du/ds" to the system and save its pointer. 00097 previous_du_ds = &(add_vector("previous_du_ds")); 00098 00099 // Add a vector to keep track of the previous nonlinear solution 00100 // at the old value of lambda. 00101 previous_u = &(add_vector("previous_u")); 00102 00103 // Add a vector to keep track of the temporary solution "y" of Ay=G_{\lambda}. 00104 y = &(add_vector("y")); 00105 00106 // Add a vector to keep track of the "old value" of "y" which is the solution of Ay=G_{\lambda}. 00107 y_old = &(add_vector("y_old")); 00108 00109 // Add a vector to keep track of the temporary solution "z" of Az=-G. 00110 z = &(add_vector("z")); 00111 00112 // Add a vector to keep track of the Newton update during the constrained PDE solves. 00113 delta_u = &(add_vector("delta_u")); 00114 00115 // Call the Parent's initialization routine. 00116 Parent::init_data(); 00117 } 00118 00119 00120 00121 00122 void ContinuationSystem::solve() 00123 { 00124 // Set the Residual RHS mode, and call the normal solve routine. 00125 rhs_mode = Residual; 00126 DifferentiableSystem::solve(); 00127 } 00128 00129 00130 00131 00132 void ContinuationSystem::initialize_tangent() 00133 { 00134 // Be sure the tangent was not already initialized. 00135 libmesh_assert (!tangent_initialized); 00136 00137 // Compute delta_s_zero, the initial arclength travelled during the 00138 // first step. Here we assume that previous_u and lambda_old store 00139 // the previous solution and control parameter. You may need to 00140 // read in an old solution (or solve the non-continuation system) 00141 // first and call save_current_solution() before getting here. 00142 00143 // 1.) Compute delta_s_zero as ||u|| - ||u_old|| + ... 00144 // Compute norms of the current and previous solutions 00145 // Real norm_u = solution->l2_norm(); 00146 // Real norm_previous_u = previous_u->l2_norm(); 00147 00148 // if (!quiet) 00149 // { 00150 // libMesh::out << "norm_u=" << norm_u << std::endl; 00151 // libMesh::out << "norm_previous_u=" << norm_previous_u << std::endl; 00152 // } 00153 00154 // if (norm_u == norm_previous_u) 00155 // { 00156 // libMesh::err << "Warning, it appears u and previous_u are the " 00157 // << "same, are you sure this is correct?" 00158 // << "It's possible you forgot to set one or the other..." 00159 // << std::endl; 00160 // } 00161 00162 // Real delta_s_zero = std::sqrt( 00163 // (norm_u - norm_previous_u)*(norm_u - norm_previous_u) + 00164 // (*continuation_parameter-old_continuation_parameter)* 00165 // (*continuation_parameter-old_continuation_parameter) 00166 // ); 00167 00168 // // 2.) Compute delta_s_zero as ||u -u_old|| + ... 00169 // *delta_u = *solution; 00170 // delta_u->add(-1., *previous_u); 00171 // delta_u->close(); 00172 // Real norm_delta_u = delta_u->l2_norm(); 00173 // Real norm_u = solution->l2_norm(); 00174 // Real norm_previous_u = previous_u->l2_norm(); 00175 00176 // // Scale norm_delta_u by the bigger of either norm_u or norm_previous_u 00177 // norm_delta_u /= std::max(norm_u, norm_previous_u); 00178 00179 // if (!quiet) 00180 // { 00181 // libMesh::out << "norm_u=" << norm_u << std::endl; 00182 // libMesh::out << "norm_previous_u=" << norm_previous_u << std::endl; 00183 // //libMesh::out << "norm_delta_u=" << norm_delta_u << std::endl; 00184 // libMesh::out << "norm_delta_u/max(|u|,|u_old|)=" << norm_delta_u << std::endl; 00185 // libMesh::out << "|norm_u-norm_previous_u|=" << std::abs(norm_u - norm_previous_u) << std::endl; 00186 // } 00187 00188 // const Real dlambda = *continuation_parameter-old_continuation_parameter; 00189 00190 // if (!quiet) 00191 // libMesh::out << "dlambda=" << dlambda << std::endl; 00192 00193 // Real delta_s_zero = std::sqrt( 00194 // (norm_delta_u*norm_delta_u) + 00195 // (dlambda*dlambda) 00196 // ); 00197 00198 // if (!quiet) 00199 // libMesh::out << "delta_s_zero=" << delta_s_zero << std::endl; 00200 00201 // 1.) + 2.) 00202 // // Now approximate the initial tangent d(lambda)/ds 00203 // this->dlambda_ds = (*continuation_parameter-old_continuation_parameter) / delta_s_zero; 00204 00205 00206 // // We can also approximate the deriv. wrt s by finite differences: 00207 // // du/ds = (u1 - u0) / delta_s_zero. 00208 // // FIXME: Use delta_u from above if we decide to keep that method. 00209 // *du_ds = *solution; 00210 // du_ds->add(-1., *previous_u); 00211 // du_ds->scale(1./delta_s_zero); 00212 // du_ds->close(); 00213 00214 00215 // 3.) Treating (u-previous_u)/(lambda - lambda_old) as an approximation to du/d(lambda), 00216 // we follow the same technique as Carnes and Shadid. 00217 // const Real dlambda = *continuation_parameter-old_continuation_parameter; 00218 // libmesh_assert_greater (dlambda, 0.); 00219 00220 // // Use delta_u for temporary calculation of du/d(lambda) 00221 // *delta_u = *solution; 00222 // delta_u->add(-1., *previous_u); 00223 // delta_u->scale(1. / dlambda); 00224 // delta_u->close(); 00225 00226 // // Determine initial normalization parameter 00227 // const Real solution_size = std::max(solution->l2_norm(), previous_u->l2_norm()); 00228 // if (solution_size > 1.) 00229 // { 00230 // Theta = 1./solution_size; 00231 00232 // if (!quiet) 00233 // libMesh::out << "Setting Normalization Parameter Theta=" << Theta << std::endl; 00234 // } 00235 00236 // // Compute d(lambda)/ds 00237 // // The correct sign of d(lambda)/ds should be positive, since we assume that (lambda > lambda_old) 00238 // // but we could always double-check that as well. 00239 // Real norm_delta_u = delta_u->l2_norm(); 00240 // this->dlambda_ds = 1. / std::sqrt(1. + Theta*Theta*norm_delta_u*norm_delta_u); 00241 00242 // // Finally, compute du/ds = d(lambda)/ds * du/d(lambda) 00243 // *du_ds = *delta_u; 00244 // du_ds->scale(dlambda_ds); 00245 // du_ds->close(); 00246 00247 00248 // 4.) Use normalized arclength formula to estimate delta_s_zero 00249 // // Determine initial normalization parameter 00250 // set_Theta(); 00251 00252 // // Compute (normalized) delta_s_zero 00253 // *delta_u = *solution; 00254 // delta_u->add(-1., *previous_u); 00255 // delta_u->close(); 00256 // Real norm_delta_u = delta_u->l2_norm(); 00257 00258 // const Real dlambda = *continuation_parameter-old_continuation_parameter; 00259 00260 // if (!quiet) 00261 // libMesh::out << "dlambda=" << dlambda << std::endl; 00262 00263 // Real delta_s_zero = std::sqrt( 00264 // (Theta_LOCA*Theta_LOCA*Theta*norm_delta_u*norm_delta_u) + 00265 // (dlambda*dlambda) 00266 // ); 00267 // *du_ds = *delta_u; 00268 // du_ds->scale(1./delta_s_zero); 00269 // dlambda_ds = dlambda / delta_s_zero; 00270 00271 // if (!quiet) 00272 // { 00273 // libMesh::out << "delta_s_zero=" << delta_s_zero << std::endl; 00274 // libMesh::out << "initial d(lambda)/ds|_0 = " << dlambda_ds << std::endl; 00275 // libMesh::out << "initial ||du_ds||_0 = " << du_ds->l2_norm() << std::endl; 00276 // } 00277 00278 // // FIXME: Also store the initial finite-differenced approximation to -du/dlambda as y. 00279 // // We stick to the convention of storing negative y, since that is what we typically 00280 // // solve for anyway. 00281 // *y = *delta_u; 00282 // y->scale(-1./dlambda); 00283 // y->close(); 00284 00285 00286 00287 // 5.) Assume dlambda/ds_0 ~ 1/sqrt(2) and determine the value of Theta_LOCA which 00288 // will satisfy this criterion 00289 00290 // Initial change in parameter 00291 const Real dlambda = *continuation_parameter-old_continuation_parameter; 00292 libmesh_assert_not_equal_to (dlambda, 0.0); 00293 00294 // Ideal initial value of dlambda_ds 00295 dlambda_ds = 1. / std::sqrt(2.); 00296 if (dlambda < 0.) 00297 dlambda_ds *= -1.; 00298 00299 // This also implies the initial value of ds 00300 ds_current = dlambda / dlambda_ds; 00301 00302 if (!quiet) 00303 libMesh::out << "Setting ds_current|_0=" << ds_current << std::endl; 00304 00305 // Set y = -du/dlambda using finite difference approximation 00306 *y = *solution; 00307 y->add(-1., *previous_u); 00308 y->scale(-1./dlambda); 00309 y->close(); 00310 const Real ynorm=y->l2_norm(); 00311 00312 // Finally, set the value of du_ds to be used in the upcoming 00313 // tangent calculation. du/ds = du/dlambda * dlambda/ds 00314 *du_ds = *y; 00315 du_ds->scale(-dlambda_ds); 00316 du_ds->close(); 00317 00318 // Determine additional solution normalization parameter 00319 // (Since we just set du/ds, it will be: ||du||*||du/ds||) 00320 set_Theta(); 00321 00322 // The value of Theta_LOCA which makes dlambda_ds = 1/sqrt(2), 00323 // assuming our Theta = ||du||^2. 00324 // Theta_LOCA = std::abs(dlambda); 00325 00326 // Assuming general Theta 00327 Theta_LOCA = std::sqrt(1./Theta/ynorm/ynorm); 00328 00329 00330 if (!quiet) 00331 { 00332 libMesh::out << "Setting initial Theta_LOCA = " << Theta_LOCA << std::endl; 00333 libMesh::out << "Theta_LOCA^2*Theta = " << Theta_LOCA*Theta_LOCA*Theta << std::endl; 00334 libMesh::out << "initial d(lambda)/ds|_0 = " << dlambda_ds << std::endl; 00335 libMesh::out << "initial ||du_ds||_0 = " << du_ds->l2_norm() << std::endl; 00336 } 00337 00338 00339 00340 // OK, we estimated the tangent at point u0. 00341 // Now, to estimate the tangent at point u1, we call the solve_tangent routine. 00342 00343 // Set the flag which tells us the method has been initialized. 00344 tangent_initialized = true; 00345 00346 solve_tangent(); 00347 00348 // Advance the solution and the parameter to the next value. 00349 update_solution(); 00350 } 00351 00352 00353 00354 00355 00356 00357 // This is most of the "guts" of this class. This is where we implement 00358 // our custom Newton iterations and perform most of the solves. 00359 void ContinuationSystem::continuation_solve() 00360 { 00361 // Be sure the user has set the continuation parameter pointer 00362 if (!continuation_parameter) 00363 libmesh_error_msg("You must set the continuation_parameter pointer " \ 00364 << "to a member variable of the derived class, preferably in the " \ 00365 << "Derived class's init_data function. This is how the ContinuationSystem " \ 00366 << "updates the continuation parameter."); 00367 00368 // Use extra precision for all the numbers printed in this function. 00369 std::streamsize old_precision = libMesh::out.precision(); 00370 libMesh::out.precision(16); 00371 libMesh::out.setf(std::ios_base::scientific); 00372 00373 // We can't start solving the augmented PDE system unless the tangent 00374 // vectors have been initialized. This only needs to occur once. 00375 if (!tangent_initialized) 00376 initialize_tangent(); 00377 00378 // Save the old value of -du/dlambda. This will be used after the Newton iterations 00379 // to compute the angle between previous tangent vectors. This cosine of this angle is 00380 // 00381 // tau := abs( (du/d(lambda)_i , du/d(lambda)_{i-1}) / (||du/d(lambda)_i|| * ||du/d(lambda)_{i-1}||) ) 00382 // 00383 // The scaling factor tau (which should vary between 0 and 1) is used to shrink the step-size ds 00384 // when we are approaching a turning point. Note that it can only shrink the step size. 00385 *y_old = *y; 00386 00387 // Set pointer to underlying Newton solver 00388 if (!newton_solver) 00389 newton_solver = cast_ptr<NewtonSolver*> (this->time_solver->diff_solver().get()); 00390 00391 // A pair for catching return values from linear system solves. 00392 std::pair<unsigned int, Real> rval; 00393 00394 // Convergence flag for the entire arcstep 00395 bool arcstep_converged = false; 00396 00397 // Begin loop over arcstep reductions. 00398 for (unsigned int ns=0; ns<n_arclength_reductions; ++ns) 00399 { 00400 if (!quiet) 00401 { 00402 libMesh::out << "Current arclength stepsize, ds_current=" << ds_current << std::endl; 00403 libMesh::out << "Current parameter value, lambda=" << *continuation_parameter << std::endl; 00404 } 00405 00406 // Upon exit from the nonlinear loop, the newton_converged flag 00407 // will tell us the convergence status of Newton's method. 00408 bool newton_converged = false; 00409 00410 // The nonlinear residual before *any* nonlinear steps have been taken. 00411 Real nonlinear_residual_firststep = 0.; 00412 00413 // The nonlinear residual from the current "k" Newton step, before the Newton step 00414 Real nonlinear_residual_beforestep = 0.; 00415 00416 // The nonlinear residual from the current "k" Newton step, after the Newton step 00417 Real nonlinear_residual_afterstep = 0.; 00418 00419 // The linear solver tolerance, can be updated dynamically at each Newton step. 00420 Real current_linear_tolerance = 0.; 00421 00422 // The nonlinear loop 00423 for (newton_step=0; newton_step<newton_solver->max_nonlinear_iterations; ++newton_step) 00424 { 00425 libMesh::out << "\n === Starting Newton step " << newton_step << " ===" << std::endl; 00426 00427 // Set the linear system solver tolerance 00428 // // 1.) Set the current linear tolerance based as a multiple of the current residual of the system. 00429 // const Real residual_multiple = 1.e-4; 00430 // Real current_linear_tolerance = residual_multiple*nonlinear_residual_beforestep; 00431 00432 // // But if the current residual isn't small, don't let the solver exit with zero iterations! 00433 // if (current_linear_tolerance > 1.) 00434 // current_linear_tolerance = residual_multiple; 00435 00436 // 2.) Set the current linear tolerance based on the method based on technique of Eisenstat & Walker. 00437 if (newton_step==0) 00438 { 00439 // At first step, only try reducing the residual by a small amount 00440 current_linear_tolerance = initial_newton_tolerance;//0.01; 00441 } 00442 00443 else 00444 { 00445 // The new tolerance is based on the ratio of the most recent tolerances 00446 const Real alp=0.5*(1.+std::sqrt(5.)); 00447 const Real gam=0.9; 00448 00449 libmesh_assert_not_equal_to (nonlinear_residual_beforestep, 0.0); 00450 libmesh_assert_not_equal_to (nonlinear_residual_afterstep, 0.0); 00451 00452 current_linear_tolerance = std::min(gam*std::pow(nonlinear_residual_afterstep/nonlinear_residual_beforestep, alp), 00453 current_linear_tolerance*current_linear_tolerance 00454 ); 00455 00456 // Don't let it get ridiculously small!! 00457 if (current_linear_tolerance < 1.e-12) 00458 current_linear_tolerance = 1.e-12; 00459 } 00460 00461 if (!quiet) 00462 libMesh::out << "Using current_linear_tolerance=" << current_linear_tolerance << std::endl; 00463 00464 00465 // Assemble the residual (and Jacobian). 00466 rhs_mode = Residual; 00467 assembly(true, // Residual 00468 true); // Jacobian 00469 rhs->close(); 00470 00471 // Save the current nonlinear residual. We don't need to recompute the residual unless 00472 // this is the first step, since it was already computed as part of the convergence check 00473 // at the end of the last loop iteration. 00474 if (newton_step==0) 00475 { 00476 nonlinear_residual_beforestep = rhs->l2_norm(); 00477 00478 // Store the residual before any steps have been taken. This will *not* 00479 // be updated at each step, and can be used to see if any progress has 00480 // been made from the initial residual at later steps. 00481 nonlinear_residual_firststep = nonlinear_residual_beforestep; 00482 00483 const Real old_norm_u = solution->l2_norm(); 00484 libMesh::out << " (before step) ||R||_{L2} = " << nonlinear_residual_beforestep << std::endl; 00485 libMesh::out << " (before step) ||R||_{L2}/||u|| = " << nonlinear_residual_beforestep / old_norm_u << std::endl; 00486 00487 // In rare cases (very small arcsteps), it's possible that the residual is 00488 // already below our absolute linear tolerance. 00489 if (nonlinear_residual_beforestep < solution_tolerance) 00490 { 00491 if (!quiet) 00492 libMesh::out << "Initial guess satisfied linear tolerance, exiting with zero Newton iterations!" << std::endl; 00493 00494 // Since we go straight from here to the solve of the next tangent, we 00495 // have to close the matrix before it can be assembled again. 00496 matrix->close(); 00497 newton_converged=true; 00498 break; // out of Newton iterations, with newton_converged=true 00499 } 00500 } 00501 00502 else 00503 { 00504 nonlinear_residual_beforestep = nonlinear_residual_afterstep; 00505 } 00506 00507 00508 // Solve the linear system G_u*z = G 00509 // Initial guess? 00510 z->zero(); // It seems to be extremely important to zero z here, otherwise the solver quits early. 00511 z->close(); 00512 00513 // It's possible that we have selected the current_linear_tolerance so large that 00514 // a guess of z=zero yields a linear system residual |Az + R| small enough that the 00515 // linear solver exits in zero iterations. If this happens, we will reduce the 00516 // current_linear_tolerance until the linear solver does at least 1 iteration. 00517 do 00518 { 00519 rval = 00520 linear_solver->solve(*matrix, 00521 *z, 00522 *rhs, 00523 //1.e-12, 00524 current_linear_tolerance, 00525 newton_solver->max_linear_iterations); // max linear iterations 00526 00527 if (rval.first==0) 00528 { 00529 if (newton_step==0) 00530 { 00531 libMesh::out << "Repeating initial solve with smaller linear tolerance!" << std::endl; 00532 current_linear_tolerance *= initial_newton_tolerance; // reduce the linear tolerance to force the solver to do some work 00533 } 00534 else 00535 // We shouldn't get here ... it means the linear solver did no work on a Newton 00536 // step other than the first one. If this happens, we need to think more about our 00537 // tolerance selection. 00538 libmesh_error_msg("Linear solver did no work!"); 00539 } 00540 00541 } while (rval.first==0); 00542 00543 00544 if (!quiet) 00545 libMesh::out << " G_u*z = G solver converged at step " 00546 << rval.first 00547 << " linear tolerance = " 00548 << rval.second 00549 << "." 00550 << std::endl; 00551 00552 // Sometimes (I am not sure why) the linear solver exits after zero iterations. 00553 // Perhaps it is hitting PETSc's divergence tolerance dtol??? If this occurs, 00554 // we should break out of the Newton iteration loop because nothing further is 00555 // going to happen... Of course if the tolerance is already small enough after 00556 // zero iterations (how can this happen?!) we should not quit. 00557 if ((rval.first == 0) && (rval.second > current_linear_tolerance*nonlinear_residual_beforestep)) 00558 { 00559 if (!quiet) 00560 libMesh::out << "Linear solver exited in zero iterations!" << std::endl; 00561 00562 // Try to find out the reason for convergence/divergence 00563 linear_solver->print_converged_reason(); 00564 00565 break; // out of Newton iterations 00566 } 00567 00568 // Note: need to scale z by -1 since our code always solves Jx=R 00569 // instead of Jx=-R. 00570 z->scale(-1.); 00571 z->close(); 00572 00573 00574 00575 00576 00577 00578 // Assemble the G_Lambda vector, skip residual. 00579 rhs_mode = G_Lambda; 00580 00581 // Assemble both rhs and Jacobian 00582 assembly(true, // Residual 00583 false); // Jacobian 00584 00585 // Not sure if this is really necessary 00586 rhs->close(); 00587 const Real yrhsnorm=rhs->l2_norm(); 00588 if (yrhsnorm == 0.0) 00589 libmesh_error_msg("||G_Lambda|| = 0"); 00590 00591 // We select a tolerance for the y-system which is based on the inexact Newton 00592 // tolerance but scaled by an extra term proportional to the RHS (which is not -> 0 in this case) 00593 const Real ysystemtol=current_linear_tolerance*(nonlinear_residual_beforestep/yrhsnorm); 00594 if (!quiet) 00595 libMesh::out << "ysystemtol=" << ysystemtol << std::endl; 00596 00597 // Solve G_u*y = G_{\lambda} 00598 // FIXME: Initial guess? This is really a solve for -du/dlambda so we could try 00599 // initializing it with the latest approximation to that... du/dlambda ~ du/ds * ds/dlambda 00600 //*y = *solution; 00601 //y->add(-1., *previous_u); 00602 //y->scale(-1. / (*continuation_parameter - old_continuation_parameter)); // Be careful of divide by zero... 00603 //y->close(); 00604 00605 // const unsigned int max_attempts=1; 00606 // unsigned int attempt=0; 00607 // do 00608 // { 00609 // if (!quiet) 00610 // libMesh::out << "Trying to solve tangent system, attempt " << attempt << std::endl; 00611 00612 rval = 00613 linear_solver->solve(*matrix, 00614 *y, 00615 *rhs, 00616 //1.e-12, 00617 ysystemtol, 00618 newton_solver->max_linear_iterations); // max linear iterations 00619 00620 if (!quiet) 00621 libMesh::out << " G_u*y = G_{lambda} solver converged at step " 00622 << rval.first 00623 << ", linear tolerance = " 00624 << rval.second 00625 << "." 00626 << std::endl; 00627 00628 // Sometimes (I am not sure why) the linear solver exits after zero iterations. 00629 // Perhaps it is hitting PETSc's divergence tolerance dtol??? If this occurs, 00630 // we should break out of the Newton iteration loop because nothing further is 00631 // going to happen... 00632 if ((rval.first == 0) && (rval.second > ysystemtol)) 00633 { 00634 if (!quiet) 00635 libMesh::out << "Linear solver exited in zero iterations!" << std::endl; 00636 00637 break; // out of Newton iterations 00638 } 00639 00640 // ++attempt; 00641 // } while ((attempt<max_attempts) && (rval.first==newton_solver->max_linear_iterations)); 00642 00643 00644 00645 00646 00647 // Compute N, the residual of the arclength constraint eqn. 00648 // Note 1: N(u,lambda,s) := (u-u_{old}, du_ds) + (lambda-lambda_{old}, dlambda_ds) - _ds 00649 // We temporarily use the delta_u vector as a temporary vector for this calculation. 00650 *delta_u = *solution; 00651 delta_u->add(-1., *previous_u); 00652 00653 // First part of the arclength constraint 00654 const Number N1 = Theta_LOCA*Theta_LOCA*Theta*delta_u->dot(*du_ds); 00655 const Number N2 = ((*continuation_parameter) - old_continuation_parameter)*dlambda_ds; 00656 const Number N3 = ds_current; 00657 00658 if (!quiet) 00659 { 00660 libMesh::out << " N1=" << N1 << std::endl; 00661 libMesh::out << " N2=" << N2 << std::endl; 00662 libMesh::out << " N3=" << N3 << std::endl; 00663 } 00664 00665 // The arclength constraint value 00666 const Number N = N1+N2-N3; 00667 00668 if (!quiet) 00669 libMesh::out << " N=" << N << std::endl; 00670 00671 const Number duds_dot_z = du_ds->dot(*z); 00672 const Number duds_dot_y = du_ds->dot(*y); 00673 00674 //libMesh::out << "duds_dot_z=" << duds_dot_z << std::endl; 00675 //libMesh::out << "duds_dot_y=" << duds_dot_y << std::endl; 00676 //libMesh::out << "dlambda_ds=" << dlambda_ds << std::endl; 00677 00678 const Number delta_lambda_numerator = -(N + Theta_LOCA*Theta_LOCA*Theta*duds_dot_z); 00679 const Number delta_lambda_denominator = (dlambda_ds - Theta_LOCA*Theta_LOCA*Theta*duds_dot_y); 00680 00681 libmesh_assert_not_equal_to (delta_lambda_denominator, 0.0); 00682 00683 // Now, we are ready to compute the step delta_lambda 00684 const Number delta_lambda_comp = delta_lambda_numerator / 00685 delta_lambda_denominator; 00686 // Lambda is real-valued 00687 const Real delta_lambda = libmesh_real(delta_lambda_comp); 00688 00689 // Knowing delta_lambda, we are ready to update delta_u 00690 // delta_u = z - delta_lambda*y 00691 delta_u->zero(); 00692 delta_u->add(1., *z); 00693 delta_u->add(-delta_lambda, *y); 00694 delta_u->close(); 00695 00696 // Update the system solution and the continuation parameter. 00697 solution->add(1., *delta_u); 00698 solution->close(); 00699 *continuation_parameter += delta_lambda; 00700 00701 // Did the Newton step actually reduce the residual? 00702 rhs_mode = Residual; 00703 assembly(true, // Residual 00704 false); // Jacobian 00705 rhs->close(); 00706 nonlinear_residual_afterstep = rhs->l2_norm(); 00707 00708 00709 // In a "normal" Newton step, ||du||/||R|| > 1 since the most recent 00710 // step is where you "just were" and the current residual is where 00711 // you are now. It can occur that ||du||/||R|| < 1, but these are 00712 // likely not good cases to attempt backtracking (?). 00713 const Real norm_du_norm_R = delta_u->l2_norm() / nonlinear_residual_afterstep; 00714 if (!quiet) 00715 libMesh::out << " norm_du_norm_R=" << norm_du_norm_R << std::endl; 00716 00717 00718 // Factor to decrease the stepsize by for backtracking 00719 Real newton_stepfactor = 1.; 00720 00721 const bool attempt_backtracking = 00722 (nonlinear_residual_afterstep > solution_tolerance) 00723 && (nonlinear_residual_afterstep > nonlinear_residual_beforestep) 00724 && (n_backtrack_steps>0) 00725 && (norm_du_norm_R > 1.) 00726 ; 00727 00728 // If residual is not reduced, do Newton back tracking. 00729 if (attempt_backtracking) 00730 { 00731 if (!quiet) 00732 libMesh::out << "Newton step did not reduce residual." << std::endl; 00733 00734 // back off the previous step. 00735 solution->add(-1., *delta_u); 00736 solution->close(); 00737 *continuation_parameter -= delta_lambda; 00738 00739 // Backtracking: start cutting the Newton stepsize by halves until 00740 // the new residual is actually smaller... 00741 for (unsigned int backtrack_step=0; backtrack_step<n_backtrack_steps; ++backtrack_step) 00742 { 00743 newton_stepfactor *= 0.5; 00744 00745 if (!quiet) 00746 libMesh::out << "Shrinking step size by " << newton_stepfactor << std::endl; 00747 00748 // Take fractional step 00749 solution->add(newton_stepfactor, *delta_u); 00750 solution->close(); 00751 *continuation_parameter += newton_stepfactor*delta_lambda; 00752 00753 rhs_mode = Residual; 00754 assembly(true, // Residual 00755 false); // Jacobian 00756 rhs->close(); 00757 nonlinear_residual_afterstep = rhs->l2_norm(); 00758 00759 if (!quiet) 00760 libMesh::out << "At shrink step " 00761 << backtrack_step 00762 << ", nonlinear_residual_afterstep=" 00763 << nonlinear_residual_afterstep 00764 << std::endl; 00765 00766 if (nonlinear_residual_afterstep < nonlinear_residual_beforestep) 00767 { 00768 if (!quiet) 00769 libMesh::out << "Backtracking succeeded!" << std::endl; 00770 00771 break; // out of backtracking loop 00772 } 00773 00774 else 00775 { 00776 // Back off that step 00777 solution->add(-newton_stepfactor, *delta_u); 00778 solution->close(); 00779 *continuation_parameter -= newton_stepfactor*delta_lambda; 00780 } 00781 00782 // Save a copy of the solution from before the Newton step. 00783 //UniquePtr<NumericVector<Number> > prior_iterate = solution->clone(); 00784 } 00785 } // end if (attempte_backtracking) 00786 00787 00788 // If we tried backtracking but the residual is still not reduced, print message. 00789 if ((attempt_backtracking) && (nonlinear_residual_afterstep > nonlinear_residual_beforestep)) 00790 { 00791 //libMesh::err << "Backtracking failed." << std::endl; 00792 libMesh::out << "Backtracking failed." << std::endl; 00793 00794 // 1.) Quit, exit program. 00795 //libmesh_error_msg("Backtracking failed!"); 00796 00797 // 2.) Continue with last newton_stepfactor 00798 if (newton_step<3) 00799 { 00800 solution->add(newton_stepfactor, *delta_u); 00801 solution->close(); 00802 *continuation_parameter += newton_stepfactor*delta_lambda; 00803 if (!quiet) 00804 libMesh::out << "Backtracking could not reduce residual ... continuing anyway!" << std::endl; 00805 } 00806 00807 // 3.) Break out of Newton iteration loop with newton_converged = false, 00808 // reduce the arclength stepsize, and try again. 00809 else 00810 { 00811 break; // out of Newton iteration loop, with newton_converged=false 00812 } 00813 } 00814 00815 // Another type of convergence check: suppose the residual has not been reduced 00816 // from its initial value after half of the allowed Newton steps have occurred. 00817 // In our experience, this typically means that it isn't going to converge and 00818 // we could probably save time by dropping out of the Newton iteration loop and 00819 // trying a smaller arcstep. 00820 if (this->newton_progress_check) 00821 { 00822 if ((nonlinear_residual_afterstep > nonlinear_residual_firststep) && 00823 (newton_step+1 > static_cast<unsigned int>(0.5*newton_solver->max_nonlinear_iterations))) 00824 { 00825 libMesh::out << "Progress check failed: the current residual: " 00826 << nonlinear_residual_afterstep 00827 << ", is\n" 00828 << "larger than the initial residual, and half of the allowed\n" 00829 << "number of Newton iterations have elapsed.\n" 00830 << "Exiting Newton iterations with converged==false." << std::endl; 00831 00832 break; // out of Newton iteration loop, newton_converged = false 00833 } 00834 } 00835 00836 // Safety check: Check the current continuation parameter against user-provided min-allowable parameter value 00837 if (*continuation_parameter < min_continuation_parameter) 00838 { 00839 libMesh::out << "Continuation parameter fell below min-allowable value." << std::endl; 00840 break; // out of Newton iteration loop, newton_converged = false 00841 } 00842 00843 // Safety check: Check the current continuation parameter against user-provided max-allowable parameter value 00844 if ( (max_continuation_parameter != 0.0) && 00845 (*continuation_parameter > max_continuation_parameter) ) 00846 { 00847 libMesh::out << "Current continuation parameter value: " 00848 << *continuation_parameter 00849 << " exceeded max-allowable value." 00850 << std::endl; 00851 break; // out of Newton iteration loop, newton_converged = false 00852 } 00853 00854 00855 // Check the convergence of the parameter and the solution. If they are small 00856 // enough, we can break out of the Newton iteration loop. 00857 const Real norm_delta_u = delta_u->l2_norm(); 00858 const Real norm_u = solution->l2_norm(); 00859 libMesh::out << " delta_lambda = " << delta_lambda << std::endl; 00860 libMesh::out << " newton_stepfactor*delta_lambda = " << newton_stepfactor*delta_lambda << std::endl; 00861 libMesh::out << " lambda_current = " << *continuation_parameter << std::endl; 00862 libMesh::out << " ||delta_u|| = " << norm_delta_u << std::endl; 00863 libMesh::out << " ||delta_u||/||u|| = " << norm_delta_u / norm_u << std::endl; 00864 00865 00866 // Evaluate the residual at the current Newton iterate. We don't want to detect 00867 // convergence due to a small Newton step when the residual is still not small. 00868 rhs_mode = Residual; 00869 assembly(true, // Residual 00870 false); // Jacobian 00871 rhs->close(); 00872 const Real norm_residual = rhs->l2_norm(); 00873 libMesh::out << " ||R||_{L2} = " << norm_residual << std::endl; 00874 libMesh::out << " ||R||_{L2}/||u|| = " << norm_residual / norm_u << std::endl; 00875 00876 00877 // FIXME: The norm_delta_u tolerance (at least) should be relative. 00878 // It doesn't make sense to converge a solution whose size is ~ 10^5 to 00879 // a tolerance of 1.e-6. Oh, and we should also probably check the 00880 // (relative) size of the residual as well, instead of just the step. 00881 if ((std::abs(delta_lambda) < continuation_parameter_tolerance) && 00882 //(norm_delta_u < solution_tolerance) && // This is a *very* strict criterion we can probably skip 00883 (norm_residual < solution_tolerance)) 00884 { 00885 if (!quiet) 00886 libMesh::out << "Newton iterations converged!" << std::endl; 00887 00888 newton_converged = true; 00889 break; // out of Newton iterations 00890 } 00891 } // end nonlinear loop 00892 00893 if (!newton_converged) 00894 { 00895 libMesh::out << "Newton iterations of augmented system did not converge!" << std::endl; 00896 00897 // Reduce ds_current, recompute the solution and parameter, and continue to next 00898 // arcstep, if there is one. 00899 ds_current *= 0.5; 00900 00901 // Go back to previous solution and parameter value. 00902 *solution = *previous_u; 00903 *continuation_parameter = old_continuation_parameter; 00904 00905 // Compute new predictor with smaller ds 00906 apply_predictor(); 00907 } 00908 else 00909 { 00910 // Set step convergence and break out 00911 arcstep_converged=true; 00912 break; // out of arclength reduction loop 00913 } 00914 00915 } // end loop over arclength reductions 00916 00917 // Check for convergence of the whole arcstep. If not converged at this 00918 // point, we have no choice but to quit. 00919 if (!arcstep_converged) 00920 libmesh_error_msg("Arcstep failed to converge after max number of reductions! Exiting..."); 00921 00922 // Print converged solution control parameter and max value. 00923 libMesh::out << "lambda_current=" << *continuation_parameter << std::endl; 00924 //libMesh::out << "u_max=" << solution->max() << std::endl; 00925 00926 // Reset old stream precision and flags. 00927 libMesh::out.precision(old_precision); 00928 libMesh::out.unsetf(std::ios_base::scientific); 00929 00930 // Note: we don't want to go on to the next guess yet, since the user may 00931 // want to post-process this data. It's up to the user to call advance_arcstep() 00932 // when they are ready to go on. 00933 } 00934 00935 00936 00937 void ContinuationSystem::advance_arcstep() 00938 { 00939 // Solve for the updated tangent du1/ds, d(lambda1)/ds 00940 solve_tangent(); 00941 00942 // Advance the solution and the parameter to the next value. 00943 update_solution(); 00944 } 00945 00946 00947 00948 // This function solves the tangent system: 00949 // [ G_u G_{lambda} ][(du/ds)_new ] = [ 0 ] 00950 // [ Theta*(du/ds)_old (dlambda/ds)_old ][(dlambda/ds)_new ] [-N_s] 00951 // The solution is given by: 00952 // .) Let G_u y = G_lambda, then 00953 // .) 2nd row yields: 00954 // (dlambda/ds)_new = 1.0 / ( (dlambda/ds)_old - Theta*(du/ds)_old*y ) 00955 // .) 1st row yields 00956 // (du_ds)_new = -(dlambda/ds)_new * y 00957 void ContinuationSystem::solve_tangent() 00958 { 00959 // We shouldn't call this unless the current tangent already makes sense. 00960 libmesh_assert (tangent_initialized); 00961 00962 // Set pointer to underlying Newton solver 00963 if (!newton_solver) 00964 newton_solver = 00965 cast_ptr<NewtonSolver*> (this->time_solver->diff_solver().get()); 00966 00967 // Assemble the system matrix AND rhs, with rhs = G_{\lambda} 00968 this->rhs_mode = G_Lambda; 00969 00970 // Assemble Residual and Jacobian 00971 this->assembly(true, // Residual 00972 true); // Jacobian 00973 00974 // Not sure if this is really necessary 00975 rhs->close(); 00976 00977 // Solve G_u*y = G_{\lambda} 00978 std::pair<unsigned int, Real> rval = 00979 linear_solver->solve(*matrix, 00980 *y, 00981 *rhs, 00982 1.e-12, // relative linear tolerance 00983 2*newton_solver->max_linear_iterations); // max linear iterations 00984 00985 // FIXME: If this doesn't converge at all, the new tangent vector is 00986 // going to be really bad... 00987 00988 if (!quiet) 00989 libMesh::out << "G_u*y = G_{lambda} solver converged at step " 00990 << rval.first 00991 << " linear tolerance = " 00992 << rval.second 00993 << "." 00994 << std::endl; 00995 00996 // Save old solution and parameter tangents for possible use in higher-order 00997 // predictor schemes. 00998 previous_dlambda_ds = dlambda_ds; 00999 *previous_du_ds = *du_ds; 01000 01001 01002 // 1.) Previous, probably wrong, technique! 01003 // // Solve for the updated d(lambda)/ds 01004 // // denom = N_{lambda} - (du_ds)^t y 01005 // // = d(lambda)/ds - (du_ds)^t y 01006 // Real denom = dlambda_ds - du_ds->dot(*y); 01007 01008 // //libMesh::out << "denom=" << denom << std::endl; 01009 // libmesh_assert_not_equal_to (denom, 0.0); 01010 01011 // dlambda_ds = 1.0 / denom; 01012 01013 01014 // if (!quiet) 01015 // libMesh::out << "dlambda_ds=" << dlambda_ds << std::endl; 01016 01017 // // Compute the updated value of du/ds = -_dlambda_ds * y 01018 // du_ds->zero(); 01019 // du_ds->add(-dlambda_ds, *y); 01020 // du_ds->close(); 01021 01022 01023 // 2.) From Brian Carnes' paper... 01024 // According to Carnes, y comes from solving G_u * y = -G_{\lambda} 01025 y->scale(-1.); 01026 const Real ynorm = y->l2_norm(); 01027 dlambda_ds = 1. / std::sqrt(1. + Theta_LOCA*Theta_LOCA*Theta*ynorm*ynorm); 01028 01029 // Determine the correct sign for dlambda_ds. 01030 01031 // We will use delta_u to temporarily compute this sign. 01032 *delta_u = *solution; 01033 delta_u->add(-1., *previous_u); 01034 delta_u->close(); 01035 01036 const Real sgn_dlambda_ds = 01037 libmesh_real(Theta_LOCA*Theta_LOCA*Theta*y->dot(*delta_u) + 01038 (*continuation_parameter-old_continuation_parameter)); 01039 01040 if (sgn_dlambda_ds < 0.) 01041 { 01042 if (!quiet) 01043 libMesh::out << "dlambda_ds is negative." << std::endl; 01044 01045 dlambda_ds *= -1.; 01046 } 01047 01048 // Finally, set the new tangent vector, du/ds = dlambda/ds * y. 01049 du_ds->zero(); 01050 du_ds->add(dlambda_ds, *y); 01051 du_ds->close(); 01052 01053 if (!quiet) 01054 { 01055 libMesh::out << "d(lambda)/ds = " << dlambda_ds << std::endl; 01056 libMesh::out << "||du_ds|| = " << du_ds->l2_norm() << std::endl; 01057 } 01058 01059 // Our next solve expects y ~ -du/dlambda, so scale it back by -1 again now. 01060 y->scale(-1.); 01061 y->close(); 01062 } 01063 01064 01065 01066 void ContinuationSystem::set_Theta() 01067 { 01068 // // Use the norm of the latest solution, squared. 01069 //const Real normu = solution->l2_norm(); 01070 //libmesh_assert_not_equal_to (normu, 0.0); 01071 //Theta = 1./normu/normu; 01072 01073 // // 1.) Use the norm of du, squared 01074 // *delta_u = *solution; 01075 // delta_u->add(-1, *previous_u); 01076 // delta_u->close(); 01077 // const Real normdu = delta_u->l2_norm(); 01078 01079 // if (normdu < 1.) // don't divide by zero or make a huge scaling parameter. 01080 // Theta = 1.; 01081 // else 01082 // Theta = 1./normdu/normdu; 01083 01084 // 2.) Use 1.0, i.e. don't scale 01085 Theta=1.; 01086 01087 // 3.) Use a formula which attempts to make the "solution triangle" isosceles. 01088 // libmesh_assert_less (std::abs(dlambda_ds), 1.); 01089 01090 // *delta_u = *solution; 01091 // delta_u->add(-1, *previous_u); 01092 // delta_u->close(); 01093 // const Real normdu = delta_u->l2_norm(); 01094 01095 // Theta = std::sqrt(1. - dlambda_ds*dlambda_ds) / normdu * tau * ds; 01096 01097 01098 // // 4.) Use the norm of du and the norm of du/ds 01099 // *delta_u = *solution; 01100 // delta_u->add(-1, *previous_u); 01101 // delta_u->close(); 01102 // const Real normdu = delta_u->l2_norm(); 01103 // du_ds->close(); 01104 // const Real normduds = du_ds->l2_norm(); 01105 01106 // if (normduds < 1.e-12) 01107 // { 01108 // libMesh::out << "Setting initial Theta= 1./normdu/normdu" << std::endl; 01109 // libMesh::out << "normdu=" << normdu << std::endl; 01110 01111 // // Don't use this scaling if the solution delta is already O(1) 01112 // if (normdu > 1.) 01113 // Theta = 1./normdu/normdu; 01114 // else 01115 // Theta = 1.; 01116 // } 01117 // else 01118 // { 01119 // libMesh::out << "Setting Theta= 1./normdu/normduds" << std::endl; 01120 // libMesh::out << "normdu=" << normdu << std::endl; 01121 // libMesh::out << "normduds=" << normduds << std::endl; 01122 01123 // // Don't use this scaling if the solution delta is already O(1) 01124 // if ((normdu>1.) || (normduds>1.)) 01125 // Theta = 1./normdu/normduds; 01126 // else 01127 // Theta = 1.; 01128 // } 01129 01130 if (!quiet) 01131 libMesh::out << "Setting Normalization Parameter Theta=" << Theta << std::endl; 01132 } 01133 01134 01135 01136 void ContinuationSystem::set_Theta_LOCA() 01137 { 01138 // We also recompute the LOCA normalization parameter based on the 01139 // most recently computed value of dlambda_ds 01140 // if (!quiet) 01141 // libMesh::out << "(Theta_LOCA) dlambda_ds=" << dlambda_ds << std::endl; 01142 01143 // Formula makes no sense if |dlambda_ds| > 1 01144 libmesh_assert_less (std::abs(dlambda_ds), 1.); 01145 01146 // 1.) Attempt to implement the method in LOCA paper 01147 // const Real g = 1./std::sqrt(2.); // "desired" dlambda_ds 01148 01149 // // According to the LOCA people, we only renormalize for 01150 // // when |dlambda_ds| exceeds some pre-selected maximum (which they take to be zero, btw). 01151 // if (std::abs(dlambda_ds) > .9) 01152 // { 01153 // // Note the *= ... This is updating the previous value of Theta_LOCA 01154 // // Note: The LOCA people actually use Theta_LOCA^2 to normalize their arclength constraint. 01155 // Theta_LOCA *= std::abs( (dlambda_ds/g)*std::sqrt( (1.-g*g) / (1.-dlambda_ds*dlambda_ds) ) ); 01156 01157 // // Suggested max-allowable value for Theta_LOCA 01158 // if (Theta_LOCA > 1.e8) 01159 // { 01160 // Theta_LOCA = 1.e8; 01161 01162 // if (!quiet) 01163 // libMesh::out << "max Theta_LOCA=" << Theta_LOCA << " has been selected." << std::endl; 01164 // } 01165 // } 01166 // else 01167 // Theta_LOCA=1.0; 01168 01169 // 2.) FIXME: Should we do *= or just =? This function is of dlambda_ds is 01170 // < 1, |dlambda_ds| < 1/sqrt(2) ~~ .7071 01171 // > 1, |dlambda_ds| > 1/sqrt(2) ~~ .7071 01172 Theta_LOCA *= std::abs( dlambda_ds / std::sqrt( (1.-dlambda_ds*dlambda_ds) ) ); 01173 01174 // Suggested max-allowable value for Theta_LOCA. I've never come close 01175 // to this value in my code. 01176 if (Theta_LOCA > 1.e8) 01177 { 01178 Theta_LOCA = 1.e8; 01179 01180 if (!quiet) 01181 libMesh::out << "max Theta_LOCA=" << Theta_LOCA << " has been selected." << std::endl; 01182 } 01183 01184 // 3.) Use 1.0, i.e. don't scale 01185 //Theta_LOCA=1.0; 01186 01187 if (!quiet) 01188 libMesh::out << "Setting Theta_LOCA=" << Theta_LOCA << std::endl; 01189 } 01190 01191 01192 01193 void ContinuationSystem::update_solution() 01194 { 01195 // Set some stream formatting flags 01196 std::streamsize old_precision = libMesh::out.precision(); 01197 libMesh::out.precision(16); 01198 libMesh::out.setf(std::ios_base::scientific); 01199 01200 // We must have a tangent that makes sense before we can update the solution. 01201 libmesh_assert (tangent_initialized); 01202 01203 // Compute tau, the stepsize scaling parameter which attempts to 01204 // reduce ds when the angle between the most recent two tangent 01205 // vectors becomes large. tau is actually the (absolute value of 01206 // the) cosine of the angle between these two vectors... so if tau ~ 01207 // 0 the angle is ~ 90 degrees, while if tau ~ 1 the angle is ~ 0 01208 // degrees. 01209 y_old->close(); 01210 y->close(); 01211 const Real yoldnorm = y_old->l2_norm(); 01212 const Real ynorm = y->l2_norm(); 01213 const Number yoldy = y_old->dot(*y); 01214 const Real yold_over_y = yoldnorm/ynorm; 01215 01216 if (!quiet) 01217 { 01218 libMesh::out << "yoldnorm=" << yoldnorm << std::endl; 01219 libMesh::out << "ynorm=" << ynorm << std::endl; 01220 libMesh::out << "yoldy=" << yoldy << std::endl; 01221 libMesh::out << "yoldnorm/ynorm=" << yoldnorm/ynorm << std::endl; 01222 } 01223 01224 // Save the current value of ds before updating it 01225 previous_ds = ds_current; 01226 01227 // // 1.) Cosine method (for some reason this always predicts the angle is ~0) 01228 // // Don't try divinding by zero 01229 // if ((yoldnorm > 1.e-12) && (ynorm > 1.e-12)) 01230 // tau = std::abs(yoldy) / yoldnorm / ynorm; 01231 // else 01232 // tau = 1.; 01233 01234 // // 2.) Relative size of old and new du/dlambda method with cutoff of 0.9 01235 // if ((yold_over_y < 0.9) && (yold_over_y > 1.e-6)) 01236 // tau = yold_over_y; 01237 // else 01238 // tau = 1.; 01239 01240 // 3.) Grow (or shrink) the arclength stepsize by the ratio of du/dlambda, but do not 01241 // exceed the user-specified value of ds. 01242 if (yold_over_y > 1.e-6) 01243 { 01244 // // 1.) Scale current ds by the ratio of successive tangents. 01245 // ds_current *= yold_over_y; 01246 // if (ds_current > ds) 01247 // ds_current = ds; 01248 01249 // 2.) Technique 1 tends to shrink the step fairly well (and even if it doesn't 01250 // get very small, we still have step reduction) but it seems to grow the step 01251 // very slowly. Another possible technique is step-doubling: 01252 // if (yold_over_y > 1.) 01253 // ds_current *= 2.; 01254 // else 01255 // ds_current *= yold_over_y; 01256 01257 // 3.) Technique 2 may over-zealous when we are also using the Newton stepgrowth 01258 // factor. For technique 3 we multiply by yold_over_y unless yold_over_y > 2 01259 // in which case we use 2. 01260 // if (yold_over_y > 2.) 01261 // ds_current *= 2.; 01262 // else 01263 // ds_current *= yold_over_y; 01264 01265 // 4.) Double-or-halve. We double the arc-step if the ratio of successive tangents 01266 // is larger than 'double_threshold', halve it if it is less than 'halve_threshold' 01267 const Real double_threshold = 0.5; 01268 const Real halve_threshold = 0.5; 01269 if (yold_over_y > double_threshold) 01270 ds_current *= 2.; 01271 else if (yold_over_y < halve_threshold) 01272 ds_current *= 0.5; 01273 01274 01275 // Also possibly use the number of Newton iterations required to compute the previous 01276 // step (relative to the maximum-allowed number of Newton iterations) to grow the step. 01277 if (newton_stepgrowth_aggressiveness > 0.) 01278 { 01279 libmesh_assert(newton_solver); 01280 const unsigned int Nmax = newton_solver->max_nonlinear_iterations; 01281 01282 // // The LOCA Newton step growth technique (note: only grows step length) 01283 // const Real stepratio = static_cast<Real>(Nmax-(newton_step+1))/static_cast<Real>(Nmax-1.); 01284 // const Real newtonstep_growthfactor = 1. + newton_stepgrowth_aggressiveness*stepratio*stepratio; 01285 01286 // The "Nopt/N" method, may grow or shrink the step. Assume Nopt=Nmax/2. 01287 const Real newtonstep_growthfactor = 01288 newton_stepgrowth_aggressiveness * 0.5 * 01289 static_cast<Real>(Nmax) / static_cast<Real>(newton_step+1); 01290 01291 if (!quiet) 01292 libMesh::out << "newtonstep_growthfactor=" << newtonstep_growthfactor << std::endl; 01293 01294 ds_current *= newtonstep_growthfactor; 01295 } 01296 } 01297 01298 01299 // Don't let the stepsize get above the user's maximum-allowed stepsize. 01300 if (ds_current > ds) 01301 ds_current = ds; 01302 01303 // Check also for a minimum allowed stepsize. 01304 if (ds_current < ds_min) 01305 { 01306 libMesh::out << "Enforcing minimum-allowed arclength stepsize of " << ds_min << std::endl; 01307 ds_current = ds_min; 01308 } 01309 01310 if (!quiet) 01311 { 01312 libMesh::out << "Current step size: ds_current=" << ds_current << std::endl; 01313 } 01314 01315 // Recompute scaling factor Theta for 01316 // the current solution before updating. 01317 set_Theta(); 01318 01319 // Also, recompute the LOCA scaling factor, which attempts to 01320 // maintain a reasonable value of dlambda/ds 01321 set_Theta_LOCA(); 01322 01323 libMesh::out << "Theta*Theta_LOCA^2=" << Theta*Theta_LOCA*Theta_LOCA << std::endl; 01324 01325 // Based on the asymptotic singular behavior of du/dlambda near simple turning points, 01326 // we can compute a single parameter which may suggest that we are close to a singularity. 01327 *delta_u = *solution; 01328 delta_u->add(-1, *previous_u); 01329 delta_u->close(); 01330 const Real normdu = delta_u->l2_norm(); 01331 const Real C = (std::log (Theta_LOCA*normdu) / 01332 std::log (std::abs(*continuation_parameter-old_continuation_parameter))) - 1.0; 01333 if (!quiet) 01334 libMesh::out << "C=" << C << std::endl; 01335 01336 // Save the current value of u and lambda before updating. 01337 save_current_solution(); 01338 01339 if (!quiet) 01340 { 01341 libMesh::out << "Updating the solution with the tangent guess." << std::endl; 01342 libMesh::out << "||u_old||=" << this->solution->l2_norm() << std::endl; 01343 libMesh::out << "lambda_old=" << *continuation_parameter << std::endl; 01344 } 01345 01346 // Since we solved for the tangent vector, now we can compute an 01347 // initial guess for the new solution, and an initial guess for the 01348 // new value of lambda. 01349 apply_predictor(); 01350 01351 if (!quiet) 01352 { 01353 libMesh::out << "||u_new||=" << this->solution->l2_norm() << std::endl; 01354 libMesh::out << "lambda_new=" << *continuation_parameter << std::endl; 01355 } 01356 01357 // Unset previous stream flags 01358 libMesh::out.precision(old_precision); 01359 libMesh::out.unsetf(std::ios_base::scientific); 01360 } 01361 01362 01363 01364 01365 void ContinuationSystem::save_current_solution() 01366 { 01367 // Save the old solution vector 01368 *previous_u = *solution; 01369 01370 // Save the old value of lambda 01371 old_continuation_parameter = *continuation_parameter; 01372 } 01373 01374 01375 01376 void ContinuationSystem::apply_predictor() 01377 { 01378 if (predictor == Euler) 01379 { 01380 // 1.) Euler Predictor 01381 // Predict next the solution 01382 solution->add(ds_current, *du_ds); 01383 solution->close(); 01384 01385 // Predict next parameter value 01386 *continuation_parameter += ds_current*dlambda_ds; 01387 } 01388 01389 01390 else if (predictor == AB2) 01391 { 01392 // 2.) 2nd-order explicit AB predictor 01393 libmesh_assert_not_equal_to (previous_ds, 0.0); 01394 const Real stepratio = ds_current/previous_ds; 01395 01396 // Build up next solution value. 01397 solution->add( 0.5*ds_current*(2.+stepratio), *du_ds); 01398 solution->add(-0.5*ds_current*stepratio , *previous_du_ds); 01399 solution->close(); 01400 01401 // Next parameter value 01402 *continuation_parameter += 01403 0.5*ds_current*((2.+stepratio)*dlambda_ds - 01404 stepratio*previous_dlambda_ds); 01405 } 01406 01407 else 01408 libmesh_error_msg("Unknown predictor!"); 01409 } 01410 01411 } // namespace libMesh