Actual source code: ex37.c
slepc-3.11.2 2019-07-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, 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[] = "Computes exp(t*A)*v for an advection diffusion operator with Peclet number.\n\n"
12: "The command line options are:\n"
13: " -n <idim>, where <idim> = number of subdivisions of the mesh in each spatial direction.\n"
14: " -t <sval>, where <sval> = scalar value that multiplies the argument.\n"
15: " -peclet <sval>, where <sval> = Peclet value.\n"
16: " -steps <ival>, where <ival> = number of time steps.\n\n";
18: #include <slepcmfn.h>
20: int main(int argc,char **argv)
21: {
22: Mat A; /* problem matrix */
23: MFN mfn;
24: FN f;
25: PetscInt i,j,Istart,Iend,II,m,n=10,N,steps=5,its,totits=0,ncv,maxit;
26: PetscReal tol,norm,h,h2,peclet=0.5,epsilon=1.0,c,i1h,j1h;
27: PetscScalar t=1e-4,sone=1.0,value,upper,diag,lower;
28: Vec v;
29: PetscErrorCode ierr;
30: MFNConvergedReason reason;
32: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
34: PetscOptionsGetScalar(NULL,NULL,"-t",&t,NULL);
35: PetscOptionsGetReal(NULL,NULL,"-peclet",&peclet,NULL);
36: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
37: PetscOptionsGetInt(NULL,NULL,"-steps",&steps,NULL);
38: m = n;
39: N = m*n;
40: /* interval [0,1], homogeneous Dirichlet boundary conditions */
41: h = 1.0/(n+1.0);
42: h2 = h*h;
43: c = 2.0*epsilon*peclet/h;
44: upper = epsilon/h2+c/(2.0*h);
45: diag = 2.0*(-2.0*epsilon/h2);
46: lower = epsilon/h2-c/(2.0*h);
48: PetscPrintf(PETSC_COMM_WORLD,"\nAdvection diffusion via y=exp(%g*A), n=%D, steps=%D, Peclet=%g\n\n",(double)PetscRealPart(t),n,steps,(double)peclet);
50: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
51: Generate matrix A
52: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
53: MatCreate(PETSC_COMM_WORLD,&A);
54: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
55: MatSetFromOptions(A);
56: MatSetUp(A);
57: MatGetOwnershipRange(A,&Istart,&Iend);
58: for (II=Istart;II<Iend;II++) {
59: i = II/n; j = II-i*n;
60: if (i>0) { MatSetValue(A,II,II-n,lower,INSERT_VALUES); }
61: if (i<m-1) { MatSetValue(A,II,II+n,upper,INSERT_VALUES); }
62: if (j>0) { MatSetValue(A,II,II-1,lower,INSERT_VALUES); }
63: if (j<n-1) { MatSetValue(A,II,II+1,upper,INSERT_VALUES); }
64: MatSetValue(A,II,II,diag,INSERT_VALUES);
65: }
66: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
67: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
68: MatCreateVecs(A,NULL,&v);
70: /*
71: Set initial condition v = 256*i^2*(1-i)^2*j^2*(1-j)^2
72: */
73: for (II=Istart;II<Iend;II++) {
74: i = II/n; j = II-i*n;
75: i1h = (i+1)*h; j1h = (j+1)*h;
76: value = 256.0*i1h*i1h*(1.0-i1h)*(1.0-i1h)*(j1h*j1h)*(1.0-j1h)*(1.0-j1h);
77: VecSetValue(v,i+j*n,value,INSERT_VALUES);
78: }
79: VecAssemblyBegin(v);
80: VecAssemblyEnd(v);
82: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
83: Create the solver and set various options
84: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
85: MFNCreate(PETSC_COMM_WORLD,&mfn);
86: MFNSetOperator(mfn,A);
87: MFNGetFN(mfn,&f);
88: FNSetType(f,FNEXP);
89: FNSetScale(f,t,sone);
90: MFNSetFromOptions(mfn);
92: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
93: Solve the problem, y=exp(t*A)*v
94: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
95: for (i=0;i<steps;i++) {
96: MFNSolve(mfn,v,v);
97: MFNGetConvergedReason(mfn,&reason);
98: if (reason<0) SETERRQ(PETSC_COMM_WORLD,1,"Solver did not converge");
99: MFNGetIterationNumber(mfn,&its);
100: totits += its;
101: }
103: /*
104: Optional: Get some information from the solver and display it
105: */
106: PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %D\n",totits);
107: MFNGetDimensions(mfn,&ncv);
108: PetscPrintf(PETSC_COMM_WORLD," Subspace dimension: %D\n",ncv);
109: MFNGetTolerances(mfn,&tol,&maxit);
110: PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%D\n",(double)tol,maxit);
111: VecNorm(v,NORM_2,&norm);
112: PetscPrintf(PETSC_COMM_WORLD," Computed vector at time t=%.4g has norm %g\n",(double)PetscRealPart(t)*steps,(double)norm);
114: /*
115: Free work space
116: */
117: MFNDestroy(&mfn);
118: MatDestroy(&A);
119: VecDestroy(&v);
120: SlepcFinalize();
121: return ierr;
122: }
124: /*TEST
126: test:
127: suffix: 1
128: args: -mfn_tol 1e-6
130: TEST*/