1c28052ffSStefano Zampini static char help[] = "Poisson problem with finite elements.\n\
2c28052ffSStefano Zampini We solve the Poisson problem using a parallel unstructured mesh to discretize it.\n\
3c28052ffSStefano Zampini This example is a simplified version of ex12.c that only solves the linear problem.\n\
4c28052ffSStefano Zampini It uses discretized auxiliary fields for coefficient and right-hand side, \n\
5c28052ffSStefano Zampini supports multilevel solvers and non-conforming mesh refinement.\n\n\n";
6c28052ffSStefano Zampini
7c28052ffSStefano Zampini /*
8c28052ffSStefano Zampini Here we describe the PETSc approach to solve nonlinear problems arising from Finite Element (FE) discretizations.
9c28052ffSStefano Zampini
10c28052ffSStefano Zampini The general model requires to solve the residual equations
11c28052ffSStefano Zampini
12c28052ffSStefano Zampini residual(u) := \int_\Omega \phi \cdot f_0(u, \nabla u, a, \nabla a) + \nabla \phi : f_1(u, \nabla u, a, \nabla a) = 0
13c28052ffSStefano Zampini
14c28052ffSStefano Zampini where \phi is a test function, u is the sought FE solution, and a generically describes auxiliary data (for example material properties).
15c28052ffSStefano Zampini
16c28052ffSStefano Zampini The functions f_0 (scalar) and f_1 (vector) describe the problem, while : is the usual contraction operator for tensors, i.e. A : B = \sum_{ij} A_{ij} B_{ij}.
17c28052ffSStefano Zampini
18c28052ffSStefano Zampini The discrete residual is (with abuse of notation)
19c28052ffSStefano Zampini
20c28052ffSStefano Zampini F(u) := \sum_e E_e^T [ B^T_e W_{0,e} f_0(u_q, \nabla u_q, a_q, \nabla a_q) + D_e W_{1,e} f_1(u_q, \nabla u_q, a_q, \nabla a_q) ]
21c28052ffSStefano Zampini
22c28052ffSStefano Zampini where E are element restriction matrices (can support non-conforming meshes), W are quadrature weights, B (resp. D) are basis function (resp. derivative of basis function) matrices, and u_q,a_q are vectors sampled at quadrature points.
23c28052ffSStefano Zampini
24c28052ffSStefano Zampini Having the residual in the above form, it is straightforward to derive its Jacobian (again with abuse of notation)
25c28052ffSStefano Zampini
26c28052ffSStefano Zampini J(u) := \sum_e E_e^T [B^T_e D^T_e] W J_e [ B_e^T, D_e^T ]^T E_e,
27c28052ffSStefano Zampini
28c28052ffSStefano Zampini where J_e is the 2x2 matrix
29c28052ffSStefano Zampini
30c28052ffSStefano Zampini | \partial_u f_0, \partial_{\grad u} f_0 |
31c28052ffSStefano Zampini | \partial_u f_1, \partial_{\grad u} f_1 |
32c28052ffSStefano Zampini
33c28052ffSStefano Zampini Here we use a single-field approach, but the extension to the multi-field case is straightforward.
34c28052ffSStefano Zampini
35c28052ffSStefano Zampini To keep the presentation simple, here we are interested in solving the Poisson problem in divergence form
36c28052ffSStefano Zampini
37c28052ffSStefano Zampini - \nabla \cdot (K * \nabla u) = g
38c28052ffSStefano Zampini
39c28052ffSStefano Zampini with either u = 0 or K * \partial_n u = 0 on \partial \Omega.
40c28052ffSStefano Zampini The above problem possesses the weak form
41c28052ffSStefano Zampini
42c28052ffSStefano Zampini \int_\Omega \nabla \phi K \nabla u + \phi g = 0,
43c28052ffSStefano Zampini
44c28052ffSStefano Zampini thus we have:
45c28052ffSStefano Zampini
46c28052ffSStefano Zampini f_0 = g, f_1 = K * \nabla u, and the only non-zero term of the Jacobian is J_{11} = K
47c28052ffSStefano Zampini
48c28052ffSStefano Zampini See https://petsc.org/release/manual/fe the and the paper "Achieving High Performance with Unified Residual Evaluation" (available at https://arxiv.org/abs/1309.1204) for additional information.
49c28052ffSStefano Zampini */
50c28052ffSStefano Zampini
51c28052ffSStefano Zampini /* Include the necessary definitions */
52c28052ffSStefano Zampini #include <petscdmplex.h>
53c28052ffSStefano Zampini #include <petscsnes.h>
54c28052ffSStefano Zampini #include <petscds.h>
55c28052ffSStefano Zampini
56dd8e379bSPierre Jolivet /* The f_0 function: we read the right-hand side from the first field of the auxiliary data */
f_0(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar f0[])57c28052ffSStefano Zampini static void f_0(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
58c28052ffSStefano Zampini {
59c28052ffSStefano Zampini const PetscScalar g = a[0];
60c28052ffSStefano Zampini
61c28052ffSStefano Zampini f0[0] = g;
62c28052ffSStefano Zampini }
63c28052ffSStefano Zampini
64c28052ffSStefano Zampini /* The f_1 function: we read the conductivity tensor from the second field of the auxiliary data */
f_1(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar f1[])65c28052ffSStefano Zampini static void f_1(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
66c28052ffSStefano Zampini {
67c28052ffSStefano Zampini const PetscScalar *K = &a[1];
68c28052ffSStefano Zampini
69c28052ffSStefano Zampini for (PetscInt d1 = 0; d1 < dim; ++d1) {
70c28052ffSStefano Zampini PetscScalar v = 0;
71c28052ffSStefano Zampini for (PetscInt d2 = 0; d2 < dim; ++d2) v += K[d1 * dim + d2] * u_x[d2];
72c28052ffSStefano Zampini f1[d1] = v;
73c28052ffSStefano Zampini }
74c28052ffSStefano Zampini }
75c28052ffSStefano Zampini
76c28052ffSStefano Zampini /* The only non-zero term for the Jacobian is J_11 */
J_11(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,PetscReal u_tShift,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar J11[])77c28052ffSStefano Zampini static void J_11(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar J11[])
78c28052ffSStefano Zampini {
79c28052ffSStefano Zampini const PetscScalar *K = &a[1];
80c28052ffSStefano Zampini
81c28052ffSStefano Zampini for (PetscInt d1 = 0; d1 < dim; ++d1) {
82c28052ffSStefano Zampini for (PetscInt d2 = 0; d2 < dim; ++d2) J11[d1 * dim + d2] = K[d1 * dim + d2];
83c28052ffSStefano Zampini }
84c28052ffSStefano Zampini }
85c28052ffSStefano Zampini
86c28052ffSStefano Zampini /* The boundary conditions: Dirichlet (essential) or Neumann (natural) */
87c28052ffSStefano Zampini typedef enum {
88c28052ffSStefano Zampini BC_DIRICHLET,
89c28052ffSStefano Zampini BC_NEUMANN,
90c28052ffSStefano Zampini } bcType;
91c28052ffSStefano Zampini
92c28052ffSStefano Zampini static const char *const bcTypes[] = {"DIRICHLET", "NEUMANN", "bcType", "BC_", 0};
93c28052ffSStefano Zampini
94c28052ffSStefano Zampini /* The forcing term: constant or analytical */
95c28052ffSStefano Zampini typedef enum {
96c28052ffSStefano Zampini RHS_CONSTANT,
97c28052ffSStefano Zampini RHS_ANALYTICAL,
98c28052ffSStefano Zampini } rhsType;
99c28052ffSStefano Zampini
100c28052ffSStefano Zampini static const char *const rhsTypes[] = {"CONSTANT", "ANALYTICAL", "rhsType", "RHS_", 0};
101c28052ffSStefano Zampini
102c28052ffSStefano Zampini /* the constant case */
rhs_constant(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * g,PetscCtx ctx)103*2a8381b2SBarry Smith static PetscErrorCode rhs_constant(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *g, PetscCtx ctx)
104c28052ffSStefano Zampini {
105c28052ffSStefano Zampini *g = 1.0;
106c28052ffSStefano Zampini return PETSC_SUCCESS;
107c28052ffSStefano Zampini }
108c28052ffSStefano Zampini
109c28052ffSStefano Zampini /* the analytical case */
rhs_analytical(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * g,PetscCtx ctx)110*2a8381b2SBarry Smith static PetscErrorCode rhs_analytical(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *g, PetscCtx ctx)
111c28052ffSStefano Zampini {
112c28052ffSStefano Zampini PetscScalar v = 1.0;
113c28052ffSStefano Zampini for (PetscInt d = 0; d < dim; d++) v *= PetscSinReal(2.0 * PETSC_PI * x[d]);
114c28052ffSStefano Zampini *g = v;
115c28052ffSStefano Zampini return PETSC_SUCCESS;
116c28052ffSStefano Zampini }
117c28052ffSStefano Zampini
118c28052ffSStefano Zampini /* For the Neumann BC case we need a functional to be integrated: average -> \int_\Omega u dx */
average(PetscInt dim,PetscInt Nf,PetscInt NfAux,const PetscInt uOff[],const PetscInt uOff_x[],const PetscScalar u[],const PetscScalar u_t[],const PetscScalar u_x[],const PetscInt aOff[],const PetscInt aOff_x[],const PetscScalar a[],const PetscScalar a_t[],const PetscScalar a_x[],PetscReal t,const PetscReal x[],PetscInt numConstants,const PetscScalar constants[],PetscScalar obj[])119c28052ffSStefano Zampini static void average(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar obj[])
120c28052ffSStefano Zampini {
121c28052ffSStefano Zampini obj[0] = u[0];
122c28052ffSStefano Zampini }
123c28052ffSStefano Zampini
124c28052ffSStefano Zampini /* The conductivity coefficient term: constant, checkerboard or analytical */
125c28052ffSStefano Zampini typedef enum {
126c28052ffSStefano Zampini COEFF_CONSTANT,
127c28052ffSStefano Zampini COEFF_CHECKERBOARD,
128c28052ffSStefano Zampini COEFF_ANALYTICAL,
129c28052ffSStefano Zampini } coeffType;
130c28052ffSStefano Zampini
131c28052ffSStefano Zampini static const char *const coeffTypes[] = {"CONSTANT", "CHECKERBOARD", "ANALYTICAL", "coeffType", "COEFF_", 0};
132c28052ffSStefano Zampini
133c28052ffSStefano Zampini /* the constant coefficient case */
coefficient_constant(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * K,PetscCtx ctx)134*2a8381b2SBarry Smith static PetscErrorCode coefficient_constant(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *K, PetscCtx ctx)
135c28052ffSStefano Zampini {
136c28052ffSStefano Zampini for (PetscInt d = 0; d < dim; d++) K[d * dim + d] = 1.0;
137c28052ffSStefano Zampini return PETSC_SUCCESS;
138c28052ffSStefano Zampini }
139c28052ffSStefano Zampini
140c28052ffSStefano Zampini /* the checkerboard coefficient case: 10^2 in odd ranks, 10^-2 in even ranks */
coefficient_checkerboard(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * K,PetscCtx ctx)141*2a8381b2SBarry Smith static PetscErrorCode coefficient_checkerboard(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *K, PetscCtx ctx)
142c28052ffSStefano Zampini {
143c28052ffSStefano Zampini PetscScalar exponent = PetscGlobalRank % 2 ? 2.0 : -2.0;
144c28052ffSStefano Zampini for (PetscInt d = 0; d < dim; d++) K[d * dim + d] = PetscPowScalar(10.0, exponent);
145c28052ffSStefano Zampini return PETSC_SUCCESS;
146c28052ffSStefano Zampini }
147c28052ffSStefano Zampini
148c28052ffSStefano Zampini /* the analytical case (channels in diagonal with 4 order of magnitude in jumps) */
coefficient_analytical(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * K,PetscCtx ctx)149*2a8381b2SBarry Smith static PetscErrorCode coefficient_analytical(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *K, PetscCtx ctx)
150c28052ffSStefano Zampini {
151c28052ffSStefano Zampini for (PetscInt d = 0; d < dim; d++) {
152c28052ffSStefano Zampini const PetscReal v = PetscPowReal(10.0, 4.0 * PetscSinReal(4.0 * PETSC_PI * x[d]) * PetscCosReal(4.0 * PETSC_PI * x[d]));
153c28052ffSStefano Zampini K[d * dim + d] = v;
154c28052ffSStefano Zampini }
155c28052ffSStefano Zampini return PETSC_SUCCESS;
156c28052ffSStefano Zampini }
157c28052ffSStefano Zampini
158c28052ffSStefano Zampini /* The application context that defines our problem */
159c28052ffSStefano Zampini typedef struct {
160c28052ffSStefano Zampini bcType bc; /* type of boundary conditions */
161c28052ffSStefano Zampini rhsType rhs; /* type of right-hand side */
162c28052ffSStefano Zampini coeffType coeff; /* type of conductivity coefficient */
163c28052ffSStefano Zampini PetscInt order; /* the polynomial order for the solution */
164c28052ffSStefano Zampini PetscInt rhsOrder; /* the polynomial order for the right-hand side */
165c28052ffSStefano Zampini PetscInt coeffOrder; /* the polynomial order for the coefficient */
166c28052ffSStefano Zampini PetscBool p4est; /* if we want to use non-conforming AMR */
167c28052ffSStefano Zampini } AppCtx;
168c28052ffSStefano Zampini
169c28052ffSStefano Zampini /* Process command line options */
ProcessOptions(MPI_Comm comm,AppCtx * options)170c28052ffSStefano Zampini static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
171c28052ffSStefano Zampini {
172c28052ffSStefano Zampini PetscFunctionBeginUser;
173c28052ffSStefano Zampini options->bc = BC_DIRICHLET;
174c28052ffSStefano Zampini options->rhs = RHS_CONSTANT;
175c28052ffSStefano Zampini options->coeff = COEFF_CONSTANT;
176c28052ffSStefano Zampini options->order = 1;
177c28052ffSStefano Zampini options->rhsOrder = 1;
178c28052ffSStefano Zampini options->coeffOrder = 1;
179c28052ffSStefano Zampini options->p4est = PETSC_FALSE;
180c28052ffSStefano Zampini
181c28052ffSStefano Zampini PetscOptionsBegin(comm, "", "Poisson problem options", "DMPLEX");
182c28052ffSStefano Zampini PetscCall(PetscOptionsEnum("-bc_type", "Type of boundary condition", __FILE__, bcTypes, (PetscEnum)options->bc, (PetscEnum *)&options->bc, NULL));
183c28052ffSStefano Zampini PetscCall(PetscOptionsEnum("-rhs_type", "Type of forcing term", __FILE__, rhsTypes, (PetscEnum)options->rhs, (PetscEnum *)&options->rhs, NULL));
184c28052ffSStefano Zampini PetscCall(PetscOptionsEnum("-coefficient_type", "Type of conductivity coefficient", __FILE__, coeffTypes, (PetscEnum)options->coeff, (PetscEnum *)&options->coeff, NULL));
185c28052ffSStefano Zampini PetscCall(PetscOptionsInt("-order", "The polynomial order for the approximation of the solution", __FILE__, options->order, &options->order, NULL));
186c28052ffSStefano Zampini PetscCall(PetscOptionsInt("-rhs_order", "The polynomial order for the approximation of the right-hand side", __FILE__, options->rhsOrder, &options->rhsOrder, NULL));
187c28052ffSStefano Zampini PetscCall(PetscOptionsInt("-coefficient_order", "The polynomial order for the approximation of the coefficient", __FILE__, options->coeffOrder, &options->coeffOrder, NULL));
188c28052ffSStefano Zampini PetscCall(PetscOptionsBool("-p4est", "Use p4est to represent the mesh", __FILE__, options->p4est, &options->p4est, NULL));
189c28052ffSStefano Zampini PetscOptionsEnd();
190c28052ffSStefano Zampini
191c28052ffSStefano Zampini /* View processed options */
192c28052ffSStefano Zampini PetscCall(PetscPrintf(comm, "Simulation parameters\n"));
193c28052ffSStefano Zampini PetscCall(PetscPrintf(comm, " polynomial order: %" PetscInt_FMT "\n", options->order));
194c28052ffSStefano Zampini PetscCall(PetscPrintf(comm, " boundary conditions: %s\n", bcTypes[options->bc]));
195c28052ffSStefano Zampini PetscCall(PetscPrintf(comm, " right-hand side: %s, order %" PetscInt_FMT "\n", rhsTypes[options->rhs], options->rhsOrder));
196c28052ffSStefano Zampini PetscCall(PetscPrintf(comm, " coefficient: %s, order %" PetscInt_FMT "\n", coeffTypes[options->coeff], options->coeffOrder));
197c28052ffSStefano Zampini PetscCall(PetscPrintf(comm, " non-conforming AMR: %d\n", options->p4est));
198c28052ffSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS);
199c28052ffSStefano Zampini }
200c28052ffSStefano Zampini
201c28052ffSStefano Zampini /* Create mesh from command line options */
CreateMesh(MPI_Comm comm,AppCtx * user,DM * dm)202c28052ffSStefano Zampini static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
203c28052ffSStefano Zampini {
204c28052ffSStefano Zampini PetscFunctionBeginUser;
205c28052ffSStefano Zampini /* Create various types of meshes only with command line options */
206c28052ffSStefano Zampini PetscCall(DMCreate(comm, dm));
207c28052ffSStefano Zampini PetscCall(DMSetType(*dm, DMPLEX));
208c28052ffSStefano Zampini PetscCall(DMSetOptionsPrefix(*dm, "initial_"));
209c28052ffSStefano Zampini PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
210c28052ffSStefano Zampini PetscCall(DMSetFromOptions(*dm));
211c28052ffSStefano Zampini PetscCall(DMLocalizeCoordinates(*dm));
212c28052ffSStefano Zampini PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
213c28052ffSStefano Zampini PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_TRUE));
214c28052ffSStefano Zampini PetscCall(DMSetOptionsPrefix(*dm, "mesh_"));
215c28052ffSStefano Zampini PetscCall(DMSetFromOptions(*dm));
216c28052ffSStefano Zampini
217c28052ffSStefano Zampini /* If requested convert to a format suitable for non-conforming adaptive mesh refinement */
218c28052ffSStefano Zampini if (user->p4est) {
219c28052ffSStefano Zampini PetscInt dim;
220c28052ffSStefano Zampini DM dmConv;
221c28052ffSStefano Zampini
222c28052ffSStefano Zampini PetscCall(DMGetDimension(*dm, &dim));
223c28052ffSStefano Zampini PetscCheck(dim == 2 || dim == 3, comm, PETSC_ERR_SUP, "p4est support not for dimension %" PetscInt_FMT, dim);
224c28052ffSStefano Zampini PetscCall(DMConvert(*dm, dim == 3 ? DMP8EST : DMP4EST, &dmConv));
225c28052ffSStefano Zampini if (dmConv) {
226c28052ffSStefano Zampini PetscCall(DMDestroy(dm));
227c28052ffSStefano Zampini PetscCall(DMSetOptionsPrefix(dmConv, "mesh_"));
228c28052ffSStefano Zampini PetscCall(DMSetFromOptions(dmConv));
229c28052ffSStefano Zampini *dm = dmConv;
230c28052ffSStefano Zampini }
231c28052ffSStefano Zampini }
232c28052ffSStefano Zampini PetscCall(DMSetUp(*dm));
233c28052ffSStefano Zampini
234c28052ffSStefano Zampini /* View the mesh.
235c28052ffSStefano Zampini With a single call we can dump ASCII information, VTK data for visualization, store the mesh in HDF5 format, etc. */
236c28052ffSStefano Zampini PetscCall(PetscObjectSetName((PetscObject)*dm, "Mesh"));
237c28052ffSStefano Zampini PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
238c28052ffSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS);
239c28052ffSStefano Zampini }
240c28052ffSStefano Zampini
241c28052ffSStefano Zampini /* Setup the discrete problem */
SetupProblem(DM dm,DM fdm,AppCtx * user)242c28052ffSStefano Zampini static PetscErrorCode SetupProblem(DM dm, DM fdm, AppCtx *user)
243c28052ffSStefano Zampini {
244c28052ffSStefano Zampini DM plex, dmAux, cdm = NULL, coordDM;
245c28052ffSStefano Zampini Vec auxData, auxDataGlobal;
246c28052ffSStefano Zampini PetscDS ds;
247c28052ffSStefano Zampini DMPolytopeType ct;
248c28052ffSStefano Zampini PetscInt dim, cdim, cStart;
249c28052ffSStefano Zampini PetscFE fe, fe_rhs, fe_K;
250*2a8381b2SBarry Smith PetscErrorCode (*auxFuncs[])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar u[], PetscCtx ctx) = {NULL, NULL};
251c28052ffSStefano Zampini void *auxCtxs[] = {NULL, NULL};
252c28052ffSStefano Zampini
253c28052ffSStefano Zampini PetscFunctionBeginUser;
254c28052ffSStefano Zampini /* Create the Finite Element for the solution and pass it to the problem DM */
255c28052ffSStefano Zampini PetscCall(DMGetDimension(dm, &dim));
256c28052ffSStefano Zampini PetscCall(DMConvert(dm, DMPLEX, &plex));
257c28052ffSStefano Zampini PetscCall(DMPlexGetHeightStratum(plex, 0, &cStart, NULL));
258c28052ffSStefano Zampini PetscCall(DMPlexGetCellType(plex, cStart, &ct));
259c28052ffSStefano Zampini PetscCall(DMDestroy(&plex));
260c28052ffSStefano Zampini PetscCall(PetscFECreateLagrangeByCell(PETSC_COMM_SELF, dim, 1, ct, user->order, PETSC_DETERMINE, &fe));
261c28052ffSStefano Zampini PetscCall(PetscObjectSetName((PetscObject)fe, "potential"));
262c28052ffSStefano Zampini PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
263c28052ffSStefano Zampini PetscCall(DMCreateDS(dm));
264c28052ffSStefano Zampini
265c28052ffSStefano Zampini /* Set residual and Jacobian callbacks */
266c28052ffSStefano Zampini PetscCall(DMGetDS(dm, &ds));
267c28052ffSStefano Zampini PetscCall(PetscDSSetResidual(ds, 0, f_0, f_1));
268c28052ffSStefano Zampini PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, J_11));
269c28052ffSStefano Zampini /* Tell DMPLEX we are going to use FEM callbacks */
270c28052ffSStefano Zampini PetscCall(DMPlexSetSNESLocalFEM(dm, PETSC_FALSE, &user));
271c28052ffSStefano Zampini
272c28052ffSStefano Zampini /* Create the Finite Element for the auxiliary data */
273c28052ffSStefano Zampini PetscCall(PetscFECreateLagrangeByCell(PETSC_COMM_SELF, dim, 1, ct, user->rhsOrder, PETSC_DETERMINE, &fe_rhs));
274c28052ffSStefano Zampini PetscCall(PetscObjectSetName((PetscObject)fe_rhs, "g"));
275c28052ffSStefano Zampini PetscCall(PetscFECopyQuadrature(fe, fe_rhs));
276c28052ffSStefano Zampini PetscCall(DMGetCoordinateDM(dm, &coordDM));
277c28052ffSStefano Zampini PetscCall(DMGetDimension(coordDM, &cdim));
278c28052ffSStefano Zampini PetscCall(PetscFECreateLagrangeByCell(PETSC_COMM_SELF, dim, cdim * cdim, ct, user->coeffOrder, PETSC_DETERMINE, &fe_K));
279c28052ffSStefano Zampini PetscCall(PetscObjectSetName((PetscObject)fe_K, "K"));
280c28052ffSStefano Zampini PetscCall(PetscFECopyQuadrature(fe, fe_K));
281c28052ffSStefano Zampini PetscCall(DMClone(dm, &dmAux));
282c28052ffSStefano Zampini PetscCall(DMSetField(dmAux, 0, NULL, (PetscObject)fe_rhs));
283c28052ffSStefano Zampini PetscCall(DMSetField(dmAux, 1, NULL, (PetscObject)fe_K));
284c28052ffSStefano Zampini PetscCall(DMCreateDS(dmAux));
285c28052ffSStefano Zampini
286c28052ffSStefano Zampini /* Project the requested rhs and K to the auxiliary DM and pass it to the problem DM */
287c28052ffSStefano Zampini PetscCall(DMCreateLocalVector(dmAux, &auxData));
288c28052ffSStefano Zampini PetscCall(DMCreateGlobalVector(dmAux, &auxDataGlobal));
289c28052ffSStefano Zampini PetscCall(PetscObjectSetName((PetscObject)auxData, ""));
290c28052ffSStefano Zampini if (!fdm) {
291c28052ffSStefano Zampini switch (user->rhs) {
292c28052ffSStefano Zampini case RHS_CONSTANT:
293c28052ffSStefano Zampini auxFuncs[0] = rhs_constant;
294c28052ffSStefano Zampini auxCtxs[0] = NULL;
295c28052ffSStefano Zampini break;
296c28052ffSStefano Zampini case RHS_ANALYTICAL:
297c28052ffSStefano Zampini auxFuncs[0] = rhs_analytical;
298c28052ffSStefano Zampini auxCtxs[0] = NULL;
299c28052ffSStefano Zampini break;
300c28052ffSStefano Zampini default:
301c28052ffSStefano Zampini SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported rhs type %s", rhsTypes[user->rhs]);
302c28052ffSStefano Zampini }
303c28052ffSStefano Zampini switch (user->coeff) {
304c28052ffSStefano Zampini case COEFF_CONSTANT:
305c28052ffSStefano Zampini auxFuncs[1] = coefficient_constant;
306c28052ffSStefano Zampini auxCtxs[1] = NULL;
307c28052ffSStefano Zampini break;
308c28052ffSStefano Zampini case COEFF_CHECKERBOARD:
309c28052ffSStefano Zampini auxFuncs[1] = coefficient_checkerboard;
310c28052ffSStefano Zampini auxCtxs[1] = NULL;
311c28052ffSStefano Zampini break;
312c28052ffSStefano Zampini case COEFF_ANALYTICAL:
313c28052ffSStefano Zampini auxFuncs[1] = coefficient_analytical;
314c28052ffSStefano Zampini auxCtxs[1] = NULL;
315c28052ffSStefano Zampini break;
316c28052ffSStefano Zampini default:
317c28052ffSStefano Zampini SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Unsupported coefficient type %s", coeffTypes[user->coeff]);
318c28052ffSStefano Zampini }
319c28052ffSStefano Zampini PetscCall(DMGetDS(dmAux, &ds));
320c28052ffSStefano Zampini PetscCall(DMProjectFunction(dmAux, 0.0, auxFuncs, auxCtxs, INSERT_ALL_VALUES, auxDataGlobal));
321c28052ffSStefano Zampini if (user->bc == BC_NEUMANN) {
322c28052ffSStefano Zampini PetscScalar vals[2];
323c28052ffSStefano Zampini PetscInt rhs_id = 0;
324c28052ffSStefano Zampini IS is;
325c28052ffSStefano Zampini
326c28052ffSStefano Zampini PetscCall(PetscDSSetObjective(ds, rhs_id, average));
327c28052ffSStefano Zampini PetscCall(DMPlexComputeIntegralFEM(dmAux, auxDataGlobal, vals, NULL));
328c28052ffSStefano Zampini PetscCall(DMCreateSubDM(dmAux, 1, &rhs_id, &is, NULL));
329c28052ffSStefano Zampini PetscCall(VecISShift(auxDataGlobal, is, -vals[rhs_id]));
330c28052ffSStefano Zampini PetscCall(DMPlexComputeIntegralFEM(dmAux, auxDataGlobal, vals, NULL));
331c28052ffSStefano Zampini PetscCall(ISDestroy(&is));
332c28052ffSStefano Zampini }
333c28052ffSStefano Zampini } else {
334c28052ffSStefano Zampini Mat J;
335c28052ffSStefano Zampini Vec auxDataGlobalf, auxDataf, Jscale;
336c28052ffSStefano Zampini DM dmAuxf;
337c28052ffSStefano Zampini
338c28052ffSStefano Zampini PetscCall(DMGetAuxiliaryVec(fdm, NULL, 0, 0, &auxDataf));
339c28052ffSStefano Zampini PetscCall(VecGetDM(auxDataf, &dmAuxf));
340c28052ffSStefano Zampini PetscCall(DMSetCoarseDM(dmAuxf, dmAux));
341c28052ffSStefano Zampini PetscCall(DMCreateGlobalVector(dmAuxf, &auxDataGlobalf));
342c28052ffSStefano Zampini PetscCall(DMLocalToGlobal(dmAuxf, auxDataf, INSERT_VALUES, auxDataGlobalf));
343c28052ffSStefano Zampini PetscCall(DMCreateInterpolation(dmAux, dmAuxf, &J, &Jscale));
344c28052ffSStefano Zampini PetscCall(MatInterpolate(J, auxDataGlobalf, auxDataGlobal));
345c28052ffSStefano Zampini PetscCall(VecPointwiseMult(auxDataGlobal, auxDataGlobal, Jscale));
346c28052ffSStefano Zampini PetscCall(VecDestroy(&Jscale));
347c28052ffSStefano Zampini PetscCall(VecDestroy(&auxDataGlobalf));
348c28052ffSStefano Zampini PetscCall(MatDestroy(&J));
349c28052ffSStefano Zampini PetscCall(DMSetCoarseDM(dmAuxf, NULL));
350c28052ffSStefano Zampini }
351c28052ffSStefano Zampini /* auxiliary data must be a local vector */
352c28052ffSStefano Zampini PetscCall(DMGlobalToLocal(dmAux, auxDataGlobal, INSERT_VALUES, auxData));
353c28052ffSStefano Zampini { /* view auxiliary data */
354c28052ffSStefano Zampini PetscInt level;
355c28052ffSStefano Zampini char optionName[PETSC_MAX_PATH_LEN];
356c28052ffSStefano Zampini
357c28052ffSStefano Zampini PetscCall(DMGetRefineLevel(dm, &level));
358c28052ffSStefano Zampini PetscCall(PetscSNPrintf(optionName, sizeof(optionName), "-aux_%" PetscInt_FMT "_vec_view", level));
359c28052ffSStefano Zampini PetscCall(VecViewFromOptions(auxData, NULL, optionName));
360c28052ffSStefano Zampini }
361c28052ffSStefano Zampini PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, auxData));
362c28052ffSStefano Zampini PetscCall(VecDestroy(&auxData));
363c28052ffSStefano Zampini PetscCall(VecDestroy(&auxDataGlobal));
364c28052ffSStefano Zampini PetscCall(DMDestroy(&dmAux));
365c28052ffSStefano Zampini
366c28052ffSStefano Zampini /* Setup boundary conditions
367c28052ffSStefano Zampini Since we use homogeneous natural boundary conditions for the Neumann problem we
368c28052ffSStefano Zampini only handle the Dirichlet case here */
369c28052ffSStefano Zampini if (user->bc == BC_DIRICHLET) {
370c28052ffSStefano Zampini DMLabel label;
371c28052ffSStefano Zampini PetscInt id = 1;
372c28052ffSStefano Zampini
373c28052ffSStefano Zampini /* Label faces on the mesh boundary */
374c28052ffSStefano Zampini PetscCall(DMCreateLabel(dm, "boundary"));
375c28052ffSStefano Zampini PetscCall(DMGetLabel(dm, "boundary", &label));
376c28052ffSStefano Zampini PetscCall(DMPlexMarkBoundaryFaces(dm, 1, label));
377c28052ffSStefano Zampini PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "dirichlet", label, 1, &id, 0, 0, NULL, NULL, NULL, NULL, NULL));
378c28052ffSStefano Zampini }
379c28052ffSStefano Zampini
380c28052ffSStefano Zampini /* Iterate on coarser mesh if present */
3819d854079SStefano Zampini PetscCall(DMGetCoarseDM(dm, &cdm));
382c28052ffSStefano Zampini if (cdm) PetscCall(SetupProblem(cdm, dm, user));
383c28052ffSStefano Zampini
384c28052ffSStefano Zampini PetscCall(PetscFEDestroy(&fe));
385c28052ffSStefano Zampini PetscCall(PetscFEDestroy(&fe_rhs));
386c28052ffSStefano Zampini PetscCall(PetscFEDestroy(&fe_K));
387c28052ffSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS);
388c28052ffSStefano Zampini }
389c28052ffSStefano Zampini
390c28052ffSStefano Zampini /* We are now ready to run the simulation */
main(int argc,char ** argv)391c28052ffSStefano Zampini int main(int argc, char **argv)
392c28052ffSStefano Zampini {
393c28052ffSStefano Zampini DM dm; /* problem specification */
394c28052ffSStefano Zampini SNES snes; /* nonlinear solver */
395c28052ffSStefano Zampini KSP ksp; /* linear solver */
396c28052ffSStefano Zampini Vec u; /* solution vector */
397c28052ffSStefano Zampini AppCtx user; /* user-defined work context */
398c28052ffSStefano Zampini
399c28052ffSStefano Zampini PetscFunctionBeginUser;
400c28052ffSStefano Zampini PetscCall(PetscInitialize(&argc, &argv, NULL, help));
401c28052ffSStefano Zampini PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
402c28052ffSStefano Zampini PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
403c28052ffSStefano Zampini PetscCall(SetupProblem(dm, NULL, &user));
404c28052ffSStefano Zampini
405c28052ffSStefano Zampini PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
406c28052ffSStefano Zampini PetscCall(SNESSetDM(snes, dm));
407c28052ffSStefano Zampini PetscCall(SNESSetType(snes, SNESKSPONLY));
408c28052ffSStefano Zampini PetscCall(SNESGetKSP(snes, &ksp));
409c28052ffSStefano Zampini PetscCall(KSPSetType(ksp, KSPCG));
410c28052ffSStefano Zampini PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL));
411c28052ffSStefano Zampini PetscCall(SNESSetFromOptions(snes));
412c28052ffSStefano Zampini
413c28052ffSStefano Zampini PetscCall(DMCreateGlobalVector(dm, &u));
414c28052ffSStefano Zampini PetscCall(PetscObjectSetName((PetscObject)u, ""));
415c28052ffSStefano Zampini PetscCall(VecSetRandom(u, NULL));
416c28052ffSStefano Zampini PetscCall(SNESSolve(snes, NULL, u));
417c28052ffSStefano Zampini PetscCall(VecDestroy(&u));
418c28052ffSStefano Zampini PetscCall(SNESDestroy(&snes));
419c28052ffSStefano Zampini PetscCall(DMDestroy(&dm));
420c28052ffSStefano Zampini PetscCall(PetscFinalize());
421c28052ffSStefano Zampini return 0;
422c28052ffSStefano Zampini }
423c28052ffSStefano Zampini
424c28052ffSStefano Zampini /*TEST
425c28052ffSStefano Zampini
426c28052ffSStefano Zampini test:
427c28052ffSStefano Zampini suffix: 0
428c28052ffSStefano Zampini requires: triangle !single
429c28052ffSStefano Zampini args: -bc_type {{dirichlet neumann}separate output} -rhs_type {{constant analytical}separate output} -coefficient_type {{constant checkerboard analytical}separate output} -initial_dm_plex_simplex 1 -pc_type svd -snes_type newtonls -snes_error_if_not_converged
430c28052ffSStefano Zampini
431c28052ffSStefano Zampini test:
432c28052ffSStefano Zampini requires: !single
433c28052ffSStefano Zampini suffix: 0_quad
434c28052ffSStefano Zampini args: -bc_type {{dirichlet neumann}separate output} -rhs_type {{constant analytical}separate output} -coefficient_type {{constant checkerboard analytical}separate output} -initial_dm_plex_simplex 0 -pc_type svd -snes_type newtonls -snes_error_if_not_converged
435c28052ffSStefano Zampini
436c28052ffSStefano Zampini test:
437c28052ffSStefano Zampini suffix: 0_p4est
438c28052ffSStefano Zampini requires: p4est !single
439c28052ffSStefano Zampini args: -bc_type {{dirichlet neumann}separate output} -rhs_type {{constant analytical}separate output} -coefficient_type {{constant checkerboard analytical}separate output} -initial_dm_plex_simplex 0 -pc_type svd -snes_type newtonls -snes_error_if_not_converged
440c28052ffSStefano Zampini
4418fad524dSPierre Jolivet testset:
442c28052ffSStefano Zampini nsize: 2
443c28052ffSStefano Zampini requires: hpddm slepc !single defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES)
44443eb4678SPierre Jolivet args: -pc_type hpddm -pc_hpddm_coarse_correction balanced -pc_hpddm_coarse_mat_type aij -pc_hpddm_levels_eps_nev 1 -pc_hpddm_levels_sub_pc_type lu -ksp_monitor -initial_dm_plex_simplex 0 -petscpartitioner_type simple
4458fad524dSPierre Jolivet output_file: output/ex11_hpddm.out
446c28052ffSStefano Zampini test:
4478fad524dSPierre Jolivet suffix: hpddm
4488fad524dSPierre Jolivet args:
4498fad524dSPierre Jolivet test:
4508fad524dSPierre Jolivet suffix: hpddm_harmonic_overlap
4518fad524dSPierre Jolivet args: -pc_hpddm_harmonic_overlap 1 -pc_hpddm_has_neumann false -pc_hpddm_levels_1_pc_asm_type basic
4528fad524dSPierre Jolivet test:
453c28052ffSStefano Zampini suffix: hpddm_p4est
4548fad524dSPierre Jolivet requires: p4est
4558fad524dSPierre Jolivet args: -p4est
4568fad524dSPierre Jolivet filter: sed -e "s/non-conforming AMR: 1/non-conforming AMR: 0/g"
457c28052ffSStefano Zampini
458c28052ffSStefano Zampini test:
459c28052ffSStefano Zampini nsize: 4
460c28052ffSStefano Zampini suffix: gdsw_corner
461c28052ffSStefano Zampini requires: triangle
462c28052ffSStefano Zampini args: -pc_type mg -pc_mg_galerkin -pc_mg_adapt_interp_coarse_space gdsw -pc_mg_levels 2 -mg_levels_pc_type asm -initial_dm_plex_shape ball -initial_dm_refine 2 -mesh_dm_mat_type {{aij is}} -mg_levels_sub_pc_type lu -mg_levels_pc_asm_type basic -petscpartitioner_type simple
463c28052ffSStefano Zampini
464a0d6ca37SStefano Zampini test:
465a0d6ca37SStefano Zampini suffix: asm_seq
466a0d6ca37SStefano Zampini args: -pc_type asm -pc_asm_dm_subdomains -sub_pc_type cholesky -snes_type newtonls -snes_error_if_not_converged -initial_dm_plex_simplex 0
467a0d6ca37SStefano Zampini
468c28052ffSStefano Zampini TEST*/
469