Actual source code: svdelemen.cxx
slepc-3.15.0 2021-03-31
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2021, 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: {
24: SVD_Elemental *ctx = (SVD_Elemental*)svd->data;
25: PetscInt M,N;
28: SVDCheckStandard(svd);
29: MatGetSize(svd->A,&M,&N);
30: if (M!=N) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Not implemented for rectangular matrices");
31: svd->ncv = N;
32: if (svd->mpd!=PETSC_DEFAULT) { PetscInfo(svd,"Warning: parameter mpd ignored\n"); }
33: if (svd->max_it==PETSC_DEFAULT) svd->max_it = 1;
34: svd->leftbasis = PETSC_TRUE;
35: SVDCheckUnsupported(svd,SVD_FEATURE_STOPPING);
36: SVDAllocateSolution(svd,0);
38: /* convert matrix */
39: MatDestroy(&ctx->Ae);
40: MatConvert(svd->OP,MATELEMENTAL,MAT_INITIAL_MATRIX,&ctx->Ae);
41: return(0);
42: }
44: PetscErrorCode SVDSolve_Elemental(SVD svd)
45: {
47: SVD_Elemental *ctx = (SVD_Elemental*)svd->data;
48: Mat A = ctx->Ae,Z,Q,U,V;
49: Mat_Elemental *a = (Mat_Elemental*)A->data,*q,*z;
50: PetscInt i,rrank,ridx,erow;
53: El::DistMatrix<PetscReal,El::STAR,El::VC> sigma(*a->grid);
54: MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&Z);
55: MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&Q);
56: z = (Mat_Elemental*)Z->data;
57: q = (Mat_Elemental*)Q->data;
59: El::SVD(*a->emat,*z->emat,sigma,*q->emat);
61: for (i=0;i<svd->ncv;i++) {
62: P2RO(A,1,i,&rrank,&ridx);
63: RO2E(A,1,rrank,ridx,&erow);
64: svd->sigma[i] = sigma.Get(erow,0);
65: }
66: BVGetMat(svd->U,&U);
67: MatConvert(Z,MATDENSE,MAT_REUSE_MATRIX,&U);
68: BVRestoreMat(svd->U,&U);
69: MatDestroy(&Z);
70: BVGetMat(svd->V,&V);
71: MatConvert(Q,MATDENSE,MAT_REUSE_MATRIX,&V);
72: BVRestoreMat(svd->V,&V);
73: MatDestroy(&Q);
75: svd->nconv = svd->ncv;
76: svd->its = 1;
77: svd->reason = SVD_CONVERGED_TOL;
78: return(0);
79: }
81: PetscErrorCode SVDDestroy_Elemental(SVD svd)
82: {
86: PetscFree(svd->data);
87: return(0);
88: }
90: PetscErrorCode SVDReset_Elemental(SVD svd)
91: {
93: SVD_Elemental *ctx = (SVD_Elemental*)svd->data;
96: MatDestroy(&ctx->Ae);
97: return(0);
98: }
100: SLEPC_EXTERN PetscErrorCode SVDCreate_Elemental(SVD svd)
101: {
103: SVD_Elemental *ctx;
106: PetscNewLog(svd,&ctx);
107: svd->data = (void*)ctx;
109: svd->ops->solve = SVDSolve_Elemental;
110: svd->ops->setup = SVDSetUp_Elemental;
111: svd->ops->destroy = SVDDestroy_Elemental;
112: svd->ops->reset = SVDReset_Elemental;
113: return(0);
114: }