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