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 PetscErrorCode ierr; 16c4762a1bSJed Brown PetscBool flg; 17c4762a1bSJed Brown 18c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 19ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 20ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 21c4762a1bSJed Brown 22910cf402Sprj- /* Load AIJ matrix A */ 23589a23caSBarry Smith ierr = PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg);CHKERRQ(ierr); 24c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option"); 25c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);CHKERRQ(ierr); 26c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 27c4762a1bSJed Brown ierr = MatLoad(A,fd);CHKERRQ(ierr); 28c4762a1bSJed Brown ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr); 29c4762a1bSJed Brown 30c4762a1bSJed Brown /* Get dof, then create S and T */ 31c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-p",&p,NULL);CHKERRQ(ierr); 32c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-q",&q,NULL);CHKERRQ(ierr); 33c4762a1bSJed Brown ierr = PetscMalloc2(p*q,&S,p*q,&T);CHKERRQ(ierr); 34c4762a1bSJed Brown for (i=0; i<p*q; i++) S[i] = 0; 35c4762a1bSJed Brown 36c4762a1bSJed Brown for (i=0; i<p; i++) { 37c4762a1bSJed Brown for (j=0; j<q; j++) { 38910cf402Sprj- /* Set some random non-zero values */ 39c4762a1bSJed Brown S[i+p*j] = ((PetscReal) ((i+1)*(j+1))) / ((PetscReal) (p+q)); 40c4762a1bSJed Brown T[i+p*j] = ((PetscReal) ((p-i)+j)) / ((PetscReal) (p*q)); 41c4762a1bSJed Brown } 42c4762a1bSJed Brown } 43c4762a1bSJed Brown 44c4762a1bSJed Brown /* Test KAIJ when both S & T are not NULL */ 45c4762a1bSJed Brown 46910cf402Sprj- /* Create KAIJ matrix TA */ 47c4762a1bSJed Brown ierr = MatCreateKAIJ(A,p,q,S,T,&TA);CHKERRQ(ierr); 48c4762a1bSJed Brown ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 49c4762a1bSJed Brown ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 50c4762a1bSJed Brown 51c4762a1bSJed Brown ierr = MatConvert(TA,MATAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); 52c4762a1bSJed Brown 53910cf402Sprj- /* Test MatKAIJGetScaledIdentity() */ 54910cf402Sprj- ierr = MatKAIJGetScaledIdentity(TA,&flg);CHKERRQ(ierr); 55910cf402Sprj- if (flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in Test 1: MatKAIJGetScaledIdentity()"); 56c4762a1bSJed Brown /* Test MatMult() */ 57c4762a1bSJed Brown ierr = MatMultEqual(TA,B,10,&flg);CHKERRQ(ierr); 58c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error in Test 1: MatMult() for KAIJ matrix"); 59c4762a1bSJed Brown /* Test MatMultAdd() */ 60c4762a1bSJed Brown ierr = MatMultAddEqual(TA,B,10,&flg);CHKERRQ(ierr); 61c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error in Test 1: MatMultAdd() for KAIJ matrix"); 62c4762a1bSJed Brown 63c4762a1bSJed Brown ierr = MatDestroy(&TA);CHKERRQ(ierr); 64c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 65c4762a1bSJed Brown 66c4762a1bSJed Brown /* Test KAIJ when S is NULL */ 67c4762a1bSJed Brown 68910cf402Sprj- /* Create KAIJ matrix TA */ 69c4762a1bSJed Brown ierr = MatCreateKAIJ(A,p,q,NULL,T,&TA);CHKERRQ(ierr); 70c4762a1bSJed Brown ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 71c4762a1bSJed Brown ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 72c4762a1bSJed Brown 73c4762a1bSJed Brown ierr = MatConvert(TA,MATAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); 74c4762a1bSJed Brown 75910cf402Sprj- /* Test MatKAIJGetScaledIdentity() */ 76910cf402Sprj- ierr = MatKAIJGetScaledIdentity(TA,&flg);CHKERRQ(ierr); 77910cf402Sprj- if (flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in Test 2: MatKAIJGetScaledIdentity()"); 78c4762a1bSJed Brown /* Test MatMult() */ 79c4762a1bSJed Brown ierr = MatMultEqual(TA,B,10,&flg);CHKERRQ(ierr); 80c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error in Test 2: MatMult() for KAIJ matrix"); 81c4762a1bSJed Brown /* Test MatMultAdd() */ 82c4762a1bSJed Brown ierr = MatMultAddEqual(TA,B,10,&flg);CHKERRQ(ierr); 83c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error in Test 2: MatMultAdd() for KAIJ matrix"); 84c4762a1bSJed Brown 85c4762a1bSJed Brown ierr = MatDestroy(&TA);CHKERRQ(ierr); 86c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 87c4762a1bSJed Brown 88c4762a1bSJed Brown /* Test KAIJ when T is NULL */ 89c4762a1bSJed Brown 90910cf402Sprj- /* Create KAIJ matrix TA */ 91c4762a1bSJed Brown ierr = MatCreateKAIJ(A,p,q,S,NULL,&TA);CHKERRQ(ierr); 92c4762a1bSJed Brown ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 93c4762a1bSJed Brown ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 94c4762a1bSJed Brown 95c4762a1bSJed Brown ierr = MatConvert(TA,MATAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); 96c4762a1bSJed Brown 97910cf402Sprj- /* Test MatKAIJGetScaledIdentity() */ 98910cf402Sprj- ierr = MatKAIJGetScaledIdentity(TA,&flg);CHKERRQ(ierr); 99910cf402Sprj- if (flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in Test 3: MatKAIJGetScaledIdentity()"); 100c4762a1bSJed Brown /* Test MatMult() */ 101c4762a1bSJed Brown ierr = MatMultEqual(TA,B,10,&flg);CHKERRQ(ierr); 102c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error in Test 3: MatMult() for KAIJ matrix"); 103c4762a1bSJed Brown /* Test MatMultAdd() */ 104c4762a1bSJed Brown ierr = MatMultAddEqual(TA,B,10,&flg);CHKERRQ(ierr); 105c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error in Test 3: MatMultAdd() for KAIJ matrix"); 106c4762a1bSJed Brown 107c4762a1bSJed Brown ierr = MatDestroy(&TA);CHKERRQ(ierr); 108c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 109c4762a1bSJed Brown 110c4762a1bSJed Brown /* Test KAIJ when T is is an identity matrix */ 111c4762a1bSJed Brown 112c4762a1bSJed Brown if (p == q) { 113c4762a1bSJed Brown for (i=0; i<p; i++) { 114c4762a1bSJed Brown for (j=0; j<q; j++) { 115c4762a1bSJed Brown if (i==j) T[i+j*p] = 1.0; 116c4762a1bSJed Brown else T[i+j*p] = 0.0; 117c4762a1bSJed Brown } 118c4762a1bSJed Brown } 119c4762a1bSJed Brown 120910cf402Sprj- /* Create KAIJ matrix TA */ 121c4762a1bSJed Brown ierr = MatCreateKAIJ(A,p,q,S,T,&TA);CHKERRQ(ierr); 122c4762a1bSJed Brown ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 123c4762a1bSJed Brown ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 124c4762a1bSJed Brown 125c4762a1bSJed Brown ierr = MatConvert(TA,MATAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); 126c4762a1bSJed Brown 127910cf402Sprj- /* Test MatKAIJGetScaledIdentity() */ 128910cf402Sprj- ierr = MatKAIJGetScaledIdentity(TA,&flg);CHKERRQ(ierr); 129910cf402Sprj- if (flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in Test 4: MatKAIJGetScaledIdentity()"); 130c4762a1bSJed Brown /* Test MatMult() */ 131c4762a1bSJed Brown ierr = MatMultEqual(TA,B,10,&flg);CHKERRQ(ierr); 132c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error in Test 4: MatMult() for KAIJ matrix"); 133c4762a1bSJed Brown /* Test MatMultAdd() */ 134c4762a1bSJed Brown ierr = MatMultAddEqual(TA,B,10,&flg);CHKERRQ(ierr); 135c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error in Test 4: MatMultAdd() for KAIJ matrix"); 136c4762a1bSJed Brown 137c4762a1bSJed Brown ierr = MatDestroy(&TA);CHKERRQ(ierr); 138c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 139910cf402Sprj- 140910cf402Sprj- ierr = MatCreateKAIJ(A,p,q,NULL,T,&TA);CHKERRQ(ierr); 141910cf402Sprj- /* Test MatKAIJGetScaledIdentity() */ 142910cf402Sprj- ierr = MatKAIJGetScaledIdentity(TA,&flg);CHKERRQ(ierr); 143910cf402Sprj- if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in Test 5: MatKAIJGetScaledIdentity()"); 144910cf402Sprj- ierr = MatDestroy(&TA);CHKERRQ(ierr); 145910cf402Sprj- 146910cf402Sprj- for (i=0; i<p; i++) { 147910cf402Sprj- for (j=0; j<q; j++) { 148910cf402Sprj- if (i==j) S[i+j*p] = T[i+j*p] = 2.0; 149910cf402Sprj- else S[i+j*p] = T[i+j*p] = 0.0; 150910cf402Sprj- } 151910cf402Sprj- } 152910cf402Sprj- ierr = MatCreateKAIJ(A,p,q,S,T,&TA);CHKERRQ(ierr); 153910cf402Sprj- /* Test MatKAIJGetScaledIdentity() */ 154910cf402Sprj- ierr = MatKAIJGetScaledIdentity(TA,&flg);CHKERRQ(ierr); 155910cf402Sprj- if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Error in Test 6: MatKAIJGetScaledIdentity()"); 156910cf402Sprj- ierr = MatDestroy(&TA);CHKERRQ(ierr); 157c4762a1bSJed Brown } 158c4762a1bSJed Brown 159c4762a1bSJed Brown /* Done with all tests */ 160c4762a1bSJed Brown 161c4762a1bSJed Brown ierr = PetscFree2(S,T);CHKERRQ(ierr); 162c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 163c4762a1bSJed Brown ierr = PetscFinalize(); 164c4762a1bSJed Brown return ierr; 165c4762a1bSJed Brown } 166c4762a1bSJed Brown 167c4762a1bSJed Brown /*TEST 168c4762a1bSJed Brown 169c4762a1bSJed Brown build: 170c4762a1bSJed Brown requires: !complex 171c4762a1bSJed Brown 172c4762a1bSJed Brown test: 173*dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 174c4762a1bSJed Brown output_file: output/ex176.out 175c4762a1bSJed Brown nsize: {{1 2 3 4}} 176c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/small -p {{2 3 7}} -q {{3 7}} -viewer_binary_skip_info 177c4762a1bSJed Brown 178c4762a1bSJed Brown TEST*/ 179