163dd3a1aSKris Buschelman #define PETSCSNES_DLL 211320018SBarry Smith 37c4f633dSBarry Smith #include "private/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() 2379f36460SBarry 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) 253ec795f1SBarry Smith - -mat_fd_type - Either wp or ds (see MATMFFD_WP or MATMFFD_DS) 26ad960d00SLois Curfman McInnes 275f3c43d9SLois Curfman McInnes Notes: 285f3c43d9SLois Curfman McInnes This routine is slow and expensive, and is not currently optimized 295f3c43d9SLois Curfman McInnes to take advantage of sparsity in the problem. Although 305f3c43d9SLois Curfman McInnes SNESDefaultComputeJacobian() is not recommended for general use 315f3c43d9SLois Curfman McInnes in large-scale applications, It can be useful in checking the 325f3c43d9SLois Curfman McInnes correctness of a user-provided Jacobian. 3311320018SBarry Smith 3479f36460SBarry Smith An alternative routine that uses coloring to exploit matrix sparsity is 352d0c0e3bSBarry Smith SNESDefaultComputeJacobianColor(). 36b4fc646aSLois Curfman McInnes 3736851e7fSLois Curfman McInnes Level: intermediate 3836851e7fSLois Curfman McInnes 395f3c43d9SLois Curfman McInnes .keywords: SNES, finite differences, Jacobian 405f3c43d9SLois Curfman McInnes 4179f36460SBarry Smith .seealso: SNESSetJacobian(), SNESDefaultComputeJacobianColor(), MatCreateSNESMF() 4211320018SBarry Smith @*/ 4363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultComputeJacobian(SNES snes,Vec x1,Mat *J,Mat *B,MatStructure *flag,void *ctx) 4411320018SBarry Smith { 4588c956adSLois Curfman McInnes Vec j1a,j2a,x2; 466849ba73SBarry Smith PetscErrorCode ierr; 4786c88abdSHong Zhang PetscInt i,N,start,end,j,value,root; 4886c88abdSHong Zhang PetscScalar dx,*y,*xx,wscale; 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; 526849ba73SBarry Smith PetscErrorCode (*eval_fct)(SNES,Vec,Vec)=0; 5379f36460SBarry Smith PetscTruth assembled,use_wp = PETSC_TRUE,flg; 5479f36460SBarry Smith const char *list[2] = {"ds","wp"}; 5586c88abdSHong Zhang PetscMPIInt size; 5686c88abdSHong Zhang const PetscInt *ranges; 570521c3abSLois Curfman McInnes 583a40ed3dSBarry Smith PetscFunctionBegin; 597adad957SLisandro Dalcin ierr = PetscOptionsGetReal(((PetscObject)snes)->prefix,"-snes_test_err",&epsilon,0);CHKERRQ(ierr); 604b27c08aSLois Curfman McInnes eval_fct = SNESComputeFunction; 6123242f5aSBarry Smith 622d0c0e3bSBarry Smith ierr = PetscObjectGetComm((PetscObject)x1,&comm);CHKERRQ(ierr); 6386c88abdSHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 640b9b6f31SBarry Smith ierr = MatAssembled(*B,&assembled);CHKERRQ(ierr); 650b9b6f31SBarry Smith if (assembled) { 662d0c0e3bSBarry Smith ierr = MatZeroEntries(*B);CHKERRQ(ierr); 670b9b6f31SBarry Smith } 68aa79bc6dSLois Curfman McInnes if (!snes->nvwork) { 69aa79bc6dSLois Curfman McInnes snes->nvwork = 3; 7085385478SLisandro Dalcin ierr = VecDuplicateVecs(x1,snes->nvwork,&snes->vwork);CHKERRQ(ierr); 7185385478SLisandro Dalcin ierr = PetscLogObjectParents(snes,snes->nvwork,snes->vwork);CHKERRQ(ierr); 72aa79bc6dSLois Curfman McInnes } 7388c956adSLois Curfman McInnes j1a = snes->vwork[0]; j2a = snes->vwork[1]; x2 = snes->vwork[2]; 7423242f5aSBarry Smith 7578b31e54SBarry Smith ierr = VecGetSize(x1,&N);CHKERRQ(ierr); 7678b31e54SBarry Smith ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr); 77a86fb38aSBarry Smith ierr = (*eval_fct)(snes,x1,j1a);CHKERRQ(ierr); 78c005e166SLois Curfman McInnes 7979f36460SBarry Smith ierr = PetscOptionsEList("-mat_fd_type","Algorithm to compute difference parameter","SNESDefaultComputeJacobian",list,2,"wp",&value,&flg);CHKERRQ(ierr); 8079f36460SBarry Smith if (flg && !value) { 8179f36460SBarry Smith use_wp = PETSC_FALSE; 8279f36460SBarry Smith } 8379f36460SBarry Smith if (use_wp) { 8479f36460SBarry Smith ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr); 8579f36460SBarry Smith } 86c005e166SLois Curfman McInnes /* Compute Jacobian approximation, 1 column at a time. 8788c956adSLois Curfman McInnes x1 = current iterate, j1a = F(x1) 8888c956adSLois Curfman McInnes x2 = perturbed iterate, j2a = F(x2) 89c005e166SLois Curfman McInnes */ 9039e2f89bSBarry Smith for (i=0; i<N; i++) { 9178b31e54SBarry Smith ierr = VecCopy(x1,x2);CHKERRQ(ierr); 9223242f5aSBarry Smith if (i>= start && i<end) { 932e8a6d31SBarry Smith ierr = VecGetArray(x1,&xx);CHKERRQ(ierr); 9479f36460SBarry Smith if (use_wp) { 9579f36460SBarry Smith dx = 1.0 + unorm; 9679f36460SBarry Smith } else { 9739e2f89bSBarry Smith dx = xx[i-start]; 9879f36460SBarry Smith } 992e8a6d31SBarry Smith ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 100*6bd79f97SJed 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 } 107a86fb38aSBarry Smith ierr = (*eval_fct)(snes,x2,j2a);CHKERRQ(ierr); 10879f36460SBarry Smith ierr = VecAXPY(j2a,-1.0,j1a);CHKERRQ(ierr); 10986c88abdSHong Zhang /* Communicate scale=1/dx_i to all processors */ 11086c88abdSHong Zhang ierr = VecGetOwnershipRanges(x1,&ranges);CHKERRQ(ierr); 11186c88abdSHong Zhang root = size; 11286c88abdSHong Zhang for (j=size-1; j>-1; j--){ 11386c88abdSHong Zhang root--; 11486c88abdSHong Zhang if (i>=ranges[j]) break; 11586c88abdSHong Zhang } 11686c88abdSHong Zhang ierr = MPI_Bcast(&wscale,1,MPIU_SCALAR,root,comm);CHKERRQ(ierr); 11786c88abdSHong Zhang 11886c88abdSHong Zhang ierr = VecScale(j2a,wscale);CHKERRQ(ierr); 1192e8a6d31SBarry Smith ierr = VecNorm(j2a,NORM_INFINITY,&amax);CHKERRQ(ierr); amax *= 1.e-14; 120d142edc7SBarry Smith ierr = VecGetArray(j2a,&y);CHKERRQ(ierr); 12123242f5aSBarry Smith for (j=start; j<end; j++) { 122f6e735bbSHong Zhang if (PetscAbsScalar(y[j-start]) > amax || j == i) { 123b4fc646aSLois Curfman McInnes ierr = MatSetValues(*B,1,&j,1,&i,y+j-start,INSERT_VALUES);CHKERRQ(ierr); 12423242f5aSBarry Smith } 12523242f5aSBarry Smith } 1262e8a6d31SBarry Smith ierr = VecRestoreArray(j2a,&y);CHKERRQ(ierr); 12723242f5aSBarry Smith } 128b4fc646aSLois Curfman McInnes ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 129b4fc646aSLois Curfman McInnes ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 130f588057bSBarry Smith if (*B != *J) { 131f588057bSBarry Smith ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 132f588057bSBarry Smith ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 133f588057bSBarry Smith } 134595b82d2SBarry Smith *flag = DIFFERENT_NONZERO_PATTERN; 1353a40ed3dSBarry Smith PetscFunctionReturn(0); 13611320018SBarry Smith } 13711320018SBarry Smith 138bd45824fSLois Curfman McInnes 139