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
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5
      SLEPc - Scalable Library for Eigenvalue Problem Computations
6
      Copyright (c) 2002-2007, Universidad Politecnica de Valencia, Spain
7
 
8
      This file is part of SLEPc. See the README file for conditions of use
9
      and additional information.
10
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
545 dsic.upv.es!jroman 11
*/
1376 slepc 12
 
524 dsic.upv.es!antodo 13
#include "src/eps/epsimpl.h"   /*I "slepceps.h" I*/
14
 
15
#undef __FUNCT__  
16
#define __FUNCT__ "EPSAllocateSolution"
545 dsic.upv.es!jroman 17
/*
18
  EPSAllocateSolution - Allocate memory storage for common variables such
19
  as eigenvalues and eigenvectors.
20
*/
524 dsic.upv.es!antodo 21
PetscErrorCode EPSAllocateSolution(EPS eps)
22
{
23
  PetscErrorCode ierr;
24
 
25
  PetscFunctionBegin;
26
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
27
  if (eps->allocated_ncv != eps->ncv) {
28
    if (eps->allocated_ncv > 0) {
29
      ierr = PetscFree(eps->eigr);CHKERRQ(ierr);
30
      ierr = PetscFree(eps->eigi);CHKERRQ(ierr);
31
      ierr = PetscFree(eps->errest);CHKERRQ(ierr);
780 dsic.upv.es!jroman 32
      ierr = PetscFree(eps->errest_left);CHKERRQ(ierr);
524 dsic.upv.es!antodo 33
      ierr = VecDestroyVecs(eps->V,eps->allocated_ncv);CHKERRQ(ierr);
34
      ierr = VecDestroyVecs(eps->AV,eps->allocated_ncv);CHKERRQ(ierr);
780 dsic.upv.es!jroman 35
      if (eps->solverclass == EPS_TWO_SIDE) {
36
        ierr = VecDestroyVecs(eps->W,eps->allocated_ncv);CHKERRQ(ierr);
37
        ierr = VecDestroyVecs(eps->AW,eps->allocated_ncv);CHKERRQ(ierr);
38
      }
524 dsic.upv.es!antodo 39
    }
40
    ierr = PetscMalloc(eps->ncv*sizeof(PetscScalar),&eps->eigr);CHKERRQ(ierr);
41
    ierr = PetscMalloc(eps->ncv*sizeof(PetscScalar),&eps->eigi);CHKERRQ(ierr);
42
    ierr = PetscMalloc(eps->ncv*sizeof(PetscReal),&eps->errest);CHKERRQ(ierr);
780 dsic.upv.es!jroman 43
    ierr = PetscMalloc(eps->ncv*sizeof(PetscReal),&eps->errest_left);CHKERRQ(ierr);
1385 slepc 44
    ierr = VecDuplicateVecs(eps->vec_initial,eps->ncv,&eps->V);CHKERRQ(ierr);
45
    ierr = VecDuplicateVecs(eps->vec_initial,eps->ncv,&eps->AV);CHKERRQ(ierr);
780 dsic.upv.es!jroman 46
    if (eps->solverclass == EPS_TWO_SIDE) {
1385 slepc 47
      ierr = VecDuplicateVecs(eps->vec_initial,eps->ncv,&eps->W);CHKERRQ(ierr);
48
      ierr = VecDuplicateVecs(eps->vec_initial,eps->ncv,&eps->AW);CHKERRQ(ierr);
780 dsic.upv.es!jroman 49
    }
524 dsic.upv.es!antodo 50
    eps->allocated_ncv = eps->ncv;
51
  }
52
  PetscFunctionReturn(0);
53
}
54
 
55
#undef __FUNCT__  
56
#define __FUNCT__ "EPSFreeSolution"
545 dsic.upv.es!jroman 57
/*
58
  EPSFreeSolution - Free memory storage. This routine is related to
59
  EPSAllocateSolution().
60
*/
524 dsic.upv.es!antodo 61
PetscErrorCode EPSFreeSolution(EPS eps)
62
{
63
  PetscErrorCode ierr;
64
 
65
  PetscFunctionBegin;
66
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
67
  if (eps->allocated_ncv > 0) {
68
    ierr = PetscFree(eps->eigr);CHKERRQ(ierr);
69
    ierr = PetscFree(eps->eigi);CHKERRQ(ierr);
70
    ierr = PetscFree(eps->errest);CHKERRQ(ierr);
780 dsic.upv.es!jroman 71
    ierr = PetscFree(eps->errest_left);CHKERRQ(ierr);
524 dsic.upv.es!antodo 72
    ierr = VecDestroyVecs(eps->V,eps->allocated_ncv);CHKERRQ(ierr);
73
    ierr = VecDestroyVecs(eps->AV,eps->allocated_ncv);CHKERRQ(ierr);
780 dsic.upv.es!jroman 74
    if (eps->solverclass == EPS_TWO_SIDE) {
75
      ierr = VecDestroyVecs(eps->W,eps->allocated_ncv);CHKERRQ(ierr);
76
      ierr = VecDestroyVecs(eps->AW,eps->allocated_ncv);CHKERRQ(ierr);
77
    }
575 dsic.upv.es!antodo 78
    eps->allocated_ncv = 0;
524 dsic.upv.es!antodo 79
  }
80
  PetscFunctionReturn(0);
81
}
82
 
