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: 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 34*0df40c35SBarry Smith This routine ignores the maximum number of function evaluations set with SNESSetTolerances() and the function 35*0df40c35SBarry Smith evaluations it performs are not counted in what is returned by of SNESGetNumberFunctionEvals(). 36*0df40c35SBarry Smith 3736851e7fSLois Curfman McInnes Level: intermediate 3836851e7fSLois 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; 45*0df40c35SBarry 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; 55*0df40c35SBarry Smith DM dm; 56*0df40c35SBarry Smith DMSNES dms; 570521c3abSLois Curfman McInnes 583a40ed3dSBarry Smith PetscFunctionBegin; 59*0df40c35SBarry 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 */ 615d5a84c9SBarry Smith ierr = MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE);CHKERRQ(ierr); 629e5d0892SLisandro Dalcin ierr = PetscOptionsGetReal(((PetscObject)snes)->options,((PetscObject)snes)->prefix,"-snes_test_err",&epsilon,NULL);CHKERRQ(ierr); 6323242f5aSBarry Smith 642d0c0e3bSBarry Smith ierr = PetscObjectGetComm((PetscObject)x1,&comm);CHKERRQ(ierr); 65ffc4695bSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr); 6694ab13aaSBarry Smith ierr = MatAssembled(B,&assembled);CHKERRQ(ierr); 670b9b6f31SBarry Smith if (assembled) { 6894ab13aaSBarry Smith ierr = MatZeroEntries(B);CHKERRQ(ierr); 690b9b6f31SBarry Smith } 70aa79bc6dSLois Curfman McInnes if (!snes->nvwork) { 71b1f624c7SBarry Smith if (snes->dm) { 72b1f624c7SBarry Smith ierr = DMGetGlobalVector(snes->dm,&j1a);CHKERRQ(ierr); 73b1f624c7SBarry Smith ierr = DMGetGlobalVector(snes->dm,&j2a);CHKERRQ(ierr); 74b1f624c7SBarry Smith ierr = DMGetGlobalVector(snes->dm,&x2);CHKERRQ(ierr); 75b1f624c7SBarry Smith } else { 76aa79bc6dSLois Curfman McInnes snes->nvwork = 3; 7785385478SLisandro Dalcin ierr = VecDuplicateVecs(x1,snes->nvwork,&snes->vwork);CHKERRQ(ierr); 7885385478SLisandro Dalcin ierr = PetscLogObjectParents(snes,snes->nvwork,snes->vwork);CHKERRQ(ierr); 7988c956adSLois Curfman McInnes j1a = snes->vwork[0]; j2a = snes->vwork[1]; x2 = snes->vwork[2]; 80b1f624c7SBarry Smith } 81b1f624c7SBarry Smith } 8223242f5aSBarry Smith 8378b31e54SBarry Smith ierr = VecGetSize(x1,&N);CHKERRQ(ierr); 8478b31e54SBarry Smith ierr = VecGetOwnershipRange(x1,&start,&end);CHKERRQ(ierr); 85*0df40c35SBarry Smith ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 86*0df40c35SBarry Smith ierr = DMGetDMSNES(dm,&dms);CHKERRQ(ierr); 87*0df40c35SBarry Smith if (dms->ops->computemffunction) { 88*0df40c35SBarry Smith ierr = SNESComputeMFFunction(snes,x1,j1a);CHKERRQ(ierr); 89*0df40c35SBarry Smith } else { 9062cd874fSBarry Smith ierr = SNESComputeFunction(snes,x1,j1a);CHKERRQ(ierr); 91*0df40c35SBarry Smith } 92c005e166SLois Curfman McInnes 93302440fdSBarry Smith ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"Differencing options","SNES");CHKERRQ(ierr); 948d359177SBarry Smith ierr = PetscOptionsEList("-mat_fd_type","Algorithm to compute difference parameter","SNESComputeJacobianDefault",list,2,"wp",&value,&flg);CHKERRQ(ierr); 95302440fdSBarry Smith ierr = PetscOptionsEnd();CHKERRQ(ierr); 96f5af7f23SKarl Rupp if (flg && !value) use_wp = PETSC_FALSE; 97f5af7f23SKarl Rupp 9879f36460SBarry Smith if (use_wp) { 9979f36460SBarry Smith ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr); 10079f36460SBarry Smith } 101c005e166SLois Curfman McInnes /* Compute Jacobian approximation, 1 column at a time. 10288c956adSLois Curfman McInnes x1 = current iterate, j1a = F(x1) 10388c956adSLois Curfman McInnes x2 = perturbed iterate, j2a = F(x2) 104c005e166SLois Curfman McInnes */ 10539e2f89bSBarry Smith for (i=0; i<N; i++) { 10678b31e54SBarry Smith ierr = VecCopy(x1,x2);CHKERRQ(ierr); 10723242f5aSBarry Smith if (i>= start && i<end) { 1085edff71fSBarry Smith ierr = VecGetArrayRead(x1,&xx);CHKERRQ(ierr); 1095620d6dcSBarry Smith if (use_wp) dx = PetscSqrtReal(1.0 + unorm); 110f5af7f23SKarl Rupp else dx = xx[i-start]; 1115edff71fSBarry Smith ierr = VecRestoreArrayRead(x1,&xx);CHKERRQ(ierr); 1126bd79f97SJed Brown if (PetscAbsScalar(dx) < dx_min) dx = (PetscRealPart(dx) < 0. ? -1. : 1.) * dx_par; 11339e2f89bSBarry Smith dx *= epsilon; 11474f6f00dSLois Curfman McInnes wscale = 1.0/dx; 1152e8a6d31SBarry Smith ierr = VecSetValues(x2,1,&i,&dx,ADD_VALUES);CHKERRQ(ierr); 1166c783aadSBarry Smith } else { 117bbb6d6a8SBarry Smith wscale = 0.0; 118bbb6d6a8SBarry Smith } 119275089ecSBarry Smith ierr = VecAssemblyBegin(x2);CHKERRQ(ierr); 120275089ecSBarry Smith ierr = VecAssemblyEnd(x2);CHKERRQ(ierr); 121*0df40c35SBarry Smith if (dms->ops->computemffunction) { 122*0df40c35SBarry Smith ierr = SNESComputeMFFunction(snes,x2,j2a);CHKERRQ(ierr); 123*0df40c35SBarry Smith } else { 12462cd874fSBarry Smith ierr = SNESComputeFunction(snes,x2,j2a);CHKERRQ(ierr); 125*0df40c35SBarry Smith } 12679f36460SBarry Smith ierr = VecAXPY(j2a,-1.0,j1a);CHKERRQ(ierr); 12786c88abdSHong Zhang /* Communicate scale=1/dx_i to all processors */ 12886c88abdSHong Zhang ierr = VecGetOwnershipRanges(x1,&ranges);CHKERRQ(ierr); 12986c88abdSHong Zhang root = size; 13086c88abdSHong Zhang for (j=size-1; j>-1; j--) { 13186c88abdSHong Zhang root--; 13286c88abdSHong Zhang if (i>=ranges[j]) break; 13386c88abdSHong Zhang } 134ffc4695bSBarry Smith ierr = MPI_Bcast(&wscale,1,MPIU_SCALAR,root,comm);CHKERRMPI(ierr); 13586c88abdSHong Zhang ierr = VecScale(j2a,wscale);CHKERRQ(ierr); 1362e8a6d31SBarry Smith ierr = VecNorm(j2a,NORM_INFINITY,&amax);CHKERRQ(ierr); amax *= 1.e-14; 137d142edc7SBarry Smith ierr = VecGetArray(j2a,&y);CHKERRQ(ierr); 13823242f5aSBarry Smith for (j=start; j<end; j++) { 139f6e735bbSHong Zhang if (PetscAbsScalar(y[j-start]) > amax || j == i) { 14094ab13aaSBarry Smith ierr = MatSetValues(B,1,&j,1,&i,y+j-start,INSERT_VALUES);CHKERRQ(ierr); 14123242f5aSBarry Smith } 14223242f5aSBarry Smith } 1432e8a6d31SBarry Smith ierr = VecRestoreArray(j2a,&y);CHKERRQ(ierr); 14423242f5aSBarry Smith } 145b1f624c7SBarry Smith if (snes->dm) { 146b1f624c7SBarry Smith ierr = DMRestoreGlobalVector(snes->dm,&j1a);CHKERRQ(ierr); 147b1f624c7SBarry Smith ierr = DMRestoreGlobalVector(snes->dm,&j2a);CHKERRQ(ierr); 148b1f624c7SBarry Smith ierr = DMRestoreGlobalVector(snes->dm,&x2);CHKERRQ(ierr); 149b1f624c7SBarry Smith } 15094ab13aaSBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 15194ab13aaSBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 15294ab13aaSBarry Smith if (B != J) { 15394ab13aaSBarry Smith ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 15494ab13aaSBarry Smith ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 155f588057bSBarry Smith } 156*0df40c35SBarry Smith snes->max_funcs = max_funcs; 157*0df40c35SBarry Smith snes->nfuncs -= N; 1583a40ed3dSBarry Smith PetscFunctionReturn(0); 15911320018SBarry Smith } 16011320018SBarry Smith 161bd45824fSLois Curfman McInnes 162