1 2 static char help[] = "Newton's method to solve a two-variable system that comes from the Rosenbrock function.\n\n"; 3 4 /* 5 Include "petscsnes.h" so that we can use SNES solvers. Note that this 6 file automatically includes: 7 petscsys.h - base PETSc routines petscvec.h - vectors 8 petscmat.h - matrices 9 petscis.h - index sets petscksp.h - Krylov subspace methods 10 petscviewer.h - viewers petscpc.h - preconditioners 11 petscksp.h - linear solvers 12 */ 13 #include <petscsnes.h> 14 15 extern PetscErrorCode FormJacobian1(SNES, Vec, Mat, Mat, void *); 16 extern PetscErrorCode FormFunction1(SNES, Vec, Vec, void *); 17 18 int main(int argc, char **argv) { 19 SNES snes; /* nonlinear solver context */ 20 Vec x, r; /* solution, residual vectors */ 21 Mat J; /* Jacobian matrix */ 22 PetscInt its; 23 PetscScalar *xx; 24 SNESConvergedReason reason; 25 26 PetscFunctionBeginUser; 27 PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 28 29 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 30 Create nonlinear solver context 31 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 32 PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 33 34 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 35 Create matrix and vector data structures; set corresponding routines 36 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 37 /* 38 Create vectors for solution and nonlinear function 39 */ 40 PetscCall(VecCreate(PETSC_COMM_WORLD, &x)); 41 PetscCall(VecSetSizes(x, PETSC_DECIDE, 2)); 42 PetscCall(VecSetFromOptions(x)); 43 PetscCall(VecDuplicate(x, &r)); 44 45 /* 46 Create Jacobian matrix data structure 47 */ 48 PetscCall(MatCreate(PETSC_COMM_WORLD, &J)); 49 PetscCall(MatSetSizes(J, PETSC_DECIDE, PETSC_DECIDE, 2, 2)); 50 PetscCall(MatSetFromOptions(J)); 51 PetscCall(MatSetUp(J)); 52 53 /* 54 Set function evaluation routine and vector. 55 */ 56 PetscCall(SNESSetFunction(snes, r, FormFunction1, NULL)); 57 58 /* 59 Set Jacobian matrix data structure and Jacobian evaluation routine 60 */ 61 PetscCall(SNESSetJacobian(snes, J, J, FormJacobian1, NULL)); 62 63 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 64 Customize nonlinear solver; set runtime options 65 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 66 PetscCall(SNESSetFromOptions(snes)); 67 68 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 69 Evaluate initial guess; then solve nonlinear system 70 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 71 PetscCall(VecGetArray(x, &xx)); 72 xx[0] = -1.2; 73 xx[1] = 1.0; 74 PetscCall(VecRestoreArray(x, &xx)); 75 76 /* 77 Note: The user should initialize the vector, x, with the initial guess 78 for the nonlinear solver prior to calling SNESSolve(). In particular, 79 to employ an initial guess of zero, the user should explicitly set 80 this vector to zero by calling VecSet(). 81 */ 82 83 PetscCall(SNESSolve(snes, NULL, x)); 84 PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD)); 85 PetscCall(SNESGetIterationNumber(snes, &its)); 86 PetscCall(SNESGetConvergedReason(snes, &reason)); 87 /* 88 Some systems computes a residual that is identically zero, thus converging 89 due to FNORM_ABS, others converge due to FNORM_RELATIVE. Here, we only 90 report the reason if the iteration did not converge so that the tests are 91 reproducible. 92 */ 93 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s number of SNES iterations = %" PetscInt_FMT "\n", reason > 0 ? "CONVERGED" : SNESConvergedReasons[reason], its)); 94 95 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 96 Free work space. All PETSc objects should be destroyed when they 97 are no longer needed. 98 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 99 100 PetscCall(VecDestroy(&x)); 101 PetscCall(VecDestroy(&r)); 102 PetscCall(MatDestroy(&J)); 103 PetscCall(SNESDestroy(&snes)); 104 PetscCall(PetscFinalize()); 105 return 0; 106 } 107 /* ------------------------------------------------------------------- */ 108 /* 109 FormFunction1 - Evaluates nonlinear function, F(x). 110 111 Input Parameters: 112 . snes - the SNES context 113 . x - input vector 114 . ctx - optional user-defined context 115 116 Output Parameter: 117 . f - function vector 118 */ 119 PetscErrorCode FormFunction1(SNES snes, Vec x, Vec f, void *ctx) { 120 PetscScalar *ff; 121 const PetscScalar *xx; 122 123 /* 124 Get pointers to vector data. 125 - For default PETSc vectors, VecGetArray() returns a pointer to 126 the data array. Otherwise, the routine is implementation dependent. 127 - You MUST call VecRestoreArray() when you no longer need access to 128 the array. 129 */ 130 PetscCall(VecGetArrayRead(x, &xx)); 131 PetscCall(VecGetArray(f, &ff)); 132 133 /* Compute function */ 134 ff[0] = -2.0 + 2.0 * xx[0] + 400.0 * xx[0] * xx[0] * xx[0] - 400.0 * xx[0] * xx[1]; 135 ff[1] = -200.0 * xx[0] * xx[0] + 200.0 * xx[1]; 136 137 /* Restore vectors */ 138 PetscCall(VecRestoreArrayRead(x, &xx)); 139 PetscCall(VecRestoreArray(f, &ff)); 140 return 0; 141 } 142 /* ------------------------------------------------------------------- */ 143 /* 144 FormJacobian1 - Evaluates Jacobian matrix. 145 146 Input Parameters: 147 . snes - the SNES context 148 . x - input vector 149 . dummy - optional user-defined context (not used here) 150 151 Output Parameters: 152 . jac - Jacobian matrix 153 . B - optionally different preconditioning matrix 154 . flag - flag indicating matrix structure 155 */ 156 PetscErrorCode FormJacobian1(SNES snes, Vec x, Mat jac, Mat B, void *dummy) { 157 const PetscScalar *xx; 158 PetscScalar A[4]; 159 PetscInt idx[2] = {0, 1}; 160 161 /* 162 Get pointer to vector data 163 */ 164 PetscCall(VecGetArrayRead(x, &xx)); 165 166 /* 167 Compute Jacobian entries and insert into matrix. 168 - Since this is such a small problem, we set all entries for 169 the matrix at once. 170 */ 171 A[0] = 2.0 + 1200.0 * xx[0] * xx[0] - 400.0 * xx[1]; 172 A[1] = -400.0 * xx[0]; 173 A[2] = -400.0 * xx[0]; 174 A[3] = 200; 175 PetscCall(MatSetValues(B, 2, idx, 2, idx, A, INSERT_VALUES)); 176 177 /* 178 Restore vector 179 */ 180 PetscCall(VecRestoreArrayRead(x, &xx)); 181 182 /* 183 Assemble matrix 184 */ 185 PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 186 PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 187 if (jac != B) { 188 PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY)); 189 PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY)); 190 } 191 return 0; 192 } 193 194 /*TEST 195 196 test: 197 suffix: 1 198 args: -snes_monitor_short -snes_max_it 1000 199 requires: !single 200 201 test: 202 suffix: 2 203 args: -snes_monitor_short -snes_max_it 1000 -snes_type newtontrdc -snes_trdc_use_cauchy false 204 requires: !single 205 206 TEST*/ 207