LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine sdrvpo ( logical, dimension( * )  DOTYPE,
integer  NN,
integer, dimension( * )  NVAL,
integer  NRHS,
real  THRESH,
logical  TSTERR,
integer  NMAX,
real, dimension( * )  A,
real, dimension( * )  AFAC,
real, dimension( * )  ASAV,
real, dimension( * )  B,
real, dimension( * )  BSAV,
real, dimension( * )  X,
real, dimension( * )  XACT,
real, dimension( * )  S,
real, dimension( * )  WORK,
real, dimension( * )  RWORK,
integer, dimension( * )  IWORK,
integer  NOUT 
)

SDRVPO

SDRVPOX

Purpose:
 SDRVPO tests the driver routines SPOSV and -SVX.
Parameters
[in]DOTYPE
          DOTYPE is LOGICAL array, dimension (NTYPES)
          The matrix types to be used for testing.  Matrices of type j
          (for 1 <= j <= NTYPES) are used for testing if DOTYPE(j) =
          .TRUE.; if DOTYPE(j) = .FALSE., then type j is not used.
[in]NN
          NN is INTEGER
          The number of values of N contained in the vector NVAL.
[in]NVAL
          NVAL is INTEGER array, dimension (NN)
          The values of the matrix dimension N.
[in]NRHS
          NRHS is INTEGER
          The number of right hand side vectors to be generated for
          each linear system.
[in]THRESH
          THRESH is REAL
          The threshold value for the test ratios.  A result is
          included in the output file if RESULT >= THRESH.  To have
          every test ratio printed, use THRESH = 0.
[in]TSTERR
          TSTERR is LOGICAL
          Flag that indicates whether error exits are to be tested.
[in]NMAX
          NMAX is INTEGER
          The maximum value permitted for N, used in dimensioning the
          work arrays.
[out]A
          A is REAL array, dimension (NMAX*NMAX)
[out]AFAC
          AFAC is REAL array, dimension (NMAX*NMAX)
[out]ASAV
          ASAV is REAL array, dimension (NMAX*NMAX)
[out]B
          B is REAL array, dimension (NMAX*NRHS)
[out]BSAV
          BSAV is REAL array, dimension (NMAX*NRHS)
[out]X
          X is REAL array, dimension (NMAX*NRHS)
[out]XACT
          XACT is REAL array, dimension (NMAX*NRHS)
[out]S
          S is REAL array, dimension (NMAX)
[out]WORK
          WORK is REAL array, dimension
                      (NMAX*max(3,NRHS))
[out]RWORK
          RWORK is REAL array, dimension (NMAX+2*NRHS)
[out]IWORK
          IWORK is INTEGER array, dimension (NMAX)
[in]NOUT
          NOUT is INTEGER
          The unit number for output.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011
Purpose:
 SDRVPO tests the driver routines SPOSV, -SVX, and -SVXX.

 Note that this file is used only when the XBLAS are available,
 otherwise sdrvpo.f defines this subroutine.
Parameters
[in]DOTYPE
          DOTYPE is LOGICAL array, dimension (NTYPES)
          The matrix types to be used for testing.  Matrices of type j
          (for 1 <= j <= NTYPES) are used for testing if DOTYPE(j) =
          .TRUE.; if DOTYPE(j) = .FALSE., then type j is not used.
[in]NN
          NN is INTEGER
          The number of values of N contained in the vector NVAL.
[in]NVAL
          NVAL is INTEGER array, dimension (NN)
          The values of the matrix dimension N.
[in]NRHS
          NRHS is INTEGER
          The number of right hand side vectors to be generated for
          each linear system.
[in]THRESH
          THRESH is REAL
          The threshold value for the test ratios.  A result is
          included in the output file if RESULT >= THRESH.  To have
          every test ratio printed, use THRESH = 0.
[in]TSTERR
          TSTERR is LOGICAL
          Flag that indicates whether error exits are to be tested.
[in]NMAX
          NMAX is INTEGER
          The maximum value permitted for N, used in dimensioning the
          work arrays.
[out]A
          A is REAL array, dimension (NMAX*NMAX)
[out]AFAC
          AFAC is REAL array, dimension (NMAX*NMAX)
