LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
ztrsm.f
Go to the documentation of this file.
1*> \brief \b ZTRSM
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8* Definition:
9* ===========
10*
11* SUBROUTINE ZTRSM(SIDE,UPLO,TRANSA,DIAG,M,N,ALPHA,A,LDA,B,LDB)
12*
13* .. Scalar Arguments ..
14* COMPLEX*16 ALPHA
15* INTEGER LDA,LDB,M,N
16* CHARACTER DIAG,SIDE,TRANSA,UPLO
17* ..
18* .. Array Arguments ..
19* COMPLEX*16 A(LDA,*),B(LDB,*)
20* ..
21*
22*
23*> \par Purpose:
24* =============
25*>
26*> \verbatim
27*>
28*> ZTRSM solves one of the matrix equations
29*>
30*> op( A )*X = alpha*B, or X*op( A ) = alpha*B,
31*>
32*> where alpha is a scalar, X and B are m by n matrices, A is a unit, or
33*> non-unit, upper or lower triangular matrix and op( A ) is one of
34*>
35*> op( A ) = A or op( A ) = A**T or op( A ) = A**H.
36*>
37*> The matrix X is overwritten on B.
38*> \endverbatim
39*
40* Arguments:
41* ==========
42*
43*> \param[in] SIDE
44*> \verbatim
45*> SIDE is CHARACTER*1
46*> On entry, SIDE specifies whether op( A ) appears on the left
47*> or right of X as follows:
48*>
49*> SIDE = 'L' or 'l' op( A )*X = alpha*B.
50*>
51*> SIDE = 'R' or 'r' X*op( A ) = alpha*B.
52*> \endverbatim
53*>
54*> \param[in] UPLO
55*> \verbatim
56*> UPLO is CHARACTER*1
57*> On entry, UPLO specifies whether the matrix A is an upper or
58*> lower triangular matrix as follows:
59*>
60*> UPLO = 'U' or 'u' A is an upper triangular matrix.
61*>
62*> UPLO = 'L' or 'l' A is a lower triangular matrix.
63*> \endverbatim
64*>
65*> \param[in] TRANSA
66*> \verbatim
67*> TRANSA is CHARACTER*1
68*> On entry, TRANSA specifies the form of op( A ) to be used in
69*> the matrix multiplication as follows:
70*>
71*> TRANSA = 'N' or 'n' op( A ) = A.
72*>
73*> TRANSA = 'T' or 't' op( A ) = A**T.
74*>
75*> TRANSA = 'C' or 'c' op( A ) = A**H.
76*> \endverbatim
77*>
78*> \param[in] DIAG
79*> \verbatim
80*> DIAG is CHARACTER*1
81*> On entry, DIAG specifies whether or not A is unit triangular
82*> as follows:
83*>
84*> DIAG = 'U' or 'u' A is assumed to be unit triangular.
85*>
86*> DIAG = 'N' or 'n' A is not assumed to be unit
87*> triangular.
88*> \endverbatim
89*>
90*> \param[in] M
91*> \verbatim
92*> M is INTEGER
93*> On entry, M specifies the number of rows of B. M must be at
94*> least zero.
95*> \endverbatim
96*>
97*> \param[in] N
98*> \verbatim
99*> N is INTEGER
100*> On entry, N specifies the number of columns of B. N must be
101*> at least zero.
102*> \endverbatim
103*>
104*> \param[in] ALPHA
105*> \verbatim
106*> ALPHA is COMPLEX*16
107*> On entry, ALPHA specifies the scalar alpha. When alpha is
108*> zero then A is not referenced and B need not be set before
109*> entry.
110*> \endverbatim
111*>
112*> \param[in] A
113*> \verbatim
114*> A is COMPLEX*16 array, dimension ( LDA, k ),
115*> where k is m when SIDE = 'L' or 'l'
116*> and k is n when SIDE = 'R' or 'r'.
117*> Before entry with UPLO = 'U' or 'u', the leading k by k
118*> upper triangular part of the array A must contain the upper
119*> triangular matrix and the strictly lower triangular part of
120*> A is not referenced.
121*> Before entry with UPLO = 'L' or 'l', the leading k by k
122*> lower triangular part of the array A must contain the lower
123*> triangular matrix and the strictly upper triangular part of
124*> A is not referenced.
125*> Note that when DIAG = 'U' or 'u', the diagonal elements of
126*> A are not referenced either, but are assumed to be unity.
127*> \endverbatim
128*>
129*> \param[in] LDA
130*> \verbatim
131*> LDA is INTEGER
132*> On entry, LDA specifies the first dimension of A as declared
133*> in the calling (sub) program. When SIDE = 'L' or 'l' then
134*> LDA must be at least max( 1, m ), when SIDE = 'R' or 'r'
135*> then LDA must be at least max( 1, n ).
136*> \endverbatim
137*>
138*> \param[in,out] B
139*> \verbatim
140*> B is COMPLEX*16 array, dimension ( LDB, N )
141*> Before entry, the leading m by n part of the array B must
142*> contain the right-hand side matrix B, and on exit is
143*> overwritten by the solution matrix X.
144*> \endverbatim
145*>
146*> \param[in] LDB
147*> \verbatim
148*> LDB is INTEGER
149*> On entry, LDB specifies the first dimension of B as declared
150*> in the calling (sub) program. LDB must be at least
151*> max( 1, m ).
152*> \endverbatim
153*
154* Authors:
155* ========
156*
157*> \author Univ. of Tennessee
158*> \author Univ. of California Berkeley
159*> \author Univ. of Colorado Denver
160*> \author NAG Ltd.
161*
162*> \ingroup trsm
163*
164*> \par Further Details:
165* =====================
166*>
167*> \verbatim
168*>
169*> Level 3 Blas routine.
170*>
171*> -- Written on 8-February-1989.
172*> Jack Dongarra, Argonne National Laboratory.
173*> Iain Duff, AERE Harwell.
174*> Jeremy Du Croz, Numerical Algorithms Group Ltd.
175*> Sven Hammarling, Numerical Algorithms Group Ltd.
176*> \endverbatim
177*>
178* =====================================================================
179 SUBROUTINE ztrsm(SIDE,UPLO,TRANSA,DIAG,M,N,ALPHA,A,LDA,B,LDB)
180*
181* -- Reference BLAS level3 routine --
182* -- Reference BLAS is a software package provided by Univ. of Tennessee, --
183* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
184*
185* .. Scalar Arguments ..
186 COMPLEX*16 ALPHA
187 INTEGER LDA,LDB,M,N
188 CHARACTER DIAG,SIDE,TRANSA,UPLO
189* ..
190* .. Array Arguments ..
191 COMPLEX*16 A(LDA,*),B(LDB,*)
192* ..
193*
194* =====================================================================
195*
196* .. External Functions ..
197 LOGICAL LSAME
198 EXTERNAL lsame
199* ..
200* .. External Subroutines ..
201 EXTERNAL xerbla
202* ..
203* .. Intrinsic Functions ..
204 INTRINSIC dconjg,max
205* ..
206* .. Local Scalars ..
207 COMPLEX*16 TEMP
208 INTEGER I,INFO,J,K,NROWA
209 LOGICAL LSIDE,NOCONJ,NOUNIT,UPPER
210* ..
211* .. Parameters ..
212 COMPLEX*16 ONE
213 parameter(one= (1.0d+0,0.0d+0))
214 COMPLEX*16 ZERO
215 parameter(zero= (0.0d+0,0.0d+0))
216* ..
217*
218* Test the input parameters.
219*
220 lside = lsame(side,'L')
221 IF (lside) THEN
222 nrowa = m
223 ELSE
224 nrowa = n
225 END IF
226 noconj = lsame(transa,'T')
227 nounit = lsame(diag,'N')
228 upper = lsame(uplo,'U')
229*
230 info = 0
231 IF ((.NOT.lside) .AND. (.NOT.lsame(side,'R'))) THEN
232 info = 1
233 ELSE IF ((.NOT.upper) .AND. (.NOT.lsame(uplo,'L'))) THEN
234 info = 2
235 ELSE IF ((.NOT.lsame(transa,'N')) .AND.
236 + (.NOT.lsame(transa,'T')) .AND.
237 + (.NOT.lsame(transa,'C'))) THEN
238 info = 3
239 ELSE IF ((.NOT.lsame(diag,'U')) .AND.
240 + (.NOT.lsame(diag,'N'))) THEN
241 info = 4
242 ELSE IF (m.LT.0) THEN
243 info = 5
244 ELSE IF (n.LT.0) THEN
245 info = 6
246 ELSE IF (lda.LT.max(1,nrowa)) THEN
247 info = 9
248 ELSE IF (ldb.LT.max(1,m)) THEN
249 info = 11
250 END IF
251 IF (info.NE.0) THEN
252 CALL xerbla('ZTRSM ',info)
253 RETURN
254 END IF
255*
256* Quick return if possible.
257*
258 IF (m.EQ.0 .OR. n.EQ.0) RETURN
259*
260* And when alpha.eq.zero.
261*
262 IF (alpha.EQ.zero) THEN
263 DO 20 j = 1,n
264 DO 10 i = 1,m
265 b(i,j) = zero
266 10 CONTINUE
267 20 CONTINUE
268 RETURN
269 END IF
270*
271* Start the operations.
272*
273 IF (lside) THEN
274 IF (lsame(transa,'N')) THEN
275*
276* Form B := alpha*inv( A )*B.
277*
278 IF (upper) THEN
279 DO 60 j = 1,n
280 IF (alpha.NE.one) THEN
281 DO 30 i = 1,m
282 b(i,j) = alpha*b(i,j)
283 30 CONTINUE
284 END IF
285 DO 50 k = m,1,-1
286 IF (b(k,j).NE.zero) THEN
287 IF (nounit) b(k,j) = b(k,j)/a(k,k)
288 DO 40 i = 1,k - 1
289 b(i,j) = b(i,j) - b(k,j)*a(i,k)
290 40 CONTINUE
291 END IF
292 50 CONTINUE
293 60 CONTINUE
294 ELSE
295 DO 100 j = 1,n
296 IF (alpha.NE.one) THEN
297 DO 70 i = 1,m
298 b(i,j) = alpha*b(i,j)
299 70 CONTINUE
300 END IF
301 DO 90 k = 1,m
302 IF (b(k,j).NE.zero) THEN
303 IF (nounit) b(k,j) = b(k,j)/a(k,k)
304 DO 80 i = k + 1,m
305 b(i,j) = b(i,j) - b(k,j)*a(i,k)
306 80 CONTINUE
307 END IF
308 90 CONTINUE
309 100 CONTINUE
310 END IF
311 ELSE
312*
313* Form B := alpha*inv( A**T )*B
314* or B := alpha*inv( A**H )*B.
315*
316 IF (upper) THEN
317 DO 140 j = 1,n
318 DO 130 i = 1,m
319 temp = alpha*b(i,j)
320 IF (noconj) THEN
321 DO 110 k = 1,i - 1
322 temp = temp - a(k,i)*b(k,j)
323 110 CONTINUE
324 IF (nounit) temp = temp/a(i,i)
325 ELSE
326 DO 120 k = 1,i - 1
327 temp = temp - dconjg(a(k,i))*b(k,j)
328 120 CONTINUE
329 IF (nounit) temp = temp/dconjg(a(i,i))
330 END IF
331 b(i,j) = temp
332 130 CONTINUE
333 140 CONTINUE
334 ELSE
335 DO 180 j = 1,n
336 DO 170 i = m,1,-1
337 temp = alpha*b(i,j)
338 IF (noconj) THEN
339 DO 150 k = i + 1,m
340 temp = temp - a(k,i)*b(k,j)
341 150 CONTINUE
342 IF (nounit) temp = temp/a(i,i)
343 ELSE
344 DO 160 k = i + 1,m
345 temp = temp - dconjg(a(k,i))*b(k,j)
346 160 CONTINUE
347 IF (nounit) temp = temp/dconjg(a(i,i))
348 END IF
349 b(i,j) = temp
350 170 CONTINUE
351 180 CONTINUE
352 END IF
353 END IF
354 ELSE
355 IF (lsame(transa,'N')) THEN
356*
357* Form B := alpha*B*inv( A ).
358*
359 IF (upper) THEN
360 DO 230 j = 1,n
361 IF (alpha.NE.one) THEN
362 DO 190 i = 1,m
363 b(i,j) = alpha*b(i,j)
364 190 CONTINUE
365 END IF
366 DO 210 k = 1,j - 1
367 IF (a(k,j).NE.zero) THEN
368 DO 200 i = 1,m
369 b(i,j) = b(i,j) - a(k,j)*b(i,k)
370 200 CONTINUE
371 END IF
372 210 CONTINUE
373 IF (nounit) THEN
374 temp = one/a(j,j)
375 DO 220 i = 1,m
376 b(i,j) = temp*b(i,j)
377 220 CONTINUE
378 END IF
379 230 CONTINUE
380 ELSE
381 DO 280 j = n,1,-1
382 IF (alpha.NE.one) THEN
383 DO 240 i = 1,m
384 b(i,j) = alpha*b(i,j)
385 240 CONTINUE
386 END IF
387 DO 260 k = j + 1,n
388 IF (a(k,j).NE.zero) THEN
389 DO 250 i = 1,m
390 b(i,j) = b(i,j) - a(k,j)*b(i,k)
391 250 CONTINUE
392 END IF
393 260 CONTINUE
394 IF (nounit) THEN
395 temp = one/a(j,j)
396 DO 270 i = 1,m
397 b(i,j) = temp*b(i,j)
398 270 CONTINUE
399 END IF
400 280 CONTINUE
401 END IF
402 ELSE
403*
404* Form B := alpha*B*inv( A**T )
405* or B := alpha*B*inv( A**H ).
406*
407 IF (upper) THEN
408 DO 330 k = n,1,-1
409 IF (nounit) THEN
410 IF (noconj) THEN
411 temp = one/a(k,k)
412 ELSE
413 temp = one/dconjg(a(k,k))
414 END IF
415 DO 290 i = 1,m
416 b(i,k) = temp*b(i,k)
417 290 CONTINUE
418 END IF
419 DO 310 j = 1,k - 1
420 IF (a(j,k).NE.zero) THEN
421 IF (noconj) THEN
422 temp = a(j,k)
423 ELSE
424 temp = dconjg(a(j,k))
425 END IF
426 DO 300 i = 1,m
427 b(i,j) = b(i,j) - temp*b(i,k)
428 300 CONTINUE
429 END IF
430 310 CONTINUE
431 IF (alpha.NE.one) THEN
432 DO 320 i = 1,m
433 b(i,k) = alpha*b(i,k)
434 320 CONTINUE
435 END IF
436 330 CONTINUE
437 ELSE
438 DO 380 k = 1,n
439 IF (nounit) THEN
440 IF (noconj) THEN
441 temp = one/a(k,k)
442 ELSE
443 temp = one/dconjg(a(k,k))
444 END IF
445 DO 340 i = 1,m
446 b(i,k) = temp*b(i,k)
447 340 CONTINUE
448 END IF
449 DO 360 j = k + 1,n
450 IF (a(j,k).NE.zero) THEN
451 IF (noconj) THEN
452 temp = a(j,k)
453 ELSE
454 temp = dconjg(a(j,k))
455 END IF
456 DO 350 i = 1,m
457 b(i,j) = b(i,j) - temp*b(i,k)
458 350 CONTINUE
459 END IF
460 360 CONTINUE
461 IF (alpha.NE.one) THEN
462 DO 370 i = 1,m
463 b(i,k) = alpha*b(i,k)
464 370 CONTINUE
465 END IF
466 380 CONTINUE
467 END IF
468 END IF
469 END IF
470*
471 RETURN
472*
473* End of ZTRSM
474*
475 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine ztrsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
ZTRSM
Definition ztrsm.f:180