1c4762a1bSJed Brown static const char help[] = "Test parallel assembly of SBAIJ matrices\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown 5c4762a1bSJed Brown PetscErrorCode Assemble(MPI_Comm comm,PetscInt bs,MatType mtype) 6c4762a1bSJed Brown { 7c4762a1bSJed Brown const PetscInt rc[] = {0,1,2,3}; 8c4762a1bSJed Brown const PetscScalar vals[] = {100, 2, 3, 4, 5, 600, 7, 8, 9c4762a1bSJed Brown 9,100,11,1200,13,14,15,1600, 10c4762a1bSJed Brown 17,18,19,20,21,22,23,24, 11c4762a1bSJed Brown 25,26,27,2800,29,30,31,32, 12c4762a1bSJed Brown 33,34,35,36,37,38,39,40, 13c4762a1bSJed Brown 41,42,43,44,45,46,47,48, 14c4762a1bSJed Brown 49,50,51,52,53,54,55,56, 15c4762a1bSJed Brown 57,58,49,60,61,62,63,64}; 16c4762a1bSJed Brown Mat A; 17c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS) || defined(PETSC_HAVE_MKL_CPARDISO) 18c4762a1bSJed Brown Mat F; 1983863137SStefano Zampini MatSolverType stype = MATSOLVERPETSC; 20c4762a1bSJed Brown PetscRandom rdm; 21c4762a1bSJed Brown Vec b,x,y; 22c4762a1bSJed Brown PetscInt i,j; 23560a203cSprj- PetscReal norm2,tol=100*PETSC_SMALL; 24c4762a1bSJed Brown PetscBool issbaij; 25c4762a1bSJed Brown #endif 26c4762a1bSJed Brown PetscViewer viewer; 27c4762a1bSJed Brown 28c4762a1bSJed Brown PetscFunctionBegin; 299566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,&A)); 309566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,4*bs,4*bs)); 319566063dSJacob Faibussowitsch PetscCall(MatSetType(A,mtype)); 329566063dSJacob Faibussowitsch PetscCall(MatMPIBAIJSetPreallocation(A,bs,2,NULL,2,NULL)); 339566063dSJacob Faibussowitsch PetscCall(MatMPISBAIJSetPreallocation(A,bs,2,NULL,2,NULL)); 349566063dSJacob Faibussowitsch PetscCall(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE)); 35c4762a1bSJed Brown /* All processes contribute a global matrix */ 369566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlocked(A,4,rc,4,rc,vals,ADD_VALUES)); 379566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 389566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 399566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm,"Matrix %s(%" PetscInt_FMT ")\n",mtype,bs)); 409566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIGetStdout(comm,&viewer)); 419566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL)); 429566063dSJacob Faibussowitsch PetscCall(MatView(A,viewer)); 439566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 449566063dSJacob Faibussowitsch PetscCall(MatView(A,viewer)); 45c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS) || defined(PETSC_HAVE_MKL_CPARDISO) 469566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(mtype,MATMPISBAIJ,&issbaij)); 47c4762a1bSJed Brown if (!issbaij) { 489566063dSJacob Faibussowitsch PetscCall(MatShift(A,10)); 49c4762a1bSJed Brown } 509566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rdm)); 519566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm)); 529566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A,&x,&y)); 539566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&b)); 54c4762a1bSJed Brown for (j=0; j<2; j++) { 55c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS) 56c4762a1bSJed Brown if (j==0) stype = MATSOLVERMUMPS; 57c4762a1bSJed Brown #else 58c4762a1bSJed Brown if (j==0) continue; 59c4762a1bSJed Brown #endif 60c4762a1bSJed Brown #if defined(PETSC_HAVE_MKL_CPARDISO) 61c4762a1bSJed Brown if (j==1) stype = MATSOLVERMKL_CPARDISO; 62c4762a1bSJed Brown #else 63c4762a1bSJed Brown if (j==1) continue; 64c4762a1bSJed Brown #endif 65c4762a1bSJed Brown if (issbaij) { 669566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A,stype,MAT_FACTOR_CHOLESKY,&F)); 679566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(F,A,NULL,NULL)); 689566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(F,A,NULL)); 69c4762a1bSJed Brown } else { 709566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A,stype,MAT_FACTOR_LU,&F)); 719566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(F,A,NULL,NULL,NULL)); 729566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(F,A,NULL)); 73c4762a1bSJed Brown } 74c4762a1bSJed Brown for (i=0; i<10; i++) { 759566063dSJacob Faibussowitsch PetscCall(VecSetRandom(b,rdm)); 769566063dSJacob Faibussowitsch PetscCall(MatSolve(F,b,y)); 77c4762a1bSJed Brown /* Check the error */ 789566063dSJacob Faibussowitsch PetscCall(MatMult(A,y,x)); 799566063dSJacob Faibussowitsch PetscCall(VecAXPY(x,-1.0,b)); 809566063dSJacob Faibussowitsch PetscCall(VecNorm(x,NORM_2,&norm2)); 81c4762a1bSJed Brown if (norm2>tol) { 829566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Error:MatSolve(), norm2: %g\n",(double)norm2)); 83c4762a1bSJed Brown } 84c4762a1bSJed Brown } 859566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 86c4762a1bSJed Brown } 879566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 889566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 909566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm)); 91c4762a1bSJed Brown #endif 929566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 93c4762a1bSJed Brown PetscFunctionReturn(0); 94c4762a1bSJed Brown } 95c4762a1bSJed Brown 96c4762a1bSJed Brown int main(int argc,char *argv[]) 97c4762a1bSJed Brown { 98c4762a1bSJed Brown MPI_Comm comm; 99c4762a1bSJed Brown PetscMPIInt size; 100c4762a1bSJed Brown 1019566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,NULL,help)); 102c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 1039566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 104*be096a46SBarry Smith PetscCheck(size == 2,comm,PETSC_ERR_USER,"This example must be run with exactly two processes"); 1059566063dSJacob Faibussowitsch PetscCall(Assemble(comm,2,MATMPIBAIJ)); 1069566063dSJacob Faibussowitsch PetscCall(Assemble(comm,2,MATMPISBAIJ)); 1079566063dSJacob Faibussowitsch PetscCall(Assemble(comm,1,MATMPIBAIJ)); 1089566063dSJacob Faibussowitsch PetscCall(Assemble(comm,1,MATMPISBAIJ)); 1099566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 110b122ec5aSJacob Faibussowitsch return 0; 111c4762a1bSJed Brown } 112c4762a1bSJed Brown 113c4762a1bSJed Brown /*TEST 114c4762a1bSJed Brown 115c4762a1bSJed Brown test: 116c4762a1bSJed Brown nsize: 2 11797929ea7SJunchao Zhang args: -mat_ignore_lower_triangular 118c4762a1bSJed Brown filter: sed -e "s~mem [0-9]*~mem~g" 119c4762a1bSJed Brown 120c4762a1bSJed Brown TEST*/ 121