static char help[] = "Test MatMatMult(), MatTranspose(), MatTransposeMatMult() for Dense and Elemental matrices.\n\n"; /* Example: mpiexec -n ./ex104 -mat_type elemental */ #include int main(int argc,char **argv) { Mat A,B,C,D; PetscInt i,M=10,N=5,j,nrows,ncols,am,an,rstart,rend; PetscRandom r; PetscBool equal,Aiselemental; PetscReal fill = 1.0; IS isrows,iscols; const PetscInt *rows,*cols; PetscScalar *v,rval; #if defined(PETSC_HAVE_ELEMENTAL) PetscBool Test_MatMatMult=PETSC_TRUE; #else PetscBool Test_MatMatMult=PETSC_FALSE; #endif PetscMPIInt size; PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); PetscCall(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL)); PetscCall(PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL)); PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N)); PetscCall(MatSetType(A,MATDENSE)); PetscCall(MatSetFromOptions(A)); PetscCall(MatSetUp(A)); PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&r)); PetscCall(PetscRandomSetFromOptions(r)); /* Set local matrix entries */ PetscCall(MatGetOwnershipIS(A,&isrows,&iscols)); PetscCall(ISGetLocalSize(isrows,&nrows)); PetscCall(ISGetIndices(isrows,&rows)); PetscCall(ISGetLocalSize(iscols,&ncols)); PetscCall(ISGetIndices(iscols,&cols)); PetscCall(PetscMalloc1(nrows*ncols,&v)); for (i=0; i 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"This test requires ELEMENTAL"); #endif PetscCall(MatTranspose(A,MAT_INITIAL_MATRIX,&B)); /* B = A^T */ PetscCall(MatMatMult(B,A,MAT_INITIAL_MATRIX,fill,&C)); /* C = B*A = A^T*A */ PetscCall(MatMatMult(B,A,MAT_REUSE_MATRIX,fill,&C)); PetscCall(MatMatMultEqual(B,A,C,10,&equal)); PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"B*A*x != C*x"); /* Test MatDuplicate for matrix product */ PetscCall(MatDuplicate(C,MAT_COPY_VALUES,&D)); PetscCall(MatDestroy(&D)); PetscCall(MatDestroy(&C)); PetscCall(MatDestroy(&B)); } /* Test MatTransposeMatMult() */ if (!Aiselemental) { PetscCall(MatTransposeMatMult(A,A,MAT_INITIAL_MATRIX,fill,&D)); /* D = A^T*A */ PetscCall(MatTransposeMatMult(A,A,MAT_REUSE_MATRIX,fill,&D)); PetscCall(MatTransposeMatMultEqual(A,A,D,10,&equal)); PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"D*x != A^T*A*x"); /* Test MatDuplicate for matrix product */ PetscCall(MatDuplicate(D,MAT_COPY_VALUES,&C)); PetscCall(MatDestroy(&C)); PetscCall(MatDestroy(&D)); /* Test D*x = A^T*C*A*x, where C is in AIJ format */ PetscCall(MatGetLocalSize(A,&am,&an)); PetscCall(MatCreate(PETSC_COMM_WORLD,&C)); if (size == 1) { PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,am,am)); } else { PetscCall(MatSetSizes(C,am,am,PETSC_DECIDE,PETSC_DECIDE)); } PetscCall(MatSetFromOptions(C)); PetscCall(MatSetUp(C)); PetscCall(MatGetOwnershipRange(C,&rstart,&rend)); v[0] = 1.0; for (i=rstart; i PETSC_SMALL * scale,PetscObjectComm((PetscObject)D), PETSC_ERR_PLIB, "MatMatTransposeMult() differs between MAT_INITIAL_MATRIX and MAT_REUSE_MATRIX"); PetscCall(MatDestroy(&C)); PetscCall(MatMatTransposeMultEqual(A,B,D,10,&equal)); PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"D*x != A^T*A*x"); PetscCall(MatDestroy(&D)); PetscCall(MatDestroy(&B)); } PetscCall(MatDestroy(&A)); PetscCall(PetscFree(v)); PetscCall(PetscFinalize()); return 0; } /*TEST test: output_file: output/ex104.out test: suffix: 2 nsize: 2 output_file: output/ex104.out test: suffix: 3 nsize: 4 output_file: output/ex104.out args: -M 23 -N 31 test: suffix: 4 nsize: 4 output_file: output/ex104.out args: -M 23 -N 31 -matmattransmult_mpidense_mpidense_via cyclic test: suffix: 5 nsize: 4 output_file: output/ex104.out args: -M 23 -N 31 -matmattransmult_mpidense_mpidense_via allgatherv test: suffix: 6 args: -mat_type elemental requires: elemental output_file: output/ex104.out test: suffix: 7 nsize: 2 args: -mat_type elemental requires: elemental output_file: output/ex104.out TEST*/