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