83
#undef __FUNCT__  
84
#define __FUNCT__ "EPSAllocateSolutionContiguous"
545 dsic.upv.es!jroman 85
/*
86
  EPSAllocateSolutionContiguous - Allocate memory storage for common
87
  variables such as eigenvalues and eigenvectors. In this version, all
88
  vectors in V (and AV) share a contiguous chunk of memory. This is
89
  necessary for external packages such as Arpack.
90
*/
524 dsic.upv.es!antodo 91
PetscErrorCode EPSAllocateSolutionContiguous(EPS eps)
92
{
93
  PetscErrorCode ierr;
982 slepc 94
  int            i;
95
  PetscInt       nloc;
780 dsic.upv.es!jroman 96
  PetscScalar    *pV,*pW;
524 dsic.upv.es!antodo 97
 
98
  PetscFunctionBegin;
99
  if (eps->allocated_ncv != eps->ncv) {
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);
111
      ierr = VecGetArray(eps->AV[0],&pV);CHKERRQ(ierr);
112
      for (i=0;i<eps->allocated_ncv;i++) {
113
        ierr = VecDestroy(eps->AV[i]);CHKERRQ(ierr);
114
      }
115
      ierr = PetscFree(pV);CHKERRQ(ierr);
116
      ierr = PetscFree(eps->AV);CHKERRQ(ierr);
780 dsic.upv.es!jroman 117
      if (eps->solverclass == EPS_TWO_SIDE) {
118
        ierr = VecGetArray(eps->W[0],&pW);CHKERRQ(ierr);
119
        for (i=0;i<eps->allocated_ncv;i++) {
120
          ierr = VecDestroy(eps->W[i]);CHKERRQ(ierr);
121
        }
122
        ierr = PetscFree(pW);CHKERRQ(ierr);
123
        ierr = PetscFree(eps->W);CHKERRQ(ierr);
124
        ierr = VecGetArray(eps->AW[0],&pW);CHKERRQ(ierr);
125
        for (i=0;i<eps->allocated_ncv;i++) {
126
          ierr = VecDestroy(eps->AW[i]);CHKERRQ(ierr);
127
        }
128
        ierr = PetscFree(pW);CHKERRQ(ierr);
129
        ierr = PetscFree(eps->AW);CHKERRQ(ierr);
130
      }
524 dsic.upv.es!antodo 131
    }
132
    ierr = PetscMalloc(eps->ncv*sizeof(PetscScalar),&eps->eigr);CHKERRQ(ierr);
133
    ierr = PetscMalloc(eps->ncv*sizeof(PetscScalar),&eps->eigi);CHKERRQ(ierr);
134
    ierr = PetscMalloc(eps->ncv*sizeof(PetscReal),&eps->errest);CHKERRQ(ierr);
780 dsic.upv.es!jroman 135
    ierr = PetscMalloc(eps->ncv*sizeof(PetscReal),&eps->errest_left);CHKERRQ(ierr);
1385 slepc 136
    ierr = VecGetLocalSize(eps->vec_initial,&nloc);CHKERRQ(ierr);
524 dsic.upv.es!antodo 137
    ierr = PetscMalloc(eps->ncv*sizeof(Vec),&eps->V);CHKERRQ(ierr);
138
    ierr = PetscMalloc(eps->ncv*nloc*sizeof(PetscScalar),&pV);CHKERRQ(ierr);
139
    for (i=0;i<eps->ncv;i++) {
1422 slepc 140
      ierr = VecCreateMPIWithArray(((PetscObject)eps)->comm,nloc,PETSC_DECIDE,pV+i*nloc,&eps->V[i]);CHKERRQ(ierr);
524 dsic.upv.es!antodo 141
    }
142
    ierr = PetscMalloc(eps->ncv*sizeof(Vec),&eps->AV);CHKERRQ(ierr);
143
    ierr = PetscMalloc(eps->ncv*nloc*sizeof(PetscScalar),&pV);CHKERRQ(ierr);
