xref: /petsc/src/mat/tests/ex256.c (revision 061e922f3926be00487707c73b78dd3d40309129)
106c5243aSJose E. Roman static char help[] = "Test some operations of SeqDense matrices with an LDA larger than M.\n";
206c5243aSJose E. Roman 
306c5243aSJose E. Roman #include <petscmat.h>
406c5243aSJose E. Roman 
main(int argc,char ** argv)5*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
6*d71ae5a4SJacob Faibussowitsch {
706c5243aSJose E. Roman   Mat          A, B;
806c5243aSJose E. Roman   PetscScalar *a, *b;
906c5243aSJose E. Roman   PetscInt     n = 4, lda = 5, i;
1006c5243aSJose E. Roman 
11327415f7SBarry Smith   PetscFunctionBeginUser;
1206c5243aSJose E. Roman   PetscCall(PetscInitialize(&argc, &argv, 0, help));
1306c5243aSJose E. Roman   PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL));
1406c5243aSJose E. Roman   PetscCall(PetscOptionsGetInt(NULL, NULL, "-lda", &lda, NULL));
1506c5243aSJose E. Roman   PetscCheck(lda >= n, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "lda %" PetscInt_FMT " < n %" PetscInt_FMT, lda, n);
1606c5243aSJose E. Roman 
1706c5243aSJose E. Roman   /*
1806c5243aSJose E. Roman    * Create two identical matrices (MatDuplicate does not preserve lda)
1906c5243aSJose E. Roman    */
2006c5243aSJose E. Roman   PetscCall(PetscCalloc2(lda * n, &a, lda * n, &b));
2106c5243aSJose E. Roman   for (i = 0; i < n; i++) {
2206c5243aSJose E. Roman     a[i + i * lda] = 1.0 + 2.0 * PETSC_i;
2306c5243aSJose E. Roman     if (i > 0) a[i + (i - 1) * lda] = 3.0 - 0.5 * PETSC_i;
2406c5243aSJose E. Roman     b[i + i * lda] = 1.0 + 2.0 * PETSC_i;
2506c5243aSJose E. Roman     if (i > 0) b[i + (i - 1) * lda] = 3.0 - 0.5 * PETSC_i;
2606c5243aSJose E. Roman   }
2706c5243aSJose E. Roman   PetscCall(MatCreate(PETSC_COMM_SELF, &A));
2806c5243aSJose E. Roman   PetscCall(MatSetSizes(A, n, n, n, n));
2906c5243aSJose E. Roman   PetscCall(MatSetType(A, MATSEQDENSE));
3006c5243aSJose E. Roman   PetscCall(MatSeqDenseSetPreallocation(A, a));
3106c5243aSJose E. Roman   PetscCall(MatDenseSetLDA(A, lda));
3206c5243aSJose E. Roman 
3306c5243aSJose E. Roman   PetscCall(MatCreate(PETSC_COMM_SELF, &B));
3406c5243aSJose E. Roman   PetscCall(MatSetSizes(B, n, n, n, n));
3506c5243aSJose E. Roman   PetscCall(MatSetType(B, MATSEQDENSE));
3606c5243aSJose E. Roman   PetscCall(MatSeqDenseSetPreallocation(B, b));
3706c5243aSJose E. Roman   PetscCall(MatDenseSetLDA(B, lda));
3806c5243aSJose E. Roman 
3906c5243aSJose E. Roman   PetscCall(MatView(A, NULL));
4006c5243aSJose E. Roman   PetscCall(MatConjugate(A));
4106c5243aSJose E. Roman   PetscCall(MatView(A, NULL));
4206c5243aSJose E. Roman   PetscCall(MatRealPart(A));
4306c5243aSJose E. Roman   PetscCall(MatView(A, NULL));
4406c5243aSJose E. Roman   PetscCall(MatImaginaryPart(B));
4506c5243aSJose E. Roman   PetscCall(MatView(B, NULL));
4606c5243aSJose E. Roman 
4706c5243aSJose E. Roman   PetscCall(PetscFree2(a, b));
4806c5243aSJose E. Roman   PetscCall(MatDestroy(&A));
4906c5243aSJose E. Roman   PetscCall(MatDestroy(&B));
5006c5243aSJose E. Roman   PetscCall(PetscFinalize());
5106c5243aSJose E. Roman   return 0;
5206c5243aSJose E. Roman }
5306c5243aSJose E. Roman 
5406c5243aSJose E. Roman /*TEST
5506c5243aSJose E. Roman 
5606c5243aSJose E. Roman    build:
5706c5243aSJose E. Roman      requires: complex
5806c5243aSJose E. Roman 
5906c5243aSJose E. Roman    test:
6006c5243aSJose E. Roman 
6106c5243aSJose E. Roman TEST*/
62