1*20cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/ 2*20cf1dd8SToby Isaac #include <petscblaslapack.h> 3*20cf1dd8SToby Isaac 4*20cf1dd8SToby Isaac PetscErrorCode PetscFEDestroy_Basic(PetscFE fem) 5*20cf1dd8SToby Isaac { 6*20cf1dd8SToby Isaac PetscFE_Basic *b = (PetscFE_Basic *) fem->data; 7*20cf1dd8SToby Isaac PetscErrorCode ierr; 8*20cf1dd8SToby Isaac 9*20cf1dd8SToby Isaac PetscFunctionBegin; 10*20cf1dd8SToby Isaac ierr = PetscFree(b);CHKERRQ(ierr); 11*20cf1dd8SToby Isaac PetscFunctionReturn(0); 12*20cf1dd8SToby Isaac } 13*20cf1dd8SToby Isaac 14*20cf1dd8SToby Isaac PetscErrorCode PetscFEView_Basic_Ascii(PetscFE fe, PetscViewer viewer) 15*20cf1dd8SToby Isaac { 16*20cf1dd8SToby Isaac PetscSpace basis; 17*20cf1dd8SToby Isaac PetscDualSpace dual; 18*20cf1dd8SToby Isaac PetscQuadrature q = NULL; 19*20cf1dd8SToby Isaac PetscInt dim, Nc, Nq; 20*20cf1dd8SToby Isaac PetscViewerFormat format; 21*20cf1dd8SToby Isaac PetscErrorCode ierr; 22*20cf1dd8SToby Isaac 23*20cf1dd8SToby Isaac PetscFunctionBegin; 24*20cf1dd8SToby Isaac ierr = PetscFEGetBasisSpace(fe, &basis);CHKERRQ(ierr); 25*20cf1dd8SToby Isaac ierr = PetscFEGetDualSpace(fe, &dual);CHKERRQ(ierr); 26*20cf1dd8SToby Isaac ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr); 27*20cf1dd8SToby Isaac ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 28*20cf1dd8SToby Isaac ierr = PetscQuadratureGetData(q, &dim, NULL, &Nq, NULL, NULL);CHKERRQ(ierr); 29*20cf1dd8SToby Isaac ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr); 30*20cf1dd8SToby Isaac ierr = PetscViewerASCIIPrintf(viewer, "Basic Finite Element:\n");CHKERRQ(ierr); 31*20cf1dd8SToby Isaac ierr = PetscViewerASCIIPrintf(viewer, " dimension: %d\n", dim);CHKERRQ(ierr); 32*20cf1dd8SToby Isaac ierr = PetscViewerASCIIPrintf(viewer, " components: %d\n", Nc);CHKERRQ(ierr); 33*20cf1dd8SToby Isaac ierr = PetscViewerASCIIPrintf(viewer, " num quad points: %d\n", Nq);CHKERRQ(ierr); 34*20cf1dd8SToby Isaac if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 35*20cf1dd8SToby Isaac ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 36*20cf1dd8SToby Isaac ierr = PetscQuadratureView(q, viewer);CHKERRQ(ierr); 37*20cf1dd8SToby Isaac ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 38*20cf1dd8SToby Isaac } 39*20cf1dd8SToby Isaac ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 40*20cf1dd8SToby Isaac ierr = PetscSpaceView(basis, viewer);CHKERRQ(ierr); 41*20cf1dd8SToby Isaac ierr = PetscDualSpaceView(dual, viewer);CHKERRQ(ierr); 42*20cf1dd8SToby Isaac ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 43*20cf1dd8SToby Isaac PetscFunctionReturn(0); 44*20cf1dd8SToby Isaac } 45*20cf1dd8SToby Isaac 46*20cf1dd8SToby Isaac PetscErrorCode PetscFEView_Basic(PetscFE fe, PetscViewer viewer) 47*20cf1dd8SToby Isaac { 48*20cf1dd8SToby Isaac PetscBool iascii; 49*20cf1dd8SToby Isaac PetscErrorCode ierr; 50*20cf1dd8SToby Isaac 51*20cf1dd8SToby Isaac PetscFunctionBegin; 52*20cf1dd8SToby Isaac PetscValidHeaderSpecific(fe, PETSCFE_CLASSID, 1); 53*20cf1dd8SToby Isaac PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 54*20cf1dd8SToby Isaac ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 55*20cf1dd8SToby Isaac if (iascii) {ierr = PetscFEView_Basic_Ascii(fe, viewer);CHKERRQ(ierr);} 56*20cf1dd8SToby Isaac PetscFunctionReturn(0); 57*20cf1dd8SToby Isaac } 58*20cf1dd8SToby Isaac 59*20cf1dd8SToby Isaac /* Construct the change of basis from prime basis to nodal basis */ 60*20cf1dd8SToby Isaac PetscErrorCode PetscFESetUp_Basic(PetscFE fem) 61*20cf1dd8SToby Isaac { 62*20cf1dd8SToby Isaac PetscScalar *work, *invVscalar; 63*20cf1dd8SToby Isaac PetscBLASInt *pivots; 64*20cf1dd8SToby Isaac PetscBLASInt n, info; 65*20cf1dd8SToby Isaac PetscInt pdim, j; 66*20cf1dd8SToby Isaac PetscErrorCode ierr; 67*20cf1dd8SToby Isaac 68*20cf1dd8SToby Isaac PetscFunctionBegin; 69*20cf1dd8SToby Isaac ierr = PetscDualSpaceGetDimension(fem->dualSpace, &pdim);CHKERRQ(ierr); 70*20cf1dd8SToby Isaac ierr = PetscMalloc1(pdim*pdim,&fem->invV);CHKERRQ(ierr); 71*20cf1dd8SToby Isaac #if defined(PETSC_USE_COMPLEX) 72*20cf1dd8SToby Isaac ierr = PetscMalloc1(pdim*pdim,&invVscalar);CHKERRQ(ierr); 73*20cf1dd8SToby Isaac #else 74*20cf1dd8SToby Isaac invVscalar = fem->invV; 75*20cf1dd8SToby Isaac #endif 76*20cf1dd8SToby Isaac for (j = 0; j < pdim; ++j) { 77*20cf1dd8SToby Isaac PetscReal *Bf; 78*20cf1dd8SToby Isaac PetscQuadrature f; 79*20cf1dd8SToby Isaac const PetscReal *points, *weights; 80*20cf1dd8SToby Isaac PetscInt Nc, Nq, q, k, c; 81*20cf1dd8SToby Isaac 82*20cf1dd8SToby Isaac ierr = PetscDualSpaceGetFunctional(fem->dualSpace, j, &f);CHKERRQ(ierr); 83*20cf1dd8SToby Isaac ierr = PetscQuadratureGetData(f, NULL, &Nc, &Nq, &points, &weights);CHKERRQ(ierr); 84*20cf1dd8SToby Isaac ierr = PetscMalloc1(Nc*Nq*pdim,&Bf);CHKERRQ(ierr); 85*20cf1dd8SToby Isaac ierr = PetscSpaceEvaluate(fem->basisSpace, Nq, points, Bf, NULL, NULL);CHKERRQ(ierr); 86*20cf1dd8SToby Isaac for (k = 0; k < pdim; ++k) { 87*20cf1dd8SToby Isaac /* V_{jk} = n_j(\phi_k) = \int \phi_k(x) n_j(x) dx */ 88*20cf1dd8SToby Isaac invVscalar[j*pdim+k] = 0.0; 89*20cf1dd8SToby Isaac 90*20cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 91*20cf1dd8SToby Isaac for (c = 0; c < Nc; ++c) invVscalar[j*pdim+k] += Bf[(q*pdim + k)*Nc + c]*weights[q*Nc + c]; 92*20cf1dd8SToby Isaac } 93*20cf1dd8SToby Isaac } 94*20cf1dd8SToby Isaac ierr = PetscFree(Bf);CHKERRQ(ierr); 95*20cf1dd8SToby Isaac } 96*20cf1dd8SToby Isaac ierr = PetscMalloc2(pdim,&pivots,pdim,&work);CHKERRQ(ierr); 97*20cf1dd8SToby Isaac n = pdim; 98*20cf1dd8SToby Isaac PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&n, &n, invVscalar, &n, pivots, &info)); 99*20cf1dd8SToby Isaac PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&n, invVscalar, &n, pivots, work, &n, &info)); 100*20cf1dd8SToby Isaac #if defined(PETSC_USE_COMPLEX) 101*20cf1dd8SToby Isaac for (j = 0; j < pdim*pdim; j++) fem->invV[j] = PetscRealPart(invVscalar[j]); 102*20cf1dd8SToby Isaac ierr = PetscFree(invVscalar);CHKERRQ(ierr); 103*20cf1dd8SToby Isaac #endif 104*20cf1dd8SToby Isaac ierr = PetscFree2(pivots,work);CHKERRQ(ierr); 105*20cf1dd8SToby Isaac PetscFunctionReturn(0); 106*20cf1dd8SToby Isaac } 107*20cf1dd8SToby Isaac 108*20cf1dd8SToby Isaac PetscErrorCode PetscFEGetDimension_Basic(PetscFE fem, PetscInt *dim) 109*20cf1dd8SToby Isaac { 110*20cf1dd8SToby Isaac PetscErrorCode ierr; 111*20cf1dd8SToby Isaac 112*20cf1dd8SToby Isaac PetscFunctionBegin; 113*20cf1dd8SToby Isaac ierr = PetscDualSpaceGetDimension(fem->dualSpace, dim);CHKERRQ(ierr); 114*20cf1dd8SToby Isaac PetscFunctionReturn(0); 115*20cf1dd8SToby Isaac } 116*20cf1dd8SToby Isaac 117*20cf1dd8SToby Isaac PetscErrorCode PetscFEGetTabulation_Basic(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscReal *B, PetscReal *D, PetscReal *H) 118*20cf1dd8SToby Isaac { 119*20cf1dd8SToby Isaac DM dm; 120*20cf1dd8SToby Isaac PetscInt pdim; /* Dimension of FE space P */ 121*20cf1dd8SToby Isaac PetscInt dim; /* Spatial dimension */ 122*20cf1dd8SToby Isaac PetscInt Nc; /* Field components */ 123*20cf1dd8SToby Isaac PetscReal *tmpB, *tmpD, *tmpH; 124*20cf1dd8SToby Isaac PetscInt p, d, j, k, c; 125*20cf1dd8SToby Isaac PetscErrorCode ierr; 126*20cf1dd8SToby Isaac 127*20cf1dd8SToby Isaac PetscFunctionBegin; 128*20cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(fem->dualSpace, &dm);CHKERRQ(ierr); 129*20cf1dd8SToby Isaac ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 130*20cf1dd8SToby Isaac ierr = PetscDualSpaceGetDimension(fem->dualSpace, &pdim);CHKERRQ(ierr); 131*20cf1dd8SToby Isaac ierr = PetscFEGetNumComponents(fem, &Nc);CHKERRQ(ierr); 132*20cf1dd8SToby Isaac /* Evaluate the prime basis functions at all points */ 133*20cf1dd8SToby Isaac if (B) {ierr = DMGetWorkArray(dm, npoints*pdim*Nc, MPIU_REAL, &tmpB);CHKERRQ(ierr);} 134*20cf1dd8SToby Isaac if (D) {ierr = DMGetWorkArray(dm, npoints*pdim*Nc*dim, MPIU_REAL, &tmpD);CHKERRQ(ierr);} 135*20cf1dd8SToby Isaac if (H) {ierr = DMGetWorkArray(dm, npoints*pdim*Nc*dim*dim, MPIU_REAL, &tmpH);CHKERRQ(ierr);} 136*20cf1dd8SToby Isaac ierr = PetscSpaceEvaluate(fem->basisSpace, npoints, points, B ? tmpB : NULL, D ? tmpD : NULL, H ? tmpH : NULL);CHKERRQ(ierr); 137*20cf1dd8SToby Isaac /* Translate to the nodal basis */ 138*20cf1dd8SToby Isaac for (p = 0; p < npoints; ++p) { 139*20cf1dd8SToby Isaac if (B) { 140*20cf1dd8SToby Isaac /* Multiply by V^{-1} (pdim x pdim) */ 141*20cf1dd8SToby Isaac for (j = 0; j < pdim; ++j) { 142*20cf1dd8SToby Isaac const PetscInt i = (p*pdim + j)*Nc; 143*20cf1dd8SToby Isaac 144*20cf1dd8SToby Isaac for (c = 0; c < Nc; ++c) { 145*20cf1dd8SToby Isaac B[i+c] = 0.0; 146*20cf1dd8SToby Isaac for (k = 0; k < pdim; ++k) { 147*20cf1dd8SToby Isaac B[i+c] += fem->invV[k*pdim+j] * tmpB[(p*pdim + k)*Nc+c]; 148*20cf1dd8SToby Isaac } 149*20cf1dd8SToby Isaac } 150*20cf1dd8SToby Isaac } 151*20cf1dd8SToby Isaac } 152*20cf1dd8SToby Isaac if (D) { 153*20cf1dd8SToby Isaac /* Multiply by V^{-1} (pdim x pdim) */ 154*20cf1dd8SToby Isaac for (j = 0; j < pdim; ++j) { 155*20cf1dd8SToby Isaac for (c = 0; c < Nc; ++c) { 156*20cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 157*20cf1dd8SToby Isaac const PetscInt i = ((p*pdim + j)*Nc + c)*dim + d; 158*20cf1dd8SToby Isaac 159*20cf1dd8SToby Isaac D[i] = 0.0; 160*20cf1dd8SToby Isaac for (k = 0; k < pdim; ++k) { 161*20cf1dd8SToby Isaac D[i] += fem->invV[k*pdim+j] * tmpD[((p*pdim + k)*Nc + c)*dim + d]; 162*20cf1dd8SToby Isaac } 163*20cf1dd8SToby Isaac } 164*20cf1dd8SToby Isaac } 165*20cf1dd8SToby Isaac } 166*20cf1dd8SToby Isaac } 167*20cf1dd8SToby Isaac if (H) { 168*20cf1dd8SToby Isaac /* Multiply by V^{-1} (pdim x pdim) */ 169*20cf1dd8SToby Isaac for (j = 0; j < pdim; ++j) { 170*20cf1dd8SToby Isaac for (c = 0; c < Nc; ++c) { 171*20cf1dd8SToby Isaac for (d = 0; d < dim*dim; ++d) { 172*20cf1dd8SToby Isaac const PetscInt i = ((p*pdim + j)*Nc + c)*dim*dim + d; 173*20cf1dd8SToby Isaac 174*20cf1dd8SToby Isaac H[i] = 0.0; 175*20cf1dd8SToby Isaac for (k = 0; k < pdim; ++k) { 176*20cf1dd8SToby Isaac H[i] += fem->invV[k*pdim+j] * tmpH[((p*pdim + k)*Nc + c)*dim*dim + d]; 177*20cf1dd8SToby Isaac } 178*20cf1dd8SToby Isaac } 179*20cf1dd8SToby Isaac } 180*20cf1dd8SToby Isaac } 181*20cf1dd8SToby Isaac } 182*20cf1dd8SToby Isaac } 183*20cf1dd8SToby Isaac if (B) {ierr = DMRestoreWorkArray(dm, npoints*pdim*Nc, MPIU_REAL, &tmpB);CHKERRQ(ierr);} 184*20cf1dd8SToby Isaac if (D) {ierr = DMRestoreWorkArray(dm, npoints*pdim*Nc*dim, MPIU_REAL, &tmpD);CHKERRQ(ierr);} 185*20cf1dd8SToby Isaac if (H) {ierr = DMRestoreWorkArray(dm, npoints*pdim*Nc*dim*dim, MPIU_REAL, &tmpH);CHKERRQ(ierr);} 186*20cf1dd8SToby Isaac PetscFunctionReturn(0); 187*20cf1dd8SToby Isaac } 188*20cf1dd8SToby Isaac 189*20cf1dd8SToby Isaac PetscErrorCode PetscFEIntegrate_Basic(PetscFE fem, PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom, 190*20cf1dd8SToby Isaac const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[]) 191*20cf1dd8SToby Isaac { 192*20cf1dd8SToby Isaac const PetscInt debug = 0; 193*20cf1dd8SToby Isaac PetscPointFunc obj_func; 194*20cf1dd8SToby Isaac PetscQuadrature quad; 195*20cf1dd8SToby Isaac PetscScalar *u, *u_x, *a, *a_x, *refSpaceDer, *refSpaceDerAux; 196*20cf1dd8SToby Isaac const PetscScalar *constants; 197*20cf1dd8SToby Isaac PetscReal *x; 198*20cf1dd8SToby Isaac PetscReal **B, **D, **BAux = NULL, **DAux = NULL; 199*20cf1dd8SToby Isaac PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL; 200*20cf1dd8SToby Isaac PetscInt dim, dE, Np, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e; 201*20cf1dd8SToby Isaac PetscBool isAffine; 202*20cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 203*20cf1dd8SToby Isaac PetscInt qNc, Nq, q; 204*20cf1dd8SToby Isaac PetscErrorCode ierr; 205*20cf1dd8SToby Isaac 206*20cf1dd8SToby Isaac PetscFunctionBegin; 207*20cf1dd8SToby Isaac ierr = PetscDSGetObjective(prob, field, &obj_func);CHKERRQ(ierr); 208*20cf1dd8SToby Isaac if (!obj_func) PetscFunctionReturn(0); 209*20cf1dd8SToby Isaac ierr = PetscFEGetSpatialDimension(fem, &dim);CHKERRQ(ierr); 210*20cf1dd8SToby Isaac ierr = PetscFEGetQuadrature(fem, &quad);CHKERRQ(ierr); 211*20cf1dd8SToby Isaac ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 212*20cf1dd8SToby Isaac ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 213*20cf1dd8SToby Isaac ierr = PetscDSGetDimensions(prob, &Nb);CHKERRQ(ierr); 214*20cf1dd8SToby Isaac ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr); 215*20cf1dd8SToby Isaac ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr); 216*20cf1dd8SToby Isaac ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr); 217*20cf1dd8SToby Isaac ierr = PetscDSGetEvaluationArrays(prob, &u, NULL, &u_x);CHKERRQ(ierr); 218*20cf1dd8SToby Isaac ierr = PetscDSGetRefCoordArrays(prob, &x, &refSpaceDer);CHKERRQ(ierr); 219*20cf1dd8SToby Isaac ierr = PetscDSGetTabulation(prob, &B, &D);CHKERRQ(ierr); 220*20cf1dd8SToby Isaac ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr); 221*20cf1dd8SToby Isaac if (probAux) { 222*20cf1dd8SToby Isaac ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr); 223*20cf1dd8SToby Isaac ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr); 224*20cf1dd8SToby Isaac ierr = PetscDSGetDimensions(probAux, &NbAux);CHKERRQ(ierr); 225*20cf1dd8SToby Isaac ierr = PetscDSGetComponents(probAux, &NcAux);CHKERRQ(ierr); 226*20cf1dd8SToby Isaac ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr); 227*20cf1dd8SToby Isaac ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr); 228*20cf1dd8SToby Isaac ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);CHKERRQ(ierr); 229*20cf1dd8SToby Isaac ierr = PetscDSGetRefCoordArrays(probAux, NULL, &refSpaceDerAux);CHKERRQ(ierr); 230*20cf1dd8SToby Isaac ierr = PetscDSGetTabulation(probAux, &BAux, &DAux);CHKERRQ(ierr); 231*20cf1dd8SToby Isaac } 232*20cf1dd8SToby Isaac ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 233*20cf1dd8SToby Isaac Np = cgeom->numPoints; 234*20cf1dd8SToby Isaac dE = cgeom->dimEmbed; 235*20cf1dd8SToby Isaac isAffine = cgeom->isAffine; 236*20cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 237*20cf1dd8SToby Isaac const PetscReal *v0 = &cgeom->v[e*Np*dE]; 238*20cf1dd8SToby Isaac const PetscReal *J = &cgeom->J[e*Np*dE*dE]; 239*20cf1dd8SToby Isaac 240*20cf1dd8SToby Isaac if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 241*20cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 242*20cf1dd8SToby Isaac PetscScalar integrand; 243*20cf1dd8SToby Isaac const PetscReal *v; 244*20cf1dd8SToby Isaac const PetscReal *invJ; 245*20cf1dd8SToby Isaac PetscReal detJ; 246*20cf1dd8SToby Isaac 247*20cf1dd8SToby Isaac if (isAffine) { 248*20cf1dd8SToby Isaac CoordinatesRefToReal(dE, dim, cgeom->xi, v0, J, &quadPoints[q*dim], x); 249*20cf1dd8SToby Isaac v = x; 250*20cf1dd8SToby Isaac invJ = &cgeom->invJ[e*dE*dE]; 251*20cf1dd8SToby Isaac detJ = cgeom->detJ[e]; 252*20cf1dd8SToby Isaac } else { 253*20cf1dd8SToby Isaac v = &v0[q*dE]; 254*20cf1dd8SToby Isaac invJ = &cgeom->invJ[(e*Np+q)*dE*dE]; 255*20cf1dd8SToby Isaac detJ = cgeom->detJ[e*Np + q]; 256*20cf1dd8SToby Isaac } 257*20cf1dd8SToby Isaac if (debug > 1 && q < Np) { 258*20cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", detJ);CHKERRQ(ierr); 259*20cf1dd8SToby Isaac #ifndef PETSC_USE_COMPLEX 260*20cf1dd8SToby Isaac ierr = DMPrintCellMatrix(e, "invJ", dim, dim, invJ);CHKERRQ(ierr); 261*20cf1dd8SToby Isaac #endif 262*20cf1dd8SToby Isaac } 263*20cf1dd8SToby Isaac if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 264*20cf1dd8SToby Isaac EvaluateFieldJets(dim, Nf, Nb, Nc, q, B, D, refSpaceDer, invJ, &coefficients[cOffset], NULL, u, u_x, NULL); 265*20cf1dd8SToby Isaac if (probAux) EvaluateFieldJets(dim, NfAux, NbAux, NcAux, q, BAux, DAux, refSpaceDerAux, invJ, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL); 266*20cf1dd8SToby Isaac obj_func(dim, Nf, NfAux, uOff, uOff_x, u, NULL, u_x, aOff, aOff_x, a, NULL, a_x, 0.0, v, numConstants, constants, &integrand); 267*20cf1dd8SToby Isaac integrand *= detJ*quadWeights[q]; 268*20cf1dd8SToby Isaac integral[e*Nf+field] += integrand; 269*20cf1dd8SToby Isaac if (debug > 1) {ierr = PetscPrintf(PETSC_COMM_SELF, " int: %g %g\n", (double) PetscRealPart(integrand), (double) PetscRealPart(integral[field]));CHKERRQ(ierr);} 270*20cf1dd8SToby Isaac } 271*20cf1dd8SToby Isaac cOffset += totDim; 272*20cf1dd8SToby Isaac cOffsetAux += totDimAux; 273*20cf1dd8SToby Isaac } 274*20cf1dd8SToby Isaac PetscFunctionReturn(0); 275*20cf1dd8SToby Isaac } 276*20cf1dd8SToby Isaac 277*20cf1dd8SToby Isaac PetscErrorCode PetscFEIntegrateResidual_Basic(PetscFE fem, PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *cgeom, 278*20cf1dd8SToby Isaac const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 279*20cf1dd8SToby Isaac { 280*20cf1dd8SToby Isaac const PetscInt debug = 0; 281*20cf1dd8SToby Isaac PetscPointFunc f0_func; 282*20cf1dd8SToby Isaac PetscPointFunc f1_func; 283*20cf1dd8SToby Isaac PetscQuadrature quad; 284*20cf1dd8SToby Isaac PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *refSpaceDer, *refSpaceDerAux; 285*20cf1dd8SToby Isaac const PetscScalar *constants; 286*20cf1dd8SToby Isaac PetscReal *x; 287*20cf1dd8SToby Isaac PetscReal **B, **D, **BAux = NULL, **DAux = NULL, *BI, *DI; 288*20cf1dd8SToby Isaac PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL; 289*20cf1dd8SToby Isaac PetscInt dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NbI, NcI; 290*20cf1dd8SToby Isaac PetscInt dE, Np; 291*20cf1dd8SToby Isaac PetscBool isAffine; 292*20cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 293*20cf1dd8SToby Isaac PetscInt qNc, Nq, q; 294*20cf1dd8SToby Isaac PetscErrorCode ierr; 295*20cf1dd8SToby Isaac 296*20cf1dd8SToby Isaac PetscFunctionBegin; 297*20cf1dd8SToby Isaac ierr = PetscFEGetSpatialDimension(fem, &dim);CHKERRQ(ierr); 298*20cf1dd8SToby Isaac ierr = PetscFEGetQuadrature(fem, &quad);CHKERRQ(ierr); 299*20cf1dd8SToby Isaac ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 300*20cf1dd8SToby Isaac ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 301*20cf1dd8SToby Isaac ierr = PetscDSGetDimensions(prob, &Nb);CHKERRQ(ierr); 302*20cf1dd8SToby Isaac ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr); 303*20cf1dd8SToby Isaac ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr); 304*20cf1dd8SToby Isaac ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr); 305*20cf1dd8SToby Isaac ierr = PetscDSGetFieldOffset(prob, field, &fOffset);CHKERRQ(ierr); 306*20cf1dd8SToby Isaac ierr = PetscDSGetResidual(prob, field, &f0_func, &f1_func);CHKERRQ(ierr); 307*20cf1dd8SToby Isaac ierr = PetscDSGetEvaluationArrays(prob, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 308*20cf1dd8SToby Isaac ierr = PetscDSGetRefCoordArrays(prob, &x, &refSpaceDer);CHKERRQ(ierr); 309*20cf1dd8SToby Isaac ierr = PetscDSGetWeakFormArrays(prob, &f0, &f1, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 310*20cf1dd8SToby Isaac ierr = PetscDSGetTabulation(prob, &B, &D);CHKERRQ(ierr); 311*20cf1dd8SToby Isaac ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr); 312*20cf1dd8SToby Isaac if (probAux) { 313*20cf1dd8SToby Isaac ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr); 314*20cf1dd8SToby Isaac ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr); 315*20cf1dd8SToby Isaac ierr = PetscDSGetDimensions(probAux, &NbAux);CHKERRQ(ierr); 316*20cf1dd8SToby Isaac ierr = PetscDSGetComponents(probAux, &NcAux);CHKERRQ(ierr); 317*20cf1dd8SToby Isaac ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr); 318*20cf1dd8SToby Isaac ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr); 319*20cf1dd8SToby Isaac ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);CHKERRQ(ierr); 320*20cf1dd8SToby Isaac ierr = PetscDSGetRefCoordArrays(probAux, NULL, &refSpaceDerAux);CHKERRQ(ierr); 321*20cf1dd8SToby Isaac ierr = PetscDSGetTabulation(probAux, &BAux, &DAux);CHKERRQ(ierr); 322*20cf1dd8SToby Isaac } 323*20cf1dd8SToby Isaac NbI = Nb[field]; 324*20cf1dd8SToby Isaac NcI = Nc[field]; 325*20cf1dd8SToby Isaac BI = B[field]; 326*20cf1dd8SToby Isaac DI = D[field]; 327*20cf1dd8SToby Isaac ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 328*20cf1dd8SToby Isaac Np = cgeom->numPoints; 329*20cf1dd8SToby Isaac dE = cgeom->dimEmbed; 330*20cf1dd8SToby Isaac isAffine = cgeom->isAffine; 331*20cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 332*20cf1dd8SToby Isaac const PetscReal *v0 = &cgeom->v[e*Np*dE]; 333*20cf1dd8SToby Isaac const PetscReal *J = &cgeom->J[e*Np*dE*dE]; 334*20cf1dd8SToby Isaac 335*20cf1dd8SToby Isaac if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 336*20cf1dd8SToby Isaac ierr = PetscMemzero(f0, Nq*NcI* sizeof(PetscScalar));CHKERRQ(ierr); 337*20cf1dd8SToby Isaac ierr = PetscMemzero(f1, Nq*NcI*dim * sizeof(PetscScalar));CHKERRQ(ierr); 338*20cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 339*20cf1dd8SToby Isaac const PetscReal *v; 340*20cf1dd8SToby Isaac const PetscReal *invJ; 341*20cf1dd8SToby Isaac PetscReal detJ; 342*20cf1dd8SToby Isaac 343*20cf1dd8SToby Isaac if (isAffine) { 344*20cf1dd8SToby Isaac CoordinatesRefToReal(dE, dim, cgeom->xi, v0, J, &quadPoints[q*dim], x); 345*20cf1dd8SToby Isaac v = x; 346*20cf1dd8SToby Isaac invJ = &cgeom->invJ[e*dE*dE]; 347*20cf1dd8SToby Isaac detJ = cgeom->detJ[e]; 348*20cf1dd8SToby Isaac } else { 349*20cf1dd8SToby Isaac v = &v0[q*dE]; 350*20cf1dd8SToby Isaac invJ = &cgeom->invJ[(e*Np+q)*dE*dE]; 351*20cf1dd8SToby Isaac detJ = cgeom->detJ[e*Np + q]; 352*20cf1dd8SToby Isaac } 353*20cf1dd8SToby Isaac if (debug > 1 && q < Np) { 354*20cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", detJ);CHKERRQ(ierr); 355*20cf1dd8SToby Isaac #ifndef PETSC_USE_COMPLEX 356*20cf1dd8SToby Isaac ierr = DMPrintCellMatrix(e, "invJ", dim, dim, invJ);CHKERRQ(ierr); 357*20cf1dd8SToby Isaac #endif 358*20cf1dd8SToby Isaac } 359*20cf1dd8SToby Isaac if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 360*20cf1dd8SToby Isaac EvaluateFieldJets(dim, Nf, Nb, Nc, q, B, D, refSpaceDer, invJ, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t); 361*20cf1dd8SToby Isaac if (probAux) EvaluateFieldJets(dim, NfAux, NbAux, NcAux, q, BAux, DAux, refSpaceDerAux, invJ, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL); 362*20cf1dd8SToby Isaac if (f0_func) f0_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, v, numConstants, constants, &f0[q*NcI]); 363*20cf1dd8SToby Isaac if (f1_func) { 364*20cf1dd8SToby Isaac ierr = PetscMemzero(refSpaceDer, NcI*dim * sizeof(PetscScalar));CHKERRQ(ierr); 365*20cf1dd8SToby Isaac f1_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, v, numConstants, constants, refSpaceDer); 366*20cf1dd8SToby Isaac } 367*20cf1dd8SToby Isaac TransformF(dim, NcI, q, invJ, detJ, quadWeights, refSpaceDer, f0_func ? f0 : NULL, f1_func ? f1 : NULL); 368*20cf1dd8SToby Isaac } 369*20cf1dd8SToby Isaac UpdateElementVec(dim, Nq, NbI, NcI, BI, DI, f0, f1, &elemVec[cOffset+fOffset]); 370*20cf1dd8SToby Isaac cOffset += totDim; 371*20cf1dd8SToby Isaac cOffsetAux += totDimAux; 372*20cf1dd8SToby Isaac } 373*20cf1dd8SToby Isaac PetscFunctionReturn(0); 374*20cf1dd8SToby Isaac } 375*20cf1dd8SToby Isaac 376*20cf1dd8SToby Isaac PetscErrorCode PetscFEIntegrateBdResidual_Basic(PetscFE fem, PetscDS prob, PetscInt field, PetscInt Ne, PetscFEGeom *fgeom, 377*20cf1dd8SToby Isaac const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[]) 378*20cf1dd8SToby Isaac { 379*20cf1dd8SToby Isaac const PetscInt debug = 0; 380*20cf1dd8SToby Isaac PetscBdPointFunc f0_func; 381*20cf1dd8SToby Isaac PetscBdPointFunc f1_func; 382*20cf1dd8SToby Isaac PetscQuadrature quad; 383*20cf1dd8SToby Isaac PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *refSpaceDer, *refSpaceDerAux; 384*20cf1dd8SToby Isaac const PetscScalar *constants; 385*20cf1dd8SToby Isaac PetscReal *x; 386*20cf1dd8SToby Isaac PetscReal **B, **D, **BAux = NULL, **DAux = NULL, *BI, *DI; 387*20cf1dd8SToby Isaac PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL; 388*20cf1dd8SToby Isaac PetscInt dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NbI, NcI; 389*20cf1dd8SToby Isaac PetscBool isAffine; 390*20cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 391*20cf1dd8SToby Isaac PetscInt qNc, Nq, q, Np, dE; 392*20cf1dd8SToby Isaac PetscErrorCode ierr; 393*20cf1dd8SToby Isaac 394*20cf1dd8SToby Isaac PetscFunctionBegin; 395*20cf1dd8SToby Isaac ierr = PetscFEGetSpatialDimension(fem, &dim);CHKERRQ(ierr); 396*20cf1dd8SToby Isaac ierr = PetscFEGetFaceQuadrature(fem, &quad);CHKERRQ(ierr); 397*20cf1dd8SToby Isaac ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 398*20cf1dd8SToby Isaac ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 399*20cf1dd8SToby Isaac ierr = PetscDSGetDimensions(prob, &Nb);CHKERRQ(ierr); 400*20cf1dd8SToby Isaac ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr); 401*20cf1dd8SToby Isaac ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr); 402*20cf1dd8SToby Isaac ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr); 403*20cf1dd8SToby Isaac ierr = PetscDSGetFieldOffset(prob, field, &fOffset);CHKERRQ(ierr); 404*20cf1dd8SToby Isaac ierr = PetscDSGetBdResidual(prob, field, &f0_func, &f1_func);CHKERRQ(ierr); 405*20cf1dd8SToby Isaac if (!f0_func && !f1_func) PetscFunctionReturn(0); 406*20cf1dd8SToby Isaac ierr = PetscDSGetEvaluationArrays(prob, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 407*20cf1dd8SToby Isaac ierr = PetscDSGetRefCoordArrays(prob, &x, &refSpaceDer);CHKERRQ(ierr); 408*20cf1dd8SToby Isaac ierr = PetscDSGetWeakFormArrays(prob, &f0, &f1, NULL, NULL, NULL, NULL);CHKERRQ(ierr); 409*20cf1dd8SToby Isaac ierr = PetscDSGetFaceTabulation(prob, &B, &D);CHKERRQ(ierr); 410*20cf1dd8SToby Isaac ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr); 411*20cf1dd8SToby Isaac if (probAux) { 412*20cf1dd8SToby Isaac ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr); 413*20cf1dd8SToby Isaac ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr); 414*20cf1dd8SToby Isaac ierr = PetscDSGetDimensions(probAux, &NbAux);CHKERRQ(ierr); 415*20cf1dd8SToby Isaac ierr = PetscDSGetComponents(probAux, &NcAux);CHKERRQ(ierr); 416*20cf1dd8SToby Isaac ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr); 417*20cf1dd8SToby Isaac ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr); 418*20cf1dd8SToby Isaac ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);CHKERRQ(ierr); 419*20cf1dd8SToby Isaac ierr = PetscDSGetRefCoordArrays(probAux, NULL, &refSpaceDerAux);CHKERRQ(ierr); 420*20cf1dd8SToby Isaac ierr = PetscDSGetFaceTabulation(probAux, &BAux, &DAux);CHKERRQ(ierr); 421*20cf1dd8SToby Isaac } 422*20cf1dd8SToby Isaac NbI = Nb[field]; 423*20cf1dd8SToby Isaac NcI = Nc[field]; 424*20cf1dd8SToby Isaac BI = B[field]; 425*20cf1dd8SToby Isaac DI = D[field]; 426*20cf1dd8SToby Isaac ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 427*20cf1dd8SToby Isaac Np = fgeom->numPoints; 428*20cf1dd8SToby Isaac dE = fgeom->dimEmbed; 429*20cf1dd8SToby Isaac isAffine = fgeom->isAffine; 430*20cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 431*20cf1dd8SToby Isaac const PetscReal *v0 = &fgeom->v[e*Np*dE]; 432*20cf1dd8SToby Isaac const PetscReal *J = &fgeom->J[e*Np*dE*dE]; 433*20cf1dd8SToby Isaac const PetscInt face = fgeom->face[e][0]; 434*20cf1dd8SToby Isaac 435*20cf1dd8SToby Isaac ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 436*20cf1dd8SToby Isaac if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 437*20cf1dd8SToby Isaac ierr = PetscMemzero(f0, Nq*NcI* sizeof(PetscScalar));CHKERRQ(ierr); 438*20cf1dd8SToby Isaac ierr = PetscMemzero(f1, Nq*NcI*dim * sizeof(PetscScalar));CHKERRQ(ierr); 439*20cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 440*20cf1dd8SToby Isaac const PetscReal *v; 441*20cf1dd8SToby Isaac const PetscReal *invJ; 442*20cf1dd8SToby Isaac const PetscReal *n; 443*20cf1dd8SToby Isaac PetscReal detJ; 444*20cf1dd8SToby Isaac if (isAffine) { 445*20cf1dd8SToby Isaac CoordinatesRefToReal(dE, dim-1, fgeom->xi, v0, J, &quadPoints[q*(dim-1)], x); 446*20cf1dd8SToby Isaac v = x; 447*20cf1dd8SToby Isaac invJ = &fgeom->suppInvJ[0][e*dE*dE]; 448*20cf1dd8SToby Isaac detJ = fgeom->detJ[e]; 449*20cf1dd8SToby Isaac n = &fgeom->n[e*dE]; 450*20cf1dd8SToby Isaac } else { 451*20cf1dd8SToby Isaac v = &v0[q*dE]; 452*20cf1dd8SToby Isaac invJ = &fgeom->suppInvJ[0][(e*Np+q)*dE*dE]; 453*20cf1dd8SToby Isaac detJ = fgeom->detJ[e*Np + q]; 454*20cf1dd8SToby Isaac n = &fgeom->n[(e*Np+q)*dE]; 455*20cf1dd8SToby Isaac } 456*20cf1dd8SToby Isaac if (debug > 1 && q < Np) { 457*20cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", detJ);CHKERRQ(ierr); 458*20cf1dd8SToby Isaac #ifndef PETSC_USE_COMPLEX 459*20cf1dd8SToby Isaac ierr = DMPrintCellMatrix(e, "invJ", dim, dim, invJ);CHKERRQ(ierr); 460*20cf1dd8SToby Isaac #endif 461*20cf1dd8SToby Isaac } 462*20cf1dd8SToby Isaac if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 463*20cf1dd8SToby Isaac EvaluateFieldJets(dim, Nf, Nb, Nc, face*Nq+q, B, D, refSpaceDer, invJ, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t); 464*20cf1dd8SToby Isaac if (probAux) EvaluateFieldJets(dim, NfAux, NbAux, NcAux, face*Nq+q, BAux, DAux, refSpaceDerAux, invJ, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL); 465*20cf1dd8SToby Isaac if (f0_func) f0_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, v, n, numConstants, constants, &f0[q*NcI]); 466*20cf1dd8SToby Isaac if (f1_func) { 467*20cf1dd8SToby Isaac ierr = PetscMemzero(refSpaceDer, NcI*dim * sizeof(PetscScalar));CHKERRQ(ierr); 468*20cf1dd8SToby Isaac f1_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, v, n, numConstants, constants, refSpaceDer); 469*20cf1dd8SToby Isaac } 470*20cf1dd8SToby Isaac TransformF(dim, NcI, q, invJ, detJ, quadWeights, refSpaceDer, f0_func ? f0 : NULL, f1_func ? f1 : NULL); 471*20cf1dd8SToby Isaac } 472*20cf1dd8SToby Isaac UpdateElementVec(dim, Nq, NbI, NcI, &BI[face*Nq*NbI*NcI], &DI[face*Nq*NbI*NcI*dim], f0, f1, &elemVec[cOffset+fOffset]); 473*20cf1dd8SToby Isaac cOffset += totDim; 474*20cf1dd8SToby Isaac cOffsetAux += totDimAux; 475*20cf1dd8SToby Isaac } 476*20cf1dd8SToby Isaac PetscFunctionReturn(0); 477*20cf1dd8SToby Isaac } 478*20cf1dd8SToby Isaac 479*20cf1dd8SToby Isaac PetscErrorCode PetscFEIntegrateJacobian_Basic(PetscFE fem, PetscDS prob, PetscFEJacobianType jtype, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *geom, 480*20cf1dd8SToby Isaac const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 481*20cf1dd8SToby Isaac { 482*20cf1dd8SToby Isaac const PetscInt debug = 0; 483*20cf1dd8SToby Isaac PetscPointJac g0_func; 484*20cf1dd8SToby Isaac PetscPointJac g1_func; 485*20cf1dd8SToby Isaac PetscPointJac g2_func; 486*20cf1dd8SToby Isaac PetscPointJac g3_func; 487*20cf1dd8SToby Isaac PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 488*20cf1dd8SToby Isaac PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 489*20cf1dd8SToby Isaac PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 490*20cf1dd8SToby Isaac PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 491*20cf1dd8SToby Isaac PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 492*20cf1dd8SToby Isaac PetscQuadrature quad; 493*20cf1dd8SToby Isaac PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *refSpaceDer, *refSpaceDerAux; 494*20cf1dd8SToby Isaac const PetscScalar *constants; 495*20cf1dd8SToby Isaac PetscReal *x; 496*20cf1dd8SToby Isaac PetscReal **B, **D, **BAux = NULL, **DAux = NULL, *BI, *DI, *BJ, *DJ; 497*20cf1dd8SToby Isaac PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL; 498*20cf1dd8SToby Isaac PetscInt NbI = 0, NcI = 0, NbJ = 0, NcJ = 0; 499*20cf1dd8SToby Isaac PetscInt dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, e; 500*20cf1dd8SToby Isaac PetscInt dE, Np; 501*20cf1dd8SToby Isaac PetscBool isAffine; 502*20cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 503*20cf1dd8SToby Isaac PetscInt qNc, Nq, q; 504*20cf1dd8SToby Isaac PetscErrorCode ierr; 505*20cf1dd8SToby Isaac 506*20cf1dd8SToby Isaac PetscFunctionBegin; 507*20cf1dd8SToby Isaac ierr = PetscFEGetSpatialDimension(fem, &dim);CHKERRQ(ierr); 508*20cf1dd8SToby Isaac ierr = PetscFEGetQuadrature(fem, &quad);CHKERRQ(ierr); 509*20cf1dd8SToby Isaac ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 510*20cf1dd8SToby Isaac ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 511*20cf1dd8SToby Isaac ierr = PetscDSGetDimensions(prob, &Nb);CHKERRQ(ierr); 512*20cf1dd8SToby Isaac ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr); 513*20cf1dd8SToby Isaac ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr); 514*20cf1dd8SToby Isaac ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr); 515*20cf1dd8SToby Isaac switch(jtype) { 516*20cf1dd8SToby Isaac case PETSCFE_JACOBIAN_DYN: ierr = PetscDSGetDynamicJacobian(prob, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);break; 517*20cf1dd8SToby Isaac case PETSCFE_JACOBIAN_PRE: ierr = PetscDSGetJacobianPreconditioner(prob, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);break; 518*20cf1dd8SToby Isaac case PETSCFE_JACOBIAN: ierr = PetscDSGetJacobian(prob, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr);break; 519*20cf1dd8SToby Isaac } 520*20cf1dd8SToby Isaac if (!g0_func && !g1_func && !g2_func && !g3_func) PetscFunctionReturn(0); 521*20cf1dd8SToby Isaac ierr = PetscDSGetEvaluationArrays(prob, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 522*20cf1dd8SToby Isaac ierr = PetscDSGetRefCoordArrays(prob, &x, &refSpaceDer);CHKERRQ(ierr); 523*20cf1dd8SToby Isaac ierr = PetscDSGetWeakFormArrays(prob, NULL, NULL, &g0, &g1, &g2, &g3);CHKERRQ(ierr); 524*20cf1dd8SToby Isaac ierr = PetscDSGetTabulation(prob, &B, &D);CHKERRQ(ierr); 525*20cf1dd8SToby Isaac ierr = PetscDSGetFieldOffset(prob, fieldI, &offsetI);CHKERRQ(ierr); 526*20cf1dd8SToby Isaac ierr = PetscDSGetFieldOffset(prob, fieldJ, &offsetJ);CHKERRQ(ierr); 527*20cf1dd8SToby Isaac ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr); 528*20cf1dd8SToby Isaac if (probAux) { 529*20cf1dd8SToby Isaac ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr); 530*20cf1dd8SToby Isaac ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr); 531*20cf1dd8SToby Isaac ierr = PetscDSGetDimensions(probAux, &NbAux);CHKERRQ(ierr); 532*20cf1dd8SToby Isaac ierr = PetscDSGetComponents(probAux, &NcAux);CHKERRQ(ierr); 533*20cf1dd8SToby Isaac ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr); 534*20cf1dd8SToby Isaac ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr); 535*20cf1dd8SToby Isaac ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);CHKERRQ(ierr); 536*20cf1dd8SToby Isaac ierr = PetscDSGetRefCoordArrays(probAux, NULL, &refSpaceDerAux);CHKERRQ(ierr); 537*20cf1dd8SToby Isaac ierr = PetscDSGetTabulation(probAux, &BAux, &DAux);CHKERRQ(ierr); 538*20cf1dd8SToby Isaac } 539*20cf1dd8SToby Isaac NbI = Nb[fieldI], NbJ = Nb[fieldJ]; 540*20cf1dd8SToby Isaac NcI = Nc[fieldI], NcJ = Nc[fieldJ]; 541*20cf1dd8SToby Isaac BI = B[fieldI], BJ = B[fieldJ]; 542*20cf1dd8SToby Isaac DI = D[fieldI], DJ = D[fieldJ]; 543*20cf1dd8SToby Isaac /* Initialize here in case the function is not defined */ 544*20cf1dd8SToby Isaac ierr = PetscMemzero(g0, NcI*NcJ * sizeof(PetscScalar));CHKERRQ(ierr); 545*20cf1dd8SToby Isaac ierr = PetscMemzero(g1, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr); 546*20cf1dd8SToby Isaac ierr = PetscMemzero(g2, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr); 547*20cf1dd8SToby Isaac ierr = PetscMemzero(g3, NcI*NcJ*dim*dim * sizeof(PetscScalar));CHKERRQ(ierr); 548*20cf1dd8SToby Isaac ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 549*20cf1dd8SToby Isaac Np = geom->numPoints; 550*20cf1dd8SToby Isaac dE = geom->dimEmbed; 551*20cf1dd8SToby Isaac isAffine = geom->isAffine; 552*20cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 553*20cf1dd8SToby Isaac const PetscReal *v0 = &geom->v[e*Np*dE]; 554*20cf1dd8SToby Isaac const PetscReal *J = &geom->J[e*Np*dE*dE]; 555*20cf1dd8SToby Isaac 556*20cf1dd8SToby Isaac ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 557*20cf1dd8SToby Isaac if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 558*20cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 559*20cf1dd8SToby Isaac const PetscReal *v; 560*20cf1dd8SToby Isaac const PetscReal *invJ; 561*20cf1dd8SToby Isaac PetscReal detJ; 562*20cf1dd8SToby Isaac const PetscReal *BIq = &BI[q*NbI*NcI], *BJq = &BJ[q*NbJ*NcJ]; 563*20cf1dd8SToby Isaac const PetscReal *DIq = &DI[q*NbI*NcI*dim], *DJq = &DJ[q*NbJ*NcJ*dim]; 564*20cf1dd8SToby Isaac PetscReal w; 565*20cf1dd8SToby Isaac PetscInt f, g, fc, gc, c; 566*20cf1dd8SToby Isaac 567*20cf1dd8SToby Isaac if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 568*20cf1dd8SToby Isaac if (isAffine) { 569*20cf1dd8SToby Isaac CoordinatesRefToReal(dE, dim, geom->xi, v0, J, &quadPoints[q*dim], x); 570*20cf1dd8SToby Isaac v = x; 571*20cf1dd8SToby Isaac invJ = &geom->invJ[e*dE*dE]; 572*20cf1dd8SToby Isaac detJ = geom->detJ[e]; 573*20cf1dd8SToby Isaac } else { 574*20cf1dd8SToby Isaac v = &v0[q*dE]; 575*20cf1dd8SToby Isaac invJ = &geom->invJ[(e*Np+q)*dE*dE]; 576*20cf1dd8SToby Isaac detJ = geom->detJ[e*Np + q]; 577*20cf1dd8SToby Isaac } 578*20cf1dd8SToby Isaac w = detJ*quadWeights[q]; 579*20cf1dd8SToby Isaac EvaluateFieldJets(dim, Nf, Nb, Nc, q, B, D, refSpaceDer, invJ, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t); 580*20cf1dd8SToby Isaac if (probAux) EvaluateFieldJets(dim, NfAux, NbAux, NcAux, q, BAux, DAux, refSpaceDerAux, invJ, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL); 581*20cf1dd8SToby Isaac if (g0_func) { 582*20cf1dd8SToby Isaac ierr = PetscMemzero(g0, NcI*NcJ * sizeof(PetscScalar));CHKERRQ(ierr); 583*20cf1dd8SToby Isaac g0_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, v, numConstants, constants, g0); 584*20cf1dd8SToby Isaac for (c = 0; c < NcI*NcJ; ++c) g0[c] *= w; 585*20cf1dd8SToby Isaac } 586*20cf1dd8SToby Isaac if (g1_func) { 587*20cf1dd8SToby Isaac PetscInt d, d2; 588*20cf1dd8SToby Isaac ierr = PetscMemzero(refSpaceDer, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr); 589*20cf1dd8SToby Isaac g1_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, v, numConstants, constants, refSpaceDer); 590*20cf1dd8SToby Isaac for (fc = 0; fc < NcI; ++fc) { 591*20cf1dd8SToby Isaac for (gc = 0; gc < NcJ; ++gc) { 592*20cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 593*20cf1dd8SToby Isaac g1[(fc*NcJ+gc)*dim+d] = 0.0; 594*20cf1dd8SToby Isaac for (d2 = 0; d2 < dim; ++d2) g1[(fc*NcJ+gc)*dim+d] += invJ[d*dim+d2]*refSpaceDer[(fc*NcJ+gc)*dim+d2]; 595*20cf1dd8SToby Isaac g1[(fc*NcJ+gc)*dim+d] *= w; 596*20cf1dd8SToby Isaac } 597*20cf1dd8SToby Isaac } 598*20cf1dd8SToby Isaac } 599*20cf1dd8SToby Isaac } 600*20cf1dd8SToby Isaac if (g2_func) { 601*20cf1dd8SToby Isaac PetscInt d, d2; 602*20cf1dd8SToby Isaac ierr = PetscMemzero(refSpaceDer, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr); 603*20cf1dd8SToby Isaac g2_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, v, numConstants, constants, refSpaceDer); 604*20cf1dd8SToby Isaac for (fc = 0; fc < NcI; ++fc) { 605*20cf1dd8SToby Isaac for (gc = 0; gc < NcJ; ++gc) { 606*20cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 607*20cf1dd8SToby Isaac g2[(fc*NcJ+gc)*dim+d] = 0.0; 608*20cf1dd8SToby Isaac for (d2 = 0; d2 < dim; ++d2) g2[(fc*NcJ+gc)*dim+d] += invJ[d*dim+d2]*refSpaceDer[(fc*NcJ+gc)*dim+d2]; 609*20cf1dd8SToby Isaac g2[(fc*NcJ+gc)*dim+d] *= w; 610*20cf1dd8SToby Isaac } 611*20cf1dd8SToby Isaac } 612*20cf1dd8SToby Isaac } 613*20cf1dd8SToby Isaac } 614*20cf1dd8SToby Isaac if (g3_func) { 615*20cf1dd8SToby Isaac PetscInt d, d2, dp, d3; 616*20cf1dd8SToby Isaac ierr = PetscMemzero(refSpaceDer, NcI*NcJ*dim*dim * sizeof(PetscScalar));CHKERRQ(ierr); 617*20cf1dd8SToby Isaac g3_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, v, numConstants, constants, refSpaceDer); 618*20cf1dd8SToby Isaac for (fc = 0; fc < NcI; ++fc) { 619*20cf1dd8SToby Isaac for (gc = 0; gc < NcJ; ++gc) { 620*20cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 621*20cf1dd8SToby Isaac for (dp = 0; dp < dim; ++dp) { 622*20cf1dd8SToby Isaac g3[((fc*NcJ+gc)*dim+d)*dim+dp] = 0.0; 623*20cf1dd8SToby Isaac for (d2 = 0; d2 < dim; ++d2) { 624*20cf1dd8SToby Isaac for (d3 = 0; d3 < dim; ++d3) { 625*20cf1dd8SToby Isaac g3[((fc*NcJ+gc)*dim+d)*dim+dp] += invJ[d*dim+d2]*refSpaceDer[((fc*NcJ+gc)*dim+d2)*dim+d3]*invJ[dp*dim+d3]; 626*20cf1dd8SToby Isaac } 627*20cf1dd8SToby Isaac } 628*20cf1dd8SToby Isaac g3[((fc*NcJ+gc)*dim+d)*dim+dp] *= w; 629*20cf1dd8SToby Isaac } 630*20cf1dd8SToby Isaac } 631*20cf1dd8SToby Isaac } 632*20cf1dd8SToby Isaac } 633*20cf1dd8SToby Isaac } 634*20cf1dd8SToby Isaac 635*20cf1dd8SToby Isaac for (f = 0; f < NbI; ++f) { 636*20cf1dd8SToby Isaac for (fc = 0; fc < NcI; ++fc) { 637*20cf1dd8SToby Isaac const PetscInt fidx = f*NcI+fc; /* Test function basis index */ 638*20cf1dd8SToby Isaac const PetscInt i = offsetI+f; /* Element matrix row */ 639*20cf1dd8SToby Isaac for (g = 0; g < NbJ; ++g) { 640*20cf1dd8SToby Isaac for (gc = 0; gc < NcJ; ++gc) { 641*20cf1dd8SToby Isaac const PetscInt gidx = g*NcJ+gc; /* Trial function basis index */ 642*20cf1dd8SToby Isaac const PetscInt j = offsetJ+g; /* Element matrix column */ 643*20cf1dd8SToby Isaac const PetscInt fOff = eOffset+i*totDim+j; 644*20cf1dd8SToby Isaac PetscInt d, d2; 645*20cf1dd8SToby Isaac 646*20cf1dd8SToby Isaac elemMat[fOff] += BIq[fidx]*g0[fc*NcJ+gc]*BJq[gidx]; 647*20cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 648*20cf1dd8SToby Isaac elemMat[fOff] += BIq[fidx]*g1[(fc*NcJ+gc)*dim+d]*DJq[gidx*dim+d]; 649*20cf1dd8SToby Isaac elemMat[fOff] += DIq[fidx*dim+d]*g2[(fc*NcJ+gc)*dim+d]*BJq[gidx]; 650*20cf1dd8SToby Isaac for (d2 = 0; d2 < dim; ++d2) { 651*20cf1dd8SToby Isaac elemMat[fOff] += DIq[fidx*dim+d]*g3[((fc*NcJ+gc)*dim+d)*dim+d2]*DJq[gidx*dim+d2]; 652*20cf1dd8SToby Isaac } 653*20cf1dd8SToby Isaac } 654*20cf1dd8SToby Isaac } 655*20cf1dd8SToby Isaac } 656*20cf1dd8SToby Isaac } 657*20cf1dd8SToby Isaac } 658*20cf1dd8SToby Isaac } 659*20cf1dd8SToby Isaac if (debug > 1) { 660*20cf1dd8SToby Isaac PetscInt fc, f, gc, g; 661*20cf1dd8SToby Isaac 662*20cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);CHKERRQ(ierr); 663*20cf1dd8SToby Isaac for (fc = 0; fc < NcI; ++fc) { 664*20cf1dd8SToby Isaac for (f = 0; f < NbI; ++f) { 665*20cf1dd8SToby Isaac const PetscInt i = offsetI + f*NcI+fc; 666*20cf1dd8SToby Isaac for (gc = 0; gc < NcJ; ++gc) { 667*20cf1dd8SToby Isaac for (g = 0; g < NbJ; ++g) { 668*20cf1dd8SToby Isaac const PetscInt j = offsetJ + g*NcJ+gc; 669*20cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, " elemMat[%d,%d,%d,%d]: %g\n", f, fc, g, gc, PetscRealPart(elemMat[eOffset+i*totDim+j]));CHKERRQ(ierr); 670*20cf1dd8SToby Isaac } 671*20cf1dd8SToby Isaac } 672*20cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr); 673*20cf1dd8SToby Isaac } 674*20cf1dd8SToby Isaac } 675*20cf1dd8SToby Isaac } 676*20cf1dd8SToby Isaac cOffset += totDim; 677*20cf1dd8SToby Isaac cOffsetAux += totDimAux; 678*20cf1dd8SToby Isaac eOffset += PetscSqr(totDim); 679*20cf1dd8SToby Isaac } 680*20cf1dd8SToby Isaac PetscFunctionReturn(0); 681*20cf1dd8SToby Isaac } 682*20cf1dd8SToby Isaac 683*20cf1dd8SToby Isaac PetscErrorCode PetscFEIntegrateBdJacobian_Basic(PetscFE fem, PetscDS prob, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEGeom *fgeom, 684*20cf1dd8SToby Isaac const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[]) 685*20cf1dd8SToby Isaac { 686*20cf1dd8SToby Isaac const PetscInt debug = 0; 687*20cf1dd8SToby Isaac PetscBdPointJac g0_func; 688*20cf1dd8SToby Isaac PetscBdPointJac g1_func; 689*20cf1dd8SToby Isaac PetscBdPointJac g2_func; 690*20cf1dd8SToby Isaac PetscBdPointJac g3_func; 691*20cf1dd8SToby Isaac PetscInt cOffset = 0; /* Offset into coefficients[] for element e */ 692*20cf1dd8SToby Isaac PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */ 693*20cf1dd8SToby Isaac PetscInt eOffset = 0; /* Offset into elemMat[] for element e */ 694*20cf1dd8SToby Isaac PetscInt offsetI = 0; /* Offset into an element vector for fieldI */ 695*20cf1dd8SToby Isaac PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */ 696*20cf1dd8SToby Isaac PetscQuadrature quad; 697*20cf1dd8SToby Isaac PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *refSpaceDer, *refSpaceDerAux; 698*20cf1dd8SToby Isaac const PetscScalar *constants; 699*20cf1dd8SToby Isaac PetscReal *x; 700*20cf1dd8SToby Isaac PetscReal **B, **D, **BAux = NULL, **DAux = NULL, *BI, *DI, *BJ, *DJ; 701*20cf1dd8SToby Isaac PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL; 702*20cf1dd8SToby Isaac PetscInt NbI = 0, NcI = 0, NbJ = 0, NcJ = 0; 703*20cf1dd8SToby Isaac PetscInt dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, e; 704*20cf1dd8SToby Isaac PetscBool isAffine; 705*20cf1dd8SToby Isaac const PetscReal *quadPoints, *quadWeights; 706*20cf1dd8SToby Isaac PetscInt qNc, Nq, q, Np, dE; 707*20cf1dd8SToby Isaac PetscErrorCode ierr; 708*20cf1dd8SToby Isaac 709*20cf1dd8SToby Isaac PetscFunctionBegin; 710*20cf1dd8SToby Isaac ierr = PetscFEGetSpatialDimension(fem, &dim);CHKERRQ(ierr); 711*20cf1dd8SToby Isaac ierr = PetscFEGetFaceQuadrature(fem, &quad);CHKERRQ(ierr); 712*20cf1dd8SToby Isaac ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 713*20cf1dd8SToby Isaac ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 714*20cf1dd8SToby Isaac ierr = PetscDSGetDimensions(prob, &Nb);CHKERRQ(ierr); 715*20cf1dd8SToby Isaac ierr = PetscDSGetComponents(prob, &Nc);CHKERRQ(ierr); 716*20cf1dd8SToby Isaac ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr); 717*20cf1dd8SToby Isaac ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr); 718*20cf1dd8SToby Isaac ierr = PetscDSGetFieldOffset(prob, fieldI, &offsetI);CHKERRQ(ierr); 719*20cf1dd8SToby Isaac ierr = PetscDSGetFieldOffset(prob, fieldJ, &offsetJ);CHKERRQ(ierr); 720*20cf1dd8SToby Isaac ierr = PetscDSGetBdJacobian(prob, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);CHKERRQ(ierr); 721*20cf1dd8SToby Isaac ierr = PetscDSGetEvaluationArrays(prob, &u, coefficients_t ? &u_t : NULL, &u_x);CHKERRQ(ierr); 722*20cf1dd8SToby Isaac ierr = PetscDSGetRefCoordArrays(prob, &x, &refSpaceDer);CHKERRQ(ierr); 723*20cf1dd8SToby Isaac ierr = PetscDSGetWeakFormArrays(prob, NULL, NULL, &g0, &g1, &g2, &g3);CHKERRQ(ierr); 724*20cf1dd8SToby Isaac ierr = PetscDSGetFaceTabulation(prob, &B, &D);CHKERRQ(ierr); 725*20cf1dd8SToby Isaac ierr = PetscDSGetConstants(prob, &numConstants, &constants);CHKERRQ(ierr); 726*20cf1dd8SToby Isaac if (probAux) { 727*20cf1dd8SToby Isaac ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr); 728*20cf1dd8SToby Isaac ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr); 729*20cf1dd8SToby Isaac ierr = PetscDSGetDimensions(probAux, &NbAux);CHKERRQ(ierr); 730*20cf1dd8SToby Isaac ierr = PetscDSGetComponents(probAux, &NcAux);CHKERRQ(ierr); 731*20cf1dd8SToby Isaac ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr); 732*20cf1dd8SToby Isaac ierr = PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);CHKERRQ(ierr); 733*20cf1dd8SToby Isaac ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);CHKERRQ(ierr); 734*20cf1dd8SToby Isaac ierr = PetscDSGetRefCoordArrays(probAux, NULL, &refSpaceDerAux);CHKERRQ(ierr); 735*20cf1dd8SToby Isaac ierr = PetscDSGetFaceTabulation(probAux, &BAux, &DAux);CHKERRQ(ierr); 736*20cf1dd8SToby Isaac } 737*20cf1dd8SToby Isaac NbI = Nb[fieldI], NbJ = Nb[fieldJ]; 738*20cf1dd8SToby Isaac NcI = Nc[fieldI], NcJ = Nc[fieldJ]; 739*20cf1dd8SToby Isaac BI = B[fieldI], BJ = B[fieldJ]; 740*20cf1dd8SToby Isaac DI = D[fieldI], DJ = D[fieldJ]; 741*20cf1dd8SToby Isaac /* Initialize here in case the function is not defined */ 742*20cf1dd8SToby Isaac ierr = PetscMemzero(g0, NcI*NcJ * sizeof(PetscScalar));CHKERRQ(ierr); 743*20cf1dd8SToby Isaac ierr = PetscMemzero(g1, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr); 744*20cf1dd8SToby Isaac ierr = PetscMemzero(g2, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr); 745*20cf1dd8SToby Isaac ierr = PetscMemzero(g3, NcI*NcJ*dim*dim * sizeof(PetscScalar));CHKERRQ(ierr); 746*20cf1dd8SToby Isaac ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 747*20cf1dd8SToby Isaac Np = fgeom->numPoints; 748*20cf1dd8SToby Isaac dE = fgeom->dimEmbed; 749*20cf1dd8SToby Isaac isAffine = fgeom->isAffine; 750*20cf1dd8SToby Isaac for (e = 0; e < Ne; ++e) { 751*20cf1dd8SToby Isaac const PetscReal *v0 = &fgeom->v[e*Np*dE]; 752*20cf1dd8SToby Isaac const PetscReal *J = &fgeom->J[e*Np*dE*dE]; 753*20cf1dd8SToby Isaac const PetscInt face = fgeom->face[e][0]; 754*20cf1dd8SToby Isaac 755*20cf1dd8SToby Isaac ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 756*20cf1dd8SToby Isaac if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc); 757*20cf1dd8SToby Isaac for (q = 0; q < Nq; ++q) { 758*20cf1dd8SToby Isaac const PetscReal *BIq = &BI[(face*Nq+q)*NbI*NcI], *BJq = &BJ[(face*Nq+q)*NbJ*NcJ]; 759*20cf1dd8SToby Isaac const PetscReal *DIq = &DI[(face*Nq+q)*NbI*NcI*dim], *DJq = &DJ[(face*Nq+q)*NbJ*NcJ*dim]; 760*20cf1dd8SToby Isaac PetscReal w; 761*20cf1dd8SToby Isaac PetscInt f, g, fc, gc, c; 762*20cf1dd8SToby Isaac const PetscReal *v; 763*20cf1dd8SToby Isaac const PetscReal *invJ; 764*20cf1dd8SToby Isaac const PetscReal *n; 765*20cf1dd8SToby Isaac PetscReal detJ; 766*20cf1dd8SToby Isaac 767*20cf1dd8SToby Isaac if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);CHKERRQ(ierr);} 768*20cf1dd8SToby Isaac if (isAffine) { 769*20cf1dd8SToby Isaac CoordinatesRefToReal(dE, dim-1, fgeom->xi, v0, J, &quadPoints[q*(dim-1)], x); 770*20cf1dd8SToby Isaac v = x; 771*20cf1dd8SToby Isaac invJ = &fgeom->suppInvJ[0][e*dE*dE]; 772*20cf1dd8SToby Isaac detJ = fgeom->detJ[e]; 773*20cf1dd8SToby Isaac n = &fgeom->n[e*dE]; 774*20cf1dd8SToby Isaac } else { 775*20cf1dd8SToby Isaac v = &v0[q*dE]; 776*20cf1dd8SToby Isaac invJ = &fgeom->suppInvJ[0][(e*Np+q)*dE*dE]; 777*20cf1dd8SToby Isaac detJ = fgeom->detJ[e*Np + q]; 778*20cf1dd8SToby Isaac n = &fgeom->n[(e*Np+q)*dE]; 779*20cf1dd8SToby Isaac } 780*20cf1dd8SToby Isaac w = detJ*quadWeights[q]; 781*20cf1dd8SToby Isaac 782*20cf1dd8SToby Isaac EvaluateFieldJets(dim, Nf, Nb, Nc, face*Nq+q, B, D, refSpaceDer, invJ, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t); 783*20cf1dd8SToby Isaac if (probAux) EvaluateFieldJets(dim, NfAux, NbAux, NcAux, face*Nq+q, BAux, DAux, refSpaceDerAux, invJ, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL); 784*20cf1dd8SToby Isaac if (g0_func) { 785*20cf1dd8SToby Isaac ierr = PetscMemzero(g0, NcI*NcJ * sizeof(PetscScalar));CHKERRQ(ierr); 786*20cf1dd8SToby Isaac g0_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, v, n, numConstants, constants, g0); 787*20cf1dd8SToby Isaac for (c = 0; c < NcI*NcJ; ++c) g0[c] *= w; 788*20cf1dd8SToby Isaac } 789*20cf1dd8SToby Isaac if (g1_func) { 790*20cf1dd8SToby Isaac PetscInt d, d2; 791*20cf1dd8SToby Isaac ierr = PetscMemzero(refSpaceDer, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr); 792*20cf1dd8SToby Isaac g1_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, v, n, numConstants, constants, refSpaceDer); 793*20cf1dd8SToby Isaac for (fc = 0; fc < NcI; ++fc) { 794*20cf1dd8SToby Isaac for (gc = 0; gc < NcJ; ++gc) { 795*20cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 796*20cf1dd8SToby Isaac g1[(fc*NcJ+gc)*dim+d] = 0.0; 797*20cf1dd8SToby Isaac for (d2 = 0; d2 < dim; ++d2) g1[(fc*NcJ+gc)*dim+d] += invJ[d*dim+d2]*refSpaceDer[(fc*NcJ+gc)*dim+d2]; 798*20cf1dd8SToby Isaac g1[(fc*NcJ+gc)*dim+d] *= w; 799*20cf1dd8SToby Isaac } 800*20cf1dd8SToby Isaac } 801*20cf1dd8SToby Isaac } 802*20cf1dd8SToby Isaac } 803*20cf1dd8SToby Isaac if (g2_func) { 804*20cf1dd8SToby Isaac PetscInt d, d2; 805*20cf1dd8SToby Isaac ierr = PetscMemzero(refSpaceDer, NcI*NcJ*dim * sizeof(PetscScalar));CHKERRQ(ierr); 806*20cf1dd8SToby Isaac g2_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, v, n, numConstants, constants, refSpaceDer); 807*20cf1dd8SToby Isaac for (fc = 0; fc < NcI; ++fc) { 808*20cf1dd8SToby Isaac for (gc = 0; gc < NcJ; ++gc) { 809*20cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 810*20cf1dd8SToby Isaac g2[(fc*NcJ+gc)*dim+d] = 0.0; 811*20cf1dd8SToby Isaac for (d2 = 0; d2 < dim; ++d2) g2[(fc*NcJ+gc)*dim+d] += invJ[d*dim+d2]*refSpaceDer[(fc*NcJ+gc)*dim+d2]; 812*20cf1dd8SToby Isaac g2[(fc*NcJ+gc)*dim+d] *= w; 813*20cf1dd8SToby Isaac } 814*20cf1dd8SToby Isaac } 815*20cf1dd8SToby Isaac } 816*20cf1dd8SToby Isaac } 817*20cf1dd8SToby Isaac if (g3_func) { 818*20cf1dd8SToby Isaac PetscInt d, d2, dp, d3; 819*20cf1dd8SToby Isaac ierr = PetscMemzero(refSpaceDer, NcI*NcJ*dim*dim * sizeof(PetscScalar));CHKERRQ(ierr); 820*20cf1dd8SToby Isaac g3_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, v, n, numConstants, constants, refSpaceDer); 821*20cf1dd8SToby Isaac for (fc = 0; fc < NcI; ++fc) { 822*20cf1dd8SToby Isaac for (gc = 0; gc < NcJ; ++gc) { 823*20cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 824*20cf1dd8SToby Isaac for (dp = 0; dp < dim; ++dp) { 825*20cf1dd8SToby Isaac g3[((fc*NcJ+gc)*dim+d)*dim+dp] = 0.0; 826*20cf1dd8SToby Isaac for (d2 = 0; d2 < dim; ++d2) { 827*20cf1dd8SToby Isaac for (d3 = 0; d3 < dim; ++d3) { 828*20cf1dd8SToby Isaac g3[((fc*NcJ+gc)*dim+d)*dim+dp] += invJ[d*dim+d2]*refSpaceDer[((fc*NcJ+gc)*dim+d2)*dim+d3]*invJ[dp*dim+d3]; 829*20cf1dd8SToby Isaac } 830*20cf1dd8SToby Isaac } 831*20cf1dd8SToby Isaac g3[((fc*NcJ+gc)*dim+d)*dim+dp] *= w; 832*20cf1dd8SToby Isaac } 833*20cf1dd8SToby Isaac } 834*20cf1dd8SToby Isaac } 835*20cf1dd8SToby Isaac } 836*20cf1dd8SToby Isaac } 837*20cf1dd8SToby Isaac 838*20cf1dd8SToby Isaac for (f = 0; f < NbI; ++f) { 839*20cf1dd8SToby Isaac for (fc = 0; fc < NcI; ++fc) { 840*20cf1dd8SToby Isaac const PetscInt fidx = f*NcI+fc; /* Test function basis index */ 841*20cf1dd8SToby Isaac const PetscInt i = offsetI+f; /* Element matrix row */ 842*20cf1dd8SToby Isaac for (g = 0; g < NbJ; ++g) { 843*20cf1dd8SToby Isaac for (gc = 0; gc < NcJ; ++gc) { 844*20cf1dd8SToby Isaac const PetscInt gidx = g*NcJ+gc; /* Trial function basis index */ 845*20cf1dd8SToby Isaac const PetscInt j = offsetJ+g; /* Element matrix column */ 846*20cf1dd8SToby Isaac const PetscInt fOff = eOffset+i*totDim+j; 847*20cf1dd8SToby Isaac PetscInt d, d2; 848*20cf1dd8SToby Isaac 849*20cf1dd8SToby Isaac elemMat[fOff] += BIq[fidx]*g0[fc*NcJ+gc]*BJq[gidx]; 850*20cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 851*20cf1dd8SToby Isaac elemMat[fOff] += BIq[fidx]*g1[(fc*NcJ+gc)*dim+d]*DJq[gidx*dim+d]; 852*20cf1dd8SToby Isaac elemMat[fOff] += DIq[fidx*dim+d]*g2[(fc*NcJ+gc)*dim+d]*BJq[gidx]; 853*20cf1dd8SToby Isaac for (d2 = 0; d2 < dim; ++d2) { 854*20cf1dd8SToby Isaac elemMat[fOff] += DIq[fidx*dim+d]*g3[((fc*NcJ+gc)*dim+d)*dim+d2]*DJq[gidx*dim+d2]; 855*20cf1dd8SToby Isaac } 856*20cf1dd8SToby Isaac } 857*20cf1dd8SToby Isaac } 858*20cf1dd8SToby Isaac } 859*20cf1dd8SToby Isaac } 860*20cf1dd8SToby Isaac } 861*20cf1dd8SToby Isaac } 862*20cf1dd8SToby Isaac if (debug > 1) { 863*20cf1dd8SToby Isaac PetscInt fc, f, gc, g; 864*20cf1dd8SToby Isaac 865*20cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);CHKERRQ(ierr); 866*20cf1dd8SToby Isaac for (fc = 0; fc < NcI; ++fc) { 867*20cf1dd8SToby Isaac for (f = 0; f < NbI; ++f) { 868*20cf1dd8SToby Isaac const PetscInt i = offsetI + f*NcI+fc; 869*20cf1dd8SToby Isaac for (gc = 0; gc < NcJ; ++gc) { 870*20cf1dd8SToby Isaac for (g = 0; g < NbJ; ++g) { 871*20cf1dd8SToby Isaac const PetscInt j = offsetJ + g*NcJ+gc; 872*20cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, " elemMat[%d,%d,%d,%d]: %g\n", f, fc, g, gc, PetscRealPart(elemMat[eOffset+i*totDim+j]));CHKERRQ(ierr); 873*20cf1dd8SToby Isaac } 874*20cf1dd8SToby Isaac } 875*20cf1dd8SToby Isaac ierr = PetscPrintf(PETSC_COMM_SELF, "\n");CHKERRQ(ierr); 876*20cf1dd8SToby Isaac } 877*20cf1dd8SToby Isaac } 878*20cf1dd8SToby Isaac } 879*20cf1dd8SToby Isaac cOffset += totDim; 880*20cf1dd8SToby Isaac cOffsetAux += totDimAux; 881*20cf1dd8SToby Isaac eOffset += PetscSqr(totDim); 882*20cf1dd8SToby Isaac } 883*20cf1dd8SToby Isaac PetscFunctionReturn(0); 884*20cf1dd8SToby Isaac } 885*20cf1dd8SToby Isaac 886*20cf1dd8SToby Isaac PetscErrorCode PetscFEInitialize_Basic(PetscFE fem) 887*20cf1dd8SToby Isaac { 888*20cf1dd8SToby Isaac PetscFunctionBegin; 889*20cf1dd8SToby Isaac fem->ops->setfromoptions = NULL; 890*20cf1dd8SToby Isaac fem->ops->setup = PetscFESetUp_Basic; 891*20cf1dd8SToby Isaac fem->ops->view = PetscFEView_Basic; 892*20cf1dd8SToby Isaac fem->ops->destroy = PetscFEDestroy_Basic; 893*20cf1dd8SToby Isaac fem->ops->getdimension = PetscFEGetDimension_Basic; 894*20cf1dd8SToby Isaac fem->ops->gettabulation = PetscFEGetTabulation_Basic; 895*20cf1dd8SToby Isaac fem->ops->integrate = PetscFEIntegrate_Basic; 896*20cf1dd8SToby Isaac fem->ops->integrateresidual = PetscFEIntegrateResidual_Basic; 897*20cf1dd8SToby Isaac fem->ops->integratebdresidual = PetscFEIntegrateBdResidual_Basic; 898*20cf1dd8SToby Isaac fem->ops->integratejacobianaction = NULL/* PetscFEIntegrateJacobianAction_Basic */; 899*20cf1dd8SToby Isaac fem->ops->integratejacobian = PetscFEIntegrateJacobian_Basic; 900*20cf1dd8SToby Isaac fem->ops->integratebdjacobian = PetscFEIntegrateBdJacobian_Basic; 901*20cf1dd8SToby Isaac PetscFunctionReturn(0); 902*20cf1dd8SToby Isaac } 903*20cf1dd8SToby Isaac 904*20cf1dd8SToby Isaac /*MC 905*20cf1dd8SToby Isaac PETSCFEBASIC = "basic" - A PetscFE object that integrates with basic tiling and no vectorization 906*20cf1dd8SToby Isaac 907*20cf1dd8SToby Isaac Level: intermediate 908*20cf1dd8SToby Isaac 909*20cf1dd8SToby Isaac .seealso: PetscFEType, PetscFECreate(), PetscFESetType() 910*20cf1dd8SToby Isaac M*/ 911*20cf1dd8SToby Isaac 912*20cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreate_Basic(PetscFE fem) 913*20cf1dd8SToby Isaac { 914*20cf1dd8SToby Isaac PetscFE_Basic *b; 915*20cf1dd8SToby Isaac PetscErrorCode ierr; 916*20cf1dd8SToby Isaac 917*20cf1dd8SToby Isaac PetscFunctionBegin; 918*20cf1dd8SToby Isaac PetscValidHeaderSpecific(fem, PETSCFE_CLASSID, 1); 919*20cf1dd8SToby Isaac ierr = PetscNewLog(fem,&b);CHKERRQ(ierr); 920*20cf1dd8SToby Isaac fem->data = b; 921*20cf1dd8SToby Isaac 922*20cf1dd8SToby Isaac ierr = PetscFEInitialize_Basic(fem);CHKERRQ(ierr); 923*20cf1dd8SToby Isaac PetscFunctionReturn(0); 924*20cf1dd8SToby Isaac } 925*20cf1dd8SToby Isaac 926*20cf1dd8SToby Isaac 927