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