1 static char help[] = "Tests implementation of PetscSpace_Sum by solving the Poisson equations using a PetscSpace_Poly and a PetscSpace_Sum and checking that \ 2 solutions agree up to machine precision.\n\n"; 3 4 #include <petscdmplex.h> 5 #include <petscds.h> 6 #include <petscfe.h> 7 #include <petscsnes.h> 8 /* We are solving the system of equations: 9 * \vec{u} = -\grad{p} 10 * \div{u} = f 11 */ 12 13 /* Exact solutions for linear velocity 14 \vec{u} = \vec{x}; 15 p = -0.5*(\vec{x} \cdot \vec{x}); 16 */ 17 static PetscErrorCode linear_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) { 18 PetscInt c; 19 20 for (c = 0; c < Nc; ++c) u[c] = x[c]; 21 return 0; 22 } 23 24 static PetscErrorCode linear_p(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) { 25 PetscInt d; 26 27 u[0] = 0.; 28 for (d = 0; d < dim; ++d) u[0] += -0.5 * x[d] * x[d]; 29 return 0; 30 } 31 32 static PetscErrorCode linear_divu(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) { 33 u[0] = dim; 34 return 0; 35 } 36 37 /* fx_v are the residual functions for the equation \vec{u} = \grad{p}. f0_v is the term <v,u>.*/ 38 static void f0_v(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[]) { 39 PetscInt i; 40 41 for (i = 0; i < dim; ++i) f0[i] = u[uOff[0] + i]; 42 } 43 44 /* f1_v is the term <v,-\grad{p}> but we integrate by parts to get <\grad{v}, -p*I> */ 45 static void f1_v(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[]) { 46 PetscInt c; 47 48 for (c = 0; c < dim; ++c) { 49 PetscInt d; 50 51 for (d = 0; d < dim; ++d) f1[c * dim + d] = (c == d) ? -u[uOff[1]] : 0; 52 } 53 } 54 55 /* Residual function for enforcing \div{u} = f. */ 56 static void f0_q_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[]) { 57 PetscScalar rhs, divu = 0; 58 PetscInt i; 59 60 (void)linear_divu(dim, t, x, dim, &rhs, NULL); 61 for (i = 0; i < dim; ++i) divu += u_x[uOff_x[0] + i * dim + i]; 62 f0[0] = divu - rhs; 63 } 64 65 /* Boundary residual. Dirichlet boundary for u means u_bdy=p*n */ 66 static void f0_bd_u_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[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) { 67 PetscScalar pressure; 68 PetscInt d; 69 70 (void)linear_p(dim, t, x, dim, &pressure, NULL); 71 for (d = 0; d < dim; ++d) f0[d] = pressure * n[d]; 72 } 73 74 /* gx_yz are the jacobian functions obtained by taking the derivative of the y residual w.r.t z*/ 75 static void g0_vu(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[]) { 76 PetscInt c; 77 78 for (c = 0; c < dim; ++c) g0[c * dim + c] = 1.0; 79 } 80 81 static void g1_qu(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[]) { 82 PetscInt c; 83 84 for (c = 0; c < dim; ++c) g1[c * dim + c] = 1.0; 85 } 86 87 static void g2_vp(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[]) { 88 PetscInt c; 89 90 for (c = 0; c < dim; ++c) g2[c * dim + c] = -1.0; 91 } 92 93 typedef struct { 94 PetscInt dummy; 95 } UserCtx; 96 97 static PetscErrorCode CreateMesh(MPI_Comm comm, UserCtx *user, DM *mesh) { 98 PetscFunctionBegin; 99 PetscCall(DMCreate(comm, mesh)); 100 PetscCall(DMSetType(*mesh, DMPLEX)); 101 PetscCall(DMSetFromOptions(*mesh)); 102 PetscCall(DMSetApplicationContext(*mesh, user)); 103 PetscCall(DMViewFromOptions(*mesh, NULL, "-dm_view")); 104 PetscFunctionReturn(0); 105 } 106 107 /* Setup the system of equations that we wish to solve */ 108 static PetscErrorCode SetupProblem(DM dm, UserCtx *user) { 109 PetscDS ds; 110 DMLabel label; 111 PetscWeakForm wf; 112 const PetscInt id = 1; 113 PetscInt bd; 114 115 PetscFunctionBegin; 116 PetscCall(DMGetDS(dm, &ds)); 117 /* All of these are independent of the user's choice of solution */ 118 PetscCall(PetscDSSetResidual(ds, 0, f0_v, f1_v)); 119 PetscCall(PetscDSSetResidual(ds, 1, f0_q_linear, NULL)); 120 PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_vu, NULL, NULL, NULL)); 121 PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_vp, NULL)); 122 PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_qu, NULL, NULL)); 123 124 PetscCall(DMGetLabel(dm, "marker", &label)); 125 PetscCall(PetscDSAddBoundary(ds, DM_BC_NATURAL, "Boundary Integral", label, 1, &id, 0, 0, NULL, (void (*)(void))NULL, NULL, user, &bd)); 126 PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 127 PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, 1, 0, 0, 0, f0_bd_u_linear, 0, NULL)); 128 129 PetscCall(PetscDSSetExactSolution(ds, 0, linear_u, NULL)); 130 PetscCall(PetscDSSetExactSolution(ds, 1, linear_p, NULL)); 131 PetscFunctionReturn(0); 132 } 133 134 /* Create the finite element spaces we will use for this system */ 135 static PetscErrorCode SetupDiscretization(DM mesh, DM mesh_sum, PetscErrorCode (*setup)(DM, UserCtx *), UserCtx *user) { 136 DM cdm = mesh, cdm_sum = mesh_sum; 137 PetscFE u, divu, u_sum, divu_sum; 138 PetscInt dim; 139 PetscBool simplex; 140 141 PetscFunctionBegin; 142 PetscCall(DMGetDimension(mesh, &dim)); 143 PetscCall(DMPlexIsSimplex(mesh, &simplex)); 144 /* Create FE objects and give them names so that options can be set from 145 * command line. Each field gets 2 instances (i.e. velocity and velocity_sum)created twice so that we can compare between approaches. */ 146 PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)mesh), dim, dim, simplex, "velocity_", -1, &u)); 147 PetscCall(PetscObjectSetName((PetscObject)u, "velocity")); 148 PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)mesh_sum), dim, dim, simplex, "velocity_sum_", -1, &u_sum)); 149 PetscCall(PetscObjectSetName((PetscObject)u_sum, "velocity_sum")); 150 PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)mesh), dim, 1, simplex, "divu_", -1, &divu)); 151 PetscCall(PetscObjectSetName((PetscObject)divu, "divu")); 152 PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)mesh_sum), dim, 1, simplex, "divu_sum_", -1, &divu_sum)); 153 PetscCall(PetscObjectSetName((PetscObject)divu_sum, "divu_sum")); 154 155 PetscCall(PetscFECopyQuadrature(u, divu)); 156 PetscCall(PetscFECopyQuadrature(u_sum, divu_sum)); 157 158 /* Associate the FE objects with the mesh and setup the system */ 159 PetscCall(DMSetField(mesh, 0, NULL, (PetscObject)u)); 160 PetscCall(DMSetField(mesh, 1, NULL, (PetscObject)divu)); 161 PetscCall(DMCreateDS(mesh)); 162 PetscCall((*setup)(mesh, user)); 163 164 PetscCall(DMSetField(mesh_sum, 0, NULL, (PetscObject)u_sum)); 165 PetscCall(DMSetField(mesh_sum, 1, NULL, (PetscObject)divu_sum)); 166 PetscCall(DMCreateDS(mesh_sum)); 167 PetscCall((*setup)(mesh_sum, user)); 168 169 while (cdm) { 170 PetscCall(DMCopyDisc(mesh, cdm)); 171 PetscCall(DMGetCoarseDM(cdm, &cdm)); 172 } 173 174 while (cdm_sum) { 175 PetscCall(DMCopyDisc(mesh_sum, cdm_sum)); 176 PetscCall(DMGetCoarseDM(cdm_sum, &cdm_sum)); 177 } 178 179 /* The Mesh now owns the fields, so we can destroy the FEs created in this 180 * function */ 181 PetscCall(PetscFEDestroy(&u)); 182 PetscCall(PetscFEDestroy(&divu)); 183 PetscCall(PetscFEDestroy(&u_sum)); 184 PetscCall(PetscFEDestroy(&divu_sum)); 185 PetscCall(DMDestroy(&cdm)); 186 PetscCall(DMDestroy(&cdm_sum)); 187 PetscFunctionReturn(0); 188 } 189 190 int main(int argc, char **argv) { 191 UserCtx user; 192 DM dm, dm_sum; 193 SNES snes, snes_sum; 194 Vec u, u_sum; 195 PetscReal errNorm; 196 const PetscReal errTol = PETSC_SMALL; 197 198 PetscFunctionBeginUser; 199 PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 200 201 /* Set up a snes for the standard approach, one space with 2 components */ 202 PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 203 PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 204 PetscCall(SNESSetDM(snes, dm)); 205 206 /* Set up a snes for the sum space approach, where each subspace of the sum space represents one component */ 207 PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes_sum)); 208 PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm_sum)); 209 PetscCall(SNESSetDM(snes_sum, dm_sum)); 210 PetscCall(SetupDiscretization(dm, dm_sum, SetupProblem, &user)); 211 212 /* Set up and solve the system using standard approach. */ 213 PetscCall(DMCreateGlobalVector(dm, &u)); 214 PetscCall(VecSet(u, 0.0)); 215 PetscCall(PetscObjectSetName((PetscObject)u, "solution")); 216 PetscCall(DMPlexSetSNESLocalFEM(dm, &user, &user, &user)); 217 PetscCall(SNESSetFromOptions(snes)); 218 PetscCall(DMSNESCheckFromOptions(snes, u)); 219 PetscCall(SNESSolve(snes, NULL, u)); 220 PetscCall(SNESGetSolution(snes, &u)); 221 PetscCall(VecViewFromOptions(u, NULL, "-solution_view")); 222 223 /* Set up and solve the sum space system */ 224 PetscCall(DMCreateGlobalVector(dm_sum, &u_sum)); 225 PetscCall(VecSet(u_sum, 0.0)); 226 PetscCall(PetscObjectSetName((PetscObject)u_sum, "solution_sum")); 227 PetscCall(DMPlexSetSNESLocalFEM(dm_sum, &user, &user, &user)); 228 PetscCall(SNESSetFromOptions(snes_sum)); 229 PetscCall(DMSNESCheckFromOptions(snes_sum, u_sum)); 230 PetscCall(SNESSolve(snes_sum, NULL, u_sum)); 231 PetscCall(SNESGetSolution(snes_sum, &u_sum)); 232 PetscCall(VecViewFromOptions(u_sum, NULL, "-solution_sum_view")); 233 234 /* Check if standard solution and sum space solution match to machine precision */ 235 PetscCall(VecAXPY(u_sum, -1, u)); 236 PetscCall(VecNorm(u_sum, NORM_2, &errNorm)); 237 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Sum space provides the same solution as a regular space: %s", (errNorm < errTol) ? "true" : "false")); 238 239 /* Cleanup */ 240 PetscCall(VecDestroy(&u_sum)); 241 PetscCall(VecDestroy(&u)); 242 PetscCall(SNESDestroy(&snes_sum)); 243 PetscCall(SNESDestroy(&snes)); 244 PetscCall(DMDestroy(&dm_sum)); 245 PetscCall(DMDestroy(&dm)); 246 PetscCall(PetscFinalize()); 247 return 0; 248 } 249 250 /*TEST 251 test: 252 suffix: 2d_lagrange 253 requires: triangle 254 args: -velocity_petscspace_degree 1 \ 255 -velocity_petscspace_type poly \ 256 -velocity_petscspace_components 2\ 257 -velocity_petscdualspace_type lagrange \ 258 -divu_petscspace_degree 0 \ 259 -divu_petscspace_type poly \ 260 -divu_petscdualspace_lagrange_continuity false \ 261 -velocity_sum_petscfe_default_quadrature_order 1 \ 262 -velocity_sum_petscspace_degree 1 \ 263 -velocity_sum_petscspace_type sum \ 264 -velocity_sum_petscspace_variables 2 \ 265 -velocity_sum_petscspace_components 2 \ 266 -velocity_sum_petscspace_sum_spaces 2 \ 267 -velocity_sum_petscspace_sum_concatenate true \ 268 -velocity_sum_petscdualspace_type lagrange \ 269 -velocity_sum_sumcomp_0_petscspace_type poly \ 270 -velocity_sum_sumcomp_0_petscspace_degree 1 \ 271 -velocity_sum_sumcomp_0_petscspace_variables 2 \ 272 -velocity_sum_sumcomp_0_petscspace_components 1 \ 273 -velocity_sum_sumcomp_1_petscspace_type poly \ 274 -velocity_sum_sumcomp_1_petscspace_degree 1 \ 275 -velocity_sum_sumcomp_1_petscspace_variables 2 \ 276 -velocity_sum_sumcomp_1_petscspace_components 1 \ 277 -divu_sum_petscspace_degree 0 \ 278 -divu_sum_petscspace_type sum \ 279 -divu_sum_petscspace_variables 2 \ 280 -divu_sum_petscspace_components 1 \ 281 -divu_sum_petscspace_sum_spaces 1 \ 282 -divu_sum_petscspace_sum_concatenate true \ 283 -divu_sum_petscdualspace_lagrange_continuity false \ 284 -divu_sum_sumcomp_0_petscspace_type poly \ 285 -divu_sum_sumcomp_0_petscspace_degree 0 \ 286 -divu_sum_sumcomp_0_petscspace_variables 2 \ 287 -divu_sum_sumcomp_0_petscspace_components 1 \ 288 -dm_refine 0 \ 289 -snes_error_if_not_converged \ 290 -ksp_rtol 1e-10 \ 291 -ksp_error_if_not_converged \ 292 -pc_type fieldsplit\ 293 -pc_fieldsplit_type schur\ 294 -divu_sum_petscdualspace_lagrange_continuity false \ 295 -pc_fieldsplit_schur_precondition full 296 TEST*/ 297