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