1c4fb69feSStefano Zampini static char help[] = "Tests MATHARA\n\n"; 2c4fb69feSStefano Zampini 3c4fb69feSStefano Zampini #include <petscmat.h> 4c4fb69feSStefano Zampini #include <petscsf.h> 5c4fb69feSStefano Zampini 6c4fb69feSStefano Zampini static PetscScalar GenEntry_Symm(PetscInt sdim, PetscReal x[], PetscReal y[], void *ctx) 7c4fb69feSStefano Zampini { 8c4fb69feSStefano Zampini PetscInt d; 9c4fb69feSStefano Zampini PetscReal clength = sdim == 3 ? 0.2 : 0.1; 10c4fb69feSStefano Zampini PetscReal dist, diff = 0.0; 11c4fb69feSStefano Zampini 12c4fb69feSStefano Zampini for (d = 0; d < sdim; d++) { diff += (x[d] - y[d]) * (x[d] - y[d]); } 13c4fb69feSStefano Zampini dist = PetscSqrtReal(diff); 14c4fb69feSStefano Zampini return PetscExpReal(-dist / clength); 15c4fb69feSStefano Zampini } 16c4fb69feSStefano Zampini 17c4fb69feSStefano Zampini static PetscScalar GenEntry_Unsymm(PetscInt sdim, PetscReal x[], PetscReal y[], void *ctx) 18c4fb69feSStefano Zampini { 19c4fb69feSStefano Zampini PetscInt d; 20c4fb69feSStefano Zampini PetscReal clength = sdim == 3 ? 0.2 : 0.1; 21c4fb69feSStefano Zampini PetscReal dist, diff = 0.0, nx = 0.0, ny = 0.0; 22c4fb69feSStefano Zampini 23c4fb69feSStefano Zampini for (d = 0; d < sdim; d++) { nx += x[d]*x[d]; } 24c4fb69feSStefano Zampini for (d = 0; d < sdim; d++) { ny += y[d]*y[d]; } 25c4fb69feSStefano Zampini for (d = 0; d < sdim; d++) { diff += (x[d] - y[d]) * (x[d] - y[d]); } 26c4fb69feSStefano Zampini dist = PetscSqrtReal(diff); 27c4fb69feSStefano Zampini return nx > ny ? PetscExpReal(-dist / clength) : PetscExpReal(-dist / clength) + 1.; 28c4fb69feSStefano Zampini } 29c4fb69feSStefano Zampini 30c4fb69feSStefano Zampini int main(int argc,char **argv) 31c4fb69feSStefano Zampini { 32c4fb69feSStefano Zampini Mat A,B,C,D; 33c4fb69feSStefano Zampini Vec v,x,y,Ax,Ay,Bx,By; 34c4fb69feSStefano Zampini PetscRandom r; 35c4fb69feSStefano Zampini PetscLayout map; 36c4fb69feSStefano Zampini PetscScalar *Adata = NULL, *Cdata = NULL, scale = 1.0; 37c4fb69feSStefano Zampini PetscReal *coords,nA,nD,nB,err,nX,norms[3]; 38c4fb69feSStefano Zampini PetscInt N, n = 64, dim = 1, i, j, nrhs = 11, lda = 0, ldc = 0, nt, ntrials = 2; 39c4fb69feSStefano Zampini PetscMPIInt size,rank; 40c4fb69feSStefano Zampini PetscBool testlayout = PETSC_FALSE,flg,symm = PETSC_FALSE, Asymm = PETSC_TRUE, kernel = PETSC_TRUE; 41c4fb69feSStefano Zampini PetscBool checkexpl = PETSC_FALSE, agpu = PETSC_FALSE, bgpu = PETSC_FALSE, cgpu = PETSC_FALSE; 42c4fb69feSStefano Zampini PetscBool testtrans, testnorm, randommat = PETSC_TRUE; 43c4fb69feSStefano Zampini void (*approxnormfunc)(void); 44c4fb69feSStefano Zampini void (*Anormfunc)(void); 45c4fb69feSStefano Zampini PetscErrorCode ierr; 46c4fb69feSStefano Zampini 47c4fb69feSStefano Zampini PETSC_MPI_THREAD_REQUIRED = MPI_THREAD_MULTIPLE; 48c4fb69feSStefano Zampini ierr = PetscInitialize(&argc,&argv,(char*) 0,help);if (ierr) return ierr; 49c4fb69feSStefano Zampini ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr); 50c4fb69feSStefano Zampini ierr = PetscOptionsGetInt(NULL,NULL,"-nrhs",&nrhs,NULL);CHKERRQ(ierr); 51c4fb69feSStefano Zampini ierr = PetscOptionsGetInt(NULL,NULL,"-dim",&dim,NULL);CHKERRQ(ierr); 52c4fb69feSStefano Zampini ierr = PetscOptionsGetInt(NULL,NULL,"-lda",&lda,NULL);CHKERRQ(ierr); 53c4fb69feSStefano Zampini ierr = PetscOptionsGetInt(NULL,NULL,"-ldc",&ldc,NULL);CHKERRQ(ierr); 54c4fb69feSStefano Zampini ierr = PetscOptionsGetInt(NULL,NULL,"-matmattrials",&ntrials,NULL);CHKERRQ(ierr); 55c4fb69feSStefano Zampini ierr = PetscOptionsGetBool(NULL,NULL,"-randommat",&randommat,NULL);CHKERRQ(ierr); 56c4fb69feSStefano Zampini ierr = PetscOptionsGetBool(NULL,NULL,"-testlayout",&testlayout,NULL);CHKERRQ(ierr); 57c4fb69feSStefano Zampini ierr = PetscOptionsGetBool(NULL,NULL,"-Asymm",&Asymm,NULL);CHKERRQ(ierr); 58c4fb69feSStefano Zampini ierr = PetscOptionsGetBool(NULL,NULL,"-symm",&symm,NULL);CHKERRQ(ierr); 59c4fb69feSStefano Zampini ierr = PetscOptionsGetBool(NULL,NULL,"-kernel",&kernel,NULL);CHKERRQ(ierr); 60c4fb69feSStefano Zampini ierr = PetscOptionsGetBool(NULL,NULL,"-checkexpl",&checkexpl,NULL);CHKERRQ(ierr); 61c4fb69feSStefano Zampini ierr = PetscOptionsGetBool(NULL,NULL,"-agpu",&agpu,NULL);CHKERRQ(ierr); 62c4fb69feSStefano Zampini ierr = PetscOptionsGetBool(NULL,NULL,"-bgpu",&bgpu,NULL);CHKERRQ(ierr); 63c4fb69feSStefano Zampini ierr = PetscOptionsGetBool(NULL,NULL,"-cgpu",&cgpu,NULL);CHKERRQ(ierr); 64c4fb69feSStefano Zampini ierr = PetscOptionsGetScalar(NULL,NULL,"-scale",&scale,NULL);CHKERRQ(ierr); 65c4fb69feSStefano Zampini if (!Asymm) symm = PETSC_FALSE; 66c4fb69feSStefano Zampini 67*ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 68c4fb69feSStefano Zampini /* MatMultTranspose for nonsymmetric matrices not implemented */ 69c4fb69feSStefano Zampini testtrans = (PetscBool)(size == 1 || symm); 70c4fb69feSStefano Zampini testnorm = (PetscBool)(size == 1 || symm); 71c4fb69feSStefano Zampini 72*ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 73c4fb69feSStefano Zampini ierr = PetscLayoutCreate(PETSC_COMM_WORLD,&map);CHKERRQ(ierr); 74c4fb69feSStefano Zampini if (testlayout) { 75c4fb69feSStefano Zampini if (rank%2) n = PetscMax(2*n-5*rank,0); 76c4fb69feSStefano Zampini else n = 2*n+rank; 77c4fb69feSStefano Zampini } 78c4fb69feSStefano Zampini ierr = PetscLayoutSetLocalSize(map,n);CHKERRQ(ierr); 79c4fb69feSStefano Zampini ierr = PetscLayoutSetUp(map);CHKERRQ(ierr); 80c4fb69feSStefano Zampini ierr = PetscLayoutGetSize(map,&N);CHKERRQ(ierr); 81c4fb69feSStefano Zampini ierr = PetscLayoutDestroy(&map);CHKERRQ(ierr); 82c4fb69feSStefano Zampini 83c4fb69feSStefano Zampini ierr = PetscRandomCreate(PETSC_COMM_WORLD,&r);CHKERRQ(ierr); 84c4fb69feSStefano Zampini ierr = PetscRandomSetFromOptions(r);CHKERRQ(ierr); 85c4fb69feSStefano Zampini if (lda) { 86c4fb69feSStefano Zampini ierr = PetscMalloc1(N*(n+lda),&Adata);CHKERRQ(ierr); 87c4fb69feSStefano Zampini } 88c4fb69feSStefano Zampini ierr = MatCreateDense(PETSC_COMM_WORLD,n,n,N,N,Adata,&A);CHKERRQ(ierr); 89c4fb69feSStefano Zampini ierr = MatDenseSetLDA(A,n+lda);CHKERRQ(ierr); 90c4fb69feSStefano Zampini ierr = PetscMalloc1(n*dim,&coords);CHKERRQ(ierr); 91c4fb69feSStefano Zampini for (j = 0; j < n; j++) { 92c4fb69feSStefano Zampini for (i = 0; i < dim; i++) { 93c4fb69feSStefano Zampini PetscScalar a; 94c4fb69feSStefano Zampini 95c4fb69feSStefano Zampini ierr = PetscRandomGetValue(r,&a);CHKERRQ(ierr); 96c4fb69feSStefano Zampini coords[j*dim + i] = PetscRealPart(a); 97c4fb69feSStefano Zampini } 98c4fb69feSStefano Zampini } 99c4fb69feSStefano Zampini if (kernel || !randommat) { 100c4fb69feSStefano Zampini MatHaraKernel k = Asymm ? GenEntry_Symm : GenEntry_Unsymm; 101c4fb69feSStefano Zampini PetscInt ist,ien; 102c4fb69feSStefano Zampini PetscReal *gcoords; 103c4fb69feSStefano Zampini 104c4fb69feSStefano Zampini if (size > 1) { /* replicated coords so that we can populate the dense matrix */ 105c4fb69feSStefano Zampini PetscSF sf; 106c4fb69feSStefano Zampini MPI_Datatype dtype; 107c4fb69feSStefano Zampini 108*ffc4695bSBarry Smith ierr = MPI_Type_contiguous(dim,MPIU_REAL,&dtype);CHKERRMPI(ierr); 109*ffc4695bSBarry Smith ierr = MPI_Type_commit(&dtype);CHKERRMPI(ierr); 110c4fb69feSStefano Zampini 111c4fb69feSStefano Zampini ierr = PetscSFCreate(PETSC_COMM_WORLD,&sf);CHKERRQ(ierr); 112c4fb69feSStefano Zampini ierr = MatGetLayouts(A,&map,NULL);CHKERRQ(ierr); 113c4fb69feSStefano Zampini ierr = PetscSFSetGraphWithPattern(sf,map,PETSCSF_PATTERN_ALLGATHER);CHKERRQ(ierr); 114c4fb69feSStefano Zampini ierr = PetscMalloc1(dim*N,&gcoords);CHKERRQ(ierr); 115c4fb69feSStefano Zampini ierr = PetscSFBcastBegin(sf,dtype,coords,gcoords);CHKERRQ(ierr); 116c4fb69feSStefano Zampini ierr = PetscSFBcastEnd(sf,dtype,coords,gcoords);CHKERRQ(ierr); 117c4fb69feSStefano Zampini ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 118*ffc4695bSBarry Smith ierr = MPI_Type_free(&dtype);CHKERRMPI(ierr); 119c4fb69feSStefano Zampini } else gcoords = (PetscReal*)coords; 120c4fb69feSStefano Zampini 121c4fb69feSStefano Zampini ierr = MatGetOwnershipRange(A,&ist,&ien);CHKERRQ(ierr); 122c4fb69feSStefano Zampini for (i = ist; i < ien; i++) { 123c4fb69feSStefano Zampini for (j = 0; j < N; j++) { 124c4fb69feSStefano Zampini ierr = MatSetValue(A,i,j,(*k)(dim,gcoords + i*dim,gcoords + j*dim,NULL),INSERT_VALUES);CHKERRQ(ierr); 125c4fb69feSStefano Zampini } 126c4fb69feSStefano Zampini } 127c4fb69feSStefano Zampini ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 128c4fb69feSStefano Zampini ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 129c4fb69feSStefano Zampini if (kernel) { 130c4fb69feSStefano Zampini ierr = MatCreateHaraFromKernel(PETSC_COMM_WORLD,n,n,N,N,dim,coords,k,NULL,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,&B);CHKERRQ(ierr); 131c4fb69feSStefano Zampini } else { 132c4fb69feSStefano Zampini ierr = MatCreateHaraFromMat(A,dim,coords,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,&B);CHKERRQ(ierr); 133c4fb69feSStefano Zampini } 134c4fb69feSStefano Zampini if (gcoords != coords) { ierr = PetscFree(gcoords);CHKERRQ(ierr); } 135c4fb69feSStefano Zampini } else { 136c4fb69feSStefano Zampini ierr = MatSetRandom(A,r);CHKERRQ(ierr); 137c4fb69feSStefano Zampini if (Asymm) { 138c4fb69feSStefano Zampini ierr = MatTranspose(A,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); 139c4fb69feSStefano Zampini ierr = MatAXPY(A,1.0,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 140c4fb69feSStefano Zampini ierr = MatDestroy(&B);CHKERRQ(ierr); 141c4fb69feSStefano Zampini ierr = MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr); 142c4fb69feSStefano Zampini } 143c4fb69feSStefano Zampini ierr = MatCreateHaraFromMat(A,dim,coords,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,&B);CHKERRQ(ierr); 144c4fb69feSStefano Zampini } 145c4fb69feSStefano Zampini if (agpu) { 146c4fb69feSStefano Zampini ierr = MatConvert(A,MATDENSECUDA,MAT_INPLACE_MATRIX,&A);CHKERRQ(ierr); 147c4fb69feSStefano Zampini } 148c4fb69feSStefano Zampini ierr = MatViewFromOptions(A,NULL,"-A_view");CHKERRQ(ierr); 149c4fb69feSStefano Zampini 150c4fb69feSStefano Zampini ierr = MatSetOption(B,MAT_SYMMETRIC,symm);CHKERRQ(ierr); 151c4fb69feSStefano Zampini ierr = PetscFree(coords);CHKERRQ(ierr); 152c4fb69feSStefano Zampini 153c4fb69feSStefano Zampini /* assemble the H-matrix */ 154c4fb69feSStefano Zampini ierr = MatBindToCPU(B,(PetscBool)!bgpu);CHKERRQ(ierr); 155c4fb69feSStefano Zampini ierr = MatSetFromOptions(B);CHKERRQ(ierr); 156c4fb69feSStefano Zampini ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 157c4fb69feSStefano Zampini ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 158c4fb69feSStefano Zampini ierr = MatViewFromOptions(B,NULL,"-B_view");CHKERRQ(ierr); 159c4fb69feSStefano Zampini 160c4fb69feSStefano Zampini /* Test MatScale */ 161c4fb69feSStefano Zampini ierr = MatScale(A,scale);CHKERRQ(ierr); 162c4fb69feSStefano Zampini ierr = MatScale(B,scale);CHKERRQ(ierr); 163c4fb69feSStefano Zampini 164c4fb69feSStefano Zampini /* Test MatMult */ 165c4fb69feSStefano Zampini ierr = MatCreateVecs(A,&Ax,&Ay);CHKERRQ(ierr); 166c4fb69feSStefano Zampini ierr = MatCreateVecs(B,&Bx,&By);CHKERRQ(ierr); 167c4fb69feSStefano Zampini ierr = VecSetRandom(Ax,r);CHKERRQ(ierr); 168c4fb69feSStefano Zampini ierr = VecCopy(Ax,Bx);CHKERRQ(ierr); 169c4fb69feSStefano Zampini ierr = MatMult(A,Ax,Ay);CHKERRQ(ierr); 170c4fb69feSStefano Zampini ierr = MatMult(B,Bx,By);CHKERRQ(ierr); 171c4fb69feSStefano Zampini ierr = VecViewFromOptions(Ay,NULL,"-mult_vec_view");CHKERRQ(ierr); 172c4fb69feSStefano Zampini ierr = VecViewFromOptions(By,NULL,"-mult_vec_view");CHKERRQ(ierr); 173c4fb69feSStefano Zampini ierr = VecNorm(Ay,NORM_INFINITY,&nX);CHKERRQ(ierr); 174c4fb69feSStefano Zampini ierr = VecAXPY(Ay,-1.0,By);CHKERRQ(ierr); 175c4fb69feSStefano Zampini ierr = VecViewFromOptions(Ay,NULL,"-mult_vec_view");CHKERRQ(ierr); 176c4fb69feSStefano Zampini ierr = VecNorm(Ay,NORM_INFINITY,&err);CHKERRQ(ierr); 177c4fb69feSStefano Zampini ierr = PetscPrintf(PETSC_COMM_WORLD,"MatMult err %g\n",err/nX);CHKERRQ(ierr); 178c4fb69feSStefano Zampini ierr = VecScale(By,-1.0);CHKERRQ(ierr); 179c4fb69feSStefano Zampini ierr = MatMultAdd(B,Bx,By,By);CHKERRQ(ierr); 180c4fb69feSStefano Zampini ierr = VecNorm(By,NORM_INFINITY,&err);CHKERRQ(ierr); 181c4fb69feSStefano Zampini ierr = VecViewFromOptions(By,NULL,"-mult_vec_view");CHKERRQ(ierr); 182c4fb69feSStefano Zampini if (err > PETSC_SMALL) { 183c4fb69feSStefano Zampini ierr = PetscPrintf(PETSC_COMM_WORLD,"MatMultAdd err %g\n",err);CHKERRQ(ierr); 184c4fb69feSStefano Zampini } 185c4fb69feSStefano Zampini 186c4fb69feSStefano Zampini /* Test MatNorm */ 187c4fb69feSStefano Zampini ierr = MatNorm(A,NORM_INFINITY,&norms[0]);CHKERRQ(ierr); 188c4fb69feSStefano Zampini ierr = MatNorm(A,NORM_1,&norms[1]);CHKERRQ(ierr); 189c4fb69feSStefano Zampini norms[2] = -1.; /* NORM_2 not supported */ 190c4fb69feSStefano Zampini ierr = PetscPrintf(PETSC_COMM_WORLD,"A Matrix norms: infty=%g, norm_1=%g, norm_2=%g\n",(double)norms[0],(double)norms[1],(double)norms[2]);CHKERRQ(ierr); 191c4fb69feSStefano Zampini ierr = MatGetOperation(A,MATOP_NORM,&Anormfunc);CHKERRQ(ierr); 192c4fb69feSStefano Zampini ierr = MatGetOperation(B,MATOP_NORM,&approxnormfunc);CHKERRQ(ierr); 193c4fb69feSStefano Zampini ierr = MatSetOperation(A,MATOP_NORM,approxnormfunc);CHKERRQ(ierr); 194c4fb69feSStefano Zampini ierr = MatNorm(A,NORM_INFINITY,&norms[0]);CHKERRQ(ierr); 195c4fb69feSStefano Zampini ierr = MatNorm(A,NORM_1,&norms[1]);CHKERRQ(ierr); 196c4fb69feSStefano Zampini ierr = MatNorm(A,NORM_2,&norms[2]);CHKERRQ(ierr); 197c4fb69feSStefano Zampini ierr = PetscPrintf(PETSC_COMM_WORLD,"A Approx Matrix norms: infty=%g, norm_1=%g, norm_2=%g\n",(double)norms[0],(double)norms[1],(double)norms[2]);CHKERRQ(ierr); 198c4fb69feSStefano Zampini if (testnorm) { 199c4fb69feSStefano Zampini ierr = MatNorm(B,NORM_INFINITY,&norms[0]);CHKERRQ(ierr); 200c4fb69feSStefano Zampini ierr = MatNorm(B,NORM_1,&norms[1]);CHKERRQ(ierr); 201c4fb69feSStefano Zampini ierr = MatNorm(B,NORM_2,&norms[2]);CHKERRQ(ierr); 202c4fb69feSStefano Zampini } else { 203c4fb69feSStefano Zampini norms[0] = -1.; 204c4fb69feSStefano Zampini norms[1] = -1.; 205c4fb69feSStefano Zampini norms[2] = -1.; 206c4fb69feSStefano Zampini } 207c4fb69feSStefano Zampini ierr = PetscPrintf(PETSC_COMM_WORLD,"B Approx Matrix norms: infty=%g, norm_1=%g, norm_2=%g\n",(double)norms[0],(double)norms[1],(double)norms[2]);CHKERRQ(ierr); 208c4fb69feSStefano Zampini ierr = MatSetOperation(A,MATOP_NORM,Anormfunc);CHKERRQ(ierr); 209c4fb69feSStefano Zampini 210c4fb69feSStefano Zampini /* Test MatDuplicate */ 211c4fb69feSStefano Zampini ierr = MatDuplicate(B,MAT_COPY_VALUES,&D);CHKERRQ(ierr); 212c4fb69feSStefano Zampini ierr = MatMultEqual(B,D,10,&flg);CHKERRQ(ierr); 213c4fb69feSStefano Zampini if (!flg) { 214c4fb69feSStefano Zampini ierr = PetscPrintf(PETSC_COMM_WORLD,"MatMult error after MatDuplicate\n");CHKERRQ(ierr); 215c4fb69feSStefano Zampini } 216c4fb69feSStefano Zampini if (testtrans) { 217c4fb69feSStefano Zampini ierr = MatMultTransposeEqual(B,D,10,&flg);CHKERRQ(ierr); 218c4fb69feSStefano Zampini if (!flg) { 219c4fb69feSStefano Zampini ierr = PetscPrintf(PETSC_COMM_WORLD,"MatMultTranpose error after MatDuplicate\n");CHKERRQ(ierr); 220c4fb69feSStefano Zampini } 221c4fb69feSStefano Zampini } 222c4fb69feSStefano Zampini ierr = MatDestroy(&D);CHKERRQ(ierr); 223c4fb69feSStefano Zampini 224c4fb69feSStefano Zampini if (testtrans) { /* MatMultTranspose for nonsymmetric matrices not implemented */ 225c4fb69feSStefano Zampini ierr = VecSetRandom(Ay,r);CHKERRQ(ierr); 226c4fb69feSStefano Zampini ierr = VecCopy(Ay,By);CHKERRQ(ierr); 227c4fb69feSStefano Zampini ierr = MatMultTranspose(A,Ay,Ax);CHKERRQ(ierr); 228c4fb69feSStefano Zampini ierr = MatMultTranspose(B,By,Bx);CHKERRQ(ierr); 229c4fb69feSStefano Zampini ierr = VecViewFromOptions(Ax,NULL,"-multtrans_vec_view");CHKERRQ(ierr); 230c4fb69feSStefano Zampini ierr = VecViewFromOptions(Bx,NULL,"-multtrans_vec_view");CHKERRQ(ierr); 231c4fb69feSStefano Zampini ierr = VecNorm(Ax,NORM_INFINITY,&nX);CHKERRQ(ierr); 232c4fb69feSStefano Zampini ierr = VecAXPY(Ax,-1.0,Bx);CHKERRQ(ierr); 233c4fb69feSStefano Zampini ierr = VecViewFromOptions(Ax,NULL,"-multtrans_vec_view");CHKERRQ(ierr); 234c4fb69feSStefano Zampini ierr = VecNorm(Ax,NORM_INFINITY,&err);CHKERRQ(ierr); 235c4fb69feSStefano Zampini ierr = PetscPrintf(PETSC_COMM_WORLD,"MatMultTranspose err %g\n",err/nX);CHKERRQ(ierr); 236c4fb69feSStefano Zampini ierr = VecScale(Bx,-1.0);CHKERRQ(ierr); 237c4fb69feSStefano Zampini ierr = MatMultTransposeAdd(B,By,Bx,Bx);CHKERRQ(ierr); 238c4fb69feSStefano Zampini ierr = VecNorm(Bx,NORM_INFINITY,&err);CHKERRQ(ierr); 239c4fb69feSStefano Zampini if (err > PETSC_SMALL) { 240c4fb69feSStefano Zampini ierr = PetscPrintf(PETSC_COMM_WORLD,"MatMultTransposeAdd err %g\n",err);CHKERRQ(ierr); 241c4fb69feSStefano Zampini } 242c4fb69feSStefano Zampini } 243c4fb69feSStefano Zampini ierr = VecDestroy(&Ax);CHKERRQ(ierr); 244c4fb69feSStefano Zampini ierr = VecDestroy(&Ay);CHKERRQ(ierr); 245c4fb69feSStefano Zampini ierr = VecDestroy(&Bx);CHKERRQ(ierr); 246c4fb69feSStefano Zampini ierr = VecDestroy(&By);CHKERRQ(ierr); 247c4fb69feSStefano Zampini 248c4fb69feSStefano Zampini /* Test MatMatMult */ 249c4fb69feSStefano Zampini if (ldc) { 250c4fb69feSStefano Zampini ierr = PetscMalloc1(nrhs*(n+ldc),&Cdata);CHKERRQ(ierr); 251c4fb69feSStefano Zampini } 252c4fb69feSStefano Zampini ierr = MatCreateDense(PETSC_COMM_WORLD,n,PETSC_DECIDE,N,nrhs,Cdata,&C);CHKERRQ(ierr); 253c4fb69feSStefano Zampini ierr = MatDenseSetLDA(C,n+ldc);CHKERRQ(ierr); 254c4fb69feSStefano Zampini ierr = MatSetRandom(C,r);CHKERRQ(ierr); 255c4fb69feSStefano Zampini if (cgpu) { 256c4fb69feSStefano Zampini ierr = MatConvert(C,MATDENSECUDA,MAT_INPLACE_MATRIX,&C);CHKERRQ(ierr); 257c4fb69feSStefano Zampini } 258c4fb69feSStefano Zampini for (nt = 0; nt < ntrials; nt++) { 259c4fb69feSStefano Zampini ierr = MatMatMult(B,C,nt ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX,PETSC_DEFAULT,&D);CHKERRQ(ierr); 260c4fb69feSStefano Zampini ierr = MatViewFromOptions(D,NULL,"-bc_view");CHKERRQ(ierr); 261c4fb69feSStefano Zampini ierr = PetscObjectBaseTypeCompareAny((PetscObject)D,&flg,MATSEQDENSE,MATMPIDENSE,"");CHKERRQ(ierr); 262c4fb69feSStefano Zampini if (flg) { 263c4fb69feSStefano Zampini ierr = MatCreateVecs(B,&x,&y);CHKERRQ(ierr); 264c4fb69feSStefano Zampini ierr = MatCreateVecs(D,NULL,&v);CHKERRQ(ierr); 265c4fb69feSStefano Zampini for (i = 0; i < nrhs; i++) { 266c4fb69feSStefano Zampini ierr = MatGetColumnVector(D,v,i);CHKERRQ(ierr); 267c4fb69feSStefano Zampini ierr = MatGetColumnVector(C,x,i);CHKERRQ(ierr); 268c4fb69feSStefano Zampini ierr = MatMult(B,x,y);CHKERRQ(ierr); 269c4fb69feSStefano Zampini ierr = VecAXPY(y,-1.0,v);CHKERRQ(ierr); 270c4fb69feSStefano Zampini ierr = VecNorm(y,NORM_INFINITY,&err);CHKERRQ(ierr); 271c4fb69feSStefano Zampini if (err > PETSC_SMALL) { ierr = PetscPrintf(PETSC_COMM_WORLD,"MatMat err %D %g\n",i,err);CHKERRQ(ierr); } 272c4fb69feSStefano Zampini } 273c4fb69feSStefano Zampini ierr = VecDestroy(&y);CHKERRQ(ierr); 274c4fb69feSStefano Zampini ierr = VecDestroy(&x);CHKERRQ(ierr); 275c4fb69feSStefano Zampini ierr = VecDestroy(&v);CHKERRQ(ierr); 276c4fb69feSStefano Zampini } 277c4fb69feSStefano Zampini } 278c4fb69feSStefano Zampini ierr = MatDestroy(&D);CHKERRQ(ierr); 279c4fb69feSStefano Zampini 280c4fb69feSStefano Zampini /* Test MatTransposeMatMult */ 281c4fb69feSStefano Zampini if (testtrans) { /* MatMultTranspose for nonsymmetric matrices not implemented */ 282c4fb69feSStefano Zampini for (nt = 0; nt < ntrials; nt++) { 283c4fb69feSStefano Zampini ierr = MatTransposeMatMult(B,C,nt ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX,PETSC_DEFAULT,&D);CHKERRQ(ierr); 284c4fb69feSStefano Zampini ierr = MatViewFromOptions(D,NULL,"-btc_view");CHKERRQ(ierr); 285c4fb69feSStefano Zampini ierr = PetscObjectBaseTypeCompareAny((PetscObject)D,&flg,MATSEQDENSE,MATMPIDENSE,"");CHKERRQ(ierr); 286c4fb69feSStefano Zampini if (flg) { 287c4fb69feSStefano Zampini ierr = MatCreateVecs(B,&y,&x);CHKERRQ(ierr); 288c4fb69feSStefano Zampini ierr = MatCreateVecs(D,NULL,&v);CHKERRQ(ierr); 289c4fb69feSStefano Zampini for (i = 0; i < nrhs; i++) { 290c4fb69feSStefano Zampini ierr = MatGetColumnVector(D,v,i);CHKERRQ(ierr); 291c4fb69feSStefano Zampini ierr = MatGetColumnVector(C,x,i);CHKERRQ(ierr); 292c4fb69feSStefano Zampini ierr = MatMultTranspose(B,x,y);CHKERRQ(ierr); 293c4fb69feSStefano Zampini ierr = VecAXPY(y,-1.0,v);CHKERRQ(ierr); 294c4fb69feSStefano Zampini ierr = VecNorm(y,NORM_INFINITY,&err);CHKERRQ(ierr); 295c4fb69feSStefano Zampini if (err > PETSC_SMALL) { ierr = PetscPrintf(PETSC_COMM_WORLD,"MatTransMat err %D %g\n",i,err);CHKERRQ(ierr); } 296c4fb69feSStefano Zampini } 297c4fb69feSStefano Zampini ierr = VecDestroy(&y);CHKERRQ(ierr); 298c4fb69feSStefano Zampini ierr = VecDestroy(&x);CHKERRQ(ierr); 299c4fb69feSStefano Zampini ierr = VecDestroy(&v);CHKERRQ(ierr); 300c4fb69feSStefano Zampini } 301c4fb69feSStefano Zampini } 302c4fb69feSStefano Zampini ierr = MatDestroy(&D);CHKERRQ(ierr); 303c4fb69feSStefano Zampini } 304c4fb69feSStefano Zampini 305c4fb69feSStefano Zampini /* check explicit operator */ 306c4fb69feSStefano Zampini if (checkexpl) { 307c4fb69feSStefano Zampini Mat Be, Bet; 308c4fb69feSStefano Zampini 309c4fb69feSStefano Zampini ierr = MatComputeOperator(B,MATDENSE,&D);CHKERRQ(ierr); 310c4fb69feSStefano Zampini ierr = MatDuplicate(D,MAT_COPY_VALUES,&Be);CHKERRQ(ierr); 311c4fb69feSStefano Zampini ierr = MatNorm(D,NORM_FROBENIUS,&nB);CHKERRQ(ierr); 312c4fb69feSStefano Zampini ierr = MatViewFromOptions(D,NULL,"-expl_view");CHKERRQ(ierr); 313c4fb69feSStefano Zampini ierr = MatAXPY(D,-1.0,A,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 314c4fb69feSStefano Zampini ierr = MatViewFromOptions(D,NULL,"-diff_view");CHKERRQ(ierr); 315c4fb69feSStefano Zampini ierr = MatNorm(D,NORM_FROBENIUS,&nD);CHKERRQ(ierr); 316c4fb69feSStefano Zampini ierr = MatNorm(A,NORM_FROBENIUS,&nA);CHKERRQ(ierr); 317c4fb69feSStefano Zampini ierr = PetscPrintf(PETSC_COMM_WORLD,"Approximation error %g (%g / %g, %g)\n",nD/nA,nD,nA,nB);CHKERRQ(ierr); 318c4fb69feSStefano Zampini ierr = MatDestroy(&D);CHKERRQ(ierr); 319c4fb69feSStefano Zampini 320c4fb69feSStefano Zampini if (testtrans) { /* MatMultTranspose for nonsymmetric matrices not implemented */ 321c4fb69feSStefano Zampini ierr = MatTranspose(A,MAT_INPLACE_MATRIX,&A);CHKERRQ(ierr); 322c4fb69feSStefano Zampini ierr = MatComputeOperatorTranspose(B,MATDENSE,&D);CHKERRQ(ierr); 323c4fb69feSStefano Zampini ierr = MatDuplicate(D,MAT_COPY_VALUES,&Bet);CHKERRQ(ierr); 324c4fb69feSStefano Zampini ierr = MatNorm(D,NORM_FROBENIUS,&nB);CHKERRQ(ierr); 325c4fb69feSStefano Zampini ierr = MatViewFromOptions(D,NULL,"-expl_trans_view");CHKERRQ(ierr); 326c4fb69feSStefano Zampini ierr = MatAXPY(D,-1.0,A,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 327c4fb69feSStefano Zampini ierr = MatViewFromOptions(D,NULL,"-diff_trans_view");CHKERRQ(ierr); 328c4fb69feSStefano Zampini ierr = MatNorm(D,NORM_FROBENIUS,&nD);CHKERRQ(ierr); 329c4fb69feSStefano Zampini ierr = MatNorm(A,NORM_FROBENIUS,&nA);CHKERRQ(ierr); 330c4fb69feSStefano Zampini ierr = PetscPrintf(PETSC_COMM_WORLD,"Approximation error transpose %g (%g / %g, %g)\n",nD/nA,nD,nA,nB);CHKERRQ(ierr); 331c4fb69feSStefano Zampini ierr = MatDestroy(&D);CHKERRQ(ierr); 332c4fb69feSStefano Zampini 333c4fb69feSStefano Zampini ierr = MatTranspose(Bet,MAT_INPLACE_MATRIX,&Bet);CHKERRQ(ierr); 334c4fb69feSStefano Zampini ierr = MatAXPY(Be,-1.0,Bet,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 335c4fb69feSStefano Zampini ierr = MatViewFromOptions(Be,NULL,"-diff_expl_view");CHKERRQ(ierr); 336c4fb69feSStefano Zampini ierr = MatNorm(Be,NORM_FROBENIUS,&nB);CHKERRQ(ierr); 337c4fb69feSStefano Zampini ierr = PetscPrintf(PETSC_COMM_WORLD,"Approximation error B - (B^T)^T %g\n",nB);CHKERRQ(ierr); 338c4fb69feSStefano Zampini ierr = MatDestroy(&Be);CHKERRQ(ierr); 339c4fb69feSStefano Zampini ierr = MatDestroy(&Bet);CHKERRQ(ierr); 340c4fb69feSStefano Zampini } 341c4fb69feSStefano Zampini } 342c4fb69feSStefano Zampini ierr = MatDestroy(&A);CHKERRQ(ierr); 343c4fb69feSStefano Zampini ierr = MatDestroy(&B);CHKERRQ(ierr); 344c4fb69feSStefano Zampini ierr = MatDestroy(&C);CHKERRQ(ierr); 345c4fb69feSStefano Zampini ierr = PetscRandomDestroy(&r);CHKERRQ(ierr); 346c4fb69feSStefano Zampini ierr = PetscFree(Cdata);CHKERRQ(ierr); 347c4fb69feSStefano Zampini ierr = PetscFree(Adata);CHKERRQ(ierr); 348c4fb69feSStefano Zampini ierr = PetscFinalize(); 349c4fb69feSStefano Zampini return ierr; 350c4fb69feSStefano Zampini } 351c4fb69feSStefano Zampini 352c4fb69feSStefano Zampini /*TEST 353c4fb69feSStefano Zampini 354c4fb69feSStefano Zampini build: 355c4fb69feSStefano Zampini requires: hara 356c4fb69feSStefano Zampini 357c4fb69feSStefano Zampini #tests from kernel 358c4fb69feSStefano Zampini test: 359c4fb69feSStefano Zampini requires: hara 360c4fb69feSStefano Zampini nsize: 1 361c4fb69feSStefano Zampini suffix: 1 362c4fb69feSStefano Zampini args: -n {{17 33}} -kernel 1 -dim {{1 2 3}} -symm {{0 1}} -checkexpl -bgpu 0 363c4fb69feSStefano Zampini 364c4fb69feSStefano Zampini test: 365c4fb69feSStefano Zampini requires: hara 366c4fb69feSStefano Zampini nsize: 1 367c4fb69feSStefano Zampini suffix: 1_ld 368c4fb69feSStefano Zampini output_file: output/ex66_1.out 369c4fb69feSStefano Zampini args: -n 33 -kernel 1 -dim 1 -lda 13 -ldc 11 -symm 0 -checkexpl -bgpu 0 370c4fb69feSStefano Zampini 371c4fb69feSStefano Zampini test: 372c4fb69feSStefano Zampini requires: hara cuda 373c4fb69feSStefano Zampini nsize: 1 374c4fb69feSStefano Zampini suffix: 1_cuda 375c4fb69feSStefano Zampini output_file: output/ex66_1.out 376c4fb69feSStefano Zampini args: -n {{17 33}} -kernel 1 -dim {{1 2 3}} -symm {{0 1}} -checkexpl -bgpu 1 377c4fb69feSStefano Zampini 378c4fb69feSStefano Zampini test: 379c4fb69feSStefano Zampini requires: hara cuda 380c4fb69feSStefano Zampini nsize: 1 381c4fb69feSStefano Zampini suffix: 1_cuda_ld 382c4fb69feSStefano Zampini output_file: output/ex66_1.out 383c4fb69feSStefano Zampini args: -n 33 -kernel 1 -dim 1 -lda 13 -ldc 11 -symm 0 -checkexpl -bgpu 1 384c4fb69feSStefano Zampini 385c4fb69feSStefano Zampini test: 386c4fb69feSStefano Zampini requires: hara define(PETSC_HAVE_MPI_INIT_THREAD) 387c4fb69feSStefano Zampini nsize: 2 388c4fb69feSStefano Zampini suffix: 1_par 389c4fb69feSStefano Zampini args: -n 32 -kernel 1 -dim 1 -ldc 12 -testlayout {{0 1}} -bgpu 0 -cgpu 0 390c4fb69feSStefano Zampini 391c4fb69feSStefano Zampini test: 392c4fb69feSStefano Zampini requires: hara cuda define(PETSC_HAVE_MPI_INIT_THREAD) 393c4fb69feSStefano Zampini nsize: 2 394c4fb69feSStefano Zampini suffix: 1_par_cuda 395c4fb69feSStefano Zampini args: -n 32 -kernel 1 -dim 1 -ldc 12 -testlayout {{0 1}} -bgpu {{0 1}} -cgpu {{0 1}} 396c4fb69feSStefano Zampini output_file: output/ex66_1_par.out 397c4fb69feSStefano Zampini 398c4fb69feSStefano Zampini #tests from matrix sampling (parallel or unsymmetric not supported) 399c4fb69feSStefano Zampini test: 400c4fb69feSStefano Zampini requires: hara 401c4fb69feSStefano Zampini nsize: 1 402c4fb69feSStefano Zampini suffix: 2 403c4fb69feSStefano Zampini args: -n {{17 33}} -kernel 0 -dim 2 -symm 1 -checkexpl -bgpu 0 404c4fb69feSStefano Zampini 405c4fb69feSStefano Zampini test: 406c4fb69feSStefano Zampini requires: hara cuda 407c4fb69feSStefano Zampini nsize: 1 408c4fb69feSStefano Zampini suffix: 2_cuda 409c4fb69feSStefano Zampini output_file: output/ex66_2.out 410c4fb69feSStefano Zampini args: -n {{17 33}} -kernel 0 -dim 2 -symm 1 -checkexpl -bgpu {{0 1}} -agpu {{0 1}} 411c4fb69feSStefano Zampini 412c4fb69feSStefano Zampini TEST*/ 413