LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
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 of 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 of 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.
Date
November 2011
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 183 of file dtrsm.f.

183 *
184 * -- Reference BLAS level3 routine (version 3.4.0) --
185 * -- Reference BLAS is a software package provided by Univ. of Tennessee, --
186 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
187 * November 2011
188 *
189 * .. Scalar Arguments ..
190  DOUBLE PRECISION alpha
191  INTEGER lda,ldb,m,n
192  CHARACTER diag,side,transa,uplo
193 * ..
194 * .. Array Arguments ..
195  DOUBLE PRECISION a(lda,*),b(ldb,*)
196 * ..
197 *
198 * =====================================================================
199 *
200 * .. External Functions ..
201  LOGICAL lsame
202  EXTERNAL lsame
203 * ..
204 * .. External Subroutines ..
205  EXTERNAL xerbla
206 * ..
207 * .. Intrinsic Functions ..
208  INTRINSIC max
209 * ..
210 * .. Local Scalars ..
211  DOUBLE PRECISION temp
212  INTEGER i,info,j,k,nrowa
213  LOGICAL lside,nounit,upper
214 * ..
215 * .. Parameters ..
216  DOUBLE PRECISION one,zero
217  parameter(one=1.0d+0,zero=0.0d+0)
218 * ..
219 *
220 * Test the input parameters.
221 *
222  lside = lsame(side,'L')
223  IF (lside) THEN
224  nrowa = m
225  ELSE
226  nrowa = n
227  END IF
228  nounit = lsame(diag,'N')
229  upper = lsame(uplo,'U')
230 *
231  info = 0
232  IF ((.NOT.lside) .AND. (.NOT.lsame(side,'R'))) THEN
233  info = 1
234  ELSE IF ((.NOT.upper) .AND. (.NOT.lsame(uplo,'L'))) THEN
235  info = 2
236  ELSE IF ((.NOT.lsame(transa,'N')) .AND.
237  + (.NOT.lsame(transa,'T')) .AND.
238  + (.NOT.lsame(transa,'C'))) THEN
239  info = 3
240  ELSE IF ((.NOT.lsame(diag,'U')) .AND. (.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('DTRSM ',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 *
315  IF (upper) THEN
316  DO 130 j = 1,n
317  DO 120 i = 1,m
318  temp = alpha*b(i,j)
319  DO 110 k = 1,i - 1
320  temp = temp - a(k,i)*b(k,j)
321  110 CONTINUE
322  IF (nounit) temp = temp/a(i,i)
323  b(i,j) = temp
324  120 CONTINUE
325  130 CONTINUE
326  ELSE
327  DO 160 j = 1,n
328  DO 150 i = m,1,-1
329  temp = alpha*b(i,j)
330  DO 140 k = i + 1,m
331  temp = temp - a(k,i)*b(k,j)
332  140 CONTINUE
333  IF (nounit) temp = temp/a(i,i)
334  b(i,j) = temp
335  150 CONTINUE
336  160 CONTINUE
337  END IF
338  END IF
339  ELSE
340  IF (lsame(transa,'N')) THEN
341 *
342 * Form B := alpha*B*inv( A ).
343 *
344  IF (upper) THEN
345  DO 210 j = 1,n
346  IF (alpha.NE.one) THEN
347  DO 170 i = 1,m
348  b(i,j) = alpha*b(i,j)
349  170 CONTINUE
350  END IF
351  DO 190 k = 1,j - 1
352  IF (a(k,j).NE.zero) THEN
353  DO 180 i = 1,m
354  b(i,j) = b(i,j) - a(k,j)*b(i,k)
355  180 CONTINUE
356  END IF
357  190 CONTINUE
358  IF (nounit) THEN
359  temp = one/a(j,j)
360  DO 200 i = 1,m
361  b(i,j) = temp*b(i,j)
362  200 CONTINUE
363  END IF
364  210 CONTINUE
365  ELSE
366  DO 260 j = n,1,-1
367  IF (alpha.NE.one) THEN
368  DO 220 i = 1,m
369  b(i,j) = alpha*b(i,j)
370  220 CONTINUE
371  END IF
372  DO 240 k = j + 1,n
373  IF (a(k,j).NE.zero) THEN
374  DO 230 i = 1,m
375  b(i,j) = b(i,j) - a(k,j)*b(i,k)
376  230 CONTINUE
377  END IF
378  240 CONTINUE
379  IF (nounit) THEN
380  temp = one/a(j,j)
381  DO 250 i = 1,m
382  b(i,j) = temp*b(i,j)
383  250 CONTINUE
384  END IF
385  260 CONTINUE
386  END IF
387  ELSE
388 *
389 * Form B := alpha*B*inv( A**T ).
390 *
391  IF (upper) THEN
392  DO 310 k = n,1,-1
393  IF (nounit) THEN
394  temp = one/a(k,k)
395  DO 270 i = 1,m
396  b(i,k) = temp*b(i,k)
397  270 CONTINUE
398  END IF
399  DO 290 j = 1,k - 1
400  IF (a(j,k).NE.zero) THEN
401  temp = a(j,k)
402  DO 280 i = 1,m
403  b(i,j) = b(i,j) - temp*b(i,k)
404  280 CONTINUE
405  END IF
406  290 CONTINUE
407  IF (alpha.NE.one) THEN
408  DO 300 i = 1,m
409  b(i,k) = alpha*b(i,k)
410  300 CONTINUE
411  END IF
412  310 CONTINUE
413  ELSE
414  DO 360 k = 1,n
415  IF (nounit) THEN
416  temp = one/a(k,k)
417  DO 320 i = 1,m
418  b(i,k) = temp*b(i,k)
419  320 CONTINUE
420  END IF
421  DO 340 j = k + 1,n
422  IF (a(j,k).NE.zero) THEN
423  temp = a(j,k)
424  DO 330 i = 1,m
425  b(i,j) = b(i,j) - temp*b(i,k)
426  330 CONTINUE
427  END IF
428  340 CONTINUE
429  IF (alpha.NE.one) THEN
430  DO 350 i = 1,m
431  b(i,k) = alpha*b(i,k)
432  350 CONTINUE
433  END IF
434  360 CONTINUE
435  END IF
436  END IF
437  END IF
438 *
439  RETURN
440 *
441 * End of DTRSM .
442 *
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55

Here is the call graph for this function:

Here is the caller graph for this function: