xref: /petsc/src/mat/matfd/fdmatrix.c (revision b45d2f2cb7e031d9c0de5873eca80614ca7b863b)
1be1d678aSKris Buschelman 
2bbf0e169SBarry Smith /*
3639f9d9dSBarry Smith    This is where the abstract matrix operations are defined that are
4639f9d9dSBarry Smith   used for finite difference computations of Jacobians using coloring.
5bbf0e169SBarry Smith */
6bbf0e169SBarry Smith 
7*b45d2f2cSJed Brown #include <petsc-private/matimpl.h>        /*I "petscmat.h" I*/
8bbf0e169SBarry Smith 
94a2ae208SSatish Balay #undef __FUNCT__
103a7fca6bSBarry Smith #define __FUNCT__ "MatFDColoringSetF"
117087cfbeSBarry Smith PetscErrorCode  MatFDColoringSetF(MatFDColoring fd,Vec F)
123a7fca6bSBarry Smith {
133a7fca6bSBarry Smith   PetscFunctionBegin;
143a7fca6bSBarry Smith   fd->F = F;
153a7fca6bSBarry Smith   PetscFunctionReturn(0);
163a7fca6bSBarry Smith }
173a7fca6bSBarry Smith 
183a7fca6bSBarry Smith #undef __FUNCT__
194a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Draw_Zoom"
206849ba73SBarry Smith static PetscErrorCode MatFDColoringView_Draw_Zoom(PetscDraw draw,void *Aa)
219194eea9SBarry Smith {
229194eea9SBarry Smith   MatFDColoring  fd = (MatFDColoring)Aa;
23dfbe8321SBarry Smith   PetscErrorCode ierr;
2438baddfdSBarry Smith   PetscInt       i,j;
259194eea9SBarry Smith   PetscReal      x,y;
269194eea9SBarry Smith 
279194eea9SBarry Smith   PetscFunctionBegin;
289194eea9SBarry Smith 
299194eea9SBarry Smith   /* loop over colors  */
309194eea9SBarry Smith   for (i=0; i<fd->ncolors; i++) {
319194eea9SBarry Smith     for (j=0; j<fd->nrows[i]; j++) {
329194eea9SBarry Smith       y = fd->M - fd->rows[i][j] - fd->rstart;
339194eea9SBarry Smith       x = fd->columnsforrow[i][j];
34b0a32e0cSBarry Smith       ierr = PetscDrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr);
359194eea9SBarry Smith     }
369194eea9SBarry Smith   }
379194eea9SBarry Smith   PetscFunctionReturn(0);
389194eea9SBarry Smith }
399194eea9SBarry Smith 
404a2ae208SSatish Balay #undef __FUNCT__
414a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Draw"
426849ba73SBarry Smith static PetscErrorCode MatFDColoringView_Draw(MatFDColoring fd,PetscViewer viewer)
43005c665bSBarry Smith {
44dfbe8321SBarry Smith   PetscErrorCode ierr;
45ace3abfcSBarry Smith   PetscBool      isnull;
46b0a32e0cSBarry Smith   PetscDraw      draw;
4754d96fa3SBarry Smith   PetscReal      xr,yr,xl,yl,h,w;
48005c665bSBarry Smith 
493a40ed3dSBarry Smith   PetscFunctionBegin;
50b0a32e0cSBarry Smith   ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
51b0a32e0cSBarry Smith   ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
529194eea9SBarry Smith 
539194eea9SBarry Smith   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
54005c665bSBarry Smith 
55005c665bSBarry Smith   xr  = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0;
56005c665bSBarry Smith   xr += w;     yr += h;    xl = -w;     yl = -h;
57b0a32e0cSBarry Smith   ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
58b0a32e0cSBarry Smith   ierr = PetscDrawZoom(draw,MatFDColoringView_Draw_Zoom,fd);CHKERRQ(ierr);
599194eea9SBarry Smith   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
603a40ed3dSBarry Smith   PetscFunctionReturn(0);
61005c665bSBarry Smith }
62005c665bSBarry Smith 
634a2ae208SSatish Balay #undef __FUNCT__
644a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView"
65bbf0e169SBarry Smith /*@C
66639f9d9dSBarry Smith    MatFDColoringView - Views a finite difference coloring context.
67bbf0e169SBarry Smith 
684c49b128SBarry Smith    Collective on MatFDColoring
69fee21e36SBarry Smith 
70ef5ee4d1SLois Curfman McInnes    Input  Parameters:
71ef5ee4d1SLois Curfman McInnes +  c - the coloring context
72ef5ee4d1SLois Curfman McInnes -  viewer - visualization context
73ef5ee4d1SLois Curfman McInnes 
7415091d37SBarry Smith    Level: intermediate
7515091d37SBarry Smith 
76b4fc646aSLois Curfman McInnes    Notes:
77b4fc646aSLois Curfman McInnes    The available visualization contexts include
78b0a32e0cSBarry Smith +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
79b0a32e0cSBarry Smith .     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
80ef5ee4d1SLois Curfman McInnes         output where only the first processor opens
81ef5ee4d1SLois Curfman McInnes         the file.  All other processors send their
82ef5ee4d1SLois Curfman McInnes         data to the first processor to print.
83b0a32e0cSBarry Smith -     PETSC_VIEWER_DRAW_WORLD - graphical display of nonzero structure
84bbf0e169SBarry Smith 
859194eea9SBarry Smith    Notes:
869194eea9SBarry Smith      Since PETSc uses only a small number of basic colors (currently 33), if the coloring
87b0a32e0cSBarry Smith    involves more than 33 then some seemingly identical colors are displayed making it look
889194eea9SBarry Smith    like an illegal coloring. This is just a graphical artifact.
899194eea9SBarry Smith 
90639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
91005c665bSBarry Smith 
92b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, view
93bbf0e169SBarry Smith @*/
947087cfbeSBarry Smith PetscErrorCode  MatFDColoringView(MatFDColoring c,PetscViewer viewer)
95bbf0e169SBarry Smith {
966849ba73SBarry Smith   PetscErrorCode    ierr;
9738baddfdSBarry Smith   PetscInt          i,j;
98ace3abfcSBarry Smith   PetscBool         isdraw,iascii;
99f3ef73ceSBarry Smith   PetscViewerFormat format;
100bbf0e169SBarry Smith 
1013a40ed3dSBarry Smith   PetscFunctionBegin;
1020700a824SBarry Smith   PetscValidHeaderSpecific(c,MAT_FDCOLORING_CLASSID,1);
1033050cee2SBarry Smith   if (!viewer) {
1047adad957SLisandro Dalcin     ierr = PetscViewerASCIIGetStdout(((PetscObject)c)->comm,&viewer);CHKERRQ(ierr);
1053050cee2SBarry Smith   }
1060700a824SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
107c9780b6fSBarry Smith   PetscCheckSameComm(c,1,viewer,2);
108bbf0e169SBarry Smith 
1092692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
1102692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
1110f5bd95cSBarry Smith   if (isdraw) {
112b4fc646aSLois Curfman McInnes     ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr);
11332077d6dSBarry Smith   } else if (iascii) {
114317d6ea6SBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)c,viewer,"MatFDColoring Object");CHKERRQ(ierr);
115a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Error tolerance=%G\n",c->error_rel);CHKERRQ(ierr);
116a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Umin=%G\n",c->umin);CHKERRQ(ierr);
11777431f27SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Number of colors=%D\n",c->ncolors);CHKERRQ(ierr);
118ae09f205SBarry Smith 
119b0a32e0cSBarry Smith     ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
120fb9695e5SSatish Balay     if (format != PETSC_VIEWER_ASCII_INFO) {
121b4fc646aSLois Curfman McInnes       for (i=0; i<c->ncolors; i++) {
12277431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  Information for color %D\n",i);CHKERRQ(ierr);
12377431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    Number of columns %D\n",c->ncolumns[i]);CHKERRQ(ierr);
124b4fc646aSLois Curfman McInnes         for (j=0; j<c->ncolumns[i]; j++) {
12577431f27SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"      %D\n",c->columns[i][j]);CHKERRQ(ierr);
126639f9d9dSBarry Smith         }
12777431f27SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    Number of rows %D\n",c->nrows[i]);CHKERRQ(ierr);
128b4fc646aSLois Curfman McInnes         for (j=0; j<c->nrows[i]; j++) {
12977431f27SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"      %D %D \n",c->rows[i][j],c->columnsforrow[i][j]);CHKERRQ(ierr);
130b4fc646aSLois Curfman McInnes         }
131bbf0e169SBarry Smith       }
132bbf0e169SBarry Smith     }
133b0a32e0cSBarry Smith     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
1345cd90555SBarry Smith   } else {
135e32f2f54SBarry Smith     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for MatFDColoring",((PetscObject)viewer)->type_name);
136bbf0e169SBarry Smith   }
1373a40ed3dSBarry Smith   PetscFunctionReturn(0);
138639f9d9dSBarry Smith }
139639f9d9dSBarry Smith 
1404a2ae208SSatish Balay #undef __FUNCT__
1414a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetParameters"
142639f9d9dSBarry Smith /*@
143b4fc646aSLois Curfman McInnes    MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
144b4fc646aSLois Curfman McInnes    a Jacobian matrix using finite differences.
145639f9d9dSBarry Smith 
1463f9fe445SBarry Smith    Logically Collective on MatFDColoring
147ef5ee4d1SLois Curfman McInnes 
148ef5ee4d1SLois Curfman McInnes    The Jacobian is estimated with the differencing approximation
149ef5ee4d1SLois Curfman McInnes .vb
15065f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
151f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
152f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
153ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
154ef5ee4d1SLois Curfman McInnes .ve
155639f9d9dSBarry Smith 
156639f9d9dSBarry Smith    Input Parameters:
157ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
158639f9d9dSBarry Smith .  error_rel - relative error
159f23b5b22SLois Curfman McInnes -  umin - minimum allowable u-value magnitude
160fee21e36SBarry Smith 
16115091d37SBarry Smith    Level: advanced
16215091d37SBarry Smith 
163b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters
164b4fc646aSLois Curfman McInnes 
16595b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFromOptions()
16695b89fc3SBarry Smith 
167639f9d9dSBarry Smith @*/
1687087cfbeSBarry Smith PetscErrorCode  MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin)
169639f9d9dSBarry Smith {
1703a40ed3dSBarry Smith   PetscFunctionBegin;
1710700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
172c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(matfd,error,2);
173c5eb9154SBarry Smith   PetscValidLogicalCollectiveReal(matfd,umin,3);
174639f9d9dSBarry Smith   if (error != PETSC_DEFAULT) matfd->error_rel = error;
175639f9d9dSBarry Smith   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
1763a40ed3dSBarry Smith   PetscFunctionReturn(0);
177639f9d9dSBarry Smith }
178639f9d9dSBarry Smith 
179005c665bSBarry Smith 
180ff0cfa39SBarry Smith 
1814a2ae208SSatish Balay #undef __FUNCT__
1824a9d489dSBarry Smith #define __FUNCT__ "MatFDColoringGetFunction"
1834a9d489dSBarry Smith /*@C
1844a9d489dSBarry Smith    MatFDColoringGetFunction - Gets the function to use for computing the Jacobian.
1854a9d489dSBarry Smith 
1863f9fe445SBarry Smith    Not Collective
1874a9d489dSBarry Smith 
1884a9d489dSBarry Smith    Input Parameters:
1894a9d489dSBarry Smith .  coloring - the coloring context
1904a9d489dSBarry Smith 
1914a9d489dSBarry Smith    Output Parameters:
1924a9d489dSBarry Smith +  f - the function
1934a9d489dSBarry Smith -  fctx - the optional user-defined function context
1944a9d489dSBarry Smith 
1954a9d489dSBarry Smith    Level: intermediate
1964a9d489dSBarry Smith 
1974a9d489dSBarry Smith .keywords: Mat, Jacobian, finite differences, set, function
19895b89fc3SBarry Smith 
19995b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringSetFunction(), MatFDColoringSetFromOptions()
20095b89fc3SBarry Smith 
2014a9d489dSBarry Smith @*/
2027087cfbeSBarry Smith PetscErrorCode  MatFDColoringGetFunction(MatFDColoring matfd,PetscErrorCode (**f)(void),void **fctx)
2034a9d489dSBarry Smith {
2044a9d489dSBarry Smith   PetscFunctionBegin;
2050700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
2064a9d489dSBarry Smith   if (f) *f = matfd->f;
2074a9d489dSBarry Smith   if (fctx) *fctx = matfd->fctx;
2084a9d489dSBarry Smith   PetscFunctionReturn(0);
2094a9d489dSBarry Smith }
2104a9d489dSBarry Smith 
2114a9d489dSBarry Smith #undef __FUNCT__
2124a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFunction"
213d64ed03dSBarry Smith /*@C
214005c665bSBarry Smith    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
215005c665bSBarry Smith 
2163f9fe445SBarry Smith    Logically Collective on MatFDColoring
217fee21e36SBarry Smith 
218ef5ee4d1SLois Curfman McInnes    Input Parameters:
219ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
220ef5ee4d1SLois Curfman McInnes .  f - the function
221ef5ee4d1SLois Curfman McInnes -  fctx - the optional user-defined function context
222ef5ee4d1SLois Curfman McInnes 
2237850c7c0SBarry Smith    Calling sequence of (*f) function:
2247850c7c0SBarry Smith     For SNES:    PetscErrorCode (*f)(SNES,Vec,Vec,void*)
225ab637aeaSJed Brown     If not using SNES: PetscErrorCode (*f)(void *dummy,Vec,Vec,void*) and dummy is ignored
22615091d37SBarry Smith 
2277850c7c0SBarry Smith    Level: advanced
2287850c7c0SBarry Smith 
229ab637aeaSJed Brown    Notes: This function is usually used automatically by SNES (when one uses SNESSetJacobian() with the argument
230ab637aeaSJed Brown      SNESDefaultComputeJacobianColor()) and only needs to be used by someone computing a matrix via coloring directly by
231ab637aeaSJed Brown      calling MatFDColoringApply()
2327850c7c0SBarry Smith 
2337850c7c0SBarry Smith    Fortran Notes:
2347850c7c0SBarry Smith     In Fortran you must call MatFDColoringSetFunction() for a coloring object to
235ab637aeaSJed Brown   be used without SNES or within the SNES solvers.
236f881d145SBarry Smith 
237b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function
23895b89fc3SBarry Smith 
23995b89fc3SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringGetFunction(), MatFDColoringSetFromOptions()
24095b89fc3SBarry Smith 
241005c665bSBarry Smith @*/
2427087cfbeSBarry Smith PetscErrorCode  MatFDColoringSetFunction(MatFDColoring matfd,PetscErrorCode (*f)(void),void *fctx)
243005c665bSBarry Smith {
2443a40ed3dSBarry Smith   PetscFunctionBegin;
2450700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
246005c665bSBarry Smith   matfd->f    = f;
247005c665bSBarry Smith   matfd->fctx = fctx;
2483a40ed3dSBarry Smith   PetscFunctionReturn(0);
249005c665bSBarry Smith }
250005c665bSBarry Smith 
2514a2ae208SSatish Balay #undef __FUNCT__
2524a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringSetFromOptions"
253639f9d9dSBarry Smith /*@
254b4fc646aSLois Curfman McInnes    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
255639f9d9dSBarry Smith    the options database.
256639f9d9dSBarry Smith 
257fee21e36SBarry Smith    Collective on MatFDColoring
258fee21e36SBarry Smith 
25965f2ba5bSLois Curfman McInnes    The Jacobian, F'(u), is estimated with the differencing approximation
260ef5ee4d1SLois Curfman McInnes .vb
26165f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
262f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
263f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
264ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
265ef5ee4d1SLois Curfman McInnes .ve
266ef5ee4d1SLois Curfman McInnes 
267ef5ee4d1SLois Curfman McInnes    Input Parameter:
268ef5ee4d1SLois Curfman McInnes .  coloring - the coloring context
269ef5ee4d1SLois Curfman McInnes 
270b4fc646aSLois Curfman McInnes    Options Database Keys:
271d15ffeeaSSatish Balay +  -mat_fd_coloring_err <err> - Sets <err> (square root
272ef5ee4d1SLois Curfman McInnes            of relative error in the function)
273f23b5b22SLois Curfman McInnes .  -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
2743ec795f1SBarry Smith .  -mat_fd_type - "wp" or "ds" (see MATMFFD_WP or MATMFFD_DS)
275ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view - Activates basic viewing
276ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view_info - Activates viewing info
277ef5ee4d1SLois Curfman McInnes -  -mat_fd_coloring_view_draw - Activates drawing
278639f9d9dSBarry Smith 
27915091d37SBarry Smith     Level: intermediate
28015091d37SBarry Smith 
281b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters
282d1c39f55SBarry Smith 
283d1c39f55SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringView(), MatFDColoringSetParameters()
284d1c39f55SBarry Smith 
285639f9d9dSBarry Smith @*/
2867087cfbeSBarry Smith PetscErrorCode  MatFDColoringSetFromOptions(MatFDColoring matfd)
287639f9d9dSBarry Smith {
288dfbe8321SBarry Smith   PetscErrorCode ierr;
289ace3abfcSBarry Smith   PetscBool      flg;
290efb30889SBarry Smith   char           value[3];
2913a40ed3dSBarry Smith 
2923a40ed3dSBarry Smith   PetscFunctionBegin;
2930700a824SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_CLASSID,1);
294639f9d9dSBarry Smith 
2953194b578SJed Brown   ierr = PetscObjectOptionsBegin((PetscObject)matfd);CHKERRQ(ierr);
29687828ca2SBarry Smith     ierr = PetscOptionsReal("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr);
29787828ca2SBarry Smith     ierr = PetscOptionsReal("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr);
2981bc75a8dSBarry Smith     ierr = PetscOptionsString("-mat_fd_type","Algorithm to compute h, wp or ds","MatFDColoringCreate",matfd->htype,value,3,&flg);CHKERRQ(ierr);
299efb30889SBarry Smith     if (flg) {
300efb30889SBarry Smith       if (value[0] == 'w' && value[1] == 'p') matfd->htype = "wp";
301efb30889SBarry Smith       else if (value[0] == 'd' && value[1] == 's') matfd->htype = "ds";
302e32f2f54SBarry Smith       else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Unknown finite differencing type %s",value);
303efb30889SBarry Smith     }
304186905e3SBarry Smith     /* not used here; but so they are presented in the GUI */
305b0a32e0cSBarry Smith     ierr = PetscOptionsName("-mat_fd_coloring_view","Print entire datastructure used for Jacobian","None",0);CHKERRQ(ierr);
306b0a32e0cSBarry Smith     ierr = PetscOptionsName("-mat_fd_coloring_view_info","Print number of colors etc for Jacobian","None",0);CHKERRQ(ierr);
307b0a32e0cSBarry Smith     ierr = PetscOptionsName("-mat_fd_coloring_view_draw","Plot nonzero structure ofJacobian","None",0);CHKERRQ(ierr);
3085d973c19SBarry Smith 
3095d973c19SBarry Smith     /* process any options handlers added with PetscObjectAddOptionsHandler() */
3105d973c19SBarry Smith     ierr = PetscObjectProcessOptionsHandlers((PetscObject)matfd);CHKERRQ(ierr);
311b0a32e0cSBarry Smith   PetscOptionsEnd();CHKERRQ(ierr);
3123a40ed3dSBarry Smith   PetscFunctionReturn(0);
313005c665bSBarry Smith }
314005c665bSBarry Smith 
3154a2ae208SSatish Balay #undef __FUNCT__
3164a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringView_Private"
317dfbe8321SBarry Smith PetscErrorCode MatFDColoringView_Private(MatFDColoring fd)
318005c665bSBarry Smith {
319dfbe8321SBarry Smith   PetscErrorCode ierr;
320ace3abfcSBarry Smith   PetscBool      flg = PETSC_FALSE;
3213050cee2SBarry Smith   PetscViewer    viewer;
322005c665bSBarry Smith 
3233a40ed3dSBarry Smith   PetscFunctionBegin;
3247adad957SLisandro Dalcin   ierr = PetscViewerASCIIGetStdout(((PetscObject)fd)->comm,&viewer);CHKERRQ(ierr);
325acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-mat_fd_coloring_view",&flg,PETSC_NULL);CHKERRQ(ierr);
326005c665bSBarry Smith   if (flg) {
3273050cee2SBarry Smith     ierr = MatFDColoringView(fd,viewer);CHKERRQ(ierr);
328005c665bSBarry Smith   }
32990d69ab7SBarry Smith   flg  = PETSC_FALSE;
330acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-mat_fd_coloring_view_info",&flg,PETSC_NULL);CHKERRQ(ierr);
331ae09f205SBarry Smith   if (flg) {
3323050cee2SBarry Smith     ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr);
3333050cee2SBarry Smith     ierr = MatFDColoringView(fd,viewer);CHKERRQ(ierr);
3343050cee2SBarry Smith     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
335ae09f205SBarry Smith   }
33690d69ab7SBarry Smith   flg  = PETSC_FALSE;
337acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg,PETSC_NULL);CHKERRQ(ierr);
338005c665bSBarry Smith   if (flg) {
3397adad957SLisandro Dalcin     ierr = MatFDColoringView(fd,PETSC_VIEWER_DRAW_(((PetscObject)fd)->comm));CHKERRQ(ierr);
3407adad957SLisandro Dalcin     ierr = PetscViewerFlush(PETSC_VIEWER_DRAW_(((PetscObject)fd)->comm));CHKERRQ(ierr);
341005c665bSBarry Smith   }
3423a40ed3dSBarry Smith   PetscFunctionReturn(0);
343bbf0e169SBarry Smith }
344bbf0e169SBarry Smith 
3454a2ae208SSatish Balay #undef __FUNCT__
3464a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringCreate"
34705869f15SSatish Balay /*@
348639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
349639f9d9dSBarry Smith    computation of Jacobians.
350bbf0e169SBarry Smith 
351ef5ee4d1SLois Curfman McInnes    Collective on Mat
352ef5ee4d1SLois Curfman McInnes 
353639f9d9dSBarry Smith    Input Parameters:
354ef5ee4d1SLois Curfman McInnes +  mat - the matrix containing the nonzero structure of the Jacobian
355e727c939SJed Brown -  iscoloring - the coloring of the matrix; usually obtained with MatGetColoring() or DMCreateColoring()
356bbf0e169SBarry Smith 
357bbf0e169SBarry Smith     Output Parameter:
358639f9d9dSBarry Smith .   color - the new coloring context
359bbf0e169SBarry Smith 
36015091d37SBarry Smith     Level: intermediate
36115091d37SBarry Smith 
3626831982aSBarry Smith .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(),
363d1c39f55SBarry Smith           MatFDColoringSetFunction(), MatFDColoringSetFromOptions(), MatFDColoringApply(),
364e727c939SJed Brown           MatFDColoringView(), MatFDColoringSetParameters(), MatGetColoring(), DMCreateColoring()
365bbf0e169SBarry Smith @*/
3667087cfbeSBarry Smith PetscErrorCode  MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
367bbf0e169SBarry Smith {
368639f9d9dSBarry Smith   MatFDColoring  c;
369639f9d9dSBarry Smith   MPI_Comm       comm;
370dfbe8321SBarry Smith   PetscErrorCode ierr;
37138baddfdSBarry Smith   PetscInt       M,N;
372639f9d9dSBarry Smith 
3733a40ed3dSBarry Smith   PetscFunctionBegin;
374d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
375639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
37617186662SBarry Smith   if (M != N) SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_SUP,"Only for square matrices");
377639f9d9dSBarry Smith 
378f881d145SBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
3793194b578SJed Brown   ierr = PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_CLASSID,0,"MatFDColoring","Jacobian computation via finite differences with coloring","Mat",comm,MatFDColoringDestroy,MatFDColoringView);CHKERRQ(ierr);
380639f9d9dSBarry Smith 
381b8f8c88eSHong Zhang   c->ctype = iscoloring->ctype;
382b8f8c88eSHong Zhang 
383f830108cSBarry Smith   if (mat->ops->fdcoloringcreate) {
384f830108cSBarry Smith     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr);
38517186662SBarry Smith   } else SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_SUP,"Code not yet written for this matrix type");
386639f9d9dSBarry Smith 
387b8f8c88eSHong Zhang   ierr = MatGetVecs(mat,PETSC_NULL,&c->w1);CHKERRQ(ierr);
388b8f8c88eSHong Zhang   ierr = PetscLogObjectParent(c,c->w1);CHKERRQ(ierr);
389b8f8c88eSHong Zhang   ierr = VecDuplicate(c->w1,&c->w2);CHKERRQ(ierr);
390b8f8c88eSHong Zhang   ierr = PetscLogObjectParent(c,c->w2);CHKERRQ(ierr);
391b8f8c88eSHong Zhang 
39277d8c4bbSBarry Smith   c->error_rel         = PETSC_SQRT_MACHINE_EPSILON;
39377d8c4bbSBarry Smith   c->umin              = 100.0*PETSC_SQRT_MACHINE_EPSILON;
39449b058dcSBarry Smith   c->currentcolor      = -1;
395efb30889SBarry Smith   c->htype             = "wp";
396639f9d9dSBarry Smith 
397639f9d9dSBarry Smith   *color = c;
398d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
3993a40ed3dSBarry Smith   PetscFunctionReturn(0);
400bbf0e169SBarry Smith }
401bbf0e169SBarry Smith 
4024a2ae208SSatish Balay #undef __FUNCT__
4034a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringDestroy"
40405869f15SSatish Balay /*@
405639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
406639f9d9dSBarry Smith     via MatFDColoringCreate().
407bbf0e169SBarry Smith 
408ef5ee4d1SLois Curfman McInnes     Collective on MatFDColoring
409ef5ee4d1SLois Curfman McInnes 
410b4fc646aSLois Curfman McInnes     Input Parameter:
411639f9d9dSBarry Smith .   c - coloring context
412bbf0e169SBarry Smith 
41315091d37SBarry Smith     Level: intermediate
41415091d37SBarry Smith 
415639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
416bbf0e169SBarry Smith @*/
4176bf464f9SBarry Smith PetscErrorCode  MatFDColoringDestroy(MatFDColoring *c)
418bbf0e169SBarry Smith {
4196849ba73SBarry Smith   PetscErrorCode ierr;
42038baddfdSBarry Smith   PetscInt       i;
421bbf0e169SBarry Smith 
4223a40ed3dSBarry Smith   PetscFunctionBegin;
4236bf464f9SBarry Smith   if (!*c) PetscFunctionReturn(0);
4246bf464f9SBarry Smith   if (--((PetscObject)(*c))->refct > 0) {*c = 0; PetscFunctionReturn(0);}
425d4bb536fSBarry Smith 
4266bf464f9SBarry Smith   for (i=0; i<(*c)->ncolors; i++) {
4276bf464f9SBarry Smith     ierr = PetscFree((*c)->columns[i]);CHKERRQ(ierr);
4286bf464f9SBarry Smith     ierr = PetscFree((*c)->rows[i]);CHKERRQ(ierr);
4296bf464f9SBarry Smith     ierr = PetscFree((*c)->columnsforrow[i]);CHKERRQ(ierr);
4306bf464f9SBarry Smith     if ((*c)->vscaleforrow) {ierr = PetscFree((*c)->vscaleforrow[i]);CHKERRQ(ierr);}
431bbf0e169SBarry Smith   }
4326bf464f9SBarry Smith   ierr = PetscFree((*c)->ncolumns);CHKERRQ(ierr);
4336bf464f9SBarry Smith   ierr = PetscFree((*c)->columns);CHKERRQ(ierr);
4346bf464f9SBarry Smith   ierr = PetscFree((*c)->nrows);CHKERRQ(ierr);
4356bf464f9SBarry Smith   ierr = PetscFree((*c)->rows);CHKERRQ(ierr);
4366bf464f9SBarry Smith   ierr = PetscFree((*c)->columnsforrow);CHKERRQ(ierr);
4376bf464f9SBarry Smith   ierr = PetscFree((*c)->vscaleforrow);CHKERRQ(ierr);
4386bf464f9SBarry Smith   ierr = VecDestroy(&(*c)->vscale);CHKERRQ(ierr);
4396bf464f9SBarry Smith   ierr = VecDestroy(&(*c)->w1);CHKERRQ(ierr);
4406bf464f9SBarry Smith   ierr = VecDestroy(&(*c)->w2);CHKERRQ(ierr);
4416bf464f9SBarry Smith   ierr = VecDestroy(&(*c)->w3);CHKERRQ(ierr);
442d38fa0fbSBarry Smith   ierr = PetscHeaderDestroy(c);CHKERRQ(ierr);
4433a40ed3dSBarry Smith   PetscFunctionReturn(0);
444bbf0e169SBarry Smith }
44543a90d84SBarry Smith 
4464a2ae208SSatish Balay #undef __FUNCT__
44749b058dcSBarry Smith #define __FUNCT__ "MatFDColoringGetPerturbedColumns"
44849b058dcSBarry Smith /*@C
44949b058dcSBarry Smith     MatFDColoringGetPerturbedColumns - Returns the indices of the columns that
45049b058dcSBarry Smith       that are currently being perturbed.
45149b058dcSBarry Smith 
45249b058dcSBarry Smith     Not Collective
45349b058dcSBarry Smith 
45449b058dcSBarry Smith     Input Parameters:
45549b058dcSBarry Smith .   coloring - coloring context created with MatFDColoringCreate()
45649b058dcSBarry Smith 
45749b058dcSBarry Smith     Output Parameters:
45849b058dcSBarry Smith +   n - the number of local columns being perturbed
45949b058dcSBarry Smith -   cols - the column indices, in global numbering
46049b058dcSBarry Smith 
46149b058dcSBarry Smith    Level: intermediate
46249b058dcSBarry Smith 
46349b058dcSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringApply()
46449b058dcSBarry Smith 
46549b058dcSBarry Smith .keywords: coloring, Jacobian, finite differences
46649b058dcSBarry Smith @*/
4677087cfbeSBarry Smith PetscErrorCode  MatFDColoringGetPerturbedColumns(MatFDColoring coloring,PetscInt *n,PetscInt *cols[])
46849b058dcSBarry Smith {
46949b058dcSBarry Smith   PetscFunctionBegin;
47049b058dcSBarry Smith   if (coloring->currentcolor >= 0) {
47149b058dcSBarry Smith     *n    = coloring->ncolumns[coloring->currentcolor];
47249b058dcSBarry Smith     *cols = coloring->columns[coloring->currentcolor];
47349b058dcSBarry Smith   } else {
47449b058dcSBarry Smith     *n = 0;
47549b058dcSBarry Smith   }
47649b058dcSBarry Smith   PetscFunctionReturn(0);
47749b058dcSBarry Smith }
47849b058dcSBarry Smith 
47949b058dcSBarry Smith #undef __FUNCT__
4804a2ae208SSatish Balay #define __FUNCT__ "MatFDColoringApply"
48143a90d84SBarry Smith /*@
482e0907662SLois Curfman McInnes     MatFDColoringApply - Given a matrix for which a MatFDColoring context
483e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
48443a90d84SBarry Smith 
485fee21e36SBarry Smith     Collective on MatFDColoring
486fee21e36SBarry Smith 
487ef5ee4d1SLois Curfman McInnes     Input Parameters:
488ef5ee4d1SLois Curfman McInnes +   mat - location to store Jacobian
489ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
490ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
4917850c7c0SBarry Smith -   sctx - context required by function, if this is being used with the SNES solver then it is SNES object, otherwise it is null
492ef5ee4d1SLois Curfman McInnes 
4938bba8e72SBarry Smith     Options Database Keys:
494ebd3b9afSBarry Smith +    -mat_fd_type - "wp" or "ds"  (see MATMFFD_WP or MATMFFD_DS)
495b9722508SLois Curfman McInnes .    -mat_fd_coloring_view - Activates basic viewing or coloring
496b9722508SLois Curfman McInnes .    -mat_fd_coloring_view_draw - Activates drawing of coloring
497b9722508SLois Curfman McInnes -    -mat_fd_coloring_view_info - Activates viewing of coloring info
4988bba8e72SBarry Smith 
49915091d37SBarry Smith     Level: intermediate
50098d79826SSatish Balay 
5017850c7c0SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView(), MatFDColoringSetFunction()
50243a90d84SBarry Smith 
50343a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
50443a90d84SBarry Smith @*/
5057087cfbeSBarry Smith PetscErrorCode  MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
50643a90d84SBarry Smith {
5073acb8795SBarry Smith   PetscErrorCode ierr;
5083acb8795SBarry Smith 
5093acb8795SBarry Smith   PetscFunctionBegin;
5100700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
5110700a824SBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2);
5120700a824SBarry Smith   PetscValidHeaderSpecific(x1,VEC_CLASSID,3);
51317186662SBarry Smith   if (!coloring->f) SETERRQ(((PetscObject)J)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");
514e32f2f54SBarry Smith   if (!J->ops->fdcoloringapply) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not supported for this matrix type %s",((PetscObject)J)->type_name);
5153acb8795SBarry Smith   ierr = (*J->ops->fdcoloringapply)(J,coloring,x1,flag,sctx);CHKERRQ(ierr);
5163acb8795SBarry Smith   PetscFunctionReturn(0);
5173acb8795SBarry Smith }
5183acb8795SBarry Smith 
5193acb8795SBarry Smith #undef __FUNCT__
5203acb8795SBarry Smith #define __FUNCT__ "MatFDColoringApply_AIJ"
5217087cfbeSBarry Smith PetscErrorCode  MatFDColoringApply_AIJ(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
5223acb8795SBarry Smith {
5236849ba73SBarry Smith   PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void *))coloring->f;
5246849ba73SBarry Smith   PetscErrorCode ierr;
525b8f8c88eSHong Zhang   PetscInt       k,start,end,l,row,col,srow,**vscaleforrow,m1,m2;
526efb30889SBarry Smith   PetscScalar    dx,*y,*xx,*w3_array;
527ea709b57SSatish Balay   PetscScalar    *vscale_array;
528efb30889SBarry Smith   PetscReal      epsilon = coloring->error_rel,umin = coloring->umin,unorm;
529b8f8c88eSHong Zhang   Vec            w1=coloring->w1,w2=coloring->w2,w3;
530005c665bSBarry Smith   void           *fctx = coloring->fctx;
531ace3abfcSBarry Smith   PetscBool      flg = PETSC_FALSE;
532705d48f7SSatish Balay   PetscInt       ctype=coloring->ctype,N,col_start=0,col_end=0;
533b8f8c88eSHong Zhang   Vec            x1_tmp;
534a2e34c3dSBarry Smith 
5353a40ed3dSBarry Smith   PetscFunctionBegin;
5360700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
5370700a824SBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,2);
5380700a824SBarry Smith   PetscValidHeaderSpecific(x1,VEC_CLASSID,3);
53917186662SBarry Smith   if (!f) SETERRQ(((PetscObject)J)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must call MatFDColoringSetFunction()");
540e0907662SLois Curfman McInnes 
541d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
5424c49b128SBarry Smith   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
543acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg,PETSC_NULL);CHKERRQ(ierr);
544f1af5d2fSBarry Smith   if (flg) {
545ae15b995SBarry Smith     ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr);
546e0907662SLois Curfman McInnes   } else {
547ace3abfcSBarry Smith     PetscBool  assembled;
5480b9b6f31SBarry Smith     ierr = MatAssembled(J,&assembled);CHKERRQ(ierr);
5490b9b6f31SBarry Smith     if (assembled) {
55043a90d84SBarry Smith       ierr = MatZeroEntries(J);CHKERRQ(ierr);
551e0907662SLois Curfman McInnes     }
5520b9b6f31SBarry Smith   }
55343a90d84SBarry Smith 
554b8f8c88eSHong Zhang   x1_tmp = x1;
555dac7f5fdSBarry Smith   if (!coloring->vscale){
556b8f8c88eSHong Zhang     ierr = VecDuplicate(x1_tmp,&coloring->vscale);CHKERRQ(ierr);
557b8f8c88eSHong Zhang   }
558be2fbe1fSBarry Smith 
5593a7fca6bSBarry Smith   /*
5603a7fca6bSBarry Smith     This is a horrible, horrible, hack. See DMMGComputeJacobian_Multigrid() it inproperly sets
5613a7fca6bSBarry Smith     coloring->F for the coarser grids from the finest
5623a7fca6bSBarry Smith   */
5633a7fca6bSBarry Smith   if (coloring->F) {
5643a7fca6bSBarry Smith     ierr = VecGetLocalSize(coloring->F,&m1);CHKERRQ(ierr);
5653a7fca6bSBarry Smith     ierr = VecGetLocalSize(w1,&m2);CHKERRQ(ierr);
5663a7fca6bSBarry Smith     if (m1 != m2) {
5673a7fca6bSBarry Smith       coloring->F = 0;
5683a7fca6bSBarry Smith       }
5693a7fca6bSBarry Smith     }
5703a7fca6bSBarry Smith 
571b8f8c88eSHong Zhang   if (coloring->htype[0] == 'w') { /* tacky test; need to make systematic if we add other approaches to computing h*/
572b8f8c88eSHong Zhang     ierr = VecNorm(x1_tmp,NORM_2,&unorm);CHKERRQ(ierr);
573b8f8c88eSHong Zhang   }
574b8f8c88eSHong Zhang   ierr = VecGetOwnershipRange(w1,&start,&end);CHKERRQ(ierr); /* OwnershipRange is used by ghosted x! */
575b8f8c88eSHong Zhang 
576b8f8c88eSHong Zhang   /* Set w1 = F(x1) */
577be2fbe1fSBarry Smith   if (coloring->F) {
578be2fbe1fSBarry Smith     w1          = coloring->F; /* use already computed value of function */
579be2fbe1fSBarry Smith     coloring->F = 0;
580be2fbe1fSBarry Smith   } else {
58166f9b7ceSBarry Smith     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
582b8f8c88eSHong Zhang     ierr = (*f)(sctx,x1_tmp,w1,fctx);CHKERRQ(ierr);
58366f9b7ceSBarry Smith     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
584be2fbe1fSBarry Smith   }
58543a90d84SBarry Smith 
586b8f8c88eSHong Zhang   if (!coloring->w3) {
587b8f8c88eSHong Zhang     ierr = VecDuplicate(x1_tmp,&coloring->w3);CHKERRQ(ierr);
588b8f8c88eSHong Zhang     ierr = PetscLogObjectParent(coloring,coloring->w3);CHKERRQ(ierr);
589efb30889SBarry Smith   }
590b8f8c88eSHong Zhang   w3 = coloring->w3;
591efb30889SBarry Smith 
592b8f8c88eSHong Zhang     /* Compute all the local scale factors, including ghost points */
593b8f8c88eSHong Zhang   ierr = VecGetLocalSize(x1_tmp,&N);CHKERRQ(ierr);
594b8f8c88eSHong Zhang   ierr = VecGetArray(x1_tmp,&xx);CHKERRQ(ierr);
595b8f8c88eSHong Zhang   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
596b8f8c88eSHong Zhang   if (ctype == IS_COLORING_GHOSTED){
597b8f8c88eSHong Zhang     col_start = 0; col_end = N;
5988ee2e534SBarry Smith   } else if (ctype == IS_COLORING_GLOBAL){
599b8f8c88eSHong Zhang     xx = xx - start;
600b8f8c88eSHong Zhang     vscale_array = vscale_array - start;
601b8f8c88eSHong Zhang     col_start = start; col_end = N + start;
602b8f8c88eSHong Zhang   }
603b8f8c88eSHong Zhang   for (col=col_start; col<col_end; col++){
604b8f8c88eSHong Zhang     /* Loop over each local column, vscale[col] = 1./(epsilon*dx[col]) */
605efb30889SBarry Smith     if (coloring->htype[0] == 'w') {
606efb30889SBarry Smith       dx = 1.0 + unorm;
607efb30889SBarry Smith     } else {
6089782cf98SBarry Smith       dx  = xx[col];
609efb30889SBarry Smith     }
610d4a378daSJed Brown     if (dx == (PetscScalar)0.0) dx = 1.0;
6119782cf98SBarry Smith #if !defined(PETSC_USE_COMPLEX)
6129782cf98SBarry Smith     if (dx < umin && dx >= 0.0)      dx = umin;
6139782cf98SBarry Smith     else if (dx < 0.0 && dx > -umin) dx = -umin;
6149782cf98SBarry Smith #else
6159782cf98SBarry Smith     if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
6169782cf98SBarry Smith     else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
6179782cf98SBarry Smith #endif
6189782cf98SBarry Smith     dx               *= epsilon;
619d4a378daSJed Brown     vscale_array[col] = (PetscScalar)1.0/dx;
6209782cf98SBarry Smith   }
6218ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL)  vscale_array = vscale_array + start;
622efb30889SBarry Smith   ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
6238ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL){
62430b34957SBarry Smith     ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
62530b34957SBarry Smith     ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
626b8f8c88eSHong Zhang   }
6279782cf98SBarry Smith 
628b8f8c88eSHong Zhang   if (coloring->vscaleforrow) {
629b8f8c88eSHong Zhang     vscaleforrow = coloring->vscaleforrow;
63017186662SBarry Smith   } else SETERRQ(((PetscObject)J)->comm,PETSC_ERR_ARG_NULL,"Null Object: coloring->vscaleforrow");
631b0a32e0cSBarry Smith 
6329782cf98SBarry Smith   /*
63343a90d84SBarry Smith     Loop over each color
63443a90d84SBarry Smith   */
635b8f8c88eSHong Zhang   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
63643a90d84SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
63749b058dcSBarry Smith     coloring->currentcolor = k;
638b8f8c88eSHong Zhang     ierr = VecCopy(x1_tmp,w3);CHKERRQ(ierr);
639b8f8c88eSHong Zhang     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);
6408ee2e534SBarry Smith     if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array - start;
64143a90d84SBarry Smith     /*
642b8f8c88eSHong Zhang       Loop over each column associated with color
643b8f8c88eSHong Zhang       adding the perturbation to the vector w3.
64443a90d84SBarry Smith     */
64543a90d84SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
646b8f8c88eSHong Zhang       col = coloring->columns[k][l];    /* local column of the matrix we are probing for */
647efb30889SBarry Smith       if (coloring->htype[0] == 'w') {
648efb30889SBarry Smith         dx = 1.0 + unorm;
649efb30889SBarry Smith       } else {
65042460c72SBarry Smith         dx  = xx[col];
651efb30889SBarry Smith       }
652d4a378daSJed Brown       if (dx == (PetscScalar)0.0) dx = 1.0;
653aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
654ae09f205SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
655ae09f205SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
65643a90d84SBarry Smith #else
657329f5518SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
658329f5518SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
65943a90d84SBarry Smith #endif
66043a90d84SBarry Smith       dx            *= epsilon;
661e32f2f54SBarry Smith       if (!PetscAbsScalar(dx)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Computed 0 differencing parameter");
66242460c72SBarry Smith       w3_array[col] += dx;
66343a90d84SBarry Smith     }
66498d79826SSatish Balay     if (ctype == IS_COLORING_GLOBAL) w3_array = w3_array + start;
665b8f8c88eSHong Zhang     ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
6663b28642cSBarry Smith 
66743a90d84SBarry Smith     /*
668b8f8c88eSHong Zhang       Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations)
669b8f8c88eSHong Zhang                            w2 = F(x1 + dx) - F(x1)
67043a90d84SBarry Smith     */
67166f9b7ceSBarry Smith     ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
67243a90d84SBarry Smith     ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
67366f9b7ceSBarry Smith     ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr);
674efb30889SBarry Smith     ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr);
6759782cf98SBarry Smith 
67643a90d84SBarry Smith     /*
677e0907662SLois Curfman McInnes       Loop over rows of vector, putting results into Jacobian matrix
67843a90d84SBarry Smith     */
6793b28642cSBarry Smith     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
68043a90d84SBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
681b8f8c88eSHong Zhang       row    = coloring->rows[k][l];             /* local row index */
682b8f8c88eSHong Zhang       col    = coloring->columnsforrow[k][l];    /* global column index */
6835904e1b1SBarry Smith       y[row] *= vscale_array[vscaleforrow[k][l]];
68443a90d84SBarry Smith       srow   = row + start;
68543a90d84SBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
68643a90d84SBarry Smith     }
6873b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
688aeb7e98aSMatthew Knepley   } /* endof for each color */
6898ee2e534SBarry Smith   if (ctype == IS_COLORING_GLOBAL) xx = xx + start;
69030b34957SBarry Smith   ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
691b8f8c88eSHong Zhang   ierr = VecRestoreArray(x1_tmp,&xx);CHKERRQ(ierr);
692b8f8c88eSHong Zhang 
693b8f8c88eSHong Zhang   coloring->currentcolor = -1;
69443a90d84SBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
69543a90d84SBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
696d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
697a2e34c3dSBarry Smith 
69890d69ab7SBarry Smith   flg  = PETSC_FALSE;
699acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-mat_null_space_test",&flg,PETSC_NULL);CHKERRQ(ierr);
700a2e34c3dSBarry Smith   if (flg) {
70195902228SMatthew Knepley     ierr = MatNullSpaceTest(J->nullsp,J,PETSC_NULL);CHKERRQ(ierr);
702a2e34c3dSBarry Smith   }
703b9722508SLois Curfman McInnes   ierr = MatFDColoringView_Private(coloring);CHKERRQ(ierr);
7043a40ed3dSBarry Smith   PetscFunctionReturn(0);
70543a90d84SBarry Smith }
706