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