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