Subversion Repositories slepc-dev

Rev

Go to most recent revision | Details | Compare with Previous | Last modification | View Log | RSS feed

Rev Author Line No. Line
1302 slepc 1
/*
2
     Routines related to orthogonalization.
1303 slepc 3
     See the SLEPc Technical Report STR-1 for a detailed explanation.
1376 slepc 4
 
5
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1672 slepc 6
   SLEPc - Scalable Library for Eigenvalue Problem Computations
7
   Copyright (c) 2002-2009, Universidad Politecnica de Valencia, Spain
1376 slepc 8
 
1672 slepc 9
   This file is part of SLEPc.
10
 
11
   SLEPc is free software: you can redistribute it and/or modify it under  the
12
   terms of version 3 of the GNU Lesser General Public License as published by
13
   the Free Software Foundation.
14
 
15
   SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
16
   WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
17
   FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
18
   more details.
19
 
20
   You  should have received a copy of the GNU Lesser General  Public  License
21
   along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
1376 slepc 22
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1302 slepc 23
*/
1376 slepc 24
 
1521 slepc 25
#include "private/ipimpl.h"      /*I "slepcip.h" I*/
1302 slepc 26
#include "slepcblaslapack.h"
27
 
28
/*
1307 slepc 29
   IPOrthogonalizeMGS - Compute one step of Modified Gram-Schmidt
1302 slepc 30
*/
31
#undef __FUNCT__  
1307 slepc 32
#define __FUNCT__ "IPOrthogonalizeMGS"
1509 slepc 33
static PetscErrorCode IPOrthogonalizeMGS(IP ip,PetscInt n,PetscTruth *which,Vec *V,Vec v,PetscScalar *H)
1302 slepc 34
{
35
  PetscErrorCode ierr;
1509 slepc 36
  PetscInt       j;
1307 slepc 37
 
38
  PetscFunctionBegin;
39
  for (j=0; j<n; j++)
40
    if (!which || which[j]) {
41
      /* h_j = ( v, v_j ) */
42
      ierr = IPInnerProduct(ip,v,V[j],&H[j]);CHKERRQ(ierr);
43
      /* v <- v - h_j v_j */
44
      ierr = VecAXPY(v,-H[j],V[j]);CHKERRQ(ierr);
45
    }
46
  PetscFunctionReturn(0);
47
}
48
 
