1c4762a1bSJed Brown static char help[] = "Testing MatCreateMPIAIJSumSeqAIJ().\n\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; 8c4762a1bSJed Brown MatScalar a[1],alpha; 9c4762a1bSJed Brown PetscMPIInt size,rank; 10c4762a1bSJed Brown PetscInt m,n,i,col, prid; 11c4762a1bSJed Brown 129566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 139566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 15c4762a1bSJed Brown prid = size; 169566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-prid",&prid,NULL)); 17c4762a1bSJed Brown 18c4762a1bSJed Brown m = n = 10*size; 199566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,&A)); 209566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,PETSC_DETERMINE,PETSC_DETERMINE,m,n)); 219566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATSEQAIJ)); 229566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 23c4762a1bSJed Brown 24c4762a1bSJed Brown a[0] = rank+1; 25c4762a1bSJed Brown for (i=0; i<m-rank; i++) { 26c4762a1bSJed Brown col = i+rank; 279566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,1,&col,a,INSERT_VALUES)); 28c4762a1bSJed Brown } 299566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 309566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 31c4762a1bSJed Brown 32c4762a1bSJed Brown if (rank == prid) { 339566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"[%d] A: \n",rank)); 349566063dSJacob Faibussowitsch PetscCall(MatView(A,PETSC_VIEWER_STDOUT_SELF)); 35c4762a1bSJed Brown } 36c4762a1bSJed Brown 37c4762a1bSJed Brown /* Test MatCreateMPIAIJSumSeqAIJ */ 389566063dSJacob Faibussowitsch PetscCall(MatCreateMPIAIJSumSeqAIJ(PETSC_COMM_WORLD,A,PETSC_DECIDE,PETSC_DECIDE,MAT_INITIAL_MATRIX,&B)); 39c4762a1bSJed Brown 40c4762a1bSJed Brown /* Test MAT_REUSE_MATRIX */ 41c4762a1bSJed Brown alpha = 0.1; 42c4762a1bSJed Brown for (i=0; i<3; i++) { 439566063dSJacob Faibussowitsch PetscCall(MatScale(A,alpha)); 449566063dSJacob Faibussowitsch PetscCall(MatCreateMPIAIJSumSeqAIJ(PETSC_COMM_WORLD,A,PETSC_DECIDE,PETSC_DECIDE,MAT_REUSE_MATRIX,&B)); 45c4762a1bSJed Brown } 469566063dSJacob Faibussowitsch PetscCall(MatView(B, PETSC_VIEWER_STDOUT_WORLD)); 479566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 489566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 499566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 50b122ec5aSJacob Faibussowitsch return 0; 51c4762a1bSJed Brown } 52c4762a1bSJed Brown 53c4762a1bSJed Brown /*TEST 54c4762a1bSJed Brown 55c4762a1bSJed Brown test: 56c4762a1bSJed Brown nsize: 3 57*8cc725e6SPierre Jolivet filter: grep -v " MPI process" 58c4762a1bSJed Brown 59c4762a1bSJed Brown test: 60c4762a1bSJed Brown suffix: 2 61*8cc725e6SPierre Jolivet filter: grep -v " MPI process" 62c4762a1bSJed Brown 63c4762a1bSJed Brown TEST*/ 64