Actual source code: matstruct.c
slepc-3.22.2 2024-12-02
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: */
11: #include <slepc/private/slepcimpl.h>
13: /*@
14: MatCreateBSE - Create a matrix that can be used to define a structured eigenvalue
15: problem of type BSE (Bethe-Salpeter Equation).
17: Collective
19: Input Parameters:
20: + R - matrix for the diagonal block (resonant)
21: - C - matrix for the off-diagonal block (coupling)
23: Output Parameter:
24: . H - the resulting matrix
26: Notes:
27: The resulting matrix has the block form H = [ R C; -C^H -R^T ], where R is assumed
28: to be (complex) Hermitian and C complex symmetric. Note that this function does
29: not check these properties, so if the matrices provided by the user do not satisfy
30: them, then the solver will not behave as expected.
32: The obtained matrix can be used as an input matrix to EPS eigensolvers via
33: EPSSetOperators() for the case that the problem type is EPS_BSE. Note that the user
34: cannot just build a matrix with the required structure, it must be done via this
35: function.
37: In the current implementation, H is a MATNEST matrix, where R and C form the top
38: block row, while the bottom block row is composed of matrices of type
39: MATTRANSPOSEVIRTUAL and MATHERMITIANTRANSPOSEVIRTUAL scaled by -1.
41: Level: intermediate
43: .seealso: MatCreateNest(), EPSSetOperators(), EPSSetProblemType()
44: @*/
45: PetscErrorCode MatCreateBSE(Mat R,Mat C,Mat *H)
46: {
47: PetscInt Mr,Mc,Nr,Nc,mr,mc,nr,nc;
48: Mat block[4] = { R, C, NULL, NULL };
49: SlepcMatStruct mctx;
51: PetscFunctionBegin;
54: PetscCheckSameTypeAndComm(R,1,C,2);
55: PetscAssertPointer(H,3);
57: /* check sizes */
58: PetscCall(MatGetSize(R,&Mr,&Nr));
59: PetscCall(MatGetLocalSize(R,&mr,&nr));
60: PetscCall(MatGetSize(C,&Mc,&Nc));
61: PetscCall(MatGetLocalSize(C,&mc,&nc));
62: PetscCheck(Mc==Mr && mc==mr,PetscObjectComm((PetscObject)R),PETSC_ERR_ARG_INCOMP,"Incompatible row dimensions");
63: PetscCheck(Nc==Nr && nc==nr,PetscObjectComm((PetscObject)R),PETSC_ERR_ARG_INCOMP,"Incompatible column dimensions");
65: /* bottom block row */
66: PetscCall(MatCreateHermitianTranspose(C,&block[2]));
67: PetscCall(MatScale(block[2],-1.0));
68: PetscCall(MatCreateTranspose(R,&block[3]));
69: PetscCall(MatScale(block[3],-1.0));
71: /* create nest matrix and compose context */
72: PetscCall(MatCreateNest(PetscObjectComm((PetscObject)R),2,NULL,2,NULL,block,H));
73: PetscCall(MatDestroy(&block[2]));
74: PetscCall(MatDestroy(&block[3]));
76: PetscCall(PetscNew(&mctx));
77: mctx->cookie = SLEPC_MAT_STRUCT_BSE;
78: PetscCall(PetscObjectContainerCompose((PetscObject)*H,"SlepcMatStruct",mctx,PetscContainerUserDestroyDefault));
79: PetscFunctionReturn(PETSC_SUCCESS);
80: }