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) 109371c9d4SSatish Balay #define PetscSparseDensePlusDot_no_function(sum, r, xv, xi, nnz) \ 119371c9d4SSatish Balay { \ 1254e8760dSRichard Tran Mills if (nnz > 0) { \ 13a9ce9505SJunchao Zhang PetscInt nnz2 = nnz, rem = nnz & 0x3; \ 14a9ce9505SJunchao Zhang switch (rem) { \ 1554e8760dSRichard Tran Mills case 3: sum += *xv++ * r[*xi++]; \ 1654e8760dSRichard Tran Mills case 2: sum += *xv++ * r[*xi++]; \ 179371c9d4SSatish Balay case 1: sum += *xv++ * r[*xi++]; nnz2 -= rem; \ 18a9ce9505SJunchao Zhang } \ 199371c9d4SSatish Balay while (nnz2 > 0) { \ 209371c9d4SSatish Balay sum += xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \ 219371c9d4SSatish Balay xv += 4; \ 229371c9d4SSatish Balay xi += 4; \ 239371c9d4SSatish Balay nnz2 -= 4; \ 249371c9d4SSatish Balay } \ 259371c9d4SSatish Balay xv -= nnz; \ 269371c9d4SSatish Balay xi -= nnz; \ 27a9ce9505SJunchao Zhang } \ 28a9ce9505SJunchao Zhang } 2954e8760dSRichard Tran Mills 3054e8760dSRichard Tran Mills #elif defined(PETSC_KERNEL_USE_UNROLL_2) 319371c9d4SSatish Balay #define PetscSparseDensePlusDot_no_function(sum, r, xv, xi, nnz) \ 329371c9d4SSatish Balay { \ 3354e8760dSRichard Tran Mills PetscInt __i, __i1, __i2; \ 349371c9d4SSatish Balay for (__i = 0; __i < nnz - 1; __i += 2) { \ 359371c9d4SSatish Balay __i1 = xi[__i]; \ 369371c9d4SSatish Balay __i2 = xi[__i + 1]; \ 379371c9d4SSatish Balay sum += (xv[__i] * r[__i1] + xv[__i + 1] * r[__i2]); \ 389371c9d4SSatish Balay } \ 399371c9d4SSatish Balay if (nnz & 0x1) sum += xv[__i] * r[xi[__i]]; \ 409371c9d4SSatish Balay } 4154e8760dSRichard Tran Mills 4254e8760dSRichard Tran Mills #else 439371c9d4SSatish Balay #define PetscSparseDensePlusDot_no_function(sum, r, xv, xi, nnz) \ 449371c9d4SSatish Balay { \ 4554e8760dSRichard Tran Mills PetscInt __i; \ 469371c9d4SSatish Balay for (__i = 0; __i < nnz; __i++) sum += xv[__i] * r[xi[__i]]; \ 479371c9d4SSatish Balay } 4854e8760dSRichard Tran Mills #endif 4954e8760dSRichard Tran Mills 50018dd85eSSatish Balay #if defined(USESHORT) 51018dd85eSSatish Balay PetscErrorCode MatMult_SeqSBAIJ_1_ushort(Mat A, Vec xx, Vec zz) 52018dd85eSSatish Balay #else 53018dd85eSSatish Balay PetscErrorCode MatMult_SeqSBAIJ_1(Mat A, Vec xx, Vec zz) 54018dd85eSSatish Balay #endif 55b5b17502SBarry Smith { 56b5b17502SBarry Smith Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data; 57b5b17502SBarry Smith const PetscScalar *x; 58b5b17502SBarry Smith PetscScalar *z, x1, sum; 59b5b17502SBarry Smith const MatScalar *v; 60b5b17502SBarry Smith MatScalar vj; 61b5b17502SBarry Smith PetscInt mbs = a->mbs, i, j, nz; 62b5b17502SBarry Smith const PetscInt *ai = a->i; 63b5b17502SBarry Smith #if defined(USESHORT) 64b5b17502SBarry Smith const unsigned short *ib = a->jshort; 65b5b17502SBarry Smith unsigned short ibt; 66b5b17502SBarry Smith #else 67b5b17502SBarry Smith const PetscInt *ib = a->j; 68b5b17502SBarry Smith PetscInt ibt; 69b5b17502SBarry Smith #endif 708dca527aSShri Abhyankar PetscInt nonzerorow = 0, jmin; 71eb1ec7c1SStefano Zampini #if defined(PETSC_USE_COMPLEX) 72b94d7dedSBarry Smith const int aconj = A->hermitian == PETSC_BOOL3_TRUE; 73eb1ec7c1SStefano Zampini #else 74eb1ec7c1SStefano Zampini const int aconj = 0; 75eb1ec7c1SStefano Zampini #endif 76b5b17502SBarry Smith 77b5b17502SBarry Smith PetscFunctionBegin; 789566063dSJacob Faibussowitsch PetscCall(VecSet(zz, 0.0)); 799566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(xx, &x)); 809566063dSJacob Faibussowitsch PetscCall(VecGetArray(zz, &z)); 81b5b17502SBarry Smith 82b5b17502SBarry Smith v = a->a; 83b5b17502SBarry Smith for (i = 0; i < mbs; i++) { 84b5b17502SBarry Smith nz = ai[i + 1] - ai[i]; /* length of i_th row of A */ 8503d09022SShri Abhyankar if (!nz) continue; /* Move to the next row if the current row is empty */ 8603d09022SShri Abhyankar nonzerorow++; 878dca527aSShri Abhyankar sum = 0.0; 888dca527aSShri Abhyankar jmin = 0; 89b5b17502SBarry Smith x1 = x[i]; 908dca527aSShri Abhyankar if (ib[0] == i) { 91b5b17502SBarry Smith sum = v[0] * x1; /* diagonal term */ 928dca527aSShri Abhyankar jmin++; 938dca527aSShri Abhyankar } 94444d8c10SJed Brown PetscPrefetchBlock(ib + nz, nz, 0, PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */ 95444d8c10SJed Brown PetscPrefetchBlock(v + nz, nz, 0, PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */ 96eb1ec7c1SStefano Zampini if (aconj) { 97eb1ec7c1SStefano Zampini for (j = jmin; j < nz; j++) { 98eb1ec7c1SStefano Zampini ibt = ib[j]; 99eb1ec7c1SStefano Zampini vj = v[j]; 100eb1ec7c1SStefano Zampini z[ibt] += PetscConj(vj) * x1; /* (strict lower triangular part of A)*x */ 101eb1ec7c1SStefano Zampini sum += vj * x[ibt]; /* (strict upper triangular part of A)*x */ 102eb1ec7c1SStefano Zampini } 103eb1ec7c1SStefano Zampini } else { 1048dca527aSShri Abhyankar for (j = jmin; j < nz; j++) { 105b5b17502SBarry Smith ibt = ib[j]; 106b5b17502SBarry Smith vj = v[j]; 107b5b17502SBarry Smith z[ibt] += vj * x1; /* (strict lower triangular part of A)*x */ 108b5b17502SBarry Smith sum += vj * x[ibt]; /* (strict upper triangular part of A)*x */ 109b5b17502SBarry Smith } 110eb1ec7c1SStefano Zampini } 111b5b17502SBarry Smith z[i] += sum; 112b5b17502SBarry Smith v += nz; 113b5b17502SBarry Smith ib += nz; 114b5b17502SBarry Smith } 115b5b17502SBarry Smith 1169566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(xx, &x)); 1179566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(zz, &z)); 1189566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * (2.0 * a->nz - nonzerorow) - nonzerorow)); 119b5b17502SBarry Smith PetscFunctionReturn(0); 120b5b17502SBarry Smith } 121b5b17502SBarry Smith 122018dd85eSSatish Balay #if defined(USESHORT) 123018dd85eSSatish Balay PetscErrorCode MatSOR_SeqSBAIJ_ushort(Mat A, Vec bb, PetscReal omega, MatSORType flag, PetscReal fshift, PetscInt its, PetscInt lits, Vec xx) 124018dd85eSSatish Balay #else 125018dd85eSSatish Balay PetscErrorCode MatSOR_SeqSBAIJ(Mat A, Vec bb, PetscReal omega, MatSORType flag, PetscReal fshift, PetscInt its, PetscInt lits, Vec xx) 126018dd85eSSatish Balay #endif 127b5b17502SBarry Smith { 128b5b17502SBarry Smith Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data; 129b5b17502SBarry Smith const MatScalar *aa = a->a, *v, *v1, *aidiag; 130b5b17502SBarry Smith PetscScalar *x, *t, sum; 131b5b17502SBarry Smith const PetscScalar *b; 132b5b17502SBarry Smith MatScalar tmp; 133b5b17502SBarry Smith PetscInt m = a->mbs, bs = A->rmap->bs, j; 134b5b17502SBarry Smith const PetscInt *ai = a->i; 135b5b17502SBarry Smith #if defined(USESHORT) 136b5b17502SBarry Smith const unsigned short *aj = a->jshort, *vj, *vj1; 137b5b17502SBarry Smith #else 138b5b17502SBarry Smith const PetscInt *aj = a->j, *vj, *vj1; 139b5b17502SBarry Smith #endif 140b5b17502SBarry Smith PetscInt nz, nz1, i; 141b5b17502SBarry Smith 142b5b17502SBarry Smith PetscFunctionBegin; 143422a814eSBarry Smith if (fshift == -1.0) fshift = 0.0; /* negative fshift indicates do not error on zero diagonal; this code never errors on zero diagonal */ 144aed4548fSBarry Smith PetscCheck(!(flag & SOR_EISENSTAT), PETSC_COMM_SELF, PETSC_ERR_SUP, "No support yet for Eisenstat"); 145b5b17502SBarry Smith 146b5b17502SBarry Smith its = its * lits; 14708401ef6SPierre 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); 148b5b17502SBarry Smith 14908401ef6SPierre Jolivet PetscCheck(bs <= 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "SSOR for block size > 1 is not yet implemented"); 150b5b17502SBarry Smith 1519566063dSJacob Faibussowitsch PetscCall(VecGetArray(xx, &x)); 1529566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(bb, &b)); 153b5b17502SBarry Smith 154b5b17502SBarry Smith if (!a->idiagvalid) { 155*48a46eb9SPierre Jolivet if (!a->idiag) PetscCall(PetscMalloc1(m, &a->idiag)); 156b5b17502SBarry Smith for (i = 0; i < a->mbs; i++) a->idiag[i] = 1.0 / a->a[a->i[i]]; 157b5b17502SBarry Smith a->idiagvalid = PETSC_TRUE; 158b5b17502SBarry Smith } 159b5b17502SBarry Smith 160*48a46eb9SPierre Jolivet if (!a->sor_work) PetscCall(PetscMalloc1(m, &a->sor_work)); 16141f059aeSBarry Smith t = a->sor_work; 162b5b17502SBarry Smith 163b5b17502SBarry Smith aidiag = a->idiag; 164b5b17502SBarry Smith 165b5b17502SBarry Smith if (flag == SOR_APPLY_UPPER) { 166b5b17502SBarry Smith /* apply (U + D/omega) to the vector */ 167b5b17502SBarry Smith PetscScalar d; 168b5b17502SBarry Smith for (i = 0; i < m; i++) { 169b5b17502SBarry Smith d = fshift + aa[ai[i]]; 170b5b17502SBarry Smith nz = ai[i + 1] - ai[i] - 1; 171b5b17502SBarry Smith vj = aj + ai[i] + 1; 172b5b17502SBarry Smith v = aa + ai[i] + 1; 173b5b17502SBarry Smith sum = b[i] * d / omega; 17454e8760dSRichard Tran Mills #ifdef USESHORT 17554e8760dSRichard Tran Mills PetscSparseDensePlusDot_no_function(sum, b, v, vj, nz); 17654e8760dSRichard Tran Mills #else 177b5b17502SBarry Smith PetscSparseDensePlusDot(sum, b, v, vj, nz); 17854e8760dSRichard Tran Mills #endif 179b5b17502SBarry Smith x[i] = sum; 180b5b17502SBarry Smith } 1819566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(a->nz)); 182b5b17502SBarry Smith } 183b5b17502SBarry Smith 184b5b17502SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 185b5b17502SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) { 1869566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t, b, m)); 187b5b17502SBarry Smith 188b5b17502SBarry Smith v = aa + 1; 189b5b17502SBarry Smith vj = aj + 1; 190b5b17502SBarry Smith for (i = 0; i < m; i++) { 191b5b17502SBarry Smith nz = ai[i + 1] - ai[i] - 1; 192b5b17502SBarry Smith tmp = -(x[i] = omega * t[i] * aidiag[i]); 19326fbe8dcSKarl Rupp for (j = 0; j < nz; j++) t[vj[j]] += tmp * v[j]; 194b5b17502SBarry Smith v += nz + 1; 195b5b17502SBarry Smith vj += nz + 1; 196b5b17502SBarry Smith } 1979566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * a->nz)); 198b5b17502SBarry Smith } 199b5b17502SBarry Smith 200b5b17502SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) { 201b5b17502SBarry Smith int nz2; 202b5b17502SBarry Smith if (!(flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP)) { 203b5b17502SBarry Smith #if defined(PETSC_USE_BACKWARD_LOOP) 204b5b17502SBarry Smith v = aa + ai[m] - 1; 205b5b17502SBarry Smith vj = aj + ai[m] - 1; 206b5b17502SBarry Smith for (i = m - 1; i >= 0; i--) { 207b5b17502SBarry Smith sum = b[i]; 208b5b17502SBarry Smith nz = ai[i + 1] - ai[i] - 1; 2099371c9d4SSatish Balay { 2109371c9d4SSatish Balay PetscInt __i; 2119371c9d4SSatish Balay for (__i = 0; __i < nz; __i++) sum -= v[-__i] * x[vj[-__i]]; 2129371c9d4SSatish Balay } 213b5b17502SBarry Smith #else 214b5b17502SBarry Smith v = aa + ai[m - 1] + 1; 215b5b17502SBarry Smith vj = aj + ai[m - 1] + 1; 216b5b17502SBarry Smith nz = 0; 217b5b17502SBarry Smith for (i = m - 1; i >= 0; i--) { 218b5b17502SBarry Smith sum = b[i]; 2196c11f900SJed Brown nz2 = ai[i] - ai[PetscMax(i - 1, 0)] - 1; /* avoid referencing ai[-1], nonsense nz2 is okay on last iteration */ 22050d8bf02SJed Brown PETSC_Prefetch(v - nz2 - 1, 0, PETSC_PREFETCH_HINT_NTA); 22150d8bf02SJed Brown PETSC_Prefetch(vj - nz2 - 1, 0, PETSC_PREFETCH_HINT_NTA); 222b5b17502SBarry Smith PetscSparseDenseMinusDot(sum, x, v, vj, nz); 223b5b17502SBarry Smith nz = nz2; 224b5b17502SBarry Smith #endif 225b5b17502SBarry Smith x[i] = omega * sum * aidiag[i]; 226b5b17502SBarry Smith v -= nz + 1; 227b5b17502SBarry Smith vj -= nz + 1; 228b5b17502SBarry Smith } 2299566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * a->nz)); 230b5b17502SBarry Smith } else { 231b5b17502SBarry Smith v = aa + ai[m - 1] + 1; 232b5b17502SBarry Smith vj = aj + ai[m - 1] + 1; 233b5b17502SBarry Smith nz = 0; 234b5b17502SBarry Smith for (i = m - 1; i >= 0; i--) { 235b5b17502SBarry Smith sum = t[i]; 2366c11f900SJed Brown nz2 = ai[i] - ai[PetscMax(i - 1, 0)] - 1; /* avoid referencing ai[-1], nonsense nz2 is okay on last iteration */ 23750d8bf02SJed Brown PETSC_Prefetch(v - nz2 - 1, 0, PETSC_PREFETCH_HINT_NTA); 23850d8bf02SJed Brown PETSC_Prefetch(vj - nz2 - 1, 0, PETSC_PREFETCH_HINT_NTA); 239b5b17502SBarry Smith PetscSparseDenseMinusDot(sum, x, v, vj, nz); 240b5b17502SBarry Smith x[i] = (1 - omega) * x[i] + omega * sum * aidiag[i]; 241b5b17502SBarry Smith nz = nz2; 242b5b17502SBarry Smith v -= nz + 1; 243b5b17502SBarry Smith vj -= nz + 1; 244b5b17502SBarry Smith } 2459566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * a->nz)); 246b5b17502SBarry Smith } 247b5b17502SBarry Smith } 248b5b17502SBarry Smith its--; 249b5b17502SBarry Smith } 250b5b17502SBarry Smith 251b5b17502SBarry Smith while (its--) { 252b5b17502SBarry Smith /* 253b5b17502SBarry Smith forward sweep: 254b5b17502SBarry Smith for i=0,...,m-1: 255b5b17502SBarry Smith sum[i] = (b[i] - U(i,:)x)/d[i]; 256b5b17502SBarry Smith x[i] = (1-omega)x[i] + omega*sum[i]; 257b5b17502SBarry Smith b = b - x[i]*U^T(i,:); 258b5b17502SBarry Smith 259b5b17502SBarry Smith */ 260b5b17502SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) { 2619566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t, b, m)); 262b5b17502SBarry Smith 263b5b17502SBarry Smith for (i = 0; i < m; i++) { 2649371c9d4SSatish Balay v = aa + ai[i] + 1; 2659371c9d4SSatish Balay v1 = v; 2669371c9d4SSatish Balay vj = aj + ai[i] + 1; 2679371c9d4SSatish Balay vj1 = vj; 2689371c9d4SSatish Balay nz = ai[i + 1] - ai[i] - 1; 2699371c9d4SSatish Balay nz1 = nz; 270b5b17502SBarry Smith sum = t[i]; 271b5b17502SBarry Smith while (nz1--) sum -= (*v1++) * x[*vj1++]; 272b5b17502SBarry Smith x[i] = (1 - omega) * x[i] + omega * sum * aidiag[i]; 273b5b17502SBarry Smith while (nz--) t[*vj++] -= x[i] * (*v++); 274b5b17502SBarry Smith } 2759566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(4.0 * a->nz)); 276b5b17502SBarry Smith } 277b5b17502SBarry Smith 278b5b17502SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) { 279b5b17502SBarry Smith /* 280b5b17502SBarry Smith backward sweep: 281b5b17502SBarry Smith b = b - x[i]*U^T(i,:), i=0,...,n-2 282b5b17502SBarry Smith for i=m-1,...,0: 283b5b17502SBarry Smith sum[i] = (b[i] - U(i,:)x)/d[i]; 284b5b17502SBarry Smith x[i] = (1-omega)x[i] + omega*sum[i]; 285b5b17502SBarry Smith */ 286b5b17502SBarry Smith /* if there was a forward sweep done above then I thing the next two for loops are not needed */ 2879566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t, b, m)); 288b5b17502SBarry Smith 289b5b17502SBarry Smith for (i = 0; i < m - 1; i++) { /* update rhs */ 290b5b17502SBarry Smith v = aa + ai[i] + 1; 291b5b17502SBarry Smith vj = aj + ai[i] + 1; 292b5b17502SBarry Smith nz = ai[i + 1] - ai[i] - 1; 293b5b17502SBarry Smith while (nz--) t[*vj++] -= x[i] * (*v++); 294b5b17502SBarry Smith } 2959566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * (a->nz - m))); 296b5b17502SBarry Smith for (i = m - 1; i >= 0; i--) { 297b5b17502SBarry Smith v = aa + ai[i] + 1; 298b5b17502SBarry Smith vj = aj + ai[i] + 1; 299b5b17502SBarry Smith nz = ai[i + 1] - ai[i] - 1; 300b5b17502SBarry Smith sum = t[i]; 301b5b17502SBarry Smith while (nz--) sum -= x[*vj++] * (*v++); 302b5b17502SBarry Smith x[i] = (1 - omega) * x[i] + omega * sum * aidiag[i]; 303b5b17502SBarry Smith } 3049566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * (a->nz + m))); 305b5b17502SBarry Smith } 306b5b17502SBarry Smith } 307b5b17502SBarry Smith 3089566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xx, &x)); 3099566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(bb, &b)); 310b5b17502SBarry Smith PetscFunctionReturn(0); 311b5b17502SBarry Smith } 312