LCOV - code coverage report
Current view: top level - sys/classes/bv/interface - bvkrylov.c (source / functions) Hit Total Coverage
Test: SLEPc Lines: 97 97 100.0 %
Date: 2024-04-24 00:34:25 Functions: 2 2 100.0 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : /*
       2             :    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
       3             :    SLEPc - Scalable Library for Eigenvalue Problem Computations
       4             :    Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
       5             : 
       6             :    This file is part of SLEPc.
       7             :    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
       8             :    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
       9             : */
      10             : /*
      11             :    BV routines related to Krylov decompositions
      12             : */
      13             : 
      14             : #include <slepc/private/bvimpl.h>          /*I   "slepcbv.h"   I*/
      15             : 
      16             : /*@
      17             :    BVMatArnoldi - Computes an Arnoldi factorization associated with a matrix.
      18             : 
      19             :    Collective
      20             : 
      21             :    Input Parameters:
      22             : +  V - basis vectors context
      23             : .  A - the matrix
      24             : .  H - (optional) the upper Hessenberg matrix
      25             : .  k - number of locked columns
      26             : -  m - dimension of the Arnoldi basis, may be modified
      27             : 
      28             :    Output Parameters:
      29             : +  beta - (optional) norm of last vector before normalization
      30             : -  breakdown - (optional) flag indicating that breakdown occurred
      31             : 
      32             :    Notes:
      33             :    Computes an m-step Arnoldi factorization for matrix A. The first k columns
      34             :    are assumed to be locked and therefore they are not modified. On exit, the
      35             :    following relation is satisfied
      36             : 
      37             : $                    A * V - V * H = beta*v_m * e_m^T
      38             : 
      39             :    where the columns of V are the Arnoldi vectors (which are orthonormal), H is
      40             :    an upper Hessenberg matrix, e_m is the m-th vector of the canonical basis.
      41             :    On exit, beta contains the norm of V[m] before normalization.
      42             : 
      43             :    The breakdown flag indicates that orthogonalization failed, see
      44             :    BVOrthonormalizeColumn(). In that case, on exit m contains the index of
      45             :    the column that failed.
      46             : 
      47             :    The values of k and m are not restricted to the active columns of V.
      48             : 
      49             :    To create an Arnoldi factorization from scratch, set k=0 and make sure the
      50             :    first column contains the normalized initial vector.
      51             : 
      52             :    Level: advanced
      53             : 
      54             : .seealso: BVMatLanczos(), BVSetActiveColumns(), BVOrthonormalizeColumn()
      55             : @*/
      56        2383 : PetscErrorCode BVMatArnoldi(BV V,Mat A,Mat H,PetscInt k,PetscInt *m,PetscReal *beta,PetscBool *breakdown)
      57             : {
      58        2383 :   PetscScalar       *h;
      59        2383 :   const PetscScalar *a;
      60        2383 :   PetscInt          j,ldh,rows,cols;
      61        2383 :   PetscBool         lindep=PETSC_FALSE;
      62        2383 :   Vec               buf;
      63             : 
      64        2383 :   PetscFunctionBegin;
      65        2383 :   PetscValidHeaderSpecific(V,BV_CLASSID,1);
      66        2383 :   PetscValidHeaderSpecific(A,MAT_CLASSID,2);
      67        9532 :   PetscValidLogicalCollectiveInt(V,k,4);
      68        2383 :   PetscAssertPointer(m,5);
      69        9532 :   PetscValidLogicalCollectiveInt(V,*m,5);
      70        2383 :   PetscValidType(V,1);
      71        2383 :   BVCheckSizes(V,1);
      72        2383 :   PetscValidType(A,2);
      73        2383 :   PetscCheckSameComm(V,1,A,2);
      74             : 
      75        2383 :   PetscCheck(k>=0 && k<=V->m,PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Argument k has wrong value %" PetscInt_FMT ", should be between 0 and %" PetscInt_FMT,k,V->m);
      76        2383 :   PetscCheck(*m>0 && *m<=V->m,PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Argument m has wrong value %" PetscInt_FMT ", should be between 1 and %" PetscInt_FMT,*m,V->m);
      77        2383 :   PetscCheck(*m>k,PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Argument m should be at least equal to k+1");
      78        2383 :   if (H) {
      79        2383 :     PetscValidHeaderSpecific(H,MAT_CLASSID,3);
      80        2383 :     PetscValidType(H,3);
      81        2383 :     PetscCheckTypeName(H,MATSEQDENSE);
      82        2383 :     PetscCall(MatGetSize(H,&rows,&cols));
      83        2383 :     PetscCall(MatDenseGetLDA(H,&ldh));
      84        2383 :     PetscCheck(rows>=*m,PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_SIZ,"Matrix H has %" PetscInt_FMT " rows, should have at least %" PetscInt_FMT,rows,*m);
      85        2383 :     PetscCheck(cols>=*m,PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_SIZ,"Matrix H has %" PetscInt_FMT " columns, should have at least %" PetscInt_FMT,cols,*m);
      86             :   }
      87             : 
      88       36678 :   for (j=k;j<*m;j++) {
      89       34335 :     PetscCall(BVMatMultColumn(V,A,j));
      90       34335 :     if (PetscUnlikely(j==V->N-1)) PetscCall(BV_OrthogonalizeColumn_Safe(V,j+1,NULL,beta,&lindep)); /* safeguard in case the full basis is requested */
      91       34313 :     else PetscCall(BVOrthonormalizeColumn(V,j+1,PETSC_FALSE,beta,&lindep));
      92       34335 :     if (PetscUnlikely(lindep)) {
      93          40 :       *m = j+1;
      94          40 :       break;
      95             :     }
      96             :   }
      97        2383 :   if (breakdown) *breakdown = lindep;
      98        2383 :   if (lindep) PetscCall(PetscInfo(V,"Arnoldi finished early at m=%" PetscInt_FMT "\n",*m));
      99             : 
     100        2383 :   if (H) {
     101        2383 :     PetscCall(MatDenseGetArray(H,&h));
     102        2383 :     PetscCall(BVGetBufferVec(V,&buf));
     103        2383 :     PetscCall(VecGetArrayRead(buf,&a));
     104       34335 :     for (j=k;j<*m-1;j++) PetscCall(PetscArraycpy(h+j*ldh,a+V->nc+(j+1)*(V->nc+V->m),j+2));
     105        2383 :     PetscCall(PetscArraycpy(h+(*m-1)*ldh,a+V->nc+(*m)*(V->nc+V->m),*m));
     106        2383 :     if (ldh>*m) h[(*m)+(*m-1)*ldh] = a[V->nc+(*m)+(*m)*(V->nc+V->m)];
     107        2383 :     PetscCall(VecRestoreArrayRead(buf,&a));
     108        2383 :     PetscCall(MatDenseRestoreArray(H,&h));
     109             :   }
     110             : 
     111        2383 :   PetscCall(PetscObjectStateIncrease((PetscObject)V));
     112        2383 :   PetscFunctionReturn(PETSC_SUCCESS);
     113             : }
     114             : 
     115             : /*@C
     116             :    BVMatLanczos - Computes a Lanczos factorization associated with a matrix.
     117             : 
     118             :    Collective
     119             : 
     120             :    Input Parameters:
     121             : +  V - basis vectors context
     122             : .  A - the matrix
     123             : .  T - (optional) the tridiagonal matrix
     124             : .  k - number of locked columns
     125             : -  m - dimension of the Lanczos basis, may be modified
     126             : 
     127             :    Output Parameters:
     128             : +  beta - (optional) norm of last vector before normalization
     129             : -  breakdown - (optional) flag indicating that breakdown occurred
     130             : 
     131             :    Notes:
     132             :    Computes an m-step Lanczos factorization for matrix A, with full
     133             :    reorthogonalization. At each Lanczos step, the corresponding Lanczos
     134             :    vector is orthogonalized with respect to all previous Lanczos vectors.
     135             :    This is equivalent to computing an m-step Arnoldi factorization and
     136             :    exploting symmetry of the operator.
     137             : 
     138             :    The first k columns are assumed to be locked and therefore they are
     139             :    not modified. On exit, the following relation is satisfied
     140             : 
     141             : $                    A * V - V * T = beta*v_m * e_m^T
     142             : 
     143             :    where the columns of V are the Lanczos vectors (which are B-orthonormal),
     144             :    T is a real symmetric tridiagonal matrix, and e_m is the m-th vector of
     145             :    the canonical basis. On exit, beta contains the B-norm of V[m] before
     146             :    normalization. The T matrix is stored in a special way, its first column
     147             :    contains the diagonal elements, and its second column the off-diagonal
     148             :    ones. In complex scalars, the elements are stored as PetscReal and thus
     149             :    occupy only the first column of the Mat object. This is the same storage
     150             :    scheme used in matrix DS_MAT_T obtained with DSGetMat().
     151             : 
     152             :    The breakdown flag indicates that orthogonalization failed, see
     153             :    BVOrthonormalizeColumn(). In that case, on exit m contains the index of
     154             :    the column that failed.
     155             : 
     156             :    The values of k and m are not restricted to the active columns of V.
     157             : 
     158             :    To create a Lanczos factorization from scratch, set k=0 and make sure the
     159             :    first column contains the normalized initial vector.
     160             : 
     161             :    Level: advanced
     162             : 
     163             : .seealso: BVMatArnoldi(), BVSetActiveColumns(), BVOrthonormalizeColumn(), DSGetMat()
     164             : @*/
     165        1447 : PetscErrorCode BVMatLanczos(BV V,Mat A,Mat T,PetscInt k,PetscInt *m,PetscReal *beta,PetscBool *breakdown)
     166             : {
     167        1447 :   PetscScalar       *t;
     168        1447 :   const PetscScalar *a;
     169        1447 :   PetscReal         *alpha,*betat;
     170        1447 :   PetscInt          j,ldt,rows,cols,mincols=PetscDefined(USE_COMPLEX)?1:2;
     171        1447 :   PetscBool         lindep=PETSC_FALSE;
     172        1447 :   Vec               buf;
     173             : 
     174        1447 :   PetscFunctionBegin;
     175        1447 :   PetscValidHeaderSpecific(V,BV_CLASSID,1);
     176        1447 :   PetscValidHeaderSpecific(A,MAT_CLASSID,2);
     177        5788 :   PetscValidLogicalCollectiveInt(V,k,4);
     178        1447 :   PetscAssertPointer(m,5);
     179        5788 :   PetscValidLogicalCollectiveInt(V,*m,5);
     180        1447 :   PetscValidType(V,1);
     181        1447 :   BVCheckSizes(V,1);
     182        1447 :   PetscValidType(A,2);
     183        1447 :   PetscCheckSameComm(V,1,A,2);
     184             : 
     185        1447 :   PetscCheck(k>=0 || k<=V->m,PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Argument k has wrong value %" PetscInt_FMT ", should be between 0 and %" PetscInt_FMT,k,V->m);
     186        1447 :   PetscCheck(*m>0 || *m<=V->m,PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Argument m has wrong value %" PetscInt_FMT ", should be between 1 and %" PetscInt_FMT,*m,V->m);
     187        1447 :   PetscCheck(*m>k,PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_OUTOFRANGE,"Argument m should be at least equal to k+1");
     188        1447 :   if (T) {
     189        1447 :     PetscValidHeaderSpecific(T,MAT_CLASSID,3);
     190        1447 :     PetscValidType(T,3);
     191        1447 :     PetscCheckTypeName(T,MATSEQDENSE);
     192        1447 :     PetscCall(MatGetSize(T,&rows,&cols));
     193        1447 :     PetscCall(MatDenseGetLDA(T,&ldt));
     194        1447 :     PetscCheck(rows>=*m,PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_SIZ,"Matrix T has %" PetscInt_FMT " rows, should have at least %" PetscInt_FMT,rows,*m);
     195        1447 :     PetscCheck(cols>=mincols,PetscObjectComm((PetscObject)V),PETSC_ERR_ARG_SIZ,"Matrix T has %" PetscInt_FMT " columns, should have at least %" PetscInt_FMT,cols,mincols);
     196             :   }
     197             : 
     198       19259 :   for (j=k;j<*m;j++) {
     199       17833 :     PetscCall(BVMatMultColumn(V,A,j));
     200       17833 :     if (PetscUnlikely(j==V->N-1)) PetscCall(BV_OrthogonalizeColumn_Safe(V,j+1,NULL,beta,&lindep)); /* safeguard in case the full basis is requested */
     201       17812 :     else PetscCall(BVOrthonormalizeColumn(V,j+1,PETSC_FALSE,beta,&lindep));
     202       17833 :     if (PetscUnlikely(lindep)) {
     203          21 :       *m = j+1;
     204          21 :       break;
     205             :     }
     206             :   }
     207        1447 :   if (breakdown) *breakdown = lindep;
     208        1447 :   if (lindep) PetscCall(PetscInfo(V,"Lanczos finished early at m=%" PetscInt_FMT "\n",*m));
     209             : 
     210        1447 :   if (T) {
     211        1447 :     PetscCall(MatDenseGetArray(T,&t));
     212        1447 :     alpha = (PetscReal*)t;
     213        1447 :     betat = alpha+ldt;
     214        1447 :     PetscCall(BVGetBufferVec(V,&buf));
     215        1447 :     PetscCall(VecGetArrayRead(buf,&a));
     216       19280 :     for (j=k;j<*m;j++) {
     217       17833 :       alpha[j] = PetscRealPart(a[V->nc+j+(j+1)*(V->nc+V->m)]);
     218       17833 :       betat[j] = PetscRealPart(a[V->nc+j+1+(j+1)*(V->nc+V->m)]);
     219             :     }
     220        1447 :     PetscCall(VecRestoreArrayRead(buf,&a));
     221        1447 :     PetscCall(MatDenseRestoreArray(T,&t));
     222             :   }
     223             : 
     224        1447 :   PetscCall(PetscObjectStateIncrease((PetscObject)V));
     225        1447 :   PetscFunctionReturn(PETSC_SUCCESS);
     226             : }

Generated by: LCOV version 1.14