xref: /petsc/src/mat/matfd/fdmatrix.c (revision 29bbc08cd5461c366aba6645924e8ff42acd1de0)
1*29bbc08cSBarry Smith /*$Id: fdmatrix.c,v 1.77 2000/09/22 20:44:38 bsmith Exp bsmith $*/
2bbf0e169SBarry Smith 
3bbf0e169SBarry Smith /*
4639f9d9dSBarry Smith    This is where the abstract matrix operations are defined that are
5639f9d9dSBarry Smith   used for finite difference computations of Jacobians using coloring.
6bbf0e169SBarry Smith */
7bbf0e169SBarry Smith 
8bbf0e169SBarry Smith #include "petsc.h"
9e090d566SSatish Balay #include "src/mat/matimpl.h"        /*I "petscmat.h" I*/
10bbf0e169SBarry Smith 
115615d1e5SSatish Balay #undef __FUNC__
129194eea9SBarry Smith #define __FUNC__ /*<a name="MatFDColoringView_Draw_Zoom"></a>*/"MatFDColoringView_Draw_Zoom"
139194eea9SBarry Smith static int MatFDColoringView_Draw_Zoom(Draw draw,void *Aa)
149194eea9SBarry Smith {
159194eea9SBarry Smith   MatFDColoring fd = (MatFDColoring)Aa;
169194eea9SBarry Smith   int           ierr,i,j;
179194eea9SBarry Smith   PetscReal     x,y;
189194eea9SBarry Smith 
199194eea9SBarry Smith   PetscFunctionBegin;
209194eea9SBarry Smith 
219194eea9SBarry Smith   /* loop over colors  */
229194eea9SBarry Smith   for (i=0; i<fd->ncolors; i++) {
239194eea9SBarry Smith     for (j=0; j<fd->nrows[i]; j++) {
249194eea9SBarry Smith       y = fd->M - fd->rows[i][j] - fd->rstart;
259194eea9SBarry Smith       x = fd->columnsforrow[i][j];
269194eea9SBarry Smith       ierr = DrawRectangle(draw,x,y,x+1,y+1,i+1,i+1,i+1,i+1);CHKERRQ(ierr);
279194eea9SBarry Smith     }
289194eea9SBarry Smith   }
299194eea9SBarry Smith   PetscFunctionReturn(0);
309194eea9SBarry Smith }
319194eea9SBarry Smith 
329194eea9SBarry Smith #undef __FUNC__
339194eea9SBarry Smith #define __FUNC__ /*<a name="MatFDColoringView_Draw"></a>*/"MatFDColoringView_Draw"
34005c665bSBarry Smith static int MatFDColoringView_Draw(MatFDColoring fd,Viewer viewer)
35005c665bSBarry Smith {
369194eea9SBarry Smith   int         ierr;
37005c665bSBarry Smith   PetscTruth  isnull;
38005c665bSBarry Smith   Draw        draw;
3954d96fa3SBarry Smith   PetscReal   xr,yr,xl,yl,h,w;
40005c665bSBarry Smith 
413a40ed3dSBarry Smith   PetscFunctionBegin;
4277ed5343SBarry Smith   ierr = ViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
433a40ed3dSBarry Smith   ierr = DrawIsNull(draw,&isnull);CHKERRQ(ierr); if (isnull) PetscFunctionReturn(0);
449194eea9SBarry Smith 
459194eea9SBarry Smith   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr);
46005c665bSBarry Smith 
47005c665bSBarry Smith   xr  = fd->N; yr = fd->M; h = yr/10.0; w = xr/10.0;
48005c665bSBarry Smith   xr += w;     yr += h;    xl = -w;     yl = -h;
49005c665bSBarry Smith   ierr = DrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr);
509194eea9SBarry Smith   ierr = DrawZoom(draw,MatFDColoringView_Draw_Zoom,fd);CHKERRQ(ierr);
519194eea9SBarry Smith   ierr = PetscObjectCompose((PetscObject)fd,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr);
523a40ed3dSBarry Smith   PetscFunctionReturn(0);
53005c665bSBarry Smith }
54005c665bSBarry Smith 
55005c665bSBarry Smith #undef __FUNC__
56b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatFDColoringView"
57bbf0e169SBarry Smith /*@C
58639f9d9dSBarry Smith    MatFDColoringView - Views a finite difference coloring context.
59bbf0e169SBarry Smith 
604c49b128SBarry Smith    Collective on MatFDColoring
61fee21e36SBarry Smith 
62ef5ee4d1SLois Curfman McInnes    Input  Parameters:
63ef5ee4d1SLois Curfman McInnes +  c - the coloring context
64ef5ee4d1SLois Curfman McInnes -  viewer - visualization context
65ef5ee4d1SLois Curfman McInnes 
6615091d37SBarry Smith    Level: intermediate
6715091d37SBarry Smith 
68b4fc646aSLois Curfman McInnes    Notes:
69b4fc646aSLois Curfman McInnes    The available visualization contexts include
70ef5ee4d1SLois Curfman McInnes +     VIEWER_STDOUT_SELF - standard output (default)
71ef5ee4d1SLois Curfman McInnes .     VIEWER_STDOUT_WORLD - synchronized standard
72ef5ee4d1SLois Curfman McInnes         output where only the first processor opens
73ef5ee4d1SLois Curfman McInnes         the file.  All other processors send their
74ef5ee4d1SLois Curfman McInnes         data to the first processor to print.
75c655490fSBarry Smith -     VIEWER_DRAW_WORLD - graphical display of nonzero structure
76bbf0e169SBarry Smith 
779194eea9SBarry Smith    Notes:
789194eea9SBarry Smith      Since PETSc uses only a small number of basic colors (currently 33), if the coloring
799194eea9SBarry Smith    involves moe than 33 then some seemingly identical colors are displayed making it look
809194eea9SBarry Smith    like an illegal coloring. This is just a graphical artifact.
819194eea9SBarry Smith 
82639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
83005c665bSBarry Smith 
84b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, view
85bbf0e169SBarry Smith @*/
86b4fc646aSLois Curfman McInnes int MatFDColoringView(MatFDColoring c,Viewer viewer)
87bbf0e169SBarry Smith {
88639f9d9dSBarry Smith   int        i,j,format,ierr;
896831982aSBarry Smith   PetscTruth isdraw,isascii;
90bbf0e169SBarry Smith 
913a40ed3dSBarry Smith   PetscFunctionBegin;
92b4fc646aSLois Curfman McInnes   PetscValidHeaderSpecific(c,MAT_FDCOLORING_COOKIE);
933eda8832SBarry Smith   if (!viewer) viewer = VIEWER_STDOUT_(c->comm);
940f5bd95cSBarry Smith   PetscValidHeaderSpecific(viewer,VIEWER_COOKIE);
956831982aSBarry Smith   PetscCheckSameComm(c,viewer);
96bbf0e169SBarry Smith 
976831982aSBarry Smith   ierr  = PetscTypeCompare((PetscObject)viewer,DRAW_VIEWER,&isdraw);CHKERRQ(ierr);
986831982aSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr);
990f5bd95cSBarry Smith   if (isdraw) {
100b4fc646aSLois Curfman McInnes     ierr = MatFDColoringView_Draw(c,viewer);CHKERRQ(ierr);
1010f5bd95cSBarry Smith   } else if (isascii) {
102d132466eSBarry Smith     ierr = ViewerASCIIPrintf(viewer,"MatFDColoring Object:\n");CHKERRQ(ierr);
103d132466eSBarry Smith     ierr = ViewerASCIIPrintf(viewer,"  Error tolerance=%g\n",c->error_rel);CHKERRQ(ierr);
104d132466eSBarry Smith     ierr = ViewerASCIIPrintf(viewer,"  Umin=%g\n",c->umin);CHKERRQ(ierr);
105d132466eSBarry Smith     ierr = ViewerASCIIPrintf(viewer,"  Number of colors=%d\n",c->ncolors);CHKERRQ(ierr);
106ae09f205SBarry Smith 
107ae09f205SBarry Smith     ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr);
108ae09f205SBarry Smith     if (format != VIEWER_FORMAT_ASCII_INFO) {
109b4fc646aSLois Curfman McInnes       for (i=0; i<c->ncolors; i++) {
110d132466eSBarry Smith         ierr = ViewerASCIIPrintf(viewer,"  Information for color %d\n",i);CHKERRQ(ierr);
111d132466eSBarry Smith         ierr = ViewerASCIIPrintf(viewer,"    Number of columns %d\n",c->ncolumns[i]);CHKERRQ(ierr);
112b4fc646aSLois Curfman McInnes         for (j=0; j<c->ncolumns[i]; j++) {
113d132466eSBarry Smith           ierr = ViewerASCIIPrintf(viewer,"      %d\n",c->columns[i][j]);CHKERRQ(ierr);
114639f9d9dSBarry Smith         }
115d132466eSBarry Smith         ierr = ViewerASCIIPrintf(viewer,"    Number of rows %d\n",c->nrows[i]);CHKERRQ(ierr);
116b4fc646aSLois Curfman McInnes         for (j=0; j<c->nrows[i]; j++) {
117d132466eSBarry Smith           ierr = ViewerASCIIPrintf(viewer,"      %d %d \n",c->rows[i][j],c->columnsforrow[i][j]);CHKERRQ(ierr);
118b4fc646aSLois Curfman McInnes         }
119bbf0e169SBarry Smith       }
120bbf0e169SBarry Smith     }
1216831982aSBarry Smith     ierr = ViewerFlush(viewer);CHKERRQ(ierr);
1225cd90555SBarry Smith   } else {
123*29bbc08cSBarry Smith     SETERRQ1(1,"Viewer type %s not supported for MatFDColoring",((PetscObject)viewer)->type_name);
124bbf0e169SBarry Smith   }
1253a40ed3dSBarry Smith   PetscFunctionReturn(0);
126639f9d9dSBarry Smith }
127639f9d9dSBarry Smith 
1285615d1e5SSatish Balay #undef __FUNC__
129b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatFDColoringSetParameters"
130639f9d9dSBarry Smith /*@
131b4fc646aSLois Curfman McInnes    MatFDColoringSetParameters - Sets the parameters for the sparse approximation of
132b4fc646aSLois Curfman McInnes    a Jacobian matrix using finite differences.
133639f9d9dSBarry Smith 
134ef5ee4d1SLois Curfman McInnes    Collective on MatFDColoring
135ef5ee4d1SLois Curfman McInnes 
136ef5ee4d1SLois Curfman McInnes    The Jacobian is estimated with the differencing approximation
137ef5ee4d1SLois Curfman McInnes .vb
13865f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
139f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
140f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
141ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
142ef5ee4d1SLois Curfman McInnes .ve
143639f9d9dSBarry Smith 
144639f9d9dSBarry Smith    Input Parameters:
145ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
146639f9d9dSBarry Smith .  error_rel - relative error
147f23b5b22SLois Curfman McInnes -  umin - minimum allowable u-value magnitude
148fee21e36SBarry Smith 
14915091d37SBarry Smith    Level: advanced
15015091d37SBarry Smith 
151b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, parameters
152b4fc646aSLois Curfman McInnes 
153b4fc646aSLois Curfman McInnes .seealso: MatFDColoringCreate()
154639f9d9dSBarry Smith @*/
155329f5518SBarry Smith int MatFDColoringSetParameters(MatFDColoring matfd,PetscReal error,PetscReal umin)
156639f9d9dSBarry Smith {
1573a40ed3dSBarry Smith   PetscFunctionBegin;
158639f9d9dSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
159639f9d9dSBarry Smith 
160639f9d9dSBarry Smith   if (error != PETSC_DEFAULT) matfd->error_rel = error;
161639f9d9dSBarry Smith   if (umin != PETSC_DEFAULT)  matfd->umin      = umin;
1623a40ed3dSBarry Smith   PetscFunctionReturn(0);
163639f9d9dSBarry Smith }
164639f9d9dSBarry Smith 
1655615d1e5SSatish Balay #undef __FUNC__
166b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatFDColoringSetFrequency"
167005c665bSBarry Smith /*@
168e0907662SLois Curfman McInnes    MatFDColoringSetFrequency - Sets the frequency for computing new Jacobian
169e0907662SLois Curfman McInnes    matrices.
170005c665bSBarry Smith 
171fee21e36SBarry Smith    Collective on MatFDColoring
172fee21e36SBarry Smith 
173ef5ee4d1SLois Curfman McInnes    Input Parameters:
174ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
175ef5ee4d1SLois Curfman McInnes -  freq - frequency (default is 1)
176ef5ee4d1SLois Curfman McInnes 
17715091d37SBarry Smith    Options Database Keys:
17815091d37SBarry Smith .  -mat_fd_coloring_freq <freq>  - Sets coloring frequency
17915091d37SBarry Smith 
18015091d37SBarry Smith    Level: advanced
18115091d37SBarry Smith 
182e0907662SLois Curfman McInnes    Notes:
183e0907662SLois Curfman McInnes    Using a modified Newton strategy, where the Jacobian remains fixed for several
184e0907662SLois Curfman McInnes    iterations, can be cost effective in terms of overall nonlinear solution
185e0907662SLois Curfman McInnes    efficiency.  This parameter indicates that a new Jacobian will be computed every
186e0907662SLois Curfman McInnes    <freq> nonlinear iterations.
187e0907662SLois Curfman McInnes 
188b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, coloring, set, frequency
189ef5ee4d1SLois Curfman McInnes 
190ef5ee4d1SLois Curfman McInnes .seealso: MatFDColoringCreate(), MatFDColoringGetFrequency()
191005c665bSBarry Smith @*/
192005c665bSBarry Smith int MatFDColoringSetFrequency(MatFDColoring matfd,int freq)
193005c665bSBarry Smith {
1943a40ed3dSBarry Smith   PetscFunctionBegin;
195005c665bSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
196005c665bSBarry Smith 
197005c665bSBarry Smith   matfd->freq = freq;
1983a40ed3dSBarry Smith   PetscFunctionReturn(0);
199005c665bSBarry Smith }
200005c665bSBarry Smith 
201005c665bSBarry Smith #undef __FUNC__
202b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatFDColoringGetFrequency"
203ff0cfa39SBarry Smith /*@
204ff0cfa39SBarry Smith    MatFDColoringGetFrequency - Gets the frequency for computing new Jacobian
205ff0cfa39SBarry Smith    matrices.
206ff0cfa39SBarry Smith 
207ef5ee4d1SLois Curfman McInnes    Not Collective
208ef5ee4d1SLois Curfman McInnes 
209ff0cfa39SBarry Smith    Input Parameters:
210ff0cfa39SBarry Smith .  coloring - the coloring context
211ff0cfa39SBarry Smith 
212ff0cfa39SBarry Smith    Output Parameters:
213ff0cfa39SBarry Smith .  freq - frequency (default is 1)
214ff0cfa39SBarry Smith 
21515091d37SBarry Smith    Options Database Keys:
21615091d37SBarry Smith .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
21715091d37SBarry Smith 
21815091d37SBarry Smith    Level: advanced
21915091d37SBarry Smith 
220ff0cfa39SBarry Smith    Notes:
221ff0cfa39SBarry Smith    Using a modified Newton strategy, where the Jacobian remains fixed for several
222ff0cfa39SBarry Smith    iterations, can be cost effective in terms of overall nonlinear solution
223ff0cfa39SBarry Smith    efficiency.  This parameter indicates that a new Jacobian will be computed every
224ff0cfa39SBarry Smith    <freq> nonlinear iterations.
225ff0cfa39SBarry Smith 
226ff0cfa39SBarry Smith .keywords: Mat, finite differences, coloring, get, frequency
227ef5ee4d1SLois Curfman McInnes 
228ef5ee4d1SLois Curfman McInnes .seealso: MatFDColoringSetFrequency()
229ff0cfa39SBarry Smith @*/
230ff0cfa39SBarry Smith int MatFDColoringGetFrequency(MatFDColoring matfd,int *freq)
231ff0cfa39SBarry Smith {
2323a40ed3dSBarry Smith   PetscFunctionBegin;
233ff0cfa39SBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
234ff0cfa39SBarry Smith 
235ff0cfa39SBarry Smith   *freq = matfd->freq;
2363a40ed3dSBarry Smith   PetscFunctionReturn(0);
237ff0cfa39SBarry Smith }
238ff0cfa39SBarry Smith 
239ff0cfa39SBarry Smith #undef __FUNC__
240b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatFDColoringSetFunction"
241d64ed03dSBarry Smith /*@C
242005c665bSBarry Smith    MatFDColoringSetFunction - Sets the function to use for computing the Jacobian.
243005c665bSBarry Smith 
244fee21e36SBarry Smith    Collective on MatFDColoring
245fee21e36SBarry Smith 
246ef5ee4d1SLois Curfman McInnes    Input Parameters:
247ef5ee4d1SLois Curfman McInnes +  coloring - the coloring context
248ef5ee4d1SLois Curfman McInnes .  f - the function
249ef5ee4d1SLois Curfman McInnes -  fctx - the optional user-defined function context
250ef5ee4d1SLois Curfman McInnes 
25115091d37SBarry Smith    Level: intermediate
25215091d37SBarry Smith 
253f881d145SBarry Smith    Notes:
254f881d145SBarry Smith     In Fortran you must call MatFDColoringSetFunctionSNES() for a coloring object to
255f881d145SBarry Smith   be used with the SNES solvers and MatFDColoringSetFunctionTS() if it is to be used
256f881d145SBarry Smith   with the TS solvers.
257f881d145SBarry Smith 
258b4fc646aSLois Curfman McInnes .keywords: Mat, Jacobian, finite differences, set, function
259005c665bSBarry Smith @*/
260840b8ebdSBarry Smith int MatFDColoringSetFunction(MatFDColoring matfd,int (*f)(void),void *fctx)
261005c665bSBarry Smith {
2623a40ed3dSBarry Smith   PetscFunctionBegin;
263005c665bSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
264005c665bSBarry Smith 
265005c665bSBarry Smith   matfd->f    = f;
266005c665bSBarry Smith   matfd->fctx = fctx;
267005c665bSBarry Smith 
2683a40ed3dSBarry Smith   PetscFunctionReturn(0);
269005c665bSBarry Smith }
270005c665bSBarry Smith 
271005c665bSBarry Smith #undef __FUNC__
272b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatFDColoringSetFromOptions"
273639f9d9dSBarry Smith /*@
274b4fc646aSLois Curfman McInnes    MatFDColoringSetFromOptions - Sets coloring finite difference parameters from
275639f9d9dSBarry Smith    the options database.
276639f9d9dSBarry Smith 
277fee21e36SBarry Smith    Collective on MatFDColoring
278fee21e36SBarry Smith 
27965f2ba5bSLois Curfman McInnes    The Jacobian, F'(u), is estimated with the differencing approximation
280ef5ee4d1SLois Curfman McInnes .vb
28165f2ba5bSLois Curfman McInnes        F'(u)_{:,i} = [F(u+h*dx_{i}) - F(u)]/h where
282f23b5b22SLois Curfman McInnes        h = error_rel*u[i]                 if  abs(u[i]) > umin
283f23b5b22SLois Curfman McInnes          = +/- error_rel*umin             otherwise, with +/- determined by the sign of u[i]
284ef5ee4d1SLois Curfman McInnes        dx_{i} = (0, ... 1, .... 0)
285ef5ee4d1SLois Curfman McInnes .ve
286ef5ee4d1SLois Curfman McInnes 
287ef5ee4d1SLois Curfman McInnes    Input Parameter:
288ef5ee4d1SLois Curfman McInnes .  coloring - the coloring context
289ef5ee4d1SLois Curfman McInnes 
290b4fc646aSLois Curfman McInnes    Options Database Keys:
291d15ffeeaSSatish Balay +  -mat_fd_coloring_err <err> - Sets <err> (square root
292ef5ee4d1SLois Curfman McInnes            of relative error in the function)
293f23b5b22SLois Curfman McInnes .  -mat_fd_coloring_umin <umin> - Sets umin, the minimum allowable u-value magnitude
294ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets frequency of computing a new Jacobian
295ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view - Activates basic viewing
296ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_view_info - Activates viewing info
297ef5ee4d1SLois Curfman McInnes -  -mat_fd_coloring_view_draw - Activates drawing
298639f9d9dSBarry Smith 
29915091d37SBarry Smith     Level: intermediate
30015091d37SBarry Smith 
301b4fc646aSLois Curfman McInnes .keywords: Mat, finite differences, parameters
302639f9d9dSBarry Smith @*/
303639f9d9dSBarry Smith int MatFDColoringSetFromOptions(MatFDColoring matfd)
304639f9d9dSBarry Smith {
305186905e3SBarry Smith   int        ierr;
3063a40ed3dSBarry Smith 
3073a40ed3dSBarry Smith   PetscFunctionBegin;
308639f9d9dSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
309639f9d9dSBarry Smith 
310cd5578b5SBarry Smith   ierr = OptionsBegin(matfd->comm,matfd->prefix,"Jacobian computation via finite differences option","MatFD");CHKERRQ(ierr);
311186905e3SBarry Smith     ierr = OptionsDouble("-mat_fd_coloring_err","Square root of relative error in function","MatFDColoringSetParameters",matfd->error_rel,&matfd->error_rel,0);CHKERRQ(ierr);
312186905e3SBarry Smith     ierr = OptionsDouble("-mat_fd_coloring_umin","Minimum allowable u magnitude","MatFDColoringSetParameters",matfd->umin,&matfd->umin,0);CHKERRQ(ierr);
313186905e3SBarry Smith     ierr = OptionsInt("-mat_fd_coloring_freq","How often Jacobian is recomputed","MatFDColoringSetFrequency",matfd->freq,&matfd->freq,0);CHKERRQ(ierr);
314186905e3SBarry Smith     /* not used here; but so they are presented in the GUI */
315186905e3SBarry Smith     ierr = OptionsName("-mat_fd_coloring_view","Print entire datastructure used for Jacobian","None",0);CHKERRQ(ierr);
316186905e3SBarry Smith     ierr = OptionsName("-mat_fd_coloring_view_info","Print number of colors etc for Jacobian","None",0);CHKERRQ(ierr);
317186905e3SBarry Smith     ierr = OptionsName("-mat_fd_coloring_view_draw","Plot nonzero structure ofJacobian","None",0);CHKERRQ(ierr);
318186905e3SBarry Smith   OptionsEnd();CHKERRQ(ierr);
3193a40ed3dSBarry Smith   PetscFunctionReturn(0);
320005c665bSBarry Smith }
321005c665bSBarry Smith 
322433994e6SBarry Smith #undef __FUNC__
323b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatFDColoringView_Private"
324005c665bSBarry Smith int MatFDColoringView_Private(MatFDColoring fd)
325005c665bSBarry Smith {
326f1af5d2fSBarry Smith   int        ierr;
327f1af5d2fSBarry Smith   PetscTruth flg;
328005c665bSBarry Smith 
3293a40ed3dSBarry Smith   PetscFunctionBegin;
330005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg);CHKERRQ(ierr);
331005c665bSBarry Smith   if (flg) {
332f8590f6eSBarry Smith     ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
333005c665bSBarry Smith   }
334ae09f205SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg);CHKERRQ(ierr);
335ae09f205SBarry Smith   if (flg) {
336f8590f6eSBarry Smith     ierr = ViewerPushFormat(VIEWER_STDOUT_(fd->comm),VIEWER_FORMAT_ASCII_INFO,PETSC_NULL);CHKERRQ(ierr);
337f8590f6eSBarry Smith     ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
338f8590f6eSBarry Smith     ierr = ViewerPopFormat(VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
339ae09f205SBarry Smith   }
340005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg);CHKERRQ(ierr);
341005c665bSBarry Smith   if (flg) {
342c655490fSBarry Smith     ierr = MatFDColoringView(fd,VIEWER_DRAW_(fd->comm));CHKERRQ(ierr);
343c655490fSBarry Smith     ierr = ViewerFlush(VIEWER_DRAW_(fd->comm));CHKERRQ(ierr);
344005c665bSBarry Smith   }
3453a40ed3dSBarry Smith   PetscFunctionReturn(0);
346bbf0e169SBarry Smith }
347bbf0e169SBarry Smith 
3485615d1e5SSatish Balay #undef __FUNC__
349b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatFDColoringCreate"
350bbf0e169SBarry Smith /*@C
351639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
352639f9d9dSBarry Smith    computation of Jacobians.
353bbf0e169SBarry Smith 
354ef5ee4d1SLois Curfman McInnes    Collective on Mat
355ef5ee4d1SLois Curfman McInnes 
356639f9d9dSBarry Smith    Input Parameters:
357ef5ee4d1SLois Curfman McInnes +  mat - the matrix containing the nonzero structure of the Jacobian
358ef5ee4d1SLois Curfman McInnes -  iscoloring - the coloring of the matrix
359bbf0e169SBarry Smith 
360bbf0e169SBarry Smith     Output Parameter:
361639f9d9dSBarry Smith .   color - the new coloring context
362bbf0e169SBarry Smith 
363b4fc646aSLois Curfman McInnes     Options Database Keys:
364ef5ee4d1SLois Curfman McInnes +    -mat_fd_coloring_view - Activates basic viewing or coloring
365ef5ee4d1SLois Curfman McInnes .    -mat_fd_coloring_view_draw - Activates drawing of coloring
366ef5ee4d1SLois Curfman McInnes -    -mat_fd_coloring_view_info - Activates viewing of coloring info
367639f9d9dSBarry Smith 
36815091d37SBarry Smith     Level: intermediate
36915091d37SBarry Smith 
3706831982aSBarry Smith .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(),
3716831982aSBarry Smith           MatFDColoringSetFunction(), MatFDColoringSetFromOptions()
372bbf0e169SBarry Smith @*/
373639f9d9dSBarry Smith int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
374bbf0e169SBarry Smith {
375639f9d9dSBarry Smith   MatFDColoring c;
376639f9d9dSBarry Smith   MPI_Comm      comm;
377639f9d9dSBarry Smith   int           ierr,M,N;
378639f9d9dSBarry Smith 
3793a40ed3dSBarry Smith   PetscFunctionBegin;
38081bfdfe8SSatish Balay   ierr = PLogEventBegin(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
381639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
382*29bbc08cSBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,"Only for square matrices");
383639f9d9dSBarry Smith 
384f881d145SBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
3859194eea9SBarry Smith   PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,"MatFDColoring",comm,MatFDColoringDestroy,MatFDColoringView);
386639f9d9dSBarry Smith   PLogObjectCreate(c);
387639f9d9dSBarry Smith 
388f830108cSBarry Smith   if (mat->ops->fdcoloringcreate) {
389f830108cSBarry Smith     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr);
390639f9d9dSBarry Smith   } else {
391*29bbc08cSBarry Smith     SETERRQ(PETSC_ERR_SUP,"Code not yet written for this matrix type");
392639f9d9dSBarry Smith   }
393639f9d9dSBarry Smith 
394639f9d9dSBarry Smith   c->error_rel = 1.e-8;
395ae09f205SBarry Smith   c->umin      = 1.e-6;
396005c665bSBarry Smith   c->freq      = 1;
397005c665bSBarry Smith 
398005c665bSBarry Smith   ierr = MatFDColoringView_Private(c);CHKERRQ(ierr);
399639f9d9dSBarry Smith 
400639f9d9dSBarry Smith   *color = c;
40181bfdfe8SSatish Balay   ierr = PLogEventEnd(MAT_FDColoringCreate,mat,0,0,0);CHKERRQ(ierr);
4023a40ed3dSBarry Smith   PetscFunctionReturn(0);
403bbf0e169SBarry Smith }
404bbf0e169SBarry Smith 
4055615d1e5SSatish Balay #undef __FUNC__
406b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatFDColoringDestroy"
407bbf0e169SBarry Smith /*@C
408639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
409639f9d9dSBarry Smith     via MatFDColoringCreate().
410bbf0e169SBarry Smith 
411ef5ee4d1SLois Curfman McInnes     Collective on MatFDColoring
412ef5ee4d1SLois Curfman McInnes 
413b4fc646aSLois Curfman McInnes     Input Parameter:
414639f9d9dSBarry Smith .   c - coloring context
415bbf0e169SBarry Smith 
41615091d37SBarry Smith     Level: intermediate
41715091d37SBarry Smith 
418639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
419bbf0e169SBarry Smith @*/
420639f9d9dSBarry Smith int MatFDColoringDestroy(MatFDColoring c)
421bbf0e169SBarry Smith {
422263760aaSBarry Smith   int i,ierr;
423bbf0e169SBarry Smith 
4243a40ed3dSBarry Smith   PetscFunctionBegin;
4253a40ed3dSBarry Smith   if (--c->refct > 0) PetscFunctionReturn(0);
426d4bb536fSBarry Smith 
427639f9d9dSBarry Smith   for (i=0; i<c->ncolors; i++) {
428606d414cSSatish Balay     if (c->columns[i])         {ierr = PetscFree(c->columns[i]);CHKERRQ(ierr);}
429606d414cSSatish Balay     if (c->rows[i])            {ierr = PetscFree(c->rows[i]);CHKERRQ(ierr);}
430606d414cSSatish Balay     if (c->columnsforrow[i])   {ierr = PetscFree(c->columnsforrow[i]);CHKERRQ(ierr);}
43130b34957SBarry Smith     if (c->vscaleforrow && c->vscaleforrow[i]) {ierr = PetscFree(c->vscaleforrow[i]);CHKERRQ(ierr);}
432bbf0e169SBarry Smith   }
433606d414cSSatish Balay   ierr = PetscFree(c->ncolumns);CHKERRQ(ierr);
434606d414cSSatish Balay   ierr = PetscFree(c->columns);CHKERRQ(ierr);
435606d414cSSatish Balay   ierr = PetscFree(c->nrows);CHKERRQ(ierr);
436606d414cSSatish Balay   ierr = PetscFree(c->rows);CHKERRQ(ierr);
437606d414cSSatish Balay   ierr = PetscFree(c->columnsforrow);CHKERRQ(ierr);
43830b34957SBarry Smith   if (c->vscaleforrow) {ierr = PetscFree(c->vscaleforrow);CHKERRQ(ierr);}
4394f113abfSBarry Smith   if (c->vscale)       {ierr = VecDestroy(c->vscale);CHKERRQ(ierr);}
440005c665bSBarry Smith   if (c->w1) {
441005c665bSBarry Smith     ierr = VecDestroy(c->w1);CHKERRQ(ierr);
442005c665bSBarry Smith     ierr = VecDestroy(c->w2);CHKERRQ(ierr);
443005c665bSBarry Smith     ierr = VecDestroy(c->w3);CHKERRQ(ierr);
444005c665bSBarry Smith   }
445639f9d9dSBarry Smith   PLogObjectDestroy(c);
446639f9d9dSBarry Smith   PetscHeaderDestroy(c);
4473a40ed3dSBarry Smith   PetscFunctionReturn(0);
448bbf0e169SBarry Smith }
44943a90d84SBarry Smith 
4505615d1e5SSatish Balay #undef __FUNC__
451b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatFDColoringApply"
45243a90d84SBarry Smith /*@
453e0907662SLois Curfman McInnes     MatFDColoringApply - Given a matrix for which a MatFDColoring context
454e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
45543a90d84SBarry Smith 
456fee21e36SBarry Smith     Collective on MatFDColoring
457fee21e36SBarry Smith 
458ef5ee4d1SLois Curfman McInnes     Input Parameters:
459ef5ee4d1SLois Curfman McInnes +   mat - location to store Jacobian
460ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
461ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
462ef5ee4d1SLois Curfman McInnes -   sctx - optional context required by function (actually a SNES context)
463ef5ee4d1SLois Curfman McInnes 
4648bba8e72SBarry Smith    Options Database Keys:
465ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
4668bba8e72SBarry Smith 
46715091d37SBarry Smith    Level: intermediate
46815091d37SBarry Smith 
46943a90d84SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
47043a90d84SBarry Smith 
47143a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
47243a90d84SBarry Smith @*/
473005c665bSBarry Smith int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
47443a90d84SBarry Smith {
4755904e1b1SBarry Smith   int           (*f)(void *,Vec,Vec,void*)= (int (*)(void *,Vec,Vec,void *))coloring->f;
4765904e1b1SBarry Smith   int           k,ierr,N,start,end,l,row,col,srow,**vscaleforrow;
4775904e1b1SBarry Smith   Scalar        dx,mone = -1.0,*y,*xx,*w3_array;
4784f113abfSBarry Smith   Scalar        *vscale_array;
479329f5518SBarry Smith   PetscReal     epsilon = coloring->error_rel,umin = coloring->umin;
480005c665bSBarry Smith   Vec           w1,w2,w3;
481005c665bSBarry Smith   void          *fctx = coloring->fctx;
482f1af5d2fSBarry Smith   PetscTruth    flg;
483005c665bSBarry Smith 
484a2e34c3dSBarry Smith 
4853a40ed3dSBarry Smith   PetscFunctionBegin;
486e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(J,MAT_COOKIE);
487e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE);
488e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(x1,VEC_COOKIE);
489e0907662SLois Curfman McInnes 
49081bfdfe8SSatish Balay   ierr = PLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
491005c665bSBarry Smith   if (!coloring->w1) {
492005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr);
493005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w1);
494005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr);
495005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w2);
496005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
497005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w3);
498005c665bSBarry Smith   }
499005c665bSBarry Smith   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
50043a90d84SBarry Smith 
5014c49b128SBarry Smith   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
502f1af5d2fSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
503f1af5d2fSBarry Smith   if (flg) {
504e0907662SLois Curfman McInnes     PLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n");
505e0907662SLois Curfman McInnes   } else {
50643a90d84SBarry Smith     ierr = MatZeroEntries(J);CHKERRQ(ierr);
507e0907662SLois Curfman McInnes   }
50843a90d84SBarry Smith 
50943a90d84SBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
51043a90d84SBarry Smith   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
51143a90d84SBarry Smith   ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr);
51243a90d84SBarry Smith 
51343a90d84SBarry Smith   /*
5149782cf98SBarry Smith       Compute all the scale factors and share with other processors
5159782cf98SBarry Smith   */
5169782cf98SBarry Smith   ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start;
5174f113abfSBarry Smith   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start;
5189782cf98SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
5199782cf98SBarry Smith     /*
5209782cf98SBarry Smith        Loop over each column associated with color adding the
5219782cf98SBarry Smith        perturbation to the vector w3.
5229782cf98SBarry Smith     */
5239782cf98SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
5249782cf98SBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
5259782cf98SBarry Smith       dx  = xx[col];
5269782cf98SBarry Smith       if (dx == 0.0) dx = 1.0;
5279782cf98SBarry Smith #if !defined(PETSC_USE_COMPLEX)
5289782cf98SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
5299782cf98SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
5309782cf98SBarry Smith #else
5319782cf98SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
5329782cf98SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
5339782cf98SBarry Smith #endif
5349782cf98SBarry Smith       dx                *= epsilon;
53530b34957SBarry Smith       vscale_array[col] = 1.0/dx;
5369782cf98SBarry Smith     }
5379782cf98SBarry Smith   }
538a2e34c3dSBarry Smith   vscale_array = vscale_array + start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
53930b34957SBarry Smith   ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
54030b34957SBarry Smith   ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5419782cf98SBarry Smith 
54230b34957SBarry Smith   if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
54330b34957SBarry Smith   else                        vscaleforrow = coloring->columnsforrow;
54430b34957SBarry Smith 
54530b34957SBarry Smith   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
5469782cf98SBarry Smith   /*
54743a90d84SBarry Smith       Loop over each color
54843a90d84SBarry Smith   */
54943a90d84SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
55043a90d84SBarry Smith     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
55142460c72SBarry Smith     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start;
55243a90d84SBarry Smith     /*
55343a90d84SBarry Smith        Loop over each column associated with color adding the
55443a90d84SBarry Smith        perturbation to the vector w3.
55543a90d84SBarry Smith     */
55643a90d84SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
55743a90d84SBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
55842460c72SBarry Smith       dx  = xx[col];
559ae09f205SBarry Smith       if (dx == 0.0) dx = 1.0;
560aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
561ae09f205SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
562ae09f205SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
56343a90d84SBarry Smith #else
564329f5518SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
565329f5518SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
56643a90d84SBarry Smith #endif
56743a90d84SBarry Smith       dx            *= epsilon;
568*29bbc08cSBarry Smith       if (!PetscAbsScalar(dx)) SETERRQ(1,"Computed 0 differencing parameter");
56942460c72SBarry Smith       w3_array[col] += dx;
57043a90d84SBarry Smith     }
57142460c72SBarry Smith     w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
5723b28642cSBarry Smith 
57343a90d84SBarry Smith     /*
574e0907662SLois Curfman McInnes        Evaluate function at x1 + dx (here dx is a vector of perturbations)
57543a90d84SBarry Smith     */
576a2e34c3dSBarry Smith 
57743a90d84SBarry Smith     ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
57843a90d84SBarry Smith     ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr);
5799782cf98SBarry Smith 
58043a90d84SBarry Smith     /*
581e0907662SLois Curfman McInnes        Loop over rows of vector, putting results into Jacobian matrix
58243a90d84SBarry Smith     */
5833b28642cSBarry Smith     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
58443a90d84SBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
58543a90d84SBarry Smith       row    = coloring->rows[k][l];
58643a90d84SBarry Smith       col    = coloring->columnsforrow[k][l];
5875904e1b1SBarry Smith       y[row] *= vscale_array[vscaleforrow[k][l]];
58843a90d84SBarry Smith       srow   = row + start;
58943a90d84SBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
59043a90d84SBarry Smith     }
5913b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
59243a90d84SBarry Smith   }
59330b34957SBarry Smith   ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
59442460c72SBarry Smith   xx = xx + start; ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
59543a90d84SBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
59643a90d84SBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
59781bfdfe8SSatish Balay   ierr = PLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
598a2e34c3dSBarry Smith 
599a2e34c3dSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_null_space_test",&flg);CHKERRQ(ierr);
600a2e34c3dSBarry Smith   if (flg) {
601a2e34c3dSBarry Smith     ierr = MatNullSpaceTest(J->nullsp,J);CHKERRQ(ierr);
602a2e34c3dSBarry Smith   }
6033a40ed3dSBarry Smith   PetscFunctionReturn(0);
60443a90d84SBarry Smith }
605840b8ebdSBarry Smith 
606840b8ebdSBarry Smith #undef __FUNC__
607b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatFDColoringApplyTS"
608840b8ebdSBarry Smith /*@
609840b8ebdSBarry Smith     MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context
610840b8ebdSBarry Smith     has been created, computes the Jacobian for a function via finite differences.
611840b8ebdSBarry Smith 
612fee21e36SBarry Smith    Collective on Mat, MatFDColoring, and Vec
613fee21e36SBarry Smith 
614ef5ee4d1SLois Curfman McInnes     Input Parameters:
6153b28642cSBarry Smith +   mat - location to store Jacobian
616ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
617ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
618ef5ee4d1SLois Curfman McInnes -   sctx - optional context required by function (actually a SNES context)
619ef5ee4d1SLois Curfman McInnes 
620840b8ebdSBarry Smith    Options Database Keys:
621ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
622840b8ebdSBarry Smith 
62315091d37SBarry Smith    Level: intermediate
62415091d37SBarry Smith 
625840b8ebdSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
626840b8ebdSBarry Smith 
627840b8ebdSBarry Smith .keywords: coloring, Jacobian, finite differences
628840b8ebdSBarry Smith @*/
629329f5518SBarry Smith int MatFDColoringApplyTS(Mat J,MatFDColoring coloring,PetscReal t,Vec x1,MatStructure *flag,void *sctx)
630840b8ebdSBarry Smith {
631329f5518SBarry Smith   int           (*f)(void *,PetscReal,Vec,Vec,void*)=(int (*)(void *,PetscReal,Vec,Vec,void *))coloring->f;
6325904e1b1SBarry Smith   int           k,ierr,N,start,end,l,row,col,srow,**vscaleforrow;
6335904e1b1SBarry Smith   Scalar        dx,mone = -1.0,*y,*xx,*w3_array;
6345904e1b1SBarry Smith   Scalar        *vscale_array;
635329f5518SBarry Smith   PetscReal     epsilon = coloring->error_rel,umin = coloring->umin;
636840b8ebdSBarry Smith   Vec           w1,w2,w3;
637840b8ebdSBarry Smith   void          *fctx = coloring->fctx;
638f1af5d2fSBarry Smith   PetscTruth    flg;
639840b8ebdSBarry Smith 
6403a40ed3dSBarry Smith   PetscFunctionBegin;
641840b8ebdSBarry Smith   PetscValidHeaderSpecific(J,MAT_COOKIE);
642840b8ebdSBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE);
643840b8ebdSBarry Smith   PetscValidHeaderSpecific(x1,VEC_COOKIE);
644840b8ebdSBarry Smith 
64581bfdfe8SSatish Balay   ierr = PLogEventBegin(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
646840b8ebdSBarry Smith   if (!coloring->w1) {
647840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr);
648840b8ebdSBarry Smith     PLogObjectParent(coloring,coloring->w1);
649840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr);
650840b8ebdSBarry Smith     PLogObjectParent(coloring,coloring->w2);
651840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
652840b8ebdSBarry Smith     PLogObjectParent(coloring,coloring->w3);
653840b8ebdSBarry Smith   }
654840b8ebdSBarry Smith   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
655840b8ebdSBarry Smith 
6565904e1b1SBarry Smith   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
657f1af5d2fSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
658f1af5d2fSBarry Smith   if (flg) {
6595904e1b1SBarry Smith     PLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n");
660840b8ebdSBarry Smith   } else {
661840b8ebdSBarry Smith     ierr = MatZeroEntries(J);CHKERRQ(ierr);
662840b8ebdSBarry Smith   }
663840b8ebdSBarry Smith 
664840b8ebdSBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
665840b8ebdSBarry Smith   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
666840b8ebdSBarry Smith   ierr = (*f)(sctx,t,x1,w1,fctx);CHKERRQ(ierr);
667840b8ebdSBarry Smith 
668840b8ebdSBarry Smith   /*
6695904e1b1SBarry Smith       Compute all the scale factors and share with other processors
670840b8ebdSBarry Smith   */
6715904e1b1SBarry Smith   ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start;
6725904e1b1SBarry Smith   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start;
673840b8ebdSBarry Smith   for (k=0; k<coloring->ncolors; k++) {
674840b8ebdSBarry Smith     /*
675840b8ebdSBarry Smith        Loop over each column associated with color adding the
676840b8ebdSBarry Smith        perturbation to the vector w3.
677840b8ebdSBarry Smith     */
678840b8ebdSBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
679840b8ebdSBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
6805904e1b1SBarry Smith       dx  = xx[col];
681840b8ebdSBarry Smith       if (dx == 0.0) dx = 1.0;
682aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
683840b8ebdSBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
684840b8ebdSBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
685840b8ebdSBarry Smith #else
686329f5518SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
687329f5518SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
688840b8ebdSBarry Smith #endif
689840b8ebdSBarry Smith       dx                *= epsilon;
6905904e1b1SBarry Smith       vscale_array[col] = 1.0/dx;
691840b8ebdSBarry Smith     }
6925904e1b1SBarry Smith   }
6935904e1b1SBarry Smith   vscale_array = vscale_array - start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
6945904e1b1SBarry Smith   ierr = VecGhostUpdateBegin(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6955904e1b1SBarry Smith   ierr = VecGhostUpdateEnd(coloring->vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6965904e1b1SBarry Smith 
6975904e1b1SBarry Smith   if (coloring->vscaleforrow) vscaleforrow = coloring->vscaleforrow;
6985904e1b1SBarry Smith   else                        vscaleforrow = coloring->columnsforrow;
6995904e1b1SBarry Smith 
7005904e1b1SBarry Smith   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
7015904e1b1SBarry Smith   /*
7025904e1b1SBarry Smith       Loop over each color
7035904e1b1SBarry Smith   */
7045904e1b1SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
7055904e1b1SBarry Smith     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
7065904e1b1SBarry Smith     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start;
7075904e1b1SBarry Smith     /*
7085904e1b1SBarry Smith        Loop over each column associated with color adding the
7095904e1b1SBarry Smith        perturbation to the vector w3.
7105904e1b1SBarry Smith     */
7115904e1b1SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
7125904e1b1SBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
7135904e1b1SBarry Smith       dx  = xx[col];
7145904e1b1SBarry Smith       if (dx == 0.0) dx = 1.0;
7155904e1b1SBarry Smith #if !defined(PETSC_USE_COMPLEX)
7165904e1b1SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
7175904e1b1SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
7185904e1b1SBarry Smith #else
7195904e1b1SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
7205904e1b1SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
7215904e1b1SBarry Smith #endif
7225904e1b1SBarry Smith       dx            *= epsilon;
7235904e1b1SBarry Smith       w3_array[col] += dx;
7245904e1b1SBarry Smith     }
7255904e1b1SBarry Smith     w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
7265904e1b1SBarry Smith 
727840b8ebdSBarry Smith     /*
728840b8ebdSBarry Smith        Evaluate function at x1 + dx (here dx is a vector of perturbations)
729840b8ebdSBarry Smith     */
730840b8ebdSBarry Smith     ierr = (*f)(sctx,t,w3,w2,fctx);CHKERRQ(ierr);
731840b8ebdSBarry Smith     ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr);
7325904e1b1SBarry Smith 
733840b8ebdSBarry Smith     /*
734840b8ebdSBarry Smith        Loop over rows of vector, putting results into Jacobian matrix
735840b8ebdSBarry Smith     */
7363b28642cSBarry Smith     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
737840b8ebdSBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
738840b8ebdSBarry Smith       row    = coloring->rows[k][l];
739840b8ebdSBarry Smith       col    = coloring->columnsforrow[k][l];
7405904e1b1SBarry Smith       y[row] *= vscale_array[vscaleforrow[k][l]];
741840b8ebdSBarry Smith       srow   = row + start;
742840b8ebdSBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
743840b8ebdSBarry Smith     }
7443b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
745840b8ebdSBarry Smith   }
7465904e1b1SBarry Smith   ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
7475904e1b1SBarry Smith   xx = xx + start; ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
748840b8ebdSBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
74981bfdfe8SSatish Balay   ierr = PLogEventEnd(MAT_FDColoringApply,coloring,J,x1,0);CHKERRQ(ierr);
750840b8ebdSBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7513a40ed3dSBarry Smith   PetscFunctionReturn(0);
752840b8ebdSBarry Smith }
7533b28642cSBarry Smith 
7543b28642cSBarry Smith 
7553b28642cSBarry Smith 
756