Actual source code: unary.c

  1: /* unary.f -- translated by f2c (version of 25 March 1992  12:58:56).

  3:         This code is protected by the GNU copyright. See the file 
  4:      ilut.c in this directory for the full copyright. See below for the Author.
  5: */
 6:  #include petsc.h
  7: /* ----------------------------------------------------------------------- */
  8: static PetscErrorCode SPARSEKIT2rperm(PetscInt *nrow,PetscScalar *a,PetscInt *ja,PetscInt *ia,PetscScalar *ao,PetscInt *jao,PetscInt *iao,PetscInt *perm,PetscInt *job)
  9: {
 10:     /* System generated locals */
 11:     PetscInt i__1,i__2;

 13:     /* Local variables */
 14:     PetscInt i,j,k,ii,ko;
 15:     PetscInt values;

 17: /* -----------------------------------------------------------------------
 18:  */
 19: /* this subroutine permutes the rows of a matrix in CSR format. */
 20: /* rperm  computes B = P A  where P is a permutation matrix. */
 21: /* the permutation P is defined through the array perm: for each j, */
 22: /* perm(j) represents the destination row number of row number j. */
 23: /* Youcef Saad -- recoded Jan 28, 1991. */
 24: /* -----------------------------------------------------------------------
 25:  */
 26: /* on entry: */
 27: /* ---------- */
 28: /* n         = dimension of the matrix */
 29: /* a, ja, ia = input matrix in csr format */
 30: /* perm         = integer array of length nrow containing the permutation arrays 
 31: */
 32: /*           for the rows: perm(i) is the destination of row i in the */
 33: /*         permuted matrix. */
 34: /*         ---> a(i,j) in the original matrix becomes a(perm(i),j) */
 35: /*         in the output  matrix. */

 37: /* job        = integer indicating the work to be done: */
 38: /*                 job = 1        permute a, ja, ia into ao, jao, iao */
 39: /*                       (including the copying of real values ao and */
 40: /*                       the array iao). */
 41: /*                 job .ne. 1 :  ignore real values. */
 42: /*                     (in which case arrays a and ao are not needed nor 
 43: */
 44: /*                      used). */

 46: /* ------------ */
 47: /* on return: */
 48: /* ------------ */
 49: /* ao, jao, iao = input matrix in a, ja, ia format */
 50: /* note : */
 51: /*        if (job.ne.1)  then the arrays a and ao are not used. */
 52: /* ----------------------------------------------------------------------c
 53:  */
 54: /*           Y. Saad, May  2, 1990                                      c 
 55: */
 56: /* ----------------------------------------------------------------------c
 57:  */
 58:     /* Parameter adjustments */
 59:     --perm;
 60:     --iao;
 61:     --jao;
 62:     --ao;
 63:     --ia;
 64:     --ja;
 65:     --a;

 67:     /* Function Body */
 68:     values = *job == 1;

 70: /*     determine pointers for output matix. */

 72:     i__1 = *nrow;
 73:     for (j = 1; j <= i__1; ++j) {
 74:         i = perm[j];
 75:         iao[i + 1] = ia[j + 1] - ia[j];
 76: /* L50: */
 77:     }

 79: /* get pointers from lengths */

 81:     iao[1] = 1;
 82:     i__1 = *nrow;
 83:     for (j = 1; j <= i__1; ++j) {
 84:         iao[j + 1] += iao[j];
 85: /* L51: */
 86:     }

 88: /* copying */

 90:     i__1 = *nrow;
 91:     for (ii = 1; ii <= i__1; ++ii) {

 93: /* old row = ii  -- new row = iperm(ii) -- ko = new pointer */

 95:         ko = iao[perm[ii]];
 96:         i__2 = ia[ii + 1] - 1;
 97:         for (k = ia[ii]; k <= i__2; ++k) {
 98:             jao[ko] = ja[k];
 99:             if (values) {
100:                 ao[ko] = a[k];
101:             }
102:             ++ko;
103: /* L60: */
104:         }
105: /* L100: */
106:     }

108:     return 0;
109: /* ---------end-of-rperm -------------------------------------------------
110:  */
111: /* -----------------------------------------------------------------------
112:  */
113: } /* rperm_ */

