1*b22c5e46SPierre Jolivet static char help[] = "Test MatConvert() with a MATNEST with scaled and shifted MATTRANSPOSEVIRTUAL blocks.\n\n"; 2*b22c5e46SPierre Jolivet 3*b22c5e46SPierre Jolivet #include <petscmat.h> 4*b22c5e46SPierre Jolivet 5*b22c5e46SPierre Jolivet /* 6*b22c5e46SPierre Jolivet This example builds the matrix 7*b22c5e46SPierre Jolivet 8*b22c5e46SPierre Jolivet H = [ R C 9*b22c5e46SPierre Jolivet alpha C^H + beta I gamma R^T + delta I ], 10*b22c5e46SPierre Jolivet 11*b22c5e46SPierre Jolivet where R is Hermitian and C is complex symmetric. In particular, R and C have the 12*b22c5e46SPierre Jolivet following Toeplitz structure: 13*b22c5e46SPierre Jolivet 14*b22c5e46SPierre Jolivet R = pentadiag{a,b,c,conj(b),conj(a)} 15*b22c5e46SPierre Jolivet C = tridiag{b,d,b} 16*b22c5e46SPierre Jolivet 17*b22c5e46SPierre Jolivet where a,b,d are complex scalars, and c is real. 18*b22c5e46SPierre Jolivet */ 19*b22c5e46SPierre Jolivet 20*b22c5e46SPierre Jolivet int main(int argc, char **argv) 21*b22c5e46SPierre Jolivet { 22*b22c5e46SPierre Jolivet Mat block[4], H, R, C, M; 23*b22c5e46SPierre Jolivet PetscScalar a, b, c, d; 24*b22c5e46SPierre Jolivet PetscInt n = 13, Istart, Iend, i; 25*b22c5e46SPierre Jolivet PetscBool flg; 26*b22c5e46SPierre Jolivet 27*b22c5e46SPierre Jolivet PetscFunctionBeginUser; 28*b22c5e46SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 29*b22c5e46SPierre Jolivet 30*b22c5e46SPierre Jolivet PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL)); 31*b22c5e46SPierre Jolivet 32*b22c5e46SPierre Jolivet a = PetscCMPLX(-0.1, 0.2); 33*b22c5e46SPierre Jolivet b = PetscCMPLX(1.0, 0.5); 34*b22c5e46SPierre Jolivet c = 4.5; 35*b22c5e46SPierre Jolivet d = PetscCMPLX(2.0, 0.2); 36*b22c5e46SPierre Jolivet 37*b22c5e46SPierre Jolivet PetscCall(MatCreate(PETSC_COMM_WORLD, &R)); 38*b22c5e46SPierre Jolivet PetscCall(MatSetSizes(R, PETSC_DECIDE, PETSC_DECIDE, n, n)); 39*b22c5e46SPierre Jolivet PetscCall(MatSetFromOptions(R)); 40*b22c5e46SPierre Jolivet 41*b22c5e46SPierre Jolivet PetscCall(MatCreate(PETSC_COMM_WORLD, &C)); 42*b22c5e46SPierre Jolivet PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, n, n)); 43*b22c5e46SPierre Jolivet PetscCall(MatSetFromOptions(C)); 44*b22c5e46SPierre Jolivet 45*b22c5e46SPierre Jolivet PetscCall(MatGetOwnershipRange(R, &Istart, &Iend)); 46*b22c5e46SPierre Jolivet for (i = Istart; i < Iend; i++) { 47*b22c5e46SPierre Jolivet if (i > 1) PetscCall(MatSetValue(R, i, i - 2, a, INSERT_VALUES)); 48*b22c5e46SPierre Jolivet if (i > 0) PetscCall(MatSetValue(R, i, i - 1, b, INSERT_VALUES)); 49*b22c5e46SPierre Jolivet PetscCall(MatSetValue(R, i, i, c, INSERT_VALUES)); 50*b22c5e46SPierre Jolivet if (i < n - 1) PetscCall(MatSetValue(R, i, i + 1, PetscConj(b), INSERT_VALUES)); 51*b22c5e46SPierre Jolivet if (i < n - 2) PetscCall(MatSetValue(R, i, i + 2, PetscConj(a), INSERT_VALUES)); 52*b22c5e46SPierre Jolivet } 53*b22c5e46SPierre Jolivet 54*b22c5e46SPierre Jolivet PetscCall(MatGetOwnershipRange(C, &Istart, &Iend)); 55*b22c5e46SPierre Jolivet for (i = Istart; i < Iend; i++) { 56*b22c5e46SPierre Jolivet if (i > 0) PetscCall(MatSetValue(C, i, i - 1, b, INSERT_VALUES)); 57*b22c5e46SPierre Jolivet PetscCall(MatSetValue(C, i, i, d, INSERT_VALUES)); 58*b22c5e46SPierre Jolivet if (i < n - 1) PetscCall(MatSetValue(C, i, i + 1, b, INSERT_VALUES)); 59*b22c5e46SPierre Jolivet } 60*b22c5e46SPierre Jolivet 61*b22c5e46SPierre Jolivet PetscCall(MatAssemblyBegin(R, MAT_FINAL_ASSEMBLY)); 62*b22c5e46SPierre Jolivet PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 63*b22c5e46SPierre Jolivet PetscCall(MatAssemblyEnd(R, MAT_FINAL_ASSEMBLY)); 64*b22c5e46SPierre Jolivet PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 65*b22c5e46SPierre Jolivet 66*b22c5e46SPierre Jolivet block[0] = R; 67*b22c5e46SPierre Jolivet block[1] = C; 68*b22c5e46SPierre Jolivet 69*b22c5e46SPierre Jolivet PetscCall(MatCreateHermitianTranspose(C, &block[2])); 70*b22c5e46SPierre Jolivet PetscCall(MatScale(block[2], PetscConj(b))); 71*b22c5e46SPierre Jolivet PetscCall(MatShift(block[2], d)); 72*b22c5e46SPierre Jolivet PetscCall(MatCreateTranspose(R, &block[3])); 73*b22c5e46SPierre Jolivet PetscCall(MatScale(block[3], PetscConj(d))); 74*b22c5e46SPierre Jolivet PetscCall(MatShift(block[3], b)); 75*b22c5e46SPierre Jolivet PetscCall(MatCreateNest(PetscObjectComm((PetscObject)R), 2, NULL, 2, NULL, block, &H)); 76*b22c5e46SPierre Jolivet PetscCall(MatDestroy(&block[2])); 77*b22c5e46SPierre Jolivet PetscCall(MatDestroy(&block[3])); 78*b22c5e46SPierre Jolivet 79*b22c5e46SPierre Jolivet PetscCall(MatConvert(H, MATAIJ, MAT_INITIAL_MATRIX, &M)); 80*b22c5e46SPierre Jolivet PetscCall(MatMultEqual(H, M, 20, &flg)); 81*b22c5e46SPierre Jolivet PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "MatNest != MatAIJ"); 82*b22c5e46SPierre Jolivet 83*b22c5e46SPierre Jolivet PetscCall(MatDestroy(&R)); 84*b22c5e46SPierre Jolivet PetscCall(MatDestroy(&C)); 85*b22c5e46SPierre Jolivet PetscCall(MatDestroy(&H)); 86*b22c5e46SPierre Jolivet PetscCall(MatDestroy(&M)); 87*b22c5e46SPierre Jolivet PetscCall(PetscFinalize()); 88*b22c5e46SPierre Jolivet return 0; 89*b22c5e46SPierre Jolivet } 90*b22c5e46SPierre Jolivet 91*b22c5e46SPierre Jolivet /*TEST 92*b22c5e46SPierre Jolivet 93*b22c5e46SPierre Jolivet build: 94*b22c5e46SPierre Jolivet requires: complex 95*b22c5e46SPierre Jolivet 96*b22c5e46SPierre Jolivet test: 97*b22c5e46SPierre Jolivet output_file: output/ex109.out 98*b22c5e46SPierre Jolivet nsize: {{1 4}} 99*b22c5e46SPierre Jolivet 100*b22c5e46SPierre Jolivet TEST*/ 101