xref: /petsc/src/mat/matfd/fdmatrix.c (revision 4f113abfafedfdd5469bb1b8415f6989092d43ce)
1*4f113abfSBarry Smith /*$Id: fdmatrix.c,v 1.65 2000/05/14 03:57:55 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 {
1230f5bd95cSBarry Smith     SETERRQ1(1,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:
291ef5ee4d1SLois Curfman McInnes +  -mat_fd_coloring_error <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 {
305f1af5d2fSBarry Smith   int        ierr,freq = 1;
306329f5518SBarry Smith   PetscReal  error = PETSC_DEFAULT,umin = PETSC_DEFAULT;
307f1af5d2fSBarry Smith   PetscTruth flag;
3083a40ed3dSBarry Smith 
3093a40ed3dSBarry Smith   PetscFunctionBegin;
310639f9d9dSBarry Smith   PetscValidHeaderSpecific(matfd,MAT_FDCOLORING_COOKIE);
311639f9d9dSBarry Smith 
312f1af5d2fSBarry Smith   ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_err",&error,PETSC_NULL);CHKERRQ(ierr);
313f1af5d2fSBarry Smith   ierr = OptionsGetDouble(matfd->prefix,"-mat_fd_coloring_umin",&umin,PETSC_NULL);CHKERRQ(ierr);
314639f9d9dSBarry Smith   ierr = MatFDColoringSetParameters(matfd,error,umin);CHKERRQ(ierr);
315f1af5d2fSBarry Smith   ierr = OptionsGetInt(matfd->prefix,"-mat_fd_coloring_freq",&freq,PETSC_NULL);CHKERRQ(ierr);
316005c665bSBarry Smith   ierr = MatFDColoringSetFrequency(matfd,freq);CHKERRQ(ierr);
317005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flag);CHKERRQ(ierr);
318639f9d9dSBarry Smith   if (flag) {
319639f9d9dSBarry Smith     ierr = MatFDColoringPrintHelp(matfd);CHKERRQ(ierr);
320639f9d9dSBarry Smith   }
3213a40ed3dSBarry Smith   PetscFunctionReturn(0);
322639f9d9dSBarry Smith }
323639f9d9dSBarry Smith 
3245615d1e5SSatish Balay #undef __FUNC__
325b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatFDColoringPrintHelp"
326639f9d9dSBarry Smith /*@
327639f9d9dSBarry Smith     MatFDColoringPrintHelp - Prints help message for matrix finite difference calculations
328639f9d9dSBarry Smith     using coloring.
329639f9d9dSBarry Smith 
330ef5ee4d1SLois Curfman McInnes     Collective on MatFDColoring
331ef5ee4d1SLois Curfman McInnes 
332639f9d9dSBarry Smith     Input Parameter:
333639f9d9dSBarry Smith .   fdcoloring - the MatFDColoring context
334639f9d9dSBarry Smith 
33515091d37SBarry Smith     Level: intermediate
33615091d37SBarry Smith 
337639f9d9dSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringSetFromOptions()
338639f9d9dSBarry Smith @*/
339639f9d9dSBarry Smith int MatFDColoringPrintHelp(MatFDColoring fd)
340639f9d9dSBarry Smith {
341d132466eSBarry Smith   int ierr;
342d132466eSBarry Smith 
3433a40ed3dSBarry Smith   PetscFunctionBegin;
344639f9d9dSBarry Smith   PetscValidHeaderSpecific(fd,MAT_FDCOLORING_COOKIE);
345639f9d9dSBarry Smith 
346d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_err <err>: set sqrt rel tol in function, defaults to %g\n",fd->error_rel);CHKERRQ(ierr);
347d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_umin <umin>: see users manual, defaults to %d\n",fd->umin);CHKERRQ(ierr);
348d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_freq <freq>: frequency that Jacobian is recomputed, defaults to %d\n",fd->freq);CHKERRQ(ierr);
349d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_view\n");CHKERRQ(ierr);
350d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_view_draw\n");CHKERRQ(ierr);
351d132466eSBarry Smith   ierr = (*PetscHelpPrintf)(fd->comm,"-mat_fd_coloring_view_info\n");CHKERRQ(ierr);
3523a40ed3dSBarry Smith   PetscFunctionReturn(0);
353005c665bSBarry Smith }
354005c665bSBarry Smith 
355433994e6SBarry Smith #undef __FUNC__
356b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatFDColoringView_Private"
357005c665bSBarry Smith int MatFDColoringView_Private(MatFDColoring fd)
358005c665bSBarry Smith {
359f1af5d2fSBarry Smith   int        ierr;
360f1af5d2fSBarry Smith   PetscTruth flg;
361005c665bSBarry Smith 
3623a40ed3dSBarry Smith   PetscFunctionBegin;
363005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view",&flg);CHKERRQ(ierr);
364005c665bSBarry Smith   if (flg) {
365f8590f6eSBarry Smith     ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
366005c665bSBarry Smith   }
367ae09f205SBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_info",&flg);CHKERRQ(ierr);
368ae09f205SBarry Smith   if (flg) {
369f8590f6eSBarry Smith     ierr = ViewerPushFormat(VIEWER_STDOUT_(fd->comm),VIEWER_FORMAT_ASCII_INFO,PETSC_NULL);CHKERRQ(ierr);
370f8590f6eSBarry Smith     ierr = MatFDColoringView(fd,VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
371f8590f6eSBarry Smith     ierr = ViewerPopFormat(VIEWER_STDOUT_(fd->comm));CHKERRQ(ierr);
372ae09f205SBarry Smith   }
373005c665bSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_view_draw",&flg);CHKERRQ(ierr);
374005c665bSBarry Smith   if (flg) {
375c655490fSBarry Smith     ierr = MatFDColoringView(fd,VIEWER_DRAW_(fd->comm));CHKERRQ(ierr);
376c655490fSBarry Smith     ierr = ViewerFlush(VIEWER_DRAW_(fd->comm));CHKERRQ(ierr);
377005c665bSBarry Smith   }
3783a40ed3dSBarry Smith   PetscFunctionReturn(0);
379bbf0e169SBarry Smith }
380bbf0e169SBarry Smith 
3815615d1e5SSatish Balay #undef __FUNC__
382b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatFDColoringCreate"
383bbf0e169SBarry Smith /*@C
384639f9d9dSBarry Smith    MatFDColoringCreate - Creates a matrix coloring context for finite difference
385639f9d9dSBarry Smith    computation of Jacobians.
386bbf0e169SBarry Smith 
387ef5ee4d1SLois Curfman McInnes    Collective on Mat
388ef5ee4d1SLois Curfman McInnes 
389639f9d9dSBarry Smith    Input Parameters:
390ef5ee4d1SLois Curfman McInnes +  mat - the matrix containing the nonzero structure of the Jacobian
391ef5ee4d1SLois Curfman McInnes -  iscoloring - the coloring of the matrix
392bbf0e169SBarry Smith 
393bbf0e169SBarry Smith     Output Parameter:
394639f9d9dSBarry Smith .   color - the new coloring context
395bbf0e169SBarry Smith 
396b4fc646aSLois Curfman McInnes     Options Database Keys:
397ef5ee4d1SLois Curfman McInnes +    -mat_fd_coloring_view - Activates basic viewing or coloring
398ef5ee4d1SLois Curfman McInnes .    -mat_fd_coloring_view_draw - Activates drawing of coloring
399ef5ee4d1SLois Curfman McInnes -    -mat_fd_coloring_view_info - Activates viewing of coloring info
400639f9d9dSBarry Smith 
40115091d37SBarry Smith     Level: intermediate
40215091d37SBarry Smith 
4036831982aSBarry Smith .seealso: MatFDColoringDestroy(),SNESDefaultComputeJacobianColor(), ISColoringCreate(),
4046831982aSBarry Smith           MatFDColoringSetFunction(), MatFDColoringSetFromOptions()
405bbf0e169SBarry Smith @*/
406639f9d9dSBarry Smith int MatFDColoringCreate(Mat mat,ISColoring iscoloring,MatFDColoring *color)
407bbf0e169SBarry Smith {
408639f9d9dSBarry Smith   MatFDColoring c;
409639f9d9dSBarry Smith   MPI_Comm      comm;
410639f9d9dSBarry Smith   int           ierr,M,N;
411639f9d9dSBarry Smith 
4123a40ed3dSBarry Smith   PetscFunctionBegin;
413639f9d9dSBarry Smith   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
414e3372554SBarry Smith   if (M != N) SETERRQ(PETSC_ERR_SUP,0,"Only for square matrices");
415639f9d9dSBarry Smith 
416f881d145SBarry Smith   ierr = PetscObjectGetComm((PetscObject)mat,&comm);CHKERRQ(ierr);
4179194eea9SBarry Smith   PetscHeaderCreate(c,_p_MatFDColoring,int,MAT_FDCOLORING_COOKIE,0,"MatFDColoring",comm,MatFDColoringDestroy,MatFDColoringView);
418639f9d9dSBarry Smith   PLogObjectCreate(c);
419639f9d9dSBarry Smith 
420f830108cSBarry Smith   if (mat->ops->fdcoloringcreate) {
421f830108cSBarry Smith     ierr = (*mat->ops->fdcoloringcreate)(mat,iscoloring,c);CHKERRQ(ierr);
422639f9d9dSBarry Smith   } else {
423e3372554SBarry Smith     SETERRQ(PETSC_ERR_SUP,0,"Code not yet written for this matrix type");
424639f9d9dSBarry Smith   }
425639f9d9dSBarry Smith 
426639f9d9dSBarry Smith   c->error_rel = 1.e-8;
427ae09f205SBarry Smith   c->umin      = 1.e-6;
428005c665bSBarry Smith   c->freq      = 1;
429005c665bSBarry Smith 
430005c665bSBarry Smith   ierr = MatFDColoringView_Private(c);CHKERRQ(ierr);
431639f9d9dSBarry Smith 
432639f9d9dSBarry Smith   *color = c;
433639f9d9dSBarry Smith 
4343a40ed3dSBarry Smith   PetscFunctionReturn(0);
435bbf0e169SBarry Smith }
436bbf0e169SBarry Smith 
4375615d1e5SSatish Balay #undef __FUNC__
438b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatFDColoringDestroy"
439bbf0e169SBarry Smith /*@C
440639f9d9dSBarry Smith     MatFDColoringDestroy - Destroys a matrix coloring context that was created
441639f9d9dSBarry Smith     via MatFDColoringCreate().
442bbf0e169SBarry Smith 
443ef5ee4d1SLois Curfman McInnes     Collective on MatFDColoring
444ef5ee4d1SLois Curfman McInnes 
445b4fc646aSLois Curfman McInnes     Input Parameter:
446639f9d9dSBarry Smith .   c - coloring context
447bbf0e169SBarry Smith 
44815091d37SBarry Smith     Level: intermediate
44915091d37SBarry Smith 
450639f9d9dSBarry Smith .seealso: MatFDColoringCreate()
451bbf0e169SBarry Smith @*/
452639f9d9dSBarry Smith int MatFDColoringDestroy(MatFDColoring c)
453bbf0e169SBarry Smith {
454263760aaSBarry Smith   int i,ierr;
455bbf0e169SBarry Smith 
4563a40ed3dSBarry Smith   PetscFunctionBegin;
4573a40ed3dSBarry Smith   if (--c->refct > 0) PetscFunctionReturn(0);
458d4bb536fSBarry Smith 
459639f9d9dSBarry Smith 
460639f9d9dSBarry Smith   for (i=0; i<c->ncolors; i++) {
461606d414cSSatish Balay     if (c->columns[i])         {ierr = PetscFree(c->columns[i]);CHKERRQ(ierr);}
462606d414cSSatish Balay     if (c->rows[i])            {ierr = PetscFree(c->rows[i]);CHKERRQ(ierr);}
463606d414cSSatish Balay     if (c->columnsforrow[i])   {ierr = PetscFree(c->columnsforrow[i]);CHKERRQ(ierr);}
464*4f113abfSBarry Smith     if (c->vscaleforrow[i])    {ierr = PetscFree(c->vscaleforrow[i]);CHKERRQ(ierr);}
465bbf0e169SBarry Smith   }
466606d414cSSatish Balay   ierr = PetscFree(c->ncolumns);CHKERRQ(ierr);
467606d414cSSatish Balay   ierr = PetscFree(c->columns);CHKERRQ(ierr);
468606d414cSSatish Balay   ierr = PetscFree(c->nrows);CHKERRQ(ierr);
469606d414cSSatish Balay   ierr = PetscFree(c->rows);CHKERRQ(ierr);
470606d414cSSatish Balay   ierr = PetscFree(c->columnsforrow);CHKERRQ(ierr);
471606d414cSSatish Balay   ierr = PetscFree(c->scale);CHKERRQ(ierr);
472*4f113abfSBarry Smith   ierr = PetscFree(c->vscaleforrow);CHKERRQ(ierr);
473*4f113abfSBarry Smith   if (c->vscale) {ierr = VecDestroy(c->vscale);CHKERRQ(ierr);}
474005c665bSBarry Smith   if (c->w1) {
475005c665bSBarry Smith     ierr = VecDestroy(c->w1);CHKERRQ(ierr);
476005c665bSBarry Smith     ierr = VecDestroy(c->w2);CHKERRQ(ierr);
477005c665bSBarry Smith     ierr = VecDestroy(c->w3);CHKERRQ(ierr);
478005c665bSBarry Smith   }
479639f9d9dSBarry Smith   PLogObjectDestroy(c);
480639f9d9dSBarry Smith   PetscHeaderDestroy(c);
4813a40ed3dSBarry Smith   PetscFunctionReturn(0);
482bbf0e169SBarry Smith }
48343a90d84SBarry Smith 
484005c665bSBarry Smith 
4855615d1e5SSatish Balay #undef __FUNC__
486b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatFDColoringApply"
48743a90d84SBarry Smith /*@
488e0907662SLois Curfman McInnes     MatFDColoringApply - Given a matrix for which a MatFDColoring context
489e0907662SLois Curfman McInnes     has been created, computes the Jacobian for a function via finite differences.
49043a90d84SBarry Smith 
491fee21e36SBarry Smith     Collective on MatFDColoring
492fee21e36SBarry Smith 
493ef5ee4d1SLois Curfman McInnes     Input Parameters:
494ef5ee4d1SLois Curfman McInnes +   mat - location to store Jacobian
495ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
496ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
497ef5ee4d1SLois Curfman McInnes -   sctx - optional context required by function (actually a SNES context)
498ef5ee4d1SLois Curfman McInnes 
4998bba8e72SBarry Smith    Options Database Keys:
500ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
5018bba8e72SBarry Smith 
50215091d37SBarry Smith    Level: intermediate
50315091d37SBarry Smith 
50443a90d84SBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
50543a90d84SBarry Smith 
50643a90d84SBarry Smith .keywords: coloring, Jacobian, finite differences
50743a90d84SBarry Smith @*/
508005c665bSBarry Smith int MatFDColoringApply(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx)
50943a90d84SBarry Smith {
510f1af5d2fSBarry Smith   int           k,ierr,N,start,end,l,row,col,srow;
5119194eea9SBarry Smith   Scalar        dx,mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale,*w3_array;
512*4f113abfSBarry Smith   Scalar        *vscale_array;
513329f5518SBarry Smith   PetscReal     epsilon = coloring->error_rel,umin = coloring->umin;
51443a90d84SBarry Smith   MPI_Comm      comm = coloring->comm;
515005c665bSBarry Smith   Vec           w1,w2,w3;
516840b8ebdSBarry Smith   int           (*f)(void *,Vec,Vec,void*)= (int (*)(void *,Vec,Vec,void *))coloring->f;
517005c665bSBarry Smith   void          *fctx = coloring->fctx;
518f1af5d2fSBarry Smith   PetscTruth    flg;
519005c665bSBarry Smith 
5203a40ed3dSBarry Smith   PetscFunctionBegin;
521e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(J,MAT_COOKIE);
522e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE);
523e0907662SLois Curfman McInnes   PetscValidHeaderSpecific(x1,VEC_COOKIE);
524e0907662SLois Curfman McInnes 
525005c665bSBarry Smith   if (!coloring->w1) {
526005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr);
527005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w1);
528005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr);
529005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w2);
530005c665bSBarry Smith     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
531005c665bSBarry Smith     PLogObjectParent(coloring,coloring->w3);
532005c665bSBarry Smith   }
533005c665bSBarry Smith   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
53443a90d84SBarry Smith 
5354c49b128SBarry Smith   ierr = MatSetUnfactored(J);CHKERRQ(ierr);
536f1af5d2fSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
537f1af5d2fSBarry Smith   if (flg) {
538e0907662SLois Curfman McInnes     PLogInfo(coloring,"MatFDColoringApply: Not calling MatZeroEntries()\n");
539e0907662SLois Curfman McInnes   } else {
54043a90d84SBarry Smith     ierr = MatZeroEntries(J);CHKERRQ(ierr);
541e0907662SLois Curfman McInnes   }
54243a90d84SBarry Smith 
54343a90d84SBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
54443a90d84SBarry Smith   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
54543a90d84SBarry Smith   ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr);
54643a90d84SBarry Smith 
547549d3d68SSatish Balay   ierr = PetscMemzero(wscale,N*sizeof(Scalar));CHKERRQ(ierr);
54843a90d84SBarry Smith   /*
5499782cf98SBarry Smith       Compute all the scale factors and share with other processors
5509782cf98SBarry Smith   */
5519782cf98SBarry Smith   ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);xx = xx - start;
552*4f113abfSBarry Smith   ierr = VecGetArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);vscale_array = vscale_array - start;
5539782cf98SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
5549782cf98SBarry Smith     /*
5559782cf98SBarry Smith        Loop over each column associated with color adding the
5569782cf98SBarry Smith        perturbation to the vector w3.
5579782cf98SBarry Smith     */
5589782cf98SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
5599782cf98SBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
5609782cf98SBarry Smith       dx  = xx[col];
5619782cf98SBarry Smith       if (dx == 0.0) dx = 1.0;
5629782cf98SBarry Smith #if !defined(PETSC_USE_COMPLEX)
5639782cf98SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
5649782cf98SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
5659782cf98SBarry Smith #else
5669782cf98SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
5679782cf98SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
5689782cf98SBarry Smith #endif
5699782cf98SBarry Smith       dx            *= epsilon;
5709782cf98SBarry Smith       wscale[col]   = 1.0/dx;
5719782cf98SBarry Smith     }
5729782cf98SBarry Smith   }
573*4f113abfSBarry Smith   vscale_array = vscale_array - start;ierr = VecRestoreArray(coloring->vscale,&vscale_array);CHKERRQ(ierr);
5749782cf98SBarry Smith   ierr = MPI_Allreduce(wscale,scale,N,MPIU_SCALAR,PetscSum_Op,comm);CHKERRQ(ierr);
5759782cf98SBarry Smith 
5769782cf98SBarry Smith   /*
57743a90d84SBarry Smith       Loop over each color
57843a90d84SBarry Smith   */
57943a90d84SBarry Smith   for (k=0; k<coloring->ncolors; k++) {
58043a90d84SBarry Smith     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
58142460c72SBarry Smith     ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr);w3_array = w3_array - start;
58243a90d84SBarry Smith     /*
58343a90d84SBarry Smith        Loop over each column associated with color adding the
58443a90d84SBarry Smith        perturbation to the vector w3.
58543a90d84SBarry Smith     */
58643a90d84SBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
58743a90d84SBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
58842460c72SBarry Smith       dx  = xx[col];
589ae09f205SBarry Smith       if (dx == 0.0) dx = 1.0;
590aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
591ae09f205SBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
592ae09f205SBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
59343a90d84SBarry Smith #else
594329f5518SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
595329f5518SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
59643a90d84SBarry Smith #endif
59743a90d84SBarry Smith       dx            *= epsilon;
59842460c72SBarry Smith       w3_array[col] += dx;
59943a90d84SBarry Smith     }
60042460c72SBarry Smith     w3_array = w3_array + start; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr);
6013b28642cSBarry Smith 
60243a90d84SBarry Smith     /*
603e0907662SLois Curfman McInnes        Evaluate function at x1 + dx (here dx is a vector of perturbations)
60443a90d84SBarry Smith     */
60543a90d84SBarry Smith     ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr);
60643a90d84SBarry Smith     ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr);
6079782cf98SBarry Smith 
60843a90d84SBarry Smith     /*
609e0907662SLois Curfman McInnes        Loop over rows of vector, putting results into Jacobian matrix
61043a90d84SBarry Smith     */
6113b28642cSBarry Smith     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
61243a90d84SBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
61343a90d84SBarry Smith       row    = coloring->rows[k][l];
61443a90d84SBarry Smith       col    = coloring->columnsforrow[k][l];
61543a90d84SBarry Smith       y[row] *= scale[col];
61643a90d84SBarry Smith       srow   = row + start;
61743a90d84SBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
61843a90d84SBarry Smith     }
6193b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
62043a90d84SBarry Smith   }
62142460c72SBarry Smith   xx = xx + start; ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
62243a90d84SBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
62343a90d84SBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
6243a40ed3dSBarry Smith   PetscFunctionReturn(0);
62543a90d84SBarry Smith }
626840b8ebdSBarry Smith 
627840b8ebdSBarry Smith #undef __FUNC__
628b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"MatFDColoringApplyTS"
629840b8ebdSBarry Smith /*@
630840b8ebdSBarry Smith     MatFDColoringApplyTS - Given a matrix for which a MatFDColoring context
631840b8ebdSBarry Smith     has been created, computes the Jacobian for a function via finite differences.
632840b8ebdSBarry Smith 
633fee21e36SBarry Smith    Collective on Mat, MatFDColoring, and Vec
634fee21e36SBarry Smith 
635ef5ee4d1SLois Curfman McInnes     Input Parameters:
6363b28642cSBarry Smith +   mat - location to store Jacobian
637ef5ee4d1SLois Curfman McInnes .   coloring - coloring context created with MatFDColoringCreate()
638ef5ee4d1SLois Curfman McInnes .   x1 - location at which Jacobian is to be computed
639ef5ee4d1SLois Curfman McInnes -   sctx - optional context required by function (actually a SNES context)
640ef5ee4d1SLois Curfman McInnes 
641840b8ebdSBarry Smith    Options Database Keys:
642ef5ee4d1SLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Sets coloring frequency
643840b8ebdSBarry Smith 
64415091d37SBarry Smith    Level: intermediate
64515091d37SBarry Smith 
646840b8ebdSBarry Smith .seealso: MatFDColoringCreate(), MatFDColoringDestroy(), MatFDColoringView()
647840b8ebdSBarry Smith 
648840b8ebdSBarry Smith .keywords: coloring, Jacobian, finite differences
649840b8ebdSBarry Smith @*/
650329f5518SBarry Smith int MatFDColoringApplyTS(Mat J,MatFDColoring coloring,PetscReal t,Vec x1,MatStructure *flag,void *sctx)
651840b8ebdSBarry Smith {
652f1af5d2fSBarry Smith   int           k,ierr,N,start,end,l,row,col,srow;
653329f5518SBarry Smith   int           (*f)(void *,PetscReal,Vec,Vec,void*)= (int (*)(void *,PetscReal,Vec,Vec,void *))coloring->f;
654840b8ebdSBarry Smith   Scalar        dx,mone = -1.0,*y,*scale = coloring->scale,*xx,*wscale = coloring->wscale;
655329f5518SBarry Smith   PetscReal     epsilon = coloring->error_rel,umin = coloring->umin;
656840b8ebdSBarry Smith   MPI_Comm      comm = coloring->comm;
657840b8ebdSBarry Smith   Vec           w1,w2,w3;
658840b8ebdSBarry Smith   void          *fctx = coloring->fctx;
659f1af5d2fSBarry Smith   PetscTruth    flg;
660840b8ebdSBarry Smith 
6613a40ed3dSBarry Smith   PetscFunctionBegin;
662840b8ebdSBarry Smith   PetscValidHeaderSpecific(J,MAT_COOKIE);
663840b8ebdSBarry Smith   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_COOKIE);
664840b8ebdSBarry Smith   PetscValidHeaderSpecific(x1,VEC_COOKIE);
665840b8ebdSBarry Smith 
666840b8ebdSBarry Smith   if (!coloring->w1) {
667840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w1);CHKERRQ(ierr);
668840b8ebdSBarry Smith     PLogObjectParent(coloring,coloring->w1);
669840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w2);CHKERRQ(ierr);
670840b8ebdSBarry Smith     PLogObjectParent(coloring,coloring->w2);
671840b8ebdSBarry Smith     ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr);
672840b8ebdSBarry Smith     PLogObjectParent(coloring,coloring->w3);
673840b8ebdSBarry Smith   }
674840b8ebdSBarry Smith   w1 = coloring->w1; w2 = coloring->w2; w3 = coloring->w3;
675840b8ebdSBarry Smith 
676f1af5d2fSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-mat_fd_coloring_dont_rezero",&flg);CHKERRQ(ierr);
677f1af5d2fSBarry Smith   if (flg) {
678840b8ebdSBarry Smith     PLogInfo(coloring,"MatFDColoringApplyTS: Not calling MatZeroEntries()\n");
679840b8ebdSBarry Smith   } else {
680840b8ebdSBarry Smith     ierr = MatZeroEntries(J);CHKERRQ(ierr);
681840b8ebdSBarry Smith   }
682840b8ebdSBarry Smith 
683840b8ebdSBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
684840b8ebdSBarry Smith   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
685840b8ebdSBarry Smith   ierr = (*f)(sctx,t,x1,w1,fctx);CHKERRQ(ierr);
686840b8ebdSBarry Smith 
687549d3d68SSatish Balay   ierr = PetscMemzero(wscale,N*sizeof(Scalar));CHKERRQ(ierr);
688840b8ebdSBarry Smith   /*
689840b8ebdSBarry Smith       Loop over each color
690840b8ebdSBarry Smith   */
691840b8ebdSBarry Smith 
6923b28642cSBarry Smith   ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);
693840b8ebdSBarry Smith   for (k=0; k<coloring->ncolors; k++) {
694840b8ebdSBarry Smith     ierr = VecCopy(x1,w3);CHKERRQ(ierr);
695840b8ebdSBarry Smith     /*
696840b8ebdSBarry Smith        Loop over each column associated with color adding the
697840b8ebdSBarry Smith        perturbation to the vector w3.
698840b8ebdSBarry Smith     */
699840b8ebdSBarry Smith     for (l=0; l<coloring->ncolumns[k]; l++) {
700840b8ebdSBarry Smith       col = coloring->columns[k][l];    /* column of the matrix we are probing for */
701840b8ebdSBarry Smith       dx  = xx[col-start];
702840b8ebdSBarry Smith       if (dx == 0.0) dx = 1.0;
703aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
704840b8ebdSBarry Smith       if (dx < umin && dx >= 0.0)      dx = umin;
705840b8ebdSBarry Smith       else if (dx < 0.0 && dx > -umin) dx = -umin;
706840b8ebdSBarry Smith #else
707329f5518SBarry Smith       if (PetscAbsScalar(dx) < umin && PetscRealPart(dx) >= 0.0)     dx = umin;
708329f5518SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < umin) dx = -umin;
709840b8ebdSBarry Smith #endif
710840b8ebdSBarry Smith       dx          *= epsilon;
711840b8ebdSBarry Smith       wscale[col] = 1.0/dx;
7123b28642cSBarry Smith       ierr = VecSetValues(w3,1,&col,&dx,ADD_VALUES);CHKERRQ(ierr);
713840b8ebdSBarry Smith     }
714840b8ebdSBarry Smith     /*
715840b8ebdSBarry Smith        Evaluate function at x1 + dx (here dx is a vector of perturbations)
716840b8ebdSBarry Smith     */
717840b8ebdSBarry Smith     ierr = (*f)(sctx,t,w3,w2,fctx);CHKERRQ(ierr);
718840b8ebdSBarry Smith     ierr = VecAXPY(&mone,w1,w2);CHKERRQ(ierr);
719840b8ebdSBarry Smith     /* Communicate scale to all processors */
7206831982aSBarry Smith     ierr = MPI_Allreduce(wscale,scale,N,MPIU_SCALAR,PetscSum_Op,comm);CHKERRQ(ierr);
721840b8ebdSBarry Smith     /*
722840b8ebdSBarry Smith        Loop over rows of vector, putting results into Jacobian matrix
723840b8ebdSBarry Smith     */
7243b28642cSBarry Smith     ierr = VecGetArray(w2,&y);CHKERRQ(ierr);
725840b8ebdSBarry Smith     for (l=0; l<coloring->nrows[k]; l++) {
726840b8ebdSBarry Smith       row    = coloring->rows[k][l];
727840b8ebdSBarry Smith       col    = coloring->columnsforrow[k][l];
728840b8ebdSBarry Smith       y[row] *= scale[col];
729840b8ebdSBarry Smith       srow   = row + start;
730840b8ebdSBarry Smith       ierr   = MatSetValues(J,1,&srow,1,&col,y+row,INSERT_VALUES);CHKERRQ(ierr);
731840b8ebdSBarry Smith     }
7323b28642cSBarry Smith     ierr = VecRestoreArray(w2,&y);CHKERRQ(ierr);
733840b8ebdSBarry Smith   }
7343b28642cSBarry Smith   ierr  = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
735840b8ebdSBarry Smith   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
736840b8ebdSBarry Smith   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
7373a40ed3dSBarry Smith   PetscFunctionReturn(0);
738840b8ebdSBarry Smith }
7393b28642cSBarry Smith 
7403b28642cSBarry Smith 
7413b28642cSBarry Smith 
742