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