Actual source code: normm.c
1: #define PETSCMAT_DLL
3: #include private/matimpl.h
5: typedef struct {
6: Mat A;
7: Vec w,left,right,leftwork,rightwork;
8: PetscScalar scale;
9: } Mat_Normal;
13: PetscErrorCode MatScale_Normal(Mat inA,PetscScalar scale)
14: {
15: Mat_Normal *a = (Mat_Normal*)inA->data;
18: a->scale *= scale;
19: return(0);
20: }
24: PetscErrorCode MatDiagonalScale_Normal(Mat inA,Vec left,Vec right)
25: {
26: Mat_Normal *a = (Mat_Normal*)inA->data;
30: if (left) {
31: if (!a->left) {
32: VecDuplicate(left,&a->left);
33: VecCopy(left,a->left);
34: } else {
35: VecPointwiseMult(a->left,left,a->left);
36: }
37: }
38: if (right) {
39: if (!a->right) {
40: VecDuplicate(right,&a->right);
41: VecCopy(right,a->right);
42: } else {
43: VecPointwiseMult(a->right,right,a->right);
44: }
45: }
46: return(0);
47: }
51: PetscErrorCode MatMult_Normal(Mat N,Vec x,Vec y)
52: {
53: Mat_Normal *Na = (Mat_Normal*)N->data;
55: Vec in;
58: in = x;
59: if (Na->right) {
60: if (!Na->rightwork) {
61: VecDuplicate(Na->right,&Na->rightwork);
62: }
63: VecPointwiseMult(Na->rightwork,Na->right,in);
64: in = Na->rightwork;
65: }
66: MatMult(Na->A,in,Na->w);
67: MatMultTranspose(Na->A,Na->w,y);
68: if (Na->left) {
69: VecPointwiseMult(y,Na->left,y);
70: }
71: VecScale(y,Na->scale);
72: return(0);
73: }
74:
77: PetscErrorCode MatMultAdd_Normal(Mat N,Vec v1,Vec v2,Vec v3)
78: {
79: Mat_Normal *Na = (Mat_Normal*)N->data;
81: Vec in;
84: in = v1;
85: if (Na->right) {
86: if (!Na->rightwork) {
87: VecDuplicate(Na->right,&Na->rightwork);
88: }
89: VecPointwiseMult(Na->rightwork,Na->right,in);
90: in = Na->rightwork;
91: }
92: MatMult(Na->A,in,Na->w);
93: VecScale(Na->w,Na->scale);
94: if (Na->left) {
95: MatMultTranspose(Na->A,Na->w,v3);
96: VecPointwiseMult(v3,Na->left,v3);
97: VecAXPY(v3,1.0,v2);
98: } else {
99: MatMultTransposeAdd(Na->A,Na->w,v2,v3);
100: }
101: return(0);
102: }
106: PetscErrorCode MatMultTranspose_Normal(Mat N,Vec x,Vec y)
107: {
108: Mat_Normal *Na = (Mat_Normal*)N->data;
110: Vec in;
113: in = x;
114: if (Na->left) {
115: if (!Na->leftwork) {
116: VecDuplicate(Na->left,&Na->leftwork);
117: }
118: VecPointwiseMult(Na->leftwork,Na->left,in);
119: in = Na->leftwork;
120: }
121: MatMult(Na->A,in,Na->w);
122: MatMultTranspose(Na->A,Na->w,y);
123: if (Na->right) {
124: VecPointwiseMult(y,Na->right,y);
125: }
126: VecScale(y,Na->scale);
127: return(0);
128: }
132: PetscErrorCode MatMultTransposeAdd_Normal(Mat N,Vec v1,Vec v2,Vec v3)
133: {
134: Mat_Normal *Na = (Mat_Normal*)N->data;
136: Vec in;
139: in = v1;
140: if (Na->left) {
141: if (!Na->leftwork) {
142: VecDuplicate(Na->left,&Na->leftwork);
143: }
144: VecPointwiseMult(Na->leftwork,Na->left,in);
145: in = Na->leftwork;
146: }
147: MatMult(Na->A,in,Na->w);
148: VecScale(Na->w,Na->scale);
149: if (Na->right) {
150: MatMultTranspose(Na->A,Na->w,v3);
151: VecPointwiseMult(v3,Na->right,v3);
152: VecAXPY(v3,1.0,v2);
153: } else {
154: MatMultTransposeAdd(Na->A,Na->w,v2,v3);
155: }
156: return(0);
157: }
161: PetscErrorCode MatDestroy_Normal(Mat N)
162: {
163: Mat_Normal *Na = (Mat_Normal*)N->data;
167: if (Na->A) { MatDestroy(Na->A); }
168: if (Na->w) { VecDestroy(Na->w); }
169: if (Na->left) { VecDestroy(Na->left); }
170: if (Na->right) { VecDestroy(Na->right); }
171: if (Na->leftwork) { VecDestroy(Na->leftwork); }
172: if (Na->rightwork) { VecDestroy(Na->rightwork); }
173: PetscFree(Na);
174: return(0);
175: }
176:
177: /*
178: Slow, nonscalable version
179: */
182: PetscErrorCode MatGetDiagonal_Normal(Mat N,Vec v)
183: {
184: Mat_Normal *Na = (Mat_Normal*)N->data;
185: Mat A = Na->A;
186: PetscErrorCode ierr;
187: PetscInt i,j,rstart,rend,nnz;
188: const PetscInt *cols;
189: PetscScalar *diag,*work,*values;
190: const PetscScalar *mvalues;
193: PetscMalloc(2*A->cmap->N*sizeof(PetscScalar),&diag);
194: work = diag + A->cmap->N;
195: PetscMemzero(work,A->cmap->N*sizeof(PetscScalar));
196: MatGetOwnershipRange(A,&rstart,&rend);
197: for (i=rstart; i<rend; i++) {
198: MatGetRow(A,i,&nnz,&cols,&mvalues);
199: for (j=0; j<nnz; j++) {
200: work[cols[j]] += mvalues[j]*mvalues[j];
201: }
202: MatRestoreRow(A,i,&nnz,&cols,&mvalues);
203: }
204: MPI_Allreduce(work,diag,A->cmap->N,MPIU_SCALAR,MPI_SUM,((PetscObject)N)->comm);
205: rstart = N->cmap->rstart;
206: rend = N->cmap->rend;
207: VecGetArray(v,&values);
208: PetscMemcpy(values,diag+rstart,(rend-rstart)*sizeof(PetscScalar));
209: VecRestoreArray(v,&values);
210: PetscFree(diag);
211: VecScale(v,Na->scale);
212: return(0);
213: }
217: /*@
218: MatCreateNormal - Creates a new matrix object that behaves like A'*A.
220: Collective on Mat
222: Input Parameter:
223: . A - the (possibly rectangular) matrix
225: Output Parameter:
226: . N - the matrix that represents A'*A
228: Level: intermediate
230: Notes: The product A'*A is NOT actually formed! Rather the new matrix
231: object performs the matrix-vector product by first multiplying by
232: A and then A'
233: @*/
234: PetscErrorCode MatCreateNormal(Mat A,Mat *N)
235: {
237: PetscInt m,n;
238: Mat_Normal *Na;
241: MatGetLocalSize(A,&m,&n);
242: MatCreate(((PetscObject)A)->comm,N);
243: MatSetSizes(*N,n,n,PETSC_DECIDE,PETSC_DECIDE);
244: PetscObjectChangeTypeName((PetscObject)*N,MATNORMAL);
245:
246: PetscNewLog(*N,Mat_Normal,&Na);
247: (*N)->data = (void*) Na;
248: PetscObjectReference((PetscObject)A);
249: Na->A = A;
250: Na->scale = 1.0;
252: VecCreateMPI(((PetscObject)A)->comm,m,PETSC_DECIDE,&Na->w);
253: (*N)->ops->destroy = MatDestroy_Normal;
254: (*N)->ops->mult = MatMult_Normal;
255: (*N)->ops->multtranspose = MatMultTranspose_Normal;
256: (*N)->ops->multtransposeadd = MatMultTransposeAdd_Normal;
257: (*N)->ops->multadd = MatMultAdd_Normal;
258: (*N)->ops->getdiagonal = MatGetDiagonal_Normal;
259: (*N)->ops->scale = MatScale_Normal;
260: (*N)->ops->diagonalscale = MatDiagonalScale_Normal;
261: (*N)->assembled = PETSC_TRUE;
262: (*N)->cmap->N = A->cmap->N;
263: (*N)->rmap->N = A->cmap->N;
264: (*N)->cmap->n = A->cmap->n;
265: (*N)->rmap->n = A->cmap->n;
266: return(0);
267: }