49
/*
50
   IPOrthogonalizeCGS - Compute |v'| (estimated), |v| and one step of CGS with only one global synchronization
51
*/
52
#undef __FUNCT__  
53
#define __FUNCT__ "IPOrthogonalizeCGS"
1509 slepc 54
PetscErrorCode IPOrthogonalizeCGS(IP ip,PetscInt n,PetscTruth *which,Vec *V,Vec v,PetscScalar *H,PetscReal *onorm,PetscReal *norm,Vec work)
1307 slepc 55
{
56
  PetscErrorCode ierr;
1509 slepc 57
  PetscInt       j;
1302 slepc 58
  PetscScalar    alpha;
59
  PetscReal      sum;
60
 
61
  PetscFunctionBegin;
1307 slepc 62
  /* h = W^* v ; alpha = (v , v) */
63
  if (which) {
64
  /* use select array */
65
    for (j=0; j<n; j++)
1329 slepc 66
      if (which[j]) { ierr = IPInnerProductBegin(ip,v,V[j],&H[j]);CHKERRQ(ierr); }
1307 slepc 67
    if (onorm || norm) { ierr = IPInnerProductBegin(ip,v,v,&alpha);CHKERRQ(ierr); }
68
    for (j=0; j<n; j++)
69
      if (which[j]) { ierr = IPInnerProductEnd(ip,v,V[j],&H[j]);CHKERRQ(ierr); }
70
    if (onorm || norm) { ierr = IPInnerProductEnd(ip,v,v,&alpha);CHKERRQ(ierr); }
71
  } else { /* merge comunications */
72
    if (onorm || norm) {
1381 slepc 73
      ierr = IPMInnerProductBegin(ip,v,n,V,H);CHKERRQ(ierr);
1307 slepc 74
      ierr = IPInnerProductBegin(ip,v,v,&alpha);CHKERRQ(ierr);
1381 slepc 75
      ierr = IPMInnerProductEnd(ip,v,n,V,H);CHKERRQ(ierr);
1307 slepc 76
      ierr = IPInnerProductEnd(ip,v,v,&alpha);CHKERRQ(ierr);
77
    } else { /* use simpler function */
1381 slepc 78
      ierr = IPMInnerProduct(ip,v,n,V,H);CHKERRQ(ierr);
1302 slepc 79
    }
1307 slepc 80
  }
1302 slepc 81
 
1307 slepc 82
  /* q = v - V h */
1329 slepc 83
  ierr = VecSet(work,0.0);CHKERRQ(ierr);
1307 slepc 84
  if (which) {
85
    for (j=0; j<n; j++)
1329 slepc 86
      if (which[j]) { ierr = VecAXPY(work,H[j],V[j]);CHKERRQ(ierr); }
1307 slepc 87
  } else {
1329 slepc 88
    ierr = VecMAXPY(work,n,H,V);CHKERRQ(ierr);  
1307 slepc 89
  }
1329 slepc 90
  ierr = VecAXPY(v,-1.0,work);CHKERRQ(ierr);
1302 slepc 91
 
1307 slepc 92
  /* compute |v| */
93
  if (onorm) *onorm = sqrt(PetscRealPart(alpha));
94
 
95
  /* compute |v'| */
96
  if (norm) {
97
    sum = 0.0;
98
    for (j=0; j<n; j++)
99
      if (!which || which[j])
100
        sum += PetscRealPart(H[j] * PetscConj(H[j]));
101
    *norm = PetscRealPart(alpha)-sum;
1329 slepc 102
    if (*norm <= 0.0) {
1307 slepc 103
      ierr = IPNorm(ip,v,norm);CHKERRQ(ierr);
104
    } else *norm = sqrt(*norm);
105
  }
106
  PetscFunctionReturn(0);
107
}
108
 
109
/*
110
   IPOrthogonalizeGS - Compute |v'|, |v| and one step of CGS or MGS
111
*/
112
#undef __FUNCT__  
113
#define __FUNCT__ "IPOrthogonalizeGS"
1509 slepc 114
static PetscErrorCode IPOrthogonalizeGS(IP ip,PetscInt n,PetscTruth *which,Vec *V,Vec v,PetscScalar *H,PetscReal *onorm,PetscReal *norm,Vec work)
1307 slepc 115
{
116
  PetscErrorCode ierr;
117
 
118
  PetscFunctionBegin;
119
  switch (ip->orthog_type) {
120
  case IP_CGS_ORTH:
1329 slepc 121
    ierr = IPOrthogonalizeCGS(ip,n,which,V,v,H,onorm,norm,work);CHKERRQ(ierr);
1302 slepc 122
    break;
123
  case IP_MGS_ORTH:
124
    /* compute |v| */
125
    if (onorm) { ierr = IPNorm(ip,v,onorm);CHKERRQ(ierr); }
1307 slepc 126
    ierr = IPOrthogonalizeMGS(ip,n,which,V,v,H);CHKERRQ(ierr);
1302 slepc 127
    /* compute |v'| */
128
    if (norm) { ierr = IPNorm(ip,v,norm);CHKERRQ(ierr); }
129
    break;
130
  default:
131
    SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown orthogonalization type");
132
  }
133
  PetscFunctionReturn(0);
134
}
135
 
