xref: /petsc/src/snes/interface/snesj.c (revision d142edc73a0a8d3dcd9542dc26565ed3dfd55de5)
1*d142edc7SBarry Smith /*$Id: snesj.c,v 1.74 2001/08/21 21:03:49 bsmith Exp bsmith $*/
211320018SBarry Smith 
3e090d566SSatish Balay #include "src/snes/snesimpl.h"    /*I  "petscsnes.h"  I*/
411320018SBarry Smith 
54a2ae208SSatish Balay #undef __FUNCT__
64a2ae208SSatish Balay #define __FUNCT__ "SNESDefaultComputeJacobian"
74b828684SBarry Smith /*@C
884a7bf8dSLois Curfman McInnes    SNESDefaultComputeJacobian - Computes the Jacobian using finite differences.
911320018SBarry Smith 
10fee21e36SBarry Smith    Collective on SNES
11fee21e36SBarry Smith 
12c7afd0dbSLois Curfman McInnes    Input Parameters:
13c7afd0dbSLois Curfman McInnes +  x1 - compute Jacobian at this point
14c7afd0dbSLois Curfman McInnes -  ctx - application's function context, as set with SNESSetFunction()
15c7afd0dbSLois Curfman McInnes 
16c7afd0dbSLois Curfman McInnes    Output Parameters:
17c7afd0dbSLois Curfman McInnes +  J - Jacobian matrix (not altered in this routine)
18c7afd0dbSLois Curfman McInnes .  B - newly computed Jacobian matrix to use with preconditioner (generally the same as J)
19c7afd0dbSLois Curfman McInnes -  flag - flag indicating whether the matrix sparsity structure has changed
20c7afd0dbSLois Curfman McInnes 
21ad960d00SLois Curfman McInnes    Options Database Key:
2208f75eefSBarry Smith +  -snes_fd - Activates SNESDefaultComputeJacobian()
2377d8c4bbSBarry Smith -  -snes_test_err - Square root of function error tolerance, default square root of machine
2477d8c4bbSBarry Smith                     epsilon (1.e-8 in double, 3.e-4 in single)
25ad960d00SLois Curfman McInnes 
265f3c43d9SLois Curfman McInnes    Notes:
275f3c43d9SLois Curfman McInnes    This routine is slow and expensive, and is not currently optimized
285f3c43d9SLois Curfman McInnes    to take advantage of sparsity in the problem.  Although
295f3c43d9SLois Curfman McInnes    SNESDefaultComputeJacobian() is not recommended for general use
305f3c43d9SLois Curfman McInnes    in large-scale applications, It can be useful in checking the
315f3c43d9SLois Curfman McInnes    correctness of a user-provided Jacobian.
3211320018SBarry Smith 
33b4fc646aSLois Curfman McInnes    An alternative routine that uses coloring to explot matrix sparsity is
342d0c0e3bSBarry Smith    SNESDefaultComputeJacobianColor().
35b4fc646aSLois Curfman McInnes 
3636851e7fSLois Curfman McInnes    Level: intermediate
3736851e7fSLois Curfman McInnes 
385f3c43d9SLois Curfman McInnes .keywords: SNES, finite differences, Jacobian
395f3c43d9SLois Curfman McInnes 
402d0c0e3bSBarry Smith .seealso: SNESSetJacobian(), SNESDefaultComputeJacobianColor()
4111320018SBarry Smith @*/
422d0c0e3bSBarry Smith int SNESDefaultComputeJacobian(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
4311320018SBarry Smith {
4488c956adSLois Curfman McInnes   Vec         j1a,j2a,x2;
4523242f5aSBarry Smith   int         i,ierr,N,start,end,j;
46ea709b57SSatish Balay   PetscScalar dx,mone = -1.0,*y,scale,*xx,wscale;
4777d8c4bbSBarry Smith   PetscReal   amax,epsilon = PETSC_SQRT_MACHINE_EPSILON;
48329f5518SBarry Smith   PetscReal   dx_min = 1.e-16,dx_par = 1.e-1;
49bbb6d6a8SBarry Smith   MPI_Comm    comm;
50a305c92eSSatish Balay   int         (*eval_fct)(SNES,Vec,Vec)=0;
510521c3abSLois Curfman McInnes 
523a40ed3dSBarry Smith   PetscFunctionBegin;
5387828ca2SBarry Smith   ierr = PetscOptionsGetReal(snes->prefix,"-snes_test_err",&epsilon,0);CHKERRQ(ierr);
543a40ed3dSBarry Smith   if (snes->method_class == SNES_NONLINEAR_EQUATIONS) eval_fct = SNESComputeFunction;
553a40ed3dSBarry Smith   else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) eval_fct = SNESComputeGradient;
5629bbc08cSBarry Smith   else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid method class");
5723242f5aSBarry Smith 
582d0c0e3bSBarry Smith   ierr = PetscObjectGetComm((PetscObject)x1,&comm);CHKERRQ(ierr);
592d0c0e3bSBarry Smith   ierr = MatZeroEntries(*B);CHKERRQ(ierr);
60aa79bc6dSLois Curfman McInnes   if (!snes->nvwork) {
61aa79bc6dSLois Curfman McInnes     ierr = VecDuplicateVecs(x1,3,&snes->vwork);CHKERRQ(ierr);
62aa79bc6dSLois Curfman McInnes     snes->nvwork = 3;
63b0a32e0cSBarry Smith     PetscLogObjectParents(snes,3,snes->vwork);
64aa79bc6dSLois Curfman McInnes   }
6588c956adSLois Curfman McInnes   j1a = snes->vwork[0]; j2a = snes->vwork[1]; x2 = snes->vwork[2];
6623242f5aSBarry Smith 
6778b31e54SBarry Smith   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
6878b31e54SBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
69a86fb38aSBarry Smith   ierr = (*eval_fct)(snes,x1,j1a);CHKERRQ(ierr);
70c005e166SLois Curfman McInnes 
71c005e166SLois Curfman McInnes   /* Compute Jacobian approximation, 1 column at a time.
7288c956adSLois Curfman McInnes       x1 = current iterate, j1a = F(x1)
7388c956adSLois Curfman McInnes       x2 = perturbed iterate, j2a = F(x2)
74c005e166SLois Curfman McInnes    */
7539e2f89bSBarry Smith   for (i=0; i<N; i++) {
7678b31e54SBarry Smith     ierr = VecCopy(x1,x2);CHKERRQ(ierr);
7723242f5aSBarry Smith     if (i>= start && i<end) {
782e8a6d31SBarry Smith       ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);
7939e2f89bSBarry Smith       dx = xx[i-start];
802e8a6d31SBarry Smith       ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
81aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
8254d73c9fSLois Curfman McInnes       if (dx < dx_min && dx >= 0.0) dx = dx_par;
8354d73c9fSLois Curfman McInnes       else if (dx < 0.0 && dx > -dx_min) dx = -dx_par;
8419a167f6SBarry Smith #else
85329f5518SBarry Smith       if (PetscAbsScalar(dx) < dx_min && PetscRealPart(dx) >= 0.0) dx = dx_par;
86329f5518SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < dx_min) dx = -dx_par;
8719a167f6SBarry Smith #endif
8839e2f89bSBarry Smith       dx *= epsilon;
8974f6f00dSLois Curfman McInnes       wscale = 1.0/dx;
902e8a6d31SBarry Smith       ierr = VecSetValues(x2,1,&i,&dx,ADD_VALUES);CHKERRQ(ierr);
916c783aadSBarry Smith     } else {
92bbb6d6a8SBarry Smith       wscale = 0.0;
93bbb6d6a8SBarry Smith     }
94a86fb38aSBarry Smith     ierr = (*eval_fct)(snes,x2,j2a);CHKERRQ(ierr);
9588c956adSLois Curfman McInnes     ierr = VecAXPY(&mone,j1a,j2a);CHKERRQ(ierr);
96c005e166SLois Curfman McInnes     /* Communicate scale to all processors */
976831982aSBarry Smith     ierr = MPI_Allreduce(&wscale,&scale,1,MPIU_SCALAR,PetscSum_Op,comm);CHKERRQ(ierr);
982e8a6d31SBarry Smith     ierr = VecScale(&scale,j2a);CHKERRQ(ierr);
992e8a6d31SBarry Smith     ierr = VecNorm(j2a,NORM_INFINITY,&amax);CHKERRQ(ierr); amax *= 1.e-14;
100*d142edc7SBarry Smith     ierr = VecGetArray(j2a,&y);CHKERRQ(ierr);
10123242f5aSBarry Smith     for (j=start; j<end; j++) {
102cddf8d76SBarry Smith       if (PetscAbsScalar(y[j-start]) > amax) {
103b4fc646aSLois Curfman McInnes         ierr = MatSetValues(*B,1,&j,1,&i,y+j-start,INSERT_VALUES);CHKERRQ(ierr);
10423242f5aSBarry Smith       }
10523242f5aSBarry Smith     }
1062e8a6d31SBarry Smith     ierr = VecRestoreArray(j2a,&y);CHKERRQ(ierr);
10723242f5aSBarry Smith   }
108b4fc646aSLois Curfman McInnes   ierr  = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
109b4fc646aSLois Curfman McInnes   ierr  = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
110595b82d2SBarry Smith   *flag =  DIFFERENT_NONZERO_PATTERN;
1113a40ed3dSBarry Smith   PetscFunctionReturn(0);
11211320018SBarry Smith }
11311320018SBarry Smith 
1144a2ae208SSatish Balay #undef __FUNCT__
1154a2ae208SSatish Balay #define __FUNCT__ "SNESDefaultComputeHessian"
1160521c3abSLois Curfman McInnes /*@C
11784a7bf8dSLois Curfman McInnes    SNESDefaultComputeHessian - Computes the Hessian using finite differences.
1180521c3abSLois Curfman McInnes 
119fee21e36SBarry Smith    Collective on SNES
120fee21e36SBarry Smith 
121c7afd0dbSLois Curfman McInnes    Input Parameters:
122c7afd0dbSLois Curfman McInnes +  x1 - compute Hessian at this point
123c7afd0dbSLois Curfman McInnes -  ctx - application's gradient context, as set with SNESSetGradient()
124c7afd0dbSLois Curfman McInnes 
125c7afd0dbSLois Curfman McInnes    Output Parameters:
126c7afd0dbSLois Curfman McInnes +  J - Hessian matrix (not altered in this routine)
127c7afd0dbSLois Curfman McInnes .  B - newly computed Hessian matrix to use with preconditioner (generally the same as J)
128c7afd0dbSLois Curfman McInnes -  flag - flag indicating whether the matrix sparsity structure has changed
129c7afd0dbSLois Curfman McInnes 
1300521c3abSLois Curfman McInnes    Options Database Key:
131c7afd0dbSLois Curfman McInnes $  -snes_fd - Activates SNESDefaultComputeHessian()
1320521c3abSLois Curfman McInnes 
13315091d37SBarry Smith 
13415091d37SBarry Smith    Level: intermediate
13515091d37SBarry Smith 
1360521c3abSLois Curfman McInnes    Notes:
1370521c3abSLois Curfman McInnes    This routine is slow and expensive, and is not currently optimized
1380521c3abSLois Curfman McInnes    to take advantage of sparsity in the problem.  Although
1390521c3abSLois Curfman McInnes    SNESDefaultComputeHessian() is not recommended for general use
1400521c3abSLois Curfman McInnes    in large-scale applications, It can be useful in checking the
1410521c3abSLois Curfman McInnes    correctness of a user-provided Hessian.
1420521c3abSLois Curfman McInnes 
1430521c3abSLois Curfman McInnes .keywords: SNES, finite differences, Hessian
1440521c3abSLois Curfman McInnes 
1459984d6fbSBarry Smith .seealso: SNESSetHessian()
1460521c3abSLois Curfman McInnes @*/
147329f5518SBarry Smith int SNESDefaultComputeHessian(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
1480521c3abSLois Curfman McInnes {
1493a40ed3dSBarry Smith   int ierr;
1503a40ed3dSBarry Smith 
1513a40ed3dSBarry Smith   PetscFunctionBegin;
1523a40ed3dSBarry Smith   ierr = SNESDefaultComputeJacobian(snes,x1,J,B,flag,ctx);CHKERRQ(ierr);
1533a40ed3dSBarry Smith   PetscFunctionReturn(0);
1540521c3abSLois Curfman McInnes }
155bd45824fSLois Curfman McInnes 
1564a2ae208SSatish Balay #undef __FUNCT__
1574a2ae208SSatish Balay #define __FUNCT__ "SNESDefaultComputeHessianColor"
158bd45824fSLois Curfman McInnes /*@C
159bd45824fSLois Curfman McInnes    SNESDefaultComputeHessianColor - Computes the Hessian using colored finite differences.
160bd45824fSLois Curfman McInnes 
161bd45824fSLois Curfman McInnes    Collective on SNES
162bd45824fSLois Curfman McInnes 
163bd45824fSLois Curfman McInnes    Input Parameters:
164bd45824fSLois Curfman McInnes +  x1 - compute Hessian at this point
165bd45824fSLois Curfman McInnes -  ctx - application's gradient context, as set with SNESSetGradient()
166bd45824fSLois Curfman McInnes 
167bd45824fSLois Curfman McInnes    Output Parameters:
168bd45824fSLois Curfman McInnes +  J - Hessian matrix (not altered in this routine)
169bd45824fSLois Curfman McInnes .  B - newly computed Hessian matrix to use with preconditioner (generally the same as J)
170bd45824fSLois Curfman McInnes -  flag - flag indicating whether the matrix sparsity structure has changed
171bd45824fSLois Curfman McInnes 
172bd45824fSLois Curfman McInnes     Options Database Keys:
173bd45824fSLois Curfman McInnes .  -mat_fd_coloring_freq <freq> - Activates SNESDefaultComputeJacobianColor()
174bd45824fSLois Curfman McInnes 
175bd45824fSLois Curfman McInnes    Level: intermediate
176bd45824fSLois Curfman McInnes 
177bd45824fSLois Curfman McInnes  .keywords: SNES, finite differences, Hessian, coloring, sparse
178bd45824fSLois Curfman McInnes 
179bd45824fSLois Curfman McInnes .seealso: SNESSetHessian()
180bd45824fSLois Curfman McInnes @*/
181329f5518SBarry Smith int SNESDefaultComputeHessianColor(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
182bd45824fSLois Curfman McInnes {
183bd45824fSLois Curfman McInnes   int ierr;
184bd45824fSLois Curfman McInnes 
185bd45824fSLois Curfman McInnes   PetscFunctionBegin;
186bd45824fSLois Curfman McInnes   ierr = SNESDefaultComputeJacobianColor(snes,x1,J,B,flag,ctx);CHKERRQ(ierr);
187bd45824fSLois Curfman McInnes   PetscFunctionReturn(0);
188bd45824fSLois Curfman McInnes }
189bd45824fSLois Curfman McInnes 
190