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 PetscBool puremf = PETSC_FALSE; 33c4762a1bSJed Brown 349566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 359566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 369566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-variant",&user.variant)); 37c4762a1bSJed Brown h = 1.0/(n-1); 38c4762a1bSJed Brown 39c4762a1bSJed Brown /* Set up data structures */ 409566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&x)); 419566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)x,"Approximate Solution")); 429566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&r)); 439566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&F)); 449566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&U)); 459566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)U,"Exact Solution")); 46c4762a1bSJed Brown 47a5b23f4aSJose E. Roman /* create explicit matrix preconditioner */ 489566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF,n,n,3,NULL,&B)); 49c4762a1bSJed Brown 50c4762a1bSJed Brown /* Store right-hand-side of PDE and exact solution */ 51c4762a1bSJed Brown for (i=0; i<n; i++) { 52c4762a1bSJed Brown v = 6.0*xp + PetscPowScalar(xp+1.e-12,6.0); /* +1.e-12 is to prevent 0^6 */ 539566063dSJacob Faibussowitsch PetscCall(VecSetValues(F,1,&i,&v,INSERT_VALUES)); 54c4762a1bSJed Brown v = xp*xp*xp; 559566063dSJacob Faibussowitsch PetscCall(VecSetValues(U,1,&i,&v,INSERT_VALUES)); 56c4762a1bSJed Brown xp += h; 57c4762a1bSJed Brown } 58c4762a1bSJed Brown 59c4762a1bSJed Brown /* Create nonlinear solver */ 609566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD,&snes)); 619566063dSJacob Faibussowitsch PetscCall(SNESSetType(snes,type)); 62c4762a1bSJed Brown 63c4762a1bSJed Brown /* Set various routines and options */ 649566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes,r,FormFunction,F)); 65c4762a1bSJed Brown if (user.variant) { 66c4762a1bSJed Brown /* this approach is not normally needed, one should use the MatCreateSNESMF() below usually */ 679566063dSJacob Faibussowitsch PetscCall(MatCreateMFFD(PETSC_COMM_WORLD,n,n,n,n,&J)); 689566063dSJacob Faibussowitsch PetscCall(MatMFFDSetFunction(J,(PetscErrorCode (*)(void*, Vec, Vec))SNESComputeFunction,snes)); 699566063dSJacob Faibussowitsch PetscCall(MatMFFDSetFunctioni(J,FormFunctioni)); 70c4762a1bSJed Brown /* Use the matrix free operator for both the Jacobian used to define the linear system and used to define the preconditioner */ 71c4762a1bSJed Brown /* This tests MatGetDiagonal() for MATMFFD */ 729566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-puremf",&puremf)); 73c4762a1bSJed Brown } else { 74c4762a1bSJed Brown /* create matrix free matrix for Jacobian */ 759566063dSJacob Faibussowitsch PetscCall(MatCreateSNESMF(snes,&J)); 76c4762a1bSJed Brown /* demonstrates differencing a different function than FormFunction() to apply a matrix operator */ 77c4762a1bSJed Brown /* note we use the same context for this function as FormFunction, the F vector */ 789566063dSJacob Faibussowitsch PetscCall(MatMFFDSetFunction(J,OtherFunctionForDifferencing,F)); 79c4762a1bSJed Brown } 80c4762a1bSJed Brown 81c4762a1bSJed Brown /* Set various routines and options */ 829566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes,J,puremf ? J : B,puremf ? FormJacobianNoMatrix : FormJacobian,&user)); 839566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 84c4762a1bSJed Brown 85c4762a1bSJed Brown /* Solve nonlinear system */ 869566063dSJacob Faibussowitsch PetscCall(FormInitialGuess(snes,x)); 879566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes,NULL,x)); 889566063dSJacob Faibussowitsch PetscCall(SNESGetIterationNumber(snes,&its)); 89*63a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"number of SNES iterations = %" PetscInt_FMT "\n\n",its)); 90c4762a1bSJed Brown 91c4762a1bSJed Brown /* Free data structures */ 929566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); PetscCall(VecDestroy(&r)); 939566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); PetscCall(VecDestroy(&F)); 949566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); PetscCall(MatDestroy(&B)); 959566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 969566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 97b122ec5aSJacob Faibussowitsch return 0; 98c4762a1bSJed Brown } 99c4762a1bSJed Brown /* -------------------- Evaluate Function F(x) --------------------- */ 100c4762a1bSJed Brown 101c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec x,Vec f,void *dummy) 102c4762a1bSJed Brown { 103c4762a1bSJed Brown const PetscScalar *xx,*FF; 104c4762a1bSJed Brown PetscScalar *ff,d; 105c4762a1bSJed Brown PetscInt i,n; 106c4762a1bSJed Brown 1079566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x,&xx)); 1089566063dSJacob Faibussowitsch PetscCall(VecGetArray(f,&ff)); 1099566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead((Vec) dummy,&FF)); 1109566063dSJacob Faibussowitsch PetscCall(VecGetSize(x,&n)); 111c4762a1bSJed Brown d = (PetscReal)(n - 1); d = d*d; 112c4762a1bSJed Brown ff[0] = xx[0]; 113c4762a1bSJed 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]; 114c4762a1bSJed Brown ff[n-1] = xx[n-1] - 1.0; 1159566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x,&xx)); 1169566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(f,&ff)); 1179566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead((Vec)dummy,&FF)); 118c4762a1bSJed Brown return 0; 119c4762a1bSJed Brown } 120c4762a1bSJed Brown 121c4762a1bSJed Brown PetscErrorCode FormFunctioni(void *dummy,PetscInt i,Vec x,PetscScalar *s) 122c4762a1bSJed Brown { 123c4762a1bSJed Brown const PetscScalar *xx,*FF; 124c4762a1bSJed Brown PetscScalar d; 125c4762a1bSJed Brown PetscInt n; 126c4762a1bSJed Brown SNES snes = (SNES) dummy; 127c4762a1bSJed Brown Vec F; 128c4762a1bSJed Brown 1299566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(snes,NULL,NULL,(void**)&F)); 1309566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x,&xx)); 1319566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(F,&FF)); 1329566063dSJacob Faibussowitsch PetscCall(VecGetSize(x,&n)); 133c4762a1bSJed Brown d = (PetscReal)(n - 1); d = d*d; 134c4762a1bSJed Brown if (i == 0) { 135c4762a1bSJed Brown *s = xx[0]; 136c4762a1bSJed Brown } else if (i == n-1) { 137c4762a1bSJed Brown *s = xx[n-1] - 1.0; 138c4762a1bSJed Brown } else { 139c4762a1bSJed Brown *s = d*(xx[i-1] - 2.0*xx[i] + xx[i+1]) + xx[i]*xx[i] - FF[i]; 140c4762a1bSJed Brown } 1419566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x,&xx)); 1429566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(F,&FF)); 143c4762a1bSJed Brown return 0; 144c4762a1bSJed Brown } 145c4762a1bSJed Brown 146c4762a1bSJed Brown /* 147c4762a1bSJed Brown 148c4762a1bSJed Brown Example function that when differenced produces the same matrix free Jacobian as FormFunction() 149c4762a1bSJed Brown this is provided to show how a user can provide a different function 150c4762a1bSJed Brown */ 151c4762a1bSJed Brown PetscErrorCode OtherFunctionForDifferencing(void *dummy,Vec x,Vec f) 152c4762a1bSJed Brown { 153c4762a1bSJed Brown 1549566063dSJacob Faibussowitsch PetscCall(FormFunction(NULL,x,f,dummy)); 1559566063dSJacob Faibussowitsch PetscCall(VecShift(f,1.0)); 156c4762a1bSJed Brown return 0; 157c4762a1bSJed Brown } 158c4762a1bSJed Brown 159c4762a1bSJed Brown /* -------------------- Form initial approximation ----------------- */ 160c4762a1bSJed Brown 161c4762a1bSJed Brown PetscErrorCode FormInitialGuess(SNES snes,Vec x) 162c4762a1bSJed Brown { 163c4762a1bSJed Brown PetscScalar pfive = .50; 1649566063dSJacob Faibussowitsch PetscCall(VecSet(x,pfive)); 165c4762a1bSJed Brown return 0; 166c4762a1bSJed Brown } 167c4762a1bSJed Brown /* -------------------- Evaluate Jacobian F'(x) -------------------- */ 168c4762a1bSJed Brown /* Evaluates a matrix that is used to precondition the matrix-free 169a5b23f4aSJose E. Roman jacobian. In this case, the explicit preconditioner matrix is 170c4762a1bSJed Brown also EXACTLY the Jacobian. In general, it would be some lower 171c4762a1bSJed Brown order, simplified apprioximation */ 172c4762a1bSJed Brown 173c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES snes,Vec x,Mat jac,Mat B,void *dummy) 174c4762a1bSJed Brown { 175c4762a1bSJed Brown const PetscScalar *xx; 176c4762a1bSJed Brown PetscScalar A[3],d; 177c4762a1bSJed Brown PetscInt i,n,j[3]; 178c4762a1bSJed Brown AppCtx *user = (AppCtx*) dummy; 179c4762a1bSJed Brown 1809566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x,&xx)); 1819566063dSJacob Faibussowitsch PetscCall(VecGetSize(x,&n)); 182c4762a1bSJed Brown d = (PetscReal)(n - 1); d = d*d; 183c4762a1bSJed Brown 184c4762a1bSJed Brown i = 0; A[0] = 1.0; 1859566063dSJacob Faibussowitsch PetscCall(MatSetValues(B,1,&i,1,&i,&A[0],INSERT_VALUES)); 186c4762a1bSJed Brown for (i=1; i<n-1; i++) { 187c4762a1bSJed Brown j[0] = i - 1; j[1] = i; j[2] = i + 1; 188c4762a1bSJed Brown A[0] = d; A[1] = -2.0*d + 2.0*xx[i]; A[2] = d; 1899566063dSJacob Faibussowitsch PetscCall(MatSetValues(B,1,&i,3,j,A,INSERT_VALUES)); 190c4762a1bSJed Brown } 191c4762a1bSJed Brown i = n-1; A[0] = 1.0; 1929566063dSJacob Faibussowitsch PetscCall(MatSetValues(B,1,&i,1,&i,&A[0],INSERT_VALUES)); 1939566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 1949566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 1959566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x,&xx)); 196c4762a1bSJed Brown 197c4762a1bSJed Brown if (user->variant) { 1989566063dSJacob Faibussowitsch PetscCall(MatMFFDSetBase(jac,x,NULL)); 199c4762a1bSJed Brown } 2009566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); 2019566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); 202c4762a1bSJed Brown return 0; 203c4762a1bSJed Brown } 204c4762a1bSJed Brown 205c4762a1bSJed Brown PetscErrorCode FormJacobianNoMatrix(SNES snes,Vec x,Mat jac,Mat B,void *dummy) 206c4762a1bSJed Brown { 207c4762a1bSJed Brown AppCtx *user = (AppCtx*) dummy; 208c4762a1bSJed Brown 209c4762a1bSJed Brown if (user->variant) { 2109566063dSJacob Faibussowitsch PetscCall(MatMFFDSetBase(jac,x,NULL)); 211c4762a1bSJed Brown } 2129566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); 2139566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); 214c4762a1bSJed Brown return 0; 215c4762a1bSJed Brown } 216c4762a1bSJed Brown 217c4762a1bSJed Brown /* -------------------- User-defined monitor ----------------------- */ 218c4762a1bSJed Brown 219c4762a1bSJed Brown PetscErrorCode Monitor(SNES snes,PetscInt its,PetscReal fnorm,void *dummy) 220c4762a1bSJed Brown { 221c4762a1bSJed Brown MonitorCtx *monP = (MonitorCtx*) dummy; 222c4762a1bSJed Brown Vec x; 223c4762a1bSJed Brown MPI_Comm comm; 224c4762a1bSJed Brown 2259566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)snes,&comm)); 226*63a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,stdout,"iter = %" PetscInt_FMT ", SNES Function norm %g \n",its,(double)fnorm)); 2279566063dSJacob Faibussowitsch PetscCall(SNESGetSolution(snes,&x)); 2289566063dSJacob Faibussowitsch PetscCall(VecView(x,monP->viewer)); 229c4762a1bSJed Brown return 0; 230c4762a1bSJed Brown } 231c4762a1bSJed Brown 232c4762a1bSJed Brown /*TEST 233c4762a1bSJed Brown 234c4762a1bSJed Brown test: 235c4762a1bSJed Brown args: -ksp_gmres_cgs_refinement_type refine_always -snes_monitor_short 236c4762a1bSJed Brown 237c4762a1bSJed Brown test: 238c4762a1bSJed Brown suffix: 2 239c4762a1bSJed Brown args: -variant -ksp_gmres_cgs_refinement_type refine_always -snes_monitor_short 240c4762a1bSJed Brown output_file: output/ex7_1.out 241c4762a1bSJed Brown 242c4762a1bSJed Brown # uses AIJ matrix to define diagonal matrix for Jacobian preconditioning 243c4762a1bSJed Brown test: 244c4762a1bSJed Brown suffix: 3 245c4762a1bSJed Brown args: -variant -pc_type jacobi -snes_view -ksp_monitor 246c4762a1bSJed Brown 247c4762a1bSJed Brown # uses MATMFFD matrix to define diagonal matrix for Jacobian preconditioning 248c4762a1bSJed Brown test: 249c4762a1bSJed Brown suffix: 4 250c4762a1bSJed Brown args: -variant -pc_type jacobi -puremf -snes_view -ksp_monitor 251c4762a1bSJed Brown 252c4762a1bSJed Brown TEST*/ 253