111320018SBarry Smith 2af0996ceSBarry Smith #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/ 3*b1f624c7SBarry 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 365f3c43d9SLois Curfman McInnes .keywords: SNES, finite differences, Jacobian 375f3c43d9SLois Curfman McInnes 388d359177SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobianDefaultColor(), MatCreateSNESMF() 3911320018SBarry Smith @*/ 40d1e9a80fSBarry Smith PetscErrorCode SNESComputeJacobianDefault(SNES snes,Vec x1,Mat J,Mat B,void *ctx) 4111320018SBarry Smith { 4288c956adSLois Curfman McInnes Vec j1a,j2a,x2; 436849ba73SBarry Smith PetscErrorCode ierr; 4486c88abdSHong Zhang PetscInt i,N,start,end,j,value,root; 455edff71fSBarry Smith PetscScalar dx,*y,wscale; 465edff71fSBarry Smith const PetscScalar *xx; 4777d8c4bbSBarry Smith PetscReal amax,epsilon = PETSC_SQRT_MACHINE_EPSILON; 4879f36460SBarry Smith PetscReal dx_min = 1.e-16,dx_par = 1.e-1,unorm; 49bbb6d6a8SBarry Smith MPI_Comm comm; 50ace3abfcSBarry Smith PetscBool assembled,use_wp = PETSC_TRUE,flg; 5179f36460SBarry Smith const char *list[2] = {"ds","wp"}; 5286c88abdSHong Zhang PetscMPIInt size; 5386c88abdSHong Zhang const PetscInt *ranges; 540521c3abSLois Curfman McInnes 553a40ed3dSBarry Smith PetscFunctionBegin; 565d5a84c9SBarry Smith /* Since this Jacobian will possibly have "extra" nonzero locations just turn off errors for these locations */ 575d5a84c9SBarry Smith ierr = MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE);CHKERRQ(ierr); 58c5929fdfSBarry Smith ierr = PetscOptionsGetReal(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_test_err",&epsilon,0);CHKERRQ(ierr); 5923242f5aSBarry Smith 602d0c0e3bSBarry Smith ierr = PetscObjectGetComm((PetscObject)x1,&comm);CHKERRQ(ierr); 6186c88abdSHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 6294ab13aaSBarry Smith ierr = MatAssembled(B,&assembled);CHKERRQ(ierr); 630b9b6f31SBarry Smith if (assembled) { 6494ab13aaSBarry Smith ierr = MatZeroEntries(B);CHKERRQ(ierr); 650b9b6f31SBarry Smith } 66aa79bc6dSLois Curfman McInnes if (!snes->nvwork) { 67*b1f624c7SBarry Smith if (snes->dm) { 68*b1f624c7SBarry Smith ierr = DMGetGlobalVector(snes->dm,&j1a);CHKERRQ(ierr); 69*b1f624c7SBarry Smith ierr = DMGetGlobalVector(snes->dm,&j2a);CHKERRQ(ierr); 70*b1f624c7SBarry Smith ierr = DMGetGlobalVector(snes->dm,&x2);CHKERRQ(ierr); 71*b1f624c7SBarry Smith } else { 72aa79bc6dSLois Curfman McInnes snes->nvwork = 3; 7385385478SLisandro Dalcin ierr = VecDuplicateVecs(x1,snes->nvwork,&snes->vwork);CHKERRQ(ierr); 7485385478SLisandro Dalcin ierr = PetscLogObjectParents(snes,snes->nvwork,snes->vwork);CHKERRQ(ierr); 7588c956adSLois Curfman McInnes j1a = snes->vwork[0]; j2a = snes->vwork[1]; x2 = snes->vwork[2]; 76*b1f624c7SBarry Smith } 77*b1f624c7SBarry Smith } 7823242f5aSBarry Smith 7978b31e54SBarry Smith ierr = VecGetSize(x1,&N);CHKERRQ(ierr); 8078b31e54SBarry Smith ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr); 8162cd874fSBarry Smith ierr = SNESComputeFunction(snes,x1,j1a);CHKERRQ(ierr); 82c005e166SLois Curfman McInnes 83302440fdSBarry Smith ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"Differencing options","SNES");CHKERRQ(ierr); 848d359177SBarry Smith ierr = PetscOptionsEList("-mat_fd_type","Algorithm to compute difference parameter","SNESComputeJacobianDefault",list,2,"wp",&value,&flg);CHKERRQ(ierr); 85302440fdSBarry Smith ierr = PetscOptionsEnd();CHKERRQ(ierr); 86f5af7f23SKarl Rupp if (flg && !value) use_wp = PETSC_FALSE; 87f5af7f23SKarl Rupp 8879f36460SBarry Smith if (use_wp) { 8979f36460SBarry Smith ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr); 9079f36460SBarry Smith } 91c005e166SLois Curfman McInnes /* Compute Jacobian approximation, 1 column at a time. 9288c956adSLois Curfman McInnes x1 = current iterate, j1a = F(x1) 9388c956adSLois Curfman McInnes x2 = perturbed iterate, j2a = F(x2) 94c005e166SLois Curfman McInnes */ 9539e2f89bSBarry Smith for (i=0; i<N; i++) { 9678b31e54SBarry Smith ierr = VecCopy(x1,x2);CHKERRQ(ierr); 9723242f5aSBarry Smith if (i>= start && i<end) { 985edff71fSBarry Smith ierr = VecGetArrayRead(x1,&xx);CHKERRQ(ierr); 995620d6dcSBarry Smith if (use_wp) dx = PetscSqrtReal(1.0 + unorm); 100f5af7f23SKarl Rupp else dx = xx[i-start]; 1015edff71fSBarry Smith ierr = VecRestoreArrayRead(x1,&xx);CHKERRQ(ierr); 1026bd79f97SJed Brown if (PetscAbsScalar(dx) < dx_min) dx = (PetscRealPart(dx) < 0. ? -1. : 1.) * dx_par; 10339e2f89bSBarry Smith dx *= epsilon; 10474f6f00dSLois Curfman McInnes wscale = 1.0/dx; 1052e8a6d31SBarry Smith ierr = VecSetValues(x2,1,&i,&dx,ADD_VALUES);CHKERRQ(ierr); 1066c783aadSBarry Smith } else { 107bbb6d6a8SBarry Smith wscale = 0.0; 108bbb6d6a8SBarry Smith } 109275089ecSBarry Smith ierr = VecAssemblyBegin(x2);CHKERRQ(ierr); 110275089ecSBarry Smith ierr = VecAssemblyEnd(x2);CHKERRQ(ierr); 11162cd874fSBarry Smith ierr = SNESComputeFunction(snes,x2,j2a);CHKERRQ(ierr); 11279f36460SBarry Smith ierr = VecAXPY(j2a,-1.0,j1a);CHKERRQ(ierr); 11386c88abdSHong Zhang /* Communicate scale=1/dx_i to all processors */ 11486c88abdSHong Zhang ierr = VecGetOwnershipRanges(x1,&ranges);CHKERRQ(ierr); 11586c88abdSHong Zhang root = size; 11686c88abdSHong Zhang for (j=size-1; j>-1; j--) { 11786c88abdSHong Zhang root--; 11886c88abdSHong Zhang if (i>=ranges[j]) break; 11986c88abdSHong Zhang } 12086c88abdSHong Zhang ierr = MPI_Bcast(&wscale,1,MPIU_SCALAR,root,comm);CHKERRQ(ierr); 12186c88abdSHong Zhang 12286c88abdSHong Zhang ierr = VecScale(j2a,wscale);CHKERRQ(ierr); 1232e8a6d31SBarry Smith ierr = VecNorm(j2a,NORM_INFINITY,&amax);CHKERRQ(ierr); amax *= 1.e-14; 124d142edc7SBarry Smith ierr = VecGetArray(j2a,&y);CHKERRQ(ierr); 12523242f5aSBarry Smith for (j=start; j<end; j++) { 126f6e735bbSHong Zhang if (PetscAbsScalar(y[j-start]) > amax || j == i) { 12794ab13aaSBarry Smith ierr = MatSetValues(B,1,&j,1,&i,y+j-start,INSERT_VALUES);CHKERRQ(ierr); 12823242f5aSBarry Smith } 12923242f5aSBarry Smith } 1302e8a6d31SBarry Smith ierr = VecRestoreArray(j2a,&y);CHKERRQ(ierr); 13123242f5aSBarry Smith } 132*b1f624c7SBarry Smith if (snes->dm) { 133*b1f624c7SBarry Smith ierr = DMRestoreGlobalVector(snes->dm,&j1a);CHKERRQ(ierr); 134*b1f624c7SBarry Smith ierr = DMRestoreGlobalVector(snes->dm,&j2a);CHKERRQ(ierr); 135*b1f624c7SBarry Smith ierr = DMRestoreGlobalVector(snes->dm,&x2);CHKERRQ(ierr); 136*b1f624c7SBarry Smith } 13794ab13aaSBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13894ab13aaSBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13994ab13aaSBarry Smith if (B != J) { 14094ab13aaSBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 14194ab13aaSBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 142f588057bSBarry Smith } 1433a40ed3dSBarry Smith PetscFunctionReturn(0); 14411320018SBarry Smith } 14511320018SBarry Smith 146bd45824fSLois Curfman McInnes 147