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