PetscErrorCode DSView_NHEP(DS ds,PetscViewer viewer) { PetscErrorCode ierr; PetscFunctionBegin; ierr = DSViewMat(ds,viewer,DS_MAT_A);CHKERRQ(ierr); if (ds->state>DS_STATE_INTERMEDIATE) { ierr = DSViewMat(ds,viewer,DS_MAT_Q);CHKERRQ(ierr); } if (ds->mat[DS_MAT_X]) { ierr = DSViewMat(ds,viewer,DS_MAT_X);CHKERRQ(ierr); } if (ds->mat[DS_MAT_Y]) { ierr = DSViewMat(ds,viewer,DS_MAT_Y);CHKERRQ(ierr); } PetscFunctionReturn(0); }
int main(int argc,char **argv) { PetscErrorCode ierr; DS ds; PetscScalar *A; PetscInt i,j,n=10,ld; PetscViewer viewer; PetscBool verbose; SlepcInitialize(&argc,&argv,(char*)0,help); ierr = PetscOptionsGetInt(NULL,"-n",&n,NULL);CHKERRQ(ierr); ierr = PetscPrintf(PETSC_COMM_WORLD,"Compute symmetric matrix exponential - dimension %D.\n",n);CHKERRQ(ierr); ierr = PetscOptionsHasName(NULL,"-verbose",&verbose);CHKERRQ(ierr); /* Create DS object */ ierr = DSCreate(PETSC_COMM_WORLD,&ds);CHKERRQ(ierr); ierr = DSSetType(ds,DSHEP);CHKERRQ(ierr); ierr = DSSetFromOptions(ds);CHKERRQ(ierr); ld = n+2; /* test leading dimension larger than n */ ierr = DSAllocate(ds,ld);CHKERRQ(ierr); ierr = DSSetDimensions(ds,n,0,0,0);CHKERRQ(ierr); /* Set up viewer */ ierr = PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);CHKERRQ(ierr); ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);CHKERRQ(ierr); ierr = DSView(ds,viewer);CHKERRQ(ierr); ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); if (verbose) { ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);CHKERRQ(ierr); } /* Fill with a symmetric Toeplitz matrix */ ierr = DSGetArray(ds,DS_MAT_A,&A);CHKERRQ(ierr); for (i=0;i<n;i++) A[i+i*ld]=2.0; for (j=1;j<3;j++) { for (i=0;i<n-j;i++) { A[i+(i+j)*ld]=1.0; A[(i+j)+i*ld]=1.0; } } ierr = DSRestoreArray(ds,DS_MAT_A,&A);CHKERRQ(ierr); ierr = DSSetState(ds,DS_STATE_RAW);CHKERRQ(ierr); if (verbose) { ierr = PetscPrintf(PETSC_COMM_WORLD,"Matrix A - - - - - - - -\n");CHKERRQ(ierr); ierr = DSView(ds,viewer);CHKERRQ(ierr); } /* Compute matrix exponential */ ierr = DSComputeFunction(ds,SLEPC_FUNCTION_EXP);CHKERRQ(ierr); if (verbose) { ierr = PetscPrintf(PETSC_COMM_WORLD,"Computed f(A) - - - - - - -\n");CHKERRQ(ierr); ierr = DSViewMat(ds,viewer,DS_MAT_F);CHKERRQ(ierr); } ierr = DSDestroy(&ds);CHKERRQ(ierr); ierr = SlepcFinalize(); return 0; }