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 PetscErrorCode ierr; 237cfd0b8cSBarry Smith 247cfd0b8cSBarry Smith PetscFunctionBegin; 257cfd0b8cSBarry Smith if (left) { 267cfd0b8cSBarry Smith if (!a->left) { 277cfd0b8cSBarry Smith ierr = VecDuplicate(left,&a->left);CHKERRQ(ierr); 287cfd0b8cSBarry Smith ierr = VecCopy(left,a->left);CHKERRQ(ierr); 297cfd0b8cSBarry Smith } else { 307cfd0b8cSBarry Smith ierr = VecPointwiseMult(a->left,left,a->left);CHKERRQ(ierr); 317cfd0b8cSBarry Smith } 327cfd0b8cSBarry Smith } 337cfd0b8cSBarry Smith if (right) { 347cfd0b8cSBarry Smith if (!a->right) { 357cfd0b8cSBarry Smith ierr = VecDuplicate(right,&a->right);CHKERRQ(ierr); 367cfd0b8cSBarry Smith ierr = VecCopy(right,a->right);CHKERRQ(ierr); 377cfd0b8cSBarry Smith } else { 387cfd0b8cSBarry Smith ierr = VecPointwiseMult(a->right,right,a->right);CHKERRQ(ierr); 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 PetscErrorCode ierr; 49fa80d070SPierre Jolivet 50fa80d070SPierre Jolivet PetscFunctionBegin; 51*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(ov < 0,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified"); 52fa80d070SPierre Jolivet ierr = MatProductCreate(a->A,a->A,NULL,&pattern);CHKERRQ(ierr); 53fa80d070SPierre Jolivet ierr = MatProductSetType(pattern,MATPRODUCT_AtB);CHKERRQ(ierr); 54fa80d070SPierre Jolivet ierr = MatProductSetFromOptions(pattern);CHKERRQ(ierr); 55fa80d070SPierre Jolivet ierr = MatProductSymbolic(pattern);CHKERRQ(ierr); 56fa80d070SPierre Jolivet ierr = MatIncreaseOverlap(pattern,is_max,is,ov);CHKERRQ(ierr); 57fa80d070SPierre Jolivet ierr = MatDestroy(&pattern);CHKERRQ(ierr); 58fa80d070SPierre Jolivet PetscFunctionReturn(0); 59fa80d070SPierre Jolivet } 60fa80d070SPierre Jolivet 61fa80d070SPierre Jolivet PetscErrorCode MatCreateSubMatrices_Normal(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[]) 62fa80d070SPierre Jolivet { 63fa80d070SPierre Jolivet Mat_Normal *a = (Mat_Normal*)mat->data; 64fa80d070SPierre Jolivet Mat B = a->A, *suba; 65fa80d070SPierre Jolivet IS *row; 66fa80d070SPierre Jolivet PetscInt M; 67fa80d070SPierre Jolivet PetscErrorCode ierr; 68fa80d070SPierre Jolivet 69fa80d070SPierre Jolivet PetscFunctionBegin; 70*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(a->left || a->right || irow != icol,PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Not implemented"); 71fa80d070SPierre Jolivet if (scall != MAT_REUSE_MATRIX) { 72fa80d070SPierre Jolivet ierr = PetscCalloc1(n,submat);CHKERRQ(ierr); 73fa80d070SPierre Jolivet } 74fa80d070SPierre Jolivet ierr = MatGetSize(B,&M,NULL);CHKERRQ(ierr); 75fa80d070SPierre Jolivet ierr = PetscMalloc1(n,&row);CHKERRQ(ierr); 76fa80d070SPierre Jolivet ierr = ISCreateStride(PETSC_COMM_SELF,M,0,1,&row[0]);CHKERRQ(ierr); 77fa80d070SPierre Jolivet ierr = ISSetIdentity(row[0]);CHKERRQ(ierr); 78fa80d070SPierre Jolivet for (M = 1; M < n; ++M) row[M] = row[0]; 79fa80d070SPierre Jolivet ierr = MatCreateSubMatrices(B,n,row,icol,MAT_INITIAL_MATRIX,&suba);CHKERRQ(ierr); 80fa80d070SPierre Jolivet for (M = 0; M < n; ++M) { 81fa80d070SPierre Jolivet ierr = MatCreateNormal(suba[M],*submat+M);CHKERRQ(ierr); 82fa80d070SPierre Jolivet ((Mat_Normal*)(*submat)[M]->data)->scale = a->scale; 83fa80d070SPierre Jolivet } 84fa80d070SPierre Jolivet ierr = ISDestroy(&row[0]);CHKERRQ(ierr); 85fa80d070SPierre Jolivet ierr = PetscFree(row);CHKERRQ(ierr); 86fa80d070SPierre Jolivet ierr = MatDestroySubMatrices(n,&suba);CHKERRQ(ierr); 87fa80d070SPierre Jolivet PetscFunctionReturn(0); 88fa80d070SPierre Jolivet } 89fa80d070SPierre Jolivet 90fa80d070SPierre Jolivet PetscErrorCode MatPermute_Normal(Mat A,IS rowp,IS colp,Mat *B) 91fa80d070SPierre Jolivet { 92fa80d070SPierre Jolivet Mat_Normal *a = (Mat_Normal*)A->data; 93fa80d070SPierre Jolivet Mat C,Aa = a->A; 94fa80d070SPierre Jolivet IS row; 95fa80d070SPierre Jolivet PetscErrorCode ierr; 96fa80d070SPierre Jolivet 97fa80d070SPierre Jolivet PetscFunctionBegin; 98*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(rowp != colp,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Row permutation and column permutation must be the same"); 99fa80d070SPierre Jolivet ierr = ISCreateStride(PetscObjectComm((PetscObject)Aa),Aa->rmap->n,Aa->rmap->rstart,1,&row);CHKERRQ(ierr); 100fa80d070SPierre Jolivet ierr = ISSetIdentity(row);CHKERRQ(ierr); 101fa80d070SPierre Jolivet ierr = MatPermute(Aa,row,colp,&C);CHKERRQ(ierr); 102fa80d070SPierre Jolivet ierr = ISDestroy(&row);CHKERRQ(ierr); 103fa80d070SPierre Jolivet ierr = MatCreateNormal(C,B);CHKERRQ(ierr); 104fa80d070SPierre Jolivet ierr = MatDestroy(&C);CHKERRQ(ierr); 105fa80d070SPierre Jolivet PetscFunctionReturn(0); 106fa80d070SPierre Jolivet } 107fa80d070SPierre Jolivet 108fa80d070SPierre Jolivet PetscErrorCode MatDuplicate_Normal(Mat A, MatDuplicateOption op, Mat *B) 109fa80d070SPierre Jolivet { 110fa80d070SPierre Jolivet Mat_Normal *a = (Mat_Normal*)A->data; 111fa80d070SPierre Jolivet Mat C; 112fa80d070SPierre Jolivet PetscErrorCode ierr; 113fa80d070SPierre Jolivet 114fa80d070SPierre Jolivet PetscFunctionBegin; 115*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(a->left || a->right,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Not implemented"); 116fa80d070SPierre Jolivet ierr = MatDuplicate(a->A,op,&C);CHKERRQ(ierr); 117fa80d070SPierre Jolivet ierr = MatCreateNormal(C,B);CHKERRQ(ierr); 118fa80d070SPierre Jolivet ierr = MatDestroy(&C);CHKERRQ(ierr); 119fa80d070SPierre Jolivet if (op == MAT_COPY_VALUES) ((Mat_Normal*)(*B)->data)->scale = a->scale; 120fa80d070SPierre Jolivet PetscFunctionReturn(0); 121fa80d070SPierre Jolivet } 122fa80d070SPierre Jolivet 123fa80d070SPierre Jolivet PetscErrorCode MatCopy_Normal(Mat A,Mat B,MatStructure str) 124fa80d070SPierre Jolivet { 125fa80d070SPierre Jolivet Mat_Normal *a = (Mat_Normal*)A->data,*b = (Mat_Normal*)B->data; 126fa80d070SPierre Jolivet PetscErrorCode ierr; 127fa80d070SPierre Jolivet 128fa80d070SPierre Jolivet PetscFunctionBegin; 129*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(a->left || a->right,PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Not implemented"); 130fa80d070SPierre Jolivet ierr = MatCopy(a->A,b->A,str);CHKERRQ(ierr); 131fa80d070SPierre Jolivet b->scale = a->scale; 132fa80d070SPierre Jolivet ierr = VecDestroy(&b->left);CHKERRQ(ierr); 133fa80d070SPierre Jolivet ierr = VecDestroy(&b->right);CHKERRQ(ierr); 134fa80d070SPierre Jolivet ierr = VecDestroy(&b->leftwork);CHKERRQ(ierr); 135fa80d070SPierre Jolivet ierr = VecDestroy(&b->rightwork);CHKERRQ(ierr); 136fa80d070SPierre Jolivet PetscFunctionReturn(0); 137fa80d070SPierre Jolivet } 138fa80d070SPierre Jolivet 139dfbe8321SBarry Smith PetscErrorCode MatMult_Normal(Mat N,Vec x,Vec y) 140c8a8475eSBarry Smith { 141c8a8475eSBarry Smith Mat_Normal *Na = (Mat_Normal*)N->data; 142dfbe8321SBarry Smith PetscErrorCode ierr; 1437cfd0b8cSBarry Smith Vec in; 144c8a8475eSBarry Smith 145c8a8475eSBarry Smith PetscFunctionBegin; 1467cfd0b8cSBarry Smith in = x; 1477cfd0b8cSBarry Smith if (Na->right) { 1487cfd0b8cSBarry Smith if (!Na->rightwork) { 1497cfd0b8cSBarry Smith ierr = VecDuplicate(Na->right,&Na->rightwork);CHKERRQ(ierr); 1507cfd0b8cSBarry Smith } 1517cfd0b8cSBarry Smith ierr = VecPointwiseMult(Na->rightwork,Na->right,in);CHKERRQ(ierr); 1527cfd0b8cSBarry Smith in = Na->rightwork; 1537cfd0b8cSBarry Smith } 1547cfd0b8cSBarry Smith ierr = MatMult(Na->A,in,Na->w);CHKERRQ(ierr); 155c8a8475eSBarry Smith ierr = MatMultTranspose(Na->A,Na->w,y);CHKERRQ(ierr); 1567cfd0b8cSBarry Smith if (Na->left) { 1577cfd0b8cSBarry Smith ierr = VecPointwiseMult(y,Na->left,y);CHKERRQ(ierr); 1587cfd0b8cSBarry Smith } 159b20f02adSBarry Smith ierr = VecScale(y,Na->scale);CHKERRQ(ierr); 160c8a8475eSBarry Smith PetscFunctionReturn(0); 161c8a8475eSBarry Smith } 162c8a8475eSBarry Smith 163db090513SMatthew Knepley PetscErrorCode MatMultAdd_Normal(Mat N,Vec v1,Vec v2,Vec v3) 164db090513SMatthew Knepley { 165db090513SMatthew Knepley Mat_Normal *Na = (Mat_Normal*)N->data; 166db090513SMatthew Knepley PetscErrorCode ierr; 1677cfd0b8cSBarry Smith Vec in; 168db090513SMatthew Knepley 169db090513SMatthew Knepley PetscFunctionBegin; 1707cfd0b8cSBarry Smith in = v1; 1717cfd0b8cSBarry Smith if (Na->right) { 1727cfd0b8cSBarry Smith if (!Na->rightwork) { 1737cfd0b8cSBarry Smith ierr = VecDuplicate(Na->right,&Na->rightwork);CHKERRQ(ierr); 1747cfd0b8cSBarry Smith } 1757cfd0b8cSBarry Smith ierr = VecPointwiseMult(Na->rightwork,Na->right,in);CHKERRQ(ierr); 1767cfd0b8cSBarry Smith in = Na->rightwork; 1777cfd0b8cSBarry Smith } 1787cfd0b8cSBarry Smith ierr = MatMult(Na->A,in,Na->w);CHKERRQ(ierr); 179b20f02adSBarry Smith ierr = VecScale(Na->w,Na->scale);CHKERRQ(ierr); 1807cfd0b8cSBarry Smith if (Na->left) { 1817cfd0b8cSBarry Smith ierr = MatMultTranspose(Na->A,Na->w,v3);CHKERRQ(ierr); 1827cfd0b8cSBarry Smith ierr = VecPointwiseMult(v3,Na->left,v3);CHKERRQ(ierr); 1837cfd0b8cSBarry Smith ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 1847cfd0b8cSBarry Smith } else { 185b20f02adSBarry Smith ierr = MatMultTransposeAdd(Na->A,Na->w,v2,v3);CHKERRQ(ierr); 1867cfd0b8cSBarry Smith } 187b20f02adSBarry Smith PetscFunctionReturn(0); 188b20f02adSBarry Smith } 189b20f02adSBarry Smith 190b20f02adSBarry Smith PetscErrorCode MatMultTranspose_Normal(Mat N,Vec x,Vec y) 191b20f02adSBarry Smith { 192b20f02adSBarry Smith Mat_Normal *Na = (Mat_Normal*)N->data; 193b20f02adSBarry Smith PetscErrorCode ierr; 1947cfd0b8cSBarry Smith Vec in; 195b20f02adSBarry Smith 196b20f02adSBarry Smith PetscFunctionBegin; 1977cfd0b8cSBarry Smith in = x; 1987cfd0b8cSBarry Smith if (Na->left) { 1997cfd0b8cSBarry Smith if (!Na->leftwork) { 2007cfd0b8cSBarry Smith ierr = VecDuplicate(Na->left,&Na->leftwork);CHKERRQ(ierr); 2017cfd0b8cSBarry Smith } 2027cfd0b8cSBarry Smith ierr = VecPointwiseMult(Na->leftwork,Na->left,in);CHKERRQ(ierr); 2037cfd0b8cSBarry Smith in = Na->leftwork; 2047cfd0b8cSBarry Smith } 2057cfd0b8cSBarry Smith ierr = MatMult(Na->A,in,Na->w);CHKERRQ(ierr); 206b20f02adSBarry Smith ierr = MatMultTranspose(Na->A,Na->w,y);CHKERRQ(ierr); 2077cfd0b8cSBarry Smith if (Na->right) { 2087cfd0b8cSBarry Smith ierr = VecPointwiseMult(y,Na->right,y);CHKERRQ(ierr); 2097cfd0b8cSBarry Smith } 210b20f02adSBarry Smith ierr = VecScale(y,Na->scale);CHKERRQ(ierr); 211b20f02adSBarry Smith PetscFunctionReturn(0); 212b20f02adSBarry Smith } 213b20f02adSBarry Smith 214b20f02adSBarry Smith PetscErrorCode MatMultTransposeAdd_Normal(Mat N,Vec v1,Vec v2,Vec v3) 215b20f02adSBarry Smith { 216b20f02adSBarry Smith Mat_Normal *Na = (Mat_Normal*)N->data; 217b20f02adSBarry Smith PetscErrorCode ierr; 2187cfd0b8cSBarry Smith Vec in; 219b20f02adSBarry Smith 220b20f02adSBarry Smith PetscFunctionBegin; 2217cfd0b8cSBarry Smith in = v1; 2227cfd0b8cSBarry Smith if (Na->left) { 2237cfd0b8cSBarry Smith if (!Na->leftwork) { 2247cfd0b8cSBarry Smith ierr = VecDuplicate(Na->left,&Na->leftwork);CHKERRQ(ierr); 2257cfd0b8cSBarry Smith } 2267cfd0b8cSBarry Smith ierr = VecPointwiseMult(Na->leftwork,Na->left,in);CHKERRQ(ierr); 2277cfd0b8cSBarry Smith in = Na->leftwork; 2287cfd0b8cSBarry Smith } 2297cfd0b8cSBarry Smith ierr = MatMult(Na->A,in,Na->w);CHKERRQ(ierr); 230b20f02adSBarry Smith ierr = VecScale(Na->w,Na->scale);CHKERRQ(ierr); 2317cfd0b8cSBarry Smith if (Na->right) { 2327cfd0b8cSBarry Smith ierr = MatMultTranspose(Na->A,Na->w,v3);CHKERRQ(ierr); 2337cfd0b8cSBarry Smith ierr = VecPointwiseMult(v3,Na->right,v3);CHKERRQ(ierr); 2347cfd0b8cSBarry Smith ierr = VecAXPY(v3,1.0,v2);CHKERRQ(ierr); 2357cfd0b8cSBarry Smith } else { 236db090513SMatthew Knepley ierr = MatMultTransposeAdd(Na->A,Na->w,v2,v3);CHKERRQ(ierr); 2377cfd0b8cSBarry Smith } 238db090513SMatthew Knepley PetscFunctionReturn(0); 239db090513SMatthew Knepley } 240db090513SMatthew Knepley 241dfbe8321SBarry Smith PetscErrorCode MatDestroy_Normal(Mat N) 242c8a8475eSBarry Smith { 243c8a8475eSBarry Smith Mat_Normal *Na = (Mat_Normal*)N->data; 244dfbe8321SBarry Smith PetscErrorCode ierr; 245c8a8475eSBarry Smith 246c8a8475eSBarry Smith PetscFunctionBegin; 2476bf464f9SBarry Smith ierr = MatDestroy(&Na->A);CHKERRQ(ierr); 2486bf464f9SBarry Smith ierr = VecDestroy(&Na->w);CHKERRQ(ierr); 2496bf464f9SBarry Smith ierr = VecDestroy(&Na->left);CHKERRQ(ierr); 2506bf464f9SBarry Smith ierr = VecDestroy(&Na->right);CHKERRQ(ierr); 2516bf464f9SBarry Smith ierr = VecDestroy(&Na->leftwork);CHKERRQ(ierr); 2526bf464f9SBarry Smith ierr = VecDestroy(&Na->rightwork);CHKERRQ(ierr); 253bf0cc555SLisandro Dalcin ierr = PetscFree(N->data);CHKERRQ(ierr); 254fa80d070SPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)N,"MatNormalGetMat_C",NULL);CHKERRQ(ierr); 255fa80d070SPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)N,"MatConvert_normal_seqaij_C",NULL);CHKERRQ(ierr); 256fa80d070SPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)N,"MatConvert_normal_mpiaij_C",NULL);CHKERRQ(ierr); 257fa80d070SPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)N,"MatProductSetFromOptions_normal_seqdense_C",NULL);CHKERRQ(ierr); 258fa80d070SPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)N,"MatProductSetFromOptions_normal_mpidense_C",NULL);CHKERRQ(ierr); 259fa80d070SPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)N,"MatProductSetFromOptions_normal_dense_C",NULL);CHKERRQ(ierr); 260c8a8475eSBarry Smith PetscFunctionReturn(0); 261c8a8475eSBarry Smith } 262c8a8475eSBarry Smith 26317a6fd94SBarry Smith /* 26417a6fd94SBarry Smith Slow, nonscalable version 26517a6fd94SBarry Smith */ 266dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_Normal(Mat N,Vec v) 26717a6fd94SBarry Smith { 26817a6fd94SBarry Smith Mat_Normal *Na = (Mat_Normal*)N->data; 26917a6fd94SBarry Smith Mat A = Na->A; 270dfbe8321SBarry Smith PetscErrorCode ierr; 271b24ad042SBarry Smith PetscInt i,j,rstart,rend,nnz; 272b24ad042SBarry Smith const PetscInt *cols; 27317a6fd94SBarry Smith PetscScalar *diag,*work,*values; 274b3cc6726SBarry Smith const PetscScalar *mvalues; 27517a6fd94SBarry Smith 27617a6fd94SBarry Smith PetscFunctionBegin; 277dcca6d9dSJed Brown ierr = PetscMalloc2(A->cmap->N,&diag,A->cmap->N,&work);CHKERRQ(ierr); 278580bdb30SBarry Smith ierr = PetscArrayzero(work,A->cmap->N);CHKERRQ(ierr); 27917a6fd94SBarry Smith ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 28017a6fd94SBarry Smith for (i=rstart; i<rend; i++) { 281b3cc6726SBarry Smith ierr = MatGetRow(A,i,&nnz,&cols,&mvalues);CHKERRQ(ierr); 28217a6fd94SBarry Smith for (j=0; j<nnz; j++) { 283b3cc6726SBarry Smith work[cols[j]] += mvalues[j]*mvalues[j]; 28417a6fd94SBarry Smith } 285b3cc6726SBarry Smith ierr = MatRestoreRow(A,i,&nnz,&cols,&mvalues);CHKERRQ(ierr); 28617a6fd94SBarry Smith } 287820f2d46SBarry Smith ierr = MPIU_Allreduce(work,diag,A->cmap->N,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)N));CHKERRMPI(ierr); 288d0f46423SBarry Smith rstart = N->cmap->rstart; 289d0f46423SBarry Smith rend = N->cmap->rend; 29017a6fd94SBarry Smith ierr = VecGetArray(v,&values);CHKERRQ(ierr); 291580bdb30SBarry Smith ierr = PetscArraycpy(values,diag+rstart,rend-rstart);CHKERRQ(ierr); 29217a6fd94SBarry Smith ierr = VecRestoreArray(v,&values);CHKERRQ(ierr); 29374ed9c26SBarry Smith ierr = PetscFree2(diag,work);CHKERRQ(ierr); 294b20f02adSBarry Smith ierr = VecScale(v,Na->scale);CHKERRQ(ierr); 29517a6fd94SBarry Smith PetscFunctionReturn(0); 29617a6fd94SBarry Smith } 297c8a8475eSBarry Smith 298fa80d070SPierre Jolivet PetscErrorCode MatNormalGetMat_Normal(Mat A,Mat *M) 299fa80d070SPierre Jolivet { 300fa80d070SPierre Jolivet Mat_Normal *Aa = (Mat_Normal*)A->data; 301fa80d070SPierre Jolivet 302fa80d070SPierre Jolivet PetscFunctionBegin; 303fa80d070SPierre Jolivet *M = Aa->A; 304fa80d070SPierre Jolivet PetscFunctionReturn(0); 305fa80d070SPierre Jolivet } 306fa80d070SPierre Jolivet 307fa80d070SPierre Jolivet /*@ 308fa80d070SPierre Jolivet MatNormalGetMat - Gets the Mat object stored inside a MATNORMAL 309fa80d070SPierre Jolivet 310fa80d070SPierre Jolivet Logically collective on Mat 311fa80d070SPierre Jolivet 312fa80d070SPierre Jolivet Input Parameter: 313fa80d070SPierre Jolivet . A - the MATNORMAL matrix 314fa80d070SPierre Jolivet 315fa80d070SPierre Jolivet Output Parameter: 316fa80d070SPierre Jolivet . M - the matrix object stored inside A 317fa80d070SPierre Jolivet 318fa80d070SPierre Jolivet Level: intermediate 319fa80d070SPierre Jolivet 320fa80d070SPierre Jolivet .seealso: MatCreateNormal() 321fa80d070SPierre Jolivet 322fa80d070SPierre Jolivet @*/ 323fa80d070SPierre Jolivet PetscErrorCode MatNormalGetMat(Mat A,Mat *M) 324fa80d070SPierre Jolivet { 325fa80d070SPierre Jolivet PetscErrorCode ierr; 326fa80d070SPierre Jolivet 327fa80d070SPierre Jolivet PetscFunctionBegin; 328fa80d070SPierre Jolivet PetscValidHeaderSpecific(A,MAT_CLASSID,1); 329fa80d070SPierre Jolivet PetscValidType(A,1); 330fa80d070SPierre Jolivet PetscValidPointer(M,2); 331fa80d070SPierre Jolivet ierr = PetscUseMethod(A,"MatNormalGetMat_C",(Mat,Mat*),(A,M));CHKERRQ(ierr); 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 PetscErrorCode ierr; 341fa80d070SPierre Jolivet 342fa80d070SPierre Jolivet PetscFunctionBegin; 343fa80d070SPierre Jolivet ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 344fa80d070SPierre Jolivet ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 345fa80d070SPierre Jolivet if (reuse == MAT_REUSE_MATRIX) { 346fa80d070SPierre Jolivet B = *newmat; 347fa80d070SPierre Jolivet ierr = MatProductReplaceMats(Aa->A,Aa->A,NULL,B);CHKERRQ(ierr); 348fa80d070SPierre Jolivet } else { 349fa80d070SPierre Jolivet ierr = MatProductCreate(Aa->A,Aa->A,NULL,&B);CHKERRQ(ierr); 350fa80d070SPierre Jolivet ierr = MatProductSetType(B,MATPRODUCT_AtB);CHKERRQ(ierr); 351fa80d070SPierre Jolivet ierr = MatProductSetFromOptions(B);CHKERRQ(ierr); 352fa80d070SPierre Jolivet ierr = MatProductSymbolic(B);CHKERRQ(ierr); 353fa80d070SPierre Jolivet ierr = MatSetOption(B,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr); 354fa80d070SPierre Jolivet } 355fa80d070SPierre Jolivet ierr = MatProductNumeric(B);CHKERRQ(ierr); 356fa80d070SPierre Jolivet if (reuse == MAT_INPLACE_MATRIX) { 357fa80d070SPierre Jolivet ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr); 358fa80d070SPierre Jolivet } else if (reuse == MAT_INITIAL_MATRIX) *newmat = B; 359fa80d070SPierre Jolivet ierr = MatConvert(*newmat,MATAIJ,MAT_INPLACE_MATRIX,newmat);CHKERRQ(ierr); 360fa80d070SPierre Jolivet PetscFunctionReturn(0); 361fa80d070SPierre Jolivet } 362fa80d070SPierre Jolivet 363fa80d070SPierre Jolivet typedef struct { 364fa80d070SPierre Jolivet Mat work[2]; 365fa80d070SPierre Jolivet } Normal_Dense; 366fa80d070SPierre Jolivet 367fa80d070SPierre Jolivet PetscErrorCode MatProductNumeric_Normal_Dense(Mat C) 368fa80d070SPierre Jolivet { 369fa80d070SPierre Jolivet Mat A,B; 370fa80d070SPierre Jolivet Normal_Dense *contents; 371fa80d070SPierre Jolivet Mat_Normal *a; 372fa80d070SPierre Jolivet PetscScalar *array; 373fa80d070SPierre Jolivet PetscErrorCode ierr; 374fa80d070SPierre Jolivet 375fa80d070SPierre Jolivet PetscFunctionBegin; 376fa80d070SPierre Jolivet MatCheckProduct(C,3); 377fa80d070SPierre Jolivet A = C->product->A; 378fa80d070SPierre Jolivet a = (Mat_Normal*)A->data; 379fa80d070SPierre Jolivet B = C->product->B; 380fa80d070SPierre Jolivet contents = (Normal_Dense*)C->product->data; 381*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!contents,PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data empty"); 382fa80d070SPierre Jolivet if (a->right) { 383fa80d070SPierre Jolivet ierr = MatCopy(B,C,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 384fa80d070SPierre Jolivet ierr = MatDiagonalScale(C,a->right,NULL);CHKERRQ(ierr); 385fa80d070SPierre Jolivet } 386fa80d070SPierre Jolivet ierr = MatProductNumeric(contents->work[0]);CHKERRQ(ierr); 387fa80d070SPierre Jolivet ierr = MatDenseGetArrayWrite(C,&array);CHKERRQ(ierr); 388fa80d070SPierre Jolivet ierr = MatDensePlaceArray(contents->work[1],array);CHKERRQ(ierr); 389fa80d070SPierre Jolivet ierr = MatProductNumeric(contents->work[1]);CHKERRQ(ierr); 390fa80d070SPierre Jolivet ierr = MatDenseRestoreArrayWrite(C,&array);CHKERRQ(ierr); 391fa80d070SPierre Jolivet ierr = MatDenseResetArray(contents->work[1]);CHKERRQ(ierr); 392fa80d070SPierre Jolivet ierr = MatSetOption(C,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 393fa80d070SPierre Jolivet ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 394fa80d070SPierre Jolivet ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 395fa80d070SPierre Jolivet ierr = MatScale(C,a->scale);CHKERRQ(ierr); 396fa80d070SPierre Jolivet PetscFunctionReturn(0); 397fa80d070SPierre Jolivet } 398fa80d070SPierre Jolivet 399fa80d070SPierre Jolivet PetscErrorCode MatNormal_DenseDestroy(void *ctx) 400fa80d070SPierre Jolivet { 401fa80d070SPierre Jolivet Normal_Dense *contents = (Normal_Dense*)ctx; 402fa80d070SPierre Jolivet PetscErrorCode ierr; 403fa80d070SPierre Jolivet 404fa80d070SPierre Jolivet PetscFunctionBegin; 405fa80d070SPierre Jolivet ierr = MatDestroy(contents->work);CHKERRQ(ierr); 406fa80d070SPierre Jolivet ierr = MatDestroy(contents->work+1);CHKERRQ(ierr); 407fa80d070SPierre Jolivet ierr = PetscFree(contents);CHKERRQ(ierr); 408fa80d070SPierre Jolivet PetscFunctionReturn(0); 409fa80d070SPierre Jolivet } 410fa80d070SPierre Jolivet 411fa80d070SPierre Jolivet PetscErrorCode MatProductSymbolic_Normal_Dense(Mat C) 412fa80d070SPierre Jolivet { 413fa80d070SPierre Jolivet Mat A,B; 414fa80d070SPierre Jolivet Normal_Dense *contents = NULL; 415fa80d070SPierre Jolivet Mat_Normal *a; 416fa80d070SPierre Jolivet PetscScalar *array; 417fa80d070SPierre Jolivet PetscInt n,N,m,M; 418fa80d070SPierre Jolivet PetscErrorCode ierr; 419fa80d070SPierre Jolivet 420fa80d070SPierre Jolivet PetscFunctionBegin; 421fa80d070SPierre Jolivet MatCheckProduct(C,4); 422*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(C->product->data,PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data not empty"); 423fa80d070SPierre Jolivet A = C->product->A; 424fa80d070SPierre Jolivet a = (Mat_Normal*)A->data; 425*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(a->left,PetscObjectComm((PetscObject)C),PETSC_ERR_SUP,"Not implemented"); 426fa80d070SPierre Jolivet B = C->product->B; 427fa80d070SPierre Jolivet ierr = MatGetLocalSize(C,&m,&n);CHKERRQ(ierr); 428fa80d070SPierre Jolivet ierr = MatGetSize(C,&M,&N);CHKERRQ(ierr); 429fa80d070SPierre Jolivet if (m == PETSC_DECIDE || n == PETSC_DECIDE || M == PETSC_DECIDE || N == PETSC_DECIDE) { 430fa80d070SPierre Jolivet ierr = MatGetLocalSize(B,NULL,&n);CHKERRQ(ierr); 431fa80d070SPierre Jolivet ierr = MatGetSize(B,NULL,&N);CHKERRQ(ierr); 432fa80d070SPierre Jolivet ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr); 433fa80d070SPierre Jolivet ierr = MatGetSize(A,&M,NULL);CHKERRQ(ierr); 434fa80d070SPierre Jolivet ierr = MatSetSizes(C,m,n,M,N);CHKERRQ(ierr); 435fa80d070SPierre Jolivet } 436fa80d070SPierre Jolivet ierr = MatSetType(C,((PetscObject)B)->type_name);CHKERRQ(ierr); 437fa80d070SPierre Jolivet ierr = MatSetUp(C);CHKERRQ(ierr); 438fa80d070SPierre Jolivet ierr = PetscNew(&contents);CHKERRQ(ierr); 439fa80d070SPierre Jolivet C->product->data = contents; 440fa80d070SPierre Jolivet C->product->destroy = MatNormal_DenseDestroy; 441fa80d070SPierre Jolivet if (a->right) { 442fa80d070SPierre Jolivet ierr = MatProductCreate(a->A,C,NULL,contents->work);CHKERRQ(ierr); 443fa80d070SPierre Jolivet } else { 444fa80d070SPierre Jolivet ierr = MatProductCreate(a->A,B,NULL,contents->work);CHKERRQ(ierr); 445fa80d070SPierre Jolivet } 446fa80d070SPierre Jolivet ierr = MatProductSetType(contents->work[0],MATPRODUCT_AB);CHKERRQ(ierr); 447fa80d070SPierre Jolivet ierr = MatProductSetFromOptions(contents->work[0]);CHKERRQ(ierr); 448fa80d070SPierre Jolivet ierr = MatProductSymbolic(contents->work[0]);CHKERRQ(ierr); 449fa80d070SPierre Jolivet ierr = MatProductCreate(a->A,contents->work[0],NULL,contents->work+1);CHKERRQ(ierr); 450fa80d070SPierre Jolivet ierr = MatProductSetType(contents->work[1],MATPRODUCT_AtB);CHKERRQ(ierr); 451fa80d070SPierre Jolivet ierr = MatProductSetFromOptions(contents->work[1]);CHKERRQ(ierr); 452fa80d070SPierre Jolivet ierr = MatProductSymbolic(contents->work[1]);CHKERRQ(ierr); 453fa80d070SPierre Jolivet ierr = MatDenseGetArrayWrite(C,&array);CHKERRQ(ierr); 454fa80d070SPierre Jolivet ierr = MatSeqDenseSetPreallocation(contents->work[1],array);CHKERRQ(ierr); 455fa80d070SPierre Jolivet ierr = MatMPIDenseSetPreallocation(contents->work[1],array);CHKERRQ(ierr); 456fa80d070SPierre Jolivet ierr = MatDenseRestoreArrayWrite(C,&array);CHKERRQ(ierr); 457fa80d070SPierre Jolivet C->ops->productnumeric = MatProductNumeric_Normal_Dense; 458fa80d070SPierre Jolivet PetscFunctionReturn(0); 459fa80d070SPierre Jolivet } 460fa80d070SPierre Jolivet 461fa80d070SPierre Jolivet PetscErrorCode MatProductSetFromOptions_Normal_Dense_AB(Mat C) 462fa80d070SPierre Jolivet { 463fa80d070SPierre Jolivet PetscFunctionBegin; 464fa80d070SPierre Jolivet C->ops->productsymbolic = MatProductSymbolic_Normal_Dense; 465fa80d070SPierre Jolivet PetscFunctionReturn(0); 466fa80d070SPierre Jolivet } 467fa80d070SPierre Jolivet 468fa80d070SPierre Jolivet PetscErrorCode MatProductSetFromOptions_Normal_Dense(Mat C) 469fa80d070SPierre Jolivet { 470fa80d070SPierre Jolivet Mat_Product *product = C->product; 471fa80d070SPierre Jolivet PetscErrorCode ierr; 472fa80d070SPierre Jolivet 473fa80d070SPierre Jolivet PetscFunctionBegin; 474fa80d070SPierre Jolivet if (product->type == MATPRODUCT_AB) { 475fa80d070SPierre Jolivet ierr = MatProductSetFromOptions_Normal_Dense_AB(C);CHKERRQ(ierr); 476fa80d070SPierre Jolivet } 477fa80d070SPierre Jolivet PetscFunctionReturn(0); 478fa80d070SPierre Jolivet } 479fa80d070SPierre Jolivet 480c8a8475eSBarry Smith /*@ 481c8a8475eSBarry Smith MatCreateNormal - Creates a new matrix object that behaves like A'*A. 482c8a8475eSBarry Smith 483c8a8475eSBarry Smith Collective on Mat 484c8a8475eSBarry Smith 485c8a8475eSBarry Smith Input Parameter: 486c8a8475eSBarry Smith . A - the (possibly rectangular) matrix 487c8a8475eSBarry Smith 488c8a8475eSBarry Smith Output Parameter: 489c8a8475eSBarry Smith . N - the matrix that represents A'*A 490c8a8475eSBarry Smith 491c30e7cdbSSatish Balay Level: intermediate 492c30e7cdbSSatish Balay 49395452b02SPatrick Sanan Notes: 49495452b02SPatrick Sanan The product A'*A is NOT actually formed! Rather the new matrix 495c8a8475eSBarry Smith object performs the matrix-vector product by first multiplying by 496c8a8475eSBarry Smith A and then A' 497c8a8475eSBarry Smith @*/ 4987087cfbeSBarry Smith PetscErrorCode MatCreateNormal(Mat A,Mat *N) 499c8a8475eSBarry Smith { 500dfbe8321SBarry Smith PetscErrorCode ierr; 5019ca0e1b6SStefano Zampini PetscInt n,nn; 502c8a8475eSBarry Smith Mat_Normal *Na; 5039ca0e1b6SStefano Zampini VecType vtype; 504c8a8475eSBarry Smith 505c8a8475eSBarry Smith PetscFunctionBegin; 5069ca0e1b6SStefano Zampini ierr = MatGetSize(A,NULL,&nn);CHKERRQ(ierr); 5079ca0e1b6SStefano Zampini ierr = MatGetLocalSize(A,NULL,&n);CHKERRQ(ierr); 508ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),N);CHKERRQ(ierr); 5099ca0e1b6SStefano Zampini ierr = MatSetSizes(*N,n,n,nn,nn);CHKERRQ(ierr); 510c8a8475eSBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)*N,MATNORMAL);CHKERRQ(ierr); 511fa80d070SPierre Jolivet ierr = PetscLayoutReference(A->cmap,&(*N)->rmap);CHKERRQ(ierr); 512fa80d070SPierre Jolivet ierr = PetscLayoutReference(A->cmap,&(*N)->cmap);CHKERRQ(ierr); 513c8a8475eSBarry Smith 514b00a9115SJed Brown ierr = PetscNewLog(*N,&Na);CHKERRQ(ierr); 51538f2d2fdSLisandro Dalcin (*N)->data = (void*) Na; 516c8a8475eSBarry Smith ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 517c3122656SLisandro Dalcin Na->A = A; 518b20f02adSBarry Smith Na->scale = 1.0; 519c8a8475eSBarry Smith 5209ca0e1b6SStefano Zampini ierr = MatCreateVecs(A,NULL,&Na->w);CHKERRQ(ierr); 5212205254eSKarl Rupp 522c8a8475eSBarry Smith (*N)->ops->destroy = MatDestroy_Normal; 523c8a8475eSBarry Smith (*N)->ops->mult = MatMult_Normal; 524b20f02adSBarry Smith (*N)->ops->multtranspose = MatMultTranspose_Normal; 525b20f02adSBarry Smith (*N)->ops->multtransposeadd = MatMultTransposeAdd_Normal; 526db090513SMatthew Knepley (*N)->ops->multadd = MatMultAdd_Normal; 52717a6fd94SBarry Smith (*N)->ops->getdiagonal = MatGetDiagonal_Normal; 52869ef1043SBarry Smith (*N)->ops->scale = MatScale_Normal; 52969ef1043SBarry Smith (*N)->ops->diagonalscale = MatDiagonalScale_Normal; 530fa80d070SPierre Jolivet (*N)->ops->increaseoverlap = MatIncreaseOverlap_Normal; 531fa80d070SPierre Jolivet (*N)->ops->createsubmatrices = MatCreateSubMatrices_Normal; 532fa80d070SPierre Jolivet (*N)->ops->permute = MatPermute_Normal; 533fa80d070SPierre Jolivet (*N)->ops->duplicate = MatDuplicate_Normal; 534fa80d070SPierre Jolivet (*N)->ops->copy = MatCopy_Normal; 535c8a8475eSBarry Smith (*N)->assembled = PETSC_TRUE; 53648331cefSBarry Smith (*N)->preallocated = PETSC_TRUE; 5379ca0e1b6SStefano Zampini 538fa80d070SPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)(*N),"MatNormalGetMat_C",MatNormalGetMat_Normal);CHKERRQ(ierr); 539fa80d070SPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)(*N),"MatConvert_normal_seqaij_C",MatConvert_Normal_AIJ);CHKERRQ(ierr); 540fa80d070SPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)(*N),"MatConvert_normal_mpiaij_C",MatConvert_Normal_AIJ);CHKERRQ(ierr); 541fa80d070SPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)(*N),"MatProductSetFromOptions_normal_seqdense_C",MatProductSetFromOptions_Normal_Dense);CHKERRQ(ierr); 542fa80d070SPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)(*N),"MatProductSetFromOptions_normal_mpidense_C",MatProductSetFromOptions_Normal_Dense);CHKERRQ(ierr); 543fa80d070SPierre Jolivet ierr = PetscObjectComposeFunction((PetscObject)(*N),"MatProductSetFromOptions_normal_dense_C",MatProductSetFromOptions_Normal_Dense);CHKERRQ(ierr); 54412ab1b64SPierre Jolivet ierr = MatSetOption(*N,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr); 5459ca0e1b6SStefano Zampini ierr = MatGetVecType(A,&vtype);CHKERRQ(ierr); 5469ca0e1b6SStefano Zampini ierr = MatSetVecType(*N,vtype);CHKERRQ(ierr); 5479ca0e1b6SStefano Zampini #if defined(PETSC_HAVE_DEVICE) 5489ca0e1b6SStefano Zampini ierr = MatBindToCPU(*N,A->boundtocpu);CHKERRQ(ierr); 5499ca0e1b6SStefano Zampini #endif 550c8a8475eSBarry Smith PetscFunctionReturn(0); 551c8a8475eSBarry Smith } 552