xref: /petsc/src/mat/tests/ex210.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
1c4762a1bSJed Brown static char help[] = "Test MatCreateNest with block sizes.\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown int main(int argc, char **argv)
6c4762a1bSJed Brown {
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   PetscErrorCode         ierr;
12c4762a1bSJed Brown 
13c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL, help); if (ierr) return ierr;
14*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
152c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size > 1,PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"Only coded for one process");
16*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD, &A));
17*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 1, 2));
18*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetBlockSizes(A, 1, 2));
19*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A,MATAIJ));
20*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD, 2, 1, indices,PETSC_COPY_VALUES, &cmap));
21*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD, 1, 1, indices,PETSC_COPY_VALUES, &rmap));
22*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetLocalToGlobalMapping(A, rmap, cmap));
23*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD, &B));
24*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(B, PETSC_DECIDE, PETSC_DECIDE, 1, 1));
25*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetBlockSizes(A, 1, 1));
26*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(B,MATAIJ));
27*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetLocalToGlobalMapping(B, rmap, rmap));
28*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISLocalToGlobalMappingDestroy(&rmap));
29*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISLocalToGlobalMappingDestroy(&cmap));
30*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(A));
31*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(B));
32c4762a1bSJed Brown   mats[0] = A;
33c4762a1bSJed Brown   mats[1] = B;
34*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateNest(PETSC_COMM_WORLD, 1, NULL, 2, NULL, mats,&C));
35*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(C));
36*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
37*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
38*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(C, NULL));
39*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
40*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
41*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
42c4762a1bSJed Brown   ierr = PetscFinalize();
43c4762a1bSJed Brown   return ierr;
44c4762a1bSJed Brown }
45c4762a1bSJed Brown 
46c4762a1bSJed Brown /*TEST
47c4762a1bSJed Brown 
48c4762a1bSJed Brown     test:
49c4762a1bSJed Brown 
50c4762a1bSJed Brown TEST*/
51