1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[] = "Tests various routines in MatKAIJ format.\n"; 3*c4762a1bSJed Brown 4*c4762a1bSJed Brown #include <petscmat.h> 5*c4762a1bSJed Brown #define IMAX 15 6*c4762a1bSJed Brown 7*c4762a1bSJed Brown int main(int argc,char **args) 8*c4762a1bSJed Brown { 9*c4762a1bSJed Brown Mat A,B,TA; 10*c4762a1bSJed Brown PetscScalar *S,*T; 11*c4762a1bSJed Brown PetscViewer fd; 12*c4762a1bSJed Brown char file[PETSC_MAX_PATH_LEN]; 13*c4762a1bSJed Brown PetscInt m,n,M,N,p=1,q=1,i,j; 14*c4762a1bSJed Brown PetscMPIInt rank,size; 15*c4762a1bSJed Brown PetscErrorCode ierr; 16*c4762a1bSJed Brown PetscBool flg; 17*c4762a1bSJed Brown 18*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 19*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 20*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 21*c4762a1bSJed Brown 22*c4762a1bSJed Brown /* Load aij matrix A */ 23*c4762a1bSJed Brown ierr = PetscOptionsGetString(NULL,NULL,"-f",file,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 24*c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option"); 25*c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);CHKERRQ(ierr); 26*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 27*c4762a1bSJed Brown ierr = MatLoad(A,fd);CHKERRQ(ierr); 28*c4762a1bSJed Brown ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr); 29*c4762a1bSJed Brown 30*c4762a1bSJed Brown /* Get dof, then create S and T */ 31*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-p",&p,NULL);CHKERRQ(ierr); 32*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-q",&q,NULL);CHKERRQ(ierr); 33*c4762a1bSJed Brown ierr = PetscMalloc2(p*q,&S,p*q,&T);CHKERRQ(ierr); 34*c4762a1bSJed Brown for (i=0; i<p*q; i++) S[i] = 0; 35*c4762a1bSJed Brown 36*c4762a1bSJed Brown for (i=0; i<p; i++) { 37*c4762a1bSJed Brown for (j = 0; j<q; j++) { 38*c4762a1bSJed Brown /* set some random non-zero values */ 39*c4762a1bSJed Brown S[i+p*j] = ((PetscReal) ((i+1)*(j+1))) / ((PetscReal) (p+q)); 40*c4762a1bSJed Brown T[i+p*j] = ((PetscReal) ((p-i)+j)) / ((PetscReal) (p*q)); 41*c4762a1bSJed Brown } 42*c4762a1bSJed Brown } 43*c4762a1bSJed Brown 44*c4762a1bSJed Brown /* Test KAIJ when both S & T are not NULL */ 45*c4762a1bSJed Brown 46*c4762a1bSJed Brown /* create kaij matrix TA */ 47*c4762a1bSJed Brown ierr = MatCreateKAIJ(A,p,q,S,T,&TA);CHKERRQ(ierr); 48*c4762a1bSJed Brown ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 49*c4762a1bSJed Brown ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 50*c4762a1bSJed Brown 51*c4762a1bSJed Brown ierr = MatConvert(TA,MATAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); 52*c4762a1bSJed Brown 53*c4762a1bSJed Brown /* Test MatMult() */ 54*c4762a1bSJed Brown ierr = MatMultEqual(TA,B,10,&flg);CHKERRQ(ierr); 55*c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error in Test 1: MatMult() for KAIJ matrix"); 56*c4762a1bSJed Brown /* Test MatMultAdd() */ 57*c4762a1bSJed Brown ierr = MatMultAddEqual(TA,B,10,&flg);CHKERRQ(ierr); 58*c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error in Test 1: MatMultAdd() for KAIJ matrix"); 59*c4762a1bSJed Brown 60*c4762a1bSJed Brown ierr = MatDestroy(&TA);CHKERRQ(ierr); 61*c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 62*c4762a1bSJed Brown 63*c4762a1bSJed Brown /* Test KAIJ when S is NULL */ 64*c4762a1bSJed Brown 65*c4762a1bSJed Brown /* create kaij matrix TA */ 66*c4762a1bSJed Brown ierr = MatCreateKAIJ(A,p,q,NULL,T,&TA);CHKERRQ(ierr); 67*c4762a1bSJed Brown ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 68*c4762a1bSJed Brown ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 69*c4762a1bSJed Brown 70*c4762a1bSJed Brown ierr = MatConvert(TA,MATAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); 71*c4762a1bSJed Brown 72*c4762a1bSJed Brown /* Test MatMult() */ 73*c4762a1bSJed Brown ierr = MatMultEqual(TA,B,10,&flg);CHKERRQ(ierr); 74*c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error in Test 2: MatMult() for KAIJ matrix"); 75*c4762a1bSJed Brown /* Test MatMultAdd() */ 76*c4762a1bSJed Brown ierr = MatMultAddEqual(TA,B,10,&flg);CHKERRQ(ierr); 77*c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error in Test 2: MatMultAdd() for KAIJ matrix"); 78*c4762a1bSJed Brown 79*c4762a1bSJed Brown ierr = MatDestroy(&TA);CHKERRQ(ierr); 80*c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 81*c4762a1bSJed Brown 82*c4762a1bSJed Brown /* Test KAIJ when T is NULL */ 83*c4762a1bSJed Brown 84*c4762a1bSJed Brown /* create kaij matrix TA */ 85*c4762a1bSJed Brown ierr = MatCreateKAIJ(A,p,q,S,NULL,&TA);CHKERRQ(ierr); 86*c4762a1bSJed Brown ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 87*c4762a1bSJed Brown ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 88*c4762a1bSJed Brown 89*c4762a1bSJed Brown ierr = MatConvert(TA,MATAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); 90*c4762a1bSJed Brown 91*c4762a1bSJed Brown /* Test MatMult() */ 92*c4762a1bSJed Brown ierr = MatMultEqual(TA,B,10,&flg);CHKERRQ(ierr); 93*c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error in Test 3: MatMult() for KAIJ matrix"); 94*c4762a1bSJed Brown /* Test MatMultAdd() */ 95*c4762a1bSJed Brown ierr = MatMultAddEqual(TA,B,10,&flg);CHKERRQ(ierr); 96*c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error in Test 3: MatMultAdd() for KAIJ matrix"); 97*c4762a1bSJed Brown 98*c4762a1bSJed Brown ierr = MatDestroy(&TA);CHKERRQ(ierr); 99*c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 100*c4762a1bSJed Brown 101*c4762a1bSJed Brown /* Test KAIJ when T is is an identity matrix */ 102*c4762a1bSJed Brown 103*c4762a1bSJed Brown if (p == q) { 104*c4762a1bSJed Brown for (i=0; i<p; i++) { 105*c4762a1bSJed Brown for (j=0; j<q; j++) { 106*c4762a1bSJed Brown if (i==j) T[i+j*p] = 1.0; 107*c4762a1bSJed Brown else T[i+j*p] = 0.0; 108*c4762a1bSJed Brown } 109*c4762a1bSJed Brown } 110*c4762a1bSJed Brown 111*c4762a1bSJed Brown /* create kaij matrix TA */ 112*c4762a1bSJed Brown ierr = MatCreateKAIJ(A,p,q,S,T,&TA);CHKERRQ(ierr); 113*c4762a1bSJed Brown ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 114*c4762a1bSJed Brown ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 115*c4762a1bSJed Brown 116*c4762a1bSJed Brown ierr = MatConvert(TA,MATAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); 117*c4762a1bSJed Brown 118*c4762a1bSJed Brown /* Test MatMult() */ 119*c4762a1bSJed Brown ierr = MatMultEqual(TA,B,10,&flg);CHKERRQ(ierr); 120*c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error in Test 4: MatMult() for KAIJ matrix"); 121*c4762a1bSJed Brown /* Test MatMultAdd() */ 122*c4762a1bSJed Brown ierr = MatMultAddEqual(TA,B,10,&flg);CHKERRQ(ierr); 123*c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error in Test 4: MatMultAdd() for KAIJ matrix"); 124*c4762a1bSJed Brown 125*c4762a1bSJed Brown ierr = MatDestroy(&TA);CHKERRQ(ierr); 126*c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 127*c4762a1bSJed Brown } 128*c4762a1bSJed Brown 129*c4762a1bSJed Brown /* Done with all tests */ 130*c4762a1bSJed Brown 131*c4762a1bSJed Brown ierr = PetscFree2(S,T);CHKERRQ(ierr); 132*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 133*c4762a1bSJed Brown ierr = PetscFinalize(); 134*c4762a1bSJed Brown return ierr; 135*c4762a1bSJed Brown } 136*c4762a1bSJed Brown 137*c4762a1bSJed Brown /*TEST 138*c4762a1bSJed Brown 139*c4762a1bSJed Brown build: 140*c4762a1bSJed Brown requires: !complex 141*c4762a1bSJed Brown 142*c4762a1bSJed Brown test: 143*c4762a1bSJed Brown requires: datafilespath !complex double !define(PETSC_USE_64BIT_INDICES) 144*c4762a1bSJed Brown output_file: output/ex176.out 145*c4762a1bSJed Brown nsize: {{1 2 3 4}} 146*c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/small -p {{2 3 7}} -q {{3 7}} -viewer_binary_skip_info 147*c4762a1bSJed Brown 148*c4762a1bSJed Brown TEST*/ 149