[out]ASAV
          ASAV is REAL array, dimension (NMAX*NMAX)
[out]B
          B is REAL array, dimension (NMAX*NRHS)
[out]BSAV
          BSAV is REAL array, dimension (NMAX*NRHS)
[out]X
          X is REAL array, dimension (NMAX*NRHS)
[out]XACT
          XACT is REAL array, dimension (NMAX*NRHS)
[out]S
          S is REAL array, dimension (NMAX)
[out]WORK
          WORK is REAL array, dimension
                      (NMAX*max(3,NRHS))
[out]RWORK
          RWORK is REAL array, dimension (NMAX+2*NRHS)
[out]IWORK
          IWORK is INTEGER array, dimension (NMAX)
[in]NOUT
          NOUT is INTEGER
          The unit number for output.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2013

Definition at line 166 of file sdrvpo.f.

166 *
167 * -- LAPACK test routine (version 3.4.0) --
168 * -- LAPACK is a software package provided by Univ. of Tennessee, --
169 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
170 * November 2011
171 *
172 * .. Scalar Arguments ..
173  LOGICAL tsterr
174  INTEGER nmax, nn, nout, nrhs
175  REAL thresh
176 * ..
177 * .. Array Arguments ..
178  LOGICAL dotype( * )
179  INTEGER iwork( * ), nval( * )
180  REAL a( * ), afac( * ), asav( * ), b( * ),
181  $ bsav( * ), rwork( * ), s( * ), work( * ),
182  $ x( * ), xact( * )
183 * ..
184 *
185 * =====================================================================
186 *
187 * .. Parameters ..
188  REAL one, zero
189  parameter ( one = 1.0e+0, zero = 0.0e+0 )
190  INTEGER ntypes
191  parameter ( ntypes = 9 )
192  INTEGER ntests
193  parameter ( ntests = 6 )
194 * ..
195 * .. Local Scalars ..
196  LOGICAL equil, nofact, prefac, zerot
197  CHARACTER dist, equed, fact, TYPE, uplo, xtype
198  CHARACTER*3 path
199  INTEGER i, iequed, ifact, imat, in, info, ioff, iuplo,
200  $ izero, k, k1, kl, ku, lda, mode, n, nb, nbmin,
201  $ nerrs, nfact, nfail, nimat, nrun, nt
202  REAL ainvnm, amax, anorm, cndnum, rcond, rcondc,
203  $ roldc, scond
204 * ..
205 * .. Local Arrays ..
206  CHARACTER equeds( 2 ), facts( 3 ), uplos( 2 )
207  INTEGER iseed( 4 ), iseedy( 4 )
208  REAL result( ntests )
209 * ..
210 * .. External Functions ..
211  LOGICAL lsame
212  REAL sget06, slansy
213  EXTERNAL lsame, sget06, slansy
214 * ..
215 * .. External Subroutines ..
216  EXTERNAL aladhd, alaerh, alasvm, serrvx, sget04, slacpy,
219  $ spotri, xlaenv
220 * ..
221 * .. Intrinsic Functions ..
222  INTRINSIC max
223 * ..
224 * .. Scalars in Common ..
225  LOGICAL lerr, ok
226  CHARACTER*32 srnamt
227  INTEGER infot, nunit
228 * ..
229 * .. Common blocks ..
230  COMMON / infoc / infot, nunit, ok, lerr
231  COMMON / srnamc / srnamt
232 * ..
233 * .. Data statements ..
234  DATA iseedy / 1988, 1989, 1990, 1991 /
235  DATA uplos / 'U', 'L' /
236  DATA facts / 'F', 'N', 'E' /
237  DATA equeds / 'N', 'Y' /
238 * ..
239 * .. Executable Statements ..
240 *
241 * Initialize constants and the random number seed.
242 *
243  path( 1: 1 ) = 'Single precision'
244  path( 2: 3 ) = 'PO'
245  nrun = 0
246  nfail = 0
247  nerrs = 0
248  DO 10 i = 1, 4
249  iseed( i ) = iseedy( i )
250  10 CONTINUE
251 *
252 * Test the error exits
253 *
254  IF( tsterr )
255  $ CALL serrvx( path, nout )
256  infot = 0
257 *
258 * Set the block size and minimum block size for testing.
259 *
260  nb = 1
261  nbmin = 2
262  CALL xlaenv( 1, nb )
263  CALL xlaenv( 2, nbmin )
264 *
265 * Do for each value of N in NVAL
266 *
267  DO 130 in = 1, nn
268  n = nval( in )
269  lda = max( n, 1 )
270  xtype = 'N'
271  nimat = ntypes
272  IF( n.LE.0 )
273  $ nimat = 1
274 *
275  DO 120 imat = 1, nimat
276 *
277 * Do the tests only if DOTYPE( IMAT ) is true.
278 *
279  IF( .NOT.dotype( imat ) )
280  $ GO TO 120
281 *
282 * Skip types 3, 4, or 5 if the matrix size is too small.
283 *
284  zerot = imat.GE.3 .AND. imat.LE.5
285  IF( zerot .AND. n.LT.imat-2 )
286  $ GO TO 120
287 *
288 * Do first for UPLO = 'U', then for UPLO = 'L'
289 *
290  DO 110 iuplo = 1, 2
291  uplo = uplos( iuplo )
292 *
293 * Set up parameters with SLATB4 and generate a test matrix
294 * with SLATMS.
295 *
296  CALL slatb4( path, imat, n, n, TYPE, kl, ku, anorm, mode,
297  $ cndnum, dist )
298 *
299  srnamt = 'SLATMS'
300  CALL slatms( n, n, dist, iseed, TYPE, rwork, mode,
301  $ cndnum, anorm, kl, ku, uplo, a, lda, work,
302  $ info )
303 *
304 * Check error code from SLATMS.
305 *
306  IF( info.NE.0 ) THEN
307  CALL alaerh( path, 'SLATMS', info, 0, uplo, n, n, -1,
308  $ -1, -1, imat, nfail, nerrs, nout )
309  GO TO 110
310  END IF
311 *
312 * For types 3-5, zero one row and column of the matrix to
313 * test that INFO is returned correctly.
314 *
315  IF( zerot ) THEN
316  IF( imat.EQ.3 ) THEN
317  izero = 1
318  ELSE IF( imat.EQ.4 ) THEN
319  izero = n
320  ELSE
321  izero = n / 2 + 1
322  END IF
323  ioff = ( izero-1 )*lda
324 *
325 * Set row and column IZERO of A to 0.
326 *
327  IF( iuplo.EQ.1 ) THEN
328  DO 20 i = 1, izero - 1
329  a( ioff+i ) = zero
330  20 CONTINUE
331  ioff = ioff + izero
332  DO 30 i = izero, n
333  a( ioff ) = zero
334  ioff = ioff + lda
335  30 CONTINUE
336  ELSE
337  ioff = izero
338  DO 40 i = 1, izero - 1
339  a( ioff ) = zero
340  ioff = ioff + lda
341  40 CONTINUE
342  ioff = ioff - izero
343  DO 50 i = izero, n
344  a( ioff+i ) = zero
345  50 CONTINUE
346  END IF
347  ELSE
348  izero = 0
349  END IF
350 *
351 * Save a copy of the matrix A in ASAV.
352 *
353  CALL slacpy( uplo, n, n, a, lda, asav, lda )
354 *
355  DO 100 iequed = 1, 2
356  equed = equeds( iequed )
357  IF( iequed.EQ.1 ) THEN
358  nfact = 3
359  ELSE
360  nfact = 1
361  END IF
362 *
363  DO 90 ifact = 1, nfact
364  fact = facts( ifact )
365  prefac = lsame( fact, 'F' )
366  nofact = lsame( fact, 'N' )
367  equil = lsame( fact, 'E' )
368 *
369  IF( zerot ) THEN
370  IF( prefac )
371  $ GO TO 90
372  rcondc = zero
373 *
374  ELSE IF( .NOT.lsame( fact, 'N' ) ) THEN
375 *
376 * Compute the condition number for comparison with
377 * the value returned by SPOSVX (FACT = 'N' reuses
378 * the condition number from the previous iteration
379 * with FACT = 'F').
380 *
381  CALL slacpy( uplo, n, n, asav, lda, afac, lda )
382  IF( equil .OR. iequed.GT.1 ) THEN
383 *
384 * Compute row and column scale factors to
385 * equilibrate the matrix A.
386 *
387  CALL spoequ( n, afac, lda, s, scond, amax,
388  $ info )
389  IF( info.EQ.0 .AND. n.GT.0 ) THEN
390  IF( iequed.GT.1 )
391  $ scond = zero
392 *
393 * Equilibrate the matrix.
394 *
395  CALL slaqsy( uplo, n, afac, lda, s, scond,
396  $ amax, equed )
397  END IF
398  END IF
399 *
400 * Save the condition number of the
401 * non-equilibrated system for use in SGET04.
402 *
403  IF( equil )
404  $ roldc = rcondc
405 *
406 * Compute the 1-norm of A.
407 *
408  anorm = slansy( '1', uplo, n, afac, lda, rwork )
409 *
410 * Factor the matrix A.
411 *
412  CALL spotrf( uplo, n, afac, lda, info )
413 *
414 * Form the inverse of A.
415 *
416  CALL slacpy( uplo, n, n, afac, lda, a, lda )
417  CALL spotri( uplo, n, a, lda, info )
418 *
419 * Compute the 1-norm condition number of A.
420 *
421  ainvnm = slansy( '1', uplo, n, a, lda, rwork )
422  IF( anorm.LE.zero .OR. ainvnm.LE.zero ) THEN
423  rcondc = one
424  ELSE
425  rcondc = ( one / anorm ) / ainvnm
426  END IF
427  END IF
428 *
429 * Restore the matrix A.
430 *
431  CALL slacpy( uplo, n, n, asav, lda, a, lda )
432 *
433 * Form an exact solution and set the right hand side.
434 *
435  srnamt = 'SLARHS'
436  CALL slarhs( path, xtype, uplo, ' ', n, n, kl, ku,
437  $ nrhs, a, lda, xact, lda, b, lda,
438  $ iseed, info )
439  xtype = 'C'
440  CALL slacpy( 'Full', n, nrhs, b, lda, bsav, lda )
441 *
442  IF( nofact ) THEN
443 *
444 * --- Test SPOSV ---
445 *
446 * Compute the L*L' or U'*U factorization of the
447 * matrix and solve the system.
448 *
449  CALL slacpy( uplo, n, n, a, lda, afac, lda )
450  CALL slacpy( 'Full', n, nrhs, b, lda, x, lda )
451 *
452  srnamt = 'SPOSV '
453  CALL sposv( uplo, n, nrhs, afac, lda, x, lda,
454  $ info )
455 *
456 * Check error code from SPOSV .
457 *
458  IF( info.NE.izero ) THEN
459  CALL alaerh( path, 'SPOSV ', info, izero,
460  $ uplo, n, n, -1, -1, nrhs, imat,
461  $ nfail, nerrs, nout )
462  GO TO 70
463  ELSE IF( info.NE.0 ) THEN
464  GO TO 70
465  END IF
466 *
467 * Reconstruct matrix from factors and compute
468 * residual.
469 *
470  CALL spot01( uplo, n, a, lda, afac, lda, rwork,
471  $ result( 1 ) )
472 *
473 * Compute residual of the computed solution.
474 *
475  CALL slacpy( 'Full', n, nrhs, b, lda, work,
476  $ lda )
477  CALL spot02( uplo, n, nrhs, a, lda, x, lda,
478  $ work, lda, rwork, result( 2 ) )
479 *
480 * Check solution from generated exact solution.
481 *
482  CALL sget04( n, nrhs, x, lda, xact, lda, rcondc,
483  $ result( 3 ) )
484  nt = 3
485 *
486 * Print information about the tests that did not
487 * pass the threshold.
488 *
489  DO 60 k = 1, nt
490  IF( result( k ).GE.thresh ) THEN
491  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
492  $ CALL aladhd( nout, path )
493  WRITE( nout, fmt = 9999 )'SPOSV ', uplo,
494  $ n, imat, k, result( k )
495  nfail = nfail + 1
496  END IF
497  60 CONTINUE
498  nrun = nrun + nt
499  70 CONTINUE
500  END IF
501 *
502 * --- Test SPOSVX ---
503 *
504  IF( .NOT.prefac )
505  $ CALL slaset( uplo, n, n, zero, zero, afac, lda )
506  CALL slaset( 'Full', n, nrhs, zero, zero, x, lda )
507  IF( iequed.GT.1 .AND. n.GT.0 ) THEN
508 *
509 * Equilibrate the matrix if FACT='F' and
510 * EQUED='Y'.
511 *
512  CALL slaqsy( uplo, n, a, lda, s, scond, amax,
513  $ equed )
514  END IF
515 *
516 * Solve the system and compute the condition number
517 * and error bounds using SPOSVX.
518 *
519  srnamt = 'SPOSVX'
520  CALL sposvx( fact, uplo, n, nrhs, a, lda, afac,
521  $ lda, equed, s, b, lda, x, lda, rcond,
522  $ rwork, rwork( nrhs+1 ), work, iwork,
523  $ info )
524 *
525 * Check the error code from SPOSVX.
526 *
527  IF( info.NE.izero ) THEN
528  CALL alaerh( path, 'SPOSVX', info, izero,
529  $ fact // uplo, n, n, -1, -1, nrhs,
530  $ imat, nfail, nerrs, nout )
531  GO TO 90
532  END IF
533 *
534  IF( info.EQ.0 ) THEN
535  IF( .NOT.prefac ) THEN
536 *
537 * Reconstruct matrix from factors and compute
538 * residual.
539 *
540  CALL spot01( uplo, n, a, lda, afac, lda,
541  $ rwork( 2*nrhs+1 ), result( 1 ) )
542  k1 = 1
543  ELSE
544  k1 = 2
545  END IF
546 *
547 * Compute residual of the computed solution.
548 *
549  CALL slacpy( 'Full', n, nrhs, bsav, lda, work,
550  $ lda )
551  CALL spot02( uplo, n, nrhs, asav, lda, x, lda,
552  $ work, lda, rwork( 2*nrhs+1 ),
553  $ result( 2 ) )
554 *
555 * Check solution from generated exact solution.
556 *
557  IF( nofact .OR. ( prefac .AND. lsame( equed,
558  $ 'N' ) ) ) THEN
559  CALL sget04( n, nrhs, x, lda, xact, lda,
560  $ rcondc, result( 3 ) )
561  ELSE
562  CALL sget04( n, nrhs, x, lda, xact, lda,
563  $ roldc, result( 3 ) )
564  END IF
565 *
566 * Check the error bounds from iterative
567 * refinement.
568 *
569  CALL spot05( uplo, n, nrhs, asav, lda, b, lda,
570  $ x, lda, xact, lda, rwork,
571  $ rwork( nrhs+1 ), result( 4 ) )
572  ELSE
573  k1 = 6
574  END IF
575 *
576 * Compare RCOND from SPOSVX with the computed value
577 * in RCONDC.
578 *
579  result( 6 ) = sget06( rcond, rcondc )
580 *
581 * Print information about the tests that did not pass
582 * the threshold.
583 *
584  DO 80 k = k1, 6
585  IF( result( k ).GE.thresh ) THEN
586  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
587  $ CALL aladhd( nout, path )
588  IF( prefac ) THEN
589  WRITE( nout, fmt = 9997 )'SPOSVX', fact,
590  $ uplo, n, equed, imat, k, result( k )
591  ELSE
592  WRITE( nout, fmt = 9998 )'SPOSVX', fact,
593  $ uplo, n, imat, k, result( k )
594  END IF
595  nfail = nfail + 1
596  END IF
597  80 CONTINUE
598  nrun = nrun + 7 - k1
599  90 CONTINUE
600  100 CONTINUE
601  110 CONTINUE
602  120 CONTINUE
603  130 CONTINUE
604 *
605 * Print a summary of the results.
606 *
607  CALL alasvm( path, nout, nfail, nrun, nerrs )
608 *
609  9999 FORMAT( 1x, a, ', UPLO=''', a1, ''', N =', i5, ', type ', i1,
610  $ ', test(', i1, ')=', g12.5 )
611  9998 FORMAT( 1x, a, ', FACT=''', a1, ''', UPLO=''', a1, ''', N=', i5,
612  $ ', type ', i1, ', test(', i1, ')=', g12.5 )
613  9997 FORMAT( 1x, a, ', FACT=''', a1, ''', UPLO=''', a1, ''', N=', i5,
614  $ ', EQUED=''', a1, ''', type ', i1, ', test(', i1, ') =',
615  $ g12.5 )
616  RETURN
617 *
618 * End of SDRVPO
619 *
subroutine alasvm(TYPE, NOUT, NFAIL, NRUN, NERRS)
ALASVM
Definition: alasvm.f:75
subroutine alaerh(PATH, SUBNAM, INFO, INFOE, OPTS, M, N, KL, KU, N5, IMAT, NFAIL, NERRS, NOUT)
ALAERH
Definition: alaerh.f:149
subroutine slatb4(PATH, IMAT, M, N, TYPE, KL, KU, ANORM, MODE, CNDNUM, DIST)
SLATB4
Definition: slatb4.f:122
subroutine slarhs(PATH, XTYPE, UPLO, TRANS, M, N, KL, KU, NRHS, A, LDA, X, LDX, B, LDB, ISEED, INFO)
SLARHS
Definition: slarhs.f:206
subroutine spot01(UPLO, N, A, LDA, AFAC, LDAFAC, RWORK, RESID)
SPOT01
Definition: spot01.f:106
real function sget06(RCOND, RCONDC)
SGET06
Definition: sget06.f:57
subroutine spot05(UPLO, N, NRHS, A, LDA, B, LDB, X, LDX, XACT, LDXACT, FERR, BERR, RESLTS)
SPOT05
Definition: spot05.f:166
subroutine xlaenv(ISPEC, NVALUE)
XLAENV
Definition: xlaenv.f:83
subroutine slacpy(UPLO, M, N, A, LDA, B, LDB)
SLACPY copies all or part of one two-dimensional array to another.
Definition: slacpy.f:105
subroutine slatms(M, N, DIST, ISEED, SYM, D, MODE, COND, DMAX, KL, KU, PACK, A, LDA, WORK, INFO)
SLATMS
Definition: slatms.f:323
subroutine spotrf(UPLO, N, A, LDA, INFO)
SPOTRF
Definition: spotrf.f:109
subroutine slaset(UPLO, M, N, ALPHA, BETA, A, LDA)
SLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values...
Definition: slaset.f:112
subroutine slaqsy(UPLO, N, A, LDA, S, SCOND, AMAX, EQUED)
SLAQSY scales a symmetric/Hermitian matrix, using scaling factors computed by spoequ.
Definition: slaqsy.f:135
subroutine sposv(UPLO, N, NRHS, A, LDA, B, LDB, INFO)
SPOSV computes the solution to system of linear equations A * X = B for PO matrices ...
Definition: sposv.f:132
subroutine sget04(N, NRHS, X, LDX, XACT, LDXACT, RCOND, RESID)
SGET04
Definition: sget04.f:104
subroutine aladhd(IOUNIT, PATH)
ALADHD
Definition: aladhd.f:80
subroutine spoequ(N, A, LDA, S, SCOND, AMAX, INFO)
SPOEQU
Definition: spoequ.f:114
subroutine serrvx(PATH, NUNIT)
SERRVX
Definition: serrvx.f:57
subroutine spot02(UPLO, N, NRHS, A, LDA, X, LDX, B, LDB, RWORK, RESID)
SPOT02
Definition: spot02.f:129
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine spotri(UPLO, N, A, LDA, INFO)
SPOTRI
Definition: spotri.f:97
subroutine sposvx(FACT, UPLO, N, NRHS, A, LDA, AF, LDAF, EQUED, S, B, LDB, X, LDX, RCOND, FERR, BERR, WORK, IWORK, INFO)
SPOSVX computes the solution to system of linear equations A * X = B for PO matrices ...
Definition: sposvx.f:309
real function slansy(NORM, UPLO, N, A, LDA, WORK)
SLANSY returns the value of the 1-norm, or the Frobenius norm, or the infinity norm, or the element of largest absolute value of a real symmetric matrix.
Definition: slansy.f:124

Here is the call graph for this function:

Here is the caller graph for this function: