1 2 #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/ 3 #include <petscdm.h> /*I "petscdm.h" I*/ 4 5 /* 6 MatFDColoringSetFunction() takes a function with four arguments, we want to use SNESComputeFunction() 7 since it logs function computation information. 8 */ 9 static PetscErrorCode SNESComputeFunctionCtx(SNES snes,Vec x,Vec f,void *ctx) 10 { 11 return SNESComputeFunction(snes,x,f); 12 } 13 14 /*@C 15 SNESComputeJacobianDefaultColor - Computes the Jacobian using 16 finite differences and coloring to exploit matrix sparsity. 17 18 Collective on SNES 19 20 Input Parameters: 21 + snes - nonlinear solver object 22 . x1 - location at which to evaluate Jacobian 23 - ctx - MatFDColoring context or NULL 24 25 Output Parameters: 26 + J - Jacobian matrix (not altered in this routine) 27 - B - newly computed Jacobian matrix to use with preconditioner (generally the same as J) 28 29 Level: intermediate 30 31 Options Database Key: 32 + -snes_fd_color_use_mat - use a matrix coloring from the explicit matrix nonzero pattern instead of from the DM providing the matrix 33 . -snes_fd_color - Activates SNESComputeJacobianDefaultColor() in SNESSetFromOptions() 34 . -mat_fd_coloring_err <err> - Sets <err> (square root of relative error in the function) 35 . -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude 36 - -mat_fd_type - Either wp or ds (see MATMFFD_WP or MATMFFD_DS) 37 38 Notes: If the coloring is not provided through the context, this will first try to get the 39 coloring from the DM. If the DM type has no coloring routine, then it will try to 40 get the coloring from the matrix. This requires that the matrix have nonzero entries 41 precomputed. 42 43 .keywords: SNES, finite differences, Jacobian, coloring, sparse 44 45 .seealso: SNESSetJacobian(), SNESTestJacobian(), SNESComputeJacobianDefault() 46 MatFDColoringCreate(), MatFDColoringSetFunction() 47 48 @*/ 49 50 PetscErrorCode SNESComputeJacobianDefaultColor(SNES snes,Vec x1,Mat J,Mat B,void *ctx) 51 { 52 MatFDColoring color = (MatFDColoring)ctx; 53 PetscErrorCode ierr; 54 DM dm; 55 MatColoring mc; 56 ISColoring iscoloring; 57 PetscBool hascolor; 58 PetscBool solvec,matcolor = PETSC_FALSE; 59 60 PetscFunctionBegin; 61 if (color) PetscValidHeaderSpecific(color,MAT_FDCOLORING_CLASSID,6); 62 else {ierr = PetscObjectQuery((PetscObject)B,"SNESMatFDColoring",(PetscObject*)&color);CHKERRQ(ierr);} 63 64 if (!color) { 65 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 66 ierr = DMHasColoring(dm,&hascolor);CHKERRQ(ierr); 67 matcolor = PETSC_FALSE; 68 ierr = PetscOptionsGetBool(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_fd_color_use_mat",&matcolor,NULL);CHKERRQ(ierr); 69 if (hascolor && !matcolor) { 70 ierr = DMCreateColoring(dm,IS_COLORING_GLOBAL,&iscoloring);CHKERRQ(ierr); 71 ierr = MatFDColoringCreate(B,iscoloring,&color);CHKERRQ(ierr); 72 ierr = MatFDColoringSetFunction(color,(PetscErrorCode (*)(void))SNESComputeFunctionCtx,NULL);CHKERRQ(ierr); 73 ierr = MatFDColoringSetFromOptions(color);CHKERRQ(ierr); 74 ierr = MatFDColoringSetUp(B,iscoloring,color);CHKERRQ(ierr); 75 ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr); 76 } else { 77 ierr = MatColoringCreate(B,&mc);CHKERRQ(ierr); 78 ierr = MatColoringSetDistance(mc,2);CHKERRQ(ierr); 79 ierr = MatColoringSetType(mc,MATCOLORINGSL);CHKERRQ(ierr); 80 ierr = MatColoringSetFromOptions(mc);CHKERRQ(ierr); 81 ierr = MatColoringApply(mc,&iscoloring);CHKERRQ(ierr); 82 ierr = MatColoringDestroy(&mc);CHKERRQ(ierr); 83 ierr = MatFDColoringCreate(B,iscoloring,&color);CHKERRQ(ierr); 84 ierr = MatFDColoringSetFunction(color,(PetscErrorCode (*)(void))SNESComputeFunctionCtx,NULL);CHKERRQ(ierr); 85 ierr = MatFDColoringSetFromOptions(color);CHKERRQ(ierr); 86 ierr = MatFDColoringSetUp(B,iscoloring,color);CHKERRQ(ierr); 87 ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr); 88 } 89 ierr = PetscObjectCompose((PetscObject)B,"SNESMatFDColoring",(PetscObject)color);CHKERRQ(ierr); 90 ierr = PetscObjectDereference((PetscObject)color);CHKERRQ(ierr); 91 } 92 93 /* F is only usable if there is no RHS on the SNES and the full solution corresponds to x1 */ 94 ierr = VecEqual(x1,snes->vec_sol,&solvec);CHKERRQ(ierr); 95 if (!snes->vec_rhs && solvec) { 96 Vec F; 97 ierr = SNESGetFunction(snes,&F,NULL,NULL);CHKERRQ(ierr); 98 ierr = MatFDColoringSetF(color,F);CHKERRQ(ierr); 99 } 100 ierr = MatFDColoringApply(B,color,x1,snes);CHKERRQ(ierr); 101 if (J != B) { 102 ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 103 ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 104 } 105 PetscFunctionReturn(0); 106 } 107