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 PetscBool flg; 18c4762a1bSJed Brown 19*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&M,NULL)); 215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL)); 22c4762a1bSJed Brown 23c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 24c4762a1bSJed Brown Create the sparse matrix 25c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 265f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 275f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N)); 285f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOptionsPrefix(A,"A_")); 295f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 305f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJSetPreallocation(A,5,NULL)); 315f80ce2aSJacob Faibussowitsch CHKERRQ(MatMPIAIJSetPreallocation(A,5,NULL,5,NULL)); 325f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(A)); 335f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetRandom(A,NULL)); 34c4762a1bSJed Brown 35c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 36c4762a1bSJed Brown Create the dense matrices 37c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 385f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&U)); 395f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(U,PETSC_DECIDE,PETSC_DECIDE,M,3)); 405f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(U,MATDENSE)); 415f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOptionsPrefix(U,"U_")); 425f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(U)); 435f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(U)); 445f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetRandom(U,NULL)); 45c4762a1bSJed Brown 465f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&V)); 475f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(V,PETSC_DECIDE,PETSC_DECIDE,N,3)); 485f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(V,MATDENSE)); 495f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOptionsPrefix(V,"V_")); 505f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(V)); 515f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(V)); 525f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetRandom(V,NULL)); 53c4762a1bSJed Brown 54c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 55c4762a1bSJed Brown Create a vector to hold the diagonal of C 56d751fb92SStefano Zampini A sequential vector can be created as well on each process 57d751fb92SStefano Zampini It is user responsibility to ensure the data in the vector 58d751fb92SStefano Zampini is consistent across processors 59c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 605f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-use_c",&flg)); 61c4762a1bSJed Brown if (flg) { 625f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPI(PETSC_COMM_WORLD,PETSC_DECIDE,3,&c)); 635f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(c,NULL)); 64c4762a1bSJed Brown } 65c4762a1bSJed Brown 66c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 67c4762a1bSJed Brown Create low rank correction matrix 68c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 695f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-low_rank",&flg)); 70c4762a1bSJed Brown if (flg) { 71c4762a1bSJed Brown /* create a low-rank matrix, with no A-matrix */ 725f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateLRC(NULL,U,c,V,&LR)); 735f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 74c4762a1bSJed Brown } else { 755f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateLRC(A,U,c,V,&LR)); 76c4762a1bSJed Brown } 77c4762a1bSJed Brown 78c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 79d751fb92SStefano Zampini Create the low rank correction matrix explicitly to check for 80d751fb92SStefano Zampini correctness 81d751fb92SStefano Zampini - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 825f80ce2aSJacob Faibussowitsch CHKERRQ(MatHermitianTranspose(V,MAT_INITIAL_MATRIX,&X)); 835f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(X,c,NULL)); 845f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(U,X,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&LRe)); 855f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&X)); 86d751fb92SStefano Zampini if (A) { 875f80ce2aSJacob Faibussowitsch CHKERRQ(MatAYPX(LRe,1.0,A,DIFFERENT_NONZERO_PATTERN)); 88d751fb92SStefano Zampini } 89d751fb92SStefano Zampini 90d751fb92SStefano Zampini /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 91c4762a1bSJed Brown Create test vectors 92c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 935f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(LR,&x,&b)); 945f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(x,NULL)); 955f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(LR,x,b)); 965f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(LR,b,x)); 975f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 985f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&b)); 99d751fb92SStefano Zampini 100d751fb92SStefano Zampini /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 101d751fb92SStefano Zampini Check correctness 102d751fb92SStefano Zampini - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1035f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(LR,LRe,10,&flg)); 1045f80ce2aSJacob Faibussowitsch if (!flg) CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Error in MatMult\n")); 105d751fb92SStefano Zampini #if !defined(PETSC_USE_COMPLEX) 1065f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultHermitianTransposeEqual(LR,LRe,10,&flg)); 1075f80ce2aSJacob Faibussowitsch if (!flg) CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Error in MatMultTranspose\n")); 108d751fb92SStefano Zampini #endif 109d751fb92SStefano Zampini 1105f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 1115f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&LRe)); 1125f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&U)); 1135f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&V)); 1145f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&c)); 1155f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&LR)); 116c4762a1bSJed Brown 117c4762a1bSJed Brown /* 118c4762a1bSJed Brown Always call PetscFinalize() before exiting a program. This routine 119c4762a1bSJed Brown - finalizes the PETSc libraries as well as MPI 120c4762a1bSJed Brown - provides summary and diagnostic information if certain runtime 121c4762a1bSJed Brown options are chosen (e.g., -log_view). 122c4762a1bSJed Brown */ 123*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 124*b122ec5aSJacob Faibussowitsch return 0; 125c4762a1bSJed Brown } 126c4762a1bSJed Brown 127c4762a1bSJed Brown /*TEST 128c4762a1bSJed Brown 129d751fb92SStefano Zampini testset: 130d751fb92SStefano Zampini output_file: output/ex102_1.out 131d751fb92SStefano Zampini nsize: {{1 2}} 132d751fb92SStefano Zampini args: -low_rank {{0 1}} -use_c {{0 1}} 133c4762a1bSJed Brown test: 134d751fb92SStefano Zampini suffix: standard 135c4762a1bSJed Brown test: 136d751fb92SStefano Zampini suffix: cuda 137d751fb92SStefano Zampini requires: cuda 138d751fb92SStefano Zampini args: -A_mat_type aijcusparse -U_mat_type densecuda -V_mat_type densecuda 139c4762a1bSJed Brown 140c4762a1bSJed Brown TEST*/ 141