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