Actual source code: ex54.c

petsc-3.8.3 2017-12-09
Report Typos and Errors
  1: /*

  3:      Tests MPIDENSE matrix operations MatMultTranspose() with processes with no rows or columns
  4: */

  6:  #include <petscksp.h>

  8: PetscErrorCode fill(Mat m, Vec v)
  9: {
 10:   PetscInt       idxn[3] = {0, 1, 2};
 11:   PetscInt       localRows = 0;
 12:   PetscMPIInt    rank,size;

 16:   MPI_Comm_rank(MPI_COMM_WORLD, &rank);
 17:   MPI_Comm_size(MPI_COMM_WORLD, &size);

 19:   if (rank == 1 || rank == 2) localRows = 4;
 20:   if (size == 1) localRows = 8;
 21:   MatSetSizes(m, localRows, PETSC_DECIDE, PETSC_DECIDE, 3);
 22:   VecSetSizes(v, localRows, PETSC_DECIDE);

 24:   MatSetFromOptions(m);
 25:   VecSetFromOptions(v);
 26:   MatSetUp(m);

 28:   if (size == 1) {
 29:     PetscInt    idxm1[4] = {0, 1, 2, 3};
 30:     PetscScalar values1[12] = {1, 0, 0, 1, 0, 1, 1, 1, 0, 1, 1, 1};
 31:     PetscInt    idxm2[4] = {4, 5, 6, 7};
 32:     PetscScalar values2[12] = {1, 2, 0, 1, 2, 1, 1, 3, 0, 1, 3, 1};

 34:     MatSetValues(m, 4, idxm1, 3, idxn, values1, INSERT_VALUES);
 35:     VecSetValue(v, 0, 1.1, INSERT_VALUES); VecSetValue(v, 1, 2.5, INSERT_VALUES);
 36:     VecSetValue(v, 2, 3, INSERT_VALUES); VecSetValue(v, 3, 4, INSERT_VALUES);

 38:     MatSetValues(m, 4, idxm2, 3, idxn, values2, INSERT_VALUES);
 39:     VecSetValue(v, 4, 5, INSERT_VALUES); VecSetValue(v, 5, 6, INSERT_VALUES);
 40:     VecSetValue(v, 6, 7, INSERT_VALUES); VecSetValue(v, 7, 8, INSERT_VALUES);
 41:   } else if (rank == 1) {
 42:     PetscInt    idxm[4] = {0, 1, 2, 3};
 43:     PetscScalar values[12] = {1, 0, 0, 1, 0, 1, 1, 1, 0, 1, 1, 1};

 45:     MatSetValues(m, 4, idxm, 3, idxn, values, INSERT_VALUES);
 46:     VecSetValue(v, 0, 1.1, INSERT_VALUES); VecSetValue(v, 1, 2.5, INSERT_VALUES);
 47:     VecSetValue(v, 2, 3, INSERT_VALUES); VecSetValue(v, 3, 4, INSERT_VALUES);
 48:   } else if (rank == 2) {
 49:     PetscInt    idxm[4] = {4, 5, 6, 7};
 50:     PetscScalar values[12] = {1, 2, 0, 1, 2, 1, 1, 3, 0, 1, 3, 1};

 52:     MatSetValues(m, 4, idxm, 3, idxn, values, INSERT_VALUES);
 53:     VecSetValue(v, 4, 5, INSERT_VALUES); VecSetValue(v, 5, 6, INSERT_VALUES);
 54:     VecSetValue(v, 6, 7, INSERT_VALUES); VecSetValue(v, 7, 8, INSERT_VALUES);
 55:   }
 56:   MatAssemblyBegin(m, MAT_FINAL_ASSEMBLY);
 57:   MatAssemblyEnd(m, MAT_FINAL_ASSEMBLY);
 58:   VecAssemblyBegin(v);
 59:   VecAssemblyEnd(v);
 60:   return(0);
 61: }


 64: int main(int argc, char** argv)
 65: {
 66:   Mat            Q;
 67:   Vec            v, a;
 68:   KSP            QRsolver;
 69:   PC             pc;

 72:   PetscInitialize(&argc, &argv, NULL, NULL);if (ierr) return ierr;

 74:   VecCreate(PETSC_COMM_WORLD, &v);
 75:   MatCreate(PETSC_COMM_WORLD, &Q);
 76:   MatSetType(Q, MATDENSE);
 77:   fill(Q, v);

 79:   MatCreateVecs(Q, &a, NULL);
 80:   KSPCreate(PETSC_COMM_WORLD, &QRsolver);
 81:   KSPGetPC(QRsolver, &pc);
 82:   PCSetType(pc, PCNONE);
 83:   KSPSetType(QRsolver, KSPLSQR);
 84:   KSPSetFromOptions(QRsolver);
 85:   KSPSetOperators(QRsolver, Q, Q);
 86:   MatViewFromOptions(Q, NULL, "-sys_view");
 87:   VecViewFromOptions(a, NULL, "-rhs_view");
 88:   KSPSolve(QRsolver, v, a);
 89:   KSPDestroy(&QRsolver);
 90:   VecDestroy(&a);
 91:   VecDestroy(&v);
 92:   MatDestroy(&Q);

 94:   PetscFinalize();
 95:   return ierr;
 96: }