Subversion Repositories slepc-dev

Rev

Details | Compare with Previous | Last modification | View Log | RSS feed

Rev Author Line No. Line
545 dsic.upv.es!jroman 1
/*
2
      EPS routines related to memory management.
1376 slepc 3
 
4
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1672 slepc 5
   SLEPc - Scalable Library for Eigenvalue Problem Computations
6
   Copyright (c) 2002-2009, Universidad Politecnica de Valencia, Spain
1376 slepc 7
 
1672 slepc 8
   This file is part of SLEPc.
9
 
10
   SLEPc is free software: you can redistribute it and/or modify it under  the
11
   terms of version 3 of the GNU Lesser General Public License as published by
12
   the Free Software Foundation.
13
 
14
   SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
15
   WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
16
   FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
17
   more details.
18
 
19
   You  should have received a copy of the GNU Lesser General  Public  License
20
   along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
1376 slepc 21
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
545 dsic.upv.es!jroman 22
*/
1376 slepc 23
 
1521 slepc 24
#include "private/epsimpl.h"   /*I "slepceps.h" I*/
524 dsic.upv.es!antodo 25
 
26
#undef __FUNCT__  
27
#define __FUNCT__ "EPSAllocateSolution"
545 dsic.upv.es!jroman 28
/*
29
  EPSAllocateSolution - Allocate memory storage for common variables such
1596 slepc 30
  as eigenvalues and eigenvectors. In this version, all
31
  vectors in V (and W) share a contiguous chunk of memory.
545 dsic.upv.es!jroman 32
*/
524 dsic.upv.es!antodo 33
PetscErrorCode EPSAllocateSolution(EPS eps)
34
{
35
  PetscErrorCode ierr;
1596 slepc 36
  PetscInt       i,nloc;
37
  PetscScalar    *pV,*pW;
524 dsic.upv.es!antodo 38
 
39
  PetscFunctionBegin;
40
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
41
  if (eps->allocated_ncv != eps->ncv) {
42
    if (eps->allocated_ncv > 0) {
43
      ierr = PetscFree(eps->eigr);CHKERRQ(ierr);
44
      ierr = PetscFree(eps->eigi);CHKERRQ(ierr);
45
      ierr = PetscFree(eps->errest);CHKERRQ(ierr);
780 dsic.upv.es!jroman 46
      ierr = PetscFree(eps->errest_left);CHKERRQ(ierr);
524 dsic.upv.es!antodo 47
      ierr = VecGetArray(eps->V[0],&pV);CHKERRQ(ierr);
48
      for (i=0;i<eps->allocated_ncv;i++) {
49
        ierr = VecDestroy(eps->V[i]);CHKERRQ(ierr);
50
      }
51
      ierr = PetscFree(pV);CHKERRQ(ierr);
52
      ierr = PetscFree(eps->V);CHKERRQ(ierr);
1596 slepc 53
      if (eps->AV) { ierr = VecDestroyVecs(eps->AV,eps->allocated_ncv);CHKERRQ(ierr); }
780 dsic.upv.es!jroman 54
      if (eps->solverclass == EPS_TWO_SIDE) {
55
        ierr = VecGetArray(eps->W[0],&pW);CHKERRQ(ierr);
56
        for (i=0;i<eps->allocated_ncv;i++) {
57
          ierr = VecDestroy(eps->W[i]);CHKERRQ(ierr);
58
        }
59
        ierr = PetscFree(pW);CHKERRQ(ierr);
60
        ierr = PetscFree(eps->W);CHKERRQ(ierr);
61
      }
524 dsic.upv.es!antodo 62
    }
63
    ierr = PetscMalloc(eps->ncv*sizeof(PetscScalar),&eps->eigr);CHKERRQ(ierr);
64
    ierr = PetscMalloc(eps->ncv*sizeof(PetscScalar),&eps->eigi);CHKERRQ(ierr);
65
    ierr = PetscMalloc(eps->ncv*sizeof(PetscReal),&eps->errest);CHKERRQ(ierr);
780 dsic.upv.es!jroman 66
    ierr = PetscMalloc(eps->ncv*sizeof(PetscReal),&eps->errest_left);CHKERRQ(ierr);
1385 slepc 67
    ierr = VecGetLocalSize(eps->vec_initial,&nloc);CHKERRQ(ierr);
524 dsic.upv.es!antodo 68
    ierr = PetscMalloc(eps->ncv*sizeof(Vec),&eps->V);CHKERRQ(ierr);
69
    ierr = PetscMalloc(eps->ncv*nloc*sizeof(PetscScalar),&pV);CHKERRQ(ierr);
70
    for (i=0;i<eps->ncv;i++) {
1422 slepc 71
      ierr = VecCreateMPIWithArray(((PetscObject)eps)->comm,nloc,PETSC_DECIDE,pV+i*nloc,&eps->V[i]);CHKERRQ(ierr);
524 dsic.upv.es!antodo 72
    }
1596 slepc 73
    eps->AV = PETSC_NULL;
780 dsic.upv.es!jroman 74
    if (eps->solverclass == EPS_TWO_SIDE) {
75
      ierr = PetscMalloc(eps->ncv*sizeof(Vec),&eps->W);CHKERRQ(ierr);
76
      ierr = PetscMalloc(eps->ncv*nloc*sizeof(PetscScalar),&pW);CHKERRQ(ierr);
77
      for (i=0;i<eps->ncv;i++) {
1422 slepc 78
        ierr = VecCreateMPIWithArray(((PetscObject)eps)->comm,nloc,PETSC_DECIDE,pW+i*nloc,&eps->W[i]);CHKERRQ(ierr);
780 dsic.upv.es!jroman 79
      }
80
    }
524 dsic.upv.es!antodo 81
    eps->allocated_ncv = eps->ncv;
82
  }
83
  PetscFunctionReturn(0);
84
}
85
 
