LAPACK  3.10.0
LAPACK: Linear Algebra PACKage

◆ dtrsm()

subroutine dtrsm ( character  SIDE,
character  UPLO,
character  TRANSA,
character  DIAG,
integer  M,
integer  N,
double precision  ALPHA,
double precision, dimension(lda,*)  A,
integer  LDA,
double precision, dimension(ldb,*)  B,
integer  LDB 
)

DTRSM

Purpose:
 DTRSM  solves one of the matrix equations

    op( A )*X = alpha*B,   or   X*op( A ) = alpha*B,

 where alpha is a scalar, X and B are m by n matrices, A is a unit, or
 non-unit,  upper or lower triangular matrix  and  op( A )  is one  of

    op( A ) = A   or   op( A ) = A**T.

 The matrix X is overwritten on B.
Parameters
[in]SIDE
          SIDE is CHARACTER*1
           On entry, SIDE specifies whether op( A ) appears on the left
           or right of X as follows:

              SIDE = 'L' or 'l'   op( A )*X = alpha*B.

              SIDE = 'R' or 'r'   X*op( A ) = alpha*B.
[in]UPLO
          UPLO is CHARACTER*1
           On entry, UPLO specifies whether the matrix A is an upper or
           lower triangular matrix as follows:

              UPLO = 'U' or 'u'   A is an upper triangular matrix.

              UPLO = 'L' or 'l'   A is a lower triangular matrix.
[in]TRANSA
          TRANSA is CHARACTER*1
           On entry, TRANSA specifies the form of op( A ) to be used in
           the matrix multiplication as follows:

              TRANSA = 'N' or 'n'   op( A ) = A.

              TRANSA = 'T' or 't'   op( A ) = A**T.

              TRANSA = 'C' or 'c'   op( A ) = A**T.
[in]DIAG
          DIAG is CHARACTER*1
           On entry, DIAG specifies whether or not A is unit triangular
           as follows:

              DIAG = 'U' or 'u'   A is assumed to be unit triangular.

              DIAG = 'N' or 'n'   A is not assumed to be unit
                                  triangular.
[in]M
          M is INTEGER
           On entry, M specifies the number of rows of B. M must be at
           least zero.
[in]N
          N is INTEGER
           On entry, N specifies the number of columns of B.  N must be
           at least zero.
[in]ALPHA
          ALPHA is DOUBLE PRECISION.
           On entry,  ALPHA specifies the scalar  alpha. When  alpha is
           zero then  A is not referenced and  B need not be set before
           entry.
[in]A
          A is DOUBLE PRECISION array, dimension ( LDA, k ),
           where k is m when SIDE = 'L' or 'l'
             and k is n when SIDE = 'R' or 'r'.
           Before entry  with  UPLO = 'U' or 'u',  the  leading  k by k
           upper triangular part of the array  A must contain the upper
           triangular matrix  and the strictly lower triangular part of
           A is not referenced.
           Before entry  with  UPLO = 'L' or 'l',  the  leading  k by k
           lower triangular part of the array  A must contain the lower
           triangular matrix  and the strictly upper triangular part of
           A is not referenced.
           Note that when  DIAG = 'U' or 'u',  the diagonal elements of
           A  are not referenced either,  but are assumed to be  unity.
[in]LDA
          LDA is INTEGER
           On entry, LDA specifies the first dimension of A as declared
           in the calling (sub) program.  When  SIDE = 'L' or 'l'  then
           LDA  must be at least  max( 1, m ),  when  SIDE = 'R' or 'r'
           then LDA must be at least max( 1, n ).
[in,out]B
          B is DOUBLE PRECISION array, dimension ( LDB, N )
           Before entry,  the leading  m by n part of the array  B must
           contain  the  right-hand  side  matrix  B,  and  on exit  is
           overwritten by the solution matrix  X.
[in]LDB
          LDB is INTEGER
           On entry, LDB specifies the first dimension of B as declared
           in  the  calling  (sub)  program.   LDB  must  be  at  least
           max( 1, m ).
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Further Details:
  Level 3 Blas routine.


  -- Written on 8-February-1989.
     Jack Dongarra, Argonne National Laboratory.
     Iain Duff, AERE Harwell.
     Jeremy Du Croz, Numerical Algorithms Group Ltd.
     Sven Hammarling, Numerical Algorithms Group Ltd.

Definition at line 180 of file dtrsm.f.

