1c8a8475eSBarry Smith 2af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 3c8a8475eSBarry Smith 4c8a8475eSBarry Smith typedef struct { 57cfd0b8cSBarry Smith Mat A; 6*a48507d3SJose E. Roman Mat D; /* local submatrix for diagonal part */ 77cfd0b8cSBarry Smith Vec w,left,right,leftwork,rightwork; 8b20f02adSBarry Smith PetscScalar scale; 9c8a8475eSBarry Smith } Mat_Normal; 10c8a8475eSBarry Smith 1169ef1043SBarry Smith PetscErrorCode MatScale_Normal(Mat inA,PetscScalar scale) 1269ef1043SBarry Smith { 1369ef1043SBarry Smith Mat_Normal *a = (Mat_Normal*)inA->data; 1469ef1043SBarry Smith 1569ef1043SBarry Smith PetscFunctionBegin; 1669ef1043SBarry Smith a->scale *= scale; 1769ef1043SBarry Smith PetscFunctionReturn(0); 1869ef1043SBarry Smith } 1969ef1043SBarry Smith 207cfd0b8cSBarry Smith PetscErrorCode MatDiagonalScale_Normal(Mat inA,Vec left,Vec right) 217cfd0b8cSBarry Smith { 227cfd0b8cSBarry Smith Mat_Normal *a = (Mat_Normal*)inA->data; 237cfd0b8cSBarry Smith 247cfd0b8cSBarry Smith PetscFunctionBegin; 257cfd0b8cSBarry Smith if (left) { 267cfd0b8cSBarry Smith if (!a->left) { 279566063dSJacob Faibussowitsch PetscCall(VecDuplicate(left,&a->left)); 289566063dSJacob Faibussowitsch PetscCall(VecCopy(left,a->left)); 297cfd0b8cSBarry Smith } else { 309566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(a->left,left,a->left)); 317cfd0b8cSBarry Smith } 327cfd0b8cSBarry Smith } 337cfd0b8cSBarry Smith if (right) { 347cfd0b8cSBarry Smith if (!a->right) { 359566063dSJacob Faibussowitsch PetscCall(VecDuplicate(right,&a->right)); 369566063dSJacob Faibussowitsch PetscCall(VecCopy(right,a->right)); 377cfd0b8cSBarry Smith } else { 389566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(a->right,right,a->right)); 397cfd0b8cSBarry Smith } 407cfd0b8cSBarry Smith } 417cfd0b8cSBarry Smith PetscFunctionReturn(0); 427cfd0b8cSBarry Smith } 437cfd0b8cSBarry Smith 44fa80d070SPierre Jolivet PetscErrorCode MatIncreaseOverlap_Normal(Mat A,PetscInt is_max,IS is[],PetscInt ov) 45fa80d070SPierre Jolivet { 46fa80d070SPierre Jolivet Mat_Normal *a = (Mat_Normal*)A->data; 47fa80d070SPierre Jolivet Mat pattern; 48fa80d070SPierre Jolivet 49fa80d070SPierre Jolivet PetscFunctionBegin; 5008401ef6SPierre Jolivet PetscCheck(ov >= 0,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified"); 519566063dSJacob Faibussowitsch PetscCall(MatProductCreate(a->A,a->A,NULL,&pattern)); 529566063dSJacob Faibussowitsch PetscCall(MatProductSetType(pattern,MATPRODUCT_AtB)); 539566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(pattern)); 549566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(pattern)); 559566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(pattern,is_max,is,ov)); 569566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pattern)); 57fa80d070SPierre Jolivet PetscFunctionReturn(0); 58fa80d070SPierre Jolivet } 59fa80d070SPierre Jolivet 60fa80d070SPierre Jolivet PetscErrorCode MatCreateSubMatrices_Normal(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[]) 61fa80d070SPierre Jolivet { 62fa80d070SPierre Jolivet Mat_Normal *a = (Mat_Normal*)mat->data; 63fa80d070SPierre Jolivet Mat B = a->A, *suba; 64fa80d070SPierre Jolivet IS *row; 65fa80d070SPierre Jolivet PetscInt M; 66fa80d070SPierre Jolivet 67fa80d070SPierre Jolivet PetscFunctionBegin; 68aed4548fSBarry Smith PetscCheck(!a->left && !a->right && irow == icol,PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Not implemented"); 69fa80d070SPierre Jolivet if (scall != MAT_REUSE_MATRIX) { 709566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(n,submat)); 71fa80d070SPierre Jolivet } 729566063dSJacob Faibussowitsch PetscCall(MatGetSize(B,&M,NULL)); 739566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&row)); 749566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,M,0,1,&row[0])); 759566063dSJacob Faibussowitsch PetscCall(ISSetIdentity(row[0])); 76fa80d070SPierre Jolivet for (M = 1; M < n; ++M) row[M] = row[0]; 779566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(B,n,row,icol,MAT_INITIAL_MATRIX,&suba)); 78fa80d070SPierre Jolivet for (M = 0; M < n; ++M) { 799566063dSJacob Faibussowitsch PetscCall(MatCreateNormal(suba[M],*submat+M)); 80fa80d070SPierre Jolivet ((Mat_Normal*)(*submat)[M]->data)->scale = a->scale; 81fa80d070SPierre Jolivet } 829566063dSJacob Faibussowitsch PetscCall(ISDestroy(&row[0])); 839566063dSJacob Faibussowitsch PetscCall(PetscFree(row)); 849566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(n,&suba)); 85fa80d070SPierre Jolivet PetscFunctionReturn(0); 86fa80d070SPierre Jolivet } 87fa80d070SPierre Jolivet 88fa80d070SPierre Jolivet PetscErrorCode MatPermute_Normal(Mat A,IS rowp,IS colp,Mat *B) 89fa80d070SPierre Jolivet { 90fa80d070SPierre Jolivet Mat_Normal *a = (Mat_Normal*)A->data; 91fa80d070SPierre Jolivet Mat C,Aa = a->A; 92fa80d070SPierre Jolivet IS row; 93fa80d070SPierre Jolivet 94fa80d070SPierre Jolivet PetscFunctionBegin; 9508401ef6SPierre Jolivet PetscCheck(rowp == colp,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Row permutation and column permutation must be the same"); 969566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PetscObjectComm((PetscObject)Aa),Aa->rmap->n,Aa->rmap->rstart,1,&row)); 979566063dSJacob Faibussowitsch PetscCall(ISSetIdentity(row)); 989566063dSJacob Faibussowitsch PetscCall(MatPermute(Aa,row,colp,&C)); 999566063dSJacob Faibussowitsch PetscCall(ISDestroy(&row)); 1009566063dSJacob Faibussowitsch PetscCall(MatCreateNormal(C,B)); 1019566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 102fa80d070SPierre Jolivet PetscFunctionReturn(0); 103fa80d070SPierre Jolivet } 104fa80d070SPierre Jolivet 105fa80d070SPierre Jolivet PetscErrorCode MatDuplicate_Normal(Mat A, MatDuplicateOption op, Mat *B) 106fa80d070SPierre Jolivet { 107fa80d070SPierre Jolivet Mat_Normal *a = (Mat_Normal*)A->data; 108fa80d070SPierre Jolivet Mat C; 109fa80d070SPierre Jolivet 110fa80d070SPierre Jolivet PetscFunctionBegin; 11108401ef6SPierre Jolivet PetscCheck(!a->left && !a->right,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Not implemented"); 1129566063dSJacob Faibussowitsch PetscCall(MatDuplicate(a->A,op,&C)); 1139566063dSJacob Faibussowitsch PetscCall(MatCreateNormal(C,B)); 1149566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 115fa80d070SPierre Jolivet if (op == MAT_COPY_VALUES) ((Mat_Normal*)(*B)->data)->scale = a->scale; 116fa80d070SPierre Jolivet PetscFunctionReturn(0); 117fa80d070SPierre Jolivet } 118fa80d070SPierre Jolivet 119fa80d070SPierre Jolivet PetscErrorCode MatCopy_Normal(Mat A,Mat B,MatStructure str) 120fa80d070SPierre Jolivet { 121fa80d070SPierre Jolivet Mat_Normal *a = (Mat_Normal*)A->data,*b = (Mat_Normal*)B->data; 122fa80d070SPierre Jolivet 123fa80d070SPierre Jolivet PetscFunctionBegin; 12408401ef6SPierre Jolivet PetscCheck(!a->left && !a->right,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Not implemented"); 1259566063dSJacob Faibussowitsch PetscCall(MatCopy(a->A,b->A,str)); 126fa80d070SPierre Jolivet b->scale = a->scale; 1279566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b->left)); 1289566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b->right)); 1299566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b->leftwork)); 1309566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b->rightwork)); 131fa80d070SPierre Jolivet PetscFunctionReturn(0); 132fa80d070SPierre Jolivet } 133fa80d070SPierre Jolivet 134dfbe8321SBarry Smith PetscErrorCode MatMult_Normal(Mat N,Vec x,Vec y) 135c8a8475eSBarry Smith { 136c8a8475eSBarry Smith Mat_Normal *Na = (Mat_Normal*)N->data; 1377cfd0b8cSBarry Smith Vec in; 138c8a8475eSBarry Smith 139c8a8475eSBarry Smith PetscFunctionBegin; 1407cfd0b8cSBarry Smith in = x; 1417cfd0b8cSBarry Smith if (Na->right) { 1427cfd0b8cSBarry Smith if (!Na->rightwork) { 1439566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Na->right,&Na->rightwork)); 1447cfd0b8cSBarry Smith } 1459566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Na->rightwork,Na->right,in)); 1467cfd0b8cSBarry Smith in = Na->rightwork; 1477cfd0b8cSBarry Smith } 1489566063dSJacob Faibussowitsch PetscCall(MatMult(Na->A,in,Na->w)); 1499566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(Na->A,Na->w,y)); 1507cfd0b8cSBarry Smith if (Na->left) { 1519566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(y,Na->left,y)); 1527cfd0b8cSBarry Smith } 1539566063dSJacob Faibussowitsch PetscCall(VecScale(y,Na->scale)); 154c8a8475eSBarry Smith PetscFunctionReturn(0); 155c8a8475eSBarry Smith } 156c8a8475eSBarry Smith 157db090513SMatthew Knepley PetscErrorCode MatMultAdd_Normal(Mat N,Vec v1,Vec v2,Vec v3) 158db090513SMatthew Knepley { 159db090513SMatthew Knepley Mat_Normal *Na = (Mat_Normal*)N->data; 1607cfd0b8cSBarry Smith Vec in; 161db090513SMatthew Knepley 162db090513SMatthew Knepley PetscFunctionBegin; 1637cfd0b8cSBarry Smith in = v1; 1647cfd0b8cSBarry Smith if (Na->right) { 1657cfd0b8cSBarry Smith if (!Na->rightwork) { 1669566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Na->right,&Na->rightwork)); 1677cfd0b8cSBarry Smith } 1689566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Na->rightwork,Na->right,in)); 1697cfd0b8cSBarry Smith in = Na->rightwork; 1707cfd0b8cSBarry Smith } 1719566063dSJacob Faibussowitsch PetscCall(MatMult(Na->A,in,Na->w)); 1729566063dSJacob Faibussowitsch PetscCall(VecScale(Na->w,Na->scale)); 1737cfd0b8cSBarry Smith if (Na->left) { 1749566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(Na->A,Na->w,v3)); 1759566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(v3,Na->left,v3)); 1769566063dSJacob Faibussowitsch PetscCall(VecAXPY(v3,1.0,v2)); 1777cfd0b8cSBarry Smith } else { 1789566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(Na->A,Na->w,v2,v3)); 1797cfd0b8cSBarry Smith } 180b20f02adSBarry Smith PetscFunctionReturn(0); 181b20f02adSBarry Smith } 182b20f02adSBarry Smith 183b20f02adSBarry Smith PetscErrorCode MatMultTranspose_Normal(Mat N,Vec x,Vec y) 184b20f02adSBarry Smith { 185b20f02adSBarry Smith Mat_Normal *Na = (Mat_Normal*)N->data; 1867cfd0b8cSBarry Smith Vec in; 187b20f02adSBarry Smith 188b20f02adSBarry Smith PetscFunctionBegin; 1897cfd0b8cSBarry Smith in = x; 1907cfd0b8cSBarry Smith if (Na->left) { 1917cfd0b8cSBarry Smith if (!Na->leftwork) { 1929566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Na->left,&Na->leftwork)); 1937cfd0b8cSBarry Smith } 1949566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Na->leftwork,Na->left,in)); 1957cfd0b8cSBarry Smith in = Na->leftwork; 1967cfd0b8cSBarry Smith } 1979566063dSJacob Faibussowitsch PetscCall(MatMult(Na->A,in,Na->w)); 1989566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(Na->A,Na->w,y)); 1997cfd0b8cSBarry Smith if (Na->right) { 2009566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(y,Na->right,y)); 2017cfd0b8cSBarry Smith } 2029566063dSJacob Faibussowitsch PetscCall(VecScale(y,Na->scale)); 203b20f02adSBarry Smith PetscFunctionReturn(0); 204b20f02adSBarry Smith } 205b20f02adSBarry Smith 206b20f02adSBarry Smith PetscErrorCode MatMultTransposeAdd_Normal(Mat N,Vec v1,Vec v2,Vec v3) 207b20f02adSBarry Smith { 208b20f02adSBarry Smith Mat_Normal *Na = (Mat_Normal*)N->data; 2097cfd0b8cSBarry Smith Vec in; 210b20f02adSBarry Smith 211b20f02adSBarry Smith PetscFunctionBegin; 2127cfd0b8cSBarry Smith in = v1; 2137cfd0b8cSBarry Smith if (Na->left) { 2147cfd0b8cSBarry Smith if (!Na->leftwork) { 2159566063dSJacob Faibussowitsch PetscCall(VecDuplicate(Na->left,&Na->leftwork)); 2167cfd0b8cSBarry Smith } 2179566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Na->leftwork,Na->left,in)); 2187cfd0b8cSBarry Smith in = Na->leftwork; 2197cfd0b8cSBarry Smith } 2209566063dSJacob Faibussowitsch PetscCall(MatMult(Na->A,in,Na->w)); 2219566063dSJacob Faibussowitsch PetscCall(VecScale(Na->w,Na->scale)); 2227cfd0b8cSBarry Smith if (Na->right) { 2239566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(Na->A,Na->w,v3)); 2249566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(v3,Na->right,v3)); 2259566063dSJacob Faibussowitsch PetscCall(VecAXPY(v3,1.0,v2)); 2267cfd0b8cSBarry Smith } else { 2279566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(Na->A,Na->w,v2,v3)); 2287cfd0b8cSBarry Smith } 229db090513SMatthew Knepley PetscFunctionReturn(0); 230db090513SMatthew Knepley } 231db090513SMatthew Knepley 232dfbe8321SBarry Smith PetscErrorCode MatDestroy_Normal(Mat N) 233c8a8475eSBarry Smith { 234c8a8475eSBarry Smith Mat_Normal *Na = (Mat_Normal*)N->data; 235c8a8475eSBarry Smith 236c8a8475eSBarry Smith PetscFunctionBegin; 2379566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Na->A)); 238*a48507d3SJose E. Roman PetscCall(MatDestroy(&Na->D)); 2399566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->w)); 2409566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->left)); 2419566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->right)); 2429566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->leftwork)); 2439566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Na->rightwork)); 2449566063dSJacob Faibussowitsch PetscCall(PetscFree(N->data)); 2459566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)N,"MatNormalGetMat_C",NULL)); 2469566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)N,"MatConvert_normal_seqaij_C",NULL)); 2479566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)N,"MatConvert_normal_mpiaij_C",NULL)); 2489566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)N,"MatProductSetFromOptions_normal_seqdense_C",NULL)); 2499566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)N,"MatProductSetFromOptions_normal_mpidense_C",NULL)); 2509566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)N,"MatProductSetFromOptions_normal_dense_C",NULL)); 251c8a8475eSBarry Smith PetscFunctionReturn(0); 252c8a8475eSBarry Smith } 253c8a8475eSBarry Smith 25417a6fd94SBarry Smith /* 25517a6fd94SBarry Smith Slow, nonscalable version 25617a6fd94SBarry Smith */ 257dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_Normal(Mat N,Vec v) 25817a6fd94SBarry Smith { 25917a6fd94SBarry Smith Mat_Normal *Na = (Mat_Normal*)N->data; 26017a6fd94SBarry Smith Mat A = Na->A; 261b24ad042SBarry Smith PetscInt i,j,rstart,rend,nnz; 262b24ad042SBarry Smith const PetscInt *cols; 26317a6fd94SBarry Smith PetscScalar *diag,*work,*values; 264b3cc6726SBarry Smith const PetscScalar *mvalues; 26517a6fd94SBarry Smith 26617a6fd94SBarry Smith PetscFunctionBegin; 2679566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(A->cmap->N,&diag,A->cmap->N,&work)); 2689566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(work,A->cmap->N)); 2699566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&rstart,&rend)); 27017a6fd94SBarry Smith for (i=rstart; i<rend; i++) { 2719566063dSJacob Faibussowitsch PetscCall(MatGetRow(A,i,&nnz,&cols,&mvalues)); 27217a6fd94SBarry Smith for (j=0; j<nnz; j++) { 273b3cc6726SBarry Smith work[cols[j]] += mvalues[j]*mvalues[j]; 27417a6fd94SBarry Smith } 2759566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A,i,&nnz,&cols,&mvalues)); 27617a6fd94SBarry Smith } 2771c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(work,diag,A->cmap->N,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)N))); 278d0f46423SBarry Smith rstart = N->cmap->rstart; 279d0f46423SBarry Smith rend = N->cmap->rend; 2809566063dSJacob Faibussowitsch PetscCall(VecGetArray(v,&values)); 2819566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(values,diag+rstart,rend-rstart)); 2829566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(v,&values)); 2839566063dSJacob Faibussowitsch PetscCall(PetscFree2(diag,work)); 2849566063dSJacob Faibussowitsch PetscCall(VecScale(v,Na->scale)); 28517a6fd94SBarry Smith PetscFunctionReturn(0); 28617a6fd94SBarry Smith } 287c8a8475eSBarry Smith 288*a48507d3SJose E. Roman PetscErrorCode MatGetDiagonalBlock_Normal(Mat N,Mat *D) 289*a48507d3SJose E. Roman { 290*a48507d3SJose E. Roman Mat_Normal *Na = (Mat_Normal*)N->data; 291*a48507d3SJose E. Roman Mat M,A = Na->A; 292*a48507d3SJose E. Roman 293*a48507d3SJose E. Roman PetscFunctionBegin; 294*a48507d3SJose E. Roman PetscCall(MatGetDiagonalBlock(A,&M)); 295*a48507d3SJose E. Roman PetscCall(MatCreateNormal(M,&Na->D)); 296*a48507d3SJose E. Roman *D = Na->D; 297*a48507d3SJose E. Roman PetscFunctionReturn(0); 298*a48507d3SJose E. Roman } 299*a48507d3SJose E. Roman 300fa80d070SPierre Jolivet PetscErrorCode MatNormalGetMat_Normal(Mat A,Mat *M) 301fa80d070SPierre Jolivet { 302fa80d070SPierre Jolivet Mat_Normal *Aa = (Mat_Normal*)A->data; 303fa80d070SPierre Jolivet 304fa80d070SPierre Jolivet PetscFunctionBegin; 305fa80d070SPierre Jolivet *M = Aa->A; 306fa80d070SPierre Jolivet PetscFunctionReturn(0); 307fa80d070SPierre Jolivet } 308fa80d070SPierre Jolivet 309fa80d070SPierre Jolivet /*@ 310fa80d070SPierre Jolivet MatNormalGetMat - Gets the Mat object stored inside a MATNORMAL 311fa80d070SPierre Jolivet 312fa80d070SPierre Jolivet Logically collective on Mat 313fa80d070SPierre Jolivet 314fa80d070SPierre Jolivet Input Parameter: 315fa80d070SPierre Jolivet . A - the MATNORMAL matrix 316fa80d070SPierre Jolivet 317fa80d070SPierre Jolivet Output Parameter: 318fa80d070SPierre Jolivet . M - the matrix object stored inside A 319fa80d070SPierre Jolivet 320fa80d070SPierre Jolivet Level: intermediate 321fa80d070SPierre Jolivet 322db781477SPatrick Sanan .seealso: `MatCreateNormal()` 323fa80d070SPierre Jolivet 324fa80d070SPierre Jolivet @*/ 325fa80d070SPierre Jolivet PetscErrorCode MatNormalGetMat(Mat A,Mat *M) 326fa80d070SPierre Jolivet { 327fa80d070SPierre Jolivet PetscFunctionBegin; 328fa80d070SPierre Jolivet PetscValidHeaderSpecific(A,MAT_CLASSID,1); 329fa80d070SPierre Jolivet PetscValidType(A,1); 330fa80d070SPierre Jolivet PetscValidPointer(M,2); 331cac4c232SBarry Smith PetscUseMethod(A,"MatNormalGetMat_C",(Mat,Mat*),(A,M)); 332fa80d070SPierre Jolivet PetscFunctionReturn(0); 333fa80d070SPierre Jolivet } 334fa80d070SPierre Jolivet 335fa80d070SPierre Jolivet PetscErrorCode MatConvert_Normal_AIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat) 336fa80d070SPierre Jolivet { 337fa80d070SPierre Jolivet Mat_Normal *Aa = (Mat_Normal*)A->data; 338fa80d070SPierre Jolivet Mat B; 339fa80d070SPierre Jolivet PetscInt m,n,M,N; 340fa80d070SPierre Jolivet 341fa80d070SPierre Jolivet PetscFunctionBegin; 3429566063dSJacob Faibussowitsch PetscCall(MatGetSize(A,&M,&N)); 3439566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A,&m,&n)); 344fa80d070SPierre Jolivet if (reuse == MAT_REUSE_MATRIX) { 345fa80d070SPierre Jolivet B = *newmat; 3469566063dSJacob Faibussowitsch PetscCall(MatProductReplaceMats(Aa->A,Aa->A,NULL,B)); 347fa80d070SPierre Jolivet } else { 3489566063dSJacob Faibussowitsch PetscCall(MatProductCreate(Aa->A,Aa->A,NULL,&B)); 3499566063dSJacob Faibussowitsch PetscCall(MatProductSetType(B,MATPRODUCT_AtB)); 3509566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(B)); 3519566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(B)); 3529566063dSJacob Faibussowitsch PetscCall(MatSetOption(B,MAT_SYMMETRIC,PETSC_TRUE)); 353fa80d070SPierre Jolivet } 3549566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(B)); 355fa80d070SPierre Jolivet if (reuse == MAT_INPLACE_MATRIX) { 3569566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(A,&B)); 357fa80d070SPierre Jolivet } else if (reuse == MAT_INITIAL_MATRIX) *newmat = B; 3589566063dSJacob Faibussowitsch PetscCall(MatConvert(*newmat,MATAIJ,MAT_INPLACE_MATRIX,newmat)); 359fa80d070SPierre Jolivet PetscFunctionReturn(0); 360fa80d070SPierre Jolivet } 361fa80d070SPierre Jolivet 362fa80d070SPierre Jolivet typedef struct { 363fa80d070SPierre Jolivet Mat work[2]; 364fa80d070SPierre Jolivet } Normal_Dense; 365fa80d070SPierre Jolivet 366fa80d070SPierre Jolivet PetscErrorCode MatProductNumeric_Normal_Dense(Mat C) 367fa80d070SPierre Jolivet { 368fa80d070SPierre Jolivet Mat A,B; 369fa80d070SPierre Jolivet Normal_Dense *contents; 370fa80d070SPierre Jolivet Mat_Normal *a; 371fa80d070SPierre Jolivet PetscScalar *array; 372fa80d070SPierre Jolivet 373fa80d070SPierre Jolivet PetscFunctionBegin; 374fa80d070SPierre Jolivet MatCheckProduct(C,3); 375fa80d070SPierre Jolivet A = C->product->A; 376fa80d070SPierre Jolivet a = (Mat_Normal*)A->data; 377fa80d070SPierre Jolivet B = C->product->B; 378fa80d070SPierre Jolivet contents = (Normal_Dense*)C->product->data; 37928b400f6SJacob Faibussowitsch PetscCheck(contents,PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data empty"); 380fa80d070SPierre Jolivet if (a->right) { 3819566063dSJacob Faibussowitsch PetscCall(MatCopy(B,C,SAME_NONZERO_PATTERN)); 3829566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(C,a->right,NULL)); 383fa80d070SPierre Jolivet } 3849566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(contents->work[0])); 3859566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(C,&array)); 3869566063dSJacob Faibussowitsch PetscCall(MatDensePlaceArray(contents->work[1],array)); 3879566063dSJacob Faibussowitsch PetscCall(MatProductNumeric(contents->work[1])); 3889566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite(C,&array)); 3899566063dSJacob Faibussowitsch PetscCall(MatDenseResetArray(contents->work[1])); 3909566063dSJacob Faibussowitsch PetscCall(MatSetOption(C,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE)); 3919566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 3929566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 3939566063dSJacob Faibussowitsch PetscCall(MatScale(C,a->scale)); 394fa80d070SPierre Jolivet PetscFunctionReturn(0); 395fa80d070SPierre Jolivet } 396fa80d070SPierre Jolivet 397fa80d070SPierre Jolivet PetscErrorCode MatNormal_DenseDestroy(void *ctx) 398fa80d070SPierre Jolivet { 399fa80d070SPierre Jolivet Normal_Dense *contents = (Normal_Dense*)ctx; 400fa80d070SPierre Jolivet 401fa80d070SPierre Jolivet PetscFunctionBegin; 4029566063dSJacob Faibussowitsch PetscCall(MatDestroy(contents->work)); 4039566063dSJacob Faibussowitsch PetscCall(MatDestroy(contents->work+1)); 4049566063dSJacob Faibussowitsch PetscCall(PetscFree(contents)); 405fa80d070SPierre Jolivet PetscFunctionReturn(0); 406fa80d070SPierre Jolivet } 407fa80d070SPierre Jolivet 408fa80d070SPierre Jolivet PetscErrorCode MatProductSymbolic_Normal_Dense(Mat C) 409fa80d070SPierre Jolivet { 410fa80d070SPierre Jolivet Mat A,B; 411fa80d070SPierre Jolivet Normal_Dense *contents = NULL; 412fa80d070SPierre Jolivet Mat_Normal *a; 413fa80d070SPierre Jolivet PetscScalar *array; 414fa80d070SPierre Jolivet PetscInt n,N,m,M; 415fa80d070SPierre Jolivet 416fa80d070SPierre Jolivet PetscFunctionBegin; 417fa80d070SPierre Jolivet MatCheckProduct(C,4); 41828b400f6SJacob Faibussowitsch PetscCheck(!C->product->data,PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data not empty"); 419fa80d070SPierre Jolivet A = C->product->A; 420fa80d070SPierre Jolivet a = (Mat_Normal*)A->data; 42128b400f6SJacob Faibussowitsch PetscCheck(!a->left,PetscObjectComm((PetscObject)C),PETSC_ERR_SUP,"Not implemented"); 422fa80d070SPierre Jolivet B = C->product->B; 4239566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(C,&m,&n)); 4249566063dSJacob Faibussowitsch PetscCall(MatGetSize(C,&M,&N)); 425fa80d070SPierre Jolivet if (m == PETSC_DECIDE || n == PETSC_DECIDE || M == PETSC_DECIDE || N == PETSC_DECIDE) { 4269566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(B,NULL,&n)); 4279566063dSJacob Faibussowitsch PetscCall(MatGetSize(B,NULL,&N)); 4289566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A,&m,NULL)); 4299566063dSJacob Faibussowitsch PetscCall(MatGetSize(A,&M,NULL)); 4309566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C,m,n,M,N)); 431fa80d070SPierre Jolivet } 4329566063dSJacob Faibussowitsch PetscCall(MatSetType(C,((PetscObject)B)->type_name)); 4339566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 4349566063dSJacob Faibussowitsch PetscCall(PetscNew(&contents)); 435fa80d070SPierre Jolivet C->product->data = contents; 436fa80d070SPierre Jolivet C->product->destroy = MatNormal_DenseDestroy; 437fa80d070SPierre Jolivet if (a->right) { 4389566063dSJacob Faibussowitsch PetscCall(MatProductCreate(a->A,C,NULL,contents->work)); 439fa80d070SPierre Jolivet } else { 4409566063dSJacob Faibussowitsch PetscCall(MatProductCreate(a->A,B,NULL,contents->work)); 441fa80d070SPierre Jolivet } 4429566063dSJacob Faibussowitsch PetscCall(MatProductSetType(contents->work[0],MATPRODUCT_AB)); 4439566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(contents->work[0])); 4449566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(contents->work[0])); 4459566063dSJacob Faibussowitsch PetscCall(MatProductCreate(a->A,contents->work[0],NULL,contents->work+1)); 4469566063dSJacob Faibussowitsch PetscCall(MatProductSetType(contents->work[1],MATPRODUCT_AtB)); 4479566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(contents->work[1])); 4489566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(contents->work[1])); 4499566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(C,&array)); 4509566063dSJacob Faibussowitsch PetscCall(MatSeqDenseSetPreallocation(contents->work[1],array)); 4519566063dSJacob Faibussowitsch PetscCall(MatMPIDenseSetPreallocation(contents->work[1],array)); 4529566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayWrite(C,&array)); 453fa80d070SPierre Jolivet C->ops->productnumeric = MatProductNumeric_Normal_Dense; 454fa80d070SPierre Jolivet PetscFunctionReturn(0); 455fa80d070SPierre Jolivet } 456fa80d070SPierre Jolivet 457fa80d070SPierre Jolivet PetscErrorCode MatProductSetFromOptions_Normal_Dense_AB(Mat C) 458fa80d070SPierre Jolivet { 459fa80d070SPierre Jolivet PetscFunctionBegin; 460fa80d070SPierre Jolivet C->ops->productsymbolic = MatProductSymbolic_Normal_Dense; 461fa80d070SPierre Jolivet PetscFunctionReturn(0); 462fa80d070SPierre Jolivet } 463fa80d070SPierre Jolivet 464fa80d070SPierre Jolivet PetscErrorCode MatProductSetFromOptions_Normal_Dense(Mat C) 465fa80d070SPierre Jolivet { 466fa80d070SPierre Jolivet Mat_Product *product = C->product; 467fa80d070SPierre Jolivet 468fa80d070SPierre Jolivet PetscFunctionBegin; 469fa80d070SPierre Jolivet if (product->type == MATPRODUCT_AB) { 4709566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions_Normal_Dense_AB(C)); 471fa80d070SPierre Jolivet } 472fa80d070SPierre Jolivet PetscFunctionReturn(0); 473fa80d070SPierre Jolivet } 474fa80d070SPierre Jolivet 475c8a8475eSBarry Smith /*@ 476c8a8475eSBarry Smith MatCreateNormal - Creates a new matrix object that behaves like A'*A. 477c8a8475eSBarry Smith 478c8a8475eSBarry Smith Collective on Mat 479c8a8475eSBarry Smith 480c8a8475eSBarry Smith Input Parameter: 481c8a8475eSBarry Smith . A - the (possibly rectangular) matrix 482c8a8475eSBarry Smith 483c8a8475eSBarry Smith Output Parameter: 484c8a8475eSBarry Smith . N - the matrix that represents A'*A 485c8a8475eSBarry Smith 486c30e7cdbSSatish Balay Level: intermediate 487c30e7cdbSSatish Balay 48895452b02SPatrick Sanan Notes: 48995452b02SPatrick Sanan The product A'*A is NOT actually formed! Rather the new matrix 490c8a8475eSBarry Smith object performs the matrix-vector product by first multiplying by 491c8a8475eSBarry Smith A and then A' 492c8a8475eSBarry Smith @*/ 4937087cfbeSBarry Smith PetscErrorCode MatCreateNormal(Mat A,Mat *N) 494c8a8475eSBarry Smith { 4959ca0e1b6SStefano Zampini PetscInt n,nn; 496c8a8475eSBarry Smith Mat_Normal *Na; 4979ca0e1b6SStefano Zampini VecType vtype; 498c8a8475eSBarry Smith 499c8a8475eSBarry Smith PetscFunctionBegin; 5009566063dSJacob Faibussowitsch PetscCall(MatGetSize(A,NULL,&nn)); 5019566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A,NULL,&n)); 5029566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A),N)); 5039566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*N,n,n,nn,nn)); 5049566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)*N,MATNORMAL)); 5059566063dSJacob Faibussowitsch PetscCall(PetscLayoutReference(A->cmap,&(*N)->rmap)); 5069566063dSJacob Faibussowitsch PetscCall(PetscLayoutReference(A->cmap,&(*N)->cmap)); 507c8a8475eSBarry Smith 5089566063dSJacob Faibussowitsch PetscCall(PetscNewLog(*N,&Na)); 50938f2d2fdSLisandro Dalcin (*N)->data = (void*) Na; 5109566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)A)); 511c3122656SLisandro Dalcin Na->A = A; 512b20f02adSBarry Smith Na->scale = 1.0; 513c8a8475eSBarry Smith 5149566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A,NULL,&Na->w)); 5152205254eSKarl Rupp 516c8a8475eSBarry Smith (*N)->ops->destroy = MatDestroy_Normal; 517c8a8475eSBarry Smith (*N)->ops->mult = MatMult_Normal; 518b20f02adSBarry Smith (*N)->ops->multtranspose = MatMultTranspose_Normal; 519b20f02adSBarry Smith (*N)->ops->multtransposeadd = MatMultTransposeAdd_Normal; 520db090513SMatthew Knepley (*N)->ops->multadd = MatMultAdd_Normal; 52117a6fd94SBarry Smith (*N)->ops->getdiagonal = MatGetDiagonal_Normal; 522*a48507d3SJose E. Roman (*N)->ops->getdiagonalblock = MatGetDiagonalBlock_Normal; 52369ef1043SBarry Smith (*N)->ops->scale = MatScale_Normal; 52469ef1043SBarry Smith (*N)->ops->diagonalscale = MatDiagonalScale_Normal; 525fa80d070SPierre Jolivet (*N)->ops->increaseoverlap = MatIncreaseOverlap_Normal; 526fa80d070SPierre Jolivet (*N)->ops->createsubmatrices = MatCreateSubMatrices_Normal; 527fa80d070SPierre Jolivet (*N)->ops->permute = MatPermute_Normal; 528fa80d070SPierre Jolivet (*N)->ops->duplicate = MatDuplicate_Normal; 529fa80d070SPierre Jolivet (*N)->ops->copy = MatCopy_Normal; 530c8a8475eSBarry Smith (*N)->assembled = PETSC_TRUE; 53148331cefSBarry Smith (*N)->preallocated = PETSC_TRUE; 5329ca0e1b6SStefano Zampini 5339566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)(*N),"MatNormalGetMat_C",MatNormalGetMat_Normal)); 5349566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)(*N),"MatConvert_normal_seqaij_C",MatConvert_Normal_AIJ)); 5359566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)(*N),"MatConvert_normal_mpiaij_C",MatConvert_Normal_AIJ)); 5369566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)(*N),"MatProductSetFromOptions_normal_seqdense_C",MatProductSetFromOptions_Normal_Dense)); 5379566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)(*N),"MatProductSetFromOptions_normal_mpidense_C",MatProductSetFromOptions_Normal_Dense)); 5389566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)(*N),"MatProductSetFromOptions_normal_dense_C",MatProductSetFromOptions_Normal_Dense)); 5399566063dSJacob Faibussowitsch PetscCall(MatSetOption(*N,MAT_SYMMETRIC,PETSC_TRUE)); 5409566063dSJacob Faibussowitsch PetscCall(MatGetVecType(A,&vtype)); 5419566063dSJacob Faibussowitsch PetscCall(MatSetVecType(*N,vtype)); 5429ca0e1b6SStefano Zampini #if defined(PETSC_HAVE_DEVICE) 5439566063dSJacob Faibussowitsch PetscCall(MatBindToCPU(*N,A->boundtocpu)); 5449ca0e1b6SStefano Zampini #endif 545c8a8475eSBarry Smith PetscFunctionReturn(0); 546c8a8475eSBarry Smith } 547