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