xref: /petsc/src/snes/interface/snesj.c (revision f588057b4e26d9feb8696312515b7a30eda4868c)
163dd3a1aSKris Buschelman #define PETSCSNES_DLL
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 @*/
4263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultComputeJacobian(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx)
4311320018SBarry Smith {
4488c956adSLois Curfman McInnes   Vec            j1a,j2a,x2;
456849ba73SBarry Smith   PetscErrorCode ierr;
460b9b6f31SBarry Smith   PetscInt       i,N,start,end,j;
47ea709b57SSatish Balay   PetscScalar    dx,mone = -1.0,*y,scale,*xx,wscale;
4877d8c4bbSBarry Smith   PetscReal      amax,epsilon = PETSC_SQRT_MACHINE_EPSILON;
49329f5518SBarry Smith   PetscReal      dx_min = 1.e-16,dx_par = 1.e-1;
50bbb6d6a8SBarry Smith   MPI_Comm       comm;
516849ba73SBarry Smith   PetscErrorCode (*eval_fct)(SNES,Vec,Vec)=0;
520b9b6f31SBarry Smith   PetscTruth     assembled;
530521c3abSLois Curfman McInnes 
543a40ed3dSBarry Smith   PetscFunctionBegin;
5587828ca2SBarry Smith   ierr = PetscOptionsGetReal(snes->prefix,"-snes_test_err",&epsilon,0);CHKERRQ(ierr);
564b27c08aSLois Curfman McInnes   eval_fct = SNESComputeFunction;
5723242f5aSBarry Smith 
582d0c0e3bSBarry Smith   ierr = PetscObjectGetComm((PetscObject)x1,&comm);CHKERRQ(ierr);
590b9b6f31SBarry Smith   ierr = MatAssembled(*B,&assembled);CHKERRQ(ierr);
600b9b6f31SBarry Smith   if (assembled) {
612d0c0e3bSBarry Smith     ierr = MatZeroEntries(*B);CHKERRQ(ierr);
620b9b6f31SBarry Smith   }
63aa79bc6dSLois Curfman McInnes   if (!snes->nvwork) {
64aa79bc6dSLois Curfman McInnes     ierr = VecDuplicateVecs(x1,3,&snes->vwork);CHKERRQ(ierr);
65aa79bc6dSLois Curfman McInnes     snes->nvwork = 3;
66efee365bSSatish Balay     ierr = PetscLogObjectParents(snes,3,snes->vwork);CHKERRQ(ierr);
67aa79bc6dSLois Curfman McInnes   }
6888c956adSLois Curfman McInnes   j1a = snes->vwork[0]; j2a = snes->vwork[1]; x2 = snes->vwork[2];
6923242f5aSBarry Smith 
7078b31e54SBarry Smith   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
7178b31e54SBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
72a86fb38aSBarry Smith   ierr = (*eval_fct)(snes,x1,j1a);CHKERRQ(ierr);
73c005e166SLois Curfman McInnes 
74c005e166SLois Curfman McInnes   /* Compute Jacobian approximation, 1 column at a time.
7588c956adSLois Curfman McInnes       x1 = current iterate, j1a = F(x1)
7688c956adSLois Curfman McInnes       x2 = perturbed iterate, j2a = F(x2)
77c005e166SLois Curfman McInnes    */
7839e2f89bSBarry Smith   for (i=0; i<N; i++) {
7978b31e54SBarry Smith     ierr = VecCopy(x1,x2);CHKERRQ(ierr);
8023242f5aSBarry Smith     if (i>= start && i<end) {
812e8a6d31SBarry Smith       ierr = VecGetArray(x1,&xx);CHKERRQ(ierr);
8239e2f89bSBarry Smith       dx = xx[i-start];
832e8a6d31SBarry Smith       ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr);
84aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX)
8554d73c9fSLois Curfman McInnes       if (dx < dx_min && dx >= 0.0) dx = dx_par;
8654d73c9fSLois Curfman McInnes       else if (dx < 0.0 && dx > -dx_min) dx = -dx_par;
8719a167f6SBarry Smith #else
88329f5518SBarry Smith       if (PetscAbsScalar(dx) < dx_min && PetscRealPart(dx) >= 0.0) dx = dx_par;
89329f5518SBarry Smith       else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < dx_min) dx = -dx_par;
9019a167f6SBarry Smith #endif
9139e2f89bSBarry Smith       dx *= epsilon;
9274f6f00dSLois Curfman McInnes       wscale = 1.0/dx;
932e8a6d31SBarry Smith       ierr = VecSetValues(x2,1,&i,&dx,ADD_VALUES);CHKERRQ(ierr);
946c783aadSBarry Smith     } else {
95bbb6d6a8SBarry Smith       wscale = 0.0;
96bbb6d6a8SBarry Smith     }
97a86fb38aSBarry Smith     ierr = (*eval_fct)(snes,x2,j2a);CHKERRQ(ierr);
9888c956adSLois Curfman McInnes     ierr = VecAXPY(&mone,j1a,j2a);CHKERRQ(ierr);
99c005e166SLois Curfman McInnes     /* Communicate scale to all processors */
1006831982aSBarry Smith     ierr = MPI_Allreduce(&wscale,&scale,1,MPIU_SCALAR,PetscSum_Op,comm);CHKERRQ(ierr);
1012e8a6d31SBarry Smith     ierr = VecScale(&scale,j2a);CHKERRQ(ierr);
1022e8a6d31SBarry Smith     ierr = VecNorm(j2a,NORM_INFINITY,&amax);CHKERRQ(ierr); amax *= 1.e-14;
103d142edc7SBarry Smith     ierr = VecGetArray(j2a,&y);CHKERRQ(ierr);
10423242f5aSBarry Smith     for (j=start; j<end; j++) {
105cddf8d76SBarry Smith       if (PetscAbsScalar(y[j-start]) > amax) {
106b4fc646aSLois Curfman McInnes         ierr = MatSetValues(*B,1,&j,1,&i,y+j-start,INSERT_VALUES);CHKERRQ(ierr);
10723242f5aSBarry Smith       }
10823242f5aSBarry Smith     }
1092e8a6d31SBarry Smith     ierr = VecRestoreArray(j2a,&y);CHKERRQ(ierr);
11023242f5aSBarry Smith   }
111b4fc646aSLois Curfman McInnes   ierr  = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
112b4fc646aSLois Curfman McInnes   ierr  = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
113*f588057bSBarry Smith   if (*B != *J) {
114*f588057bSBarry Smith     ierr  = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
115*f588057bSBarry Smith     ierr  = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
116*f588057bSBarry Smith   }
117595b82d2SBarry Smith   *flag =  DIFFERENT_NONZERO_PATTERN;
1183a40ed3dSBarry Smith   PetscFunctionReturn(0);
11911320018SBarry Smith }
12011320018SBarry Smith 
121bd45824fSLois Curfman McInnes 
122