1d21efd2eSMatthew G. Knepley static const char help[] = "Tests for determining whether a new finite element works";
2d21efd2eSMatthew G. Knepley
3d21efd2eSMatthew G. Knepley /*
4d21efd2eSMatthew G. Knepley Use -interpolation_view and -l2_projection_view to look at the interpolants.
5d21efd2eSMatthew G. Knepley */
6d21efd2eSMatthew G. Knepley
7d21efd2eSMatthew G. Knepley #include <petscdmplex.h>
8d21efd2eSMatthew G. Knepley #include <petscfe.h>
9d21efd2eSMatthew G. Knepley #include <petscds.h>
10d21efd2eSMatthew G. Knepley #include <petscsnes.h>
11d21efd2eSMatthew G. Knepley
constant(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[])12d71ae5a4SJacob Faibussowitsch static void constant(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[])
13d71ae5a4SJacob Faibussowitsch {
14d21efd2eSMatthew G. Knepley const PetscInt Nc = uOff[1] - uOff[0];
15d21efd2eSMatthew G. Knepley for (PetscInt c = 0; c < Nc; ++c) f0[c] += 5.;
16d21efd2eSMatthew G. Knepley }
17d21efd2eSMatthew G. Knepley
linear(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[])18d71ae5a4SJacob Faibussowitsch static void linear(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[])
19d71ae5a4SJacob Faibussowitsch {
20d21efd2eSMatthew G. Knepley const PetscInt Nc = uOff[1] - uOff[0];
21d21efd2eSMatthew G. Knepley for (PetscInt c = 0; c < Nc; ++c) f0[c] += 5. * x[c];
22d21efd2eSMatthew G. Knepley }
23d21efd2eSMatthew G. Knepley
quadratic(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[])24d71ae5a4SJacob Faibussowitsch static void quadratic(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[])
25d71ae5a4SJacob Faibussowitsch {
26d21efd2eSMatthew G. Knepley const PetscInt Nc = uOff[1] - uOff[0];
27d21efd2eSMatthew G. Knepley for (PetscInt c = 0; c < Nc; ++c) f0[c] += 5. * x[c] * x[c];
28d21efd2eSMatthew G. Knepley }
29d21efd2eSMatthew G. Knepley
trig(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[])30d71ae5a4SJacob Faibussowitsch static void trig(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[])
31d71ae5a4SJacob Faibussowitsch {
32d21efd2eSMatthew G. Knepley const PetscInt Nc = uOff[1] - uOff[0];
33d21efd2eSMatthew G. Knepley for (PetscInt c = 0; c < Nc; ++c) f0[c] += PetscCosReal(2. * PETSC_PI * x[c]);
34d21efd2eSMatthew G. Knepley }
35d21efd2eSMatthew G. Knepley
36e239af90SMatthew G. Knepley /*
37e239af90SMatthew G. Knepley The prime basis for the Wheeler-Yotov-Xue prism.
38e239af90SMatthew G. Knepley */
prime(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[])39d71ae5a4SJacob Faibussowitsch static void prime(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[])
40d71ae5a4SJacob Faibussowitsch {
41e239af90SMatthew G. Knepley PetscReal x = X[0], y = X[1], z = X[2], b = 1 + x + y + z;
42e239af90SMatthew G. Knepley f0[0] += b + 2.0 * x * z + 2.0 * y * z + x * y + x * x;
43e239af90SMatthew G. Knepley f0[1] += b + 2.0 * x * z + 2.0 * y * z + x * y + y * y;
44e239af90SMatthew G. Knepley f0[2] += b - 3.0 * x * z - 3.0 * y * z - 2.0 * z * z;
45e239af90SMatthew G. Knepley }
46e239af90SMatthew G. Knepley
47e239af90SMatthew G. Knepley static const char *names[] = {"constant", "linear", "quadratic", "trig", "prime"};
482192575eSBarry Smith static PetscPointFn *functions[] = {constant, linear, quadratic, trig, prime};
49d21efd2eSMatthew G. Knepley
50d21efd2eSMatthew G. Knepley typedef struct {
512192575eSBarry Smith PetscPointFn *exactSol;
52e239af90SMatthew G. Knepley PetscReal shear, flatten;
53d21efd2eSMatthew G. Knepley } AppCtx;
54d21efd2eSMatthew G. Knepley
ProcessOptions(MPI_Comm comm,AppCtx * options)55d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
56d71ae5a4SJacob Faibussowitsch {
57d21efd2eSMatthew G. Knepley char name[PETSC_MAX_PATH_LEN] = "constant";
58dd39110bSPierre Jolivet PetscInt Nfunc = PETSC_STATIC_ARRAY_LENGTH(names), i;
59d21efd2eSMatthew G. Knepley
60d21efd2eSMatthew G. Knepley PetscFunctionBeginUser;
61d21efd2eSMatthew G. Knepley options->exactSol = NULL;
62e239af90SMatthew G. Knepley options->shear = 0.;
63e239af90SMatthew G. Knepley options->flatten = 1.;
64d21efd2eSMatthew G. Knepley
65d0609cedSBarry Smith PetscOptionsBegin(comm, "", "FE Test Options", "PETSCFE");
669566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-func", "Function to project into space", "", name, name, PETSC_MAX_PATH_LEN, NULL));
67f4f49eeaSPierre Jolivet PetscCall(PetscOptionsReal("-shear", "Factor by which to shear along the x-direction", "", options->shear, &options->shear, NULL));
68f4f49eeaSPierre Jolivet PetscCall(PetscOptionsReal("-flatten", "Factor by which to flatten", "", options->flatten, &options->flatten, NULL));
69d0609cedSBarry Smith PetscOptionsEnd();
70d21efd2eSMatthew G. Knepley
71d21efd2eSMatthew G. Knepley for (i = 0; i < Nfunc; ++i) {
72d21efd2eSMatthew G. Knepley PetscBool flg;
73d21efd2eSMatthew G. Knepley
749566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(name, names[i], &flg));
759371c9d4SSatish Balay if (flg) {
769371c9d4SSatish Balay options->exactSol = functions[i];
779371c9d4SSatish Balay break;
789371c9d4SSatish Balay }
79d21efd2eSMatthew G. Knepley }
80e239af90SMatthew G. Knepley PetscCheck(options->exactSol, comm, PETSC_ERR_ARG_WRONG, "Invalid test function %s", name);
813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
82d21efd2eSMatthew G. Knepley }
83d21efd2eSMatthew G. Knepley
84d21efd2eSMatthew G. Knepley /* The exact solution is the negative of the f0 contribution */
exactSolution(PetscInt dim,PetscReal time,const PetscReal x[],PetscInt Nc,PetscScalar * u,PetscCtx ctx)85*2a8381b2SBarry Smith static PetscErrorCode exactSolution(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx)
86d71ae5a4SJacob Faibussowitsch {
87d21efd2eSMatthew G. Knepley AppCtx *user = (AppCtx *)ctx;
88e239af90SMatthew G. Knepley PetscInt uOff[2] = {0, Nc};
89d21efd2eSMatthew G. Knepley
90d21efd2eSMatthew G. Knepley user->exactSol(dim, 1, 0, uOff, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, time, x, 0, NULL, u);
91d21efd2eSMatthew G. Knepley for (PetscInt c = 0; c < Nc; ++c) u[c] *= -1.;
923ba16761SJacob Faibussowitsch return PETSC_SUCCESS;
93d21efd2eSMatthew G. Knepley }
94d21efd2eSMatthew G. Knepley
f0(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[])95d71ae5a4SJacob Faibussowitsch static void f0(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[])
96d71ae5a4SJacob Faibussowitsch {
97d21efd2eSMatthew G. Knepley const PetscInt Nc = uOff[1] - uOff[0];
98d21efd2eSMatthew G. Knepley for (PetscInt c = 0; c < Nc; ++c) f0[c] += u[c];
99d21efd2eSMatthew G. Knepley }
100d21efd2eSMatthew G. Knepley
g0(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[])101d71ae5a4SJacob Faibussowitsch static void g0(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[])
102d71ae5a4SJacob Faibussowitsch {
103d21efd2eSMatthew G. Knepley const PetscInt Nc = uOff[1] - uOff[0];
104d21efd2eSMatthew G. Knepley for (PetscInt c = 0; c < Nc; ++c) g0[c * Nc + c] = 1.0;
105d21efd2eSMatthew G. Knepley }
106d21efd2eSMatthew G. Knepley
SetupProblem(DM dm,AppCtx * user)107d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupProblem(DM dm, AppCtx *user)
108d71ae5a4SJacob Faibussowitsch {
109d21efd2eSMatthew G. Knepley PetscDS ds;
110d21efd2eSMatthew G. Knepley PetscWeakForm wf;
111d21efd2eSMatthew G. Knepley
112d21efd2eSMatthew G. Knepley PetscFunctionBeginUser;
1139566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds));
1149566063dSJacob Faibussowitsch PetscCall(PetscDSGetWeakForm(ds, &wf));
1159566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, NULL, 0, 0, 0, 0, f0, 0, NULL));
1169566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexResidual(wf, NULL, 0, 0, 0, 1, user->exactSol, 0, NULL));
1179566063dSJacob Faibussowitsch PetscCall(PetscWeakFormSetIndexJacobian(wf, NULL, 0, 0, 0, 0, 0, g0, 0, NULL, 0, NULL, 0, NULL));
1189566063dSJacob Faibussowitsch PetscCall(PetscDSSetExactSolution(ds, 0, exactSolution, user));
1193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
120d21efd2eSMatthew G. Knepley }
121d21efd2eSMatthew G. Knepley
SetupDiscretization(DM dm,const char name[],AppCtx * user)122d71ae5a4SJacob Faibussowitsch static PetscErrorCode SetupDiscretization(DM dm, const char name[], AppCtx *user)
123d71ae5a4SJacob Faibussowitsch {
124d21efd2eSMatthew G. Knepley DM cdm = dm;
125d21efd2eSMatthew G. Knepley PetscFE fe;
126d21efd2eSMatthew G. Knepley char prefix[PETSC_MAX_PATH_LEN];
127d21efd2eSMatthew G. Knepley
128d21efd2eSMatthew G. Knepley PetscFunctionBeginUser;
12963a3b9bcSJacob Faibussowitsch if (name) PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%s_", name));
1309566063dSJacob Faibussowitsch PetscCall(DMCreateFEDefault(dm, 1, name ? prefix : NULL, -1, &fe));
1319566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe, name ? name : "Solution"));
132d21efd2eSMatthew G. Knepley /* Set discretization and boundary conditions for each mesh */
1339566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
1349566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm));
1359566063dSJacob Faibussowitsch PetscCall(SetupProblem(dm, user));
136d21efd2eSMatthew G. Knepley while (cdm) {
1379566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm));
1389566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm));
139d21efd2eSMatthew G. Knepley }
1409566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe));
1413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
142d21efd2eSMatthew G. Knepley }
143d21efd2eSMatthew G. Knepley
144d21efd2eSMatthew G. Knepley /* This test tells us whether the given function is contained in the approximation space */
CheckInterpolation(DM dm,AppCtx * user)145d71ae5a4SJacob Faibussowitsch static PetscErrorCode CheckInterpolation(DM dm, AppCtx *user)
146d71ae5a4SJacob Faibussowitsch {
1478434afd1SBarry Smith PetscSimplePointFn *exactSol[1];
148d21efd2eSMatthew G. Knepley void *exactCtx[1];
149d21efd2eSMatthew G. Knepley PetscDS ds;
150d21efd2eSMatthew G. Knepley Vec u;
151d21efd2eSMatthew G. Knepley PetscReal error, tol = PETSC_SMALL;
152d21efd2eSMatthew G. Knepley MPI_Comm comm;
153d21efd2eSMatthew G. Knepley
154d21efd2eSMatthew G. Knepley PetscFunctionBeginUser;
1559566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
1569566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds));
1579566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &u));
1589566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 0, &exactSol[0], &exactCtx[0]));
1599566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, 0.0, exactSol, exactCtx, INSERT_ALL_VALUES, u));
1609566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-interpolation_view"));
1619566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(dm, 0.0, exactSol, exactCtx, u, &error));
1629566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &u));
1639566063dSJacob Faibussowitsch if (error > tol) PetscCall(PetscPrintf(comm, "Interpolation tests FAIL at tolerance %g error %g\n", (double)tol, (double)error));
1649566063dSJacob Faibussowitsch else PetscCall(PetscPrintf(comm, "Interpolation tests pass at tolerance %g\n", (double)tol));
1653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
166d21efd2eSMatthew G. Knepley }
167d21efd2eSMatthew G. Knepley
168d21efd2eSMatthew G. Knepley /* This test tells us whether the element is unisolvent (the mass matrix has full rank), and what rate of convergence we achieve */
CheckL2Projection(DM dm,AppCtx * user)169d71ae5a4SJacob Faibussowitsch static PetscErrorCode CheckL2Projection(DM dm, AppCtx *user)
170d71ae5a4SJacob Faibussowitsch {
1718434afd1SBarry Smith PetscSimplePointFn *exactSol[1];
172d21efd2eSMatthew G. Knepley void *exactCtx[1];
173d21efd2eSMatthew G. Knepley SNES snes;
174d21efd2eSMatthew G. Knepley PetscDS ds;
175d21efd2eSMatthew G. Knepley Vec u;
176d21efd2eSMatthew G. Knepley PetscReal error, tol = PETSC_SMALL;
177d21efd2eSMatthew G. Knepley MPI_Comm comm;
178d21efd2eSMatthew G. Knepley
179d21efd2eSMatthew G. Knepley PetscFunctionBeginUser;
1809566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
1819566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &ds));
1829566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &u));
1839566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, 0, &exactSol[0], &exactCtx[0]));
1849566063dSJacob Faibussowitsch PetscCall(SNESCreate(comm, &snes));
1859566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, dm));
1869566063dSJacob Faibussowitsch PetscCall(VecSet(u, 0.0));
1879566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)u, "solution"));
1886493148fSStefano Zampini PetscCall(DMPlexSetSNESLocalFEM(dm, PETSC_FALSE, user));
1899566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes));
1909566063dSJacob Faibussowitsch PetscCall(DMSNESCheckFromOptions(snes, u));
1919566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, u));
1929566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes));
1939566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-l2_projection_view"));
1949566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(dm, 0.0, exactSol, exactCtx, u, &error));
1959566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &u));
1969566063dSJacob Faibussowitsch if (error > tol) PetscCall(PetscPrintf(comm, "L2 projection tests FAIL at tolerance %g error %g\n", (double)tol, (double)error));
1979566063dSJacob Faibussowitsch else PetscCall(PetscPrintf(comm, "L2 projection tests pass at tolerance %g\n", (double)tol));
1983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
199d21efd2eSMatthew G. Knepley }
200d21efd2eSMatthew G. Knepley
201e239af90SMatthew G. Knepley /* Distorts the mesh by shearing in the x-direction and flattening, factors provided in the options. */
DistortMesh(DM dm,AppCtx * user)202d71ae5a4SJacob Faibussowitsch static PetscErrorCode DistortMesh(DM dm, AppCtx *user)
203d71ae5a4SJacob Faibussowitsch {
204e239af90SMatthew G. Knepley Vec coordinates;
205e239af90SMatthew G. Knepley PetscScalar *ca;
206e239af90SMatthew G. Knepley PetscInt dE, n, i;
207e239af90SMatthew G. Knepley
208e239af90SMatthew G. Knepley PetscFunctionBeginUser;
2099566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &dE));
2109566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(dm, &coordinates));
2119566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(coordinates, &n));
2129566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordinates, &ca));
213e239af90SMatthew G. Knepley for (i = 0; i < (n / dE); ++i) {
214e239af90SMatthew G. Knepley ca[i * dE + 0] += user->shear * ca[i * dE + 0];
215e239af90SMatthew G. Knepley ca[i * dE + 1] *= user->flatten;
216e239af90SMatthew G. Knepley }
2179566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordinates, &ca));
2183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
219e239af90SMatthew G. Knepley }
220e239af90SMatthew G. Knepley
main(int argc,char ** argv)221d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
222d71ae5a4SJacob Faibussowitsch {
223d21efd2eSMatthew G. Knepley DM dm;
224d21efd2eSMatthew G. Knepley AppCtx user;
225d21efd2eSMatthew G. Knepley PetscMPIInt size;
226d21efd2eSMatthew G. Knepley
227327415f7SBarry Smith PetscFunctionBeginUser;
2289566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
2299566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
230bd158744SPierre Jolivet PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only.");
2319566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
2329566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
2339566063dSJacob Faibussowitsch PetscCall(DMSetType(dm, DMPLEX));
2349566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm));
2359566063dSJacob Faibussowitsch PetscCall(DistortMesh(dm, &user));
2369566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
2379566063dSJacob Faibussowitsch PetscCall(SetupDiscretization(dm, NULL, &user));
238d21efd2eSMatthew G. Knepley
2399566063dSJacob Faibussowitsch PetscCall(CheckInterpolation(dm, &user));
2409566063dSJacob Faibussowitsch PetscCall(CheckL2Projection(dm, &user));
241d21efd2eSMatthew G. Knepley
2429566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm));
2439566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
244b122ec5aSJacob Faibussowitsch return 0;
245d21efd2eSMatthew G. Knepley }
246d21efd2eSMatthew G. Knepley
247d21efd2eSMatthew G. Knepley /*TEST
248d21efd2eSMatthew G. Knepley
249d21efd2eSMatthew G. Knepley testset:
250d21efd2eSMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -petscspace_degree 1\
251d21efd2eSMatthew G. Knepley -snes_error_if_not_converged -ksp_error_if_not_converged -pc_type lu
252d21efd2eSMatthew G. Knepley
253d21efd2eSMatthew G. Knepley test:
254d21efd2eSMatthew G. Knepley suffix: p1_0
255d21efd2eSMatthew G. Knepley args: -func {{constant linear}}
256d21efd2eSMatthew G. Knepley
257d21efd2eSMatthew G. Knepley # Using -dm_refine 2 -convest_num_refine 4 gives convergence rate 2.0
258d21efd2eSMatthew G. Knepley test:
259d21efd2eSMatthew G. Knepley suffix: p1_1
260d21efd2eSMatthew G. Knepley args: -func {{quadratic trig}} \
261ea5adb04SJed Brown -snes_convergence_estimate -convest_num_refine 2 -dm_refine 1
262d21efd2eSMatthew G. Knepley
263d21efd2eSMatthew G. Knepley testset:
264e239af90SMatthew G. Knepley requires: !complex double
265d21efd2eSMatthew G. Knepley args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism \
266d21efd2eSMatthew G. Knepley -petscspace_type sum \
267d21efd2eSMatthew G. Knepley -petscspace_variables 3 \
268d21efd2eSMatthew G. Knepley -petscspace_components 3 \
269d21efd2eSMatthew G. Knepley -petscspace_sum_spaces 2 \
270d21efd2eSMatthew G. Knepley -petscspace_sum_concatenate false \
271d21efd2eSMatthew G. Knepley -sumcomp_0_petscspace_variables 3 \
272d21efd2eSMatthew G. Knepley -sumcomp_0_petscspace_components 3 \
273d21efd2eSMatthew G. Knepley -sumcomp_0_petscspace_degree 1 \
274d21efd2eSMatthew G. Knepley -sumcomp_1_petscspace_variables 3 \
275d21efd2eSMatthew G. Knepley -sumcomp_1_petscspace_components 3 \
276d21efd2eSMatthew G. Knepley -sumcomp_1_petscspace_type wxy \
277d21efd2eSMatthew G. Knepley -petscdualspace_form_degree 0 \
278d21efd2eSMatthew G. Knepley -petscdualspace_order 1 \
279d21efd2eSMatthew G. Knepley -petscdualspace_components 3 \
280d21efd2eSMatthew G. Knepley -snes_error_if_not_converged -ksp_error_if_not_converged -pc_type lu
281d21efd2eSMatthew G. Knepley
282d21efd2eSMatthew G. Knepley test:
283d21efd2eSMatthew G. Knepley suffix: wxy_0
284d21efd2eSMatthew G. Knepley args: -func constant
285d21efd2eSMatthew G. Knepley
286d21efd2eSMatthew G. Knepley test:
287d21efd2eSMatthew G. Knepley suffix: wxy_1
288d21efd2eSMatthew G. Knepley args: -func linear
289d21efd2eSMatthew G. Knepley
290e239af90SMatthew G. Knepley test:
291e239af90SMatthew G. Knepley suffix: wxy_2
292e239af90SMatthew G. Knepley args: -func prime
293e239af90SMatthew G. Knepley
294e239af90SMatthew G. Knepley test:
295e239af90SMatthew G. Knepley suffix: wxy_3
296e239af90SMatthew G. Knepley args: -func linear -shear 1 -flatten 1e-5
297e239af90SMatthew G. Knepley
298d21efd2eSMatthew G. Knepley TEST*/
299