144
    for (i=0;i<eps->ncv;i++) {
1422 slepc 145
      ierr = VecCreateMPIWithArray(((PetscObject)eps)->comm,nloc,PETSC_DECIDE,pV+i*nloc,&eps->AV[i]);CHKERRQ(ierr);
524 dsic.upv.es!antodo 146
    }
780 dsic.upv.es!jroman 147
    if (eps->solverclass == EPS_TWO_SIDE) {
148
      ierr = PetscMalloc(eps->ncv*sizeof(Vec),&eps->W);CHKERRQ(ierr);
149
      ierr = PetscMalloc(eps->ncv*nloc*sizeof(PetscScalar),&pW);CHKERRQ(ierr);
150
      for (i=0;i<eps->ncv;i++) {
1422 slepc 151
        ierr = VecCreateMPIWithArray(((PetscObject)eps)->comm,nloc,PETSC_DECIDE,pW+i*nloc,&eps->W[i]);CHKERRQ(ierr);
780 dsic.upv.es!jroman 152
      }
153
      ierr = PetscMalloc(eps->ncv*sizeof(Vec),&eps->AW);CHKERRQ(ierr);
154
      ierr = PetscMalloc(eps->ncv*nloc*sizeof(PetscScalar),&pW);CHKERRQ(ierr);
155
      for (i=0;i<eps->ncv;i++) {
1422 slepc 156
        ierr = VecCreateMPIWithArray(((PetscObject)eps)->comm,nloc,PETSC_DECIDE,pW+i*nloc,&eps->AW[i]);CHKERRQ(ierr);
780 dsic.upv.es!jroman 157
      }
158
    }
524 dsic.upv.es!antodo 159
    eps->allocated_ncv = eps->ncv;
160
  }
161
  PetscFunctionReturn(0);
162
}
163
 
164
#undef __FUNCT__  
165
#define __FUNCT__ "EPSFreeSolutionContiguous"
545 dsic.upv.es!jroman 166
/*
167
  EPSFreeSolution - Free memory storage. This routine is related to
168
  EPSAllocateSolutionContiguous().
169
*/
524 dsic.upv.es!antodo 170
PetscErrorCode EPSFreeSolutionContiguous(EPS eps)
171
{
172
  PetscErrorCode ierr;
173
  int            i;
780 dsic.upv.es!jroman 174
  PetscScalar    *pV,*pW;
524 dsic.upv.es!antodo 175
 
176
  PetscFunctionBegin;
177
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
178
  if (eps->allocated_ncv > 0) {
179
    ierr = PetscFree(eps->eigr);CHKERRQ(ierr);
180
    ierr = PetscFree(eps->eigi);CHKERRQ(ierr);
181
    ierr = PetscFree(eps->errest);CHKERRQ(ierr);
780 dsic.upv.es!jroman 182
    ierr = PetscFree(eps->errest_left);CHKERRQ(ierr);
524 dsic.upv.es!antodo 183
    ierr = VecGetArray(eps->V[0],&pV);CHKERRQ(ierr);
184
    for (i=0;i<eps->allocated_ncv;i++) {
185
      ierr = VecDestroy(eps->V[i]);CHKERRQ(ierr);
186
    }
187
    ierr = PetscFree(pV);CHKERRQ(ierr);
188
    ierr = PetscFree(eps->V);CHKERRQ(ierr);
189
    ierr = VecGetArray(eps->AV[0],&pV);CHKERRQ(ierr);
190
    for (i=0;i<eps->allocated_ncv;i++) {
191
      ierr = VecDestroy(eps->AV[i]);CHKERRQ(ierr);
192
    }
193
    ierr = PetscFree(pV);CHKERRQ(ierr);
194
    ierr = PetscFree(eps->AV);CHKERRQ(ierr);
780 dsic.upv.es!jroman 195
    if (eps->solverclass == EPS_TWO_SIDE) {
196
      ierr = VecGetArray(eps->W[0],&pW);CHKERRQ(ierr);
197
      for (i=0;i<eps->allocated_ncv;i++) {
198
        ierr = VecDestroy(eps->W[i]);CHKERRQ(ierr);
199
      }
200
      ierr = PetscFree(pW);CHKERRQ(ierr);
201
      ierr = PetscFree(eps->W);CHKERRQ(ierr);
202
      ierr = VecGetArray(eps->AW[0],&pW);CHKERRQ(ierr);
203
      for (i=0;i<eps->allocated_ncv;i++) {
204
        ierr = VecDestroy(eps->AW[i]);CHKERRQ(ierr);
205
      }
206
      ierr = PetscFree(pW);CHKERRQ(ierr);
207
      ierr = PetscFree(eps->AW);CHKERRQ(ierr);
208
    }
575 dsic.upv.es!antodo 209
    eps->allocated_ncv = 0;
524 dsic.upv.es!antodo 210
  }
211
  PetscFunctionReturn(0);
212
}