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