1b5b17502SBarry Smith 2b5b17502SBarry Smith /* 3b5b17502SBarry Smith This is included by sbaij.c to generate unsigned short and regular versions of these two functions 4b5b17502SBarry Smith */ 554e8760dSRichard Tran Mills 654e8760dSRichard Tran Mills /* We cut-and-past below from aij.h to make a "no_function" version of PetscSparseDensePlusDot(). 754e8760dSRichard Tran Mills * This is necessary because the USESHORT case cannot use the inlined functions that may be employed. */ 854e8760dSRichard Tran Mills 954e8760dSRichard Tran Mills #if defined(PETSC_KERNEL_USE_UNROLL_4) 1054e8760dSRichard Tran Mills #define PetscSparseDensePlusDot_no_function(sum,r,xv,xi,nnz) { \ 1154e8760dSRichard Tran Mills if (nnz > 0) { \ 12a9ce9505SJunchao Zhang PetscInt nnz2=nnz,rem=nnz&0x3; \ 13a9ce9505SJunchao Zhang switch (rem) { \ 1454e8760dSRichard Tran Mills case 3: sum += *xv++ *r[*xi++]; \ 1554e8760dSRichard Tran Mills case 2: sum += *xv++ *r[*xi++]; \ 1654e8760dSRichard Tran Mills case 1: sum += *xv++ *r[*xi++]; \ 17a9ce9505SJunchao Zhang nnz2 -= rem;} \ 18a9ce9505SJunchao Zhang while (nnz2 > 0) { \ 1954e8760dSRichard Tran Mills sum += xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + \ 2054e8760dSRichard Tran Mills xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \ 21a9ce9505SJunchao Zhang xv += 4; xi += 4; nnz2 -= 4; \ 22a9ce9505SJunchao Zhang } \ 23a9ce9505SJunchao Zhang xv -= nnz; xi -= nnz; \ 24a9ce9505SJunchao Zhang } \ 25a9ce9505SJunchao Zhang } 2654e8760dSRichard Tran Mills 2754e8760dSRichard Tran Mills #elif defined(PETSC_KERNEL_USE_UNROLL_2) 2854e8760dSRichard Tran Mills #define PetscSparseDensePlusDot_no_function(sum,r,xv,xi,nnz) { \ 2954e8760dSRichard Tran Mills PetscInt __i,__i1,__i2; \ 3054e8760dSRichard Tran Mills for (__i=0; __i<nnz-1; __i+=2) {__i1 = xi[__i]; __i2=xi[__i+1]; \ 3154e8760dSRichard Tran Mills sum += (xv[__i]*r[__i1] + xv[__i+1]*r[__i2]);} \ 3254e8760dSRichard Tran Mills if (nnz & 0x1) sum += xv[__i] * r[xi[__i]];} 3354e8760dSRichard Tran Mills 3454e8760dSRichard Tran Mills #else 3554e8760dSRichard Tran Mills #define PetscSparseDensePlusDot_no_function(sum,r,xv,xi,nnz) { \ 3654e8760dSRichard Tran Mills PetscInt __i; \ 3754e8760dSRichard Tran Mills for (__i=0; __i<nnz; __i++) sum += xv[__i] * r[xi[__i]];} 3854e8760dSRichard Tran Mills #endif 3954e8760dSRichard Tran Mills 40018dd85eSSatish Balay #if defined(USESHORT) 41018dd85eSSatish Balay PetscErrorCode MatMult_SeqSBAIJ_1_ushort(Mat A,Vec xx,Vec zz) 42018dd85eSSatish Balay #else 43018dd85eSSatish Balay PetscErrorCode MatMult_SeqSBAIJ_1(Mat A,Vec xx,Vec zz) 44018dd85eSSatish Balay #endif 45b5b17502SBarry Smith { 46b5b17502SBarry Smith Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 47b5b17502SBarry Smith const PetscScalar *x; 48b5b17502SBarry Smith PetscScalar *z,x1,sum; 49b5b17502SBarry Smith const MatScalar *v; 50b5b17502SBarry Smith MatScalar vj; 51b5b17502SBarry Smith PetscErrorCode ierr; 52b5b17502SBarry Smith PetscInt mbs=a->mbs,i,j,nz; 53b5b17502SBarry Smith const PetscInt *ai=a->i; 54b5b17502SBarry Smith #if defined(USESHORT) 55b5b17502SBarry Smith const unsigned short *ib=a->jshort; 56b5b17502SBarry Smith unsigned short ibt; 57b5b17502SBarry Smith #else 58b5b17502SBarry Smith const PetscInt *ib=a->j; 59b5b17502SBarry Smith PetscInt ibt; 60b5b17502SBarry Smith #endif 618dca527aSShri Abhyankar PetscInt nonzerorow=0,jmin; 62eb1ec7c1SStefano Zampini #if defined(PETSC_USE_COMPLEX) 63eb1ec7c1SStefano Zampini const int aconj = A->hermitian; 64eb1ec7c1SStefano Zampini #else 65eb1ec7c1SStefano Zampini const int aconj = 0; 66eb1ec7c1SStefano Zampini #endif 67b5b17502SBarry Smith 68b5b17502SBarry Smith PetscFunctionBegin; 69b5b17502SBarry Smith ierr = VecSet(zz,0.0);CHKERRQ(ierr); 703649974fSBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 71b5b17502SBarry Smith ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 72b5b17502SBarry Smith 73b5b17502SBarry Smith v = a->a; 74b5b17502SBarry Smith for (i=0; i<mbs; i++) { 75b5b17502SBarry Smith nz = ai[i+1] - ai[i]; /* length of i_th row of A */ 7603d09022SShri Abhyankar if (!nz) continue; /* Move to the next row if the current row is empty */ 7703d09022SShri Abhyankar nonzerorow++; 788dca527aSShri Abhyankar sum = 0.0; 798dca527aSShri Abhyankar jmin = 0; 80b5b17502SBarry Smith x1 = x[i]; 818dca527aSShri Abhyankar if (ib[0] == i) { 82b5b17502SBarry Smith sum = v[0]*x1; /* diagonal term */ 838dca527aSShri Abhyankar jmin++; 848dca527aSShri Abhyankar } 85444d8c10SJed Brown PetscPrefetchBlock(ib+nz,nz,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */ 86444d8c10SJed Brown PetscPrefetchBlock(v+nz,nz,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */ 87eb1ec7c1SStefano Zampini if (aconj) { 88eb1ec7c1SStefano Zampini for (j=jmin; j<nz; j++) { 89eb1ec7c1SStefano Zampini ibt = ib[j]; 90eb1ec7c1SStefano Zampini vj = v[j]; 91eb1ec7c1SStefano Zampini z[ibt] += PetscConj(vj) * x1; /* (strict lower triangular part of A)*x */ 92eb1ec7c1SStefano Zampini sum += vj * x[ibt]; /* (strict upper triangular part of A)*x */ 93eb1ec7c1SStefano Zampini } 94eb1ec7c1SStefano Zampini } else { 958dca527aSShri Abhyankar for (j=jmin; j<nz; j++) { 96b5b17502SBarry Smith ibt = ib[j]; 97b5b17502SBarry Smith vj = v[j]; 98b5b17502SBarry Smith z[ibt] += vj * x1; /* (strict lower triangular part of A)*x */ 99b5b17502SBarry Smith sum += vj * x[ibt]; /* (strict upper triangular part of A)*x */ 100b5b17502SBarry Smith } 101eb1ec7c1SStefano Zampini } 102b5b17502SBarry Smith z[i] += sum; 103b5b17502SBarry Smith v += nz; 104b5b17502SBarry Smith ib += nz; 105b5b17502SBarry Smith } 106b5b17502SBarry Smith 1073649974fSBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 108b5b17502SBarry Smith ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 10903d09022SShri Abhyankar ierr = PetscLogFlops(2.0*(2.0*a->nz - nonzerorow) - nonzerorow);CHKERRQ(ierr); 110b5b17502SBarry Smith PetscFunctionReturn(0); 111b5b17502SBarry Smith } 112b5b17502SBarry Smith 113018dd85eSSatish Balay #if defined(USESHORT) 114018dd85eSSatish Balay PetscErrorCode MatSOR_SeqSBAIJ_ushort(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx) 115018dd85eSSatish Balay #else 116018dd85eSSatish Balay PetscErrorCode MatSOR_SeqSBAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx) 117018dd85eSSatish Balay #endif 118b5b17502SBarry Smith { 119b5b17502SBarry Smith Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 120b5b17502SBarry Smith const MatScalar *aa=a->a,*v,*v1,*aidiag; 121b5b17502SBarry Smith PetscScalar *x,*t,sum; 122b5b17502SBarry Smith const PetscScalar *b; 123b5b17502SBarry Smith MatScalar tmp; 124b5b17502SBarry Smith PetscErrorCode ierr; 125b5b17502SBarry Smith PetscInt m =a->mbs,bs=A->rmap->bs,j; 126b5b17502SBarry Smith const PetscInt *ai=a->i; 127b5b17502SBarry Smith #if defined(USESHORT) 128b5b17502SBarry Smith const unsigned short *aj=a->jshort,*vj,*vj1; 129b5b17502SBarry Smith #else 130b5b17502SBarry Smith const PetscInt *aj=a->j,*vj,*vj1; 131b5b17502SBarry Smith #endif 132b5b17502SBarry Smith PetscInt nz,nz1,i; 133b5b17502SBarry Smith 134b5b17502SBarry Smith PetscFunctionBegin; 135422a814eSBarry Smith if (fshift == -1.0) fshift = 0.0; /* negative fshift indicates do not error on zero diagonal; this code never errors on zero diagonal */ 136*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(flag & SOR_EISENSTAT,PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for Eisenstat"); 137b5b17502SBarry Smith 138b5b17502SBarry Smith its = its*lits; 139*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(its <= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %" PetscInt_FMT " and local its %" PetscInt_FMT " both positive",its,lits); 140b5b17502SBarry Smith 141*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(bs > 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented"); 142b5b17502SBarry Smith 143b5b17502SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 1443649974fSBarry Smith ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr); 145b5b17502SBarry Smith 146b5b17502SBarry Smith if (!a->idiagvalid) { 147b5b17502SBarry Smith if (!a->idiag) { 148785e854fSJed Brown ierr = PetscMalloc1(m,&a->idiag);CHKERRQ(ierr); 149b5b17502SBarry Smith } 150b5b17502SBarry Smith for (i=0; i<a->mbs; i++) a->idiag[i] = 1.0/a->a[a->i[i]]; 151b5b17502SBarry Smith a->idiagvalid = PETSC_TRUE; 152b5b17502SBarry Smith } 153b5b17502SBarry Smith 15441f059aeSBarry Smith if (!a->sor_work) { 155785e854fSJed Brown ierr = PetscMalloc1(m,&a->sor_work);CHKERRQ(ierr); 156b5b17502SBarry Smith } 15741f059aeSBarry Smith t = a->sor_work; 158b5b17502SBarry Smith 159b5b17502SBarry Smith aidiag = a->idiag; 160b5b17502SBarry Smith 161b5b17502SBarry Smith if (flag == SOR_APPLY_UPPER) { 162b5b17502SBarry Smith /* apply (U + D/omega) to the vector */ 163b5b17502SBarry Smith PetscScalar d; 164b5b17502SBarry Smith for (i=0; i<m; i++) { 165b5b17502SBarry Smith d = fshift + aa[ai[i]]; 166b5b17502SBarry Smith nz = ai[i+1] - ai[i] - 1; 167b5b17502SBarry Smith vj = aj + ai[i] + 1; 168b5b17502SBarry Smith v = aa + ai[i] + 1; 169b5b17502SBarry Smith sum = b[i]*d/omega; 17054e8760dSRichard Tran Mills #ifdef USESHORT 17154e8760dSRichard Tran Mills PetscSparseDensePlusDot_no_function(sum,b,v,vj,nz); 17254e8760dSRichard Tran Mills #else 173b5b17502SBarry Smith PetscSparseDensePlusDot(sum,b,v,vj,nz); 17454e8760dSRichard Tran Mills #endif 175b5b17502SBarry Smith x[i] = sum; 176b5b17502SBarry Smith } 177b5b17502SBarry Smith ierr = PetscLogFlops(a->nz);CHKERRQ(ierr); 178b5b17502SBarry Smith } 179b5b17502SBarry Smith 180b5b17502SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 181b5b17502SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) { 182580bdb30SBarry Smith ierr = PetscArraycpy(t,b,m);CHKERRQ(ierr); 183b5b17502SBarry Smith 184b5b17502SBarry Smith v = aa + 1; 185b5b17502SBarry Smith vj = aj + 1; 186b5b17502SBarry Smith for (i=0; i<m; i++) { 187b5b17502SBarry Smith nz = ai[i+1] - ai[i] - 1; 188b5b17502SBarry Smith tmp = -(x[i] = omega*t[i]*aidiag[i]); 18926fbe8dcSKarl Rupp for (j=0; j<nz; j++) t[vj[j]] += tmp*v[j]; 190b5b17502SBarry Smith v += nz + 1; 191b5b17502SBarry Smith vj += nz + 1; 192b5b17502SBarry Smith } 193ca0c957dSBarry Smith ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr); 194b5b17502SBarry Smith } 195b5b17502SBarry Smith 196b5b17502SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) { 197b5b17502SBarry Smith int nz2; 198b5b17502SBarry Smith if (!(flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP)) { 199b5b17502SBarry Smith #if defined(PETSC_USE_BACKWARD_LOOP) 200b5b17502SBarry Smith v = aa + ai[m] - 1; 201b5b17502SBarry Smith vj = aj + ai[m] - 1; 202b5b17502SBarry Smith for (i=m-1; i>=0; i--) { 203b5b17502SBarry Smith sum = b[i]; 204b5b17502SBarry Smith nz = ai[i+1] - ai[i] - 1; 205b5b17502SBarry Smith {PetscInt __i;for (__i=0; __i<nz; __i++) sum -= v[-__i] * x[vj[-__i]];} 206b5b17502SBarry Smith #else 207b5b17502SBarry Smith v = aa + ai[m-1] + 1; 208b5b17502SBarry Smith vj = aj + ai[m-1] + 1; 209b5b17502SBarry Smith nz = 0; 210b5b17502SBarry Smith for (i=m-1; i>=0; i--) { 211b5b17502SBarry Smith sum = b[i]; 2126c11f900SJed Brown nz2 = ai[i] - ai[PetscMax(i-1,0)] - 1; /* avoid referencing ai[-1], nonsense nz2 is okay on last iteration */ 21350d8bf02SJed Brown PETSC_Prefetch(v-nz2-1,0,PETSC_PREFETCH_HINT_NTA); 21450d8bf02SJed Brown PETSC_Prefetch(vj-nz2-1,0,PETSC_PREFETCH_HINT_NTA); 215b5b17502SBarry Smith PetscSparseDenseMinusDot(sum,x,v,vj,nz); 216b5b17502SBarry Smith nz = nz2; 217b5b17502SBarry Smith #endif 218b5b17502SBarry Smith x[i] = omega*sum*aidiag[i]; 219b5b17502SBarry Smith v -= nz + 1; 220b5b17502SBarry Smith vj -= nz + 1; 221b5b17502SBarry Smith } 222ca0c957dSBarry Smith ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr); 223b5b17502SBarry Smith } else { 224b5b17502SBarry Smith v = aa + ai[m-1] + 1; 225b5b17502SBarry Smith vj = aj + ai[m-1] + 1; 226b5b17502SBarry Smith nz = 0; 227b5b17502SBarry Smith for (i=m-1; i>=0; i--) { 228b5b17502SBarry Smith sum = t[i]; 2296c11f900SJed Brown nz2 = ai[i] - ai[PetscMax(i-1,0)] - 1; /* avoid referencing ai[-1], nonsense nz2 is okay on last iteration */ 23050d8bf02SJed Brown PETSC_Prefetch(v-nz2-1,0,PETSC_PREFETCH_HINT_NTA); 23150d8bf02SJed Brown PETSC_Prefetch(vj-nz2-1,0,PETSC_PREFETCH_HINT_NTA); 232b5b17502SBarry Smith PetscSparseDenseMinusDot(sum,x,v,vj,nz); 233b5b17502SBarry Smith x[i] = (1-omega)*x[i] + omega*sum*aidiag[i]; 234b5b17502SBarry Smith nz = nz2; 235b5b17502SBarry Smith v -= nz + 1; 236b5b17502SBarry Smith vj -= nz + 1; 237b5b17502SBarry Smith } 238ca0c957dSBarry Smith ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr); 239b5b17502SBarry Smith } 240b5b17502SBarry Smith } 241b5b17502SBarry Smith its--; 242b5b17502SBarry Smith } 243b5b17502SBarry Smith 244b5b17502SBarry Smith while (its--) { 245b5b17502SBarry Smith /* 246b5b17502SBarry Smith forward sweep: 247b5b17502SBarry Smith for i=0,...,m-1: 248b5b17502SBarry Smith sum[i] = (b[i] - U(i,:)x)/d[i]; 249b5b17502SBarry Smith x[i] = (1-omega)x[i] + omega*sum[i]; 250b5b17502SBarry Smith b = b - x[i]*U^T(i,:); 251b5b17502SBarry Smith 252b5b17502SBarry Smith */ 253b5b17502SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) { 254580bdb30SBarry Smith ierr = PetscArraycpy(t,b,m);CHKERRQ(ierr); 255b5b17502SBarry Smith 256b5b17502SBarry Smith for (i=0; i<m; i++) { 257b5b17502SBarry Smith v = aa + ai[i] + 1; v1=v; 258b5b17502SBarry Smith vj = aj + ai[i] + 1; vj1=vj; 259b5b17502SBarry Smith nz = ai[i+1] - ai[i] - 1; nz1=nz; 260b5b17502SBarry Smith sum = t[i]; 261b5b17502SBarry Smith while (nz1--) sum -= (*v1++)*x[*vj1++]; 262b5b17502SBarry Smith x[i] = (1-omega)*x[i] + omega*sum*aidiag[i]; 263b5b17502SBarry Smith while (nz--) t[*vj++] -= x[i]*(*v++); 264b5b17502SBarry Smith } 265dcca0238SPatrick Lacasse ierr = PetscLogFlops(4.0*a->nz);CHKERRQ(ierr); 266b5b17502SBarry Smith } 267b5b17502SBarry Smith 268b5b17502SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) { 269b5b17502SBarry Smith /* 270b5b17502SBarry Smith backward sweep: 271b5b17502SBarry Smith b = b - x[i]*U^T(i,:), i=0,...,n-2 272b5b17502SBarry Smith for i=m-1,...,0: 273b5b17502SBarry Smith sum[i] = (b[i] - U(i,:)x)/d[i]; 274b5b17502SBarry Smith x[i] = (1-omega)x[i] + omega*sum[i]; 275b5b17502SBarry Smith */ 276b5b17502SBarry Smith /* if there was a forward sweep done above then I thing the next two for loops are not needed */ 277580bdb30SBarry Smith ierr = PetscArraycpy(t,b,m);CHKERRQ(ierr); 278b5b17502SBarry Smith 279b5b17502SBarry Smith for (i=0; i<m-1; i++) { /* update rhs */ 280b5b17502SBarry Smith v = aa + ai[i] + 1; 281b5b17502SBarry Smith vj = aj + ai[i] + 1; 282b5b17502SBarry Smith nz = ai[i+1] - ai[i] - 1; 283b5b17502SBarry Smith while (nz--) t[*vj++] -= x[i]*(*v++); 284b5b17502SBarry Smith } 285dcca0238SPatrick Lacasse ierr = PetscLogFlops(2.0*(a->nz - m));CHKERRQ(ierr); 286b5b17502SBarry Smith for (i=m-1; i>=0; i--) { 287b5b17502SBarry Smith v = aa + ai[i] + 1; 288b5b17502SBarry Smith vj = aj + ai[i] + 1; 289b5b17502SBarry Smith nz = ai[i+1] - ai[i] - 1; 290b5b17502SBarry Smith sum = t[i]; 291b5b17502SBarry Smith while (nz--) sum -= x[*vj++]*(*v++); 292b5b17502SBarry Smith x[i] = (1-omega)*x[i] + omega*sum*aidiag[i]; 293b5b17502SBarry Smith } 294dcca0238SPatrick Lacasse ierr = PetscLogFlops(2.0*(a->nz + m));CHKERRQ(ierr); 295b5b17502SBarry Smith } 296b5b17502SBarry Smith } 297b5b17502SBarry Smith 298b5b17502SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2993649974fSBarry Smith ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr); 300b5b17502SBarry Smith PetscFunctionReturn(0); 301b5b17502SBarry Smith } 302