1c4762a1bSJed Brown static char help[] = "Example of inverting a block diagonal matrix.\n" 2c4762a1bSJed Brown "\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown int main(int argc, char **args) 7c4762a1bSJed Brown { 8c4762a1bSJed Brown Mat A,A_inv; 9c4762a1bSJed Brown PetscMPIInt rank,size; 10c4762a1bSJed Brown PetscInt M,m,bs,rstart,rend,j,x,y; 11c4762a1bSJed Brown PetscInt* dnnz; 12c4762a1bSJed Brown PetscScalar *v; 13c4762a1bSJed Brown Vec X, Y; 14c4762a1bSJed Brown PetscReal norm; 15c4762a1bSJed Brown 169566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 179566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 189566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 19c4762a1bSJed Brown 20*d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"ex184","Mat"); 21c4762a1bSJed Brown M=8; 229566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-mat_size",&M,NULL)); 23c4762a1bSJed Brown bs=3; 249566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-mat_block_size",&bs,NULL)); 25*d0609cedSBarry Smith PetscOptionsEnd(); 26c4762a1bSJed Brown 279566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 289566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 299566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M*bs,M*bs)); 309566063dSJacob Faibussowitsch PetscCall(MatSetBlockSize(A,bs)); 319566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 329566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A,&m,NULL)); 339566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m/bs,&dnnz)); 34c4762a1bSJed Brown for (j = 0; j < m/bs; j++) { 35c4762a1bSJed Brown dnnz[j] = 1; 36c4762a1bSJed Brown } 379566063dSJacob Faibussowitsch PetscCall(MatXAIJSetPreallocation(A,bs,dnnz,NULL,NULL,NULL)); 389566063dSJacob Faibussowitsch PetscCall(PetscFree(dnnz)); 39c4762a1bSJed Brown 409566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs*bs,&v)); 419566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&rstart,&rend)); 42c4762a1bSJed Brown for (j = rstart/bs; j < rend/bs; j++) { 43c4762a1bSJed Brown for (x = 0; x < bs; x++) { 44c4762a1bSJed Brown for (y = 0; y < bs; y++) { 45c4762a1bSJed Brown if (x == y) { 46c4762a1bSJed Brown v[y+bs*x] = 2*bs; 47c4762a1bSJed Brown } else { 48c4762a1bSJed Brown v[y+bs*x] = -1 * (x < y) - 2 * (x > y); 49c4762a1bSJed Brown } 50c4762a1bSJed Brown } 51c4762a1bSJed Brown } 529566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlocked(A,1,&j,1,&j,v,INSERT_VALUES)); 53c4762a1bSJed Brown } 549566063dSJacob Faibussowitsch PetscCall(PetscFree(v)); 559566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 569566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 57c4762a1bSJed Brown 58c4762a1bSJed Brown /* check that A = inv(inv(A)) */ 599566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A_inv)); 609566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A_inv)); 619566063dSJacob Faibussowitsch PetscCall(MatInvertBlockDiagonalMat(A,A_inv)); 62c4762a1bSJed Brown 63c4762a1bSJed Brown /* Test A_inv * A on a random vector */ 649566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &X, &Y)); 659566063dSJacob Faibussowitsch PetscCall(VecSetRandom(X, NULL)); 669566063dSJacob Faibussowitsch PetscCall(MatMult(A, X, Y)); 679566063dSJacob Faibussowitsch PetscCall(VecScale(X, -1)); 689566063dSJacob Faibussowitsch PetscCall(MatMultAdd(A_inv, Y, X, X)); 699566063dSJacob Faibussowitsch PetscCall(VecNorm(X, NORM_MAX, &norm)); 70c4762a1bSJed Brown if (norm > PETSC_SMALL) { 719566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Norm of error exceeds tolerance.\nInverse of block diagonal A\n")); 729566063dSJacob Faibussowitsch PetscCall(MatView(A_inv,PETSC_VIEWER_STDOUT_WORLD)); 73c4762a1bSJed Brown } 74c4762a1bSJed Brown 759566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 769566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A_inv)); 779566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X)); 789566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Y)); 79c4762a1bSJed Brown 809566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 81b122ec5aSJacob Faibussowitsch return 0; 82c4762a1bSJed Brown } 83c4762a1bSJed Brown 84c4762a1bSJed Brown /*TEST 85c4762a1bSJed Brown test: 86c4762a1bSJed Brown suffix: seqaij 87c4762a1bSJed Brown args: -mat_type seqaij -mat_size 12 -mat_block_size 3 88c4762a1bSJed Brown nsize: 1 89c4762a1bSJed Brown test: 90c4762a1bSJed Brown suffix: seqbaij 91c4762a1bSJed Brown args: -mat_type seqbaij -mat_size 12 -mat_block_size 3 92c4762a1bSJed Brown nsize: 1 93c4762a1bSJed Brown test: 94c4762a1bSJed Brown suffix: mpiaij 95c4762a1bSJed Brown args: -mat_type mpiaij -mat_size 12 -mat_block_size 3 96c4762a1bSJed Brown nsize: 2 97c4762a1bSJed Brown test: 98c4762a1bSJed Brown suffix: mpibaij 99c4762a1bSJed Brown args: -mat_type mpibaij -mat_size 12 -mat_block_size 3 100c4762a1bSJed Brown nsize: 2 101c4762a1bSJed Brown TEST*/ 102