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 PetscInt mbs=a->mbs,i,j,nz; 52b5b17502SBarry Smith const PetscInt *ai=a->i; 53b5b17502SBarry Smith #if defined(USESHORT) 54b5b17502SBarry Smith const unsigned short *ib=a->jshort; 55b5b17502SBarry Smith unsigned short ibt; 56b5b17502SBarry Smith #else 57b5b17502SBarry Smith const PetscInt *ib=a->j; 58b5b17502SBarry Smith PetscInt ibt; 59b5b17502SBarry Smith #endif 608dca527aSShri Abhyankar PetscInt nonzerorow=0,jmin; 61eb1ec7c1SStefano Zampini #if defined(PETSC_USE_COMPLEX) 62eb1ec7c1SStefano Zampini const int aconj = A->hermitian; 63eb1ec7c1SStefano Zampini #else 64eb1ec7c1SStefano Zampini const int aconj = 0; 65eb1ec7c1SStefano Zampini #endif 66b5b17502SBarry Smith 67b5b17502SBarry Smith PetscFunctionBegin; 689566063dSJacob Faibussowitsch PetscCall(VecSet(zz,0.0)); 699566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(xx,&x)); 709566063dSJacob Faibussowitsch PetscCall(VecGetArray(zz,&z)); 71b5b17502SBarry Smith 72b5b17502SBarry Smith v = a->a; 73b5b17502SBarry Smith for (i=0; i<mbs; i++) { 74b5b17502SBarry Smith nz = ai[i+1] - ai[i]; /* length of i_th row of A */ 7503d09022SShri Abhyankar if (!nz) continue; /* Move to the next row if the current row is empty */ 7603d09022SShri Abhyankar nonzerorow++; 778dca527aSShri Abhyankar sum = 0.0; 788dca527aSShri Abhyankar jmin = 0; 79b5b17502SBarry Smith x1 = x[i]; 808dca527aSShri Abhyankar if (ib[0] == i) { 81b5b17502SBarry Smith sum = v[0]*x1; /* diagonal term */ 828dca527aSShri Abhyankar jmin++; 838dca527aSShri Abhyankar } 84444d8c10SJed Brown PetscPrefetchBlock(ib+nz,nz,0,PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */ 85444d8c10SJed Brown PetscPrefetchBlock(v+nz,nz,0,PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */ 86eb1ec7c1SStefano Zampini if (aconj) { 87eb1ec7c1SStefano Zampini for (j=jmin; j<nz; j++) { 88eb1ec7c1SStefano Zampini ibt = ib[j]; 89eb1ec7c1SStefano Zampini vj = v[j]; 90eb1ec7c1SStefano Zampini z[ibt] += PetscConj(vj) * x1; /* (strict lower triangular part of A)*x */ 91eb1ec7c1SStefano Zampini sum += vj * x[ibt]; /* (strict upper triangular part of A)*x */ 92eb1ec7c1SStefano Zampini } 93eb1ec7c1SStefano Zampini } else { 948dca527aSShri Abhyankar for (j=jmin; j<nz; j++) { 95b5b17502SBarry Smith ibt = ib[j]; 96b5b17502SBarry Smith vj = v[j]; 97b5b17502SBarry Smith z[ibt] += vj * x1; /* (strict lower triangular part of A)*x */ 98b5b17502SBarry Smith sum += vj * x[ibt]; /* (strict upper triangular part of A)*x */ 99b5b17502SBarry Smith } 100eb1ec7c1SStefano Zampini } 101b5b17502SBarry Smith z[i] += sum; 102b5b17502SBarry Smith v += nz; 103b5b17502SBarry Smith ib += nz; 104b5b17502SBarry Smith } 105b5b17502SBarry Smith 1069566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(xx,&x)); 1079566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(zz,&z)); 1089566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0*(2.0*a->nz - nonzerorow) - nonzerorow)); 109b5b17502SBarry Smith PetscFunctionReturn(0); 110b5b17502SBarry Smith } 111b5b17502SBarry Smith 112018dd85eSSatish Balay #if defined(USESHORT) 113018dd85eSSatish Balay PetscErrorCode MatSOR_SeqSBAIJ_ushort(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx) 114018dd85eSSatish Balay #else 115018dd85eSSatish Balay PetscErrorCode MatSOR_SeqSBAIJ(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx) 116018dd85eSSatish Balay #endif 117b5b17502SBarry Smith { 118b5b17502SBarry Smith Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 119b5b17502SBarry Smith const MatScalar *aa=a->a,*v,*v1,*aidiag; 120b5b17502SBarry Smith PetscScalar *x,*t,sum; 121b5b17502SBarry Smith const PetscScalar *b; 122b5b17502SBarry Smith MatScalar tmp; 123b5b17502SBarry Smith PetscInt m =a->mbs,bs=A->rmap->bs,j; 124b5b17502SBarry Smith const PetscInt *ai=a->i; 125b5b17502SBarry Smith #if defined(USESHORT) 126b5b17502SBarry Smith const unsigned short *aj=a->jshort,*vj,*vj1; 127b5b17502SBarry Smith #else 128b5b17502SBarry Smith const PetscInt *aj=a->j,*vj,*vj1; 129b5b17502SBarry Smith #endif 130b5b17502SBarry Smith PetscInt nz,nz1,i; 131b5b17502SBarry Smith 132b5b17502SBarry Smith PetscFunctionBegin; 133422a814eSBarry Smith if (fshift == -1.0) fshift = 0.0; /* negative fshift indicates do not error on zero diagonal; this code never errors on zero diagonal */ 1342c71b3e2SJacob Faibussowitsch PetscCheckFalse(flag & SOR_EISENSTAT,PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for Eisenstat"); 135b5b17502SBarry Smith 136b5b17502SBarry Smith its = its*lits; 137*08401ef6SPierre Jolivet PetscCheck(its > 0,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %" PetscInt_FMT " and local its %" PetscInt_FMT " both positive",its,lits); 138b5b17502SBarry Smith 139*08401ef6SPierre Jolivet PetscCheck(bs <= 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented"); 140b5b17502SBarry Smith 1419566063dSJacob Faibussowitsch PetscCall(VecGetArray(xx,&x)); 1429566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(bb,&b)); 143b5b17502SBarry Smith 144b5b17502SBarry Smith if (!a->idiagvalid) { 145b5b17502SBarry Smith if (!a->idiag) { 1469566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m,&a->idiag)); 147b5b17502SBarry Smith } 148b5b17502SBarry Smith for (i=0; i<a->mbs; i++) a->idiag[i] = 1.0/a->a[a->i[i]]; 149b5b17502SBarry Smith a->idiagvalid = PETSC_TRUE; 150b5b17502SBarry Smith } 151b5b17502SBarry Smith 15241f059aeSBarry Smith if (!a->sor_work) { 1539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m,&a->sor_work)); 154b5b17502SBarry Smith } 15541f059aeSBarry Smith t = a->sor_work; 156b5b17502SBarry Smith 157b5b17502SBarry Smith aidiag = a->idiag; 158b5b17502SBarry Smith 159b5b17502SBarry Smith if (flag == SOR_APPLY_UPPER) { 160b5b17502SBarry Smith /* apply (U + D/omega) to the vector */ 161b5b17502SBarry Smith PetscScalar d; 162b5b17502SBarry Smith for (i=0; i<m; i++) { 163b5b17502SBarry Smith d = fshift + aa[ai[i]]; 164b5b17502SBarry Smith nz = ai[i+1] - ai[i] - 1; 165b5b17502SBarry Smith vj = aj + ai[i] + 1; 166b5b17502SBarry Smith v = aa + ai[i] + 1; 167b5b17502SBarry Smith sum = b[i]*d/omega; 16854e8760dSRichard Tran Mills #ifdef USESHORT 16954e8760dSRichard Tran Mills PetscSparseDensePlusDot_no_function(sum,b,v,vj,nz); 17054e8760dSRichard Tran Mills #else 171b5b17502SBarry Smith PetscSparseDensePlusDot(sum,b,v,vj,nz); 17254e8760dSRichard Tran Mills #endif 173b5b17502SBarry Smith x[i] = sum; 174b5b17502SBarry Smith } 1759566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(a->nz)); 176b5b17502SBarry Smith } 177b5b17502SBarry Smith 178b5b17502SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 179b5b17502SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) { 1809566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t,b,m)); 181b5b17502SBarry Smith 182b5b17502SBarry Smith v = aa + 1; 183b5b17502SBarry Smith vj = aj + 1; 184b5b17502SBarry Smith for (i=0; i<m; i++) { 185b5b17502SBarry Smith nz = ai[i+1] - ai[i] - 1; 186b5b17502SBarry Smith tmp = -(x[i] = omega*t[i]*aidiag[i]); 18726fbe8dcSKarl Rupp for (j=0; j<nz; j++) t[vj[j]] += tmp*v[j]; 188b5b17502SBarry Smith v += nz + 1; 189b5b17502SBarry Smith vj += nz + 1; 190b5b17502SBarry Smith } 1919566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0*a->nz)); 192b5b17502SBarry Smith } 193b5b17502SBarry Smith 194b5b17502SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) { 195b5b17502SBarry Smith int nz2; 196b5b17502SBarry Smith if (!(flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP)) { 197b5b17502SBarry Smith #if defined(PETSC_USE_BACKWARD_LOOP) 198b5b17502SBarry Smith v = aa + ai[m] - 1; 199b5b17502SBarry Smith vj = aj + ai[m] - 1; 200b5b17502SBarry Smith for (i=m-1; i>=0; i--) { 201b5b17502SBarry Smith sum = b[i]; 202b5b17502SBarry Smith nz = ai[i+1] - ai[i] - 1; 203b5b17502SBarry Smith {PetscInt __i;for (__i=0; __i<nz; __i++) sum -= v[-__i] * x[vj[-__i]];} 204b5b17502SBarry Smith #else 205b5b17502SBarry Smith v = aa + ai[m-1] + 1; 206b5b17502SBarry Smith vj = aj + ai[m-1] + 1; 207b5b17502SBarry Smith nz = 0; 208b5b17502SBarry Smith for (i=m-1; i>=0; i--) { 209b5b17502SBarry Smith sum = b[i]; 2106c11f900SJed Brown nz2 = ai[i] - ai[PetscMax(i-1,0)] - 1; /* avoid referencing ai[-1], nonsense nz2 is okay on last iteration */ 21150d8bf02SJed Brown PETSC_Prefetch(v-nz2-1,0,PETSC_PREFETCH_HINT_NTA); 21250d8bf02SJed Brown PETSC_Prefetch(vj-nz2-1,0,PETSC_PREFETCH_HINT_NTA); 213b5b17502SBarry Smith PetscSparseDenseMinusDot(sum,x,v,vj,nz); 214b5b17502SBarry Smith nz = nz2; 215b5b17502SBarry Smith #endif 216b5b17502SBarry Smith x[i] = omega*sum*aidiag[i]; 217b5b17502SBarry Smith v -= nz + 1; 218b5b17502SBarry Smith vj -= nz + 1; 219b5b17502SBarry Smith } 2209566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0*a->nz)); 221b5b17502SBarry Smith } else { 222b5b17502SBarry Smith v = aa + ai[m-1] + 1; 223b5b17502SBarry Smith vj = aj + ai[m-1] + 1; 224b5b17502SBarry Smith nz = 0; 225b5b17502SBarry Smith for (i=m-1; i>=0; i--) { 226b5b17502SBarry Smith sum = t[i]; 2276c11f900SJed Brown nz2 = ai[i] - ai[PetscMax(i-1,0)] - 1; /* avoid referencing ai[-1], nonsense nz2 is okay on last iteration */ 22850d8bf02SJed Brown PETSC_Prefetch(v-nz2-1,0,PETSC_PREFETCH_HINT_NTA); 22950d8bf02SJed Brown PETSC_Prefetch(vj-nz2-1,0,PETSC_PREFETCH_HINT_NTA); 230b5b17502SBarry Smith PetscSparseDenseMinusDot(sum,x,v,vj,nz); 231b5b17502SBarry Smith x[i] = (1-omega)*x[i] + omega*sum*aidiag[i]; 232b5b17502SBarry Smith nz = nz2; 233b5b17502SBarry Smith v -= nz + 1; 234b5b17502SBarry Smith vj -= nz + 1; 235b5b17502SBarry Smith } 2369566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0*a->nz)); 237b5b17502SBarry Smith } 238b5b17502SBarry Smith } 239b5b17502SBarry Smith its--; 240b5b17502SBarry Smith } 241b5b17502SBarry Smith 242b5b17502SBarry Smith while (its--) { 243b5b17502SBarry Smith /* 244b5b17502SBarry Smith forward sweep: 245b5b17502SBarry Smith for i=0,...,m-1: 246b5b17502SBarry Smith sum[i] = (b[i] - U(i,:)x)/d[i]; 247b5b17502SBarry Smith x[i] = (1-omega)x[i] + omega*sum[i]; 248b5b17502SBarry Smith b = b - x[i]*U^T(i,:); 249b5b17502SBarry Smith 250b5b17502SBarry Smith */ 251b5b17502SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) { 2529566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t,b,m)); 253b5b17502SBarry Smith 254b5b17502SBarry Smith for (i=0; i<m; i++) { 255b5b17502SBarry Smith v = aa + ai[i] + 1; v1=v; 256b5b17502SBarry Smith vj = aj + ai[i] + 1; vj1=vj; 257b5b17502SBarry Smith nz = ai[i+1] - ai[i] - 1; nz1=nz; 258b5b17502SBarry Smith sum = t[i]; 259b5b17502SBarry Smith while (nz1--) sum -= (*v1++)*x[*vj1++]; 260b5b17502SBarry Smith x[i] = (1-omega)*x[i] + omega*sum*aidiag[i]; 261b5b17502SBarry Smith while (nz--) t[*vj++] -= x[i]*(*v++); 262b5b17502SBarry Smith } 2639566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(4.0*a->nz)); 264b5b17502SBarry Smith } 265b5b17502SBarry Smith 266b5b17502SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) { 267b5b17502SBarry Smith /* 268b5b17502SBarry Smith backward sweep: 269b5b17502SBarry Smith b = b - x[i]*U^T(i,:), i=0,...,n-2 270b5b17502SBarry Smith for i=m-1,...,0: 271b5b17502SBarry Smith sum[i] = (b[i] - U(i,:)x)/d[i]; 272b5b17502SBarry Smith x[i] = (1-omega)x[i] + omega*sum[i]; 273b5b17502SBarry Smith */ 274b5b17502SBarry Smith /* if there was a forward sweep done above then I thing the next two for loops are not needed */ 2759566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t,b,m)); 276b5b17502SBarry Smith 277b5b17502SBarry Smith for (i=0; i<m-1; i++) { /* update rhs */ 278b5b17502SBarry Smith v = aa + ai[i] + 1; 279b5b17502SBarry Smith vj = aj + ai[i] + 1; 280b5b17502SBarry Smith nz = ai[i+1] - ai[i] - 1; 281b5b17502SBarry Smith while (nz--) t[*vj++] -= x[i]*(*v++); 282b5b17502SBarry Smith } 2839566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0*(a->nz - m))); 284b5b17502SBarry Smith for (i=m-1; i>=0; i--) { 285b5b17502SBarry Smith v = aa + ai[i] + 1; 286b5b17502SBarry Smith vj = aj + ai[i] + 1; 287b5b17502SBarry Smith nz = ai[i+1] - ai[i] - 1; 288b5b17502SBarry Smith sum = t[i]; 289b5b17502SBarry Smith while (nz--) sum -= x[*vj++]*(*v++); 290b5b17502SBarry Smith x[i] = (1-omega)*x[i] + omega*sum*aidiag[i]; 291b5b17502SBarry Smith } 2929566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0*(a->nz + m))); 293b5b17502SBarry Smith } 294b5b17502SBarry Smith } 295b5b17502SBarry Smith 2969566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xx,&x)); 2979566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(bb,&b)); 298b5b17502SBarry Smith PetscFunctionReturn(0); 299b5b17502SBarry Smith } 300