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