1c4762a1bSJed Brown static char help[] = "Test MatCreateNest with block sizes.\n";
2c4762a1bSJed Brown
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown
main(int argc,char ** argv)5*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
6*d71ae5a4SJacob Faibussowitsch {
7c4762a1bSJed Brown Mat A, B, C, mats[2];
8c4762a1bSJed Brown ISLocalToGlobalMapping cmap, rmap;
9c4762a1bSJed Brown const PetscInt indices[1] = {0};
10c4762a1bSJed Brown PetscMPIInt size;
11c4762a1bSJed Brown
12327415f7SBarry Smith PetscFunctionBeginUser;
139566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
15be096a46SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only coded for one process");
169566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
179566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 1, 2));
189566063dSJacob Faibussowitsch PetscCall(MatSetBlockSizes(A, 1, 2));
199566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATAIJ));
209566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD, 2, 1, indices, PETSC_COPY_VALUES, &cmap));
219566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD, 1, 1, indices, PETSC_COPY_VALUES, &rmap));
229566063dSJacob Faibussowitsch PetscCall(MatSetLocalToGlobalMapping(A, rmap, cmap));
239566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &B));
249566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, PETSC_DECIDE, PETSC_DECIDE, 1, 1));
259566063dSJacob Faibussowitsch PetscCall(MatSetBlockSizes(A, 1, 1));
269566063dSJacob Faibussowitsch PetscCall(MatSetType(B, MATAIJ));
279566063dSJacob Faibussowitsch PetscCall(MatSetLocalToGlobalMapping(B, rmap, rmap));
289566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&rmap));
299566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&cmap));
309566063dSJacob Faibussowitsch PetscCall(MatSetUp(A));
319566063dSJacob Faibussowitsch PetscCall(MatSetUp(B));
32c4762a1bSJed Brown mats[0] = A;
33c4762a1bSJed Brown mats[1] = B;
349566063dSJacob Faibussowitsch PetscCall(MatCreateNest(PETSC_COMM_WORLD, 1, NULL, 2, NULL, mats, &C));
359566063dSJacob Faibussowitsch PetscCall(MatSetUp(C));
369566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
379566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
389566063dSJacob Faibussowitsch PetscCall(MatView(C, NULL));
399566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C));
409566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B));
419566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
429566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
43b122ec5aSJacob Faibussowitsch return 0;
44c4762a1bSJed Brown }
45c4762a1bSJed Brown
46c4762a1bSJed Brown /*TEST
47c4762a1bSJed Brown
48c4762a1bSJed Brown test:
49c4762a1bSJed Brown
50c4762a1bSJed Brown TEST*/
51