xref: /petsc/src/mat/impls/normal/normm.c (revision fa80d070c7523742fd9c5276039a7d5c181a5908)
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 
44*fa80d070SPierre Jolivet PetscErrorCode MatIncreaseOverlap_Normal(Mat A,PetscInt is_max,IS is[],PetscInt ov)
45*fa80d070SPierre Jolivet {
46*fa80d070SPierre Jolivet   Mat_Normal     *a = (Mat_Normal*)A->data;
47*fa80d070SPierre Jolivet   Mat            pattern;
48*fa80d070SPierre Jolivet   PetscErrorCode ierr;
49*fa80d070SPierre Jolivet 
50*fa80d070SPierre Jolivet   PetscFunctionBegin;
51*fa80d070SPierre Jolivet   if (ov < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap specified");
52*fa80d070SPierre Jolivet   ierr = MatProductCreate(a->A,a->A,NULL,&pattern);CHKERRQ(ierr);
53*fa80d070SPierre Jolivet   ierr = MatProductSetType(pattern,MATPRODUCT_AtB);CHKERRQ(ierr);
54*fa80d070SPierre Jolivet   ierr = MatProductSetFromOptions(pattern);CHKERRQ(ierr);
55*fa80d070SPierre Jolivet   ierr = MatProductSymbolic(pattern);CHKERRQ(ierr);
56*fa80d070SPierre Jolivet   ierr = MatIncreaseOverlap(pattern,is_max,is,ov);CHKERRQ(ierr);
57*fa80d070SPierre Jolivet   ierr = MatDestroy(&pattern);CHKERRQ(ierr);
58*fa80d070SPierre Jolivet   PetscFunctionReturn(0);
59*fa80d070SPierre Jolivet }
60*fa80d070SPierre Jolivet 
61*fa80d070SPierre Jolivet PetscErrorCode MatCreateSubMatrices_Normal(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[])
62*fa80d070SPierre Jolivet {
63*fa80d070SPierre Jolivet   Mat_Normal     *a = (Mat_Normal*)mat->data;
64*fa80d070SPierre Jolivet   Mat            B = a->A, *suba;
65*fa80d070SPierre Jolivet   IS             *row;
66*fa80d070SPierre Jolivet   PetscInt       M;
67*fa80d070SPierre Jolivet   PetscErrorCode ierr;
68*fa80d070SPierre Jolivet 
69*fa80d070SPierre Jolivet   PetscFunctionBegin;
70*fa80d070SPierre Jolivet   if (a->left || a->right || irow != icol) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Not implemented");
71*fa80d070SPierre Jolivet   if (scall != MAT_REUSE_MATRIX) {
72*fa80d070SPierre Jolivet     ierr = PetscCalloc1(n,submat);CHKERRQ(ierr);
73*fa80d070SPierre Jolivet   }
74*fa80d070SPierre Jolivet   ierr = MatGetSize(B,&M,NULL);CHKERRQ(ierr);
75*fa80d070SPierre Jolivet   ierr = PetscMalloc1(n,&row);CHKERRQ(ierr);
76*fa80d070SPierre Jolivet   ierr = ISCreateStride(PETSC_COMM_SELF,M,0,1,&row[0]);CHKERRQ(ierr);
77*fa80d070SPierre Jolivet   ierr = ISSetIdentity(row[0]);CHKERRQ(ierr);
78*fa80d070SPierre Jolivet   for (M = 1; M < n; ++M) row[M] = row[0];
79*fa80d070SPierre Jolivet   ierr = MatCreateSubMatrices(B,n,row,icol,MAT_INITIAL_MATRIX,&suba);CHKERRQ(ierr);
80*fa80d070SPierre Jolivet   for (M = 0; M < n; ++M) {
81*fa80d070SPierre Jolivet     ierr = MatCreateNormal(suba[M],*submat+M);CHKERRQ(ierr);
82*fa80d070SPierre Jolivet     ((Mat_Normal*)(*submat)[M]->data)->scale = a->scale;
83*fa80d070SPierre Jolivet   }
84*fa80d070SPierre Jolivet   ierr = ISDestroy(&row[0]);CHKERRQ(ierr);
85*fa80d070SPierre Jolivet   ierr = PetscFree(row);CHKERRQ(ierr);
86*fa80d070SPierre Jolivet   ierr = MatDestroySubMatrices(n,&suba);CHKERRQ(ierr);
87*fa80d070SPierre Jolivet   PetscFunctionReturn(0);
88*fa80d070SPierre Jolivet }
89*fa80d070SPierre Jolivet 
90*fa80d070SPierre Jolivet PetscErrorCode MatPermute_Normal(Mat A,IS rowp,IS colp,Mat *B)
91*fa80d070SPierre Jolivet {
92*fa80d070SPierre Jolivet   Mat_Normal     *a = (Mat_Normal*)A->data;
93*fa80d070SPierre Jolivet   Mat            C,Aa = a->A;
94*fa80d070SPierre Jolivet   IS             row;
95*fa80d070SPierre Jolivet   PetscErrorCode ierr;
96*fa80d070SPierre Jolivet 
97*fa80d070SPierre Jolivet   PetscFunctionBegin;
98*fa80d070SPierre Jolivet   if (rowp != colp) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Row permutation and column permutation must be the same");
99*fa80d070SPierre Jolivet   ierr = ISCreateStride(PetscObjectComm((PetscObject)Aa),Aa->rmap->n,Aa->rmap->rstart,1,&row);CHKERRQ(ierr);
100*fa80d070SPierre Jolivet   ierr = ISSetIdentity(row);CHKERRQ(ierr);
101*fa80d070SPierre Jolivet   ierr = MatPermute(Aa,row,colp,&C);CHKERRQ(ierr);
102*fa80d070SPierre Jolivet   ierr = ISDestroy(&row);CHKERRQ(ierr);
103*fa80d070SPierre Jolivet   ierr = MatCreateNormal(C,B);CHKERRQ(ierr);
104*fa80d070SPierre Jolivet   ierr = MatDestroy(&C);CHKERRQ(ierr);
105*fa80d070SPierre Jolivet   PetscFunctionReturn(0);
106*fa80d070SPierre Jolivet }
107*fa80d070SPierre Jolivet 
108*fa80d070SPierre Jolivet PetscErrorCode MatDuplicate_Normal(Mat A, MatDuplicateOption op, Mat *B)
109*fa80d070SPierre Jolivet {
110*fa80d070SPierre Jolivet   Mat_Normal     *a = (Mat_Normal*)A->data;
111*fa80d070SPierre Jolivet   Mat            C;
112*fa80d070SPierre Jolivet   PetscErrorCode ierr;
113*fa80d070SPierre Jolivet 
114*fa80d070SPierre Jolivet   PetscFunctionBegin;
115*fa80d070SPierre Jolivet   if (a->left || a->right) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Not implemented");
116*fa80d070SPierre Jolivet   ierr = MatDuplicate(a->A,op,&C);CHKERRQ(ierr);
117*fa80d070SPierre Jolivet   ierr = MatCreateNormal(C,B);CHKERRQ(ierr);
118*fa80d070SPierre Jolivet   ierr = MatDestroy(&C);CHKERRQ(ierr);
119*fa80d070SPierre Jolivet   if (op == MAT_COPY_VALUES) ((Mat_Normal*)(*B)->data)->scale = a->scale;
120*fa80d070SPierre Jolivet   PetscFunctionReturn(0);
121*fa80d070SPierre Jolivet }
122*fa80d070SPierre Jolivet 
123*fa80d070SPierre Jolivet PetscErrorCode MatCopy_Normal(Mat A,Mat B,MatStructure str)
124*fa80d070SPierre Jolivet {
125*fa80d070SPierre Jolivet   Mat_Normal     *a = (Mat_Normal*)A->data,*b = (Mat_Normal*)B->data;
126*fa80d070SPierre Jolivet   PetscErrorCode ierr;
127*fa80d070SPierre Jolivet 
128*fa80d070SPierre Jolivet   PetscFunctionBegin;
129*fa80d070SPierre Jolivet   if (a->left || a->right) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Not implemented");
130*fa80d070SPierre Jolivet   ierr = MatCopy(a->A,b->A,str);CHKERRQ(ierr);
131*fa80d070SPierre Jolivet   b->scale = a->scale;
132*fa80d070SPierre Jolivet   ierr = VecDestroy(&b->left);CHKERRQ(ierr);
133*fa80d070SPierre Jolivet   ierr = VecDestroy(&b->right);CHKERRQ(ierr);
134*fa80d070SPierre Jolivet   ierr = VecDestroy(&b->leftwork);CHKERRQ(ierr);
135*fa80d070SPierre Jolivet   ierr = VecDestroy(&b->rightwork);CHKERRQ(ierr);
136*fa80d070SPierre Jolivet   PetscFunctionReturn(0);
137*fa80d070SPierre Jolivet }
138*fa80d070SPierre 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);
254*fa80d070SPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)N,"MatNormalGetMat_C",NULL);CHKERRQ(ierr);
255*fa80d070SPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)N,"MatConvert_normal_seqaij_C",NULL);CHKERRQ(ierr);
256*fa80d070SPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)N,"MatConvert_normal_mpiaij_C",NULL);CHKERRQ(ierr);
257*fa80d070SPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)N,"MatProductSetFromOptions_normal_seqdense_C",NULL);CHKERRQ(ierr);
258*fa80d070SPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)N,"MatProductSetFromOptions_normal_mpidense_C",NULL);CHKERRQ(ierr);
259*fa80d070SPierre 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 
298*fa80d070SPierre Jolivet PetscErrorCode MatNormalGetMat_Normal(Mat A,Mat *M)
299*fa80d070SPierre Jolivet {
300*fa80d070SPierre Jolivet   Mat_Normal *Aa = (Mat_Normal*)A->data;
301*fa80d070SPierre Jolivet 
302*fa80d070SPierre Jolivet   PetscFunctionBegin;
303*fa80d070SPierre Jolivet   *M = Aa->A;
304*fa80d070SPierre Jolivet   PetscFunctionReturn(0);
305*fa80d070SPierre Jolivet }
306*fa80d070SPierre Jolivet 
307*fa80d070SPierre Jolivet /*@
308*fa80d070SPierre Jolivet       MatNormalGetMat - Gets the Mat object stored inside a MATNORMAL
309*fa80d070SPierre Jolivet 
310*fa80d070SPierre Jolivet    Logically collective on Mat
311*fa80d070SPierre Jolivet 
312*fa80d070SPierre Jolivet    Input Parameter:
313*fa80d070SPierre Jolivet .   A  - the MATNORMAL matrix
314*fa80d070SPierre Jolivet 
315*fa80d070SPierre Jolivet    Output Parameter:
316*fa80d070SPierre Jolivet .   M - the matrix object stored inside A
317*fa80d070SPierre Jolivet 
318*fa80d070SPierre Jolivet    Level: intermediate
319*fa80d070SPierre Jolivet 
320*fa80d070SPierre Jolivet .seealso: MatCreateNormal()
321*fa80d070SPierre Jolivet 
322*fa80d070SPierre Jolivet @*/
323*fa80d070SPierre Jolivet PetscErrorCode MatNormalGetMat(Mat A,Mat *M)
324*fa80d070SPierre Jolivet {
325*fa80d070SPierre Jolivet   PetscErrorCode ierr;
326*fa80d070SPierre Jolivet 
327*fa80d070SPierre Jolivet   PetscFunctionBegin;
328*fa80d070SPierre Jolivet   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
329*fa80d070SPierre Jolivet   PetscValidType(A,1);
330*fa80d070SPierre Jolivet   PetscValidPointer(M,2);
331*fa80d070SPierre Jolivet   ierr = PetscUseMethod(A,"MatNormalGetMat_C",(Mat,Mat*),(A,M));CHKERRQ(ierr);
332*fa80d070SPierre Jolivet   PetscFunctionReturn(0);
333*fa80d070SPierre Jolivet }
334*fa80d070SPierre Jolivet 
335*fa80d070SPierre Jolivet PetscErrorCode MatConvert_Normal_AIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
336*fa80d070SPierre Jolivet {
337*fa80d070SPierre Jolivet   Mat_Normal     *Aa = (Mat_Normal*)A->data;
338*fa80d070SPierre Jolivet   Mat            B;
339*fa80d070SPierre Jolivet   PetscInt       m,n,M,N;
340*fa80d070SPierre Jolivet   PetscErrorCode ierr;
341*fa80d070SPierre Jolivet 
342*fa80d070SPierre Jolivet   PetscFunctionBegin;
343*fa80d070SPierre Jolivet   ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr);
344*fa80d070SPierre Jolivet   ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr);
345*fa80d070SPierre Jolivet   if (reuse == MAT_REUSE_MATRIX) {
346*fa80d070SPierre Jolivet     B = *newmat;
347*fa80d070SPierre Jolivet     ierr = MatProductReplaceMats(Aa->A,Aa->A,NULL,B);CHKERRQ(ierr);
348*fa80d070SPierre Jolivet   } else {
349*fa80d070SPierre Jolivet     ierr = MatProductCreate(Aa->A,Aa->A,NULL,&B);CHKERRQ(ierr);
350*fa80d070SPierre Jolivet     ierr = MatProductSetType(B,MATPRODUCT_AtB);CHKERRQ(ierr);
351*fa80d070SPierre Jolivet     ierr = MatProductSetFromOptions(B);CHKERRQ(ierr);
352*fa80d070SPierre Jolivet     ierr = MatProductSymbolic(B);CHKERRQ(ierr);
353*fa80d070SPierre Jolivet     ierr = MatSetOption(B,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr);
354*fa80d070SPierre Jolivet   }
355*fa80d070SPierre Jolivet   ierr = MatProductNumeric(B);CHKERRQ(ierr);
356*fa80d070SPierre Jolivet   if (reuse == MAT_INPLACE_MATRIX) {
357*fa80d070SPierre Jolivet     ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr);
358*fa80d070SPierre Jolivet   } else if (reuse == MAT_INITIAL_MATRIX) *newmat = B;
359*fa80d070SPierre Jolivet   ierr = MatConvert(*newmat,MATAIJ,MAT_INPLACE_MATRIX,newmat);CHKERRQ(ierr);
360*fa80d070SPierre Jolivet   PetscFunctionReturn(0);
361*fa80d070SPierre Jolivet }
362*fa80d070SPierre Jolivet 
363*fa80d070SPierre Jolivet typedef struct {
364*fa80d070SPierre Jolivet   Mat work[2];
365*fa80d070SPierre Jolivet } Normal_Dense;
366*fa80d070SPierre Jolivet 
367*fa80d070SPierre Jolivet PetscErrorCode MatProductNumeric_Normal_Dense(Mat C)
368*fa80d070SPierre Jolivet {
369*fa80d070SPierre Jolivet   Mat            A,B;
370*fa80d070SPierre Jolivet   Normal_Dense   *contents;
371*fa80d070SPierre Jolivet   Mat_Normal     *a;
372*fa80d070SPierre Jolivet   PetscScalar    *array;
373*fa80d070SPierre Jolivet   PetscErrorCode ierr;
374*fa80d070SPierre Jolivet 
375*fa80d070SPierre Jolivet   PetscFunctionBegin;
376*fa80d070SPierre Jolivet   MatCheckProduct(C,3);
377*fa80d070SPierre Jolivet   A = C->product->A;
378*fa80d070SPierre Jolivet   a = (Mat_Normal*)A->data;
379*fa80d070SPierre Jolivet   B = C->product->B;
380*fa80d070SPierre Jolivet   contents = (Normal_Dense*)C->product->data;
381*fa80d070SPierre Jolivet   if (!contents) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data empty");
382*fa80d070SPierre Jolivet   if (a->right) {
383*fa80d070SPierre Jolivet     ierr = MatCopy(B,C,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
384*fa80d070SPierre Jolivet     ierr = MatDiagonalScale(C,a->right,NULL);CHKERRQ(ierr);
385*fa80d070SPierre Jolivet   }
386*fa80d070SPierre Jolivet   ierr = MatProductNumeric(contents->work[0]);CHKERRQ(ierr);
387*fa80d070SPierre Jolivet   ierr = MatDenseGetArrayWrite(C,&array);CHKERRQ(ierr);
388*fa80d070SPierre Jolivet   ierr = MatDensePlaceArray(contents->work[1],array);CHKERRQ(ierr);
389*fa80d070SPierre Jolivet   ierr = MatProductNumeric(contents->work[1]);CHKERRQ(ierr);
390*fa80d070SPierre Jolivet   ierr = MatDenseRestoreArrayWrite(C,&array);CHKERRQ(ierr);
391*fa80d070SPierre Jolivet   ierr = MatDenseResetArray(contents->work[1]);CHKERRQ(ierr);
392*fa80d070SPierre Jolivet   ierr = MatSetOption(C,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
393*fa80d070SPierre Jolivet   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
394*fa80d070SPierre Jolivet   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
395*fa80d070SPierre Jolivet   ierr = MatScale(C,a->scale);CHKERRQ(ierr);
396*fa80d070SPierre Jolivet   PetscFunctionReturn(0);
397*fa80d070SPierre Jolivet }
398*fa80d070SPierre Jolivet 
399*fa80d070SPierre Jolivet PetscErrorCode MatNormal_DenseDestroy(void *ctx)
400*fa80d070SPierre Jolivet {
401*fa80d070SPierre Jolivet   Normal_Dense   *contents = (Normal_Dense*)ctx;
402*fa80d070SPierre Jolivet   PetscErrorCode ierr;
403*fa80d070SPierre Jolivet 
404*fa80d070SPierre Jolivet   PetscFunctionBegin;
405*fa80d070SPierre Jolivet   ierr = MatDestroy(contents->work);CHKERRQ(ierr);
406*fa80d070SPierre Jolivet   ierr = MatDestroy(contents->work+1);CHKERRQ(ierr);
407*fa80d070SPierre Jolivet   ierr = PetscFree(contents);CHKERRQ(ierr);
408*fa80d070SPierre Jolivet   PetscFunctionReturn(0);
409*fa80d070SPierre Jolivet }
410*fa80d070SPierre Jolivet 
411*fa80d070SPierre Jolivet PetscErrorCode MatProductSymbolic_Normal_Dense(Mat C)
412*fa80d070SPierre Jolivet {
413*fa80d070SPierre Jolivet   Mat            A,B;
414*fa80d070SPierre Jolivet   Normal_Dense   *contents = NULL;
415*fa80d070SPierre Jolivet   Mat_Normal     *a;
416*fa80d070SPierre Jolivet   PetscScalar    *array;
417*fa80d070SPierre Jolivet   PetscInt       n,N,m,M;
418*fa80d070SPierre Jolivet   PetscErrorCode ierr;
419*fa80d070SPierre Jolivet 
420*fa80d070SPierre Jolivet   PetscFunctionBegin;
421*fa80d070SPierre Jolivet   MatCheckProduct(C,4);
422*fa80d070SPierre Jolivet   if (C->product->data) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data not empty");
423*fa80d070SPierre Jolivet   A = C->product->A;
424*fa80d070SPierre Jolivet   a = (Mat_Normal*)A->data;
425*fa80d070SPierre Jolivet   if (a->left) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_SUP,"Not implemented");
426*fa80d070SPierre Jolivet   B = C->product->B;
427*fa80d070SPierre Jolivet   ierr = MatGetLocalSize(C,&m,&n);CHKERRQ(ierr);
428*fa80d070SPierre Jolivet   ierr = MatGetSize(C,&M,&N);CHKERRQ(ierr);
429*fa80d070SPierre Jolivet   if (m == PETSC_DECIDE || n == PETSC_DECIDE || M == PETSC_DECIDE || N == PETSC_DECIDE) {
430*fa80d070SPierre Jolivet     ierr = MatGetLocalSize(B,NULL,&n);CHKERRQ(ierr);
431*fa80d070SPierre Jolivet     ierr = MatGetSize(B,NULL,&N);CHKERRQ(ierr);
432*fa80d070SPierre Jolivet     ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr);
433*fa80d070SPierre Jolivet     ierr = MatGetSize(A,&M,NULL);CHKERRQ(ierr);
434*fa80d070SPierre Jolivet     ierr = MatSetSizes(C,m,n,M,N);CHKERRQ(ierr);
435*fa80d070SPierre Jolivet   }
436*fa80d070SPierre Jolivet   ierr = MatSetType(C,((PetscObject)B)->type_name);CHKERRQ(ierr);
437*fa80d070SPierre Jolivet   ierr = MatSetUp(C);CHKERRQ(ierr);
438*fa80d070SPierre Jolivet   ierr = PetscNew(&contents);CHKERRQ(ierr);
439*fa80d070SPierre Jolivet   C->product->data = contents;
440*fa80d070SPierre Jolivet   C->product->destroy = MatNormal_DenseDestroy;
441*fa80d070SPierre Jolivet   if (a->right) {
442*fa80d070SPierre Jolivet     ierr = MatProductCreate(a->A,C,NULL,contents->work);CHKERRQ(ierr);
443*fa80d070SPierre Jolivet   } else {
444*fa80d070SPierre Jolivet     ierr = MatProductCreate(a->A,B,NULL,contents->work);CHKERRQ(ierr);
445*fa80d070SPierre Jolivet   }
446*fa80d070SPierre Jolivet   ierr = MatProductSetType(contents->work[0],MATPRODUCT_AB);CHKERRQ(ierr);
447*fa80d070SPierre Jolivet   ierr = MatProductSetFromOptions(contents->work[0]);CHKERRQ(ierr);
448*fa80d070SPierre Jolivet   ierr = MatProductSymbolic(contents->work[0]);CHKERRQ(ierr);
449*fa80d070SPierre Jolivet   ierr = MatProductCreate(a->A,contents->work[0],NULL,contents->work+1);CHKERRQ(ierr);
450*fa80d070SPierre Jolivet   ierr = MatProductSetType(contents->work[1],MATPRODUCT_AtB);CHKERRQ(ierr);
451*fa80d070SPierre Jolivet   ierr = MatProductSetFromOptions(contents->work[1]);CHKERRQ(ierr);
452*fa80d070SPierre Jolivet   ierr = MatProductSymbolic(contents->work[1]);CHKERRQ(ierr);
453*fa80d070SPierre Jolivet   ierr = MatDenseGetArrayWrite(C,&array);CHKERRQ(ierr);
454*fa80d070SPierre Jolivet   ierr = MatSeqDenseSetPreallocation(contents->work[1],array);CHKERRQ(ierr);
455*fa80d070SPierre Jolivet   ierr = MatMPIDenseSetPreallocation(contents->work[1],array);CHKERRQ(ierr);
456*fa80d070SPierre Jolivet   ierr = MatDenseRestoreArrayWrite(C,&array);CHKERRQ(ierr);
457*fa80d070SPierre Jolivet   C->ops->productnumeric = MatProductNumeric_Normal_Dense;
458*fa80d070SPierre Jolivet   PetscFunctionReturn(0);
459*fa80d070SPierre Jolivet }
460*fa80d070SPierre Jolivet 
461*fa80d070SPierre Jolivet PetscErrorCode MatProductSetFromOptions_Normal_Dense_AB(Mat C)
462*fa80d070SPierre Jolivet {
463*fa80d070SPierre Jolivet   PetscFunctionBegin;
464*fa80d070SPierre Jolivet   C->ops->productsymbolic = MatProductSymbolic_Normal_Dense;
465*fa80d070SPierre Jolivet   PetscFunctionReturn(0);
466*fa80d070SPierre Jolivet }
467*fa80d070SPierre Jolivet 
468*fa80d070SPierre Jolivet PetscErrorCode MatProductSetFromOptions_Normal_Dense(Mat C)
469*fa80d070SPierre Jolivet {
470*fa80d070SPierre Jolivet   Mat_Product    *product = C->product;
471*fa80d070SPierre Jolivet   PetscErrorCode ierr;
472*fa80d070SPierre Jolivet 
473*fa80d070SPierre Jolivet   PetscFunctionBegin;
474*fa80d070SPierre Jolivet   if (product->type == MATPRODUCT_AB) {
475*fa80d070SPierre Jolivet     ierr = MatProductSetFromOptions_Normal_Dense_AB(C);CHKERRQ(ierr);
476*fa80d070SPierre Jolivet   }
477*fa80d070SPierre Jolivet   PetscFunctionReturn(0);
478*fa80d070SPierre Jolivet }
479*fa80d070SPierre 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);
511*fa80d070SPierre Jolivet   ierr = PetscLayoutReference(A->cmap,&(*N)->rmap);CHKERRQ(ierr);
512*fa80d070SPierre 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;
530*fa80d070SPierre Jolivet   (*N)->ops->increaseoverlap   = MatIncreaseOverlap_Normal;
531*fa80d070SPierre Jolivet   (*N)->ops->createsubmatrices = MatCreateSubMatrices_Normal;
532*fa80d070SPierre Jolivet   (*N)->ops->permute           = MatPermute_Normal;
533*fa80d070SPierre Jolivet   (*N)->ops->duplicate         = MatDuplicate_Normal;
534*fa80d070SPierre Jolivet   (*N)->ops->copy              = MatCopy_Normal;
535c8a8475eSBarry Smith   (*N)->assembled              = PETSC_TRUE;
53648331cefSBarry Smith   (*N)->preallocated           = PETSC_TRUE;
5379ca0e1b6SStefano Zampini 
538*fa80d070SPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)(*N),"MatNormalGetMat_C",MatNormalGetMat_Normal);CHKERRQ(ierr);
539*fa80d070SPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)(*N),"MatConvert_normal_seqaij_C",MatConvert_Normal_AIJ);CHKERRQ(ierr);
540*fa80d070SPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)(*N),"MatConvert_normal_mpiaij_C",MatConvert_Normal_AIJ);CHKERRQ(ierr);
541*fa80d070SPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)(*N),"MatProductSetFromOptions_normal_seqdense_C",MatProductSetFromOptions_Normal_Dense);CHKERRQ(ierr);
542*fa80d070SPierre Jolivet   ierr = PetscObjectComposeFunction((PetscObject)(*N),"MatProductSetFromOptions_normal_mpidense_C",MatProductSetFromOptions_Normal_Dense);CHKERRQ(ierr);
543*fa80d070SPierre 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