Actual source code: qeplin.c

slepc-3.7.3 2016-09-29
Report Typos and Errors
  1: /*

  3:    Various types of linearization for quadratic eigenvalue problem.

  5:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  6:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  7:    Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain

  9:    This file is part of SLEPc.

 11:    SLEPc is free software: you can redistribute it and/or modify it under  the
 12:    terms of version 3 of the GNU Lesser General Public License as published by
 13:    the Free Software Foundation.

 15:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 16:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 17:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 18:    more details.

 20:    You  should have received a copy of the GNU Lesser General  Public  License
 21:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 22:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 23: */

 25: #include <slepc/private/pepimpl.h>
 26:  #include linearp.h

 28: /*
 29:     Given the quadratic problem (l^2*M + l*C + K)*x = 0 the linearization is
 30:     A*z = l*B*z for z = [  x  ] and A,B defined as follows:
 31:                         [ l*x ]

 33:             N1:
 34:                      A = [  0   I ]     B = [ I  0 ]
 35:                          [ -K  -C ]         [ 0  M ]

 37:             N2:
 38:                      A = [ -K   0 ]     B = [ C  M ]
 39:                          [  0   I ]         [ I  0 ]

 41:             S1:
 42:                      A = [  0  -K ]     B = [-K  0 ]
 43:                          [ -K  -C ]         [ 0  M ]

 45:             S2:
 46:                      A = [ -K   0 ]     B = [ C  M ]
 47:                          [  0   M ]         [ M  0 ]

 49:             H1:
 50:                      A = [  K   0 ]     B = [ 0  K ]
 51:                          [  C   K ]         [-M  0 ]

 53:             H2:
 54:                      A = [  0  -K ]     B = [ M  C ]
 55:                          [  M   0 ]         [ 0  M ]
 56:  */

 58: /* - - - N1 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 62: PetscErrorCode MatCreateExplicit_Linear_N1A(MPI_Comm comm,PEP_LINEAR *ctx,Mat *A)
 63: {
 65:   PetscInt       M,N,m,n,i,Istart,Iend;
 66:   Mat            Id;

 69:   MatGetSize(ctx->M,&M,&N);
 70:   MatGetLocalSize(ctx->M,&m,&n);
 71:   MatCreate(PetscObjectComm((PetscObject)ctx->M),&Id);
 72:   MatSetSizes(Id,m,n,M,N);
 73:   MatSetFromOptions(Id);
 74:   MatSetUp(Id);
 75:   MatGetOwnershipRange(Id,&Istart,&Iend);
 76:   for (i=Istart;i<Iend;i++) {
 77:     MatSetValue(Id,i,i,1.0,INSERT_VALUES);
 78:   }
 79:   MatAssemblyBegin(Id,MAT_FINAL_ASSEMBLY);
 80:   MatAssemblyEnd(Id,MAT_FINAL_ASSEMBLY);
 81:   SlepcMatTile(0.0,Id,1.0,Id,-ctx->dsfactor,ctx->K,-ctx->sfactor*ctx->dsfactor,ctx->C,A);
 82:   MatDestroy(&Id);
 83:   return(0);
 84: }

 88: PetscErrorCode MatCreateExplicit_Linear_N1B(MPI_Comm comm,PEP_LINEAR *ctx,Mat *B)
 89: {
 91:   PetscInt       M,N,m,n,i,Istart,Iend;
 92:   Mat            Id;

 95:   MatGetSize(ctx->M,&M,&N);
 96:   MatGetLocalSize(ctx->M,&m,&n);
 97:   MatCreate(PetscObjectComm((PetscObject)ctx->M),&Id);
 98:   MatSetSizes(Id,m,n,M,N);
 99:   MatSetFromOptions(Id);
100:   MatSetUp(Id);
101:   MatGetOwnershipRange(Id,&Istart,&Iend);
102:   for (i=Istart;i<Iend;i++) {
103:     MatSetValue(Id,i,i,1.0,INSERT_VALUES);
104:   }
105:   MatAssemblyBegin(Id,MAT_FINAL_ASSEMBLY);
106:   MatAssemblyEnd(Id,MAT_FINAL_ASSEMBLY);
107:   SlepcMatTile(1.0,Id,0.0,Id,0.0,Id,ctx->sfactor*ctx->sfactor*ctx->dsfactor,ctx->M,B);
108:   MatDestroy(&Id);
109:   return(0);
110: }

112: /* - - - N2 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

116: PetscErrorCode MatCreateExplicit_Linear_N2A(MPI_Comm comm,PEP_LINEAR *ctx,Mat *A)
117: {
119:   PetscInt       M,N,m,n,i,Istart,Iend;
120:   Mat            Id;

123:   MatGetSize(ctx->M,&M,&N);
124:   MatGetLocalSize(ctx->M,&m,&n);
125:   MatCreate(PetscObjectComm((PetscObject)ctx->M),&Id);
126:   MatSetSizes(Id,m,n,M,N);
127:   MatSetFromOptions(Id);
128:   MatSetUp(Id);
129:   MatGetOwnershipRange(Id,&Istart,&Iend);
130:   for (i=Istart;i<Iend;i++) {
131:     MatSetValue(Id,i,i,1.0,INSERT_VALUES);
132:   }
133:   MatAssemblyBegin(Id,MAT_FINAL_ASSEMBLY);
134:   MatAssemblyEnd(Id,MAT_FINAL_ASSEMBLY);
135:   SlepcMatTile(-1.0,ctx->K,0.0,Id,0.0,Id,1.0,Id,A);
136:   MatDestroy(&Id);
137:   return(0);
138: }

142: PetscErrorCode MatCreateExplicit_Linear_N2B(MPI_Comm comm,PEP_LINEAR *ctx,Mat *B)
143: {
145:   PetscInt       M,N,m,n,i,Istart,Iend;
146:   Mat            Id;

149:   MatGetSize(ctx->M,&M,&N);
150:   MatGetLocalSize(ctx->M,&m,&n);
151:   MatCreate(PetscObjectComm((PetscObject)ctx->M),&Id);
152:   MatSetSizes(Id,m,n,M,N);
153:   MatSetFromOptions(Id);
154:   MatSetUp(Id);
155:   MatGetOwnershipRange(Id,&Istart,&Iend);
156:   for (i=Istart;i<Iend;i++) {
157:     MatSetValue(Id,i,i,1.0,INSERT_VALUES);
158:   }
159:   MatAssemblyBegin(Id,MAT_FINAL_ASSEMBLY);
160:   MatAssemblyEnd(Id,MAT_FINAL_ASSEMBLY);
161:   SlepcMatTile(ctx->sfactor,ctx->C,ctx->sfactor*ctx->sfactor,ctx->M,1.0,Id,0.0,Id,B);
162:   MatDestroy(&Id);
163:   return(0);
164: }

166: /* - - - S1 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

170: PetscErrorCode MatCreateExplicit_Linear_S1A(MPI_Comm comm,PEP_LINEAR *ctx,Mat *A)
171: {

175:   SlepcMatTile(0.0,ctx->K,-1.0,ctx->K,-1.0,ctx->K,-ctx->sfactor,ctx->C,A);
176:   return(0);
177: }

181: PetscErrorCode MatCreateExplicit_Linear_S1B(MPI_Comm comm,PEP_LINEAR *ctx,Mat *B)
182: {

186:   SlepcMatTile(-1.0,ctx->K,0.0,ctx->M,0.0,ctx->M,ctx->sfactor*ctx->sfactor,ctx->M,B);
187:   return(0);
188: }

190: /* - - - S2 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

194: PetscErrorCode MatCreateExplicit_Linear_S2A(MPI_Comm comm,PEP_LINEAR *ctx,Mat *A)
195: {

199:   SlepcMatTile(-1.0,ctx->K,0.0,ctx->M,0.0,ctx->M,ctx->sfactor*ctx->sfactor,ctx->M,A);
200:   return(0);
201: }

205: PetscErrorCode MatCreateExplicit_Linear_S2B(MPI_Comm comm,PEP_LINEAR *ctx,Mat *B)
206: {

210:   SlepcMatTile(ctx->sfactor,ctx->C,ctx->sfactor*ctx->sfactor,ctx->M,ctx->sfactor*ctx->sfactor,ctx->M,0.0,ctx->M,B);
211:   return(0);
212: }

214: /* - - - H1 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

218: PetscErrorCode MatCreateExplicit_Linear_H1A(MPI_Comm comm,PEP_LINEAR *ctx,Mat *A)
219: {

223:   SlepcMatTile(1.0,ctx->K,0.0,ctx->K,ctx->sfactor,ctx->C,1.0,ctx->K,A);
224:   return(0);
225: }

229: PetscErrorCode MatCreateExplicit_Linear_H1B(MPI_Comm comm,PEP_LINEAR *ctx,Mat *B)
230: {

234:   SlepcMatTile(0.0,ctx->K,1.0,ctx->K,-ctx->sfactor*ctx->sfactor,ctx->M,0.0,ctx->K,B);
235:   return(0);
236: }

238: /* - - - H2 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

242: PetscErrorCode MatCreateExplicit_Linear_H2A(MPI_Comm comm,PEP_LINEAR *ctx,Mat *A)
243: {

247:   SlepcMatTile(0.0,ctx->K,-1.0,ctx->K,ctx->sfactor*ctx->sfactor,ctx->M,0.0,ctx->K,A);
248:   return(0);
249: }

253: PetscErrorCode MatCreateExplicit_Linear_H2B(MPI_Comm comm,PEP_LINEAR *ctx,Mat *B)
254: {

258:   SlepcMatTile(ctx->sfactor*ctx->sfactor,ctx->M,ctx->sfactor,ctx->C,0.0,ctx->C,ctx->sfactor*ctx->sfactor,ctx->M,B);
259:   return(0);
260: }