163dd3a1aSKris Buschelman #define PETSCSNES_DLL 211320018SBarry Smith 3b9147fbbSdalcinl #include "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) 25*3ec795f1SBarry 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; 4779f36460SBarry Smith PetscInt i,N,start,end,j,value; 4879f36460SBarry Smith PetscScalar dx,*y,scale,*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"}; 550521c3abSLois Curfman McInnes 563a40ed3dSBarry Smith PetscFunctionBegin; 5787828ca2SBarry Smith ierr = PetscOptionsGetReal(snes->prefix,"-snes_test_err",&epsilon,0);CHKERRQ(ierr); 584b27c08aSLois Curfman McInnes eval_fct = SNESComputeFunction; 5923242f5aSBarry Smith 602d0c0e3bSBarry Smith ierr = PetscObjectGetComm((PetscObject)x1,&comm);CHKERRQ(ierr); 610b9b6f31SBarry Smith ierr = MatAssembled(*B,&assembled);CHKERRQ(ierr); 620b9b6f31SBarry Smith if (assembled) { 632d0c0e3bSBarry Smith ierr = MatZeroEntries(*B);CHKERRQ(ierr); 640b9b6f31SBarry Smith } 65aa79bc6dSLois Curfman McInnes if (!snes->nvwork) { 66aa79bc6dSLois Curfman McInnes ierr = VecDuplicateVecs(x1,3,&snes->vwork);CHKERRQ(ierr); 67aa79bc6dSLois Curfman McInnes snes->nvwork = 3; 68efee365bSSatish Balay ierr = PetscLogObjectParents(snes,3,snes->vwork);CHKERRQ(ierr); 69aa79bc6dSLois Curfman McInnes } 7088c956adSLois Curfman McInnes j1a = snes->vwork[0]; j2a = snes->vwork[1]; x2 = snes->vwork[2]; 7123242f5aSBarry Smith 7278b31e54SBarry Smith ierr = VecGetSize(x1,&N);CHKERRQ(ierr); 7378b31e54SBarry Smith ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr); 74a86fb38aSBarry Smith ierr = (*eval_fct)(snes,x1,j1a);CHKERRQ(ierr); 75c005e166SLois Curfman McInnes 7679f36460SBarry Smith ierr = PetscOptionsEList("-mat_fd_type","Algorithm to compute difference parameter","SNESDefaultComputeJacobian",list,2,"wp",&value,&flg);CHKERRQ(ierr); 7779f36460SBarry Smith if (flg && !value) { 7879f36460SBarry Smith use_wp = PETSC_FALSE; 7979f36460SBarry Smith } 8079f36460SBarry Smith if (use_wp) { 8179f36460SBarry Smith ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr); 8279f36460SBarry Smith } 83c005e166SLois Curfman McInnes /* Compute Jacobian approximation, 1 column at a time. 8488c956adSLois Curfman McInnes x1 = current iterate, j1a = F(x1) 8588c956adSLois Curfman McInnes x2 = perturbed iterate, j2a = F(x2) 86c005e166SLois Curfman McInnes */ 8739e2f89bSBarry Smith for (i=0; i<N; i++) { 8878b31e54SBarry Smith ierr = VecCopy(x1,x2);CHKERRQ(ierr); 8923242f5aSBarry Smith if (i>= start && i<end) { 902e8a6d31SBarry Smith ierr = VecGetArray(x1,&xx);CHKERRQ(ierr); 9179f36460SBarry Smith if (use_wp) { 9279f36460SBarry Smith dx = 1.0 + unorm; 9379f36460SBarry Smith } else { 9439e2f89bSBarry Smith dx = xx[i-start]; 9579f36460SBarry Smith } 962e8a6d31SBarry Smith ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 97aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX) 9854d73c9fSLois Curfman McInnes if (dx < dx_min && dx >= 0.0) dx = dx_par; 9954d73c9fSLois Curfman McInnes else if (dx < 0.0 && dx > -dx_min) dx = -dx_par; 10019a167f6SBarry Smith #else 101329f5518SBarry Smith if (PetscAbsScalar(dx) < dx_min && PetscRealPart(dx) >= 0.0) dx = dx_par; 102329f5518SBarry Smith else if (PetscRealPart(dx) < 0.0 && PetscAbsScalar(dx) < dx_min) dx = -dx_par; 10319a167f6SBarry Smith #endif 10439e2f89bSBarry Smith dx *= epsilon; 10574f6f00dSLois Curfman McInnes wscale = 1.0/dx; 1062e8a6d31SBarry Smith ierr = VecSetValues(x2,1,&i,&dx,ADD_VALUES);CHKERRQ(ierr); 1076c783aadSBarry Smith } else { 108bbb6d6a8SBarry Smith wscale = 0.0; 109bbb6d6a8SBarry Smith } 110a86fb38aSBarry Smith ierr = (*eval_fct)(snes,x2,j2a);CHKERRQ(ierr); 11179f36460SBarry Smith ierr = VecAXPY(j2a,-1.0,j1a);CHKERRQ(ierr); 112c005e166SLois Curfman McInnes /* Communicate scale to all processors */ 1136831982aSBarry Smith ierr = MPI_Allreduce(&wscale,&scale,1,MPIU_SCALAR,PetscSum_Op,comm);CHKERRQ(ierr); 1142dcb1b2aSMatthew Knepley ierr = VecScale(j2a,scale);CHKERRQ(ierr); 1152e8a6d31SBarry Smith ierr = VecNorm(j2a,NORM_INFINITY,&amax);CHKERRQ(ierr); amax *= 1.e-14; 116d142edc7SBarry Smith ierr = VecGetArray(j2a,&y);CHKERRQ(ierr); 11723242f5aSBarry Smith for (j=start; j<end; j++) { 118cddf8d76SBarry Smith if (PetscAbsScalar(y[j-start]) > amax) { 119b4fc646aSLois Curfman McInnes ierr = MatSetValues(*B,1,&j,1,&i,y+j-start,INSERT_VALUES);CHKERRQ(ierr); 12023242f5aSBarry Smith } 12123242f5aSBarry Smith } 1222e8a6d31SBarry Smith ierr = VecRestoreArray(j2a,&y);CHKERRQ(ierr); 12323242f5aSBarry Smith } 124b4fc646aSLois Curfman McInnes ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 125b4fc646aSLois Curfman McInnes ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 126f588057bSBarry Smith if (*B != *J) { 127f588057bSBarry Smith ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 128f588057bSBarry Smith ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 129f588057bSBarry Smith } 130595b82d2SBarry Smith *flag = DIFFERENT_NONZERO_PATTERN; 1313a40ed3dSBarry Smith PetscFunctionReturn(0); 13211320018SBarry Smith } 13311320018SBarry Smith 134bd45824fSLois Curfman McInnes 135