1 #define PETSCMAT_DLL 2 3 #include "src/mat/matimpl.h" /*I "petscmat.h" I*/ 4 5 typedef struct { 6 Mat A; 7 Vec w; 8 } Mat_Normal; 9 10 #undef __FUNCT__ 11 #define __FUNCT__ "MatMult_Normal" 12 PetscErrorCode MatMult_Normal(Mat N,Vec x,Vec y) 13 { 14 Mat_Normal *Na = (Mat_Normal*)N->data; 15 PetscErrorCode ierr; 16 17 PetscFunctionBegin; 18 ierr = MatMult(Na->A,x,Na->w);CHKERRQ(ierr); 19 ierr = MatMultTranspose(Na->A,Na->w,y);CHKERRQ(ierr); 20 PetscFunctionReturn(0); 21 } 22 23 #undef __FUNCT__ 24 #define __FUNCT__ "MatMultAdd_Normal" 25 PetscErrorCode MatMultAdd_Normal(Mat N,Vec v1,Vec v2,Vec v3) 26 { 27 Mat_Normal *Na = (Mat_Normal*)N->data; 28 PetscErrorCode ierr; 29 30 PetscFunctionBegin; 31 ierr = MatMult(Na->A,v1,Na->w);CHKERRQ(ierr); 32 ierr = MatMultTransposeAdd(Na->A,Na->w,v2,v3);CHKERRQ(ierr); 33 PetscFunctionReturn(0); 34 } 35 36 #undef __FUNCT__ 37 #define __FUNCT__ "MatDestroy_Normal" 38 PetscErrorCode MatDestroy_Normal(Mat N) 39 { 40 Mat_Normal *Na = (Mat_Normal*)N->data; 41 PetscErrorCode ierr; 42 43 PetscFunctionBegin; 44 ierr = PetscObjectDereference((PetscObject)Na->A);CHKERRQ(ierr); 45 ierr = VecDestroy(Na->w);CHKERRQ(ierr); 46 ierr = PetscFree(Na);CHKERRQ(ierr); 47 PetscFunctionReturn(0); 48 } 49 50 /* 51 Slow, nonscalable version 52 */ 53 #undef __FUNCT__ 54 #define __FUNCT__ "MatGetDiagonal_Normal" 55 PetscErrorCode MatGetDiagonal_Normal(Mat N,Vec v) 56 { 57 Mat_Normal *Na = (Mat_Normal*)N->data; 58 Mat A = Na->A; 59 PetscErrorCode ierr; 60 PetscInt i,j,rstart,rend,nnz; 61 const PetscInt *cols; 62 PetscScalar *diag,*work,*values; 63 const PetscScalar *mvalues; 64 PetscMap cmap; 65 66 PetscFunctionBegin; 67 ierr = PetscMalloc(2*A->N*sizeof(PetscScalar),&diag);CHKERRQ(ierr); 68 work = diag + A->N; 69 ierr = PetscMemzero(work,A->N*sizeof(PetscScalar));CHKERRQ(ierr); 70 ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 71 for (i=rstart; i<rend; i++) { 72 ierr = MatGetRow(A,i,&nnz,&cols,&mvalues);CHKERRQ(ierr); 73 for (j=0; j<nnz; j++) { 74 work[cols[j]] += mvalues[j]*mvalues[j]; 75 } 76 ierr = MatRestoreRow(A,i,&nnz,&cols,&mvalues);CHKERRQ(ierr); 77 } 78 ierr = MPI_Allreduce(work,diag,A->N,MPIU_SCALAR,MPI_SUM,N->comm);CHKERRQ(ierr); 79 ierr = MatGetPetscMaps(A,PETSC_NULL,&cmap);CHKERRQ(ierr); 80 ierr = PetscMapGetLocalRange(cmap,&rstart,&rend);CHKERRQ(ierr); 81 ierr = VecGetArray(v,&values);CHKERRQ(ierr); 82 ierr = PetscMemcpy(values,diag+rstart,(rend-rstart)*sizeof(PetscScalar));CHKERRQ(ierr); 83 ierr = VecRestoreArray(v,&values);CHKERRQ(ierr); 84 ierr = PetscFree(diag);CHKERRQ(ierr); 85 PetscFunctionReturn(0); 86 } 87 88 #undef __FUNCT__ 89 #define __FUNCT__ "MatCreateNormal" 90 /*@ 91 MatCreateNormal - Creates a new matrix object that behaves like A'*A. 92 93 Collective on Mat 94 95 Input Parameter: 96 . A - the (possibly rectangular) matrix 97 98 Output Parameter: 99 . N - the matrix that represents A'*A 100 101 Level: intermediate 102 103 Notes: The product A'*A is NOT actually formed! Rather the new matrix 104 object performs the matrix-vector product by first multiplying by 105 A and then A' 106 @*/ 107 PetscErrorCode PETSCMAT_DLLEXPORT MatCreateNormal(Mat A,Mat *N) 108 { 109 PetscErrorCode ierr; 110 PetscInt m,n; 111 Mat_Normal *Na; 112 113 PetscFunctionBegin; 114 ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 115 ierr = MatCreate(A->comm,N);CHKERRQ(ierr); 116 ierr = MatSetSizes(*N,n,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 117 ierr = PetscObjectChangeTypeName((PetscObject)*N,MATNORMAL);CHKERRQ(ierr); 118 119 ierr = PetscNew(Mat_Normal,&Na);CHKERRQ(ierr); 120 Na->A = A; 121 ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 122 (*N)->data = (void*) Na; 123 124 ierr = VecCreateMPI(A->comm,m,PETSC_DECIDE,&Na->w);CHKERRQ(ierr); 125 (*N)->ops->destroy = MatDestroy_Normal; 126 (*N)->ops->mult = MatMult_Normal; 127 (*N)->ops->multadd = MatMultAdd_Normal; 128 (*N)->ops->getdiagonal = MatGetDiagonal_Normal; 129 (*N)->assembled = PETSC_TRUE; 130 (*N)->N = A->N; 131 (*N)->M = A->N; 132 (*N)->n = A->n; 133 (*N)->m = A->n; 134 PetscFunctionReturn(0); 135 } 136 137