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