16aad120cSJose E. Roman static char help[] = "testing SeqDense matrices with an LDA (leading dimension of the user-allocated array) larger than M.\n"; 2c4762a1bSJed Brown /* 3c4762a1bSJed Brown * Example code testing SeqDense matrices with an LDA (leading dimension 46aad120cSJose E. Roman * of the user-allocated array) larger than M. 5c4762a1bSJed Brown * This example tests the functionality of MatSetValues(), MatMult(), 6c4762a1bSJed Brown * and MatMultTranspose(). 7c4762a1bSJed Brown */ 8c4762a1bSJed Brown 9c4762a1bSJed Brown #include <petscmat.h> 10c4762a1bSJed Brown 11*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 12*d71ae5a4SJacob Faibussowitsch { 13c4762a1bSJed Brown Mat A, A11, A12, A21, A22; 14c4762a1bSJed Brown Vec X, X1, X2, Y, Z, Z1, Z2; 15c4762a1bSJed Brown PetscScalar *a, *b, *x, *y, *z, v, one = 1; 16c4762a1bSJed Brown PetscReal nrm; 17c4762a1bSJed Brown PetscInt size = 8, size1 = 6, size2 = 2, i, j; 18c4762a1bSJed Brown PetscRandom rnd; 19c4762a1bSJed Brown 20327415f7SBarry Smith PetscFunctionBeginUser; 219566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, 0, help)); 229566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rnd)); 23c4762a1bSJed Brown 24c4762a1bSJed Brown /* 25c4762a1bSJed Brown * Create matrix and three vectors: these are all normal 26c4762a1bSJed Brown */ 279566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size * size, &a)); 289566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size * size, &b)); 29c4762a1bSJed Brown for (i = 0; i < size; i++) { 30c4762a1bSJed Brown for (j = 0; j < size; j++) { 319566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rnd, &a[i + j * size])); 32c4762a1bSJed Brown b[i + j * size] = a[i + j * size]; 33c4762a1bSJed Brown } 34c4762a1bSJed Brown } 359566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &A)); 369566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, size, size, size, size)); 379566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATSEQDENSE)); 389566063dSJacob Faibussowitsch PetscCall(MatSeqDenseSetPreallocation(A, a)); 399566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 409566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 41c4762a1bSJed Brown 429566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &x)); 43ad540459SPierre Jolivet for (i = 0; i < size; i++) x[i] = one; 449566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, size, x, &X)); 459566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(X)); 469566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(X)); 47c4762a1bSJed Brown 489566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &y)); 499566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, size, y, &Y)); 509566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(Y)); 519566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(Y)); 52c4762a1bSJed Brown 539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &z)); 549566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, size, z, &Z)); 559566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(Z)); 569566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(Z)); 57c4762a1bSJed Brown 58c4762a1bSJed Brown /* 59c4762a1bSJed Brown * Now create submatrices and subvectors 60c4762a1bSJed Brown */ 619566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &A11)); 629566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A11, size1, size1, size1, size1)); 639566063dSJacob Faibussowitsch PetscCall(MatSetType(A11, MATSEQDENSE)); 649566063dSJacob Faibussowitsch PetscCall(MatSeqDenseSetPreallocation(A11, b)); 659566063dSJacob Faibussowitsch PetscCall(MatDenseSetLDA(A11, size)); 669566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A11, MAT_FINAL_ASSEMBLY)); 679566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A11, MAT_FINAL_ASSEMBLY)); 68c4762a1bSJed Brown 699566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &A12)); 709566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A12, size1, size2, size1, size2)); 719566063dSJacob Faibussowitsch PetscCall(MatSetType(A12, MATSEQDENSE)); 729566063dSJacob Faibussowitsch PetscCall(MatSeqDenseSetPreallocation(A12, b + size1 * size)); 739566063dSJacob Faibussowitsch PetscCall(MatDenseSetLDA(A12, size)); 749566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A12, MAT_FINAL_ASSEMBLY)); 759566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A12, MAT_FINAL_ASSEMBLY)); 76c4762a1bSJed Brown 779566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &A21)); 789566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A21, size2, size1, size2, size1)); 799566063dSJacob Faibussowitsch PetscCall(MatSetType(A21, MATSEQDENSE)); 809566063dSJacob Faibussowitsch PetscCall(MatSeqDenseSetPreallocation(A21, b + size1)); 819566063dSJacob Faibussowitsch PetscCall(MatDenseSetLDA(A21, size)); 829566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A21, MAT_FINAL_ASSEMBLY)); 839566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A21, MAT_FINAL_ASSEMBLY)); 84c4762a1bSJed Brown 859566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF, &A22)); 869566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A22, size2, size2, size2, size2)); 879566063dSJacob Faibussowitsch PetscCall(MatSetType(A22, MATSEQDENSE)); 889566063dSJacob Faibussowitsch PetscCall(MatSeqDenseSetPreallocation(A22, b + size1 * size + size1)); 899566063dSJacob Faibussowitsch PetscCall(MatDenseSetLDA(A22, size)); 909566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A22, MAT_FINAL_ASSEMBLY)); 919566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A22, MAT_FINAL_ASSEMBLY)); 92c4762a1bSJed Brown 939566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, size1, x, &X1)); 949566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, size2, x + size1, &X2)); 959566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, size1, z, &Z1)); 969566063dSJacob Faibussowitsch PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, 1, size2, z + size1, &Z2)); 97c4762a1bSJed Brown 98c4762a1bSJed Brown /* 99c4762a1bSJed Brown * Now multiple matrix times input in two ways; 100c4762a1bSJed Brown * compare the results 101c4762a1bSJed Brown */ 1029566063dSJacob Faibussowitsch PetscCall(MatMult(A, X, Y)); 1039566063dSJacob Faibussowitsch PetscCall(MatMult(A11, X1, Z1)); 1049566063dSJacob Faibussowitsch PetscCall(MatMultAdd(A12, X2, Z1, Z1)); 1059566063dSJacob Faibussowitsch PetscCall(MatMult(A22, X2, Z2)); 1069566063dSJacob Faibussowitsch PetscCall(MatMultAdd(A21, X1, Z2, Z2)); 1079566063dSJacob Faibussowitsch PetscCall(VecAXPY(Z, -1.0, Y)); 1089566063dSJacob Faibussowitsch PetscCall(VecNorm(Z, NORM_2, &nrm)); 10948a46eb9SPierre Jolivet if (nrm > 100.0 * PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test1; error norm=%g\n", (double)nrm)); 110c4762a1bSJed Brown 111c4762a1bSJed Brown /* 112c4762a1bSJed Brown * Next test: change both matrices 113c4762a1bSJed Brown */ 1149566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rnd, &v)); 115c4762a1bSJed Brown i = 1; 116c4762a1bSJed Brown j = size - 2; 1179566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 1, &j, &v, INSERT_VALUES)); 118c4762a1bSJed Brown j -= size1; 1199566063dSJacob Faibussowitsch PetscCall(MatSetValues(A12, 1, &i, 1, &j, &v, INSERT_VALUES)); 1209566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rnd, &v)); 121c4762a1bSJed Brown i = j = size1 + 1; 1229566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 1, &j, &v, INSERT_VALUES)); 123c4762a1bSJed Brown i = j = 1; 1249566063dSJacob Faibussowitsch PetscCall(MatSetValues(A22, 1, &i, 1, &j, &v, INSERT_VALUES)); 1259566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1269566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 1279566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A12, MAT_FINAL_ASSEMBLY)); 1289566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A12, MAT_FINAL_ASSEMBLY)); 1299566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A22, MAT_FINAL_ASSEMBLY)); 1309566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A22, MAT_FINAL_ASSEMBLY)); 131c4762a1bSJed Brown 1329566063dSJacob Faibussowitsch PetscCall(MatMult(A, X, Y)); 1339566063dSJacob Faibussowitsch PetscCall(MatMult(A11, X1, Z1)); 1349566063dSJacob Faibussowitsch PetscCall(MatMultAdd(A12, X2, Z1, Z1)); 1359566063dSJacob Faibussowitsch PetscCall(MatMult(A22, X2, Z2)); 1369566063dSJacob Faibussowitsch PetscCall(MatMultAdd(A21, X1, Z2, Z2)); 1379566063dSJacob Faibussowitsch PetscCall(VecAXPY(Z, -1.0, Y)); 1389566063dSJacob Faibussowitsch PetscCall(VecNorm(Z, NORM_2, &nrm)); 13948a46eb9SPierre Jolivet if (nrm > 100.0 * PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test2; error norm=%g\n", (double)nrm)); 140c4762a1bSJed Brown 141c4762a1bSJed Brown /* 142c4762a1bSJed Brown * Transpose product 143c4762a1bSJed Brown */ 1449566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A, X, Y)); 1459566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A11, X1, Z1)); 1469566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(A21, X2, Z1, Z1)); 1479566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A22, X2, Z2)); 1489566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(A12, X1, Z2, Z2)); 1499566063dSJacob Faibussowitsch PetscCall(VecAXPY(Z, -1.0, Y)); 1509566063dSJacob Faibussowitsch PetscCall(VecNorm(Z, NORM_2, &nrm)); 15148a46eb9SPierre Jolivet if (nrm > 100.0 * PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Test3; error norm=%g\n", (double)nrm)); 1529566063dSJacob Faibussowitsch PetscCall(PetscFree(a)); 1539566063dSJacob Faibussowitsch PetscCall(PetscFree(b)); 1549566063dSJacob Faibussowitsch PetscCall(PetscFree(x)); 1559566063dSJacob Faibussowitsch PetscCall(PetscFree(y)); 1569566063dSJacob Faibussowitsch PetscCall(PetscFree(z)); 1579566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1589566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A11)); 1599566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A12)); 1609566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A21)); 1619566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A22)); 162c4762a1bSJed Brown 1639566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X)); 1649566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y)); 1659566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Z)); 166c4762a1bSJed Brown 1679566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X1)); 1689566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X2)); 1699566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Z1)); 1709566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Z2)); 1719566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rnd)); 1729566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 173b122ec5aSJacob Faibussowitsch return 0; 174c4762a1bSJed Brown } 175c4762a1bSJed Brown 176c4762a1bSJed Brown /*TEST 177c4762a1bSJed Brown 178c4762a1bSJed Brown test: 179c4762a1bSJed Brown 180c4762a1bSJed Brown TEST*/ 181