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