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