1c4762a1bSJed Brown static char help[] = "testing SeqDense matrices with an LDA (leading dimension of the user-allocated arrray) larger than M.\n"; 2c4762a1bSJed Brown /* 3c4762a1bSJed Brown * Example code testing SeqDense matrices with an LDA (leading dimension 4c4762a1bSJed Brown * of the user-allocated arrray) larger than M. 5c4762a1bSJed Brown * This example tests the functionality of MatSetValues(), MatMult(), 6c4762a1bSJed Brown * and MatMultTranspose(). 7c4762a1bSJed Brown */ 8c4762a1bSJed Brown 9c4762a1bSJed Brown #include <petscmat.h> 10c4762a1bSJed Brown 11c4762a1bSJed Brown int main(int argc,char **argv) 12c4762a1bSJed Brown { 13c4762a1bSJed Brown Mat A,A11,A12,A21,A22; 14c4762a1bSJed Brown Vec X,X1,X2,Y,Z,Z1,Z2; 15c4762a1bSJed Brown PetscScalar *a,*b,*x,*y,*z,v,one=1; 16c4762a1bSJed Brown PetscReal nrm; 17c4762a1bSJed Brown PetscErrorCode ierr; 18c4762a1bSJed Brown PetscInt size=8,size1=6,size2=2, i,j; 19c4762a1bSJed Brown PetscRandom rnd; 20c4762a1bSJed Brown 21c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,0,help);if (ierr) return ierr; 22*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_SELF,&rnd)); 23c4762a1bSJed Brown 24c4762a1bSJed Brown /* 25c4762a1bSJed Brown * Create matrix and three vectors: these are all normal 26c4762a1bSJed Brown */ 27*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(size*size,&a)); 28*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(size*size,&b)); 29c4762a1bSJed Brown for (i=0; i<size; i++) { 30c4762a1bSJed Brown for (j=0; j<size; j++) { 31*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValue(rnd,&a[i+j*size])); 32c4762a1bSJed Brown b[i+j*size] = a[i+j*size]; 33c4762a1bSJed Brown } 34c4762a1bSJed Brown } 35*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_SELF,&A)); 36*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,size,size,size,size)); 37*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATSEQDENSE)); 38*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqDenseSetPreallocation(A,a)); 39*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 40*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 41c4762a1bSJed Brown 42*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(size,&x)); 43c4762a1bSJed Brown for (i=0; i<size; i++) { 44c4762a1bSJed Brown x[i] = one; 45c4762a1bSJed Brown } 46*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeqWithArray(PETSC_COMM_SELF,1,size,x,&X)); 47*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyBegin(X)); 48*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyEnd(X)); 49c4762a1bSJed Brown 50*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(size,&y)); 51*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeqWithArray(PETSC_COMM_SELF,1,size,y,&Y)); 52*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyBegin(Y)); 53*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyEnd(Y)); 54c4762a1bSJed Brown 55*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(size,&z)); 56*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeqWithArray(PETSC_COMM_SELF,1,size,z,&Z)); 57*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyBegin(Z)); 58*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyEnd(Z)); 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* 61c4762a1bSJed Brown * Now create submatrices and subvectors 62c4762a1bSJed Brown */ 63*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_SELF,&A11)); 64*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A11,size1,size1,size1,size1)); 65*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A11,MATSEQDENSE)); 66*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqDenseSetPreallocation(A11,b)); 67*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseSetLDA(A11,size)); 68*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A11,MAT_FINAL_ASSEMBLY)); 69*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A11,MAT_FINAL_ASSEMBLY)); 70c4762a1bSJed Brown 71*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_SELF,&A12)); 72*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A12,size1,size2,size1,size2)); 73*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A12,MATSEQDENSE)); 74*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqDenseSetPreallocation(A12,b+size1*size)); 75*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseSetLDA(A12,size)); 76*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A12,MAT_FINAL_ASSEMBLY)); 77*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A12,MAT_FINAL_ASSEMBLY)); 78c4762a1bSJed Brown 79*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_SELF,&A21)); 80*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A21,size2,size1,size2,size1)); 81*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A21,MATSEQDENSE)); 82*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqDenseSetPreallocation(A21,b+size1)); 83*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseSetLDA(A21,size)); 84*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A21,MAT_FINAL_ASSEMBLY)); 85*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A21,MAT_FINAL_ASSEMBLY)); 86c4762a1bSJed Brown 87*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_SELF,&A22)); 88*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A22,size2,size2,size2,size2)); 89*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A22,MATSEQDENSE)); 90*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqDenseSetPreallocation(A22,b+size1*size+size1)); 91*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseSetLDA(A22,size)); 92*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A22,MAT_FINAL_ASSEMBLY)); 93*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A22,MAT_FINAL_ASSEMBLY)); 94c4762a1bSJed Brown 95*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeqWithArray(PETSC_COMM_SELF,1,size1,x,&X1)); 96*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeqWithArray(PETSC_COMM_SELF,1,size2,x+size1,&X2)); 97*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeqWithArray(PETSC_COMM_SELF,1,size1,z,&Z1)); 98*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeqWithArray(PETSC_COMM_SELF,1,size2,z+size1,&Z2)); 99c4762a1bSJed Brown 100c4762a1bSJed Brown /* 101c4762a1bSJed Brown * Now multiple matrix times input in two ways; 102c4762a1bSJed Brown * compare the results 103c4762a1bSJed Brown */ 104*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,X,Y)); 105*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A11,X1,Z1)); 106*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAdd(A12,X2,Z1,Z1)); 107*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A22,X2,Z2)); 108*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAdd(A21,X1,Z2,Z2)); 109*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(Z,-1.0,Y)); 110*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(Z,NORM_2,&nrm)); 111c4762a1bSJed Brown if (nrm > 100.0*PETSC_MACHINE_EPSILON) { 112*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Test1; error norm=%g\n",(double)nrm)); 113c4762a1bSJed Brown } 114c4762a1bSJed Brown 115c4762a1bSJed Brown /* 116c4762a1bSJed Brown * Next test: change both matrices 117c4762a1bSJed Brown */ 118*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValue(rnd,&v)); 119c4762a1bSJed Brown i = 1; 120c4762a1bSJed Brown j = size-2; 121*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES)); 122c4762a1bSJed Brown j -= size1; 123*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A12,1,&i,1,&j,&v,INSERT_VALUES)); 124*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValue(rnd,&v)); 125c4762a1bSJed Brown i = j = size1+1; 126*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES)); 127c4762a1bSJed Brown i =j = 1; 128*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A22,1,&i,1,&j,&v,INSERT_VALUES)); 129*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 130*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 131*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A12,MAT_FINAL_ASSEMBLY)); 132*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A12,MAT_FINAL_ASSEMBLY)); 133*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A22,MAT_FINAL_ASSEMBLY)); 134*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A22,MAT_FINAL_ASSEMBLY)); 135c4762a1bSJed Brown 136*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,X,Y)); 137*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A11,X1,Z1)); 138*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAdd(A12,X2,Z1,Z1)); 139*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A22,X2,Z2)); 140*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAdd(A21,X1,Z2,Z2)); 141*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(Z,-1.0,Y)); 142*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(Z,NORM_2,&nrm)); 143c4762a1bSJed Brown if (nrm > 100.0*PETSC_MACHINE_EPSILON) { 144*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Test2; error norm=%g\n",(double)nrm)); 145c4762a1bSJed Brown } 146c4762a1bSJed Brown 147c4762a1bSJed Brown /* 148c4762a1bSJed Brown * Transpose product 149c4762a1bSJed Brown */ 150*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(A,X,Y)); 151*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(A11,X1,Z1)); 152*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTransposeAdd(A21,X2,Z1,Z1)); 153*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(A22,X2,Z2)); 154*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTransposeAdd(A12,X1,Z2,Z2)); 155*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(Z,-1.0,Y)); 156*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(Z,NORM_2,&nrm)); 157c4762a1bSJed Brown if (nrm > 100.0*PETSC_MACHINE_EPSILON) { 158*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Test3; error norm=%g\n",(double)nrm)); 159c4762a1bSJed Brown } 160*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(a)); 161*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(b)); 162*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(x)); 163*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(y)); 164*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(z)); 165*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 166*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A11)); 167*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A12)); 168*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A21)); 169*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A22)); 170c4762a1bSJed Brown 171*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&X)); 172*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&Y)); 173*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&Z)); 174c4762a1bSJed Brown 175*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&X1)); 176*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&X2)); 177*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&Z1)); 178*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&Z2)); 179*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rnd)); 180c4762a1bSJed Brown ierr = PetscFinalize(); 181c4762a1bSJed Brown return ierr; 182c4762a1bSJed Brown } 183c4762a1bSJed Brown 184c4762a1bSJed Brown /*TEST 185c4762a1bSJed Brown 186c4762a1bSJed Brown test: 187c4762a1bSJed Brown 188c4762a1bSJed Brown TEST*/ 189