1*6280154eSStefano Zampini #include <petscmat.h> 2*6280154eSStefano Zampini 3*6280154eSStefano Zampini static char help[] = "Tests MatMatMult with MAT_REUSE_MATRIX and already allocated dense result.\n\n"; 4*6280154eSStefano Zampini 5*6280154eSStefano Zampini static PetscErrorCode CheckLocal(Mat A, Mat B, PetscScalar *a, PetscScalar *b) 6*6280154eSStefano Zampini { 7*6280154eSStefano Zampini PetscErrorCode ierr; 8*6280154eSStefano Zampini PetscBool wA = PETSC_FALSE, wB = PETSC_FALSE; 9*6280154eSStefano Zampini 10*6280154eSStefano Zampini PetscFunctionBegin; 11*6280154eSStefano Zampini if (a) { 12*6280154eSStefano Zampini const PetscScalar *Aa; 13*6280154eSStefano Zampini ierr = MatDenseGetArrayRead(A,&Aa);CHKERRQ(ierr); 14*6280154eSStefano Zampini wA = (PetscBool)(a != Aa); 15*6280154eSStefano Zampini ierr = MatDenseRestoreArrayRead(A,&Aa);CHKERRQ(ierr); 16*6280154eSStefano Zampini } 17*6280154eSStefano Zampini if (b) { 18*6280154eSStefano Zampini const PetscScalar *Bb; 19*6280154eSStefano Zampini ierr = MatDenseGetArrayRead(B,&Bb);CHKERRQ(ierr); 20*6280154eSStefano Zampini wB = (PetscBool)(b != Bb); 21*6280154eSStefano Zampini ierr = MatDenseRestoreArrayRead(B,&Bb);CHKERRQ(ierr); 22*6280154eSStefano Zampini } 23*6280154eSStefano Zampini if (wA || wB) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong array in A Mat %d, Wrong array in B Mat %d",wA,wB); 24*6280154eSStefano Zampini PetscFunctionReturn(0); 25*6280154eSStefano Zampini } 26*6280154eSStefano Zampini 27*6280154eSStefano Zampini int main(int argc,char **args) 28*6280154eSStefano Zampini { 29*6280154eSStefano Zampini Mat X,B,A,T; 30*6280154eSStefano Zampini PetscInt m,n,M = 10,N = 10,K = 5, ldx = 0, ldb = 0; 31*6280154eSStefano Zampini const char *deft = MATAIJ; 32*6280154eSStefano Zampini char mattype[256]; 33*6280154eSStefano Zampini PetscBool flg,symm = PETSC_FALSE,testtt = PETSC_TRUE, testnest = PETSC_TRUE, testtranspose = PETSC_TRUE, testcircular = PETSC_FALSE, local = PETSC_TRUE; 34*6280154eSStefano Zampini PetscScalar *dataX = NULL,*dataB = NULL; 35*6280154eSStefano Zampini PetscScalar *aX,*aB; 36*6280154eSStefano Zampini PetscErrorCode ierr; 37*6280154eSStefano Zampini 38*6280154eSStefano Zampini ierr = PetscInitialize(&argc,&args,NULL,help);if (ierr) return ierr; 39*6280154eSStefano Zampini ierr = PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL);CHKERRQ(ierr); 40*6280154eSStefano Zampini ierr = PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL);CHKERRQ(ierr); 41*6280154eSStefano Zampini ierr = PetscOptionsGetInt(NULL,NULL,"-K",&K,NULL);CHKERRQ(ierr); 42*6280154eSStefano Zampini ierr = PetscOptionsGetBool(NULL,NULL,"-symm",&symm,NULL);CHKERRQ(ierr); 43*6280154eSStefano Zampini ierr = PetscOptionsGetBool(NULL,NULL,"-local",&local,NULL);CHKERRQ(ierr); 44*6280154eSStefano Zampini ierr = PetscOptionsGetInt(NULL,NULL,"-ldx",&ldx,NULL);CHKERRQ(ierr); 45*6280154eSStefano Zampini ierr = PetscOptionsGetInt(NULL,NULL,"-ldb",&ldb,NULL);CHKERRQ(ierr); 46*6280154eSStefano Zampini ierr = PetscOptionsGetBool(NULL,NULL,"-testtranspose",&testtranspose,NULL);CHKERRQ(ierr); 47*6280154eSStefano Zampini ierr = PetscOptionsGetBool(NULL,NULL,"-testnest",&testnest,NULL);CHKERRQ(ierr); 48*6280154eSStefano Zampini ierr = PetscOptionsGetBool(NULL,NULL,"-testtt",&testtt,NULL);CHKERRQ(ierr); 49*6280154eSStefano Zampini ierr = PetscOptionsGetBool(NULL,NULL,"-testcircular",&testcircular,NULL);CHKERRQ(ierr); 50*6280154eSStefano Zampini #if defined(PETSC_USE_COMPLEX) 51*6280154eSStefano Zampini testtranspose = PETSC_FALSE; 52*6280154eSStefano Zampini testtt = PETSC_FALSE; 53*6280154eSStefano Zampini #endif 54*6280154eSStefano Zampini ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 55*6280154eSStefano Zampini ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr); 56*6280154eSStefano Zampini ierr = MatSetType(A,MATAIJ);CHKERRQ(ierr); 57*6280154eSStefano Zampini ierr = MatSetUp(A);CHKERRQ(ierr); 58*6280154eSStefano Zampini ierr = MatSetRandom(A,NULL);CHKERRQ(ierr); 59*6280154eSStefano Zampini if (M==N && symm) { 60*6280154eSStefano Zampini Mat AT; 61*6280154eSStefano Zampini 62*6280154eSStefano Zampini ierr = MatHermitianTranspose(A,MAT_INITIAL_MATRIX,&AT);CHKERRQ(ierr); 63*6280154eSStefano Zampini ierr = MatAXPY(A,1.0,AT,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 64*6280154eSStefano Zampini ierr = MatDestroy(&AT);CHKERRQ(ierr); 65*6280154eSStefano Zampini ierr = MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);CHKERRQ(ierr); 66*6280154eSStefano Zampini } 67*6280154eSStefano Zampini ierr = MatViewFromOptions(A,NULL,"-A_init_view");CHKERRQ(ierr); 68*6280154eSStefano Zampini ierr = PetscOptionsBegin(PETSC_COMM_WORLD,"","","");CHKERRQ(ierr); 69*6280154eSStefano Zampini ierr = PetscOptionsFList("-mat_type","Matrix type","MatSetType",MatList,deft,mattype,256,&flg);CHKERRQ(ierr); 70*6280154eSStefano Zampini ierr = PetscOptionsEnd();CHKERRQ(ierr); 71*6280154eSStefano Zampini if (flg) { 72*6280154eSStefano Zampini Mat A2; 73*6280154eSStefano Zampini 74*6280154eSStefano Zampini ierr = MatConvert(A,mattype,MAT_INITIAL_MATRIX,&A2);CHKERRQ(ierr); 75*6280154eSStefano Zampini ierr = MatMultEqual(A,A2,10,&flg);CHKERRQ(ierr); 76*6280154eSStefano Zampini if (!flg) { 77*6280154eSStefano Zampini Mat AE,A2E; 78*6280154eSStefano Zampini 79*6280154eSStefano Zampini ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with convert\n");CHKERRQ(ierr); 80*6280154eSStefano Zampini ierr = MatComputeOperator(A,MATDENSE,&AE);CHKERRQ(ierr); 81*6280154eSStefano Zampini ierr = MatComputeOperator(A2,MATDENSE,&A2E);CHKERRQ(ierr); 82*6280154eSStefano Zampini ierr = MatView(AE,NULL);CHKERRQ(ierr); 83*6280154eSStefano Zampini ierr = MatView(A2E,NULL);CHKERRQ(ierr); 84*6280154eSStefano Zampini ierr = MatAXPY(A2E,-1.0,A,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 85*6280154eSStefano Zampini ierr = MatView(A2E,NULL);CHKERRQ(ierr); 86*6280154eSStefano Zampini ierr = MatDestroy(&A2E);CHKERRQ(ierr); 87*6280154eSStefano Zampini ierr = MatDestroy(&AE);CHKERRQ(ierr); 88*6280154eSStefano Zampini } 89*6280154eSStefano Zampini ierr = MatDestroy(&A);CHKERRQ(ierr); 90*6280154eSStefano Zampini A = A2; 91*6280154eSStefano Zampini } 92*6280154eSStefano Zampini ierr = MatViewFromOptions(A,NULL,"-A_view");CHKERRQ(ierr); 93*6280154eSStefano Zampini 94*6280154eSStefano Zampini ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 95*6280154eSStefano Zampini if (local) { 96*6280154eSStefano Zampini ierr = PetscMalloc1(PetscMax((m+ldx)*K,1),&dataX);CHKERRQ(ierr); 97*6280154eSStefano Zampini ierr = PetscMalloc1(PetscMax((n+ldb)*K,1),&dataB);CHKERRQ(ierr); 98*6280154eSStefano Zampini } 99*6280154eSStefano Zampini ierr = MatCreateDense(PETSC_COMM_WORLD,n,PETSC_DECIDE,N,K,dataB,&B);CHKERRQ(ierr); 100*6280154eSStefano Zampini ierr = MatCreateDense(PETSC_COMM_WORLD,m,PETSC_DECIDE,M,K,dataX,&X);CHKERRQ(ierr); 101*6280154eSStefano Zampini 102*6280154eSStefano Zampini /* store pointer to dense data for testing */ 103*6280154eSStefano Zampini ierr = MatDenseGetArrayRead(B,(const PetscScalar**)&dataB);CHKERRQ(ierr); 104*6280154eSStefano Zampini ierr = MatDenseGetArrayRead(X,(const PetscScalar**)&dataX);CHKERRQ(ierr); 105*6280154eSStefano Zampini aX = dataX; 106*6280154eSStefano Zampini aB = dataB; 107*6280154eSStefano Zampini ierr = MatDenseRestoreArrayRead(B,(const PetscScalar**)&dataB);CHKERRQ(ierr); 108*6280154eSStefano Zampini ierr = MatDenseRestoreArrayRead(X,(const PetscScalar**)&dataX);CHKERRQ(ierr); 109*6280154eSStefano Zampini if (local) { 110*6280154eSStefano Zampini dataX = aX; 111*6280154eSStefano Zampini dataB = aB; 112*6280154eSStefano Zampini } 113*6280154eSStefano Zampini ierr = MatSetRandom(B,NULL);CHKERRQ(ierr); 114*6280154eSStefano Zampini 115*6280154eSStefano Zampini /* Test reusing a previously allocated dense buffer */ 116*6280154eSStefano Zampini ierr = MatMatMult(A,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&X);CHKERRQ(ierr); 117*6280154eSStefano Zampini ierr = CheckLocal(B,X,aB,aX);CHKERRQ(ierr); 118*6280154eSStefano Zampini ierr = MatMatMultEqual(A,B,X,10,&flg);CHKERRQ(ierr); 119*6280154eSStefano Zampini if (!flg) { 120*6280154eSStefano Zampini ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with reusage\n");CHKERRQ(ierr); 121*6280154eSStefano Zampini ierr = MatMatMult(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T);CHKERRQ(ierr); 122*6280154eSStefano Zampini ierr = MatAXPY(T,-1.0,X,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 123*6280154eSStefano Zampini ierr = MatView(T,NULL);CHKERRQ(ierr); 124*6280154eSStefano Zampini ierr = MatDestroy(&T);CHKERRQ(ierr); 125*6280154eSStefano Zampini } 126*6280154eSStefano Zampini 127*6280154eSStefano Zampini if (testnest) { /* test with MatNest */ 128*6280154eSStefano Zampini Mat NA; 129*6280154eSStefano Zampini 130*6280154eSStefano Zampini ierr = MatCreateNest(PETSC_COMM_WORLD,1,NULL,1,NULL,&A,&NA);CHKERRQ(ierr); 131*6280154eSStefano Zampini ierr = MatViewFromOptions(NA,NULL,"-NA_view");CHKERRQ(ierr); 132*6280154eSStefano Zampini ierr = MatMatMult(NA,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&X);CHKERRQ(ierr); 133*6280154eSStefano Zampini ierr = CheckLocal(B,X,aB,aX);CHKERRQ(ierr); 134*6280154eSStefano Zampini ierr = MatMatMultEqual(NA,B,X,10,&flg);CHKERRQ(ierr); 135*6280154eSStefano Zampini if (!flg) { 136*6280154eSStefano Zampini ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with Nest\n");CHKERRQ(ierr); 137*6280154eSStefano Zampini ierr = MatMatMult(NA,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T);CHKERRQ(ierr); 138*6280154eSStefano Zampini ierr = MatAXPY(T,-1.0,X,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 139*6280154eSStefano Zampini ierr = MatView(T,NULL);CHKERRQ(ierr); 140*6280154eSStefano Zampini ierr = MatDestroy(&T);CHKERRQ(ierr); 141*6280154eSStefano Zampini } 142*6280154eSStefano Zampini ierr = MatDestroy(&NA);CHKERRQ(ierr); 143*6280154eSStefano Zampini } 144*6280154eSStefano Zampini 145*6280154eSStefano Zampini if (testtranspose) { /* test with Transpose */ 146*6280154eSStefano Zampini Mat TA; 147*6280154eSStefano Zampini 148*6280154eSStefano Zampini ierr = MatCreateHermitianTranspose(A,&TA);CHKERRQ(ierr); 149*6280154eSStefano Zampini ierr = MatMatMult(TA,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&B);CHKERRQ(ierr); 150*6280154eSStefano Zampini ierr = CheckLocal(B,X,aB,aX);CHKERRQ(ierr); 151*6280154eSStefano Zampini ierr = MatMatMultEqual(TA,X,B,10,&flg);CHKERRQ(ierr); 152*6280154eSStefano Zampini if (!flg) { 153*6280154eSStefano Zampini ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with Transpose\n");CHKERRQ(ierr); 154*6280154eSStefano Zampini ierr = MatMatMult(TA,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T);CHKERRQ(ierr); 155*6280154eSStefano Zampini ierr = MatAXPY(T,-1.0,X,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 156*6280154eSStefano Zampini ierr = MatView(T,NULL);CHKERRQ(ierr); 157*6280154eSStefano Zampini ierr = MatDestroy(&T);CHKERRQ(ierr); 158*6280154eSStefano Zampini } 159*6280154eSStefano Zampini ierr = MatDestroy(&TA);CHKERRQ(ierr); 160*6280154eSStefano Zampini } 161*6280154eSStefano Zampini 162*6280154eSStefano Zampini if (testtt) { /* test with Transpose(Transpose) */ 163*6280154eSStefano Zampini Mat TA, TTA; 164*6280154eSStefano Zampini 165*6280154eSStefano Zampini ierr = MatCreateHermitianTranspose(A,&TA);CHKERRQ(ierr); 166*6280154eSStefano Zampini ierr = MatCreateHermitianTranspose(TA,&TTA);CHKERRQ(ierr); 167*6280154eSStefano Zampini ierr = MatMatMult(TTA,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&X);CHKERRQ(ierr); 168*6280154eSStefano Zampini ierr = CheckLocal(B,X,aB,aX);CHKERRQ(ierr); 169*6280154eSStefano Zampini ierr = MatMatMultEqual(TTA,B,X,10,&flg);CHKERRQ(ierr); 170*6280154eSStefano Zampini if (!flg) { 171*6280154eSStefano Zampini ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with Transpose(Transpose)\n");CHKERRQ(ierr); 172*6280154eSStefano Zampini ierr = MatMatMult(TTA,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T);CHKERRQ(ierr); 173*6280154eSStefano Zampini ierr = MatAXPY(T,-1.0,X,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 174*6280154eSStefano Zampini ierr = MatView(T,NULL);CHKERRQ(ierr); 175*6280154eSStefano Zampini ierr = MatDestroy(&T);CHKERRQ(ierr); 176*6280154eSStefano Zampini } 177*6280154eSStefano Zampini ierr = MatDestroy(&TA);CHKERRQ(ierr); 178*6280154eSStefano Zampini ierr = MatDestroy(&TTA);CHKERRQ(ierr); 179*6280154eSStefano Zampini } 180*6280154eSStefano Zampini 181*6280154eSStefano Zampini if (testcircular) { /* test circular */ 182*6280154eSStefano Zampini Mat AB; 183*6280154eSStefano Zampini 184*6280154eSStefano Zampini ierr = MatMatMult(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&AB);CHKERRQ(ierr); 185*6280154eSStefano Zampini ierr = MatMatMult(A,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&X);CHKERRQ(ierr); 186*6280154eSStefano Zampini ierr = CheckLocal(B,X,aB,aX);CHKERRQ(ierr); 187*6280154eSStefano Zampini if (M == N && N == K) { 188*6280154eSStefano Zampini ierr = MatMatMult(A,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&B);CHKERRQ(ierr); 189*6280154eSStefano Zampini } else { 190*6280154eSStefano Zampini ierr = MatTransposeMatMult(A,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&B);CHKERRQ(ierr); 191*6280154eSStefano Zampini } 192*6280154eSStefano Zampini ierr = CheckLocal(B,X,aB,aX);CHKERRQ(ierr); 193*6280154eSStefano Zampini ierr = MatDestroy(&AB);CHKERRQ(ierr); 194*6280154eSStefano Zampini } 195*6280154eSStefano Zampini ierr = MatDestroy(&X);CHKERRQ(ierr); 196*6280154eSStefano Zampini ierr = MatDestroy(&B);CHKERRQ(ierr); 197*6280154eSStefano Zampini ierr = MatDestroy(&A);CHKERRQ(ierr); 198*6280154eSStefano Zampini ierr = PetscFree(dataX);CHKERRQ(ierr); 199*6280154eSStefano Zampini ierr = PetscFree(dataB);CHKERRQ(ierr); 200*6280154eSStefano Zampini ierr = PetscFinalize(); 201*6280154eSStefano Zampini return ierr; 202*6280154eSStefano Zampini } 203*6280154eSStefano Zampini 204*6280154eSStefano Zampini /*TEST 205*6280154eSStefano Zampini 206*6280154eSStefano Zampini test: 207*6280154eSStefano Zampini suffix: 1 208*6280154eSStefano Zampini args: -local {{0 1}} 209*6280154eSStefano Zampini 210*6280154eSStefano Zampini test: 211*6280154eSStefano Zampini output_file: output/ex70_1.out 212*6280154eSStefano Zampini nsize: 2 213*6280154eSStefano Zampini suffix: 1_par 214*6280154eSStefano Zampini args: -testtranspose 0 -local {{0 1}} 215*6280154eSStefano Zampini 216*6280154eSStefano Zampini test: 217*6280154eSStefano Zampini TODO: MPIAIJ x MPIDENSE broken for MatTransposeMatMult 218*6280154eSStefano Zampini output_file: output/ex70_1.out 219*6280154eSStefano Zampini nsize: 2 220*6280154eSStefano Zampini suffix: 1_par_broken 221*6280154eSStefano Zampini args: -testtranspose -local {{0 1}} 222*6280154eSStefano Zampini 223*6280154eSStefano Zampini test: 224*6280154eSStefano Zampini output_file: output/ex70_1.out 225*6280154eSStefano Zampini suffix: 2 226*6280154eSStefano Zampini nsize: 1 227*6280154eSStefano Zampini args: -M {{7 11}} -N {{12 9}} -K {{1 3}} -local {{0 1}} 228*6280154eSStefano Zampini 229*6280154eSStefano Zampini test: 230*6280154eSStefano Zampini output_file: output/ex70_1.out 231*6280154eSStefano Zampini suffix: 2_par 232*6280154eSStefano Zampini nsize: 2 233*6280154eSStefano Zampini args: -M {{7 11}} -N {{12 9}} -K {{1 3}} -local {{0 1}} -testcircular 0 234*6280154eSStefano Zampini 235*6280154eSStefano Zampini test: 236*6280154eSStefano Zampini TODO: MatTransposeMatMultSymbolic_SeqAIJ_SeqDense plays with the destroy routine 237*6280154eSStefano Zampini output_file: output/ex70_1.out 238*6280154eSStefano Zampini suffix: 2_broken 239*6280154eSStefano Zampini nsize: 2 240*6280154eSStefano Zampini args: -M {{7 11}} -N {{12 9}} -K {{1 3}} -local {{0 1}} -testcircular 1 241*6280154eSStefano Zampini 242*6280154eSStefano Zampini test: 243*6280154eSStefano Zampini output_file: output/ex70_1.out 244*6280154eSStefano Zampini suffix: 3 245*6280154eSStefano Zampini nsize: 1 246*6280154eSStefano Zampini args: -M 13 -N 13 -K {{1 3}} -local {{0 1}} -mat_type sbaij -symm -testtranspose 0 247*6280154eSStefano Zampini 248*6280154eSStefano Zampini test: 249*6280154eSStefano Zampini output_file: output/ex70_1.out 250*6280154eSStefano Zampini suffix: 4 251*6280154eSStefano Zampini nsize: 1 252*6280154eSStefano Zampini args: -M 3 -N 3 -K 3 -local {{0 1}} -testcircular 253*6280154eSStefano Zampini 254*6280154eSStefano Zampini test: 255*6280154eSStefano Zampini output_file: output/ex70_1.out 256*6280154eSStefano Zampini suffix: 5 257*6280154eSStefano Zampini nsize: {{2 4}} 258*6280154eSStefano Zampini args: -M 3 -N 3 -K 3 -local 1 -testcircular -testtranspose 0 259*6280154eSStefano Zampini 260*6280154eSStefano Zampini test: 261*6280154eSStefano Zampini TODO: MatCreateDense broken with some processors not having local rows 262*6280154eSStefano Zampini output_file: output/ex70_1.out 263*6280154eSStefano Zampini suffix: 5_broken 264*6280154eSStefano Zampini nsize: {{2 4}} 265*6280154eSStefano Zampini args: -M 3 -N 3 -K 3 -local 0 -testcircular -testtranspose 0 266*6280154eSStefano Zampini 267*6280154eSStefano Zampini test: 268*6280154eSStefano Zampini output_file: output/ex70_1.out 269*6280154eSStefano Zampini suffix: 6 270*6280154eSStefano Zampini nsize: 1 271*6280154eSStefano Zampini args: -M {{1 3}} -N {{2 5}} -K {{1 2}} -local {{0 1}} -testcircular -testtranspose 0 272*6280154eSStefano Zampini 273*6280154eSStefano Zampini test: 274*6280154eSStefano Zampini TODO: MatTransposeMatMultSymbolic_SeqAIJ_SeqDense plays with the destroy routine 275*6280154eSStefano Zampini output_file: output/ex70_1.out 276*6280154eSStefano Zampini suffix: 6_broken 277*6280154eSStefano Zampini nsize: 1 278*6280154eSStefano Zampini args: -M {{1 3}} -N {{2 5}} -K {{1 2}} -local {{0 1}} -testcircular -testtranspose 279*6280154eSStefano Zampini 280*6280154eSStefano Zampini test: 281*6280154eSStefano Zampini output_file: output/ex70_1.out 282*6280154eSStefano Zampini suffix: 7 283*6280154eSStefano Zampini nsize: 1 284*6280154eSStefano Zampini args: -M 13 -N 13 -K {{1 3}} -local {{0 1}} -mat_type dense -testnest 0 -testcircular 285*6280154eSStefano Zampini 286*6280154eSStefano Zampini test: 287*6280154eSStefano Zampini TODO: NEST x DENSE with dense nested matrices seems broken in this case 288*6280154eSStefano Zampini output_file: output/ex70_1.out 289*6280154eSStefano Zampini suffix: 7_broken 290*6280154eSStefano Zampini nsize: 1 291*6280154eSStefano Zampini args: -M 13 -N 13 -K {{1 3}} -local {{0 1}} -mat_type dense -testnest -testcircular 292*6280154eSStefano Zampini TEST*/ 293