1c8a8475eSBarry Smith 2af0996ceSBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 3c8a8475eSBarry Smith 4c8a8475eSBarry Smith typedef struct { 57cfd0b8cSBarry Smith Mat A; 67cfd0b8cSBarry Smith Vec w,left,right,leftwork,rightwork; 7b20f02adSBarry Smith PetscScalar scale; 8c8a8475eSBarry Smith } Mat_Normal; 9c8a8475eSBarry Smith 1069ef1043SBarry Smith PetscErrorCode MatScale_Normal(Mat inA,PetscScalar scale) 1169ef1043SBarry Smith { 1269ef1043SBarry Smith Mat_Normal *a = (Mat_Normal*)inA->data; 1369ef1043SBarry Smith 1469ef1043SBarry Smith PetscFunctionBegin; 1569ef1043SBarry Smith a->scale *= scale; 1669ef1043SBarry Smith PetscFunctionReturn(0); 1769ef1043SBarry Smith } 1869ef1043SBarry Smith 197cfd0b8cSBarry Smith PetscErrorCode MatDiagonalScale_Normal(Mat inA,Vec left,Vec right) 207cfd0b8cSBarry Smith { 217cfd0b8cSBarry Smith Mat_Normal *a = (Mat_Normal*)inA->data; 227cfd0b8cSBarry Smith 237cfd0b8cSBarry Smith PetscFunctionBegin; 247cfd0b8cSBarry Smith if (left) { 257cfd0b8cSBarry Smith if (!a->left) { 265f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(left,&a->left)); 275f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(left,a->left)); 287cfd0b8cSBarry Smith } else { 295f80ce2aSJacob Faibussowitsch CHKERRQ(VecPointwiseMult(a->left,left,a->left)); 307cfd0b8cSBarry Smith } 317cfd0b8cSBarry Smith } 327cfd0b8cSBarry Smith if (right) { 337cfd0b8cSBarry Smith if (!a->right) { 345f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(right,&a->right)); 355f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(right,a->right)); 367cfd0b8cSBarry Smith } else { 375f80ce2aSJacob Faibussowitsch CHKERRQ(VecPointwiseMult(a->right,right,a->right)); 387cfd0b8cSBarry Smith } 397cfd0b8cSBarry Smith } 407cfd0b8cSBarry Smith PetscFunctionReturn(0); 417cfd0b8cSBarry Smith } 427cfd0b8cSBarry Smith 43fa80d070SPierre Jolivet PetscErrorCode MatIncreaseOverlap_Normal(Mat A,PetscInt is_max,IS is[],PetscInt ov) 44fa80d070SPierre Jolivet { 45fa80d070SPierre Jolivet Mat_Normal *a = (Mat_Normal*)A->data; 46fa80d070SPierre Jolivet Mat pattern; 47fa80d070SPierre Jolivet 48fa80d070SPierre Jolivet PetscFunctionBegin; 492c71b3e2SJacob Faibussowitsch PetscCheckFalse(ov < 0,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified"); 505f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductCreate(a->A,a->A,NULL,&pattern)); 515f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSetType(pattern,MATPRODUCT_AtB)); 525f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSetFromOptions(pattern)); 535f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSymbolic(pattern)); 545f80ce2aSJacob Faibussowitsch CHKERRQ(MatIncreaseOverlap(pattern,is_max,is,ov)); 555f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&pattern)); 56fa80d070SPierre Jolivet PetscFunctionReturn(0); 57fa80d070SPierre Jolivet } 58fa80d070SPierre Jolivet 59fa80d070SPierre Jolivet PetscErrorCode MatCreateSubMatrices_Normal(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[]) 60fa80d070SPierre Jolivet { 61fa80d070SPierre Jolivet Mat_Normal *a = (Mat_Normal*)mat->data; 62fa80d070SPierre Jolivet Mat B = a->A, *suba; 63fa80d070SPierre Jolivet IS *row; 64fa80d070SPierre Jolivet PetscInt M; 65fa80d070SPierre Jolivet 66fa80d070SPierre Jolivet PetscFunctionBegin; 672c71b3e2SJacob Faibussowitsch PetscCheckFalse(a->left || a->right || irow != icol,PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Not implemented"); 68fa80d070SPierre Jolivet if (scall != MAT_REUSE_MATRIX) { 695f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc1(n,submat)); 70fa80d070SPierre Jolivet } 715f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(B,&M,NULL)); 725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(n,&row)); 735f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateStride(PETSC_COMM_SELF,M,0,1,&row[0])); 745f80ce2aSJacob Faibussowitsch CHKERRQ(ISSetIdentity(row[0])); 75fa80d070SPierre Jolivet for (M = 1; M < n; ++M) row[M] = row[0]; 765f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSubMatrices(B,n,row,icol,MAT_INITIAL_MATRIX,&suba)); 77fa80d070SPierre Jolivet for (M = 0; M < n; ++M) { 785f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateNormal(suba[M],*submat+M)); 79fa80d070SPierre Jolivet ((Mat_Normal*)(*submat)[M]->data)->scale = a->scale; 80fa80d070SPierre Jolivet } 815f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&row[0])); 825f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(row)); 835f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroySubMatrices(n,&suba)); 84fa80d070SPierre Jolivet PetscFunctionReturn(0); 85fa80d070SPierre Jolivet } 86fa80d070SPierre Jolivet 87fa80d070SPierre Jolivet PetscErrorCode MatPermute_Normal(Mat A,IS rowp,IS colp,Mat *B) 88fa80d070SPierre Jolivet { 89fa80d070SPierre Jolivet Mat_Normal *a = (Mat_Normal*)A->data; 90fa80d070SPierre Jolivet Mat C,Aa = a->A; 91fa80d070SPierre Jolivet IS row; 92fa80d070SPierre Jolivet 93fa80d070SPierre Jolivet PetscFunctionBegin; 942c71b3e2SJacob Faibussowitsch PetscCheckFalse(rowp != colp,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Row permutation and column permutation must be the same"); 955f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateStride(PetscObjectComm((PetscObject)Aa),Aa->rmap->n,Aa->rmap->rstart,1,&row)); 965f80ce2aSJacob Faibussowitsch CHKERRQ(ISSetIdentity(row)); 975f80ce2aSJacob Faibussowitsch CHKERRQ(MatPermute(Aa,row,colp,&C)); 985f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&row)); 995f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateNormal(C,B)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 101fa80d070SPierre Jolivet PetscFunctionReturn(0); 102fa80d070SPierre Jolivet } 103fa80d070SPierre Jolivet 104fa80d070SPierre Jolivet PetscErrorCode MatDuplicate_Normal(Mat A, MatDuplicateOption op, Mat *B) 105fa80d070SPierre Jolivet { 106fa80d070SPierre Jolivet Mat_Normal *a = (Mat_Normal*)A->data; 107fa80d070SPierre Jolivet Mat C; 108fa80d070SPierre Jolivet 109fa80d070SPierre Jolivet PetscFunctionBegin; 1102c71b3e2SJacob Faibussowitsch PetscCheckFalse(a->left || a->right,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Not implemented"); 1115f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(a->A,op,&C)); 1125f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateNormal(C,B)); 1135f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 114fa80d070SPierre Jolivet if (op == MAT_COPY_VALUES) ((Mat_Normal*)(*B)->data)->scale = a->scale; 115fa80d070SPierre Jolivet PetscFunctionReturn(0); 116fa80d070SPierre Jolivet } 117fa80d070SPierre Jolivet 118fa80d070SPierre Jolivet PetscErrorCode MatCopy_Normal(Mat A,Mat B,MatStructure str) 119fa80d070SPierre Jolivet { 120fa80d070SPierre Jolivet Mat_Normal *a = (Mat_Normal*)A->data,*b = (Mat_Normal*)B->data; 121fa80d070SPierre Jolivet 122fa80d070SPierre Jolivet PetscFunctionBegin; 1232c71b3e2SJacob Faibussowitsch PetscCheckFalse(a->left || a->right,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Not implemented"); 1245f80ce2aSJacob Faibussowitsch CHKERRQ(MatCopy(a->A,b->A,str)); 125fa80d070SPierre Jolivet b->scale = a->scale; 1265f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&b->left)); 1275f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&b->right)); 1285f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&b->leftwork)); 1295f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&b->rightwork)); 130fa80d070SPierre Jolivet PetscFunctionReturn(0); 131fa80d070SPierre Jolivet } 132fa80d070SPierre Jolivet 133dfbe8321SBarry Smith PetscErrorCode MatMult_Normal(Mat N,Vec x,Vec y) 134c8a8475eSBarry Smith { 135c8a8475eSBarry Smith Mat_Normal *Na = (Mat_Normal*)N->data; 1367cfd0b8cSBarry Smith Vec in; 137c8a8475eSBarry Smith 138c8a8475eSBarry Smith PetscFunctionBegin; 1397cfd0b8cSBarry Smith in = x; 1407cfd0b8cSBarry Smith if (Na->right) { 1417cfd0b8cSBarry Smith if (!Na->rightwork) { 1425f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(Na->right,&Na->rightwork)); 1437cfd0b8cSBarry Smith } 1445f80ce2aSJacob Faibussowitsch CHKERRQ(VecPointwiseMult(Na->rightwork,Na->right,in)); 1457cfd0b8cSBarry Smith in = Na->rightwork; 1467cfd0b8cSBarry Smith } 1475f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(Na->A,in,Na->w)); 1485f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(Na->A,Na->w,y)); 1497cfd0b8cSBarry Smith if (Na->left) { 1505f80ce2aSJacob Faibussowitsch CHKERRQ(VecPointwiseMult(y,Na->left,y)); 1517cfd0b8cSBarry Smith } 1525f80ce2aSJacob Faibussowitsch CHKERRQ(VecScale(y,Na->scale)); 153c8a8475eSBarry Smith PetscFunctionReturn(0); 154c8a8475eSBarry Smith } 155c8a8475eSBarry Smith 156db090513SMatthew Knepley PetscErrorCode MatMultAdd_Normal(Mat N,Vec v1,Vec v2,Vec v3) 157db090513SMatthew Knepley { 158db090513SMatthew Knepley Mat_Normal *Na = (Mat_Normal*)N->data; 1597cfd0b8cSBarry Smith Vec in; 160db090513SMatthew Knepley 161db090513SMatthew Knepley PetscFunctionBegin; 1627cfd0b8cSBarry Smith in = v1; 1637cfd0b8cSBarry Smith if (Na->right) { 1647cfd0b8cSBarry Smith if (!Na->rightwork) { 1655f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(Na->right,&Na->rightwork)); 1667cfd0b8cSBarry Smith } 1675f80ce2aSJacob Faibussowitsch CHKERRQ(VecPointwiseMult(Na->rightwork,Na->right,in)); 1687cfd0b8cSBarry Smith in = Na->rightwork; 1697cfd0b8cSBarry Smith } 1705f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(Na->A,in,Na->w)); 1715f80ce2aSJacob Faibussowitsch CHKERRQ(VecScale(Na->w,Na->scale)); 1727cfd0b8cSBarry Smith if (Na->left) { 1735f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(Na->A,Na->w,v3)); 1745f80ce2aSJacob Faibussowitsch CHKERRQ(VecPointwiseMult(v3,Na->left,v3)); 1755f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(v3,1.0,v2)); 1767cfd0b8cSBarry Smith } else { 1775f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTransposeAdd(Na->A,Na->w,v2,v3)); 1787cfd0b8cSBarry Smith } 179b20f02adSBarry Smith PetscFunctionReturn(0); 180b20f02adSBarry Smith } 181b20f02adSBarry Smith 182b20f02adSBarry Smith PetscErrorCode MatMultTranspose_Normal(Mat N,Vec x,Vec y) 183b20f02adSBarry Smith { 184b20f02adSBarry Smith Mat_Normal *Na = (Mat_Normal*)N->data; 1857cfd0b8cSBarry Smith Vec in; 186b20f02adSBarry Smith 187b20f02adSBarry Smith PetscFunctionBegin; 1887cfd0b8cSBarry Smith in = x; 1897cfd0b8cSBarry Smith if (Na->left) { 1907cfd0b8cSBarry Smith if (!Na->leftwork) { 1915f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(Na->left,&Na->leftwork)); 1927cfd0b8cSBarry Smith } 1935f80ce2aSJacob Faibussowitsch CHKERRQ(VecPointwiseMult(Na->leftwork,Na->left,in)); 1947cfd0b8cSBarry Smith in = Na->leftwork; 1957cfd0b8cSBarry Smith } 1965f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(Na->A,in,Na->w)); 1975f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(Na->A,Na->w,y)); 1987cfd0b8cSBarry Smith if (Na->right) { 1995f80ce2aSJacob Faibussowitsch CHKERRQ(VecPointwiseMult(y,Na->right,y)); 2007cfd0b8cSBarry Smith } 2015f80ce2aSJacob Faibussowitsch CHKERRQ(VecScale(y,Na->scale)); 202b20f02adSBarry Smith PetscFunctionReturn(0); 203b20f02adSBarry Smith } 204b20f02adSBarry Smith 205b20f02adSBarry Smith PetscErrorCode MatMultTransposeAdd_Normal(Mat N,Vec v1,Vec v2,Vec v3) 206b20f02adSBarry Smith { 207b20f02adSBarry Smith Mat_Normal *Na = (Mat_Normal*)N->data; 2087cfd0b8cSBarry Smith Vec in; 209b20f02adSBarry Smith 210b20f02adSBarry Smith PetscFunctionBegin; 2117cfd0b8cSBarry Smith in = v1; 2127cfd0b8cSBarry Smith if (Na->left) { 2137cfd0b8cSBarry Smith if (!Na->leftwork) { 2145f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(Na->left,&Na->leftwork)); 2157cfd0b8cSBarry Smith } 2165f80ce2aSJacob Faibussowitsch CHKERRQ(VecPointwiseMult(Na->leftwork,Na->left,in)); 2177cfd0b8cSBarry Smith in = Na->leftwork; 2187cfd0b8cSBarry Smith } 2195f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(Na->A,in,Na->w)); 2205f80ce2aSJacob Faibussowitsch CHKERRQ(VecScale(Na->w,Na->scale)); 2217cfd0b8cSBarry Smith if (Na->right) { 2225f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(Na->A,Na->w,v3)); 2235f80ce2aSJacob Faibussowitsch CHKERRQ(VecPointwiseMult(v3,Na->right,v3)); 2245f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(v3,1.0,v2)); 2257cfd0b8cSBarry Smith } else { 2265f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTransposeAdd(Na->A,Na->w,v2,v3)); 2277cfd0b8cSBarry Smith } 228db090513SMatthew Knepley PetscFunctionReturn(0); 229db090513SMatthew Knepley } 230db090513SMatthew Knepley 231dfbe8321SBarry Smith PetscErrorCode MatDestroy_Normal(Mat N) 232c8a8475eSBarry Smith { 233c8a8475eSBarry Smith Mat_Normal *Na = (Mat_Normal*)N->data; 234c8a8475eSBarry Smith 235c8a8475eSBarry Smith PetscFunctionBegin; 2365f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Na->A)); 2375f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&Na->w)); 2385f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&Na->left)); 2395f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&Na->right)); 2405f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&Na->leftwork)); 2415f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&Na->rightwork)); 2425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(N->data)); 2435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)N,"MatNormalGetMat_C",NULL)); 2445f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)N,"MatConvert_normal_seqaij_C",NULL)); 2455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)N,"MatConvert_normal_mpiaij_C",NULL)); 2465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)N,"MatProductSetFromOptions_normal_seqdense_C",NULL)); 2475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)N,"MatProductSetFromOptions_normal_mpidense_C",NULL)); 2485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)N,"MatProductSetFromOptions_normal_dense_C",NULL)); 249c8a8475eSBarry Smith PetscFunctionReturn(0); 250c8a8475eSBarry Smith } 251c8a8475eSBarry Smith 25217a6fd94SBarry Smith /* 25317a6fd94SBarry Smith Slow, nonscalable version 25417a6fd94SBarry Smith */ 255dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_Normal(Mat N,Vec v) 25617a6fd94SBarry Smith { 25717a6fd94SBarry Smith Mat_Normal *Na = (Mat_Normal*)N->data; 25817a6fd94SBarry Smith Mat A = Na->A; 259b24ad042SBarry Smith PetscInt i,j,rstart,rend,nnz; 260b24ad042SBarry Smith const PetscInt *cols; 26117a6fd94SBarry Smith PetscScalar *diag,*work,*values; 262b3cc6726SBarry Smith const PetscScalar *mvalues; 26317a6fd94SBarry Smith 26417a6fd94SBarry Smith PetscFunctionBegin; 2655f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(A->cmap->N,&diag,A->cmap->N,&work)); 2665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArrayzero(work,A->cmap->N)); 2675f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(A,&rstart,&rend)); 26817a6fd94SBarry Smith for (i=rstart; i<rend; i++) { 2695f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRow(A,i,&nnz,&cols,&mvalues)); 27017a6fd94SBarry Smith for (j=0; j<nnz; j++) { 271b3cc6726SBarry Smith work[cols[j]] += mvalues[j]*mvalues[j]; 27217a6fd94SBarry Smith } 2735f80ce2aSJacob Faibussowitsch CHKERRQ(MatRestoreRow(A,i,&nnz,&cols,&mvalues)); 27417a6fd94SBarry Smith } 2755f80ce2aSJacob Faibussowitsch CHKERRMPI(MPIU_Allreduce(work,diag,A->cmap->N,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)N))); 276d0f46423SBarry Smith rstart = N->cmap->rstart; 277d0f46423SBarry Smith rend = N->cmap->rend; 2785f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(v,&values)); 2795f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(values,diag+rstart,rend-rstart)); 2805f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(v,&values)); 2815f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(diag,work)); 2825f80ce2aSJacob Faibussowitsch CHKERRQ(VecScale(v,Na->scale)); 28317a6fd94SBarry Smith PetscFunctionReturn(0); 28417a6fd94SBarry Smith } 285c8a8475eSBarry Smith 286fa80d070SPierre Jolivet PetscErrorCode MatNormalGetMat_Normal(Mat A,Mat *M) 287fa80d070SPierre Jolivet { 288fa80d070SPierre Jolivet Mat_Normal *Aa = (Mat_Normal*)A->data; 289fa80d070SPierre Jolivet 290fa80d070SPierre Jolivet PetscFunctionBegin; 291fa80d070SPierre Jolivet *M = Aa->A; 292fa80d070SPierre Jolivet PetscFunctionReturn(0); 293fa80d070SPierre Jolivet } 294fa80d070SPierre Jolivet 295fa80d070SPierre Jolivet /*@ 296fa80d070SPierre Jolivet MatNormalGetMat - Gets the Mat object stored inside a MATNORMAL 297fa80d070SPierre Jolivet 298fa80d070SPierre Jolivet Logically collective on Mat 299fa80d070SPierre Jolivet 300fa80d070SPierre Jolivet Input Parameter: 301fa80d070SPierre Jolivet . A - the MATNORMAL matrix 302fa80d070SPierre Jolivet 303fa80d070SPierre Jolivet Output Parameter: 304fa80d070SPierre Jolivet . M - the matrix object stored inside A 305fa80d070SPierre Jolivet 306fa80d070SPierre Jolivet Level: intermediate 307fa80d070SPierre Jolivet 308fa80d070SPierre Jolivet .seealso: MatCreateNormal() 309fa80d070SPierre Jolivet 310fa80d070SPierre Jolivet @*/ 311fa80d070SPierre Jolivet PetscErrorCode MatNormalGetMat(Mat A,Mat *M) 312fa80d070SPierre Jolivet { 313fa80d070SPierre Jolivet PetscFunctionBegin; 314fa80d070SPierre Jolivet PetscValidHeaderSpecific(A,MAT_CLASSID,1); 315fa80d070SPierre Jolivet PetscValidType(A,1); 316fa80d070SPierre Jolivet PetscValidPointer(M,2); 3175f80ce2aSJacob Faibussowitsch CHKERRQ(PetscUseMethod(A,"MatNormalGetMat_C",(Mat,Mat*),(A,M))); 318fa80d070SPierre Jolivet PetscFunctionReturn(0); 319fa80d070SPierre Jolivet } 320fa80d070SPierre Jolivet 321fa80d070SPierre Jolivet PetscErrorCode MatConvert_Normal_AIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat) 322fa80d070SPierre Jolivet { 323fa80d070SPierre Jolivet Mat_Normal *Aa = (Mat_Normal*)A->data; 324fa80d070SPierre Jolivet Mat B; 325fa80d070SPierre Jolivet PetscInt m,n,M,N; 326fa80d070SPierre Jolivet 327fa80d070SPierre Jolivet PetscFunctionBegin; 3285f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(A,&M,&N)); 3295f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(A,&m,&n)); 330fa80d070SPierre Jolivet if (reuse == MAT_REUSE_MATRIX) { 331fa80d070SPierre Jolivet B = *newmat; 3325f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductReplaceMats(Aa->A,Aa->A,NULL,B)); 333fa80d070SPierre Jolivet } else { 3345f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductCreate(Aa->A,Aa->A,NULL,&B)); 3355f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSetType(B,MATPRODUCT_AtB)); 3365f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSetFromOptions(B)); 3375f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSymbolic(B)); 3385f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(B,MAT_SYMMETRIC,PETSC_TRUE)); 339fa80d070SPierre Jolivet } 3405f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductNumeric(B)); 341fa80d070SPierre Jolivet if (reuse == MAT_INPLACE_MATRIX) { 3425f80ce2aSJacob Faibussowitsch CHKERRQ(MatHeaderReplace(A,&B)); 343fa80d070SPierre Jolivet } else if (reuse == MAT_INITIAL_MATRIX) *newmat = B; 3445f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(*newmat,MATAIJ,MAT_INPLACE_MATRIX,newmat)); 345fa80d070SPierre Jolivet PetscFunctionReturn(0); 346fa80d070SPierre Jolivet } 347fa80d070SPierre Jolivet 348fa80d070SPierre Jolivet typedef struct { 349fa80d070SPierre Jolivet Mat work[2]; 350fa80d070SPierre Jolivet } Normal_Dense; 351fa80d070SPierre Jolivet 352fa80d070SPierre Jolivet PetscErrorCode MatProductNumeric_Normal_Dense(Mat C) 353fa80d070SPierre Jolivet { 354fa80d070SPierre Jolivet Mat A,B; 355fa80d070SPierre Jolivet Normal_Dense *contents; 356fa80d070SPierre Jolivet Mat_Normal *a; 357fa80d070SPierre Jolivet PetscScalar *array; 358fa80d070SPierre Jolivet 359fa80d070SPierre Jolivet PetscFunctionBegin; 360fa80d070SPierre Jolivet MatCheckProduct(C,3); 361fa80d070SPierre Jolivet A = C->product->A; 362fa80d070SPierre Jolivet a = (Mat_Normal*)A->data; 363fa80d070SPierre Jolivet B = C->product->B; 364fa80d070SPierre Jolivet contents = (Normal_Dense*)C->product->data; 365*28b400f6SJacob Faibussowitsch PetscCheck(contents,PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data empty"); 366fa80d070SPierre Jolivet if (a->right) { 3675f80ce2aSJacob Faibussowitsch CHKERRQ(MatCopy(B,C,SAME_NONZERO_PATTERN)); 3685f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(C,a->right,NULL)); 369fa80d070SPierre Jolivet } 3705f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductNumeric(contents->work[0])); 3715f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetArrayWrite(C,&array)); 3725f80ce2aSJacob Faibussowitsch CHKERRQ(MatDensePlaceArray(contents->work[1],array)); 3735f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductNumeric(contents->work[1])); 3745f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreArrayWrite(C,&array)); 3755f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseResetArray(contents->work[1])); 3765f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(C,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE)); 3775f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 3785f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 3795f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(C,a->scale)); 380fa80d070SPierre Jolivet PetscFunctionReturn(0); 381fa80d070SPierre Jolivet } 382fa80d070SPierre Jolivet 383fa80d070SPierre Jolivet PetscErrorCode MatNormal_DenseDestroy(void *ctx) 384fa80d070SPierre Jolivet { 385fa80d070SPierre Jolivet Normal_Dense *contents = (Normal_Dense*)ctx; 386fa80d070SPierre Jolivet 387fa80d070SPierre Jolivet PetscFunctionBegin; 3885f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(contents->work)); 3895f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(contents->work+1)); 3905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(contents)); 391fa80d070SPierre Jolivet PetscFunctionReturn(0); 392fa80d070SPierre Jolivet } 393fa80d070SPierre Jolivet 394fa80d070SPierre Jolivet PetscErrorCode MatProductSymbolic_Normal_Dense(Mat C) 395fa80d070SPierre Jolivet { 396fa80d070SPierre Jolivet Mat A,B; 397fa80d070SPierre Jolivet Normal_Dense *contents = NULL; 398fa80d070SPierre Jolivet Mat_Normal *a; 399fa80d070SPierre Jolivet PetscScalar *array; 400fa80d070SPierre Jolivet PetscInt n,N,m,M; 401fa80d070SPierre Jolivet 402fa80d070SPierre Jolivet PetscFunctionBegin; 403fa80d070SPierre Jolivet MatCheckProduct(C,4); 404*28b400f6SJacob Faibussowitsch PetscCheck(!C->product->data,PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data not empty"); 405fa80d070SPierre Jolivet A = C->product->A; 406fa80d070SPierre Jolivet a = (Mat_Normal*)A->data; 407*28b400f6SJacob Faibussowitsch PetscCheck(!a->left,PetscObjectComm((PetscObject)C),PETSC_ERR_SUP,"Not implemented"); 408fa80d070SPierre Jolivet B = C->product->B; 4095f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(C,&m,&n)); 4105f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(C,&M,&N)); 411fa80d070SPierre Jolivet if (m == PETSC_DECIDE || n == PETSC_DECIDE || M == PETSC_DECIDE || N == PETSC_DECIDE) { 4125f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(B,NULL,&n)); 4135f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(B,NULL,&N)); 4145f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(A,&m,NULL)); 4155f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(A,&M,NULL)); 4165f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(C,m,n,M,N)); 417fa80d070SPierre Jolivet } 4185f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(C,((PetscObject)B)->type_name)); 4195f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(C)); 4205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNew(&contents)); 421fa80d070SPierre Jolivet C->product->data = contents; 422fa80d070SPierre Jolivet C->product->destroy = MatNormal_DenseDestroy; 423fa80d070SPierre Jolivet if (a->right) { 4245f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductCreate(a->A,C,NULL,contents->work)); 425fa80d070SPierre Jolivet } else { 4265f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductCreate(a->A,B,NULL,contents->work)); 427fa80d070SPierre Jolivet } 4285f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSetType(contents->work[0],MATPRODUCT_AB)); 4295f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSetFromOptions(contents->work[0])); 4305f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSymbolic(contents->work[0])); 4315f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductCreate(a->A,contents->work[0],NULL,contents->work+1)); 4325f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSetType(contents->work[1],MATPRODUCT_AtB)); 4335f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSetFromOptions(contents->work[1])); 4345f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSymbolic(contents->work[1])); 4355f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetArrayWrite(C,&array)); 4365f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqDenseSetPreallocation(contents->work[1],array)); 4375f80ce2aSJacob Faibussowitsch CHKERRQ(MatMPIDenseSetPreallocation(contents->work[1],array)); 4385f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreArrayWrite(C,&array)); 439fa80d070SPierre Jolivet C->ops->productnumeric = MatProductNumeric_Normal_Dense; 440fa80d070SPierre Jolivet PetscFunctionReturn(0); 441fa80d070SPierre Jolivet } 442fa80d070SPierre Jolivet 443fa80d070SPierre Jolivet PetscErrorCode MatProductSetFromOptions_Normal_Dense_AB(Mat C) 444fa80d070SPierre Jolivet { 445fa80d070SPierre Jolivet PetscFunctionBegin; 446fa80d070SPierre Jolivet C->ops->productsymbolic = MatProductSymbolic_Normal_Dense; 447fa80d070SPierre Jolivet PetscFunctionReturn(0); 448fa80d070SPierre Jolivet } 449fa80d070SPierre Jolivet 450fa80d070SPierre Jolivet PetscErrorCode MatProductSetFromOptions_Normal_Dense(Mat C) 451fa80d070SPierre Jolivet { 452fa80d070SPierre Jolivet Mat_Product *product = C->product; 453fa80d070SPierre Jolivet 454fa80d070SPierre Jolivet PetscFunctionBegin; 455fa80d070SPierre Jolivet if (product->type == MATPRODUCT_AB) { 4565f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductSetFromOptions_Normal_Dense_AB(C)); 457fa80d070SPierre Jolivet } 458fa80d070SPierre Jolivet PetscFunctionReturn(0); 459fa80d070SPierre Jolivet } 460fa80d070SPierre Jolivet 461c8a8475eSBarry Smith /*@ 462c8a8475eSBarry Smith MatCreateNormal - Creates a new matrix object that behaves like A'*A. 463c8a8475eSBarry Smith 464c8a8475eSBarry Smith Collective on Mat 465c8a8475eSBarry Smith 466c8a8475eSBarry Smith Input Parameter: 467c8a8475eSBarry Smith . A - the (possibly rectangular) matrix 468c8a8475eSBarry Smith 469c8a8475eSBarry Smith Output Parameter: 470c8a8475eSBarry Smith . N - the matrix that represents A'*A 471c8a8475eSBarry Smith 472c30e7cdbSSatish Balay Level: intermediate 473c30e7cdbSSatish Balay 47495452b02SPatrick Sanan Notes: 47595452b02SPatrick Sanan The product A'*A is NOT actually formed! Rather the new matrix 476c8a8475eSBarry Smith object performs the matrix-vector product by first multiplying by 477c8a8475eSBarry Smith A and then A' 478c8a8475eSBarry Smith @*/ 4797087cfbeSBarry Smith PetscErrorCode MatCreateNormal(Mat A,Mat *N) 480c8a8475eSBarry Smith { 4819ca0e1b6SStefano Zampini PetscInt n,nn; 482c8a8475eSBarry Smith Mat_Normal *Na; 4839ca0e1b6SStefano Zampini VecType vtype; 484c8a8475eSBarry Smith 485c8a8475eSBarry Smith PetscFunctionBegin; 4865f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(A,NULL,&nn)); 4875f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(A,NULL,&n)); 4885f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PetscObjectComm((PetscObject)A),N)); 4895f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(*N,n,n,nn,nn)); 4905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectChangeTypeName((PetscObject)*N,MATNORMAL)); 4915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLayoutReference(A->cmap,&(*N)->rmap)); 4925f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLayoutReference(A->cmap,&(*N)->cmap)); 493c8a8475eSBarry Smith 4945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNewLog(*N,&Na)); 49538f2d2fdSLisandro Dalcin (*N)->data = (void*) Na; 4965f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)A)); 497c3122656SLisandro Dalcin Na->A = A; 498b20f02adSBarry Smith Na->scale = 1.0; 499c8a8475eSBarry Smith 5005f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(A,NULL,&Na->w)); 5012205254eSKarl Rupp 502c8a8475eSBarry Smith (*N)->ops->destroy = MatDestroy_Normal; 503c8a8475eSBarry Smith (*N)->ops->mult = MatMult_Normal; 504b20f02adSBarry Smith (*N)->ops->multtranspose = MatMultTranspose_Normal; 505b20f02adSBarry Smith (*N)->ops->multtransposeadd = MatMultTransposeAdd_Normal; 506db090513SMatthew Knepley (*N)->ops->multadd = MatMultAdd_Normal; 50717a6fd94SBarry Smith (*N)->ops->getdiagonal = MatGetDiagonal_Normal; 50869ef1043SBarry Smith (*N)->ops->scale = MatScale_Normal; 50969ef1043SBarry Smith (*N)->ops->diagonalscale = MatDiagonalScale_Normal; 510fa80d070SPierre Jolivet (*N)->ops->increaseoverlap = MatIncreaseOverlap_Normal; 511fa80d070SPierre Jolivet (*N)->ops->createsubmatrices = MatCreateSubMatrices_Normal; 512fa80d070SPierre Jolivet (*N)->ops->permute = MatPermute_Normal; 513fa80d070SPierre Jolivet (*N)->ops->duplicate = MatDuplicate_Normal; 514fa80d070SPierre Jolivet (*N)->ops->copy = MatCopy_Normal; 515c8a8475eSBarry Smith (*N)->assembled = PETSC_TRUE; 51648331cefSBarry Smith (*N)->preallocated = PETSC_TRUE; 5179ca0e1b6SStefano Zampini 5185f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)(*N),"MatNormalGetMat_C",MatNormalGetMat_Normal)); 5195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)(*N),"MatConvert_normal_seqaij_C",MatConvert_Normal_AIJ)); 5205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)(*N),"MatConvert_normal_mpiaij_C",MatConvert_Normal_AIJ)); 5215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)(*N),"MatProductSetFromOptions_normal_seqdense_C",MatProductSetFromOptions_Normal_Dense)); 5225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)(*N),"MatProductSetFromOptions_normal_mpidense_C",MatProductSetFromOptions_Normal_Dense)); 5235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)(*N),"MatProductSetFromOptions_normal_dense_C",MatProductSetFromOptions_Normal_Dense)); 5245f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(*N,MAT_SYMMETRIC,PETSC_TRUE)); 5255f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetVecType(A,&vtype)); 5265f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetVecType(*N,vtype)); 5279ca0e1b6SStefano Zampini #if defined(PETSC_HAVE_DEVICE) 5285f80ce2aSJacob Faibussowitsch CHKERRQ(MatBindToCPU(*N,A->boundtocpu)); 5299ca0e1b6SStefano Zampini #endif 530c8a8475eSBarry Smith PetscFunctionReturn(0); 531c8a8475eSBarry Smith } 532