115: /* ----------------------------------------------------------------------- */
116: static PetscErrorCode SPARSEKIT2cperm(PetscInt *nrow,PetscScalar * a,PetscInt * ja,PetscInt * ia,PetscScalar * ao,PetscInt * jao,PetscInt * iao,PetscInt * perm,PetscInt * job)
117: {
118:     /* System generated locals */
119:     PetscInt i__1;

121:     /* Local variables */
122:     PetscInt i,k,nnz;

124: /* -----------------------------------------------------------------------
125:  */
126: /* this subroutine permutes the columns of a matrix a, ja, ia. */
127: /* the result is written in the output matrix  ao, jao, iao. */
128: /* cperm computes B = A P, where  P is a permutation matrix */
129: /* that maps column j into column perm(j), i.e., on return */
130: /*      a(i,j) becomes a(i,perm(j)) in new matrix */
131: /* Y. Saad, May 2, 1990 / modified Jan. 28, 1991. */
132: /* -----------------------------------------------------------------------
133:  */
134: /* on entry: */
135: /* ---------- */
136: /* nrow         = row dimension of the matrix */

138: /* a, ja, ia = input matrix in csr format. */

140: /* perm        = integer array of length ncol (number of columns of A */
141: /*         containing the permutation array  the columns: */
142: /*         a(i,j) in the original matrix becomes a(i,perm(j)) */
143: /*         in the output matrix. */

145: /* job        = integer indicating the work to be done: */
146: /*                 job = 1        permute a, ja, ia into ao, jao, iao */
147: /*                       (including the copying of real values ao and */
148: /*                       the array iao). */
149: /*                 job .ne. 1 :  ignore real values ao and ignore iao. */

151: /* ------------ */
152: /* on return: */
153: /* ------------ */
154: /* ao, jao, iao = input matrix in a, ja, ia format (array ao not needed) 
155: */

157: /* Notes: */
158: /* ------- */
159: /* 1. if job=1 then ao, iao are not used. */
160: /* 2. This routine is in place: ja, jao can be the same. */
161: /* 3. If the matrix is initially sorted (by increasing column number) */
162: /*    then ao,jao,iao  may not be on return. */

164: /* ----------------------------------------------------------------------c
165:  */
166: /* local parameters: */

168:     /* Parameter adjustments */
169:     --perm;
170:     --iao;
171:     --jao;
172:     --ao;
173:     --ia;
174:     --ja;
175:     --a;

177:     /* Function Body */
178:     nnz = ia[*nrow + 1] - 1;
179:     i__1 = nnz;
180:     for (k = 1; k <= i__1; ++k) {
181:         jao[k] = perm[ja[k]];
182: /* L100: */
183:     }

185: /*     done with ja array. return if no need to touch values. */

187:     if (*job != 1) {
188:         return 0;
189:     }

191: /* else get new pointers -- and copy values too. */

193:     i__1 = *nrow + 1;
194:     for (i = 1; i <= i__1; ++i) {
195:         iao[i] = ia[i];
196: /* L1: */
197:     }

199:     i__1 = nnz;
200:     for (k = 1; k <= i__1; ++k) {
201:         ao[k] = a[k];
202: /* L2: */
203:     }

205:     return 0;
206: /* ---------end-of-cperm--------------------------------------------------
207:  */
208: /* -----------------------------------------------------------------------
209:  */
210: } /* cperm_ */

212: /* ----------------------------------------------------------------------- */
213: PetscErrorCode SPARSEKIT2dperm(PetscInt *nrow,PetscScalar *a,PetscInt *ja,PetscInt *ia,PetscScalar *ao,PetscInt *jao,PetscInt *iao,PetscInt *perm,PetscInt *qperm,PetscInt *job)
214: {
215:     PetscInt locjob;

217: /* -----------------------------------------------------------------------
218:  */
219: /* This routine permutes the rows and columns of a matrix stored in CSR */

221: /* format. i.e., it computes P A Q, where P, Q are permutation matrices. 
222: */
223: /* P maps row i into row perm(i) and Q maps column j into column qperm(j):
224:  */
225: /*      a(i,j)    becomes   a(perm(i),qperm(j)) in new matrix */
226: /* In the particular case where Q is the transpose of P (symmetric */
227: /* permutation of A) then qperm is not needed. */
228: /* note that qperm should be of length ncol (number of columns) but this 
229: */
230: /* is not checked. */
231: /* -----------------------------------------------------------------------
232:  */
233: /* Y. Saad, Sep. 21 1989 / recoded Jan. 28 1991. */
234: /* -----------------------------------------------------------------------
235:  */
236: /* on entry: */
237: /* ---------- */
238: /* n         = dimension of the matrix */
239: /* a, ja, */
240: /*    ia = input matrix in a, ja, ia format */
241: /* perm         = integer array of length n containing the permutation arrays */
242: /*           for the rows: perm(i) is the destination of row i in the */
243: /*         permuted matrix -- also the destination of column i in case */
244: /*         permutation is symmetric (job .le. 2) */

246: /* qperm        = same thing for the columns. This should be provided only */
247: /*         if job=3 or job=4, i.e., only in the case of a nonsymmetric */
248: /*           permutation of rows and columns. Otherwise qperm is a dummy */

250: /* job        = integer indicating the work to be done: */
251: /* * job = 1,2 permutation is symmetric  Ao :== P * A * transp(P) */
252: /*                 job = 1        permute a, ja, ia into ao, jao, iao */
253: /*                 job = 2 permute matrix ignoring real values. */
254: /* * job = 3,4 permutation is non-symmetric  Ao :== P * A * Q */
255: /*                 job = 3        permute a, ja, ia into ao, jao, iao */
256: /*                 job = 4 permute matrix ignoring real values. */

258: /* on return: */
259: /* ----------- */
260: /* ao, jao, iao = input matrix in a, ja, ia format */

262: /* in case job .eq. 2 or job .eq. 4, a and ao are never referred to */
263: /* and can be dummy arguments. */
264: /* Notes: */
265: /* ------- */
266: /*  1) algorithm is in place */
267: /*  2) column indices may not be sorted on return even  though they may be
268:  */
269: /*     on entry. */
270: /* ----------------------------------------------------------------------c
271:  */
272: /* local variables */

274: /*     locjob indicates whether or not real values must be copied. */

276:     /* Parameter adjustments */
277:     --qperm;
278:     --perm;
279:     --iao;
280:     --jao;
281:     --ao;
282:     --ia;
283:     --ja;
284:     --a;

286:     /* Function Body */
287:     locjob = *job % 2;

289: /* permute rows first */

291:     SPARSEKIT2rperm(nrow, &a[1], &ja[1], &ia[1], &ao[1], &jao[1], &iao[1], &perm[1], &locjob);

293: /* then permute columns */

295:     locjob = 0;

297:     if (*job <= 2) {
298:         SPARSEKIT2cperm(nrow, &ao[1], &jao[1], &iao[1], &ao[1], &jao[1], &iao[1], &perm[1], &locjob);
299:     } else {
300:         SPARSEKIT2cperm(nrow, &ao[1], &jao[1], &iao[1], &ao[1], &jao[1], &iao[1], &qperm[1], &locjob);
301:     }

303:     return 0;
304: /* -------end-of-dperm----------------------------------------------------
305:  */
306: /* -----------------------------------------------------------------------
307:  */
308: } /* dperm_ */

