1af0996ceSBarry Smith #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/ 2f2bd0052SStefano Zampini #include <petsc/private/vecimpl.h> /* for Vec->ops->setvalues */ 3b1f624c7SBarry Smith #include <petscdm.h> 411320018SBarry Smith 55d83a8b1SBarry Smith /*@ 68d359177SBarry Smith SNESComputeJacobianDefault - Computes the Jacobian using finite differences. 711320018SBarry Smith 8c3339decSBarry Smith Collective 9fee21e36SBarry Smith 10c7afd0dbSLois Curfman McInnes Input Parameters: 11f6dfbefdSBarry Smith + snes - the `SNES` context 126b867d5aSJose E. Roman . x1 - compute Jacobian at this point 13f6dfbefdSBarry Smith - 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) 17dc4c0fb0SBarry Smith - B - newly computed Jacobian matrix to use with preconditioner (generally the same as `J`) 18c7afd0dbSLois Curfman McInnes 19f6dfbefdSBarry Smith Options Database Keys: 20f6dfbefdSBarry Smith + -snes_fd - Activates `SNESComputeJacobianDefault()` 21dc4c0fb0SBarry Smith . -snes_fd_coloring - Activates a faster computation that uses a graph coloring of the matrix 2279f36460SBarry Smith . -snes_test_err - Square root of function error tolerance, default square root of machine 2377d8c4bbSBarry Smith epsilon (1.e-8 in double, 3.e-4 in single) 24f6dfbefdSBarry Smith - -mat_fd_type - Either wp or ds (see `MATMFFD_WP` or `MATMFFD_DS`) 25ad960d00SLois Curfman McInnes 26dc4c0fb0SBarry Smith Level: intermediate 27dc4c0fb0SBarry Smith 285f3c43d9SLois Curfman McInnes Notes: 295f3c43d9SLois Curfman McInnes This routine is slow and expensive, and is not currently optimized 305f3c43d9SLois Curfman McInnes to take advantage of sparsity in the problem. Although 31f6dfbefdSBarry Smith `SNESComputeJacobianDefault()` is not recommended for general use 325f3c43d9SLois Curfman McInnes in large-scale applications, It can be useful in checking the 335f3c43d9SLois Curfman McInnes correctness of a user-provided Jacobian. 3411320018SBarry Smith 3579f36460SBarry Smith An alternative routine that uses coloring to exploit matrix sparsity is 36f6dfbefdSBarry Smith `SNESComputeJacobianDefaultColor()`. 37b4fc646aSLois Curfman McInnes 38f6dfbefdSBarry Smith This routine ignores the maximum number of function evaluations set with `SNESSetTolerances()` and the function 39f6dfbefdSBarry Smith evaluations it performs are not counted in what is returned by of `SNESGetNumberFunctionEvals()`. 40f6dfbefdSBarry Smith 41f6dfbefdSBarry Smith This function can be provided to `SNESSetJacobian()` along with a dense matrix to hold the Jacobian 420df40c35SBarry Smith 43420bcc1bSBarry Smith Developer Note: 44420bcc1bSBarry Smith The function has a poorly chosen name since it does not mention the use of finite differences 45420bcc1bSBarry Smith 46420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSetJacobian()`, `SNESComputeJacobianDefaultColor()`, `MatCreateSNESMF()` 4711320018SBarry Smith @*/ 48*2a8381b2SBarry Smith PetscErrorCode SNESComputeJacobianDefault(SNES snes, Vec x1, Mat J, Mat B, PetscCtx ctx) 49d71ae5a4SJacob Faibussowitsch { 5088c956adSLois Curfman McInnes Vec j1a, j2a, x2; 516497c311SBarry Smith PetscInt i, N, start, end, j, value, max_funcs = snes->max_funcs; 525edff71fSBarry Smith PetscScalar dx, *y, wscale; 535edff71fSBarry Smith const PetscScalar *xx; 5477d8c4bbSBarry Smith PetscReal amax, epsilon = PETSC_SQRT_MACHINE_EPSILON; 5579f36460SBarry Smith PetscReal dx_min = 1.e-16, dx_par = 1.e-1, unorm; 56bbb6d6a8SBarry Smith MPI_Comm comm; 57ace3abfcSBarry Smith PetscBool assembled, use_wp = PETSC_TRUE, flg; 5879f36460SBarry Smith const char *list[2] = {"ds", "wp"}; 596497c311SBarry Smith PetscMPIInt size, root; 6086c88abdSHong Zhang const PetscInt *ranges; 610df40c35SBarry Smith DM dm; 620df40c35SBarry Smith DMSNES dms; 630521c3abSLois Curfman McInnes 643a40ed3dSBarry Smith PetscFunctionBegin; 651690c2aeSBarry Smith snes->max_funcs = PETSC_INT_MAX; 665d5a84c9SBarry Smith /* Since this Jacobian will possibly have "extra" nonzero locations just turn off errors for these locations */ 679566063dSJacob Faibussowitsch PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 689566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(((PetscObject)snes)->options, ((PetscObject)snes)->prefix, "-snes_test_err", &epsilon, NULL)); 6923242f5aSBarry Smith 709566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)x1, &comm)); 719566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 729566063dSJacob Faibussowitsch PetscCall(MatAssembled(B, &assembled)); 731baa6e33SBarry Smith if (assembled) PetscCall(MatZeroEntries(B)); 74aa79bc6dSLois Curfman McInnes if (!snes->nvwork) { 75b1f624c7SBarry Smith if (snes->dm) { 769566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(snes->dm, &j1a)); 779566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(snes->dm, &j2a)); 789566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(snes->dm, &x2)); 79b1f624c7SBarry Smith } else { 80aa79bc6dSLois Curfman McInnes snes->nvwork = 3; 819566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(x1, snes->nvwork, &snes->vwork)); 829371c9d4SSatish Balay j1a = snes->vwork[0]; 839371c9d4SSatish Balay j2a = snes->vwork[1]; 849371c9d4SSatish Balay x2 = snes->vwork[2]; 85b1f624c7SBarry Smith } 86b1f624c7SBarry Smith } 8723242f5aSBarry Smith 889566063dSJacob Faibussowitsch PetscCall(VecGetSize(x1, &N)); 899566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(x1, &start, &end)); 909566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 919566063dSJacob Faibussowitsch PetscCall(DMGetDMSNES(dm, &dms)); 920df40c35SBarry Smith if (dms->ops->computemffunction) { 939566063dSJacob Faibussowitsch PetscCall(SNESComputeMFFunction(snes, x1, j1a)); 940df40c35SBarry Smith } else { 9576c63389SBarry Smith PetscCall(SNESComputeFunction(snes, x1, j1a)); /* does not handle use of SNESSetFunctionDomainError() correctly */ 960df40c35SBarry Smith } 97c005e166SLois Curfman McInnes 98d0609cedSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)snes), ((PetscObject)snes)->prefix, "Differencing options", "SNES"); 999566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-mat_fd_type", "Algorithm to compute difference parameter", "SNESComputeJacobianDefault", list, 2, "wp", &value, &flg)); 100d0609cedSBarry Smith PetscOptionsEnd(); 101f5af7f23SKarl Rupp if (flg && !value) use_wp = PETSC_FALSE; 102f5af7f23SKarl Rupp 10348a46eb9SPierre Jolivet if (use_wp) PetscCall(VecNorm(x1, NORM_2, &unorm)); 104c005e166SLois Curfman McInnes /* Compute Jacobian approximation, 1 column at a time. 10588c956adSLois Curfman McInnes x1 = current iterate, j1a = F(x1) 10688c956adSLois Curfman McInnes x2 = perturbed iterate, j2a = F(x2) 107c005e166SLois Curfman McInnes */ 10839e2f89bSBarry Smith for (i = 0; i < N; i++) { 1099566063dSJacob Faibussowitsch PetscCall(VecCopy(x1, x2)); 11023242f5aSBarry Smith if (i >= start && i < end) { 1119566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x1, &xx)); 1125620d6dcSBarry Smith if (use_wp) dx = PetscSqrtReal(1.0 + unorm); 113f5af7f23SKarl Rupp else dx = xx[i - start]; 1149566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x1, &xx)); 1156bd79f97SJed Brown if (PetscAbsScalar(dx) < dx_min) dx = (PetscRealPart(dx) < 0. ? -1. : 1.) * dx_par; 11639e2f89bSBarry Smith dx *= epsilon; 11774f6f00dSLois Curfman McInnes wscale = 1.0 / dx; 118f2bd0052SStefano Zampini if (x2->ops->setvalues) PetscCall(VecSetValues(x2, 1, &i, &dx, ADD_VALUES)); 119f2bd0052SStefano Zampini else { 120f2bd0052SStefano Zampini PetscCall(VecGetArray(x2, &y)); 121f2bd0052SStefano Zampini y[i - start] += dx; 122f2bd0052SStefano Zampini PetscCall(VecRestoreArray(x2, &y)); 123f2bd0052SStefano Zampini } 1246c783aadSBarry Smith } else { 125bbb6d6a8SBarry Smith wscale = 0.0; 126bbb6d6a8SBarry Smith } 1279566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(x2)); 1289566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(x2)); 1290df40c35SBarry Smith if (dms->ops->computemffunction) { 1309566063dSJacob Faibussowitsch PetscCall(SNESComputeMFFunction(snes, x2, j2a)); 1310df40c35SBarry Smith } else { 13276c63389SBarry Smith PetscCall(SNESComputeFunction(snes, x2, j2a)); /* does not handle use of SNESSetFunctionDomainError() correctly */ 1330df40c35SBarry Smith } 1349566063dSJacob Faibussowitsch PetscCall(VecAXPY(j2a, -1.0, j1a)); 13586c88abdSHong Zhang /* Communicate scale=1/dx_i to all processors */ 1369566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRanges(x1, &ranges)); 13786c88abdSHong Zhang root = size; 13886c88abdSHong Zhang for (j = size - 1; j > -1; j--) { 13986c88abdSHong Zhang root--; 14086c88abdSHong Zhang if (i >= ranges[j]) break; 14186c88abdSHong Zhang } 1429566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&wscale, 1, MPIU_SCALAR, root, comm)); 1439566063dSJacob Faibussowitsch PetscCall(VecScale(j2a, wscale)); 1449371c9d4SSatish Balay PetscCall(VecNorm(j2a, NORM_INFINITY, &amax)); 1459371c9d4SSatish Balay amax *= 1.e-14; 1469566063dSJacob Faibussowitsch PetscCall(VecGetArray(j2a, &y)); 14723242f5aSBarry Smith for (j = start; j < end; j++) { 14848a46eb9SPierre Jolivet if (PetscAbsScalar(y[j - start]) > amax || j == i) PetscCall(MatSetValues(B, 1, &j, 1, &i, y + j - start, INSERT_VALUES)); 14923242f5aSBarry Smith } 1509566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(j2a, &y)); 15123242f5aSBarry Smith } 152b1f624c7SBarry Smith if (snes->dm) { 1539566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(snes->dm, &j1a)); 1549566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(snes->dm, &j2a)); 1559566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(snes->dm, &x2)); 156b1f624c7SBarry Smith } 1579566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 1589566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 15994ab13aaSBarry Smith if (B != J) { 1609566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 1619566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 162f588057bSBarry Smith } 1630df40c35SBarry Smith snes->max_funcs = max_funcs; 1640df40c35SBarry Smith snes->nfuncs -= N; 1653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16611320018SBarry Smith } 167