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