1 #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/ 2 #include <petscblaslapack.h> 3 4 static PetscErrorCode PetscFEDestroy_Basic(PetscFE fem) 5 { 6 PetscFE_Basic *b = (PetscFE_Basic *)fem->data; 7 8 PetscFunctionBegin; 9 PetscCall(PetscFree(b)); 10 PetscFunctionReturn(PETSC_SUCCESS); 11 } 12 13 static PetscErrorCode PetscFEView_Basic_Ascii(PetscFE fe, PetscViewer v) 14 { 15 PetscInt dim, Nc; 16 PetscSpace basis = NULL; 17 PetscDualSpace dual = NULL; 18 PetscQuadrature quad = NULL; 19 20 PetscFunctionBegin; 21 PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 22 PetscCall(PetscFEGetNumComponents(fe, &Nc)); 23 PetscCall(PetscFEGetBasisSpace(fe, &basis)); 24 PetscCall(PetscFEGetDualSpace(fe, &dual)); 25 PetscCall(PetscFEGetQuadrature(fe, &quad)); 26 PetscCall(PetscViewerASCIIPushTab(v)); 27 PetscCall(PetscViewerASCIIPrintf(v, "Basic Finite Element in %" PetscInt_FMT " dimensions with %" PetscInt_FMT " components\n", dim, Nc)); 28 if (basis) PetscCall(PetscSpaceView(basis, v)); 29 if (dual) PetscCall(PetscDualSpaceView(dual, v)); 30 if (quad) PetscCall(PetscQuadratureView(quad, v)); 31 PetscCall(PetscViewerASCIIPopTab(v)); 32 PetscFunctionReturn(PETSC_SUCCESS); 33 } 34 35 static PetscErrorCode PetscFEView_Basic(PetscFE fe, PetscViewer v) 36 { 37 PetscBool iascii; 38 39 PetscFunctionBegin; 40 PetscCall(PetscObjectTypeCompare((PetscObject)v, PETSCVIEWERASCII, &iascii)); 41 if (iascii) PetscCall(PetscFEView_Basic_Ascii(fe, v)); 42 PetscFunctionReturn(PETSC_SUCCESS); 43 } 44 45 /* Construct the change of basis from prime basis to nodal basis */ 46 PETSC_INTERN PetscErrorCode PetscFESetUp_Basic(PetscFE fem) 47 { 48 PetscReal *work; 49 PetscBLASInt *pivots; 50 PetscBLASInt n, info; 51 PetscInt pdim, j; 52 53 PetscFunctionBegin; 54 PetscCall(PetscDualSpaceGetDimension(fem->dualSpace, &pdim)); 55 PetscCall(PetscMalloc1(pdim * pdim, &fem->invV)); 56 for (j = 0; j < pdim; ++j) { 57 PetscReal *Bf; 58 PetscQuadrature f; 59 const PetscReal *points, *weights; 60 PetscInt Nc, Nq, q, k, c; 61 62 PetscCall(PetscDualSpaceGetFunctional(fem->dualSpace, j, &f)); 63 PetscCall(PetscQuadratureGetData(f, NULL, &Nc, &Nq, &points, &weights)); 64 PetscCall(PetscMalloc1(Nc * Nq * pdim, &Bf)); 65 PetscCall(PetscSpaceEvaluate(fem->basisSpace, Nq, points, Bf, NULL, NULL)); 66 for (k = 0; k < pdim; ++k) { 67 /* V_{jk} = n_j(\phi_k) = \int \phi_k(x) n_j(x) dx */ 68 fem->invV[j * pdim + k] = 0.0; 69 70 for (q = 0; q < Nq; ++q) { 71 for (c = 0; c < Nc; ++c) fem->invV[j * pdim + k] += Bf[(q * pdim + k) * Nc + c] * weights[q * Nc + c]; 72 } 73 } 74 PetscCall(PetscFree(Bf)); 75 } 76 77 PetscCall(PetscMalloc2(pdim, &pivots, pdim, &work)); 78 PetscCall(PetscBLASIntCast(pdim, &n)); 79 PetscCallBLAS("LAPACKgetrf", LAPACKREALgetrf_(&n, &n, fem->invV, &n, pivots, &info)); 80 PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error returned from LAPACKgetrf %" PetscBLASInt_FMT, info); 81 PetscCallBLAS("LAPACKgetri", LAPACKREALgetri_(&n, fem->invV, &n, pivots, work, &n, &info)); 82 PetscCheck(!info, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error returned from LAPACKgetri %" PetscBLASInt_FMT, info); 83 PetscCall(PetscFree2(pivots, work)); 84 PetscFunctionReturn(PETSC_SUCCESS); 85 } 86 87 PetscErrorCode PetscFEGetDimension_Basic(PetscFE fem, PetscInt *dim) 88 { 89 PetscFunctionBegin; 90 PetscCall(PetscDualSpaceGetDimension(fem->dualSpace, dim)); 91 PetscFunctionReturn(PETSC_SUCCESS); 92 } 93 94 /* Tensor contraction on the middle index, 95 * C[m,n,p] = A[m,k,p] * B[k,n] 96 * where all matrices use C-style ordering. 97 */ 98 static PetscErrorCode TensorContract_Private(PetscInt m, PetscInt n, PetscInt p, PetscInt k, const PetscReal *A, const PetscReal *B, PetscReal *C) 99 { 100 PetscInt i; 101 102 PetscFunctionBegin; 103 PetscCheck(n && p, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Empty tensor is not allowed %" PetscInt_FMT " %" PetscInt_FMT, n, p); 104 for (i = 0; i < m; i++) { 105 PetscBLASInt n_, p_, k_, lda, ldb, ldc; 106 PetscReal one = 1, zero = 0; 107 /* Taking contiguous submatrices, we wish to comput c[n,p] = a[k,p] * B[k,n] 108 * or, in Fortran ordering, c(p,n) = a(p,k) * B(n,k) 109 */ 110 PetscCall(PetscBLASIntCast(n, &n_)); 111 PetscCall(PetscBLASIntCast(p, &p_)); 112 PetscCall(PetscBLASIntCast(k, &k_)); 113 lda = p_; 114 ldb = n_; 115 ldc = p_; 116 PetscCallBLAS("BLASgemm", BLASREALgemm_("N", "T", &p_, &n_, &k_, &one, A + i * k * p, &lda, B, &ldb, &zero, C + i * n * p, &ldc)); 117 } 118 PetscCall(PetscLogFlops(2. * m * n * p * k)); 119 PetscFunctionReturn(PETSC_SUCCESS); 120 } 121 122 PETSC_INTERN PetscErrorCode PetscFEComputeTabulation_Basic(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation T) 123 { 124 DM dm; 125 PetscInt pdim; /* Dimension of FE space P */ 126 PetscInt dim; /* Spatial dimension */ 127 PetscInt Nc; /* Field components */ 128 PetscReal *B = K >= 0 ? T->T[0] : NULL; 129 PetscReal *D = K >= 1 ? T->T[1] : NULL; 130 PetscReal *H = K >= 2 ? T->T[2] : NULL; 131 PetscReal *tmpB = NULL, *tmpD = NULL, *tmpH = NULL; 132 133 PetscFunctionBegin; 134 PetscCall(PetscDualSpaceGetDM(fem->dualSpace, &dm)); 135 PetscCall(DMGetDimension(dm, &dim)); 136 PetscCall(PetscDualSpaceGetDimension(fem->dualSpace, &pdim)); 137 PetscCall(PetscFEGetNumComponents(fem, &Nc)); 138 /* Evaluate the prime basis functions at all points */ 139 if (K >= 0) PetscCall(DMGetWorkArray(dm, npoints * pdim * Nc, MPIU_REAL, &tmpB)); 140 if (K >= 1) PetscCall(DMGetWorkArray(dm, npoints * pdim * Nc * dim, MPIU_REAL, &tmpD)); 141 if (K >= 2) PetscCall(DMGetWorkArray(dm, npoints * pdim * Nc * dim * dim, MPIU_REAL, &tmpH)); 142 PetscCall(PetscSpaceEvaluate(fem->basisSpace, npoints, points, tmpB, tmpD, tmpH)); 143 /* Translate from prime to nodal basis */ 144 if (B) { 145 /* B[npoints, nodes, Nc] = tmpB[npoints, prime, Nc] * invV[prime, nodes] */ 146 PetscCall(TensorContract_Private(npoints, pdim, Nc, pdim, tmpB, fem->invV, B)); 147 } 148 if (D && dim) { 149 /* D[npoints, nodes, Nc, dim] = tmpD[npoints, prime, Nc, dim] * invV[prime, nodes] */ 150 PetscCall(TensorContract_Private(npoints, pdim, Nc * dim, pdim, tmpD, fem->invV, D)); 151 } 152 if (H && dim) { 153 /* H[npoints, nodes, Nc, dim, dim] = tmpH[npoints, prime, Nc, dim, dim] * invV[prime, nodes] */ 154 PetscCall(TensorContract_Private(npoints, pdim, Nc * dim * dim, pdim, tmpH, fem->invV, H)); 155 } 156 if (K >= 0) PetscCall(DMRestoreWorkArray(dm, npoints * pdim * Nc, MPIU_REAL, &tmpB)); 157 if (K >= 1) PetscCall(DMRestoreWorkArray(dm, npoints * pdim * Nc * dim, MPIU_REAL, &tmpD)); 158 if (K >= 2) PetscCall(DMRestoreWorkArray(dm, npoints * pdim * Nc * dim * dim, MPIU_REAL, &tmpH)); 159 PetscFunctionReturn(PETSC_SUCCESS); 160 } 161 162 PETSC_INTERN PetscErrorCode PetscFEIntegrate_Basic(PetscDS ds, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom, const PetscScalar coefficients[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscScalar integral[]) 163 { 164 const PetscInt debug = ds->printIntegrate; 165 PetscFE fe; 166 PetscPointFunc obj_func; 167 PetscQuadrature quad; 168 PetscTabulation *T, *TAux = NULL; 169 PetscScalar *u, *u_x, *a, *a_x; 170 const PetscScalar *constants; 171 PetscReal *x, cellScale; 172 PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 173 PetscInt dim, dE, Np, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e; 174 PetscBool isAffine; 175 const PetscReal *quadPoints, *quadWeights; 176 PetscInt qNc, Nq, q; 177 178 PetscFunctionBegin; 179 PetscCall(PetscDSGetObjective(ds, field, &obj_func)); 180 if (!obj_func) PetscFunctionReturn(PETSC_SUCCESS); 181 PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe)); 182 PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 183 cellScale = (PetscReal)PetscPowInt(2, dim); 184 PetscCall(PetscFEGetQuadrature(fe, &quad)); 185 PetscCall(PetscDSGetNumFields(ds, &Nf)); 186 PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 187 PetscCall(PetscDSGetComponentOffsets(ds, &uOff)); 188 PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x)); 189 PetscCall(PetscDSGetTabulation(ds, &T)); 190 PetscCall(PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x)); 191 PetscCall(PetscDSGetWorkspace(ds, &x, NULL, NULL, NULL, NULL)); 192 PetscCall(PetscDSSetIntegrationParameters(ds, field, PETSC_DETERMINE)); 193 PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 194 if (dsAux) { 195 PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 196 PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 197 PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 198 PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 199 PetscCall(PetscDSGetTabulation(dsAux, &TAux)); 200 PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 201 PetscCheck(T[0]->Np == TAux[0]->Np, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %" PetscInt_FMT " != %" PetscInt_FMT " number of auxiliary tabulation points", T[0]->Np, TAux[0]->Np); 202 } 203 PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights)); 204 PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 205 Np = cgeom->numPoints; 206 dE = cgeom->dimEmbed; 207 isAffine = cgeom->isAffine; 208 for (e = 0; e < Ne; ++e) { 209 PetscFEGeom fegeom; 210 211 fegeom.dim = cgeom->dim; 212 fegeom.dimEmbed = cgeom->dimEmbed; 213 fegeom.xi = NULL; 214 if (isAffine) { 215 fegeom.v = x; 216 fegeom.xi = cgeom->xi; 217 fegeom.J = &cgeom->J[e * Np * dE * dE]; 218 fegeom.invJ = &cgeom->invJ[e * Np * dE * dE]; 219 fegeom.detJ = &cgeom->detJ[e * Np]; 220 } 221 for (q = 0; q < Nq; ++q) { 222 PetscScalar integrand = 0.; 223 PetscReal w; 224 225 if (isAffine) { 226 CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * dim], x); 227 } else { 228 fegeom.v = &cgeom->v[(e * Np + q) * dE]; 229 fegeom.J = &cgeom->J[(e * Np + q) * dE * dE]; 230 fegeom.invJ = &cgeom->invJ[(e * Np + q) * dE * dE]; 231 fegeom.detJ = &cgeom->detJ[e * Np + q]; 232 } 233 PetscCall(PetscDSSetCellParameters(ds, fegeom.detJ[0] * cellScale)); 234 w = fegeom.detJ[0] * quadWeights[q]; 235 if (debug > 1 && q < Np) { 236 PetscCall(PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", (double)fegeom.detJ[0])); 237 #if !defined(PETSC_USE_COMPLEX) 238 PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ)); 239 #endif 240 } 241 if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, " quad point %" PetscInt_FMT "\n", q)); 242 PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], NULL, u, u_x, NULL)); 243 if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL)); 244 obj_func(dim, Nf, NfAux, uOff, uOff_x, u, NULL, u_x, aOff, aOff_x, a, NULL, a_x, 0.0, fegeom.v, numConstants, constants, &integrand); 245 integrand *= w; 246 integral[e * Nf + field] += integrand; 247 } 248 if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, " Element Field %" PetscInt_FMT " integral: %g\n", Nf, (double)PetscRealPart(integral[e * Nf + field]))); 249 cOffset += totDim; 250 cOffsetAux += totDimAux; 251 } 252 PetscFunctionReturn(PETSC_SUCCESS); 253 } 254 255 PETSC_INTERN PetscErrorCode PetscFEIntegrateBd_Basic(PetscDS ds, PetscInt field, PetscBdPointFunc obj_func, PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscScalar integral[]) 256 { 257 const PetscInt debug = ds->printIntegrate; 258 PetscFE fe; 259 PetscQuadrature quad; 260 PetscTabulation *Tf, *TfAux = NULL; 261 PetscScalar *u, *u_x, *a, *a_x, *basisReal, *basisDerReal; 262 const PetscScalar *constants; 263 PetscReal *x, cellScale; 264 PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 265 PetscBool isAffine, auxOnBd; 266 const PetscReal *quadPoints, *quadWeights; 267 PetscInt qNc, Nq, q, Np, dE; 268 PetscInt dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e; 269 270 PetscFunctionBegin; 271 if (!obj_func) PetscFunctionReturn(PETSC_SUCCESS); 272 PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe)); 273 PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 274 cellScale = (PetscReal)PetscPowInt(2, dim); 275 PetscCall(PetscFEGetFaceQuadrature(fe, &quad)); 276 PetscCall(PetscDSGetNumFields(ds, &Nf)); 277 PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 278 PetscCall(PetscDSGetComponentOffsets(ds, &uOff)); 279 PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x)); 280 PetscCall(PetscDSGetEvaluationArrays(ds, &u, NULL, &u_x)); 281 PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL)); 282 PetscCall(PetscDSGetFaceTabulation(ds, &Tf)); 283 PetscCall(PetscDSSetIntegrationParameters(ds, field, PETSC_DETERMINE)); 284 PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 285 if (dsAux) { 286 PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux)); 287 PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 288 PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 289 PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 290 PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 291 PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 292 auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE; 293 if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux)); 294 else PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux)); 295 PetscCheck(Tf[0]->Np == TfAux[0]->Np, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %" PetscInt_FMT " != %" PetscInt_FMT " number of auxiliary tabulation points", Tf[0]->Np, TfAux[0]->Np); 296 } 297 PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights)); 298 PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 299 if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Field: %" PetscInt_FMT " Nface: %" PetscInt_FMT " Nq: %" PetscInt_FMT "\n", field, Ne, Nq)); 300 Np = fgeom->numPoints; 301 dE = fgeom->dimEmbed; 302 isAffine = fgeom->isAffine; 303 for (e = 0; e < Ne; ++e) { 304 PetscFEGeom fegeom, cgeom; 305 const PetscInt face = fgeom->face[e][0]; /* Local face number in cell */ 306 fegeom.n = NULL; 307 fegeom.v = NULL; 308 fegeom.xi = NULL; 309 fegeom.J = NULL; 310 fegeom.invJ = NULL; 311 fegeom.detJ = NULL; 312 fegeom.dim = fgeom->dim; 313 fegeom.dimEmbed = fgeom->dimEmbed; 314 cgeom.dim = fgeom->dim; 315 cgeom.dimEmbed = fgeom->dimEmbed; 316 if (isAffine) { 317 fegeom.v = x; 318 fegeom.xi = fgeom->xi; 319 fegeom.J = &fgeom->J[e * Np * dE * dE]; 320 fegeom.invJ = &fgeom->invJ[e * Np * dE * dE]; 321 fegeom.detJ = &fgeom->detJ[e * Np]; 322 fegeom.n = &fgeom->n[e * Np * dE]; 323 324 cgeom.J = &fgeom->suppJ[0][e * Np * dE * dE]; 325 cgeom.invJ = &fgeom->suppInvJ[0][e * Np * dE * dE]; 326 cgeom.detJ = &fgeom->suppDetJ[0][e * Np]; 327 } 328 for (q = 0; q < Nq; ++q) { 329 PetscScalar integrand = 0.; 330 PetscReal w; 331 332 if (isAffine) { 333 CoordinatesRefToReal(dE, dim - 1, fegeom.xi, &fgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * (dim - 1)], x); 334 } else { 335 fegeom.v = &fgeom->v[(e * Np + q) * dE]; 336 fegeom.J = &fgeom->J[(e * Np + q) * dE * dE]; 337 fegeom.invJ = &fgeom->invJ[(e * Np + q) * dE * dE]; 338 fegeom.detJ = &fgeom->detJ[e * Np + q]; 339 fegeom.n = &fgeom->n[(e * Np + q) * dE]; 340 341 cgeom.J = &fgeom->suppJ[0][(e * Np + q) * dE * dE]; 342 cgeom.invJ = &fgeom->suppInvJ[0][(e * Np + q) * dE * dE]; 343 cgeom.detJ = &fgeom->suppDetJ[0][e * Np + q]; 344 } 345 PetscCall(PetscDSSetCellParameters(ds, fegeom.detJ[0] * cellScale)); 346 w = fegeom.detJ[0] * quadWeights[q]; 347 if (debug > 1 && q < Np) { 348 PetscCall(PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", (double)fegeom.detJ[0])); 349 #ifndef PETSC_USE_COMPLEX 350 PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ)); 351 #endif 352 } 353 if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, " quad point %" PetscInt_FMT "\n", q)); 354 if (debug > 3) { 355 PetscCall(PetscPrintf(PETSC_COMM_SELF, " x_q (")); 356 for (PetscInt d = 0; d < dE; ++d) { 357 if (d) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", ")); 358 PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)fegeom.v[d])); 359 } 360 PetscCall(PetscPrintf(PETSC_COMM_SELF, ")\n")); 361 PetscCall(PetscPrintf(PETSC_COMM_SELF, " n_q (")); 362 for (PetscInt d = 0; d < dE; ++d) { 363 if (d) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", ")); 364 PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)fegeom.n[d])); 365 } 366 PetscCall(PetscPrintf(PETSC_COMM_SELF, ")\n")); 367 for (PetscInt f = 0; f < Nf; ++f) { 368 PetscCall(PetscPrintf(PETSC_COMM_SELF, " u_%" PetscInt_FMT " (", f)); 369 for (PetscInt c = 0; c < uOff[f + 1] - uOff[f]; ++c) { 370 if (c) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", ")); 371 PetscCall(PetscPrintf(PETSC_COMM_SELF, "%g", (double)PetscRealPart(u[uOff[f] + c]))); 372 } 373 PetscCall(PetscPrintf(PETSC_COMM_SELF, ")\n")); 374 } 375 } 376 PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], NULL, u, u_x, NULL)); 377 if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL)); 378 obj_func(dim, Nf, NfAux, uOff, uOff_x, u, NULL, u_x, aOff, aOff_x, a, NULL, a_x, 0.0, fegeom.v, fegeom.n, numConstants, constants, &integrand); 379 integrand *= w; 380 integral[e * Nf + field] += integrand; 381 if (debug > 1) PetscCall(PetscPrintf(PETSC_COMM_SELF, " int: %g tot: %g\n", (double)PetscRealPart(integrand), (double)PetscRealPart(integral[e * Nf + field]))); 382 } 383 cOffset += totDim; 384 cOffsetAux += totDimAux; 385 } 386 PetscFunctionReturn(PETSC_SUCCESS); 387 } 388 389 PetscErrorCode PetscFEIntegrateResidual_Basic(PetscDS ds, PetscFormKey key, PetscInt Ne, PetscFEGeom *cgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 390 { 391 const PetscInt debug = ds->printIntegrate; 392 const PetscInt field = key.field; 393 PetscFE fe; 394 PetscWeakForm wf; 395 PetscInt n0, n1, i; 396 PetscPointFunc *f0_func, *f1_func; 397 PetscQuadrature quad; 398 PetscTabulation *T, *TAux = NULL; 399 PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal; 400 const PetscScalar *constants; 401 PetscReal *x, cellScale; 402 PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 403 PetscInt dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e; 404 const PetscReal *quadPoints, *quadWeights; 405 PetscInt qdim, qNc, Nq, q, dE; 406 407 PetscFunctionBegin; 408 PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe)); 409 PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 410 cellScale = (PetscReal)PetscPowInt(2, dim); 411 PetscCall(PetscFEGetQuadrature(fe, &quad)); 412 PetscCall(PetscDSGetNumFields(ds, &Nf)); 413 PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 414 PetscCall(PetscDSGetComponentOffsets(ds, &uOff)); 415 PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x)); 416 PetscCall(PetscDSGetFieldOffset(ds, field, &fOffset)); 417 PetscCall(PetscDSGetWeakForm(ds, &wf)); 418 PetscCall(PetscWeakFormGetResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func)); 419 if (!n0 && !n1) PetscFunctionReturn(PETSC_SUCCESS); 420 PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x)); 421 PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL)); 422 PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL)); 423 PetscCall(PetscDSGetTabulation(ds, &T)); 424 PetscCall(PetscDSSetIntegrationParameters(ds, field, PETSC_DETERMINE)); 425 PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 426 if (dsAux) { 427 PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 428 PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 429 PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 430 PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 431 PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 432 PetscCall(PetscDSGetTabulation(dsAux, &TAux)); 433 PetscCheck(T[0]->Np == TAux[0]->Np, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %" PetscInt_FMT " != %" PetscInt_FMT " number of auxiliary tabulation points", T[0]->Np, TAux[0]->Np); 434 } 435 PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights)); 436 PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 437 dE = cgeom->dimEmbed; 438 PetscCheck(cgeom->dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", cgeom->dim, qdim); 439 for (e = 0; e < Ne; ++e) { 440 PetscFEGeom fegeom; 441 442 fegeom.v = x; /* workspace */ 443 PetscCall(PetscArrayzero(f0, Nq * T[field]->Nc)); 444 PetscCall(PetscArrayzero(f1, Nq * T[field]->Nc * dE)); 445 for (q = 0; q < Nq; ++q) { 446 PetscReal w; 447 PetscInt c, d; 448 449 PetscCall(PetscFEGeomGetPoint(cgeom, e, q, &quadPoints[q * cgeom->dim], &fegeom)); 450 PetscCall(PetscDSSetCellParameters(ds, fegeom.detJ[0] * cellScale)); 451 w = fegeom.detJ[0] * quadWeights[q]; 452 if (debug > 1 && q < cgeom->numPoints) { 453 PetscCall(PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", (double)fegeom.detJ[0])); 454 #if !defined(PETSC_USE_COMPLEX) 455 PetscCall(DMPrintCellMatrix(e, "invJ", dE, dE, fegeom.invJ)); 456 #endif 457 } 458 PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], PetscSafePointerPlusOffset(coefficients_t, cOffset), u, u_x, u_t)); 459 if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL)); 460 for (i = 0; i < n0; ++i) f0_func[i](dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, numConstants, constants, &f0[q * T[field]->Nc]); 461 for (c = 0; c < T[field]->Nc; ++c) f0[q * T[field]->Nc + c] *= w; 462 for (i = 0; i < n1; ++i) f1_func[i](dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, numConstants, constants, &f1[q * T[field]->Nc * dE]); 463 for (c = 0; c < T[field]->Nc; ++c) 464 for (d = 0; d < dE; ++d) f1[(q * T[field]->Nc + c) * dE + d] *= w; 465 if (debug) { 466 // LCOV_EXCL_START 467 PetscCall(PetscPrintf(PETSC_COMM_SELF, " quad point %" PetscInt_FMT " wt %g x:", q, (double)quadWeights[q])); 468 for (c = 0; c < dE; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)fegeom.v[c])); 469 PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 470 if (debug > 2) { 471 PetscCall(PetscPrintf(PETSC_COMM_SELF, " field %" PetscInt_FMT ":", field)); 472 for (c = 0; c < T[field]->Nc; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(u[uOff[field] + c]))); 473 PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 474 PetscCall(PetscPrintf(PETSC_COMM_SELF, " field der %" PetscInt_FMT ":", field)); 475 for (c = 0; c < T[field]->Nc * dE; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(u_x[uOff[field] + c]))); 476 PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 477 PetscCall(PetscPrintf(PETSC_COMM_SELF, " resid %" PetscInt_FMT ":", field)); 478 for (c = 0; c < T[field]->Nc; ++c) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(f0[q * T[field]->Nc + c]))); 479 PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 480 PetscCall(PetscPrintf(PETSC_COMM_SELF, " res der %" PetscInt_FMT ":", field)); 481 for (c = 0; c < T[field]->Nc; ++c) { 482 for (d = 0; d < dE; ++d) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %g", (double)PetscRealPart(f1[(q * T[field]->Nc + c) * dE + d]))); 483 } 484 PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 485 } 486 // LCOV_EXCL_STOP 487 } 488 } 489 PetscCall(PetscFEUpdateElementVec_Internal(fe, T[field], 0, basisReal, basisDerReal, e, cgeom, f0, f1, &elemVec[cOffset + fOffset])); 490 cOffset += totDim; 491 cOffsetAux += totDimAux; 492 } 493 PetscFunctionReturn(PETSC_SUCCESS); 494 } 495 496 PetscErrorCode PetscFEIntegrateBdResidual_Basic(PetscDS ds, PetscWeakForm wf, PetscFormKey key, PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 497 { 498 const PetscInt debug = ds->printIntegrate; 499 const PetscInt field = key.field; 500 PetscFE fe; 501 PetscInt n0, n1, i; 502 PetscBdPointFunc *f0_func, *f1_func; 503 PetscQuadrature quad; 504 PetscTabulation *Tf, *TfAux = NULL; 505 PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal; 506 const PetscScalar *constants; 507 PetscReal *x, cellScale; 508 PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 509 PetscInt dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NcI; 510 PetscBool auxOnBd = PETSC_FALSE; 511 const PetscReal *quadPoints, *quadWeights; 512 PetscInt qdim, qNc, Nq, q, dE; 513 514 PetscFunctionBegin; 515 PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe)); 516 PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 517 cellScale = (PetscReal)PetscPowInt(2, dim); 518 PetscCall(PetscFEGetFaceQuadrature(fe, &quad)); 519 PetscCall(PetscDSGetNumFields(ds, &Nf)); 520 PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 521 PetscCall(PetscDSGetComponentOffsets(ds, &uOff)); 522 PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x)); 523 PetscCall(PetscDSGetFieldOffset(ds, field, &fOffset)); 524 PetscCall(PetscWeakFormGetBdResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func)); 525 if (!n0 && !n1) PetscFunctionReturn(PETSC_SUCCESS); 526 PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x)); 527 PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL)); 528 PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL)); 529 PetscCall(PetscDSGetFaceTabulation(ds, &Tf)); 530 PetscCall(PetscDSSetIntegrationParameters(ds, field, PETSC_DETERMINE)); 531 PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 532 if (dsAux) { 533 PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux)); 534 PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 535 PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 536 PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 537 PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 538 PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 539 auxOnBd = dimAux < dim ? PETSC_TRUE : PETSC_FALSE; 540 if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux)); 541 else PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux)); 542 PetscCheck(Tf[0]->Np == TfAux[0]->Np, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %" PetscInt_FMT " != %" PetscInt_FMT " number of auxiliary tabulation points", Tf[0]->Np, TfAux[0]->Np); 543 } 544 NcI = Tf[field]->Nc; 545 PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights)); 546 PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 547 dE = fgeom->dimEmbed; 548 /* TODO FIX THIS */ 549 fgeom->dim = dim - 1; 550 PetscCheck(fgeom->dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", fgeom->dim, qdim); 551 for (e = 0; e < Ne; ++e) { 552 PetscFEGeom fegeom, cgeom; 553 const PetscInt face = fgeom->face[e][0]; 554 555 fegeom.v = x; /* Workspace */ 556 PetscCall(PetscArrayzero(f0, Nq * NcI)); 557 PetscCall(PetscArrayzero(f1, Nq * NcI * dE)); 558 for (q = 0; q < Nq; ++q) { 559 PetscReal w; 560 PetscInt c, d; 561 562 PetscCall(PetscFEGeomGetPoint(fgeom, e, q, &quadPoints[q * fgeom->dim], &fegeom)); 563 PetscCall(PetscFEGeomGetCellPoint(fgeom, e, q, &cgeom)); 564 PetscCall(PetscDSSetCellParameters(ds, fegeom.detJ[0] * cellScale)); 565 w = fegeom.detJ[0] * quadWeights[q]; 566 if (debug > 1) { 567 if ((fgeom->isAffine && q == 0) || (!fgeom->isAffine)) { 568 PetscCall(PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", (double)fegeom.detJ[0])); 569 #if !defined(PETSC_USE_COMPLEX) 570 PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ)); 571 PetscCall(DMPrintCellVector(e, "n", dim, fegeom.n)); 572 #endif 573 } 574 } 575 PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], PetscSafePointerPlusOffset(coefficients_t, cOffset), u, u_x, u_t)); 576 if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, auxOnBd ? 0 : face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL)); 577 for (i = 0; i < n0; ++i) f0_func[i](dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, fegeom.n, numConstants, constants, &f0[q * NcI]); 578 for (c = 0; c < NcI; ++c) f0[q * NcI + c] *= w; 579 for (i = 0; i < n1; ++i) f1_func[i](dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, fegeom.n, numConstants, constants, &f1[q * NcI * dE]); 580 for (c = 0; c < NcI; ++c) 581 for (d = 0; d < dE; ++d) f1[(q * NcI + c) * dE + d] *= w; 582 if (debug) { 583 PetscCall(PetscPrintf(PETSC_COMM_SELF, " elem %" PetscInt_FMT " quad point %" PetscInt_FMT "\n", e, q)); 584 for (c = 0; c < NcI; ++c) { 585 if (n0) PetscCall(PetscPrintf(PETSC_COMM_SELF, " f0[%" PetscInt_FMT "] %g\n", c, (double)PetscRealPart(f0[q * NcI + c]))); 586 if (n1) { 587 for (d = 0; d < dim; ++d) PetscCall(PetscPrintf(PETSC_COMM_SELF, " f1[%" PetscInt_FMT ",%" PetscInt_FMT "] %g", c, d, (double)PetscRealPart(f1[(q * NcI + c) * dim + d]))); 588 PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 589 } 590 } 591 } 592 } 593 PetscCall(PetscFEUpdateElementVec_Internal(fe, Tf[field], face, basisReal, basisDerReal, e, fgeom, f0, f1, &elemVec[cOffset + fOffset])); 594 cOffset += totDim; 595 cOffsetAux += totDimAux; 596 } 597 PetscFunctionReturn(PETSC_SUCCESS); 598 } 599 600 /* 601 BdIntegral: Operates completely in the embedding dimension. The trick is to have special "face quadrature" so we only integrate over the face, but 602 all transforms operate in the full space and are square. 603 604 HybridIntegral: The discretization is lower dimensional. That means the transforms are non-square. 605 1) DMPlexGetCellFields() retrieves from the hybrid cell, so it gets fields from both faces 606 2) We need to assume that the orientation is 0 for both 607 3) TODO We need to use a non-square Jacobian for the derivative maps, meaning the embedding dimension has to go to EvaluateFieldJets() and UpdateElementVec() 608 */ 609 PETSC_INTERN PetscErrorCode PetscFEIntegrateHybridResidual_Basic(PetscDS ds, PetscDS dsIn, PetscFormKey key, PetscInt s, PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 610 { 611 const PetscInt debug = ds->printIntegrate; 612 const PetscInt field = key.field; 613 PetscFE fe; 614 PetscWeakForm wf; 615 PetscInt n0, n1, i; 616 PetscBdPointFunc *f0_func, *f1_func; 617 PetscQuadrature quad; 618 DMPolytopeType ct; 619 PetscTabulation *Tf, *TfIn, *TfAux = NULL; 620 PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal; 621 const PetscScalar *constants; 622 PetscReal *x; 623 PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 624 PetscInt dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimIn, totDimAux = 0, cOffset = 0, cOffsetIn = 0, cOffsetAux = 0, fOffset, e, NcI, NcS; 625 PetscBool isCohesiveField, auxOnBd = PETSC_FALSE; 626 const PetscReal *quadPoints, *quadWeights; 627 PetscInt qdim, qNc, Nq, q, dE; 628 629 PetscFunctionBegin; 630 /* Hybrid discretization is posed directly on faces */ 631 PetscCall(PetscDSGetDiscretization(ds, field, (PetscObject *)&fe)); 632 PetscCall(PetscFEGetSpatialDimension(fe, &dim)); 633 PetscCall(PetscFEGetQuadrature(fe, &quad)); 634 PetscCall(PetscDSGetNumFields(ds, &Nf)); 635 PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 636 PetscCall(PetscDSGetTotalDimension(dsIn, &totDimIn)); 637 PetscCall(PetscDSGetComponentOffsetsCohesive(dsIn, 0, &uOff)); // Change 0 to s for one-sided offsets 638 PetscCall(PetscDSGetComponentDerivativeOffsetsCohesive(dsIn, s, &uOff_x)); 639 PetscCall(PetscDSGetFieldOffsetCohesive(ds, field, &fOffset)); 640 PetscCall(PetscDSGetWeakForm(ds, &wf)); 641 PetscCall(PetscWeakFormGetBdResidual(wf, key.label, key.value, key.field, key.part, &n0, &f0_func, &n1, &f1_func)); 642 if (!n0 && !n1) PetscFunctionReturn(PETSC_SUCCESS); 643 PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x)); 644 PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL)); 645 PetscCall(PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL)); 646 /* NOTE This is a bulk tabulation because the DS is a face discretization */ 647 PetscCall(PetscDSGetTabulation(ds, &Tf)); 648 PetscCall(PetscDSGetFaceTabulation(dsIn, &TfIn)); 649 PetscCall(PetscDSSetIntegrationParameters(ds, field, PETSC_DETERMINE)); 650 PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 651 if (dsAux) { 652 PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux)); 653 PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 654 PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 655 PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 656 PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 657 PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 658 auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE; 659 if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TfAux)); 660 else PetscCall(PetscDSGetFaceTabulation(dsAux, &TfAux)); 661 PetscCheck(Tf[0]->Np == TfAux[0]->Np, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %" PetscInt_FMT " != %" PetscInt_FMT " number of auxiliary tabulation points", Tf[0]->Np, TfAux[0]->Np); 662 } 663 PetscCall(PetscDSGetCohesive(ds, field, &isCohesiveField)); 664 NcI = Tf[field]->Nc; 665 NcS = NcI; 666 if (!isCohesiveField && s == 2) { 667 // If we are integrating over a cohesive cell (s = 2) for a non-cohesive fields, we use both sides 668 NcS *= 2; 669 } 670 PetscCall(PetscQuadratureGetData(quad, &qdim, &qNc, &Nq, &quadPoints, &quadWeights)); 671 PetscCall(PetscQuadratureGetCellType(quad, &ct)); 672 PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 673 dE = fgeom->dimEmbed; 674 PetscCheck(fgeom->dim == qdim, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "FEGeom dim %" PetscInt_FMT " != %" PetscInt_FMT " quadrature dim", fgeom->dim, qdim); 675 for (e = 0; e < Ne; ++e) { 676 PetscFEGeom fegeom; 677 const PetscInt face[2] = {fgeom->face[e * 2 + 0][0], fgeom->face[e * 2 + 1][2]}; 678 const PetscInt ornt[2] = {fgeom->face[e * 2 + 0][1], fgeom->face[e * 2 + 1][3]}; 679 const PetscInt cornt[2] = {fgeom->face[e * 2 + 0][3], fgeom->face[e * 2 + 1][1]}; 680 681 fegeom.v = x; /* Workspace */ 682 PetscCall(PetscArrayzero(f0, Nq * NcS)); 683 PetscCall(PetscArrayzero(f1, Nq * NcS * dE)); 684 for (q = 0; q < Nq; ++q) { 685 PetscInt qpt[2]; 686 PetscReal w; 687 PetscInt c, d; 688 689 PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, cornt[0], ornt[0]), field, q, &qpt[0])); 690 PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, ornt[1], cornt[1]), field, q, &qpt[1])); 691 PetscCall(PetscFEGeomGetPoint(fgeom, e * 2, q, &quadPoints[q * fgeom->dim], &fegeom)); 692 w = fegeom.detJ[0] * quadWeights[q]; 693 if (debug > 1 && q < fgeom->numPoints) { 694 PetscCall(PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", (double)fegeom.detJ[0])); 695 #if !defined(PETSC_USE_COMPLEX) 696 PetscCall(DMPrintCellMatrix(e, "invJ", dim, dE, fegeom.invJ)); 697 #endif 698 } 699 if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, " quad point %" PetscInt_FMT " weight %g detJ %g\n", q, (double)quadWeights[q], (double)fegeom.detJ[0])); 700 /* TODO Is this cell or face quadrature, meaning should we use 'q' or 'face*Nq+q' */ 701 PetscCall(PetscFEEvaluateFieldJets_Hybrid_Internal(ds, Nf, 0, q, Tf, face, qpt, TfIn, &fegeom, &coefficients[cOffsetIn], PetscSafePointerPlusOffset(coefficients_t, cOffsetIn), u, u_x, u_t)); 702 if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, auxOnBd ? 0 : face[s], auxOnBd ? q : qpt[s], TfAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL)); 703 for (i = 0; i < n0; ++i) f0_func[i](dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, fegeom.n, numConstants, constants, &f0[q * NcS]); 704 for (c = 0; c < NcS; ++c) f0[q * NcS + c] *= w; 705 for (i = 0; i < n1; ++i) f1_func[i](dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, fegeom.v, fegeom.n, numConstants, constants, &f1[q * NcS * dE]); 706 for (c = 0; c < NcS; ++c) 707 for (d = 0; d < dE; ++d) f1[(q * NcS + c) * dE + d] *= w; 708 } 709 if (isCohesiveField) { 710 PetscCall(PetscFEUpdateElementVec_Internal(fe, Tf[field], 0, basisReal, basisDerReal, e, fgeom, f0, f1, &elemVec[cOffset + fOffset])); 711 } else { 712 PetscCall(PetscFEUpdateElementVec_Hybrid_Internal(fe, Tf[field], 0, s, basisReal, basisDerReal, fgeom, f0, f1, &elemVec[cOffset + fOffset])); 713 } 714 cOffset += totDim; 715 cOffsetIn += totDimIn; 716 cOffsetAux += totDimAux; 717 } 718 PetscFunctionReturn(PETSC_SUCCESS); 719 } 720 721 PetscErrorCode PetscFEIntegrateJacobian_Basic(PetscDS ds, PetscFEJacobianType jtype, PetscFormKey key, PetscInt Ne, PetscFEGeom *cgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 722 { 723 const PetscInt debug = ds->printIntegrate; 724 PetscFE feI, feJ; 725 PetscWeakForm wf; 726 PetscPointJac *g0_func, *g1_func, *g2_func, *g3_func; 727 PetscInt n0, n1, n2, n3, i; 728 PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 729 PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 730 PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 731 PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 732 PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 733 PetscQuadrature quad; 734 PetscTabulation *T, *TAux = NULL; 735 PetscScalar *g0 = NULL, *g1 = NULL, *g2 = NULL, *g3 = NULL, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal; 736 const PetscScalar *constants; 737 PetscReal *x, cellScale; 738 PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 739 PetscInt NcI = 0, NcJ = 0; 740 PetscInt dim, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e; 741 PetscInt dE, Np; 742 PetscBool isAffine; 743 const PetscReal *quadPoints, *quadWeights; 744 PetscInt qNc, Nq, q; 745 746 PetscFunctionBegin; 747 PetscCall(PetscDSGetNumFields(ds, &Nf)); 748 fieldI = key.field / Nf; 749 fieldJ = key.field % Nf; 750 PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&feI)); 751 PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *)&feJ)); 752 PetscCall(PetscFEGetSpatialDimension(feI, &dim)); 753 cellScale = (PetscReal)PetscPowInt(2, dim); 754 PetscCall(PetscFEGetQuadrature(feI, &quad)); 755 PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 756 PetscCall(PetscDSGetComponentOffsets(ds, &uOff)); 757 PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x)); 758 PetscCall(PetscDSGetWeakForm(ds, &wf)); 759 switch (jtype) { 760 case PETSCFE_JACOBIAN_DYN: 761 PetscCall(PetscWeakFormGetDynamicJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); 762 break; 763 case PETSCFE_JACOBIAN_PRE: 764 PetscCall(PetscWeakFormGetJacobianPreconditioner(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); 765 break; 766 case PETSCFE_JACOBIAN: 767 PetscCall(PetscWeakFormGetJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); 768 break; 769 } 770 if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(PETSC_SUCCESS); 771 PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x)); 772 PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal)); 773 PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, n0 ? &g0 : NULL, n1 ? &g1 : NULL, n2 ? &g2 : NULL, n3 ? &g3 : NULL)); 774 775 PetscCall(PetscDSGetTabulation(ds, &T)); 776 PetscCall(PetscDSGetFieldOffset(ds, fieldI, &offsetI)); 777 PetscCall(PetscDSGetFieldOffset(ds, fieldJ, &offsetJ)); 778 PetscCall(PetscDSSetIntegrationParameters(ds, fieldI, fieldJ)); 779 PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 780 if (dsAux) { 781 PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 782 PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 783 PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 784 PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 785 PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 786 PetscCall(PetscDSGetTabulation(dsAux, &TAux)); 787 PetscCheck(T[0]->Np == TAux[0]->Np, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %" PetscInt_FMT " != %" PetscInt_FMT " number of auxiliary tabulation points", T[0]->Np, TAux[0]->Np); 788 } 789 NcI = T[fieldI]->Nc; 790 NcJ = T[fieldJ]->Nc; 791 Np = cgeom->numPoints; 792 dE = cgeom->dimEmbed; 793 isAffine = cgeom->isAffine; 794 PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights)); 795 PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 796 797 for (e = 0; e < Ne; ++e) { 798 PetscFEGeom fegeom; 799 800 fegeom.dim = cgeom->dim; 801 fegeom.dimEmbed = cgeom->dimEmbed; 802 fegeom.xi = NULL; 803 if (isAffine) { 804 fegeom.v = x; 805 fegeom.xi = cgeom->xi; 806 fegeom.J = &cgeom->J[e * Np * dE * dE]; 807 fegeom.invJ = &cgeom->invJ[e * Np * dE * dE]; 808 fegeom.detJ = &cgeom->detJ[e * Np]; 809 } 810 for (q = 0; q < Nq; ++q) { 811 PetscReal w; 812 PetscInt c; 813 814 if (isAffine) { 815 CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * dim], x); 816 } else { 817 fegeom.v = &cgeom->v[(e * Np + q) * dE]; 818 fegeom.J = &cgeom->J[(e * Np + q) * dE * dE]; 819 fegeom.invJ = &cgeom->invJ[(e * Np + q) * dE * dE]; 820 fegeom.detJ = &cgeom->detJ[e * Np + q]; 821 } 822 PetscCall(PetscDSSetCellParameters(ds, fegeom.detJ[0] * cellScale)); 823 if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, " quad point %" PetscInt_FMT " weight %g detJ %g\n", q, (double)quadWeights[q], (double)fegeom.detJ[0])); 824 w = fegeom.detJ[0] * quadWeights[q]; 825 if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], PetscSafePointerPlusOffset(coefficients_t, cOffset), u, u_x, u_t)); 826 if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL)); 827 if (n0) { 828 PetscCall(PetscArrayzero(g0, NcI * NcJ)); 829 for (i = 0; i < n0; ++i) g0_func[i](dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, numConstants, constants, g0); 830 for (c = 0; c < NcI * NcJ; ++c) g0[c] *= w; 831 } 832 if (n1) { 833 PetscCall(PetscArrayzero(g1, NcI * NcJ * dE)); 834 for (i = 0; i < n1; ++i) g1_func[i](dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, numConstants, constants, g1); 835 for (c = 0; c < NcI * NcJ * dE; ++c) g1[c] *= w; 836 } 837 if (n2) { 838 PetscCall(PetscArrayzero(g2, NcI * NcJ * dE)); 839 for (i = 0; i < n2; ++i) g2_func[i](dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, numConstants, constants, g2); 840 for (c = 0; c < NcI * NcJ * dE; ++c) g2[c] *= w; 841 } 842 if (n3) { 843 PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE)); 844 for (i = 0; i < n3; ++i) g3_func[i](dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, numConstants, constants, g3); 845 for (c = 0; c < NcI * NcJ * dE * dE; ++c) g3[c] *= w; 846 } 847 848 PetscCall(PetscFEUpdateElementMat_Internal(feI, feJ, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, totDim, offsetI, offsetJ, elemMat + eOffset)); 849 } 850 if (debug > 1) { 851 PetscInt f, g; 852 853 PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %" PetscInt_FMT " and %" PetscInt_FMT "\n", fieldI, fieldJ)); 854 for (f = 0; f < T[fieldI]->Nb; ++f) { 855 const PetscInt i = offsetI + f; 856 for (g = 0; g < T[fieldJ]->Nb; ++g) { 857 const PetscInt j = offsetJ + g; 858 PetscCall(PetscPrintf(PETSC_COMM_SELF, " elemMat[%" PetscInt_FMT ", %" PetscInt_FMT "]: %g\n", f, g, (double)PetscRealPart(elemMat[eOffset + i * totDim + j]))); 859 } 860 PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 861 } 862 } 863 cOffset += totDim; 864 cOffsetAux += totDimAux; 865 eOffset += PetscSqr(totDim); 866 } 867 PetscFunctionReturn(PETSC_SUCCESS); 868 } 869 870 PETSC_INTERN PetscErrorCode PetscFEIntegrateBdJacobian_Basic(PetscDS ds, PetscWeakForm wf, PetscFEJacobianType jtype, PetscFormKey key, PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 871 { 872 const PetscInt debug = ds->printIntegrate; 873 PetscFE feI, feJ; 874 PetscBdPointJac *g0_func, *g1_func, *g2_func, *g3_func; 875 PetscInt n0, n1, n2, n3, i; 876 PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 877 PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 878 PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 879 PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 880 PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 881 PetscQuadrature quad; 882 PetscTabulation *T, *TAux = NULL; 883 PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal; 884 const PetscScalar *constants; 885 PetscReal *x, cellScale; 886 PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 887 PetscInt NcI = 0, NcJ = 0; 888 PetscInt dim, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e; 889 PetscBool isAffine; 890 const PetscReal *quadPoints, *quadWeights; 891 PetscInt qNc, Nq, q, Np, dE; 892 893 PetscFunctionBegin; 894 PetscCall(PetscDSGetNumFields(ds, &Nf)); 895 fieldI = key.field / Nf; 896 fieldJ = key.field % Nf; 897 PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&feI)); 898 PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *)&feJ)); 899 PetscCall(PetscFEGetSpatialDimension(feI, &dim)); 900 cellScale = (PetscReal)PetscPowInt(2, dim); 901 PetscCall(PetscFEGetFaceQuadrature(feI, &quad)); 902 PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 903 PetscCall(PetscDSGetComponentOffsets(ds, &uOff)); 904 PetscCall(PetscDSGetComponentDerivativeOffsets(ds, &uOff_x)); 905 PetscCall(PetscDSGetFieldOffset(ds, fieldI, &offsetI)); 906 PetscCall(PetscDSGetFieldOffset(ds, fieldJ, &offsetJ)); 907 switch (jtype) { 908 case PETSCFE_JACOBIAN_PRE: 909 PetscCall(PetscWeakFormGetBdJacobianPreconditioner(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); 910 break; 911 case PETSCFE_JACOBIAN: 912 PetscCall(PetscWeakFormGetBdJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); 913 break; 914 case PETSCFE_JACOBIAN_DYN: 915 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "PETSCFE_JACOBIAN_DYN is not supported for PetscFEIntegrateBdJacobian()"); 916 } 917 if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(PETSC_SUCCESS); 918 PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x)); 919 PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal)); 920 PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3)); 921 PetscCall(PetscDSGetFaceTabulation(ds, &T)); 922 PetscCall(PetscDSSetIntegrationParameters(ds, fieldI, fieldJ)); 923 PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 924 if (dsAux) { 925 PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 926 PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 927 PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 928 PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 929 PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 930 PetscCall(PetscDSGetFaceTabulation(dsAux, &TAux)); 931 } 932 NcI = T[fieldI]->Nc, NcJ = T[fieldJ]->Nc; 933 Np = fgeom->numPoints; 934 dE = fgeom->dimEmbed; 935 isAffine = fgeom->isAffine; 936 /* Initialize here in case the function is not defined */ 937 PetscCall(PetscArrayzero(g0, NcI * NcJ)); 938 PetscCall(PetscArrayzero(g1, NcI * NcJ * dE)); 939 PetscCall(PetscArrayzero(g2, NcI * NcJ * dE)); 940 PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE)); 941 PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights)); 942 PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 943 for (e = 0; e < Ne; ++e) { 944 PetscFEGeom fegeom, cgeom; 945 const PetscInt face = fgeom->face[e][0]; 946 fegeom.n = NULL; 947 fegeom.v = NULL; 948 fegeom.xi = NULL; 949 fegeom.J = NULL; 950 fegeom.detJ = NULL; 951 fegeom.dim = fgeom->dim; 952 fegeom.dimEmbed = fgeom->dimEmbed; 953 cgeom.dim = fgeom->dim; 954 cgeom.dimEmbed = fgeom->dimEmbed; 955 if (isAffine) { 956 fegeom.v = x; 957 fegeom.xi = fgeom->xi; 958 fegeom.J = &fgeom->J[e * Np * dE * dE]; 959 fegeom.invJ = &fgeom->invJ[e * Np * dE * dE]; 960 fegeom.detJ = &fgeom->detJ[e * Np]; 961 fegeom.n = &fgeom->n[e * Np * dE]; 962 963 cgeom.J = &fgeom->suppJ[0][e * Np * dE * dE]; 964 cgeom.invJ = &fgeom->suppInvJ[0][e * Np * dE * dE]; 965 cgeom.detJ = &fgeom->suppDetJ[0][e * Np]; 966 } 967 for (q = 0; q < Nq; ++q) { 968 PetscReal w; 969 PetscInt c; 970 971 if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, " quad point %" PetscInt_FMT "\n", q)); 972 if (isAffine) { 973 CoordinatesRefToReal(dE, dim - 1, fegeom.xi, &fgeom->v[e * Np * dE], fegeom.J, &quadPoints[q * (dim - 1)], x); 974 } else { 975 fegeom.v = &fgeom->v[(e * Np + q) * dE]; 976 fegeom.J = &fgeom->J[(e * Np + q) * dE * dE]; 977 fegeom.invJ = &fgeom->invJ[(e * Np + q) * dE * dE]; 978 fegeom.detJ = &fgeom->detJ[e * Np + q]; 979 fegeom.n = &fgeom->n[(e * Np + q) * dE]; 980 981 cgeom.J = &fgeom->suppJ[0][(e * Np + q) * dE * dE]; 982 cgeom.invJ = &fgeom->suppInvJ[0][(e * Np + q) * dE * dE]; 983 cgeom.detJ = &fgeom->suppDetJ[0][e * Np + q]; 984 } 985 PetscCall(PetscDSSetCellParameters(ds, fegeom.detJ[0] * cellScale)); 986 w = fegeom.detJ[0] * quadWeights[q]; 987 if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, T, &cgeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t)); 988 if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL)); 989 if (n0) { 990 PetscCall(PetscArrayzero(g0, NcI * NcJ)); 991 for (i = 0; i < n0; ++i) g0_func[i](dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, fegeom.n, numConstants, constants, g0); 992 for (c = 0; c < NcI * NcJ; ++c) g0[c] *= w; 993 } 994 if (n1) { 995 PetscCall(PetscArrayzero(g1, NcI * NcJ * dE)); 996 for (i = 0; i < n1; ++i) g1_func[i](dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, fegeom.n, numConstants, constants, g1); 997 for (c = 0; c < NcI * NcJ * dim; ++c) g1[c] *= w; 998 } 999 if (n2) { 1000 PetscCall(PetscArrayzero(g2, NcI * NcJ * dE)); 1001 for (i = 0; i < n2; ++i) g2_func[i](dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, fegeom.n, numConstants, constants, g2); 1002 for (c = 0; c < NcI * NcJ * dim; ++c) g2[c] *= w; 1003 } 1004 if (n3) { 1005 PetscCall(PetscArrayzero(g3, NcI * NcJ * dE * dE)); 1006 for (i = 0; i < n3; ++i) g3_func[i](dE, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, fegeom.n, numConstants, constants, g3); 1007 for (c = 0; c < NcI * NcJ * dim * dim; ++c) g3[c] *= w; 1008 } 1009 1010 PetscCall(PetscFEUpdateElementMat_Internal(feI, feJ, face, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &cgeom, g0, g1, g2, g3, totDim, offsetI, offsetJ, elemMat + eOffset)); 1011 } 1012 if (debug > 1) { 1013 PetscInt fc, f, gc, g; 1014 1015 PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %" PetscInt_FMT " and %" PetscInt_FMT "\n", fieldI, fieldJ)); 1016 for (fc = 0; fc < T[fieldI]->Nc; ++fc) { 1017 for (f = 0; f < T[fieldI]->Nb; ++f) { 1018 const PetscInt i = offsetI + f * T[fieldI]->Nc + fc; 1019 for (gc = 0; gc < T[fieldJ]->Nc; ++gc) { 1020 for (g = 0; g < T[fieldJ]->Nb; ++g) { 1021 const PetscInt j = offsetJ + g * T[fieldJ]->Nc + gc; 1022 PetscCall(PetscPrintf(PETSC_COMM_SELF, " elemMat[%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT "]: %g\n", f, fc, g, gc, (double)PetscRealPart(elemMat[eOffset + i * totDim + j]))); 1023 } 1024 } 1025 PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 1026 } 1027 } 1028 } 1029 cOffset += totDim; 1030 cOffsetAux += totDimAux; 1031 eOffset += PetscSqr(totDim); 1032 } 1033 PetscFunctionReturn(PETSC_SUCCESS); 1034 } 1035 1036 PETSC_INTERN PetscErrorCode PetscFEIntegrateHybridJacobian_Basic(PetscDS ds, PetscDS dsIn, PetscFEJacobianType jtype, PetscFormKey key, PetscInt s, PetscInt Ne, PetscFEGeom *fgeom, const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 1037 { 1038 const PetscInt debug = ds->printIntegrate; 1039 PetscFE feI, feJ; 1040 PetscWeakForm wf; 1041 PetscBdPointJac *g0_func, *g1_func, *g2_func, *g3_func; 1042 PetscInt n0, n1, n2, n3, i; 1043 PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 1044 PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 1045 PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 1046 PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 1047 PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 1048 PetscQuadrature quad; 1049 DMPolytopeType ct; 1050 PetscTabulation *T, *TfIn, *TAux = NULL; 1051 PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal; 1052 const PetscScalar *constants; 1053 PetscReal *x; 1054 PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 1055 PetscInt NcI = 0, NcJ = 0, NcS, NcT; 1056 PetscInt dim, dimAux, numConstants, Nf, fieldI, fieldJ, NfAux = 0, totDim, totDimAux = 0, e; 1057 PetscBool isCohesiveFieldI, isCohesiveFieldJ, auxOnBd = PETSC_FALSE; 1058 const PetscReal *quadPoints, *quadWeights; 1059 PetscInt qNc, Nq, q; 1060 1061 PetscFunctionBegin; 1062 PetscCall(PetscDSGetNumFields(ds, &Nf)); 1063 fieldI = key.field / Nf; 1064 fieldJ = key.field % Nf; 1065 /* Hybrid discretization is posed directly on faces */ 1066 PetscCall(PetscDSGetDiscretization(ds, fieldI, (PetscObject *)&feI)); 1067 PetscCall(PetscDSGetDiscretization(ds, fieldJ, (PetscObject *)&feJ)); 1068 PetscCall(PetscFEGetSpatialDimension(feI, &dim)); 1069 PetscCall(PetscFEGetQuadrature(feI, &quad)); 1070 PetscCall(PetscDSGetTotalDimension(ds, &totDim)); 1071 PetscCall(PetscDSGetComponentOffsetsCohesive(ds, 0, &uOff)); // Change 0 to s for one-sided offsets 1072 PetscCall(PetscDSGetComponentDerivativeOffsetsCohesive(ds, s, &uOff_x)); 1073 PetscCall(PetscDSGetWeakForm(ds, &wf)); 1074 switch (jtype) { 1075 case PETSCFE_JACOBIAN_PRE: 1076 PetscCall(PetscWeakFormGetBdJacobianPreconditioner(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); 1077 break; 1078 case PETSCFE_JACOBIAN: 1079 PetscCall(PetscWeakFormGetBdJacobian(wf, key.label, key.value, fieldI, fieldJ, key.part, &n0, &g0_func, &n1, &g1_func, &n2, &g2_func, &n3, &g3_func)); 1080 break; 1081 case PETSCFE_JACOBIAN_DYN: 1082 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No boundary hybrid Jacobians :)"); 1083 } 1084 if (!n0 && !n1 && !n2 && !n3) PetscFunctionReturn(PETSC_SUCCESS); 1085 PetscCall(PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x)); 1086 PetscCall(PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal)); 1087 PetscCall(PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3)); 1088 PetscCall(PetscDSGetTabulation(ds, &T)); 1089 PetscCall(PetscDSGetFaceTabulation(dsIn, &TfIn)); 1090 PetscCall(PetscDSGetFieldOffsetCohesive(ds, fieldI, &offsetI)); 1091 PetscCall(PetscDSGetFieldOffsetCohesive(ds, fieldJ, &offsetJ)); 1092 PetscCall(PetscDSSetIntegrationParameters(ds, fieldI, fieldJ)); 1093 PetscCall(PetscDSGetConstants(ds, &numConstants, &constants)); 1094 if (dsAux) { 1095 PetscCall(PetscDSGetSpatialDimension(dsAux, &dimAux)); 1096 PetscCall(PetscDSGetNumFields(dsAux, &NfAux)); 1097 PetscCall(PetscDSGetTotalDimension(dsAux, &totDimAux)); 1098 PetscCall(PetscDSGetComponentOffsets(dsAux, &aOff)); 1099 PetscCall(PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x)); 1100 PetscCall(PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x)); 1101 auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE; 1102 if (auxOnBd) PetscCall(PetscDSGetTabulation(dsAux, &TAux)); 1103 else PetscCall(PetscDSGetFaceTabulation(dsAux, &TAux)); 1104 PetscCheck(T[0]->Np == TAux[0]->Np, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of tabulation points %" PetscInt_FMT " != %" PetscInt_FMT " number of auxiliary tabulation points", T[0]->Np, TAux[0]->Np); 1105 } 1106 PetscCall(PetscDSGetCohesive(ds, fieldI, &isCohesiveFieldI)); 1107 PetscCall(PetscDSGetCohesive(ds, fieldJ, &isCohesiveFieldJ)); 1108 NcI = T[fieldI]->Nc; 1109 NcJ = T[fieldJ]->Nc; 1110 NcS = isCohesiveFieldI ? NcI : 2 * NcI; 1111 NcT = isCohesiveFieldJ ? NcJ : 2 * NcJ; 1112 if (!isCohesiveFieldI && s == 2) { 1113 // If we are integrating over a cohesive cell (s = 2) for a non-cohesive fields, we use both sides 1114 NcS *= 2; 1115 } 1116 if (!isCohesiveFieldJ && s == 2) { 1117 // If we are integrating over a cohesive cell (s = 2) for a non-cohesive fields, we use both sides 1118 NcT *= 2; 1119 } 1120 // The derivatives are constrained to be along the cell, so there are dim, not dE, components, even though 1121 // the coordinates are in dE dimensions 1122 PetscCall(PetscArrayzero(g0, NcS * NcT)); 1123 PetscCall(PetscArrayzero(g1, NcS * NcT * dim)); 1124 PetscCall(PetscArrayzero(g2, NcS * NcT * dim)); 1125 PetscCall(PetscArrayzero(g3, NcS * NcT * dim * dim)); 1126 PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights)); 1127 PetscCall(PetscQuadratureGetCellType(quad, &ct)); 1128 PetscCheck(qNc == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %" PetscInt_FMT " components", qNc); 1129 for (e = 0; e < Ne; ++e) { 1130 PetscFEGeom fegeom; 1131 const PetscInt face[2] = {fgeom->face[e * 2 + 0][0], fgeom->face[e * 2 + 1][2]}; 1132 const PetscInt ornt[2] = {fgeom->face[e * 2 + 0][1], fgeom->face[e * 2 + 1][3]}; 1133 const PetscInt cornt[2] = {fgeom->face[e * 2 + 0][3], fgeom->face[e * 2 + 1][1]}; 1134 1135 fegeom.v = x; /* Workspace */ 1136 for (q = 0; q < Nq; ++q) { 1137 PetscInt qpt[2]; 1138 PetscReal w; 1139 PetscInt c; 1140 1141 PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, cornt[0], ornt[0]), fieldI, q, &qpt[0])); 1142 PetscCall(PetscDSPermuteQuadPoint(ds, DMPolytopeTypeComposeOrientationInv(ct, ornt[1], cornt[1]), fieldI, q, &qpt[1])); 1143 PetscCall(PetscFEGeomGetPoint(fgeom, e * 2, q, &quadPoints[q * fgeom->dim], &fegeom)); 1144 w = fegeom.detJ[0] * quadWeights[q]; 1145 if (debug > 1 && q < fgeom->numPoints) { 1146 PetscCall(PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", (double)fegeom.detJ[0])); 1147 #if !defined(PETSC_USE_COMPLEX) 1148 PetscCall(DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ)); 1149 #endif 1150 } 1151 if (debug) PetscCall(PetscPrintf(PETSC_COMM_SELF, " quad point %" PetscInt_FMT "\n", q)); 1152 if (coefficients) PetscCall(PetscFEEvaluateFieldJets_Hybrid_Internal(ds, Nf, 0, q, T, face, qpt, TfIn, &fegeom, &coefficients[cOffset], PetscSafePointerPlusOffset(coefficients_t, cOffset), u, u_x, u_t)); 1153 if (dsAux) PetscCall(PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, auxOnBd ? 0 : face[s], auxOnBd ? q : qpt[s], TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL)); 1154 if (n0) { 1155 PetscCall(PetscArrayzero(g0, NcS * NcT)); 1156 for (i = 0; i < n0; ++i) g0_func[i](dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, fegeom.n, numConstants, constants, g0); 1157 for (c = 0; c < NcS * NcT; ++c) g0[c] *= w; 1158 } 1159 if (n1) { 1160 PetscCall(PetscArrayzero(g1, NcS * NcT * dim)); 1161 for (i = 0; i < n1; ++i) g1_func[i](dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, fegeom.n, numConstants, constants, g1); 1162 for (c = 0; c < NcS * NcT * dim; ++c) g1[c] *= w; 1163 } 1164 if (n2) { 1165 PetscCall(PetscArrayzero(g2, NcS * NcT * dim)); 1166 for (i = 0; i < n2; ++i) g2_func[i](dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, fegeom.n, numConstants, constants, g2); 1167 for (c = 0; c < NcS * NcT * dim; ++c) g2[c] *= w; 1168 } 1169 if (n3) { 1170 PetscCall(PetscArrayzero(g3, NcS * NcT * dim * dim)); 1171 for (i = 0; i < n3; ++i) g3_func[i](dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, fegeom.v, fegeom.n, numConstants, constants, g3); 1172 for (c = 0; c < NcS * NcT * dim * dim; ++c) g3[c] *= w; 1173 } 1174 1175 if (isCohesiveFieldI) { 1176 if (isCohesiveFieldJ) { 1177 PetscCall(PetscFEUpdateElementMat_Internal(feI, feJ, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, totDim, offsetI, offsetJ, elemMat + eOffset)); 1178 } else { 1179 PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 0, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat)); 1180 PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 1, 1, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, &g0[NcI * NcJ], &g1[NcI * NcJ * dim], &g2[NcI * NcJ * dim], &g3[NcI * NcJ * dim * dim], eOffset, totDim, offsetI, offsetJ, elemMat)); 1181 } 1182 } else { 1183 if (s == 2) { 1184 if (isCohesiveFieldJ) { 1185 PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 0, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat)); 1186 PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 1, 1, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, &g0[NcI * NcJ], &g1[NcI * NcJ * dim], &g2[NcI * NcJ * dim], &g3[NcI * NcJ * dim * dim], eOffset, totDim, offsetI, offsetJ, elemMat)); 1187 } else { 1188 PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 0, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat)); 1189 PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 0, 1, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, &g0[NcI * NcJ], &g1[NcI * NcJ * dim], &g2[NcI * NcJ * dim], &g3[NcI * NcJ * dim * dim], eOffset, totDim, offsetI, offsetJ, elemMat)); 1190 PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 1, 0, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, &g0[NcI * NcJ * 2], &g1[NcI * NcJ * dim * 2], &g2[NcI * NcJ * dim * 2], &g3[NcI * NcJ * dim * dim * 2], eOffset, totDim, offsetI, offsetJ, elemMat)); 1191 PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, 1, 1, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, &g0[NcI * NcJ * 3], &g1[NcI * NcJ * dim * 3], &g2[NcI * NcJ * dim * 3], &g3[NcI * NcJ * dim * dim * 3], eOffset, totDim, offsetI, offsetJ, elemMat)); 1192 } 1193 } else 1194 PetscCall(PetscFEUpdateElementMat_Hybrid_Internal(feI, isCohesiveFieldI, feJ, isCohesiveFieldJ, 0, s, s, q, T[fieldI], basisReal, basisDerReal, T[fieldJ], testReal, testDerReal, &fegeom, g0, g1, g2, g3, eOffset, totDim, offsetI, offsetJ, elemMat)); 1195 } 1196 } 1197 if (debug > 1) { 1198 const PetscInt fS = 0 + (isCohesiveFieldI ? 0 : (s == 2 ? 0 : s * T[fieldI]->Nb)); 1199 const PetscInt fE = T[fieldI]->Nb + (isCohesiveFieldI ? 0 : (s == 2 ? T[fieldI]->Nb : s * T[fieldI]->Nb)); 1200 const PetscInt gS = 0 + (isCohesiveFieldJ ? 0 : (s == 2 ? 0 : s * T[fieldJ]->Nb)); 1201 const PetscInt gE = T[fieldJ]->Nb + (isCohesiveFieldJ ? 0 : (s == 2 ? T[fieldJ]->Nb : s * T[fieldJ]->Nb)); 1202 PetscInt f, g; 1203 1204 PetscCall(PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %" PetscInt_FMT " and %" PetscInt_FMT " s %s totDim %" PetscInt_FMT " offsets (%" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT ")\n", fieldI, fieldJ, s ? (s > 1 ? "Coh" : "Pos") : "Neg", totDim, eOffset, offsetI, offsetJ)); 1205 for (f = fS; f < fE; ++f) { 1206 const PetscInt i = offsetI + f; 1207 for (g = gS; g < gE; ++g) { 1208 const PetscInt j = offsetJ + g; 1209 PetscCheck(i < totDim && j < totDim, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Fuck up %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT, f, i, g, j); 1210 PetscCall(PetscPrintf(PETSC_COMM_SELF, " elemMat[%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT "]: %g\n", f / NcI, f % NcI, g / NcJ, g % NcJ, (double)PetscRealPart(elemMat[eOffset + i * totDim + j]))); 1211 } 1212 PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n")); 1213 } 1214 } 1215 cOffset += totDim; 1216 cOffsetAux += totDimAux; 1217 eOffset += PetscSqr(totDim); 1218 } 1219 PetscFunctionReturn(PETSC_SUCCESS); 1220 } 1221 1222 static PetscErrorCode PetscFEInitialize_Basic(PetscFE fem) 1223 { 1224 PetscFunctionBegin; 1225 fem->ops->setfromoptions = NULL; 1226 fem->ops->setup = PetscFESetUp_Basic; 1227 fem->ops->view = PetscFEView_Basic; 1228 fem->ops->destroy = PetscFEDestroy_Basic; 1229 fem->ops->getdimension = PetscFEGetDimension_Basic; 1230 fem->ops->computetabulation = PetscFEComputeTabulation_Basic; 1231 fem->ops->integrate = PetscFEIntegrate_Basic; 1232 fem->ops->integratebd = PetscFEIntegrateBd_Basic; 1233 fem->ops->integrateresidual = PetscFEIntegrateResidual_Basic; 1234 fem->ops->integratebdresidual = PetscFEIntegrateBdResidual_Basic; 1235 fem->ops->integratehybridresidual = PetscFEIntegrateHybridResidual_Basic; 1236 fem->ops->integratejacobianaction = NULL /* PetscFEIntegrateJacobianAction_Basic */; 1237 fem->ops->integratejacobian = PetscFEIntegrateJacobian_Basic; 1238 fem->ops->integratebdjacobian = PetscFEIntegrateBdJacobian_Basic; 1239 fem->ops->integratehybridjacobian = PetscFEIntegrateHybridJacobian_Basic; 1240 PetscFunctionReturn(PETSC_SUCCESS); 1241 } 1242 1243 /*MC 1244 PETSCFEBASIC = "basic" - A `PetscFE` object that integrates with basic tiling and no vectorization 1245 1246 Level: intermediate 1247 1248 .seealso: `PetscFE`, `PetscFEType`, `PetscFECreate()`, `PetscFESetType()` 1249 M*/ 1250 1251 PETSC_EXTERN PetscErrorCode PetscFECreate_Basic(PetscFE fem) 1252 { 1253 PetscFE_Basic *b; 1254 1255 PetscFunctionBegin; 1256 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1257 PetscCall(PetscNew(&b)); 1258 fem->data = b; 1259 1260 PetscCall(PetscFEInitialize_Basic(fem)); 1261 PetscFunctionReturn(PETSC_SUCCESS); 1262 } 1263