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 singular value solver: "cyclic"
12 :
13 : Method: Uses a Hermitian eigensolver for H(A) = [ 0 A ; A^T 0 ]
14 : */
15 :
16 : #include <slepc/private/svdimpl.h> /*I "slepcsvd.h" I*/
17 : #include <slepc/private/bvimpl.h>
18 : #include "cyclic.h"
19 :
20 40873 : static PetscErrorCode MatMult_Cyclic(Mat B,Vec x,Vec y)
21 : {
22 40873 : SVD_CYCLIC_SHELL *ctx;
23 40873 : const PetscScalar *px;
24 40873 : PetscScalar *py;
25 40873 : PetscInt m;
26 :
27 40873 : PetscFunctionBegin;
28 40873 : PetscCall(MatShellGetContext(B,&ctx));
29 40873 : PetscCall(MatGetLocalSize(ctx->A,&m,NULL));
30 40873 : PetscCall(VecGetArrayRead(x,&px));
31 40873 : PetscCall(VecGetArrayWrite(y,&py));
32 40873 : PetscCall(VecPlaceArray(ctx->x1,px));
33 40873 : PetscCall(VecPlaceArray(ctx->x2,px+m));
34 40873 : PetscCall(VecPlaceArray(ctx->y1,py));
35 40873 : PetscCall(VecPlaceArray(ctx->y2,py+m));
36 40873 : PetscCall(MatMult(ctx->A,ctx->x2,ctx->y1));
37 40873 : PetscCall(MatMult(ctx->AT,ctx->x1,ctx->y2));
38 40873 : PetscCall(VecResetArray(ctx->x1));
39 40873 : PetscCall(VecResetArray(ctx->x2));
40 40873 : PetscCall(VecResetArray(ctx->y1));
41 40873 : PetscCall(VecResetArray(ctx->y2));
42 40873 : PetscCall(VecRestoreArrayRead(x,&px));
43 40873 : PetscCall(VecRestoreArrayWrite(y,&py));
44 40873 : PetscFunctionReturn(PETSC_SUCCESS);
45 : }
46 :
47 0 : static PetscErrorCode MatGetDiagonal_Cyclic(Mat B,Vec diag)
48 : {
49 0 : PetscFunctionBegin;
50 0 : PetscCall(VecSet(diag,0.0));
51 0 : PetscFunctionReturn(PETSC_SUCCESS);
52 : }
53 :
54 29 : static PetscErrorCode MatDestroy_Cyclic(Mat B)
55 : {
56 29 : SVD_CYCLIC_SHELL *ctx;
57 :
58 29 : PetscFunctionBegin;
59 29 : PetscCall(MatShellGetContext(B,&ctx));
60 29 : PetscCall(VecDestroy(&ctx->x1));
61 29 : PetscCall(VecDestroy(&ctx->x2));
62 29 : PetscCall(VecDestroy(&ctx->y1));
63 29 : PetscCall(VecDestroy(&ctx->y2));
64 29 : if (ctx->misaligned) {
65 0 : PetscCall(VecDestroy(&ctx->wx2));
66 0 : PetscCall(VecDestroy(&ctx->wy2));
67 : }
68 29 : PetscCall(PetscFree(ctx));
69 29 : PetscFunctionReturn(PETSC_SUCCESS);
70 : }
71 :
72 : /*
73 : Builds cyclic matrix C = | 0 A |
74 : | AT 0 |
75 : */
76 46 : static PetscErrorCode SVDCyclicGetCyclicMat(SVD svd,Mat A,Mat AT,Mat *C)
77 : {
78 46 : SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
79 46 : SVD_CYCLIC_SHELL *ctx;
80 46 : PetscInt i,M,N,m,n,Istart,Iend;
81 46 : VecType vtype;
82 46 : Mat Zm,Zn;
83 : #if defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP)
84 : PetscBool gpu;
85 : const PetscInt *ranges;
86 : PetscMPIInt size;
87 : #endif
88 :
89 46 : PetscFunctionBegin;
90 46 : PetscCall(MatGetSize(A,&M,&N));
91 46 : PetscCall(MatGetLocalSize(A,&m,&n));
92 :
93 46 : if (cyclic->explicitmatrix) {
94 17 : PetscCheck(svd->expltrans,PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Cannot use explicit cyclic matrix with implicit transpose");
95 17 : PetscCall(MatCreate(PetscObjectComm((PetscObject)svd),&Zm));
96 17 : PetscCall(MatSetSizes(Zm,m,m,M,M));
97 17 : PetscCall(MatSetFromOptions(Zm));
98 17 : PetscCall(MatGetOwnershipRange(Zm,&Istart,&Iend));
99 2794 : for (i=Istart;i<Iend;i++) PetscCall(MatSetValue(Zm,i,i,0.0,INSERT_VALUES));
100 17 : PetscCall(MatAssemblyBegin(Zm,MAT_FINAL_ASSEMBLY));
101 17 : PetscCall(MatAssemblyEnd(Zm,MAT_FINAL_ASSEMBLY));
102 17 : PetscCall(MatCreate(PetscObjectComm((PetscObject)svd),&Zn));
103 17 : PetscCall(MatSetSizes(Zn,n,n,N,N));
104 17 : PetscCall(MatSetFromOptions(Zn));
105 17 : PetscCall(MatGetOwnershipRange(Zn,&Istart,&Iend));
106 1360 : for (i=Istart;i<Iend;i++) PetscCall(MatSetValue(Zn,i,i,0.0,INSERT_VALUES));
107 17 : PetscCall(MatAssemblyBegin(Zn,MAT_FINAL_ASSEMBLY));
108 17 : PetscCall(MatAssemblyEnd(Zn,MAT_FINAL_ASSEMBLY));
109 17 : PetscCall(MatCreateTile(1.0,Zm,1.0,A,1.0,AT,1.0,Zn,C));
110 17 : PetscCall(MatDestroy(&Zm));
111 17 : PetscCall(MatDestroy(&Zn));
112 : } else {
113 29 : PetscCall(PetscNew(&ctx));
114 29 : ctx->A = A;
115 29 : ctx->AT = AT;
116 29 : ctx->swapped = svd->swapped;
117 29 : PetscCall(MatCreateVecsEmpty(A,&ctx->x2,&ctx->x1));
118 29 : PetscCall(MatCreateVecsEmpty(A,&ctx->y2,&ctx->y1));
119 29 : PetscCall(MatCreateShell(PetscObjectComm((PetscObject)svd),m+n,m+n,M+N,M+N,ctx,C));
120 29 : PetscCall(MatShellSetOperation(*C,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_Cyclic));
121 29 : PetscCall(MatShellSetOperation(*C,MATOP_DESTROY,(void(*)(void))MatDestroy_Cyclic));
122 : #if defined(PETSC_HAVE_CUDA)
123 : PetscCall(PetscObjectTypeCompareAny((PetscObject)(svd->swapped?AT:A),&gpu,MATSEQAIJCUSPARSE,MATMPIAIJCUSPARSE,""));
124 : if (gpu) PetscCall(MatShellSetOperation(*C,MATOP_MULT,(void(*)(void))MatMult_Cyclic_CUDA));
125 : else
126 : #elif defined(PETSC_HAVE_HIP)
127 : PetscCall(PetscObjectTypeCompareAny((PetscObject)(svd->swapped?AT:A),&gpu,MATSEQAIJHIPSPARSE,MATMPIAIJHIPSPARSE,""));
128 : if (gpu) PetscCall(MatShellSetOperation(*C,MATOP_MULT,(void(*)(void))MatMult_Cyclic_HIP));
129 : else
130 : #endif
131 29 : PetscCall(MatShellSetOperation(*C,MATOP_MULT,(void(*)(void))MatMult_Cyclic));
132 29 : PetscCall(MatGetVecType(A,&vtype));
133 29 : PetscCall(MatSetVecType(*C,vtype));
134 : #if defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP)
135 : if (gpu) {
136 : /* check alignment of bottom block */
137 : PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)ctx->x1),&size));
138 : PetscCall(VecGetOwnershipRanges(ctx->x1,&ranges));
139 : for (i=0;i<size;i++) {
140 : ctx->misaligned = (((ranges[i+1]-ranges[i])*sizeof(PetscScalar))%16)? PETSC_TRUE: PETSC_FALSE;
141 : if (ctx->misaligned) break;
142 : }
143 : if (ctx->misaligned) { /* create work vectors for MatMult */
144 : PetscCall(VecDuplicate(ctx->x2,&ctx->wx2));
145 : PetscCall(VecDuplicate(ctx->y2,&ctx->wy2));
146 : }
147 : }
148 : #endif
149 : }
150 46 : PetscFunctionReturn(PETSC_SUCCESS);
151 : }
152 :
153 7649 : static PetscErrorCode MatMult_ECross(Mat B,Vec x,Vec y)
154 : {
155 7649 : SVD_CYCLIC_SHELL *ctx;
156 7649 : const PetscScalar *px;
157 7649 : PetscScalar *py;
158 7649 : PetscInt mn,m,n;
159 :
160 7649 : PetscFunctionBegin;
161 7649 : PetscCall(MatShellGetContext(B,&ctx));
162 7649 : PetscCall(MatGetLocalSize(ctx->A,NULL,&n));
163 7649 : PetscCall(VecGetLocalSize(y,&mn));
164 7649 : m = mn-n;
165 7649 : PetscCall(VecGetArrayRead(x,&px));
166 7649 : PetscCall(VecGetArrayWrite(y,&py));
167 7649 : PetscCall(VecPlaceArray(ctx->x1,px));
168 7649 : PetscCall(VecPlaceArray(ctx->x2,px+m));
169 7649 : PetscCall(VecPlaceArray(ctx->y1,py));
170 7649 : PetscCall(VecPlaceArray(ctx->y2,py+m));
171 7649 : PetscCall(VecCopy(ctx->x1,ctx->y1));
172 7649 : PetscCall(MatMult(ctx->A,ctx->x2,ctx->w));
173 7649 : PetscCall(MatMult(ctx->AT,ctx->w,ctx->y2));
174 7649 : PetscCall(VecResetArray(ctx->x1));
175 7649 : PetscCall(VecResetArray(ctx->x2));
176 7649 : PetscCall(VecResetArray(ctx->y1));
177 7649 : PetscCall(VecResetArray(ctx->y2));
178 7649 : PetscCall(VecRestoreArrayRead(x,&px));
179 7649 : PetscCall(VecRestoreArrayWrite(y,&py));
180 7649 : PetscFunctionReturn(PETSC_SUCCESS);
181 : }
182 :
183 6 : static PetscErrorCode MatGetDiagonal_ECross(Mat B,Vec d)
184 : {
185 6 : SVD_CYCLIC_SHELL *ctx;
186 6 : PetscScalar *pd;
187 6 : PetscMPIInt len;
188 6 : PetscInt mn,m,n,N,i,j,start,end,ncols;
189 6 : PetscScalar *work1,*work2,*diag;
190 6 : const PetscInt *cols;
191 6 : const PetscScalar *vals;
192 :
193 6 : PetscFunctionBegin;
194 6 : PetscCall(MatShellGetContext(B,&ctx));
195 6 : PetscCall(MatGetLocalSize(ctx->A,NULL,&n));
196 6 : PetscCall(VecGetLocalSize(d,&mn));
197 6 : m = mn-n;
198 6 : PetscCall(VecGetArrayWrite(d,&pd));
199 6 : PetscCall(VecPlaceArray(ctx->y1,pd));
200 6 : PetscCall(VecSet(ctx->y1,1.0));
201 6 : PetscCall(VecResetArray(ctx->y1));
202 6 : PetscCall(VecPlaceArray(ctx->y2,pd+m));
203 6 : if (!ctx->diag) {
204 : /* compute diagonal from rows and store in ctx->diag */
205 6 : PetscCall(VecDuplicate(ctx->y2,&ctx->diag));
206 6 : PetscCall(MatGetSize(ctx->A,NULL,&N));
207 6 : PetscCall(PetscCalloc2(N,&work1,N,&work2));
208 6 : if (ctx->swapped) {
209 0 : PetscCall(MatGetOwnershipRange(ctx->AT,&start,&end));
210 0 : for (i=start;i<end;i++) {
211 0 : PetscCall(MatGetRow(ctx->AT,i,&ncols,NULL,&vals));
212 0 : for (j=0;j<ncols;j++) work1[i] += vals[j]*vals[j];
213 0 : PetscCall(MatRestoreRow(ctx->AT,i,&ncols,NULL,&vals));
214 : }
215 : } else {
216 6 : PetscCall(MatGetOwnershipRange(ctx->A,&start,&end));
217 137 : for (i=start;i<end;i++) {
218 131 : PetscCall(MatGetRow(ctx->A,i,&ncols,&cols,&vals));
219 894 : for (j=0;j<ncols;j++) work1[cols[j]] += vals[j]*vals[j];
220 131 : PetscCall(MatRestoreRow(ctx->A,i,&ncols,&cols,&vals));
221 : }
222 : }
223 6 : PetscCall(PetscMPIIntCast(N,&len));
224 18 : PetscCallMPI(MPIU_Allreduce(work1,work2,len,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)B)));
225 6 : PetscCall(VecGetOwnershipRange(ctx->diag,&start,&end));
226 6 : PetscCall(VecGetArrayWrite(ctx->diag,&diag));
227 121 : for (i=start;i<end;i++) diag[i-start] = work2[i];
228 6 : PetscCall(VecRestoreArrayWrite(ctx->diag,&diag));
229 6 : PetscCall(PetscFree2(work1,work2));
230 : }
231 6 : PetscCall(VecCopy(ctx->diag,ctx->y2));
232 6 : PetscCall(VecResetArray(ctx->y2));
233 6 : PetscCall(VecRestoreArrayWrite(d,&pd));
234 6 : PetscFunctionReturn(PETSC_SUCCESS);
235 : }
236 :
237 6 : static PetscErrorCode MatDestroy_ECross(Mat B)
238 : {
239 6 : SVD_CYCLIC_SHELL *ctx;
240 :
241 6 : PetscFunctionBegin;
242 6 : PetscCall(MatShellGetContext(B,&ctx));
243 6 : PetscCall(VecDestroy(&ctx->x1));
244 6 : PetscCall(VecDestroy(&ctx->x2));
245 6 : PetscCall(VecDestroy(&ctx->y1));
246 6 : PetscCall(VecDestroy(&ctx->y2));
247 6 : PetscCall(VecDestroy(&ctx->diag));
248 6 : PetscCall(VecDestroy(&ctx->w));
249 6 : if (ctx->misaligned) {
250 0 : PetscCall(VecDestroy(&ctx->wx2));
251 0 : PetscCall(VecDestroy(&ctx->wy2));
252 : }
253 6 : PetscCall(PetscFree(ctx));
254 6 : PetscFunctionReturn(PETSC_SUCCESS);
255 : }
256 :
257 : /*
258 : Builds extended cross product matrix C = | I_m 0 |
259 : | 0 AT*A |
260 : t is an auxiliary Vec used to take the dimensions of the upper block
261 : */
262 13 : static PetscErrorCode SVDCyclicGetECrossMat(SVD svd,Mat A,Mat AT,Mat *C,Vec t)
263 : {
264 13 : SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
265 13 : SVD_CYCLIC_SHELL *ctx;
266 13 : PetscInt i,M,N,m,n,Istart,Iend;
267 13 : VecType vtype;
268 13 : Mat Id,Zm,Zn,ATA;
269 : #if defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP)
270 : PetscBool gpu;
271 : const PetscInt *ranges;
272 : PetscMPIInt size;
273 : #endif
274 :
275 13 : PetscFunctionBegin;
276 13 : PetscCall(MatGetSize(A,NULL,&N));
277 13 : PetscCall(MatGetLocalSize(A,NULL,&n));
278 13 : PetscCall(VecGetSize(t,&M));
279 13 : PetscCall(VecGetLocalSize(t,&m));
280 :
281 13 : if (cyclic->explicitmatrix) {
282 7 : PetscCheck(svd->expltrans,PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Cannot use explicit cyclic matrix with implicit transpose");
283 7 : PetscCall(MatCreateConstantDiagonal(PetscObjectComm((PetscObject)svd),m,m,M,M,1.0,&Id));
284 7 : PetscCall(MatCreate(PetscObjectComm((PetscObject)svd),&Zm));
285 7 : PetscCall(MatSetSizes(Zm,m,n,M,N));
286 7 : PetscCall(MatSetFromOptions(Zm));
287 7 : PetscCall(MatGetOwnershipRange(Zm,&Istart,&Iend));
288 355 : for (i=Istart;i<Iend;i++) {
289 348 : if (i<N) PetscCall(MatSetValue(Zm,i,i,0.0,INSERT_VALUES));
290 : }
291 7 : PetscCall(MatAssemblyBegin(Zm,MAT_FINAL_ASSEMBLY));
292 7 : PetscCall(MatAssemblyEnd(Zm,MAT_FINAL_ASSEMBLY));
293 7 : PetscCall(MatCreate(PetscObjectComm((PetscObject)svd),&Zn));
294 7 : PetscCall(MatSetSizes(Zn,n,m,N,M));
295 7 : PetscCall(MatSetFromOptions(Zn));
296 7 : PetscCall(MatGetOwnershipRange(Zn,&Istart,&Iend));
297 364 : for (i=Istart;i<Iend;i++) {
298 357 : if (i<m) PetscCall(MatSetValue(Zn,i,i,0.0,INSERT_VALUES));
299 : }
300 7 : PetscCall(MatAssemblyBegin(Zn,MAT_FINAL_ASSEMBLY));
301 7 : PetscCall(MatAssemblyEnd(Zn,MAT_FINAL_ASSEMBLY));
302 7 : PetscCall(MatProductCreate(AT,A,NULL,&ATA));
303 7 : PetscCall(MatProductSetType(ATA,MATPRODUCT_AB));
304 7 : PetscCall(MatProductSetFromOptions(ATA));
305 7 : PetscCall(MatProductSymbolic(ATA));
306 7 : PetscCall(MatProductNumeric(ATA));
307 7 : PetscCall(MatCreateTile(1.0,Id,1.0,Zm,1.0,Zn,1.0,ATA,C));
308 7 : PetscCall(MatDestroy(&Id));
309 7 : PetscCall(MatDestroy(&Zm));
310 7 : PetscCall(MatDestroy(&Zn));
311 7 : PetscCall(MatDestroy(&ATA));
312 : } else {
313 6 : PetscCall(PetscNew(&ctx));
314 6 : ctx->A = A;
315 6 : ctx->AT = AT;
316 6 : ctx->swapped = svd->swapped;
317 6 : PetscCall(VecDuplicateEmpty(t,&ctx->x1));
318 6 : PetscCall(VecDuplicateEmpty(t,&ctx->y1));
319 6 : PetscCall(MatCreateVecsEmpty(A,&ctx->x2,NULL));
320 6 : PetscCall(MatCreateVecsEmpty(A,&ctx->y2,NULL));
321 6 : PetscCall(MatCreateVecs(A,NULL,&ctx->w));
322 6 : PetscCall(MatCreateShell(PetscObjectComm((PetscObject)svd),m+n,m+n,M+N,M+N,ctx,C));
323 6 : PetscCall(MatShellSetOperation(*C,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_ECross));
324 6 : PetscCall(MatShellSetOperation(*C,MATOP_DESTROY,(void(*)(void))MatDestroy_ECross));
325 : #if defined(PETSC_HAVE_CUDA)
326 : PetscCall(PetscObjectTypeCompareAny((PetscObject)(svd->swapped?AT:A),&gpu,MATSEQAIJCUSPARSE,MATMPIAIJCUSPARSE,""));
327 : if (gpu) PetscCall(MatShellSetOperation(*C,MATOP_MULT,(void(*)(void))MatMult_ECross_CUDA));
328 : else
329 : #elif defined(PETSC_HAVE_HIP)
330 : PetscCall(PetscObjectTypeCompareAny((PetscObject)(svd->swapped?AT:A),&gpu,MATSEQAIJHIPSPARSE,MATMPIAIJHIPSPARSE,""));
331 : if (gpu) PetscCall(MatShellSetOperation(*C,MATOP_MULT,(void(*)(void))MatMult_ECross_HIP));
332 : else
333 : #endif
334 6 : PetscCall(MatShellSetOperation(*C,MATOP_MULT,(void(*)(void))MatMult_ECross));
335 6 : PetscCall(MatGetVecType(A,&vtype));
336 6 : PetscCall(MatSetVecType(*C,vtype));
337 : #if defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP)
338 : if (gpu) {
339 : /* check alignment of bottom block */
340 : PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)ctx->x1),&size));
341 : PetscCall(VecGetOwnershipRanges(ctx->x1,&ranges));
342 : for (i=0;i<size;i++) {
343 : ctx->misaligned = (((ranges[i+1]-ranges[i])*sizeof(PetscScalar))%16)? PETSC_TRUE: PETSC_FALSE;
344 : if (ctx->misaligned) break;
345 : }
346 : if (ctx->misaligned) { /* create work vectors for MatMult */
347 : PetscCall(VecDuplicate(ctx->x2,&ctx->wx2));
348 : PetscCall(VecDuplicate(ctx->y2,&ctx->wy2));
349 : }
350 : }
351 : #endif
352 : }
353 13 : PetscFunctionReturn(PETSC_SUCCESS);
354 : }
355 :
356 : /* Convergence test relative to the norm of R (used in GSVD only) */
357 175 : static PetscErrorCode EPSConv_Cyclic(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
358 : {
359 175 : SVD svd = (SVD)ctx;
360 :
361 175 : PetscFunctionBegin;
362 175 : *errest = res/PetscMax(svd->nrma,svd->nrmb);
363 175 : PetscFunctionReturn(PETSC_SUCCESS);
364 : }
365 :
366 46 : static PetscErrorCode SVDSetUp_Cyclic(SVD svd)
367 : {
368 46 : SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
369 46 : PetscInt M,N,m,n,p,k,i,isl,offset,nev,ncv,mpd,maxit;
370 46 : PetscReal tol;
371 46 : const PetscScalar *isa,*oa;
372 46 : PetscScalar *va;
373 46 : EPSProblemType ptype;
374 46 : PetscBool trackall,issinv;
375 46 : Vec v,t;
376 46 : ST st;
377 46 : Mat Omega;
378 46 : MatType Atype;
379 :
380 46 : PetscFunctionBegin;
381 46 : PetscCall(MatGetSize(svd->A,&M,&N));
382 46 : PetscCall(MatGetLocalSize(svd->A,&m,&n));
383 46 : if (!cyclic->eps) PetscCall(SVDCyclicGetEPS(svd,&cyclic->eps));
384 46 : PetscCall(MatDestroy(&cyclic->C));
385 46 : PetscCall(MatDestroy(&cyclic->D));
386 46 : if (svd->isgeneralized) {
387 13 : if (svd->which==SVD_SMALLEST) { /* alternative pencil */
388 3 : PetscCall(MatCreateVecs(svd->B,NULL,&t));
389 3 : PetscCall(SVDCyclicGetCyclicMat(svd,svd->B,svd->BT,&cyclic->C));
390 3 : PetscCall(SVDCyclicGetECrossMat(svd,svd->A,svd->AT,&cyclic->D,t));
391 : } else {
392 10 : PetscCall(MatCreateVecs(svd->A,NULL,&t));
393 10 : PetscCall(SVDCyclicGetCyclicMat(svd,svd->A,svd->AT,&cyclic->C));
394 10 : PetscCall(SVDCyclicGetECrossMat(svd,svd->B,svd->BT,&cyclic->D,t));
395 : }
396 13 : PetscCall(VecDestroy(&t));
397 13 : PetscCall(EPSSetOperators(cyclic->eps,cyclic->C,cyclic->D));
398 13 : PetscCall(EPSGetProblemType(cyclic->eps,&ptype));
399 13 : if (!ptype) PetscCall(EPSSetProblemType(cyclic->eps,EPS_GHEP));
400 33 : } else if (svd->ishyperbolic) {
401 8 : PetscCall(SVDCyclicGetCyclicMat(svd,svd->A,svd->AT,&cyclic->C));
402 8 : PetscCall(MatCreateVecs(cyclic->C,&v,NULL));
403 8 : PetscCall(VecSet(v,1.0));
404 8 : PetscCall(VecGetArrayRead(svd->omega,&oa));
405 8 : PetscCall(VecGetArray(v,&va));
406 8 : if (svd->swapped) PetscCall(PetscArraycpy(va+m,oa,n));
407 6 : else PetscCall(PetscArraycpy(va,oa,m));
408 8 : PetscCall(VecRestoreArrayRead(svd->omega,&oa));
409 8 : PetscCall(VecRestoreArray(v,&va));
410 8 : PetscCall(MatGetType(svd->OP,&Atype));
411 8 : PetscCall(MatCreate(PetscObjectComm((PetscObject)svd),&Omega));
412 8 : PetscCall(MatSetSizes(Omega,m+n,m+n,M+N,M+N));
413 8 : PetscCall(MatSetType(Omega,Atype));
414 8 : PetscCall(MatDiagonalSet(Omega,v,INSERT_VALUES));
415 8 : PetscCall(EPSSetOperators(cyclic->eps,cyclic->C,Omega));
416 8 : PetscCall(EPSSetProblemType(cyclic->eps,EPS_GHIEP));
417 8 : PetscCall(MatDestroy(&Omega));
418 8 : PetscCall(VecDestroy(&v));
419 : } else {
420 25 : PetscCall(SVDCyclicGetCyclicMat(svd,svd->A,svd->AT,&cyclic->C));
421 25 : PetscCall(EPSSetOperators(cyclic->eps,cyclic->C,NULL));
422 25 : PetscCall(EPSSetProblemType(cyclic->eps,EPS_HEP));
423 : }
424 46 : if (!cyclic->usereps) {
425 45 : if (svd->which == SVD_LARGEST) {
426 39 : PetscCall(EPSGetST(cyclic->eps,&st));
427 39 : PetscCall(PetscObjectTypeCompare((PetscObject)st,STSINVERT,&issinv));
428 39 : if (issinv) PetscCall(EPSSetWhichEigenpairs(cyclic->eps,EPS_TARGET_MAGNITUDE));
429 38 : else if (svd->ishyperbolic) PetscCall(EPSSetWhichEigenpairs(cyclic->eps,EPS_LARGEST_MAGNITUDE));
430 31 : else PetscCall(EPSSetWhichEigenpairs(cyclic->eps,EPS_LARGEST_REAL));
431 : } else {
432 6 : if (svd->isgeneralized) { /* computes sigma^{-1} via alternative pencil */
433 3 : PetscCall(EPSSetWhichEigenpairs(cyclic->eps,EPS_LARGEST_REAL));
434 : } else {
435 3 : if (svd->ishyperbolic) PetscCall(EPSSetWhichEigenpairs(cyclic->eps,EPS_TARGET_MAGNITUDE));
436 2 : else PetscCall(EPSSetEigenvalueComparison(cyclic->eps,SlepcCompareSmallestPosReal,NULL));
437 3 : PetscCall(EPSSetTarget(cyclic->eps,0.0));
438 : }
439 : }
440 45 : PetscCall(EPSGetDimensions(cyclic->eps,&nev,&ncv,&mpd));
441 45 : PetscCheck(nev==1 || nev>=2*svd->nsv,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_WRONG,"The number of requested eigenvalues %" PetscInt_FMT " must be at least 2*%" PetscInt_FMT,nev,svd->nsv);
442 45 : nev = PetscMax(nev,2*svd->nsv);
443 45 : if (ncv==PETSC_DETERMINE && svd->ncv!=PETSC_DETERMINE) ncv = PetscMax(3*svd->nsv,svd->ncv);
444 45 : if (mpd==PETSC_DETERMINE && svd->mpd!=PETSC_DETERMINE) mpd = svd->mpd;
445 45 : PetscCall(EPSSetDimensions(cyclic->eps,nev,ncv,mpd));
446 45 : PetscCall(EPSGetTolerances(cyclic->eps,&tol,&maxit));
447 51 : if (tol==(PetscReal)PETSC_DETERMINE) tol = svd->tol==(PetscReal)PETSC_DETERMINE? SLEPC_DEFAULT_TOL/10.0: svd->tol;
448 45 : if (maxit==PETSC_DETERMINE) maxit = svd->max_it;
449 45 : PetscCall(EPSSetTolerances(cyclic->eps,tol,maxit));
450 45 : switch (svd->conv) {
451 3 : case SVD_CONV_ABS:
452 3 : PetscCall(EPSSetConvergenceTest(cyclic->eps,EPS_CONV_ABS));break;
453 29 : case SVD_CONV_REL:
454 29 : PetscCall(EPSSetConvergenceTest(cyclic->eps,EPS_CONV_REL));break;
455 13 : case SVD_CONV_NORM:
456 13 : if (svd->isgeneralized) {
457 13 : if (!svd->nrma) PetscCall(MatNorm(svd->OP,NORM_INFINITY,&svd->nrma));
458 13 : if (!svd->nrmb) PetscCall(MatNorm(svd->OPb,NORM_INFINITY,&svd->nrmb));
459 13 : PetscCall(EPSSetConvergenceTestFunction(cyclic->eps,EPSConv_Cyclic,svd,NULL));
460 : } else {
461 0 : PetscCall(EPSSetConvergenceTest(cyclic->eps,EPS_CONV_NORM));break;
462 : }
463 : break;
464 0 : case SVD_CONV_MAXIT:
465 0 : SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Maxit convergence test not supported in this solver");
466 0 : case SVD_CONV_USER:
467 0 : SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"User-defined convergence test not supported in this solver");
468 : }
469 1 : }
470 46 : SVDCheckUnsupported(svd,SVD_FEATURE_STOPPING);
471 : /* Transfer the trackall option from svd to eps */
472 46 : PetscCall(SVDGetTrackAll(svd,&trackall));
473 46 : PetscCall(EPSSetTrackAll(cyclic->eps,trackall));
474 : /* Transfer the initial subspace from svd to eps */
475 46 : if (svd->nini<0 || svd->ninil<0) {
476 8 : for (i=0;i<-PetscMin(svd->nini,svd->ninil);i++) {
477 4 : PetscCall(MatCreateVecs(cyclic->C,&v,NULL));
478 4 : PetscCall(VecGetArrayWrite(v,&va));
479 4 : if (svd->isgeneralized) PetscCall(MatGetLocalSize(svd->B,&p,NULL));
480 4 : k = (svd->isgeneralized && svd->which==SVD_SMALLEST)? p: m; /* size of upper block row */
481 4 : if (i<-svd->ninil) {
482 4 : PetscCall(VecGetArrayRead(svd->ISL[i],&isa));
483 4 : if (svd->isgeneralized) {
484 1 : PetscCall(VecGetLocalSize(svd->ISL[i],&isl));
485 1 : PetscCheck(isl==m+p,PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Size mismatch for left initial vector");
486 1 : offset = (svd->which==SVD_SMALLEST)? m: 0;
487 1 : PetscCall(PetscArraycpy(va,isa+offset,k));
488 : } else {
489 3 : PetscCall(VecGetLocalSize(svd->ISL[i],&isl));
490 3 : PetscCheck(isl==k,PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Size mismatch for left initial vector");
491 3 : PetscCall(PetscArraycpy(va,isa,k));
492 : }
493 4 : PetscCall(VecRestoreArrayRead(svd->IS[i],&isa));
494 0 : } else PetscCall(PetscArrayzero(&va,k));
495 4 : if (i<-svd->nini) {
496 4 : PetscCall(VecGetLocalSize(svd->IS[i],&isl));
497 4 : PetscCheck(isl==n,PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Size mismatch for right initial vector");
498 4 : PetscCall(VecGetArrayRead(svd->IS[i],&isa));
499 4 : PetscCall(PetscArraycpy(va+k,isa,n));
500 4 : PetscCall(VecRestoreArrayRead(svd->IS[i],&isa));
501 0 : } else PetscCall(PetscArrayzero(va+k,n));
502 4 : PetscCall(VecRestoreArrayWrite(v,&va));
503 4 : PetscCall(VecDestroy(&svd->IS[i]));
504 4 : svd->IS[i] = v;
505 : }
506 4 : svd->nini = PetscMin(svd->nini,svd->ninil);
507 4 : PetscCall(EPSSetInitialSpace(cyclic->eps,-svd->nini,svd->IS));
508 4 : PetscCall(SlepcBasisDestroy_Private(&svd->nini,&svd->IS));
509 4 : PetscCall(SlepcBasisDestroy_Private(&svd->ninil,&svd->ISL));
510 : }
511 46 : PetscCall(EPSSetUp(cyclic->eps));
512 46 : PetscCall(EPSGetDimensions(cyclic->eps,NULL,&svd->ncv,&svd->mpd));
513 46 : svd->ncv = PetscMin(svd->ncv,PetscMin(M,N));
514 46 : PetscCall(EPSGetTolerances(cyclic->eps,NULL,&svd->max_it));
515 46 : if (svd->tol==(PetscReal)PETSC_DETERMINE) svd->tol = SLEPC_DEFAULT_TOL;
516 :
517 46 : svd->leftbasis = PETSC_TRUE;
518 46 : PetscCall(SVDAllocateSolution(svd,0));
519 46 : PetscFunctionReturn(PETSC_SUCCESS);
520 : }
521 :
522 84714 : static PetscErrorCode SVDCyclicCheckEigenvalue(SVD svd,PetscScalar er,PetscScalar ei,PetscReal *sigma,PetscBool *isreal)
523 : {
524 84714 : PetscFunctionBegin;
525 84714 : if (svd->ishyperbolic && PetscDefined(USE_COMPLEX) && PetscAbsReal(PetscImaginaryPart(er))>10*PetscAbsReal(PetscRealPart(er))) {
526 : *sigma = PetscImaginaryPart(er);
527 : if (isreal) *isreal = PETSC_FALSE;
528 84714 : } else if (svd->ishyperbolic && !PetscDefined(USE_COMPLEX) && PetscAbsScalar(ei)>10*PetscAbsScalar(er)) {
529 4956 : *sigma = PetscRealPart(ei);
530 4956 : if (isreal) *isreal = PETSC_FALSE;
531 : } else {
532 79758 : *sigma = PetscRealPart(er);
533 79758 : if (isreal) *isreal = PETSC_TRUE;
534 : }
535 84714 : PetscFunctionReturn(PETSC_SUCCESS);
536 : }
537 :
538 46 : static PetscErrorCode SVDSolve_Cyclic(SVD svd)
539 : {
540 46 : SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
541 46 : PetscInt i,j,nconv;
542 46 : PetscScalar er,ei;
543 46 : PetscReal sigma;
544 :
545 46 : PetscFunctionBegin;
546 46 : PetscCall(EPSSolve(cyclic->eps));
547 46 : PetscCall(EPSGetConverged(cyclic->eps,&nconv));
548 46 : PetscCall(EPSGetIterationNumber(cyclic->eps,&svd->its));
549 46 : PetscCall(EPSGetConvergedReason(cyclic->eps,(EPSConvergedReason*)&svd->reason));
550 743 : for (i=0,j=0;i<nconv;i++) {
551 697 : PetscCall(EPSGetEigenvalue(cyclic->eps,i,&er,&ei));
552 697 : PetscCall(SVDCyclicCheckEigenvalue(svd,er,ei,&sigma,NULL));
553 697 : if (sigma>0.0) {
554 469 : if (svd->isgeneralized && svd->which==SVD_SMALLEST) svd->sigma[j] = 1.0/sigma;
555 453 : else svd->sigma[j] = sigma;
556 469 : j++;
557 : }
558 : }
559 46 : svd->nconv = j;
560 46 : PetscFunctionReturn(PETSC_SUCCESS);
561 : }
562 :
563 21 : static PetscErrorCode SVDComputeVectors_Cyclic_Standard(SVD svd)
564 : {
565 21 : SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
566 21 : PetscInt i,j,m,nconv;
567 21 : PetscScalar er,ei;
568 21 : PetscReal sigma;
569 21 : const PetscScalar *px;
570 21 : Vec x,x1,x2;
571 :
572 21 : PetscFunctionBegin;
573 21 : PetscCall(MatCreateVecs(cyclic->C,&x,NULL));
574 21 : PetscCall(MatGetLocalSize(svd->A,&m,NULL));
575 21 : PetscCall(MatCreateVecsEmpty(svd->A,&x2,&x1));
576 21 : PetscCall(EPSGetConverged(cyclic->eps,&nconv));
577 186 : for (i=0,j=0;i<nconv;i++) {
578 165 : PetscCall(EPSGetEigenpair(cyclic->eps,i,&er,&ei,x,NULL));
579 165 : PetscCall(SVDCyclicCheckEigenvalue(svd,er,ei,&sigma,NULL));
580 165 : if (sigma<0.0) continue;
581 153 : PetscCall(VecGetArrayRead(x,&px));
582 153 : PetscCall(VecPlaceArray(x1,px));
583 153 : PetscCall(VecPlaceArray(x2,px+m));
584 153 : PetscCall(BVInsertVec(svd->U,j,x1));
585 153 : PetscCall(BVScaleColumn(svd->U,j,PETSC_SQRT2));
586 153 : PetscCall(BVInsertVec(svd->V,j,x2));
587 153 : PetscCall(BVScaleColumn(svd->V,j,PETSC_SQRT2));
588 153 : PetscCall(VecResetArray(x1));
589 153 : PetscCall(VecResetArray(x2));
590 153 : PetscCall(VecRestoreArrayRead(x,&px));
591 153 : j++;
592 : }
593 21 : PetscCall(VecDestroy(&x));
594 21 : PetscCall(VecDestroy(&x1));
595 21 : PetscCall(VecDestroy(&x2));
596 21 : PetscFunctionReturn(PETSC_SUCCESS);
597 : }
598 :
599 13 : static PetscErrorCode SVDComputeVectors_Cyclic_Generalized(SVD svd)
600 : {
601 13 : SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
602 13 : PetscInt i,j,m,p,nconv;
603 13 : PetscScalar *dst,er,ei;
604 13 : PetscReal sigma;
605 13 : const PetscScalar *src,*px;
606 13 : Vec u,v,x,x1,x2,uv;
607 :
608 13 : PetscFunctionBegin;
609 13 : PetscCall(MatGetLocalSize(svd->A,&m,NULL));
610 13 : PetscCall(MatGetLocalSize(svd->B,&p,NULL));
611 13 : PetscCall(MatCreateVecs(cyclic->C,&x,NULL));
612 13 : if (svd->which==SVD_SMALLEST) PetscCall(MatCreateVecsEmpty(svd->B,&x1,&x2));
613 10 : else PetscCall(MatCreateVecsEmpty(svd->A,&x2,&x1));
614 13 : PetscCall(MatCreateVecs(svd->A,NULL,&u));
615 13 : PetscCall(MatCreateVecs(svd->B,NULL,&v));
616 13 : PetscCall(EPSGetConverged(cyclic->eps,&nconv));
617 100 : for (i=0,j=0;i<nconv;i++) {
618 87 : PetscCall(EPSGetEigenpair(cyclic->eps,i,&er,&ei,x,NULL));
619 87 : PetscCall(SVDCyclicCheckEigenvalue(svd,er,ei,&sigma,NULL));
620 87 : if (sigma<0.0) continue;
621 87 : if (svd->which==SVD_SMALLEST) {
622 : /* evec_i = 1/sqrt(2)*[ v_i; w_i ], w_i = x_i/c_i */
623 16 : PetscCall(VecGetArrayRead(x,&px));
624 16 : PetscCall(VecPlaceArray(x2,px));
625 16 : PetscCall(VecPlaceArray(x1,px+p));
626 16 : PetscCall(VecCopy(x2,v));
627 16 : PetscCall(VecScale(v,PETSC_SQRT2)); /* v_i = sqrt(2)*evec_i_1 */
628 16 : PetscCall(VecScale(x1,PETSC_SQRT2)); /* w_i = sqrt(2)*evec_i_2 */
629 16 : PetscCall(MatMult(svd->A,x1,u)); /* A*w_i = u_i */
630 16 : PetscCall(VecScale(x1,1.0/PetscSqrtScalar(1.0+sigma*sigma))); /* x_i = w_i*c_i */
631 16 : PetscCall(BVInsertVec(svd->V,j,x1));
632 16 : PetscCall(VecResetArray(x2));
633 16 : PetscCall(VecResetArray(x1));
634 16 : PetscCall(VecRestoreArrayRead(x,&px));
635 : } else {
636 : /* evec_i = 1/sqrt(2)*[ u_i; w_i ], w_i = x_i/s_i */
637 71 : PetscCall(VecGetArrayRead(x,&px));
638 71 : PetscCall(VecPlaceArray(x1,px));
639 71 : PetscCall(VecPlaceArray(x2,px+m));
640 71 : PetscCall(VecCopy(x1,u));
641 71 : PetscCall(VecScale(u,PETSC_SQRT2)); /* u_i = sqrt(2)*evec_i_1 */
642 71 : PetscCall(VecScale(x2,PETSC_SQRT2)); /* w_i = sqrt(2)*evec_i_2 */
643 71 : PetscCall(MatMult(svd->B,x2,v)); /* B*w_i = v_i */
644 71 : PetscCall(VecScale(x2,1.0/PetscSqrtScalar(1.0+sigma*sigma))); /* x_i = w_i*s_i */
645 71 : PetscCall(BVInsertVec(svd->V,j,x2));
646 71 : PetscCall(VecResetArray(x1));
647 71 : PetscCall(VecResetArray(x2));
648 71 : PetscCall(VecRestoreArrayRead(x,&px));
649 : }
650 : /* copy [u;v] to U[j] */
651 87 : PetscCall(BVGetColumn(svd->U,j,&uv));
652 87 : PetscCall(VecGetArrayWrite(uv,&dst));
653 87 : PetscCall(VecGetArrayRead(u,&src));
654 87 : PetscCall(PetscArraycpy(dst,src,m));
655 87 : PetscCall(VecRestoreArrayRead(u,&src));
656 87 : PetscCall(VecGetArrayRead(v,&src));
657 87 : PetscCall(PetscArraycpy(dst+m,src,p));
658 87 : PetscCall(VecRestoreArrayRead(v,&src));
659 87 : PetscCall(VecRestoreArrayWrite(uv,&dst));
660 87 : PetscCall(BVRestoreColumn(svd->U,j,&uv));
661 87 : j++;
662 : }
663 13 : PetscCall(VecDestroy(&x));
664 13 : PetscCall(VecDestroy(&x1));
665 13 : PetscCall(VecDestroy(&x2));
666 13 : PetscCall(VecDestroy(&u));
667 13 : PetscCall(VecDestroy(&v));
668 13 : PetscFunctionReturn(PETSC_SUCCESS);
669 : }
670 :
671 : #if defined(PETSC_USE_COMPLEX)
672 : /* VecMaxAbs: returns the entry of x that has max(abs(x(i))), using w as a workspace vector */
673 : static PetscErrorCode VecMaxAbs(Vec x,Vec w,PetscScalar *v)
674 : {
675 : PetscMPIInt size,rank,root;
676 : const PetscScalar *xx;
677 : const PetscInt *ranges;
678 : PetscReal val;
679 : PetscInt p;
680 :
681 : PetscFunctionBegin;
682 : PetscCall(VecCopy(x,w));
683 : PetscCall(VecAbs(w));
684 : PetscCall(VecMax(w,&p,&val));
685 : PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)x),&size));
686 : PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)x),&rank));
687 : PetscCall(VecGetOwnershipRanges(x,&ranges));
688 : for (root=0;root<size;root++) if (p>=ranges[root] && p<ranges[root+1]) break;
689 : if (rank==root) {
690 : PetscCall(VecGetArrayRead(x,&xx));
691 : *v = xx[p-ranges[root]];
692 : PetscCall(VecRestoreArrayRead(x,&xx));
693 : }
694 : PetscCallMPI(MPI_Bcast(v,1,MPIU_SCALAR,root,PetscObjectComm((PetscObject)x)));
695 : PetscFunctionReturn(PETSC_SUCCESS);
696 : }
697 : #endif
698 :
699 8 : static PetscErrorCode SVDComputeVectors_Cyclic_Hyperbolic(SVD svd)
700 : {
701 8 : SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
702 8 : PetscInt i,j,m,n,nconv;
703 8 : PetscScalar er,ei;
704 8 : PetscReal sigma,nrm;
705 8 : PetscBool isreal;
706 8 : const PetscScalar *px;
707 8 : Vec u,x,xi=NULL,x1,x2,x1i=NULL,x2i;
708 8 : BV U=NULL,V=NULL;
709 : #if !defined(PETSC_USE_COMPLEX)
710 8 : const PetscScalar *pxi;
711 8 : PetscReal nrmr,nrmi;
712 : #else
713 : PetscScalar alpha;
714 : #endif
715 :
716 8 : PetscFunctionBegin;
717 8 : PetscCall(MatCreateVecs(cyclic->C,&x,svd->ishyperbolic?&xi:NULL));
718 8 : PetscCall(MatGetLocalSize(svd->A,&m,NULL));
719 8 : PetscCall(MatCreateVecsEmpty(svd->OP,&x2,&x1));
720 : #if defined(PETSC_USE_COMPLEX)
721 : PetscCall(MatCreateVecs(svd->OP,&x2i,&x1i));
722 : #else
723 8 : PetscCall(MatCreateVecsEmpty(svd->OP,&x2i,&x1i));
724 : #endif
725 : /* set-up Omega-normalization of U */
726 8 : U = svd->swapped? svd->V: svd->U;
727 8 : V = svd->swapped? svd->U: svd->V;
728 8 : PetscCall(BVGetSizes(U,&n,NULL,NULL));
729 8 : PetscCall(BV_SetMatrixDiagonal(U,svd->omega,svd->A));
730 8 : PetscCall(EPSGetConverged(cyclic->eps,&nconv));
731 441 : for (i=0,j=0;i<nconv;i++) {
732 433 : PetscCall(EPSGetEigenpair(cyclic->eps,i,&er,&ei,x,xi));
733 433 : PetscCall(SVDCyclicCheckEigenvalue(svd,er,ei,&sigma,&isreal));
734 433 : if (sigma<0.0) continue;
735 217 : PetscCall(VecGetArrayRead(x,&px));
736 217 : if (svd->swapped) {
737 12 : PetscCall(VecPlaceArray(x2,px));
738 12 : PetscCall(VecPlaceArray(x1,px+m));
739 : } else {
740 205 : PetscCall(VecPlaceArray(x1,px));
741 205 : PetscCall(VecPlaceArray(x2,px+n));
742 : }
743 : #if defined(PETSC_USE_COMPLEX)
744 : PetscCall(BVInsertVec(U,j,x1));
745 : PetscCall(BVInsertVec(V,j,x2));
746 : if (!isreal) {
747 : PetscCall(VecMaxAbs(x1,x1i,&alpha));
748 : PetscCall(BVScaleColumn(U,j,PetscAbsScalar(alpha)/alpha));
749 : PetscCall(BVScaleColumn(V,j,PetscAbsScalar(alpha)/(alpha*PETSC_i)));
750 : }
751 : #else
752 217 : PetscCall(VecGetArrayRead(xi,&pxi));
753 217 : if (svd->swapped) {
754 12 : PetscCall(VecPlaceArray(x2i,pxi));
755 12 : PetscCall(VecPlaceArray(x1i,pxi+m));
756 : } else {
757 205 : PetscCall(VecPlaceArray(x1i,pxi));
758 205 : PetscCall(VecPlaceArray(x2i,pxi+n));
759 : }
760 217 : PetscCall(VecNorm(x2,NORM_2,&nrmr));
761 217 : PetscCall(VecNorm(x2i,NORM_2,&nrmi));
762 217 : if (nrmi>nrmr) {
763 3 : if (isreal) {
764 0 : PetscCall(BVInsertVec(U,j,x1i));
765 0 : PetscCall(BVInsertVec(V,j,x2i));
766 : } else {
767 3 : PetscCall(BVInsertVec(U,j,x1));
768 3 : PetscCall(BVInsertVec(V,j,x2i));
769 : }
770 : } else {
771 214 : if (isreal) {
772 213 : PetscCall(BVInsertVec(U,j,x1));
773 213 : PetscCall(BVInsertVec(V,j,x2));
774 : } else {
775 1 : PetscCall(BVInsertVec(U,j,x1i));
776 1 : PetscCall(BVScaleColumn(U,j,-1.0));
777 1 : PetscCall(BVInsertVec(V,j,x2));
778 : }
779 : }
780 217 : PetscCall(VecResetArray(x1i));
781 217 : PetscCall(VecResetArray(x2i));
782 217 : PetscCall(VecRestoreArrayRead(xi,&pxi));
783 : #endif
784 217 : PetscCall(VecResetArray(x1));
785 217 : PetscCall(VecResetArray(x2));
786 217 : PetscCall(VecRestoreArrayRead(x,&px));
787 217 : PetscCall(BVGetColumn(U,j,&u));
788 217 : PetscCall(VecPointwiseMult(u,u,svd->omega));
789 217 : PetscCall(BVRestoreColumn(U,j,&u));
790 217 : PetscCall(BVNormColumn(U,j,NORM_2,&nrm));
791 217 : PetscCall(BVScaleColumn(U,j,1.0/PetscAbs(nrm)));
792 217 : svd->sign[j] = PetscSign(nrm);
793 217 : PetscCall(BVNormColumn(V,j,NORM_2,&nrm));
794 217 : PetscCall(BVScaleColumn(V,j,1.0/nrm));
795 217 : j++;
796 : }
797 8 : PetscCall(VecDestroy(&x));
798 8 : PetscCall(VecDestroy(&x1));
799 8 : PetscCall(VecDestroy(&x2));
800 8 : PetscCall(VecDestroy(&xi));
801 8 : PetscCall(VecDestroy(&x1i));
802 8 : PetscCall(VecDestroy(&x2i));
803 8 : PetscFunctionReturn(PETSC_SUCCESS);
804 : }
805 :
806 42 : static PetscErrorCode SVDComputeVectors_Cyclic(SVD svd)
807 : {
808 42 : PetscFunctionBegin;
809 42 : switch (svd->problem_type) {
810 21 : case SVD_STANDARD:
811 21 : PetscCall(SVDComputeVectors_Cyclic_Standard(svd));
812 : break;
813 13 : case SVD_GENERALIZED:
814 13 : PetscCall(SVDComputeVectors_Cyclic_Generalized(svd));
815 : break;
816 8 : case SVD_HYPERBOLIC:
817 8 : PetscCall(SVDComputeVectors_Cyclic_Hyperbolic(svd));
818 : break;
819 0 : default:
820 0 : SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_WRONG,"Unknown singular value problem type");
821 : }
822 42 : PetscFunctionReturn(PETSC_SUCCESS);
823 : }
824 :
825 2883 : static PetscErrorCode EPSMonitor_Cyclic(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
826 : {
827 2883 : PetscInt i,j;
828 2883 : SVD svd = (SVD)ctx;
829 2883 : PetscScalar er,ei;
830 2883 : PetscReal sigma;
831 2883 : ST st;
832 :
833 2883 : PetscFunctionBegin;
834 2883 : nconv = 0;
835 2883 : PetscCall(EPSGetST(eps,&st));
836 86215 : for (i=0,j=0;i<PetscMin(nest,svd->ncv);i++) {
837 83332 : er = eigr[i]; ei = eigi[i];
838 83332 : PetscCall(STBackTransform(st,1,&er,&ei));
839 83332 : PetscCall(SVDCyclicCheckEigenvalue(svd,er,ei,&sigma,NULL));
840 83332 : if (sigma>0.0) {
841 36176 : svd->sigma[j] = sigma;
842 36176 : svd->errest[j] = errest[i];
843 36176 : if (errest[i] && errest[i] < svd->tol) nconv++;
844 36176 : j++;
845 : }
846 : }
847 2883 : nest = j;
848 2883 : PetscCall(SVDMonitor(svd,its,nconv,svd->sigma,svd->errest,nest));
849 2883 : PetscFunctionReturn(PETSC_SUCCESS);
850 : }
851 :
852 34 : static PetscErrorCode SVDSetFromOptions_Cyclic(SVD svd,PetscOptionItems *PetscOptionsObject)
853 : {
854 34 : PetscBool set,val;
855 34 : SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
856 34 : ST st;
857 :
858 34 : PetscFunctionBegin;
859 34 : PetscOptionsHeadBegin(PetscOptionsObject,"SVD Cyclic Options");
860 :
861 34 : PetscCall(PetscOptionsBool("-svd_cyclic_explicitmatrix","Use cyclic explicit matrix","SVDCyclicSetExplicitMatrix",cyclic->explicitmatrix,&val,&set));
862 34 : if (set) PetscCall(SVDCyclicSetExplicitMatrix(svd,val));
863 :
864 34 : PetscOptionsHeadEnd();
865 :
866 34 : if (!cyclic->eps) PetscCall(SVDCyclicGetEPS(svd,&cyclic->eps));
867 34 : if (!cyclic->explicitmatrix && !cyclic->usereps) {
868 : /* use as default an ST with shell matrix and Jacobi */
869 20 : PetscCall(EPSGetST(cyclic->eps,&st));
870 20 : PetscCall(STSetMatMode(st,ST_MATMODE_SHELL));
871 : }
872 34 : PetscCall(EPSSetFromOptions(cyclic->eps));
873 34 : PetscFunctionReturn(PETSC_SUCCESS);
874 : }
875 :
876 22 : static PetscErrorCode SVDCyclicSetExplicitMatrix_Cyclic(SVD svd,PetscBool explicitmat)
877 : {
878 22 : SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
879 :
880 22 : PetscFunctionBegin;
881 22 : if (cyclic->explicitmatrix != explicitmat) {
882 14 : cyclic->explicitmatrix = explicitmat;
883 14 : svd->state = SVD_STATE_INITIAL;
884 : }
885 22 : PetscFunctionReturn(PETSC_SUCCESS);
886 : }
887 :
888 : /*@
889 : SVDCyclicSetExplicitMatrix - Indicate if the eigensolver operator
890 : H(A) = [ 0 A ; A^T 0 ] must be computed explicitly.
891 :
892 : Logically Collective
893 :
894 : Input Parameters:
895 : + svd - singular value solver
896 : - explicitmat - boolean flag indicating if H(A) is built explicitly
897 :
898 : Options Database Key:
899 : . -svd_cyclic_explicitmatrix <boolean> - Indicates the boolean flag
900 :
901 : Level: advanced
902 :
903 : .seealso: SVDCyclicGetExplicitMatrix()
904 : @*/
905 22 : PetscErrorCode SVDCyclicSetExplicitMatrix(SVD svd,PetscBool explicitmat)
906 : {
907 22 : PetscFunctionBegin;
908 22 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
909 66 : PetscValidLogicalCollectiveBool(svd,explicitmat,2);
910 22 : PetscTryMethod(svd,"SVDCyclicSetExplicitMatrix_C",(SVD,PetscBool),(svd,explicitmat));
911 22 : PetscFunctionReturn(PETSC_SUCCESS);
912 : }
913 :
914 1 : static PetscErrorCode SVDCyclicGetExplicitMatrix_Cyclic(SVD svd,PetscBool *explicitmat)
915 : {
916 1 : SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
917 :
918 1 : PetscFunctionBegin;
919 1 : *explicitmat = cyclic->explicitmatrix;
920 1 : PetscFunctionReturn(PETSC_SUCCESS);
921 : }
922 :
923 : /*@
924 : SVDCyclicGetExplicitMatrix - Returns the flag indicating if H(A) is built explicitly.
925 :
926 : Not Collective
927 :
928 : Input Parameter:
929 : . svd - singular value solver
930 :
931 : Output Parameter:
932 : . explicitmat - the mode flag
933 :
934 : Level: advanced
935 :
936 : .seealso: SVDCyclicSetExplicitMatrix()
937 : @*/
938 1 : PetscErrorCode SVDCyclicGetExplicitMatrix(SVD svd,PetscBool *explicitmat)
939 : {
940 1 : PetscFunctionBegin;
941 1 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
942 1 : PetscAssertPointer(explicitmat,2);
943 1 : PetscUseMethod(svd,"SVDCyclicGetExplicitMatrix_C",(SVD,PetscBool*),(svd,explicitmat));
944 1 : PetscFunctionReturn(PETSC_SUCCESS);
945 : }
946 :
947 1 : static PetscErrorCode SVDCyclicSetEPS_Cyclic(SVD svd,EPS eps)
948 : {
949 1 : SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
950 :
951 1 : PetscFunctionBegin;
952 1 : PetscCall(PetscObjectReference((PetscObject)eps));
953 1 : PetscCall(EPSDestroy(&cyclic->eps));
954 1 : cyclic->eps = eps;
955 1 : cyclic->usereps = PETSC_TRUE;
956 1 : svd->state = SVD_STATE_INITIAL;
957 1 : PetscFunctionReturn(PETSC_SUCCESS);
958 : }
959 :
960 : /*@
961 : SVDCyclicSetEPS - Associate an eigensolver object (EPS) to the
962 : singular value solver.
963 :
964 : Collective
965 :
966 : Input Parameters:
967 : + svd - singular value solver
968 : - eps - the eigensolver object
969 :
970 : Level: advanced
971 :
972 : .seealso: SVDCyclicGetEPS()
973 : @*/
974 1 : PetscErrorCode SVDCyclicSetEPS(SVD svd,EPS eps)
975 : {
976 1 : PetscFunctionBegin;
977 1 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
978 1 : PetscValidHeaderSpecific(eps,EPS_CLASSID,2);
979 1 : PetscCheckSameComm(svd,1,eps,2);
980 1 : PetscTryMethod(svd,"SVDCyclicSetEPS_C",(SVD,EPS),(svd,eps));
981 1 : PetscFunctionReturn(PETSC_SUCCESS);
982 : }
983 :
984 35 : static PetscErrorCode SVDCyclicGetEPS_Cyclic(SVD svd,EPS *eps)
985 : {
986 35 : SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
987 :
988 35 : PetscFunctionBegin;
989 35 : if (!cyclic->eps) {
990 35 : PetscCall(EPSCreate(PetscObjectComm((PetscObject)svd),&cyclic->eps));
991 35 : PetscCall(PetscObjectIncrementTabLevel((PetscObject)cyclic->eps,(PetscObject)svd,1));
992 35 : PetscCall(EPSSetOptionsPrefix(cyclic->eps,((PetscObject)svd)->prefix));
993 35 : PetscCall(EPSAppendOptionsPrefix(cyclic->eps,"svd_cyclic_"));
994 35 : PetscCall(PetscObjectSetOptions((PetscObject)cyclic->eps,((PetscObject)svd)->options));
995 35 : PetscCall(EPSSetWhichEigenpairs(cyclic->eps,EPS_LARGEST_REAL));
996 35 : PetscCall(EPSMonitorSet(cyclic->eps,EPSMonitor_Cyclic,svd,NULL));
997 : }
998 35 : *eps = cyclic->eps;
999 35 : PetscFunctionReturn(PETSC_SUCCESS);
1000 : }
1001 :
1002 : /*@
1003 : SVDCyclicGetEPS - Retrieve the eigensolver object (EPS) associated
1004 : to the singular value solver.
1005 :
1006 : Collective
1007 :
1008 : Input Parameter:
1009 : . svd - singular value solver
1010 :
1011 : Output Parameter:
1012 : . eps - the eigensolver object
1013 :
1014 : Level: advanced
1015 :
1016 : .seealso: SVDCyclicSetEPS()
1017 : @*/
1018 35 : PetscErrorCode SVDCyclicGetEPS(SVD svd,EPS *eps)
1019 : {
1020 35 : PetscFunctionBegin;
1021 35 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
1022 35 : PetscAssertPointer(eps,2);
1023 35 : PetscUseMethod(svd,"SVDCyclicGetEPS_C",(SVD,EPS*),(svd,eps));
1024 35 : PetscFunctionReturn(PETSC_SUCCESS);
1025 : }
1026 :
1027 1 : static PetscErrorCode SVDView_Cyclic(SVD svd,PetscViewer viewer)
1028 : {
1029 1 : SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
1030 1 : PetscBool isascii;
1031 :
1032 1 : PetscFunctionBegin;
1033 1 : PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
1034 1 : if (isascii) {
1035 1 : if (!cyclic->eps) PetscCall(SVDCyclicGetEPS(svd,&cyclic->eps));
1036 1 : PetscCall(PetscViewerASCIIPrintf(viewer," %s matrix\n",cyclic->explicitmatrix?"explicit":"implicit"));
1037 1 : PetscCall(PetscViewerASCIIPushTab(viewer));
1038 1 : PetscCall(EPSView(cyclic->eps,viewer));
1039 1 : PetscCall(PetscViewerASCIIPopTab(viewer));
1040 : }
1041 1 : PetscFunctionReturn(PETSC_SUCCESS);
1042 : }
1043 :
1044 37 : static PetscErrorCode SVDReset_Cyclic(SVD svd)
1045 : {
1046 37 : SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
1047 :
1048 37 : PetscFunctionBegin;
1049 37 : PetscCall(EPSReset(cyclic->eps));
1050 37 : PetscCall(MatDestroy(&cyclic->C));
1051 37 : PetscCall(MatDestroy(&cyclic->D));
1052 37 : PetscFunctionReturn(PETSC_SUCCESS);
1053 : }
1054 :
1055 36 : static PetscErrorCode SVDDestroy_Cyclic(SVD svd)
1056 : {
1057 36 : SVD_CYCLIC *cyclic = (SVD_CYCLIC*)svd->data;
1058 :
1059 36 : PetscFunctionBegin;
1060 36 : PetscCall(EPSDestroy(&cyclic->eps));
1061 36 : PetscCall(PetscFree(svd->data));
1062 36 : PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicSetEPS_C",NULL));
1063 36 : PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicGetEPS_C",NULL));
1064 36 : PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicSetExplicitMatrix_C",NULL));
1065 36 : PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicGetExplicitMatrix_C",NULL));
1066 36 : PetscFunctionReturn(PETSC_SUCCESS);
1067 : }
1068 :
1069 36 : SLEPC_EXTERN PetscErrorCode SVDCreate_Cyclic(SVD svd)
1070 : {
1071 36 : SVD_CYCLIC *cyclic;
1072 :
1073 36 : PetscFunctionBegin;
1074 36 : PetscCall(PetscNew(&cyclic));
1075 36 : svd->data = (void*)cyclic;
1076 36 : svd->ops->solve = SVDSolve_Cyclic;
1077 36 : svd->ops->solveg = SVDSolve_Cyclic;
1078 36 : svd->ops->solveh = SVDSolve_Cyclic;
1079 36 : svd->ops->setup = SVDSetUp_Cyclic;
1080 36 : svd->ops->setfromoptions = SVDSetFromOptions_Cyclic;
1081 36 : svd->ops->destroy = SVDDestroy_Cyclic;
1082 36 : svd->ops->reset = SVDReset_Cyclic;
1083 36 : svd->ops->view = SVDView_Cyclic;
1084 36 : svd->ops->computevectors = SVDComputeVectors_Cyclic;
1085 36 : PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicSetEPS_C",SVDCyclicSetEPS_Cyclic));
1086 36 : PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicGetEPS_C",SVDCyclicGetEPS_Cyclic));
1087 36 : PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicSetExplicitMatrix_C",SVDCyclicSetExplicitMatrix_Cyclic));
1088 36 : PetscCall(PetscObjectComposeFunction((PetscObject)svd,"SVDCyclicGetExplicitMatrix_C",SVDCyclicGetExplicitMatrix_Cyclic));
1089 36 : PetscFunctionReturn(PETSC_SUCCESS);
1090 : }
|