136
#undef __FUNCT__  
137
#define __FUNCT__ "IPOrthogonalize"
138
/*@
139
   IPOrthogonalize - Orthogonalize a vector with respect to a set of vectors.
140
 
141
   Collective on IP
142
 
143
   Input Parameters:
1303 slepc 144
+  ip    - the inner product (IP) context
1302 slepc 145
.  n      - number of columns of V
146
.  which  - logical array indicating columns of V to be used
1381 slepc 147
.  V      - set of vectors
1538 slepc 148
.  work   - workspace vector
149
-  swork  - workspace array
1302 slepc 150
 
151
   Input/Output Parameter:
152
.  v      - (input) vector to be orthogonalized and (output) result of
153
            orthogonalization
154
 
155
   Output Parameter:
156
+  H      - coefficients computed during orthogonalization
157
.  norm   - norm of the vector after being orthogonalized
158
-  lindep - flag indicating that refinement did not improve the quality
159
            of orthogonalization
160
 
161
   Notes:
162
   This function applies an orthogonal projector to project vector v onto the
163
   orthogonal complement of the span of the columns of V.
164
 
165
   On exit, v0 = [V v]*H, where v0 is the original vector v.
166
 
167
   This routine does not normalize the resulting vector.
168
 
169
   Level: developer
170
 
171
.seealso: IPSetOrthogonalization(), IPBiOrthogonalize()
172
@*/
1538 slepc 173
PetscErrorCode IPOrthogonalize(IP ip,PetscInt n,PetscTruth *which,Vec *V,Vec v,PetscScalar *H,PetscReal *norm,PetscTruth *lindep,Vec work,PetscScalar* swork)
1302 slepc 174
{
175
  PetscErrorCode ierr;
176
  PetscScalar    lh[100],*h,lc[100],*c;
1345 slepc 177
  PetscTruth     allocatedh = PETSC_FALSE,allocatedc = PETSC_FALSE,allocatedw = PETSC_FALSE;
1302 slepc 178
  PetscReal      onrm,nrm;
1509 slepc 179
  PetscInt       j,k;
1302 slepc 180
  PetscFunctionBegin;
181
  if (n==0) {
182
    if (norm) { ierr = IPNorm(ip,v,norm);CHKERRQ(ierr); }
183
    if (lindep) *lindep = PETSC_FALSE;
184
    PetscFunctionReturn(0);
185
  }
186
  ierr = PetscLogEventBegin(IP_Orthogonalize,ip,0,0,0);CHKERRQ(ierr);
187
 
1345 slepc 188
  /* allocate H, c and work if needed */
1302 slepc 189
  if (!H) {
190
    if (n<=100) h = lh;
191
    else {
192
      ierr = PetscMalloc(n*sizeof(PetscScalar),&h);CHKERRQ(ierr);
193
      allocatedh = PETSC_TRUE;
194
    }
195
  } else h = H;
1538 slepc 196
  if (!swork) {
197
    if (ip->orthog_ref != IP_ORTH_REFINE_NEVER) {
198
      if (n<=100) c = lc;
199
      else {
200
        ierr = PetscMalloc(n*sizeof(PetscScalar),&c);CHKERRQ(ierr);
201
        allocatedc = PETSC_TRUE;
202
      }
1302 slepc 203
    }
1538 slepc 204
  } else c = swork;
1436 slepc 205
  if (!work && ip->orthog_type == IP_CGS_ORTH) {
1345 slepc 206
    ierr = VecDuplicate(v,&work);CHKERRQ(ierr);
1348 slepc 207
    allocatedw = PETSC_TRUE;
1345 slepc 208
  }
1302 slepc 209
 
210
  /* orthogonalize and compute onorm */
211
  switch (ip->orthog_ref) {
212
 
213
  case IP_ORTH_REFINE_NEVER:
1329 slepc 214
    ierr = IPOrthogonalizeGS(ip,n,which,V,v,h,PETSC_NULL,PETSC_NULL,work);CHKERRQ(ierr);
1302 slepc 215
    /* compute |v| */
216
    if (norm) { ierr = IPNorm(ip,v,norm);CHKERRQ(ierr); }
217
    /* linear dependence check does not work without refinement */
218
    if (lindep) *lindep = PETSC_FALSE;
219
    break;
220
 
221
  case IP_ORTH_REFINE_ALWAYS:
1329 slepc 222
    ierr = IPOrthogonalizeGS(ip,n,which,V,v,h,PETSC_NULL,PETSC_NULL,work);CHKERRQ(ierr);
1302 slepc 223
    if (lindep) {
1329 slepc 224
      ierr = IPOrthogonalizeGS(ip,n,which,V,v,c,&onrm,&nrm,work);CHKERRQ(ierr);
1302 slepc 225
      if (norm) *norm = nrm;
226
      if (nrm < ip->orthog_eta * onrm) *lindep = PETSC_TRUE;
227
      else *lindep = PETSC_FALSE;
228
    } else {
1329 slepc 229
      ierr = IPOrthogonalizeGS(ip,n,which,V,v,c,PETSC_NULL,norm,work);CHKERRQ(ierr);
1302 slepc 230
    }
231
    for (j=0;j<n;j++)
232
      if (!which || which[j]) h[j] += c[j];
233
    break;
234
 
235
  case IP_ORTH_REFINE_IFNEEDED:
1329 slepc 236
    ierr = IPOrthogonalizeGS(ip,n,which,V,v,h,&onrm,&nrm,work);CHKERRQ(ierr);
1302 slepc 237
    /* ||q|| < eta ||h|| */
238
    k = 1;
239
    while (k<3 && nrm < ip->orthog_eta * onrm) {
240
      k++;
241
      switch (ip->orthog_type) {
242
      case IP_CGS_ORTH:
1329 slepc 243
        ierr = IPOrthogonalizeGS(ip,n,which,V,v,c,&onrm,&nrm,work);CHKERRQ(ierr);
1302 slepc 244
        break;
245
      case IP_MGS_ORTH:
246
        onrm = nrm;
1329 slepc 247
        ierr = IPOrthogonalizeGS(ip,n,which,V,v,c,PETSC_NULL,&nrm,work);CHKERRQ(ierr);
1302 slepc 248
        break;
249
      default:
250
        SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown orthogonalization type");
251
      }        
252
      for (j=0;j<n;j++)
253
        if (!which || which[j]) h[j] += c[j];
254
    }
255
    if (norm) *norm = nrm;
256
    if (lindep) {
257
      if (nrm < ip->orthog_eta * onrm) *lindep = PETSC_TRUE;
258
      else *lindep = PETSC_FALSE;
259
    }
260
 
261
    break;
262
  default:
263
    SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown orthogonalization refinement");
264
  }
265
 
266
  /* free work space */
267
  if (allocatedc) { ierr = PetscFree(c);CHKERRQ(ierr); }
268
  if (allocatedh) { ierr = PetscFree(h);CHKERRQ(ierr); }
1345 slepc 269
  if (allocatedw) { ierr = VecDestroy(work);CHKERRQ(ierr); }
1302 slepc 270
 
271
  ierr = PetscLogEventEnd(IP_Orthogonalize,ip,0,0,0);CHKERRQ(ierr);
272
  PetscFunctionReturn(0);
273
}
1345 slepc 274
 