310: /* ----------------------------------------------------------------------- */
311: PetscErrorCode SPARSEKIT2msrcsr(PetscInt *n,PetscScalar * a,PetscInt * ja,PetscScalar * ao,PetscInt * jao,PetscInt * iao,PetscScalar * wk,PetscInt * iwk)
312: {
313:     /* System generated locals */
314:     PetscInt i__1, i__2;

316:     /* Local variables */
317:     PetscInt iptr;
318:     PetscInt added;
319:     PetscInt i, j, k, idiag, ii;

321: /* -----------------------------------------------------------------------
322:  */
323: /*       Modified - Sparse Row  to   Compressed Sparse Row */

325: /* -----------------------------------------------------------------------
326:  */
327: /* converts a compressed matrix using a separated diagonal */
328: /* (modified sparse row format) in the Compressed Sparse Row */
329: /* format. */
330: /* does not check for zero elements in the diagonal. */


333: /* on entry : */
334: /* --------- */
335: /* n          = row dimension of matrix */
336: /* a, ja      = sparse matrix in msr sparse storage format */
337: /*              see routine csrmsr for details on data structure */

339: /* on return : */
340: /* ----------- */

342: /* ao,jao,iao = output matrix in csr format. */

344: /* work arrays: */
345: /* ------------ */
346: /* wk       = real work array of length n */
347: /* iwk      = integer work array of length n+1 */

349: /* notes: */
350: /*   The original version of this was NOT in place, but has */
351: /*   been modified by adding the vector iwk to be in place. */
352: /*   The original version had ja instead of iwk everywhere in */
353: /*   loop 500.  Modified  Sun 29 May 1994 by R. Bramley (Indiana). */

355: /* -----------------------------------------------------------------------
356:  */
357:     /* Parameter adjustments */
358:     --iwk;
359:     --wk;
360:     --iao;
361:     --jao;
362:     --ao;
363:     --ja;
364:     --a;

366:     /* Function Body */
367:     i__1 = *n;
368:     for (i = 1; i <= i__1; ++i) {
369:         wk[i] = a[i];
370:         iwk[i] = ja[i];
371: /* L1: */
372:     }
373:     iwk[*n + 1] = ja[*n + 1];
374:     iao[1] = 1;
375:     iptr = 1;
376: /* --------- */
377:     i__1 = *n;
378:     for (ii = 1; ii <= i__1; ++ii) {
379:         added = 0;
380:         idiag = iptr + (iwk[ii + 1] - iwk[ii]);
381:         i__2 = iwk[ii + 1] - 1;
382:         for (k = iwk[ii]; k <= i__2; ++k) {
383:             j = ja[k];
384:             if (j < ii) {
385:                 ao[iptr] = a[k];
386:                 jao[iptr] = j;
387:                 ++iptr;
388:             } else if (added) {
389:                 ao[iptr] = a[k];
390:                 jao[iptr] = j;
391:                 ++iptr;
392:             } else {
393: /* add diag element - only reserve a position for it. */
394:                 idiag = iptr;
395:                 ++iptr;
396:                 added = 1;
397: /*     then other element */
398:                 ao[iptr] = a[k];
399:                 jao[iptr] = j;
400:                 ++iptr;
401:             }
402: /* L100: */
403:         }
404:         ao[idiag] = wk[ii];
405:         jao[idiag] = ii;
406:         if (! added) {
407:             ++iptr;
408:         }
409:         iao[ii + 1] = iptr;
410: /* L500: */
411:     }
412:     return 0;
413: } /* msrcsr_ */