1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests repeated use of assembly for matrices.\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown int main(int argc,char **args) 7c4762a1bSJed Brown { 8c4762a1bSJed Brown Mat C; 9c4762a1bSJed Brown PetscInt i,j,m = 5,n = 2,Ii,J; 10c4762a1bSJed Brown PetscMPIInt rank,size; 11c4762a1bSJed Brown PetscScalar v; 12c4762a1bSJed Brown 13*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 14*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 15*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 16c4762a1bSJed Brown n = 2*size; 17c4762a1bSJed Brown 18c4762a1bSJed Brown /* create the matrix for the five point stencil, YET AGAIN*/ 19*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&C)); 20*9566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n)); 21*9566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 22*9566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 23c4762a1bSJed Brown for (i=0; i<m; i++) { 24c4762a1bSJed Brown for (j=2*rank; j<2*rank+2; j++) { 25c4762a1bSJed Brown v = -1.0; Ii = j + n*i; 26*9566063dSJacob Faibussowitsch if (i>0) {J = Ii - n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 27*9566063dSJacob Faibussowitsch if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 28*9566063dSJacob Faibussowitsch if (j>0) {J = Ii - 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 29*9566063dSJacob Faibussowitsch if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 30*9566063dSJacob Faibussowitsch v = 4.0; PetscCall(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 31c4762a1bSJed Brown } 32c4762a1bSJed Brown } 33*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 34*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 35c4762a1bSJed Brown for (i=0; i<m; i++) { 36c4762a1bSJed Brown for (j=2*rank; j<2*rank+2; j++) { 37c4762a1bSJed Brown v = 1.0; Ii = j + n*i; 38*9566063dSJacob Faibussowitsch if (i>0) {J = Ii - n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 39*9566063dSJacob Faibussowitsch if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 40*9566063dSJacob Faibussowitsch if (j>0) {J = Ii - 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 41*9566063dSJacob Faibussowitsch if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES));} 42*9566063dSJacob Faibussowitsch v = -4.0; PetscCall(MatSetValues(C,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 43c4762a1bSJed Brown } 44c4762a1bSJed Brown } 45*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 46*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 47c4762a1bSJed Brown 48*9566063dSJacob Faibussowitsch PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 49c4762a1bSJed Brown 50*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 51*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 52b122ec5aSJacob Faibussowitsch return 0; 53c4762a1bSJed Brown } 54c4762a1bSJed Brown 55c4762a1bSJed Brown /*TEST 56c4762a1bSJed Brown 57c4762a1bSJed Brown test: 58c4762a1bSJed Brown 59c4762a1bSJed Brown TEST*/ 60