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 12*327415f7SBarry Smith PetscFunctionBeginUser; 139566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 159566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 16c4762a1bSJed Brown prid = size; 179566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-prid",&prid,NULL)); 18c4762a1bSJed Brown 19c4762a1bSJed Brown m = n = 10*size; 209566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,&A)); 219566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,PETSC_DETERMINE,PETSC_DETERMINE,m,n)); 229566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATSEQAIJ)); 239566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 24c4762a1bSJed Brown 25c4762a1bSJed Brown a[0] = rank+1; 26c4762a1bSJed Brown for (i=0; i<m-rank; i++) { 27c4762a1bSJed Brown col = i+rank; 289566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,1,&col,a,INSERT_VALUES)); 29c4762a1bSJed Brown } 309566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 319566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 32c4762a1bSJed Brown 33c4762a1bSJed Brown if (rank == prid) { 349566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"[%d] A: \n",rank)); 359566063dSJacob Faibussowitsch PetscCall(MatView(A,PETSC_VIEWER_STDOUT_SELF)); 36c4762a1bSJed Brown } 37c4762a1bSJed Brown 38c4762a1bSJed Brown /* Test MatCreateMPIAIJSumSeqAIJ */ 399566063dSJacob Faibussowitsch PetscCall(MatCreateMPIAIJSumSeqAIJ(PETSC_COMM_WORLD,A,PETSC_DECIDE,PETSC_DECIDE,MAT_INITIAL_MATRIX,&B)); 40c4762a1bSJed Brown 41c4762a1bSJed Brown /* Test MAT_REUSE_MATRIX */ 42c4762a1bSJed Brown alpha = 0.1; 43c4762a1bSJed Brown for (i=0; i<3; i++) { 449566063dSJacob Faibussowitsch PetscCall(MatScale(A,alpha)); 459566063dSJacob Faibussowitsch PetscCall(MatCreateMPIAIJSumSeqAIJ(PETSC_COMM_WORLD,A,PETSC_DECIDE,PETSC_DECIDE,MAT_REUSE_MATRIX,&B)); 46c4762a1bSJed Brown } 479566063dSJacob Faibussowitsch PetscCall(MatView(B, PETSC_VIEWER_STDOUT_WORLD)); 489566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 499566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 509566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 51b122ec5aSJacob Faibussowitsch return 0; 52c4762a1bSJed Brown } 53c4762a1bSJed Brown 54c4762a1bSJed Brown /*TEST 55c4762a1bSJed Brown 56c4762a1bSJed Brown test: 57c4762a1bSJed Brown nsize: 3 588cc725e6SPierre Jolivet filter: grep -v " MPI process" 59c4762a1bSJed Brown 60c4762a1bSJed Brown test: 61c4762a1bSJed Brown suffix: 2 628cc725e6SPierre Jolivet filter: grep -v " MPI process" 63c4762a1bSJed Brown 64c4762a1bSJed Brown TEST*/ 65