1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests various routines in MatKAIJ format.\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown #define IMAX 15 6c4762a1bSJed Brown 7c4762a1bSJed Brown int main(int argc,char **args) 8c4762a1bSJed Brown { 9c4762a1bSJed Brown Mat A,B,TA; 10c4762a1bSJed Brown PetscScalar *S,*T; 11c4762a1bSJed Brown PetscViewer fd; 12c4762a1bSJed Brown char file[PETSC_MAX_PATH_LEN]; 13c4762a1bSJed Brown PetscInt m,n,M,N,p=1,q=1,i,j; 14c4762a1bSJed Brown PetscMPIInt rank,size; 15c4762a1bSJed Brown PetscBool flg; 16c4762a1bSJed Brown 17*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 185f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 195f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 20c4762a1bSJed Brown 21910cf402Sprj- /* Load AIJ matrix A */ 225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg)); 2328b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option"); 245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd)); 255f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 265f80ce2aSJacob Faibussowitsch CHKERRQ(MatLoad(A,fd)); 275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&fd)); 28c4762a1bSJed Brown 29c4762a1bSJed Brown /* Get dof, then create S and T */ 305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-p",&p,NULL)); 315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-q",&q,NULL)); 325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(p*q,&S,p*q,&T)); 33c4762a1bSJed Brown for (i=0; i<p*q; i++) S[i] = 0; 34c4762a1bSJed Brown 35c4762a1bSJed Brown for (i=0; i<p; i++) { 36c4762a1bSJed Brown for (j=0; j<q; j++) { 37910cf402Sprj- /* Set some random non-zero values */ 38c4762a1bSJed Brown S[i+p*j] = ((PetscReal) ((i+1)*(j+1))) / ((PetscReal) (p+q)); 39c4762a1bSJed Brown T[i+p*j] = ((PetscReal) ((p-i)+j)) / ((PetscReal) (p*q)); 40c4762a1bSJed Brown } 41c4762a1bSJed Brown } 42c4762a1bSJed Brown 43c4762a1bSJed Brown /* Test KAIJ when both S & T are not NULL */ 44c4762a1bSJed Brown 45910cf402Sprj- /* Create KAIJ matrix TA */ 465f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateKAIJ(A,p,q,S,T,&TA)); 475f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(A,&m,&n)); 485f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(A,&M,&N)); 49c4762a1bSJed Brown 505f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(TA,MATAIJ,MAT_INITIAL_MATRIX,&B)); 51c4762a1bSJed Brown 52910cf402Sprj- /* Test MatKAIJGetScaledIdentity() */ 535f80ce2aSJacob Faibussowitsch CHKERRQ(MatKAIJGetScaledIdentity(TA,&flg)); 5428b400f6SJacob Faibussowitsch PetscCheck(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in Test 1: MatKAIJGetScaledIdentity()"); 55c4762a1bSJed Brown /* Test MatMult() */ 565f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(TA,B,10,&flg)); 5728b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error in Test 1: MatMult() for KAIJ matrix"); 58c4762a1bSJed Brown /* Test MatMultAdd() */ 595f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAddEqual(TA,B,10,&flg)); 6028b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error in Test 1: MatMultAdd() for KAIJ matrix"); 61c4762a1bSJed Brown 625f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&TA)); 635f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 64c4762a1bSJed Brown 65c4762a1bSJed Brown /* Test KAIJ when S is NULL */ 66c4762a1bSJed Brown 67910cf402Sprj- /* Create KAIJ matrix TA */ 685f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateKAIJ(A,p,q,NULL,T,&TA)); 695f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(A,&m,&n)); 705f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(A,&M,&N)); 71c4762a1bSJed Brown 725f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(TA,MATAIJ,MAT_INITIAL_MATRIX,&B)); 73c4762a1bSJed Brown 74910cf402Sprj- /* Test MatKAIJGetScaledIdentity() */ 755f80ce2aSJacob Faibussowitsch CHKERRQ(MatKAIJGetScaledIdentity(TA,&flg)); 7628b400f6SJacob Faibussowitsch PetscCheck(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in Test 2: MatKAIJGetScaledIdentity()"); 77c4762a1bSJed Brown /* Test MatMult() */ 785f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(TA,B,10,&flg)); 7928b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error in Test 2: MatMult() for KAIJ matrix"); 80c4762a1bSJed Brown /* Test MatMultAdd() */ 815f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAddEqual(TA,B,10,&flg)); 8228b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error in Test 2: MatMultAdd() for KAIJ matrix"); 83c4762a1bSJed Brown 845f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&TA)); 855f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 86c4762a1bSJed Brown 87c4762a1bSJed Brown /* Test KAIJ when T is NULL */ 88c4762a1bSJed Brown 89910cf402Sprj- /* Create KAIJ matrix TA */ 905f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateKAIJ(A,p,q,S,NULL,&TA)); 915f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(A,&m,&n)); 925f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(A,&M,&N)); 93c4762a1bSJed Brown 945f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(TA,MATAIJ,MAT_INITIAL_MATRIX,&B)); 95c4762a1bSJed Brown 96910cf402Sprj- /* Test MatKAIJGetScaledIdentity() */ 975f80ce2aSJacob Faibussowitsch CHKERRQ(MatKAIJGetScaledIdentity(TA,&flg)); 9828b400f6SJacob Faibussowitsch PetscCheck(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in Test 3: MatKAIJGetScaledIdentity()"); 99c4762a1bSJed Brown /* Test MatMult() */ 1005f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(TA,B,10,&flg)); 10128b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error in Test 3: MatMult() for KAIJ matrix"); 102c4762a1bSJed Brown /* Test MatMultAdd() */ 1035f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAddEqual(TA,B,10,&flg)); 10428b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error in Test 3: MatMultAdd() for KAIJ matrix"); 105c4762a1bSJed Brown 1065f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&TA)); 1075f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 108c4762a1bSJed Brown 109c4762a1bSJed Brown /* Test KAIJ when T is is an identity matrix */ 110c4762a1bSJed Brown 111c4762a1bSJed Brown if (p == q) { 112c4762a1bSJed Brown for (i=0; i<p; i++) { 113c4762a1bSJed Brown for (j=0; j<q; j++) { 114c4762a1bSJed Brown if (i==j) T[i+j*p] = 1.0; 115c4762a1bSJed Brown else T[i+j*p] = 0.0; 116c4762a1bSJed Brown } 117c4762a1bSJed Brown } 118c4762a1bSJed Brown 119910cf402Sprj- /* Create KAIJ matrix TA */ 1205f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateKAIJ(A,p,q,S,T,&TA)); 1215f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(A,&m,&n)); 1225f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(A,&M,&N)); 123c4762a1bSJed Brown 1245f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(TA,MATAIJ,MAT_INITIAL_MATRIX,&B)); 125c4762a1bSJed Brown 126910cf402Sprj- /* Test MatKAIJGetScaledIdentity() */ 1275f80ce2aSJacob Faibussowitsch CHKERRQ(MatKAIJGetScaledIdentity(TA,&flg)); 12828b400f6SJacob Faibussowitsch PetscCheck(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in Test 4: MatKAIJGetScaledIdentity()"); 129c4762a1bSJed Brown /* Test MatMult() */ 1305f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(TA,B,10,&flg)); 13128b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error in Test 4: MatMult() for KAIJ matrix"); 132c4762a1bSJed Brown /* Test MatMultAdd() */ 1335f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAddEqual(TA,B,10,&flg)); 13428b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error in Test 4: MatMultAdd() for KAIJ matrix"); 135c4762a1bSJed Brown 1365f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&TA)); 1375f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 138910cf402Sprj- 1395f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateKAIJ(A,p,q,NULL,T,&TA)); 140910cf402Sprj- /* Test MatKAIJGetScaledIdentity() */ 1415f80ce2aSJacob Faibussowitsch CHKERRQ(MatKAIJGetScaledIdentity(TA,&flg)); 14228b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in Test 5: MatKAIJGetScaledIdentity()"); 1435f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&TA)); 144910cf402Sprj- 145910cf402Sprj- for (i=0; i<p; i++) { 146910cf402Sprj- for (j=0; j<q; j++) { 147910cf402Sprj- if (i==j) S[i+j*p] = T[i+j*p] = 2.0; 148910cf402Sprj- else S[i+j*p] = T[i+j*p] = 0.0; 149910cf402Sprj- } 150910cf402Sprj- } 1515f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateKAIJ(A,p,q,S,T,&TA)); 152910cf402Sprj- /* Test MatKAIJGetScaledIdentity() */ 1535f80ce2aSJacob Faibussowitsch CHKERRQ(MatKAIJGetScaledIdentity(TA,&flg)); 15428b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in Test 6: MatKAIJGetScaledIdentity()"); 1555f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&TA)); 156c4762a1bSJed Brown } 157c4762a1bSJed Brown 158c4762a1bSJed Brown /* Done with all tests */ 159c4762a1bSJed Brown 1605f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(S,T)); 1615f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 162*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 163*b122ec5aSJacob Faibussowitsch return 0; 164c4762a1bSJed Brown } 165c4762a1bSJed Brown 166c4762a1bSJed Brown /*TEST 167c4762a1bSJed Brown 168c4762a1bSJed Brown build: 169c4762a1bSJed Brown requires: !complex 170c4762a1bSJed Brown 171c4762a1bSJed Brown test: 172dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 173c4762a1bSJed Brown output_file: output/ex176.out 174c4762a1bSJed Brown nsize: {{1 2 3 4}} 175c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/small -p {{2 3 7}} -q {{3 7}} -viewer_binary_skip_info 176c4762a1bSJed Brown 177c4762a1bSJed Brown TEST*/ 178