Actual source code: qeplin.c
slepc-3.16.1 2021-11-17
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: Various types of linearization for quadratic eigenvalue problem
12: */
14: #include <slepc/private/pepimpl.h>
15: #include "linear.h"
17: /*
18: Given the quadratic problem (l^2*M + l*C + K)*x = 0 the linearization is
19: A*z = l*B*z for z = [ x ] and A,B defined as follows:
20: [ l*x ]
22: N:
23: A = [-bK aI ] B = [ aI+bC bM ]
24: [-aK -aC+bI ] [ bI aM ]
26: S:
27: A = [ bK aK ] B = [ aK-bC -bM ]
28: [ aK aC-bM ] [-bM -aM ]
30: H:
31: A = [ aK -bK ] B = [ bM aK+bC ]
32: [ aC+bM aK ] [-aM bM ]
33: */
35: /* - - - N - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
37: PetscErrorCode MatCreateExplicit_Linear_NA(MPI_Comm comm,PEP_LINEAR *ctx,Mat *A)
38: {
40: PetscInt M,N,m,n,i,Istart,Iend;
41: Mat Id,T=NULL;
42: PetscReal a=ctx->alpha,b=ctx->beta;
43: PetscScalar scalt=1.0;
46: MatGetSize(ctx->M,&M,&N);
47: MatGetLocalSize(ctx->M,&m,&n);
48: MatCreate(PetscObjectComm((PetscObject)ctx->M),&Id);
49: MatSetSizes(Id,m,n,M,N);
50: MatSetFromOptions(Id);
51: MatSetUp(Id);
52: MatGetOwnershipRange(Id,&Istart,&Iend);
53: for (i=Istart;i<Iend;i++) {
54: MatSetValue(Id,i,i,1.0,INSERT_VALUES);
55: }
56: MatAssemblyBegin(Id,MAT_FINAL_ASSEMBLY);
57: MatAssemblyEnd(Id,MAT_FINAL_ASSEMBLY);
58: if (a!=0.0 && b!=0.0) {
59: MatDuplicate(ctx->C,MAT_COPY_VALUES,&T);
60: MatScale(T,-a*ctx->dsfactor*ctx->sfactor);
61: MatShift(T,b);
62: } else {
63: if (a==0.0) { T = Id; scalt = b; }
64: else { T = ctx->C; scalt = -a*ctx->dsfactor*ctx->sfactor; }
65: }
66: MatCreateTile(-b*ctx->dsfactor,ctx->K,a,Id,-ctx->dsfactor*a,ctx->K,scalt,T,A);
67: MatDestroy(&Id);
68: if (a!=0.0 && b!=0.0) {
69: MatDestroy(&T);
70: }
71: return(0);
72: }
74: PetscErrorCode MatCreateExplicit_Linear_NB(MPI_Comm comm,PEP_LINEAR *ctx,Mat *B)
75: {
77: PetscInt M,N,m,n,i,Istart,Iend;
78: Mat Id,T=NULL;
79: PetscReal a=ctx->alpha,b=ctx->beta;
80: PetscScalar scalt=1.0;
83: MatGetSize(ctx->M,&M,&N);
84: MatGetLocalSize(ctx->M,&m,&n);
85: MatCreate(PetscObjectComm((PetscObject)ctx->M),&Id);
86: MatSetSizes(Id,m,n,M,N);
87: MatSetFromOptions(Id);
88: MatSetUp(Id);
89: MatGetOwnershipRange(Id,&Istart,&Iend);
90: for (i=Istart;i<Iend;i++) {
91: MatSetValue(Id,i,i,1.0,INSERT_VALUES);
92: }
93: MatAssemblyBegin(Id,MAT_FINAL_ASSEMBLY);
94: MatAssemblyEnd(Id,MAT_FINAL_ASSEMBLY);
95: if (a!=0.0 && b!=0.0) {
96: MatDuplicate(ctx->C,MAT_COPY_VALUES,&T);
97: MatScale(T,b*ctx->dsfactor*ctx->sfactor);
98: MatShift(T,a);
99: } else {
100: if (b==0.0) { T = Id; scalt = a; }
101: else { T = ctx->C; scalt = b*ctx->dsfactor*ctx->sfactor; }
102: }
103: MatCreateTile(scalt,T,b*ctx->dsfactor*ctx->sfactor*ctx->sfactor,ctx->M,b,Id,a*ctx->sfactor*ctx->sfactor*ctx->dsfactor,ctx->M,B);
104: MatDestroy(&Id);
105: if (a!=0.0 && b!=0.0) {
106: MatDestroy(&T);
107: }
108: return(0);
109: }
111: /* - - - S - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
113: PetscErrorCode MatCreateExplicit_Linear_SA(MPI_Comm comm,PEP_LINEAR *ctx,Mat *A)
114: {
116: Mat T=NULL;
117: PetscScalar scalt=1.0;
118: PetscReal a=ctx->alpha,b=ctx->beta;
121: if (a!=0.0 && b!=0.0) {
122: MatDuplicate(ctx->C,MAT_COPY_VALUES,&T);
123: MatScale(T,a*ctx->dsfactor*ctx->sfactor);
124: MatAXPY(T,-b*ctx->dsfactor*ctx->sfactor*ctx->sfactor,ctx->M,UNKNOWN_NONZERO_PATTERN);
125: } else {
126: if (a==0.0) { T = ctx->M; scalt = -b*ctx->dsfactor*ctx->sfactor*ctx->sfactor; }
127: else { T = ctx->C; scalt = a*ctx->dsfactor*ctx->sfactor; }
128: }
129: MatCreateTile(b*ctx->dsfactor,ctx->K,a*ctx->dsfactor,ctx->K,a*ctx->dsfactor,ctx->K,scalt,T,A);
130: if (a!=0.0 && b!=0.0) {
131: MatDestroy(&T);
132: }
133: return(0);
134: }
136: PetscErrorCode MatCreateExplicit_Linear_SB(MPI_Comm comm,PEP_LINEAR *ctx,Mat *B)
137: {
139: Mat T=NULL;
140: PetscScalar scalt=1.0;
141: PetscReal a=ctx->alpha,b=ctx->beta;
144: if (a!=0.0 && b!=0.0) {
145: MatDuplicate(ctx->C,MAT_COPY_VALUES,&T);
146: MatScale(T,-b*ctx->dsfactor*ctx->sfactor);
147: MatAXPY(T,a*ctx->dsfactor,ctx->K,UNKNOWN_NONZERO_PATTERN);
148: } else {
149: if (b==0.0) { T = ctx->K; scalt = a*ctx->dsfactor; }
150: else { T = ctx->C; scalt = -b*ctx->dsfactor*ctx->sfactor; }
151: }
152: MatCreateTile(scalt,T,-b*ctx->dsfactor*ctx->sfactor*ctx->sfactor,ctx->M,-b*ctx->dsfactor*ctx->sfactor*ctx->sfactor,ctx->M,-a*ctx->dsfactor*ctx->sfactor*ctx->sfactor,ctx->M,B);
153: if (a!=0.0 && b!=0.0) {
154: MatDestroy(&T);
155: }
156: return(0);
157: }
159: /* - - - H - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
161: PetscErrorCode MatCreateExplicit_Linear_HA(MPI_Comm comm,PEP_LINEAR *ctx,Mat *A)
162: {
164: Mat T=NULL;
165: PetscScalar scalt=1.0;
166: PetscReal a=ctx->alpha,b=ctx->beta;
169: if (a!=0.0 && b!=0.0) {
170: MatDuplicate(ctx->C,MAT_COPY_VALUES,&T);
171: MatScale(T,a*ctx->dsfactor*ctx->sfactor);
172: MatAXPY(T,b*ctx->dsfactor*ctx->sfactor*ctx->sfactor,ctx->M,UNKNOWN_NONZERO_PATTERN);
173: } else {
174: if (a==0.0) { T = ctx->M; scalt = b*ctx->dsfactor*ctx->sfactor*ctx->sfactor; }
175: else { T = ctx->C; scalt = a*ctx->dsfactor*ctx->sfactor; }
176: }
177: MatCreateTile(a*ctx->dsfactor,ctx->K,-b*ctx->dsfactor,ctx->K,scalt,T,a*ctx->dsfactor,ctx->K,A);
178: if (a!=0.0 && b!=0.0) {
179: MatDestroy(&T);
180: }
181: return(0);
182: }
184: PetscErrorCode MatCreateExplicit_Linear_HB(MPI_Comm comm,PEP_LINEAR *ctx,Mat *B)
185: {
187: Mat T=NULL;
188: PetscScalar scalt=1.0;
189: PetscReal a=ctx->alpha,b=ctx->beta;
192: if (a!=0.0 && b!=0.0) {
193: MatDuplicate(ctx->C,MAT_COPY_VALUES,&T);
194: MatScale(T,b*ctx->dsfactor*ctx->sfactor);
195: MatAXPY(T,a*ctx->dsfactor,ctx->K,UNKNOWN_NONZERO_PATTERN);
196: } else {
197: if (b==0.0) { T = ctx->K; scalt = a*ctx->dsfactor; }
198: else { T = ctx->C; scalt = b*ctx->dsfactor*ctx->sfactor; }
199: }
200: MatCreateTile(b*ctx->dsfactor*ctx->sfactor*ctx->sfactor,ctx->M,scalt,T,-a*ctx->dsfactor*ctx->sfactor*ctx->sfactor,ctx->M,b*ctx->dsfactor*ctx->sfactor*ctx->sfactor,ctx->M,B);
201: if (a!=0.0 && b!=0.0) {
202: MatDestroy(&T);
203: }
204: return(0);
205: }