111320018SBarry Smith 2b45d2f2cSJed Brown #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/ 311320018SBarry Smith 44a2ae208SSatish Balay #undef __FUNCT__ 58d359177SBarry Smith #define __FUNCT__ "SNESComputeJacobianDefault" 64b828684SBarry Smith /*@C 78d359177SBarry Smith SNESComputeJacobianDefault - 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) 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 3536851e7fSLois Curfman McInnes Level: intermediate 3636851e7fSLois Curfman McInnes 375f3c43d9SLois Curfman McInnes .keywords: SNES, finite differences, Jacobian 385f3c43d9SLois Curfman McInnes 398d359177SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobianDefaultColor(), MatCreateSNESMF() 4011320018SBarry Smith @*/ 41d1e9a80fSBarry Smith PetscErrorCode SNESComputeJacobianDefault(SNES snes,Vec x1,Mat J,Mat B,void *ctx) 4211320018SBarry Smith { 4388c956adSLois Curfman McInnes Vec j1a,j2a,x2; 446849ba73SBarry Smith PetscErrorCode ierr; 4586c88abdSHong Zhang PetscInt i,N,start,end,j,value,root; 465edff71fSBarry Smith PetscScalar dx,*y,wscale; 475edff71fSBarry Smith const PetscScalar *xx; 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); 6394ab13aaSBarry Smith ierr = MatAssembled(B,&assembled);CHKERRQ(ierr); 640b9b6f31SBarry Smith if (assembled) { 6594ab13aaSBarry Smith ierr = MatZeroEntries(B);CHKERRQ(ierr); 660b9b6f31SBarry Smith } 67aa79bc6dSLois Curfman McInnes if (!snes->nvwork) { 68aa79bc6dSLois Curfman McInnes snes->nvwork = 3; 69f5af7f23SKarl Rupp 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 79e55864a3SBarry Smith ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"Differencing options","SNES"); 808d359177SBarry Smith ierr = PetscOptionsEList("-mat_fd_type","Algorithm to compute difference parameter","SNESComputeJacobianDefault",list,2,"wp",&value,&flg);CHKERRQ(ierr); 81e55864a3SBarry Smith ierr = PetscOptionsEnd(); 82f5af7f23SKarl Rupp if (flg && !value) use_wp = PETSC_FALSE; 83f5af7f23SKarl Rupp 8479f36460SBarry Smith if (use_wp) { 8579f36460SBarry Smith ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr); 8679f36460SBarry Smith } 87c005e166SLois Curfman McInnes /* Compute Jacobian approximation, 1 column at a time. 8888c956adSLois Curfman McInnes x1 = current iterate, j1a = F(x1) 8988c956adSLois Curfman McInnes x2 = perturbed iterate, j2a = F(x2) 90c005e166SLois Curfman McInnes */ 9139e2f89bSBarry Smith for (i=0; i<N; i++) { 9278b31e54SBarry Smith ierr = VecCopy(x1,x2);CHKERRQ(ierr); 9323242f5aSBarry Smith if (i>= start && i<end) { 945edff71fSBarry Smith ierr = VecGetArrayRead(x1,&xx);CHKERRQ(ierr); 95*5620d6dcSBarry Smith if (use_wp) dx = PetscSqrtReal(1.0 + unorm); 96f5af7f23SKarl Rupp else dx = xx[i-start]; 975edff71fSBarry Smith ierr = VecRestoreArrayRead(x1,&xx);CHKERRQ(ierr); 986bd79f97SJed Brown if (PetscAbsScalar(dx) < dx_min) dx = (PetscRealPart(dx) < 0. ? -1. : 1.) * dx_par; 9939e2f89bSBarry Smith dx *= epsilon; 10074f6f00dSLois Curfman McInnes wscale = 1.0/dx; 1012e8a6d31SBarry Smith ierr = VecSetValues(x2,1,&i,&dx,ADD_VALUES);CHKERRQ(ierr); 1026c783aadSBarry Smith } else { 103bbb6d6a8SBarry Smith wscale = 0.0; 104bbb6d6a8SBarry Smith } 105275089ecSBarry Smith ierr = VecAssemblyBegin(x2);CHKERRQ(ierr); 106275089ecSBarry Smith ierr = VecAssemblyEnd(x2);CHKERRQ(ierr); 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) { 12394ab13aaSBarry Smith 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 } 12894ab13aaSBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 12994ab13aaSBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13094ab13aaSBarry Smith if (B != J) { 13194ab13aaSBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13294ab13aaSBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 133f588057bSBarry Smith } 1343a40ed3dSBarry Smith PetscFunctionReturn(0); 13511320018SBarry Smith } 13611320018SBarry Smith 137bd45824fSLois Curfman McInnes 138