Actual source code: normm.c

 2:  #include src/mat/matimpl.h

  4: typedef struct {
  5:   Mat A;
  6:   Vec w;
  7: } Mat_Normal;

 11: PetscErrorCode MatMult_Normal(Mat N,Vec x,Vec y)
 12: {
 13:   Mat_Normal     *Na = (Mat_Normal*)N->data;

 17:   MatMult(Na->A,x,Na->w);
 18:   MatMultTranspose(Na->A,Na->w,y);
 19:   return(0);
 20: }

 24: PetscErrorCode MatDestroy_Normal(Mat N)
 25: {
 26:   Mat_Normal     *Na = (Mat_Normal*)N->data;

 30:   PetscObjectDereference((PetscObject)Na->A);
 31:   VecDestroy(Na->w);
 32:   PetscFree(Na);
 33:   return(0);
 34: }
 35: 
 36: /*
 37:       Slow, nonscalable version
 38: */
 41: PetscErrorCode MatGetDiagonal_Normal(Mat N,Vec v)
 42: {
 43:   Mat_Normal        *Na = (Mat_Normal*)N->data;
 44:   Mat               A = Na->A;
 45:   PetscErrorCode    ierr;
 46:   PetscInt          i,j,rstart,rend,nnz;
 47:   const PetscInt    *cols;
 48:   PetscScalar       *diag,*work,*values;
 49:   const PetscScalar *mvalues;
 50:   PetscMap          cmap;

 53:   PetscMalloc(2*A->N*sizeof(PetscScalar),&diag);
 54:   work = diag + A->N;
 55:   PetscMemzero(work,A->N*sizeof(PetscScalar));
 56:   MatGetOwnershipRange(A,&rstart,&rend);
 57:   for (i=rstart; i<rend; i++) {
 58:     MatGetRow(A,i,&nnz,&cols,&mvalues);
 59:     for (j=0; j<nnz; j++) {
 60:       work[cols[j]] += mvalues[j]*mvalues[j];
 61:     }
 62:     MatRestoreRow(A,i,&nnz,&cols,&mvalues);
 63:   }
 64:   MPI_Allreduce(work,diag,A->N,MPIU_SCALAR,MPI_SUM,N->comm);
 65:   MatGetPetscMaps(A,PETSC_NULL,&cmap);
 66:   PetscMapGetLocalRange(cmap,&rstart,&rend);
 67:   VecGetArray(v,&values);
 68:   PetscMemcpy(values,diag+rstart,(rend-rstart)*sizeof(PetscScalar));
 69:   VecRestoreArray(v,&values);
 70:   PetscFree(diag);
 71:   return(0);
 72: }

 76: /*@
 77:       MatCreateNormal - Creates a new matrix object that behaves like A'*A.

 79:    Collective on Mat

 81:    Input Parameter:
 82: .   A  - the (possibly rectangular) matrix

 84:    Output Parameter:
 85: .   N - the matrix that represents A'*A

 87:    Level: intermediate

 89:    Notes: The product A'*A is NOT actually formed! Rather the new matrix
 90:           object performs the matrix-vector product by first multiplying by
 91:           A and then A'
 92: @*/
 93: PetscErrorCode MatCreateNormal(Mat A,Mat *N)
 94: {
 96:   PetscInt       m,n;
 97:   Mat_Normal     *Na;

100:   MatGetLocalSize(A,&m,&n);
101:   MatCreate(A->comm,n,n,PETSC_DECIDE,PETSC_DECIDE,N);
102:   PetscObjectChangeTypeName((PetscObject)*N,MATNORMAL);
103: 
104:   PetscNew(Mat_Normal,&Na);
105:   Na->A     = A;
106:   PetscObjectReference((PetscObject)A);
107:   (*N)->data = (void*) Na;

109:   VecCreateMPI(A->comm,m,PETSC_DECIDE,&Na->w);
110:   (*N)->ops->destroy     = MatDestroy_Normal;
111:   (*N)->ops->mult        = MatMult_Normal;
112:   (*N)->ops->getdiagonal = MatGetDiagonal_Normal;
113:   (*N)->assembled        = PETSC_TRUE;
114:   (*N)->N                = A->N;
115:   (*N)->M                = A->N;
116:   (*N)->n                = A->n;
117:   (*N)->m                = A->n;
118:   return(0);
119: }