xref: /petsc/src/snes/interface/snesj.c (revision 9e5d08929599755f0fc0004c3575ee121dac326a)
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:
11c7afd0dbSLois Curfman McInnes +  x1 - compute Jacobian at this point
12c7afd0dbSLois Curfman McInnes -  ctx - application's function context, as set with SNESSetFunction()
13c7afd0dbSLois Curfman McInnes 
14c7afd0dbSLois Curfman McInnes    Output Parameters:
15c7afd0dbSLois Curfman McInnes +  J - Jacobian matrix (not altered in this routine)
16d1e9a80fSBarry Smith -  B - newly computed Jacobian matrix to use with preconditioner (generally the same as J)
17c7afd0dbSLois Curfman McInnes 
18ad960d00SLois Curfman McInnes    Options Database Key:
198d359177SBarry Smith +  -snes_fd - Activates SNESComputeJacobianDefault()
2079f36460SBarry Smith .  -snes_test_err - Square root of function error tolerance, default square root of machine
2177d8c4bbSBarry Smith                     epsilon (1.e-8 in double, 3.e-4 in single)
223ec795f1SBarry Smith -  -mat_fd_type - Either wp or ds (see MATMFFD_WP or MATMFFD_DS)
23ad960d00SLois Curfman McInnes 
245f3c43d9SLois Curfman McInnes    Notes:
255f3c43d9SLois Curfman McInnes    This routine is slow and expensive, and is not currently optimized
265f3c43d9SLois Curfman McInnes    to take advantage of sparsity in the problem.  Although
278d359177SBarry Smith    SNESComputeJacobianDefault() is not recommended for general use
285f3c43d9SLois Curfman McInnes    in large-scale applications, It can be useful in checking the
295f3c43d9SLois Curfman McInnes    correctness of a user-provided Jacobian.
3011320018SBarry Smith 
3179f36460SBarry Smith    An alternative routine that uses coloring to exploit matrix sparsity is
328d359177SBarry Smith    SNESComputeJacobianDefaultColor().
33b4fc646aSLois Curfman McInnes 
3436851e7fSLois Curfman McInnes    Level: intermediate
3536851e7fSLois Curfman McInnes 
368d359177SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobianDefaultColor(), MatCreateSNESMF()
3711320018SBarry Smith @*/
38d1e9a80fSBarry Smith PetscErrorCode  SNESComputeJacobianDefault(SNES snes,Vec x1,Mat J,Mat B,void *ctx)
3911320018SBarry Smith {
4088c956adSLois Curfman McInnes   Vec               j1a,j2a,x2;
416849ba73SBarry Smith   PetscErrorCode    ierr;
4286c88abdSHong Zhang   PetscInt          i,N,start,end,j,value,root;
435edff71fSBarry Smith   PetscScalar       dx,*y,wscale;
445edff71fSBarry Smith   const PetscScalar *xx;
4577d8c4bbSBarry Smith   PetscReal         amax,epsilon = PETSC_SQRT_MACHINE_EPSILON;
4679f36460SBarry Smith   PetscReal         dx_min = 1.e-16,dx_par = 1.e-1,unorm;
47bbb6d6a8SBarry Smith   MPI_Comm          comm;
48ace3abfcSBarry Smith   PetscBool         assembled,use_wp = PETSC_TRUE,flg;
4979f36460SBarry Smith   const char        *list[2] = {"ds","wp"};
5086c88abdSHong Zhang   PetscMPIInt       size;
5186c88abdSHong Zhang   const PetscInt    *ranges;
520521c3abSLois Curfman McInnes 
533a40ed3dSBarry Smith   PetscFunctionBegin;
545d5a84c9SBarry Smith   /* Since this Jacobian will possibly have "extra" nonzero locations just turn off errors for these locations */
555d5a84c9SBarry Smith   ierr = MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE);CHKERRQ(ierr);
56*9e5d0892SLisandro Dalcin   ierr = PetscOptionsGetReal(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_test_err",&epsilon,NULL);CHKERRQ(ierr);
5723242f5aSBarry Smith 
582d0c0e3bSBarry Smith   ierr = PetscObjectGetComm((PetscObject)x1,&comm);CHKERRQ(ierr);
5986c88abdSHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
6094ab13aaSBarry Smith   ierr = MatAssembled(B,&assembled);CHKERRQ(ierr);
610b9b6f31SBarry Smith   if (assembled) {
6294ab13aaSBarry Smith     ierr = MatZeroEntries(B);CHKERRQ(ierr);
630b9b6f31SBarry Smith   }
64aa79bc6dSLois Curfman McInnes   if (!snes->nvwork) {
65b1f624c7SBarry Smith     if (snes->dm) {
66b1f624c7SBarry Smith       ierr = DMGetGlobalVector(snes->dm,&j1a);CHKERRQ(ierr);
67b1f624c7SBarry Smith       ierr = DMGetGlobalVector(snes->dm,&j2a);CHKERRQ(ierr);
68b1f624c7SBarry Smith       ierr = DMGetGlobalVector(snes->dm,&x2);CHKERRQ(ierr);
69b1f624c7SBarry Smith     } else {
70aa79bc6dSLois Curfman McInnes       snes->nvwork = 3;
7185385478SLisandro Dalcin       ierr = VecDuplicateVecs(x1,snes->nvwork,&snes->vwork);CHKERRQ(ierr);
7285385478SLisandro Dalcin       ierr = PetscLogObjectParents(snes,snes->nvwork,snes->vwork);CHKERRQ(ierr);
7388c956adSLois Curfman McInnes       j1a = snes->vwork[0]; j2a = snes->vwork[1]; x2 = snes->vwork[2];
74b1f624c7SBarry Smith     }
75b1f624c7SBarry Smith   }
7623242f5aSBarry Smith 
7778b31e54SBarry Smith   ierr = VecGetSize(x1,&N);CHKERRQ(ierr);
7878b31e54SBarry Smith   ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr);
7962cd874fSBarry Smith   ierr = SNESComputeFunction(snes,x1,j1a);CHKERRQ(ierr);
80c005e166SLois Curfman McInnes 
81302440fdSBarry Smith   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"Differencing options","SNES");CHKERRQ(ierr);
828d359177SBarry Smith   ierr = PetscOptionsEList("-mat_fd_type","Algorithm to compute difference parameter","SNESComputeJacobianDefault",list,2,"wp",&value,&flg);CHKERRQ(ierr);
83302440fdSBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
84f5af7f23SKarl Rupp   if (flg && !value) use_wp = PETSC_FALSE;
85f5af7f23SKarl Rupp 
8679f36460SBarry Smith   if (use_wp) {
8779f36460SBarry Smith     ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr);
8879f36460SBarry Smith   }
89c005e166SLois Curfman McInnes   /* Compute Jacobian approximation, 1 column at a time.
9088c956adSLois Curfman McInnes       x1 = current iterate, j1a = F(x1)
9188c956adSLois Curfman McInnes       x2 = perturbed iterate, j2a = F(x2)
92c005e166SLois Curfman McInnes    */
9339e2f89bSBarry Smith   for (i=0; i<N; i++) {
9478b31e54SBarry Smith     ierr = VecCopy(x1,x2);CHKERRQ(ierr);
9523242f5aSBarry Smith     if (i>= start && i<end) {
965edff71fSBarry Smith       ierr = VecGetArrayRead(x1,&xx);CHKERRQ(ierr);
975620d6dcSBarry Smith       if (use_wp) dx = PetscSqrtReal(1.0 + unorm);
98f5af7f23SKarl Rupp       else        dx = xx[i-start];
995edff71fSBarry Smith       ierr = VecRestoreArrayRead(x1,&xx);CHKERRQ(ierr);
1006bd79f97SJed Brown       if (PetscAbsScalar(dx) < dx_min) dx = (PetscRealPart(dx) < 0. ? -1. : 1.) * dx_par;
10139e2f89bSBarry Smith       dx    *= epsilon;
10274f6f00dSLois Curfman McInnes       wscale = 1.0/dx;
1032e8a6d31SBarry Smith       ierr   = VecSetValues(x2,1,&i,&dx,ADD_VALUES);CHKERRQ(ierr);
1046c783aadSBarry Smith     } else {
105bbb6d6a8SBarry Smith       wscale = 0.0;
106bbb6d6a8SBarry Smith     }
107275089ecSBarry Smith     ierr = VecAssemblyBegin(x2);CHKERRQ(ierr);
108275089ecSBarry Smith     ierr = VecAssemblyEnd(x2);CHKERRQ(ierr);
10962cd874fSBarry Smith     ierr = SNESComputeFunction(snes,x2,j2a);CHKERRQ(ierr);
11079f36460SBarry Smith     ierr = VecAXPY(j2a,-1.0,j1a);CHKERRQ(ierr);
11186c88abdSHong Zhang     /* Communicate scale=1/dx_i to all processors */
11286c88abdSHong Zhang     ierr = VecGetOwnershipRanges(x1,&ranges);CHKERRQ(ierr);
11386c88abdSHong Zhang     root = size;
11486c88abdSHong Zhang     for (j=size-1; j>-1; j--) {
11586c88abdSHong Zhang       root--;
11686c88abdSHong Zhang       if (i>=ranges[j]) break;
11786c88abdSHong Zhang     }
11886c88abdSHong Zhang     ierr = MPI_Bcast(&wscale,1,MPIU_SCALAR,root,comm);CHKERRQ(ierr);
11986c88abdSHong Zhang 
12086c88abdSHong Zhang     ierr = VecScale(j2a,wscale);CHKERRQ(ierr);
1212e8a6d31SBarry Smith     ierr = VecNorm(j2a,NORM_INFINITY,&amax);CHKERRQ(ierr); amax *= 1.e-14;
122d142edc7SBarry Smith     ierr = VecGetArray(j2a,&y);CHKERRQ(ierr);
12323242f5aSBarry Smith     for (j=start; j<end; j++) {
124f6e735bbSHong Zhang       if (PetscAbsScalar(y[j-start]) > amax || j == i) {
12594ab13aaSBarry Smith         ierr = MatSetValues(B,1,&j,1,&i,y+j-start,INSERT_VALUES);CHKERRQ(ierr);
12623242f5aSBarry Smith       }
12723242f5aSBarry Smith     }
1282e8a6d31SBarry Smith     ierr = VecRestoreArray(j2a,&y);CHKERRQ(ierr);
12923242f5aSBarry Smith   }
130b1f624c7SBarry Smith   if (snes->dm) {
131b1f624c7SBarry Smith     ierr = DMRestoreGlobalVector(snes->dm,&j1a);CHKERRQ(ierr);
132b1f624c7SBarry Smith     ierr = DMRestoreGlobalVector(snes->dm,&j2a);CHKERRQ(ierr);
133b1f624c7SBarry Smith     ierr = DMRestoreGlobalVector(snes->dm,&x2);CHKERRQ(ierr);
134b1f624c7SBarry Smith   }
13594ab13aaSBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13694ab13aaSBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13794ab13aaSBarry Smith   if (B != J) {
13894ab13aaSBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
13994ab13aaSBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
140f588057bSBarry Smith   }
1413a40ed3dSBarry Smith   PetscFunctionReturn(0);
14211320018SBarry Smith }
14311320018SBarry Smith 
144bd45824fSLois Curfman McInnes 
145