86
#undef __FUNCT__  
1596 slepc 87
#define __FUNCT__ "EPSFreeSolution"
545 dsic.upv.es!jroman 88
/*
89
  EPSFreeSolution - Free memory storage. This routine is related to
1596 slepc 90
  EPSAllocateSolution().
545 dsic.upv.es!jroman 91
*/
1596 slepc 92
PetscErrorCode EPSFreeSolution(EPS eps)
524 dsic.upv.es!antodo 93
{
94
  PetscErrorCode ierr;
1509 slepc 95
  PetscInt       i;
780 dsic.upv.es!jroman 96
  PetscScalar    *pV,*pW;
524 dsic.upv.es!antodo 97
 
98
  PetscFunctionBegin;
99
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
100
  if (eps->allocated_ncv > 0) {
101
    ierr = PetscFree(eps->eigr);CHKERRQ(ierr);
102
    ierr = PetscFree(eps->eigi);CHKERRQ(ierr);
103
    ierr = PetscFree(eps->errest);CHKERRQ(ierr);
780 dsic.upv.es!jroman 104
    ierr = PetscFree(eps->errest_left);CHKERRQ(ierr);
524 dsic.upv.es!antodo 105
    ierr = VecGetArray(eps->V[0],&pV);CHKERRQ(ierr);
106
    for (i=0;i<eps->allocated_ncv;i++) {
107
      ierr = VecDestroy(eps->V[i]);CHKERRQ(ierr);
108
    }
109
    ierr = PetscFree(pV);CHKERRQ(ierr);
110
    ierr = PetscFree(eps->V);CHKERRQ(ierr);
1596 slepc 111
    if (eps->AV) { ierr = VecDestroyVecs(eps->AV,eps->allocated_ncv);CHKERRQ(ierr); }
780 dsic.upv.es!jroman 112
    if (eps->solverclass == EPS_TWO_SIDE) {
113
      ierr = VecGetArray(eps->W[0],&pW);CHKERRQ(ierr);
114
      for (i=0;i<eps->allocated_ncv;i++) {
115
        ierr = VecDestroy(eps->W[i]);CHKERRQ(ierr);
116
      }
117
      ierr = PetscFree(pW);CHKERRQ(ierr);
118
      ierr = PetscFree(eps->W);CHKERRQ(ierr);
119
    }
575 dsic.upv.es!antodo 120
    eps->allocated_ncv = 0;
524 dsic.upv.es!antodo 121
  }
122
  PetscFunctionReturn(0);
123
}