1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Solves u`` + u^{2} = f with Newton-like methods. Using\n\ 3c4762a1bSJed Brown matrix-free techniques with user-provided explicit preconditioner matrix.\n\n"; 4c4762a1bSJed Brown 5c4762a1bSJed Brown #include <petscsnes.h> 6c4762a1bSJed Brown 7c4762a1bSJed Brown extern PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*); 8c4762a1bSJed Brown extern PetscErrorCode FormJacobianNoMatrix(SNES,Vec,Mat,Mat,void*); 9c4762a1bSJed Brown extern PetscErrorCode FormFunction(SNES,Vec,Vec,void*); 10c4762a1bSJed Brown extern PetscErrorCode FormFunctioni(void *,PetscInt,Vec,PetscScalar *); 11c4762a1bSJed Brown extern PetscErrorCode OtherFunctionForDifferencing(void*,Vec,Vec); 12c4762a1bSJed Brown extern PetscErrorCode FormInitialGuess(SNES,Vec); 13c4762a1bSJed Brown extern PetscErrorCode Monitor(SNES,PetscInt,PetscReal,void*); 14c4762a1bSJed Brown 15c4762a1bSJed Brown typedef struct { 16c4762a1bSJed Brown PetscViewer viewer; 17c4762a1bSJed Brown } MonitorCtx; 18c4762a1bSJed Brown 19c4762a1bSJed Brown typedef struct { 20c4762a1bSJed Brown PetscBool variant; 21c4762a1bSJed Brown } AppCtx; 22c4762a1bSJed Brown 23c4762a1bSJed Brown int main(int argc,char **argv) 24c4762a1bSJed Brown { 25c4762a1bSJed Brown SNES snes; /* SNES context */ 26c4762a1bSJed Brown SNESType type = SNESNEWTONLS; /* default nonlinear solution method */ 27c4762a1bSJed Brown Vec x,r,F,U; /* vectors */ 28c4762a1bSJed Brown Mat J,B; /* Jacobian matrix-free, explicit preconditioner */ 29c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 30c4762a1bSJed Brown PetscScalar h,xp = 0.0,v; 31c4762a1bSJed Brown PetscInt its,n = 5,i; 32c4762a1bSJed Brown PetscErrorCode ierr; 33c4762a1bSJed Brown PetscBool puremf = PETSC_FALSE; 34c4762a1bSJed Brown 35c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 36*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 37*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-variant",&user.variant)); 38c4762a1bSJed Brown h = 1.0/(n-1); 39c4762a1bSJed Brown 40c4762a1bSJed Brown /* Set up data structures */ 41*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&x)); 42*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)x,"Approximate Solution")); 43*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&r)); 44*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&F)); 45*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&U)); 46*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject)U,"Exact Solution")); 47c4762a1bSJed Brown 48a5b23f4aSJose E. Roman /* create explicit matrix preconditioner */ 49*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_SELF,n,n,3,NULL,&B)); 50c4762a1bSJed Brown 51c4762a1bSJed Brown /* Store right-hand-side of PDE and exact solution */ 52c4762a1bSJed Brown for (i=0; i<n; i++) { 53c4762a1bSJed Brown v = 6.0*xp + PetscPowScalar(xp+1.e-12,6.0); /* +1.e-12 is to prevent 0^6 */ 54*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValues(F,1,&i,&v,INSERT_VALUES)); 55c4762a1bSJed Brown v = xp*xp*xp; 56*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValues(U,1,&i,&v,INSERT_VALUES)); 57c4762a1bSJed Brown xp += h; 58c4762a1bSJed Brown } 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* Create nonlinear solver */ 61*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESCreate(PETSC_COMM_WORLD,&snes)); 62*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetType(snes,type)); 63c4762a1bSJed Brown 64c4762a1bSJed Brown /* Set various routines and options */ 65*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFunction(snes,r,FormFunction,F)); 66c4762a1bSJed Brown if (user.variant) { 67c4762a1bSJed Brown /* this approach is not normally needed, one should use the MatCreateSNESMF() below usually */ 68*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateMFFD(PETSC_COMM_WORLD,n,n,n,n,&J)); 69*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMFFDSetFunction(J,(PetscErrorCode (*)(void*, Vec, Vec))SNESComputeFunction,snes)); 70*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMFFDSetFunctioni(J,FormFunctioni)); 71c4762a1bSJed Brown /* Use the matrix free operator for both the Jacobian used to define the linear system and used to define the preconditioner */ 72c4762a1bSJed Brown /* This tests MatGetDiagonal() for MATMFFD */ 73*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-puremf",&puremf)); 74c4762a1bSJed Brown } else { 75c4762a1bSJed Brown /* create matrix free matrix for Jacobian */ 76*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSNESMF(snes,&J)); 77c4762a1bSJed Brown /* demonstrates differencing a different function than FormFunction() to apply a matrix operator */ 78c4762a1bSJed Brown /* note we use the same context for this function as FormFunction, the F vector */ 79*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMFFDSetFunction(J,OtherFunctionForDifferencing,F)); 80c4762a1bSJed Brown } 81c4762a1bSJed Brown 82c4762a1bSJed Brown /* Set various routines and options */ 83*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetJacobian(snes,J,puremf ? J : B,puremf ? FormJacobianNoMatrix : FormJacobian,&user)); 84*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFromOptions(snes)); 85c4762a1bSJed Brown 86c4762a1bSJed Brown /* Solve nonlinear system */ 87*5f80ce2aSJacob Faibussowitsch CHKERRQ(FormInitialGuess(snes,x)); 88*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSolve(snes,NULL,x)); 89*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetIterationNumber(snes,&its)); 90*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"number of SNES iterations = %D\n\n",its)); 91c4762a1bSJed Brown 92c4762a1bSJed Brown /* Free data structures */ 93*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); CHKERRQ(VecDestroy(&r)); 94*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&U)); CHKERRQ(VecDestroy(&F)); 95*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&J)); CHKERRQ(MatDestroy(&B)); 96*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESDestroy(&snes)); 97c4762a1bSJed Brown ierr = PetscFinalize(); 98c4762a1bSJed Brown return ierr; 99c4762a1bSJed Brown } 100c4762a1bSJed Brown /* -------------------- Evaluate Function F(x) --------------------- */ 101c4762a1bSJed Brown 102c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec x,Vec f,void *dummy) 103c4762a1bSJed Brown { 104c4762a1bSJed Brown const PetscScalar *xx,*FF; 105c4762a1bSJed Brown PetscScalar *ff,d; 106c4762a1bSJed Brown PetscInt i,n; 107c4762a1bSJed Brown 108*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(x,&xx)); 109*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(f,&ff)); 110*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead((Vec) dummy,&FF)); 111*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(x,&n)); 112c4762a1bSJed Brown d = (PetscReal)(n - 1); d = d*d; 113c4762a1bSJed Brown ff[0] = xx[0]; 114c4762a1bSJed Brown for (i=1; i<n-1; i++) ff[i] = d*(xx[i-1] - 2.0*xx[i] + xx[i+1]) + xx[i]*xx[i] - FF[i]; 115c4762a1bSJed Brown ff[n-1] = xx[n-1] - 1.0; 116*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(x,&xx)); 117*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(f,&ff)); 118*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead((Vec)dummy,&FF)); 119c4762a1bSJed Brown return 0; 120c4762a1bSJed Brown } 121c4762a1bSJed Brown 122c4762a1bSJed Brown PetscErrorCode FormFunctioni(void *dummy,PetscInt i,Vec x,PetscScalar *s) 123c4762a1bSJed Brown { 124c4762a1bSJed Brown const PetscScalar *xx,*FF; 125c4762a1bSJed Brown PetscScalar d; 126c4762a1bSJed Brown PetscInt n; 127c4762a1bSJed Brown SNES snes = (SNES) dummy; 128c4762a1bSJed Brown Vec F; 129c4762a1bSJed Brown 130*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetFunction(snes,NULL,NULL,(void**)&F)); 131*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(x,&xx)); 132*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(F,&FF)); 133*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(x,&n)); 134c4762a1bSJed Brown d = (PetscReal)(n - 1); d = d*d; 135c4762a1bSJed Brown if (i == 0) { 136c4762a1bSJed Brown *s = xx[0]; 137c4762a1bSJed Brown } else if (i == n-1) { 138c4762a1bSJed Brown *s = xx[n-1] - 1.0; 139c4762a1bSJed Brown } else { 140c4762a1bSJed Brown *s = d*(xx[i-1] - 2.0*xx[i] + xx[i+1]) + xx[i]*xx[i] - FF[i]; 141c4762a1bSJed Brown } 142*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(x,&xx)); 143*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(F,&FF)); 144c4762a1bSJed Brown return 0; 145c4762a1bSJed Brown } 146c4762a1bSJed Brown 147c4762a1bSJed Brown /* 148c4762a1bSJed Brown 149c4762a1bSJed Brown Example function that when differenced produces the same matrix free Jacobian as FormFunction() 150c4762a1bSJed Brown this is provided to show how a user can provide a different function 151c4762a1bSJed Brown */ 152c4762a1bSJed Brown PetscErrorCode OtherFunctionForDifferencing(void *dummy,Vec x,Vec f) 153c4762a1bSJed Brown { 154c4762a1bSJed Brown 155*5f80ce2aSJacob Faibussowitsch CHKERRQ(FormFunction(NULL,x,f,dummy)); 156*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecShift(f,1.0)); 157c4762a1bSJed Brown return 0; 158c4762a1bSJed Brown } 159c4762a1bSJed Brown 160c4762a1bSJed Brown /* -------------------- Form initial approximation ----------------- */ 161c4762a1bSJed Brown 162c4762a1bSJed Brown PetscErrorCode FormInitialGuess(SNES snes,Vec x) 163c4762a1bSJed Brown { 164c4762a1bSJed Brown PetscScalar pfive = .50; 165*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(x,pfive)); 166c4762a1bSJed Brown return 0; 167c4762a1bSJed Brown } 168c4762a1bSJed Brown /* -------------------- Evaluate Jacobian F'(x) -------------------- */ 169c4762a1bSJed Brown /* Evaluates a matrix that is used to precondition the matrix-free 170a5b23f4aSJose E. Roman jacobian. In this case, the explicit preconditioner matrix is 171c4762a1bSJed Brown also EXACTLY the Jacobian. In general, it would be some lower 172c4762a1bSJed Brown order, simplified apprioximation */ 173c4762a1bSJed Brown 174c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES snes,Vec x,Mat jac,Mat B,void *dummy) 175c4762a1bSJed Brown { 176c4762a1bSJed Brown const PetscScalar *xx; 177c4762a1bSJed Brown PetscScalar A[3],d; 178c4762a1bSJed Brown PetscInt i,n,j[3]; 179c4762a1bSJed Brown AppCtx *user = (AppCtx*) dummy; 180c4762a1bSJed Brown 181*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(x,&xx)); 182*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(x,&n)); 183c4762a1bSJed Brown d = (PetscReal)(n - 1); d = d*d; 184c4762a1bSJed Brown 185c4762a1bSJed Brown i = 0; A[0] = 1.0; 186*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(B,1,&i,1,&i,&A[0],INSERT_VALUES)); 187c4762a1bSJed Brown for (i=1; i<n-1; i++) { 188c4762a1bSJed Brown j[0] = i - 1; j[1] = i; j[2] = i + 1; 189c4762a1bSJed Brown A[0] = d; A[1] = -2.0*d + 2.0*xx[i]; A[2] = d; 190*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(B,1,&i,3,j,A,INSERT_VALUES)); 191c4762a1bSJed Brown } 192c4762a1bSJed Brown i = n-1; A[0] = 1.0; 193*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(B,1,&i,1,&i,&A[0],INSERT_VALUES)); 194*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 195*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 196*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(x,&xx)); 197c4762a1bSJed Brown 198c4762a1bSJed Brown if (user->variant) { 199*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMFFDSetBase(jac,x,NULL)); 200c4762a1bSJed Brown } 201*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); 202*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); 203c4762a1bSJed Brown return 0; 204c4762a1bSJed Brown } 205c4762a1bSJed Brown 206c4762a1bSJed Brown PetscErrorCode FormJacobianNoMatrix(SNES snes,Vec x,Mat jac,Mat B,void *dummy) 207c4762a1bSJed Brown { 208c4762a1bSJed Brown AppCtx *user = (AppCtx*) dummy; 209c4762a1bSJed Brown 210c4762a1bSJed Brown if (user->variant) { 211*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMFFDSetBase(jac,x,NULL)); 212c4762a1bSJed Brown } 213*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); 214*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); 215c4762a1bSJed Brown return 0; 216c4762a1bSJed Brown } 217c4762a1bSJed Brown 218c4762a1bSJed Brown /* -------------------- User-defined monitor ----------------------- */ 219c4762a1bSJed Brown 220c4762a1bSJed Brown PetscErrorCode Monitor(SNES snes,PetscInt its,PetscReal fnorm,void *dummy) 221c4762a1bSJed Brown { 222c4762a1bSJed Brown MonitorCtx *monP = (MonitorCtx*) dummy; 223c4762a1bSJed Brown Vec x; 224c4762a1bSJed Brown MPI_Comm comm; 225c4762a1bSJed Brown 226*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject)snes,&comm)); 227*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFPrintf(comm,stdout,"iter = %D, SNES Function norm %g \n",its,(double)fnorm)); 228*5f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetSolution(snes,&x)); 229*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(x,monP->viewer)); 230c4762a1bSJed Brown return 0; 231c4762a1bSJed Brown } 232c4762a1bSJed Brown 233c4762a1bSJed Brown /*TEST 234c4762a1bSJed Brown 235c4762a1bSJed Brown test: 236c4762a1bSJed Brown args: -ksp_gmres_cgs_refinement_type refine_always -snes_monitor_short 237c4762a1bSJed Brown 238c4762a1bSJed Brown test: 239c4762a1bSJed Brown suffix: 2 240c4762a1bSJed Brown args: -variant -ksp_gmres_cgs_refinement_type refine_always -snes_monitor_short 241c4762a1bSJed Brown output_file: output/ex7_1.out 242c4762a1bSJed Brown 243c4762a1bSJed Brown # uses AIJ matrix to define diagonal matrix for Jacobian preconditioning 244c4762a1bSJed Brown test: 245c4762a1bSJed Brown suffix: 3 246c4762a1bSJed Brown args: -variant -pc_type jacobi -snes_view -ksp_monitor 247c4762a1bSJed Brown 248c4762a1bSJed Brown # uses MATMFFD matrix to define diagonal matrix for Jacobian preconditioning 249c4762a1bSJed Brown test: 250c4762a1bSJed Brown suffix: 4 251c4762a1bSJed Brown args: -variant -pc_type jacobi -puremf -snes_view -ksp_monitor 252c4762a1bSJed Brown 253c4762a1bSJed Brown TEST*/ 254