1 /*$Id: mmaij.c,v 1.59 2001/08/07 03:02:49 balay Exp $*/ 2 3 /* 4 Support for the parallel AIJ matrix vector multiply 5 */ 6 #include "src/mat/impls/aij/mpi/mpiaij.h" 7 8 #undef __FUNCT__ 9 #define __FUNCT__ "MatSetUpMultiply_MPIAIJ" 10 int MatSetUpMultiply_MPIAIJ(Mat mat) 11 { 12 Mat_MPIAIJ *aij = (Mat_MPIAIJ*)mat->data; 13 Mat_SeqAIJ *B = (Mat_SeqAIJ*)(aij->B->data); 14 int i,j,*aj = B->j,ierr,ec = 0,*garray; 15 IS from,to; 16 Vec gvec; 17 #if defined (PETSC_USE_CTABLE) 18 PetscTable gid1_lid1; 19 PetscTablePosition tpos; 20 int gid,lid; 21 #else 22 int N = mat->N,*indices; 23 24 #endif 25 26 PetscFunctionBegin; 27 28 #if defined (PETSC_USE_CTABLE) 29 /* use a table - Mark Adams (this has not been tested with "shift") */ 30 ierr = PetscTableCreate(aij->B->m,&gid1_lid1);CHKERRQ(ierr); 31 for (i=0; i<aij->B->m; i++) { 32 for (j=0; j<B->ilen[i]; j++) { 33 int data,gid1 = aj[B->i[i] + j] + 1; 34 ierr = PetscTableFind(gid1_lid1,gid1,&data);CHKERRQ(ierr); 35 if (!data) { 36 /* one based table */ 37 ierr = PetscTableAdd(gid1_lid1,gid1,++ec);CHKERRQ(ierr); 38 } 39 } 40 } 41 /* form array of columns we need */ 42 ierr = PetscMalloc((ec+1)*sizeof(int),&garray);CHKERRQ(ierr); 43 ierr = PetscTableGetHeadPosition(gid1_lid1,&tpos);CHKERRQ(ierr); 44 while (tpos) { 45 ierr = PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);CHKERRQ(ierr); 46 gid--; 47 lid--; 48 garray[lid] = gid; 49 } 50 ierr = PetscSortInt(ec,garray);CHKERRQ(ierr); /* sort, and rebuild */ 51 ierr = PetscTableRemoveAll(gid1_lid1);CHKERRQ(ierr); 52 for (i=0; i<ec; i++) { 53 ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1);CHKERRQ(ierr); 54 } 55 /* compact out the extra columns in B */ 56 for (i=0; i<aij->B->m; i++) { 57 for (j=0; j<B->ilen[i]; j++) { 58 int gid1 = aj[B->i[i] + j] + 1; 59 ierr = PetscTableFind(gid1_lid1,gid1,&lid);CHKERRQ(ierr); 60 lid --; 61 aj[B->i[i] + j] = lid; 62 } 63 } 64 aij->B->n = aij->B->N = ec; 65 ierr = PetscTableDelete(gid1_lid1);CHKERRQ(ierr); 66 /* Mark Adams */ 67 #else 68 /* For the first stab we make an array as long as the number of columns */ 69 /* mark those columns that are in aij->B */ 70 ierr = PetscMalloc((N+1)*sizeof(int),&indices);CHKERRQ(ierr); 71 ierr = PetscMemzero(indices,N*sizeof(int));CHKERRQ(ierr); 72 for (i=0; i<aij->B->m; i++) { 73 for (j=0; j<B->ilen[i]; j++) { 74 if (!indices[aj[B->i[i] + j] ]) ec++; 75 indices[aj[B->i[i] + j] ] = 1; 76 } 77 } 78 79 /* form array of columns we need */ 80 ierr = PetscMalloc((ec+1)*sizeof(int),&garray);CHKERRQ(ierr); 81 ec = 0; 82 for (i=0; i<N; i++) { 83 if (indices[i]) garray[ec++] = i; 84 } 85 86 /* make indices now point into garray */ 87 for (i=0; i<ec; i++) { 88 indices[garray[i]] = i; 89 } 90 91 /* compact out the extra columns in B */ 92 for (i=0; i<aij->B->m; i++) { 93 for (j=0; j<B->ilen[i]; j++) { 94 aj[B->i[i] + j] = indices[aj[B->i[i] + j]]; 95 } 96 } 97 aij->B->n = aij->B->N = ec; 98 ierr = PetscFree(indices);CHKERRQ(ierr); 99 #endif 100 /* create local vector that is used to scatter into */ 101 ierr = VecCreateSeq(PETSC_COMM_SELF,ec,&aij->lvec);CHKERRQ(ierr); 102 103 /* create two temporary Index sets for build scatter gather */ 104 ierr = ISCreateGeneral(mat->comm,ec,garray,&from);CHKERRQ(ierr); 105 ierr = ISCreateStride(PETSC_COMM_SELF,ec,0,1,&to);CHKERRQ(ierr); 106 107 /* create temporary global vector to generate scatter context */ 108 /* this is inefficient, but otherwise we must do either 109 1) save garray until the first actual scatter when the vector is known or 110 2) have another way of generating a scatter context without a vector.*/ 111 ierr = VecCreateMPI(mat->comm,mat->n,mat->N,&gvec);CHKERRQ(ierr); 112 113 /* generate the scatter context */ 114 ierr = VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx);CHKERRQ(ierr); 115 PetscLogObjectParent(mat,aij->Mvctx); 116 PetscLogObjectParent(mat,aij->lvec); 117 PetscLogObjectParent(mat,from); 118 PetscLogObjectParent(mat,to); 119 aij->garray = garray; 120 PetscLogObjectMemory(mat,(ec+1)*sizeof(int)); 121 ierr = ISDestroy(from);CHKERRQ(ierr); 122 ierr = ISDestroy(to);CHKERRQ(ierr); 123 ierr = VecDestroy(gvec);CHKERRQ(ierr); 124 PetscFunctionReturn(0); 125 } 126 127 128 #undef __FUNCT__ 129 #define __FUNCT__ "DisAssemble_MPIAIJ" 130 /* 131 Takes the local part of an already assembled MPIAIJ matrix 132 and disassembles it. This is to allow new nonzeros into the matrix 133 that require more communication in the matrix vector multiply. 134 Thus certain data-structures must be rebuilt. 135 136 Kind of slow! But that's what application programmers get when 137 they are sloppy. 138 */ 139 int DisAssemble_MPIAIJ(Mat A) 140 { 141 Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 142 Mat B = aij->B,Bnew; 143 Mat_SeqAIJ *Baij = (Mat_SeqAIJ*)B->data; 144 int ierr,i,j,m = B->m,n = A->N,col,ct = 0,*garray = aij->garray; 145 int *nz,ec; 146 PetscScalar v; 147 148 PetscFunctionBegin; 149 /* free stuff related to matrix-vec multiply */ 150 ierr = VecGetSize(aij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */ 151 ierr = VecDestroy(aij->lvec);CHKERRQ(ierr); aij->lvec = 0; 152 ierr = VecScatterDestroy(aij->Mvctx);CHKERRQ(ierr); aij->Mvctx = 0; 153 if (aij->colmap) { 154 #if defined (PETSC_USE_CTABLE) 155 ierr = PetscTableDelete(aij->colmap);CHKERRQ(ierr); 156 aij->colmap = 0; 157 #else 158 ierr = PetscFree(aij->colmap);CHKERRQ(ierr); 159 aij->colmap = 0; 160 PetscLogObjectMemory(A,-aij->B->n*sizeof(int)); 161 #endif 162 } 163 164 /* make sure that B is assembled so we can access its values */ 165 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 166 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 167 168 /* invent new B and copy stuff over */ 169 ierr = PetscMalloc((m+1)*sizeof(int),&nz);CHKERRQ(ierr); 170 for (i=0; i<m; i++) { 171 nz[i] = Baij->i[i+1] - Baij->i[i]; 172 } 173 ierr = MatCreate(PETSC_COMM_SELF,m,n,m,n,&Bnew);CHKERRQ(ierr); 174 ierr = MatSetType(Bnew,B->type_name);CHKERRQ(ierr); 175 ierr = MatSeqAIJSetPreallocation(Bnew,0,nz);CHKERRQ(ierr); 176 ierr = PetscFree(nz);CHKERRQ(ierr); 177 for (i=0; i<m; i++) { 178 for (j=Baij->i[i]; j<Baij->i[i+1]; j++) { 179 col = garray[Baij->j[ct]]; 180 v = Baij->a[ct++]; 181 ierr = MatSetValues(Bnew,1,&i,1,&col,&v,B->insertmode);CHKERRQ(ierr); 182 } 183 } 184 ierr = PetscFree(aij->garray);CHKERRQ(ierr); 185 aij->garray = 0; 186 PetscLogObjectMemory(A,-ec*sizeof(int)); 187 ierr = MatDestroy(B);CHKERRQ(ierr); 188 PetscLogObjectParent(A,Bnew); 189 aij->B = Bnew; 190 A->was_assembled = PETSC_FALSE; 191 PetscFunctionReturn(0); 192 } 193 194 /* ugly stuff added for Glenn someday we should fix this up */ 195 196 static int *auglyrmapd = 0,*auglyrmapo = 0; /* mapping from the local ordering to the "diagonal" and "off-diagonal" 197 parts of the local matrix */ 198 static Vec auglydd = 0,auglyoo = 0; /* work vectors used to scale the two parts of the local matrix */ 199 200 201 #undef __FUNCT__ 202 #define __FUNCT__ "MatMPIAIJDiagonalScaleLocalSetUp" 203 int MatMPIAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale) 204 { 205 Mat_MPIAIJ *ina = (Mat_MPIAIJ*) inA->data; /*access private part of matrix */ 206 int ierr,i,n,nt,cstart,cend,no,*garray = ina->garray,*lindices; 207 int *r_rmapd,*r_rmapo; 208 209 PetscFunctionBegin; 210 ierr = MatGetOwnershipRange(inA,&cstart,&cend);CHKERRQ(ierr); 211 ierr = MatGetSize(ina->A,PETSC_NULL,&n);CHKERRQ(ierr); 212 ierr = PetscMalloc((inA->mapping->n+1)*sizeof(int),&r_rmapd);CHKERRQ(ierr); 213 ierr = PetscMemzero(r_rmapd,inA->mapping->n*sizeof(int));CHKERRQ(ierr); 214 nt = 0; 215 for (i=0; i<inA->mapping->n; i++) { 216 if (inA->mapping->indices[i] >= cstart && inA->mapping->indices[i] < cend) { 217 nt++; 218 r_rmapd[i] = inA->mapping->indices[i] + 1; 219 } 220 } 221 if (nt != n) SETERRQ2(1,"Hmm nt %d n %d",nt,n); 222 ierr = PetscMalloc((n+1)*sizeof(int),&auglyrmapd);CHKERRQ(ierr); 223 for (i=0; i<inA->mapping->n; i++) { 224 if (r_rmapd[i]){ 225 auglyrmapd[(r_rmapd[i]-1)-cstart] = i; 226 } 227 } 228 ierr = PetscFree(r_rmapd);CHKERRQ(ierr); 229 ierr = VecCreateSeq(PETSC_COMM_SELF,n,&auglydd);CHKERRQ(ierr); 230 231 ierr = PetscMalloc((inA->N+1)*sizeof(int),&lindices);CHKERRQ(ierr); 232 ierr = PetscMemzero(lindices,inA->N*sizeof(int));CHKERRQ(ierr); 233 for (i=0; i<ina->B->n; i++) { 234 lindices[garray[i]] = i+1; 235 } 236 no = inA->mapping->n - nt; 237 ierr = PetscMalloc((inA->mapping->n+1)*sizeof(int),&r_rmapo);CHKERRQ(ierr); 238 ierr = PetscMemzero(r_rmapo,inA->mapping->n*sizeof(int));CHKERRQ(ierr); 239 nt = 0; 240 for (i=0; i<inA->mapping->n; i++) { 241 if (lindices[inA->mapping->indices[i]]) { 242 nt++; 243 r_rmapo[i] = lindices[inA->mapping->indices[i]]; 244 } 245 } 246 if (nt > no) SETERRQ2(1,"Hmm nt %d no %d",nt,n); 247 ierr = PetscFree(lindices);CHKERRQ(ierr); 248 ierr = PetscMalloc((nt+1)*sizeof(int),&auglyrmapo);CHKERRQ(ierr); 249 for (i=0; i<inA->mapping->n; i++) { 250 if (r_rmapo[i]){ 251 auglyrmapo[(r_rmapo[i]-1)] = i; 252 } 253 } 254 ierr = PetscFree(r_rmapo);CHKERRQ(ierr); 255 ierr = VecCreateSeq(PETSC_COMM_SELF,nt,&auglyoo);CHKERRQ(ierr); 256 257 PetscFunctionReturn(0); 258 } 259 260 #undef __FUNCT__ 261 #define __FUNCT__ "MatMPIAIJDiagonalScaleLocal" 262 int MatMPIAIJDiagonalScaleLocal(Mat A,Vec scale) 263 { 264 /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */ 265 int ierr,(*f)(Mat,Vec); 266 267 PetscFunctionBegin; 268 ierr = PetscObjectQueryFunction((PetscObject)A,"MatDiagonalScaleLocal_C",(void (**)(void))&f);CHKERRQ(ierr); 269 if (f) { 270 ierr = (*f)(A,scale);CHKERRQ(ierr); 271 } 272 PetscFunctionReturn(0); 273 } 274 275 EXTERN_C_BEGIN 276 #undef __FUNCT__ 277 #define __FUNCT__ "MatDiagonalScaleLocal_MPIAIJ" 278 int MatDiagonalScaleLocal_MPIAIJ(Mat A,Vec scale) 279 { 280 Mat_MPIAIJ *a = (Mat_MPIAIJ*) A->data; /*access private part of matrix */ 281 int ierr,n,i; 282 PetscScalar *d,*o,*s; 283 284 PetscFunctionBegin; 285 if (!auglyrmapd) { 286 ierr = MatMPIAIJDiagonalScaleLocalSetUp(A,scale);CHKERRQ(ierr); 287 } 288 289 ierr = VecGetArray(scale,&s);CHKERRQ(ierr); 290 291 ierr = VecGetLocalSize(auglydd,&n);CHKERRQ(ierr); 292 ierr = VecGetArray(auglydd,&d);CHKERRQ(ierr); 293 for (i=0; i<n; i++) { 294 d[i] = s[auglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */ 295 } 296 ierr = VecRestoreArray(auglydd,&d);CHKERRQ(ierr); 297 /* column scale "diagonal" portion of local matrix */ 298 ierr = MatDiagonalScale(a->A,PETSC_NULL,auglydd);CHKERRQ(ierr); 299 300 ierr = VecGetLocalSize(auglyoo,&n);CHKERRQ(ierr); 301 ierr = VecGetArray(auglyoo,&o);CHKERRQ(ierr); 302 for (i=0; i<n; i++) { 303 o[i] = s[auglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */ 304 } 305 ierr = VecRestoreArray(scale,&s);CHKERRQ(ierr); 306 ierr = VecRestoreArray(auglyoo,&o);CHKERRQ(ierr); 307 /* column scale "off-diagonal" portion of local matrix */ 308 ierr = MatDiagonalScale(a->B,PETSC_NULL,auglyoo);CHKERRQ(ierr); 309 310 PetscFunctionReturn(0); 311 } 312 EXTERN_C_END 313 314 315