147bfdf3fSMatthew G. Knepley static char help[] = "Poisson Problem in 2d and 3d with simplicial finite elements in both primal and mixed form.\n\ 247bfdf3fSMatthew G. Knepley We solve the Poisson problem in a rectangular\n\ 347bfdf3fSMatthew G. Knepley domain, using a parallel unstructured mesh (DMPLEX) to discretize it.\n\ 447bfdf3fSMatthew G. Knepley This example solves mixed form equation to get the flux field to calculate flux norm. We then use that for adaptive mesh refinement. \n\n\n"; 547bfdf3fSMatthew G. Knepley 647bfdf3fSMatthew G. Knepley /* 747bfdf3fSMatthew G. Knepley The primal (or original) Poisson problem, in the strong form, is given by, 847bfdf3fSMatthew G. Knepley 947bfdf3fSMatthew G. Knepley \begin{align} 1047bfdf3fSMatthew G. Knepley - \nabla \cdot ( \nabla u ) = f 1147bfdf3fSMatthew G. Knepley \end{align} 1247bfdf3fSMatthew G. Knepley where $u$ is potential. 1347bfdf3fSMatthew G. Knepley 1447bfdf3fSMatthew G. Knepley The weak form of this equation is 1547bfdf3fSMatthew G. Knepley 1647bfdf3fSMatthew G. Knepley \begin{align} 1747bfdf3fSMatthew G. Knepley < \nabla v, \nabla u > - < v, \nabla u \cdot \hat{n} >_\Gamma - < v, f > = 0 1847bfdf3fSMatthew G. Knepley \end{align} 1947bfdf3fSMatthew G. Knepley 2047bfdf3fSMatthew G. Knepley The mixed Poisson problem, in the strong form, is given by, 2147bfdf3fSMatthew G. Knepley 2247bfdf3fSMatthew G. Knepley \begin{align} 2347bfdf3fSMatthew G. Knepley q - \nabla u &= 0 \\ 2447bfdf3fSMatthew G. Knepley - \nabla \cdot q &= f 2547bfdf3fSMatthew G. Knepley \end{align} 2647bfdf3fSMatthew G. Knepley where $u$ is the potential and $q$ is the flux. 2747bfdf3fSMatthew G. Knepley 2847bfdf3fSMatthew G. Knepley The weak form of this equation is 2947bfdf3fSMatthew G. Knepley 3047bfdf3fSMatthew G. Knepley \begin{align} 3147bfdf3fSMatthew G. Knepley < t, q > + < \nabla \cdot t, u > - < t \cdot \hat{n}, u >_\Gamma &= 0 \\ 3247bfdf3fSMatthew G. Knepley <v, \nabla \cdot q> - < v, f > &= 0 3347bfdf3fSMatthew G. Knepley \end{align} 3447bfdf3fSMatthew G. Knepley 3547bfdf3fSMatthew G. Knepley We solve both primal and mixed problem and calculate the error in the flux norm, namely || e || = || q^m_h - \nabla u^p_h ||. Here superscript 'm' represents field from mixed form and 'p' represents field from the primal form. 3647bfdf3fSMatthew G. Knepley 3747bfdf3fSMatthew G. Knepley The following boundary conditions are prescribed. 3847bfdf3fSMatthew G. Knepley 3947bfdf3fSMatthew G. Knepley Primal problem: 4047bfdf3fSMatthew G. Knepley \begin{align} 4147bfdf3fSMatthew G. Knepley u = u_0 on \Gamma_D 4247bfdf3fSMatthew G. Knepley \nabla u \cdot \hat{n} = g on \Gamma_N 4347bfdf3fSMatthew G. Knepley \end{align} 4447bfdf3fSMatthew G. Knepley 4547bfdf3fSMatthew G. Knepley Mixed problem: 4647bfdf3fSMatthew G. Knepley \begin{align} 4747bfdf3fSMatthew G. Knepley u = u_0 on \Gamma_D 4847bfdf3fSMatthew G. Knepley q \cdot \hat{n} = g on \Gamma_N 4947bfdf3fSMatthew G. Knepley \end{align} 5047bfdf3fSMatthew G. Knepley __________\Gamma_D_____________ 5147bfdf3fSMatthew G. Knepley | | 5247bfdf3fSMatthew G. Knepley | | 5347bfdf3fSMatthew G. Knepley | | 5447bfdf3fSMatthew G. Knepley \Gamma_N \Gamma_N 5547bfdf3fSMatthew G. Knepley | | 5647bfdf3fSMatthew G. Knepley | | 5747bfdf3fSMatthew G. Knepley | | 5847bfdf3fSMatthew G. Knepley |_________\Gamma_D_____________| 5947bfdf3fSMatthew G. Knepley 6047bfdf3fSMatthew G. Knepley To visualize the automated adaptation 6147bfdf3fSMatthew G. Knepley 6247bfdf3fSMatthew G. Knepley -dm_adapt_pre_view draw -dm_adapt_view draw -draw_pause -1 -geometry 0,0,1024,1024 6347bfdf3fSMatthew G. Knepley 6447bfdf3fSMatthew G. Knepley and to compare with a naice gradient estimator use 6547bfdf3fSMatthew G. Knepley 6647bfdf3fSMatthew G. Knepley -adaptor_type gradient 6747bfdf3fSMatthew G. Knepley 6847bfdf3fSMatthew G. Knepley To see a sequence of adaptations 6947bfdf3fSMatthew G. Knepley 7047bfdf3fSMatthew G. Knepley -snes_adapt_sequence 8 -adaptor_monitor_error draw::draw_lg 7147bfdf3fSMatthew G. Knepley -dm_adapt_pre_view draw -dm_adapt_iter_view draw -dm_adapt_view draw -draw_pause 1 -geometry 0,0,1024,1024 7247bfdf3fSMatthew G. Knepley 7347bfdf3fSMatthew G. Knepley To get a better view of the by-hand process, use 7447bfdf3fSMatthew G. Knepley 7547bfdf3fSMatthew G. Knepley -dm_view hdf5:${PWD}/mesh.h5 7647bfdf3fSMatthew G. Knepley -primal_sol_vec_view hdf5:${PWD}/mesh.h5::append 7747bfdf3fSMatthew G. Knepley -flux_error_vec_view hdf5:${PWD}/mesh.h5::append 7847bfdf3fSMatthew G. Knepley -exact_error_vec_view hdf5:${PWD}/mesh.h5::append 7947bfdf3fSMatthew G. Knepley -refine_vec_view hdf5:${PWD}/mesh.h5::append 8047bfdf3fSMatthew G. Knepley -adapt_dm_view draw -draw_pause -1 8147bfdf3fSMatthew G. Knepley 8247bfdf3fSMatthew G. Knepley This is also possible with the automated path 8347bfdf3fSMatthew G. Knepley 8447bfdf3fSMatthew G. Knepley -dm_view hdf5:${PWD}/mesh.h5 8547bfdf3fSMatthew G. Knepley -adapt_primal_sol_vec_view hdf5:${PWD}/mesh.h5::append 8647bfdf3fSMatthew G. Knepley -adapt_error_vec_view hdf5:${PWD}/mesh.h5::append 8747bfdf3fSMatthew G. Knepley -adapt_vec_view hdf5:${PWD}/mesh.h5::append 8847bfdf3fSMatthew G. Knepley */ 8947bfdf3fSMatthew G. Knepley 9047bfdf3fSMatthew G. Knepley #include <petsc/private/petscfeimpl.h> 9147bfdf3fSMatthew G. Knepley #include <petscdmplex.h> 9247bfdf3fSMatthew G. Knepley #include <petscsnes.h> 9347bfdf3fSMatthew G. Knepley #include <petscdmadaptor.h> 9447bfdf3fSMatthew G. Knepley #include <petscds.h> 9547bfdf3fSMatthew G. Knepley #include <petscviewerhdf5.h> 9647bfdf3fSMatthew G. Knepley #include <petscbag.h> 9747bfdf3fSMatthew G. Knepley 9847bfdf3fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode SetupMixed(DMAdaptor, DM); 9947bfdf3fSMatthew G. Knepley 10047bfdf3fSMatthew G. Knepley typedef enum { 10147bfdf3fSMatthew G. Knepley SOL_QUADRATIC, 10247bfdf3fSMatthew G. Knepley SOL_TRIG, 10347bfdf3fSMatthew G. Knepley SOL_SENSOR, 10447bfdf3fSMatthew G. Knepley SOL_UNKNOWN, 10547bfdf3fSMatthew G. Knepley NUM_SOL_TYPE 10647bfdf3fSMatthew G. Knepley } SolType; 10747bfdf3fSMatthew G. Knepley const char *SolTypeNames[NUM_SOL_TYPE + 4] = {"quadratic", "trig", "sensor", "unknown", "SolType", "SOL_", NULL}; 10847bfdf3fSMatthew G. Knepley 10947bfdf3fSMatthew G. Knepley typedef struct { 11047bfdf3fSMatthew G. Knepley PetscBag param; 11147bfdf3fSMatthew G. Knepley SolType solType; 11247bfdf3fSMatthew G. Knepley PetscBool byHand; 113a5cee669SMatthew G. Knepley PetscInt numAdapt; 11447bfdf3fSMatthew G. Knepley } AppCtx; 11547bfdf3fSMatthew G. Knepley 11647bfdf3fSMatthew G. Knepley typedef struct { 11747bfdf3fSMatthew G. Knepley PetscReal k; 11847bfdf3fSMatthew G. Knepley } Parameter; 11947bfdf3fSMatthew G. Knepley 12047bfdf3fSMatthew G. Knepley /* Exact solution: u = x^2 + y^2 */ 12147bfdf3fSMatthew G. Knepley static PetscErrorCode quadratic_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 12247bfdf3fSMatthew G. Knepley { 12347bfdf3fSMatthew G. Knepley PetscInt d; 12447bfdf3fSMatthew G. Knepley 12547bfdf3fSMatthew G. Knepley u[0] = 0.0; 12647bfdf3fSMatthew G. Knepley for (d = 0; d < dim; ++d) u[0] += x[d] * x[d]; 12747bfdf3fSMatthew G. Knepley return PETSC_SUCCESS; 12847bfdf3fSMatthew G. Knepley } 12947bfdf3fSMatthew G. Knepley /* Exact solution: q = (2x, 2y) */ 13047bfdf3fSMatthew G. Knepley static PetscErrorCode quadratic_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 13147bfdf3fSMatthew G. Knepley { 13247bfdf3fSMatthew G. Knepley PetscInt c; 13347bfdf3fSMatthew G. Knepley for (c = 0; c < Nc; ++c) u[c] = 2.0 * x[c]; 13447bfdf3fSMatthew G. Knepley return PETSC_SUCCESS; 13547bfdf3fSMatthew G. Knepley } 13647bfdf3fSMatthew G. Knepley 13747bfdf3fSMatthew G. Knepley /* Exact solution: u = sin( n \pi x ) * sin( n \pi y ) */ 13847bfdf3fSMatthew G. Knepley static PetscErrorCode trig_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 13947bfdf3fSMatthew G. Knepley { 14047bfdf3fSMatthew G. Knepley const PetscReal n = 5.5; 14147bfdf3fSMatthew G. Knepley 14247bfdf3fSMatthew G. Knepley u[0] = 1.0; 14347bfdf3fSMatthew G. Knepley for (PetscInt d = 0; d < dim; ++d) u[0] *= PetscSinReal(n * PETSC_PI * x[d]); 14447bfdf3fSMatthew G. Knepley return PETSC_SUCCESS; 14547bfdf3fSMatthew G. Knepley } 14647bfdf3fSMatthew G. Knepley static PetscErrorCode trig_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 14747bfdf3fSMatthew G. Knepley { 14847bfdf3fSMatthew G. Knepley const PetscReal n = 5.5; 14947bfdf3fSMatthew G. Knepley 15047bfdf3fSMatthew G. Knepley for (PetscInt c = 0; c < Nc; c++) u[c] = n * PETSC_PI * PetscCosReal(n * PETSC_PI * x[c]) * PetscSinReal(n * PETSC_PI * x[Nc - c - 1]); 15147bfdf3fSMatthew G. Knepley return PETSC_SUCCESS; 15247bfdf3fSMatthew G. Knepley } 15347bfdf3fSMatthew G. Knepley 15447bfdf3fSMatthew G. Knepley /* 15547bfdf3fSMatthew G. Knepley Classic hyperbolic sensor function for testing multi-scale anisotropic mesh adaptation: 15647bfdf3fSMatthew G. Knepley 15747bfdf3fSMatthew G. Knepley f:[-1, 1]x[-1, 1] \to R, 15847bfdf3fSMatthew G. Knepley f(x, y) = sin(50xy)/100 if |xy| > 2\pi/50 else sin(50xy) 15947bfdf3fSMatthew G. Knepley 16047bfdf3fSMatthew G. Knepley (mapped to have domain [0,1] x [0,1] in this case). 16147bfdf3fSMatthew G. Knepley */ 16247bfdf3fSMatthew G. Knepley static PetscErrorCode sensor_u(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx) 16347bfdf3fSMatthew G. Knepley { 16447bfdf3fSMatthew G. Knepley const PetscReal xref = 2. * x[0] - 1.; 16547bfdf3fSMatthew G. Knepley const PetscReal yref = 2. * x[1] - 1.; 16647bfdf3fSMatthew G. Knepley const PetscReal xy = xref * yref; 16747bfdf3fSMatthew G. Knepley 16847bfdf3fSMatthew G. Knepley u[0] = PetscSinReal(50. * xy); 16947bfdf3fSMatthew G. Knepley if (PetscAbsReal(xy) > 2. * PETSC_PI / 50.) u[0] *= 0.01; 17047bfdf3fSMatthew G. Knepley 17147bfdf3fSMatthew G. Knepley return PETSC_SUCCESS; 17247bfdf3fSMatthew G. Knepley } 17347bfdf3fSMatthew G. Knepley 17447bfdf3fSMatthew G. Knepley /* Flux is (cos(50xy) * 50y/100, cos(50xy) * 50x/100) if |xy| > 2\pi/50 else (cos(50xy) * 50y, cos(50xy) * 50x) */ 17547bfdf3fSMatthew G. Knepley static PetscErrorCode sensor_q(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx) 17647bfdf3fSMatthew G. Knepley { 17747bfdf3fSMatthew G. Knepley const PetscReal xref = 2. * x[0] - 1.; 17847bfdf3fSMatthew G. Knepley const PetscReal yref = 2. * x[1] - 1.; 17947bfdf3fSMatthew G. Knepley const PetscReal xy = xref * yref; 18047bfdf3fSMatthew G. Knepley 18147bfdf3fSMatthew G. Knepley u[0] = 50. * yref * PetscCosReal(50. * xy) * 2.0; 18247bfdf3fSMatthew G. Knepley u[1] = 50. * xref * PetscCosReal(50. * xy) * 2.0; 18347bfdf3fSMatthew G. Knepley if (PetscAbsReal(xy) > 2. * PETSC_PI / 50.) { 18447bfdf3fSMatthew G. Knepley u[0] *= 0.01; 18547bfdf3fSMatthew G. Knepley u[1] *= 0.01; 18647bfdf3fSMatthew G. Knepley } 18747bfdf3fSMatthew G. Knepley return PETSC_SUCCESS; 18847bfdf3fSMatthew G. Knepley } 18947bfdf3fSMatthew G. Knepley 19047bfdf3fSMatthew G. Knepley /* We set up residuals and Jacobians for the primal problem. */ 19147bfdf3fSMatthew G. Knepley static void f0_quadratic_primal(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[]) 19247bfdf3fSMatthew G. Knepley { 19347bfdf3fSMatthew G. Knepley f0[0] = 4.0; 19447bfdf3fSMatthew G. Knepley } 19547bfdf3fSMatthew G. Knepley 19647bfdf3fSMatthew G. Knepley static void f0_trig_primal(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[]) 19747bfdf3fSMatthew G. Knepley { 19847bfdf3fSMatthew G. Knepley const PetscReal n = 5.5; 19947bfdf3fSMatthew G. Knepley 20047bfdf3fSMatthew G. Knepley f0[0] = -2.0 * PetscSqr(n * PETSC_PI) * PetscSinReal(n * PETSC_PI * x[0]) * PetscSinReal(n * PETSC_PI * x[1]); 20147bfdf3fSMatthew G. Knepley } 20247bfdf3fSMatthew G. Knepley 20347bfdf3fSMatthew G. Knepley static void f0_sensor_primal(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[]) 20447bfdf3fSMatthew G. Knepley { 20547bfdf3fSMatthew G. Knepley const PetscReal xref = 2. * x[0] - 1.; 20647bfdf3fSMatthew G. Knepley const PetscReal yref = 2. * x[1] - 1.; 20747bfdf3fSMatthew G. Knepley const PetscReal xy = xref * yref; 20847bfdf3fSMatthew G. Knepley 20947bfdf3fSMatthew G. Knepley f0[0] = -2500.0 * PetscSinReal(50. * xy) * (xref * xref + yref * yref) * 4.0; 21047bfdf3fSMatthew G. Knepley if (PetscAbsReal(xy) > 2. * PETSC_PI / 50.) f0[0] *= 0.01; 21147bfdf3fSMatthew G. Knepley } 21247bfdf3fSMatthew G. Knepley 21347bfdf3fSMatthew G. Knepley static void f1_primal(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[]) 21447bfdf3fSMatthew G. Knepley { 21547bfdf3fSMatthew G. Knepley const PetscReal k = PetscRealPart(constants[0]); 21647bfdf3fSMatthew G. Knepley 21747bfdf3fSMatthew G. Knepley for (PetscInt d = 0; d < dim; ++d) f1[d] = k * u_x[uOff_x[0] + d]; 21847bfdf3fSMatthew G. Knepley } 21947bfdf3fSMatthew G. Knepley 22047bfdf3fSMatthew G. Knepley static void f0_quadratic_bd_primal(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[]) 22147bfdf3fSMatthew G. Knepley { 22247bfdf3fSMatthew G. Knepley const PetscReal k = PetscRealPart(constants[0]); 22347bfdf3fSMatthew G. Knepley PetscScalar flux; 22447bfdf3fSMatthew G. Knepley 22547bfdf3fSMatthew G. Knepley PetscCallAbort(PETSC_COMM_SELF, quadratic_q(dim, t, x, dim, &flux, NULL)); 22647bfdf3fSMatthew G. Knepley for (PetscInt d = 0; d < dim; ++d) f0[d] = -k * flux * n[d]; 22747bfdf3fSMatthew G. Knepley } 22847bfdf3fSMatthew G. Knepley 22947bfdf3fSMatthew G. Knepley static void f0_trig_bd_primal(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[]) 23047bfdf3fSMatthew G. Knepley { 23147bfdf3fSMatthew G. Knepley const PetscReal k = PetscRealPart(constants[0]); 23247bfdf3fSMatthew G. Knepley PetscScalar flux; 23347bfdf3fSMatthew G. Knepley 23447bfdf3fSMatthew G. Knepley PetscCallAbort(PETSC_COMM_SELF, trig_q(dim, t, x, dim, &flux, NULL)); 23547bfdf3fSMatthew G. Knepley for (PetscInt d = 0; d < dim; ++d) f0[d] = -k * flux * n[d]; 23647bfdf3fSMatthew G. Knepley } 23747bfdf3fSMatthew G. Knepley 23847bfdf3fSMatthew G. Knepley static void f0_sensor_bd_primal(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[]) 23947bfdf3fSMatthew G. Knepley { 24047bfdf3fSMatthew G. Knepley const PetscReal k = PetscRealPart(constants[0]); 241*f5c5fea7SStefano Zampini PetscScalar flux[2]; 24247bfdf3fSMatthew G. Knepley 243*f5c5fea7SStefano Zampini PetscCallAbort(PETSC_COMM_SELF, sensor_q(dim, t, x, dim, flux, NULL)); 244*f5c5fea7SStefano Zampini for (PetscInt d = 0; d < dim; ++d) f0[d] = -k * flux[d] * n[d]; 24547bfdf3fSMatthew G. Knepley } 24647bfdf3fSMatthew G. Knepley 24747bfdf3fSMatthew G. Knepley static void g3_primal_uu(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 g3[]) 24847bfdf3fSMatthew G. Knepley { 24947bfdf3fSMatthew G. Knepley const PetscReal k = PetscRealPart(constants[0]); 25047bfdf3fSMatthew G. Knepley 25147bfdf3fSMatthew G. Knepley for (PetscInt d = 0; d < dim; ++d) g3[d * dim + d] = k; 25247bfdf3fSMatthew G. Knepley } 25347bfdf3fSMatthew G. Knepley 25447bfdf3fSMatthew G. Knepley /* Now we set up the residuals and Jacobians mixed problem. */ 25547bfdf3fSMatthew G. Knepley static void f0_mixed_quadratic_u(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[]) 25647bfdf3fSMatthew G. Knepley { 25747bfdf3fSMatthew G. Knepley f0[0] = 4.0; 25847bfdf3fSMatthew G. Knepley for (PetscInt d = 0; d < dim; ++d) f0[0] += -u_x[uOff_x[0] + d * dim + d]; 25947bfdf3fSMatthew G. Knepley } 26047bfdf3fSMatthew G. Knepley static void f0_mixed_trig_u(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[]) 26147bfdf3fSMatthew G. Knepley { 26247bfdf3fSMatthew G. Knepley PetscReal n = 5.5; 26347bfdf3fSMatthew G. Knepley 26447bfdf3fSMatthew G. Knepley f0[0] = -2.0 * PetscSqr(n * PETSC_PI) * PetscSinReal(n * PETSC_PI * x[0]) * PetscSinReal(n * PETSC_PI * x[1]); 26547bfdf3fSMatthew G. Knepley for (PetscInt d = 0; d < dim; ++d) f0[0] += -u_x[uOff_x[0] + d * dim + d]; 26647bfdf3fSMatthew G. Knepley } 26747bfdf3fSMatthew G. Knepley static void f0_mixed_sensor_u(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[]) 26847bfdf3fSMatthew G. Knepley { 26947bfdf3fSMatthew G. Knepley const PetscReal xref = 2. * x[0] - 1.; 27047bfdf3fSMatthew G. Knepley const PetscReal yref = 2. * x[1] - 1.; 27147bfdf3fSMatthew G. Knepley const PetscReal xy = xref * yref; 27247bfdf3fSMatthew G. Knepley 27347bfdf3fSMatthew G. Knepley f0[0] = -2500.0 * PetscSinReal(50. * xy) * (xref * xref + yref * yref) * 4.0; 27447bfdf3fSMatthew G. Knepley if (PetscAbsReal(xy) > 2. * PETSC_PI / 50.) f0[0] *= 0.01; 27547bfdf3fSMatthew G. Knepley for (PetscInt d = 0; d < dim; ++d) f0[0] += -u_x[uOff_x[0] + d * dim + d]; 27647bfdf3fSMatthew G. Knepley } 27747bfdf3fSMatthew G. Knepley 27847bfdf3fSMatthew G. Knepley static void f0_mixed_q(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[]) 27947bfdf3fSMatthew G. Knepley { 28047bfdf3fSMatthew G. Knepley for (PetscInt d = 0; d < dim; d++) f0[d] = u[uOff[0] + d]; 28147bfdf3fSMatthew G. Knepley } 28247bfdf3fSMatthew G. Knepley 28347bfdf3fSMatthew G. Knepley static void f1_mixed_q(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[]) 28447bfdf3fSMatthew G. Knepley { 28547bfdf3fSMatthew G. Knepley const PetscReal k = PetscRealPart(constants[0]); 28647bfdf3fSMatthew G. Knepley 28747bfdf3fSMatthew G. Knepley for (PetscInt d = 0; d < dim; d++) f1[d * dim + d] = k * u[uOff[1]]; 28847bfdf3fSMatthew G. Knepley } 28947bfdf3fSMatthew G. Knepley 29047bfdf3fSMatthew G. Knepley static void f0_quadratic_bd_mixed_q(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[]) 29147bfdf3fSMatthew G. Knepley { 29247bfdf3fSMatthew G. Knepley const PetscReal k = PetscRealPart(constants[0]); 29347bfdf3fSMatthew G. Knepley PetscScalar potential; 29447bfdf3fSMatthew G. Knepley 29547bfdf3fSMatthew G. Knepley PetscCallAbort(PETSC_COMM_SELF, quadratic_u(dim, t, x, dim, &potential, NULL)); 29647bfdf3fSMatthew G. Knepley for (PetscInt d = 0; d < dim; ++d) f0[d] = -k * potential * n[d]; 29747bfdf3fSMatthew G. Knepley } 29847bfdf3fSMatthew G. Knepley 29947bfdf3fSMatthew G. Knepley static void f0_trig_bd_mixed_q(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[]) 30047bfdf3fSMatthew G. Knepley { 30147bfdf3fSMatthew G. Knepley const PetscReal k = PetscRealPart(constants[0]); 30247bfdf3fSMatthew G. Knepley PetscScalar potential; 30347bfdf3fSMatthew G. Knepley 30447bfdf3fSMatthew G. Knepley PetscCallAbort(PETSC_COMM_SELF, trig_u(dim, t, x, dim, &potential, NULL)); 30547bfdf3fSMatthew G. Knepley for (PetscInt d = 0; d < dim; ++d) f0[d * dim + d] = -k * potential * n[d]; 30647bfdf3fSMatthew G. Knepley } 30747bfdf3fSMatthew G. Knepley 30847bfdf3fSMatthew G. Knepley static void f0_sensor_bd_mixed_q(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[]) 30947bfdf3fSMatthew G. Knepley { 31047bfdf3fSMatthew G. Knepley const PetscReal k = PetscRealPart(constants[0]); 31147bfdf3fSMatthew G. Knepley PetscScalar potential; 31247bfdf3fSMatthew G. Knepley 31347bfdf3fSMatthew G. Knepley PetscCallAbort(PETSC_COMM_SELF, sensor_u(dim, t, x, dim, &potential, NULL)); 31447bfdf3fSMatthew G. Knepley for (PetscInt d = 0; d < dim; ++d) f0[d * dim + d] = -k * potential * n[d]; 31547bfdf3fSMatthew G. Knepley } 31647bfdf3fSMatthew G. Knepley 31747bfdf3fSMatthew G. Knepley /* <v, \nabla\cdot q> */ 31847bfdf3fSMatthew G. Knepley static void g1_mixed_uq(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[]) 31947bfdf3fSMatthew G. Knepley { 32047bfdf3fSMatthew G. Knepley for (PetscInt d = 0; d < dim; d++) g1[d * dim + d] = -1.0; 32147bfdf3fSMatthew G. Knepley } 32247bfdf3fSMatthew G. Knepley 32347bfdf3fSMatthew G. Knepley /* < t, q> */ 32447bfdf3fSMatthew G. Knepley static void g0_mixed_qq(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[]) 32547bfdf3fSMatthew G. Knepley { 32647bfdf3fSMatthew G. Knepley for (PetscInt d = 0; d < dim; d++) g0[d * dim + d] = 1.0; 32747bfdf3fSMatthew G. Knepley } 32847bfdf3fSMatthew G. Knepley 32947bfdf3fSMatthew G. Knepley /* <\nabla\cdot t, u> = <\nabla t, Iu> */ 33047bfdf3fSMatthew G. Knepley static void g2_mixed_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 g2[]) 33147bfdf3fSMatthew G. Knepley { 33247bfdf3fSMatthew G. Knepley const PetscReal k = PetscRealPart(constants[0]); 33347bfdf3fSMatthew G. Knepley 33447bfdf3fSMatthew G. Knepley for (PetscInt d = 0; d < dim; d++) g2[d * dim + d] = k; 33547bfdf3fSMatthew G. Knepley } 33647bfdf3fSMatthew G. Knepley 33747bfdf3fSMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *user) 33847bfdf3fSMatthew G. Knepley { 33947bfdf3fSMatthew G. Knepley PetscFunctionBeginUser; 34047bfdf3fSMatthew G. Knepley PetscOptionsBegin(comm, "", "Flux norm error in primal poisson problem Options", "DMPLEX"); 34147bfdf3fSMatthew G. Knepley user->byHand = PETSC_TRUE; 342a5cee669SMatthew G. Knepley user->numAdapt = 1; 34347bfdf3fSMatthew G. Knepley user->solType = SOL_QUADRATIC; 34447bfdf3fSMatthew G. Knepley 34547bfdf3fSMatthew G. Knepley PetscCall(PetscOptionsGetBool(NULL, NULL, "-by_hand", &user->byHand, NULL)); 346a5cee669SMatthew G. Knepley PetscCall(PetscOptionsGetInt(NULL, NULL, "-num_adapt", &user->numAdapt, NULL)); 34747bfdf3fSMatthew G. Knepley PetscCall(PetscOptionsEnum("-sol_type", "Type of exact solution", "ex27.c", SolTypeNames, (PetscEnum)user->solType, (PetscEnum *)&user->solType, NULL)); 34847bfdf3fSMatthew G. Knepley PetscOptionsEnd(); 34947bfdf3fSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 35047bfdf3fSMatthew G. Knepley } 35147bfdf3fSMatthew G. Knepley 35247bfdf3fSMatthew G. Knepley static PetscErrorCode SetupParameters(PetscBag bag, AppCtx *user) 35347bfdf3fSMatthew G. Knepley { 35447bfdf3fSMatthew G. Knepley Parameter *param; 35547bfdf3fSMatthew G. Knepley 35647bfdf3fSMatthew G. Knepley PetscFunctionBeginUser; 35747bfdf3fSMatthew G. Knepley PetscCall(PetscBagGetData(bag, (void **)¶m)); 35847bfdf3fSMatthew G. Knepley PetscCall(PetscBagSetName(bag, "par", "Poisson parameters")); 35947bfdf3fSMatthew G. Knepley PetscCall(PetscBagRegisterReal(bag, ¶m->k, 1.0, "k", "Thermal conductivity")); 36047bfdf3fSMatthew G. Knepley PetscCall(PetscBagSetFromOptions(bag)); 36147bfdf3fSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 36247bfdf3fSMatthew G. Knepley } 36347bfdf3fSMatthew G. Knepley 36447bfdf3fSMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 36547bfdf3fSMatthew G. Knepley { 36647bfdf3fSMatthew G. Knepley PetscFunctionBeginUser; 36747bfdf3fSMatthew G. Knepley PetscCall(DMCreate(comm, dm)); 36847bfdf3fSMatthew G. Knepley PetscCall(DMSetType(*dm, DMPLEX)); 36947bfdf3fSMatthew G. Knepley PetscCall(DMSetFromOptions(*dm)); 37047bfdf3fSMatthew G. Knepley PetscCall(DMSetApplicationContext(*dm, &user)); 37147bfdf3fSMatthew G. Knepley PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 37247bfdf3fSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 37347bfdf3fSMatthew G. Knepley } 37447bfdf3fSMatthew G. Knepley 37547bfdf3fSMatthew G. Knepley static PetscErrorCode SetupPrimalProblem(DM dm, AppCtx *user) 37647bfdf3fSMatthew G. Knepley { 37747bfdf3fSMatthew G. Knepley PetscDS ds; 37847bfdf3fSMatthew G. Knepley DMLabel label; 37947bfdf3fSMatthew G. Knepley PetscInt id, bd; 38047bfdf3fSMatthew G. Knepley Parameter *param; 38147bfdf3fSMatthew G. Knepley PetscWeakForm wf; 38247bfdf3fSMatthew G. Knepley 38347bfdf3fSMatthew G. Knepley PetscFunctionBeginUser; 38447bfdf3fSMatthew G. Knepley PetscCall(DMGetDS(dm, &ds)); 38547bfdf3fSMatthew G. Knepley PetscCall(DMGetLabel(dm, "marker", &label)); 38647bfdf3fSMatthew G. Knepley 38747bfdf3fSMatthew G. Knepley PetscCall(PetscDSSetJacobian(ds, 0, 0, NULL, NULL, NULL, g3_primal_uu)); 38847bfdf3fSMatthew G. Knepley 38947bfdf3fSMatthew G. Knepley switch (user->solType) { 39047bfdf3fSMatthew G. Knepley case SOL_QUADRATIC: 39147bfdf3fSMatthew G. Knepley PetscCall(PetscDSSetResidual(ds, 0, f0_quadratic_primal, f1_primal)); 39247bfdf3fSMatthew G. Knepley 39347bfdf3fSMatthew G. Knepley id = 1; 39447bfdf3fSMatthew G. Knepley PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "bottom wall primal potential", label, 1, &id, 0, 0, NULL, (void (*)(void))quadratic_u, NULL, user, NULL)); 39547bfdf3fSMatthew G. Knepley 39647bfdf3fSMatthew G. Knepley id = 2; 39747bfdf3fSMatthew G. Knepley PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "right wall flux", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd)); 39847bfdf3fSMatthew G. Knepley PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 39947bfdf3fSMatthew G. Knepley PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_quadratic_bd_primal, 0, NULL)); 40047bfdf3fSMatthew G. Knepley 40147bfdf3fSMatthew G. Knepley id = 3; 40247bfdf3fSMatthew G. Knepley PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "top wall primal potential", label, 1, &id, 0, 0, NULL, (void (*)(void))quadratic_u, NULL, user, NULL)); 40347bfdf3fSMatthew G. Knepley 40447bfdf3fSMatthew G. Knepley id = 4; 40547bfdf3fSMatthew G. Knepley PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "left wall flux", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd)); 40647bfdf3fSMatthew G. Knepley PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 40747bfdf3fSMatthew G. Knepley PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_quadratic_bd_primal, 0, NULL)); 40847bfdf3fSMatthew G. Knepley 40947bfdf3fSMatthew G. Knepley PetscCall(PetscDSSetExactSolution(ds, 0, quadratic_u, user)); 41047bfdf3fSMatthew G. Knepley break; 41147bfdf3fSMatthew G. Knepley case SOL_TRIG: 41247bfdf3fSMatthew G. Knepley PetscCall(PetscDSSetResidual(ds, 0, f0_trig_primal, f1_primal)); 41347bfdf3fSMatthew G. Knepley 41447bfdf3fSMatthew G. Knepley id = 1; 41547bfdf3fSMatthew G. Knepley PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "bottom wall primal potential", label, 1, &id, 0, 0, NULL, (void (*)(void))trig_u, NULL, user, NULL)); 41647bfdf3fSMatthew G. Knepley 41747bfdf3fSMatthew G. Knepley id = 2; 41847bfdf3fSMatthew G. Knepley PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "right wall flux", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd)); 41947bfdf3fSMatthew G. Knepley PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 42047bfdf3fSMatthew G. Knepley PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_trig_bd_primal, 0, NULL)); 42147bfdf3fSMatthew G. Knepley 42247bfdf3fSMatthew G. Knepley id = 3; 42347bfdf3fSMatthew G. Knepley PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "top wall primal potential", label, 1, &id, 0, 0, NULL, (void (*)(void))trig_u, NULL, user, NULL)); 42447bfdf3fSMatthew G. Knepley 42547bfdf3fSMatthew G. Knepley id = 4; 42647bfdf3fSMatthew G. Knepley PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "left wall flux", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd)); 42747bfdf3fSMatthew G. Knepley PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 42847bfdf3fSMatthew G. Knepley PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_trig_bd_primal, 0, NULL)); 42947bfdf3fSMatthew G. Knepley 43047bfdf3fSMatthew G. Knepley PetscCall(PetscDSSetExactSolution(ds, 0, trig_u, user)); 43147bfdf3fSMatthew G. Knepley break; 43247bfdf3fSMatthew G. Knepley case SOL_SENSOR: 43347bfdf3fSMatthew G. Knepley PetscCall(PetscDSSetResidual(ds, 0, f0_sensor_primal, f1_primal)); 43447bfdf3fSMatthew G. Knepley 43547bfdf3fSMatthew G. Knepley id = 1; 43647bfdf3fSMatthew G. Knepley PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "bottom wall primal potential", label, 1, &id, 0, 0, NULL, (void (*)(void))sensor_u, NULL, user, NULL)); 43747bfdf3fSMatthew G. Knepley 43847bfdf3fSMatthew G. Knepley id = 2; 43947bfdf3fSMatthew G. Knepley PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "right wall flux", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd)); 44047bfdf3fSMatthew G. Knepley PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 44147bfdf3fSMatthew G. Knepley PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_sensor_bd_primal, 0, NULL)); 44247bfdf3fSMatthew G. Knepley 44347bfdf3fSMatthew G. Knepley id = 3; 44447bfdf3fSMatthew G. Knepley PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "top wall primal potential", label, 1, &id, 0, 0, NULL, (void (*)(void))sensor_u, NULL, user, NULL)); 44547bfdf3fSMatthew G. Knepley 44647bfdf3fSMatthew G. Knepley id = 4; 44747bfdf3fSMatthew G. Knepley PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "left wall flux", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd)); 44847bfdf3fSMatthew G. Knepley PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 44947bfdf3fSMatthew G. Knepley PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_sensor_bd_primal, 0, NULL)); 45047bfdf3fSMatthew G. Knepley 45147bfdf3fSMatthew G. Knepley PetscCall(PetscDSSetExactSolution(ds, 0, sensor_u, user)); 45247bfdf3fSMatthew G. Knepley break; 45347bfdf3fSMatthew G. Knepley default: 45447bfdf3fSMatthew G. Knepley SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid exact solution type %s", SolTypeNames[PetscMin(user->solType, SOL_UNKNOWN)]); 45547bfdf3fSMatthew G. Knepley } 45647bfdf3fSMatthew G. Knepley 45747bfdf3fSMatthew G. Knepley /* Setup constants */ 45847bfdf3fSMatthew G. Knepley { 45947bfdf3fSMatthew G. Knepley PetscCall(PetscBagGetData(user->param, (void **)¶m)); 46047bfdf3fSMatthew G. Knepley PetscScalar constants[1]; 46147bfdf3fSMatthew G. Knepley 46247bfdf3fSMatthew G. Knepley constants[0] = param->k; 46347bfdf3fSMatthew G. Knepley 46447bfdf3fSMatthew G. Knepley PetscCall(PetscDSSetConstants(ds, 1, constants)); 46547bfdf3fSMatthew G. Knepley } 46647bfdf3fSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 46747bfdf3fSMatthew G. Knepley } 46847bfdf3fSMatthew G. Knepley 46947bfdf3fSMatthew G. Knepley static PetscErrorCode SetupPrimalDiscretization(DM dm, AppCtx *user) 47047bfdf3fSMatthew G. Knepley { 47147bfdf3fSMatthew G. Knepley DM cdm = dm; 47247bfdf3fSMatthew G. Knepley PetscFE fe[1]; 47347bfdf3fSMatthew G. Knepley DMPolytopeType ct; 47447bfdf3fSMatthew G. Knepley PetscInt dim, cStart; 47547bfdf3fSMatthew G. Knepley MPI_Comm comm; 47647bfdf3fSMatthew G. Knepley 47747bfdf3fSMatthew G. Knepley PetscFunctionBeginUser; 47847bfdf3fSMatthew G. Knepley PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 47947bfdf3fSMatthew G. Knepley PetscCall(DMGetDimension(dm, &dim)); 48047bfdf3fSMatthew G. Knepley 48147bfdf3fSMatthew G. Knepley /* Create finite element */ 48247bfdf3fSMatthew G. Knepley PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 48347bfdf3fSMatthew G. Knepley PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 48447bfdf3fSMatthew G. Knepley PetscCall(PetscFECreateByCell(comm, dim, 1, ct, "primal_potential_", PETSC_DEFAULT, &fe[0])); 48547bfdf3fSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)fe[0], "primal potential")); 48647bfdf3fSMatthew G. Knepley 48747bfdf3fSMatthew G. Knepley /* Set discretization and boundary conditions for each mesh */ 48847bfdf3fSMatthew G. Knepley PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe[0])); 48947bfdf3fSMatthew G. Knepley PetscCall(DMCreateDS(dm)); 49047bfdf3fSMatthew G. Knepley while (cdm) { 49147bfdf3fSMatthew G. Knepley PetscCall(DMCopyDisc(dm, cdm)); 49247bfdf3fSMatthew G. Knepley PetscCall(DMGetCoarseDM(cdm, &cdm)); 49347bfdf3fSMatthew G. Knepley } 49447bfdf3fSMatthew G. Knepley 49547bfdf3fSMatthew G. Knepley PetscCall(PetscFEDestroy(&fe[0])); 49647bfdf3fSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 49747bfdf3fSMatthew G. Knepley } 49847bfdf3fSMatthew G. Knepley 49947bfdf3fSMatthew G. Knepley static PetscErrorCode SetupMixedProblem(DM dm, AppCtx *user) 50047bfdf3fSMatthew G. Knepley { 50147bfdf3fSMatthew G. Knepley PetscDS ds; 50247bfdf3fSMatthew G. Knepley DMLabel label; 50347bfdf3fSMatthew G. Knepley PetscInt id, bd; 50447bfdf3fSMatthew G. Knepley Parameter *param; 50547bfdf3fSMatthew G. Knepley PetscWeakForm wf; 50647bfdf3fSMatthew G. Knepley 50747bfdf3fSMatthew G. Knepley PetscFunctionBeginUser; 50847bfdf3fSMatthew G. Knepley PetscCall(DMGetDS(dm, &ds)); 50947bfdf3fSMatthew G. Knepley PetscCall(DMGetLabel(dm, "marker", &label)); 51047bfdf3fSMatthew G. Knepley 51147bfdf3fSMatthew G. Knepley /* Residual terms */ 51247bfdf3fSMatthew G. Knepley PetscCall(PetscDSSetResidual(ds, 0, f0_mixed_q, f1_mixed_q)); 51347bfdf3fSMatthew G. Knepley 51447bfdf3fSMatthew G. Knepley PetscCall(PetscDSSetJacobian(ds, 0, 0, g0_mixed_qq, NULL, NULL, NULL)); 51547bfdf3fSMatthew G. Knepley PetscCall(PetscDSSetJacobian(ds, 0, 1, NULL, NULL, g2_mixed_qu, NULL)); 51647bfdf3fSMatthew G. Knepley 51747bfdf3fSMatthew G. Knepley PetscCall(PetscDSSetJacobian(ds, 1, 0, NULL, g1_mixed_uq, NULL, NULL)); 51847bfdf3fSMatthew G. Knepley 51947bfdf3fSMatthew G. Knepley switch (user->solType) { 52047bfdf3fSMatthew G. Knepley case SOL_QUADRATIC: 52147bfdf3fSMatthew G. Knepley PetscCall(PetscDSSetResidual(ds, 1, f0_mixed_quadratic_u, NULL)); 52247bfdf3fSMatthew G. Knepley 52347bfdf3fSMatthew G. Knepley id = 1; 52447bfdf3fSMatthew G. Knepley PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "Dirichlet Bd Integral bottom wall", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd)); 52547bfdf3fSMatthew G. Knepley PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 52647bfdf3fSMatthew G. Knepley PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_quadratic_bd_mixed_q, 0, NULL)); 52747bfdf3fSMatthew G. Knepley 52847bfdf3fSMatthew G. Knepley id = 2; 52947bfdf3fSMatthew G. Knepley PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "right wall flux", label, 1, &id, 0, 0, NULL, (void (*)(void))quadratic_q, NULL, user, NULL)); 53047bfdf3fSMatthew G. Knepley 53147bfdf3fSMatthew G. Knepley id = 4; 53247bfdf3fSMatthew G. Knepley PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "left wall flux", label, 1, &id, 0, 0, NULL, (void (*)(void))quadratic_q, NULL, user, NULL)); 53347bfdf3fSMatthew G. Knepley 53447bfdf3fSMatthew G. Knepley id = 3; 53547bfdf3fSMatthew G. Knepley PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "Dirichlet Bd Integral top wall", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd)); 53647bfdf3fSMatthew G. Knepley PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 53747bfdf3fSMatthew G. Knepley PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_quadratic_bd_mixed_q, 0, NULL)); 53847bfdf3fSMatthew G. Knepley 53947bfdf3fSMatthew G. Knepley PetscCall(PetscDSSetExactSolution(ds, 0, quadratic_q, user)); 54047bfdf3fSMatthew G. Knepley PetscCall(PetscDSSetExactSolution(ds, 1, quadratic_u, user)); 54147bfdf3fSMatthew G. Knepley break; 54247bfdf3fSMatthew G. Knepley case SOL_TRIG: 54347bfdf3fSMatthew G. Knepley PetscCall(PetscDSSetResidual(ds, 1, f0_mixed_trig_u, NULL)); 54447bfdf3fSMatthew G. Knepley 54547bfdf3fSMatthew G. Knepley id = 1; 54647bfdf3fSMatthew G. Knepley PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "Dirichlet Bd Integral bottom wall", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd)); 54747bfdf3fSMatthew G. Knepley PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 54847bfdf3fSMatthew G. Knepley PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_trig_bd_mixed_q, 0, NULL)); 54947bfdf3fSMatthew G. Knepley 55047bfdf3fSMatthew G. Knepley id = 2; 55147bfdf3fSMatthew G. Knepley PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "right wall flux", label, 1, &id, 1, 0, NULL, (void (*)(void))trig_q, NULL, user, NULL)); 55247bfdf3fSMatthew G. Knepley 55347bfdf3fSMatthew G. Knepley id = 4; 55447bfdf3fSMatthew G. Knepley PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "left wall flux", label, 1, &id, 1, 0, NULL, (void (*)(void))trig_q, NULL, user, NULL)); 55547bfdf3fSMatthew G. Knepley 55647bfdf3fSMatthew G. Knepley id = 3; 55747bfdf3fSMatthew G. Knepley PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "Dirichlet Bd Integral top wall", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd)); 55847bfdf3fSMatthew G. Knepley PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 55947bfdf3fSMatthew G. Knepley PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_trig_bd_mixed_q, 0, NULL)); 56047bfdf3fSMatthew G. Knepley 56147bfdf3fSMatthew G. Knepley PetscCall(PetscDSSetExactSolution(ds, 0, trig_q, user)); 56247bfdf3fSMatthew G. Knepley PetscCall(PetscDSSetExactSolution(ds, 1, trig_u, user)); 56347bfdf3fSMatthew G. Knepley break; 56447bfdf3fSMatthew G. Knepley case SOL_SENSOR: 56547bfdf3fSMatthew G. Knepley PetscCall(PetscDSSetResidual(ds, 1, f0_mixed_sensor_u, NULL)); 56647bfdf3fSMatthew G. Knepley 56747bfdf3fSMatthew G. Knepley id = 1; 56847bfdf3fSMatthew G. Knepley PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "Dirichlet Bd Integral bottom wall", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd)); 56947bfdf3fSMatthew G. Knepley PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 57047bfdf3fSMatthew G. Knepley PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_sensor_bd_mixed_q, 0, NULL)); 57147bfdf3fSMatthew G. Knepley 57247bfdf3fSMatthew G. Knepley id = 2; 57347bfdf3fSMatthew G. Knepley PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "right wall flux", label, 1, &id, 1, 0, NULL, (void (*)(void))sensor_q, NULL, user, NULL)); 57447bfdf3fSMatthew G. Knepley 57547bfdf3fSMatthew G. Knepley id = 4; 57647bfdf3fSMatthew G. Knepley PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "left wall flux", label, 1, &id, 1, 0, NULL, (void (*)(void))sensor_q, NULL, user, NULL)); 57747bfdf3fSMatthew G. Knepley 57847bfdf3fSMatthew G. Knepley id = 3; 57947bfdf3fSMatthew G. Knepley PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "Dirichlet Bd Integral top wall", label, 1, &id, 0, 0, NULL, NULL, NULL, user, &bd)); 58047bfdf3fSMatthew G. Knepley PetscCall(PetscDSGetBoundary(ds, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 58147bfdf3fSMatthew G. Knepley PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, id, 0, 0, 0, f0_sensor_bd_mixed_q, 0, NULL)); 58247bfdf3fSMatthew G. Knepley 58347bfdf3fSMatthew G. Knepley PetscCall(PetscDSSetExactSolution(ds, 0, sensor_q, user)); 58447bfdf3fSMatthew G. Knepley PetscCall(PetscDSSetExactSolution(ds, 1, sensor_u, user)); 58547bfdf3fSMatthew G. Knepley break; 58647bfdf3fSMatthew G. Knepley default: 58747bfdf3fSMatthew G. Knepley SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid exact solution type %s", SolTypeNames[PetscMin(user->solType, SOL_UNKNOWN)]); 58847bfdf3fSMatthew G. Knepley } 58947bfdf3fSMatthew G. Knepley 59047bfdf3fSMatthew G. Knepley /* Setup constants */ 59147bfdf3fSMatthew G. Knepley { 59247bfdf3fSMatthew G. Knepley PetscCall(PetscBagGetData(user->param, (void **)¶m)); 59347bfdf3fSMatthew G. Knepley PetscScalar constants[1]; 59447bfdf3fSMatthew G. Knepley 59547bfdf3fSMatthew G. Knepley constants[0] = param->k; 59647bfdf3fSMatthew G. Knepley 59747bfdf3fSMatthew G. Knepley PetscCall(PetscDSSetConstants(ds, 1, constants)); 59847bfdf3fSMatthew G. Knepley } 59947bfdf3fSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 60047bfdf3fSMatthew G. Knepley } 60147bfdf3fSMatthew G. Knepley 60247bfdf3fSMatthew G. Knepley static PetscErrorCode SetupMixedDiscretization(DM dm, AppCtx *user) 60347bfdf3fSMatthew G. Knepley { 60447bfdf3fSMatthew G. Knepley DM cdm = dm; 60547bfdf3fSMatthew G. Knepley PetscFE fe[2]; 60647bfdf3fSMatthew G. Knepley DMPolytopeType ct; 60747bfdf3fSMatthew G. Knepley PetscInt dim, cStart; 60847bfdf3fSMatthew G. Knepley MPI_Comm comm; 60947bfdf3fSMatthew G. Knepley 61047bfdf3fSMatthew G. Knepley PetscFunctionBeginUser; 61147bfdf3fSMatthew G. Knepley PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 61247bfdf3fSMatthew G. Knepley PetscCall(DMGetDimension(dm, &dim)); 61347bfdf3fSMatthew G. Knepley 61447bfdf3fSMatthew G. Knepley /* Create finite element */ 61547bfdf3fSMatthew G. Knepley PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL)); 61647bfdf3fSMatthew G. Knepley PetscCall(DMPlexGetCellType(dm, cStart, &ct)); 61747bfdf3fSMatthew G. Knepley PetscCall(PetscFECreateByCell(comm, dim, dim, ct, "mixed_flux_", PETSC_DEFAULT, &fe[0])); 61847bfdf3fSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)fe[0], "mixed flux")); 61947bfdf3fSMatthew G. Knepley /* NOTE: Set the same quadrature order as the primal problem here or use the command line option. */ 62047bfdf3fSMatthew G. Knepley 62147bfdf3fSMatthew G. Knepley PetscCall(PetscFECreateByCell(comm, dim, 1, ct, "mixed_potential_", PETSC_DEFAULT, &fe[1])); 62247bfdf3fSMatthew G. Knepley PetscCall(PetscFECopyQuadrature(fe[0], fe[1])); 62347bfdf3fSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)fe[1], "mixed potential")); 62447bfdf3fSMatthew G. Knepley 62547bfdf3fSMatthew G. Knepley /* Set discretization and boundary conditions for each mesh */ 62647bfdf3fSMatthew G. Knepley PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe[0])); 62747bfdf3fSMatthew G. Knepley PetscCall(DMSetField(dm, 1, NULL, (PetscObject)fe[1])); 62847bfdf3fSMatthew G. Knepley PetscCall(DMCreateDS(dm)); 62947bfdf3fSMatthew G. Knepley while (cdm) { 63047bfdf3fSMatthew G. Knepley PetscCall(DMCopyDisc(dm, cdm)); 63147bfdf3fSMatthew G. Knepley PetscCall(DMGetCoarseDM(cdm, &cdm)); 63247bfdf3fSMatthew G. Knepley } 63347bfdf3fSMatthew G. Knepley PetscCall(PetscFEDestroy(&fe[0])); 63447bfdf3fSMatthew G. Knepley PetscCall(PetscFEDestroy(&fe[1])); 63547bfdf3fSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 63647bfdf3fSMatthew G. Knepley } 63747bfdf3fSMatthew G. Knepley 63847bfdf3fSMatthew G. Knepley PetscErrorCode SetupMixed(DMAdaptor adaptor, DM mdm) 63947bfdf3fSMatthew G. Knepley { 64047bfdf3fSMatthew G. Knepley AppCtx *ctx; 64147bfdf3fSMatthew G. Knepley 64247bfdf3fSMatthew G. Knepley PetscFunctionBeginUser; 64347bfdf3fSMatthew G. Knepley PetscCall(DMGetApplicationContext(mdm, (void **)&ctx)); 64447bfdf3fSMatthew G. Knepley PetscCall(SetupMixedDiscretization(mdm, ctx)); 64547bfdf3fSMatthew G. Knepley PetscCall(SetupMixedProblem(mdm, ctx)); 64647bfdf3fSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 64747bfdf3fSMatthew G. Knepley } 64847bfdf3fSMatthew G. Knepley 64947bfdf3fSMatthew G. Knepley int main(int argc, char **argv) 65047bfdf3fSMatthew G. Knepley { 65147bfdf3fSMatthew G. Knepley DM dm, mdm; /* problem specification */ 65247bfdf3fSMatthew G. Knepley SNES snes, msnes; /* nonlinear solvers */ 65347bfdf3fSMatthew G. Knepley Vec u, mu; /* solution vectors */ 65447bfdf3fSMatthew G. Knepley Vec fluxError, exactError; /* Element wise error vector */ 65547bfdf3fSMatthew G. Knepley PetscReal fluxNorm, exactNorm; /* Flux error norm */ 65647bfdf3fSMatthew G. Knepley AppCtx user; /* user-defined work context */ 65747bfdf3fSMatthew G. Knepley 65847bfdf3fSMatthew G. Knepley PetscFunctionBeginUser; 65947bfdf3fSMatthew G. Knepley PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 66047bfdf3fSMatthew G. Knepley PetscCall(PetscBagCreate(PETSC_COMM_WORLD, sizeof(Parameter), &user.param)); 66147bfdf3fSMatthew G. Knepley PetscCall(SetupParameters(user.param, &user)); 66247bfdf3fSMatthew G. Knepley PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 66347bfdf3fSMatthew G. Knepley 66447bfdf3fSMatthew G. Knepley // Set up and solve primal problem 66547bfdf3fSMatthew G. Knepley PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 66647bfdf3fSMatthew G. Knepley PetscCall(DMSetApplicationContext(dm, &user)); 66747bfdf3fSMatthew G. Knepley PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 66847bfdf3fSMatthew G. Knepley PetscCall(SNESSetDM(snes, dm)); 66947bfdf3fSMatthew G. Knepley 67047bfdf3fSMatthew G. Knepley PetscCall(SetupPrimalDiscretization(dm, &user)); 67147bfdf3fSMatthew G. Knepley PetscCall(SetupPrimalProblem(dm, &user)); 67247bfdf3fSMatthew G. Knepley PetscCall(DMCreateGlobalVector(dm, &u)); 67347bfdf3fSMatthew G. Knepley PetscCall(VecSet(u, 0.0)); 67447bfdf3fSMatthew G. Knepley PetscCall(DMPlexSetSNESLocalFEM(dm, PETSC_FALSE, &user)); 67547bfdf3fSMatthew G. Knepley PetscCall(SNESSetFromOptions(snes)); 67647bfdf3fSMatthew G. Knepley PetscCall(DMSNESCheckFromOptions(snes, u)); 677a5cee669SMatthew G. Knepley 678a5cee669SMatthew G. Knepley for (PetscInt a = 0; a < user.numAdapt; ++a) { 679a5cee669SMatthew G. Knepley if (a > 0) { 680a5cee669SMatthew G. Knepley char prefix[16]; 681a5cee669SMatthew G. Knepley 682a5cee669SMatthew G. Knepley PetscCall(PetscSNPrintf(prefix, 16, "a%d_", (int)a)); 683a5cee669SMatthew G. Knepley PetscCall(SNESSetOptionsPrefix(snes, prefix)); 684a5cee669SMatthew G. Knepley } 68547bfdf3fSMatthew G. Knepley PetscCall(SNESSolve(snes, NULL, u)); 68647bfdf3fSMatthew G. Knepley 68747bfdf3fSMatthew G. Knepley // Needed if you allow SNES to refine 68847bfdf3fSMatthew G. Knepley PetscCall(SNESGetSolution(snes, &u)); 68947bfdf3fSMatthew G. Knepley PetscCall(VecGetDM(u, &dm)); 690a5cee669SMatthew G. Knepley } 69147bfdf3fSMatthew G. Knepley 69247bfdf3fSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)u, "Primal Solution ")); 69347bfdf3fSMatthew G. Knepley PetscCall(VecViewFromOptions(u, NULL, "-primal_sol_vec_view")); 69447bfdf3fSMatthew G. Knepley 69547bfdf3fSMatthew G. Knepley if (user.byHand) { 69647bfdf3fSMatthew G. Knepley // Set up and solve mixed problem 69747bfdf3fSMatthew G. Knepley PetscCall(DMClone(dm, &mdm)); 69847bfdf3fSMatthew G. Knepley PetscCall(SNESCreate(PETSC_COMM_WORLD, &msnes)); 69947bfdf3fSMatthew G. Knepley PetscCall(SNESSetOptionsPrefix(msnes, "mixed_")); 70047bfdf3fSMatthew G. Knepley PetscCall(SNESSetDM(msnes, mdm)); 70147bfdf3fSMatthew G. Knepley 70247bfdf3fSMatthew G. Knepley PetscCall(SetupMixedDiscretization(mdm, &user)); 70347bfdf3fSMatthew G. Knepley PetscCall(SetupMixedProblem(mdm, &user)); 70447bfdf3fSMatthew G. Knepley PetscCall(DMCreateGlobalVector(mdm, &mu)); 70547bfdf3fSMatthew G. Knepley PetscCall(VecSet(mu, 0.0)); 70647bfdf3fSMatthew G. Knepley PetscCall(DMPlexSetSNESLocalFEM(mdm, PETSC_FALSE, &user)); 70747bfdf3fSMatthew G. Knepley PetscCall(SNESSetFromOptions(msnes)); 70847bfdf3fSMatthew G. Knepley 70947bfdf3fSMatthew G. Knepley PetscCall(DMSNESCheckFromOptions(msnes, mu)); 71047bfdf3fSMatthew G. Knepley PetscCall(SNESSolve(msnes, NULL, mu)); 71147bfdf3fSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)mu, "Mixed Solution ")); 71247bfdf3fSMatthew G. Knepley PetscCall(VecViewFromOptions(mu, NULL, "-mixed_sol_vec_view")); 71347bfdf3fSMatthew G. Knepley 71447bfdf3fSMatthew G. Knepley // Create the error space of piecewise constants 71547bfdf3fSMatthew G. Knepley DM dmErr; 71647bfdf3fSMatthew G. Knepley PetscFE feErr; 71747bfdf3fSMatthew G. Knepley DMPolytopeType ct; 71847bfdf3fSMatthew G. Knepley PetscInt dim, cStart; 71947bfdf3fSMatthew G. Knepley 72047bfdf3fSMatthew G. Knepley PetscCall(DMClone(dm, &dmErr)); 72147bfdf3fSMatthew G. Knepley PetscCall(DMGetDimension(dmErr, &dim)); 72247bfdf3fSMatthew G. Knepley PetscCall(DMPlexGetHeightStratum(dmErr, 0, &cStart, NULL)); 72347bfdf3fSMatthew G. Knepley PetscCall(DMPlexGetCellType(dmErr, cStart, &ct)); 72447bfdf3fSMatthew G. Knepley PetscCall(PetscFECreateLagrangeByCell(PETSC_COMM_SELF, dim, 1, ct, 0, PETSC_DEFAULT, &feErr)); 72547bfdf3fSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)feErr, "Error")); 72647bfdf3fSMatthew G. Knepley PetscCall(DMSetField(dmErr, 0, NULL, (PetscObject)feErr)); 72747bfdf3fSMatthew G. Knepley PetscCall(PetscFEDestroy(&feErr)); 72847bfdf3fSMatthew G. Knepley PetscCall(DMCreateDS(dmErr)); 72947bfdf3fSMatthew G. Knepley PetscCall(DMViewFromOptions(dmErr, NULL, "-dmerr_view")); 73047bfdf3fSMatthew G. Knepley 73147bfdf3fSMatthew G. Knepley // Compute the flux norm 73247bfdf3fSMatthew G. Knepley PetscCall(DMGetGlobalVector(dmErr, &fluxError)); 73347bfdf3fSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)fluxError, "Flux Error")); 73447bfdf3fSMatthew G. Knepley PetscCall(DMGetGlobalVector(dmErr, &exactError)); 73547bfdf3fSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)exactError, "Analytical Error")); 73647bfdf3fSMatthew G. Knepley PetscCall(DMPlexComputeL2FluxDiffVec(u, 0, mu, 0, fluxError)); 73747bfdf3fSMatthew G. Knepley { 73847bfdf3fSMatthew G. Knepley PetscDS ds; 73947bfdf3fSMatthew G. Knepley PetscSimplePointFn *func[2] = {NULL, NULL}; 74047bfdf3fSMatthew G. Knepley void *ctx[2] = {NULL, NULL}; 74147bfdf3fSMatthew G. Knepley 74247bfdf3fSMatthew G. Knepley PetscCall(DMGetDS(mdm, &ds)); 74347bfdf3fSMatthew G. Knepley PetscCall(PetscDSGetExactSolution(ds, 0, &func[0], &ctx[0])); 74447bfdf3fSMatthew G. Knepley PetscCall(DMPlexComputeL2DiffVec(mdm, 0.0, func, ctx, mu, exactError)); 74547bfdf3fSMatthew G. Knepley } 74647bfdf3fSMatthew G. Knepley PetscCall(VecNorm(fluxError, NORM_2, &fluxNorm)); 74747bfdf3fSMatthew G. Knepley PetscCall(VecNorm(exactError, NORM_2, &exactNorm)); 74847bfdf3fSMatthew G. Knepley PetscCall(PetscPrintf(PETSC_COMM_SELF, "Flux error norm = %g\t Exact flux error norm = %g\n", (double)fluxNorm, (double)exactNorm)); 74947bfdf3fSMatthew G. Knepley PetscCall(VecViewFromOptions(fluxError, NULL, "-flux_error_vec_view")); 75047bfdf3fSMatthew G. Knepley PetscCall(VecViewFromOptions(exactError, NULL, "-exact_error_vec_view")); 75147bfdf3fSMatthew G. Knepley 75247bfdf3fSMatthew G. Knepley // Adaptive refinement based on calculated error 75347bfdf3fSMatthew G. Knepley DM rdm; 75447bfdf3fSMatthew G. Knepley VecTagger refineTag; 75547bfdf3fSMatthew G. Knepley DMLabel adaptLabel; 75647bfdf3fSMatthew G. Knepley IS refineIS; 75747bfdf3fSMatthew G. Knepley Vec ref; 75847bfdf3fSMatthew G. Knepley 75947bfdf3fSMatthew G. Knepley PetscCall(DMLabelCreate(PETSC_COMM_WORLD, "adapt", &adaptLabel)); 76047bfdf3fSMatthew G. Knepley PetscCall(DMLabelSetDefaultValue(adaptLabel, DM_ADAPT_COARSEN)); 76147bfdf3fSMatthew G. Knepley PetscCall(VecTaggerCreate(PETSC_COMM_WORLD, &refineTag)); 76247bfdf3fSMatthew G. Knepley PetscCall(VecTaggerSetFromOptions(refineTag)); 76347bfdf3fSMatthew G. Knepley PetscCall(VecTaggerSetUp(refineTag)); 76447bfdf3fSMatthew G. Knepley PetscCall(PetscObjectViewFromOptions((PetscObject)refineTag, NULL, "-tag_view")); 76547bfdf3fSMatthew G. Knepley 76647bfdf3fSMatthew G. Knepley PetscCall(VecTaggerComputeIS(refineTag, fluxError, &refineIS, NULL)); 76747bfdf3fSMatthew G. Knepley PetscCall(VecTaggerDestroy(&refineTag)); 76847bfdf3fSMatthew G. Knepley PetscCall(DMLabelSetStratumIS(adaptLabel, DM_ADAPT_REFINE, refineIS)); 76947bfdf3fSMatthew G. Knepley PetscCall(ISViewFromOptions(refineIS, NULL, "-refine_is_view")); 77047bfdf3fSMatthew G. Knepley PetscCall(ISDestroy(&refineIS)); 77147bfdf3fSMatthew G. Knepley 77247bfdf3fSMatthew G. Knepley PetscCall(DMPlexCreateLabelField(dm, adaptLabel, &ref)); 77347bfdf3fSMatthew G. Knepley PetscCall(VecViewFromOptions(ref, NULL, "-refine_vec_view")); 77447bfdf3fSMatthew G. Knepley PetscCall(VecDestroy(&ref)); 77547bfdf3fSMatthew G. Knepley 77647bfdf3fSMatthew G. Knepley // Mark adaptation phase with prefix ref_ 77747bfdf3fSMatthew G. Knepley PetscCall(PetscObjectSetOptionsPrefix((PetscObject)dm, "adapt_")); 77847bfdf3fSMatthew G. Knepley PetscCall(DMAdaptLabel(dm, adaptLabel, &rdm)); 77947bfdf3fSMatthew G. Knepley PetscCall(PetscObjectSetOptionsPrefix((PetscObject)dm, NULL)); 78047bfdf3fSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)rdm, "Adaptively Refined DM")); 78147bfdf3fSMatthew G. Knepley PetscCall(DMViewFromOptions(rdm, NULL, "-adapt_dm_view")); 78247bfdf3fSMatthew G. Knepley PetscCall(DMDestroy(&rdm)); 78347bfdf3fSMatthew G. Knepley PetscCall(DMLabelDestroy(&adaptLabel)); 78447bfdf3fSMatthew G. Knepley 78547bfdf3fSMatthew G. Knepley // Destroy the error structures 78647bfdf3fSMatthew G. Knepley PetscCall(DMRestoreGlobalVector(dmErr, &fluxError)); 78747bfdf3fSMatthew G. Knepley PetscCall(DMRestoreGlobalVector(dmErr, &exactError)); 78847bfdf3fSMatthew G. Knepley PetscCall(DMDestroy(&dmErr)); 78947bfdf3fSMatthew G. Knepley 79047bfdf3fSMatthew G. Knepley // Destroy the mixed structures 79147bfdf3fSMatthew G. Knepley PetscCall(VecDestroy(&mu)); 79247bfdf3fSMatthew G. Knepley PetscCall(DMDestroy(&mdm)); 79347bfdf3fSMatthew G. Knepley PetscCall(SNESDestroy(&msnes)); 79447bfdf3fSMatthew G. Knepley } 79547bfdf3fSMatthew G. Knepley 79647bfdf3fSMatthew G. Knepley // Destroy the primal structures 79747bfdf3fSMatthew G. Knepley PetscCall(VecDestroy(&u)); 79847bfdf3fSMatthew G. Knepley PetscCall(DMDestroy(&dm)); 79947bfdf3fSMatthew G. Knepley PetscCall(SNESDestroy(&snes)); 80047bfdf3fSMatthew G. Knepley PetscCall(PetscBagDestroy(&user.param)); 80147bfdf3fSMatthew G. Knepley PetscCall(PetscFinalize()); 80247bfdf3fSMatthew G. Knepley return 0; 80347bfdf3fSMatthew G. Knepley } 80447bfdf3fSMatthew G. Knepley 80547bfdf3fSMatthew G. Knepley /*TEST 80647bfdf3fSMatthew G. Knepley 80747bfdf3fSMatthew G. Knepley # Tests using the explicit code above 80847bfdf3fSMatthew G. Knepley testset: 80947bfdf3fSMatthew G. Knepley suffix: 2d_p2_rt0p0_byhand 81047bfdf3fSMatthew G. Knepley requires: triangle 81147bfdf3fSMatthew G. Knepley args: -dm_plex_box_faces 1,1 -dm_plex_box_lower 0,0 -dm_plex_box_upper 1,1 -dm_refine 3 \ 81247bfdf3fSMatthew G. Knepley -primal_potential_petscspace_degree 2 \ 81347bfdf3fSMatthew G. Knepley -mixed_potential_petscdualspace_lagrange_continuity 0 \ 81447bfdf3fSMatthew G. Knepley -mixed_flux_petscspace_type ptrimmed \ 81547bfdf3fSMatthew G. Knepley -mixed_flux_petscspace_components 2 \ 81647bfdf3fSMatthew G. Knepley -mixed_flux_petscspace_ptrimmed_form_degree -1 \ 81747bfdf3fSMatthew G. Knepley -mixed_flux_petscdualspace_order 1 \ 81847bfdf3fSMatthew G. Knepley -mixed_flux_petscdualspace_form_degree -1 \ 81947bfdf3fSMatthew G. Knepley -mixed_flux_petscdualspace_lagrange_trimmed true \ 82047bfdf3fSMatthew G. Knepley -mixed_flux_petscfe_default_quadrature_order 2 \ 82147bfdf3fSMatthew G. Knepley -vec_tagger_type cdf -vec_tagger_box 0.9,1.0 \ 82247bfdf3fSMatthew G. Knepley -tag_view \ 82347bfdf3fSMatthew G. Knepley -adapt_dm_adaptor cellrefiner -adapt_dm_plex_transform_type refine_sbr \ 82447bfdf3fSMatthew G. Knepley -dmsnes_check 0.001 -mixed_dmsnes_check 0.001 -pc_type jacobi -mixed_pc_type jacobi 82547bfdf3fSMatthew G. Knepley test: 82647bfdf3fSMatthew G. Knepley suffix: quadratic 82747bfdf3fSMatthew G. Knepley args: -sol_type quadratic 82847bfdf3fSMatthew G. Knepley test: 82947bfdf3fSMatthew G. Knepley suffix: trig 83047bfdf3fSMatthew G. Knepley args: -sol_type trig 83147bfdf3fSMatthew G. Knepley test: 83247bfdf3fSMatthew G. Knepley suffix: sensor 83347bfdf3fSMatthew G. Knepley args: -sol_type sensor 83447bfdf3fSMatthew G. Knepley 83547bfdf3fSMatthew G. Knepley # Tests using the embedded adaptor in SNES 83647bfdf3fSMatthew G. Knepley testset: 83747bfdf3fSMatthew G. Knepley suffix: 2d_p2_rt0p0 83847bfdf3fSMatthew G. Knepley requires: triangle defined(PETSC_HAVE_EXECUTABLE_EXPORT) 83947bfdf3fSMatthew G. Knepley args: -dm_plex_box_faces 1,1 -dm_plex_box_lower 0,0 -dm_plex_box_upper 1,1 -dm_refine 3 \ 84047bfdf3fSMatthew G. Knepley -primal_potential_petscspace_degree 2 \ 84147bfdf3fSMatthew G. Knepley -mixed_potential_petscdualspace_lagrange_continuity 0 \ 84247bfdf3fSMatthew G. Knepley -mixed_flux_petscspace_type ptrimmed \ 84347bfdf3fSMatthew G. Knepley -mixed_flux_petscspace_components 2 \ 84447bfdf3fSMatthew G. Knepley -mixed_flux_petscspace_ptrimmed_form_degree -1 \ 84547bfdf3fSMatthew G. Knepley -mixed_flux_petscdualspace_order 1 \ 84647bfdf3fSMatthew G. Knepley -mixed_flux_petscdualspace_form_degree -1 \ 84747bfdf3fSMatthew G. Knepley -mixed_flux_petscdualspace_lagrange_trimmed true \ 84847bfdf3fSMatthew G. Knepley -mixed_flux_petscfe_default_quadrature_order 2 \ 84947bfdf3fSMatthew G. Knepley -by_hand 0 \ 85047bfdf3fSMatthew G. Knepley -refine_vec_tagger_type cdf -refine_vec_tagger_box 0.9,1.0 \ 85147bfdf3fSMatthew G. Knepley -snes_adapt_view \ 85247bfdf3fSMatthew G. Knepley -adapt_dm_adaptor cellrefiner -adapt_dm_plex_transform_type refine_sbr \ 85347bfdf3fSMatthew G. Knepley -adaptor_criterion label -adaptor_type flux -adaptor_mixed_setup_function SetupMixed \ 85447bfdf3fSMatthew G. Knepley -snes_adapt_sequence 1 -pc_type jacobi -mixed_pc_type jacobi 85547bfdf3fSMatthew G. Knepley test: 85647bfdf3fSMatthew G. Knepley suffix: quadratic 85747bfdf3fSMatthew G. Knepley args: -sol_type quadratic -adaptor_monitor_error 85847bfdf3fSMatthew G. Knepley test: 85947bfdf3fSMatthew G. Knepley suffix: trig 86047bfdf3fSMatthew G. Knepley args: -sol_type trig -adaptor_monitor_error 86147bfdf3fSMatthew G. Knepley test: 86247bfdf3fSMatthew G. Knepley suffix: sensor 86347bfdf3fSMatthew G. Knepley args: -sol_type sensor 86447bfdf3fSMatthew G. Knepley 865a5cee669SMatthew G. Knepley # Tests using multiple adaptor loops 866a5cee669SMatthew G. Knepley testset: 867a5cee669SMatthew G. Knepley suffix: 2d_p2_rt0p0_a2 868a5cee669SMatthew G. Knepley requires: triangle defined(PETSC_HAVE_EXECUTABLE_EXPORT) 869a5cee669SMatthew G. Knepley args: -dm_plex_box_faces 1,1 -dm_plex_box_lower 0,0 -dm_plex_box_upper 1,1 -dm_refine 3 \ 870a5cee669SMatthew G. Knepley -primal_potential_petscspace_degree 2 \ 871a5cee669SMatthew G. Knepley -mixed_potential_petscdualspace_lagrange_continuity 0 \ 872a5cee669SMatthew G. Knepley -mixed_flux_petscspace_type ptrimmed \ 873a5cee669SMatthew G. Knepley -mixed_flux_petscspace_components 2 \ 874a5cee669SMatthew G. Knepley -mixed_flux_petscspace_ptrimmed_form_degree -1 \ 875a5cee669SMatthew G. Knepley -mixed_flux_petscdualspace_order 1 \ 876a5cee669SMatthew G. Knepley -mixed_flux_petscdualspace_form_degree -1 \ 877a5cee669SMatthew G. Knepley -mixed_flux_petscdualspace_lagrange_trimmed true \ 878a5cee669SMatthew G. Knepley -mixed_flux_petscfe_default_quadrature_order 2 \ 879a5cee669SMatthew G. Knepley -by_hand 0 \ 880a5cee669SMatthew G. Knepley -num_adapt 2 \ 881a5cee669SMatthew G. Knepley -refine_vec_tagger_type cdf -refine_vec_tagger_box 0.9,1.0 \ 882a5cee669SMatthew G. Knepley -snes_adapt_view \ 883a5cee669SMatthew G. Knepley -adapt_dm_adaptor cellrefiner -adapt_dm_plex_transform_type refine_sbr \ 884a5cee669SMatthew G. Knepley -adaptor_criterion label -adaptor_type gradient -adaptor_mixed_setup_function SetupMixed \ 885a5cee669SMatthew G. Knepley -snes_adapt_sequence 2 -pc_type jacobi \ 886a5cee669SMatthew G. Knepley -a1_refine_vec_tagger_type cdf -a1_refine_vec_tagger_box 0.9,1.0 \ 887a5cee669SMatthew G. Knepley -a1_snes_adapt_view \ 888a5cee669SMatthew G. Knepley -a1_adaptor_criterion label -a1_adaptor_type flux -a1_adaptor_mixed_setup_function SetupMixed \ 889a5cee669SMatthew G. Knepley -a1_snes_adapt_sequence 1 -a1_pc_type jacobi -a1_mixed_pc_type jacobi 890a5cee669SMatthew G. Knepley test: 891a5cee669SMatthew G. Knepley suffix: sensor 892a5cee669SMatthew G. Knepley args: -sol_type sensor 893a5cee669SMatthew G. Knepley 89447bfdf3fSMatthew G. Knepley TEST*/ 895