1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests reusing MPI parallel matrices and MatGetValues().\n\ 3c4762a1bSJed Brown To test the parallel matrix assembly, this example intentionally lays out\n\ 4c4762a1bSJed Brown the matrix across processors differently from the way it is assembled.\n\ 5c4762a1bSJed Brown This example uses bilinear elements on the unit square. Input arguments are:\n\ 6c4762a1bSJed Brown -m <size> : problem size\n\n"; 7c4762a1bSJed Brown 8c4762a1bSJed Brown #include <petscmat.h> 9c4762a1bSJed Brown 10c4762a1bSJed Brown int FormElementStiffness(PetscReal H,PetscScalar *Ke) 11c4762a1bSJed Brown { 12c4762a1bSJed Brown PetscFunctionBegin; 13c4762a1bSJed Brown Ke[0] = H/6.0; Ke[1] = -.125*H; Ke[2] = H/12.0; Ke[3] = -.125*H; 14c4762a1bSJed Brown Ke[4] = -.125*H; Ke[5] = H/6.0; Ke[6] = -.125*H; Ke[7] = H/12.0; 15c4762a1bSJed Brown Ke[8] = H/12.0; Ke[9] = -.125*H; Ke[10] = H/6.0; Ke[11] = -.125*H; 16c4762a1bSJed Brown Ke[12] = -.125*H; Ke[13] = H/12.0; Ke[14] = -.125*H; Ke[15] = H/6.0; 17c4762a1bSJed Brown PetscFunctionReturn(0); 18c4762a1bSJed Brown } 19c4762a1bSJed Brown 20c4762a1bSJed Brown int main(int argc,char **args) 21c4762a1bSJed Brown { 22c4762a1bSJed Brown Mat C; 23c4762a1bSJed Brown Vec u,b; 24c4762a1bSJed Brown PetscMPIInt size,rank; 25c4762a1bSJed Brown PetscInt i,m = 5,N,start,end,M,idx[4]; 26c4762a1bSJed Brown PetscInt j,nrsub,ncsub,*rsub,*csub,mystart,myend; 27c4762a1bSJed Brown PetscBool flg; 28c4762a1bSJed Brown PetscScalar one = 1.0,Ke[16],*vals; 29c4762a1bSJed Brown PetscReal h,norm; 30c4762a1bSJed Brown 31*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 32*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 33c4762a1bSJed Brown 34c4762a1bSJed Brown N = (m+1)*(m+1); /* dimension of matrix */ 35c4762a1bSJed Brown M = m*m; /* number of elements */ 36c4762a1bSJed Brown h = 1.0/m; /* mesh width */ 37*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 38*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 39c4762a1bSJed Brown 40c4762a1bSJed Brown /* Create stiffness matrix */ 41*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&C)); 42*9566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N)); 43*9566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 44*9566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 45c4762a1bSJed Brown 46c4762a1bSJed Brown start = rank*(M/size) + ((M%size) < rank ? (M%size) : rank); 47c4762a1bSJed Brown end = start + M/size + ((M%size) > rank); 48c4762a1bSJed Brown 49c4762a1bSJed Brown /* Form the element stiffness for the Laplacian */ 50*9566063dSJacob Faibussowitsch PetscCall(FormElementStiffness(h*h,Ke)); 51c4762a1bSJed Brown for (i=start; i<end; i++) { 52c4762a1bSJed Brown /* location of lower left corner of element */ 53c4762a1bSJed Brown /* node numbers for the four corners of element */ 54c4762a1bSJed Brown idx[0] = (m+1)*(i/m) + (i % m); 55c4762a1bSJed Brown idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1; 56*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(C,4,idx,4,idx,Ke,ADD_VALUES)); 57c4762a1bSJed Brown } 58*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 59*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 60c4762a1bSJed Brown 61c4762a1bSJed Brown /* Assemble the matrix again */ 62*9566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(C)); 63c4762a1bSJed Brown 64c4762a1bSJed Brown for (i=start; i<end; i++) { 65c4762a1bSJed Brown /* location of lower left corner of element */ 66c4762a1bSJed Brown /* node numbers for the four corners of element */ 67c4762a1bSJed Brown idx[0] = (m+1)*(i/m) + (i % m); 68c4762a1bSJed Brown idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1; 69*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(C,4,idx,4,idx,Ke,ADD_VALUES)); 70c4762a1bSJed Brown } 71*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 72*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 73c4762a1bSJed Brown 74c4762a1bSJed Brown /* Create test vectors */ 75*9566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&u)); 76*9566063dSJacob Faibussowitsch PetscCall(VecSetSizes(u,PETSC_DECIDE,N)); 77*9566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(u)); 78*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u,&b)); 79*9566063dSJacob Faibussowitsch PetscCall(VecSet(u,one)); 80c4762a1bSJed Brown 81c4762a1bSJed Brown /* Check error */ 82*9566063dSJacob Faibussowitsch PetscCall(MatMult(C,u,b)); 83*9566063dSJacob Faibussowitsch PetscCall(VecNorm(b,NORM_2,&norm)); 84c4762a1bSJed Brown if (norm > PETSC_SQRT_MACHINE_EPSILON) { 85*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Norm of error b %g should be near 0\n",(double)norm)); 86c4762a1bSJed Brown } 87c4762a1bSJed Brown 88c4762a1bSJed Brown /* Now test MatGetValues() */ 89*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-get_values",&flg)); 90c4762a1bSJed Brown if (flg) { 91*9566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(C,&mystart,&myend)); 92c4762a1bSJed Brown nrsub = myend - mystart; ncsub = 4; 93*9566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrsub*ncsub,&vals)); 94*9566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrsub,&rsub)); 95*9566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ncsub,&csub)); 96c4762a1bSJed Brown for (i=myend-1; i>=mystart; i--) rsub[myend-i-1] = i; 97c4762a1bSJed Brown for (i=0; i<ncsub; i++) csub[i] = 2*(ncsub-i) + mystart; 98*9566063dSJacob Faibussowitsch PetscCall(MatGetValues(C,nrsub,rsub,ncsub,csub,vals)); 99*9566063dSJacob Faibussowitsch PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 100*9566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"processor number %d: start=%" PetscInt_FMT ", end=%" PetscInt_FMT ", mystart=%" PetscInt_FMT ", myend=%" PetscInt_FMT "\n",rank,start,end,mystart,myend)); 101c4762a1bSJed Brown for (i=0; i<nrsub; i++) { 102c4762a1bSJed Brown for (j=0; j<ncsub; j++) { 103c4762a1bSJed Brown if (PetscImaginaryPart(vals[i*ncsub+j]) != 0.0) { 104*9566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD," C[%" PetscInt_FMT ", %" PetscInt_FMT "] = %g + %g i\n",rsub[i],csub[j],(double)PetscRealPart(vals[i*ncsub+j]),(double)PetscImaginaryPart(vals[i*ncsub+j]))); 105c4762a1bSJed Brown } else { 106*9566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD," C[%" PetscInt_FMT ", %" PetscInt_FMT "] = %g\n",rsub[i],csub[j],(double)PetscRealPart(vals[i*ncsub+j]))); 107c4762a1bSJed Brown } 108c4762a1bSJed Brown } 109c4762a1bSJed Brown } 110*9566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT)); 111*9566063dSJacob Faibussowitsch PetscCall(PetscFree(rsub)); 112*9566063dSJacob Faibussowitsch PetscCall(PetscFree(csub)); 113*9566063dSJacob Faibussowitsch PetscCall(PetscFree(vals)); 114c4762a1bSJed Brown } 115c4762a1bSJed Brown 116c4762a1bSJed Brown /* Free data structures */ 117*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 118*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 119*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 120*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 121b122ec5aSJacob Faibussowitsch return 0; 122c4762a1bSJed Brown } 123c4762a1bSJed Brown 124c4762a1bSJed Brown /*TEST 125c4762a1bSJed Brown 126c4762a1bSJed Brown test: 127c4762a1bSJed Brown nsize: 4 128c4762a1bSJed Brown 129c4762a1bSJed Brown TEST*/ 130