Line data Source code
1 : /*
2 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3 : SLEPc - Scalable Library for Eigenvalue Problem Computations
4 : Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
5 :
6 : This file is part of SLEPc.
7 : SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9 : */
10 : /*
11 : SLEPc eigensolver: "davidson"
12 :
13 : Step: test for restarting, updateV, restartV
14 : */
15 :
16 : #include "davidson.h"
17 :
18 : typedef struct {
19 : PetscInt min_size_V; /* restart with this number of eigenvectors */
20 : PetscInt plusk; /* at restart, save plusk vectors from last iteration */
21 : PetscInt mpd; /* max size of the searching subspace */
22 : void *old_updateV_data; /* old updateV data */
23 : PetscErrorCode (*old_isRestarting)(dvdDashboard*,PetscBool*); /* old isRestarting */
24 : Mat oldU; /* previous projected right igenvectors */
25 : Mat oldV; /* previous projected left eigenvectors */
26 : PetscInt size_oldU; /* size of oldU */
27 : PetscBool allResiduals; /* if computing all the residuals */
28 : } dvdManagV_basic;
29 :
30 95 : static PetscErrorCode dvd_updateV_start(dvdDashboard *d)
31 : {
32 95 : dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
33 95 : PetscInt i;
34 :
35 95 : PetscFunctionBegin;
36 1920 : for (i=0;i<d->eps->ncv;i++) d->eigi[i] = 0.0;
37 95 : d->nR = d->real_nR;
38 1920 : for (i=0;i<d->eps->ncv;i++) d->nR[i] = 1.0;
39 95 : d->nX = d->real_nX;
40 1920 : for (i=0;i<d->eps->ncv;i++) d->errest[i] = 1.0;
41 95 : data->size_oldU = 0;
42 95 : d->nconv = 0;
43 95 : d->npreconv = 0;
44 95 : d->V_tra_s = d->V_tra_e = d->V_new_s = d->V_new_e = 0;
45 95 : d->size_D = 0;
46 95 : PetscFunctionReturn(PETSC_SUCCESS);
47 : }
48 :
49 7851 : static PetscErrorCode dvd_isrestarting_fullV(dvdDashboard *d,PetscBool *r)
50 : {
51 7851 : PetscInt l,k;
52 7851 : PetscBool restart;
53 7851 : dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
54 :
55 7851 : PetscFunctionBegin;
56 7851 : PetscCall(BVGetActiveColumns(d->eps->V,&l,&k));
57 7851 : restart = (k+2 > d->eps->ncv)? PETSC_TRUE: PETSC_FALSE;
58 :
59 : /* Check old isRestarting function */
60 7851 : if (PetscUnlikely(!restart && data->old_isRestarting)) PetscCall(data->old_isRestarting(d,&restart));
61 7851 : *r = restart;
62 7851 : PetscFunctionReturn(PETSC_SUCCESS);
63 : }
64 :
65 95 : static PetscErrorCode dvd_managementV_basic_d(dvdDashboard *d)
66 : {
67 95 : dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
68 :
69 95 : PetscFunctionBegin;
70 : /* Restore changes in dvdDashboard */
71 95 : d->updateV_data = data->old_updateV_data;
72 :
73 : /* Free local data */
74 95 : PetscCall(MatDestroy(&data->oldU));
75 95 : PetscCall(MatDestroy(&data->oldV));
76 95 : PetscCall(PetscFree(d->real_nR));
77 95 : PetscCall(PetscFree(d->real_nX));
78 95 : PetscCall(PetscFree(data));
79 95 : PetscFunctionReturn(PETSC_SUCCESS);
80 : }
81 :
82 308 : static PetscErrorCode dvd_updateV_conv_gen(dvdDashboard *d)
83 : {
84 308 : dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
85 308 : PetscInt npreconv,cMT,cMTX,lV,kV,nV;
86 308 : Mat Z,Z0,Q,Q0;
87 308 : PetscBool t;
88 : #if !defined(PETSC_USE_COMPLEX)
89 308 : PetscInt i;
90 : #endif
91 :
92 308 : PetscFunctionBegin;
93 308 : npreconv = d->npreconv;
94 : /* Constrains the converged pairs to nev */
95 : #if !defined(PETSC_USE_COMPLEX)
96 : /* Tries to maintain together conjugate eigenpairs */
97 616 : for (i=0; (i + (d->eigi[i]!=0.0?1:0) < npreconv) && (d->nconv + i < d->nev); i+= (d->eigi[i]!=0.0?2:1));
98 308 : npreconv = i;
99 : #else
100 : npreconv = PetscMax(PetscMin(d->nev-d->nconv,npreconv),0);
101 : #endif
102 : /* For GHEP without B-ortho, converge all of the requested pairs at once */
103 308 : PetscCall(PetscObjectTypeCompare((PetscObject)d->eps->ds,DSGHEP,&t));
104 308 : if (t && d->nconv+npreconv<d->nev) npreconv = 0;
105 : /* Quick exit */
106 308 : if (npreconv == 0) PetscFunctionReturn(PETSC_SUCCESS);
107 :
108 308 : PetscCall(BVGetActiveColumns(d->eps->V,&lV,&kV));
109 308 : nV = kV - lV;
110 308 : cMT = nV - npreconv;
111 : /* Harmonics restarts with right eigenvectors, and other with the left ones.
112 : If the problem is standard or hermitian, left and right vectors are the same */
113 308 : if (!(d->W||DVD_IS(d->sEP,DVD_EP_STD)||DVD_IS(d->sEP,DVD_EP_HERMITIAN))) {
114 : /* ps.Q <- [ps.Q(0:npreconv-1) ps.Z(npreconv:size_H-1)] */
115 0 : PetscCall(DSGetMat(d->eps->ds,DS_MAT_Q,&Q));
116 0 : PetscCall(DSGetMat(d->eps->ds,DS_MAT_Z,&Z));
117 0 : PetscCall(MatDenseGetSubMatrix(Q,0,npreconv,nV,npreconv+cMT,&Q0));
118 0 : PetscCall(MatDenseGetSubMatrix(Z,0,npreconv,nV,npreconv+cMT,&Z0));
119 0 : PetscCall(MatCopy(Z0,Q0,SAME_NONZERO_PATTERN));
120 0 : PetscCall(MatDenseRestoreSubMatrix(Q,&Q0));
121 0 : PetscCall(MatDenseRestoreSubMatrix(Z,&Z0));
122 0 : PetscCall(DSRestoreMat(d->eps->ds,DS_MAT_Q,&Q));
123 0 : PetscCall(DSRestoreMat(d->eps->ds,DS_MAT_Z,&Z));
124 : }
125 308 : if (DVD_IS(d->sEP,DVD_EP_INDEFINITE)) PetscCall(DSPseudoOrthogonalize(d->eps->ds,DS_MAT_Q,nV,d->nBds,&cMTX,d->nBds));
126 308 : else PetscCall(DSOrthogonalize(d->eps->ds,DS_MAT_Q,nV,&cMTX));
127 308 : cMT = cMTX - npreconv;
128 :
129 308 : if (d->W) {
130 55 : PetscCall(DSOrthogonalize(d->eps->ds,DS_MAT_Z,nV,&cMTX));
131 55 : cMT = PetscMin(cMT,cMTX - npreconv);
132 : }
133 :
134 : /* Lock the converged pairs */
135 308 : d->eigr+= npreconv;
136 : #if !defined(PETSC_USE_COMPLEX)
137 308 : if (d->eigi) d->eigi+= npreconv;
138 : #endif
139 308 : d->nconv+= npreconv;
140 308 : d->errest+= npreconv;
141 : /* Notify the changes in V and update the other subspaces */
142 308 : d->V_tra_s = npreconv; d->V_tra_e = nV;
143 308 : d->V_new_s = cMT; d->V_new_e = d->V_new_s;
144 : /* Remove oldU */
145 308 : data->size_oldU = 0;
146 :
147 308 : d->npreconv-= npreconv;
148 308 : PetscFunctionReturn(PETSC_SUCCESS);
149 : }
150 :
151 1129 : static PetscErrorCode dvd_updateV_restart_gen(dvdDashboard *d)
152 : {
153 1129 : dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
154 1129 : PetscInt lV,kV,nV,size_plusk,size_X,cMTX,cMTY,max_restart_size;
155 1129 : Mat Q,Q0,Z,Z0,U,V;
156 :
157 1129 : PetscFunctionBegin;
158 : /* Select size_X desired pairs from V */
159 : /* The restarted basis should:
160 : - have at least one spot to add a new direction;
161 : - keep converged vectors, npreconv;
162 : - keep at least 1 oldU direction if possible.
163 : */
164 1129 : PetscCall(BVGetActiveColumns(d->eps->V,&lV,&kV));
165 1129 : nV = kV - lV;
166 1129 : max_restart_size = PetscMax(0,PetscMin(d->eps->mpd - 1,d->eps->ncv - lV - 2));
167 1129 : size_X = PetscMin(PetscMin(data->min_size_V+d->npreconv,max_restart_size - (max_restart_size - d->npreconv > 1 && data->plusk > 0 && data->size_oldU > 0 ? 1 : 0)), nV);
168 :
169 : /* Add plusk eigenvectors from the previous iteration */
170 1129 : size_plusk = PetscMax(0,PetscMin(PetscMin(PetscMin(data->plusk,data->size_oldU),max_restart_size - size_X),nV - size_X));
171 :
172 1129 : d->size_MT = nV;
173 : /* ps.Q <- orth([pX(0:size_X-1) [oldU(0:size_plusk-1); 0] ]) */
174 : /* Harmonics restarts with right eigenvectors, and other with the left ones.
175 : If the problem is standard or hermitian, left and right vectors are the same */
176 1129 : if (!(d->W||DVD_IS(d->sEP,DVD_EP_STD)||DVD_IS(d->sEP,DVD_EP_HERMITIAN))) {
177 0 : PetscCall(DSGetMat(d->eps->ds,DS_MAT_Q,&Q));
178 0 : PetscCall(DSGetMat(d->eps->ds,DS_MAT_Z,&Z));
179 0 : PetscCall(MatDenseGetSubMatrix(Q,0,nV,0,size_X,&Q0));
180 0 : PetscCall(MatDenseGetSubMatrix(Z,0,nV,0,size_X,&Z0));
181 0 : PetscCall(MatCopy(Z0,Q0,SAME_NONZERO_PATTERN));
182 0 : PetscCall(MatDenseRestoreSubMatrix(Q,&Q0));
183 0 : PetscCall(MatDenseRestoreSubMatrix(Z,&Z0));
184 0 : PetscCall(DSRestoreMat(d->eps->ds,DS_MAT_Q,&Q));
185 0 : PetscCall(DSRestoreMat(d->eps->ds,DS_MAT_Z,&Z));
186 : }
187 1129 : PetscCheck(size_plusk<=0 || !DVD_IS(d->sEP,DVD_EP_INDEFINITE),PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported plusk>0 in indefinite eigenvalue problems");
188 1129 : if (size_plusk > 0) {
189 1112 : PetscCall(DSGetMat(d->eps->ds,DS_MAT_Q,&Q));
190 1112 : PetscCall(MatDenseGetSubMatrix(Q,0,nV,size_X,size_X+size_plusk,&Q0));
191 1112 : PetscCall(MatDenseGetSubMatrix(data->oldU,0,nV,0,size_plusk,&U));
192 1112 : PetscCall(MatCopy(U,Q0,SAME_NONZERO_PATTERN));
193 1112 : PetscCall(MatDenseRestoreSubMatrix(Q,&Q0));
194 1112 : PetscCall(MatDenseRestoreSubMatrix(data->oldU,&U));
195 1112 : PetscCall(DSRestoreMat(d->eps->ds,DS_MAT_Q,&Q));
196 : }
197 1129 : if (DVD_IS(d->sEP,DVD_EP_INDEFINITE)) PetscCall(DSPseudoOrthogonalize(d->eps->ds,DS_MAT_Q,size_X,d->nBds,&cMTX,d->nBds));
198 1129 : else PetscCall(DSOrthogonalize(d->eps->ds,DS_MAT_Q,size_X+size_plusk,&cMTX));
199 :
200 1129 : if (d->W && size_plusk > 0) {
201 : /* ps.Z <- orth([ps.Z(0:size_X-1) [oldV(0:size_plusk-1); 0] ]) */
202 122 : PetscCall(DSGetMat(d->eps->ds,DS_MAT_Z,&Z));
203 122 : PetscCall(MatDenseGetSubMatrix(Z,0,nV,size_X,size_X+size_plusk,&Z0));
204 122 : PetscCall(MatDenseGetSubMatrix(data->oldV,0,nV,0,size_plusk,&V));
205 122 : PetscCall(MatCopy(V,Z0,SAME_NONZERO_PATTERN));
206 122 : PetscCall(MatDenseRestoreSubMatrix(Z,&Z0));
207 122 : PetscCall(MatDenseRestoreSubMatrix(data->oldV,&V));
208 122 : PetscCall(DSRestoreMat(d->eps->ds,DS_MAT_Z,&Z));
209 122 : PetscCall(DSOrthogonalize(d->eps->ds,DS_MAT_Z,size_X+size_plusk,&cMTY));
210 122 : cMTX = PetscMin(cMTX, cMTY);
211 : }
212 1129 : PetscAssert(cMTX<=size_X+size_plusk,PETSC_COMM_SELF,PETSC_ERR_SUP,"Invalid number of columns to restart");
213 :
214 : /* Notify the changes in V and update the other subspaces */
215 1129 : d->V_tra_s = 0; d->V_tra_e = cMTX;
216 1129 : d->V_new_s = d->V_tra_e; d->V_new_e = d->V_new_s;
217 :
218 : /* Remove oldU */
219 1129 : data->size_oldU = 0;
220 :
221 : /* Remove npreconv */
222 1129 : d->npreconv = 0;
223 1129 : PetscFunctionReturn(PETSC_SUCCESS);
224 : }
225 :
226 13136 : static PetscErrorCode dvd_updateV_testConv(dvdDashboard *d,PetscInt s,PetscInt pre,PetscInt e,PetscInt *nConv)
227 : {
228 13136 : PetscInt i,j,b;
229 13136 : PetscReal norm;
230 13136 : PetscBool conv, c;
231 13136 : dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
232 :
233 13136 : PetscFunctionBegin;
234 13136 : if (nConv) *nConv = s;
235 22311 : for (i=s,conv=PETSC_TRUE;(conv || data->allResiduals) && (i < e);i+=b) {
236 : #if !defined(PETSC_USE_COMPLEX)
237 9175 : b = d->eigi[i]!=0.0?2:1;
238 : #else
239 : b = 1;
240 : #endif
241 9175 : if (i+b-1 >= pre) PetscCall(d->calcpairs_residual(d,i,i+b));
242 : /* Test the Schur vector */
243 18352 : for (j=0,c=PETSC_TRUE;j<b && c;j++) {
244 9177 : norm = d->nR[i+j]/d->nX[i+j];
245 9177 : c = d->testConv(d,d->eigr[i+j],d->eigi[i+j],norm,&d->errest[i+j]);
246 : }
247 9175 : if (conv && c) { if (nConv) *nConv = i+b; }
248 : else conv = PETSC_FALSE;
249 : }
250 13136 : pre = PetscMax(pre,i);
251 :
252 : #if !defined(PETSC_USE_COMPLEX)
253 : /* Enforce converged conjugate complex eigenpairs */
254 13136 : if (nConv) {
255 7030 : for (j=0;j<*nConv;j++) if (d->eigi[j] != 0.0) j++;
256 6722 : if (j>*nConv) (*nConv)--;
257 : }
258 : #endif
259 13282 : for (i=pre;i<e;i++) d->errest[i] = d->nR[i] = 1.0;
260 13136 : PetscFunctionReturn(PETSC_SUCCESS);
261 : }
262 :
263 6722 : static PetscErrorCode dvd_updateV_update_gen(dvdDashboard *d)
264 : {
265 6722 : dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
266 6722 : PetscInt size_D,s,lV,kV,nV;
267 6722 : Mat Q,Q0,Z,Z0,U,V;
268 :
269 6722 : PetscFunctionBegin;
270 : /* Select the desired pairs */
271 6722 : PetscCall(BVGetActiveColumns(d->eps->V,&lV,&kV));
272 6722 : nV = kV - lV;
273 6722 : size_D = PetscMin(PetscMin(PetscMin(d->bs,nV),d->eps->ncv-nV),nV);
274 6722 : if (size_D == 0) PetscFunctionReturn(PETSC_SUCCESS);
275 :
276 : /* Fill V with D */
277 6722 : PetscCall(d->improveX(d,d->npreconv,d->npreconv+size_D,&size_D));
278 :
279 : /* If D is empty, exit */
280 6722 : d->size_D = size_D;
281 6722 : if (size_D == 0) PetscFunctionReturn(PETSC_SUCCESS);
282 :
283 : /* Get the residual of all pairs */
284 : #if !defined(PETSC_USE_COMPLEX)
285 6414 : s = (d->eigi[0]!=0.0)? 2: 1;
286 : #else
287 : s = 1;
288 : #endif
289 6414 : PetscCall(BVGetActiveColumns(d->eps->V,&lV,&kV));
290 6414 : nV = kV - lV;
291 6414 : PetscCall(dvd_updateV_testConv(d,s,s,data->allResiduals?nV:size_D,NULL));
292 :
293 : /* Notify the changes in V */
294 6414 : d->V_tra_s = 0; d->V_tra_e = 0;
295 6414 : d->V_new_s = nV; d->V_new_e = nV+size_D;
296 :
297 : /* Save the projected eigenvectors */
298 6414 : if (data->plusk > 0) {
299 6370 : PetscCall(MatZeroEntries(data->oldU));
300 6370 : data->size_oldU = nV;
301 6370 : PetscCall(DSGetMat(d->eps->ds,DS_MAT_Q,&Q));
302 6370 : PetscCall(MatDenseGetSubMatrix(Q,0,nV,0,nV,&Q0));
303 6370 : PetscCall(MatDenseGetSubMatrix(data->oldU,0,nV,0,nV,&U));
304 6370 : PetscCall(MatCopy(Q0,U,SAME_NONZERO_PATTERN));
305 6370 : PetscCall(MatDenseRestoreSubMatrix(Q,&Q0));
306 6370 : PetscCall(MatDenseRestoreSubMatrix(data->oldU,&U));
307 6370 : PetscCall(DSRestoreMat(d->eps->ds,DS_MAT_Q,&Q));
308 6370 : if (d->W) {
309 874 : PetscCall(MatZeroEntries(data->oldV));
310 874 : PetscCall(DSGetMat(d->eps->ds,DS_MAT_Z,&Z));
311 874 : PetscCall(MatDenseGetSubMatrix(Z,0,nV,0,nV,&Z0));
312 874 : PetscCall(MatDenseGetSubMatrix(data->oldV,0,nV,0,nV,&V));
313 874 : PetscCall(MatCopy(Z0,V,SAME_NONZERO_PATTERN));
314 874 : PetscCall(MatDenseRestoreSubMatrix(Z,&Z0));
315 874 : PetscCall(MatDenseRestoreSubMatrix(data->oldV,&V));
316 874 : PetscCall(DSRestoreMat(d->eps->ds,DS_MAT_Z,&Z));
317 : }
318 : }
319 6414 : PetscFunctionReturn(PETSC_SUCCESS);
320 : }
321 :
322 7851 : static PetscErrorCode dvd_updateV_extrapol(dvdDashboard *d)
323 : {
324 7851 : dvdManagV_basic *data = (dvdManagV_basic*)d->updateV_data;
325 7851 : PetscInt i;
326 7851 : PetscBool restart,t;
327 :
328 7851 : PetscFunctionBegin;
329 : /* TODO: restrict select pairs to each case */
330 7851 : PetscCall(d->calcpairs_selectPairs(d, data->min_size_V+d->npreconv));
331 :
332 : /* If the subspaces doesn't need restart, add new vector */
333 7851 : PetscCall(d->isRestarting(d,&restart));
334 7851 : if (!restart) {
335 6722 : d->size_D = 0;
336 6722 : PetscCall(dvd_updateV_update_gen(d));
337 :
338 : /* If no vector were converged, exit */
339 : /* For GHEP without B-ortho, converge all of the requested pairs at once */
340 6722 : PetscCall(PetscObjectTypeCompare((PetscObject)d->eps->ds,DSGHEP,&t));
341 6722 : if (d->nconv+d->npreconv < d->nev && (t || d->npreconv == 0)) PetscFunctionReturn(PETSC_SUCCESS);
342 : }
343 :
344 : /* If some eigenpairs were converged, lock them */
345 1437 : if (d->npreconv > 0) {
346 308 : i = d->npreconv;
347 308 : PetscCall(dvd_updateV_conv_gen(d));
348 :
349 : /* If some eigenpair was locked, exit */
350 308 : if (i > d->npreconv) PetscFunctionReturn(PETSC_SUCCESS);
351 : }
352 :
353 : /* Else, a restarting is performed */
354 1129 : PetscCall(dvd_updateV_restart_gen(d));
355 1129 : PetscFunctionReturn(PETSC_SUCCESS);
356 : }
357 :
358 285 : PetscErrorCode dvd_managementV_basic(dvdDashboard *d,dvdBlackboard *b,PetscInt bs,PetscInt mpd,PetscInt min_size_V,PetscInt plusk,PetscBool harm,PetscBool allResiduals)
359 : {
360 285 : dvdManagV_basic *data;
361 : #if !defined(PETSC_USE_COMPLEX)
362 285 : PetscBool her_probl,std_probl;
363 : #endif
364 :
365 285 : PetscFunctionBegin;
366 : /* Setting configuration constrains */
367 : #if !defined(PETSC_USE_COMPLEX)
368 : /* if the last converged eigenvalue is complex its conjugate pair is also
369 : converged */
370 285 : her_probl = DVD_IS(d->sEP,DVD_EP_HERMITIAN)? PETSC_TRUE: PETSC_FALSE;
371 285 : std_probl = DVD_IS(d->sEP,DVD_EP_STD)? PETSC_TRUE: PETSC_FALSE;
372 285 : b->max_size_X = PetscMax(b->max_size_X,bs+((her_probl && std_probl)?0:1));
373 : #else
374 : b->max_size_X = PetscMax(b->max_size_X,bs);
375 : #endif
376 :
377 285 : b->max_size_V = PetscMax(b->max_size_V,mpd);
378 285 : min_size_V = PetscMin(min_size_V,mpd-bs);
379 285 : b->size_V = PetscMax(b->size_V,b->max_size_V+b->max_size_P+b->max_nev);
380 285 : b->max_size_oldX = plusk;
381 :
382 : /* Setup the step */
383 285 : if (b->state >= DVD_STATE_CONF) {
384 95 : PetscCall(PetscNew(&data));
385 95 : data->mpd = b->max_size_V;
386 95 : data->min_size_V = min_size_V;
387 95 : d->bs = bs;
388 95 : data->plusk = plusk;
389 95 : data->allResiduals = allResiduals;
390 :
391 95 : d->eigr = d->eps->eigr;
392 95 : d->eigi = d->eps->eigi;
393 95 : d->errest = d->eps->errest;
394 95 : PetscCall(PetscMalloc1(d->eps->ncv,&d->real_nR));
395 95 : PetscCall(PetscMalloc1(d->eps->ncv,&d->real_nX));
396 95 : if (plusk > 0) PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,d->eps->ncv,d->eps->ncv,NULL,&data->oldU));
397 2 : else data->oldU = NULL;
398 95 : if (harm && plusk>0) PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,d->eps->ncv,d->eps->ncv,NULL,&data->oldV));
399 74 : else data->oldV = NULL;
400 :
401 95 : data->old_updateV_data = d->updateV_data;
402 95 : d->updateV_data = data;
403 95 : data->old_isRestarting = d->isRestarting;
404 95 : d->isRestarting = dvd_isrestarting_fullV;
405 95 : d->updateV = dvd_updateV_extrapol;
406 95 : d->preTestConv = dvd_updateV_testConv;
407 95 : PetscCall(EPSDavidsonFLAdd(&d->startList,dvd_updateV_start));
408 95 : PetscCall(EPSDavidsonFLAdd(&d->destroyList,dvd_managementV_basic_d));
409 : }
410 285 : PetscFunctionReturn(PETSC_SUCCESS);
411 : }
|