00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030 #include "SimpleNewtonNonlinearSolver.hpp"
00031 #include <iostream>
00032 #include <cassert>
00033 #include "Exception.hpp"
00034
00035 SimpleNewtonNonlinearSolver::SimpleNewtonNonlinearSolver(double linearSolverRelativeTolerance)
00036 : mLinearSolverRelativeTolerance(linearSolverRelativeTolerance),
00037 mTolerance(1e-5),
00038 mWriteStats(false)
00039 {
00040 mTestDampingValues.push_back(-0.1);
00041 mTestDampingValues.push_back(0.05);
00042 for (unsigned i=1; i<=12; i++)
00043 {
00044 double val = double(i)/10;
00045 mTestDampingValues.push_back(val);
00046 }
00047 }
00048
00049 SimpleNewtonNonlinearSolver::~SimpleNewtonNonlinearSolver()
00050 {}
00051
00052 Vec SimpleNewtonNonlinearSolver::Solve(PetscErrorCode (*pComputeResidual)(SNES,Vec,Vec,void*),
00053 PetscErrorCode (*pComputeJacobian)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),
00054 Vec initialGuess,
00055 unsigned fill,
00056 void* pContext)
00057 {
00058 PetscInt size;
00059 VecGetSize(initialGuess, &size);
00060
00061 Vec current_solution;
00062 VecDuplicate(initialGuess, ¤t_solution);
00063 VecCopy(initialGuess, current_solution);
00064
00065 LinearSystem linear_system(current_solution, fill);
00066
00067 (*pComputeResidual)(NULL, current_solution, linear_system.rGetRhsVector(), pContext);
00068
00069
00070 double residual_norm;
00071 VecNorm(linear_system.rGetRhsVector(), NORM_2, &residual_norm);
00072 double scaled_residual_norm = residual_norm/size;
00073
00074 if (mWriteStats)
00075 {
00076 std::cout << "Newton's method:\n Initial ||residual||/N = " << scaled_residual_norm
00077 << "\n Attempting to solve to tolerance " << mTolerance << "..\n";
00078 }
00079
00080 double old_scaled_residual_norm;
00081 unsigned counter = 0;
00082 while (scaled_residual_norm > mTolerance)
00083 {
00084 counter++;
00085
00086
00087 old_scaled_residual_norm = scaled_residual_norm;
00088
00089
00090 (*pComputeJacobian)(NULL, current_solution, &(linear_system.rGetLhsMatrix()), NULL, NULL, pContext);
00091
00092 Vec negative_update = linear_system.Solve();
00093
00094
00095 Vec test_vec;
00096 VecDuplicate(initialGuess, &test_vec);
00097
00098 double best_damping_factor = 1.0;
00099 double best_scaled_residual = 1e20;
00100
00101
00102
00103 for (unsigned i=0; i<mTestDampingValues.size(); i++)
00104 {
00105
00106 PetscVecTools::WAXPY(test_vec,-mTestDampingValues[i],negative_update,current_solution);
00107
00108
00109 linear_system.ZeroLinearSystem();
00110 (*pComputeResidual)(NULL, test_vec, linear_system.rGetRhsVector(), pContext);
00111 VecNorm(linear_system.rGetRhsVector(), NORM_2, &residual_norm);
00112 scaled_residual_norm = residual_norm/size;
00113
00114 if (scaled_residual_norm < best_scaled_residual)
00115 {
00116 best_scaled_residual = scaled_residual_norm;
00117 best_damping_factor = mTestDampingValues[i];
00118 }
00119 }
00120 VecDestroy(test_vec);
00121
00122
00123
00124 if (best_scaled_residual > old_scaled_residual_norm)
00125 {
00126
00127 VecDestroy(current_solution);
00128 VecDestroy(negative_update);
00129
00130 std::stringstream error_message;
00131 error_message << "Iteration " << counter << ", unable to find damping factor such that residual decreases in update direction";
00132 EXCEPTION(error_message.str());
00133 }
00134
00135 if (mWriteStats)
00136 {
00137 std::cout << " Best damping factor = " << best_damping_factor << "\n";
00138 }
00139
00140
00141
00142 PetscVecTools::AddScaledVector(current_solution, negative_update, -best_damping_factor);
00143 scaled_residual_norm = best_scaled_residual;
00144 VecDestroy(negative_update);
00145
00146
00147 linear_system.ZeroLinearSystem();
00148 (*pComputeResidual)(NULL, current_solution, linear_system.rGetRhsVector(), pContext);
00149
00150 if (mWriteStats)
00151 {
00152 std::cout << " Iteration " << counter << ": ||residual||/N = " << scaled_residual_norm << "\n";
00153 }
00154 }
00155
00156 if (mWriteStats)
00157 {
00158 std::cout << " ..solved!\n\n";
00159 }
00160
00161 return current_solution;
00162 }
00163
00164
00165 void SimpleNewtonNonlinearSolver::SetTolerance(double tolerance)
00166 {
00167 assert(tolerance > 0);
00168 mTolerance = tolerance;
00169 }
00170