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 */ 5b5b17502SBarry Smith #undef __FUNCT__ 6018dd85eSSatish Balay #if defined(USESHORT) 7*b583852cSJed Brown #define __FUNCT__ "MatMult_SeqSBAIJ_1_Hermitian_ushort" 8018dd85eSSatish Balay PetscErrorCode MatMult_SeqSBAIJ_1_Hermitian_ushort(Mat A,Vec xx,Vec zz) 9018dd85eSSatish Balay #else 10*b583852cSJed Brown #define __FUNCT__ "MatMult_SeqSBAIJ_1_Hermitian" 11018dd85eSSatish Balay PetscErrorCode MatMult_SeqSBAIJ_1_Hermitian(Mat A,Vec xx,Vec zz) 12018dd85eSSatish Balay #endif 13eeffb40dSHong Zhang { 14eeffb40dSHong Zhang Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 15eeffb40dSHong Zhang const PetscScalar *x; 16eeffb40dSHong Zhang PetscScalar *z,x1,sum; 17eeffb40dSHong Zhang const MatScalar *v; 18eeffb40dSHong Zhang MatScalar vj; 19eeffb40dSHong Zhang PetscErrorCode ierr; 20eeffb40dSHong Zhang PetscInt mbs=a->mbs,i,j,nz; 21eeffb40dSHong Zhang const PetscInt *ai=a->i; 22eeffb40dSHong Zhang #if defined(USESHORT) 23eeffb40dSHong Zhang const unsigned short *ib=a->jshort; 24eeffb40dSHong Zhang unsigned short ibt; 25eeffb40dSHong Zhang #else 26eeffb40dSHong Zhang const PetscInt *ib=a->j; 27eeffb40dSHong Zhang PetscInt ibt; 28eeffb40dSHong Zhang #endif 29eeffb40dSHong Zhang 30eeffb40dSHong Zhang PetscFunctionBegin; 31eeffb40dSHong Zhang ierr = VecSet(zz,0.0);CHKERRQ(ierr); 323649974fSBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 33eeffb40dSHong Zhang ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 34eeffb40dSHong Zhang 35eeffb40dSHong Zhang v = a->a; 36eeffb40dSHong Zhang for (i=0; i<mbs; i++) { 37eeffb40dSHong Zhang nz = ai[i+1] - ai[i]; /* length of i_th row of A */ 38eeffb40dSHong Zhang x1 = x[i]; 39eeffb40dSHong Zhang sum = v[0]*x1; /* diagonal term */ 40eeffb40dSHong Zhang for (j=1; j<nz; j++) { 41a6f0575fSBarry Smith ibt = ib[j]; 42eeffb40dSHong Zhang vj = v[j]; 43eeffb40dSHong Zhang sum += vj * x[ibt]; /* (strict upper triangular part of A)*x */ 44a6f0575fSBarry Smith z[ibt] += PetscConj(v[j]) * x1; /* (strict lower triangular part of A)*x */ 45eeffb40dSHong Zhang } 46eeffb40dSHong Zhang z[i] += sum; 47eeffb40dSHong Zhang v += nz; 48eeffb40dSHong Zhang ib += nz; 49eeffb40dSHong Zhang } 50eeffb40dSHong Zhang 513649974fSBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 52eeffb40dSHong Zhang ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 53eeffb40dSHong Zhang ierr = PetscLogFlops(2.0*(2.0*a->nz - mbs) - mbs);CHKERRQ(ierr); 54eeffb40dSHong Zhang PetscFunctionReturn(0); 55eeffb40dSHong Zhang } 56eeffb40dSHong Zhang 57eeffb40dSHong Zhang #undef __FUNCT__ 58018dd85eSSatish Balay #if defined(USESHORT) 59*b583852cSJed Brown #define __FUNCT__ "MatMult_SeqSBAIJ_1_ushort" 60018dd85eSSatish Balay PetscErrorCode MatMult_SeqSBAIJ_1_ushort(Mat A,Vec xx,Vec zz) 61018dd85eSSatish Balay #else 62*b583852cSJed Brown #define __FUNCT__ "MatMult_SeqSBAIJ_1" 63018dd85eSSatish Balay PetscErrorCode MatMult_SeqSBAIJ_1(Mat A,Vec xx,Vec zz) 64018dd85eSSatish Balay #endif 65b5b17502SBarry Smith { 66b5b17502SBarry Smith Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 67b5b17502SBarry Smith const PetscScalar *x; 68b5b17502SBarry Smith PetscScalar *z,x1,sum; 69b5b17502SBarry Smith const MatScalar *v; 70b5b17502SBarry Smith MatScalar vj; 71b5b17502SBarry Smith PetscErrorCode ierr; 72b5b17502SBarry Smith PetscInt mbs=a->mbs,i,j,nz; 73b5b17502SBarry Smith const PetscInt *ai=a->i; 74b5b17502SBarry Smith #if defined(USESHORT) 75b5b17502SBarry Smith const unsigned short *ib=a->jshort; 76b5b17502SBarry Smith unsigned short ibt; 77b5b17502SBarry Smith #else 78b5b17502SBarry Smith const PetscInt *ib=a->j; 79b5b17502SBarry Smith PetscInt ibt; 80b5b17502SBarry Smith #endif 81b5b17502SBarry Smith 82b5b17502SBarry Smith PetscFunctionBegin; 83b5b17502SBarry Smith ierr = VecSet(zz,0.0);CHKERRQ(ierr); 843649974fSBarry Smith ierr = VecGetArrayRead(xx,&x);CHKERRQ(ierr); 85b5b17502SBarry Smith ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 86b5b17502SBarry Smith 87b5b17502SBarry Smith v = a->a; 88b5b17502SBarry Smith for (i=0; i<mbs; i++) { 89b5b17502SBarry Smith nz = ai[i+1] - ai[i]; /* length of i_th row of A */ 90b5b17502SBarry Smith x1 = x[i]; 91b5b17502SBarry Smith sum = v[0]*x1; /* diagonal term */ 92b5b17502SBarry Smith for (j=1; j<nz; j++) { 93b5b17502SBarry Smith ibt = ib[j]; 94b5b17502SBarry Smith vj = v[j]; 95b5b17502SBarry Smith z[ibt] += vj * x1; /* (strict lower triangular part of A)*x */ 96b5b17502SBarry Smith sum += vj * x[ibt]; /* (strict upper triangular part of A)*x */ 97b5b17502SBarry Smith } 98b5b17502SBarry Smith z[i] += sum; 99b5b17502SBarry Smith v += nz; 100b5b17502SBarry Smith ib += nz; 101b5b17502SBarry Smith } 102b5b17502SBarry Smith 1033649974fSBarry Smith ierr = VecRestoreArrayRead(xx,&x);CHKERRQ(ierr); 104b5b17502SBarry Smith ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 105b5b17502SBarry Smith ierr = PetscLogFlops(2.0*(2.0*a->nz - mbs) - mbs);CHKERRQ(ierr); 106b5b17502SBarry Smith PetscFunctionReturn(0); 107b5b17502SBarry Smith } 108b5b17502SBarry Smith 109b5b17502SBarry Smith #undef __FUNCT__ 110018dd85eSSatish Balay #if defined(USESHORT) 111*b583852cSJed Brown #define __FUNCT__ "MatSOR_SeqSBAIJ_ushort" 112018dd85eSSatish Balay PetscErrorCode MatSOR_SeqSBAIJ_ushort(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec xx) 113018dd85eSSatish Balay #else 114*b583852cSJed Brown #define __FUNCT__ "MatSOR_SeqSBAIJ" 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 PetscErrorCode ierr; 124b5b17502SBarry Smith PetscInt m=a->mbs,bs=A->rmap->bs,j; 125b5b17502SBarry Smith const PetscInt *ai=a->i; 126b5b17502SBarry Smith #if defined(USESHORT) 127b5b17502SBarry Smith const unsigned short *aj=a->jshort,*vj,*vj1; 128b5b17502SBarry Smith #else 129b5b17502SBarry Smith const PetscInt *aj=a->j,*vj,*vj1; 130b5b17502SBarry Smith #endif 131b5b17502SBarry Smith PetscInt nz,nz1,i; 132b5b17502SBarry Smith 133b5b17502SBarry Smith PetscFunctionBegin; 134e32f2f54SBarry Smith if (flag & SOR_EISENSTAT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support yet for Eisenstat"); 135b5b17502SBarry Smith 136b5b17502SBarry Smith its = its*lits; 137e32f2f54SBarry Smith if (its <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits); 138b5b17502SBarry Smith 139e32f2f54SBarry Smith if (bs > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"SSOR for block size > 1 is not yet implemented"); 140b5b17502SBarry Smith 141b5b17502SBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 1423649974fSBarry Smith ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr); 143b5b17502SBarry Smith 144b5b17502SBarry Smith if (!a->idiagvalid) { 145b5b17502SBarry Smith if (!a->idiag) { 146b5b17502SBarry Smith ierr = PetscMalloc(m*sizeof(PetscScalar),&a->idiag);CHKERRQ(ierr); 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) { 15341f059aeSBarry Smith ierr = PetscMalloc(m*sizeof(PetscScalar),&a->sor_work);CHKERRQ(ierr); 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; 168b5b17502SBarry Smith PetscSparseDensePlusDot(sum,b,v,vj,nz); 169b5b17502SBarry Smith x[i] = sum; 170b5b17502SBarry Smith } 171b5b17502SBarry Smith ierr = PetscLogFlops(a->nz);CHKERRQ(ierr); 172b5b17502SBarry Smith } 173b5b17502SBarry Smith 174b5b17502SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 175b5b17502SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 176b5b17502SBarry Smith ierr = PetscMemcpy(t,b,m*sizeof(PetscScalar));CHKERRQ(ierr); 177b5b17502SBarry Smith 178b5b17502SBarry Smith v = aa + 1; 179b5b17502SBarry Smith vj = aj + 1; 180b5b17502SBarry Smith for (i=0; i<m; i++){ 181b5b17502SBarry Smith nz = ai[i+1] - ai[i] - 1; 182b5b17502SBarry Smith tmp = - (x[i] = omega*t[i]*aidiag[i]); 183b5b17502SBarry Smith for (j=0; j<nz; j++) { 184b5b17502SBarry Smith t[vj[j]] += tmp*v[j]; 185b5b17502SBarry Smith } 186b5b17502SBarry Smith v += nz + 1; 187b5b17502SBarry Smith vj += nz + 1; 188b5b17502SBarry Smith } 189b5b17502SBarry Smith ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr); 190b5b17502SBarry Smith } 191b5b17502SBarry Smith 192b5b17502SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 193b5b17502SBarry Smith int nz2; 194b5b17502SBarry Smith if (!(flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP)){ 195b5b17502SBarry Smith #if defined(PETSC_USE_BACKWARD_LOOP) 196b5b17502SBarry Smith v = aa + ai[m] - 1; 197b5b17502SBarry Smith vj = aj + ai[m] - 1; 198b5b17502SBarry Smith for (i=m-1; i>=0; i--){ 199b5b17502SBarry Smith sum = b[i]; 200b5b17502SBarry Smith nz = ai[i+1] - ai[i] - 1; 201b5b17502SBarry Smith {PetscInt __i;for(__i=0;__i<nz;__i++) sum -= v[-__i] * x[vj[-__i]];} 202b5b17502SBarry Smith #else 203b5b17502SBarry Smith v = aa + ai[m-1] + 1; 204b5b17502SBarry Smith vj = aj + ai[m-1] + 1; 205b5b17502SBarry Smith nz = 0; 206b5b17502SBarry Smith for (i=m-1; i>=0; i--){ 207b5b17502SBarry Smith sum = b[i]; 208b5b17502SBarry Smith nz2 = ai[i] - ai[i-1] - 1; 20950d8bf02SJed Brown PETSC_Prefetch(v-nz2-1,0,PETSC_PREFETCH_HINT_NTA); 21050d8bf02SJed Brown PETSC_Prefetch(vj-nz2-1,0,PETSC_PREFETCH_HINT_NTA); 211b5b17502SBarry Smith PetscSparseDenseMinusDot(sum,x,v,vj,nz); 212b5b17502SBarry Smith nz = nz2; 213b5b17502SBarry Smith #endif 214b5b17502SBarry Smith x[i] = omega*sum*aidiag[i]; 215b5b17502SBarry Smith v -= nz + 1; 216b5b17502SBarry Smith vj -= nz + 1; 217b5b17502SBarry Smith } 218b5b17502SBarry Smith ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr); 219b5b17502SBarry Smith } else { 220b5b17502SBarry Smith v = aa + ai[m-1] + 1; 221b5b17502SBarry Smith vj = aj + ai[m-1] + 1; 222b5b17502SBarry Smith nz = 0; 223b5b17502SBarry Smith for (i=m-1; i>=0; i--){ 224b5b17502SBarry Smith sum = t[i]; 225b5b17502SBarry Smith nz2 = ai[i] - ai[i-1] - 1; 22650d8bf02SJed Brown PETSC_Prefetch(v-nz2-1,0,PETSC_PREFETCH_HINT_NTA); 22750d8bf02SJed Brown PETSC_Prefetch(vj-nz2-1,0,PETSC_PREFETCH_HINT_NTA); 228b5b17502SBarry Smith PetscSparseDenseMinusDot(sum,x,v,vj,nz); 229b5b17502SBarry Smith x[i] = (1-omega)*x[i] + omega*sum*aidiag[i]; 230b5b17502SBarry Smith nz = nz2; 231b5b17502SBarry Smith v -= nz + 1; 232b5b17502SBarry Smith vj -= nz + 1; 233b5b17502SBarry Smith } 234b5b17502SBarry Smith ierr = PetscLogFlops(2*a->nz);CHKERRQ(ierr); 235b5b17502SBarry Smith } 236b5b17502SBarry Smith } 237b5b17502SBarry Smith its--; 238b5b17502SBarry Smith } 239b5b17502SBarry Smith 240b5b17502SBarry Smith while (its--) { 241b5b17502SBarry Smith /* 242b5b17502SBarry Smith forward sweep: 243b5b17502SBarry Smith for i=0,...,m-1: 244b5b17502SBarry Smith sum[i] = (b[i] - U(i,:)x )/d[i]; 245b5b17502SBarry Smith x[i] = (1-omega)x[i] + omega*sum[i]; 246b5b17502SBarry Smith b = b - x[i]*U^T(i,:); 247b5b17502SBarry Smith 248b5b17502SBarry Smith */ 249b5b17502SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 250b5b17502SBarry Smith ierr = PetscMemcpy(t,b,m*sizeof(PetscScalar));CHKERRQ(ierr); 251b5b17502SBarry Smith 252b5b17502SBarry Smith for (i=0; i<m; i++){ 253b5b17502SBarry Smith v = aa + ai[i] + 1; v1=v; 254b5b17502SBarry Smith vj = aj + ai[i] + 1; vj1=vj; 255b5b17502SBarry Smith nz = ai[i+1] - ai[i] - 1; nz1=nz; 256b5b17502SBarry Smith sum = t[i]; 257b5b17502SBarry Smith ierr = PetscLogFlops(4.0*nz-2);CHKERRQ(ierr); 258b5b17502SBarry Smith while (nz1--) sum -= (*v1++)*x[*vj1++]; 259b5b17502SBarry Smith x[i] = (1-omega)*x[i] + omega*sum*aidiag[i]; 260b5b17502SBarry Smith while (nz--) t[*vj++] -= x[i]*(*v++); 261b5b17502SBarry Smith } 262b5b17502SBarry Smith } 263b5b17502SBarry Smith 264b5b17502SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 265b5b17502SBarry Smith /* 266b5b17502SBarry Smith backward sweep: 267b5b17502SBarry Smith b = b - x[i]*U^T(i,:), i=0,...,n-2 268b5b17502SBarry Smith for i=m-1,...,0: 269b5b17502SBarry Smith sum[i] = (b[i] - U(i,:)x )/d[i]; 270b5b17502SBarry Smith x[i] = (1-omega)x[i] + omega*sum[i]; 271b5b17502SBarry Smith */ 272b5b17502SBarry Smith /* if there was a forward sweep done above then I thing the next two for loops are not needed */ 273b5b17502SBarry Smith ierr = PetscMemcpy(t,b,m*sizeof(PetscScalar));CHKERRQ(ierr); 274b5b17502SBarry Smith 275b5b17502SBarry Smith for (i=0; i<m-1; i++){ /* update rhs */ 276b5b17502SBarry Smith v = aa + ai[i] + 1; 277b5b17502SBarry Smith vj = aj + ai[i] + 1; 278b5b17502SBarry Smith nz = ai[i+1] - ai[i] - 1; 279b5b17502SBarry Smith ierr = PetscLogFlops(2.0*nz-1);CHKERRQ(ierr); 280b5b17502SBarry Smith while (nz--) t[*vj++] -= x[i]*(*v++); 281b5b17502SBarry Smith } 282b5b17502SBarry Smith for (i=m-1; i>=0; i--){ 283b5b17502SBarry Smith v = aa + ai[i] + 1; 284b5b17502SBarry Smith vj = aj + ai[i] + 1; 285b5b17502SBarry Smith nz = ai[i+1] - ai[i] - 1; 286b5b17502SBarry Smith ierr = PetscLogFlops(2.0*nz-1);CHKERRQ(ierr); 287b5b17502SBarry Smith sum = t[i]; 288b5b17502SBarry Smith while (nz--) sum -= x[*vj++]*(*v++); 289b5b17502SBarry Smith x[i] = (1-omega)*x[i] + omega*sum*aidiag[i]; 290b5b17502SBarry Smith } 291b5b17502SBarry Smith } 292b5b17502SBarry Smith } 293b5b17502SBarry Smith 294b5b17502SBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2953649974fSBarry Smith ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr); 296b5b17502SBarry Smith PetscFunctionReturn(0); 297b5b17502SBarry Smith } 298