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 19*2a8381b2SBarry Smith PetscErrorCode ComputeFunctionLinear(SNES snes, Vec x, Vec f, PetscCtx ctx) 20d71ae5a4SJacob Faibussowitsch { 21c4762a1bSJed Brown Mat A = (Mat)ctx; 22c4762a1bSJed Brown 23c4762a1bSJed Brown PetscFunctionBeginUser; 249566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, f)); 253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 26c4762a1bSJed Brown } 27c4762a1bSJed Brown 28*2a8381b2SBarry Smith PetscErrorCode ComputeJacobianLinear(SNES snes, Vec x, Mat A, Mat J, PetscCtx ctx) 29d71ae5a4SJacob Faibussowitsch { 30c4762a1bSJed Brown PetscFunctionBeginUser; 313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 32c4762a1bSJed Brown } 33c4762a1bSJed Brown 34d71ae5a4SJacob Faibussowitsch PetscErrorCode ConstructProblem1(Mat A, Vec b) 35d71ae5a4SJacob Faibussowitsch { 36c4762a1bSJed Brown PetscInt rStart, rEnd, row; 37c4762a1bSJed Brown 38c4762a1bSJed Brown PetscFunctionBeginUser; 399566063dSJacob Faibussowitsch PetscCall(VecSet(b, -3.0)); 409566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rStart, &rEnd)); 41c4762a1bSJed Brown for (row = rStart; row < rEnd; ++row) { 42c4762a1bSJed Brown PetscScalar val = 1.0; 43c4762a1bSJed Brown 449566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, &row, &val, INSERT_VALUES)); 45c4762a1bSJed Brown } 469566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 479566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 49c4762a1bSJed Brown } 50c4762a1bSJed Brown 51d71ae5a4SJacob Faibussowitsch PetscErrorCode CheckProblem1(Mat A, Vec b, Vec u) 52d71ae5a4SJacob Faibussowitsch { 53c4762a1bSJed Brown Vec errorVec; 54c4762a1bSJed Brown PetscReal norm, error; 55c4762a1bSJed Brown 56c4762a1bSJed Brown PetscFunctionBeginUser; 579566063dSJacob Faibussowitsch PetscCall(VecDuplicate(b, &errorVec)); 589566063dSJacob Faibussowitsch PetscCall(VecWAXPY(errorVec, -1.0, b, u)); 599566063dSJacob Faibussowitsch PetscCall(VecNorm(errorVec, NORM_2, &error)); 609566063dSJacob Faibussowitsch PetscCall(VecNorm(b, NORM_2, &norm)); 6163a3b9bcSJacob Faibussowitsch PetscCheck(error / norm <= 1000. * PETSC_MACHINE_EPSILON, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Relative error %g is too large", (double)(error / norm)); 629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&errorVec)); 633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 64c4762a1bSJed Brown } 65c4762a1bSJed Brown 66d71ae5a4SJacob Faibussowitsch PetscErrorCode ConstructProblem2(Mat A, Vec b) 67d71ae5a4SJacob Faibussowitsch { 68c4762a1bSJed Brown PetscInt N = 10, constraintSize = 4; 69c4762a1bSJed Brown PetscInt row; 70c4762a1bSJed Brown 71c4762a1bSJed Brown PetscFunctionBeginUser; 729566063dSJacob Faibussowitsch PetscCall(VecSet(b, -3.0)); 73c4762a1bSJed Brown for (row = 0; row < constraintSize; ++row) { 74c4762a1bSJed Brown PetscScalar vals[2] = {1.0, 1.0}; 75c4762a1bSJed Brown PetscInt cols[2]; 76c4762a1bSJed Brown 779371c9d4SSatish Balay cols[0] = row; 789371c9d4SSatish Balay cols[1] = row + N - constraintSize; 799566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 2, cols, vals, INSERT_VALUES)); 80c4762a1bSJed Brown } 81c4762a1bSJed Brown for (row = constraintSize; row < N - constraintSize; ++row) { 82c4762a1bSJed Brown PetscScalar val = 1.0; 83c4762a1bSJed Brown 849566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, &row, &val, INSERT_VALUES)); 85c4762a1bSJed Brown } 86c4762a1bSJed Brown for (row = N - constraintSize; row < N; ++row) { 87c4762a1bSJed Brown PetscInt col = row - (N - constraintSize); 88c4762a1bSJed Brown PetscScalar val = 1.0; 89c4762a1bSJed Brown 909566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, &col, &val, INSERT_VALUES)); 91c4762a1bSJed Brown } 929566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 939566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 95c4762a1bSJed Brown } 96c4762a1bSJed Brown 97d71ae5a4SJacob Faibussowitsch PetscErrorCode CheckProblem2(Mat A, Vec b, Vec u) 98d71ae5a4SJacob Faibussowitsch { 99c4762a1bSJed Brown PetscInt N = 10, constraintSize = 4, r; 100c4762a1bSJed Brown PetscReal norm, error; 101c4762a1bSJed Brown const PetscScalar *uArray, *bArray; 102c4762a1bSJed Brown 103c4762a1bSJed Brown PetscFunctionBeginUser; 1049566063dSJacob Faibussowitsch PetscCall(VecNorm(b, NORM_2, &norm)); 1059566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(u, &uArray)); 1069566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(b, &bArray)); 107c4762a1bSJed Brown error = 0.0; 108c4762a1bSJed Brown for (r = 0; r < constraintSize; ++r) error += PetscRealPart(PetscSqr(uArray[r] - bArray[r + N - constraintSize])); 109c4762a1bSJed Brown 11063a3b9bcSJacob Faibussowitsch PetscCheck(error / norm <= 10000 * PETSC_MACHINE_EPSILON, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Relative error %g is too large", (double)(error / norm)); 111c4762a1bSJed Brown error = 0.0; 112c4762a1bSJed Brown for (r = constraintSize; r < N - constraintSize; ++r) error += PetscRealPart(PetscSqr(uArray[r] - bArray[r])); 113c4762a1bSJed Brown 11463a3b9bcSJacob Faibussowitsch PetscCheck(error / norm <= 10000 * PETSC_MACHINE_EPSILON, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Relative error %g is too large", (double)(error / norm)); 115c4762a1bSJed Brown error = 0.0; 116c4762a1bSJed Brown for (r = N - constraintSize; r < N; ++r) error += PetscRealPart(PetscSqr(uArray[r] - (bArray[r - (N - constraintSize)] - bArray[r]))); 117c4762a1bSJed Brown 11863a3b9bcSJacob Faibussowitsch PetscCheck(error / norm <= 10000 * PETSC_MACHINE_EPSILON, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Relative error %g is too large", (double)(error / norm)); 1199566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(u, &uArray)); 1209566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(b, &bArray)); 1213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 122c4762a1bSJed Brown } 123c4762a1bSJed Brown 124d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 125d71ae5a4SJacob Faibussowitsch { 126c4762a1bSJed Brown MPI_Comm comm; 127c4762a1bSJed Brown SNES snes; /* nonlinear solver */ 128c4762a1bSJed Brown Vec u, r, b; /* solution, residual, and rhs vectors */ 129c4762a1bSJed Brown Mat A, J; /* Jacobian matrix */ 130c4762a1bSJed Brown PetscInt problem = 1, N = 10; 131c4762a1bSJed Brown 132327415f7SBarry Smith PetscFunctionBeginUser; 1339566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 134c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 1359566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-problem", &problem, NULL)); 1369566063dSJacob Faibussowitsch PetscCall(VecCreate(comm, &u)); 1379566063dSJacob Faibussowitsch PetscCall(VecSetSizes(u, PETSC_DETERMINE, N)); 1389566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(u)); 1399566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &r)); 1409566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &b)); 141c4762a1bSJed Brown 1429566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, &A)); 1439566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DETERMINE, PETSC_DETERMINE, N, N)); 1449566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 1459566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(A, 5, NULL)); 146c4762a1bSJed Brown J = A; 147c4762a1bSJed Brown 148c4762a1bSJed Brown switch (problem) { 149d71ae5a4SJacob Faibussowitsch case 1: 150d71ae5a4SJacob Faibussowitsch PetscCall(ConstructProblem1(A, b)); 151d71ae5a4SJacob Faibussowitsch break; 152d71ae5a4SJacob Faibussowitsch case 2: 153d71ae5a4SJacob Faibussowitsch PetscCall(ConstructProblem2(A, b)); 154d71ae5a4SJacob Faibussowitsch break; 155d71ae5a4SJacob Faibussowitsch default: 156d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid problem number %" PetscInt_FMT, problem); 157c4762a1bSJed Brown } 158c4762a1bSJed Brown 1599566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 1609566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, A, J, ComputeJacobianLinear, NULL)); 1619566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, r, ComputeFunctionLinear, A)); 1629566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 163c4762a1bSJed Brown 1649566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, b, u)); 1659566063dSJacob Faibussowitsch PetscCall(VecView(u, NULL)); 166c4762a1bSJed Brown 167c4762a1bSJed Brown switch (problem) { 168d71ae5a4SJacob Faibussowitsch case 1: 169d71ae5a4SJacob Faibussowitsch PetscCall(CheckProblem1(A, b, u)); 170d71ae5a4SJacob Faibussowitsch break; 171d71ae5a4SJacob Faibussowitsch case 2: 172d71ae5a4SJacob Faibussowitsch PetscCall(CheckProblem2(A, b, u)); 173d71ae5a4SJacob Faibussowitsch break; 174d71ae5a4SJacob Faibussowitsch default: 175d71ae5a4SJacob Faibussowitsch SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid problem number %" PetscInt_FMT, problem); 176c4762a1bSJed Brown } 177c4762a1bSJed Brown 17848a46eb9SPierre Jolivet if (A != J) PetscCall(MatDestroy(&A)); 1799566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 1809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 1819566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 1829566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 1839566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 1849566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 185b122ec5aSJacob Faibussowitsch return 0; 186c4762a1bSJed Brown } 187c4762a1bSJed Brown 188c4762a1bSJed Brown /*TEST 189c4762a1bSJed Brown 190c4762a1bSJed Brown test: 191c4762a1bSJed Brown args: -snes_monitor 192c4762a1bSJed Brown 193c4762a1bSJed Brown test: 194c4762a1bSJed Brown suffix: 2 195c4762a1bSJed Brown args: -problem 2 -pc_type jacobi -snes_monitor 196c4762a1bSJed Brown 197c4762a1bSJed Brown TEST*/ 198