twostep_time_solver.C
Go to the documentation of this file.00001 // The libMesh Finite Element Library. 00002 // Copyright (C) 2002-2012 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner 00003 00004 // This library is free software; you can redistribute it and/or 00005 // modify it under the terms of the GNU Lesser General Public 00006 // License as published by the Free Software Foundation; either 00007 // version 2.1 of the License, or (at your option) any later version. 00008 00009 // This library is distributed in the hope that it will be useful, 00010 // but WITHOUT ANY WARRANTY; without even the implied warranty of 00011 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU 00012 // Lesser General Public License for more details. 00013 00014 // You should have received a copy of the GNU Lesser General Public 00015 // License along with this library; if not, write to the Free Software 00016 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA 00017 00018 00019 #include "libmesh/twostep_time_solver.h" 00020 #include "libmesh/diff_system.h" 00021 #include "libmesh/euler_solver.h" 00022 #include "libmesh/numeric_vector.h" 00023 00024 namespace libMesh 00025 { 00026 00027 00028 00029 TwostepTimeSolver::TwostepTimeSolver (sys_type& s) 00030 : AdaptiveTimeSolver(s) 00031 00032 { 00033 // We start with a reasonable time solver: implicit Euler 00034 core_time_solver.reset(new EulerSolver(s)); 00035 } 00036 00037 00038 00039 TwostepTimeSolver::~TwostepTimeSolver () 00040 { 00041 } 00042 00043 00044 00045 void TwostepTimeSolver::solve() 00046 { 00047 libmesh_assert(core_time_solver.get()); 00048 00049 // The core_time_solver will handle any first_solve actions 00050 first_solve = false; 00051 00052 // We may have to repeat timesteps entirely if our error is bad 00053 // enough 00054 bool max_tolerance_met = false; 00055 00056 // Calculating error values each time 00057 Real single_norm(0.), double_norm(0.), error_norm(0.), 00058 relative_error(0.); 00059 00060 while (!max_tolerance_met) 00061 { 00062 // If we've been asked to reduce deltat if necessary, make sure 00063 // the core timesolver does so 00064 core_time_solver->reduce_deltat_on_diffsolver_failure = 00065 this->reduce_deltat_on_diffsolver_failure; 00066 00067 if (!quiet) 00068 { 00069 libMesh::out << "\n === Computing adaptive timestep === " 00070 << std::endl; 00071 } 00072 00073 // Use the double-length timestep first (so the 00074 // old_nonlinear_solution won't have to change) 00075 core_time_solver->solve(); 00076 00077 // Save a copy of the double-length nonlinear solution 00078 // and the old nonlinear solution 00079 AutoPtr<NumericVector<Number> > double_solution = 00080 _system.solution->clone(); 00081 AutoPtr<NumericVector<Number> > old_solution = 00082 _system.get_vector("_old_nonlinear_solution").clone(); 00083 00084 double_norm = calculate_norm(_system, *double_solution); 00085 if (!quiet) 00086 { 00087 libMesh::out << "Double norm = " << double_norm << std::endl; 00088 } 00089 00090 // Then reset the initial guess for our single-length calcs 00091 *(_system.solution) = _system.get_vector("_old_nonlinear_solution"); 00092 00093 // Call two single-length timesteps 00094 // Be sure that the core_time_solver does not change the 00095 // timestep here. (This is unlikely because it just succeeded 00096 // with a timestep twice as large!) 00097 // FIXME: even if diffsolver failure is unlikely, we ought to 00098 // do *something* if it happens 00099 core_time_solver->reduce_deltat_on_diffsolver_failure = 0; 00100 00101 Real old_time = _system.time; 00102 Real old_deltat = _system.deltat; 00103 _system.deltat *= 0.5; 00104 core_time_solver->solve(); 00105 core_time_solver->advance_timestep(); 00106 core_time_solver->solve(); 00107 00108 single_norm = calculate_norm(_system, *_system.solution); 00109 if (!quiet) 00110 { 00111 libMesh::out << "Single norm = " << single_norm << std::endl; 00112 } 00113 00114 // Reset the core_time_solver's reduce_deltat... value. 00115 core_time_solver->reduce_deltat_on_diffsolver_failure = 00116 this->reduce_deltat_on_diffsolver_failure; 00117 00118 // But then back off just in case our advance_timestep() isn't 00119 // called. 00120 // FIXME: this probably doesn't work with multistep methods 00121 _system.get_vector("_old_nonlinear_solution") = *old_solution; 00122 _system.time = old_time; 00123 _system.deltat = old_deltat; 00124 00125 // Find the relative error 00126 *double_solution -= *(_system.solution); 00127 error_norm = calculate_norm(_system, *double_solution); 00128 relative_error = error_norm / _system.deltat / 00129 std::max(double_norm, single_norm); 00130 00131 // If the relative error makes no sense, we're done 00132 if (!double_norm && !single_norm) 00133 return; 00134 00135 if (!quiet) 00136 { 00137 libMesh::out << "Error norm = " << error_norm << std::endl; 00138 libMesh::out << "Local relative error = " 00139 << (error_norm / 00140 std::max(double_norm, single_norm)) 00141 << std::endl; 00142 libMesh::out << "Global relative error = " 00143 << (error_norm / _system.deltat / 00144 std::max(double_norm, single_norm)) 00145 << std::endl; 00146 libMesh::out << "old delta t = " << _system.deltat << std::endl; 00147 } 00148 00149 // If we haven't met our upper error tolerance, we'll have to 00150 // repeat this timestep entirely 00151 if (this->upper_tolerance && relative_error > this->upper_tolerance) 00152 { 00153 // Reset the initial guess for our next try 00154 *(_system.solution) = 00155 _system.get_vector("_old_nonlinear_solution"); 00156 00157 // Chop delta t in half 00158 _system.deltat /= 2.; 00159 00160 if (!quiet) 00161 { 00162 libMesh::out << "Failed to meet upper error tolerance" 00163 << std::endl; 00164 libMesh::out << "Retrying with delta t = " 00165 << _system.deltat << std::endl; 00166 } 00167 } 00168 else 00169 max_tolerance_met = true; 00170 } 00171 00172 00173 // Otherwise, compare the relative error to the tolerance 00174 // and adjust deltat 00175 last_deltat = _system.deltat; 00176 00177 const Real global_shrink_or_growth_factor = 00178 std::pow(this->target_tolerance / relative_error, 00179 static_cast<Real>(1. / core_time_solver->error_order())); 00180 00181 const Real local_shrink_or_growth_factor = 00182 std::pow(this->target_tolerance / 00183 (error_norm/std::max(double_norm, single_norm)), 00184 static_cast<Real>(1. / (core_time_solver->error_order()+1.))); 00185 00186 if (!quiet) 00187 { 00188 libMesh::out << "The global growth/shrink factor is: " 00189 << global_shrink_or_growth_factor << std::endl; 00190 libMesh::out << "The local growth/shrink factor is: " 00191 << local_shrink_or_growth_factor << std::endl; 00192 } 00193 00194 // The local s.o.g. factor is based on the expected **local** 00195 // truncation error for the timestepping method, the global 00196 // s.o.g. factor is based on the method's **global** truncation 00197 // error. You can shrink/grow the timestep to attempt to satisfy 00198 // either a global or local time-discretization error tolerance. 00199 00200 Real shrink_or_growth_factor = 00201 this->global_tolerance ? global_shrink_or_growth_factor : 00202 local_shrink_or_growth_factor; 00203 00204 if (this->max_growth && this->max_growth < shrink_or_growth_factor) 00205 { 00206 if (!quiet && this->global_tolerance) 00207 { 00208 libMesh::out << "delta t is constrained by max_growth" << std::endl; 00209 } 00210 shrink_or_growth_factor = this->max_growth; 00211 } 00212 00213 _system.deltat *= shrink_or_growth_factor; 00214 00215 // Restrict deltat to max-allowable value if necessary 00216 if ((this->max_deltat != 0.0) && (_system.deltat > this->max_deltat)) 00217 { 00218 if (!quiet) 00219 { 00220 libMesh::out << "delta t is constrained by maximum-allowable delta t." 00221 << std::endl; 00222 } 00223 _system.deltat = this->max_deltat; 00224 } 00225 00226 // Restrict deltat to min-allowable value if necessary 00227 if ((this->min_deltat != 0.0) && (_system.deltat < this->min_deltat)) 00228 { 00229 if (!quiet) 00230 { 00231 libMesh::out << "delta t is constrained by minimum-allowable delta t." 00232 << std::endl; 00233 } 00234 _system.deltat = this->min_deltat; 00235 } 00236 00237 if (!quiet) 00238 { 00239 libMesh::out << "new delta t = " << _system.deltat << std::endl; 00240 } 00241 } 00242 00243 } // namespace libMesh
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:49 UTC
Hosted By: