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