1557cf195SMatthew G. Knepley static char help[] = "Test adaptive interpolation of functions of a given polynomial order\n\n"; 2557cf195SMatthew G. Knepley 3557cf195SMatthew G. Knepley #include <petscdmplex.h> 4557cf195SMatthew G. Knepley #include <petscsnes.h> 5557cf195SMatthew G. Knepley 6557cf195SMatthew G. Knepley /* 7557cf195SMatthew G. Knepley What properties does the adapted interpolator have? 8557cf195SMatthew G. Knepley 9557cf195SMatthew G. Knepley 1) If we adapt to quadratics, we can get lower interpolation error for quadratics (than local interpolation) when using a linear basis 10557cf195SMatthew G. Knepley 11557cf195SMatthew G. Knepley $ ./ex8 -dm_refine 2 -petscspace_degree 1 -qorder 1 -dim 2 -porder 2 -K 2 -num_comp 1 -use_poly 1 12557cf195SMatthew G. Knepley Function tests FAIL for order 2 at tolerance 1e-10 error 0.00273757 13557cf195SMatthew G. Knepley Function tests FAIL for order 2 derivatives at tolerance 1e-10 error 0.0721688 14557cf195SMatthew G. Knepley Interpolation tests FAIL for order 2 at tolerance 1e-10 error 0.00284555 15557cf195SMatthew G. Knepley Interpolation tests FAIL for order 2 derivatives at tolerance 1e-10 error 0.0721688 16557cf195SMatthew G. Knepley Adapting interpolator using polynomials 17557cf195SMatthew G. Knepley The number of input vectors 4 < 7 the maximum number of column entries 18557cf195SMatthew G. Knepley Interpolation poly tests FAIL for order 2 at tolerance 1e-10 error 0.00659864 19557cf195SMatthew G. Knepley Interpolation poly tests FAIL for order 2 derivatives at tolerance 1e-10 error 0.0836582 20557cf195SMatthew G. Knepley Interpolation trig (0, 0) tests FAIL for order 1 at tolerance 1e-10 error 0.476194 21557cf195SMatthew G. Knepley Interpolation trig (0, 0) tests FAIL for order 1 derivatives at tolerance 1e-10 error 2.22144 22557cf195SMatthew G. Knepley Interpolation trig (0, 1) tests FAIL for order 1 at tolerance 1e-10 error 1.39768 23557cf195SMatthew G. Knepley Interpolation trig (0, 1) tests FAIL for order 1 derivatives at tolerance 1e-10 error 2.22144 24557cf195SMatthew G. Knepley Interpolation trig (1, 0) tests FAIL for order 2 at tolerance 1e-10 error 1.07315 25557cf195SMatthew G. Knepley Interpolation trig (1, 0) tests FAIL for order 2 derivatives at tolerance 1e-10 error 4.55403 26557cf195SMatthew G. Knepley Interpolation trig (1, 1) tests FAIL for order 2 at tolerance 1e-10 error 1.07315 27557cf195SMatthew G. Knepley Interpolation trig (1, 1) tests FAIL for order 2 derivatives at tolerance 1e-10 error 4.55403 28557cf195SMatthew G. Knepley 29557cf195SMatthew G. Knepley $ ./ex8 -dm_refine 2 -petscspace_degree 1 -qorder 1 -dim 2 -porder 2 -K 3 -num_comp 1 -use_poly 1 30557cf195SMatthew G. Knepley Function tests FAIL for order 2 at tolerance 1e-10 error 0.00273757 31557cf195SMatthew G. Knepley Function tests FAIL for order 2 derivatives at tolerance 1e-10 error 0.0721688 32557cf195SMatthew G. Knepley Interpolation tests FAIL for order 2 at tolerance 1e-10 error 0.00284555 33557cf195SMatthew G. Knepley Interpolation tests FAIL for order 2 derivatives at tolerance 1e-10 error 0.0721688 34557cf195SMatthew G. Knepley Adapting interpolator using polynomials 35557cf195SMatthew G. Knepley The number of input vectors 6 < 7 the maximum number of column entries 36557cf195SMatthew G. Knepley Interpolation poly tests FAIL for order 2 at tolerance 1e-10 error 0.00194055 37557cf195SMatthew G. Knepley Interpolation poly tests FAIL for order 2 derivatives at tolerance 1e-10 error 0.0525591 38557cf195SMatthew G. Knepley Interpolation trig (0, 0) tests FAIL for order 1 at tolerance 1e-10 error 0.476255 39557cf195SMatthew G. Knepley Interpolation trig (0, 0) tests FAIL for order 1 derivatives at tolerance 1e-10 error 2.22132 40557cf195SMatthew G. Knepley Interpolation trig (0, 1) tests FAIL for order 1 at tolerance 1e-10 error 1.39785 41557cf195SMatthew G. Knepley Interpolation trig (0, 1) tests FAIL for order 1 derivatives at tolerance 1e-10 error 2.22119 42557cf195SMatthew G. Knepley Interpolation trig (1, 0) tests FAIL for order 2 at tolerance 1e-10 error 1.0727 43557cf195SMatthew G. Knepley Interpolation trig (1, 0) tests FAIL for order 2 derivatives at tolerance 1e-10 error 4.55364 44557cf195SMatthew G. Knepley Interpolation trig (1, 1) tests FAIL for order 2 at tolerance 1e-10 error 1.0727 45557cf195SMatthew G. Knepley Interpolation trig (1, 1) tests FAIL for order 2 derivatives at tolerance 1e-10 error 4.55364 46557cf195SMatthew G. Knepley Interpolation trig (2, 0) tests FAIL for order 3 at tolerance 1e-10 error 0.705258 47557cf195SMatthew G. Knepley Interpolation trig (2, 0) tests FAIL for order 3 derivatives at tolerance 1e-10 error 6.82037 48557cf195SMatthew G. Knepley Interpolation trig (2, 1) tests FAIL for order 3 at tolerance 1e-10 error 0.705258 49557cf195SMatthew G. Knepley Interpolation trig (2, 1) tests FAIL for order 3 derivatives at tolerance 1e-10 error 6.82037 50557cf195SMatthew G. Knepley 51557cf195SMatthew G. Knepley 2) We can more accurately capture low harmonics 52557cf195SMatthew G. Knepley 53557cf195SMatthew G. Knepley If we adapt polynomials, we can be exact 54557cf195SMatthew G. Knepley 55557cf195SMatthew G. Knepley $ ./ex8 -dm_refine 2 -petscspace_degree 1 -qorder 1 -dim 2 -porder 1 -K 2 -num_comp 1 -use_poly 1 56557cf195SMatthew G. Knepley Function tests pass for order 1 at tolerance 1e-10 57557cf195SMatthew G. Knepley Function tests pass for order 1 derivatives at tolerance 1e-10 58557cf195SMatthew G. Knepley Interpolation tests pass for order 1 at tolerance 1e-10 59557cf195SMatthew G. Knepley Interpolation tests pass for order 1 derivatives at tolerance 1e-10 60557cf195SMatthew G. Knepley Adapting interpolator using polynomials 61557cf195SMatthew G. Knepley The number of input vectors 4 < 7 the maximum number of column entries 62557cf195SMatthew G. Knepley Interpolation poly tests pass for order 1 at tolerance 1e-10 63557cf195SMatthew G. Knepley Interpolation poly tests pass for order 1 derivatives at tolerance 1e-10 64557cf195SMatthew G. Knepley Interpolation trig (0, 0) tests FAIL for order 1 at tolerance 1e-10 error 0.476194 65557cf195SMatthew G. Knepley Interpolation trig (0, 0) tests FAIL for order 1 derivatives at tolerance 1e-10 error 2.22144 66557cf195SMatthew G. Knepley Interpolation trig (0, 1) tests FAIL for order 1 at tolerance 1e-10 error 1.39768 67557cf195SMatthew G. Knepley Interpolation trig (0, 1) tests FAIL for order 1 derivatives at tolerance 1e-10 error 2.22144 68557cf195SMatthew G. Knepley Interpolation trig (1, 0) tests FAIL for order 2 at tolerance 1e-10 error 1.07315 69557cf195SMatthew G. Knepley Interpolation trig (1, 0) tests FAIL for order 2 derivatives at tolerance 1e-10 error 4.55403 70557cf195SMatthew G. Knepley Interpolation trig (1, 1) tests FAIL for order 2 at tolerance 1e-10 error 1.07315 71557cf195SMatthew G. Knepley Interpolation trig (1, 1) tests FAIL for order 2 derivatives at tolerance 1e-10 error 4.55403 72557cf195SMatthew G. Knepley 73557cf195SMatthew G. Knepley and least for small K, 74557cf195SMatthew G. Knepley 75557cf195SMatthew G. Knepley $ ./ex8 -dm_refine 2 -petscspace_degree 1 -qorder 1 -dim 2 -porder 1 -K 4 -num_comp 1 -use_poly 1 76557cf195SMatthew G. Knepley Function tests pass for order 1 at tolerance 1e-10 77557cf195SMatthew G. Knepley Function tests pass for order 1 derivatives at tolerance 1e-10 78557cf195SMatthew G. Knepley Interpolation tests pass for order 1 at tolerance 1e-10 79557cf195SMatthew G. Knepley Interpolation tests pass for order 1 derivatives at tolerance 1e-10 80557cf195SMatthew G. Knepley Adapting interpolator using polynomials 81557cf195SMatthew G. Knepley Interpolation poly tests FAIL for order 1 at tolerance 1e-10 error 0.0015351 82557cf195SMatthew G. Knepley Interpolation poly tests FAIL for order 1 derivatives at tolerance 1e-10 error 0.0427369 83557cf195SMatthew G. Knepley Interpolation trig (0, 0) tests FAIL for order 1 at tolerance 1e-10 error 0.476359 84557cf195SMatthew G. Knepley Interpolation trig (0, 0) tests FAIL for order 1 derivatives at tolerance 1e-10 error 2.22115 85557cf195SMatthew G. Knepley Interpolation trig (0, 1) tests FAIL for order 1 at tolerance 1e-10 error 1.3981 86557cf195SMatthew G. Knepley Interpolation trig (0, 1) tests FAIL for order 1 derivatives at tolerance 1e-10 error 2.22087 87557cf195SMatthew G. Knepley Interpolation trig (1, 0) tests FAIL for order 2 at tolerance 1e-10 error 1.07228 88557cf195SMatthew G. Knepley Interpolation trig (1, 0) tests FAIL for order 2 derivatives at tolerance 1e-10 error 4.55238 89557cf195SMatthew G. Knepley Interpolation trig (1, 1) tests FAIL for order 2 at tolerance 1e-10 error 1.07228 90557cf195SMatthew G. Knepley Interpolation trig (1, 1) tests FAIL for order 2 derivatives at tolerance 1e-10 error 4.55238 91557cf195SMatthew G. Knepley Interpolation trig (2, 0) tests FAIL for order 3 at tolerance 1e-10 error 0.704947 92557cf195SMatthew G. Knepley Interpolation trig (2, 0) tests FAIL for order 3 derivatives at tolerance 1e-10 error 6.82254 93557cf195SMatthew G. Knepley Interpolation trig (2, 1) tests FAIL for order 3 at tolerance 1e-10 error 0.704948 94557cf195SMatthew G. Knepley Interpolation trig (2, 1) tests FAIL for order 3 derivatives at tolerance 1e-10 error 6.82254 95557cf195SMatthew G. Knepley Interpolation trig (3, 0) tests FAIL for order 4 at tolerance 1e-10 error 0.893279 96557cf195SMatthew G. Knepley Interpolation trig (3, 0) tests FAIL for order 4 derivatives at tolerance 1e-10 error 8.93718 97557cf195SMatthew G. Knepley Interpolation trig (3, 1) tests FAIL for order 4 at tolerance 1e-10 error 0.89328 98557cf195SMatthew G. Knepley Interpolation trig (3, 1) tests FAIL for order 4 derivatives at tolerance 1e-10 error 8.93717 99557cf195SMatthew G. Knepley 100557cf195SMatthew G. Knepley but adapting to harmonics gives alright polynomials errors and much better harmonics errors. 101557cf195SMatthew G. Knepley 102557cf195SMatthew G. Knepley $ ./ex8 -dm_refine 2 -petscspace_degree 1 -qorder 1 -dim 2 -porder 1 -K 4 -num_comp 1 -use_poly 0 103557cf195SMatthew G. Knepley Function tests pass for order 1 at tolerance 1e-10 104557cf195SMatthew G. Knepley Function tests pass for order 1 derivatives at tolerance 1e-10 105557cf195SMatthew G. Knepley Interpolation tests pass for order 1 at tolerance 1e-10 106557cf195SMatthew G. Knepley Interpolation tests pass for order 1 derivatives at tolerance 1e-10 107557cf195SMatthew G. Knepley Adapting interpolator using harmonics 108557cf195SMatthew G. Knepley Interpolation poly tests FAIL for order 1 at tolerance 1e-10 error 0.0720606 109557cf195SMatthew G. Knepley Interpolation poly tests FAIL for order 1 derivatives at tolerance 1e-10 error 1.97779 110557cf195SMatthew G. Knepley Interpolation trig (0, 0) tests FAIL for order 1 at tolerance 1e-10 error 0.0398055 111557cf195SMatthew G. Knepley Interpolation trig (0, 0) tests FAIL for order 1 derivatives at tolerance 1e-10 error 0.995963 112557cf195SMatthew G. Knepley Interpolation trig (0, 1) tests FAIL for order 1 at tolerance 1e-10 error 0.0398051 113557cf195SMatthew G. Knepley Interpolation trig (0, 1) tests FAIL for order 1 derivatives at tolerance 1e-10 error 0.995964 114557cf195SMatthew G. Knepley Interpolation trig (1, 0) tests FAIL for order 2 at tolerance 1e-10 error 0.0238441 115557cf195SMatthew G. Knepley Interpolation trig (1, 0) tests FAIL for order 2 derivatives at tolerance 1e-10 error 0.888611 116557cf195SMatthew G. Knepley Interpolation trig (1, 1) tests FAIL for order 2 at tolerance 1e-10 error 0.0238346 117557cf195SMatthew G. Knepley Interpolation trig (1, 1) tests FAIL for order 2 derivatives at tolerance 1e-10 error 0.888612 118557cf195SMatthew G. Knepley Interpolation trig (2, 0) tests FAIL for order 3 at tolerance 1e-10 error 0.0537968 119557cf195SMatthew G. Knepley Interpolation trig (2, 0) tests FAIL for order 3 derivatives at tolerance 1e-10 error 1.57665 120557cf195SMatthew G. Knepley Interpolation trig (2, 1) tests FAIL for order 3 at tolerance 1e-10 error 0.0537779 121557cf195SMatthew G. Knepley Interpolation trig (2, 1) tests FAIL for order 3 derivatives at tolerance 1e-10 error 1.57666 122557cf195SMatthew G. Knepley Interpolation trig (3, 0) tests FAIL for order 4 at tolerance 1e-10 error 0.0775838 123557cf195SMatthew G. Knepley Interpolation trig (3, 0) tests FAIL for order 4 derivatives at tolerance 1e-10 error 2.36926 124557cf195SMatthew G. Knepley Interpolation trig (3, 1) tests FAIL for order 4 at tolerance 1e-10 error 0.0775464 125557cf195SMatthew G. Knepley Interpolation trig (3, 1) tests FAIL for order 4 derivatives at tolerance 1e-10 error 2.36929 126557cf195SMatthew G. Knepley */ 127557cf195SMatthew G. Knepley 128557cf195SMatthew G. Knepley typedef struct { 129557cf195SMatthew G. Knepley /* Element definition */ 130557cf195SMatthew G. Knepley PetscInt qorder; /* Order of the quadrature */ 131557cf195SMatthew G. Knepley PetscInt Nc; /* Number of field components */ 132557cf195SMatthew G. Knepley /* Testing space */ 133557cf195SMatthew G. Knepley PetscInt porder; /* Order of polynomials to test */ 134557cf195SMatthew G. Knepley PetscReal constants[3]; /* Constant values for each dimension */ 135557cf195SMatthew G. Knepley PetscInt m; /* The frequency of sinusoids to use */ 136557cf195SMatthew G. Knepley PetscInt dir; /* The direction of sinusoids to use */ 137557cf195SMatthew G. Knepley /* Adaptation */ 138557cf195SMatthew G. Knepley PetscInt K; /* Number of coarse modes used for optimization */ 139557cf195SMatthew G. Knepley PetscBool usePoly; /* Use polynomials, or harmonics, to adapt interpolator */ 140557cf195SMatthew G. Knepley } AppCtx; 141557cf195SMatthew G. Knepley 142*9371c9d4SSatish Balay typedef enum { 143*9371c9d4SSatish Balay INTERPOLATION, 144*9371c9d4SSatish Balay RESTRICTION, 145*9371c9d4SSatish Balay INJECTION 146*9371c9d4SSatish Balay } InterpType; 147557cf195SMatthew G. Knepley 148557cf195SMatthew G. Knepley /* u = 1 */ 149*9371c9d4SSatish Balay PetscErrorCode constant(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, void *ctx) { 150557cf195SMatthew G. Knepley AppCtx *user = (AppCtx *)ctx; 151557cf195SMatthew G. Knepley PetscInt d = user->dir; 152557cf195SMatthew G. Knepley 153557cf195SMatthew G. Knepley if (Nc > 1) { 154557cf195SMatthew G. Knepley for (d = 0; d < Nc; ++d) u[d] = user->constants[d]; 155557cf195SMatthew G. Knepley } else { 156557cf195SMatthew G. Knepley u[0] = user->constants[d]; 157557cf195SMatthew G. Knepley } 158557cf195SMatthew G. Knepley return 0; 159557cf195SMatthew G. Knepley } 160*9371c9d4SSatish Balay PetscErrorCode constantDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nc, PetscScalar *u, void *ctx) { 161557cf195SMatthew G. Knepley AppCtx *user = (AppCtx *)ctx; 162557cf195SMatthew G. Knepley PetscInt d = user->dir; 163557cf195SMatthew G. Knepley 164557cf195SMatthew G. Knepley if (Nc > 1) { 165557cf195SMatthew G. Knepley for (d = 0; d < Nc; ++d) u[d] = 0.0; 166557cf195SMatthew G. Knepley } else { 167557cf195SMatthew G. Knepley u[0] = user->constants[d]; 168557cf195SMatthew G. Knepley } 169557cf195SMatthew G. Knepley return 0; 170557cf195SMatthew G. Knepley } 171557cf195SMatthew G. Knepley 172557cf195SMatthew G. Knepley /* u = x */ 173*9371c9d4SSatish Balay PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, void *ctx) { 174557cf195SMatthew G. Knepley AppCtx *user = (AppCtx *)ctx; 175557cf195SMatthew G. Knepley PetscInt d = user->dir; 176557cf195SMatthew G. Knepley 177557cf195SMatthew G. Knepley if (Nc > 1) { 178557cf195SMatthew G. Knepley for (d = 0; d < Nc; ++d) u[d] = coords[d]; 179557cf195SMatthew G. Knepley } else { 180557cf195SMatthew G. Knepley u[0] = coords[d]; 181557cf195SMatthew G. Knepley } 182557cf195SMatthew G. Knepley return 0; 183557cf195SMatthew G. Knepley } 184*9371c9d4SSatish Balay PetscErrorCode linearDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nc, PetscScalar *u, void *ctx) { 185557cf195SMatthew G. Knepley AppCtx *user = (AppCtx *)ctx; 186557cf195SMatthew G. Knepley PetscInt d = user->dir; 187557cf195SMatthew G. Knepley 188557cf195SMatthew G. Knepley if (Nc > 1) { 189557cf195SMatthew G. Knepley PetscInt e; 190557cf195SMatthew G. Knepley for (d = 0; d < Nc; ++d) { 191557cf195SMatthew G. Knepley u[d] = 0.0; 192557cf195SMatthew G. Knepley for (e = 0; e < dim; ++e) u[d] += (d == e ? 1.0 : 0.0) * n[e]; 193557cf195SMatthew G. Knepley } 194557cf195SMatthew G. Knepley } else { 195557cf195SMatthew G. Knepley u[0] = n[d]; 196557cf195SMatthew G. Knepley } 197557cf195SMatthew G. Knepley return 0; 198557cf195SMatthew G. Knepley } 199557cf195SMatthew G. Knepley 200557cf195SMatthew G. Knepley /* u = x^2 or u = (x^2, xy) or u = (xy, yz, zx) */ 201*9371c9d4SSatish Balay PetscErrorCode quadratic(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, void *ctx) { 202557cf195SMatthew G. Knepley AppCtx *user = (AppCtx *)ctx; 203557cf195SMatthew G. Knepley PetscInt d = user->dir; 204557cf195SMatthew G. Knepley 205557cf195SMatthew G. Knepley if (Nc > 1) { 206*9371c9d4SSatish Balay if (Nc > 2) { 207*9371c9d4SSatish Balay u[0] = coords[0] * coords[1]; 208*9371c9d4SSatish Balay u[1] = coords[1] * coords[2]; 209*9371c9d4SSatish Balay u[2] = coords[2] * coords[0]; 210*9371c9d4SSatish Balay } else { 211*9371c9d4SSatish Balay u[0] = coords[0] * coords[0]; 212*9371c9d4SSatish Balay u[1] = coords[0] * coords[1]; 213*9371c9d4SSatish Balay } 214557cf195SMatthew G. Knepley } else { 215557cf195SMatthew G. Knepley u[0] = coords[d] * coords[d]; 216557cf195SMatthew G. Knepley } 217557cf195SMatthew G. Knepley return 0; 218557cf195SMatthew G. Knepley } 219*9371c9d4SSatish Balay PetscErrorCode quadraticDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nc, PetscScalar *u, void *ctx) { 220557cf195SMatthew G. Knepley AppCtx *user = (AppCtx *)ctx; 221557cf195SMatthew G. Knepley PetscInt d = user->dir; 222557cf195SMatthew G. Knepley 223557cf195SMatthew G. Knepley if (Nc > 1) { 224*9371c9d4SSatish Balay if (Nc > 2) { 225*9371c9d4SSatish Balay u[0] = coords[1] * n[0] + coords[0] * n[1]; 226*9371c9d4SSatish Balay u[1] = coords[2] * n[1] + coords[1] * n[2]; 227*9371c9d4SSatish Balay u[2] = coords[2] * n[0] + coords[0] * n[2]; 228*9371c9d4SSatish Balay } else { 229*9371c9d4SSatish Balay u[0] = 2.0 * coords[0] * n[0]; 230*9371c9d4SSatish Balay u[1] = coords[1] * n[0] + coords[0] * n[1]; 231*9371c9d4SSatish Balay } 232557cf195SMatthew G. Knepley } else { 233557cf195SMatthew G. Knepley u[0] = 2.0 * coords[d] * n[d]; 234557cf195SMatthew G. Knepley } 235557cf195SMatthew G. Knepley return 0; 236557cf195SMatthew G. Knepley } 237557cf195SMatthew G. Knepley 238557cf195SMatthew G. Knepley /* u = x^3 or u = (x^3, x^2y) or u = (x^2y, y^2z, z^2x) */ 239*9371c9d4SSatish Balay PetscErrorCode cubic(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, void *ctx) { 240557cf195SMatthew G. Knepley AppCtx *user = (AppCtx *)ctx; 241557cf195SMatthew G. Knepley PetscInt d = user->dir; 242557cf195SMatthew G. Knepley 243557cf195SMatthew G. Knepley if (Nc > 1) { 244*9371c9d4SSatish Balay if (Nc > 2) { 245*9371c9d4SSatish Balay u[0] = coords[0] * coords[0] * coords[1]; 246*9371c9d4SSatish Balay u[1] = coords[1] * coords[1] * coords[2]; 247*9371c9d4SSatish Balay u[2] = coords[2] * coords[2] * coords[0]; 248*9371c9d4SSatish Balay } else { 249*9371c9d4SSatish Balay u[0] = coords[0] * coords[0] * coords[0]; 250*9371c9d4SSatish Balay u[1] = coords[0] * coords[0] * coords[1]; 251*9371c9d4SSatish Balay } 252557cf195SMatthew G. Knepley } else { 253557cf195SMatthew G. Knepley u[0] = coords[d] * coords[d] * coords[d]; 254557cf195SMatthew G. Knepley } 255557cf195SMatthew G. Knepley return 0; 256557cf195SMatthew G. Knepley } 257*9371c9d4SSatish Balay PetscErrorCode cubicDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nc, PetscScalar *u, void *ctx) { 258557cf195SMatthew G. Knepley AppCtx *user = (AppCtx *)ctx; 259557cf195SMatthew G. Knepley PetscInt d = user->dir; 260557cf195SMatthew G. Knepley 261557cf195SMatthew G. Knepley if (Nc > 1) { 262*9371c9d4SSatish Balay if (Nc > 2) { 263*9371c9d4SSatish Balay u[0] = 2.0 * coords[0] * coords[1] * n[0] + coords[0] * coords[0] * n[1]; 264*9371c9d4SSatish Balay u[1] = 2.0 * coords[1] * coords[2] * n[1] + coords[1] * coords[1] * n[2]; 265*9371c9d4SSatish Balay u[2] = 2.0 * coords[2] * coords[0] * n[2] + coords[2] * coords[2] * n[0]; 266*9371c9d4SSatish Balay } else { 267*9371c9d4SSatish Balay u[0] = 3.0 * coords[0] * coords[0] * n[0]; 268*9371c9d4SSatish Balay u[1] = 2.0 * coords[0] * coords[1] * n[0] + coords[0] * coords[0] * n[1]; 269*9371c9d4SSatish Balay } 270557cf195SMatthew G. Knepley } else { 271557cf195SMatthew G. Knepley u[0] = 3.0 * coords[d] * coords[d] * n[d]; 272557cf195SMatthew G. Knepley } 273557cf195SMatthew G. Knepley return 0; 274557cf195SMatthew G. Knepley } 275557cf195SMatthew G. Knepley 276557cf195SMatthew G. Knepley /* u = x^4 or u = (x^4, x^2y^2) or u = (x^2y^2, y^2z^2, z^2x^2) */ 277*9371c9d4SSatish Balay PetscErrorCode quartic(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, void *ctx) { 278557cf195SMatthew G. Knepley AppCtx *user = (AppCtx *)ctx; 279557cf195SMatthew G. Knepley PetscInt d = user->dir; 280557cf195SMatthew G. Knepley 281557cf195SMatthew G. Knepley if (Nc > 1) { 282*9371c9d4SSatish Balay if (Nc > 2) { 283*9371c9d4SSatish Balay u[0] = coords[0] * coords[0] * coords[1] * coords[1]; 284*9371c9d4SSatish Balay u[1] = coords[1] * coords[1] * coords[2] * coords[2]; 285*9371c9d4SSatish Balay u[2] = coords[2] * coords[2] * coords[0] * coords[0]; 286*9371c9d4SSatish Balay } else { 287*9371c9d4SSatish Balay u[0] = coords[0] * coords[0] * coords[0] * coords[0]; 288*9371c9d4SSatish Balay u[1] = coords[0] * coords[0] * coords[1] * coords[1]; 289*9371c9d4SSatish Balay } 290557cf195SMatthew G. Knepley } else { 291557cf195SMatthew G. Knepley u[0] = coords[d] * coords[d] * coords[d] * coords[d]; 292557cf195SMatthew G. Knepley } 293557cf195SMatthew G. Knepley return 0; 294557cf195SMatthew G. Knepley } 295*9371c9d4SSatish Balay PetscErrorCode quarticDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nc, PetscScalar *u, void *ctx) { 296557cf195SMatthew G. Knepley AppCtx *user = (AppCtx *)ctx; 297557cf195SMatthew G. Knepley PetscInt d = user->dir; 298557cf195SMatthew G. Knepley 299557cf195SMatthew G. Knepley if (Nc > 1) { 300*9371c9d4SSatish Balay if (Nc > 2) { 301*9371c9d4SSatish Balay u[0] = 2.0 * coords[0] * coords[1] * coords[1] * n[0] + 2.0 * coords[0] * coords[0] * coords[1] * n[1]; 302557cf195SMatthew G. Knepley u[1] = 2.0 * coords[1] * coords[2] * coords[2] * n[1] + 2.0 * coords[1] * coords[1] * coords[2] * n[2]; 303*9371c9d4SSatish Balay u[2] = 2.0 * coords[2] * coords[0] * coords[0] * n[2] + 2.0 * coords[2] * coords[2] * coords[0] * n[0]; 304*9371c9d4SSatish Balay } else { 305*9371c9d4SSatish Balay u[0] = 4.0 * coords[0] * coords[0] * coords[0] * n[0]; 306*9371c9d4SSatish Balay u[1] = 2.0 * coords[0] * coords[1] * coords[1] * n[0] + 2.0 * coords[0] * coords[0] * coords[1] * n[1]; 307*9371c9d4SSatish Balay } 308557cf195SMatthew G. Knepley } else { 309557cf195SMatthew G. Knepley u[0] = 4.0 * coords[d] * coords[d] * coords[d] * n[d]; 310557cf195SMatthew G. Knepley } 311557cf195SMatthew G. Knepley return 0; 312557cf195SMatthew G. Knepley } 313557cf195SMatthew G. Knepley 314*9371c9d4SSatish Balay PetscErrorCode mytanh(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, void *ctx) { 315557cf195SMatthew G. Knepley AppCtx *user = (AppCtx *)ctx; 316557cf195SMatthew G. Knepley PetscInt d = user->dir; 317557cf195SMatthew G. Knepley 318557cf195SMatthew G. Knepley if (Nc > 1) { 319557cf195SMatthew G. Knepley for (d = 0; d < Nc; ++d) u[d] = PetscTanhReal(coords[d] - 0.5); 320557cf195SMatthew G. Knepley } else { 321557cf195SMatthew G. Knepley u[0] = PetscTanhReal(coords[d] - 0.5); 322557cf195SMatthew G. Knepley } 323557cf195SMatthew G. Knepley return 0; 324557cf195SMatthew G. Knepley } 325*9371c9d4SSatish Balay PetscErrorCode mytanhDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nc, PetscScalar *u, void *ctx) { 326557cf195SMatthew G. Knepley AppCtx *user = (AppCtx *)ctx; 327557cf195SMatthew G. Knepley PetscInt d = user->dir; 328557cf195SMatthew G. Knepley 329557cf195SMatthew G. Knepley if (Nc > 1) { 330557cf195SMatthew G. Knepley for (d = 0; d < Nc; ++d) u[d] = 1.0 / PetscSqr(PetscCoshReal(coords[d] - 0.5)) * n[d]; 331557cf195SMatthew G. Knepley } else { 332557cf195SMatthew G. Knepley u[0] = 1.0 / PetscSqr(PetscCoshReal(coords[d] - 0.5)) * n[d]; 333557cf195SMatthew G. Knepley } 334557cf195SMatthew G. Knepley return 0; 335557cf195SMatthew G. Knepley } 336557cf195SMatthew G. Knepley 337*9371c9d4SSatish Balay PetscErrorCode trig(PetscInt dim, PetscReal time, const PetscReal coords[], PetscInt Nc, PetscScalar *u, void *ctx) { 338557cf195SMatthew G. Knepley AppCtx *user = (AppCtx *)ctx; 339557cf195SMatthew G. Knepley PetscInt m = user->m, d = user->dir; 340557cf195SMatthew G. Knepley 341557cf195SMatthew G. Knepley if (Nc > 1) { 342557cf195SMatthew G. Knepley for (d = 0; d < Nc; ++d) u[d] = PetscSinReal(PETSC_PI * m * coords[d]); 343557cf195SMatthew G. Knepley } else { 344557cf195SMatthew G. Knepley u[0] = PetscSinReal(PETSC_PI * m * coords[d]); 345557cf195SMatthew G. Knepley } 346557cf195SMatthew G. Knepley return 0; 347557cf195SMatthew G. Knepley } 348*9371c9d4SSatish Balay PetscErrorCode trigDer(PetscInt dim, PetscReal time, const PetscReal coords[], const PetscReal n[], PetscInt Nc, PetscScalar *u, void *ctx) { 349557cf195SMatthew G. Knepley AppCtx *user = (AppCtx *)ctx; 350557cf195SMatthew G. Knepley PetscInt m = user->m, d = user->dir; 351557cf195SMatthew G. Knepley 352557cf195SMatthew G. Knepley if (Nc > 1) { 353557cf195SMatthew G. Knepley for (d = 0; d < Nc; ++d) u[d] = PETSC_PI * m * PetscCosReal(PETSC_PI * m * coords[d]) * n[d]; 354557cf195SMatthew G. Knepley } else { 355557cf195SMatthew G. Knepley u[0] = PETSC_PI * m * PetscCosReal(PETSC_PI * m * coords[d]) * n[d]; 356557cf195SMatthew G. Knepley } 357557cf195SMatthew G. Knepley return 0; 358557cf195SMatthew G. Knepley } 359557cf195SMatthew G. Knepley 360*9371c9d4SSatish Balay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) { 361557cf195SMatthew G. Knepley PetscFunctionBeginUser; 362557cf195SMatthew G. Knepley options->qorder = 0; 363557cf195SMatthew G. Knepley options->Nc = PETSC_DEFAULT; 364557cf195SMatthew G. Knepley options->porder = 0; 365557cf195SMatthew G. Knepley options->m = 1; 366557cf195SMatthew G. Knepley options->dir = 0; 367557cf195SMatthew G. Knepley options->K = 0; 368557cf195SMatthew G. Knepley options->usePoly = PETSC_TRUE; 369557cf195SMatthew G. Knepley 370d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Projection Test Options", "DMPlex"); 3719566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-qorder", "The quadrature order", "ex8.c", options->qorder, &options->qorder, NULL)); 3729566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-num_comp", "The number of field components", "ex8.c", options->Nc, &options->Nc, NULL)); 3739566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-porder", "The order of polynomials to test", "ex8.c", options->porder, &options->porder, NULL)); 3749566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-K", "The number of coarse modes used in optimization", "ex8.c", options->K, &options->K, NULL)); 3759566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-use_poly", "Use polynomials (or harmonics) to adapt interpolator", "ex8.c", options->usePoly, &options->usePoly, NULL)); 376d0609cedSBarry Smith PetscOptionsEnd(); 377557cf195SMatthew G. Knepley PetscFunctionReturn(0); 378557cf195SMatthew G. Knepley } 379557cf195SMatthew G. Knepley 380*9371c9d4SSatish Balay static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) { 381557cf195SMatthew G. Knepley PetscFunctionBeginUser; 3829566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 3839566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 3849566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 3859566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 386557cf195SMatthew G. Knepley PetscFunctionReturn(0); 387557cf195SMatthew G. Knepley } 388557cf195SMatthew G. Knepley 389557cf195SMatthew G. Knepley /* Setup functions to approximate */ 390*9371c9d4SSatish Balay static PetscErrorCode SetupFunctions(DM dm, PetscBool usePoly, PetscInt order, PetscInt dir, PetscErrorCode (**exactFuncs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), PetscErrorCode (**exactFuncDers)(PetscInt, PetscReal, const PetscReal[], const PetscReal[], PetscInt, PetscScalar *, void *), AppCtx *user) { 391557cf195SMatthew G. Knepley PetscInt dim; 392557cf195SMatthew G. Knepley 393557cf195SMatthew G. Knepley PetscFunctionBeginUser; 394557cf195SMatthew G. Knepley user->dir = dir; 395557cf195SMatthew G. Knepley if (usePoly) { 396557cf195SMatthew G. Knepley switch (order) { 397557cf195SMatthew G. Knepley case 0: 398557cf195SMatthew G. Knepley exactFuncs[0] = constant; 399557cf195SMatthew G. Knepley exactFuncDers[0] = constantDer; 400557cf195SMatthew G. Knepley break; 401557cf195SMatthew G. Knepley case 1: 402557cf195SMatthew G. Knepley exactFuncs[0] = linear; 403557cf195SMatthew G. Knepley exactFuncDers[0] = linearDer; 404557cf195SMatthew G. Knepley break; 405557cf195SMatthew G. Knepley case 2: 406557cf195SMatthew G. Knepley exactFuncs[0] = quadratic; 407557cf195SMatthew G. Knepley exactFuncDers[0] = quadraticDer; 408557cf195SMatthew G. Knepley break; 409557cf195SMatthew G. Knepley case 3: 410557cf195SMatthew G. Knepley exactFuncs[0] = cubic; 411557cf195SMatthew G. Knepley exactFuncDers[0] = cubicDer; 412557cf195SMatthew G. Knepley break; 413557cf195SMatthew G. Knepley case 4: 414557cf195SMatthew G. Knepley exactFuncs[0] = quartic; 415557cf195SMatthew G. Knepley exactFuncDers[0] = quarticDer; 416557cf195SMatthew G. Knepley break; 417*9371c9d4SSatish Balay default: PetscCall(DMGetDimension(dm, &dim)); SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Could not determine functions to test for dimension %" PetscInt_FMT " order %" PetscInt_FMT, dim, order); 418557cf195SMatthew G. Knepley } 419557cf195SMatthew G. Knepley } else { 420557cf195SMatthew G. Knepley user->m = order; 421557cf195SMatthew G. Knepley exactFuncs[0] = trig; 422557cf195SMatthew G. Knepley exactFuncDers[0] = trigDer; 423557cf195SMatthew G. Knepley } 424557cf195SMatthew G. Knepley PetscFunctionReturn(0); 425557cf195SMatthew G. Knepley } 426557cf195SMatthew G. Knepley 427*9371c9d4SSatish Balay static PetscErrorCode ComputeError(DM dm, PetscErrorCode (**exactFuncs)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), PetscErrorCode (**exactFuncDers)(PetscInt, PetscReal, const PetscReal[], const PetscReal[], PetscInt, PetscScalar *, void *), void **exactCtxs, PetscReal *error, PetscReal *errorDer, AppCtx *user) { 428557cf195SMatthew G. Knepley Vec u; 429557cf195SMatthew G. Knepley PetscReal n[3] = {1.0, 1.0, 1.0}; 430557cf195SMatthew G. Knepley 431557cf195SMatthew G. Knepley PetscFunctionBeginUser; 4329566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &u)); 433557cf195SMatthew G. Knepley /* Project function into FE function space */ 4349566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, 0.0, exactFuncs, exactCtxs, INSERT_ALL_VALUES, u)); 4359566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-projection_view")); 436557cf195SMatthew G. Knepley /* Compare approximation to exact in L_2 */ 4379566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(dm, 0.0, exactFuncs, exactCtxs, u, error)); 4389566063dSJacob Faibussowitsch PetscCall(DMComputeL2GradientDiff(dm, 0.0, exactFuncDers, exactCtxs, u, n, errorDer)); 4399566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &u)); 440557cf195SMatthew G. Knepley PetscFunctionReturn(0); 441557cf195SMatthew G. Knepley } 442557cf195SMatthew G. Knepley 443*9371c9d4SSatish Balay static PetscErrorCode CheckFunctions(DM dm, PetscInt order, AppCtx *user) { 444557cf195SMatthew G. Knepley PetscErrorCode (*exactFuncs[1])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 445557cf195SMatthew G. Knepley PetscErrorCode (*exactFuncDers[1])(PetscInt dim, PetscReal time, const PetscReal x[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx); 446557cf195SMatthew G. Knepley void *exactCtxs[3]; 447557cf195SMatthew G. Knepley MPI_Comm comm; 448557cf195SMatthew G. Knepley PetscReal error, errorDer, tol = PETSC_SMALL; 449557cf195SMatthew G. Knepley 450557cf195SMatthew G. Knepley PetscFunctionBeginUser; 451557cf195SMatthew G. Knepley exactCtxs[0] = user; 452557cf195SMatthew G. Knepley exactCtxs[1] = user; 453557cf195SMatthew G. Knepley exactCtxs[2] = user; 454557cf195SMatthew G. Knepley user->constants[0] = 1.0; 455557cf195SMatthew G. Knepley user->constants[1] = 2.0; 456557cf195SMatthew G. Knepley user->constants[2] = 3.0; 4579566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 4589566063dSJacob Faibussowitsch PetscCall(SetupFunctions(dm, PETSC_TRUE, order, 0, exactFuncs, exactFuncDers, user)); 4599566063dSJacob Faibussowitsch PetscCall(ComputeError(dm, exactFuncs, exactFuncDers, exactCtxs, &error, &errorDer, user)); 460557cf195SMatthew G. Knepley /* Report result */ 46163a3b9bcSJacob Faibussowitsch if (error > tol) PetscCall(PetscPrintf(comm, "Function tests FAIL for order %" PetscInt_FMT " at tolerance %g error %g\n", order, (double)tol, (double)error)); 46263a3b9bcSJacob Faibussowitsch else PetscCall(PetscPrintf(comm, "Function tests pass for order %" PetscInt_FMT " at tolerance %g\n", order, (double)tol)); 46363a3b9bcSJacob Faibussowitsch if (errorDer > tol) PetscCall(PetscPrintf(comm, "Function tests FAIL for order %" PetscInt_FMT " derivatives at tolerance %g error %g\n", order, (double)tol, (double)errorDer)); 46463a3b9bcSJacob Faibussowitsch else PetscCall(PetscPrintf(comm, "Function tests pass for order %" PetscInt_FMT " derivatives at tolerance %g\n", order, (double)tol)); 465557cf195SMatthew G. Knepley PetscFunctionReturn(0); 466557cf195SMatthew G. Knepley } 467557cf195SMatthew G. Knepley 468557cf195SMatthew G. Knepley /* Compare approximation to exact in L_2 */ 469*9371c9d4SSatish Balay static PetscErrorCode CheckTransferError(DM fdm, PetscBool usePoly, PetscInt order, PetscInt dir, const char *testname, Vec fu, AppCtx *user) { 470557cf195SMatthew G. Knepley PetscErrorCode (*exactFuncs[1])(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx); 471557cf195SMatthew G. Knepley PetscErrorCode (*exactFuncDers[1])(PetscInt dim, PetscReal time, const PetscReal x[], const PetscReal n[], PetscInt Nf, PetscScalar *u, void *ctx); 472557cf195SMatthew G. Knepley PetscReal n[3] = {1.0, 1.0, 1.0}; 473557cf195SMatthew G. Knepley void *exactCtxs[3]; 474557cf195SMatthew G. Knepley MPI_Comm comm; 475557cf195SMatthew G. Knepley PetscReal error, errorDer, tol = PETSC_SMALL; 476557cf195SMatthew G. Knepley 477557cf195SMatthew G. Knepley PetscFunctionBeginUser; 478557cf195SMatthew G. Knepley exactCtxs[0] = user; 479557cf195SMatthew G. Knepley exactCtxs[1] = user; 480557cf195SMatthew G. Knepley exactCtxs[2] = user; 481557cf195SMatthew G. Knepley user->constants[0] = 1.0; 482557cf195SMatthew G. Knepley user->constants[1] = 2.0; 483557cf195SMatthew G. Knepley user->constants[2] = 3.0; 4849566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)fdm, &comm)); 4859566063dSJacob Faibussowitsch PetscCall(SetupFunctions(fdm, usePoly, order, dir, exactFuncs, exactFuncDers, user)); 4869566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(fdm, 0.0, exactFuncs, exactCtxs, fu, &error)); 4879566063dSJacob Faibussowitsch PetscCall(DMComputeL2GradientDiff(fdm, 0.0, exactFuncDers, exactCtxs, fu, n, &errorDer)); 488557cf195SMatthew G. Knepley /* Report result */ 48963a3b9bcSJacob Faibussowitsch if (error > tol) PetscCall(PetscPrintf(comm, "%s tests FAIL for order %" PetscInt_FMT " at tolerance %g error %g\n", testname, order, (double)tol, (double)error)); 49063a3b9bcSJacob Faibussowitsch else PetscCall(PetscPrintf(comm, "%s tests pass for order %" PetscInt_FMT " at tolerance %g\n", testname, order, (double)tol)); 49163a3b9bcSJacob Faibussowitsch if (errorDer > tol) PetscCall(PetscPrintf(comm, "%s tests FAIL for order %" PetscInt_FMT " derivatives at tolerance %g error %g\n", testname, order, (double)tol, (double)errorDer)); 49263a3b9bcSJacob Faibussowitsch else PetscCall(PetscPrintf(comm, "%s tests pass for order %" PetscInt_FMT " derivatives at tolerance %g\n", testname, order, (double)tol)); 493557cf195SMatthew G. Knepley PetscFunctionReturn(0); 494557cf195SMatthew G. Knepley } 495557cf195SMatthew G. Knepley 496*9371c9d4SSatish Balay static PetscErrorCode CheckTransfer(DM dm, InterpType inType, PetscInt order, AppCtx *user) { 497557cf195SMatthew G. Knepley PetscErrorCode (*exactFuncs[1])(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx); 498557cf195SMatthew G. Knepley PetscErrorCode (*exactFuncDers[1])(PetscInt, PetscReal, const PetscReal x[], const PetscReal n[], PetscInt, PetscScalar *u, void *ctx); 499d6837840SMatthew G. Knepley void *exactCtxs[3]; 500d6837840SMatthew G. Knepley DM rdm = NULL, idm = NULL, fdm = NULL; 501557cf195SMatthew G. Knepley Mat Interp, InterpAdapt = NULL; 502d6837840SMatthew G. Knepley Vec iu, fu, scaling = NULL; 503557cf195SMatthew G. Knepley MPI_Comm comm; 504d6837840SMatthew G. Knepley const char *testname = "Unknown"; 505557cf195SMatthew G. Knepley char checkname[PETSC_MAX_PATH_LEN]; 506557cf195SMatthew G. Knepley 507557cf195SMatthew G. Knepley PetscFunctionBeginUser; 508d6837840SMatthew G. Knepley exactCtxs[0] = exactCtxs[1] = exactCtxs[2] = user; 5099566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 5109566063dSJacob Faibussowitsch PetscCall(DMRefine(dm, comm, &rdm)); 5119566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(rdm, NULL, "-ref_dm_view")); 5129566063dSJacob Faibussowitsch PetscCall(DMSetCoarseDM(rdm, dm)); 5139566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, rdm)); 514557cf195SMatthew G. Knepley switch (inType) { 515557cf195SMatthew G. Knepley case INTERPOLATION: 516557cf195SMatthew G. Knepley testname = "Interpolation"; 517557cf195SMatthew G. Knepley idm = dm; 518557cf195SMatthew G. Knepley fdm = rdm; 519557cf195SMatthew G. Knepley break; 520557cf195SMatthew G. Knepley case RESTRICTION: 521557cf195SMatthew G. Knepley testname = "Restriction"; 522557cf195SMatthew G. Knepley idm = rdm; 523557cf195SMatthew G. Knepley fdm = dm; 524557cf195SMatthew G. Knepley break; 525557cf195SMatthew G. Knepley case INJECTION: 526557cf195SMatthew G. Knepley testname = "Injection"; 527557cf195SMatthew G. Knepley idm = rdm; 528557cf195SMatthew G. Knepley fdm = dm; 529557cf195SMatthew G. Knepley break; 530557cf195SMatthew G. Knepley } 5319566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(idm, &iu)); 5329566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(fdm, &fu)); 5339566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(dm, user)); 5349566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(rdm, user)); 535557cf195SMatthew G. Knepley /* Project function into initial FE function space */ 5369566063dSJacob Faibussowitsch PetscCall(SetupFunctions(dm, PETSC_TRUE, order, 0, exactFuncs, exactFuncDers, user)); 5379566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(idm, 0.0, exactFuncs, exactCtxs, INSERT_ALL_VALUES, iu)); 538557cf195SMatthew G. Knepley /* Interpolate function into final FE function space */ 539557cf195SMatthew G. Knepley switch (inType) { 540557cf195SMatthew G. Knepley case INTERPOLATION: 5419566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(dm, rdm, &Interp, &scaling)); 5429566063dSJacob Faibussowitsch PetscCall(MatInterpolate(Interp, iu, fu)); 543557cf195SMatthew G. Knepley break; 544557cf195SMatthew G. Knepley case RESTRICTION: 5459566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(dm, rdm, &Interp, &scaling)); 5469566063dSJacob Faibussowitsch PetscCall(MatRestrict(Interp, iu, fu)); 5479566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(fu, scaling, fu)); 548557cf195SMatthew G. Knepley break; 549557cf195SMatthew G. Knepley case INJECTION: 5509566063dSJacob Faibussowitsch PetscCall(DMCreateInjection(dm, rdm, &Interp)); 5519566063dSJacob Faibussowitsch PetscCall(MatRestrict(Interp, iu, fu)); 552557cf195SMatthew G. Knepley break; 553557cf195SMatthew G. Knepley } 5549566063dSJacob Faibussowitsch PetscCall(CheckTransferError(fdm, PETSC_TRUE, order, 0, testname, fu, user)); 555557cf195SMatthew G. Knepley if (user->K && (inType == INTERPOLATION)) { 556557cf195SMatthew G. Knepley KSP smoother; 5572b3cbbdaSStefano Zampini Mat A, iVM, fVM; 5582b3cbbdaSStefano Zampini Vec iV, fV; 5592b3cbbdaSStefano Zampini PetscInt k, dim, d, im, fm; 560557cf195SMatthew G. Knepley 5619566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, " Adapting interpolator using %s\n", user->usePoly ? "polynomials" : "harmonics")); 5629566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 563557cf195SMatthew G. Knepley /* Project coarse modes into initial and final FE function space */ 5642b3cbbdaSStefano Zampini PetscCall(DMGetGlobalVector(idm, &iV)); 5652b3cbbdaSStefano Zampini PetscCall(DMGetGlobalVector(fdm, &fV)); 5662b3cbbdaSStefano Zampini PetscCall(VecGetLocalSize(iV, &im)); 5672b3cbbdaSStefano Zampini PetscCall(VecGetLocalSize(fV, &fm)); 5682b3cbbdaSStefano Zampini PetscCall(MatCreateDense(PetscObjectComm((PetscObject)dm), im, PETSC_DECIDE, PETSC_DECIDE, user->K * dim, NULL, &iVM)); 5692b3cbbdaSStefano Zampini PetscCall(MatCreateDense(PetscObjectComm((PetscObject)dm), fm, PETSC_DECIDE, PETSC_DECIDE, user->K * dim, NULL, &fVM)); 5702b3cbbdaSStefano Zampini PetscCall(DMRestoreGlobalVector(idm, &iV)); 5712b3cbbdaSStefano Zampini PetscCall(DMRestoreGlobalVector(fdm, &fV)); 572557cf195SMatthew G. Knepley for (k = 0; k < user->K; ++k) { 573557cf195SMatthew G. Knepley for (d = 0; d < dim; ++d) { 5742b3cbbdaSStefano Zampini PetscCall(MatDenseGetColumnVecWrite(iVM, k * dim + d, &iV)); 5752b3cbbdaSStefano Zampini PetscCall(MatDenseGetColumnVecWrite(fVM, k * dim + d, &fV)); 5769566063dSJacob Faibussowitsch PetscCall(SetupFunctions(idm, user->usePoly, user->usePoly ? k : k + 1, d, exactFuncs, exactFuncDers, user)); 5772b3cbbdaSStefano Zampini PetscCall(DMProjectFunction(idm, 0.0, exactFuncs, exactCtxs, INSERT_ALL_VALUES, iV)); 5782b3cbbdaSStefano Zampini PetscCall(DMProjectFunction(fdm, 0.0, exactFuncs, exactCtxs, INSERT_ALL_VALUES, fV)); 5792b3cbbdaSStefano Zampini PetscCall(MatDenseRestoreColumnVecWrite(iVM, k * dim + d, &iV)); 5802b3cbbdaSStefano Zampini PetscCall(MatDenseRestoreColumnVecWrite(fVM, k * dim + d, &fV)); 581557cf195SMatthew G. Knepley } 582557cf195SMatthew G. Knepley } 5832b3cbbdaSStefano Zampini 584557cf195SMatthew G. Knepley /* Adapt interpolator */ 5859566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(rdm, &A)); 5869566063dSJacob Faibussowitsch PetscCall(MatShift(A, 1.0)); 5879566063dSJacob Faibussowitsch PetscCall(KSPCreate(comm, &smoother)); 5889566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(smoother)); 5899566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(smoother, A, A)); 5902b3cbbdaSStefano Zampini PetscCall(DMAdaptInterpolator(dm, rdm, Interp, smoother, fVM, iVM, &InterpAdapt, user)); 591557cf195SMatthew G. Knepley /* Interpolate function into final FE function space */ 5929566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(checkname, PETSC_MAX_PATH_LEN, " %s poly", testname)); 5939566063dSJacob Faibussowitsch PetscCall(MatInterpolate(InterpAdapt, iu, fu)); 5949566063dSJacob Faibussowitsch PetscCall(CheckTransferError(fdm, PETSC_TRUE, order, 0, checkname, fu, user)); 595557cf195SMatthew G. Knepley for (k = 0; k < user->K; ++k) { 596557cf195SMatthew G. Knepley for (d = 0; d < dim; ++d) { 59763a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(checkname, PETSC_MAX_PATH_LEN, " %s trig (%" PetscInt_FMT ", %" PetscInt_FMT ")", testname, k, d)); 5982b3cbbdaSStefano Zampini PetscCall(MatDenseGetColumnVecRead(iVM, k * dim + d, &iV)); 5992b3cbbdaSStefano Zampini PetscCall(MatDenseGetColumnVecWrite(fVM, k * dim + d, &fV)); 6002b3cbbdaSStefano Zampini PetscCall(MatInterpolate(InterpAdapt, iV, fV)); 6012b3cbbdaSStefano Zampini PetscCall(CheckTransferError(fdm, PETSC_FALSE, k + 1, d, checkname, fV, user)); 6022b3cbbdaSStefano Zampini PetscCall(MatDenseRestoreColumnVecRead(iVM, k * dim + d, &iV)); 6032b3cbbdaSStefano Zampini PetscCall(MatDenseRestoreColumnVecWrite(fVM, k * dim + d, &fV)); 604557cf195SMatthew G. Knepley } 605557cf195SMatthew G. Knepley } 606557cf195SMatthew G. Knepley /* Cleanup */ 6079566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&smoother)); 6089566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 6099566063dSJacob Faibussowitsch PetscCall(MatDestroy(&InterpAdapt)); 6102b3cbbdaSStefano Zampini PetscCall(MatDestroy(&iVM)); 6112b3cbbdaSStefano Zampini PetscCall(MatDestroy(&fVM)); 612557cf195SMatthew G. Knepley } 6139566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(idm, &iu)); 6149566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(fdm, &fu)); 6159566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Interp)); 6169566063dSJacob Faibussowitsch PetscCall(VecDestroy(&scaling)); 6179566063dSJacob Faibussowitsch PetscCall(DMDestroy(&rdm)); 618557cf195SMatthew G. Knepley PetscFunctionReturn(0); 619557cf195SMatthew G. Knepley } 620557cf195SMatthew G. Knepley 621*9371c9d4SSatish Balay int main(int argc, char **argv) { 622557cf195SMatthew G. Knepley DM dm; 62330602db0SMatthew G. Knepley PetscFE fe; 62430602db0SMatthew G. Knepley AppCtx user; 62530602db0SMatthew G. Knepley PetscInt dim; 62630602db0SMatthew G. Knepley PetscBool simplex; 627557cf195SMatthew G. Knepley 628327415f7SBarry Smith PetscFunctionBeginUser; 6299566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 6309566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 6319566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 63230602db0SMatthew G. Knepley 6339566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 6349566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex)); 6359566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_WORLD, dim, user.Nc < 0 ? dim : user.Nc, simplex, NULL, user.qorder, &fe)); 6369566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe)); 6379566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 6389566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 63930602db0SMatthew G. Knepley 6409566063dSJacob Faibussowitsch PetscCall(CheckFunctions(dm, user.porder, &user)); 6419566063dSJacob Faibussowitsch PetscCall(CheckTransfer(dm, INTERPOLATION, user.porder, &user)); 6429566063dSJacob Faibussowitsch PetscCall(CheckTransfer(dm, INJECTION, user.porder, &user)); 6439566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 6449566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 645b122ec5aSJacob Faibussowitsch return 0; 646557cf195SMatthew G. Knepley } 647557cf195SMatthew G. Knepley 648557cf195SMatthew G. Knepley /*TEST 649557cf195SMatthew G. Knepley 650557cf195SMatthew G. Knepley # TODO dim 3 will not work until I get composite elements in 3D (see plexrefine.c:34) 651557cf195SMatthew G. Knepley # 2D/3D P_1 on a simplex 652557cf195SMatthew G. Knepley test: 653557cf195SMatthew G. Knepley suffix: p1 654557cf195SMatthew G. Knepley requires: triangle ctetgen 65530602db0SMatthew G. Knepley args: -dm_plex_dim {{2}separate output} -dm_plex_box_faces 2,2,2 -petscspace_degree 1 -num_comp 1 -qorder 1 -porder {{1}separate output} 656557cf195SMatthew G. Knepley test: 657557cf195SMatthew G. Knepley suffix: p1_pragmatic 658557cf195SMatthew G. Knepley requires: triangle ctetgen pragmatic 65930602db0SMatthew G. Knepley args: -dm_plex_dim {{2}separate output} -dm_plex_box_faces 2,2,2 -petscspace_degree 1 -qorder 1 -dm_plex_hash_location -porder {{1 2}separate output} 660557cf195SMatthew G. Knepley test: 661557cf195SMatthew G. Knepley suffix: p1_adapt 662557cf195SMatthew G. Knepley requires: triangle ctetgen 66330602db0SMatthew G. Knepley args: -dm_plex_dim {{2}separate output} -dm_plex_box_faces 2,2,2 -dm_refine 3 -petscspace_degree 1 -qorder 1 -porder {{1 2}separate output} 664557cf195SMatthew G. Knepley 665557cf195SMatthew G. Knepley # TODO dim 3 will not work until I get composite elements in 3D (see plexrefine.c:34) 666557cf195SMatthew G. Knepley # 2D/3D P_2 on a simplex 667557cf195SMatthew G. Knepley test: 668557cf195SMatthew G. Knepley suffix: p2 669557cf195SMatthew G. Knepley requires: triangle ctetgen 67030602db0SMatthew G. Knepley args: -dm_plex_dim {{2}separate output} -dm_plex_box_faces 2,2,2 -petscspace_degree 2 -qorder 2 -porder {{1 2 3}separate output} 671557cf195SMatthew G. Knepley test: 672557cf195SMatthew G. Knepley suffix: p2_pragmatic 673557cf195SMatthew G. Knepley requires: triangle ctetgen pragmatic 67430602db0SMatthew G. Knepley args: -dm_plex_dim {{2}separate output} -dm_plex_box_faces 2,2,2 -petscspace_degree 2 -qorder 2 -dm_plex_hash_location -porder {{1 2 3}separate output} 675557cf195SMatthew G. Knepley 676557cf195SMatthew G. Knepley # TODO dim 3 will not work until I get composite elements in 3D (see plexrefine.c:34) 677557cf195SMatthew G. Knepley # TODO This is broken. Check ex3 which worked 678557cf195SMatthew G. Knepley # 2D/3D P_3 on a simplex 679557cf195SMatthew G. Knepley test: 680d6837840SMatthew G. Knepley TODO: gll Lagrange nodes break this 681557cf195SMatthew G. Knepley suffix: p3 682557cf195SMatthew G. Knepley requires: triangle ctetgen !single 68330602db0SMatthew G. Knepley args: -dm_plex_dim {{2}separate output} -dm_plex_box_faces 2,2,2 -petscspace_degree 3 -qorder 3 -porder {{1 2 3 4}separate output} 684557cf195SMatthew G. Knepley test: 685d6837840SMatthew G. Knepley TODO: gll Lagrange nodes break this 686557cf195SMatthew G. Knepley suffix: p3_pragmatic 687557cf195SMatthew G. Knepley requires: triangle ctetgen pragmatic !single 68830602db0SMatthew G. Knepley args: -dm_plex_dim {{2}separate output} -dm_plex_box_faces 2,2,2 -petscspace_degree 3 -qorder 3 -dm_plex_hash_location -porder {{1 2 3 4}separate output} 689557cf195SMatthew G. Knepley 690557cf195SMatthew G. Knepley # 2D/3D Q_1 on a tensor cell 691557cf195SMatthew G. Knepley test: 692557cf195SMatthew G. Knepley suffix: q1 69330602db0SMatthew G. Knepley args: -dm_plex_dim {{2 3}separate output} -dm_plex_box_faces 2,2,2 -dm_plex_simplex 0 -petscspace_degree 1 -qorder 1 -porder {{1 2}separate output} 694557cf195SMatthew G. Knepley 695557cf195SMatthew G. Knepley # 2D/3D Q_2 on a tensor cell 696557cf195SMatthew G. Knepley test: 697557cf195SMatthew G. Knepley suffix: q2 69899a192c5SJunchao Zhang requires: !single 69930602db0SMatthew G. Knepley args: -dm_plex_dim {{2 3}separate output} -dm_plex_box_faces 2,2,2 -dm_plex_simplex 0 -petscspace_degree 2 -qorder 2 -porder {{1 2 3}separate output} 700557cf195SMatthew G. Knepley 701557cf195SMatthew G. Knepley # 2D/3D Q_3 on a tensor cell 702557cf195SMatthew G. Knepley test: 703d6837840SMatthew G. Knepley TODO: gll Lagrange nodes break this 704557cf195SMatthew G. Knepley suffix: q3 70599a192c5SJunchao Zhang requires: !single 70630602db0SMatthew G. Knepley args: -dm_plex_dim {{2 3}separate output} -dm_plex_box_faces 2,2,2 -dm_plex_simplex 0 -petscspace_degree 3 -qorder 3 -porder {{1 2 3 4}separate output} 707557cf195SMatthew G. Knepley 708557cf195SMatthew G. Knepley # 2D/3D P_1disc on a triangle/quadrilateral 709557cf195SMatthew G. Knepley # TODO Missing injection functional for simplices 710557cf195SMatthew G. Knepley test: 711557cf195SMatthew G. Knepley suffix: p1d 712557cf195SMatthew G. Knepley requires: triangle ctetgen 71330602db0SMatthew G. Knepley args: -dm_plex_dim {{2}separate output} -dm_plex_box_faces 2,2,2 -dm_plex_simplex {{0}separate output} -petscspace_degree 1 -petscdualspace_lagrange_continuity 0 -qorder 1 -porder {{1 2}separate output} 714557cf195SMatthew G. Knepley 715557cf195SMatthew G. Knepley TEST*/ 716