18a9c020eSBarry Smith 28a9c020eSBarry Smith static char help[] = "Tests MatInvertVariableBlockEnvelope()\n\n"; 38a9c020eSBarry Smith 48a9c020eSBarry Smith #include <petscmat.h> 58a9c020eSBarry Smith extern PetscErrorCode MatIsDiagonal(Mat); 68a9c020eSBarry Smith extern PetscErrorCode BuildMatrix(const PetscInt*, PetscInt,const PetscInt *,Mat*); 78a9c020eSBarry Smith 88a9c020eSBarry Smith int main(int argc,char **argv) 98a9c020eSBarry Smith { 108a9c020eSBarry Smith Mat A,C,D,F; 118a9c020eSBarry Smith PetscInt i,j,rows[2],*parts,cnt, N = 21,nblocks,*blocksizes; 128a9c020eSBarry Smith PetscScalar values[2][2]; 138a9c020eSBarry Smith PetscReal rand; 148a9c020eSBarry Smith PetscRandom rctx; 158a9c020eSBarry Smith PetscMPIInt size; 168a9c020eSBarry Smith 178a9c020eSBarry Smith PetscCall(PetscInitialize(&argc,&argv,(char*) 0,help)); 188a9c020eSBarry Smith PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_DENSE)); 198a9c020eSBarry Smith 208a9c020eSBarry Smith PetscCall(MatCreate(PETSC_COMM_WORLD,&C)); 218a9c020eSBarry Smith PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,6,18)); 228a9c020eSBarry Smith PetscCall(MatSetFromOptions(C)); 238a9c020eSBarry Smith PetscCall(MatSetUp(C)); 248a9c020eSBarry Smith values[0][0] = 2; values[0][1] = 1; 258a9c020eSBarry Smith values[1][0] = 1; values[1][1] = 2; 268a9c020eSBarry Smith for (i=0;i<3;i++){ 278a9c020eSBarry Smith rows[0] = 2*i; rows[1] = 2*i + 1; 288a9c020eSBarry Smith PetscCall(MatSetValues(C,2,rows,2,rows,(PetscScalar*)values,INSERT_VALUES)); 298a9c020eSBarry Smith } 308a9c020eSBarry Smith PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 318a9c020eSBarry Smith PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 328a9c020eSBarry Smith PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 338a9c020eSBarry Smith 348a9c020eSBarry Smith PetscCall(MatMatTransposeMult(C,C,MAT_INITIAL_MATRIX,PETSC_DETERMINE,&A)); 358a9c020eSBarry Smith PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 368a9c020eSBarry Smith 378a9c020eSBarry Smith PetscCall(MatInvertVariableBlockEnvelope(A,MAT_INITIAL_MATRIX,&D)); 388a9c020eSBarry Smith PetscCall(MatView(D,PETSC_VIEWER_STDOUT_WORLD)); 398a9c020eSBarry Smith 408a9c020eSBarry Smith PetscCall(MatMatMult(A,D,MAT_INITIAL_MATRIX,1.0,&F)); 418a9c020eSBarry Smith PetscCall(MatView(F,PETSC_VIEWER_STDOUT_WORLD)); 428a9c020eSBarry Smith PetscCall(MatIsDiagonal(F)); 438a9c020eSBarry Smith 448a9c020eSBarry Smith PetscCall(MatDestroy(&A)); 458a9c020eSBarry Smith PetscCall(MatDestroy(&D)); 468a9c020eSBarry Smith PetscCall(MatDestroy(&C)); 478a9c020eSBarry Smith PetscCall(MatDestroy(&F)); 488a9c020eSBarry Smith 498a9c020eSBarry Smith PetscCall(PetscRandomCreate(PETSC_COMM_SELF,&rctx)); 508a9c020eSBarry Smith PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 518a9c020eSBarry Smith PetscCall(PetscMalloc1(size,&parts)); 528a9c020eSBarry Smith 538a9c020eSBarry Smith for (j=0; j<128; j++) { 548a9c020eSBarry Smith cnt = 0; 558a9c020eSBarry Smith for (i=0; i<size-1; i++) { 568a9c020eSBarry Smith PetscCall(PetscRandomGetValueReal(rctx,&rand)); 578a9c020eSBarry Smith parts[i] = (PetscInt) N*rand; 588a9c020eSBarry Smith parts[i] = PetscMin(parts[i],N-cnt); 598a9c020eSBarry Smith cnt += parts[i]; 608a9c020eSBarry Smith } 618a9c020eSBarry Smith parts[size-1] = N - cnt; 628a9c020eSBarry Smith 638a9c020eSBarry Smith PetscCall(PetscRandomGetValueReal(rctx,&rand)); 648a9c020eSBarry Smith nblocks = rand*10; 658a9c020eSBarry Smith nblocks = PetscMax(nblocks,2); 668a9c020eSBarry Smith cnt = 0; 678a9c020eSBarry Smith PetscCall(PetscMalloc1(nblocks,&blocksizes)); 688a9c020eSBarry Smith for (i=0; i<nblocks-1; i++) { 698a9c020eSBarry Smith PetscCall(PetscRandomGetValueReal(rctx,&rand)); 708a9c020eSBarry Smith blocksizes[i] = PetscMax(1,(PetscInt) N*rand); 718a9c020eSBarry Smith blocksizes[i] = PetscMin(blocksizes[i],N-cnt); 728a9c020eSBarry Smith cnt += blocksizes[i]; 738a9c020eSBarry Smith if (cnt == N) { 748a9c020eSBarry Smith nblocks = i + 1; 758a9c020eSBarry Smith break; 768a9c020eSBarry Smith } 778a9c020eSBarry Smith } 788a9c020eSBarry Smith if (cnt < N) { 798a9c020eSBarry Smith blocksizes[nblocks-1] = N - cnt; 808a9c020eSBarry Smith } 818a9c020eSBarry Smith 828a9c020eSBarry Smith PetscCall(BuildMatrix(parts,nblocks,blocksizes,&A)); 838a9c020eSBarry Smith PetscCall(PetscFree(blocksizes)); 848a9c020eSBarry Smith 858a9c020eSBarry Smith PetscCall(MatInvertVariableBlockEnvelope(A,MAT_INITIAL_MATRIX,&D)); 868a9c020eSBarry Smith 878a9c020eSBarry Smith PetscCall(MatMatMult(A,D,MAT_INITIAL_MATRIX,1.0,&F)); 888a9c020eSBarry Smith PetscCall(MatIsDiagonal(F)); 898a9c020eSBarry Smith 908a9c020eSBarry Smith PetscCall(MatDestroy(&A)); 918a9c020eSBarry Smith PetscCall(MatDestroy(&D)); 928a9c020eSBarry Smith PetscCall(MatDestroy(&F)); 938a9c020eSBarry Smith } 948a9c020eSBarry Smith PetscCall(PetscFree(parts)); 958a9c020eSBarry Smith PetscCall(PetscRandomDestroy(&rctx)); 968a9c020eSBarry Smith 978a9c020eSBarry Smith PetscCall(PetscFinalize()); 988a9c020eSBarry Smith return 0; 998a9c020eSBarry Smith } 1008a9c020eSBarry Smith 1018a9c020eSBarry Smith PetscErrorCode MatIsDiagonal(Mat A) 1028a9c020eSBarry Smith { 1038a9c020eSBarry Smith PetscInt ncols,i,j,rstart,rend; 1048a9c020eSBarry Smith const PetscInt *cols; 1058a9c020eSBarry Smith const PetscScalar *vals; 1068a9c020eSBarry Smith PetscBool founddiag; 1078a9c020eSBarry Smith 1088a9c020eSBarry Smith PetscFunctionBeginUser; 1098a9c020eSBarry Smith PetscCall(MatGetOwnershipRange(A,&rstart,&rend)); 1108a9c020eSBarry Smith for (i=rstart; i<rend; i++) { 1118a9c020eSBarry Smith founddiag = PETSC_FALSE; 1128a9c020eSBarry Smith PetscCall(MatGetRow(A,i,&ncols,&cols,&vals)); 1138a9c020eSBarry Smith for (j=0; j<ncols; j++) { 1148a9c020eSBarry Smith if (cols[j] == i) { 1158a9c020eSBarry 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])); 1168a9c020eSBarry Smith founddiag = PETSC_TRUE; 1178a9c020eSBarry Smith } else { 1188a9c020eSBarry 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]); 1198a9c020eSBarry Smith } 1208a9c020eSBarry Smith } 1218a9c020eSBarry Smith PetscCall(MatRestoreRow(A,i,&ncols,&cols,&vals)); 1228a9c020eSBarry Smith PetscCheck(founddiag,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Row %" PetscInt_FMT " does not have diagonal entrie",i); 1238a9c020eSBarry Smith } 1248a9c020eSBarry Smith PetscFunctionReturn(0); 1258a9c020eSBarry Smith } 1268a9c020eSBarry Smith 1278a9c020eSBarry Smith /* 1288a9c020eSBarry Smith All processes receive all the block information 1298a9c020eSBarry Smith */ 1308a9c020eSBarry Smith PetscErrorCode BuildMatrix(const PetscInt *parts, PetscInt nblocks, const PetscInt *blocksizes,Mat *A) 1318a9c020eSBarry Smith { 1328a9c020eSBarry Smith PetscInt i,cnt = 0; 1338a9c020eSBarry Smith PetscMPIInt rank; 1348a9c020eSBarry Smith 1358a9c020eSBarry Smith PetscFunctionBeginUser; 1368a9c020eSBarry Smith PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 1378a9c020eSBarry Smith PetscCall(MatCreateAIJ(PETSC_COMM_WORLD,parts[rank],parts[rank],PETSC_DETERMINE,PETSC_DETERMINE,0,NULL,0,NULL,A)); 1388a9c020eSBarry Smith PetscCall(MatSetOption(*A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 139*c5853193SPierre Jolivet if (rank == 0) { 1408a9c020eSBarry Smith for (i=0; i<nblocks; i++) { 1418a9c020eSBarry Smith PetscCall(MatSetValue(*A,cnt,cnt+blocksizes[i]-1,1.0,INSERT_VALUES)); 1428a9c020eSBarry Smith PetscCall(MatSetValue(*A,cnt+blocksizes[i]-1,cnt,1.0,INSERT_VALUES)); 1438a9c020eSBarry Smith cnt += blocksizes[i]; 1448a9c020eSBarry Smith } 1458a9c020eSBarry Smith } 1468a9c020eSBarry Smith PetscCall(MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY)); 1478a9c020eSBarry Smith PetscCall(MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY)); 1488a9c020eSBarry Smith PetscCall(MatShift(*A,10)); 1498a9c020eSBarry Smith PetscFunctionReturn(0); 1508a9c020eSBarry Smith } 1518a9c020eSBarry Smith 1528a9c020eSBarry Smith /*TEST 1538a9c020eSBarry Smith 1548a9c020eSBarry Smith test: 1558a9c020eSBarry Smith 1568a9c020eSBarry Smith test: 1578a9c020eSBarry Smith suffix: 2 1588a9c020eSBarry Smith nsize: 2 1598a9c020eSBarry Smith 1608a9c020eSBarry Smith test: 1618a9c020eSBarry Smith suffix: 5 1628a9c020eSBarry Smith nsize: 5 1638a9c020eSBarry Smith 1648a9c020eSBarry Smith TEST*/ 165