111320018SBarry Smith 2af0996ceSBarry Smith #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/ 3*f2bd0052SStefano Zampini #include <petsc/private/vecimpl.h> /* for Vec->ops->setvalues */ 4b1f624c7SBarry Smith #include <petscdm.h> 511320018SBarry Smith 64b828684SBarry Smith /*@C 78d359177SBarry Smith SNESComputeJacobianDefault - Computes the Jacobian using finite differences. 811320018SBarry Smith 9c3339decSBarry Smith Collective 10fee21e36SBarry Smith 11c7afd0dbSLois Curfman McInnes Input Parameters: 12f6dfbefdSBarry Smith + snes - the `SNES` context 136b867d5aSJose E. Roman . x1 - compute Jacobian at this point 14f6dfbefdSBarry Smith - ctx - application's function context, as set with `SNESSetFunction()` 15c7afd0dbSLois Curfman McInnes 16c7afd0dbSLois Curfman McInnes Output Parameters: 17c7afd0dbSLois Curfman McInnes + J - Jacobian matrix (not altered in this routine) 18dc4c0fb0SBarry Smith - B - newly computed Jacobian matrix to use with preconditioner (generally the same as `J`) 19c7afd0dbSLois Curfman McInnes 20f6dfbefdSBarry Smith Options Database Keys: 21f6dfbefdSBarry Smith + -snes_fd - Activates `SNESComputeJacobianDefault()` 22dc4c0fb0SBarry Smith . -snes_fd_coloring - Activates a faster computation that uses a graph coloring of the matrix 2379f36460SBarry Smith . -snes_test_err - Square root of function error tolerance, default square root of machine 2477d8c4bbSBarry Smith epsilon (1.e-8 in double, 3.e-4 in single) 25f6dfbefdSBarry Smith - -mat_fd_type - Either wp or ds (see `MATMFFD_WP` or `MATMFFD_DS`) 26ad960d00SLois Curfman McInnes 27dc4c0fb0SBarry Smith Level: intermediate 28dc4c0fb0SBarry Smith 295f3c43d9SLois Curfman McInnes Notes: 305f3c43d9SLois Curfman McInnes This routine is slow and expensive, and is not currently optimized 315f3c43d9SLois Curfman McInnes to take advantage of sparsity in the problem. Although 32f6dfbefdSBarry Smith `SNESComputeJacobianDefault()` is not recommended for general use 335f3c43d9SLois Curfman McInnes in large-scale applications, It can be useful in checking the 345f3c43d9SLois Curfman McInnes correctness of a user-provided Jacobian. 3511320018SBarry Smith 3679f36460SBarry Smith An alternative routine that uses coloring to exploit matrix sparsity is 37f6dfbefdSBarry Smith `SNESComputeJacobianDefaultColor()`. 38b4fc646aSLois Curfman McInnes 39f6dfbefdSBarry Smith This routine ignores the maximum number of function evaluations set with `SNESSetTolerances()` and the function 40f6dfbefdSBarry Smith evaluations it performs are not counted in what is returned by of `SNESGetNumberFunctionEvals()`. 41f6dfbefdSBarry Smith 42f6dfbefdSBarry Smith This function can be provided to `SNESSetJacobian()` along with a dense matrix to hold the Jacobian 430df40c35SBarry Smith 44f6dfbefdSBarry Smith .seealso: `SNES`, `SNESSetJacobian()`, `SNESSetJacobian()`, `SNESComputeJacobianDefaultColor()`, `MatCreateSNESMF()` 4511320018SBarry Smith @*/ 46d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESComputeJacobianDefault(SNES snes, Vec x1, Mat J, Mat B, void *ctx) 47d71ae5a4SJacob Faibussowitsch { 4888c956adSLois Curfman McInnes Vec j1a, j2a, x2; 490df40c35SBarry Smith PetscInt i, N, start, end, j, value, root, max_funcs = snes->max_funcs; 505edff71fSBarry Smith PetscScalar dx, *y, wscale; 515edff71fSBarry Smith const PetscScalar *xx; 5277d8c4bbSBarry Smith PetscReal amax, epsilon = PETSC_SQRT_MACHINE_EPSILON; 5379f36460SBarry Smith PetscReal dx_min = 1.e-16, dx_par = 1.e-1, unorm; 54bbb6d6a8SBarry Smith MPI_Comm comm; 55ace3abfcSBarry Smith PetscBool assembled, use_wp = PETSC_TRUE, flg; 5679f36460SBarry Smith const char *list[2] = {"ds", "wp"}; 5786c88abdSHong Zhang PetscMPIInt size; 5886c88abdSHong Zhang const PetscInt *ranges; 590df40c35SBarry Smith DM dm; 600df40c35SBarry Smith DMSNES dms; 610521c3abSLois Curfman McInnes 623a40ed3dSBarry Smith PetscFunctionBegin; 630df40c35SBarry Smith snes->max_funcs = PETSC_MAX_INT; 645d5a84c9SBarry Smith /* Since this Jacobian will possibly have "extra" nonzero locations just turn off errors for these locations */ 659566063dSJacob Faibussowitsch PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 669566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(((PetscObject)snes)->options, ((PetscObject)snes)->prefix, "-snes_test_err", &epsilon, NULL)); 6723242f5aSBarry Smith 689566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)x1, &comm)); 699566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 709566063dSJacob Faibussowitsch PetscCall(MatAssembled(B, &assembled)); 711baa6e33SBarry Smith if (assembled) PetscCall(MatZeroEntries(B)); 72aa79bc6dSLois Curfman McInnes if (!snes->nvwork) { 73b1f624c7SBarry Smith if (snes->dm) { 749566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(snes->dm, &j1a)); 759566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(snes->dm, &j2a)); 769566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(snes->dm, &x2)); 77b1f624c7SBarry Smith } else { 78aa79bc6dSLois Curfman McInnes snes->nvwork = 3; 799566063dSJacob Faibussowitsch PetscCall(VecDuplicateVecs(x1, snes->nvwork, &snes->vwork)); 809371c9d4SSatish Balay j1a = snes->vwork[0]; 819371c9d4SSatish Balay j2a = snes->vwork[1]; 829371c9d4SSatish Balay x2 = snes->vwork[2]; 83b1f624c7SBarry Smith } 84b1f624c7SBarry Smith } 8523242f5aSBarry Smith 869566063dSJacob Faibussowitsch PetscCall(VecGetSize(x1, &N)); 879566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(x1, &start, &end)); 889566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 899566063dSJacob Faibussowitsch PetscCall(DMGetDMSNES(dm, &dms)); 900df40c35SBarry Smith if (dms->ops->computemffunction) { 919566063dSJacob Faibussowitsch PetscCall(SNESComputeMFFunction(snes, x1, j1a)); 920df40c35SBarry Smith } else { 939566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, x1, j1a)); 940df40c35SBarry Smith } 95c005e166SLois Curfman McInnes 96d0609cedSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)snes), ((PetscObject)snes)->prefix, "Differencing options", "SNES"); 979566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-mat_fd_type", "Algorithm to compute difference parameter", "SNESComputeJacobianDefault", list, 2, "wp", &value, &flg)); 98d0609cedSBarry Smith PetscOptionsEnd(); 99f5af7f23SKarl Rupp if (flg && !value) use_wp = PETSC_FALSE; 100f5af7f23SKarl Rupp 10148a46eb9SPierre Jolivet if (use_wp) PetscCall(VecNorm(x1, NORM_2, &unorm)); 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++) { 1079566063dSJacob Faibussowitsch PetscCall(VecCopy(x1, x2)); 10823242f5aSBarry Smith if (i >= start && i < end) { 1099566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x1, &xx)); 1105620d6dcSBarry Smith if (use_wp) dx = PetscSqrtReal(1.0 + unorm); 111f5af7f23SKarl Rupp else dx = xx[i - start]; 1129566063dSJacob 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*f2bd0052SStefano Zampini if (x2->ops->setvalues) PetscCall(VecSetValues(x2, 1, &i, &dx, ADD_VALUES)); 117*f2bd0052SStefano Zampini else { 118*f2bd0052SStefano Zampini PetscCall(VecGetArray(x2, &y)); 119*f2bd0052SStefano Zampini y[i - start] += dx; 120*f2bd0052SStefano Zampini PetscCall(VecRestoreArray(x2, &y)); 121*f2bd0052SStefano Zampini } 1226c783aadSBarry Smith } else { 123bbb6d6a8SBarry Smith wscale = 0.0; 124bbb6d6a8SBarry Smith } 1259566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(x2)); 1269566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(x2)); 1270df40c35SBarry Smith if (dms->ops->computemffunction) { 1289566063dSJacob Faibussowitsch PetscCall(SNESComputeMFFunction(snes, x2, j2a)); 1290df40c35SBarry Smith } else { 1309566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, x2, j2a)); 1310df40c35SBarry Smith } 1329566063dSJacob Faibussowitsch PetscCall(VecAXPY(j2a, -1.0, j1a)); 13386c88abdSHong Zhang /* Communicate scale=1/dx_i to all processors */ 1349566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRanges(x1, &ranges)); 13586c88abdSHong Zhang root = size; 13686c88abdSHong Zhang for (j = size - 1; j > -1; j--) { 13786c88abdSHong Zhang root--; 13886c88abdSHong Zhang if (i >= ranges[j]) break; 13986c88abdSHong Zhang } 1409566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&wscale, 1, MPIU_SCALAR, root, comm)); 1419566063dSJacob Faibussowitsch PetscCall(VecScale(j2a, wscale)); 1429371c9d4SSatish Balay PetscCall(VecNorm(j2a, NORM_INFINITY, &amax)); 1439371c9d4SSatish Balay amax *= 1.e-14; 1449566063dSJacob Faibussowitsch PetscCall(VecGetArray(j2a, &y)); 14523242f5aSBarry Smith for (j = start; j < end; j++) { 14648a46eb9SPierre Jolivet if (PetscAbsScalar(y[j - start]) > amax || j == i) PetscCall(MatSetValues(B, 1, &j, 1, &i, y + j - start, INSERT_VALUES)); 14723242f5aSBarry Smith } 1489566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(j2a, &y)); 14923242f5aSBarry Smith } 150b1f624c7SBarry Smith if (snes->dm) { 1519566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(snes->dm, &j1a)); 1529566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(snes->dm, &j2a)); 1539566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(snes->dm, &x2)); 154b1f624c7SBarry Smith } 1559566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 1569566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 15794ab13aaSBarry Smith if (B != J) { 1589566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 1599566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 160f588057bSBarry Smith } 1610df40c35SBarry Smith snes->max_funcs = max_funcs; 1620df40c35SBarry Smith snes->nfuncs -= N; 1633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16411320018SBarry Smith } 165