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