181 *
182 * -- Reference BLAS level3 routine --
183 * -- Reference BLAS is a software package provided by Univ. of Tennessee, --
184 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
185 *
186 * .. Scalar Arguments ..
187  DOUBLE PRECISION ALPHA
188  INTEGER LDA,LDB,M,N
189  CHARACTER DIAG,SIDE,TRANSA,UPLO
190 * ..
191 * .. Array Arguments ..
192  DOUBLE PRECISION A(LDA,*),B(LDB,*)
193 * ..
194 *
195 * =====================================================================
196 *
197 * .. External Functions ..
198  LOGICAL LSAME
199  EXTERNAL lsame
200 * ..
201 * .. External Subroutines ..
202  EXTERNAL xerbla
203 * ..
204 * .. Intrinsic Functions ..
205  INTRINSIC max
206 * ..
207 * .. Local Scalars ..
208  DOUBLE PRECISION TEMP
209  INTEGER I,INFO,J,K,NROWA
210  LOGICAL LSIDE,NOUNIT,UPPER
211 * ..
212 * .. Parameters ..
213  DOUBLE PRECISION ONE,ZERO
214  parameter(one=1.0d+0,zero=0.0d+0)
215 * ..
216 *
217 * Test the input parameters.
218 *
219  lside = lsame(side,'L')
220  IF (lside) THEN
221  nrowa = m
222  ELSE
223  nrowa = n
224  END IF
225  nounit = lsame(diag,'N')
226  upper = lsame(uplo,'U')
227 *
228  info = 0
229  IF ((.NOT.lside) .AND. (.NOT.lsame(side,'R'))) THEN
230  info = 1
231  ELSE IF ((.NOT.upper) .AND. (.NOT.lsame(uplo,'L'))) THEN
232  info = 2
233  ELSE IF ((.NOT.lsame(transa,'N')) .AND.
234  + (.NOT.lsame(transa,'T')) .AND.
235  + (.NOT.lsame(transa,'C'))) THEN
236  info = 3
237  ELSE IF ((.NOT.lsame(diag,'U')) .AND. (.NOT.lsame(diag,'N'))) THEN
238  info = 4
239  ELSE IF (m.LT.0) THEN
240  info = 5
241  ELSE IF (n.LT.0) THEN
242  info = 6
243  ELSE IF (lda.LT.max(1,nrowa)) THEN
244  info = 9
245  ELSE IF (ldb.LT.max(1,m)) THEN
246  info = 11
247  END IF
248  IF (info.NE.0) THEN
249  CALL xerbla('DTRSM ',info)
250  RETURN
251  END IF
252 *
253 * Quick return if possible.
254 *
255  IF (m.EQ.0 .OR. n.EQ.0) RETURN
256 *
257 * And when alpha.eq.zero.
258 *
259  IF (alpha.EQ.zero) THEN
260  DO 20 j = 1,n
261  DO 10 i = 1,m
262  b(i,j) = zero
263  10 CONTINUE
264  20 CONTINUE
265  RETURN
266  END IF
267 *
268 * Start the operations.
269 *
270  IF (lside) THEN
271  IF (lsame(transa,'N')) THEN
272 *
273 * Form B := alpha*inv( A )*B.
274 *
275  IF (upper) THEN
276  DO 60 j = 1,n
277  IF (alpha.NE.one) THEN
278  DO 30 i = 1,m
279  b(i,j) = alpha*b(i,j)
280  30 CONTINUE
281  END IF
282  DO 50 k = m,1,-1
283  IF (b(k,j).NE.zero) THEN
284  IF (nounit) b(k,j) = b(k,j)/a(k,k)
285  DO 40 i = 1,k - 1
286  b(i,j) = b(i,j) - b(k,j)*a(i,k)
287  40 CONTINUE
288  END IF
289  50 CONTINUE
290  60 CONTINUE
291  ELSE
292  DO 100 j = 1,n
293  IF (alpha.NE.one) THEN
294  DO 70 i = 1,m
295  b(i,j) = alpha*b(i,j)
296  70 CONTINUE
297  END IF
298  DO 90 k = 1,m
299  IF (b(k,j).NE.zero) THEN
300  IF (nounit) b(k,j) = b(k,j)/a(k,k)
301  DO 80 i = k + 1,m
302  b(i,j) = b(i,j) - b(k,j)*a(i,k)
303  80 CONTINUE
304  END IF
305  90 CONTINUE
306  100 CONTINUE
307  END IF
308  ELSE
309 *
310 * Form B := alpha*inv( A**T )*B.
311 *
312  IF (upper) THEN
313  DO 130 j = 1,n
314  DO 120 i = 1,m
315  temp = alpha*b(i,j)
316  DO 110 k = 1,i - 1
317  temp = temp - a(k,i)*b(k,j)
318  110 CONTINUE
319  IF (nounit) temp = temp/a(i,i)
320  b(i,j) = temp
321  120 CONTINUE
322  130 CONTINUE
323  ELSE
324  DO 160 j = 1,n
325  DO 150 i = m,1,-1
326  temp = alpha*b(i,j)
327  DO 140 k = i + 1,m
328  temp = temp - a(k,i)*b(k,j)
329  140 CONTINUE
330  IF (nounit) temp = temp/a(i,i)
331  b(i,j) = temp
332  150 CONTINUE
333  160 CONTINUE
334  END IF
335  END IF
336  ELSE
337  IF (lsame(transa,'N')) THEN
338 *
339 * Form B := alpha*B*inv( A ).
340 *
341  IF (upper) THEN
342  DO 210 j = 1,n
343  IF (alpha.NE.one) THEN
344  DO 170 i = 1,m
345  b(i,j) = alpha*b(i,j)
346  170 CONTINUE
347  END IF
348  DO 190 k = 1,j - 1
349  IF (a(k,j).NE.zero) THEN
350  DO 180 i = 1,m
351  b(i,j) = b(i,j) - a(k,j)*b(i,k)
352  180 CONTINUE
353  END IF
354  190 CONTINUE
355  IF (nounit) THEN
356  temp = one/a(j,j)
357  DO 200 i = 1,m
358  b(i,j) = temp*b(i,j)
359  200 CONTINUE
360  END IF
361  210 CONTINUE
362  ELSE
363  DO 260 j = n,1,-1
364  IF (alpha.NE.one) THEN
365  DO 220 i = 1,m
366  b(i,j) = alpha*b(i,j)
367  220 CONTINUE
368  END IF
369  DO 240 k = j + 1,n
370  IF (a(k,j).NE.zero) THEN
371  DO 230 i = 1,m
372  b(i,j) = b(i,j) - a(k,j)*b(i,k)
373  230 CONTINUE
374  END IF
375  240 CONTINUE
376  IF (nounit) THEN
377  temp = one/a(j,j)
378  DO 250 i = 1,m
379  b(i,j) = temp*b(i,j)
380  250 CONTINUE
381  END IF
382  260 CONTINUE
383  END IF
384  ELSE
385 *
386 * Form B := alpha*B*inv( A**T ).
387 *
388  IF (upper) THEN
389  DO 310 k = n,1,-1
390  IF (nounit) THEN
391  temp = one/a(k,k)
392  DO 270 i = 1,m
393  b(i,k) = temp*b(i,k)
394  270 CONTINUE
395  END IF
396  DO 290 j = 1,k - 1
397  IF (a(j,k).NE.zero) THEN
398  temp = a(j,k)
399  DO 280 i = 1,m
400  b(i,j) = b(i,j) - temp*b(i,k)
401  280 CONTINUE
402  END IF
403  290 CONTINUE
404  IF (alpha.NE.one) THEN
405  DO 300 i = 1,m
406  b(i,k) = alpha*b(i,k)
407  300 CONTINUE
408  END IF
409  310 CONTINUE
410  ELSE
411  DO 360 k = 1,n
412  IF (nounit) THEN
413  temp = one/a(k,k)
414  DO 320 i = 1,m
415  b(i,k) = temp*b(i,k)
416  320 CONTINUE
417  END IF
418  DO 340 j = k + 1,n
419  IF (a(j,k).NE.zero) THEN
420  temp = a(j,k)
421  DO 330 i = 1,m
422  b(i,j) = b(i,j) - temp*b(i,k)
423  330 CONTINUE
424  END IF
425  340 CONTINUE
426  IF (alpha.NE.one) THEN
427  DO 350 i = 1,m
428  b(i,k) = alpha*b(i,k)
429  350 CONTINUE
430  END IF
431  360 CONTINUE
432  END IF
433  END IF
434  END IF
435 *
436  RETURN
437 *
438 * End of DTRSM
439 *
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:60
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:53
Here is the call graph for this function:
Here is the caller graph for this function: