xref: /petsc/src/snes/interface/snesj.c (revision 6b867d5ac32ed0c728f185df9d084acdf26f70bf)
111320018SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/snesimpl.h>    /*I  "petscsnes.h"  I*/
3b1f624c7SBarry Smith #include <petscdm.h>
411320018SBarry Smith 
54b828684SBarry Smith /*@C
68d359177SBarry Smith    SNESComputeJacobianDefault - Computes the Jacobian using finite differences.
711320018SBarry Smith 
8fee21e36SBarry Smith    Collective on SNES
9fee21e36SBarry Smith 
10c7afd0dbSLois Curfman McInnes    Input Parameters:
11*6b867d5aSJose E. Roman +  snes - the SNES context
12*6b867d5aSJose E. Roman .  x1 - compute Jacobian at this point
13c7afd0dbSLois Curfman McInnes -  ctx - application's function context, as set with SNESSetFunction()
14c7afd0dbSLois Curfman McInnes 
15c7afd0dbSLois Curfman McInnes    Output Parameters:
16c7afd0dbSLois Curfman McInnes +  J - Jacobian matrix (not altered in this routine)
17d1e9a80fSBarry Smith -  B - newly computed Jacobian matrix to use with preconditioner (generally the same as J)
18c7afd0dbSLois Curfman McInnes 
19ad960d00SLois Curfman McInnes    Options Database Key:
208d359177SBarry Smith +  -snes_fd - Activates SNESComputeJacobianDefault()
2179f36460SBarry Smith .  -snes_test_err - Square root of function error tolerance, default square root of machine
2277d8c4bbSBarry Smith                     epsilon (1.e-8 in double, 3.e-4 in single)
233ec795f1SBarry Smith -  -mat_fd_type - Either wp or ds (see MATMFFD_WP or MATMFFD_DS)
24ad960d00SLois Curfman McInnes 
255f3c43d9SLois Curfman McInnes    Notes:
265f3c43d9SLois Curfman McInnes    This routine is slow and expensive, and is not currently optimized
275f3c43d9SLois Curfman McInnes    to take advantage of sparsity in the problem.  Although
288d359177SBarry Smith    SNESComputeJacobianDefault() is not recommended for general use
295f3c43d9SLois Curfman McInnes    in large-scale applications, It can be useful in checking the
305f3c43d9SLois Curfman McInnes    correctness of a user-provided Jacobian.
3111320018SBarry Smith 
3279f36460SBarry Smith    An alternative routine that uses coloring to exploit matrix sparsity is
338d359177SBarry Smith    SNESComputeJacobianDefaultColor().
34b4fc646aSLois Curfman McInnes 
350df40c35SBarry Smith    This routine ignores the maximum number of function evaluations set with SNESSetTolerances() and the function
360df40c35SBarry Smith    evaluations it performs are not counted in what is returned by of SNESGetNumberFunctionEvals().
370df40c35SBarry Smith 
3836851e7fSLois Curfman McInnes    Level: intermediate
3936851e7fSLois Curfman McInnes 
408d359177SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobianDefaultColor(), MatCreateSNESMF()
4111320018SBarry Smith @*/
42d1e9a80fSBarry Smith PetscErrorCode  SNESComputeJacobianDefault(SNES snes,Vec x1,Mat J,Mat B,void *ctx)
4311320018SBarry Smith {
4488c956adSLois Curfman McInnes   Vec               j1a,j2a,x2;
456849ba73SBarry Smith   PetscErrorCode    ierr;
460df40c35SBarry Smith   PetscInt          i,N,start,end,j,value,root,max_funcs = snes->max_funcs;
475edff71fSBarry Smith   PetscScalar       dx,*y,wscale;
485edff71fSBarry Smith   const PetscScalar *xx;
4977d8c4bbSBarry Smith   PetscReal         amax,epsilon = PETSC_SQRT_MACHINE_EPSILON;
5079f36460SBarry Smith   PetscReal         dx_min = 1.e-16,dx_par = 1.e-1,unorm;
51bbb6d6a8SBarry Smith   MPI_Comm          comm;
52ace3abfcSBarry Smith   PetscBool         assembled,use_wp = PETSC_TRUE,flg;
5379f36460SBarry Smith   const char        *list[2] = {"ds","wp"};
5486c88abdSHong Zhang   PetscMPIInt       size;
5586c88abdSHong Zhang   const PetscInt    *ranges;
560df40c35SBarry Smith   DM                dm;
570df40c35SBarry Smith   DMSNES            dms;
580521c3abSLois Curfman McInnes 
593a40ed3dSBarry Smith   PetscFunctionBegin;
600df40c35SBarry Smith   snes->max_funcs = PETSC_MAX_INT;
615d5a84c9SBarry Smith   /* Since this Jacobian will possibly have "extra" nonzero locations just turn off errors for these locations */
625d5a84c9SBarry Smith   ierr = MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE);CHKERRQ(ierr);
639e5d0892SLisandro Dalcin   ierr = PetscOptionsGetReal(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_test_err",&epsilon,NULL);CHKERRQ(ierr);
6423242f5aSBarry Smith 
652d0c0e3bSBarry Smith   ierr = PetscObjectGetComm((PetscObject)x1,&comm);CHKERRQ(ierr);
66ffc4695bSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr);
6794ab13aaSBarry Smith   ierr = MatAssembled(B,&assembled);CHKERRQ(ierr);
680b9b6f31SBarry Smith   if (assembled) {
6994ab13aaSBarry Smith     ierr = MatZeroEntries(B);CHKERRQ(ierr);
700b9b6f31SBarry Smith   }
71aa79bc6dSLois Curfman McInnes   if (!snes->nvwork) {
72b1f624c7SBarry Smith     if (snes->dm) {
73b1f624c7SBarry Smith       ierr = DMGetGlobalVector(snes->dm,&j1a);CHKERRQ(ierr);
74b1f624c7SBarry Smith       ierr = DMGetGlobalVector(snes->dm,&j2a);CHKERRQ(ierr);
75b1f624c7SBarry Smith       ierr = DMGetGlobalVector(snes->dm,&x2);CHKERRQ(ierr);
76b1f624c7SBarry Smith     } else {
77aa79bc6dSLois Curfman McInnes       snes->nvwork = 3;
7885385478SLisandro Dalcin       ierr = VecDuplicateVecs(x1,snes->nvwork,&snes->vwork);CHKERRQ(ierr);
7985385478SLisandro Dalcin       ierr = PetscLogObjectParents(snes,snes->nvwork,snes->vwork);CHKERRQ(ierr);
8088c956adSLois Curfman McInnes       j1a = snes->vwork[0]; j2a = snes->vwork[1]; x2 = snes->vwork[2];
81b1f624c7SBarry Smith     }
82b1f624c7SBarry Smith   }
8323242f5aSBarry Smith 
8478b31e54SBarry Smith   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
8578b31e54SBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
860df40c35SBarry Smith   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
870df40c35SBarry Smith   ierr = DMGetDMSNES(dm,&dms);CHKERRQ(ierr);
880df40c35SBarry Smith   if (dms->ops->computemffunction) {
890df40c35SBarry Smith     ierr = SNESComputeMFFunction(snes,x1,j1a);CHKERRQ(ierr);
900df40c35SBarry Smith   } else {
9162cd874fSBarry Smith     ierr = SNESComputeFunction(snes,x1,j1a);CHKERRQ(ierr);
920df40c35SBarry Smith   }
93c005e166SLois Curfman McInnes 
94302440fdSBarry Smith   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"Differencing options","SNES");CHKERRQ(ierr);
958d359177SBarry Smith   ierr = PetscOptionsEList("-mat_fd_type","Algorithm to compute difference parameter","SNESComputeJacobianDefault",list,2,"wp",&value,&flg);CHKERRQ(ierr);
96302440fdSBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
97f5af7f23SKarl Rupp   if (flg && !value) use_wp = PETSC_FALSE;
98f5af7f23SKarl Rupp 
9979f36460SBarry Smith   if (use_wp) {
10079f36460SBarry Smith     ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr);
10179f36460SBarry Smith   }
102c005e166SLois Curfman McInnes   /* Compute Jacobian approximation, 1 column at a time.
10388c956adSLois Curfman McInnes       x1 = current iterate, j1a = F(x1)
10488c956adSLois Curfman McInnes       x2 = perturbed iterate, j2a = F(x2)
105c005e166SLois Curfman McInnes    */
10639e2f89bSBarry Smith   for (i=0; i<N; i++) {
10778b31e54SBarry Smith     ierr = VecCopy(x1,x2);CHKERRQ(ierr);
10823242f5aSBarry Smith     if (i>= start && i<end) {
1095edff71fSBarry Smith       ierr = VecGetArrayRead(x1,&xx);CHKERRQ(ierr);
1105620d6dcSBarry Smith       if (use_wp) dx = PetscSqrtReal(1.0 + unorm);
111f5af7f23SKarl Rupp       else        dx = xx[i-start];
1125edff71fSBarry Smith       ierr = VecRestoreArrayRead(x1,&xx);CHKERRQ(ierr);
1136bd79f97SJed Brown       if (PetscAbsScalar(dx) < dx_min) dx = (PetscRealPart(dx) < 0. ? -1. : 1.) * dx_par;
11439e2f89bSBarry Smith       dx    *= epsilon;
11574f6f00dSLois Curfman McInnes       wscale = 1.0/dx;
1162e8a6d31SBarry Smith       ierr   = VecSetValues(x2,1,&i,&dx,ADD_VALUES);CHKERRQ(ierr);
1176c783aadSBarry Smith     } else {
118bbb6d6a8SBarry Smith       wscale = 0.0;
119bbb6d6a8SBarry Smith     }
120275089ecSBarry Smith     ierr = VecAssemblyBegin(x2);CHKERRQ(ierr);
121275089ecSBarry Smith     ierr = VecAssemblyEnd(x2);CHKERRQ(ierr);
1220df40c35SBarry Smith     if (dms->ops->computemffunction) {
1230df40c35SBarry Smith       ierr = SNESComputeMFFunction(snes,x2,j2a);CHKERRQ(ierr);
1240df40c35SBarry Smith     } else {
12562cd874fSBarry Smith       ierr = SNESComputeFunction(snes,x2,j2a);CHKERRQ(ierr);
1260df40c35SBarry Smith     }
12779f36460SBarry Smith     ierr = VecAXPY(j2a,-1.0,j1a);CHKERRQ(ierr);
12886c88abdSHong Zhang     /* Communicate scale=1/dx_i to all processors */
12986c88abdSHong Zhang     ierr = VecGetOwnershipRanges(x1,&ranges);CHKERRQ(ierr);
13086c88abdSHong Zhang     root = size;
13186c88abdSHong Zhang     for (j=size-1; j>-1; j--) {
13286c88abdSHong Zhang       root--;
13386c88abdSHong Zhang       if (i>=ranges[j]) break;
13486c88abdSHong Zhang     }
135ffc4695bSBarry Smith     ierr = MPI_Bcast(&wscale,1,MPIU_SCALAR,root,comm);CHKERRMPI(ierr);
13686c88abdSHong Zhang     ierr = VecScale(j2a,wscale);CHKERRQ(ierr);
1372e8a6d31SBarry Smith     ierr = VecNorm(j2a,NORM_INFINITY,&amax);CHKERRQ(ierr); amax *= 1.e-14;
138d142edc7SBarry Smith     ierr = VecGetArray(j2a,&y);CHKERRQ(ierr);
13923242f5aSBarry Smith     for (j=start; j<end; j++) {
140f6e735bbSHong Zhang       if (PetscAbsScalar(y[j-start]) > amax || j == i) {
14194ab13aaSBarry Smith         ierr = MatSetValues(B,1,&j,1,&i,y+j-start,INSERT_VALUES);CHKERRQ(ierr);
14223242f5aSBarry Smith       }
14323242f5aSBarry Smith     }
1442e8a6d31SBarry Smith     ierr = VecRestoreArray(j2a,&y);CHKERRQ(ierr);
14523242f5aSBarry Smith   }
146b1f624c7SBarry Smith   if (snes->dm) {
147b1f624c7SBarry Smith     ierr = DMRestoreGlobalVector(snes->dm,&j1a);CHKERRQ(ierr);
148b1f624c7SBarry Smith     ierr = DMRestoreGlobalVector(snes->dm,&j2a);CHKERRQ(ierr);
149b1f624c7SBarry Smith     ierr = DMRestoreGlobalVector(snes->dm,&x2);CHKERRQ(ierr);
150b1f624c7SBarry Smith   }
15194ab13aaSBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15294ab13aaSBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15394ab13aaSBarry Smith   if (B != J) {
15494ab13aaSBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15594ab13aaSBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
156f588057bSBarry Smith   }
1570df40c35SBarry Smith   snes->max_funcs = max_funcs;
1580df40c35SBarry Smith   snes->nfuncs    -= N;
1593a40ed3dSBarry Smith   PetscFunctionReturn(0);
16011320018SBarry Smith }
16111320018SBarry Smith 
162