Actual source code: svdelemen.cxx
slepc-3.18.0 2022-10-01
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
10: /*
11: This file implements a wrapper to the Elemental SVD solver
12: */
14: #include <slepc/private/svdimpl.h>
15: #include <petsc/private/petscelemental.h>
17: typedef struct {
18: Mat Ae; /* converted matrix */
19: } SVD_Elemental;
21: PetscErrorCode SVDSetUp_Elemental(SVD svd)
22: {
23: SVD_Elemental *ctx = (SVD_Elemental*)svd->data;
24: PetscInt M,N;
26: SVDCheckStandard(svd);
27: SVDCheckDefinite(svd);
28: MatGetSize(svd->A,&M,&N);
30: svd->ncv = N;
31: if (svd->mpd!=PETSC_DEFAULT) PetscInfo(svd,"Warning: parameter mpd ignored\n");
32: if (svd->max_it==PETSC_DEFAULT) svd->max_it = 1;
33: svd->leftbasis = PETSC_TRUE;
34: SVDCheckUnsupported(svd,SVD_FEATURE_STOPPING);
35: SVDAllocateSolution(svd,0);
37: /* convert matrix */
38: MatDestroy(&ctx->Ae);
39: MatConvert(svd->OP,MATELEMENTAL,MAT_INITIAL_MATRIX,&ctx->Ae);
40: return 0;
41: }
43: PetscErrorCode SVDSolve_Elemental(SVD svd)
44: {
45: SVD_Elemental *ctx = (SVD_Elemental*)svd->data;
46: Mat A = ctx->Ae,Z,Q,U,V;
47: Mat_Elemental *a = (Mat_Elemental*)A->data,*q,*z;
48: PetscInt i,rrank,ridx,erow;
50: El::DistMatrix<PetscReal,El::STAR,El::VC> sigma(*a->grid);
51: MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&Z);
52: MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&Q);
53: z = (Mat_Elemental*)Z->data;
54: q = (Mat_Elemental*)Q->data;
56: El::SVD(*a->emat,*z->emat,sigma,*q->emat);
58: for (i=0;i<svd->ncv;i++) {
59: P2RO(A,1,i,&rrank,&ridx);
60: RO2E(A,1,rrank,ridx,&erow);
61: svd->sigma[i] = sigma.Get(erow,0);
62: }
63: BVGetMat(svd->U,&U);
64: MatConvert(Z,MATDENSE,MAT_REUSE_MATRIX,&U);
65: BVRestoreMat(svd->U,&U);
66: MatDestroy(&Z);
67: BVGetMat(svd->V,&V);
68: MatConvert(Q,MATDENSE,MAT_REUSE_MATRIX,&V);
69: BVRestoreMat(svd->V,&V);
70: MatDestroy(&Q);
72: svd->nconv = svd->ncv;
73: svd->its = 1;
74: svd->reason = SVD_CONVERGED_TOL;
75: return 0;
76: }
78: PetscErrorCode SVDDestroy_Elemental(SVD svd)
79: {
80: PetscFree(svd->data);
81: return 0;
82: }
84: PetscErrorCode SVDReset_Elemental(SVD svd)
85: {
86: SVD_Elemental *ctx = (SVD_Elemental*)svd->data;
88: MatDestroy(&ctx->Ae);
89: return 0;
90: }
92: SLEPC_EXTERN PetscErrorCode SVDCreate_Elemental(SVD svd)
93: {
94: SVD_Elemental *ctx;
96: PetscNew(&ctx);
97: svd->data = (void*)ctx;
99: svd->ops->solve = SVDSolve_Elemental;
100: svd->ops->setup = SVDSetUp_Elemental;
101: svd->ops->destroy = SVDDestroy_Elemental;
102: svd->ops->reset = SVDReset_Elemental;
103: return 0;
104: }