Actual source code: svdelemen.cxx
slepc-3.17.0 2022-03-31
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: MatGetSize(svd->A,&M,&N);
29: svd->ncv = N;
30: if (svd->mpd!=PETSC_DEFAULT) PetscInfo(svd,"Warning: parameter mpd ignored\n");
31: if (svd->max_it==PETSC_DEFAULT) svd->max_it = 1;
32: svd->leftbasis = PETSC_TRUE;
33: SVDCheckUnsupported(svd,SVD_FEATURE_STOPPING);
34: SVDAllocateSolution(svd,0);
36: /* convert matrix */
37: MatDestroy(&ctx->Ae);
38: MatConvert(svd->OP,MATELEMENTAL,MAT_INITIAL_MATRIX,&ctx->Ae);
39: PetscFunctionReturn(0);
40: }
42: PetscErrorCode SVDSolve_Elemental(SVD svd)
43: {
44: SVD_Elemental *ctx = (SVD_Elemental*)svd->data;
45: Mat A = ctx->Ae,Z,Q,U,V;
46: Mat_Elemental *a = (Mat_Elemental*)A->data,*q,*z;
47: PetscInt i,rrank,ridx,erow;
49: El::DistMatrix<PetscReal,El::STAR,El::VC> sigma(*a->grid);
50: MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&Z);
51: MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&Q);
52: z = (Mat_Elemental*)Z->data;
53: q = (Mat_Elemental*)Q->data;
55: El::SVD(*a->emat,*z->emat,sigma,*q->emat);
57: for (i=0;i<svd->ncv;i++) {
58: P2RO(A,1,i,&rrank,&ridx);
59: RO2E(A,1,rrank,ridx,&erow);
60: svd->sigma[i] = sigma.Get(erow,0);
61: }
62: BVGetMat(svd->U,&U);
63: MatConvert(Z,MATDENSE,MAT_REUSE_MATRIX,&U);
64: BVRestoreMat(svd->U,&U);
65: MatDestroy(&Z);
66: BVGetMat(svd->V,&V);
67: MatConvert(Q,MATDENSE,MAT_REUSE_MATRIX,&V);
68: BVRestoreMat(svd->V,&V);
69: MatDestroy(&Q);
71: svd->nconv = svd->ncv;
72: svd->its = 1;
73: svd->reason = SVD_CONVERGED_TOL;
74: PetscFunctionReturn(0);
75: }
77: PetscErrorCode SVDDestroy_Elemental(SVD svd)
78: {
79: PetscFree(svd->data);
80: PetscFunctionReturn(0);
81: }
83: PetscErrorCode SVDReset_Elemental(SVD svd)
84: {
85: SVD_Elemental *ctx = (SVD_Elemental*)svd->data;
87: MatDestroy(&ctx->Ae);
88: PetscFunctionReturn(0);
89: }
91: SLEPC_EXTERN PetscErrorCode SVDCreate_Elemental(SVD svd)
92: {
93: SVD_Elemental *ctx;
95: PetscNewLog(svd,&ctx);
96: svd->data = (void*)ctx;
98: svd->ops->solve = SVDSolve_Elemental;
99: svd->ops->setup = SVDSetUp_Elemental;
100: svd->ops->destroy = SVDDestroy_Elemental;
101: svd->ops->reset = SVDReset_Elemental;
102: PetscFunctionReturn(0);
103: }