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