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 199371c9d4SSatish Balay PetscErrorCode ComputeFunctionLinear(SNES snes, Vec x, Vec f, void *ctx) { 20c4762a1bSJed Brown Mat A = (Mat)ctx; 21c4762a1bSJed Brown 22c4762a1bSJed Brown PetscFunctionBeginUser; 239566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, f)); 24c4762a1bSJed Brown PetscFunctionReturn(0); 25c4762a1bSJed Brown } 26c4762a1bSJed Brown 279371c9d4SSatish Balay PetscErrorCode ComputeJacobianLinear(SNES snes, Vec x, Mat A, Mat J, void *ctx) { 28c4762a1bSJed Brown PetscFunctionBeginUser; 29c4762a1bSJed Brown PetscFunctionReturn(0); 30c4762a1bSJed Brown } 31c4762a1bSJed Brown 329371c9d4SSatish Balay PetscErrorCode ConstructProblem1(Mat A, Vec b) { 33c4762a1bSJed Brown PetscInt rStart, rEnd, row; 34c4762a1bSJed Brown 35c4762a1bSJed Brown PetscFunctionBeginUser; 369566063dSJacob Faibussowitsch PetscCall(VecSet(b, -3.0)); 379566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rStart, &rEnd)); 38c4762a1bSJed Brown for (row = rStart; row < rEnd; ++row) { 39c4762a1bSJed Brown PetscScalar val = 1.0; 40c4762a1bSJed Brown 419566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, &row, &val, INSERT_VALUES)); 42c4762a1bSJed Brown } 439566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 449566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 45c4762a1bSJed Brown PetscFunctionReturn(0); 46c4762a1bSJed Brown } 47c4762a1bSJed Brown 489371c9d4SSatish Balay PetscErrorCode CheckProblem1(Mat A, Vec b, Vec u) { 49c4762a1bSJed Brown Vec errorVec; 50c4762a1bSJed Brown PetscReal norm, error; 51c4762a1bSJed Brown 52c4762a1bSJed Brown PetscFunctionBeginUser; 539566063dSJacob Faibussowitsch PetscCall(VecDuplicate(b, &errorVec)); 549566063dSJacob Faibussowitsch PetscCall(VecWAXPY(errorVec, -1.0, b, u)); 559566063dSJacob Faibussowitsch PetscCall(VecNorm(errorVec, NORM_2, &error)); 569566063dSJacob Faibussowitsch PetscCall(VecNorm(b, NORM_2, &norm)); 5763a3b9bcSJacob Faibussowitsch PetscCheck(error / norm <= 1000. * PETSC_MACHINE_EPSILON, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Relative error %g is too large", (double)(error / norm)); 589566063dSJacob Faibussowitsch PetscCall(VecDestroy(&errorVec)); 59c4762a1bSJed Brown PetscFunctionReturn(0); 60c4762a1bSJed Brown } 61c4762a1bSJed Brown 629371c9d4SSatish Balay PetscErrorCode ConstructProblem2(Mat A, Vec b) { 63c4762a1bSJed Brown PetscInt N = 10, constraintSize = 4; 64c4762a1bSJed Brown PetscInt row; 65c4762a1bSJed Brown 66c4762a1bSJed Brown PetscFunctionBeginUser; 679566063dSJacob Faibussowitsch PetscCall(VecSet(b, -3.0)); 68c4762a1bSJed Brown for (row = 0; row < constraintSize; ++row) { 69c4762a1bSJed Brown PetscScalar vals[2] = {1.0, 1.0}; 70c4762a1bSJed Brown PetscInt cols[2]; 71c4762a1bSJed Brown 729371c9d4SSatish Balay cols[0] = row; 739371c9d4SSatish Balay cols[1] = row + N - constraintSize; 749566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 2, cols, vals, INSERT_VALUES)); 75c4762a1bSJed Brown } 76c4762a1bSJed Brown for (row = constraintSize; row < N - constraintSize; ++row) { 77c4762a1bSJed Brown PetscScalar val = 1.0; 78c4762a1bSJed Brown 799566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, &row, &val, INSERT_VALUES)); 80c4762a1bSJed Brown } 81c4762a1bSJed Brown for (row = N - constraintSize; row < N; ++row) { 82c4762a1bSJed Brown PetscInt col = row - (N - constraintSize); 83c4762a1bSJed Brown PetscScalar val = 1.0; 84c4762a1bSJed Brown 859566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, &col, &val, INSERT_VALUES)); 86c4762a1bSJed Brown } 879566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 889566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 89c4762a1bSJed Brown PetscFunctionReturn(0); 90c4762a1bSJed Brown } 91c4762a1bSJed Brown 929371c9d4SSatish Balay PetscErrorCode CheckProblem2(Mat A, Vec b, Vec u) { 93c4762a1bSJed Brown PetscInt N = 10, constraintSize = 4, r; 94c4762a1bSJed Brown PetscReal norm, error; 95c4762a1bSJed Brown const PetscScalar *uArray, *bArray; 96c4762a1bSJed Brown 97c4762a1bSJed Brown PetscFunctionBeginUser; 989566063dSJacob Faibussowitsch PetscCall(VecNorm(b, NORM_2, &norm)); 999566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(u, &uArray)); 1009566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(b, &bArray)); 101c4762a1bSJed Brown error = 0.0; 102c4762a1bSJed Brown for (r = 0; r < constraintSize; ++r) error += PetscRealPart(PetscSqr(uArray[r] - bArray[r + N - constraintSize])); 103c4762a1bSJed Brown 10463a3b9bcSJacob Faibussowitsch PetscCheck(error / norm <= 10000 * PETSC_MACHINE_EPSILON, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Relative error %g is too large", (double)(error / norm)); 105c4762a1bSJed Brown error = 0.0; 106c4762a1bSJed Brown for (r = constraintSize; r < N - constraintSize; ++r) error += PetscRealPart(PetscSqr(uArray[r] - bArray[r])); 107c4762a1bSJed Brown 10863a3b9bcSJacob Faibussowitsch PetscCheck(error / norm <= 10000 * PETSC_MACHINE_EPSILON, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Relative error %g is too large", (double)(error / norm)); 109c4762a1bSJed Brown error = 0.0; 110c4762a1bSJed Brown for (r = N - constraintSize; r < N; ++r) error += PetscRealPart(PetscSqr(uArray[r] - (bArray[r - (N - constraintSize)] - bArray[r]))); 111c4762a1bSJed Brown 11263a3b9bcSJacob Faibussowitsch PetscCheck(error / norm <= 10000 * PETSC_MACHINE_EPSILON, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Relative error %g is too large", (double)(error / norm)); 1139566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(u, &uArray)); 1149566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(b, &bArray)); 115c4762a1bSJed Brown PetscFunctionReturn(0); 116c4762a1bSJed Brown } 117c4762a1bSJed Brown 1189371c9d4SSatish Balay int main(int argc, char **argv) { 119c4762a1bSJed Brown MPI_Comm comm; 120c4762a1bSJed Brown SNES snes; /* nonlinear solver */ 121c4762a1bSJed Brown Vec u, r, b; /* solution, residual, and rhs vectors */ 122c4762a1bSJed Brown Mat A, J; /* Jacobian matrix */ 123c4762a1bSJed Brown PetscInt problem = 1, N = 10; 124c4762a1bSJed Brown 125327415f7SBarry Smith PetscFunctionBeginUser; 1269566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 127c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 1289566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-problem", &problem, NULL)); 1299566063dSJacob Faibussowitsch PetscCall(VecCreate(comm, &u)); 1309566063dSJacob Faibussowitsch PetscCall(VecSetSizes(u, PETSC_DETERMINE, N)); 1319566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(u)); 1329566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &r)); 1339566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &b)); 134c4762a1bSJed Brown 1359566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &A)); 1369566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DETERMINE, PETSC_DETERMINE, N, N)); 1379566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 1389566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(A, 5, NULL)); 139c4762a1bSJed Brown J = A; 140c4762a1bSJed Brown 141c4762a1bSJed Brown switch (problem) { 1429371c9d4SSatish Balay case 1: PetscCall(ConstructProblem1(A, b)); break; 1439371c9d4SSatish Balay case 2: PetscCall(ConstructProblem2(A, b)); break; 1449371c9d4SSatish Balay default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid problem number %" PetscInt_FMT, problem); 145c4762a1bSJed Brown } 146c4762a1bSJed Brown 1479566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 1489566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, A, J, ComputeJacobianLinear, NULL)); 1499566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, r, ComputeFunctionLinear, A)); 1509566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 151c4762a1bSJed Brown 1529566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, b, u)); 1539566063dSJacob Faibussowitsch PetscCall(VecView(u, NULL)); 154c4762a1bSJed Brown 155c4762a1bSJed Brown switch (problem) { 1569371c9d4SSatish Balay case 1: PetscCall(CheckProblem1(A, b, u)); break; 1579371c9d4SSatish Balay case 2: PetscCall(CheckProblem2(A, b, u)); break; 1589371c9d4SSatish Balay default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid problem number %" PetscInt_FMT, problem); 159c4762a1bSJed Brown } 160c4762a1bSJed Brown 161*48a46eb9SPierre Jolivet if (A != J) PetscCall(MatDestroy(&A)); 1629566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 1639566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 1649566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 1659566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 1669566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 1679566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 168b122ec5aSJacob Faibussowitsch return 0; 169c4762a1bSJed Brown } 170c4762a1bSJed Brown 171c4762a1bSJed Brown /*TEST 172c4762a1bSJed Brown 173c4762a1bSJed Brown test: 174c4762a1bSJed Brown args: -snes_monitor 175c4762a1bSJed Brown 176c4762a1bSJed Brown test: 177c4762a1bSJed Brown suffix: 2 178c4762a1bSJed Brown args: -problem 2 -pc_type jacobi -snes_monitor 179c4762a1bSJed Brown 180c4762a1bSJed Brown TEST*/ 181