Actual source code: test18.c
slepc-3.17.0 2022-03-31
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: static char help[] = "Test DSSynchronize() on a NHEP.\n\n";
13: #include <slepcds.h>
15: PetscErrorCode CheckArray(PetscScalar *A,const char *label,PetscInt k)
16: {
17: PetscInt j;
18: PetscMPIInt p,size,rank;
19: PetscScalar dif,*buf;
20: PetscReal error;
23: MPI_Comm_size(PETSC_COMM_WORLD,&size);
24: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
25: if (rank) MPI_Send(A,k,MPIU_SCALAR,0,111,PETSC_COMM_WORLD);
26: else {
27: PetscMalloc1(k,&buf);
28: for (p=1;p<size;p++) {
29: MPI_Recv(buf,k,MPIU_SCALAR,p,111,PETSC_COMM_WORLD,MPI_STATUS_IGNORE);
30: dif = 0.0;
31: for (j=0;j<k;j++) dif += A[j]-buf[j];
32: error = PetscAbsScalar(dif);
33: if (error>10*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD," Array %s differs in proc %d: %g\n",label,(int)p,(double)error);
34: }
35: PetscFree(buf);
36: }
37: PetscFunctionReturn(0);
38: }
40: int main(int argc,char **argv)
41: {
42: DS ds;
43: SlepcSC sc;
44: PetscScalar *A,*Q,*wr,*wi;
45: PetscReal re,im;
46: PetscInt i,j,n=10;
47: PetscMPIInt size;
49: SlepcInitialize(&argc,&argv,(char*)0,help);
50: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
51: PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type NHEP - dimension %" PetscInt_FMT ".\n",n);
53: /* Create DS object */
54: DSCreate(PETSC_COMM_WORLD,&ds);
55: DSSetType(ds,DSNHEP);
56: DSSetFromOptions(ds);
57: DSAllocate(ds,n);
58: DSSetDimensions(ds,n,0,0);
60: /* Fill with Grcar matrix */
61: DSGetArray(ds,DS_MAT_A,&A);
62: for (i=1;i<n;i++) A[i+(i-1)*n]=-1.0;
63: for (j=0;j<4;j++) {
64: for (i=0;i<n-j;i++) A[i+(i+j)*n]=1.0;
65: }
66: DSRestoreArray(ds,DS_MAT_A,&A);
67: DSSetState(ds,DS_STATE_INTERMEDIATE);
69: /* Solve */
70: PetscMalloc2(n,&wr,n,&wi);
71: DSGetSlepcSC(ds,&sc);
72: sc->comparison = SlepcCompareLargestMagnitude;
73: sc->comparisonctx = NULL;
74: sc->map = NULL;
75: sc->mapobj = NULL;
76: DSSolve(ds,wr,wi);
77: DSSort(ds,wr,wi,NULL,NULL,NULL);
79: /* Print eigenvalues */
80: PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n");
81: for (i=0;i<n;i++) {
82: #if defined(PETSC_USE_COMPLEX)
83: re = PetscRealPart(wr[i]);
84: im = PetscImaginaryPart(wr[i]);
85: #else
86: re = wr[i];
87: im = wi[i];
88: #endif
89: if (PetscAbs(im)<1e-10) PetscPrintf(PETSC_COMM_WORLD," %.5f\n",(double)re);
90: else PetscPrintf(PETSC_COMM_WORLD," %.5f%+.5fi\n",(double)re,(double)im);
91: }
93: /* Synchronize data and check */
94: DSSynchronize(ds,wr,wi);
95: MPI_Comm_size(PETSC_COMM_WORLD,&size);
96: if (size>1) {
97: CheckArray(wr,"wr",n);
98: #if !defined(PETSC_USE_COMPLEX)
99: CheckArray(wi,"wi",n);
100: #endif
101: DSGetArray(ds,DS_MAT_A,&A);
102: CheckArray(A,"A",n*n);
103: DSRestoreArray(ds,DS_MAT_A,&A);
104: DSGetArray(ds,DS_MAT_Q,&Q);
105: CheckArray(Q,"Q",n*n);
106: DSRestoreArray(ds,DS_MAT_Q,&Q);
107: }
109: PetscFree2(wr,wi);
110: DSDestroy(&ds);
111: SlepcFinalize();
112: return 0;
113: }
115: /*TEST
117: testset:
118: nsize: {{1 2 3}}
119: filter: sed -e "s/[+-]\([0-9]\.[0-9]*i\)/+-\\1/" | sed -e "s/1.58254/1.58255/" | sed -e "s/1.75989/1.75988/"
120: output_file: output/test18_1.out
121: test:
122: suffix: 1
123: args: -ds_parallel redundant
124: test:
125: suffix: 2
126: args: -ds_parallel synchronized
128: TEST*/