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 40db781477SPatrick Sanan .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; 450df40c35SBarry Smith PetscInt i,N,start,end,j,value,root,max_funcs = snes->max_funcs; 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; 51ace3abfcSBarry Smith PetscBool assembled,use_wp = PETSC_TRUE,flg; 5279f36460SBarry Smith const char *list[2] = {"ds","wp"}; 5386c88abdSHong Zhang PetscMPIInt size; 5486c88abdSHong Zhang const PetscInt *ranges; 550df40c35SBarry Smith DM dm; 560df40c35SBarry Smith DMSNES dms; 570521c3abSLois Curfman McInnes 583a40ed3dSBarry Smith PetscFunctionBegin; 590df40c35SBarry Smith snes->max_funcs = PETSC_MAX_INT; 605d5a84c9SBarry Smith /* Since this Jacobian will possibly have "extra" nonzero locations just turn off errors for these locations */ 619566063dSJacob Faibussowitsch PetscCall(MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 629566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_test_err",&epsilon,NULL)); 6323242f5aSBarry Smith 649566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)x1,&comm)); 659566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 669566063dSJacob Faibussowitsch PetscCall(MatAssembled(B,&assembled)); 67*1baa6e33SBarry Smith if (assembled) PetscCall(MatZeroEntries(B)); 68aa79bc6dSLois Curfman McInnes if (!snes->nvwork) { 69b1f624c7SBarry Smith if (snes->dm) { 709566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(snes->dm,&j1a)); 719566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(snes->dm,&j2a)); 729566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(snes->dm,&x2)); 73b1f624c7SBarry Smith } else { 74aa79bc6dSLois Curfman McInnes snes->nvwork = 3; 759566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(x1,snes->nvwork,&snes->vwork)); 769566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParents(snes,snes->nvwork,snes->vwork)); 7788c956adSLois Curfman McInnes j1a = snes->vwork[0]; j2a = snes->vwork[1]; x2 = snes->vwork[2]; 78b1f624c7SBarry Smith } 79b1f624c7SBarry Smith } 8023242f5aSBarry Smith 819566063dSJacob Faibussowitsch PetscCall(VecGetSize(x1,&N)); 829566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(x1,&start,&end)); 839566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes,&dm)); 849566063dSJacob Faibussowitsch PetscCall(DMGetDMSNES(dm,&dms)); 850df40c35SBarry Smith if (dms->ops->computemffunction) { 869566063dSJacob Faibussowitsch PetscCall(SNESComputeMFFunction(snes,x1,j1a)); 870df40c35SBarry Smith } else { 889566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes,x1,j1a)); 890df40c35SBarry Smith } 90c005e166SLois Curfman McInnes 91d0609cedSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"Differencing options","SNES"); 929566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-mat_fd_type","Algorithm to compute difference parameter","SNESComputeJacobianDefault",list,2,"wp",&value,&flg)); 93d0609cedSBarry Smith PetscOptionsEnd(); 94f5af7f23SKarl Rupp if (flg && !value) use_wp = PETSC_FALSE; 95f5af7f23SKarl Rupp 9679f36460SBarry Smith if (use_wp) { 979566063dSJacob Faibussowitsch PetscCall(VecNorm(x1,NORM_2,&unorm)); 9879f36460SBarry Smith } 99c005e166SLois Curfman McInnes /* Compute Jacobian approximation, 1 column at a time. 10088c956adSLois Curfman McInnes x1 = current iterate, j1a = F(x1) 10188c956adSLois Curfman McInnes x2 = perturbed iterate, j2a = F(x2) 102c005e166SLois Curfman McInnes */ 10339e2f89bSBarry Smith for (i=0; i<N; i++) { 1049566063dSJacob Faibussowitsch PetscCall(VecCopy(x1,x2)); 10523242f5aSBarry Smith if (i>= start && i<end) { 1069566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x1,&xx)); 1075620d6dcSBarry Smith if (use_wp) dx = PetscSqrtReal(1.0 + unorm); 108f5af7f23SKarl Rupp else dx = xx[i-start]; 1099566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x1,&xx)); 1106bd79f97SJed Brown if (PetscAbsScalar(dx) < dx_min) dx = (PetscRealPart(dx) < 0. ? -1. : 1.) * dx_par; 11139e2f89bSBarry Smith dx *= epsilon; 11274f6f00dSLois Curfman McInnes wscale = 1.0/dx; 1139566063dSJacob Faibussowitsch PetscCall(VecSetValues(x2,1,&i,&dx,ADD_VALUES)); 1146c783aadSBarry Smith } else { 115bbb6d6a8SBarry Smith wscale = 0.0; 116bbb6d6a8SBarry Smith } 1179566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(x2)); 1189566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(x2)); 1190df40c35SBarry Smith if (dms->ops->computemffunction) { 1209566063dSJacob Faibussowitsch PetscCall(SNESComputeMFFunction(snes,x2,j2a)); 1210df40c35SBarry Smith } else { 1229566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes,x2,j2a)); 1230df40c35SBarry Smith } 1249566063dSJacob Faibussowitsch PetscCall(VecAXPY(j2a,-1.0,j1a)); 12586c88abdSHong Zhang /* Communicate scale=1/dx_i to all processors */ 1269566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRanges(x1,&ranges)); 12786c88abdSHong Zhang root = size; 12886c88abdSHong Zhang for (j=size-1; j>-1; j--) { 12986c88abdSHong Zhang root--; 13086c88abdSHong Zhang if (i>=ranges[j]) break; 13186c88abdSHong Zhang } 1329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&wscale,1,MPIU_SCALAR,root,comm)); 1339566063dSJacob Faibussowitsch PetscCall(VecScale(j2a,wscale)); 1349566063dSJacob Faibussowitsch PetscCall(VecNorm(j2a,NORM_INFINITY,&amax)); amax *= 1.e-14; 1359566063dSJacob Faibussowitsch PetscCall(VecGetArray(j2a,&y)); 13623242f5aSBarry Smith for (j=start; j<end; j++) { 137f6e735bbSHong Zhang if (PetscAbsScalar(y[j-start]) > amax || j == i) { 1389566063dSJacob Faibussowitsch PetscCall(MatSetValues(B,1,&j,1,&i,y+j-start,INSERT_VALUES)); 13923242f5aSBarry Smith } 14023242f5aSBarry Smith } 1419566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(j2a,&y)); 14223242f5aSBarry Smith } 143b1f624c7SBarry Smith if (snes->dm) { 1449566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(snes->dm,&j1a)); 1459566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(snes->dm,&j2a)); 1469566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(snes->dm,&x2)); 147b1f624c7SBarry Smith } 1489566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 1499566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 15094ab13aaSBarry Smith if (B != J) { 1519566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY)); 1529566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY)); 153f588057bSBarry Smith } 1540df40c35SBarry Smith snes->max_funcs = max_funcs; 1550df40c35SBarry Smith snes->nfuncs -= N; 1563a40ed3dSBarry Smith PetscFunctionReturn(0); 15711320018SBarry Smith } 158