SCALAPACK 2.2.2
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
pzgbsv.f
Go to the documentation of this file.
1 SUBROUTINE pzgbsv( N, BWL, BWU, NRHS, A, JA, DESCA, IPIV, B, IB,
2 $ DESCB, WORK, LWORK, INFO )
3*
4*
5*
6* -- ScaLAPACK routine (version 1.7) --
7* University of Tennessee, Knoxville, Oak Ridge National Laboratory,
8* and University of California, Berkeley.
9* November 15, 1997
10*
11* .. Scalar Arguments ..
12 INTEGER BWL, BWU, IB, INFO, JA, LWORK, N, NRHS
13* ..
14* .. Array Arguments ..
15 INTEGER DESCA( * ), DESCB( * ), IPIV( * )
16 COMPLEX*16 A( * ), B( * ), WORK( * )
17* ..
18*
19*
20* Purpose
21* =======
22*
23* PZGBSV solves a system of linear equations
24*
25* A(1:N, JA:JA+N-1) * X = B(IB:IB+N-1, 1:NRHS)
26*
27* where A(1:N, JA:JA+N-1) is an N-by-N complex
28* banded distributed
29* matrix with bandwidth BWL, BWU.
30*
31* Gaussian elimination with pivoting
32* is used to factor a reordering
33* of the matrix into P L U.
34*
35* See PZGBTRF and PZGBTRS for details.
36*
37* =====================================================================
38*
39* Arguments
40* =========
41*
42*
43* N (global input) INTEGER
44* The number of rows and columns to be operated on, i.e. the
45* order of the distributed submatrix A(1:N, JA:JA+N-1). N >= 0.
46*
47* BWL (global input) INTEGER
48* Number of subdiagonals. 0 <= BWL <= N-1
49*
50* BWU (global input) INTEGER
51* Number of superdiagonals. 0 <= BWU <= N-1
52*
53* NRHS (global input) INTEGER
54* The number of right hand sides, i.e., the number of columns
55* of the distributed submatrix B(IB:IB+N-1, 1:NRHS).
56* NRHS >= 0.
57*
58* A (local input/local output) COMPLEX*16 pointer into
59* local memory to an array with first dimension
60* LLD_A >=(2*bwl+2*bwu+1) (stored in DESCA).
61* On entry, this array contains the local pieces of the
62* This local portion is stored in the packed banded format
63* used in LAPACK. Please see the Notes below and the
64* ScaLAPACK manual for more detail on the format of
65* distributed matrices.
66* On exit, this array contains information containing details
67* of the factorization.
68* Note that permutations are performed on the matrix, so that
69* the factors returned are different from those returned
70* by LAPACK.
71*
72* JA (global input) INTEGER
73* The index in the global array A that points to the start of
74* the matrix to be operated on (which may be either all of A
75* or a submatrix of A).
76*
77* DESCA (global and local input) INTEGER array of dimension DLEN.
78* if 1D type (DTYPE_A=501), DLEN >= 7;
79* if 2D type (DTYPE_A=1), DLEN >= 9 .
80* The array descriptor for the distributed matrix A.
81* Contains information of mapping of A to memory. Please
82* see NOTES below for full description and options.
83*
84* IPIV (local output) INTEGER array, dimension >= DESCA( NB ).
85* Pivot indices for local factorizations.
86* Users *should not* alter the contents between
87* factorization and solve.
88*
89* B (local input/local output) COMPLEX*16 pointer into
90* local memory to an array of local lead dimension lld_b>=NB.
91* On entry, this array contains the
92* the local pieces of the right hand sides
93* B(IB:IB+N-1, 1:NRHS).
94* On exit, this contains the local piece of the solutions
95* distributed matrix X.
96*
97* IB (global input) INTEGER
98* The row index in the global array B that points to the first
99* row of the matrix to be operated on (which may be either
100* all of B or a submatrix of B).
101*
102* DESCB (global and local input) INTEGER array of dimension DLEN.
103* if 1D type (DTYPE_B=502), DLEN >=7;
104* if 2D type (DTYPE_B=1), DLEN >= 9.
105* The array descriptor for the distributed matrix B.
106* Contains information of mapping of B to memory. Please
107* see NOTES below for full description and options.
108*
109* WORK (local workspace/local output)
110* COMPLEX*16 temporary workspace. This space may
111* be overwritten in between calls to routines. WORK must be
112* the size given in LWORK.
113* On exit, WORK( 1 ) contains the minimal LWORK.
114*
115* LWORK (local input or global input) INTEGER
116* Size of user-input workspace WORK.
117* If LWORK is too small, the minimal acceptable size will be
118* returned in WORK(1) and an error code is returned. LWORK>=
119* (NB+bwu)*(bwl+bwu)+6*(bwl+bwu)*(bwl+2*bwu)
120* +max(NRHS*(NB+2*bwl+4*bwu), 1)
121*
122* INFO (global output) INTEGER
123* = 0: successful exit
124* < 0: If the i-th argument is an array and the j-entry had
125* an illegal value, then INFO = -(i*100+j), if the i-th
126* argument is a scalar and had an illegal value, then
127* INFO = -i.
128* > 0: If INFO = K<=NPROCS, the submatrix stored on processor
129* INFO and factored locally was not
130* nonsingular, and
131* the factorization was not completed.
132* If INFO = K>NPROCS, the submatrix stored on processor
133* INFO-NPROCS representing interactions with other
134* processors was not
135* nonsingular,
136* and the factorization was not completed.
137*
138* =====================================================================
139*
140*
141* Restrictions
142* ============
143*
144* The following are restrictions on the input parameters. Some of these
145* are temporary and will be removed in future releases, while others
146* may reflect fundamental technical limitations.
147*
148* Non-cyclic restriction: VERY IMPORTANT!
149* P*NB>= mod(JA-1,NB)+N.
150* The mapping for matrices must be blocked, reflecting the nature
151* of the divide and conquer algorithm as a task-parallel algorithm.
152* This formula in words is: no processor may have more than one
153* chunk of the matrix.
154*
155* Blocksize cannot be too small:
156* If the matrix spans more than one processor, the following
157* restriction on NB, the size of each block on each processor,
158* must hold:
159* NB >= (BWL+BWU)+1
160* The bulk of parallel computation is done on the matrix of size
161* O(NB) on each processor. If this is too small, divide and conquer
162* is a poor choice of algorithm.
163*
164* Submatrix reference:
165* JA = IB
166* Alignment restriction that prevents unnecessary communication.
167*
168*
169* =====================================================================
170*
171*
172* Notes
173* =====
174*
175* If the factorization routine and the solve routine are to be called
176* separately (to solve various sets of righthand sides using the same
177* coefficient matrix), the auxiliary space AF *must not be altered*
178* between calls to the factorization routine and the solve routine.
179*
180* The best algorithm for solving banded and tridiagonal linear systems
181* depends on a variety of parameters, especially the bandwidth.
182* Currently, only algorithms designed for the case N/P >> bw are
183* implemented. These go by many names, including Divide and Conquer,
184* Partitioning, domain decomposition-type, etc.
185*
186* Algorithm description: Divide and Conquer
187*
188* The Divide and Conqer algorithm assumes the matrix is narrowly
189* banded compared with the number of equations. In this situation,
190* it is best to distribute the input matrix A one-dimensionally,
191* with columns atomic and rows divided amongst the processes.
192* The basic algorithm divides the banded matrix up into
193* P pieces with one stored on each processor,
194* and then proceeds in 2 phases for the factorization or 3 for the
195* solution of a linear system.
196* 1) Local Phase:
197* The individual pieces are factored independently and in
198* parallel. These factors are applied to the matrix creating
199* fillin, which is stored in a non-inspectable way in auxiliary
200* space AF. Mathematically, this is equivalent to reordering
201* the matrix A as P A P^T and then factoring the principal
202* leading submatrix of size equal to the sum of the sizes of
203* the matrices factored on each processor. The factors of
204* these submatrices overwrite the corresponding parts of A
205* in memory.
206* 2) Reduced System Phase:
207* A small (max(bwl,bwu)* (P-1)) system is formed representing
208* interaction of the larger blocks, and is stored (as are its
209* factors) in the space AF. A parallel Block Cyclic Reduction
210* algorithm is used. For a linear system, a parallel front solve
211* followed by an analagous backsolve, both using the structure
212* of the factored matrix, are performed.
213* 3) Backsubsitution Phase:
214* For a linear system, a local backsubstitution is performed on
215* each processor in parallel.
216*
217*
218* Descriptors
219* ===========
220*
221* Descriptors now have *types* and differ from ScaLAPACK 1.0.
222*
223* Note: banded codes can use either the old two dimensional
224* or new one-dimensional descriptors, though the processor grid in
225* both cases *must be one-dimensional*. We describe both types below.
226*
227* Each global data object is described by an associated description
228* vector. This vector stores the information required to establish
229* the mapping between an object element and its corresponding process
230* and memory location.
231*
232* Let A be a generic term for any 2D block cyclicly distributed array.
233* Such a global array has an associated description vector DESCA.
234* In the following comments, the character _ should be read as
235* "of the global array".
236*
237* NOTATION STORED IN EXPLANATION
238* --------------- -------------- --------------------------------------
239* DTYPE_A(global) DESCA( DTYPE_ )The descriptor type. In this case,
240* DTYPE_A = 1.
241* CTXT_A (global) DESCA( CTXT_ ) The BLACS context handle, indicating
242* the BLACS process grid A is distribu-
243* ted over. The context itself is glo-
244* bal, but the handle (the integer
245* value) may vary.
246* M_A (global) DESCA( M_ ) The number of rows in the global
247* array A.
248* N_A (global) DESCA( N_ ) The number of columns in the global
249* array A.
250* MB_A (global) DESCA( MB_ ) The blocking factor used to distribute
251* the rows of the array.
252* NB_A (global) DESCA( NB_ ) The blocking factor used to distribute
253* the columns of the array.
254* RSRC_A (global) DESCA( RSRC_ ) The process row over which the first
255* row of the array A is distributed.
256* CSRC_A (global) DESCA( CSRC_ ) The process column over which the
257* first column of the array A is
258* distributed.
259* LLD_A (local) DESCA( LLD_ ) The leading dimension of the local
260* array. LLD_A >= MAX(1,LOCr(M_A)).
261*
262* Let K be the number of rows or columns of a distributed matrix,
263* and assume that its process grid has dimension p x q.
264* LOCr( K ) denotes the number of elements of K that a process
265* would receive if K were distributed over the p processes of its
266* process column.
267* Similarly, LOCc( K ) denotes the number of elements of K that a
268* process would receive if K were distributed over the q processes of
269* its process row.
270* The values of LOCr() and LOCc() may be determined via a call to the
271* ScaLAPACK tool function, NUMROC:
272* LOCr( M ) = NUMROC( M, MB_A, MYROW, RSRC_A, NPROW ),
273* LOCc( N ) = NUMROC( N, NB_A, MYCOL, CSRC_A, NPCOL ).
274* An upper bound for these quantities may be computed by:
275* LOCr( M ) <= ceil( ceil(M/MB_A)/NPROW )*MB_A
276* LOCc( N ) <= ceil( ceil(N/NB_A)/NPCOL )*NB_A
277*
278*
279* One-dimensional descriptors:
280*
281* One-dimensional descriptors are a new addition to ScaLAPACK since
282* version 1.0. They simplify and shorten the descriptor for 1D
283* arrays.
284*
285* Since ScaLAPACK supports two-dimensional arrays as the fundamental
286* object, we allow 1D arrays to be distributed either over the
287* first dimension of the array (as if the grid were P-by-1) or the
288* 2nd dimension (as if the grid were 1-by-P). This choice is
289* indicated by the descriptor type (501 or 502)
290* as described below.
291*
292* IMPORTANT NOTE: the actual BLACS grid represented by the
293* CTXT entry in the descriptor may be *either* P-by-1 or 1-by-P
294* irrespective of which one-dimensional descriptor type
295* (501 or 502) is input.
296* This routine will interpret the grid properly either way.
297* ScaLAPACK routines *do not support intercontext operations* so that
298* the grid passed to a single ScaLAPACK routine *must be the same*
299* for all array descriptors passed to that routine.
300*
301* NOTE: In all cases where 1D descriptors are used, 2D descriptors
302* may also be used, since a one-dimensional array is a special case
303* of a two-dimensional array with one dimension of size unity.
304* The two-dimensional array used in this case *must* be of the
305* proper orientation:
306* If the appropriate one-dimensional descriptor is DTYPEA=501
307* (1 by P type), then the two dimensional descriptor must
308* have a CTXT value that refers to a 1 by P BLACS grid;
309* If the appropriate one-dimensional descriptor is DTYPEA=502
310* (P by 1 type), then the two dimensional descriptor must
311* have a CTXT value that refers to a P by 1 BLACS grid.
312*
313*
314* Summary of allowed descriptors, types, and BLACS grids:
315* DTYPE 501 502 1 1
316* BLACS grid 1xP or Px1 1xP or Px1 1xP Px1
317* -----------------------------------------------------
318* A OK NO OK NO
319* B NO OK NO OK
320*
321* Note that a consequence of this chart is that it is not possible
322* for *both* DTYPE_A and DTYPE_B to be 2D_type(1), as these lead
323* to opposite requirements for the orientation of the BLACS grid,
324* and as noted before, the *same* BLACS context must be used in
325* all descriptors in a single ScaLAPACK subroutine call.
326*
327* Let A be a generic term for any 1D block cyclicly distributed array.
328* Such a global array has an associated description vector DESCA.
329* In the following comments, the character _ should be read as
330* "of the global array".
331*
332* NOTATION STORED IN EXPLANATION
333* --------------- ---------- ------------------------------------------
334* DTYPE_A(global) DESCA( 1 ) The descriptor type. For 1D grids,
335* TYPE_A = 501: 1-by-P grid.
336* TYPE_A = 502: P-by-1 grid.
337* CTXT_A (global) DESCA( 2 ) The BLACS context handle, indicating
338* the BLACS process grid A is distribu-
339* ted over. The context itself is glo-
340* bal, but the handle (the integer
341* value) may vary.
342* N_A (global) DESCA( 3 ) The size of the array dimension being
343* distributed.
344* NB_A (global) DESCA( 4 ) The blocking factor used to distribute
345* the distributed dimension of the array.
346* SRC_A (global) DESCA( 5 ) The process row or column over which the
347* first row or column of the array
348* is distributed.
349* LLD_A (local) DESCA( 6 ) The leading dimension of the local array
350* storing the local blocks of the distri-
351* buted array A. Minimum value of LLD_A
352* depends on TYPE_A.
353* TYPE_A = 501: LLD_A >=
354* size of undistributed dimension, 1.
355* TYPE_A = 502: LLD_A >=NB_A, 1.
356* Reserved DESCA( 7 ) Reserved for future use.
357*
358*
359*
360* =====================================================================
361*
362* Code Developer: Andrew J. Cleary, University of Tennessee.
363* Current address: Lawrence Livermore National Labs.
364* This version released: August, 2001.
365*
366* =====================================================================
367*
368* ..
369* .. Parameters ..
370 DOUBLE PRECISION ONE, ZERO
371 parameter( one = 1.0d+0 )
372 parameter( zero = 0.0d+0 )
373 COMPLEX*16 CONE, CZERO
374 parameter( cone = ( 1.0d+0, 0.0d+0 ) )
375 parameter( czero = ( 0.0d+0, 0.0d+0 ) )
376 INTEGER INT_ONE
377 parameter( int_one = 1 )
378 INTEGER DESCMULT, BIGNUM
379 parameter(descmult = 100, bignum = descmult * descmult)
380 INTEGER BLOCK_CYCLIC_2D, CSRC_, CTXT_, DLEN_, DTYPE_,
381 $ lld_, mb_, m_, nb_, n_, rsrc_
382 parameter( block_cyclic_2d = 1, dlen_ = 9, dtype_ = 1,
383 $ ctxt_ = 2, m_ = 3, n_ = 4, mb_ = 5, nb_ = 6,
384 $ rsrc_ = 7, csrc_ = 8, lld_ = 9 )
385* ..
386* .. Local Scalars ..
387 INTEGER ICTXT, MYCOL, MYROW, NB, NPCOL, NPROW,
388 $ ws_factor
389* ..
390* .. External Subroutines ..
391 EXTERNAL pxerbla, pzgbtrf, pzgbtrs
392* ..
393* .. Executable Statements ..
394*
395* Note: to avoid duplication, most error checking is not performed
396* in this routine and is left to routines
397* PZGBTRF and PZGBTRS.
398*
399* Begin main code
400*
401 info = 0
402*
403* Get block size to calculate workspace requirements
404*
405 IF( desca( dtype_ ) .EQ. block_cyclic_2d ) THEN
406 nb = desca( nb_ )
407 ictxt = desca( ctxt_ )
408 ELSEIF( desca( dtype_ ) .EQ. 501 ) THEN
409 nb = desca( 4 )
410 ictxt = desca( 2 )
411 ELSE
412 info = -( 6*100 + dtype_ )
413 CALL pxerbla( ictxt,
414 $ 'PZGBSV',
415 $ -info )
416 RETURN
417 ENDIF
418*
419 CALL blacs_gridinfo( ictxt, nprow, npcol, myrow, mycol )
420*
421*
422* Size needed for AF in factorization
423*
424 ws_factor = (nb+bwu)*(bwl+bwu)+6*(bwl+bwu)*(bwl+2*bwu)
425*
426* Factor the matrix
427*
428 CALL pzgbtrf( n, bwl, bwu, a, ja, desca, ipiv, work,
429 $ min( lwork, ws_factor ), work( 1+ws_factor ),
430 $ lwork-ws_factor, info )
431*
432* Check info for error conditions
433*
434 IF( info.NE.0 ) THEN
435 IF( info .LT. 0 ) THEN
436 CALL pxerbla( ictxt, 'PZGBSV', -info )
437 ENDIF
438 RETURN
439 END IF
440*
441* Solve the system using the factorization
442*
443 CALL pzgbtrs( 'N', n, bwl, bwu, nrhs, a, ja, desca, ipiv, b, ib,
444 $ descb, work, min( lwork, ws_factor ),
445 $ work( 1+ws_factor), lwork-ws_factor, info )
446*
447* Check info for error conditions
448*
449 IF( info.NE.0 ) THEN
450 CALL pxerbla( ictxt, 'PZGBSV', -info )
451 RETURN
452 END IF
453*
454 RETURN
455*
456* End of PZGBSV
457*
458 END
#define min(A, B)
Definition pcgemr.c:181
subroutine pxerbla(ictxt, srname, info)
Definition pxerbla.f:2
subroutine pzgbsv(n, bwl, bwu, nrhs, a, ja, desca, ipiv, b, ib, descb, work, lwork, info)
Definition pzgbsv.f:3
subroutine pzgbtrf(n, bwl, bwu, a, ja, desca, ipiv, af, laf, work, lwork, info)
Definition pzgbtrf.f:3
subroutine pzgbtrs(trans, n, bwl, bwu, nrhs, a, ja, desca, ipiv, b, ib, descb, af, laf, work, lwork, info)
Definition pzgbtrs.f:3