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