Actual source code: svdelemen.cxx
slepc-3.14.0 2020-09-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2020, 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: SVDMatGetSize(svd,&M,&N);
29: if (M!=N) SETERRQ(PetscObjectComm((PetscObject)svd),PETSC_ERR_SUP,"Not implemented for rectangular matrices");
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: {
46: SVD_Elemental *ctx = (SVD_Elemental*)svd->data;
47: Mat A = ctx->Ae,Z,Q,U,V;
48: Mat_Elemental *a = (Mat_Elemental*)A->data,*q,*z;
49: PetscInt i,rrank,ridx,erow;
52: El::DistMatrix<PetscReal,El::STAR,El::VC> sigma(*a->grid);
53: MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&Z);
54: MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&Q);
55: z = (Mat_Elemental*)Z->data;
56: q = (Mat_Elemental*)Q->data;
58: El::SVD(*a->emat,*z->emat,sigma,*q->emat);
60: for (i=0;i<svd->ncv;i++) {
61: P2RO(A,1,i,&rrank,&ridx);
62: RO2E(A,1,rrank,ridx,&erow);
63: svd->sigma[i] = sigma.Get(erow,0);
64: }
65: BVGetMat(svd->U,&U);
66: MatConvert(Z,MATDENSE,MAT_REUSE_MATRIX,&U);
67: BVRestoreMat(svd->U,&U);
68: MatDestroy(&Z);
69: BVGetMat(svd->V,&V);
70: MatConvert(Q,MATDENSE,MAT_REUSE_MATRIX,&V);
71: BVRestoreMat(svd->V,&V);
72: MatDestroy(&Q);
74: svd->nconv = svd->ncv;
75: svd->its = 1;
76: svd->reason = SVD_CONVERGED_TOL;
77: return(0);
78: }
80: PetscErrorCode SVDDestroy_Elemental(SVD svd)
81: {
85: PetscFree(svd->data);
86: return(0);
87: }
89: PetscErrorCode SVDReset_Elemental(SVD svd)
90: {
92: SVD_Elemental *ctx = (SVD_Elemental*)svd->data;
95: MatDestroy(&ctx->Ae);
96: return(0);
97: }
99: SLEPC_EXTERN PetscErrorCode SVDCreate_Elemental(SVD svd)
100: {
102: SVD_Elemental *ctx;
105: PetscNewLog(svd,&ctx);
106: svd->data = (void*)ctx;
108: svd->ops->solve = SVDSolve_Elemental;
109: svd->ops->setup = SVDSetUp_Elemental;
110: svd->ops->destroy = SVDDestroy_Elemental;
111: svd->ops->reset = SVDReset_Elemental;
112: return(0);
113: }