xref: /petsc/src/snes/tests/ex68.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
1c4762a1bSJed Brown static char help[] = "Test problems for Schur complement solvers.\n\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscsnes.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown /*
6c4762a1bSJed Brown Test 1:
7c4762a1bSJed Brown   I u = b
8c4762a1bSJed Brown 
9c4762a1bSJed Brown   solution: u = b
10c4762a1bSJed Brown 
11c4762a1bSJed Brown Test 2:
12c4762a1bSJed Brown   / I 0 I \  / u_1 \   / b_1 \
13c4762a1bSJed Brown   | 0 I 0 | |  u_2 | = | b_2 |
14c4762a1bSJed Brown   \ I 0 0 /  \ u_3 /   \ b_3 /
15c4762a1bSJed Brown 
16c4762a1bSJed Brown   solution: u_1 = b_3, u_2 = b_2, u_3 = b_1 - b_3
17c4762a1bSJed Brown */
18c4762a1bSJed Brown 
19c4762a1bSJed Brown PetscErrorCode ComputeFunctionLinear(SNES snes, Vec x, Vec f, void *ctx)
20c4762a1bSJed Brown {
21c4762a1bSJed Brown   Mat            A = (Mat) ctx;
22c4762a1bSJed Brown   PetscErrorCode ierr;
23c4762a1bSJed Brown 
24c4762a1bSJed Brown   PetscFunctionBeginUser;
25c4762a1bSJed Brown   ierr = MatMult(A, x, f);CHKERRQ(ierr);
26c4762a1bSJed Brown   PetscFunctionReturn(0);
27c4762a1bSJed Brown }
28c4762a1bSJed Brown 
29c4762a1bSJed Brown PetscErrorCode ComputeJacobianLinear(SNES snes, Vec x, Mat A, Mat J, void *ctx)
30c4762a1bSJed Brown {
31c4762a1bSJed Brown   PetscFunctionBeginUser;
32c4762a1bSJed Brown   PetscFunctionReturn(0);
33c4762a1bSJed Brown }
34c4762a1bSJed Brown 
35c4762a1bSJed Brown PetscErrorCode ConstructProblem1(Mat A, Vec b)
36c4762a1bSJed Brown {
37c4762a1bSJed Brown   PetscInt       rStart, rEnd, row;
38c4762a1bSJed Brown   PetscErrorCode ierr;
39c4762a1bSJed Brown 
40c4762a1bSJed Brown   PetscFunctionBeginUser;
41c4762a1bSJed Brown   ierr = VecSet(b, -3.0);CHKERRQ(ierr);
42c4762a1bSJed Brown   ierr = MatGetOwnershipRange(A, &rStart, &rEnd);CHKERRQ(ierr);
43c4762a1bSJed Brown   for (row = rStart; row < rEnd; ++row) {
44c4762a1bSJed Brown     PetscScalar val = 1.0;
45c4762a1bSJed Brown 
46c4762a1bSJed Brown     ierr = MatSetValues(A, 1, &row, 1, &row, &val, INSERT_VALUES);CHKERRQ(ierr);
47c4762a1bSJed Brown   }
48c4762a1bSJed Brown   ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
49c4762a1bSJed Brown   ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
50c4762a1bSJed Brown   PetscFunctionReturn(0);
51c4762a1bSJed Brown }
52c4762a1bSJed Brown 
53c4762a1bSJed Brown PetscErrorCode CheckProblem1(Mat A, Vec b, Vec u)
54c4762a1bSJed Brown {
55c4762a1bSJed Brown   Vec            errorVec;
56c4762a1bSJed Brown   PetscReal      norm, error;
57c4762a1bSJed Brown   PetscErrorCode ierr;
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   PetscFunctionBeginUser;
60c4762a1bSJed Brown   ierr = VecDuplicate(b, &errorVec);CHKERRQ(ierr);
61c4762a1bSJed Brown   ierr = VecWAXPY(errorVec, -1.0, b, u);CHKERRQ(ierr);
62c4762a1bSJed Brown   ierr = VecNorm(errorVec, NORM_2, &error);CHKERRQ(ierr);
63c4762a1bSJed Brown   ierr = VecNorm(b, NORM_2, &norm);CHKERRQ(ierr);
64*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(error/norm > 1000.*PETSC_MACHINE_EPSILON,PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Relative error %g is too large", error/norm);
65c4762a1bSJed Brown   ierr = VecDestroy(&errorVec);CHKERRQ(ierr);
66c4762a1bSJed Brown   PetscFunctionReturn(0);
67c4762a1bSJed Brown }
68c4762a1bSJed Brown 
69c4762a1bSJed Brown PetscErrorCode ConstructProblem2(Mat A, Vec b)
70c4762a1bSJed Brown {
71c4762a1bSJed Brown   PetscInt       N = 10, constraintSize = 4;
72c4762a1bSJed Brown   PetscInt       row;
73c4762a1bSJed Brown   PetscErrorCode ierr;
74c4762a1bSJed Brown 
75c4762a1bSJed Brown   PetscFunctionBeginUser;
76c4762a1bSJed Brown   ierr = VecSet(b, -3.0);CHKERRQ(ierr);
77c4762a1bSJed Brown   for (row = 0; row < constraintSize; ++row) {
78c4762a1bSJed Brown     PetscScalar vals[2] = {1.0, 1.0};
79c4762a1bSJed Brown     PetscInt    cols[2];
80c4762a1bSJed Brown 
81c4762a1bSJed Brown     cols[0] = row; cols[1] = row + N - constraintSize;
82c4762a1bSJed Brown     ierr    = MatSetValues(A, 1, &row, 2, cols, vals, INSERT_VALUES);CHKERRQ(ierr);
83c4762a1bSJed Brown   }
84c4762a1bSJed Brown   for (row = constraintSize; row < N - constraintSize; ++row) {
85c4762a1bSJed Brown     PetscScalar val = 1.0;
86c4762a1bSJed Brown 
87c4762a1bSJed Brown     ierr = MatSetValues(A, 1, &row, 1, &row, &val, INSERT_VALUES);CHKERRQ(ierr);
88c4762a1bSJed Brown   }
89c4762a1bSJed Brown   for (row = N - constraintSize; row < N; ++row) {
90c4762a1bSJed Brown     PetscInt    col = row - (N - constraintSize);
91c4762a1bSJed Brown     PetscScalar val = 1.0;
92c4762a1bSJed Brown 
93c4762a1bSJed Brown     ierr = MatSetValues(A, 1, &row, 1, &col, &val, INSERT_VALUES);CHKERRQ(ierr);
94c4762a1bSJed Brown   }
95c4762a1bSJed Brown   ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
96c4762a1bSJed Brown   ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
97c4762a1bSJed Brown   PetscFunctionReturn(0);
98c4762a1bSJed Brown }
99c4762a1bSJed Brown 
100c4762a1bSJed Brown PetscErrorCode CheckProblem2(Mat A, Vec b, Vec u)
101c4762a1bSJed Brown {
102c4762a1bSJed Brown   PetscInt          N = 10, constraintSize = 4, r;
103c4762a1bSJed Brown   PetscReal         norm, error;
104c4762a1bSJed Brown   const PetscScalar *uArray, *bArray;
105c4762a1bSJed Brown   PetscErrorCode    ierr;
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   PetscFunctionBeginUser;
108c4762a1bSJed Brown   ierr  = VecNorm(b, NORM_2, &norm);CHKERRQ(ierr);
109c4762a1bSJed Brown   ierr  = VecGetArrayRead(u, &uArray);CHKERRQ(ierr);
110c4762a1bSJed Brown   ierr  = VecGetArrayRead(b, &bArray);CHKERRQ(ierr);
111c4762a1bSJed Brown   error = 0.0;
112c4762a1bSJed Brown   for (r = 0; r < constraintSize; ++r) error += PetscRealPart(PetscSqr(uArray[r] - bArray[r + N-constraintSize]));
113c4762a1bSJed Brown 
114*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(error/norm > 10000*PETSC_MACHINE_EPSILON,PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Relative error %g is too large", error/norm);
115c4762a1bSJed Brown   error = 0.0;
116c4762a1bSJed Brown   for (r = constraintSize; r < N - constraintSize; ++r) error += PetscRealPart(PetscSqr(uArray[r] - bArray[r]));
117c4762a1bSJed Brown 
118*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(error/norm > 10000*PETSC_MACHINE_EPSILON,PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Relative error %g is too large", error/norm);
119c4762a1bSJed Brown   error = 0.0;
120c4762a1bSJed Brown   for (r = N - constraintSize; r < N; ++r) error += PetscRealPart(PetscSqr(uArray[r] - (bArray[r - (N-constraintSize)] - bArray[r])));
121c4762a1bSJed Brown 
122*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(error/norm > 10000*PETSC_MACHINE_EPSILON,PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Relative error %g is too large", error/norm);
123c4762a1bSJed Brown   ierr = VecRestoreArrayRead(u, &uArray);CHKERRQ(ierr);
124c4762a1bSJed Brown   ierr = VecRestoreArrayRead(b, &bArray);CHKERRQ(ierr);
125c4762a1bSJed Brown   PetscFunctionReturn(0);
126c4762a1bSJed Brown }
127c4762a1bSJed Brown 
128c4762a1bSJed Brown int main(int argc, char **argv)
129c4762a1bSJed Brown {
130c4762a1bSJed Brown   MPI_Comm       comm;
131c4762a1bSJed Brown   SNES           snes;                 /* nonlinear solver */
132c4762a1bSJed Brown   Vec            u,r,b;                /* solution, residual, and rhs vectors */
133c4762a1bSJed Brown   Mat            A,J;                  /* Jacobian matrix */
134c4762a1bSJed Brown   PetscInt       problem = 1, N = 10;
135c4762a1bSJed Brown   PetscErrorCode ierr;
136c4762a1bSJed Brown 
137c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
138c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
139c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL, "-problem", &problem, NULL);CHKERRQ(ierr);
140c4762a1bSJed Brown   ierr = VecCreate(comm, &u);CHKERRQ(ierr);
141c4762a1bSJed Brown   ierr = VecSetSizes(u, PETSC_DETERMINE, N);CHKERRQ(ierr);
142c4762a1bSJed Brown   ierr = VecSetFromOptions(u);CHKERRQ(ierr);
143c4762a1bSJed Brown   ierr = VecDuplicate(u, &r);CHKERRQ(ierr);
144c4762a1bSJed Brown   ierr = VecDuplicate(u, &b);CHKERRQ(ierr);
145c4762a1bSJed Brown 
146c4762a1bSJed Brown   ierr = MatCreate(comm, &A);CHKERRQ(ierr);
147c4762a1bSJed Brown   ierr = MatSetSizes(A, PETSC_DETERMINE, PETSC_DETERMINE, N, N);CHKERRQ(ierr);
148c4762a1bSJed Brown   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
149c4762a1bSJed Brown   ierr = MatSeqAIJSetPreallocation(A, 5, NULL);CHKERRQ(ierr);
150c4762a1bSJed Brown   J    = A;
151c4762a1bSJed Brown 
152c4762a1bSJed Brown   switch (problem) {
153c4762a1bSJed Brown   case 1:
154c4762a1bSJed Brown     ierr = ConstructProblem1(A, b);CHKERRQ(ierr);
155c4762a1bSJed Brown     break;
156c4762a1bSJed Brown   case 2:
157c4762a1bSJed Brown     ierr = ConstructProblem2(A, b);CHKERRQ(ierr);
158c4762a1bSJed Brown     break;
159c4762a1bSJed Brown   default:
16098921bdaSJacob Faibussowitsch     SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid problem number %d", problem);
161c4762a1bSJed Brown   }
162c4762a1bSJed Brown 
163c4762a1bSJed Brown   ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr);
164c4762a1bSJed Brown   ierr = SNESSetJacobian(snes, A, J, ComputeJacobianLinear, NULL);CHKERRQ(ierr);
165c4762a1bSJed Brown   ierr = SNESSetFunction(snes, r, ComputeFunctionLinear, A);CHKERRQ(ierr);
166c4762a1bSJed Brown   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
167c4762a1bSJed Brown 
168c4762a1bSJed Brown   ierr = SNESSolve(snes, b, u);CHKERRQ(ierr);
169c4762a1bSJed Brown   ierr = VecView(u, NULL);CHKERRQ(ierr);
170c4762a1bSJed Brown 
171c4762a1bSJed Brown   switch (problem) {
172c4762a1bSJed Brown   case 1:
173c4762a1bSJed Brown     ierr = CheckProblem1(A, b, u);CHKERRQ(ierr);
174c4762a1bSJed Brown     break;
175c4762a1bSJed Brown   case 2:
176c4762a1bSJed Brown     ierr = CheckProblem2(A, b, u);CHKERRQ(ierr);
177c4762a1bSJed Brown     break;
178c4762a1bSJed Brown   default:
17998921bdaSJacob Faibussowitsch     SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid problem number %d", problem);
180c4762a1bSJed Brown   }
181c4762a1bSJed Brown 
182c4762a1bSJed Brown   if (A != J) {
183c4762a1bSJed Brown     ierr = MatDestroy(&A);CHKERRQ(ierr);
184c4762a1bSJed Brown   }
185c4762a1bSJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
186c4762a1bSJed Brown   ierr = VecDestroy(&u);CHKERRQ(ierr);
187c4762a1bSJed Brown   ierr = VecDestroy(&r);CHKERRQ(ierr);
188c4762a1bSJed Brown   ierr = VecDestroy(&b);CHKERRQ(ierr);
189c4762a1bSJed Brown   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
190c4762a1bSJed Brown   ierr = PetscFinalize();
191c4762a1bSJed Brown   return ierr;
192c4762a1bSJed Brown }
193c4762a1bSJed Brown 
194c4762a1bSJed Brown /*TEST
195c4762a1bSJed Brown 
196c4762a1bSJed Brown    test:
197c4762a1bSJed Brown      args: -snes_monitor
198c4762a1bSJed Brown 
199c4762a1bSJed Brown    test:
200c4762a1bSJed Brown      suffix: 2
201c4762a1bSJed Brown      args: -problem 2 -pc_type jacobi -snes_monitor
202c4762a1bSJed Brown 
203c4762a1bSJed Brown TEST*/
204