xref: /petsc/src/snes/tutorials/ex31.c (revision fcf85c8c01162926ce932c17797a2648ce6ee273)
1*fcf85c8cSAdelekeBankole 
2*fcf85c8cSAdelekeBankole static char help[] = "A Chebyshev spectral method for the compressible Blasius boundary layer equations.\n\n";
3*fcf85c8cSAdelekeBankole 
4*fcf85c8cSAdelekeBankole /*
5*fcf85c8cSAdelekeBankole    Include "petscsnes.h" so that we can use SNES solvers.  Note that this
6*fcf85c8cSAdelekeBankole    file automatically includes:
7*fcf85c8cSAdelekeBankole      petscsys.h       - base PETSc routines   petscvec.h - vectors
8*fcf85c8cSAdelekeBankole      petscmat.h - matrices
9*fcf85c8cSAdelekeBankole      petscis.h     - index sets            petscksp.h - Krylov subspace methods
10*fcf85c8cSAdelekeBankole      petscviewer.h - viewers               petscpc.h  - preconditioners
11*fcf85c8cSAdelekeBankole      petscksp.h   - linear solvers
12*fcf85c8cSAdelekeBankole    Include "petscdt.h" so that we can have support for use of Quadrature formulas
13*fcf85c8cSAdelekeBankole */
14*fcf85c8cSAdelekeBankole /*F
15*fcf85c8cSAdelekeBankole This examples solves the compressible Blasius boundary layer equations
16*fcf85c8cSAdelekeBankole 2(\rho\muf'')' + ff'' = 0
17*fcf85c8cSAdelekeBankole (\rho\muh')' + Prfh' + Pr(\gamma-1)Ma^{2}\rho\muf''^{2} = 0
18*fcf85c8cSAdelekeBankole following Howarth-Dorodnitsyn transformation with boundary conditions
19*fcf85c8cSAdelekeBankole f(0) = f'(0) = 0, f'(\infty) = 1, h(\infty) = 1, h = \theta(0). Where \theta = T/T_{\infty}
20*fcf85c8cSAdelekeBankole Note: density (\rho) and viscosity (\mu) are treated as constants in this example
21*fcf85c8cSAdelekeBankole F*/
22*fcf85c8cSAdelekeBankole #include <petscsnes.h>
23*fcf85c8cSAdelekeBankole #include <petscdt.h>
24*fcf85c8cSAdelekeBankole 
25*fcf85c8cSAdelekeBankole /*
26*fcf85c8cSAdelekeBankole    User-defined routines
27*fcf85c8cSAdelekeBankole */
28*fcf85c8cSAdelekeBankole 
29*fcf85c8cSAdelekeBankole extern PetscErrorCode FormFunction(SNES,Vec,Vec,void*);
30*fcf85c8cSAdelekeBankole 
31*fcf85c8cSAdelekeBankole typedef struct {
32*fcf85c8cSAdelekeBankole    PetscReal Ma, Pr, h_0;
33*fcf85c8cSAdelekeBankole    PetscInt N;
34*fcf85c8cSAdelekeBankole    PetscReal dx_deta;
35*fcf85c8cSAdelekeBankole    PetscReal *x;
36*fcf85c8cSAdelekeBankole    PetscReal gamma;
37*fcf85c8cSAdelekeBankole } Blasius;
38*fcf85c8cSAdelekeBankole 
39*fcf85c8cSAdelekeBankole int main(int argc,char **argv)
40*fcf85c8cSAdelekeBankole {
41*fcf85c8cSAdelekeBankole   SNES           snes;         /* nonlinear solver context */
42*fcf85c8cSAdelekeBankole   Vec            x,r;          /* solution, residual vectors */
43*fcf85c8cSAdelekeBankole   PetscMPIInt    size;
44*fcf85c8cSAdelekeBankole   Blasius        *blasius;
45*fcf85c8cSAdelekeBankole   PetscReal      L, *weight;   /* L is size of the domain */
46*fcf85c8cSAdelekeBankole 
47*fcf85c8cSAdelekeBankole   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
48*fcf85c8cSAdelekeBankole   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
49*fcf85c8cSAdelekeBankole   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Example is only for sequential runs");
50*fcf85c8cSAdelekeBankole 
51*fcf85c8cSAdelekeBankole   // Read command-line arguments
52*fcf85c8cSAdelekeBankole   PetscCall(PetscCalloc1(1, &blasius));
53*fcf85c8cSAdelekeBankole   blasius->Ma = 2;            /* Mach number */
54*fcf85c8cSAdelekeBankole   blasius->Pr = 0.7;          /* Prandtl number */
55*fcf85c8cSAdelekeBankole   blasius->h_0 = 2.;          /* relative temperature at the wall */
56*fcf85c8cSAdelekeBankole   blasius->N = 10;            /* Number of Chebyshev terms */
57*fcf85c8cSAdelekeBankole   blasius->gamma = 1.4;       /* specific heat ratio */
58*fcf85c8cSAdelekeBankole   L = 5;
59*fcf85c8cSAdelekeBankole   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Compressible Blasius boundary layer equations", "");
60*fcf85c8cSAdelekeBankole   PetscCall(PetscOptionsReal("-mach", "Mach number at freestream", "", blasius->Ma, &blasius->Ma, NULL));
61*fcf85c8cSAdelekeBankole   PetscCall(PetscOptionsReal("-prandtl", "Prandtl number", "", blasius->Pr, &blasius->Pr, NULL));
62*fcf85c8cSAdelekeBankole   PetscCall(PetscOptionsReal("-h_0", "Relative enthalpy at wall", "", blasius->h_0, &blasius->h_0, NULL));
63*fcf85c8cSAdelekeBankole   PetscCall(PetscOptionsReal("-gamma", "Ratio of specific heats", "", blasius->gamma, &blasius->gamma, NULL));
64*fcf85c8cSAdelekeBankole   PetscCall(PetscOptionsInt("-N", "Number of Chebyshev terms for f", "", blasius->N, &blasius->N, NULL));
65*fcf85c8cSAdelekeBankole   PetscCall(PetscOptionsReal("-L", "Extent of the domain", "", L, &L, NULL));
66*fcf85c8cSAdelekeBankole   PetscOptionsEnd();
67*fcf85c8cSAdelekeBankole   blasius->dx_deta = 2 / L;   /* this helps to map [-1,1] to [0,L] */
68*fcf85c8cSAdelekeBankole   PetscCall(PetscMalloc2(blasius->N-3, &blasius->x, blasius->N-3, &weight));
69*fcf85c8cSAdelekeBankole   PetscCall(PetscDTGaussQuadrature(blasius->N-3, -1., 1., blasius->x, weight));
70*fcf85c8cSAdelekeBankole 
71*fcf85c8cSAdelekeBankole   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
72*fcf85c8cSAdelekeBankole      Create nonlinear solver context
73*fcf85c8cSAdelekeBankole      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
74*fcf85c8cSAdelekeBankole   PetscCall(SNESCreate(PETSC_COMM_WORLD,&snes));
75*fcf85c8cSAdelekeBankole 
76*fcf85c8cSAdelekeBankole   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
77*fcf85c8cSAdelekeBankole      Create matrix and vector data structures; set corresponding routines
78*fcf85c8cSAdelekeBankole      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
79*fcf85c8cSAdelekeBankole   /*
80*fcf85c8cSAdelekeBankole      Create vectors for solution and nonlinear function
81*fcf85c8cSAdelekeBankole   */
82*fcf85c8cSAdelekeBankole   PetscCall(VecCreate(PETSC_COMM_WORLD,&x));
83*fcf85c8cSAdelekeBankole   PetscCall(VecSetSizes(x,PETSC_DECIDE,2*blasius->N-1));
84*fcf85c8cSAdelekeBankole   PetscCall(VecSetFromOptions(x));
85*fcf85c8cSAdelekeBankole   PetscCall(VecDuplicate(x,&r));
86*fcf85c8cSAdelekeBankole 
87*fcf85c8cSAdelekeBankole    /*
88*fcf85c8cSAdelekeBankole       Set function evaluation routine and vector.
89*fcf85c8cSAdelekeBankole    */
90*fcf85c8cSAdelekeBankole   PetscCall(SNESSetFunction(snes,r,FormFunction,blasius));
91*fcf85c8cSAdelekeBankole   {
92*fcf85c8cSAdelekeBankole      KSP ksp;
93*fcf85c8cSAdelekeBankole      PC pc;
94*fcf85c8cSAdelekeBankole      SNESGetKSP(snes, &ksp);
95*fcf85c8cSAdelekeBankole      KSPSetType(ksp, KSPPREONLY);
96*fcf85c8cSAdelekeBankole      KSPGetPC(ksp, &pc);
97*fcf85c8cSAdelekeBankole      PCSetType(pc, PCLU);
98*fcf85c8cSAdelekeBankole   }
99*fcf85c8cSAdelekeBankole   /*
100*fcf85c8cSAdelekeBankole      Set SNES/KSP/KSP/PC runtime options, e.g.,
101*fcf85c8cSAdelekeBankole          -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
102*fcf85c8cSAdelekeBankole      These options will override those specified above as long as
103*fcf85c8cSAdelekeBankole      SNESSetFromOptions() is called _after_ any other customization
104*fcf85c8cSAdelekeBankole      routines.
105*fcf85c8cSAdelekeBankole   */
106*fcf85c8cSAdelekeBankole   PetscCall(SNESSetFromOptions(snes));
107*fcf85c8cSAdelekeBankole 
108*fcf85c8cSAdelekeBankole   PetscCall(SNESSolve(snes,NULL,x));
109*fcf85c8cSAdelekeBankole   //PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD));
110*fcf85c8cSAdelekeBankole 
111*fcf85c8cSAdelekeBankole   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
112*fcf85c8cSAdelekeBankole      Free work space.  All PETSc objects should be destroyed when they
113*fcf85c8cSAdelekeBankole      are no longer needed.
114*fcf85c8cSAdelekeBankole    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
115*fcf85c8cSAdelekeBankole 
116*fcf85c8cSAdelekeBankole   PetscCall(PetscFree2(blasius->x, weight));
117*fcf85c8cSAdelekeBankole   PetscCall(PetscFree(blasius));
118*fcf85c8cSAdelekeBankole   PetscCall(VecDestroy(&x));
119*fcf85c8cSAdelekeBankole   PetscCall(VecDestroy(&r));
120*fcf85c8cSAdelekeBankole   PetscCall(SNESDestroy(&snes));
121*fcf85c8cSAdelekeBankole   PetscCall(PetscFinalize());
122*fcf85c8cSAdelekeBankole   return 0;
123*fcf85c8cSAdelekeBankole }
124*fcf85c8cSAdelekeBankole 
125*fcf85c8cSAdelekeBankole /*-------------------------------------------------------------------------------
126*fcf85c8cSAdelekeBankole    Helper function to evaluate Chebyshev polynomials with a set of coefficients
127*fcf85c8cSAdelekeBankole    with all their derivatives represented as a recurrence table
128*fcf85c8cSAdelekeBankole -------------------------------------------------------------------------------*/
129*fcf85c8cSAdelekeBankole static void ChebyshevEval(PetscInt N, const PetscScalar *Tf, PetscReal x, PetscReal dx_deta, PetscScalar *f){
130*fcf85c8cSAdelekeBankole    PetscScalar table[4][3] = {
131*fcf85c8cSAdelekeBankole       {1, x, 2*x*x - 1}, {0, 1, 4*x}, {0, 0, 4}, {0, 0, 0} /* Chebyshev polynomials T_0, T_1, T_2 of the first kind in (-1,1)  */
132*fcf85c8cSAdelekeBankole    };
133*fcf85c8cSAdelekeBankole    for (int i=0; i<4; i++) {
134*fcf85c8cSAdelekeBankole       f[i] = table[i][0] * Tf[0] + table[i][1] * Tf[1] + table[i][2] * Tf[2];  /* i-th derivative of f */
135*fcf85c8cSAdelekeBankole    }
136*fcf85c8cSAdelekeBankole    for (int i=3; i<N; i++) {
137*fcf85c8cSAdelekeBankole       table[0][i%3] = 2 * x * table[0][(i-1) % 3] - table[0][(i-2)%3];  /* T_n(x) = 2xT_{n-1}(x) - T_{n-2}(x) */
138*fcf85c8cSAdelekeBankole       /* Differentiate Chebyshev polynomials with the recurrence relation */
139*fcf85c8cSAdelekeBankole       for (int j=1; j<4; j++) {
140*fcf85c8cSAdelekeBankole          table[j][i%3] = i * (2 * table[j-1][(i-1) % 3] + table[j][(i-2)%3] / (i-2)); /* T'_{n}(x)/n = 2T_{n-1}(x) + T'_{n-2}(x)/n-2 */
141*fcf85c8cSAdelekeBankole       }
142*fcf85c8cSAdelekeBankole       for (int j=0; j<4; j++) {
143*fcf85c8cSAdelekeBankole          f[j] += table[j][i%3] * Tf[i];
144*fcf85c8cSAdelekeBankole       }
145*fcf85c8cSAdelekeBankole    }
146*fcf85c8cSAdelekeBankole    for (int i=1; i<4; i++) {
147*fcf85c8cSAdelekeBankole       for (int j=0; j<i; j++) f[i] *= dx_deta; /* Here happens the physics of the problem */
148*fcf85c8cSAdelekeBankole    }
149*fcf85c8cSAdelekeBankole }
150*fcf85c8cSAdelekeBankole 
151*fcf85c8cSAdelekeBankole /* ------------------------------------------------------------------- */
152*fcf85c8cSAdelekeBankole /*
153*fcf85c8cSAdelekeBankole    FormFunction - Evaluates nonlinear function, F(x).
154*fcf85c8cSAdelekeBankole 
155*fcf85c8cSAdelekeBankole    Input Parameters:
156*fcf85c8cSAdelekeBankole .  snes - the SNES context
157*fcf85c8cSAdelekeBankole .  X    - input vector
158*fcf85c8cSAdelekeBankole .  ctx  - optional user-defined context
159*fcf85c8cSAdelekeBankole 
160*fcf85c8cSAdelekeBankole    Output Parameter:
161*fcf85c8cSAdelekeBankole .  R - function vector
162*fcf85c8cSAdelekeBankole  */
163*fcf85c8cSAdelekeBankole PetscErrorCode FormFunction(SNES snes,Vec X,Vec R,void *ctx)
164*fcf85c8cSAdelekeBankole {
165*fcf85c8cSAdelekeBankole   Blasius *blasius = (Blasius *)ctx;
166*fcf85c8cSAdelekeBankole   const PetscScalar *Tf, *Th;  /* Tf and Th are Chebyshev coefficients */
167*fcf85c8cSAdelekeBankole   PetscScalar       *r, f[4], h[4];
168*fcf85c8cSAdelekeBankole   PetscInt N = blasius->N;
169*fcf85c8cSAdelekeBankole   PetscReal Ma = blasius->Ma, Pr = blasius->Pr;
170*fcf85c8cSAdelekeBankole 
171*fcf85c8cSAdelekeBankole   /*
172*fcf85c8cSAdelekeBankole    Get pointers to vector data.
173*fcf85c8cSAdelekeBankole       - For default PETSc vectors, VecGetArray() returns a pointer to
174*fcf85c8cSAdelekeBankole         the data array.  Otherwise, the routine is implementation dependent.
175*fcf85c8cSAdelekeBankole       - You MUST call VecRestoreArray() when you no longer need access to
176*fcf85c8cSAdelekeBankole         the array.
177*fcf85c8cSAdelekeBankole    */
178*fcf85c8cSAdelekeBankole   PetscCall(VecGetArrayRead(X,&Tf));
179*fcf85c8cSAdelekeBankole   Th = Tf + N;
180*fcf85c8cSAdelekeBankole   PetscCall(VecGetArray(R,&r));
181*fcf85c8cSAdelekeBankole 
182*fcf85c8cSAdelekeBankole   /* Compute function */
183*fcf85c8cSAdelekeBankole   ChebyshevEval(N, Tf, -1., blasius->dx_deta, f);
184*fcf85c8cSAdelekeBankole   r[0] = f[0];
185*fcf85c8cSAdelekeBankole   r[1] = f[1];
186*fcf85c8cSAdelekeBankole   ChebyshevEval(N, Tf, 1., blasius->dx_deta, f);
187*fcf85c8cSAdelekeBankole   r[2] = f[1]  - 1; /* Right end boundary condition */
188*fcf85c8cSAdelekeBankole   for (int i=0; i<N - 3; i++) {
189*fcf85c8cSAdelekeBankole      ChebyshevEval(N, Tf, blasius->x[i], blasius->dx_deta, f);
190*fcf85c8cSAdelekeBankole      r[3+i] = 2*f[3] + f[2] * f[0];
191*fcf85c8cSAdelekeBankole      ChebyshevEval(N-1, Th, blasius->x[i], blasius->dx_deta, h);
192*fcf85c8cSAdelekeBankole      r[N+2+i] = h[2] + Pr * f[0] * h[1] + Pr * (blasius->gamma - 1) * PetscSqr(Ma * f[2]);
193*fcf85c8cSAdelekeBankole   }
194*fcf85c8cSAdelekeBankole   ChebyshevEval(N-1, Th, -1., blasius->dx_deta, h);
195*fcf85c8cSAdelekeBankole   r[N] = h[0] - blasius->h_0; /* Left end boundary condition */
196*fcf85c8cSAdelekeBankole   ChebyshevEval(N-1, Th, 1., blasius->dx_deta, h);
197*fcf85c8cSAdelekeBankole   r[N+1] = h[0] - 1; /* Left end boundary condition */
198*fcf85c8cSAdelekeBankole 
199*fcf85c8cSAdelekeBankole   /* Restore vectors */
200*fcf85c8cSAdelekeBankole   PetscCall(VecRestoreArrayRead(X,&Tf));
201*fcf85c8cSAdelekeBankole   PetscCall(VecRestoreArray(R,&r));
202*fcf85c8cSAdelekeBankole   return 0;
203*fcf85c8cSAdelekeBankole }
204*fcf85c8cSAdelekeBankole 
205*fcf85c8cSAdelekeBankole /*TEST
206*fcf85c8cSAdelekeBankole 
207*fcf85c8cSAdelekeBankole    test:
208*fcf85c8cSAdelekeBankole       args: -snes_monitor -pc_type svd
209*fcf85c8cSAdelekeBankole       requires: !single
210*fcf85c8cSAdelekeBankole 
211*fcf85c8cSAdelekeBankole TEST*/
212