1c4762a1bSJed Brown static char help[] = "Pure advection with finite elements.\n\
2c4762a1bSJed Brown We solve the hyperbolic problem in a rectangular\n\
3c4762a1bSJed Brown domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\n\n";
4c4762a1bSJed Brown
5c4762a1bSJed Brown /*
6c4762a1bSJed Brown The continuity equation (https://en.wikipedia.org/wiki/Continuity_equation) for advection
7c4762a1bSJed Brown (https://en.wikipedia.org/wiki/Advection) of a conserved scalar quantity phi, with source q,
8c4762a1bSJed Brown
9c4762a1bSJed Brown phi_t + div (phi u) = q
10c4762a1bSJed Brown
11c4762a1bSJed Brown if used with a solenoidal velocity field u (div u = 0) is given by
12c4762a1bSJed Brown
13c4762a1bSJed Brown phi_t + u . grad phi = q
14c4762a1bSJed Brown
15c4762a1bSJed Brown For a vector quantity a, we likewise have
16c4762a1bSJed Brown
17c4762a1bSJed Brown a_t + u . grad a = q
18c4762a1bSJed Brown */
19c4762a1bSJed Brown
20c4762a1bSJed Brown /*
21c4762a1bSJed Brown r1: 8 SOR
22c4762a1bSJed Brown r2: 1128 SOR
23c4762a1bSJed Brown r3: > 10000 SOR
24c4762a1bSJed Brown
25c4762a1bSJed Brown SOR is completely unreliable as a smoother, use Jacobi
26c4762a1bSJed Brown r1: 8 MG
27c4762a1bSJed Brown r2:
28c4762a1bSJed Brown */
29c4762a1bSJed Brown
30c4762a1bSJed Brown #include <petscdmplex.h>
31c4762a1bSJed Brown #include <petscts.h>
32c4762a1bSJed Brown #include <petscds.h>
33c4762a1bSJed Brown
349371c9d4SSatish Balay typedef enum {
359371c9d4SSatish Balay PRIMITIVE,
369371c9d4SSatish Balay INT_BY_PARTS
379371c9d4SSatish Balay } WeakFormType;
38c4762a1bSJed Brown
39c4762a1bSJed Brown typedef struct {
40c4762a1bSJed Brown WeakFormType formType;
41c4762a1bSJed Brown } AppCtx;
42c4762a1bSJed Brown
43c4762a1bSJed Brown /* MMS1:
44c4762a1bSJed Brown
45c4762a1bSJed Brown 2D:
46c4762a1bSJed Brown u = <1, 1>
47c4762a1bSJed Brown phi = x + y - 2t
48c4762a1bSJed Brown
49c4762a1bSJed Brown phi_t + u . grad phi = -2 + <1, 1> . <1, 1> = 0
50c4762a1bSJed Brown
51c4762a1bSJed Brown 3D:
52c4762a1bSJed Brown u = <1, 1, 1>
53c4762a1bSJed Brown phi = x + y + z - 3t
54c4762a1bSJed Brown
55c4762a1bSJed Brown phi_t + u . grad phi = -3 + <1, 1, 1> . <1, 1, 1> = 0
56c4762a1bSJed Brown */
57c4762a1bSJed Brown
analytic_phi(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)58*2a8381b2SBarry Smith static PetscErrorCode analytic_phi(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
59d71ae5a4SJacob Faibussowitsch {
60c4762a1bSJed Brown PetscInt d;
61c4762a1bSJed Brown
62c4762a1bSJed Brown *u = -dim * time;
63c4762a1bSJed Brown for (d = 0; d < dim; ++d) *u += x[d];
643ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
65c4762a1bSJed Brown }
66c4762a1bSJed Brown
velocity(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nf,PetscScalar * u,PetscCtx ctx)67*2a8381b2SBarry Smith static PetscErrorCode velocity(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, PetscCtx ctx)
68d71ae5a4SJacob Faibussowitsch {
69c4762a1bSJed Brown PetscInt d;
70c4762a1bSJed Brown for (d = 0; d < dim; ++d) u[d] = 1.0;
713ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
72c4762a1bSJed Brown }
73c4762a1bSJed Brown
74c4762a1bSJed Brown /* <psi, phi_t> + <psi, u . grad phi> */
f0_prim_phi(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[])75d71ae5a4SJacob Faibussowitsch static void f0_prim_phi(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[])
76d71ae5a4SJacob Faibussowitsch {
77c4762a1bSJed Brown PetscInt d;
78c4762a1bSJed Brown
79c4762a1bSJed Brown f0[0] = u_t[0];
80c4762a1bSJed Brown for (d = 0; d < dim; ++d) f0[0] += a[d] * u_x[d];
81c4762a1bSJed Brown }
82c4762a1bSJed Brown
83c4762a1bSJed Brown /* <psi, phi_t> */
f0_ibp_phi(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[])84d71ae5a4SJacob Faibussowitsch static void f0_ibp_phi(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[])
85d71ae5a4SJacob Faibussowitsch {
86c4762a1bSJed Brown f0[0] = u_t[0];
87c4762a1bSJed Brown }
88c4762a1bSJed Brown
89c4762a1bSJed Brown /* <grad psi, u phi> */
f1_ibp_phi(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[])90d71ae5a4SJacob Faibussowitsch static void f1_ibp_phi(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[])
91d71ae5a4SJacob Faibussowitsch {
92c4762a1bSJed Brown PetscInt d;
93c4762a1bSJed Brown for (d = 0; d < dim; ++d) f1[d] = a[d] * u[0];
94c4762a1bSJed Brown }
95c4762a1bSJed Brown
96c4762a1bSJed Brown /* <psi, phi_t> */
g0_prim_phi(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 g0[])97d71ae5a4SJacob Faibussowitsch static void g0_prim_phi(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 g0[])
98d71ae5a4SJacob Faibussowitsch {
99c4762a1bSJed Brown g0[0] = u_tShift * 1.0;
100c4762a1bSJed Brown }
101c4762a1bSJed Brown
102c4762a1bSJed Brown /* <psi, u . grad phi> */
g1_prim_phi(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 g1[])103d71ae5a4SJacob Faibussowitsch static void g1_prim_phi(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 g1[])
104d71ae5a4SJacob Faibussowitsch {
105c4762a1bSJed Brown PetscInt d;
106c4762a1bSJed Brown for (d = 0; d < dim; ++d) g1[d] = a[d];
107c4762a1bSJed Brown }
108c4762a1bSJed Brown
109c4762a1bSJed Brown /* <grad psi, u phi> */
g2_ibp_phi(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 g2[])110d71ae5a4SJacob Faibussowitsch static void g2_ibp_phi(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 g2[])
111d71ae5a4SJacob Faibussowitsch {
112c4762a1bSJed Brown PetscInt d;
113c4762a1bSJed Brown for (d = 0; d < dim; ++d) g2[d] = a[d];
114c4762a1bSJed Brown }
115c4762a1bSJed Brown
ProcessOptions(MPI_Comm comm,AppCtx * options)116d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
117d71ae5a4SJacob Faibussowitsch {
118c4762a1bSJed Brown const char *formTypes[2] = {"primitive", "int_by_parts"};
119c4762a1bSJed Brown PetscInt form;
120c4762a1bSJed Brown
121c4762a1bSJed Brown PetscFunctionBeginUser;
122c4762a1bSJed Brown options->formType = PRIMITIVE;
123d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Advection Equation Options", "DMPLEX");
124c4762a1bSJed Brown form = options->formType;
1259566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-form_type", "The weak form type", "ex47.c", formTypes, 2, formTypes[options->formType], &form, NULL));
126c4762a1bSJed Brown options->formType = (WeakFormType)form;
127d0609cedSBarry Smith PetscOptionsEnd();
1283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
129c4762a1bSJed Brown }
130c4762a1bSJed Brown
CreateMesh(MPI_Comm comm,DM * dm,AppCtx * ctx)131d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm, AppCtx *ctx)
132d71ae5a4SJacob Faibussowitsch {
133c4762a1bSJed Brown PetscFunctionBeginUser;
1349566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm));
1359566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX));
1369566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm));
1379566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
1383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
139c4762a1bSJed Brown }
140c4762a1bSJed Brown
SetupProblem(DM dm,AppCtx * ctx)141d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupProblem(DM dm, AppCtx *ctx)
142d71ae5a4SJacob Faibussowitsch {
14345480ffeSMatthew G. Knepley PetscDS ds;
14445480ffeSMatthew G. Knepley DMLabel label;
145c4762a1bSJed Brown const PetscInt id = 1;
146c4762a1bSJed Brown
147c4762a1bSJed Brown PetscFunctionBeginUser;
1489566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds));
149c4762a1bSJed Brown switch (ctx->formType) {
150c4762a1bSJed Brown case PRIMITIVE:
1519566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_prim_phi, NULL));
1529566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_prim_phi, g1_prim_phi, NULL, NULL));
153c4762a1bSJed Brown break;
154c4762a1bSJed Brown case INT_BY_PARTS:
1559566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(ds, 0, f0_ibp_phi, f1_ibp_phi));
1569566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_prim_phi, NULL, g2_ibp_phi, NULL));
157c4762a1bSJed Brown break;
158c4762a1bSJed Brown }
1599566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, analytic_phi, ctx));
1609566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label));
16157d50842SBarry Smith PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (PetscVoidFn *)analytic_phi, NULL, ctx, NULL));
1623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
163c4762a1bSJed Brown }
164c4762a1bSJed Brown
SetupVelocity(DM dm,DM dmAux,AppCtx * user)165d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupVelocity(DM dm, DM dmAux, AppCtx *user)
166d71ae5a4SJacob Faibussowitsch {
1678434afd1SBarry Smith PetscSimplePointFn *funcs[1] = {velocity};
168c4762a1bSJed Brown Vec v;
169c4762a1bSJed Brown
170c4762a1bSJed Brown PetscFunctionBeginUser;
1719566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(dmAux, &v));
1729566063dSJacob Faibussowitsch PetscCall(DMProjectFunctionLocal(dmAux, 0.0, funcs, NULL, INSERT_ALL_VALUES, v));
1739566063dSJacob Faibussowitsch PetscCall(DMSetAuxiliaryVec(dm, NULL, 0, 0, v));
1749566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v));
1753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
176c4762a1bSJed Brown }
177c4762a1bSJed Brown
SetupAuxDM(DM dm,PetscFE feAux,AppCtx * user)178d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupAuxDM(DM dm, PetscFE feAux, AppCtx *user)
179d71ae5a4SJacob Faibussowitsch {
180c4762a1bSJed Brown DM dmAux, coordDM;
181c4762a1bSJed Brown
1827510d9b0SBarry Smith PetscFunctionBeginUser;
183c4762a1bSJed Brown /* MUST call DMGetCoordinateDM() in order to get p4est setup if present */
1849566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &coordDM));
1853ba16761SJacob Faibussowitsch if (!feAux) PetscFunctionReturn(PETSC_SUCCESS);
1869566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &dmAux));
1879566063dSJacob Faibussowitsch PetscCall(DMSetCoordinateDM(dmAux, coordDM));
1889566063dSJacob Faibussowitsch PetscCall(DMSetField(dmAux, 0, NULL, (PetscObject)feAux));
1899566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dmAux));
1909566063dSJacob Faibussowitsch PetscCall(SetupVelocity(dm, dmAux, user));
1919566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmAux));
1923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
193c4762a1bSJed Brown }
194c4762a1bSJed Brown
SetupDiscretization(DM dm,AppCtx * ctx)195d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupDiscretization(DM dm, AppCtx *ctx)
196d71ae5a4SJacob Faibussowitsch {
197c4762a1bSJed Brown DM cdm = dm;
198c4762a1bSJed Brown PetscFE fe, feAux;
199c4762a1bSJed Brown MPI_Comm comm;
20030602db0SMatthew G. Knepley PetscInt dim;
20130602db0SMatthew G. Knepley PetscBool simplex;
202c4762a1bSJed Brown
203c4762a1bSJed Brown PetscFunctionBeginUser;
2049566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim));
2059566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex));
2069566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
2079566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, 1, simplex, "phi_", -1, &fe));
2089566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe, "phi"));
2099566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(comm, dim, dim, simplex, "vel_", -1, &feAux));
2109566063dSJacob Faibussowitsch PetscCall(PetscFECopyQuadrature(fe, feAux));
2119566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
2129566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm));
2139566063dSJacob Faibussowitsch PetscCall(SetupProblem(dm, ctx));
214c4762a1bSJed Brown while (cdm) {
2159566063dSJacob Faibussowitsch PetscCall(SetupAuxDM(cdm, feAux, ctx));
2169566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm));
2179566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm));
218c4762a1bSJed Brown }
2199566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe));
2209566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&feAux));
2213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
222c4762a1bSJed Brown }
223c4762a1bSJed Brown
MonitorError(KSP ksp,PetscInt it,PetscReal rnorm,PetscCtx ctx)224*2a8381b2SBarry Smith static PetscErrorCode MonitorError(KSP ksp, PetscInt it, PetscReal rnorm, PetscCtx ctx)
225d71ae5a4SJacob Faibussowitsch {
226c4762a1bSJed Brown DM dm;
22730602db0SMatthew G. Knepley PetscDS ds;
2288434afd1SBarry Smith PetscSimplePointFn *func[1];
22930602db0SMatthew G. Knepley void *ctxs[1];
230c4762a1bSJed Brown Vec u, r, error;
231c4762a1bSJed Brown PetscReal time = 0.5, res;
232c4762a1bSJed Brown
233c4762a1bSJed Brown PetscFunctionBeginUser;
2349566063dSJacob Faibussowitsch PetscCall(KSPGetDM(ksp, &dm));
2359566063dSJacob Faibussowitsch PetscCall(DMSetOutputSequenceNumber(dm, it, time));
236c4762a1bSJed Brown /* Calculate residual */
2379566063dSJacob Faibussowitsch PetscCall(KSPBuildResidual(ksp, NULL, NULL, &r));
2389566063dSJacob Faibussowitsch PetscCall(VecNorm(r, NORM_2, &res));
2399566063dSJacob Faibussowitsch PetscCall(DMSetOutputSequenceNumber(dm, it, res));
2409566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)r, "residual"));
2419566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(r, NULL, "-res_vec_view"));
2429566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r));
243c4762a1bSJed Brown /* Calculate error */
2449566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds));
2459566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 0, &func[0], &ctxs[0]));
2469566063dSJacob Faibussowitsch PetscCall(KSPBuildSolution(ksp, NULL, &u));
2479566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &error));
2489566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, time, func, ctxs, INSERT_ALL_VALUES, error));
2499566063dSJacob Faibussowitsch PetscCall(VecAXPY(error, -1.0, u));
2509566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)error, "error"));
2519566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(error, NULL, "-err_vec_view"));
2529566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &error));
2533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
254c4762a1bSJed Brown }
255c4762a1bSJed Brown
MyTSMonitorError(TS ts,PetscInt step,PetscReal crtime,Vec u,PetscCtx ctx)256*2a8381b2SBarry Smith static PetscErrorCode MyTSMonitorError(TS ts, PetscInt step, PetscReal crtime, Vec u, PetscCtx ctx)
257d71ae5a4SJacob Faibussowitsch {
258c4762a1bSJed Brown DM dm;
25930602db0SMatthew G. Knepley PetscDS ds;
2608434afd1SBarry Smith PetscSimplePointFn *func[1];
26130602db0SMatthew G. Knepley void *ctxs[1];
262c4762a1bSJed Brown PetscReal error;
263c4762a1bSJed Brown
264c4762a1bSJed Brown PetscFunctionBeginUser;
2659566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm));
2669566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds));
2679566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 0, &func[0], &ctxs[0]));
2689566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(dm, crtime, func, ctxs, u, &error));
2699566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Timestep: %04d time = %-8.4g \t L_2 Error: %2.5g\n", (int)step, (double)crtime, (double)error));
2703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
271c4762a1bSJed Brown }
272c4762a1bSJed Brown
main(int argc,char ** argv)273d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
274d71ae5a4SJacob Faibussowitsch {
275c4762a1bSJed Brown AppCtx ctx;
276c4762a1bSJed Brown DM dm;
277c4762a1bSJed Brown TS ts;
278c4762a1bSJed Brown Vec u, r;
279c4762a1bSJed Brown PetscReal t = 0.0;
280c4762a1bSJed Brown
281327415f7SBarry Smith PetscFunctionBeginUser;
2829566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
2839566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx));
2849566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm, &ctx));
2859566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(dm, &ctx));
2869566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, &ctx));
287c4762a1bSJed Brown
2889566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &u));
2899566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)u, "phi"));
2909566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &r));
291c4762a1bSJed Brown
2929566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
2939566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts, MyTSMonitorError, &ctx, NULL));
2949566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, dm));
2959566063dSJacob Faibussowitsch PetscCall(DMTSSetBoundaryLocal(dm, DMPlexTSComputeBoundary, &ctx));
2969566063dSJacob Faibussowitsch PetscCall(DMTSSetIFunctionLocal(dm, DMPlexTSComputeIFunctionFEM, &ctx));
2979566063dSJacob Faibussowitsch PetscCall(DMTSSetIJacobianLocal(dm, DMPlexTSComputeIJacobianFEM, &ctx));
2989566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
2999566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts));
300c4762a1bSJed Brown
30130602db0SMatthew G. Knepley {
30230602db0SMatthew G. Knepley PetscDS ds;
3038434afd1SBarry Smith PetscSimplePointFn *func[1];
30430602db0SMatthew G. Knepley void *ctxs[1];
30530602db0SMatthew G. Knepley
3069566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds));
3079566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 0, &func[0], &ctxs[0]));
3089566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, t, func, ctxs, INSERT_ALL_VALUES, u));
30930602db0SMatthew G. Knepley }
310c4762a1bSJed Brown {
311c4762a1bSJed Brown SNES snes;
312c4762a1bSJed Brown KSP ksp;
313c4762a1bSJed Brown
3149566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes));
3159566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp));
3169566063dSJacob Faibussowitsch PetscCall(KSPMonitorSet(ksp, MonitorError, &ctx, NULL));
317c4762a1bSJed Brown }
3189566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, u));
3199566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-sol_vec_view"));
320c4762a1bSJed Brown
3219566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u));
3229566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r));
3239566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts));
3249566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm));
3259566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
326b122ec5aSJacob Faibussowitsch return 0;
327c4762a1bSJed Brown }
328c4762a1bSJed Brown
329c4762a1bSJed Brown /*TEST
330c4762a1bSJed Brown
331c4762a1bSJed Brown # Full solves
332c4762a1bSJed Brown test:
333c4762a1bSJed Brown suffix: 2d_p1p1_r1
334c4762a1bSJed Brown requires: triangle
335188af4bfSBarry Smith args: -dm_refine 1 -phi_petscspace_degree 1 -vel_petscspace_degree 1 -ts_type beuler -ts_max_steps 10 -ts_time_step 0.1 -pc_type lu -snes_monitor_short -snes_converged_reason -ts_monitor
336c4762a1bSJed Brown
337c4762a1bSJed Brown test:
338c4762a1bSJed Brown suffix: 2d_p1p1_sor_r1
339c4762a1bSJed Brown requires: triangle !single
340188af4bfSBarry Smith args: -dm_refine 1 -phi_petscspace_degree 1 -vel_petscspace_degree 1 -ts_type beuler -ts_max_steps 10 -ts_time_step 0.1 -ksp_rtol 1.0e-9 -pc_type sor -snes_monitor_short -snes_converged_reason -ksp_monitor_short -ts_monitor
341c4762a1bSJed Brown
342c4762a1bSJed Brown test:
343c4762a1bSJed Brown suffix: 2d_p1p1_mg_r1
344c4762a1bSJed Brown requires: triangle !single
345188af4bfSBarry Smith args: -dm_refine_hierarchy 1 -phi_petscspace_degree 1 -vel_petscspace_degree 1 -ts_type beuler -ts_max_steps 10 -ts_time_step 0.1 -ksp_type fgmres -ksp_rtol 1.0e-9 -pc_type mg -pc_mg_levels 2 -snes_monitor_short -snes_converged_reason -snes_view -ksp_monitor_true_residual -ts_monitor
346c4762a1bSJed Brown
347c4762a1bSJed Brown TEST*/
348