275
#undef __FUNCT__  
276
#define __FUNCT__ "IPQRDecomposition"
277
/*@
278
   IPQRDecomposition - Compute the QR factorization of a set of vectors.
279
 
280
   Collective on IP
281
 
282
   Input Parameters:
283
+  ip - the eigenproblem solver context
284
.  V - set of vectors
285
.  m - starting column
286
.  n - ending column
1381 slepc 287
.  ldr - leading dimension of R
288
-  work - workspace vector
1345 slepc 289
 
290
   Output Parameter:
291
.  R  - triangular matrix of QR factorization
292
 
293
   Notes:
294
   This routine orthonormalizes the columns of V so that V'*V=I up to
295
   working precision. It assumes that the first m columns (from 0 to m-1)
296
   are already orthonormal. The coefficients of the orthogonalization are
297
   stored in R so that V*R is equal to the original V.
298
 
299
   In some cases, this routine makes V B-orthonormal, that is, V'*B*V=I.
300
 
301
   If one of the vectors is linearly dependent wrt the rest (at working
302
   precision) then it is replaced by a random vector.
303
 
304
   Level: developer
305
 
306
.seealso: IPOrthogonalize(), IPNorm(), IPInnerProduct().
307
@*/
1509 slepc 308
PetscErrorCode IPQRDecomposition(IP ip,Vec *V,PetscInt m,PetscInt n,PetscScalar *R,PetscInt ldr,Vec work)
1345 slepc 309
{
310
  PetscErrorCode ierr;
1509 slepc 311
  PetscInt       k;
1345 slepc 312
  PetscReal      norm;
313
  PetscTruth     lindep;
314
 
315
  PetscFunctionBegin;
316
 
317
  for (k=m; k<n; k++) {
318
 
319
    /* orthogonalize v_k with respect to v_0, ..., v_{k-1} */
1538 slepc 320
    if (R) { ierr = IPOrthogonalize(ip,k,PETSC_NULL,V,V[k],&R[0+ldr*k],&norm,&lindep,work,PETSC_NULL);CHKERRQ(ierr); }
321
    else   { ierr = IPOrthogonalize(ip,k,PETSC_NULL,V,V[k],PETSC_NULL,&norm,&lindep,work,PETSC_NULL);CHKERRQ(ierr); }
1345 slepc 322
 
323
    /* normalize v_k: r_{k,k} = ||v_k||_2; v_k = v_k/r_{k,k} */
324
    if (norm==0.0 || lindep) {
325
      PetscInfo(ip,"Linearly dependent vector found, generating a new random vector\n");
326
      ierr = SlepcVecSetRandom(V[k]);CHKERRQ(ierr);
327
      ierr = IPNorm(ip,V[k],&norm);CHKERRQ(ierr);
328
    }
329
    ierr = VecScale(V[k],1.0/norm);CHKERRQ(ierr);
330
    if (R) R[k+ldr*k] = norm;
331
 
332
  }
333
 
334
  PetscFunctionReturn(0);
335
}
336
 
