191627157SBarry Smith 2b45d2f2cSJed Brown #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/ 31e25c274SJed Brown #include <petscdm.h> /*I "petscdm.h" I*/ 491627157SBarry Smith 54a2ae208SSatish Balay #undef __FUNCT__ 68d359177SBarry Smith #define __FUNCT__ "SNESComputeJacobianDefaultColor" 791627157SBarry Smith /*@C 88d359177SBarry Smith SNESComputeJacobianDefaultColor - Computes the Jacobian using 9b4fc646aSLois Curfman McInnes finite differences and coloring to exploit matrix sparsity. 1091627157SBarry Smith 11fee21e36SBarry Smith Collective on SNES 12fee21e36SBarry Smith 13c7afd0dbSLois Curfman McInnes Input Parameters: 14c7afd0dbSLois Curfman McInnes + snes - nonlinear solver object 15c7afd0dbSLois Curfman McInnes . x1 - location at which to evaluate Jacobian 167fb41025SPeter Brune - ctx - MatFDColoring context or NULL 17c7afd0dbSLois Curfman McInnes 18c7afd0dbSLois Curfman McInnes Output Parameters: 19c7afd0dbSLois Curfman McInnes + J - Jacobian matrix (not altered in this routine) 2056744491SBarry Smith - B - newly computed Jacobian matrix to use with preconditioner (generally the same as J) 21c7afd0dbSLois Curfman McInnes 2236851e7fSLois Curfman McInnes Level: intermediate 2336851e7fSLois Curfman McInnes 24*5620d6dcSBarry Smith Options Database Key: 25*5620d6dcSBarry Smith + -snes_fd_color_use_mat - use a matrix coloring from the explicit matrix nonzero pattern instead of from the DM providing the matrix 26*5620d6dcSBarry Smith . -snes_fd_color - Activates SNESComputeJacobianDefaultColor() in SNESSetFromOptions() 27*5620d6dcSBarry Smith . -mat_fd_coloring_err <err> - Sets <err> (square root of relative error in the function) 28*5620d6dcSBarry Smith . -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude 29*5620d6dcSBarry Smith - -mat_fd_type - Either wp or ds (see MATMFFD_WP or MATMFFD_DS) 30*5620d6dcSBarry Smith 31*5620d6dcSBarry Smith Notes: If the coloring is not provided through the context, this will first try to get the 327fb41025SPeter Brune coloring from the DM. If the DM type has no coloring routine, then it will try to 337fb41025SPeter Brune get the coloring from the matrix. This requires that the matrix have nonzero entries 347086a01eSPeter Brune precomputed. This is discouraged, as MatColoringApply() is not parallel by default. 35b0ae01b7SPeter Brune 36b4fc646aSLois Curfman McInnes .keywords: SNES, finite differences, Jacobian, coloring, sparse 3791627157SBarry Smith 388d359177SBarry Smith .seealso: SNESSetJacobian(), SNESTestJacobian(), SNESComputeJacobianDefault() 39ab637aeaSJed Brown MatFDColoringCreate(), MatFDColoringSetFunction() 40cb5b572fSBarry Smith 4191627157SBarry Smith @*/ 4295d750ceSBarry Smith 43d1e9a80fSBarry Smith PetscErrorCode SNESComputeJacobianDefaultColor(SNES snes,Vec x1,Mat J,Mat B,void *ctx) 4491627157SBarry Smith { 457fb41025SPeter Brune MatFDColoring color = (MatFDColoring)ctx; 46dfbe8321SBarry Smith PetscErrorCode ierr; 474e269d77SPeter Brune DM dm; 484e269d77SPeter Brune PetscErrorCode (*func)(SNES,Vec,Vec,void*); 494e269d77SPeter Brune Vec F; 504e269d77SPeter Brune void *funcctx; 51335efc43SPeter Brune MatColoring mc; 524e269d77SPeter Brune ISColoring iscoloring; 53b0ae01b7SPeter Brune PetscBool hascolor; 54aa2f1b4eSJed Brown PetscBool solvec,matcolor = PETSC_FALSE; 55dff777c9SBarry Smith 563a40ed3dSBarry Smith PetscFunctionBegin; 577fb41025SPeter Brune if (color) PetscValidHeaderSpecific(color,MAT_FDCOLORING_CLASSID,6); 5894ab13aaSBarry Smith else {ierr = PetscObjectQuery((PetscObject)B,"SNESMatFDColoring",(PetscObject*)&color);CHKERRQ(ierr);} 5922d28d08SBarry Smith ierr = SNESGetFunction(snes,&F,&func,&funcctx);CHKERRQ(ierr); 604e269d77SPeter Brune if (!color) { 614e269d77SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 62b0ae01b7SPeter Brune ierr = DMHasColoring(dm,&hascolor);CHKERRQ(ierr); 63924833adSPeter Brune matcolor = PETSC_FALSE; 64c3138ed3SPeter Brune ierr = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_fd_color_use_mat",&matcolor,NULL);CHKERRQ(ierr); 65c3138ed3SPeter Brune if (hascolor && !matcolor) { 66b412c318SBarry Smith ierr = DMCreateColoring(dm,IS_COLORING_GLOBAL,&iscoloring);CHKERRQ(ierr); 6794ab13aaSBarry Smith ierr = MatFDColoringCreate(B,iscoloring,&color);CHKERRQ(ierr); 684e269d77SPeter Brune ierr = MatFDColoringSetFunction(color,(PetscErrorCode (*)(void))func,funcctx);CHKERRQ(ierr); 694e269d77SPeter Brune ierr = MatFDColoringSetFromOptions(color);CHKERRQ(ierr); 7094ab13aaSBarry Smith ierr = MatFDColoringSetUp(B,iscoloring,color);CHKERRQ(ierr); 71f86b9fbaSHong Zhang ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr); 72b0ae01b7SPeter Brune } else { 7394ab13aaSBarry Smith ierr = MatColoringCreate(B,&mc);CHKERRQ(ierr); 74335efc43SPeter Brune ierr = MatColoringSetDistance(mc,2);CHKERRQ(ierr); 75335efc43SPeter Brune ierr = MatColoringSetType(mc,MATCOLORINGSL);CHKERRQ(ierr); 76335efc43SPeter Brune ierr = MatColoringSetFromOptions(mc);CHKERRQ(ierr); 77335efc43SPeter Brune ierr = MatColoringApply(mc,&iscoloring);CHKERRQ(ierr); 78335efc43SPeter Brune ierr = MatColoringDestroy(&mc);CHKERRQ(ierr); 7994ab13aaSBarry Smith ierr = MatFDColoringCreate(B,iscoloring,&color);CHKERRQ(ierr); 80b0ae01b7SPeter Brune ierr = MatFDColoringSetFunction(color,(PetscErrorCode (*)(void))func,(void*)funcctx);CHKERRQ(ierr); 81b0ae01b7SPeter Brune ierr = MatFDColoringSetFromOptions(color);CHKERRQ(ierr); 8294ab13aaSBarry Smith ierr = MatFDColoringSetUp(B,iscoloring,color);CHKERRQ(ierr); 83f86b9fbaSHong Zhang ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr); 84b0ae01b7SPeter Brune } 8594ab13aaSBarry Smith ierr = PetscObjectCompose((PetscObject)B,"SNESMatFDColoring",(PetscObject)color);CHKERRQ(ierr); 864108615fSPeter Brune ierr = PetscObjectDereference((PetscObject)color);CHKERRQ(ierr); 874a9d489dSBarry Smith } 8895862db2SPeter Brune 89c89319d2SPeter Brune /* F is only usable if there is no RHS on the SNES and the full solution corresponds to x1 */ 90c89319d2SPeter Brune ierr = VecEqual(x1,snes->vec_sol,&solvec);CHKERRQ(ierr); 91c89319d2SPeter Brune if (!snes->vec_rhs && solvec) { 92432bb422SPeter Brune ierr = MatFDColoringSetF(color,F);CHKERRQ(ierr); 9395862db2SPeter Brune } 94d1e9a80fSBarry Smith ierr = MatFDColoringApply(B,color,x1,snes);CHKERRQ(ierr); 9594ab13aaSBarry Smith if (J != B) { 9694ab13aaSBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 9794ab13aaSBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 98194405b3SBarry Smith } 993a40ed3dSBarry Smith PetscFunctionReturn(0); 10091627157SBarry Smith } 101