16280154eSStefano Zampini #include <petscmat.h> 26280154eSStefano Zampini 36280154eSStefano Zampini static char help[] = "Tests MatMatMult with MAT_REUSE_MATRIX and already allocated dense result.\n\n"; 46280154eSStefano Zampini 56280154eSStefano Zampini static PetscErrorCode CheckLocal(Mat A, Mat B, PetscScalar *a, PetscScalar *b) 66280154eSStefano Zampini { 76280154eSStefano Zampini PetscErrorCode ierr; 86280154eSStefano Zampini PetscBool wA = PETSC_FALSE, wB = PETSC_FALSE; 96280154eSStefano Zampini 106280154eSStefano Zampini PetscFunctionBegin; 116280154eSStefano Zampini if (a) { 126280154eSStefano Zampini const PetscScalar *Aa; 136280154eSStefano Zampini ierr = MatDenseGetArrayRead(A,&Aa);CHKERRQ(ierr); 146280154eSStefano Zampini wA = (PetscBool)(a != Aa); 156280154eSStefano Zampini ierr = MatDenseRestoreArrayRead(A,&Aa);CHKERRQ(ierr); 166280154eSStefano Zampini } 176280154eSStefano Zampini if (b) { 186280154eSStefano Zampini const PetscScalar *Bb; 196280154eSStefano Zampini ierr = MatDenseGetArrayRead(B,&Bb);CHKERRQ(ierr); 206280154eSStefano Zampini wB = (PetscBool)(b != Bb); 216280154eSStefano Zampini ierr = MatDenseRestoreArrayRead(B,&Bb);CHKERRQ(ierr); 226280154eSStefano Zampini } 23456288a8SStefano Zampini if (wA || wB) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong array in first Mat? %d, Wrong array in second Mat? %d",wA,wB); 246280154eSStefano Zampini PetscFunctionReturn(0); 256280154eSStefano Zampini } 266280154eSStefano Zampini 276280154eSStefano Zampini int main(int argc,char **args) 286280154eSStefano Zampini { 29456288a8SStefano Zampini Mat X,B,A,T,T2; 30456288a8SStefano Zampini PetscInt m,n,k,M = 10,N = 10,K = 5, ldx = 3, ldb = 5; 316280154eSStefano Zampini const char *deft = MATAIJ; 326280154eSStefano Zampini char mattype[256]; 336280154eSStefano Zampini PetscBool flg,symm = PETSC_FALSE,testtt = PETSC_TRUE, testnest = PETSC_TRUE, testtranspose = PETSC_TRUE, testcircular = PETSC_FALSE, local = PETSC_TRUE; 34456288a8SStefano Zampini PetscBool xgpu = PETSC_FALSE,bgpu = PETSC_FALSE; 356280154eSStefano Zampini PetscScalar *dataX = NULL,*dataB = NULL; 366280154eSStefano Zampini PetscScalar *aX,*aB; 37456288a8SStefano Zampini PetscReal err; 386280154eSStefano Zampini PetscErrorCode ierr; 396280154eSStefano Zampini 406280154eSStefano Zampini ierr = PetscInitialize(&argc,&args,NULL,help);if (ierr) return ierr; 416280154eSStefano Zampini ierr = PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL);CHKERRQ(ierr); 426280154eSStefano Zampini ierr = PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL);CHKERRQ(ierr); 436280154eSStefano Zampini ierr = PetscOptionsGetInt(NULL,NULL,"-K",&K,NULL);CHKERRQ(ierr); 446280154eSStefano Zampini ierr = PetscOptionsGetBool(NULL,NULL,"-symm",&symm,NULL);CHKERRQ(ierr); 456280154eSStefano Zampini ierr = PetscOptionsGetBool(NULL,NULL,"-local",&local,NULL);CHKERRQ(ierr); 466280154eSStefano Zampini ierr = PetscOptionsGetInt(NULL,NULL,"-ldx",&ldx,NULL);CHKERRQ(ierr); 476280154eSStefano Zampini ierr = PetscOptionsGetInt(NULL,NULL,"-ldb",&ldb,NULL);CHKERRQ(ierr); 486280154eSStefano Zampini ierr = PetscOptionsGetBool(NULL,NULL,"-testtranspose",&testtranspose,NULL);CHKERRQ(ierr); 496280154eSStefano Zampini ierr = PetscOptionsGetBool(NULL,NULL,"-testnest",&testnest,NULL);CHKERRQ(ierr); 506280154eSStefano Zampini ierr = PetscOptionsGetBool(NULL,NULL,"-testtt",&testtt,NULL);CHKERRQ(ierr); 516280154eSStefano Zampini ierr = PetscOptionsGetBool(NULL,NULL,"-testcircular",&testcircular,NULL);CHKERRQ(ierr); 52456288a8SStefano Zampini ierr = PetscOptionsGetBool(NULL,NULL,"-xgpu",&xgpu,NULL);CHKERRQ(ierr); 53456288a8SStefano Zampini ierr = PetscOptionsGetBool(NULL,NULL,"-bgpu",&bgpu,NULL);CHKERRQ(ierr); 546280154eSStefano Zampini #if defined(PETSC_USE_COMPLEX) 556280154eSStefano Zampini testtranspose = PETSC_FALSE; 566280154eSStefano Zampini testtt = PETSC_FALSE; 576280154eSStefano Zampini #endif 586280154eSStefano Zampini ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 596280154eSStefano Zampini ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr); 606280154eSStefano Zampini ierr = MatSetType(A,MATAIJ);CHKERRQ(ierr); 616280154eSStefano Zampini ierr = MatSetUp(A);CHKERRQ(ierr); 626280154eSStefano Zampini ierr = MatSetRandom(A,NULL);CHKERRQ(ierr); 636280154eSStefano Zampini if (M==N && symm) { 646280154eSStefano Zampini Mat AT; 656280154eSStefano Zampini 666280154eSStefano Zampini ierr = MatHermitianTranspose(A,MAT_INITIAL_MATRIX,&AT);CHKERRQ(ierr); 676280154eSStefano Zampini ierr = MatAXPY(A,1.0,AT,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 686280154eSStefano Zampini ierr = MatDestroy(&AT);CHKERRQ(ierr); 696280154eSStefano Zampini ierr = MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);CHKERRQ(ierr); 706280154eSStefano Zampini } 716280154eSStefano Zampini ierr = MatViewFromOptions(A,NULL,"-A_init_view");CHKERRQ(ierr); 726280154eSStefano Zampini ierr = PetscOptionsBegin(PETSC_COMM_WORLD,"","","");CHKERRQ(ierr); 73456288a8SStefano Zampini ierr = PetscOptionsFList("-A_mat_type","Matrix type","MatSetType",MatList,deft,mattype,256,&flg);CHKERRQ(ierr); 746280154eSStefano Zampini ierr = PetscOptionsEnd();CHKERRQ(ierr); 756280154eSStefano Zampini if (flg) { 766280154eSStefano Zampini Mat A2; 776280154eSStefano Zampini 78456288a8SStefano Zampini /* MATSEQAIJCUSPARSE does not support MAT_INITIAL_MATRIX */ 79456288a8SStefano Zampini ierr = MatDuplicate(A,MAT_COPY_VALUES,&A2);CHKERRQ(ierr); 80456288a8SStefano Zampini ierr = MatConvert(A,mattype,MAT_INPLACE_MATRIX,&A);CHKERRQ(ierr); 816280154eSStefano Zampini ierr = MatMultEqual(A,A2,10,&flg);CHKERRQ(ierr); 826280154eSStefano Zampini if (!flg) { 836280154eSStefano Zampini Mat AE,A2E; 846280154eSStefano Zampini 856280154eSStefano Zampini ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with convert\n");CHKERRQ(ierr); 866280154eSStefano Zampini ierr = MatComputeOperator(A,MATDENSE,&AE);CHKERRQ(ierr); 876280154eSStefano Zampini ierr = MatComputeOperator(A2,MATDENSE,&A2E);CHKERRQ(ierr); 886280154eSStefano Zampini ierr = MatView(AE,NULL);CHKERRQ(ierr); 896280154eSStefano Zampini ierr = MatView(A2E,NULL);CHKERRQ(ierr); 906280154eSStefano Zampini ierr = MatAXPY(A2E,-1.0,A,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 916280154eSStefano Zampini ierr = MatView(A2E,NULL);CHKERRQ(ierr); 926280154eSStefano Zampini ierr = MatDestroy(&A2E);CHKERRQ(ierr); 936280154eSStefano Zampini ierr = MatDestroy(&AE);CHKERRQ(ierr); 946280154eSStefano Zampini } 95456288a8SStefano Zampini ierr = MatDestroy(&A2);CHKERRQ(ierr); 966280154eSStefano Zampini } 976280154eSStefano Zampini ierr = MatViewFromOptions(A,NULL,"-A_view");CHKERRQ(ierr); 986280154eSStefano Zampini 996280154eSStefano Zampini ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 1006280154eSStefano Zampini if (local) { 1016280154eSStefano Zampini ierr = PetscMalloc1(PetscMax((m+ldx)*K,1),&dataX);CHKERRQ(ierr); 1026280154eSStefano Zampini ierr = PetscMalloc1(PetscMax((n+ldb)*K,1),&dataB);CHKERRQ(ierr); 1036280154eSStefano Zampini } 1046280154eSStefano Zampini ierr = MatCreateDense(PETSC_COMM_WORLD,n,PETSC_DECIDE,N,K,dataB,&B);CHKERRQ(ierr); 1056280154eSStefano Zampini ierr = MatCreateDense(PETSC_COMM_WORLD,m,PETSC_DECIDE,M,K,dataX,&X);CHKERRQ(ierr); 1066280154eSStefano Zampini 1076280154eSStefano Zampini /* store pointer to dense data for testing */ 1086280154eSStefano Zampini ierr = MatDenseGetArrayRead(B,(const PetscScalar**)&dataB);CHKERRQ(ierr); 1096280154eSStefano Zampini ierr = MatDenseGetArrayRead(X,(const PetscScalar**)&dataX);CHKERRQ(ierr); 1106280154eSStefano Zampini aX = dataX; 1116280154eSStefano Zampini aB = dataB; 1126280154eSStefano Zampini ierr = MatDenseRestoreArrayRead(B,(const PetscScalar**)&dataB);CHKERRQ(ierr); 1136280154eSStefano Zampini ierr = MatDenseRestoreArrayRead(X,(const PetscScalar**)&dataX);CHKERRQ(ierr); 1146280154eSStefano Zampini if (local) { 1156280154eSStefano Zampini dataX = aX; 1166280154eSStefano Zampini dataB = aB; 1176280154eSStefano Zampini } 1186280154eSStefano Zampini ierr = MatSetRandom(B,NULL);CHKERRQ(ierr); 119456288a8SStefano Zampini /* convert to CUDA if needed */ 120456288a8SStefano Zampini if (bgpu) { 121456288a8SStefano Zampini ierr = MatConvert(B,MATDENSECUDA,MAT_INPLACE_MATRIX,&B);CHKERRQ(ierr); 122456288a8SStefano Zampini } 123456288a8SStefano Zampini if (xgpu) { 124456288a8SStefano Zampini ierr = MatAssemblyBegin(X,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 125456288a8SStefano Zampini ierr = MatAssemblyEnd(X,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 126456288a8SStefano Zampini ierr = MatConvert(X,MATDENSECUDA,MAT_INPLACE_MATRIX,&X);CHKERRQ(ierr); 127456288a8SStefano Zampini } 128456288a8SStefano Zampini ierr = CheckLocal(B,X,aB,aX);CHKERRQ(ierr); 1296280154eSStefano Zampini 1306280154eSStefano Zampini /* Test reusing a previously allocated dense buffer */ 1316280154eSStefano Zampini ierr = MatMatMult(A,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&X);CHKERRQ(ierr); 1326280154eSStefano Zampini ierr = CheckLocal(B,X,aB,aX);CHKERRQ(ierr); 1336280154eSStefano Zampini ierr = MatMatMultEqual(A,B,X,10,&flg);CHKERRQ(ierr); 1346280154eSStefano Zampini if (!flg) { 1356280154eSStefano Zampini ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with reusage\n");CHKERRQ(ierr); 1366280154eSStefano Zampini ierr = MatMatMult(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T);CHKERRQ(ierr); 1376280154eSStefano Zampini ierr = MatAXPY(T,-1.0,X,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 1386280154eSStefano Zampini ierr = MatView(T,NULL);CHKERRQ(ierr); 1396280154eSStefano Zampini ierr = MatDestroy(&T);CHKERRQ(ierr); 1406280154eSStefano Zampini } 1416280154eSStefano Zampini 142456288a8SStefano Zampini /* Test MatDenseGetColumnVec and friends */ 143456288a8SStefano Zampini ierr = MatMatMult(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T);CHKERRQ(ierr); 144456288a8SStefano Zampini ierr = MatDuplicate(T,MAT_DO_NOT_COPY_VALUES,&T2);CHKERRQ(ierr); 145456288a8SStefano Zampini for (k=0;k<K;k++) { 146456288a8SStefano Zampini Vec Xv,Tv,T2v; 147456288a8SStefano Zampini 148456288a8SStefano Zampini ierr = MatDenseGetColumnVecRead(X,k,&Xv);CHKERRQ(ierr); 149456288a8SStefano Zampini ierr = MatDenseGetColumnVec(T,k,&Tv);CHKERRQ(ierr); 150456288a8SStefano Zampini ierr = MatDenseGetColumnVecWrite(T2,k,&T2v);CHKERRQ(ierr); 151456288a8SStefano Zampini ierr = VecCopy(Xv,T2v);CHKERRQ(ierr); 152456288a8SStefano Zampini ierr = VecAXPY(Tv,-1.,Xv);CHKERRQ(ierr); 153456288a8SStefano Zampini ierr = MatDenseRestoreColumnVecRead(X,k,&Xv);CHKERRQ(ierr); 154456288a8SStefano Zampini ierr = MatDenseRestoreColumnVec(T,k,&Tv);CHKERRQ(ierr); 155456288a8SStefano Zampini ierr = MatDenseRestoreColumnVecWrite(T2,k,&T2v);CHKERRQ(ierr); 156456288a8SStefano Zampini } 157456288a8SStefano Zampini ierr = MatNorm(T,NORM_FROBENIUS,&err);CHKERRQ(ierr); 158456288a8SStefano Zampini if (err > PETSC_SMALL) { 159456288a8SStefano Zampini ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with MatDenseGetColumnVec\n");CHKERRQ(ierr); 160456288a8SStefano Zampini ierr = MatView(T,NULL);CHKERRQ(ierr); 161456288a8SStefano Zampini } 162456288a8SStefano Zampini ierr = MatAXPY(T2,-1.,X,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 163456288a8SStefano Zampini ierr = MatNorm(T2,NORM_FROBENIUS,&err);CHKERRQ(ierr); 164456288a8SStefano Zampini if (err > PETSC_SMALL) { 165456288a8SStefano Zampini ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with MatDenseGetColumnVecWrite\n");CHKERRQ(ierr); 166456288a8SStefano Zampini ierr = MatView(T2,NULL);CHKERRQ(ierr); 167456288a8SStefano Zampini } 168456288a8SStefano Zampini ierr = MatDestroy(&T);CHKERRQ(ierr); 169456288a8SStefano Zampini ierr = MatDestroy(&T2);CHKERRQ(ierr); 170456288a8SStefano Zampini 171456288a8SStefano Zampini /* Test with MatShell */ 172456288a8SStefano Zampini ierr = MatConvert(A,MATSHELL,MAT_INITIAL_MATRIX,&T2);CHKERRQ(ierr); 173456288a8SStefano Zampini ierr = MatMatMult(T2,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&X);CHKERRQ(ierr); 174456288a8SStefano Zampini ierr = CheckLocal(B,X,aB,aX);CHKERRQ(ierr); 175456288a8SStefano Zampini ierr = MatMatMultEqual(T2,B,X,10,&flg);CHKERRQ(ierr); 176456288a8SStefano Zampini if (!flg) { 177456288a8SStefano Zampini ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with reusage (MATSHELL)\n");CHKERRQ(ierr); 178456288a8SStefano Zampini ierr = MatMatMult(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T);CHKERRQ(ierr); 179456288a8SStefano Zampini ierr = MatAXPY(T,-1.0,X,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 180456288a8SStefano Zampini ierr = MatView(T,NULL);CHKERRQ(ierr); 181456288a8SStefano Zampini ierr = MatDestroy(&T);CHKERRQ(ierr); 182456288a8SStefano Zampini } 183*6718818eSStefano Zampini 184*6718818eSStefano Zampini if (testtranspose) { 185456288a8SStefano Zampini ierr = MatTransposeMatMult(T2,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&B);CHKERRQ(ierr); 186456288a8SStefano Zampini ierr = CheckLocal(B,X,aB,aX);CHKERRQ(ierr); 187456288a8SStefano Zampini ierr = MatTransposeMatMultEqual(T2,X,B,10,&flg);CHKERRQ(ierr); 188456288a8SStefano Zampini if (!flg) { 189456288a8SStefano Zampini ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with reusage (MatTranspose, MATSHELL)\n");CHKERRQ(ierr); 190*6718818eSStefano Zampini ierr = MatTransposeMatMult(A,X,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T);CHKERRQ(ierr); 191*6718818eSStefano Zampini ierr = MatAXPY(T,-1.0,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 192*6718818eSStefano Zampini ierr = MatView(T,NULL);CHKERRQ(ierr); 193*6718818eSStefano Zampini ierr = MatDestroy(&T);CHKERRQ(ierr); 194*6718818eSStefano Zampini } 195456288a8SStefano Zampini } 196456288a8SStefano Zampini ierr = MatDestroy(&T2);CHKERRQ(ierr); 197456288a8SStefano Zampini 1986280154eSStefano Zampini if (testnest) { /* test with MatNest */ 1996280154eSStefano Zampini Mat NA; 200456288a8SStefano Zampini const char *vtype; 2016280154eSStefano Zampini 2026280154eSStefano Zampini ierr = MatCreateNest(PETSC_COMM_WORLD,1,NULL,1,NULL,&A,&NA);CHKERRQ(ierr); 203456288a8SStefano Zampini /* needed to test against CUSPARSE matrices */ 204456288a8SStefano Zampini ierr = MatGetVecType(A,&vtype);CHKERRQ(ierr); 205456288a8SStefano Zampini ierr = MatSetVecType(NA,vtype);CHKERRQ(ierr); 2066280154eSStefano Zampini ierr = MatViewFromOptions(NA,NULL,"-NA_view");CHKERRQ(ierr); 2076280154eSStefano Zampini ierr = MatMatMult(NA,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&X);CHKERRQ(ierr); 2086280154eSStefano Zampini ierr = CheckLocal(B,X,aB,aX);CHKERRQ(ierr); 2096280154eSStefano Zampini ierr = MatMatMultEqual(NA,B,X,10,&flg);CHKERRQ(ierr); 2106280154eSStefano Zampini if (!flg) { 2116280154eSStefano Zampini ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with Nest\n");CHKERRQ(ierr); 2126280154eSStefano Zampini ierr = MatMatMult(NA,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T);CHKERRQ(ierr); 2136280154eSStefano Zampini ierr = MatAXPY(T,-1.0,X,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 2146280154eSStefano Zampini ierr = MatView(T,NULL);CHKERRQ(ierr); 2156280154eSStefano Zampini ierr = MatDestroy(&T);CHKERRQ(ierr); 2166280154eSStefano Zampini } 2176280154eSStefano Zampini ierr = MatDestroy(&NA);CHKERRQ(ierr); 2186280154eSStefano Zampini } 2196280154eSStefano Zampini 2206280154eSStefano Zampini if (testtranspose) { /* test with Transpose */ 2216280154eSStefano Zampini Mat TA; 2226280154eSStefano Zampini 2236280154eSStefano Zampini ierr = MatCreateHermitianTranspose(A,&TA);CHKERRQ(ierr); 2246280154eSStefano Zampini ierr = MatMatMult(TA,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&B);CHKERRQ(ierr); 2256280154eSStefano Zampini ierr = CheckLocal(B,X,aB,aX);CHKERRQ(ierr); 2266280154eSStefano Zampini ierr = MatMatMultEqual(TA,X,B,10,&flg);CHKERRQ(ierr); 2276280154eSStefano Zampini if (!flg) { 2286280154eSStefano Zampini ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with Transpose\n");CHKERRQ(ierr); 2296280154eSStefano Zampini ierr = MatMatMult(TA,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T);CHKERRQ(ierr); 2306280154eSStefano Zampini ierr = MatAXPY(T,-1.0,X,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 2316280154eSStefano Zampini ierr = MatView(T,NULL);CHKERRQ(ierr); 2326280154eSStefano Zampini ierr = MatDestroy(&T);CHKERRQ(ierr); 2336280154eSStefano Zampini } 2346280154eSStefano Zampini ierr = MatDestroy(&TA);CHKERRQ(ierr); 2356280154eSStefano Zampini } 2366280154eSStefano Zampini 2376280154eSStefano Zampini if (testtt) { /* test with Transpose(Transpose) */ 2386280154eSStefano Zampini Mat TA, TTA; 2396280154eSStefano Zampini 2406280154eSStefano Zampini ierr = MatCreateHermitianTranspose(A,&TA);CHKERRQ(ierr); 2416280154eSStefano Zampini ierr = MatCreateHermitianTranspose(TA,&TTA);CHKERRQ(ierr); 2426280154eSStefano Zampini ierr = MatMatMult(TTA,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&X);CHKERRQ(ierr); 2436280154eSStefano Zampini ierr = CheckLocal(B,X,aB,aX);CHKERRQ(ierr); 2446280154eSStefano Zampini ierr = MatMatMultEqual(TTA,B,X,10,&flg);CHKERRQ(ierr); 2456280154eSStefano Zampini if (!flg) { 2466280154eSStefano Zampini ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with Transpose(Transpose)\n");CHKERRQ(ierr); 2476280154eSStefano Zampini ierr = MatMatMult(TTA,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T);CHKERRQ(ierr); 2486280154eSStefano Zampini ierr = MatAXPY(T,-1.0,X,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 2496280154eSStefano Zampini ierr = MatView(T,NULL);CHKERRQ(ierr); 2506280154eSStefano Zampini ierr = MatDestroy(&T);CHKERRQ(ierr); 2516280154eSStefano Zampini } 2526280154eSStefano Zampini ierr = MatDestroy(&TA);CHKERRQ(ierr); 2536280154eSStefano Zampini ierr = MatDestroy(&TTA);CHKERRQ(ierr); 2546280154eSStefano Zampini } 2556280154eSStefano Zampini 2566280154eSStefano Zampini if (testcircular) { /* test circular */ 2576280154eSStefano Zampini Mat AB; 2586280154eSStefano Zampini 2596280154eSStefano Zampini ierr = MatMatMult(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&AB);CHKERRQ(ierr); 2606280154eSStefano Zampini ierr = MatMatMult(A,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&X);CHKERRQ(ierr); 2616280154eSStefano Zampini ierr = CheckLocal(B,X,aB,aX);CHKERRQ(ierr); 2626280154eSStefano Zampini if (M == N && N == K) { 2636280154eSStefano Zampini ierr = MatMatMult(A,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&B);CHKERRQ(ierr); 2646280154eSStefano Zampini } else { 2656280154eSStefano Zampini ierr = MatTransposeMatMult(A,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&B);CHKERRQ(ierr); 2666280154eSStefano Zampini } 2676280154eSStefano Zampini ierr = CheckLocal(B,X,aB,aX);CHKERRQ(ierr); 2686280154eSStefano Zampini ierr = MatDestroy(&AB);CHKERRQ(ierr); 2696280154eSStefano Zampini } 2706280154eSStefano Zampini ierr = MatDestroy(&X);CHKERRQ(ierr); 2716280154eSStefano Zampini ierr = MatDestroy(&B);CHKERRQ(ierr); 2726280154eSStefano Zampini ierr = MatDestroy(&A);CHKERRQ(ierr); 2736280154eSStefano Zampini ierr = PetscFree(dataX);CHKERRQ(ierr); 2746280154eSStefano Zampini ierr = PetscFree(dataB);CHKERRQ(ierr); 2756280154eSStefano Zampini ierr = PetscFinalize(); 2766280154eSStefano Zampini return ierr; 2776280154eSStefano Zampini } 2786280154eSStefano Zampini 2796280154eSStefano Zampini /*TEST 2806280154eSStefano Zampini 2816280154eSStefano Zampini test: 2826280154eSStefano Zampini suffix: 1 2836280154eSStefano Zampini args: -local {{0 1}} 2846280154eSStefano Zampini 2856280154eSStefano Zampini test: 2866280154eSStefano Zampini output_file: output/ex70_1.out 287456288a8SStefano Zampini requires: cuda 288456288a8SStefano Zampini suffix: 1_cuda 289456288a8SStefano Zampini args: -local {{0 1}} -xgpu {{0 1}} -bgpu {{0 1}} -A_mat_type {{seqaijcusparse seqaij}} -testnest 0 290456288a8SStefano Zampini 291456288a8SStefano Zampini test: 292456288a8SStefano Zampini TODO: VecGetSubVector seems broken with CUDA 293456288a8SStefano Zampini output_file: output/ex70_1.out 294456288a8SStefano Zampini requires: cuda 295456288a8SStefano Zampini suffix: 1_cuda_broken 296456288a8SStefano Zampini args: -local {{0 1}} -xgpu {{0 1}} -bgpu {{0 1}} -A_mat_type seqaijcusparse -testnest 297456288a8SStefano Zampini 298456288a8SStefano Zampini test: 299456288a8SStefano Zampini output_file: output/ex70_1.out 3006280154eSStefano Zampini nsize: 2 3016280154eSStefano Zampini suffix: 1_par 302*6718818eSStefano Zampini args: -testtranspose -local {{0 1}} 3036280154eSStefano Zampini 3046280154eSStefano Zampini test: 305456288a8SStefano Zampini output_file: output/ex70_1.out 306456288a8SStefano Zampini requires: cuda 307456288a8SStefano Zampini nsize: 2 308456288a8SStefano Zampini suffix: 1_par_cuda 309456288a8SStefano Zampini args: -testtranspose 0 -local {{0 1}} -xgpu {{0 1}} -bgpu {{0 1}} -A_mat_type {{mpiaijcusparse mpiaij}} -testnest 0 310456288a8SStefano Zampini 311456288a8SStefano Zampini test: 3126280154eSStefano Zampini output_file: output/ex70_1.out 3136280154eSStefano Zampini suffix: 2 3146280154eSStefano Zampini nsize: 1 3156280154eSStefano Zampini args: -M {{7 11}} -N {{12 9}} -K {{1 3}} -local {{0 1}} 3166280154eSStefano Zampini 3176280154eSStefano Zampini test: 318456288a8SStefano Zampini requires: cuda 319456288a8SStefano Zampini output_file: output/ex70_1.out 320456288a8SStefano Zampini suffix: 2_cuda 321456288a8SStefano Zampini nsize: 1 322456288a8SStefano Zampini args: -M 7 -N 9 -K 2 -local {{0 1}} -testnest 0 -A_mat_type {{seqdensecuda seqdense}} -xgpu {{0 1}} -bgpu {{0 1}} 323456288a8SStefano Zampini 324456288a8SStefano Zampini test: 3256280154eSStefano Zampini output_file: output/ex70_1.out 3266280154eSStefano Zampini suffix: 2_par 3276280154eSStefano Zampini nsize: 2 328*6718818eSStefano Zampini args: -M {{7 11}} -N {{12 9}} -K {{1 3}} -local {{0 1}} -testcircular 3296280154eSStefano Zampini 3306280154eSStefano Zampini test: 331456288a8SStefano Zampini requires: cuda 332456288a8SStefano Zampini output_file: output/ex70_1.out 333456288a8SStefano Zampini suffix: 2_par_cuda 334456288a8SStefano Zampini nsize: 2 335456288a8SStefano Zampini args: -M 11 -N 9 -K 1 -local {{0 1}} -testcircular 0 -A_mat_type mpiaijcusparse -xgpu -bgpu -testnest 0 336456288a8SStefano Zampini 337456288a8SStefano Zampini test: 3386280154eSStefano Zampini output_file: output/ex70_1.out 3396280154eSStefano Zampini suffix: 3 340456288a8SStefano Zampini nsize: {{1 3}} 341456288a8SStefano Zampini args: -M 13 -N 13 -K {{1 3}} -local {{0 1}} -A_mat_type sbaij -symm -testtranspose 0 3426280154eSStefano Zampini 3436280154eSStefano Zampini test: 3446280154eSStefano Zampini output_file: output/ex70_1.out 3456280154eSStefano Zampini suffix: 4 3466280154eSStefano Zampini nsize: 1 3476280154eSStefano Zampini args: -M 3 -N 3 -K 3 -local {{0 1}} -testcircular 3486280154eSStefano Zampini 3496280154eSStefano Zampini test: 3506280154eSStefano Zampini output_file: output/ex70_1.out 3516280154eSStefano Zampini suffix: 5 3526280154eSStefano Zampini nsize: {{2 4}} 3536280154eSStefano Zampini args: -M 3 -N 3 -K 3 -local 1 -testcircular -testtranspose 0 3546280154eSStefano Zampini 3556280154eSStefano Zampini test: 3566280154eSStefano Zampini TODO: MatCreateDense broken with some processors not having local rows 3576280154eSStefano Zampini output_file: output/ex70_1.out 3586280154eSStefano Zampini suffix: 5_broken 3596280154eSStefano Zampini nsize: {{2 4}} 3606280154eSStefano Zampini args: -M 3 -N 3 -K 3 -local 0 -testcircular -testtranspose 0 3616280154eSStefano Zampini 3626280154eSStefano Zampini test: 3636280154eSStefano Zampini output_file: output/ex70_1.out 3646280154eSStefano Zampini suffix: 6 3656280154eSStefano Zampini nsize: 1 3666280154eSStefano Zampini args: -M {{1 3}} -N {{2 5}} -K {{1 2}} -local {{0 1}} -testcircular -testtranspose 3676280154eSStefano Zampini 3686280154eSStefano Zampini test: 3696280154eSStefano Zampini output_file: output/ex70_1.out 3706280154eSStefano Zampini suffix: 7 3716280154eSStefano Zampini nsize: 1 372456288a8SStefano Zampini args: -M 13 -N 13 -K {{1 3}} -local {{0 1}} -A_mat_type dense -testnest 0 -testcircular 3736280154eSStefano Zampini 3746280154eSStefano Zampini test: 3756280154eSStefano Zampini TODO: NEST x DENSE with dense nested matrices seems broken in this case 3766280154eSStefano Zampini output_file: output/ex70_1.out 3776280154eSStefano Zampini suffix: 7_broken 3786280154eSStefano Zampini nsize: 1 379456288a8SStefano Zampini args: -M 13 -N 13 -K {{1 3}} -local {{0 1}} -A_mat_type dense -testnest -testcircular 380*6718818eSStefano Zampini 3816280154eSStefano Zampini TEST*/ 382