LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
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 of 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 of 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 *> \date November 2011
160 *
161 *> \ingroup double_blas_level3
162 *
163 *> \par Further Details:
164 * =====================
165 *>
166 *> \verbatim
167 *>
168 *> Level 3 Blas routine.
169 *>
170 *> -- Written on 8-February-1989.
171 *> Jack Dongarra, Argonne National Laboratory.
172 *> Iain Duff, AERE Harwell.
173 *> Jeremy Du Croz, Numerical Algorithms Group Ltd.
174 *> Sven Hammarling, Numerical Algorithms Group Ltd.
175 *> \endverbatim
176 *>
177 * =====================================================================
178  SUBROUTINE dtrmm(SIDE,UPLO,TRANSA,DIAG,M,N,ALPHA,A,LDA,B,LDB)
179 *
180 * -- Reference BLAS level3 routine (version 3.4.0) --
181 * -- Reference BLAS is a software package provided by Univ. of Tennessee, --
182 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
183 * November 2011
184 *
185 * .. Scalar Arguments ..
186  DOUBLE PRECISION ALPHA
187  INTEGER LDA,LDB,M,N
188  CHARACTER DIAG,SIDE,TRANSA,UPLO
189 * ..
190 * .. Array Arguments ..
191  DOUBLE PRECISION 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 max
205 * ..
206 * .. Local Scalars ..
207  DOUBLE PRECISION TEMP
208  INTEGER I,INFO,J,K,NROWA
209  LOGICAL LSIDE,NOUNIT,UPPER
210 * ..
211 * .. Parameters ..
212  DOUBLE PRECISION ONE,ZERO
213  parameter(one=1.0d+0,zero=0.0d+0)
214 * ..
215 *
216 * Test the input parameters.
217 *
218  lside = lsame(side,'L')
219  IF (lside) THEN
220  nrowa = m
221  ELSE
222  nrowa = n
223  END IF
224  nounit = lsame(diag,'N')
225  upper = lsame(uplo,'U')
226 *
227  info = 0
228  IF ((.NOT.lside) .AND. (.NOT.lsame(side,'R'))) THEN
229  info = 1
230  ELSE IF ((.NOT.upper) .AND. (.NOT.lsame(uplo,'L'))) THEN
231  info = 2
232  ELSE IF ((.NOT.lsame(transa,'N')) .AND.
233  + (.NOT.lsame(transa,'T')) .AND.
234  + (.NOT.lsame(transa,'C'))) THEN
235  info = 3
236  ELSE IF ((.NOT.lsame(diag,'U')) .AND. (.NOT.lsame(diag,'N'))) THEN
237  info = 4
238  ELSE IF (m.LT.0) THEN
239  info = 5
240  ELSE IF (n.LT.0) THEN
241  info = 6
242  ELSE IF (lda.LT.max(1,nrowa)) THEN
243  info = 9
244  ELSE IF (ldb.LT.max(1,m)) THEN
245  info = 11
246  END IF
247  IF (info.NE.0) THEN
248  CALL xerbla('DTRMM ',info)
249  RETURN
250  END IF
251 *
252 * Quick return if possible.
253 *
254  IF (m.EQ.0 .OR. n.EQ.0) RETURN
255 *
256 * And when alpha.eq.zero.
257 *
258  IF (alpha.EQ.zero) THEN
259  DO 20 j = 1,n
260  DO 10 i = 1,m
261  b(i,j) = zero
262  10 CONTINUE
263  20 CONTINUE
264  RETURN
265  END IF
266 *
267 * Start the operations.
268 *
269  IF (lside) THEN
270  IF (lsame(transa,'N')) THEN
271 *
272 * Form B := alpha*A*B.
273 *
274  IF (upper) THEN
275  DO 50 j = 1,n
276  DO 40 k = 1,m
277  IF (b(k,j).NE.zero) THEN
278  temp = alpha*b(k,j)
279  DO 30 i = 1,k - 1
280  b(i,j) = b(i,j) + temp*a(i,k)
281  30 CONTINUE
282  IF (nounit) temp = temp*a(k,k)
283  b(k,j) = temp
284  END IF
285  40 CONTINUE
286  50 CONTINUE
287  ELSE
288  DO 80 j = 1,n
289  DO 70 k = m,1,-1
290  IF (b(k,j).NE.zero) THEN
291  temp = alpha*b(k,j)
292  b(k,j) = temp
293  IF (nounit) b(k,j) = b(k,j)*a(k,k)
294  DO 60 i = k + 1,m
295  b(i,j) = b(i,j) + temp*a(i,k)
296  60 CONTINUE
297  END IF
298  70 CONTINUE
299  80 CONTINUE
300  END IF
301  ELSE
302 *
303 * Form B := alpha*A**T*B.
304 *
305  IF (upper) THEN
306  DO 110 j = 1,n
307  DO 100 i = m,1,-1
308  temp = b(i,j)
309  IF (nounit) temp = temp*a(i,i)
310  DO 90 k = 1,i - 1
311  temp = temp + a(k,i)*b(k,j)
312  90 CONTINUE
313  b(i,j) = alpha*temp
314  100 CONTINUE
315  110 CONTINUE
316  ELSE
317  DO 140 j = 1,n
318  DO 130 i = 1,m
319  temp = b(i,j)
320  IF (nounit) temp = temp*a(i,i)
321  DO 120 k = i + 1,m
322  temp = temp + a(k,i)*b(k,j)
323  120 CONTINUE
324  b(i,j) = alpha*temp
325  130 CONTINUE
326  140 CONTINUE
327  END IF
328  END IF
329  ELSE
330  IF (lsame(transa,'N')) THEN
331 *
332 * Form B := alpha*B*A.
333 *
334  IF (upper) THEN
335  DO 180 j = n,1,-1
336  temp = alpha
337  IF (nounit) temp = temp*a(j,j)
338  DO 150 i = 1,m
339  b(i,j) = temp*b(i,j)
340  150 CONTINUE
341  DO 170 k = 1,j - 1
342  IF (a(k,j).NE.zero) THEN
343  temp = alpha*a(k,j)
344  DO 160 i = 1,m
345  b(i,j) = b(i,j) + temp*b(i,k)
346  160 CONTINUE
347  END IF
348  170 CONTINUE
349  180 CONTINUE
350  ELSE
351  DO 220 j = 1,n
352  temp = alpha
353  IF (nounit) temp = temp*a(j,j)
354  DO 190 i = 1,m
355  b(i,j) = temp*b(i,j)
356  190 CONTINUE
357  DO 210 k = j + 1,n
358  IF (a(k,j).NE.zero) THEN
359  temp = alpha*a(k,j)
360  DO 200 i = 1,m
361  b(i,j) = b(i,j) + temp*b(i,k)
362  200 CONTINUE
363  END IF
364  210 CONTINUE
365  220 CONTINUE
366  END IF
367  ELSE
368 *
369 * Form B := alpha*B*A**T.
370 *
371  IF (upper) THEN
372  DO 260 k = 1,n
373  DO 240 j = 1,k - 1
374  IF (a(j,k).NE.zero) THEN
375  temp = alpha*a(j,k)
376  DO 230 i = 1,m
377  b(i,j) = b(i,j) + temp*b(i,k)
378  230 CONTINUE
379  END IF
380  240 CONTINUE
381  temp = alpha
382  IF (nounit) temp = temp*a(k,k)
383  IF (temp.NE.one) THEN
384  DO 250 i = 1,m
385  b(i,k) = temp*b(i,k)
386  250 CONTINUE
387  END IF
388  260 CONTINUE
389  ELSE
390  DO 300 k = n,1,-1
391  DO 280 j = k + 1,n
392  IF (a(j,k).NE.zero) THEN
393  temp = alpha*a(j,k)
394  DO 270 i = 1,m
395  b(i,j) = b(i,j) + temp*b(i,k)
396  270 CONTINUE
397  END IF
398  280 CONTINUE
399  temp = alpha
400  IF (nounit) temp = temp*a(k,k)
401  IF (temp.NE.one) THEN
402  DO 290 i = 1,m
403  b(i,k) = temp*b(i,k)
404  290 CONTINUE
405  END IF
406  300 CONTINUE
407  END IF
408  END IF
409  END IF
410 *
411  RETURN
412 *
413 * End of DTRMM .
414 *
415  END
subroutine dtrmm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
DTRMM
Definition: dtrmm.f:179
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62