111320018SBarry Smith 2af0996ceSBarry Smith #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/ 3b1f624c7SBarry 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: 116b867d5aSJose E. Roman + snes - the SNES context 126b867d5aSJose E. Roman . 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 350df40c35SBarry Smith This routine ignores the maximum number of function evaluations set with SNESSetTolerances() and the function 360df40c35SBarry Smith evaluations it performs are not counted in what is returned by of SNESGetNumberFunctionEvals(). 370df40c35SBarry Smith 3836851e7fSLois Curfman McInnes Level: intermediate 3936851e7fSLois Curfman McInnes 408d359177SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobianDefaultColor(), MatCreateSNESMF() 4111320018SBarry Smith @*/ 42d1e9a80fSBarry Smith PetscErrorCode SNESComputeJacobianDefault(SNES snes,Vec x1,Mat J,Mat B,void *ctx) 4311320018SBarry Smith { 4488c956adSLois Curfman McInnes Vec j1a,j2a,x2; 456849ba73SBarry Smith PetscErrorCode ierr; 460df40c35SBarry Smith PetscInt i,N,start,end,j,value,root,max_funcs = snes->max_funcs; 475edff71fSBarry Smith PetscScalar dx,*y,wscale; 485edff71fSBarry Smith const PetscScalar *xx; 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; 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; 560df40c35SBarry Smith DM dm; 570df40c35SBarry Smith DMSNES dms; 580521c3abSLois Curfman McInnes 593a40ed3dSBarry Smith PetscFunctionBegin; 600df40c35SBarry Smith snes->max_funcs = PETSC_MAX_INT; 615d5a84c9SBarry Smith /* Since this Jacobian will possibly have "extra" nonzero locations just turn off errors for these locations */ 62*9566063dSJacob Faibussowitsch PetscCall(MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 63*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_test_err",&epsilon,NULL)); 6423242f5aSBarry Smith 65*9566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)x1,&comm)); 66*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 67*9566063dSJacob Faibussowitsch PetscCall(MatAssembled(B,&assembled)); 680b9b6f31SBarry Smith if (assembled) { 69*9566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(B)); 700b9b6f31SBarry Smith } 71aa79bc6dSLois Curfman McInnes if (!snes->nvwork) { 72b1f624c7SBarry Smith if (snes->dm) { 73*9566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(snes->dm,&j1a)); 74*9566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(snes->dm,&j2a)); 75*9566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(snes->dm,&x2)); 76b1f624c7SBarry Smith } else { 77aa79bc6dSLois Curfman McInnes snes->nvwork = 3; 78*9566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(x1,snes->nvwork,&snes->vwork)); 79*9566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParents(snes,snes->nvwork,snes->vwork)); 8088c956adSLois Curfman McInnes j1a = snes->vwork[0]; j2a = snes->vwork[1]; x2 = snes->vwork[2]; 81b1f624c7SBarry Smith } 82b1f624c7SBarry Smith } 8323242f5aSBarry Smith 84*9566063dSJacob Faibussowitsch PetscCall(VecGetSize(x1,&N)); 85*9566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(x1,&start,&end)); 86*9566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes,&dm)); 87*9566063dSJacob Faibussowitsch PetscCall(DMGetDMSNES(dm,&dms)); 880df40c35SBarry Smith if (dms->ops->computemffunction) { 89*9566063dSJacob Faibussowitsch PetscCall(SNESComputeMFFunction(snes,x1,j1a)); 900df40c35SBarry Smith } else { 91*9566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes,x1,j1a)); 920df40c35SBarry Smith } 93c005e166SLois Curfman McInnes 94*9566063dSJacob Faibussowitsch ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"Differencing options","SNES");PetscCall(ierr); 95*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-mat_fd_type","Algorithm to compute difference parameter","SNESComputeJacobianDefault",list,2,"wp",&value,&flg)); 96*9566063dSJacob Faibussowitsch ierr = PetscOptionsEnd();PetscCall(ierr); 97f5af7f23SKarl Rupp if (flg && !value) use_wp = PETSC_FALSE; 98f5af7f23SKarl Rupp 9979f36460SBarry Smith if (use_wp) { 100*9566063dSJacob Faibussowitsch PetscCall(VecNorm(x1,NORM_2,&unorm)); 10179f36460SBarry Smith } 102c005e166SLois Curfman McInnes /* Compute Jacobian approximation, 1 column at a time. 10388c956adSLois Curfman McInnes x1 = current iterate, j1a = F(x1) 10488c956adSLois Curfman McInnes x2 = perturbed iterate, j2a = F(x2) 105c005e166SLois Curfman McInnes */ 10639e2f89bSBarry Smith for (i=0; i<N; i++) { 107*9566063dSJacob Faibussowitsch PetscCall(VecCopy(x1,x2)); 10823242f5aSBarry Smith if (i>= start && i<end) { 109*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x1,&xx)); 1105620d6dcSBarry Smith if (use_wp) dx = PetscSqrtReal(1.0 + unorm); 111f5af7f23SKarl Rupp else dx = xx[i-start]; 112*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x1,&xx)); 1136bd79f97SJed Brown if (PetscAbsScalar(dx) < dx_min) dx = (PetscRealPart(dx) < 0. ? -1. : 1.) * dx_par; 11439e2f89bSBarry Smith dx *= epsilon; 11574f6f00dSLois Curfman McInnes wscale = 1.0/dx; 116*9566063dSJacob Faibussowitsch PetscCall(VecSetValues(x2,1,&i,&dx,ADD_VALUES)); 1176c783aadSBarry Smith } else { 118bbb6d6a8SBarry Smith wscale = 0.0; 119bbb6d6a8SBarry Smith } 120*9566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(x2)); 121*9566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(x2)); 1220df40c35SBarry Smith if (dms->ops->computemffunction) { 123*9566063dSJacob Faibussowitsch PetscCall(SNESComputeMFFunction(snes,x2,j2a)); 1240df40c35SBarry Smith } else { 125*9566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes,x2,j2a)); 1260df40c35SBarry Smith } 127*9566063dSJacob Faibussowitsch PetscCall(VecAXPY(j2a,-1.0,j1a)); 12886c88abdSHong Zhang /* Communicate scale=1/dx_i to all processors */ 129*9566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRanges(x1,&ranges)); 13086c88abdSHong Zhang root = size; 13186c88abdSHong Zhang for (j=size-1; j>-1; j--) { 13286c88abdSHong Zhang root--; 13386c88abdSHong Zhang if (i>=ranges[j]) break; 13486c88abdSHong Zhang } 135*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&wscale,1,MPIU_SCALAR,root,comm)); 136*9566063dSJacob Faibussowitsch PetscCall(VecScale(j2a,wscale)); 137*9566063dSJacob Faibussowitsch PetscCall(VecNorm(j2a,NORM_INFINITY,&amax)); amax *= 1.e-14; 138*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(j2a,&y)); 13923242f5aSBarry Smith for (j=start; j<end; j++) { 140f6e735bbSHong Zhang if (PetscAbsScalar(y[j-start]) > amax || j == i) { 141*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(B,1,&j,1,&i,y+j-start,INSERT_VALUES)); 14223242f5aSBarry Smith } 14323242f5aSBarry Smith } 144*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(j2a,&y)); 14523242f5aSBarry Smith } 146b1f624c7SBarry Smith if (snes->dm) { 147*9566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(snes->dm,&j1a)); 148*9566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(snes->dm,&j2a)); 149*9566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(snes->dm,&x2)); 150b1f624c7SBarry Smith } 151*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 152*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 15394ab13aaSBarry Smith if (B != J) { 154*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY)); 155*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY)); 156f588057bSBarry Smith } 1570df40c35SBarry Smith snes->max_funcs = max_funcs; 1580df40c35SBarry Smith snes->nfuncs -= N; 1593a40ed3dSBarry Smith PetscFunctionReturn(0); 16011320018SBarry Smith } 161