1*8a9c020eSBarry Smith 2*8a9c020eSBarry Smith static char help[] = "Tests MatInvertVariableBlockEnvelope()\n\n"; 3*8a9c020eSBarry Smith 4*8a9c020eSBarry Smith #include <petscmat.h> 5*8a9c020eSBarry Smith extern PetscErrorCode MatIsDiagonal(Mat); 6*8a9c020eSBarry Smith extern PetscErrorCode BuildMatrix(const PetscInt*, PetscInt,const PetscInt *,Mat*); 7*8a9c020eSBarry Smith 8*8a9c020eSBarry Smith int main(int argc,char **argv) 9*8a9c020eSBarry Smith { 10*8a9c020eSBarry Smith Mat A,C,D,F; 11*8a9c020eSBarry Smith PetscInt i,j,rows[2],*parts,cnt, N = 21,nblocks,*blocksizes; 12*8a9c020eSBarry Smith PetscScalar values[2][2]; 13*8a9c020eSBarry Smith PetscReal rand; 14*8a9c020eSBarry Smith PetscRandom rctx; 15*8a9c020eSBarry Smith PetscMPIInt size; 16*8a9c020eSBarry Smith 17*8a9c020eSBarry Smith PetscCall(PetscInitialize(&argc,&argv,(char*) 0,help)); 18*8a9c020eSBarry Smith PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_DENSE)); 19*8a9c020eSBarry Smith 20*8a9c020eSBarry Smith PetscCall(MatCreate(PETSC_COMM_WORLD,&C)); 21*8a9c020eSBarry Smith PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,6,18)); 22*8a9c020eSBarry Smith PetscCall(MatSetFromOptions(C)); 23*8a9c020eSBarry Smith PetscCall(MatSetUp(C)); 24*8a9c020eSBarry Smith values[0][0] = 2; values[0][1] = 1; 25*8a9c020eSBarry Smith values[1][0] = 1; values[1][1] = 2; 26*8a9c020eSBarry Smith for (i=0;i<3;i++){ 27*8a9c020eSBarry Smith rows[0] = 2*i; rows[1] = 2*i + 1; 28*8a9c020eSBarry Smith PetscCall(MatSetValues(C,2,rows,2,rows,(PetscScalar*)values,INSERT_VALUES)); 29*8a9c020eSBarry Smith } 30*8a9c020eSBarry Smith PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 31*8a9c020eSBarry Smith PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 32*8a9c020eSBarry Smith PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 33*8a9c020eSBarry Smith 34*8a9c020eSBarry Smith PetscCall(MatMatTransposeMult(C,C,MAT_INITIAL_MATRIX,PETSC_DETERMINE,&A)); 35*8a9c020eSBarry Smith PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 36*8a9c020eSBarry Smith 37*8a9c020eSBarry Smith PetscCall(MatInvertVariableBlockEnvelope(A,MAT_INITIAL_MATRIX,&D)); 38*8a9c020eSBarry Smith PetscCall(MatView(D,PETSC_VIEWER_STDOUT_WORLD)); 39*8a9c020eSBarry Smith 40*8a9c020eSBarry Smith PetscCall(MatMatMult(A,D,MAT_INITIAL_MATRIX,1.0,&F)); 41*8a9c020eSBarry Smith PetscCall(MatView(F,PETSC_VIEWER_STDOUT_WORLD)); 42*8a9c020eSBarry Smith PetscCall(MatIsDiagonal(F)); 43*8a9c020eSBarry Smith 44*8a9c020eSBarry Smith PetscCall(MatDestroy(&A)); 45*8a9c020eSBarry Smith PetscCall(MatDestroy(&D)); 46*8a9c020eSBarry Smith PetscCall(MatDestroy(&C)); 47*8a9c020eSBarry Smith PetscCall(MatDestroy(&F)); 48*8a9c020eSBarry Smith 49*8a9c020eSBarry Smith PetscCall(PetscRandomCreate(PETSC_COMM_SELF,&rctx)); 50*8a9c020eSBarry Smith PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 51*8a9c020eSBarry Smith PetscCall(PetscMalloc1(size,&parts)); 52*8a9c020eSBarry Smith 53*8a9c020eSBarry Smith for (j=0; j<128; j++) { 54*8a9c020eSBarry Smith cnt = 0; 55*8a9c020eSBarry Smith for (i=0; i<size-1; i++) { 56*8a9c020eSBarry Smith PetscCall(PetscRandomGetValueReal(rctx,&rand)); 57*8a9c020eSBarry Smith parts[i] = (PetscInt) N*rand; 58*8a9c020eSBarry Smith parts[i] = PetscMin(parts[i],N-cnt); 59*8a9c020eSBarry Smith cnt += parts[i]; 60*8a9c020eSBarry Smith } 61*8a9c020eSBarry Smith parts[size-1] = N - cnt; 62*8a9c020eSBarry Smith 63*8a9c020eSBarry Smith PetscCall(PetscRandomGetValueReal(rctx,&rand)); 64*8a9c020eSBarry Smith nblocks = rand*10; 65*8a9c020eSBarry Smith nblocks = PetscMax(nblocks,2); 66*8a9c020eSBarry Smith cnt = 0; 67*8a9c020eSBarry Smith PetscCall(PetscMalloc1(nblocks,&blocksizes)); 68*8a9c020eSBarry Smith for (i=0; i<nblocks-1; i++) { 69*8a9c020eSBarry Smith PetscCall(PetscRandomGetValueReal(rctx,&rand)); 70*8a9c020eSBarry Smith blocksizes[i] = PetscMax(1,(PetscInt) N*rand); 71*8a9c020eSBarry Smith blocksizes[i] = PetscMin(blocksizes[i],N-cnt); 72*8a9c020eSBarry Smith cnt += blocksizes[i]; 73*8a9c020eSBarry Smith if (cnt == N) { 74*8a9c020eSBarry Smith nblocks = i + 1; 75*8a9c020eSBarry Smith break; 76*8a9c020eSBarry Smith } 77*8a9c020eSBarry Smith } 78*8a9c020eSBarry Smith if (cnt < N) { 79*8a9c020eSBarry Smith blocksizes[nblocks-1] = N - cnt; 80*8a9c020eSBarry Smith } 81*8a9c020eSBarry Smith 82*8a9c020eSBarry Smith PetscCall(BuildMatrix(parts,nblocks,blocksizes,&A)); 83*8a9c020eSBarry Smith PetscCall(PetscFree(blocksizes)); 84*8a9c020eSBarry Smith 85*8a9c020eSBarry Smith PetscCall(MatInvertVariableBlockEnvelope(A,MAT_INITIAL_MATRIX,&D)); 86*8a9c020eSBarry Smith 87*8a9c020eSBarry Smith PetscCall(MatMatMult(A,D,MAT_INITIAL_MATRIX,1.0,&F)); 88*8a9c020eSBarry Smith PetscCall(MatIsDiagonal(F)); 89*8a9c020eSBarry Smith 90*8a9c020eSBarry Smith PetscCall(MatDestroy(&A)); 91*8a9c020eSBarry Smith PetscCall(MatDestroy(&D)); 92*8a9c020eSBarry Smith PetscCall(MatDestroy(&F)); 93*8a9c020eSBarry Smith } 94*8a9c020eSBarry Smith PetscCall(PetscFree(parts)); 95*8a9c020eSBarry Smith PetscCall(PetscRandomDestroy(&rctx)); 96*8a9c020eSBarry Smith 97*8a9c020eSBarry Smith PetscCall(PetscFinalize()); 98*8a9c020eSBarry Smith return 0; 99*8a9c020eSBarry Smith } 100*8a9c020eSBarry Smith 101*8a9c020eSBarry Smith PetscErrorCode MatIsDiagonal(Mat A) 102*8a9c020eSBarry Smith { 103*8a9c020eSBarry Smith PetscInt ncols,i,j,rstart,rend; 104*8a9c020eSBarry Smith const PetscInt *cols; 105*8a9c020eSBarry Smith const PetscScalar *vals; 106*8a9c020eSBarry Smith PetscBool founddiag; 107*8a9c020eSBarry Smith 108*8a9c020eSBarry Smith PetscFunctionBeginUser; 109*8a9c020eSBarry Smith PetscCall(MatGetOwnershipRange(A,&rstart,&rend)); 110*8a9c020eSBarry Smith for (i=rstart; i<rend; i++) { 111*8a9c020eSBarry Smith founddiag = PETSC_FALSE; 112*8a9c020eSBarry Smith PetscCall(MatGetRow(A,i,&ncols,&cols,&vals)); 113*8a9c020eSBarry Smith for (j=0; j<ncols; j++) { 114*8a9c020eSBarry Smith if (cols[j] == i) { 115*8a9c020eSBarry Smith PetscCheck(PetscAbsScalar(vals[j] - 1) < PETSC_SMALL,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Row %" PetscInt_FMT " does not have 1 on the diagonal, it has %g",i,(double)PetscAbsScalar(vals[j])); 116*8a9c020eSBarry Smith founddiag = PETSC_TRUE; 117*8a9c020eSBarry Smith } else { 118*8a9c020eSBarry Smith PetscCheck(PetscAbsScalar(vals[j]) < PETSC_SMALL,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Row %" PetscInt_FMT " has off-diagonal value %g at %" PetscInt_FMT "",i,(double)PetscAbsScalar(vals[j]),cols[j]); 119*8a9c020eSBarry Smith } 120*8a9c020eSBarry Smith } 121*8a9c020eSBarry Smith PetscCall(MatRestoreRow(A,i,&ncols,&cols,&vals)); 122*8a9c020eSBarry Smith PetscCheck(founddiag,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Row %" PetscInt_FMT " does not have diagonal entrie",i); 123*8a9c020eSBarry Smith } 124*8a9c020eSBarry Smith PetscFunctionReturn(0); 125*8a9c020eSBarry Smith } 126*8a9c020eSBarry Smith 127*8a9c020eSBarry Smith /* 128*8a9c020eSBarry Smith All processes receive all the block information 129*8a9c020eSBarry Smith */ 130*8a9c020eSBarry Smith PetscErrorCode BuildMatrix(const PetscInt *parts, PetscInt nblocks, const PetscInt *blocksizes,Mat *A) 131*8a9c020eSBarry Smith { 132*8a9c020eSBarry Smith PetscInt i,cnt = 0; 133*8a9c020eSBarry Smith PetscMPIInt rank; 134*8a9c020eSBarry Smith 135*8a9c020eSBarry Smith PetscFunctionBeginUser; 136*8a9c020eSBarry Smith PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 137*8a9c020eSBarry Smith PetscCall(MatCreateAIJ(PETSC_COMM_WORLD,parts[rank],parts[rank],PETSC_DETERMINE,PETSC_DETERMINE,0,NULL,0,NULL,A)); 138*8a9c020eSBarry Smith PetscCall(MatSetOption(*A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 139*8a9c020eSBarry Smith if (!rank) { 140*8a9c020eSBarry Smith for (i=0; i<nblocks; i++) { 141*8a9c020eSBarry Smith PetscCall(MatSetValue(*A,cnt,cnt+blocksizes[i]-1,1.0,INSERT_VALUES)); 142*8a9c020eSBarry Smith PetscCall(MatSetValue(*A,cnt+blocksizes[i]-1,cnt,1.0,INSERT_VALUES)); 143*8a9c020eSBarry Smith cnt += blocksizes[i]; 144*8a9c020eSBarry Smith } 145*8a9c020eSBarry Smith } 146*8a9c020eSBarry Smith PetscCall(MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY)); 147*8a9c020eSBarry Smith PetscCall(MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY)); 148*8a9c020eSBarry Smith PetscCall(MatShift(*A,10)); 149*8a9c020eSBarry Smith PetscFunctionReturn(0); 150*8a9c020eSBarry Smith } 151*8a9c020eSBarry Smith 152*8a9c020eSBarry Smith /*TEST 153*8a9c020eSBarry Smith 154*8a9c020eSBarry Smith test: 155*8a9c020eSBarry Smith 156*8a9c020eSBarry Smith test: 157*8a9c020eSBarry Smith suffix: 2 158*8a9c020eSBarry Smith nsize: 2 159*8a9c020eSBarry Smith 160*8a9c020eSBarry Smith test: 161*8a9c020eSBarry Smith suffix: 5 162*8a9c020eSBarry Smith nsize: 5 163*8a9c020eSBarry Smith 164*8a9c020eSBarry Smith TEST*/ 165