xref: /petsc/src/snes/interface/snesj2.c (revision 924833ad85fb09d0f19f1b6dd4f66241da1ddf13)
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)
20c7afd0dbSLois Curfman McInnes .   B - newly computed Jacobian matrix to use with preconditioner (generally the same as J)
21c7afd0dbSLois Curfman McInnes -   flag - flag indicating whether the matrix sparsity structure has changed
22c7afd0dbSLois Curfman McInnes 
2336851e7fSLois Curfman McInnes     Level: intermediate
2436851e7fSLois Curfman McInnes 
257fb41025SPeter Brune .notes: If the coloring is not provided through the context, this will first try to get the
267fb41025SPeter Brune         coloring from the DM.  If the DM type has no coloring routine, then it will try to
277fb41025SPeter Brune         get the coloring from the matrix.  This requires that the matrix have nonzero entries
287086a01eSPeter Brune         precomputed.  This is discouraged, as MatColoringApply() is not parallel by default.
29b0ae01b7SPeter Brune 
30b4fc646aSLois Curfman McInnes .keywords: SNES, finite differences, Jacobian, coloring, sparse
3191627157SBarry Smith 
328d359177SBarry Smith .seealso: SNESSetJacobian(), SNESTestJacobian(), SNESComputeJacobianDefault()
33ab637aeaSJed Brown           MatFDColoringCreate(), MatFDColoringSetFunction()
34cb5b572fSBarry Smith 
3591627157SBarry Smith @*/
3695d750ceSBarry Smith 
378d359177SBarry Smith PetscErrorCode  SNESComputeJacobianDefaultColor(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
3891627157SBarry Smith {
397fb41025SPeter Brune   MatFDColoring  color = (MatFDColoring)ctx;
40dfbe8321SBarry Smith   PetscErrorCode ierr;
414e269d77SPeter Brune   DM             dm;
424e269d77SPeter Brune   PetscErrorCode (*func)(SNES,Vec,Vec,void*);
434e269d77SPeter Brune   Vec            F;
444e269d77SPeter Brune   void           *funcctx;
45335efc43SPeter Brune   MatColoring    mc;
464e269d77SPeter Brune   ISColoring     iscoloring;
47b0ae01b7SPeter Brune   PetscBool      hascolor;
48c3138ed3SPeter Brune   PetscBool      solvec,matcolor;
49dff777c9SBarry Smith 
503a40ed3dSBarry Smith   PetscFunctionBegin;
517fb41025SPeter Brune   if (color) PetscValidHeaderSpecific(color,MAT_FDCOLORING_CLASSID,6);
527fb41025SPeter Brune   else {ierr  = PetscObjectQuery((PetscObject)*B,"SNESMatFDColoring",(PetscObject*)&color);CHKERRQ(ierr);}
5392abec31SBarry Smith   *flag = SAME_NONZERO_PATTERN;
5422d28d08SBarry Smith   ierr  = SNESGetFunction(snes,&F,&func,&funcctx);CHKERRQ(ierr);
554e269d77SPeter Brune   if (!color) {
564e269d77SPeter Brune     ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
57b0ae01b7SPeter Brune     ierr = DMHasColoring(dm,&hascolor);CHKERRQ(ierr);
58*924833adSPeter Brune     matcolor = PETSC_FALSE;
59c3138ed3SPeter Brune     ierr = PetscOptionsGetBool(((PetscObject)snes)->prefix,"-snes_fd_color_use_mat",&matcolor,NULL);CHKERRQ(ierr);
60c3138ed3SPeter Brune     if (hascolor && !matcolor) {
61b412c318SBarry Smith       ierr = DMCreateColoring(dm,IS_COLORING_GLOBAL,&iscoloring);CHKERRQ(ierr);
624e269d77SPeter Brune       ierr = MatFDColoringCreate(*B,iscoloring,&color);CHKERRQ(ierr);
634e269d77SPeter Brune       ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
644e269d77SPeter Brune       ierr = MatFDColoringSetFunction(color,(PetscErrorCode (*)(void))func,funcctx);CHKERRQ(ierr);
654e269d77SPeter Brune       ierr = MatFDColoringSetFromOptions(color);CHKERRQ(ierr);
66b0ae01b7SPeter Brune     } else {
67335efc43SPeter Brune       ierr = MatColoringCreate(*B,&mc);CHKERRQ(ierr);
68335efc43SPeter Brune       ierr = MatColoringSetDistance(mc,2);CHKERRQ(ierr);
69335efc43SPeter Brune       ierr = MatColoringSetType(mc,MATCOLORINGSL);CHKERRQ(ierr);
70335efc43SPeter Brune       ierr = MatColoringSetFromOptions(mc);CHKERRQ(ierr);
71335efc43SPeter Brune       ierr = MatColoringApply(mc,&iscoloring);CHKERRQ(ierr);
72335efc43SPeter Brune       ierr = MatColoringDestroy(&mc);CHKERRQ(ierr);
73b0ae01b7SPeter Brune       ierr = MatFDColoringCreate(*B,iscoloring,&color);CHKERRQ(ierr);
740171d955SPeter Brune       ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
75b0ae01b7SPeter Brune       ierr = MatFDColoringSetFunction(color,(PetscErrorCode (*)(void))func,(void*)funcctx);CHKERRQ(ierr);
76b0ae01b7SPeter Brune       ierr = MatFDColoringSetFromOptions(color);CHKERRQ(ierr);
77b0ae01b7SPeter Brune     }
784e269d77SPeter Brune     ierr = PetscObjectCompose((PetscObject)*B,"SNESMatFDColoring",(PetscObject)color);CHKERRQ(ierr);
794108615fSPeter Brune     ierr = PetscObjectDereference((PetscObject)color);CHKERRQ(ierr);
804a9d489dSBarry Smith   }
8195862db2SPeter Brune 
82c89319d2SPeter Brune   /* F is only usable if there is no RHS on the SNES and the full solution corresponds to x1 */
83c89319d2SPeter Brune   ierr = VecEqual(x1,snes->vec_sol,&solvec);CHKERRQ(ierr);
84c89319d2SPeter Brune   if (!snes->vec_rhs && solvec) {
85432bb422SPeter Brune     ierr = MatFDColoringSetF(color,F);CHKERRQ(ierr);
8695862db2SPeter Brune   }
87005c665bSBarry Smith   ierr = MatFDColoringApply(*B,color,x1,flag,snes);CHKERRQ(ierr);
8832dfb669SBarry Smith   if (*J != *B) {
89194405b3SBarry Smith     ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
90194405b3SBarry Smith     ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
91194405b3SBarry Smith   }
923a40ed3dSBarry Smith   PetscFunctionReturn(0);
9391627157SBarry Smith }
94