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