1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MatCreateLRC()\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /*T 5c4762a1bSJed Brown Concepts: Low rank correction 6c4762a1bSJed Brown 7c4762a1bSJed Brown Processors: n 8c4762a1bSJed Brown T*/ 9c4762a1bSJed Brown 10c4762a1bSJed Brown #include <petscmat.h> 11c4762a1bSJed Brown 12c4762a1bSJed Brown int main(int argc,char **args) 13c4762a1bSJed Brown { 14c4762a1bSJed Brown Vec x,b,c=NULL; 15d751fb92SStefano Zampini Mat A,U,V,LR,X,LRe; 16d751fb92SStefano Zampini PetscInt M = 5, N = 7; 17c4762a1bSJed Brown PetscErrorCode ierr; 18c4762a1bSJed Brown PetscBool flg; 19c4762a1bSJed Brown 20c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 21*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&M,NULL)); 22*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL)); 23c4762a1bSJed Brown 24c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 25c4762a1bSJed Brown Create the sparse matrix 26c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 27*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 28*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N)); 29*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOptionsPrefix(A,"A_")); 30*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 31*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJSetPreallocation(A,5,NULL)); 32*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMPIAIJSetPreallocation(A,5,NULL,5,NULL)); 33*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(A)); 34*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetRandom(A,NULL)); 35c4762a1bSJed Brown 36c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 37c4762a1bSJed Brown Create the dense matrices 38c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 39*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&U)); 40*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(U,PETSC_DECIDE,PETSC_DECIDE,M,3)); 41*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(U,MATDENSE)); 42*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOptionsPrefix(U,"U_")); 43*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(U)); 44*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(U)); 45*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetRandom(U,NULL)); 46c4762a1bSJed Brown 47*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&V)); 48*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(V,PETSC_DECIDE,PETSC_DECIDE,N,3)); 49*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(V,MATDENSE)); 50*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOptionsPrefix(V,"V_")); 51*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(V)); 52*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(V)); 53*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetRandom(V,NULL)); 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 56c4762a1bSJed Brown Create a vector to hold the diagonal of C 57d751fb92SStefano Zampini A sequential vector can be created as well on each process 58d751fb92SStefano Zampini It is user responsibility to ensure the data in the vector 59d751fb92SStefano Zampini is consistent across processors 60c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 61*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-use_c",&flg)); 62c4762a1bSJed Brown if (flg) { 63*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPI(PETSC_COMM_WORLD,PETSC_DECIDE,3,&c)); 64*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(c,NULL)); 65c4762a1bSJed Brown } 66c4762a1bSJed Brown 67c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 68c4762a1bSJed Brown Create low rank correction matrix 69c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 70*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-low_rank",&flg)); 71c4762a1bSJed Brown if (flg) { 72c4762a1bSJed Brown /* create a low-rank matrix, with no A-matrix */ 73*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateLRC(NULL,U,c,V,&LR)); 74*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 75c4762a1bSJed Brown } else { 76*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateLRC(A,U,c,V,&LR)); 77c4762a1bSJed Brown } 78c4762a1bSJed Brown 79c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 80d751fb92SStefano Zampini Create the low rank correction matrix explicitly to check for 81d751fb92SStefano Zampini correctness 82d751fb92SStefano Zampini - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 83*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatHermitianTranspose(V,MAT_INITIAL_MATRIX,&X)); 84*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(X,c,NULL)); 85*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(U,X,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&LRe)); 86*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&X)); 87d751fb92SStefano Zampini if (A) { 88*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAYPX(LRe,1.0,A,DIFFERENT_NONZERO_PATTERN)); 89d751fb92SStefano Zampini } 90d751fb92SStefano Zampini 91d751fb92SStefano Zampini /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 92c4762a1bSJed Brown Create test vectors 93c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 94*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(LR,&x,&b)); 95*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(x,NULL)); 96*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(LR,x,b)); 97*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(LR,b,x)); 98*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 99*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&b)); 100d751fb92SStefano Zampini 101d751fb92SStefano Zampini /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 102d751fb92SStefano Zampini Check correctness 103d751fb92SStefano Zampini - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 104*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(LR,LRe,10,&flg)); 105*5f80ce2aSJacob Faibussowitsch if (!flg) CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Error in MatMult\n")); 106d751fb92SStefano Zampini #if !defined(PETSC_USE_COMPLEX) 107*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultHermitianTransposeEqual(LR,LRe,10,&flg)); 108*5f80ce2aSJacob Faibussowitsch if (!flg) CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Error in MatMultTranspose\n")); 109d751fb92SStefano Zampini #endif 110d751fb92SStefano Zampini 111*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 112*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&LRe)); 113*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&U)); 114*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&V)); 115*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&c)); 116*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&LR)); 117c4762a1bSJed Brown 118c4762a1bSJed Brown /* 119c4762a1bSJed Brown Always call PetscFinalize() before exiting a program. This routine 120c4762a1bSJed Brown - finalizes the PETSc libraries as well as MPI 121c4762a1bSJed Brown - provides summary and diagnostic information if certain runtime 122c4762a1bSJed Brown options are chosen (e.g., -log_view). 123c4762a1bSJed Brown */ 124c4762a1bSJed Brown ierr = PetscFinalize(); 125c4762a1bSJed Brown return ierr; 126c4762a1bSJed Brown } 127c4762a1bSJed Brown 128c4762a1bSJed Brown /*TEST 129c4762a1bSJed Brown 130d751fb92SStefano Zampini testset: 131d751fb92SStefano Zampini output_file: output/ex102_1.out 132d751fb92SStefano Zampini nsize: {{1 2}} 133d751fb92SStefano Zampini args: -low_rank {{0 1}} -use_c {{0 1}} 134c4762a1bSJed Brown test: 135d751fb92SStefano Zampini suffix: standard 136c4762a1bSJed Brown test: 137d751fb92SStefano Zampini suffix: cuda 138d751fb92SStefano Zampini requires: cuda 139d751fb92SStefano Zampini args: -A_mat_type aijcusparse -U_mat_type densecuda -V_mat_type densecuda 140c4762a1bSJed Brown 141c4762a1bSJed Brown TEST*/ 142