Actual source code: ex241.cxx
1: static char help[] = "Artificial test to check that snes->domainerror is being reset appropriately";
3: /* ------------------------------------------------------------------------
5: Artificial test to check that snes->domainerror is being reset appropriately
7: ------------------------------------------------------------------------- */
9: #define PETSC_SKIP_COMPLEX
10: #include <petscsnes.h>
12: typedef struct {
13: PetscReal value; /* parameter in nonlinear function */
14: } AppCtx;
16: PetscErrorCode UserFunction(SNES, Vec, Vec, void *);
17: PetscErrorCode UserJacobian(SNES, Vec, Mat, Mat, void *);
19: int main(int argc, char **argv)
20: {
21: SNES snes;
22: Vec x, r;
23: Mat J;
24: PetscInt its;
25: AppCtx user;
26: PetscMPIInt size;
28: PetscFunctionBeginUser;
29: PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
30: PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
31: PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
33: /* Allocate vectors / matrix */
34: PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
35: PetscCall(VecSetSizes(x, PETSC_DECIDE, 1));
36: PetscCall(VecSetFromOptions(x));
37: PetscCall(VecDuplicate(x, &r));
39: PetscCall(MatCreateSeqAIJ(PETSC_COMM_WORLD, 1, 1, 1, NULL, &J));
41: /* Create / set-up SNES */
42: PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
43: PetscCall(SNESSetFunction(snes, r, UserFunction, &user));
44: PetscCall(SNESSetJacobian(snes, J, J, UserJacobian, &user));
45: PetscCall(SNESSetFromOptions(snes));
47: /* Set initial guess (=1) and target value */
48: user.value = 1e-4;
50: PetscCall(VecSet(x, 1.0));
52: /* Set initial guess / solve */
53: PetscCall(SNESSolve(snes, NULL, x));
54: PetscCall(SNESGetIterationNumber(snes, &its));
55: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %" PetscInt_FMT "\n", its));
56: PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
58: /* Done */
59: PetscCall(VecDestroy(&x));
60: PetscCall(VecDestroy(&r));
61: PetscCall(MatDestroy(&J));
62: PetscCall(SNESDestroy(&snes));
63: PetscCall(PetscFinalize());
64: return 0;
65: }
67: /*
68: UserFunction - for nonlinear function x^2 - value = 0
69: */
70: PetscErrorCode UserFunction(SNES snes, Vec X, Vec F, void *ptr)
71: {
72: AppCtx *user = (AppCtx *)ptr;
73: PetscInt N, i;
74: PetscScalar *f;
75: PetscReal half;
76: const PetscScalar *x;
78: half = 0.5;
79: PetscFunctionBeginUser;
80: PetscCall(VecGetSize(X, &N));
81: PetscCall(VecGetArrayRead(X, &x));
82: PetscCall(VecGetArray(F, &f));
84: /* Calculate residual */
85: for (i = 0; i < N; ++i) {
86: /*
87: Test for domain error.
88: Artificial test is applied. With starting value 1.0, first iterate will be 0.5 + user->value/2.
89: Declare (0.5-value,0.5+value) to be infeasible.
90: In later iterations, snes->domainerror should be cleared, allowing iterations in the feasible region to be accepted.
91: */
92: if ((half - user->value) < PetscRealPart(x[i]) && PetscRealPart(x[i]) < (half + user->value)) {
93: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "DOMAIN ERROR: x=%g\n", (double)PetscRealPart(x[i])));
94: PetscCall(SNESSetFunctionDomainError(snes));
95: }
96: f[i] = x[i] * x[i] - user->value;
97: }
98: PetscCall(VecRestoreArrayRead(X, &x));
99: PetscCall(VecRestoreArray(F, &f));
100: PetscFunctionReturn(PETSC_SUCCESS);
101: }
103: /*
104: UserJacobian - for nonlinear function x^2 - value = 0
105: */
106: PetscErrorCode UserJacobian(SNES snes, Vec X, Mat J, Mat jac, void *ptr)
107: {
108: PetscInt N, i, row, col;
109: const PetscScalar *x;
110: PetscScalar v;
112: PetscFunctionBeginUser;
113: PetscCall(VecGetSize(X, &N));
114: PetscCall(VecGetArrayRead(X, &x));
116: /* Calculate Jacobian */
117: for (i = 0; i < N; ++i) {
118: row = i;
119: col = i;
120: v = 2 * x[i];
121: PetscCall(MatSetValues(jac, 1, &row, 1, &col, &v, INSERT_VALUES));
122: }
123: PetscCall(VecRestoreArrayRead(X, &x));
124: PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
125: PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
127: if (jac != J) {
128: PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
129: PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
130: }
131: PetscFunctionReturn(PETSC_SUCCESS);
132: }
134: /*TEST
136: build:
137: requires: !single !defined(PETSC_HAVE_SUN_CXX) !complex
139: test:
140: args: -snes_monitor_solution -snes_linesearch_monitor
142: TEST*/