337
/*
338
    Biorthogonalization routine using classical Gram-Schmidt with refinement.
339
 */
340
#undef __FUNCT__  
341
#define __FUNCT__ "IPCGSBiOrthogonalization"
1509 slepc 342
static PetscErrorCode IPCGSBiOrthogonalization(IP ip,PetscInt n_,Vec *V,Vec *W,Vec v,PetscScalar *H,PetscReal *hnorm,PetscReal *norm)
1345 slepc 343
{
344
#if defined(SLEPC_MISSING_LAPACK_GELQF) || defined(SLEPC_MISSING_LAPACK_ORMLQ)
345
  PetscFunctionBegin;
346
  SETERRQ(PETSC_ERR_SUP,"xGELQF - Lapack routine is unavailable.");
347
#else
348
  PetscErrorCode ierr;
1509 slepc 349
  PetscBLASInt   j,ione=1,lwork,info,n=n_;
1345 slepc 350
  PetscScalar    shh[100],*lhh,*vw,*tau,one=1.0,*work;
351
  Vec            w;
352
 
353
  PetscFunctionBegin;
354
 
355
  /* Don't allocate small arrays */
356
  if (n<=100) lhh = shh;
357
  else { ierr = PetscMalloc(n*sizeof(PetscScalar),&lhh);CHKERRQ(ierr); }
358
  ierr = PetscMalloc(n*n*sizeof(PetscScalar),&vw);CHKERRQ(ierr);
359
  ierr = VecDuplicate(v,&w);CHKERRQ(ierr);
360
 
361
  for (j=0;j<n;j++) {
1381 slepc 362
    ierr = IPMInnerProduct(ip,V[j],n,W,vw+j*n);CHKERRQ(ierr);
1345 slepc 363
  }
364
  lwork = n;
365
  ierr = PetscMalloc(n*sizeof(PetscScalar),&tau);CHKERRQ(ierr);
366
  ierr = PetscMalloc(lwork*sizeof(PetscScalar),&work);CHKERRQ(ierr);
367
  LAPACKgelqf_(&n,&n,vw,&n,tau,work,&lwork,&info);
368
  if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xGELQF %i",info);
369
 
370
  /*** First orthogonalization ***/
371
 
372
  /* h = W^* v */
373
  /* q = v - V h */
1381 slepc 374
  ierr = IPMInnerProduct(ip,v,n,W,H);CHKERRQ(ierr);
1536 slepc 375
  BLAStrsm_("L","L","N","N",&n,&ione,&one,vw,&n,H,&n);
376
  LAPACKormlq_("L","N",&n,&ione,&n,vw,&n,tau,H,&n,work,&lwork,&info);
1345 slepc 377
  if (info) SETERRQ1(PETSC_ERR_LIB,"Error in Lapack xORMLQ %i",info);
378
  ierr = VecSet(w,0.0);CHKERRQ(ierr);
379
  ierr = VecMAXPY(w,n,H,V);CHKERRQ(ierr);
380
  ierr = VecAXPY(v,-1.0,w);CHKERRQ(ierr);
381
 
382
  /* compute norm of v */
383
  if (norm) { ierr = IPNorm(ip,v,norm);CHKERRQ(ierr); }
384
 
385
  if (n>100) { ierr = PetscFree(lhh);CHKERRQ(ierr); }
386
  ierr = PetscFree(vw);CHKERRQ(ierr);
387
  ierr = PetscFree(tau);CHKERRQ(ierr);
388
  ierr = PetscFree(work);CHKERRQ(ierr);
389
  ierr = VecDestroy(w);CHKERRQ(ierr);
390
  PetscFunctionReturn(0);
391
#endif
392
}
393
 
