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