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, PetscInt field, 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 PetscFE fe; 374 PetscPointFunc f0_func; 375 PetscPointFunc f1_func; 376 PetscQuadrature quad; 377 PetscTabulation *T, *TAux = NULL; 378 PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal; 379 const PetscScalar *constants; 380 PetscReal *x; 381 PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 382 PetscInt dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e; 383 PetscBool isAffine; 384 const PetscReal *quadPoints, *quadWeights; 385 PetscInt qNc, Nq, q, Np, dE; 386 PetscErrorCode ierr; 387 388 PetscFunctionBegin; 389 ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr); 390 ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr); 391 ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 392 ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 393 ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 394 ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 395 ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 396 ierr = PetscDSGetFieldOffset(ds, field, &fOffset);CHKERRQ(ierr); 397 ierr = PetscDSGetResidual(ds, field, &f0_func, &f1_func);CHKERRQ(ierr); 398 ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 399 ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL);CHKERRQ(ierr); 400 ierr = PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 401 if (!f0_func && !f1_func) PetscFunctionReturn(0); 402 ierr = PetscDSGetTabulation(ds, &T);CHKERRQ(ierr); 403 ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 404 if (dsAux) { 405 ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 406 ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 407 ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 408 ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 409 ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 410 ierr = PetscDSGetTabulation(dsAux, &TAux);CHKERRQ(ierr); 411 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); 412 } 413 ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 414 if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 415 Np = cgeom->numPoints; 416 dE = cgeom->dimEmbed; 417 isAffine = cgeom->isAffine; 418 for (e = 0; e < Ne; ++e) { 419 PetscFEGeom fegeom; 420 421 fegeom.dim = cgeom->dim; 422 fegeom.dimEmbed = cgeom->dimEmbed; 423 if (isAffine) { 424 fegeom.v = x; 425 fegeom.xi = cgeom->xi; 426 fegeom.J = &cgeom->J[e*Np*dE*dE]; 427 fegeom.invJ = &cgeom->invJ[e*Np*dE*dE]; 428 fegeom.detJ = &cgeom->detJ[e*Np]; 429 } 430 ierr = PetscArrayzero(f0, Nq*T[field]->Nc);CHKERRQ(ierr); 431 ierr = PetscArrayzero(f1, Nq*T[field]->Nc*dE);CHKERRQ(ierr); 432 for (q = 0; q < Nq; ++q) { 433 PetscReal w; 434 PetscInt c, d; 435 436 if (isAffine) { 437 CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e*Np*dE], fegeom.J, &quadPoints[q*dim], x); 438 } else { 439 fegeom.v = &cgeom->v[(e*Np+q)*dE]; 440 fegeom.J = &cgeom->J[(e*Np+q)*dE*dE]; 441 fegeom.invJ = &cgeom->invJ[(e*Np+q)*dE*dE]; 442 fegeom.detJ = &cgeom->detJ[e*Np+q]; 443 } 444 w = fegeom.detJ[0]*quadWeights[q]; 445 if (debug > 1 && q < Np) { 446 ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr); 447 #if !defined(PETSC_USE_COMPLEX) 448 ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr); 449 #endif 450 } 451 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 452 ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr); 453 if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);} 454 if (f0_func) { 455 f0_func(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]); 456 for (c = 0; c < T[field]->Nc; ++c) f0[q*T[field]->Nc+c] *= w; 457 } 458 if (f1_func) { 459 f1_func(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 } 463 ierr = PetscFEUpdateElementVec_Internal(fe, T[field], 0, basisReal, basisDerReal, &fegeom, f0, f1, &elemVec[cOffset+fOffset]);CHKERRQ(ierr); 464 cOffset += totDim; 465 cOffsetAux += totDimAux; 466 } 467 PetscFunctionReturn(0); 468 } 469 470 PetscErrorCode PetscFEIntegrateBdResidual_Basic(PetscDS ds, PetscInt field, PetscInt Ne, PetscFEGeom *fgeom, 471 const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 472 { 473 const PetscInt debug = 0; 474 PetscFE fe; 475 PetscBdPointFunc f0_func; 476 PetscBdPointFunc 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 = PetscDSGetBdResidual(ds, field, &f0_func, &f1_func);CHKERRQ(ierr); 499 if (!f0_func && !f1_func) 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 && q < Np) { 567 ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr); 568 #if !defined(PETSC_USE_COMPLEX) 569 ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr); 570 #endif 571 } 572 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 573 ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, Tf, &cgeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr); 574 if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, auxOnBd ? 0 : face, q, TfAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);} 575 if (f0_func) { 576 f0_func(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]); 577 for (c = 0; c < NcI; ++c) f0[q*NcI+c] *= w; 578 } 579 if (f1_func) { 580 f1_func(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]); 581 for (c = 0; c < NcI; ++c) for (d = 0; d < dim; ++d) f1[(q*NcI+c)*dim+d] *= w; 582 } 583 } 584 ierr = PetscFEUpdateElementVec_Internal(fe, Tf[field], face, basisReal, basisDerReal, &cgeom, f0, f1, &elemVec[cOffset+fOffset]);CHKERRQ(ierr); 585 cOffset += totDim; 586 cOffsetAux += totDimAux; 587 } 588 PetscFunctionReturn(0); 589 } 590 591 /* 592 BdIntegral: Operates completely in the embedding dimension. The trick is to have special "face quadrature" so we only integrate over the face, but 593 all transforms operate in the full space and are square. 594 595 HybridIntegral: The discretization is lower dimensional. That means the transforms are non-square. 596 1) DMPlexGetCellFields() retrieves from the hybrid cell, so it gets fields from both faces 597 2) We need to assume that the orientation is 0 for both 598 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() 599 */ 600 static PetscErrorCode PetscFEIntegrateHybridResidual_Basic(PetscDS ds, PetscInt field, PetscInt Ne, PetscFEGeom *fgeom, 601 const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 602 { 603 const PetscInt debug = 0; 604 PetscFE fe; 605 PetscBdPointFunc f0_func; 606 PetscBdPointFunc f1_func; 607 PetscQuadrature quad; 608 PetscTabulation *Tf, *TfAux = NULL; 609 PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal; 610 const PetscScalar *constants; 611 PetscReal *x; 612 PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 613 PetscInt dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NcI, NcS; 614 PetscBool isCohesiveField, isAffine, auxOnBd = PETSC_FALSE; 615 const PetscReal *quadPoints, *quadWeights; 616 PetscInt qNc, Nq, q, Np, dE; 617 PetscErrorCode ierr; 618 619 PetscFunctionBegin; 620 /* Hybrid discretization is posed directly on faces */ 621 ierr = PetscDSGetDiscretization(ds, field, (PetscObject *) &fe);CHKERRQ(ierr); 622 ierr = PetscFEGetSpatialDimension(fe, &dim);CHKERRQ(ierr); 623 ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 624 ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 625 ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 626 ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 627 ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 628 ierr = PetscDSGetFieldOffset(ds, field, &fOffset);CHKERRQ(ierr); 629 ierr = PetscDSGetBdResidual(ds, field, &f0_func, &f1_func);CHKERRQ(ierr); 630 if (!f0_func && !f1_func) PetscFunctionReturn(0); 631 ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 632 ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, NULL, NULL);CHKERRQ(ierr); 633 ierr = PetscDSGetWeakFormArrays(ds, &f0, &f1, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 634 /* NOTE This is a bulk tabulation because the DS is a face discretization */ 635 ierr = PetscDSGetTabulation(ds, &Tf);CHKERRQ(ierr); 636 ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 637 if (dsAux) { 638 ierr = PetscDSGetSpatialDimension(dsAux, &dimAux);CHKERRQ(ierr); 639 ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 640 ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 641 ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 642 ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 643 ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 644 auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE; 645 if (auxOnBd) {ierr = PetscDSGetTabulation(dsAux, &TfAux);CHKERRQ(ierr);} 646 else {ierr = PetscDSGetFaceTabulation(dsAux, &TfAux);CHKERRQ(ierr);} 647 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); 648 } 649 isCohesiveField = field == Nf-1 ? PETSC_TRUE : PETSC_FALSE; 650 NcI = Tf[field]->Nc; 651 NcS = isCohesiveField ? NcI : 2*NcI; 652 ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 653 if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 654 Np = fgeom->numPoints; 655 dE = fgeom->dimEmbed; 656 isAffine = fgeom->isAffine; 657 for (e = 0; e < Ne; ++e) { 658 PetscFEGeom fegeom; 659 const PetscInt face = fgeom->face[e][0]; 660 661 fegeom.dim = fgeom->dim; 662 fegeom.dimEmbed = fgeom->dimEmbed; 663 if (isAffine) { 664 fegeom.v = x; 665 fegeom.xi = fgeom->xi; 666 fegeom.J = &fgeom->J[e*dE*dE]; 667 fegeom.invJ = &fgeom->invJ[e*dE*dE]; 668 fegeom.detJ = &fgeom->detJ[e]; 669 fegeom.n = &fgeom->n[e*dE]; 670 } 671 ierr = PetscArrayzero(f0, Nq*NcS);CHKERRQ(ierr); 672 ierr = PetscArrayzero(f1, Nq*NcS*dE);CHKERRQ(ierr); 673 for (q = 0; q < Nq; ++q) { 674 PetscReal w; 675 PetscInt c, d; 676 677 if (isAffine) { 678 CoordinatesRefToReal(dE, dim, fegeom.xi, &fgeom->v[e*dE], fegeom.J, &quadPoints[q*dim], x); 679 } else { 680 fegeom.v = &fgeom->v[(e*Np+q)*dE]; 681 fegeom.J = &fgeom->J[(e*Np+q)*dE*dE]; 682 fegeom.invJ = &fgeom->invJ[(e*Np+q)*dE*dE]; 683 fegeom.detJ = &fgeom->detJ[e*Np+q]; 684 fegeom.n = &fgeom->n[(e*Np+q)*dE]; 685 } 686 w = fegeom.detJ[0]*quadWeights[q]; 687 if (debug > 1 && q < Np) { 688 ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr); 689 #if !defined(PETSC_USE_COMPLEX) 690 ierr = DMPrintCellMatrix(e, "invJ", dim, dE, fegeom.invJ);CHKERRQ(ierr); 691 #endif 692 } 693 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 694 /* TODO Is this cell or face quadrature, meaning should we use 'q' or 'face*Nq+q' */ 695 ierr = PetscFEEvaluateFieldJets_Hybrid_Internal(ds, Nf, 0, q, Tf, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr); 696 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);} 697 if (f0_func) { 698 f0_func(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]); 699 for (c = 0; c < NcS; ++c) f0[q*NcS+c] *= w; 700 } 701 if (f1_func) { 702 f1_func(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]); 703 for (c = 0; c < NcS; ++c) for (d = 0; d < dim; ++d) f1[(q*NcS+c)*dim+d] *= w; 704 } 705 } 706 if (isCohesiveField) {PetscFEUpdateElementVec_Internal(fe, Tf[field], 0, basisReal, basisDerReal, &fegeom, f0, f1, &elemVec[cOffset+fOffset*2]);} 707 else {PetscFEUpdateElementVec_Hybrid_Internal(fe, Tf[field], 0, basisReal, basisDerReal, &fegeom, f0, f1, &elemVec[cOffset+fOffset*2]);} 708 cOffset += totDim; 709 cOffsetAux += totDimAux; 710 } 711 PetscFunctionReturn(0); 712 } 713 714 PetscErrorCode PetscFEIntegrateJacobian_Basic(PetscDS ds, PetscFEJacobianType jtype, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *cgeom, 715 const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 716 { 717 const PetscInt debug = 0; 718 PetscFE feI, feJ; 719 PetscPointJac g0_func, g1_func, g2_func, g3_func; 720 PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 721 PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 722 PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 723 PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 724 PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 725 PetscQuadrature quad; 726 PetscTabulation *T, *TAux = NULL; 727 PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal; 728 const PetscScalar *constants; 729 PetscReal *x; 730 PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 731 PetscInt NcI = 0, NcJ = 0; 732 PetscInt dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, e; 733 PetscInt dE, Np; 734 PetscBool isAffine; 735 const PetscReal *quadPoints, *quadWeights; 736 PetscInt qNc, Nq, q; 737 PetscErrorCode ierr; 738 739 PetscFunctionBegin; 740 ierr = PetscDSGetDiscretization(ds, fieldI, (PetscObject *) &feI);CHKERRQ(ierr); 741 ierr = PetscDSGetDiscretization(ds, fieldJ, (PetscObject *) &feJ);CHKERRQ(ierr); 742 ierr = PetscFEGetSpatialDimension(feI, &dim);CHKERRQ(ierr); 743 ierr = PetscFEGetQuadrature(feI, &quad);CHKERRQ(ierr); 744 ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 745 ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 746 ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 747 ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 748 switch(jtype) { 749 case PETSCFE_JACOBIAN_DYN: ierr = PetscDSGetDynamicJacobian(ds, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);break; 750 case PETSCFE_JACOBIAN_PRE: ierr = PetscDSGetJacobianPreconditioner(ds, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);break; 751 case PETSCFE_JACOBIAN: ierr = PetscDSGetJacobian(ds, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);break; 752 } 753 if (!g0_func && !g1_func && !g2_func && !g3_func) PetscFunctionReturn(0); 754 ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 755 ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal);CHKERRQ(ierr); 756 ierr = PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3);CHKERRQ(ierr); 757 ierr = PetscDSGetTabulation(ds, &T);CHKERRQ(ierr); 758 ierr = PetscDSGetFieldOffset(ds, fieldI, &offsetI);CHKERRQ(ierr); 759 ierr = PetscDSGetFieldOffset(ds, fieldJ, &offsetJ);CHKERRQ(ierr); 760 ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 761 if (dsAux) { 762 ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 763 ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 764 ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 765 ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 766 ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 767 ierr = PetscDSGetTabulation(dsAux, &TAux);CHKERRQ(ierr); 768 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); 769 } 770 NcI = T[fieldI]->Nc; 771 NcJ = T[fieldJ]->Nc; 772 Np = cgeom->numPoints; 773 dE = cgeom->dimEmbed; 774 isAffine = cgeom->isAffine; 775 /* Initialize here in case the function is not defined */ 776 ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr); 777 ierr = PetscArrayzero(g1, NcI*NcJ*dE);CHKERRQ(ierr); 778 ierr = PetscArrayzero(g2, NcI*NcJ*dE);CHKERRQ(ierr); 779 ierr = PetscArrayzero(g3, NcI*NcJ*dE*dE);CHKERRQ(ierr); 780 ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 781 if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 782 for (e = 0; e < Ne; ++e) { 783 PetscFEGeom fegeom; 784 785 fegeom.dim = cgeom->dim; 786 fegeom.dimEmbed = cgeom->dimEmbed; 787 if (isAffine) { 788 fegeom.v = x; 789 fegeom.xi = cgeom->xi; 790 fegeom.J = &cgeom->J[e*Np*dE*dE]; 791 fegeom.invJ = &cgeom->invJ[e*Np*dE*dE]; 792 fegeom.detJ = &cgeom->detJ[e*Np]; 793 } 794 for (q = 0; q < Nq; ++q) { 795 PetscReal w; 796 PetscInt c; 797 798 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 799 if (isAffine) { 800 CoordinatesRefToReal(dE, dim, fegeom.xi, &cgeom->v[e*Np*dE], fegeom.J, &quadPoints[q*dim], x); 801 } else { 802 fegeom.v = &cgeom->v[(e*Np+q)*dE]; 803 fegeom.J = &cgeom->J[(e*Np+q)*dE*dE]; 804 fegeom.invJ = &cgeom->invJ[(e*Np+q)*dE*dE]; 805 fegeom.detJ = &cgeom->detJ[e*Np+q]; 806 } 807 w = fegeom.detJ[0]*quadWeights[q]; 808 if (coefficients) {ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr);} 809 if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, 0, q, TAux, &fegeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);} 810 if (g0_func) { 811 ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr); 812 g0_func(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); 813 for (c = 0; c < NcI*NcJ; ++c) g0[c] *= w; 814 } 815 if (g1_func) { 816 ierr = PetscArrayzero(g1, NcI*NcJ*dE);CHKERRQ(ierr); 817 g1_func(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); 818 for (c = 0; c < NcI*NcJ*dim; ++c) g1[c] *= w; 819 } 820 if (g2_func) { 821 ierr = PetscArrayzero(g2, NcI*NcJ*dE);CHKERRQ(ierr); 822 g2_func(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); 823 for (c = 0; c < NcI*NcJ*dim; ++c) g2[c] *= w; 824 } 825 if (g3_func) { 826 ierr = PetscArrayzero(g3, NcI*NcJ*dE*dE);CHKERRQ(ierr); 827 g3_func(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); 828 for (c = 0; c < NcI*NcJ*dim*dim; ++c) g3[c] *= w; 829 } 830 831 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); 832 } 833 if (debug > 1) { 834 PetscInt fc, f, gc, g; 835 836 ierr = PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);CHKERRQ(ierr); 837 for (fc = 0; fc < T[fieldI]->Nc; ++fc) { 838 for (f = 0; f < T[fieldI]->Nb; ++f) { 839 const PetscInt i = offsetI + f*T[fieldI]->Nc+fc; 840 for (gc = 0; gc < T[fieldJ]->Nc; ++gc) { 841 for (g = 0; g < T[fieldJ]->Nb; ++g) { 842 const PetscInt j = offsetJ + g*T[fieldJ]->Nc+gc; 843 ierr = PetscPrintf(PETSC_COMM_SELF, " elemMat[%d,%d,%d,%d]: %g\n", f, fc, g, gc, PetscRealPart(elemMat[eOffset+i*totDim+j]));CHKERRQ(ierr); 844 } 845 } 846 ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr); 847 } 848 } 849 } 850 cOffset += totDim; 851 cOffsetAux += totDimAux; 852 eOffset += PetscSqr(totDim); 853 } 854 PetscFunctionReturn(0); 855 } 856 857 static PetscErrorCode PetscFEIntegrateBdJacobian_Basic(PetscDS ds, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *fgeom, 858 const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 859 { 860 const PetscInt debug = 0; 861 PetscFE feI, feJ; 862 PetscBdPointJac g0_func, g1_func, g2_func, g3_func; 863 PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 864 PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 865 PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 866 PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 867 PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 868 PetscQuadrature quad; 869 PetscTabulation *T, *TAux = NULL; 870 PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal; 871 const PetscScalar *constants; 872 PetscReal *x; 873 PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 874 PetscInt NcI = 0, NcJ = 0; 875 PetscInt dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, e; 876 PetscBool isAffine; 877 const PetscReal *quadPoints, *quadWeights; 878 PetscInt qNc, Nq, q, Np, dE; 879 PetscErrorCode ierr; 880 881 PetscFunctionBegin; 882 ierr = PetscDSGetDiscretization(ds, fieldI, (PetscObject *) &feI);CHKERRQ(ierr); 883 ierr = PetscDSGetDiscretization(ds, fieldJ, (PetscObject *) &feJ);CHKERRQ(ierr); 884 ierr = PetscFEGetSpatialDimension(feI, &dim);CHKERRQ(ierr); 885 ierr = PetscFEGetFaceQuadrature(feI, &quad);CHKERRQ(ierr); 886 ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 887 ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 888 ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 889 ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 890 ierr = PetscDSGetFieldOffset(ds, fieldI, &offsetI);CHKERRQ(ierr); 891 ierr = PetscDSGetFieldOffset(ds, fieldJ, &offsetJ);CHKERRQ(ierr); 892 ierr = PetscDSGetBdJacobian(ds, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr); 893 if (!g0_func && !g1_func && !g2_func && !g3_func) PetscFunctionReturn(0); 894 ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 895 ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal);CHKERRQ(ierr); 896 ierr = PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3);CHKERRQ(ierr); 897 ierr = PetscDSGetFaceTabulation(ds, &T);CHKERRQ(ierr); 898 ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 899 if (dsAux) { 900 ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 901 ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 902 ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 903 ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 904 ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 905 ierr = PetscDSGetFaceTabulation(dsAux, &TAux);CHKERRQ(ierr); 906 } 907 NcI = T[fieldI]->Nc, NcJ = T[fieldJ]->Nc; 908 Np = fgeom->numPoints; 909 dE = fgeom->dimEmbed; 910 isAffine = fgeom->isAffine; 911 /* Initialize here in case the function is not defined */ 912 ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr); 913 ierr = PetscArrayzero(g1, NcI*NcJ*dE);CHKERRQ(ierr); 914 ierr = PetscArrayzero(g2, NcI*NcJ*dE);CHKERRQ(ierr); 915 ierr = PetscArrayzero(g3, NcI*NcJ*dE*dE);CHKERRQ(ierr); 916 ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 917 if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 918 for (e = 0; e < Ne; ++e) { 919 PetscFEGeom fegeom, cgeom; 920 const PetscInt face = fgeom->face[e][0]; 921 fegeom.n = NULL; 922 fegeom.v = NULL; 923 fegeom.J = NULL; 924 fegeom.detJ = NULL; 925 fegeom.dim = fgeom->dim; 926 fegeom.dimEmbed = fgeom->dimEmbed; 927 cgeom.dim = fgeom->dim; 928 cgeom.dimEmbed = fgeom->dimEmbed; 929 if (isAffine) { 930 fegeom.v = x; 931 fegeom.xi = fgeom->xi; 932 fegeom.J = &fgeom->J[e*Np*dE*dE]; 933 fegeom.invJ = &fgeom->invJ[e*Np*dE*dE]; 934 fegeom.detJ = &fgeom->detJ[e*Np]; 935 fegeom.n = &fgeom->n[e*Np*dE]; 936 937 cgeom.J = &fgeom->suppJ[0][e*Np*dE*dE]; 938 cgeom.invJ = &fgeom->suppInvJ[0][e*Np*dE*dE]; 939 cgeom.detJ = &fgeom->suppDetJ[0][e*Np]; 940 } 941 for (q = 0; q < Nq; ++q) { 942 PetscReal w; 943 PetscInt c; 944 945 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 946 if (isAffine) { 947 CoordinatesRefToReal(dE, dim-1, fegeom.xi, &fgeom->v[e*Np*dE], fegeom.J, &quadPoints[q*(dim-1)], x); 948 } else { 949 fegeom.v = &fgeom->v[(e*Np+q)*dE]; 950 fegeom.J = &fgeom->J[(e*Np+q)*dE*dE]; 951 fegeom.invJ = &fgeom->invJ[(e*Np+q)*dE*dE]; 952 fegeom.detJ = &fgeom->detJ[e*Np+q]; 953 fegeom.n = &fgeom->n[(e*Np+q)*dE]; 954 955 cgeom.J = &fgeom->suppJ[0][(e*Np+q)*dE*dE]; 956 cgeom.invJ = &fgeom->suppInvJ[0][(e*Np+q)*dE*dE]; 957 cgeom.detJ = &fgeom->suppDetJ[0][e*Np+q]; 958 } 959 w = fegeom.detJ[0]*quadWeights[q]; 960 if (coefficients) {ierr = PetscFEEvaluateFieldJets_Internal(ds, Nf, face, q, T, &cgeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr);} 961 if (dsAux) {ierr = PetscFEEvaluateFieldJets_Internal(dsAux, NfAux, face, q, TAux, &cgeom, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);CHKERRQ(ierr);} 962 if (g0_func) { 963 ierr = PetscArrayzero(g0, NcI*NcJ);CHKERRQ(ierr); 964 g0_func(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); 965 for (c = 0; c < NcI*NcJ; ++c) g0[c] *= w; 966 } 967 if (g1_func) { 968 ierr = PetscArrayzero(g1, NcI*NcJ*dE);CHKERRQ(ierr); 969 g1_func(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); 970 for (c = 0; c < NcI*NcJ*dim; ++c) g1[c] *= w; 971 } 972 if (g2_func) { 973 ierr = PetscArrayzero(g2, NcI*NcJ*dE);CHKERRQ(ierr); 974 g2_func(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); 975 for (c = 0; c < NcI*NcJ*dim; ++c) g2[c] *= w; 976 } 977 if (g3_func) { 978 ierr = PetscArrayzero(g3, NcI*NcJ*dE*dE);CHKERRQ(ierr); 979 g3_func(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); 980 for (c = 0; c < NcI*NcJ*dim*dim; ++c) g3[c] *= w; 981 } 982 983 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); 984 } 985 if (debug > 1) { 986 PetscInt fc, f, gc, g; 987 988 ierr = PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);CHKERRQ(ierr); 989 for (fc = 0; fc < T[fieldI]->Nc; ++fc) { 990 for (f = 0; f < T[fieldI]->Nb; ++f) { 991 const PetscInt i = offsetI + f*T[fieldI]->Nc+fc; 992 for (gc = 0; gc < T[fieldJ]->Nc; ++gc) { 993 for (g = 0; g < T[fieldJ]->Nb; ++g) { 994 const PetscInt j = offsetJ + g*T[fieldJ]->Nc+gc; 995 ierr = PetscPrintf(PETSC_COMM_SELF, " elemMat[%d,%d,%d,%d]: %g\n", f, fc, g, gc, PetscRealPart(elemMat[eOffset+i*totDim+j]));CHKERRQ(ierr); 996 } 997 } 998 ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr); 999 } 1000 } 1001 } 1002 cOffset += totDim; 1003 cOffsetAux += totDimAux; 1004 eOffset += PetscSqr(totDim); 1005 } 1006 PetscFunctionReturn(0); 1007 } 1008 1009 PetscErrorCode PetscFEIntegrateHybridJacobian_Basic(PetscDS ds, PetscFEJacobianType jtype, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *fgeom, 1010 const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS dsAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 1011 { 1012 const PetscInt debug = 0; 1013 PetscFE feI, feJ; 1014 PetscBdPointJac g0_func, g1_func, g2_func, g3_func; 1015 PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 1016 PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 1017 PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 1018 PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 1019 PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 1020 PetscQuadrature quad; 1021 PetscTabulation *T, *TAux = NULL; 1022 PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *basisReal, *basisDerReal, *testReal, *testDerReal; 1023 const PetscScalar *constants; 1024 PetscReal *x; 1025 PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL; 1026 PetscInt NcI = 0, NcJ = 0, NcS, NcT; 1027 PetscInt dim, dimAux, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, e; 1028 PetscBool isCohesiveFieldI, isCohesiveFieldJ, isAffine, auxOnBd = PETSC_FALSE; 1029 const PetscReal *quadPoints, *quadWeights; 1030 PetscInt qNc, Nq, q, Np, dE; 1031 PetscErrorCode ierr; 1032 1033 PetscFunctionBegin; 1034 /* Hybrid discretization is posed directly on faces */ 1035 ierr = PetscDSGetDiscretization(ds, fieldI, (PetscObject *) &feI);CHKERRQ(ierr); 1036 ierr = PetscDSGetDiscretization(ds, fieldJ, (PetscObject *) &feJ);CHKERRQ(ierr); 1037 ierr = PetscFEGetSpatialDimension(feI, &dim);CHKERRQ(ierr); 1038 ierr = PetscFEGetQuadrature(feI, &quad);CHKERRQ(ierr); 1039 ierr = PetscDSGetNumFields(ds, &Nf);CHKERRQ(ierr); 1040 ierr = PetscDSGetTotalDimension(ds, &totDim);CHKERRQ(ierr); 1041 ierr = PetscDSGetComponentOffsets(ds, &uOff);CHKERRQ(ierr); 1042 ierr = PetscDSGetComponentDerivativeOffsets(ds, &uOff_x);CHKERRQ(ierr); 1043 switch(jtype) { 1044 case PETSCFE_JACOBIAN_PRE: ierr = PetscDSGetBdJacobianPreconditioner(ds, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);break; 1045 case PETSCFE_JACOBIAN: ierr = PetscDSGetBdJacobian(ds, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);break; 1046 case PETSCFE_JACOBIAN_DYN: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No boundary hybrid Jacobians :)"); 1047 } 1048 if (!g0_func && !g1_func && !g2_func && !g3_func) PetscFunctionReturn(0); 1049 ierr = PetscDSGetEvaluationArrays(ds, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 1050 ierr = PetscDSGetWorkspace(ds, &x, &basisReal, &basisDerReal, &testReal, &testDerReal);CHKERRQ(ierr); 1051 ierr = PetscDSGetWeakFormArrays(ds, NULL, NULL, &g0, &g1, &g2, &g3);CHKERRQ(ierr); 1052 ierr = PetscDSGetTabulation(ds, &T);CHKERRQ(ierr); 1053 ierr = PetscDSGetFieldOffset(ds, fieldI, &offsetI);CHKERRQ(ierr); 1054 ierr = PetscDSGetFieldOffset(ds, fieldJ, &offsetJ);CHKERRQ(ierr); 1055 ierr = PetscDSGetConstants(ds, &numConstants, &constants);CHKERRQ(ierr); 1056 if (dsAux) { 1057 ierr = PetscDSGetSpatialDimension(dsAux, &dimAux);CHKERRQ(ierr); 1058 ierr = PetscDSGetNumFields(dsAux, &NfAux);CHKERRQ(ierr); 1059 ierr = PetscDSGetTotalDimension(dsAux, &totDimAux);CHKERRQ(ierr); 1060 ierr = PetscDSGetComponentOffsets(dsAux, &aOff);CHKERRQ(ierr); 1061 ierr = PetscDSGetComponentDerivativeOffsets(dsAux, &aOff_x);CHKERRQ(ierr); 1062 ierr = PetscDSGetEvaluationArrays(dsAux, &a, NULL, &a_x);CHKERRQ(ierr); 1063 auxOnBd = dimAux == dim ? PETSC_TRUE : PETSC_FALSE; 1064 if (auxOnBd) {ierr = PetscDSGetTabulation(dsAux, &TAux);CHKERRQ(ierr);} 1065 else {ierr = PetscDSGetFaceTabulation(dsAux, &TAux);CHKERRQ(ierr);} 1066 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); 1067 } 1068 isCohesiveFieldI = fieldI == Nf-1 ? PETSC_TRUE : PETSC_FALSE; 1069 isCohesiveFieldJ = fieldJ == Nf-1 ? PETSC_TRUE : PETSC_FALSE; 1070 NcI = T[fieldI]->Nc; 1071 NcJ = T[fieldJ]->Nc; 1072 NcS = isCohesiveFieldI ? NcI : 2*NcI; 1073 NcT = isCohesiveFieldJ ? NcJ : 2*NcJ; 1074 Np = fgeom->numPoints; 1075 dE = fgeom->dimEmbed; 1076 isAffine = fgeom->isAffine; 1077 ierr = PetscArrayzero(g0, NcS*NcT);CHKERRQ(ierr); 1078 ierr = PetscArrayzero(g1, NcS*NcT*dE);CHKERRQ(ierr); 1079 ierr = PetscArrayzero(g2, NcS*NcT*dE);CHKERRQ(ierr); 1080 ierr = PetscArrayzero(g3, NcS*NcT*dE*dE);CHKERRQ(ierr); 1081 ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 1082 if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 1083 for (e = 0; e < Ne; ++e) { 1084 PetscFEGeom fegeom; 1085 const PetscInt face = fgeom->face[e][0]; 1086 1087 fegeom.dim = fgeom->dim; 1088 fegeom.dimEmbed = fgeom->dimEmbed; 1089 if (isAffine) { 1090 fegeom.v = x; 1091 fegeom.xi = fgeom->xi; 1092 fegeom.J = &fgeom->J[e*dE*dE]; 1093 fegeom.invJ = &fgeom->invJ[e*dE*dE]; 1094 fegeom.detJ = &fgeom->detJ[e]; 1095 fegeom.n = &fgeom->n[e*dE]; 1096 } 1097 for (q = 0; q < Nq; ++q) { 1098 PetscReal w; 1099 PetscInt c; 1100 1101 if (isAffine) { 1102 /* TODO Is it correct to have 'dim' here, or should it be 'dim-1'? */ 1103 CoordinatesRefToReal(dE, dim, fegeom.xi, &fgeom->v[e*dE], fegeom.J, &quadPoints[q*dim], x); 1104 } else { 1105 fegeom.v = &fegeom.v[(e*Np+q)*dE]; 1106 fegeom.J = &fgeom->J[(e*Np+q)*dE*dE]; 1107 fegeom.invJ = &fgeom->invJ[(e*Np+q)*dE*dE]; 1108 fegeom.detJ = &fgeom->detJ[e*Np+q]; 1109 fegeom.n = &fgeom->n[(e*Np+q)*dE]; 1110 } 1111 w = fegeom.detJ[0]*quadWeights[q]; 1112 if (debug > 1 && q < Np) { 1113 ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", fegeom.detJ[0]);CHKERRQ(ierr); 1114 #if !defined(PETSC_USE_COMPLEX) 1115 ierr = DMPrintCellMatrix(e, "invJ", dim, dim, fegeom.invJ);CHKERRQ(ierr); 1116 #endif 1117 } 1118 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 1119 if (coefficients) {ierr = PetscFEEvaluateFieldJets_Hybrid_Internal(ds, Nf, 0, q, T, &fegeom, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);CHKERRQ(ierr);} 1120 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);} 1121 if (g0_func) { 1122 ierr = PetscArrayzero(g0, NcS*NcT);CHKERRQ(ierr); 1123 g0_func(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); 1124 for (c = 0; c < NcS*NcT; ++c) g0[c] *= w; 1125 } 1126 if (g1_func) { 1127 ierr = PetscArrayzero(g1, NcS*NcT*dE);CHKERRQ(ierr); 1128 g1_func(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); 1129 for (c = 0; c < NcS*NcT*dE; ++c) g1[c] *= w; 1130 } 1131 if (g2_func) { 1132 ierr = PetscArrayzero(g2, NcS*NcT*dE);CHKERRQ(ierr); 1133 g2_func(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); 1134 for (c = 0; c < NcS*NcT*dE; ++c) g2[c] *= w; 1135 } 1136 if (g3_func) { 1137 ierr = PetscArrayzero(g3, NcS*NcT*dE*dE);CHKERRQ(ierr); 1138 g3_func(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); 1139 for (c = 0; c < NcS*NcT*dE*dE; ++c) g3[c] *= w; 1140 } 1141 1142 if (isCohesiveFieldI && isCohesiveFieldJ) { 1143 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); 1144 } else { 1145 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); 1146 } 1147 } 1148 if (debug > 1) { 1149 PetscInt fc, f, gc, g; 1150 1151 ierr = PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);CHKERRQ(ierr); 1152 for (fc = 0; fc < NcI; ++fc) { 1153 for (f = 0; f < T[fieldI]->Nb; ++f) { 1154 const PetscInt i = offsetI + f*NcI+fc; 1155 for (gc = 0; gc < NcJ; ++gc) { 1156 for (g = 0; g < T[fieldJ]->Nb; ++g) { 1157 const PetscInt j = offsetJ + g*NcJ+gc; 1158 ierr = PetscPrintf(PETSC_COMM_SELF, " elemMat[%d,%d,%d,%d]: %g\n", f, fc, g, gc, PetscRealPart(elemMat[eOffset+i*totDim+j]));CHKERRQ(ierr); 1159 } 1160 } 1161 ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr); 1162 } 1163 } 1164 } 1165 cOffset += totDim; 1166 cOffsetAux += totDimAux; 1167 eOffset += PetscSqr(totDim); 1168 } 1169 PetscFunctionReturn(0); 1170 } 1171 1172 static PetscErrorCode PetscFEInitialize_Basic(PetscFE fem) 1173 { 1174 PetscFunctionBegin; 1175 fem->ops->setfromoptions = NULL; 1176 fem->ops->setup = PetscFESetUp_Basic; 1177 fem->ops->view = PetscFEView_Basic; 1178 fem->ops->destroy = PetscFEDestroy_Basic; 1179 fem->ops->getdimension = PetscFEGetDimension_Basic; 1180 fem->ops->createtabulation = PetscFECreateTabulation_Basic; 1181 fem->ops->integrate = PetscFEIntegrate_Basic; 1182 fem->ops->integratebd = PetscFEIntegrateBd_Basic; 1183 fem->ops->integrateresidual = PetscFEIntegrateResidual_Basic; 1184 fem->ops->integratebdresidual = PetscFEIntegrateBdResidual_Basic; 1185 fem->ops->integratehybridresidual = PetscFEIntegrateHybridResidual_Basic; 1186 fem->ops->integratejacobianaction = NULL/* PetscFEIntegrateJacobianAction_Basic */; 1187 fem->ops->integratejacobian = PetscFEIntegrateJacobian_Basic; 1188 fem->ops->integratebdjacobian = PetscFEIntegrateBdJacobian_Basic; 1189 fem->ops->integratehybridjacobian = PetscFEIntegrateHybridJacobian_Basic; 1190 PetscFunctionReturn(0); 1191 } 1192 1193 /*MC 1194 PETSCFEBASIC = "basic" - A PetscFE object that integrates with basic tiling and no vectorization 1195 1196 Level: intermediate 1197 1198 .seealso: PetscFEType, PetscFECreate(), PetscFESetType() 1199 M*/ 1200 1201 PETSC_EXTERN PetscErrorCode PetscFECreate_Basic(PetscFE fem) 1202 { 1203 PetscFE_Basic *b; 1204 PetscErrorCode ierr; 1205 1206 PetscFunctionBegin; 1207 PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 1208 ierr = PetscNewLog(fem,&b);CHKERRQ(ierr); 1209 fem->data = b; 1210 1211 ierr = PetscFEInitialize_Basic(fem);CHKERRQ(ierr); 1212 PetscFunctionReturn(0); 1213 } 1214