1b5b17502SBarry Smith /* 2b5b17502SBarry Smith This is included by sbaij.c to generate unsigned short and regular versions of these two functions 3b5b17502SBarry Smith */ 454e8760dSRichard Tran Mills 554e8760dSRichard Tran Mills /* We cut-and-past below from aij.h to make a "no_function" version of PetscSparseDensePlusDot(). 654e8760dSRichard Tran Mills * This is necessary because the USESHORT case cannot use the inlined functions that may be employed. */ 754e8760dSRichard Tran Mills 854e8760dSRichard Tran Mills #if defined(PETSC_KERNEL_USE_UNROLL_4) 99371c9d4SSatish Balay #define PetscSparseDensePlusDot_no_function(sum, r, xv, xi, nnz) \ 10a8f51744SPierre Jolivet do { \ 1154e8760dSRichard Tran Mills if (nnz > 0) { \ 12a9ce9505SJunchao Zhang PetscInt nnz2 = nnz, rem = nnz & 0x3; \ 13a9ce9505SJunchao Zhang switch (rem) { \ 14d71ae5a4SJacob Faibussowitsch case 3: \ 15d71ae5a4SJacob Faibussowitsch sum += *xv++ * r[*xi++]; \ 16d71ae5a4SJacob Faibussowitsch case 2: \ 17d71ae5a4SJacob Faibussowitsch sum += *xv++ * r[*xi++]; \ 18d71ae5a4SJacob Faibussowitsch case 1: \ 19d71ae5a4SJacob Faibussowitsch sum += *xv++ * r[*xi++]; \ 20d71ae5a4SJacob Faibussowitsch nnz2 -= rem; \ 21a9ce9505SJunchao Zhang } \ 229371c9d4SSatish Balay while (nnz2 > 0) { \ 239371c9d4SSatish Balay sum += xv[0] * r[xi[0]] + xv[1] * r[xi[1]] + xv[2] * r[xi[2]] + xv[3] * r[xi[3]]; \ 249371c9d4SSatish Balay xv += 4; \ 259371c9d4SSatish Balay xi += 4; \ 269371c9d4SSatish Balay nnz2 -= 4; \ 279371c9d4SSatish Balay } \ 289371c9d4SSatish Balay xv -= nnz; \ 299371c9d4SSatish Balay xi -= nnz; \ 30a9ce9505SJunchao Zhang } \ 31a8f51744SPierre Jolivet } while (0) 3254e8760dSRichard Tran Mills 3354e8760dSRichard Tran Mills #elif defined(PETSC_KERNEL_USE_UNROLL_2) 349371c9d4SSatish Balay #define PetscSparseDensePlusDot_no_function(sum, r, xv, xi, nnz) \ 35a8f51744SPierre Jolivet do { \ 3654e8760dSRichard Tran Mills PetscInt __i, __i1, __i2; \ 379371c9d4SSatish Balay for (__i = 0; __i < nnz - 1; __i += 2) { \ 389371c9d4SSatish Balay __i1 = xi[__i]; \ 399371c9d4SSatish Balay __i2 = xi[__i + 1]; \ 409371c9d4SSatish Balay sum += (xv[__i] * r[__i1] + xv[__i + 1] * r[__i2]); \ 419371c9d4SSatish Balay } \ 429371c9d4SSatish Balay if (nnz & 0x1) sum += xv[__i] * r[xi[__i]]; \ 43a8f51744SPierre Jolivet } while (0) 4454e8760dSRichard Tran Mills 4554e8760dSRichard Tran Mills #else 469371c9d4SSatish Balay #define PetscSparseDensePlusDot_no_function(sum, r, xv, xi, nnz) \ 47a8f51744SPierre Jolivet do { \ 4854e8760dSRichard Tran Mills PetscInt __i; \ 499371c9d4SSatish Balay for (__i = 0; __i < nnz; __i++) sum += xv[__i] * r[xi[__i]]; \ 50a8f51744SPierre Jolivet } while (0) 5154e8760dSRichard Tran Mills #endif 5254e8760dSRichard Tran Mills 53018dd85eSSatish Balay #if defined(USESHORT) 54018dd85eSSatish Balay PetscErrorCode MatMult_SeqSBAIJ_1_ushort(Mat A, Vec xx, Vec zz) 55018dd85eSSatish Balay #else 56018dd85eSSatish Balay PetscErrorCode MatMult_SeqSBAIJ_1(Mat A, Vec xx, Vec zz) 57018dd85eSSatish Balay #endif 58b5b17502SBarry Smith { 59b5b17502SBarry Smith Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data; 60b5b17502SBarry Smith const PetscScalar *x; 61b5b17502SBarry Smith PetscScalar *z, x1, sum; 62b5b17502SBarry Smith const MatScalar *v; 63b5b17502SBarry Smith MatScalar vj; 64b5b17502SBarry Smith PetscInt mbs = a->mbs, i, j, nz; 65b5b17502SBarry Smith const PetscInt *ai = a->i; 66b5b17502SBarry Smith #if defined(USESHORT) 67b5b17502SBarry Smith const unsigned short *ib = a->jshort; 68b5b17502SBarry Smith unsigned short ibt; 69b5b17502SBarry Smith #else 70b5b17502SBarry Smith const PetscInt *ib = a->j; 71b5b17502SBarry Smith PetscInt ibt; 72b5b17502SBarry Smith #endif 738dca527aSShri Abhyankar PetscInt nonzerorow = 0, jmin; 74eb1ec7c1SStefano Zampini #if defined(PETSC_USE_COMPLEX) 75b94d7dedSBarry Smith const int aconj = A->hermitian == PETSC_BOOL3_TRUE; 76eb1ec7c1SStefano Zampini #else 77eb1ec7c1SStefano Zampini const int aconj = 0; 78eb1ec7c1SStefano Zampini #endif 79b5b17502SBarry Smith 80b5b17502SBarry Smith PetscFunctionBegin; 819566063dSJacob Faibussowitsch PetscCall(VecSet(zz, 0.0)); 829566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(xx, &x)); 839566063dSJacob Faibussowitsch PetscCall(VecGetArray(zz, &z)); 84b5b17502SBarry Smith 85b5b17502SBarry Smith v = a->a; 86b5b17502SBarry Smith for (i = 0; i < mbs; i++) { 87b5b17502SBarry Smith nz = ai[i + 1] - ai[i]; /* length of i_th row of A */ 8803d09022SShri Abhyankar if (!nz) continue; /* Move to the next row if the current row is empty */ 8903d09022SShri Abhyankar nonzerorow++; 908dca527aSShri Abhyankar sum = 0.0; 918dca527aSShri Abhyankar jmin = 0; 92b5b17502SBarry Smith x1 = x[i]; 938dca527aSShri Abhyankar if (ib[0] == i) { 94b5b17502SBarry Smith sum = v[0] * x1; /* diagonal term */ 958dca527aSShri Abhyankar jmin++; 968dca527aSShri Abhyankar } 97444d8c10SJed Brown PetscPrefetchBlock(ib + nz, nz, 0, PETSC_PREFETCH_HINT_NTA); /* Indices for the next row (assumes same size as this one) */ 98444d8c10SJed Brown PetscPrefetchBlock(v + nz, nz, 0, PETSC_PREFETCH_HINT_NTA); /* Entries for the next row */ 99eb1ec7c1SStefano Zampini if (aconj) { 100eb1ec7c1SStefano Zampini for (j = jmin; j < nz; j++) { 101eb1ec7c1SStefano Zampini ibt = ib[j]; 102eb1ec7c1SStefano Zampini vj = v[j]; 103eb1ec7c1SStefano Zampini z[ibt] += PetscConj(vj) * x1; /* (strict lower triangular part of A)*x */ 104eb1ec7c1SStefano Zampini sum += vj * x[ibt]; /* (strict upper triangular part of A)*x */ 105eb1ec7c1SStefano Zampini } 106eb1ec7c1SStefano Zampini } else { 1078dca527aSShri Abhyankar for (j = jmin; j < nz; j++) { 108b5b17502SBarry Smith ibt = ib[j]; 109b5b17502SBarry Smith vj = v[j]; 110b5b17502SBarry Smith z[ibt] += vj * x1; /* (strict lower triangular part of A)*x */ 111b5b17502SBarry Smith sum += vj * x[ibt]; /* (strict upper triangular part of A)*x */ 112b5b17502SBarry Smith } 113eb1ec7c1SStefano Zampini } 114b5b17502SBarry Smith z[i] += sum; 115b5b17502SBarry Smith v += nz; 116b5b17502SBarry Smith ib += nz; 117b5b17502SBarry Smith } 118b5b17502SBarry Smith 1199566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(xx, &x)); 1209566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(zz, &z)); 1219566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * (2.0 * a->nz - nonzerorow) - nonzerorow)); 1223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 123b5b17502SBarry Smith } 124b5b17502SBarry Smith 125018dd85eSSatish Balay #if defined(USESHORT) 126018dd85eSSatish Balay PetscErrorCode MatSOR_SeqSBAIJ_ushort(Mat A, Vec bb, PetscReal omega, MatSORType flag, PetscReal fshift, PetscInt its, PetscInt lits, Vec xx) 127018dd85eSSatish Balay #else 128018dd85eSSatish Balay PetscErrorCode MatSOR_SeqSBAIJ(Mat A, Vec bb, PetscReal omega, MatSORType flag, PetscReal fshift, PetscInt its, PetscInt lits, Vec xx) 129018dd85eSSatish Balay #endif 130b5b17502SBarry Smith { 131b5b17502SBarry Smith Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data; 132b5b17502SBarry Smith const MatScalar *aa = a->a, *v, *v1, *aidiag; 133b5b17502SBarry Smith PetscScalar *x, *t, sum; 134b5b17502SBarry Smith const PetscScalar *b; 135b5b17502SBarry Smith MatScalar tmp; 136b5b17502SBarry Smith PetscInt m = a->mbs, bs = A->rmap->bs, j; 137b5b17502SBarry Smith const PetscInt *ai = a->i; 138b5b17502SBarry Smith #if defined(USESHORT) 139b5b17502SBarry Smith const unsigned short *aj = a->jshort, *vj, *vj1; 140b5b17502SBarry Smith #else 141b5b17502SBarry Smith const PetscInt *aj = a->j, *vj, *vj1; 142b5b17502SBarry Smith #endif 143b5b17502SBarry Smith PetscInt nz, nz1, i; 144b5b17502SBarry Smith 145b5b17502SBarry Smith PetscFunctionBegin; 146422a814eSBarry Smith if (fshift == -1.0) fshift = 0.0; /* negative fshift indicates do not error on zero diagonal; this code never errors on zero diagonal */ 147aed4548fSBarry Smith PetscCheck(!(flag & SOR_EISENSTAT), PETSC_COMM_SELF, PETSC_ERR_SUP, "No support yet for Eisenstat"); 148b5b17502SBarry Smith 149b5b17502SBarry Smith its = its * lits; 15008401ef6SPierre 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); 151b5b17502SBarry Smith 15208401ef6SPierre Jolivet PetscCheck(bs <= 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "SSOR for block size > 1 is not yet implemented"); 153b5b17502SBarry Smith 1549566063dSJacob Faibussowitsch PetscCall(VecGetArray(xx, &x)); 1559566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(bb, &b)); 156b5b17502SBarry Smith 157b5b17502SBarry Smith if (!a->idiagvalid) { 15848a46eb9SPierre Jolivet if (!a->idiag) PetscCall(PetscMalloc1(m, &a->idiag)); 159b5b17502SBarry Smith for (i = 0; i < a->mbs; i++) a->idiag[i] = 1.0 / a->a[a->i[i]]; 160b5b17502SBarry Smith a->idiagvalid = PETSC_TRUE; 161b5b17502SBarry Smith } 162b5b17502SBarry Smith 16348a46eb9SPierre Jolivet if (!a->sor_work) PetscCall(PetscMalloc1(m, &a->sor_work)); 16441f059aeSBarry Smith t = a->sor_work; 165b5b17502SBarry Smith 166b5b17502SBarry Smith aidiag = a->idiag; 167b5b17502SBarry Smith 168b5b17502SBarry Smith if (flag == SOR_APPLY_UPPER) { 169b5b17502SBarry Smith /* apply (U + D/omega) to the vector */ 170b5b17502SBarry Smith PetscScalar d; 171b5b17502SBarry Smith for (i = 0; i < m; i++) { 172b5b17502SBarry Smith d = fshift + aa[ai[i]]; 173b5b17502SBarry Smith nz = ai[i + 1] - ai[i] - 1; 174b5b17502SBarry Smith vj = aj + ai[i] + 1; 175b5b17502SBarry Smith v = aa + ai[i] + 1; 176b5b17502SBarry Smith sum = b[i] * d / omega; 17754e8760dSRichard Tran Mills #ifdef USESHORT 17854e8760dSRichard Tran Mills PetscSparseDensePlusDot_no_function(sum, b, v, vj, nz); 17954e8760dSRichard Tran Mills #else 180b5b17502SBarry Smith PetscSparseDensePlusDot(sum, b, v, vj, nz); 18154e8760dSRichard Tran Mills #endif 182b5b17502SBarry Smith x[i] = sum; 183b5b17502SBarry Smith } 1849566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(a->nz)); 185b5b17502SBarry Smith } 186b5b17502SBarry Smith 187b5b17502SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 188b5b17502SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) { 1899566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t, b, m)); 190b5b17502SBarry Smith 191b5b17502SBarry Smith v = aa + 1; 192b5b17502SBarry Smith vj = aj + 1; 193b5b17502SBarry Smith for (i = 0; i < m; i++) { 194b5b17502SBarry Smith nz = ai[i + 1] - ai[i] - 1; 195b5b17502SBarry Smith tmp = -(x[i] = omega * t[i] * aidiag[i]); 19626fbe8dcSKarl Rupp for (j = 0; j < nz; j++) t[vj[j]] += tmp * v[j]; 197b5b17502SBarry Smith v += nz + 1; 198b5b17502SBarry Smith vj += nz + 1; 199b5b17502SBarry Smith } 2009566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * a->nz)); 201b5b17502SBarry Smith } 202b5b17502SBarry Smith 203b5b17502SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) { 204*6497c311SBarry Smith PetscInt nz2; 205b5b17502SBarry Smith if (!(flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP)) { 206b5b17502SBarry Smith #if defined(PETSC_USE_BACKWARD_LOOP) 207b5b17502SBarry Smith v = aa + ai[m] - 1; 208b5b17502SBarry Smith vj = aj + ai[m] - 1; 209b5b17502SBarry Smith for (i = m - 1; i >= 0; i--) { 210b5b17502SBarry Smith sum = b[i]; 211b5b17502SBarry Smith nz = ai[i + 1] - ai[i] - 1; 2129371c9d4SSatish Balay { 2139371c9d4SSatish Balay PetscInt __i; 2149371c9d4SSatish Balay for (__i = 0; __i < nz; __i++) sum -= v[-__i] * x[vj[-__i]]; 2159371c9d4SSatish Balay } 216b5b17502SBarry Smith #else 217b5b17502SBarry Smith v = aa + ai[m - 1] + 1; 218b5b17502SBarry Smith vj = aj + ai[m - 1] + 1; 219b5b17502SBarry Smith nz = 0; 220b5b17502SBarry Smith for (i = m - 1; i >= 0; i--) { 221b5b17502SBarry Smith sum = b[i]; 2226c11f900SJed Brown nz2 = ai[i] - ai[PetscMax(i - 1, 0)] - 1; /* avoid referencing ai[-1], nonsense nz2 is okay on last iteration */ 22350d8bf02SJed Brown PETSC_Prefetch(v - nz2 - 1, 0, PETSC_PREFETCH_HINT_NTA); 22450d8bf02SJed Brown PETSC_Prefetch(vj - nz2 - 1, 0, PETSC_PREFETCH_HINT_NTA); 225b5b17502SBarry Smith PetscSparseDenseMinusDot(sum, x, v, vj, nz); 226b5b17502SBarry Smith nz = nz2; 227b5b17502SBarry Smith #endif 228b5b17502SBarry Smith x[i] = omega * sum * aidiag[i]; 229b5b17502SBarry Smith v -= nz + 1; 230b5b17502SBarry Smith vj -= nz + 1; 231b5b17502SBarry Smith } 2329566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * a->nz)); 233b5b17502SBarry Smith } else { 234b5b17502SBarry Smith v = aa + ai[m - 1] + 1; 235b5b17502SBarry Smith vj = aj + ai[m - 1] + 1; 236b5b17502SBarry Smith nz = 0; 237b5b17502SBarry Smith for (i = m - 1; i >= 0; i--) { 238b5b17502SBarry Smith sum = t[i]; 2396c11f900SJed Brown nz2 = ai[i] - ai[PetscMax(i - 1, 0)] - 1; /* avoid referencing ai[-1], nonsense nz2 is okay on last iteration */ 24050d8bf02SJed Brown PETSC_Prefetch(v - nz2 - 1, 0, PETSC_PREFETCH_HINT_NTA); 24150d8bf02SJed Brown PETSC_Prefetch(vj - nz2 - 1, 0, PETSC_PREFETCH_HINT_NTA); 242b5b17502SBarry Smith PetscSparseDenseMinusDot(sum, x, v, vj, nz); 243b5b17502SBarry Smith x[i] = (1 - omega) * x[i] + omega * sum * aidiag[i]; 244b5b17502SBarry Smith nz = nz2; 245b5b17502SBarry Smith v -= nz + 1; 246b5b17502SBarry Smith vj -= nz + 1; 247b5b17502SBarry Smith } 2489566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * a->nz)); 249b5b17502SBarry Smith } 250b5b17502SBarry Smith } 251b5b17502SBarry Smith its--; 252b5b17502SBarry Smith } 253b5b17502SBarry Smith 254b5b17502SBarry Smith while (its--) { 255b5b17502SBarry Smith /* 256b5b17502SBarry Smith forward sweep: 257b5b17502SBarry Smith for i=0,...,m-1: 258b5b17502SBarry Smith sum[i] = (b[i] - U(i,:)x)/d[i]; 259b5b17502SBarry Smith x[i] = (1-omega)x[i] + omega*sum[i]; 260b5b17502SBarry Smith b = b - x[i]*U^T(i,:); 261b5b17502SBarry Smith 262b5b17502SBarry Smith */ 263b5b17502SBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) { 2649566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t, b, m)); 265b5b17502SBarry Smith 266b5b17502SBarry Smith for (i = 0; i < m; i++) { 2679371c9d4SSatish Balay v = aa + ai[i] + 1; 2689371c9d4SSatish Balay v1 = v; 2699371c9d4SSatish Balay vj = aj + ai[i] + 1; 2709371c9d4SSatish Balay vj1 = vj; 2719371c9d4SSatish Balay nz = ai[i + 1] - ai[i] - 1; 2729371c9d4SSatish Balay nz1 = nz; 273b5b17502SBarry Smith sum = t[i]; 274b5b17502SBarry Smith while (nz1--) sum -= (*v1++) * x[*vj1++]; 275b5b17502SBarry Smith x[i] = (1 - omega) * x[i] + omega * sum * aidiag[i]; 276b5b17502SBarry Smith while (nz--) t[*vj++] -= x[i] * (*v++); 277b5b17502SBarry Smith } 2789566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(4.0 * a->nz)); 279b5b17502SBarry Smith } 280b5b17502SBarry Smith 281b5b17502SBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) { 282b5b17502SBarry Smith /* 283b5b17502SBarry Smith backward sweep: 284b5b17502SBarry Smith b = b - x[i]*U^T(i,:), i=0,...,n-2 285b5b17502SBarry Smith for i=m-1,...,0: 286b5b17502SBarry Smith sum[i] = (b[i] - U(i,:)x)/d[i]; 287b5b17502SBarry Smith x[i] = (1-omega)x[i] + omega*sum[i]; 288b5b17502SBarry Smith */ 289b5b17502SBarry Smith /* if there was a forward sweep done above then I thing the next two for loops are not needed */ 2909566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(t, b, m)); 291b5b17502SBarry Smith 292b5b17502SBarry Smith for (i = 0; i < m - 1; i++) { /* update rhs */ 293b5b17502SBarry Smith v = aa + ai[i] + 1; 294b5b17502SBarry Smith vj = aj + ai[i] + 1; 295b5b17502SBarry Smith nz = ai[i + 1] - ai[i] - 1; 296b5b17502SBarry Smith while (nz--) t[*vj++] -= x[i] * (*v++); 297b5b17502SBarry Smith } 2989566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * (a->nz - m))); 299b5b17502SBarry Smith for (i = m - 1; i >= 0; i--) { 300b5b17502SBarry Smith v = aa + ai[i] + 1; 301b5b17502SBarry Smith vj = aj + ai[i] + 1; 302b5b17502SBarry Smith nz = ai[i + 1] - ai[i] - 1; 303b5b17502SBarry Smith sum = t[i]; 304b5b17502SBarry Smith while (nz--) sum -= x[*vj++] * (*v++); 305b5b17502SBarry Smith x[i] = (1 - omega) * x[i] + omega * sum * aidiag[i]; 306b5b17502SBarry Smith } 3079566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * (a->nz + m))); 308b5b17502SBarry Smith } 309b5b17502SBarry Smith } 310b5b17502SBarry Smith 3119566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xx, &x)); 3129566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(bb, &b)); 3133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 314b5b17502SBarry Smith } 315