1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests TSLINESEARCHL2 handing of Inf/Nan.\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /*T 5c4762a1bSJed Brown Concepts: SNES^basic example 6c4762a1bSJed Brown T*/ 7c4762a1bSJed Brown 8c4762a1bSJed Brown 9c4762a1bSJed Brown 10c4762a1bSJed Brown /* 11c4762a1bSJed Brown Include "petscsnes.h" so that we can use SNES solvers. Note that this 12c4762a1bSJed Brown file automatically includes: 13c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 14c4762a1bSJed Brown petscmat.h - matrices 15c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 16c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 17c4762a1bSJed Brown petscksp.h - linear solvers 18c4762a1bSJed Brown */ 19c4762a1bSJed Brown /*F 20c4762a1bSJed Brown This examples solves either 21c4762a1bSJed Brown \begin{equation} 22c4762a1bSJed Brown F\genfrac{(}{)}{0pt}{}{x_0}{x_1} = \genfrac{(}{)}{0pt}{}{\sin(3 x_0) + x_0}{x_1} 23c4762a1bSJed Brown \end{equation} 24c4762a1bSJed Brown F*/ 25c4762a1bSJed Brown #include <petscsnes.h> 26c4762a1bSJed Brown 27c4762a1bSJed Brown /* 28c4762a1bSJed Brown User-defined routines 29c4762a1bSJed Brown */ 30c4762a1bSJed Brown extern PetscErrorCode FormJacobian2(SNES,Vec,Mat,Mat,void*); 31c4762a1bSJed Brown extern PetscErrorCode FormFunction2(SNES,Vec,Vec,void*); 32c4762a1bSJed Brown extern PetscErrorCode FormObjective(SNES,Vec,PetscReal*,void*); 33c4762a1bSJed Brown 34c4762a1bSJed Brown /* 35c4762a1bSJed Brown This is a very hacking way to trigger the objective function generating an infinity at a particular count to the call FormObjective(). 36c4762a1bSJed Brown Different line searches evaluate the full step at different counts. For l2 it is the third call (infatcount == 2) while for bt it is the second call. 37c4762a1bSJed Brown */ 38c4762a1bSJed Brown PetscInt infatcount = 0; 39c4762a1bSJed Brown 40c4762a1bSJed Brown int main(int argc,char **argv) 41c4762a1bSJed Brown { 42c4762a1bSJed Brown SNES snes; /* nonlinear solver context */ 43c4762a1bSJed Brown KSP ksp; /* linear solver context */ 44c4762a1bSJed Brown PC pc; /* preconditioner context */ 45c4762a1bSJed Brown Vec x,r; /* solution, residual vectors */ 46c4762a1bSJed Brown Mat J; /* Jacobian matrix */ 47c4762a1bSJed Brown PetscErrorCode ierr; 48c4762a1bSJed Brown PetscInt its; 49c4762a1bSJed Brown PetscMPIInt size; 50c4762a1bSJed Brown PetscScalar *xx; 51c4762a1bSJed Brown PetscBool flg; 52c4762a1bSJed Brown char type[256]; 53c4762a1bSJed Brown 54c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 55*589a23caSBarry Smith ierr = PetscOptionsGetString(NULL,NULL,"-snes_linesearch_type",type,sizeof(type),&flg);CHKERRQ(ierr); 56c4762a1bSJed Brown if (flg) { 57c4762a1bSJed Brown ierr = PetscStrcmp(type,SNESLINESEARCHBT,&flg);CHKERRQ(ierr); 58c4762a1bSJed Brown if (flg) infatcount = 1; 59c4762a1bSJed Brown ierr = PetscStrcmp(type,SNESLINESEARCHL2,&flg);CHKERRQ(ierr); 60c4762a1bSJed Brown if (flg) infatcount = 2; 61c4762a1bSJed Brown } 62c4762a1bSJed Brown 63c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 64c4762a1bSJed Brown if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Example is only for sequential runs"); 65c4762a1bSJed Brown 66c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 67c4762a1bSJed Brown Create nonlinear solver context 68c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 69c4762a1bSJed Brown ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr); 70c4762a1bSJed Brown 71c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 72c4762a1bSJed Brown Create matrix and vector data structures; set corresponding routines 73c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 74c4762a1bSJed Brown /* 75c4762a1bSJed Brown Create vectors for solution and nonlinear function 76c4762a1bSJed Brown */ 77c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr); 78c4762a1bSJed Brown ierr = VecSetSizes(x,PETSC_DECIDE,2);CHKERRQ(ierr); 79c4762a1bSJed Brown ierr = VecSetFromOptions(x);CHKERRQ(ierr); 80c4762a1bSJed Brown ierr = VecDuplicate(x,&r);CHKERRQ(ierr); 81c4762a1bSJed Brown 82c4762a1bSJed Brown /* 83c4762a1bSJed Brown Create Jacobian matrix data structure 84c4762a1bSJed Brown */ 85c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&J);CHKERRQ(ierr); 86c4762a1bSJed Brown ierr = MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,2,2);CHKERRQ(ierr); 87c4762a1bSJed Brown ierr = MatSetFromOptions(J);CHKERRQ(ierr); 88c4762a1bSJed Brown ierr = MatSetUp(J);CHKERRQ(ierr); 89c4762a1bSJed Brown 90c4762a1bSJed Brown ierr = SNESSetFunction(snes,r,FormFunction2,NULL);CHKERRQ(ierr); 91c4762a1bSJed Brown ierr = SNESSetObjective(snes,FormObjective,NULL);CHKERRQ(ierr); 92c4762a1bSJed Brown ierr = SNESSetJacobian(snes,J,J,FormJacobian2,NULL);CHKERRQ(ierr); 93c4762a1bSJed Brown 94c4762a1bSJed Brown 95c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 96c4762a1bSJed Brown Customize nonlinear solver; set runtime options 97c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 98c4762a1bSJed Brown /* 99c4762a1bSJed Brown Set linear solver defaults for this problem. By extracting the 100c4762a1bSJed Brown KSP and PC contexts from the SNES context, we can then 101c4762a1bSJed Brown directly call any KSP and PC routines to set various options. 102c4762a1bSJed Brown */ 103c4762a1bSJed Brown ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 104c4762a1bSJed Brown ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 105c4762a1bSJed Brown ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 106c4762a1bSJed Brown ierr = KSPSetTolerances(ksp,1.e-4,PETSC_DEFAULT,PETSC_DEFAULT,20);CHKERRQ(ierr); 107c4762a1bSJed Brown 108c4762a1bSJed Brown /* 109c4762a1bSJed Brown Set SNES/KSP/KSP/PC runtime options, e.g., 110c4762a1bSJed Brown -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc> 111c4762a1bSJed Brown These options will override those specified above as long as 112c4762a1bSJed Brown SNESSetFromOptions() is called _after_ any other customization 113c4762a1bSJed Brown routines. 114c4762a1bSJed Brown */ 115c4762a1bSJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 116c4762a1bSJed Brown 117c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 118c4762a1bSJed Brown Evaluate initial guess; then solve nonlinear system 119c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 120c4762a1bSJed Brown ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 121c4762a1bSJed Brown xx[0] = 2.0; xx[1] = 3.0; 122c4762a1bSJed Brown ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 123c4762a1bSJed Brown 124c4762a1bSJed Brown /* 125c4762a1bSJed Brown Note: The user should initialize the vector, x, with the initial guess 126c4762a1bSJed Brown for the nonlinear solver prior to calling SNESSolve(). In particular, 127c4762a1bSJed Brown to employ an initial guess of zero, the user should explicitly set 128c4762a1bSJed Brown this vector to zero by calling VecSet(). 129c4762a1bSJed Brown */ 130c4762a1bSJed Brown 131c4762a1bSJed Brown ierr = SNESSolve(snes,NULL,x);CHKERRQ(ierr); 132c4762a1bSJed Brown ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr); 133c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Number of SNES iterations = %D\n",its);CHKERRQ(ierr); 134c4762a1bSJed Brown 135c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 136c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 137c4762a1bSJed Brown are no longer needed. 138c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 139c4762a1bSJed Brown 140c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); ierr = VecDestroy(&r);CHKERRQ(ierr); 141c4762a1bSJed Brown ierr = MatDestroy(&J);CHKERRQ(ierr); ierr = SNESDestroy(&snes);CHKERRQ(ierr); 142c4762a1bSJed Brown ierr = PetscFinalize(); 143c4762a1bSJed Brown return ierr; 144c4762a1bSJed Brown } 145c4762a1bSJed Brown 146c4762a1bSJed Brown PetscErrorCode FormObjective(SNES snes,Vec x,PetscReal *f,void *dummy) 147c4762a1bSJed Brown { 148c4762a1bSJed Brown PetscErrorCode ierr; 149c4762a1bSJed Brown Vec F; 150c4762a1bSJed Brown static PetscInt cnt = 0; 151c4762a1bSJed Brown 152c4762a1bSJed Brown if (cnt++ == infatcount) *f = INFINITY; 153c4762a1bSJed Brown else { 154c4762a1bSJed Brown ierr = VecDuplicate(x,&F);CHKERRQ(ierr); 155c4762a1bSJed Brown ierr = FormFunction2(snes,x,F,dummy);CHKERRQ(ierr); 156c4762a1bSJed Brown ierr = VecNorm(F,NORM_2,f);CHKERRQ(ierr); 157c4762a1bSJed Brown ierr = VecDestroy(&F);CHKERRQ(ierr); 158c4762a1bSJed Brown *f = (*f)*(*f); 159c4762a1bSJed Brown } 160c4762a1bSJed Brown return 0; 161c4762a1bSJed Brown } 162c4762a1bSJed Brown 163c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 164c4762a1bSJed Brown PetscErrorCode FormFunction2(SNES snes,Vec x,Vec f,void *dummy) 165c4762a1bSJed Brown { 166c4762a1bSJed Brown PetscErrorCode ierr; 167c4762a1bSJed Brown const PetscScalar *xx; 168c4762a1bSJed Brown PetscScalar *ff; 169c4762a1bSJed Brown 170c4762a1bSJed Brown /* 171c4762a1bSJed Brown Get pointers to vector data. 172c4762a1bSJed Brown - For default PETSc vectors, VecGetArray() returns a pointer to 173c4762a1bSJed Brown the data array. Otherwise, the routine is implementation dependent. 174c4762a1bSJed Brown - You MUST call VecRestoreArray() when you no longer need access to 175c4762a1bSJed Brown the array. 176c4762a1bSJed Brown */ 177c4762a1bSJed Brown ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); 178c4762a1bSJed Brown ierr = VecGetArray(f,&ff);CHKERRQ(ierr); 179c4762a1bSJed Brown 180c4762a1bSJed Brown /* 181c4762a1bSJed Brown Compute function 182c4762a1bSJed Brown */ 183c4762a1bSJed Brown ff[0] = PetscSinScalar(3.0*xx[0]) + xx[0]; 184c4762a1bSJed Brown ff[1] = xx[1]; 185c4762a1bSJed Brown 186c4762a1bSJed Brown /* 187c4762a1bSJed Brown Restore vectors 188c4762a1bSJed Brown */ 189c4762a1bSJed Brown ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr); 190c4762a1bSJed Brown ierr = VecRestoreArray(f,&ff);CHKERRQ(ierr); 191c4762a1bSJed Brown return 0; 192c4762a1bSJed Brown } 193c4762a1bSJed Brown /* ------------------------------------------------------------------- */ 194c4762a1bSJed Brown PetscErrorCode FormJacobian2(SNES snes,Vec x,Mat jac,Mat B,void *dummy) 195c4762a1bSJed Brown { 196c4762a1bSJed Brown const PetscScalar *xx; 197c4762a1bSJed Brown PetscScalar A[4]; 198c4762a1bSJed Brown PetscErrorCode ierr; 199c4762a1bSJed Brown PetscInt idx[2] = {0,1}; 200c4762a1bSJed Brown 201c4762a1bSJed Brown /* 202c4762a1bSJed Brown Get pointer to vector data 203c4762a1bSJed Brown */ 204c4762a1bSJed Brown ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); 205c4762a1bSJed Brown 206c4762a1bSJed Brown /* 207c4762a1bSJed Brown Compute Jacobian entries and insert into matrix. 208c4762a1bSJed Brown - Since this is such a small problem, we set all entries for 209c4762a1bSJed Brown the matrix at once. 210c4762a1bSJed Brown */ 211c4762a1bSJed Brown A[0] = 3.0*PetscCosScalar(3.0*xx[0]) + 1.0; A[1] = 0.0; 212c4762a1bSJed Brown A[2] = 0.0; A[3] = 1.0; 213c4762a1bSJed Brown ierr = MatSetValues(B,2,idx,2,idx,A,INSERT_VALUES);CHKERRQ(ierr); 214c4762a1bSJed Brown 215c4762a1bSJed Brown /* 216c4762a1bSJed Brown Restore vector 217c4762a1bSJed Brown */ 218c4762a1bSJed Brown ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr); 219c4762a1bSJed Brown 220c4762a1bSJed Brown /* 221c4762a1bSJed Brown Assemble matrix 222c4762a1bSJed Brown */ 223c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 224c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 225c4762a1bSJed Brown if (jac != B) { 226c4762a1bSJed Brown ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 227c4762a1bSJed Brown ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 228c4762a1bSJed Brown } 229c4762a1bSJed Brown return 0; 230c4762a1bSJed Brown } 231c4762a1bSJed Brown 232c4762a1bSJed Brown 233c4762a1bSJed Brown 234c4762a1bSJed Brown 235c4762a1bSJed Brown /*TEST 236c4762a1bSJed Brown 237c4762a1bSJed Brown build: 238c4762a1bSJed Brown requires: c99 infinity 239c4762a1bSJed Brown 240c4762a1bSJed Brown test: 241c4762a1bSJed Brown args: -snes_converged_reason -snes_linesearch_monitor -snes_linesearch_type l2 242c4762a1bSJed Brown filter: grep Inf 243c4762a1bSJed Brown 244c4762a1bSJed Brown TEST*/ 245