Actual source code: dvdinitv.c

slepc-3.6.1 2015-09-03
Report Typos and Errors
  1: /*
  2:   SLEPc eigensolver: "davidson"

  4:   Step: init subspace V

  6:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  7:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  8:    Copyright (c) 2002-2015, Universitat Politecnica de Valencia, Spain

 10:    This file is part of SLEPc.

 12:    SLEPc is free software: you can redistribute it and/or modify it under  the
 13:    terms of version 3 of the GNU Lesser General Public License as published by
 14:    the Free Software Foundation.

 16:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 17:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 18:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 19:    more details.

 21:    You  should have received a copy of the GNU Lesser General  Public  License
 22:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 23:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 24: */

 26:  #include davidson.h

 28: typedef struct {
 29:   PetscInt k;                 /* desired initial subspace size */
 30:   PetscInt user;              /* number of user initial vectors */
 31:   void     *old_initV_data;   /* old initV data */
 32: } dvdInitV;

 36: static PetscErrorCode dvd_initV_classic_0(dvdDashboard *d)
 37: {
 39:   dvdInitV       *data = (dvdInitV*)d->initV_data;
 40:   PetscInt       i,user = PetscMin(data->user,d->eps->ncv),k = PetscMin(data->k,d->eps->ncv);

 43:   /* Generate a set of random initial vectors and orthonormalize them */
 44:   for (i=user;i<k;i++) {
 45:     BVSetRandomColumn(d->eps->V,i,d->eps->rand);
 46:   }
 47:   d->V_tra_s = 0; d->V_tra_e = 0;
 48:   d->V_new_s = 0; d->V_new_e = i;

 50:   /* After that the user vectors will be destroyed */
 51:   data->user = 0;
 52:   return(0);
 53: }

 57: static PetscErrorCode dvd_initV_krylov_0(dvdDashboard *d)
 58: {
 60:   dvdInitV       *data = (dvdInitV*)d->initV_data;
 61:   PetscInt       i,user = PetscMin(data->user,d->eps->ncv),k = PetscMin(data->k,d->eps->ncv);
 62:   Vec            av,v1,v2;

 65:   /* If needed, generate a random vector for starting the arnoldi method */
 66:   if (user == 0) {
 67:     BVSetRandomColumn(d->eps->V,0,d->eps->rand);
 68:     user = 1;
 69:   }

 71:   /* Perform k steps of Arnoldi with the operator K^{-1}*(t[1]*A-t[2]*B) */
 72:   dvd_orthV(d->eps->V,0,user,d->eps->rand);
 73:   for (i=user;i<k;i++) {
 74:     /* aux <- theta[1]A*in - theta[0]*B*in */
 75:     BVGetColumn(d->eps->V,i,&v1);
 76:     BVGetColumn(d->eps->V,i-user,&v2);
 77:     BVGetColumn(d->auxBV,0,&av);
 78:     if (d->B) {
 79:       MatMult(d->A,v2,v1);
 80:       MatMult(d->B,v2,av);
 81:       VecAXPBY(av,d->target[1],-d->target[0],v1);
 82:     } else {
 83:       MatMult(d->A,v2,av);
 84:       VecAXPBY(av,-d->target[0],d->target[1],v2);
 85:     }
 86:     d->improvex_precond(d,0,av,v1);
 87:     BVRestoreColumn(d->eps->V,i,&v1);
 88:     BVRestoreColumn(d->eps->V,i-user,&v2);
 89:     BVRestoreColumn(d->auxBV,0,&av);
 90:     dvd_orthV(d->eps->V,i,i+1,d->eps->rand);
 91:   }

 93:   d->V_tra_s = 0; d->V_tra_e = 0;
 94:   d->V_new_s = 0; d->V_new_e = i;

 96:   /* After that the user vectors will be destroyed */
 97:   data->user = 0;
 98:   return(0);
 99: }

103: static PetscErrorCode dvd_initV_d(dvdDashboard *d)
104: {
106:   dvdInitV       *data = (dvdInitV*)d->initV_data;

109:   /* Restore changes in dvdDashboard */
110:   d->initV_data = data->old_initV_data;

112:   /* Free local data */
113:   PetscFree(data);
114:   return(0);
115: }

119: PetscErrorCode dvd_initV(dvdDashboard *d, dvdBlackboard *b, PetscInt k,PetscInt user, PetscBool krylov)
120: {
122:   dvdInitV       *data;

125:   /* Setting configuration constrains */
126:   b->max_size_V = PetscMax(b->max_size_V, k);

128:   /* Setup the step */
129:   if (b->state >= DVD_STATE_CONF) {
130:     PetscNewLog(d->eps,&data);
131:     data->k = k;
132:     data->user = PetscMin(k, user);
133:     data->old_initV_data = d->initV_data;
134:     d->initV_data = data;
135:     if (krylov) d->initV = dvd_initV_krylov_0;
136:     else d->initV = dvd_initV_classic_0;
137:     EPSDavidsonFLAdd(&d->destroyList,dvd_initV_d);
138:   }
139:   return(0);
140: }

144: PetscErrorCode dvd_orthV(BV V,PetscInt V_new_s,PetscInt V_new_e,PetscRandom rand)
145: {
147:   PetscInt       i,j,l,k;
148:   PetscBool      lindep;
149:   PetscReal      norm;

152:   BVGetActiveColumns(V,&l,&k);
153:   for (i=V_new_s;i<V_new_e;i++) {
154:     BVSetActiveColumns(V,0,i);
155:     for (j=0;j<3;j++) {
156:       if (j>0) {
157:         BVSetRandomColumn(V,i,rand);
158:         PetscInfo1(V,"Orthonormalization problems adding the vector %D to the searching subspace\n",i);
159:       }
160:       BVOrthogonalizeColumn(V,i,NULL,&norm,&lindep);
161:       if (!lindep && (PetscAbsReal(norm) > PETSC_SQRT_MACHINE_EPSILON)) break;
162:     }
163:     if (lindep || (PetscAbsReal(norm) < PETSC_SQRT_MACHINE_EPSILON)) SETERRQ(PetscObjectComm((PetscObject)V),1, "Error during the orthonormalization of the vectors");
164:     BVScaleColumn(V,i,1.0/norm);
165:   }
166:   BVSetActiveColumns(V,l,k);
167:   return(0);
168: }