111320018SBarry Smith 2*b45d2f2cSJed Brown #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/ 311320018SBarry Smith 44a2ae208SSatish Balay #undef __FUNCT__ 54a2ae208SSatish Balay #define __FUNCT__ "SNESDefaultComputeJacobian" 64b828684SBarry Smith /*@C 784a7bf8dSLois Curfman McInnes SNESDefaultComputeJacobian - Computes the Jacobian using finite differences. 811320018SBarry Smith 9fee21e36SBarry Smith Collective on SNES 10fee21e36SBarry Smith 11c7afd0dbSLois Curfman McInnes Input Parameters: 12c7afd0dbSLois Curfman McInnes + 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) 17c7afd0dbSLois Curfman McInnes . B - newly computed Jacobian matrix to use with preconditioner (generally the same as J) 18c7afd0dbSLois Curfman McInnes - flag - flag indicating whether the matrix sparsity structure has changed 19c7afd0dbSLois Curfman McInnes 20ad960d00SLois Curfman McInnes Options Database Key: 2108f75eefSBarry Smith + -snes_fd - Activates SNESDefaultComputeJacobian() 2279f36460SBarry Smith . -snes_test_err - Square root of function error tolerance, default square root of machine 2377d8c4bbSBarry Smith epsilon (1.e-8 in double, 3.e-4 in single) 243ec795f1SBarry Smith - -mat_fd_type - Either wp or ds (see MATMFFD_WP or MATMFFD_DS) 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 3379f36460SBarry Smith An alternative routine that uses coloring to exploit 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 4079f36460SBarry Smith .seealso: SNESSetJacobian(), SNESDefaultComputeJacobianColor(), MatCreateSNESMF() 4111320018SBarry Smith @*/ 427087cfbeSBarry Smith PetscErrorCode 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; 4686c88abdSHong Zhang PetscInt i,N,start,end,j,value,root; 4786c88abdSHong Zhang PetscScalar dx,*y,*xx,wscale; 4877d8c4bbSBarry Smith PetscReal amax,epsilon = PETSC_SQRT_MACHINE_EPSILON; 4979f36460SBarry Smith PetscReal dx_min = 1.e-16,dx_par = 1.e-1,unorm; 50bbb6d6a8SBarry Smith MPI_Comm comm; 516849ba73SBarry Smith PetscErrorCode (*eval_fct)(SNES,Vec,Vec)=0; 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; 560521c3abSLois Curfman McInnes 573a40ed3dSBarry Smith PetscFunctionBegin; 587adad957SLisandro Dalcin ierr = PetscOptionsGetReal(((PetscObject)snes)->prefix,"-snes_test_err",&epsilon,0);CHKERRQ(ierr); 594b27c08aSLois Curfman McInnes eval_fct = SNESComputeFunction; 6023242f5aSBarry Smith 612d0c0e3bSBarry Smith ierr = PetscObjectGetComm((PetscObject)x1,&comm);CHKERRQ(ierr); 6286c88abdSHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 630b9b6f31SBarry Smith ierr = MatAssembled(*B,&assembled);CHKERRQ(ierr); 640b9b6f31SBarry Smith if (assembled) { 652d0c0e3bSBarry Smith ierr = MatZeroEntries(*B);CHKERRQ(ierr); 660b9b6f31SBarry Smith } 67aa79bc6dSLois Curfman McInnes if (!snes->nvwork) { 68aa79bc6dSLois Curfman McInnes snes->nvwork = 3; 6985385478SLisandro Dalcin ierr = VecDuplicateVecs(x1,snes->nvwork,&snes->vwork);CHKERRQ(ierr); 7085385478SLisandro Dalcin ierr = PetscLogObjectParents(snes,snes->nvwork,snes->vwork);CHKERRQ(ierr); 71aa79bc6dSLois Curfman McInnes } 7288c956adSLois Curfman McInnes j1a = snes->vwork[0]; j2a = snes->vwork[1]; x2 = snes->vwork[2]; 7323242f5aSBarry Smith 7478b31e54SBarry Smith ierr = VecGetSize(x1,&N);CHKERRQ(ierr); 7578b31e54SBarry Smith ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr); 76a86fb38aSBarry Smith ierr = (*eval_fct)(snes,x1,j1a);CHKERRQ(ierr); 77c005e166SLois Curfman McInnes 7879f36460SBarry Smith ierr = PetscOptionsEList("-mat_fd_type","Algorithm to compute difference parameter","SNESDefaultComputeJacobian",list,2,"wp",&value,&flg);CHKERRQ(ierr); 7979f36460SBarry Smith if (flg && !value) { 8079f36460SBarry Smith use_wp = PETSC_FALSE; 8179f36460SBarry Smith } 8279f36460SBarry Smith if (use_wp) { 8379f36460SBarry Smith ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr); 8479f36460SBarry Smith } 85c005e166SLois Curfman McInnes /* Compute Jacobian approximation, 1 column at a time. 8688c956adSLois Curfman McInnes x1 = current iterate, j1a = F(x1) 8788c956adSLois Curfman McInnes x2 = perturbed iterate, j2a = F(x2) 88c005e166SLois Curfman McInnes */ 8939e2f89bSBarry Smith for (i=0; i<N; i++) { 9078b31e54SBarry Smith ierr = VecCopy(x1,x2);CHKERRQ(ierr); 9123242f5aSBarry Smith if (i>= start && i<end) { 922e8a6d31SBarry Smith ierr = VecGetArray(x1,&xx);CHKERRQ(ierr); 9379f36460SBarry Smith if (use_wp) { 9479f36460SBarry Smith dx = 1.0 + unorm; 9579f36460SBarry Smith } else { 9639e2f89bSBarry Smith dx = xx[i-start]; 9779f36460SBarry Smith } 982e8a6d31SBarry Smith ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); 996bd79f97SJed Brown if (PetscAbsScalar(dx) < dx_min) dx = (PetscRealPart(dx) < 0. ? -1. : 1.) * dx_par; 10039e2f89bSBarry Smith dx *= epsilon; 10174f6f00dSLois Curfman McInnes wscale = 1.0/dx; 1022e8a6d31SBarry Smith ierr = VecSetValues(x2,1,&i,&dx,ADD_VALUES);CHKERRQ(ierr); 1036c783aadSBarry Smith } else { 104bbb6d6a8SBarry Smith wscale = 0.0; 105bbb6d6a8SBarry Smith } 106275089ecSBarry Smith ierr = VecAssemblyBegin(x2);CHKERRQ(ierr); 107275089ecSBarry Smith ierr = VecAssemblyEnd(x2);CHKERRQ(ierr); 108a86fb38aSBarry Smith ierr = (*eval_fct)(snes,x2,j2a);CHKERRQ(ierr); 10979f36460SBarry Smith ierr = VecAXPY(j2a,-1.0,j1a);CHKERRQ(ierr); 11086c88abdSHong Zhang /* Communicate scale=1/dx_i to all processors */ 11186c88abdSHong Zhang ierr = VecGetOwnershipRanges(x1,&ranges);CHKERRQ(ierr); 11286c88abdSHong Zhang root = size; 11386c88abdSHong Zhang for (j=size-1; j>-1; j--){ 11486c88abdSHong Zhang root--; 11586c88abdSHong Zhang if (i>=ranges[j]) break; 11686c88abdSHong Zhang } 11786c88abdSHong Zhang ierr = MPI_Bcast(&wscale,1,MPIU_SCALAR,root,comm);CHKERRQ(ierr); 11886c88abdSHong Zhang 11986c88abdSHong Zhang ierr = VecScale(j2a,wscale);CHKERRQ(ierr); 1202e8a6d31SBarry Smith ierr = VecNorm(j2a,NORM_INFINITY,&amax);CHKERRQ(ierr); amax *= 1.e-14; 121d142edc7SBarry Smith ierr = VecGetArray(j2a,&y);CHKERRQ(ierr); 12223242f5aSBarry Smith for (j=start; j<end; j++) { 123f6e735bbSHong Zhang if (PetscAbsScalar(y[j-start]) > amax || j == i) { 124b4fc646aSLois Curfman McInnes ierr = MatSetValues(*B,1,&j,1,&i,y+j-start,INSERT_VALUES);CHKERRQ(ierr); 12523242f5aSBarry Smith } 12623242f5aSBarry Smith } 1272e8a6d31SBarry Smith ierr = VecRestoreArray(j2a,&y);CHKERRQ(ierr); 12823242f5aSBarry Smith } 129b4fc646aSLois Curfman McInnes ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 130b4fc646aSLois Curfman McInnes ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 131f588057bSBarry Smith if (*B != *J) { 132f588057bSBarry Smith ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 133f588057bSBarry Smith ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 134f588057bSBarry Smith } 135595b82d2SBarry Smith *flag = DIFFERENT_NONZERO_PATTERN; 1363a40ed3dSBarry Smith PetscFunctionReturn(0); 13711320018SBarry Smith } 13811320018SBarry Smith 139bd45824fSLois Curfman McInnes 140