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