394
#undef __FUNCT__  
395
#define __FUNCT__ "IPBiOrthogonalize"
396
/*@
397
   IPBiOrthogonalize - Bi-orthogonalize a vector with respect to a set of vectors.
398
 
399
   Collective on IP
400
 
401
   Input Parameters:
1349 slepc 402
+  ip - the inner product context
1345 slepc 403
.  n - number of columns of V
404
.  V - set of vectors
405
-  W - set of vectors
406
 
407
   Input/Output Parameter:
408
.  v - vector to be orthogonalized
409
 
410
   Output Parameter:
411
+  H  - coefficients computed during orthogonalization
412
-  norm - norm of the vector after being orthogonalized
413
 
414
   Notes:
415
   This function applies an oblique projector to project vector v onto the
416
   span of the columns of V along the orthogonal complement of the column
417
   space of W.
418
 
419
   On exit, v0 = [V v]*H, where v0 is the original vector v.
420
 
421
   This routine does not normalize the resulting vector.
422
 
423
   Level: developer
424
 
425
.seealso: IPSetOrthogonalization(), IPOrthogonalize()
426
@*/
1509 slepc 427
PetscErrorCode IPBiOrthogonalize(IP ip,PetscInt n,Vec *V,Vec *W,Vec v,PetscScalar *H,PetscReal *norm)
1345 slepc 428
{
429
  PetscErrorCode ierr;
430
  PetscScalar    lh[100],*h;
431
  PetscTruth     allocated = PETSC_FALSE;
432
  PetscReal      lhnrm,*hnrm,lnrm,*nrm;
433
  PetscFunctionBegin;
434
  if (n==0) {
435
    if (norm) { ierr = IPNorm(ip,v,norm);CHKERRQ(ierr); }
436
  } else {
437
    ierr = PetscLogEventBegin(IP_Orthogonalize,ip,0,0,0);CHKERRQ(ierr);
438
 
439
    /* allocate H if needed */
440
    if (!H) {
441
      if (n<=100) h = lh;
442
      else {
443
        ierr = PetscMalloc(n*sizeof(PetscScalar),&h);CHKERRQ(ierr);
444
        allocated = PETSC_TRUE;
445
      }
446
    } else h = H;
447
 
448
    /* retrieve hnrm and nrm for linear dependence check or conditional refinement */
449
    if (ip->orthog_ref == IP_ORTH_REFINE_IFNEEDED) {
450
      hnrm = &lhnrm;
451
      if (norm) nrm = norm;
452
      else nrm = &lnrm;
453
    } else {
454
      hnrm = PETSC_NULL;
455
      nrm = norm;
456
    }
457
 
458
    switch (ip->orthog_type) {
459
      case IP_CGS_ORTH:
460
        ierr = IPCGSBiOrthogonalization(ip,n,V,W,v,h,hnrm,nrm);CHKERRQ(ierr);
461
        break;
462
      default:
463
        SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown orthogonalization type");
464
    }
465
 
466
    if (allocated) { ierr = PetscFree(h);CHKERRQ(ierr); }
467
 
468
    ierr = PetscLogEventEnd(IP_Orthogonalize,ip,0,0,0);CHKERRQ(ierr);
469
  }
470
  PetscFunctionReturn(0);
471
}
472