Actual source code: dtfe.c
petsc-3.9.4 2018-09-11
1: /* Basis Jet Tabulation
3: We would like to tabulate the nodal basis functions and derivatives at a set of points, usually quadrature points. We
4: follow here the derviation in http://www.math.ttu.edu/~kirby/papers/fiat-toms-2004.pdf. The nodal basis $\psi_i$ can
5: be expressed in terms of a prime basis $\phi_i$ which can be stably evaluated. In PETSc, we will use the Legendre basis
6: as a prime basis.
8: \psi_i = \sum_k \alpha_{ki} \phi_k
10: Our nodal basis is defined in terms of the dual basis $n_j$
12: n_j \cdot \psi_i = \delta_{ji}
14: and we may act on the first equation to obtain
16: n_j \cdot \psi_i = \sum_k \alpha_{ki} n_j \cdot \phi_k
17: \delta_{ji} = \sum_k \alpha_{ki} V_{jk}
18: I = V \alpha
20: so the coefficients of the nodal basis in the prime basis are
22: \alpha = V^{-1}
24: We will define the dual basis vectors $n_j$ using a quadrature rule.
26: Right now, we will just use the polynomial spaces P^k. I know some elements use the space of symmetric polynomials
27: (I think Nedelec), but we will neglect this for now. Constraints in the space, e.g. Arnold-Winther elements, can
28: be implemented exactly as in FIAT using functionals $L_j$.
30: I will have to count the degrees correctly for the Legendre product when we are on simplices.
32: We will have three objects:
33: - Space, P: this just need point evaluation I think
34: - Dual Space, P'+K: This looks like a set of functionals that can act on members of P, each n is defined by a Q
35: - FEM: This keeps {P, P', Q}
36: */
37: #include <petsc/private/petscfeimpl.h>
38: #include <petsc/private/dtimpl.h>
39: #include <petsc/private/dmpleximpl.h>
40: #include <petscdmshell.h>
41: #include <petscdmplex.h>
42: #include <petscblaslapack.h>
44: PetscBool FEcite = PETSC_FALSE;
45: const char FECitation[] = "@article{kirby2004,\n"
46: " title = {Algorithm 839: FIAT, a New Paradigm for Computing Finite Element Basis Functions},\n"
47: " journal = {ACM Transactions on Mathematical Software},\n"
48: " author = {Robert C. Kirby},\n"
49: " volume = {30},\n"
50: " number = {4},\n"
51: " pages = {502--516},\n"
52: " doi = {10.1145/1039813.1039820},\n"
53: " year = {2004}\n}\n";
55: PetscClassId PETSCSPACE_CLASSID = 0;
57: PetscFunctionList PetscSpaceList = NULL;
58: PetscBool PetscSpaceRegisterAllCalled = PETSC_FALSE;
60: /*@C
61: PetscSpaceRegister - Adds a new PetscSpace implementation
63: Not Collective
65: Input Parameters:
66: + name - The name of a new user-defined creation routine
67: - create_func - The creation routine itself
69: Notes:
70: PetscSpaceRegister() may be called multiple times to add several user-defined PetscSpaces
72: Sample usage:
73: .vb
74: PetscSpaceRegister("my_space", MyPetscSpaceCreate);
75: .ve
77: Then, your PetscSpace type can be chosen with the procedural interface via
78: .vb
79: PetscSpaceCreate(MPI_Comm, PetscSpace *);
80: PetscSpaceSetType(PetscSpace, "my_space");
81: .ve
82: or at runtime via the option
83: .vb
84: -petscspace_type my_space
85: .ve
87: Level: advanced
89: .keywords: PetscSpace, register
90: .seealso: PetscSpaceRegisterAll(), PetscSpaceRegisterDestroy()
92: @*/
93: PetscErrorCode PetscSpaceRegister(const char sname[], PetscErrorCode (*function)(PetscSpace))
94: {
98: PetscFunctionListAdd(&PetscSpaceList, sname, function);
99: return(0);
100: }
102: /*@C
103: PetscSpaceSetType - Builds a particular PetscSpace
105: Collective on PetscSpace
107: Input Parameters:
108: + sp - The PetscSpace object
109: - name - The kind of space
111: Options Database Key:
112: . -petscspace_type <type> - Sets the PetscSpace type; use -help for a list of available types
114: Level: intermediate
116: .keywords: PetscSpace, set, type
117: .seealso: PetscSpaceGetType(), PetscSpaceCreate()
118: @*/
119: PetscErrorCode PetscSpaceSetType(PetscSpace sp, PetscSpaceType name)
120: {
121: PetscErrorCode (*r)(PetscSpace);
122: PetscBool match;
127: PetscObjectTypeCompare((PetscObject) sp, name, &match);
128: if (match) return(0);
130: PetscSpaceRegisterAll();
131: PetscFunctionListFind(PetscSpaceList, name, &r);
132: if (!r) SETERRQ1(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscSpace type: %s", name);
134: if (sp->ops->destroy) {
135: (*sp->ops->destroy)(sp);
136: sp->ops->destroy = NULL;
137: }
138: (*r)(sp);
139: PetscObjectChangeTypeName((PetscObject) sp, name);
140: return(0);
141: }
143: /*@C
144: PetscSpaceGetType - Gets the PetscSpace type name (as a string) from the object.
146: Not Collective
148: Input Parameter:
149: . sp - The PetscSpace
151: Output Parameter:
152: . name - The PetscSpace type name
154: Level: intermediate
156: .keywords: PetscSpace, get, type, name
157: .seealso: PetscSpaceSetType(), PetscSpaceCreate()
158: @*/
159: PetscErrorCode PetscSpaceGetType(PetscSpace sp, PetscSpaceType *name)
160: {
166: if (!PetscSpaceRegisterAllCalled) {
167: PetscSpaceRegisterAll();
168: }
169: *name = ((PetscObject) sp)->type_name;
170: return(0);
171: }
173: /*@C
174: PetscSpaceView - Views a PetscSpace
176: Collective on PetscSpace
178: Input Parameter:
179: + sp - the PetscSpace object to view
180: - v - the viewer
182: Level: developer
184: .seealso PetscSpaceDestroy()
185: @*/
186: PetscErrorCode PetscSpaceView(PetscSpace sp, PetscViewer v)
187: {
192: if (!v) {PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) sp), &v);}
193: if (sp->ops->view) {(*sp->ops->view)(sp, v);}
194: return(0);
195: }
197: /*@
198: PetscSpaceSetFromOptions - sets parameters in a PetscSpace from the options database
200: Collective on PetscSpace
202: Input Parameter:
203: . sp - the PetscSpace object to set options for
205: Options Database:
206: . -petscspace_order the approximation order of the space
208: Level: developer
210: .seealso PetscSpaceView()
211: @*/
212: PetscErrorCode PetscSpaceSetFromOptions(PetscSpace sp)
213: {
214: const char *defaultType;
215: char name[256];
216: PetscBool flg;
221: if (!((PetscObject) sp)->type_name) {
222: defaultType = PETSCSPACEPOLYNOMIAL;
223: } else {
224: defaultType = ((PetscObject) sp)->type_name;
225: }
226: if (!PetscSpaceRegisterAllCalled) {PetscSpaceRegisterAll();}
228: PetscObjectOptionsBegin((PetscObject) sp);
229: PetscOptionsFList("-petscspace_type", "Linear space", "PetscSpaceSetType", PetscSpaceList, defaultType, name, 256, &flg);
230: if (flg) {
231: PetscSpaceSetType(sp, name);
232: } else if (!((PetscObject) sp)->type_name) {
233: PetscSpaceSetType(sp, defaultType);
234: }
235: PetscOptionsInt("-petscspace_order", "The approximation order", "PetscSpaceSetOrder", sp->order, &sp->order, NULL);
236: PetscOptionsInt("-petscspace_components", "The number of components", "PetscSpaceSetNumComponents", sp->Nc, &sp->Nc, NULL);
237: if (sp->ops->setfromoptions) {
238: (*sp->ops->setfromoptions)(PetscOptionsObject,sp);
239: }
240: /* process any options handlers added with PetscObjectAddOptionsHandler() */
241: PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) sp);
242: PetscOptionsEnd();
243: PetscSpaceViewFromOptions(sp, NULL, "-petscspace_view");
244: return(0);
245: }
247: /*@C
248: PetscSpaceSetUp - Construct data structures for the PetscSpace
250: Collective on PetscSpace
252: Input Parameter:
253: . sp - the PetscSpace object to setup
255: Level: developer
257: .seealso PetscSpaceView(), PetscSpaceDestroy()
258: @*/
259: PetscErrorCode PetscSpaceSetUp(PetscSpace sp)
260: {
265: if (sp->ops->setup) {(*sp->ops->setup)(sp);}
266: return(0);
267: }
269: /*@
270: PetscSpaceDestroy - Destroys a PetscSpace object
272: Collective on PetscSpace
274: Input Parameter:
275: . sp - the PetscSpace object to destroy
277: Level: developer
279: .seealso PetscSpaceView()
280: @*/
281: PetscErrorCode PetscSpaceDestroy(PetscSpace *sp)
282: {
286: if (!*sp) return(0);
289: if (--((PetscObject)(*sp))->refct > 0) {*sp = 0; return(0);}
290: ((PetscObject) (*sp))->refct = 0;
291: DMDestroy(&(*sp)->dm);
293: (*(*sp)->ops->destroy)(*sp);
294: PetscHeaderDestroy(sp);
295: return(0);
296: }
298: /*@
299: PetscSpaceCreate - Creates an empty PetscSpace object. The type can then be set with PetscSpaceSetType().
301: Collective on MPI_Comm
303: Input Parameter:
304: . comm - The communicator for the PetscSpace object
306: Output Parameter:
307: . sp - The PetscSpace object
309: Level: beginner
311: .seealso: PetscSpaceSetType(), PETSCSPACEPOLYNOMIAL
312: @*/
313: PetscErrorCode PetscSpaceCreate(MPI_Comm comm, PetscSpace *sp)
314: {
315: PetscSpace s;
320: PetscCitationsRegister(FECitation,&FEcite);
321: *sp = NULL;
322: PetscFEInitializePackage();
324: PetscHeaderCreate(s, PETSCSPACE_CLASSID, "PetscSpace", "Linear Space", "PetscSpace", comm, PetscSpaceDestroy, PetscSpaceView);
326: s->order = 0;
327: s->Nc = 1;
328: DMShellCreate(comm, &s->dm);
329: PetscSpaceSetType(s, PETSCSPACEPOLYNOMIAL);
331: *sp = s;
332: return(0);
333: }
335: /*@
336: PetscSpaceGetDimension - Return the dimension of this space, i.e. the number of basis vectors
338: Input Parameter:
339: . sp - The PetscSpace
341: Output Parameter:
342: . dim - The dimension
344: Level: intermediate
346: .seealso: PetscSpaceGetOrder(), PetscSpaceCreate(), PetscSpace
347: @*/
348: PetscErrorCode PetscSpaceGetDimension(PetscSpace sp, PetscInt *dim)
349: {
355: *dim = 0;
356: if (sp->ops->getdimension) {(*sp->ops->getdimension)(sp, dim);}
357: return(0);
358: }
360: /*@
361: PetscSpaceGetOrder - Return the order of approximation for this space
363: Input Parameter:
364: . sp - The PetscSpace
366: Output Parameter:
367: . order - The approximation order
369: Level: intermediate
371: .seealso: PetscSpaceSetOrder(), PetscSpaceGetDimension(), PetscSpaceCreate(), PetscSpace
372: @*/
373: PetscErrorCode PetscSpaceGetOrder(PetscSpace sp, PetscInt *order)
374: {
378: *order = sp->order;
379: return(0);
380: }
382: /*@
383: PetscSpaceSetOrder - Set the order of approximation for this space
385: Input Parameters:
386: + sp - The PetscSpace
387: - order - The approximation order
389: Level: intermediate
391: .seealso: PetscSpaceGetOrder(), PetscSpaceCreate(), PetscSpace
392: @*/
393: PetscErrorCode PetscSpaceSetOrder(PetscSpace sp, PetscInt order)
394: {
397: sp->order = order;
398: return(0);
399: }
401: /*@
402: PetscSpaceGetNumComponents - Return the number of components for this space
404: Input Parameter:
405: . sp - The PetscSpace
407: Output Parameter:
408: . Nc - The number of components
410: Note: A vector space, for example, will have d components, where d is the spatial dimension
412: Level: intermediate
414: .seealso: PetscSpaceSetNumComponents(), PetscSpaceGetDimension(), PetscSpaceCreate(), PetscSpace
415: @*/
416: PetscErrorCode PetscSpaceGetNumComponents(PetscSpace sp, PetscInt *Nc)
417: {
421: *Nc = sp->Nc;
422: return(0);
423: }
425: /*@
426: PetscSpaceSetNumComponents - Set the number of components for this space
428: Input Parameters:
429: + sp - The PetscSpace
430: - order - The number of components
432: Level: intermediate
434: .seealso: PetscSpaceGetNumComponents(), PetscSpaceCreate(), PetscSpace
435: @*/
436: PetscErrorCode PetscSpaceSetNumComponents(PetscSpace sp, PetscInt Nc)
437: {
440: sp->Nc = Nc;
441: return(0);
442: }
444: /*@C
445: PetscSpaceEvaluate - Evaluate the basis functions and their derivatives (jet) at each point
447: Input Parameters:
448: + sp - The PetscSpace
449: . npoints - The number of evaluation points, in reference coordinates
450: - points - The point coordinates
452: Output Parameters:
453: + B - The function evaluations in a npoints x nfuncs array
454: . D - The derivative evaluations in a npoints x nfuncs x dim array
455: - H - The second derivative evaluations in a npoints x nfuncs x dim x dim array
457: Note: Above nfuncs is the dimension of the space, and dim is the spatial dimension. The coordinates are given
458: on the reference cell, not in real space.
460: Level: advanced
462: .seealso: PetscFEGetTabulation(), PetscFEGetDefaultTabulation(), PetscSpaceCreate()
463: @*/
464: PetscErrorCode PetscSpaceEvaluate(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[])
465: {
469: if (!npoints) return(0);
475: if (sp->ops->evaluate) {(*sp->ops->evaluate)(sp, npoints, points, B, D, H);}
476: return(0);
477: }
479: /*@
480: PetscSpaceGetHeightSubspace - Get the subset of the primal space basis that is supported on a mesh point of a given height.
482: If the space is not defined on mesh points of the given height (e.g. if the space is discontinuous and
483: pointwise values are not defined on the element boundaries), or if the implementation of PetscSpace does not
484: support extracting subspaces, then NULL is returned.
486: This does not increment the reference count on the returned space, and the user should not destroy it.
488: Not collective
490: Input Parameters:
491: + sp - the PetscSpace object
492: - height - the height of the mesh point for which the subspace is desired
494: Output Parameter:
495: . subsp - the subspace
497: Level: advanced
499: .seealso: PetscDualSpaceGetHeightSubspace(), PetscSpace
500: @*/
501: PetscErrorCode PetscSpaceGetHeightSubspace(PetscSpace sp, PetscInt height, PetscSpace *subsp)
502: {
508: *subsp = NULL;
509: if (sp->ops->getheightsubspace) {
510: (*sp->ops->getheightsubspace)(sp, height, subsp);
511: }
512: return(0);
513: }
515: PetscErrorCode PetscSpaceSetFromOptions_Polynomial(PetscOptionItems *PetscOptionsObject,PetscSpace sp)
516: {
517: PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
518: PetscErrorCode ierr;
521: PetscOptionsHead(PetscOptionsObject,"PetscSpace polynomial options");
522: PetscOptionsInt("-petscspace_poly_num_variables", "The number of different variables, e.g. x and y", "PetscSpacePolynomialSetNumVariables", poly->numVariables, &poly->numVariables, NULL);
523: PetscOptionsBool("-petscspace_poly_sym", "Use only symmetric polynomials", "PetscSpacePolynomialSetSymmetric", poly->symmetric, &poly->symmetric, NULL);
524: PetscOptionsBool("-petscspace_poly_tensor", "Use the tensor product polynomials", "PetscSpacePolynomialSetTensor", poly->tensor, &poly->tensor, NULL);
525: PetscOptionsTail();
526: return(0);
527: }
529: static PetscErrorCode PetscSpacePolynomialView_Ascii(PetscSpace sp, PetscViewer viewer)
530: {
531: PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
532: PetscErrorCode ierr;
535: if (sp->Nc > 1) {
536: if (poly->tensor) {PetscViewerASCIIPrintf(viewer, "Tensor polynomial space in %D variables of degree %D with %D components\n", poly->numVariables, sp->order, sp->Nc);}
537: else {PetscViewerASCIIPrintf(viewer, "Polynomial space in %D variables of degree %D with %D components\n", poly->numVariables, sp->order, sp->Nc);}
538: } else {
539: if (poly->tensor) {PetscViewerASCIIPrintf(viewer, "Tensor polynomial space in %d variables of degree %d\n", poly->numVariables, sp->order);}
540: else {PetscViewerASCIIPrintf(viewer, "Polynomial space in %d variables of degree %d\n", poly->numVariables, sp->order);}
541: }
542: return(0);
543: }
545: PetscErrorCode PetscSpaceView_Polynomial(PetscSpace sp, PetscViewer viewer)
546: {
547: PetscBool iascii;
553: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
554: if (iascii) {PetscSpacePolynomialView_Ascii(sp, viewer);}
555: return(0);
556: }
558: PetscErrorCode PetscSpaceSetUp_Polynomial(PetscSpace sp)
559: {
560: PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
561: PetscInt ndegree = sp->order+1;
562: PetscInt deg;
563: PetscErrorCode ierr;
566: PetscMalloc1(ndegree, &poly->degrees);
567: for (deg = 0; deg < ndegree; ++deg) poly->degrees[deg] = deg;
568: return(0);
569: }
571: PetscErrorCode PetscSpaceDestroy_Polynomial(PetscSpace sp)
572: {
573: PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
574: PetscErrorCode ierr;
577: PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialGetTensor_C", NULL);
578: PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialSetTensor_C", NULL);
579: PetscFree(poly->degrees);
580: if (poly->subspaces) {
581: PetscInt d;
583: for (d = 0; d < poly->numVariables; ++d) {
584: PetscSpaceDestroy(&poly->subspaces[d]);
585: }
586: }
587: PetscFree(poly->subspaces);
588: PetscFree(poly);
589: return(0);
590: }
592: /* We treat the space as a tensor product of scalar polynomial spaces, so the dimension is multiplied by Nc */
593: PetscErrorCode PetscSpaceGetDimension_Polynomial(PetscSpace sp, PetscInt *dim)
594: {
595: PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
596: PetscInt deg = sp->order;
597: PetscInt n = poly->numVariables, i;
598: PetscReal D = 1.0;
601: if (poly->tensor) {
602: *dim = 1;
603: for (i = 0; i < n; ++i) *dim *= (deg+1);
604: } else {
605: for (i = 1; i <= n; ++i) {
606: D *= ((PetscReal) (deg+i))/i;
607: }
608: *dim = (PetscInt) (D + 0.5);
609: }
610: *dim *= sp->Nc;
611: return(0);
612: }
614: /*
615: LatticePoint_Internal - Returns all tuples of size 'len' with nonnegative integers that sum up to 'sum'.
617: Input Parameters:
618: + len - The length of the tuple
619: . sum - The sum of all entries in the tuple
620: - ind - The current multi-index of the tuple, initialized to the 0 tuple
622: Output Parameter:
623: + ind - The multi-index of the tuple, -1 indicates the iteration has terminated
624: . tup - A tuple of len integers addig to sum
626: Level: developer
628: .seealso:
629: */
630: static PetscErrorCode LatticePoint_Internal(PetscInt len, PetscInt sum, PetscInt ind[], PetscInt tup[])
631: {
632: PetscInt i;
636: if (len == 1) {
637: ind[0] = -1;
638: tup[0] = sum;
639: } else if (sum == 0) {
640: for (i = 0; i < len; ++i) {ind[0] = -1; tup[i] = 0;}
641: } else {
642: tup[0] = sum - ind[0];
643: LatticePoint_Internal(len-1, ind[0], &ind[1], &tup[1]);
644: if (ind[1] < 0) {
645: if (ind[0] == sum) {ind[0] = -1;}
646: else {ind[1] = 0; ++ind[0];}
647: }
648: }
649: return(0);
650: }
652: /*
653: LatticePointLexicographic_Internal - Returns all tuples of size 'len' with nonnegative integers that sum up to at most 'max'.
654: Ordering is lexicographic with lowest index as least significant in ordering.
655: e.g. for len == 2 and max == 2, this will return, in order, {0,0}, {1,0}, {2,0}, {0,1}, {1,1}, {2,0}.
657: Input Parameters:
658: + len - The length of the tuple
659: . max - The maximum sum
660: - tup - A tuple of length len+1: tup[len] > 0 indicates a stopping condition
662: Output Parameter:
663: . tup - A tuple of len integers whos sum is at most 'max'
664: */
665: static PetscErrorCode LatticePointLexicographic_Internal(PetscInt len, PetscInt max, PetscInt tup[])
666: {
668: while (len--) {
669: max -= tup[len];
670: if (!max) {
671: tup[len] = 0;
672: break;
673: }
674: }
675: tup[++len]++;
676: return(0);
677: }
679: /*
680: TensorPoint_Internal - Returns all tuples of size 'len' with nonnegative integers that are less than 'max'.
682: Input Parameters:
683: + len - The length of the tuple
684: . max - The max for all entries in the tuple
685: - ind - The current multi-index of the tuple, initialized to the 0 tuple
687: Output Parameter:
688: + ind - The multi-index of the tuple, -1 indicates the iteration has terminated
689: . tup - A tuple of len integers less than max
691: Level: developer
693: .seealso:
694: */
695: static PetscErrorCode TensorPoint_Internal(PetscInt len, PetscInt max, PetscInt ind[], PetscInt tup[])
696: {
697: PetscInt i;
701: if (len == 1) {
702: tup[0] = ind[0]++;
703: ind[0] = ind[0] >= max ? -1 : ind[0];
704: } else if (max == 0) {
705: for (i = 0; i < len; ++i) {ind[0] = -1; tup[i] = 0;}
706: } else {
707: tup[0] = ind[0];
708: TensorPoint_Internal(len-1, max, &ind[1], &tup[1]);
709: if (ind[1] < 0) {
710: ind[1] = 0;
711: if (ind[0] == max-1) {ind[0] = -1;}
712: else {++ind[0];}
713: }
714: }
715: return(0);
716: }
718: /*
719: TensorPointLexicographic_Internal - Returns all tuples of size 'len' with nonnegative integers that are all less than or equal to 'max'.
720: Ordering is lexicographic with lowest index as least significant in ordering.
721: e.g. for len == 2 and max == 2, this will return, in order, {0,0}, {1,0}, {2,0}, {0,1}, {1,1}, {2,1}, {0,2}, {1,2}, {2,2}.
723: Input Parameters:
724: + len - The length of the tuple
725: . max - The maximum value
726: - tup - A tuple of length len+1: tup[len] > 0 indicates a stopping condition
728: Output Parameter:
729: . tup - A tuple of len integers whos sum is at most 'max'
730: */
731: static PetscErrorCode TensorPointLexicographic_Internal(PetscInt len, PetscInt max, PetscInt tup[])
732: {
733: PetscInt i;
736: for (i = 0; i < len; i++) {
737: if (tup[i] < max) {
738: break;
739: } else {
740: tup[i] = 0;
741: }
742: }
743: tup[i]++;
744: return(0);
745: }
747: /*
748: p in [0, npoints), i in [0, pdim), c in [0, Nc)
750: B[p][i][c] = B[p][i_scalar][c][c]
751: */
752: PetscErrorCode PetscSpaceEvaluate_Polynomial(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[])
753: {
754: PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
755: DM dm = sp->dm;
756: PetscInt Nc = sp->Nc;
757: PetscInt ndegree = sp->order+1;
758: PetscInt *degrees = poly->degrees;
759: PetscInt dim = poly->numVariables;
760: PetscReal *lpoints, *tmp, *LB, *LD, *LH;
761: PetscInt *ind, *tup;
762: PetscInt c, pdim, d, der, i, p, deg, o;
763: PetscErrorCode ierr;
766: PetscSpaceGetDimension(sp, &pdim);
767: pdim /= Nc;
768: DMGetWorkArray(dm, npoints, MPIU_REAL, &lpoints);
769: DMGetWorkArray(dm, npoints*ndegree*3, MPIU_REAL, &tmp);
770: if (B) {DMGetWorkArray(dm, npoints*dim*ndegree, MPIU_REAL, &LB);}
771: if (D) {DMGetWorkArray(dm, npoints*dim*ndegree, MPIU_REAL, &LD);}
772: if (H) {DMGetWorkArray(dm, npoints*dim*ndegree, MPIU_REAL, &LH);}
773: for (d = 0; d < dim; ++d) {
774: for (p = 0; p < npoints; ++p) {
775: lpoints[p] = points[p*dim+d];
776: }
777: PetscDTLegendreEval(npoints, lpoints, ndegree, degrees, tmp, &tmp[1*npoints*ndegree], &tmp[2*npoints*ndegree]);
778: /* LB, LD, LH (ndegree * dim x npoints) */
779: for (deg = 0; deg < ndegree; ++deg) {
780: for (p = 0; p < npoints; ++p) {
781: if (B) LB[(deg*dim + d)*npoints + p] = tmp[(0*npoints + p)*ndegree+deg];
782: if (D) LD[(deg*dim + d)*npoints + p] = tmp[(1*npoints + p)*ndegree+deg];
783: if (H) LH[(deg*dim + d)*npoints + p] = tmp[(2*npoints + p)*ndegree+deg];
784: }
785: }
786: }
787: /* Multiply by A (pdim x ndegree * dim) */
788: PetscMalloc2(dim,&ind,dim,&tup);
789: if (B) {
790: /* B (npoints x pdim x Nc) */
791: PetscMemzero(B, npoints*pdim*Nc*Nc * sizeof(PetscReal));
792: if (poly->tensor) {
793: i = 0;
794: PetscMemzero(ind, dim * sizeof(PetscInt));
795: while (ind[0] >= 0) {
796: TensorPoint_Internal(dim, sp->order+1, ind, tup);
797: for (p = 0; p < npoints; ++p) {
798: B[(p*pdim + i)*Nc*Nc] = 1.0;
799: for (d = 0; d < dim; ++d) {
800: B[(p*pdim + i)*Nc*Nc] *= LB[(tup[d]*dim + d)*npoints + p];
801: }
802: }
803: ++i;
804: }
805: } else {
806: i = 0;
807: for (o = 0; o <= sp->order; ++o) {
808: PetscMemzero(ind, dim * sizeof(PetscInt));
809: while (ind[0] >= 0) {
810: LatticePoint_Internal(dim, o, ind, tup);
811: for (p = 0; p < npoints; ++p) {
812: B[(p*pdim + i)*Nc*Nc] = 1.0;
813: for (d = 0; d < dim; ++d) {
814: B[(p*pdim + i)*Nc*Nc] *= LB[(tup[d]*dim + d)*npoints + p];
815: }
816: }
817: ++i;
818: }
819: }
820: }
821: /* Make direct sum basis for multicomponent space */
822: for (p = 0; p < npoints; ++p) {
823: for (i = 0; i < pdim; ++i) {
824: for (c = 1; c < Nc; ++c) {
825: B[(p*pdim*Nc + i*Nc + c)*Nc + c] = B[(p*pdim + i)*Nc*Nc];
826: }
827: }
828: }
829: }
830: if (D) {
831: /* D (npoints x pdim x Nc x dim) */
832: PetscMemzero(D, npoints*pdim*Nc*Nc*dim * sizeof(PetscReal));
833: if (poly->tensor) {
834: i = 0;
835: PetscMemzero(ind, dim * sizeof(PetscInt));
836: while (ind[0] >= 0) {
837: TensorPoint_Internal(dim, sp->order+1, ind, tup);
838: for (p = 0; p < npoints; ++p) {
839: for (der = 0; der < dim; ++der) {
840: D[(p*pdim + i)*Nc*Nc*dim + der] = 1.0;
841: for (d = 0; d < dim; ++d) {
842: if (d == der) {
843: D[(p*pdim + i)*Nc*Nc*dim + der] *= LD[(tup[d]*dim + d)*npoints + p];
844: } else {
845: D[(p*pdim + i)*Nc*Nc*dim + der] *= LB[(tup[d]*dim + d)*npoints + p];
846: }
847: }
848: }
849: }
850: ++i;
851: }
852: } else {
853: i = 0;
854: for (o = 0; o <= sp->order; ++o) {
855: PetscMemzero(ind, dim * sizeof(PetscInt));
856: while (ind[0] >= 0) {
857: LatticePoint_Internal(dim, o, ind, tup);
858: for (p = 0; p < npoints; ++p) {
859: for (der = 0; der < dim; ++der) {
860: D[(p*pdim + i)*Nc*Nc*dim + der] = 1.0;
861: for (d = 0; d < dim; ++d) {
862: if (d == der) {
863: D[(p*pdim + i)*Nc*Nc*dim + der] *= LD[(tup[d]*dim + d)*npoints + p];
864: } else {
865: D[(p*pdim + i)*Nc*Nc*dim + der] *= LB[(tup[d]*dim + d)*npoints + p];
866: }
867: }
868: }
869: }
870: ++i;
871: }
872: }
873: }
874: /* Make direct sum basis for multicomponent space */
875: for (p = 0; p < npoints; ++p) {
876: for (i = 0; i < pdim; ++i) {
877: for (c = 1; c < Nc; ++c) {
878: for (d = 0; d < dim; ++d) {
879: D[((p*pdim*Nc + i*Nc + c)*Nc + c)*dim + d] = D[(p*pdim + i)*Nc*Nc*dim + d];
880: }
881: }
882: }
883: }
884: }
885: if (H) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Too lazy to code second derivatives");
886: PetscFree2(ind,tup);
887: if (B) {DMRestoreWorkArray(dm, npoints*dim*ndegree, MPIU_REAL, &LB);}
888: if (D) {DMRestoreWorkArray(dm, npoints*dim*ndegree, MPIU_REAL, &LD);}
889: if (H) {DMRestoreWorkArray(dm, npoints*dim*ndegree, MPIU_REAL, &LH);}
890: DMRestoreWorkArray(dm, npoints*ndegree*3, MPIU_REAL, &tmp);
891: DMRestoreWorkArray(dm, npoints, MPIU_REAL, &lpoints);
892: return(0);
893: }
895: /*@
896: PetscSpacePolynomialSetTensor - Set whether a function space is a space of tensor polynomials (the space is spanned
897: by polynomials whose degree in each variabl is bounded by the given order), as opposed to polynomials (the space is
898: spanned by polynomials whose total degree---summing over all variables---is bounded by the given order).
900: Input Parameters:
901: + sp - the function space object
902: - tensor - PETSC_TRUE for a tensor polynomial space, PETSC_FALSE for a polynomial space
904: Level: beginner
906: .seealso: PetscSpacePolynomialGetTensor(), PetscSpaceSetOrder(), PetscSpacePolynomialSetNumVariables()
907: @*/
908: PetscErrorCode PetscSpacePolynomialSetTensor(PetscSpace sp, PetscBool tensor)
909: {
914: PetscTryMethod(sp,"PetscSpacePolynomialSetTensor_C",(PetscSpace,PetscBool),(sp,tensor));
915: return(0);
916: }
918: /*@
919: PetscSpacePolynomialGetTensor - Get whether a function space is a space of tensor polynomials (the space is spanned
920: by polynomials whose degree in each variabl is bounded by the given order), as opposed to polynomials (the space is
921: spanned by polynomials whose total degree---summing over all variables---is bounded by the given order).
923: Input Parameters:
924: . sp - the function space object
926: Output Parameters:
927: . tensor - PETSC_TRUE for a tensor polynomial space, PETSC_FALSE for a polynomial space
929: Level: beginner
931: .seealso: PetscSpacePolynomialSetTensor(), PetscSpaceSetOrder(), PetscSpacePolynomialSetNumVariables()
932: @*/
933: PetscErrorCode PetscSpacePolynomialGetTensor(PetscSpace sp, PetscBool *tensor)
934: {
940: PetscTryMethod(sp,"PetscSpacePolynomialGetTensor_C",(PetscSpace,PetscBool*),(sp,tensor));
941: return(0);
942: }
944: static PetscErrorCode PetscSpacePolynomialSetTensor_Polynomial(PetscSpace sp, PetscBool tensor)
945: {
946: PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
949: poly->tensor = tensor;
950: return(0);
951: }
953: static PetscErrorCode PetscSpacePolynomialGetTensor_Polynomial(PetscSpace sp, PetscBool *tensor)
954: {
955: PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
960: *tensor = poly->tensor;
961: return(0);
962: }
964: static PetscErrorCode PetscSpaceGetHeightSubspace_Polynomial(PetscSpace sp, PetscInt height, PetscSpace *subsp)
965: {
966: PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
967: PetscInt Nc, dim, order;
968: PetscBool tensor;
969: PetscErrorCode ierr;
972: PetscSpaceGetNumComponents(sp, &Nc);
973: PetscSpacePolynomialGetNumVariables(sp, &dim);
974: PetscSpaceGetOrder(sp, &order);
975: PetscSpacePolynomialGetTensor(sp, &tensor);
976: if (height > dim || height < 0) {SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Asked for space at height %D for dimension %D space", height, dim);}
977: if (!poly->subspaces) {PetscCalloc1(dim, &poly->subspaces);}
978: if (height <= dim) {
979: if (!poly->subspaces[height-1]) {
980: PetscSpace sub;
982: PetscSpaceCreate(PetscObjectComm((PetscObject) sp), &sub);
983: PetscSpaceSetNumComponents(sub, Nc);
984: PetscSpaceSetOrder(sub, order);
985: PetscSpaceSetType(sub, PETSCSPACEPOLYNOMIAL);
986: PetscSpacePolynomialSetNumVariables(sub, dim-height);
987: PetscSpacePolynomialSetTensor(sub, tensor);
988: PetscSpaceSetUp(sub);
989: poly->subspaces[height-1] = sub;
990: }
991: *subsp = poly->subspaces[height-1];
992: } else {
993: *subsp = NULL;
994: }
995: return(0);
996: }
998: PetscErrorCode PetscSpaceInitialize_Polynomial(PetscSpace sp)
999: {
1003: sp->ops->setfromoptions = PetscSpaceSetFromOptions_Polynomial;
1004: sp->ops->setup = PetscSpaceSetUp_Polynomial;
1005: sp->ops->view = PetscSpaceView_Polynomial;
1006: sp->ops->destroy = PetscSpaceDestroy_Polynomial;
1007: sp->ops->getdimension = PetscSpaceGetDimension_Polynomial;
1008: sp->ops->evaluate = PetscSpaceEvaluate_Polynomial;
1009: sp->ops->getheightsubspace = PetscSpaceGetHeightSubspace_Polynomial;
1010: PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialGetTensor_C", PetscSpacePolynomialGetTensor_Polynomial);
1011: PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialSetTensor_C", PetscSpacePolynomialSetTensor_Polynomial);
1012: return(0);
1013: }
1015: /*MC
1016: PETSCSPACEPOLYNOMIAL = "poly" - A PetscSpace object that encapsulates a polynomial space, e.g. P1 is the space of
1017: linear polynomials. The space is replicated for each component.
1019: Level: intermediate
1021: .seealso: PetscSpaceType, PetscSpaceCreate(), PetscSpaceSetType()
1022: M*/
1024: PETSC_EXTERN PetscErrorCode PetscSpaceCreate_Polynomial(PetscSpace sp)
1025: {
1026: PetscSpace_Poly *poly;
1027: PetscErrorCode ierr;
1031: PetscNewLog(sp,&poly);
1032: sp->data = poly;
1034: poly->numVariables = 0;
1035: poly->symmetric = PETSC_FALSE;
1036: poly->tensor = PETSC_FALSE;
1037: poly->degrees = NULL;
1038: poly->subspaces = NULL;
1040: PetscSpaceInitialize_Polynomial(sp);
1041: return(0);
1042: }
1044: PetscErrorCode PetscSpacePolynomialSetSymmetric(PetscSpace sp, PetscBool sym)
1045: {
1046: PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
1050: poly->symmetric = sym;
1051: return(0);
1052: }
1054: PetscErrorCode PetscSpacePolynomialGetSymmetric(PetscSpace sp, PetscBool *sym)
1055: {
1056: PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
1061: *sym = poly->symmetric;
1062: return(0);
1063: }
1065: PetscErrorCode PetscSpacePolynomialSetNumVariables(PetscSpace sp, PetscInt n)
1066: {
1067: PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
1071: poly->numVariables = n;
1072: return(0);
1073: }
1075: PetscErrorCode PetscSpacePolynomialGetNumVariables(PetscSpace sp, PetscInt *n)
1076: {
1077: PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
1082: *n = poly->numVariables;
1083: return(0);
1084: }
1086: PetscErrorCode PetscSpaceSetFromOptions_Point(PetscOptionItems *PetscOptionsObject,PetscSpace sp)
1087: {
1088: PetscSpace_Point *pt = (PetscSpace_Point *) sp->data;
1089: PetscErrorCode ierr;
1092: PetscOptionsHead(PetscOptionsObject,"PetscSpace Point options");
1093: PetscOptionsInt("-petscspace_point_num_variables", "The number of different variables, e.g. x and y", "PetscSpacePointSetNumVariables", pt->numVariables, &pt->numVariables, NULL);
1094: PetscOptionsTail();
1095: return(0);
1096: }
1098: PetscErrorCode PetscSpacePointView_Ascii(PetscSpace sp, PetscViewer viewer)
1099: {
1100: PetscSpace_Point *pt = (PetscSpace_Point *) sp->data;
1101: PetscViewerFormat format;
1102: PetscErrorCode ierr;
1105: PetscViewerGetFormat(viewer, &format);
1106: if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
1107: PetscViewerASCIIPrintf(viewer, "Point space in dimension %d:\n", pt->numVariables);
1108: PetscViewerASCIIPushTab(viewer);
1109: PetscQuadratureView(pt->quad, viewer);
1110: PetscViewerASCIIPopTab(viewer);
1111: } else {
1112: PetscViewerASCIIPrintf(viewer, "Point space in dimension %d on %d points\n", pt->numVariables, pt->quad->numPoints);
1113: }
1114: return(0);
1115: }
1117: PetscErrorCode PetscSpaceView_Point(PetscSpace sp, PetscViewer viewer)
1118: {
1119: PetscBool iascii;
1125: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
1126: if (iascii) {PetscSpacePointView_Ascii(sp, viewer);}
1127: return(0);
1128: }
1130: PetscErrorCode PetscSpaceSetUp_Point(PetscSpace sp)
1131: {
1132: PetscSpace_Point *pt = (PetscSpace_Point *) sp->data;
1133: PetscErrorCode ierr;
1136: if (!pt->quad->points && sp->order >= 0) {
1137: PetscQuadratureDestroy(&pt->quad);
1138: PetscDTGaussJacobiQuadrature(pt->numVariables, sp->Nc, PetscMax(sp->order + 1, 1), -1.0, 1.0, &pt->quad);
1139: }
1140: return(0);
1141: }
1143: PetscErrorCode PetscSpaceDestroy_Point(PetscSpace sp)
1144: {
1145: PetscSpace_Point *pt = (PetscSpace_Point *) sp->data;
1146: PetscErrorCode ierr;
1149: PetscQuadratureDestroy(&pt->quad);
1150: PetscFree(pt);
1151: return(0);
1152: }
1154: PetscErrorCode PetscSpaceGetDimension_Point(PetscSpace sp, PetscInt *dim)
1155: {
1156: PetscSpace_Point *pt = (PetscSpace_Point *) sp->data;
1159: *dim = pt->quad->numPoints;
1160: return(0);
1161: }
1163: PetscErrorCode PetscSpaceEvaluate_Point(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[])
1164: {
1165: PetscSpace_Point *pt = (PetscSpace_Point *) sp->data;
1166: PetscInt dim = pt->numVariables, pdim = pt->quad->numPoints, d, p, i, c;
1167: PetscErrorCode ierr;
1170: if (npoints != pt->quad->numPoints) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot evaluate Point space on %d points != %d size", npoints, pt->quad->numPoints);
1171: PetscMemzero(B, npoints*pdim * sizeof(PetscReal));
1172: for (p = 0; p < npoints; ++p) {
1173: for (i = 0; i < pdim; ++i) {
1174: for (d = 0; d < dim; ++d) {
1175: if (PetscAbsReal(points[p*dim+d] - pt->quad->points[p*dim+d]) > 1.0e-10) break;
1176: }
1177: if (d >= dim) {B[p*pdim+i] = 1.0; break;}
1178: }
1179: }
1180: /* Replicate for other components */
1181: for (c = 1; c < sp->Nc; ++c) {
1182: for (p = 0; p < npoints; ++p) {
1183: for (i = 0; i < pdim; ++i) {
1184: B[(c*npoints + p)*pdim + i] = B[p*pdim + i];
1185: }
1186: }
1187: }
1188: if (D) {PetscMemzero(D, npoints*pdim*dim * sizeof(PetscReal));}
1189: if (H) {PetscMemzero(H, npoints*pdim*dim*dim * sizeof(PetscReal));}
1190: return(0);
1191: }
1193: PetscErrorCode PetscSpaceInitialize_Point(PetscSpace sp)
1194: {
1196: sp->ops->setfromoptions = PetscSpaceSetFromOptions_Point;
1197: sp->ops->setup = PetscSpaceSetUp_Point;
1198: sp->ops->view = PetscSpaceView_Point;
1199: sp->ops->destroy = PetscSpaceDestroy_Point;
1200: sp->ops->getdimension = PetscSpaceGetDimension_Point;
1201: sp->ops->evaluate = PetscSpaceEvaluate_Point;
1202: return(0);
1203: }
1205: /*MC
1206: PETSCSPACEPOINT = "point" - A PetscSpace object that encapsulates functions defined on a set of quadrature points.
1208: Level: intermediate
1210: .seealso: PetscSpaceType, PetscSpaceCreate(), PetscSpaceSetType()
1211: M*/
1213: PETSC_EXTERN PetscErrorCode PetscSpaceCreate_Point(PetscSpace sp)
1214: {
1215: PetscSpace_Point *pt;
1216: PetscErrorCode ierr;
1220: PetscNewLog(sp,&pt);
1221: sp->data = pt;
1223: pt->numVariables = 0;
1224: PetscQuadratureCreate(PETSC_COMM_SELF, &pt->quad);
1225: PetscQuadratureSetData(pt->quad, 0, 1, 0, NULL, NULL);
1227: PetscSpaceInitialize_Point(sp);
1228: return(0);
1229: }
1231: /*@
1232: PetscSpacePointSetPoints - Sets the evaluation points for the space to coincide with the points of a quadrature rule
1234: Logically collective
1236: Input Parameters:
1237: + sp - The PetscSpace
1238: - q - The PetscQuadrature defining the points
1240: Level: intermediate
1242: .keywords: PetscSpacePoint
1243: .seealso: PetscSpaceCreate(), PetscSpaceSetType()
1244: @*/
1245: PetscErrorCode PetscSpacePointSetPoints(PetscSpace sp, PetscQuadrature q)
1246: {
1247: PetscSpace_Point *pt = (PetscSpace_Point *) sp->data;
1248: PetscErrorCode ierr;
1253: PetscQuadratureDestroy(&pt->quad);
1254: PetscQuadratureDuplicate(q, &pt->quad);
1255: return(0);
1256: }
1258: /*@
1259: PetscSpacePointGetPoints - Gets the evaluation points for the space as the points of a quadrature rule
1261: Logically collective
1263: Input Parameter:
1264: . sp - The PetscSpace
1266: Output Parameter:
1267: . q - The PetscQuadrature defining the points
1269: Level: intermediate
1271: .keywords: PetscSpacePoint
1272: .seealso: PetscSpaceCreate(), PetscSpaceSetType()
1273: @*/
1274: PetscErrorCode PetscSpacePointGetPoints(PetscSpace sp, PetscQuadrature *q)
1275: {
1276: PetscSpace_Point *pt = (PetscSpace_Point *) sp->data;
1281: *q = pt->quad;
1282: return(0);
1283: }
1286: PetscClassId PETSCDUALSPACE_CLASSID = 0;
1288: PetscFunctionList PetscDualSpaceList = NULL;
1289: PetscBool PetscDualSpaceRegisterAllCalled = PETSC_FALSE;
1291: /*@C
1292: PetscDualSpaceRegister - Adds a new PetscDualSpace implementation
1294: Not Collective
1296: Input Parameters:
1297: + name - The name of a new user-defined creation routine
1298: - create_func - The creation routine itself
1300: Notes:
1301: PetscDualSpaceRegister() may be called multiple times to add several user-defined PetscDualSpaces
1303: Sample usage:
1304: .vb
1305: PetscDualSpaceRegister("my_space", MyPetscDualSpaceCreate);
1306: .ve
1308: Then, your PetscDualSpace type can be chosen with the procedural interface via
1309: .vb
1310: PetscDualSpaceCreate(MPI_Comm, PetscDualSpace *);
1311: PetscDualSpaceSetType(PetscDualSpace, "my_dual_space");
1312: .ve
1313: or at runtime via the option
1314: .vb
1315: -petscdualspace_type my_dual_space
1316: .ve
1318: Level: advanced
1320: .keywords: PetscDualSpace, register
1321: .seealso: PetscDualSpaceRegisterAll(), PetscDualSpaceRegisterDestroy()
1323: @*/
1324: PetscErrorCode PetscDualSpaceRegister(const char sname[], PetscErrorCode (*function)(PetscDualSpace))
1325: {
1329: PetscFunctionListAdd(&PetscDualSpaceList, sname, function);
1330: return(0);
1331: }
1333: /*@C
1334: PetscDualSpaceSetType - Builds a particular PetscDualSpace
1336: Collective on PetscDualSpace
1338: Input Parameters:
1339: + sp - The PetscDualSpace object
1340: - name - The kind of space
1342: Options Database Key:
1343: . -petscdualspace_type <type> - Sets the PetscDualSpace type; use -help for a list of available types
1345: Level: intermediate
1347: .keywords: PetscDualSpace, set, type
1348: .seealso: PetscDualSpaceGetType(), PetscDualSpaceCreate()
1349: @*/
1350: PetscErrorCode PetscDualSpaceSetType(PetscDualSpace sp, PetscDualSpaceType name)
1351: {
1352: PetscErrorCode (*r)(PetscDualSpace);
1353: PetscBool match;
1358: PetscObjectTypeCompare((PetscObject) sp, name, &match);
1359: if (match) return(0);
1361: if (!PetscDualSpaceRegisterAllCalled) {PetscDualSpaceRegisterAll();}
1362: PetscFunctionListFind(PetscDualSpaceList, name, &r);
1363: if (!r) SETERRQ1(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscDualSpace type: %s", name);
1365: if (sp->ops->destroy) {
1366: (*sp->ops->destroy)(sp);
1367: sp->ops->destroy = NULL;
1368: }
1369: (*r)(sp);
1370: PetscObjectChangeTypeName((PetscObject) sp, name);
1371: return(0);
1372: }
1374: /*@C
1375: PetscDualSpaceGetType - Gets the PetscDualSpace type name (as a string) from the object.
1377: Not Collective
1379: Input Parameter:
1380: . sp - The PetscDualSpace
1382: Output Parameter:
1383: . name - The PetscDualSpace type name
1385: Level: intermediate
1387: .keywords: PetscDualSpace, get, type, name
1388: .seealso: PetscDualSpaceSetType(), PetscDualSpaceCreate()
1389: @*/
1390: PetscErrorCode PetscDualSpaceGetType(PetscDualSpace sp, PetscDualSpaceType *name)
1391: {
1397: if (!PetscDualSpaceRegisterAllCalled) {
1398: PetscDualSpaceRegisterAll();
1399: }
1400: *name = ((PetscObject) sp)->type_name;
1401: return(0);
1402: }
1404: /*@
1405: PetscDualSpaceView - Views a PetscDualSpace
1407: Collective on PetscDualSpace
1409: Input Parameter:
1410: + sp - the PetscDualSpace object to view
1411: - v - the viewer
1413: Level: developer
1415: .seealso PetscDualSpaceDestroy()
1416: @*/
1417: PetscErrorCode PetscDualSpaceView(PetscDualSpace sp, PetscViewer v)
1418: {
1423: if (!v) {PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) sp), &v);}
1424: if (sp->ops->view) {(*sp->ops->view)(sp, v);}
1425: return(0);
1426: }
1428: /*@
1429: PetscDualSpaceSetFromOptions - sets parameters in a PetscDualSpace from the options database
1431: Collective on PetscDualSpace
1433: Input Parameter:
1434: . sp - the PetscDualSpace object to set options for
1436: Options Database:
1437: . -petscspace_order the approximation order of the space
1439: Level: developer
1441: .seealso PetscDualSpaceView()
1442: @*/
1443: PetscErrorCode PetscDualSpaceSetFromOptions(PetscDualSpace sp)
1444: {
1445: const char *defaultType;
1446: char name[256];
1447: PetscBool flg;
1452: if (!((PetscObject) sp)->type_name) {
1453: defaultType = PETSCDUALSPACELAGRANGE;
1454: } else {
1455: defaultType = ((PetscObject) sp)->type_name;
1456: }
1457: if (!PetscSpaceRegisterAllCalled) {PetscSpaceRegisterAll();}
1459: PetscObjectOptionsBegin((PetscObject) sp);
1460: PetscOptionsFList("-petscdualspace_type", "Dual space", "PetscDualSpaceSetType", PetscDualSpaceList, defaultType, name, 256, &flg);
1461: if (flg) {
1462: PetscDualSpaceSetType(sp, name);
1463: } else if (!((PetscObject) sp)->type_name) {
1464: PetscDualSpaceSetType(sp, defaultType);
1465: }
1466: PetscOptionsInt("-petscdualspace_order", "The approximation order", "PetscDualSpaceSetOrder", sp->order, &sp->order, NULL);
1467: PetscOptionsInt("-petscdualspace_components", "The number of components", "PetscDualSpaceSetNumComponents", sp->Nc, &sp->Nc, NULL);
1468: if (sp->ops->setfromoptions) {
1469: (*sp->ops->setfromoptions)(PetscOptionsObject,sp);
1470: }
1471: /* process any options handlers added with PetscObjectAddOptionsHandler() */
1472: PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) sp);
1473: PetscOptionsEnd();
1474: PetscDualSpaceViewFromOptions(sp, NULL, "-petscdualspace_view");
1475: return(0);
1476: }
1478: /*@
1479: PetscDualSpaceSetUp - Construct a basis for the PetscDualSpace
1481: Collective on PetscDualSpace
1483: Input Parameter:
1484: . sp - the PetscDualSpace object to setup
1486: Level: developer
1488: .seealso PetscDualSpaceView(), PetscDualSpaceDestroy()
1489: @*/
1490: PetscErrorCode PetscDualSpaceSetUp(PetscDualSpace sp)
1491: {
1496: if (sp->setupcalled) return(0);
1497: sp->setupcalled = PETSC_TRUE;
1498: if (sp->ops->setup) {(*sp->ops->setup)(sp);}
1499: return(0);
1500: }
1502: /*@
1503: PetscDualSpaceDestroy - Destroys a PetscDualSpace object
1505: Collective on PetscDualSpace
1507: Input Parameter:
1508: . sp - the PetscDualSpace object to destroy
1510: Level: developer
1512: .seealso PetscDualSpaceView()
1513: @*/
1514: PetscErrorCode PetscDualSpaceDestroy(PetscDualSpace *sp)
1515: {
1516: PetscInt dim, f;
1520: if (!*sp) return(0);
1523: if (--((PetscObject)(*sp))->refct > 0) {*sp = 0; return(0);}
1524: ((PetscObject) (*sp))->refct = 0;
1526: PetscDualSpaceGetDimension(*sp, &dim);
1527: for (f = 0; f < dim; ++f) {
1528: PetscQuadratureDestroy(&(*sp)->functional[f]);
1529: }
1530: PetscFree((*sp)->functional);
1531: DMDestroy(&(*sp)->dm);
1533: if ((*sp)->ops->destroy) {(*(*sp)->ops->destroy)(*sp);}
1534: PetscHeaderDestroy(sp);
1535: return(0);
1536: }
1538: /*@
1539: PetscDualSpaceCreate - Creates an empty PetscDualSpace object. The type can then be set with PetscDualSpaceSetType().
1541: Collective on MPI_Comm
1543: Input Parameter:
1544: . comm - The communicator for the PetscDualSpace object
1546: Output Parameter:
1547: . sp - The PetscDualSpace object
1549: Level: beginner
1551: .seealso: PetscDualSpaceSetType(), PETSCDUALSPACELAGRANGE
1552: @*/
1553: PetscErrorCode PetscDualSpaceCreate(MPI_Comm comm, PetscDualSpace *sp)
1554: {
1555: PetscDualSpace s;
1560: PetscCitationsRegister(FECitation,&FEcite);
1561: *sp = NULL;
1562: PetscFEInitializePackage();
1564: PetscHeaderCreate(s, PETSCDUALSPACE_CLASSID, "PetscDualSpace", "Dual Space", "PetscDualSpace", comm, PetscDualSpaceDestroy, PetscDualSpaceView);
1566: s->order = 0;
1567: s->Nc = 1;
1568: s->setupcalled = PETSC_FALSE;
1570: *sp = s;
1571: return(0);
1572: }
1574: /*@
1575: PetscDualSpaceDuplicate - Creates a duplicate PetscDualSpace object, however it is not setup.
1577: Collective on PetscDualSpace
1579: Input Parameter:
1580: . sp - The original PetscDualSpace
1582: Output Parameter:
1583: . spNew - The duplicate PetscDualSpace
1585: Level: beginner
1587: .seealso: PetscDualSpaceCreate(), PetscDualSpaceSetType()
1588: @*/
1589: PetscErrorCode PetscDualSpaceDuplicate(PetscDualSpace sp, PetscDualSpace *spNew)
1590: {
1596: (*sp->ops->duplicate)(sp, spNew);
1597: return(0);
1598: }
1600: /*@
1601: PetscDualSpaceGetDM - Get the DM representing the reference cell
1603: Not collective
1605: Input Parameter:
1606: . sp - The PetscDualSpace
1608: Output Parameter:
1609: . dm - The reference cell
1611: Level: intermediate
1613: .seealso: PetscDualSpaceSetDM(), PetscDualSpaceCreate()
1614: @*/
1615: PetscErrorCode PetscDualSpaceGetDM(PetscDualSpace sp, DM *dm)
1616: {
1620: *dm = sp->dm;
1621: return(0);
1622: }
1624: /*@
1625: PetscDualSpaceSetDM - Get the DM representing the reference cell
1627: Not collective
1629: Input Parameters:
1630: + sp - The PetscDualSpace
1631: - dm - The reference cell
1633: Level: intermediate
1635: .seealso: PetscDualSpaceGetDM(), PetscDualSpaceCreate()
1636: @*/
1637: PetscErrorCode PetscDualSpaceSetDM(PetscDualSpace sp, DM dm)
1638: {
1644: DMDestroy(&sp->dm);
1645: PetscObjectReference((PetscObject) dm);
1646: sp->dm = dm;
1647: return(0);
1648: }
1650: /*@
1651: PetscDualSpaceGetOrder - Get the order of the dual space
1653: Not collective
1655: Input Parameter:
1656: . sp - The PetscDualSpace
1658: Output Parameter:
1659: . order - The order
1661: Level: intermediate
1663: .seealso: PetscDualSpaceSetOrder(), PetscDualSpaceCreate()
1664: @*/
1665: PetscErrorCode PetscDualSpaceGetOrder(PetscDualSpace sp, PetscInt *order)
1666: {
1670: *order = sp->order;
1671: return(0);
1672: }
1674: /*@
1675: PetscDualSpaceSetOrder - Set the order of the dual space
1677: Not collective
1679: Input Parameters:
1680: + sp - The PetscDualSpace
1681: - order - The order
1683: Level: intermediate
1685: .seealso: PetscDualSpaceGetOrder(), PetscDualSpaceCreate()
1686: @*/
1687: PetscErrorCode PetscDualSpaceSetOrder(PetscDualSpace sp, PetscInt order)
1688: {
1691: sp->order = order;
1692: return(0);
1693: }
1695: /*@
1696: PetscDualSpaceGetNumComponents - Return the number of components for this space
1698: Input Parameter:
1699: . sp - The PetscDualSpace
1701: Output Parameter:
1702: . Nc - The number of components
1704: Note: A vector space, for example, will have d components, where d is the spatial dimension
1706: Level: intermediate
1708: .seealso: PetscDualSpaceSetNumComponents(), PetscDualSpaceGetDimension(), PetscDualSpaceCreate(), PetscDualSpace
1709: @*/
1710: PetscErrorCode PetscDualSpaceGetNumComponents(PetscDualSpace sp, PetscInt *Nc)
1711: {
1715: *Nc = sp->Nc;
1716: return(0);
1717: }
1719: /*@
1720: PetscDualSpaceSetNumComponents - Set the number of components for this space
1722: Input Parameters:
1723: + sp - The PetscDualSpace
1724: - order - The number of components
1726: Level: intermediate
1728: .seealso: PetscDualSpaceGetNumComponents(), PetscDualSpaceCreate(), PetscDualSpace
1729: @*/
1730: PetscErrorCode PetscDualSpaceSetNumComponents(PetscDualSpace sp, PetscInt Nc)
1731: {
1734: sp->Nc = Nc;
1735: return(0);
1736: }
1738: /*@
1739: PetscDualSpaceLagrangeGetTensor - Get the tensor nature of the dual space
1741: Not collective
1743: Input Parameter:
1744: . sp - The PetscDualSpace
1746: Output Parameter:
1747: . tensor - Whether the dual space has tensor layout (vs. simplicial)
1749: Level: intermediate
1751: .seealso: PetscDualSpaceLagrangeSetTensor(), PetscDualSpaceCreate()
1752: @*/
1753: PetscErrorCode PetscDualSpaceLagrangeGetTensor(PetscDualSpace sp, PetscBool *tensor)
1754: {
1760: PetscTryMethod(sp,"PetscDualSpaceLagrangeGetTensor_C",(PetscDualSpace,PetscBool *),(sp,tensor));
1761: return(0);
1762: }
1764: /*@
1765: PetscDualSpaceLagrangeSetTensor - Set the tensor nature of the dual space
1767: Not collective
1769: Input Parameters:
1770: + sp - The PetscDualSpace
1771: - tensor - Whether the dual space has tensor layout (vs. simplicial)
1773: Level: intermediate
1775: .seealso: PetscDualSpaceLagrangeGetTensor(), PetscDualSpaceCreate()
1776: @*/
1777: PetscErrorCode PetscDualSpaceLagrangeSetTensor(PetscDualSpace sp, PetscBool tensor)
1778: {
1783: PetscTryMethod(sp,"PetscDualSpaceLagrangeSetTensor_C",(PetscDualSpace,PetscBool),(sp,tensor));
1784: return(0);
1785: }
1787: /*@
1788: PetscDualSpaceGetFunctional - Get the i-th basis functional in the dual space
1790: Not collective
1792: Input Parameters:
1793: + sp - The PetscDualSpace
1794: - i - The basis number
1796: Output Parameter:
1797: . functional - The basis functional
1799: Level: intermediate
1801: .seealso: PetscDualSpaceGetDimension(), PetscDualSpaceCreate()
1802: @*/
1803: PetscErrorCode PetscDualSpaceGetFunctional(PetscDualSpace sp, PetscInt i, PetscQuadrature *functional)
1804: {
1805: PetscInt dim;
1811: PetscDualSpaceGetDimension(sp, &dim);
1812: if ((i < 0) || (i >= dim)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Functional index %d must be in [0, %d)", i, dim);
1813: *functional = sp->functional[i];
1814: return(0);
1815: }
1817: /*@
1818: PetscDualSpaceGetDimension - Get the dimension of the dual space, i.e. the number of basis functionals
1820: Not collective
1822: Input Parameter:
1823: . sp - The PetscDualSpace
1825: Output Parameter:
1826: . dim - The dimension
1828: Level: intermediate
1830: .seealso: PetscDualSpaceGetFunctional(), PetscDualSpaceCreate()
1831: @*/
1832: PetscErrorCode PetscDualSpaceGetDimension(PetscDualSpace sp, PetscInt *dim)
1833: {
1839: *dim = 0;
1840: if (sp->ops->getdimension) {(*sp->ops->getdimension)(sp, dim);}
1841: return(0);
1842: }
1844: /*@C
1845: PetscDualSpaceGetNumDof - Get the number of degrees of freedom for each spatial (topological) dimension
1847: Not collective
1849: Input Parameter:
1850: . sp - The PetscDualSpace
1852: Output Parameter:
1853: . numDof - An array of length dim+1 which holds the number of dofs for each dimension
1855: Level: intermediate
1857: .seealso: PetscDualSpaceGetFunctional(), PetscDualSpaceCreate()
1858: @*/
1859: PetscErrorCode PetscDualSpaceGetNumDof(PetscDualSpace sp, const PetscInt **numDof)
1860: {
1866: (*sp->ops->getnumdof)(sp, numDof);
1867: if (!*numDof) SETERRQ(PetscObjectComm((PetscObject) sp), PETSC_ERR_LIB, "Empty numDof[] returned from dual space implementation");
1868: return(0);
1869: }
1871: /*@
1872: PetscDualSpaceCreateReferenceCell - Create a DMPLEX with the appropriate FEM reference cell
1874: Collective on PetscDualSpace
1876: Input Parameters:
1877: + sp - The PetscDualSpace
1878: . dim - The spatial dimension
1879: - simplex - Flag for simplex, otherwise use a tensor-product cell
1881: Output Parameter:
1882: . refdm - The reference cell
1884: Level: advanced
1886: .keywords: PetscDualSpace, reference cell
1887: .seealso: PetscDualSpaceCreate(), DMPLEX
1888: @*/
1889: PetscErrorCode PetscDualSpaceCreateReferenceCell(PetscDualSpace sp, PetscInt dim, PetscBool simplex, DM *refdm)
1890: {
1894: DMPlexCreateReferenceCell(PetscObjectComm((PetscObject) sp), dim, simplex, refdm);
1895: return(0);
1896: }
1898: /*@C
1899: PetscDualSpaceApply - Apply a functional from the dual space basis to an input function
1901: Input Parameters:
1902: + sp - The PetscDualSpace object
1903: . f - The basis functional index
1904: . time - The time
1905: . cgeom - A context with geometric information for this cell, we use v0 (the initial vertex) and J (the Jacobian)
1906: . numComp - The number of components for the function
1907: . func - The input function
1908: - ctx - A context for the function
1910: Output Parameter:
1911: . value - numComp output values
1913: Note: The calling sequence for the callback func is given by:
1915: $ func(PetscInt dim, PetscReal time, const PetscReal x[],
1916: $ PetscInt numComponents, PetscScalar values[], void *ctx)
1918: Level: developer
1920: .seealso: PetscDualSpaceCreate()
1921: @*/
1922: PetscErrorCode PetscDualSpaceApply(PetscDualSpace sp, PetscInt f, PetscReal time, PetscFECellGeom *cgeom, PetscInt numComp, PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void *ctx, PetscScalar *value)
1923: {
1930: (*sp->ops->apply)(sp, f, time, cgeom, numComp, func, ctx, value);
1931: return(0);
1932: }
1934: /*@C
1935: PetscDualSpaceApplyDefault - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional.
1937: Input Parameters:
1938: + sp - The PetscDualSpace object
1939: . f - The basis functional index
1940: . time - The time
1941: . cgeom - A context with geometric information for this cell, we use v0 (the initial vertex) and J (the Jacobian)
1942: . Nc - The number of components for the function
1943: . func - The input function
1944: - ctx - A context for the function
1946: Output Parameter:
1947: . value - The output value
1949: Note: The calling sequence for the callback func is given by:
1951: $ func(PetscInt dim, PetscReal time, const PetscReal x[],
1952: $ PetscInt numComponents, PetscScalar values[], void *ctx)
1954: and the idea is to evaluate the functional as an integral
1956: $ n(f) = int dx n(x) . f(x)
1958: where both n and f have Nc components.
1960: Level: developer
1962: .seealso: PetscDualSpaceCreate()
1963: @*/
1964: PetscErrorCode PetscDualSpaceApplyDefault(PetscDualSpace sp, PetscInt f, PetscReal time, PetscFECellGeom *cgeom, PetscInt Nc, PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void *ctx, PetscScalar *value)
1965: {
1966: DM dm;
1967: PetscQuadrature n;
1968: const PetscReal *points, *weights;
1969: PetscReal x[3];
1970: PetscScalar *val;
1971: PetscInt dim, qNc, c, Nq, q;
1972: PetscErrorCode ierr;
1977: PetscDualSpaceGetDM(sp, &dm);
1978: PetscDualSpaceGetFunctional(sp, f, &n);
1979: PetscQuadratureGetData(n, &dim, &qNc, &Nq, &points, &weights);
1980: if (dim != cgeom->dim) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature spatial dimension %D != cell geometry dimension %D", dim, cgeom->dim);
1981: if (qNc != Nc) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc);
1982: DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val);
1983: *value = 0.0;
1984: for (q = 0; q < Nq; ++q) {
1985: CoordinatesRefToReal(cgeom->dimEmbed, dim, cgeom->v0, cgeom->J, &points[q*dim], x);
1986: (*func)(cgeom->dimEmbed, time, x, Nc, val, ctx);
1987: for (c = 0; c < Nc; ++c) {
1988: *value += val[c]*weights[q*Nc+c];
1989: }
1990: }
1991: DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val);
1992: return(0);
1993: }
1995: /*@C
1996: PetscDualSpaceApplyFVM - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional at the cell centroid.
1998: Input Parameters:
1999: + sp - The PetscDualSpace object
2000: . f - The basis functional index
2001: . time - The time
2002: . cgeom - A context with geometric information for this cell, we currently just use the centroid
2003: . Nc - The number of components for the function
2004: . func - The input function
2005: - ctx - A context for the function
2007: Output Parameter:
2008: . value - The output value (scalar)
2010: Note: The calling sequence for the callback func is given by:
2012: $ func(PetscInt dim, PetscReal time, const PetscReal x[],
2013: $ PetscInt numComponents, PetscScalar values[], void *ctx)
2015: and the idea is to evaluate the functional as an integral
2017: $ n(f) = int dx n(x) . f(x)
2019: where both n and f have Nc components.
2021: Level: developer
2023: .seealso: PetscDualSpaceCreate()
2024: @*/
2025: PetscErrorCode PetscDualSpaceApplyFVM(PetscDualSpace sp, PetscInt f, PetscReal time, PetscFVCellGeom *cgeom, PetscInt Nc, PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void *ctx, PetscScalar *value)
2026: {
2027: DM dm;
2028: PetscQuadrature n;
2029: const PetscReal *points, *weights;
2030: PetscScalar *val;
2031: PetscInt dimEmbed, qNc, c, Nq, q;
2032: PetscErrorCode ierr;
2037: PetscDualSpaceGetDM(sp, &dm);
2038: DMGetCoordinateDim(dm, &dimEmbed);
2039: PetscDualSpaceGetFunctional(sp, f, &n);
2040: PetscQuadratureGetData(n, NULL, &qNc, &Nq, &points, &weights);
2041: if (qNc != Nc) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc);
2042: DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val);
2043: *value = 0.;
2044: for (q = 0; q < Nq; ++q) {
2045: (*func)(dimEmbed, time, cgeom->centroid, Nc, val, ctx);
2046: for (c = 0; c < Nc; ++c) {
2047: *value += val[c]*weights[q*Nc+c];
2048: }
2049: }
2050: DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val);
2051: return(0);
2052: }
2054: /*@
2055: PetscDualSpaceGetHeightSubspace - Get the subset of the dual space basis that is supported on a mesh point of a given height.
2057: If the dual space is not defined on mesh points of the given height (e.g. if the space is discontinuous and
2058: pointwise values are not defined on the element boundaries), or if the implementation of PetscDualSpace does not
2059: support extracting subspaces, then NULL is returned.
2061: This does not increment the reference count on the returned dual space, and the user should not destroy it.
2063: Not collective
2065: Input Parameters:
2066: + sp - the PetscDualSpace object
2067: - height - the height of the mesh point for which the subspace is desired
2069: Output Parameter:
2070: . subsp - the subspace
2072: Level: advanced
2074: .seealso: PetscSpaceGetHeightSubspace(), PetscDualSpace
2075: @*/
2076: PetscErrorCode PetscDualSpaceGetHeightSubspace(PetscDualSpace sp, PetscInt height, PetscDualSpace *subsp)
2077: {
2083: *subsp = NULL;
2084: if (sp->ops->getheightsubspace) {
2085: (*sp->ops->getheightsubspace)(sp, height, subsp);
2086: }
2087: return(0);
2088: }
2090: static PetscErrorCode PetscDualSpaceLagrangeGetTensor_Lagrange(PetscDualSpace sp, PetscBool *tensor)
2091: {
2092: PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data;
2095: *tensor = lag->tensorSpace;
2096: return(0);
2097: }
2099: static PetscErrorCode PetscDualSpaceLagrangeSetTensor_Lagrange(PetscDualSpace sp, PetscBool tensor)
2100: {
2101: PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *)sp->data;
2104: lag->tensorSpace = tensor;
2105: return(0);
2106: }
2108: #define BaryIndex(perEdge,a,b,c) (((b)*(2*perEdge+1-(b)))/2)+(c)
2110: #define CartIndex(perEdge,a,b) (perEdge*(a)+b)
2112: static PetscErrorCode PetscDualSpaceGetSymmetries_Lagrange(PetscDualSpace sp, const PetscInt ****perms, const PetscScalar ****flips)
2113: {
2115: PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data;
2116: PetscInt dim, order, p, Nc;
2117: PetscErrorCode ierr;
2120: PetscDualSpaceGetOrder(sp,&order);
2121: PetscDualSpaceGetNumComponents(sp,&Nc);
2122: DMGetDimension(sp->dm,&dim);
2123: if (!dim || !lag->continuous || order < 3) return(0);
2124: if (dim > 3) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Lagrange symmetries not implemented for dim = %D > 3",dim);
2125: if (!lag->symmetries) { /* store symmetries */
2126: PetscDualSpace hsp;
2127: DM K;
2128: PetscInt numPoints = 1, d;
2129: PetscInt numFaces;
2130: PetscInt ***symmetries;
2131: const PetscInt ***hsymmetries;
2133: if (lag->simplexCell) {
2134: numFaces = 1 + dim;
2135: for (d = 0; d < dim; d++) numPoints = numPoints * 2 + 1;
2136: }
2137: else {
2138: numPoints = PetscPowInt(3,dim);
2139: numFaces = 2 * dim;
2140: }
2141: PetscCalloc1(numPoints,&symmetries);
2142: if (0 < dim && dim < 3) { /* compute self symmetries */
2143: PetscInt **cellSymmetries;
2145: lag->numSelfSym = 2 * numFaces;
2146: lag->selfSymOff = numFaces;
2147: PetscCalloc1(2*numFaces,&cellSymmetries);
2148: /* we want to be able to index symmetries directly with the orientations, which range from [-numFaces,numFaces) */
2149: symmetries[0] = &cellSymmetries[numFaces];
2150: if (dim == 1) {
2151: PetscInt dofPerEdge = order - 1;
2153: if (dofPerEdge > 1) {
2154: PetscInt i, j, *reverse;
2156: PetscMalloc1(dofPerEdge*Nc,&reverse);
2157: for (i = 0; i < dofPerEdge; i++) {
2158: for (j = 0; j < Nc; j++) {
2159: reverse[i*Nc + j] = Nc * (dofPerEdge - 1 - i) + j;
2160: }
2161: }
2162: symmetries[0][-2] = reverse;
2164: /* yes, this is redundant, but it makes it easier to cleanup if I don't have to worry about what not to free */
2165: PetscMalloc1(dofPerEdge*Nc,&reverse);
2166: for (i = 0; i < dofPerEdge; i++) {
2167: for (j = 0; j < Nc; j++) {
2168: reverse[i*Nc + j] = Nc * (dofPerEdge - 1 - i) + j;
2169: }
2170: }
2171: symmetries[0][1] = reverse;
2172: }
2173: } else {
2174: PetscInt dofPerEdge = lag->simplexCell ? (order - 2) : (order - 1), s;
2175: PetscInt dofPerFace;
2177: if (dofPerEdge > 1) {
2178: for (s = -numFaces; s < numFaces; s++) {
2179: PetscInt *sym, i, j, k, l;
2181: if (!s) continue;
2182: if (lag->simplexCell) {
2183: dofPerFace = (dofPerEdge * (dofPerEdge + 1))/2;
2184: PetscMalloc1(Nc*dofPerFace,&sym);
2185: for (j = 0, l = 0; j < dofPerEdge; j++) {
2186: for (k = 0; k < dofPerEdge - j; k++, l++) {
2187: i = dofPerEdge - 1 - j - k;
2188: switch (s) {
2189: case -3:
2190: sym[Nc*l] = BaryIndex(dofPerEdge,i,k,j);
2191: break;
2192: case -2:
2193: sym[Nc*l] = BaryIndex(dofPerEdge,j,i,k);
2194: break;
2195: case -1:
2196: sym[Nc*l] = BaryIndex(dofPerEdge,k,j,i);
2197: break;
2198: case 1:
2199: sym[Nc*l] = BaryIndex(dofPerEdge,k,i,j);
2200: break;
2201: case 2:
2202: sym[Nc*l] = BaryIndex(dofPerEdge,j,k,i);
2203: break;
2204: }
2205: }
2206: }
2207: } else {
2208: dofPerFace = dofPerEdge * dofPerEdge;
2209: PetscMalloc1(Nc*dofPerFace,&sym);
2210: for (j = 0, l = 0; j < dofPerEdge; j++) {
2211: for (k = 0; k < dofPerEdge; k++, l++) {
2212: switch (s) {
2213: case -4:
2214: sym[Nc*l] = CartIndex(dofPerEdge,k,j);
2215: break;
2216: case -3:
2217: sym[Nc*l] = CartIndex(dofPerEdge,(dofPerEdge - 1 - j),k);
2218: break;
2219: case -2:
2220: sym[Nc*l] = CartIndex(dofPerEdge,(dofPerEdge - 1 - k),(dofPerEdge - 1 - j));
2221: break;
2222: case -1:
2223: sym[Nc*l] = CartIndex(dofPerEdge,j,(dofPerEdge - 1 - k));
2224: break;
2225: case 1:
2226: sym[Nc*l] = CartIndex(dofPerEdge,(dofPerEdge - 1 - k),j);
2227: break;
2228: case 2:
2229: sym[Nc*l] = CartIndex(dofPerEdge,(dofPerEdge - 1 - j),(dofPerEdge - 1 - k));
2230: break;
2231: case 3:
2232: sym[Nc*l] = CartIndex(dofPerEdge,k,(dofPerEdge - 1 - j));
2233: break;
2234: }
2235: }
2236: }
2237: }
2238: for (i = 0; i < dofPerFace; i++) {
2239: sym[Nc*i] *= Nc;
2240: for (j = 1; j < Nc; j++) {
2241: sym[Nc*i+j] = sym[Nc*i] + j;
2242: }
2243: }
2244: symmetries[0][s] = sym;
2245: }
2246: }
2247: }
2248: }
2249: PetscDualSpaceGetHeightSubspace(sp,1,&hsp);
2250: PetscDualSpaceGetSymmetries(hsp,&hsymmetries,NULL);
2251: if (hsymmetries) {
2252: PetscBool *seen;
2253: const PetscInt *cone;
2254: PetscInt KclosureSize, *Kclosure = NULL;
2256: PetscDualSpaceGetDM(sp,&K);
2257: PetscCalloc1(numPoints,&seen);
2258: DMPlexGetCone(K,0,&cone);
2259: DMPlexGetTransitiveClosure(K,0,PETSC_TRUE,&KclosureSize,&Kclosure);
2260: for (p = 0; p < numFaces; p++) {
2261: PetscInt closureSize, *closure = NULL, q;
2263: DMPlexGetTransitiveClosure(K,cone[p],PETSC_TRUE,&closureSize,&closure);
2264: for (q = 0; q < closureSize; q++) {
2265: PetscInt point = closure[2*q], r;
2267: if(!seen[point]) {
2268: for (r = 0; r < KclosureSize; r++) {
2269: if (Kclosure[2 * r] == point) break;
2270: }
2271: seen[point] = PETSC_TRUE;
2272: symmetries[r] = (PetscInt **) hsymmetries[q];
2273: }
2274: }
2275: DMPlexRestoreTransitiveClosure(K,cone[p],PETSC_TRUE,&closureSize,&closure);
2276: }
2277: DMPlexRestoreTransitiveClosure(K,0,PETSC_TRUE,&KclosureSize,&Kclosure);
2278: PetscFree(seen);
2279: }
2280: lag->symmetries = symmetries;
2281: }
2282: if (perms) *perms = (const PetscInt ***) lag->symmetries;
2283: return(0);
2284: }
2286: /*@C
2287: PetscDualSpaceGetSymmetries - Returns a description of the symmetries of this basis
2289: Not collective
2291: Input Parameter:
2292: . sp - the PetscDualSpace object
2294: Output Parameters:
2295: + perms - Permutations of the local degrees of freedom, parameterized by the point orientation
2296: - flips - Sign reversal of the local degrees of freedom, parameterized by the point orientation
2298: Note: The permutation and flip arrays are organized in the following way
2299: $ perms[p][ornt][dof # on point] = new local dof #
2300: $ flips[p][ornt][dof # on point] = reversal or not
2302: Level: developer
2304: .seealso: PetscDualSpaceSetSymmetries()
2305: @*/
2306: PetscErrorCode PetscDualSpaceGetSymmetries(PetscDualSpace sp, const PetscInt ****perms, const PetscScalar ****flips)
2307: {
2312: if (perms) {
2314: *perms = NULL;
2315: }
2316: if (flips) {
2318: *flips = NULL;
2319: }
2320: if (sp->ops->getsymmetries) {
2321: (sp->ops->getsymmetries)(sp,perms,flips);
2322: }
2323: return(0);
2324: }
2326: static PetscErrorCode PetscDualSpaceLagrangeView_Ascii(PetscDualSpace sp, PetscViewer viewer)
2327: {
2328: PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data;
2329: PetscErrorCode ierr;
2332: PetscViewerASCIIPrintf(viewer, "%s %sLagrange dual space of order %D", lag->continuous ? "Continuous" : "Discontinuous", lag->tensorSpace ? "Tensor " : "", sp->order, sp->Nc);
2333: if (sp->Nc > 1) {PetscViewerASCIIPrintf(viewer, " with %D components", sp->Nc);}
2334: PetscViewerASCIIPrintf(viewer, "\n");
2335: return(0);
2336: }
2338: PetscErrorCode PetscDualSpaceView_Lagrange(PetscDualSpace sp, PetscViewer viewer)
2339: {
2340: PetscBool iascii;
2346: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
2347: if (iascii) {PetscDualSpaceLagrangeView_Ascii(sp, viewer);}
2348: return(0);
2349: }
2351: static PetscErrorCode PetscDualSpaceGetDimension_SingleCell_Lagrange(PetscDualSpace sp, PetscInt order, PetscInt *dim)
2352: {
2353: PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data;
2354: PetscReal D = 1.0;
2355: PetscInt n, i;
2356: PetscErrorCode ierr;
2359: *dim = -1; /* Ensure that the compiler knows *dim is set. */
2360: DMGetDimension(sp->dm, &n);
2361: if (!lag->tensorSpace) {
2362: for (i = 1; i <= n; ++i) {
2363: D *= ((PetscReal) (order+i))/i;
2364: }
2365: *dim = (PetscInt) (D + 0.5);
2366: } else {
2367: *dim = 1;
2368: for (i = 0; i < n; ++i) *dim *= (order+1);
2369: }
2370: *dim *= sp->Nc;
2371: return(0);
2372: }
2374: static PetscErrorCode PetscDualSpaceCreateHeightSubspace_Lagrange(PetscDualSpace sp, PetscInt height, PetscDualSpace *bdsp)
2375: {
2376: PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data;
2377: PetscBool continuous, tensor;
2378: PetscInt order;
2379: PetscErrorCode ierr;
2384: PetscDualSpaceLagrangeGetContinuity(sp,&continuous);
2385: PetscDualSpaceGetOrder(sp,&order);
2386: if (height == 0) {
2387: PetscObjectReference((PetscObject)sp);
2388: *bdsp = sp;
2389: } else if (continuous == PETSC_FALSE || !order) {
2390: *bdsp = NULL;
2391: } else {
2392: DM dm, K;
2393: PetscInt dim;
2395: PetscDualSpaceGetDM(sp,&dm);
2396: DMGetDimension(dm,&dim);
2397: if (height > dim || height < 0) {SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Asked for dual space at height %d for dimension %d reference element\n",height,dim);}
2398: PetscDualSpaceDuplicate(sp,bdsp);
2399: PetscDualSpaceCreateReferenceCell(*bdsp, dim-height, lag->simplexCell, &K);
2400: PetscDualSpaceSetDM(*bdsp, K);
2401: DMDestroy(&K);
2402: PetscDualSpaceLagrangeGetTensor(sp,&tensor);
2403: PetscDualSpaceLagrangeSetTensor(*bdsp,tensor);
2404: PetscDualSpaceSetUp(*bdsp);
2405: }
2406: return(0);
2407: }
2409: PetscErrorCode PetscDualSpaceSetUp_Lagrange(PetscDualSpace sp)
2410: {
2411: PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data;
2412: DM dm = sp->dm;
2413: PetscInt order = sp->order;
2414: PetscInt Nc = sp->Nc;
2415: PetscBool continuous;
2416: PetscSection csection;
2417: Vec coordinates;
2418: PetscReal *qpoints, *qweights;
2419: PetscInt depth, dim, pdimMax, pStart, pEnd, p, *pStratStart, *pStratEnd, coneSize, d, f = 0, c;
2420: PetscBool simplex, tensorSpace;
2421: PetscErrorCode ierr;
2424: /* Classify element type */
2425: if (!order) lag->continuous = PETSC_FALSE;
2426: continuous = lag->continuous;
2427: DMGetDimension(dm, &dim);
2428: DMPlexGetDepth(dm, &depth);
2429: DMPlexGetChart(dm, &pStart, &pEnd);
2430: PetscCalloc1(dim+1, &lag->numDof);
2431: PetscMalloc2(depth+1,&pStratStart,depth+1,&pStratEnd);
2432: for (d = 0; d <= depth; ++d) {DMPlexGetDepthStratum(dm, d, &pStratStart[d], &pStratEnd[d]);}
2433: DMPlexGetConeSize(dm, pStratStart[depth], &coneSize);
2434: DMGetCoordinateSection(dm, &csection);
2435: DMGetCoordinatesLocal(dm, &coordinates);
2436: if (depth == 1) {
2437: if (coneSize == dim+1) simplex = PETSC_TRUE;
2438: else if (coneSize == 1 << dim) simplex = PETSC_FALSE;
2439: else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only support simplices and tensor product cells");
2440: } else if (depth == dim) {
2441: if (coneSize == dim+1) simplex = PETSC_TRUE;
2442: else if (coneSize == 2 * dim) simplex = PETSC_FALSE;
2443: else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only support simplices and tensor product cells");
2444: } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only support cell-vertex meshes or interpolated meshes");
2445: lag->simplexCell = simplex;
2446: if (dim > 1 && continuous && lag->simplexCell == lag->tensorSpace) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "Mismatching simplex/tensor cells and spaces only allowed for discontinuous elements");
2447: tensorSpace = lag->tensorSpace;
2448: lag->height = 0;
2449: lag->subspaces = NULL;
2450: if (continuous && sp->order > 0 && dim > 0) {
2451: PetscInt i;
2453: lag->height = dim;
2454: PetscMalloc1(dim,&lag->subspaces);
2455: PetscDualSpaceCreateHeightSubspace_Lagrange(sp,1,&lag->subspaces[0]);
2456: PetscDualSpaceSetUp(lag->subspaces[0]);
2457: for (i = 1; i < dim; i++) {
2458: PetscDualSpaceGetHeightSubspace(lag->subspaces[i-1],1,&lag->subspaces[i]);
2459: PetscObjectReference((PetscObject)(lag->subspaces[i]));
2460: }
2461: }
2462: PetscDualSpaceGetDimension_SingleCell_Lagrange(sp, sp->order, &pdimMax);
2463: pdimMax *= (pStratEnd[depth] - pStratStart[depth]);
2464: PetscMalloc1(pdimMax, &sp->functional);
2465: if (!dim) {
2466: for (c = 0; c < Nc; ++c) {
2467: PetscQuadratureCreate(PETSC_COMM_SELF, &sp->functional[f]);
2468: PetscCalloc1(Nc, &qweights);
2469: PetscQuadratureSetOrder(sp->functional[f], 0);
2470: PetscQuadratureSetData(sp->functional[f], 0, Nc, 1, NULL, qweights);
2471: qweights[c] = 1.0;
2472: ++f;
2473: lag->numDof[0]++;
2474: }
2475: } else {
2476: PetscInt *tup;
2477: PetscReal *v0, *hv0, *J, *invJ, detJ, hdetJ;
2478: PetscSection section;
2480: PetscSectionCreate(PETSC_COMM_SELF,§ion);
2481: PetscSectionSetChart(section,pStart,pEnd);
2482: PetscCalloc5(dim+1,&tup,dim,&v0,dim,&hv0,dim*dim,&J,dim*dim,&invJ);
2483: for (p = pStart; p < pEnd; p++) {
2484: PetscInt pointDim, d, nFunc = 0;
2485: PetscDualSpace hsp;
2487: DMPlexComputeCellGeometryFEM(dm, p, NULL, v0, J, invJ, &detJ);
2488: for (d = 0; d < depth; d++) {if (p >= pStratStart[d] && p < pStratEnd[d]) break;}
2489: pointDim = (depth == 1 && d == 1) ? dim : d;
2490: hsp = ((pointDim < dim) && lag->subspaces) ? lag->subspaces[dim - pointDim - 1] : NULL;
2491: if (hsp) {
2492: PetscDualSpace_Lag *hlag = (PetscDualSpace_Lag *) hsp->data;
2493: DM hdm;
2495: PetscDualSpaceGetDM(hsp,&hdm);
2496: DMPlexComputeCellGeometryFEM(hdm, 0, NULL, hv0, NULL, NULL, &hdetJ);
2497: nFunc = lag->numDof[pointDim] = hlag->numDof[pointDim];
2498: }
2499: if (pointDim == dim) {
2500: /* Cells, create for self */
2501: PetscInt orderEff = continuous ? (!tensorSpace ? order-1-dim : order-2) : order;
2502: PetscReal denom = continuous ? order : (!tensorSpace ? order+1+dim : order+2);
2503: PetscReal numer = (!simplex || !tensorSpace) ? 2. : (2./dim);
2504: PetscReal dx = numer/denom;
2505: PetscInt cdim, d, d2;
2507: if (orderEff < 0) continue;
2508: PetscDualSpaceGetDimension_SingleCell_Lagrange(sp, orderEff, &cdim);
2509: PetscMemzero(tup,(dim+1)*sizeof(PetscInt));
2510: if (!tensorSpace) {
2511: while (!tup[dim]) {
2512: for (c = 0; c < Nc; ++c) {
2513: PetscQuadratureCreate(PETSC_COMM_SELF, &sp->functional[f]);
2514: PetscMalloc1(dim, &qpoints);
2515: PetscCalloc1(Nc, &qweights);
2516: PetscQuadratureSetOrder(sp->functional[f], 0);
2517: PetscQuadratureSetData(sp->functional[f], dim, Nc, 1, qpoints, qweights);
2518: for (d = 0; d < dim; ++d) {
2519: qpoints[d] = v0[d];
2520: for (d2 = 0; d2 < dim; ++d2) qpoints[d] += J[d*dim+d2]*((tup[d2]+1)*dx);
2521: }
2522: qweights[c] = 1.0;
2523: ++f;
2524: }
2525: LatticePointLexicographic_Internal(dim, orderEff, tup);
2526: }
2527: } else {
2528: while (!tup[dim]) {
2529: for (c = 0; c < Nc; ++c) {
2530: PetscQuadratureCreate(PETSC_COMM_SELF, &sp->functional[f]);
2531: PetscMalloc1(dim, &qpoints);
2532: PetscCalloc1(Nc, &qweights);
2533: PetscQuadratureSetOrder(sp->functional[f], 0);
2534: PetscQuadratureSetData(sp->functional[f], dim, Nc, 1, qpoints, qweights);
2535: for (d = 0; d < dim; ++d) {
2536: qpoints[d] = v0[d];
2537: for (d2 = 0; d2 < dim; ++d2) qpoints[d] += J[d*dim+d2]*((tup[d2]+1)*dx);
2538: }
2539: qweights[c] = 1.0;
2540: ++f;
2541: }
2542: TensorPointLexicographic_Internal(dim, orderEff, tup);
2543: }
2544: }
2545: lag->numDof[dim] = cdim;
2546: } else { /* transform functionals from subspaces */
2547: PetscInt q;
2549: for (q = 0; q < nFunc; q++, f++) {
2550: PetscQuadrature fn;
2551: PetscInt fdim, Nc, c, nPoints, i;
2552: const PetscReal *points;
2553: const PetscReal *weights;
2554: PetscReal *qpoints;
2555: PetscReal *qweights;
2557: PetscDualSpaceGetFunctional(hsp, q, &fn);
2558: PetscQuadratureGetData(fn,&fdim,&Nc,&nPoints,&points,&weights);
2559: if (fdim != pointDim) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Expected height dual space dim %D, got %D",pointDim,fdim);
2560: PetscMalloc1(nPoints * dim, &qpoints);
2561: PetscCalloc1(nPoints * Nc, &qweights);
2562: for (i = 0; i < nPoints; i++) {
2563: PetscInt j, k;
2564: PetscReal *qp = &qpoints[i * dim];
2566: for (c = 0; c < Nc; ++c) qweights[i*Nc+c] = weights[i*Nc+c];
2567: for (j = 0; j < dim; ++j) qp[j] = v0[j];
2568: for (j = 0; j < dim; ++j) {
2569: for (k = 0; k < pointDim; k++) qp[j] += J[dim * j + k] * (points[pointDim * i + k] - hv0[k]);
2570: }
2571: }
2572: PetscQuadratureCreate(PETSC_COMM_SELF, &sp->functional[f]);
2573: PetscQuadratureSetOrder(sp->functional[f],0);
2574: PetscQuadratureSetData(sp->functional[f],dim,Nc,nPoints,qpoints,qweights);
2575: }
2576: }
2577: PetscSectionSetDof(section,p,lag->numDof[pointDim]);
2578: }
2579: PetscFree5(tup,v0,hv0,J,invJ);
2580: PetscSectionSetUp(section);
2581: { /* reorder to closure order */
2582: PetscInt *key, count;
2583: PetscQuadrature *reorder = NULL;
2585: PetscCalloc1(f,&key);
2586: PetscMalloc1(f*sp->Nc,&reorder);
2588: for (p = pStratStart[depth], count = 0; p < pStratEnd[depth]; p++) {
2589: PetscInt *closure = NULL, closureSize, c;
2591: DMPlexGetTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);
2592: for (c = 0; c < closureSize; c++) {
2593: PetscInt point = closure[2 * c], dof, off, i;
2595: PetscSectionGetDof(section,point,&dof);
2596: PetscSectionGetOffset(section,point,&off);
2597: for (i = 0; i < dof; i++) {
2598: PetscInt fi = i + off;
2599: if (!key[fi]) {
2600: key[fi] = 1;
2601: reorder[count++] = sp->functional[fi];
2602: }
2603: }
2604: }
2605: DMPlexRestoreTransitiveClosure(dm,p,PETSC_TRUE,&closureSize,&closure);
2606: }
2607: PetscFree(sp->functional);
2608: sp->functional = reorder;
2609: PetscFree(key);
2610: }
2611: PetscSectionDestroy(§ion);
2612: }
2613: if (pStratEnd[depth] == 1 && f != pdimMax) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of dual basis vectors %d not equal to dimension %d", f, pdimMax);
2614: PetscFree2(pStratStart, pStratEnd);
2615: if (f > pdimMax) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of dual basis vectors %d is greater than dimension %d", f, pdimMax);
2616: return(0);
2617: }
2619: PetscErrorCode PetscDualSpaceDestroy_Lagrange(PetscDualSpace sp)
2620: {
2621: PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data;
2622: PetscInt i;
2623: PetscErrorCode ierr;
2626: if (lag->symmetries) {
2627: PetscInt **selfSyms = lag->symmetries[0];
2629: if (selfSyms) {
2630: PetscInt i, **allocated = &selfSyms[-lag->selfSymOff];
2632: for (i = 0; i < lag->numSelfSym; i++) {
2633: PetscFree(allocated[i]);
2634: }
2635: PetscFree(allocated);
2636: }
2637: PetscFree(lag->symmetries);
2638: }
2639: for (i = 0; i < lag->height; i++) {
2640: PetscDualSpaceDestroy(&lag->subspaces[i]);
2641: }
2642: PetscFree(lag->subspaces);
2643: PetscFree(lag->numDof);
2644: PetscFree(lag);
2645: PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetContinuity_C", NULL);
2646: PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetContinuity_C", NULL);
2647: PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetTensor_C", NULL);
2648: PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetTensor_C", NULL);
2649: return(0);
2650: }
2652: PetscErrorCode PetscDualSpaceDuplicate_Lagrange(PetscDualSpace sp, PetscDualSpace *spNew)
2653: {
2654: PetscInt order, Nc;
2655: PetscBool cont, tensor;
2659: PetscDualSpaceCreate(PetscObjectComm((PetscObject) sp), spNew);
2660: PetscDualSpaceSetType(*spNew, PETSCDUALSPACELAGRANGE);
2661: PetscDualSpaceGetOrder(sp, &order);
2662: PetscDualSpaceSetOrder(*spNew, order);
2663: PetscDualSpaceGetNumComponents(sp, &Nc);
2664: PetscDualSpaceSetNumComponents(*spNew, Nc);
2665: PetscDualSpaceLagrangeGetContinuity(sp, &cont);
2666: PetscDualSpaceLagrangeSetContinuity(*spNew, cont);
2667: PetscDualSpaceLagrangeGetTensor(sp, &tensor);
2668: PetscDualSpaceLagrangeSetTensor(*spNew, tensor);
2669: return(0);
2670: }
2672: PetscErrorCode PetscDualSpaceSetFromOptions_Lagrange(PetscOptionItems *PetscOptionsObject,PetscDualSpace sp)
2673: {
2674: PetscBool continuous, tensor, flg;
2678: PetscDualSpaceLagrangeGetContinuity(sp, &continuous);
2679: PetscDualSpaceLagrangeGetTensor(sp, &tensor);
2680: PetscOptionsHead(PetscOptionsObject,"PetscDualSpace Lagrange Options");
2681: PetscOptionsBool("-petscdualspace_lagrange_continuity", "Flag for continuous element", "PetscDualSpaceLagrangeSetContinuity", continuous, &continuous, &flg);
2682: if (flg) {PetscDualSpaceLagrangeSetContinuity(sp, continuous);}
2683: PetscOptionsBool("-petscdualspace_lagrange_tensor", "Flag for tensor dual space", "PetscDualSpaceLagrangeSetContinuity", tensor, &tensor, &flg);
2684: if (flg) {PetscDualSpaceLagrangeSetTensor(sp, tensor);}
2685: PetscOptionsTail();
2686: return(0);
2687: }
2689: PetscErrorCode PetscDualSpaceGetDimension_Lagrange(PetscDualSpace sp, PetscInt *dim)
2690: {
2691: DM K;
2692: const PetscInt *numDof;
2693: PetscInt spatialDim, Nc, size = 0, d;
2694: PetscErrorCode ierr;
2697: PetscDualSpaceGetDM(sp, &K);
2698: PetscDualSpaceGetNumDof(sp, &numDof);
2699: DMGetDimension(K, &spatialDim);
2700: DMPlexGetHeightStratum(K, 0, NULL, &Nc);
2701: if (Nc == 1) {PetscDualSpaceGetDimension_SingleCell_Lagrange(sp, sp->order, dim); return(0);}
2702: for (d = 0; d <= spatialDim; ++d) {
2703: PetscInt pStart, pEnd;
2705: DMPlexGetDepthStratum(K, d, &pStart, &pEnd);
2706: size += (pEnd-pStart)*numDof[d];
2707: }
2708: *dim = size;
2709: return(0);
2710: }
2712: PetscErrorCode PetscDualSpaceGetNumDof_Lagrange(PetscDualSpace sp, const PetscInt **numDof)
2713: {
2714: PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data;
2717: *numDof = lag->numDof;
2718: return(0);
2719: }
2721: static PetscErrorCode PetscDualSpaceLagrangeGetContinuity_Lagrange(PetscDualSpace sp, PetscBool *continuous)
2722: {
2723: PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data;
2728: *continuous = lag->continuous;
2729: return(0);
2730: }
2732: static PetscErrorCode PetscDualSpaceLagrangeSetContinuity_Lagrange(PetscDualSpace sp, PetscBool continuous)
2733: {
2734: PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data;
2738: lag->continuous = continuous;
2739: return(0);
2740: }
2742: /*@
2743: PetscDualSpaceLagrangeGetContinuity - Retrieves the flag for element continuity
2745: Not Collective
2747: Input Parameter:
2748: . sp - the PetscDualSpace
2750: Output Parameter:
2751: . continuous - flag for element continuity
2753: Level: intermediate
2755: .keywords: PetscDualSpace, Lagrange, continuous, discontinuous
2756: .seealso: PetscDualSpaceLagrangeSetContinuity()
2757: @*/
2758: PetscErrorCode PetscDualSpaceLagrangeGetContinuity(PetscDualSpace sp, PetscBool *continuous)
2759: {
2765: PetscUseMethod(sp, "PetscDualSpaceLagrangeGetContinuity_C", (PetscDualSpace,PetscBool*),(sp,continuous));
2766: return(0);
2767: }
2769: /*@
2770: PetscDualSpaceLagrangeSetContinuity - Indicate whether the element is continuous
2772: Logically Collective on PetscDualSpace
2774: Input Parameters:
2775: + sp - the PetscDualSpace
2776: - continuous - flag for element continuity
2778: Options Database:
2779: . -petscdualspace_lagrange_continuity <bool>
2781: Level: intermediate
2783: .keywords: PetscDualSpace, Lagrange, continuous, discontinuous
2784: .seealso: PetscDualSpaceLagrangeGetContinuity()
2785: @*/
2786: PetscErrorCode PetscDualSpaceLagrangeSetContinuity(PetscDualSpace sp, PetscBool continuous)
2787: {
2793: PetscTryMethod(sp, "PetscDualSpaceLagrangeSetContinuity_C", (PetscDualSpace,PetscBool),(sp,continuous));
2794: return(0);
2795: }
2797: PetscErrorCode PetscDualSpaceGetHeightSubspace_Lagrange(PetscDualSpace sp, PetscInt height, PetscDualSpace *bdsp)
2798: {
2799: PetscDualSpace_Lag *lag = (PetscDualSpace_Lag *) sp->data;
2800: PetscErrorCode ierr;
2805: if (height == 0) {
2806: *bdsp = sp;
2807: }
2808: else {
2809: DM dm;
2810: PetscInt dim;
2812: PetscDualSpaceGetDM(sp,&dm);
2813: DMGetDimension(dm,&dim);
2814: if (height > dim || height < 0) {SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Asked for dual space at height %d for dimension %d reference element\n",height,dim);}
2815: if (height <= lag->height) {
2816: *bdsp = lag->subspaces[height-1];
2817: }
2818: else {
2819: *bdsp = NULL;
2820: }
2821: }
2822: return(0);
2823: }
2825: PetscErrorCode PetscDualSpaceInitialize_Lagrange(PetscDualSpace sp)
2826: {
2828: sp->ops->setfromoptions = PetscDualSpaceSetFromOptions_Lagrange;
2829: sp->ops->setup = PetscDualSpaceSetUp_Lagrange;
2830: sp->ops->view = PetscDualSpaceView_Lagrange;
2831: sp->ops->destroy = PetscDualSpaceDestroy_Lagrange;
2832: sp->ops->duplicate = PetscDualSpaceDuplicate_Lagrange;
2833: sp->ops->getdimension = PetscDualSpaceGetDimension_Lagrange;
2834: sp->ops->getnumdof = PetscDualSpaceGetNumDof_Lagrange;
2835: sp->ops->getheightsubspace = PetscDualSpaceGetHeightSubspace_Lagrange;
2836: sp->ops->getsymmetries = PetscDualSpaceGetSymmetries_Lagrange;
2837: sp->ops->apply = PetscDualSpaceApplyDefault;
2838: return(0);
2839: }
2841: /*MC
2842: PETSCDUALSPACELAGRANGE = "lagrange" - A PetscDualSpace object that encapsulates a dual space of pointwise evaluation functionals
2844: Level: intermediate
2846: .seealso: PetscDualSpaceType, PetscDualSpaceCreate(), PetscDualSpaceSetType()
2847: M*/
2849: PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate_Lagrange(PetscDualSpace sp)
2850: {
2851: PetscDualSpace_Lag *lag;
2852: PetscErrorCode ierr;
2856: PetscNewLog(sp,&lag);
2857: sp->data = lag;
2859: lag->numDof = NULL;
2860: lag->simplexCell = PETSC_TRUE;
2861: lag->tensorSpace = PETSC_FALSE;
2862: lag->continuous = PETSC_TRUE;
2864: PetscDualSpaceInitialize_Lagrange(sp);
2865: PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetContinuity_C", PetscDualSpaceLagrangeGetContinuity_Lagrange);
2866: PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetContinuity_C", PetscDualSpaceLagrangeSetContinuity_Lagrange);
2867: PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeGetTensor_C", PetscDualSpaceLagrangeGetTensor_Lagrange);
2868: PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceLagrangeSetTensor_C", PetscDualSpaceLagrangeSetTensor_Lagrange);
2869: return(0);
2870: }
2872: PetscErrorCode PetscDualSpaceSetUp_Simple(PetscDualSpace sp)
2873: {
2874: PetscDualSpace_Simple *s = (PetscDualSpace_Simple *) sp->data;
2875: DM dm = sp->dm;
2876: PetscInt dim;
2877: PetscErrorCode ierr;
2880: DMGetDimension(dm, &dim);
2881: PetscCalloc1(dim+1, &s->numDof);
2882: return(0);
2883: }
2885: PetscErrorCode PetscDualSpaceDestroy_Simple(PetscDualSpace sp)
2886: {
2887: PetscDualSpace_Simple *s = (PetscDualSpace_Simple *) sp->data;
2888: PetscErrorCode ierr;
2891: PetscFree(s->numDof);
2892: PetscFree(s);
2893: PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceSimpleSetDimension_C", NULL);
2894: PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceSimpleSetFunctional_C", NULL);
2895: return(0);
2896: }
2898: PetscErrorCode PetscDualSpaceDuplicate_Simple(PetscDualSpace sp, PetscDualSpace *spNew)
2899: {
2900: PetscInt dim, d, Nc;
2904: PetscDualSpaceCreate(PetscObjectComm((PetscObject) sp), spNew);
2905: PetscDualSpaceSetType(*spNew, PETSCDUALSPACESIMPLE);
2906: PetscDualSpaceGetNumComponents(sp, &Nc);
2907: PetscDualSpaceSetNumComponents(sp, Nc);
2908: PetscDualSpaceGetDimension(sp, &dim);
2909: PetscDualSpaceSimpleSetDimension(*spNew, dim);
2910: for (d = 0; d < dim; ++d) {
2911: PetscQuadrature q;
2913: PetscDualSpaceGetFunctional(sp, d, &q);
2914: PetscDualSpaceSimpleSetFunctional(*spNew, d, q);
2915: }
2916: return(0);
2917: }
2919: PetscErrorCode PetscDualSpaceSetFromOptions_Simple(PetscOptionItems *PetscOptionsObject,PetscDualSpace sp)
2920: {
2922: return(0);
2923: }
2925: PetscErrorCode PetscDualSpaceGetDimension_Simple(PetscDualSpace sp, PetscInt *dim)
2926: {
2927: PetscDualSpace_Simple *s = (PetscDualSpace_Simple *) sp->data;
2930: *dim = s->dim;
2931: return(0);
2932: }
2934: PetscErrorCode PetscDualSpaceSimpleSetDimension_Simple(PetscDualSpace sp, const PetscInt dim)
2935: {
2936: PetscDualSpace_Simple *s = (PetscDualSpace_Simple *) sp->data;
2937: DM dm;
2938: PetscInt spatialDim, f;
2939: PetscErrorCode ierr;
2942: for (f = 0; f < s->dim; ++f) {PetscQuadratureDestroy(&sp->functional[f]);}
2943: PetscFree(sp->functional);
2944: s->dim = dim;
2945: PetscCalloc1(s->dim, &sp->functional);
2946: PetscFree(s->numDof);
2947: PetscDualSpaceGetDM(sp, &dm);
2948: DMGetCoordinateDim(dm, &spatialDim);
2949: PetscCalloc1(spatialDim+1, &s->numDof);
2950: s->numDof[spatialDim] = dim;
2951: return(0);
2952: }
2954: PetscErrorCode PetscDualSpaceGetNumDof_Simple(PetscDualSpace sp, const PetscInt **numDof)
2955: {
2956: PetscDualSpace_Simple *s = (PetscDualSpace_Simple *) sp->data;
2959: *numDof = s->numDof;
2960: return(0);
2961: }
2963: PetscErrorCode PetscDualSpaceSimpleSetFunctional_Simple(PetscDualSpace sp, PetscInt f, PetscQuadrature q)
2964: {
2965: PetscDualSpace_Simple *s = (PetscDualSpace_Simple *) sp->data;
2966: PetscReal *weights;
2967: PetscInt Nc, c, Nq, p;
2968: PetscErrorCode ierr;
2971: if ((f < 0) || (f >= s->dim)) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_OUTOFRANGE, "Basis index %d not in [0, %d)", f, s->dim);
2972: PetscQuadratureDuplicate(q, &sp->functional[f]);
2973: /* Reweight so that it has unit volume: Do we want to do this for Nc > 1? */
2974: PetscQuadratureGetData(sp->functional[f], NULL, &Nc, &Nq, NULL, (const PetscReal **) &weights);
2975: for (c = 0; c < Nc; ++c) {
2976: PetscReal vol = 0.0;
2978: for (p = 0; p < Nq; ++p) vol += weights[p*Nc+c];
2979: for (p = 0; p < Nq; ++p) weights[p*Nc+c] /= (vol == 0.0 ? 1.0 : vol);
2980: }
2981: return(0);
2982: }
2984: /*@
2985: PetscDualSpaceSimpleSetDimension - Set the number of functionals in the dual space basis
2987: Logically Collective on PetscDualSpace
2989: Input Parameters:
2990: + sp - the PetscDualSpace
2991: - dim - the basis dimension
2993: Level: intermediate
2995: .keywords: PetscDualSpace, dimension
2996: .seealso: PetscDualSpaceSimpleSetFunctional()
2997: @*/
2998: PetscErrorCode PetscDualSpaceSimpleSetDimension(PetscDualSpace sp, PetscInt dim)
2999: {
3005: PetscTryMethod(sp, "PetscDualSpaceSimpleSetDimension_C", (PetscDualSpace,PetscInt),(sp,dim));
3006: return(0);
3007: }
3009: /*@
3010: PetscDualSpaceSimpleSetFunctional - Set the given basis element for this dual space
3012: Not Collective
3014: Input Parameters:
3015: + sp - the PetscDualSpace
3016: . f - the basis index
3017: - q - the basis functional
3019: Level: intermediate
3021: Note: The quadrature will be reweighted so that it has unit volume.
3023: .keywords: PetscDualSpace, functional
3024: .seealso: PetscDualSpaceSimpleSetDimension()
3025: @*/
3026: PetscErrorCode PetscDualSpaceSimpleSetFunctional(PetscDualSpace sp, PetscInt func, PetscQuadrature q)
3027: {
3032: PetscTryMethod(sp, "PetscDualSpaceSimpleSetFunctional_C", (PetscDualSpace,PetscInt,PetscQuadrature),(sp,func,q));
3033: return(0);
3034: }
3036: PetscErrorCode PetscDualSpaceInitialize_Simple(PetscDualSpace sp)
3037: {
3039: sp->ops->setfromoptions = PetscDualSpaceSetFromOptions_Simple;
3040: sp->ops->setup = PetscDualSpaceSetUp_Simple;
3041: sp->ops->view = NULL;
3042: sp->ops->destroy = PetscDualSpaceDestroy_Simple;
3043: sp->ops->duplicate = PetscDualSpaceDuplicate_Simple;
3044: sp->ops->getdimension = PetscDualSpaceGetDimension_Simple;
3045: sp->ops->getnumdof = PetscDualSpaceGetNumDof_Simple;
3046: sp->ops->getheightsubspace = NULL;
3047: sp->ops->getsymmetries = NULL;
3048: sp->ops->apply = PetscDualSpaceApplyDefault;
3049: return(0);
3050: }
3052: /*MC
3053: PETSCDUALSPACESIMPLE = "simple" - A PetscDualSpace object that encapsulates a dual space of arbitrary functionals
3055: Level: intermediate
3057: .seealso: PetscDualSpaceType, PetscDualSpaceCreate(), PetscDualSpaceSetType()
3058: M*/
3060: PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate_Simple(PetscDualSpace sp)
3061: {
3062: PetscDualSpace_Simple *s;
3063: PetscErrorCode ierr;
3067: PetscNewLog(sp,&s);
3068: sp->data = s;
3070: s->dim = 0;
3071: s->numDof = NULL;
3073: PetscDualSpaceInitialize_Simple(sp);
3074: PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceSimpleSetDimension_C", PetscDualSpaceSimpleSetDimension_Simple);
3075: PetscObjectComposeFunction((PetscObject) sp, "PetscDualSpaceSimpleSetFunctional_C", PetscDualSpaceSimpleSetFunctional_Simple);
3076: return(0);
3077: }
3080: PetscClassId PETSCFE_CLASSID = 0;
3082: PetscFunctionList PetscFEList = NULL;
3083: PetscBool PetscFERegisterAllCalled = PETSC_FALSE;
3085: /*@C
3086: PetscFERegister - Adds a new PetscFE implementation
3088: Not Collective
3090: Input Parameters:
3091: + name - The name of a new user-defined creation routine
3092: - create_func - The creation routine itself
3094: Notes:
3095: PetscFERegister() may be called multiple times to add several user-defined PetscFEs
3097: Sample usage:
3098: .vb
3099: PetscFERegister("my_fe", MyPetscFECreate);
3100: .ve
3102: Then, your PetscFE type can be chosen with the procedural interface via
3103: .vb
3104: PetscFECreate(MPI_Comm, PetscFE *);
3105: PetscFESetType(PetscFE, "my_fe");
3106: .ve
3107: or at runtime via the option
3108: .vb
3109: -petscfe_type my_fe
3110: .ve
3112: Level: advanced
3114: .keywords: PetscFE, register
3115: .seealso: PetscFERegisterAll(), PetscFERegisterDestroy()
3117: @*/
3118: PetscErrorCode PetscFERegister(const char sname[], PetscErrorCode (*function)(PetscFE))
3119: {
3123: PetscFunctionListAdd(&PetscFEList, sname, function);
3124: return(0);
3125: }
3127: /*@C
3128: PetscFESetType - Builds a particular PetscFE
3130: Collective on PetscFE
3132: Input Parameters:
3133: + fem - The PetscFE object
3134: - name - The kind of FEM space
3136: Options Database Key:
3137: . -petscfe_type <type> - Sets the PetscFE type; use -help for a list of available types
3139: Level: intermediate
3141: .keywords: PetscFE, set, type
3142: .seealso: PetscFEGetType(), PetscFECreate()
3143: @*/
3144: PetscErrorCode PetscFESetType(PetscFE fem, PetscFEType name)
3145: {
3146: PetscErrorCode (*r)(PetscFE);
3147: PetscBool match;
3152: PetscObjectTypeCompare((PetscObject) fem, name, &match);
3153: if (match) return(0);
3155: if (!PetscFERegisterAllCalled) {PetscFERegisterAll();}
3156: PetscFunctionListFind(PetscFEList, name, &r);
3157: if (!r) SETERRQ1(PetscObjectComm((PetscObject) fem), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscFE type: %s", name);
3159: if (fem->ops->destroy) {
3160: (*fem->ops->destroy)(fem);
3161: fem->ops->destroy = NULL;
3162: }
3163: (*r)(fem);
3164: PetscObjectChangeTypeName((PetscObject) fem, name);
3165: return(0);
3166: }
3168: /*@C
3169: PetscFEGetType - Gets the PetscFE type name (as a string) from the object.
3171: Not Collective
3173: Input Parameter:
3174: . fem - The PetscFE
3176: Output Parameter:
3177: . name - The PetscFE type name
3179: Level: intermediate
3181: .keywords: PetscFE, get, type, name
3182: .seealso: PetscFESetType(), PetscFECreate()
3183: @*/
3184: PetscErrorCode PetscFEGetType(PetscFE fem, PetscFEType *name)
3185: {
3191: if (!PetscFERegisterAllCalled) {
3192: PetscFERegisterAll();
3193: }
3194: *name = ((PetscObject) fem)->type_name;
3195: return(0);
3196: }
3198: /*@C
3199: PetscFEView - Views a PetscFE
3201: Collective on PetscFE
3203: Input Parameter:
3204: + fem - the PetscFE object to view
3205: - v - the viewer
3207: Level: developer
3209: .seealso PetscFEDestroy()
3210: @*/
3211: PetscErrorCode PetscFEView(PetscFE fem, PetscViewer v)
3212: {
3217: if (!v) {PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) fem), &v);}
3218: if (fem->ops->view) {(*fem->ops->view)(fem, v);}
3219: return(0);
3220: }
3222: /*@
3223: PetscFESetFromOptions - sets parameters in a PetscFE from the options database
3225: Collective on PetscFE
3227: Input Parameter:
3228: . fem - the PetscFE object to set options for
3230: Options Database:
3231: . -petscfe_num_blocks the number of cell blocks to integrate concurrently
3232: . -petscfe_num_batches the number of cell batches to integrate serially
3234: Level: developer
3236: .seealso PetscFEView()
3237: @*/
3238: PetscErrorCode PetscFESetFromOptions(PetscFE fem)
3239: {
3240: const char *defaultType;
3241: char name[256];
3242: PetscBool flg;
3247: if (!((PetscObject) fem)->type_name) {
3248: defaultType = PETSCFEBASIC;
3249: } else {
3250: defaultType = ((PetscObject) fem)->type_name;
3251: }
3252: if (!PetscFERegisterAllCalled) {PetscFERegisterAll();}
3254: PetscObjectOptionsBegin((PetscObject) fem);
3255: PetscOptionsFList("-petscfe_type", "Finite element space", "PetscFESetType", PetscFEList, defaultType, name, 256, &flg);
3256: if (flg) {
3257: PetscFESetType(fem, name);
3258: } else if (!((PetscObject) fem)->type_name) {
3259: PetscFESetType(fem, defaultType);
3260: }
3261: PetscOptionsInt("-petscfe_num_blocks", "The number of cell blocks to integrate concurrently", "PetscSpaceSetTileSizes", fem->numBlocks, &fem->numBlocks, NULL);
3262: PetscOptionsInt("-petscfe_num_batches", "The number of cell batches to integrate serially", "PetscSpaceSetTileSizes", fem->numBatches, &fem->numBatches, NULL);
3263: if (fem->ops->setfromoptions) {
3264: (*fem->ops->setfromoptions)(PetscOptionsObject,fem);
3265: }
3266: /* process any options handlers added with PetscObjectAddOptionsHandler() */
3267: PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) fem);
3268: PetscOptionsEnd();
3269: PetscFEViewFromOptions(fem, NULL, "-petscfe_view");
3270: return(0);
3271: }
3273: /*@C
3274: PetscFESetUp - Construct data structures for the PetscFE
3276: Collective on PetscFE
3278: Input Parameter:
3279: . fem - the PetscFE object to setup
3281: Level: developer
3283: .seealso PetscFEView(), PetscFEDestroy()
3284: @*/
3285: PetscErrorCode PetscFESetUp(PetscFE fem)
3286: {
3291: if (fem->ops->setup) {(*fem->ops->setup)(fem);}
3292: return(0);
3293: }
3295: /*@
3296: PetscFEDestroy - Destroys a PetscFE object
3298: Collective on PetscFE
3300: Input Parameter:
3301: . fem - the PetscFE object to destroy
3303: Level: developer
3305: .seealso PetscFEView()
3306: @*/
3307: PetscErrorCode PetscFEDestroy(PetscFE *fem)
3308: {
3312: if (!*fem) return(0);
3315: if (--((PetscObject)(*fem))->refct > 0) {*fem = 0; return(0);}
3316: ((PetscObject) (*fem))->refct = 0;
3318: if ((*fem)->subspaces) {
3319: PetscInt dim, d;
3321: PetscDualSpaceGetDimension((*fem)->dualSpace, &dim);
3322: for (d = 0; d < dim; ++d) {PetscFEDestroy(&(*fem)->subspaces[d]);}
3323: }
3324: PetscFree((*fem)->subspaces);
3325: PetscFree((*fem)->invV);
3326: PetscFERestoreTabulation((*fem), 0, NULL, &(*fem)->B, &(*fem)->D, NULL /*&(*fem)->H*/);
3327: PetscFERestoreTabulation((*fem), 0, NULL, &(*fem)->Bf, &(*fem)->Df, NULL /*&(*fem)->Hf*/);
3328: PetscFERestoreTabulation((*fem), 0, NULL, &(*fem)->F, NULL, NULL);
3329: PetscSpaceDestroy(&(*fem)->basisSpace);
3330: PetscDualSpaceDestroy(&(*fem)->dualSpace);
3331: PetscQuadratureDestroy(&(*fem)->quadrature);
3332: PetscQuadratureDestroy(&(*fem)->faceQuadrature);
3334: if ((*fem)->ops->destroy) {(*(*fem)->ops->destroy)(*fem);}
3335: PetscHeaderDestroy(fem);
3336: return(0);
3337: }
3339: /*@
3340: PetscFECreate - Creates an empty PetscFE object. The type can then be set with PetscFESetType().
3342: Collective on MPI_Comm
3344: Input Parameter:
3345: . comm - The communicator for the PetscFE object
3347: Output Parameter:
3348: . fem - The PetscFE object
3350: Level: beginner
3352: .seealso: PetscFESetType(), PETSCFEGALERKIN
3353: @*/
3354: PetscErrorCode PetscFECreate(MPI_Comm comm, PetscFE *fem)
3355: {
3356: PetscFE f;
3361: PetscCitationsRegister(FECitation,&FEcite);
3362: *fem = NULL;
3363: PetscFEInitializePackage();
3365: PetscHeaderCreate(f, PETSCFE_CLASSID, "PetscFE", "Finite Element", "PetscFE", comm, PetscFEDestroy, PetscFEView);
3367: f->basisSpace = NULL;
3368: f->dualSpace = NULL;
3369: f->numComponents = 1;
3370: f->subspaces = NULL;
3371: f->invV = NULL;
3372: f->B = NULL;
3373: f->D = NULL;
3374: f->H = NULL;
3375: f->Bf = NULL;
3376: f->Df = NULL;
3377: f->Hf = NULL;
3378: PetscMemzero(&f->quadrature, sizeof(PetscQuadrature));
3379: PetscMemzero(&f->faceQuadrature, sizeof(PetscQuadrature));
3380: f->blockSize = 0;
3381: f->numBlocks = 1;
3382: f->batchSize = 0;
3383: f->numBatches = 1;
3385: *fem = f;
3386: return(0);
3387: }
3389: /*@
3390: PetscFEGetSpatialDimension - Returns the spatial dimension of the element
3392: Not collective
3394: Input Parameter:
3395: . fem - The PetscFE object
3397: Output Parameter:
3398: . dim - The spatial dimension
3400: Level: intermediate
3402: .seealso: PetscFECreate()
3403: @*/
3404: PetscErrorCode PetscFEGetSpatialDimension(PetscFE fem, PetscInt *dim)
3405: {
3406: DM dm;
3412: PetscDualSpaceGetDM(fem->dualSpace, &dm);
3413: DMGetDimension(dm, dim);
3414: return(0);
3415: }
3417: /*@
3418: PetscFESetNumComponents - Sets the number of components in the element
3420: Not collective
3422: Input Parameters:
3423: + fem - The PetscFE object
3424: - comp - The number of field components
3426: Level: intermediate
3428: .seealso: PetscFECreate()
3429: @*/
3430: PetscErrorCode PetscFESetNumComponents(PetscFE fem, PetscInt comp)
3431: {
3434: fem->numComponents = comp;
3435: return(0);
3436: }
3438: /*@
3439: PetscFEGetNumComponents - Returns the number of components in the element
3441: Not collective
3443: Input Parameter:
3444: . fem - The PetscFE object
3446: Output Parameter:
3447: . comp - The number of field components
3449: Level: intermediate
3451: .seealso: PetscFECreate()
3452: @*/
3453: PetscErrorCode PetscFEGetNumComponents(PetscFE fem, PetscInt *comp)
3454: {
3458: *comp = fem->numComponents;
3459: return(0);
3460: }
3462: /*@
3463: PetscFESetTileSizes - Sets the tile sizes for evaluation
3465: Not collective
3467: Input Parameters:
3468: + fem - The PetscFE object
3469: . blockSize - The number of elements in a block
3470: . numBlocks - The number of blocks in a batch
3471: . batchSize - The number of elements in a batch
3472: - numBatches - The number of batches in a chunk
3474: Level: intermediate
3476: .seealso: PetscFECreate()
3477: @*/
3478: PetscErrorCode PetscFESetTileSizes(PetscFE fem, PetscInt blockSize, PetscInt numBlocks, PetscInt batchSize, PetscInt numBatches)
3479: {
3482: fem->blockSize = blockSize;
3483: fem->numBlocks = numBlocks;
3484: fem->batchSize = batchSize;
3485: fem->numBatches = numBatches;
3486: return(0);
3487: }
3489: /*@
3490: PetscFEGetTileSizes - Returns the tile sizes for evaluation
3492: Not collective
3494: Input Parameter:
3495: . fem - The PetscFE object
3497: Output Parameters:
3498: + blockSize - The number of elements in a block
3499: . numBlocks - The number of blocks in a batch
3500: . batchSize - The number of elements in a batch
3501: - numBatches - The number of batches in a chunk
3503: Level: intermediate
3505: .seealso: PetscFECreate()
3506: @*/
3507: PetscErrorCode PetscFEGetTileSizes(PetscFE fem, PetscInt *blockSize, PetscInt *numBlocks, PetscInt *batchSize, PetscInt *numBatches)
3508: {
3515: if (blockSize) *blockSize = fem->blockSize;
3516: if (numBlocks) *numBlocks = fem->numBlocks;
3517: if (batchSize) *batchSize = fem->batchSize;
3518: if (numBatches) *numBatches = fem->numBatches;
3519: return(0);
3520: }
3522: /*@
3523: PetscFEGetBasisSpace - Returns the PetscSpace used for approximation of the solution
3525: Not collective
3527: Input Parameter:
3528: . fem - The PetscFE object
3530: Output Parameter:
3531: . sp - The PetscSpace object
3533: Level: intermediate
3535: .seealso: PetscFECreate()
3536: @*/
3537: PetscErrorCode PetscFEGetBasisSpace(PetscFE fem, PetscSpace *sp)
3538: {
3542: *sp = fem->basisSpace;
3543: return(0);
3544: }
3546: /*@
3547: PetscFESetBasisSpace - Sets the PetscSpace used for approximation of the solution
3549: Not collective
3551: Input Parameters:
3552: + fem - The PetscFE object
3553: - sp - The PetscSpace object
3555: Level: intermediate
3557: .seealso: PetscFECreate()
3558: @*/
3559: PetscErrorCode PetscFESetBasisSpace(PetscFE fem, PetscSpace sp)
3560: {
3566: PetscSpaceDestroy(&fem->basisSpace);
3567: fem->basisSpace = sp;
3568: PetscObjectReference((PetscObject) fem->basisSpace);
3569: return(0);
3570: }
3572: /*@
3573: PetscFEGetDualSpace - Returns the PetscDualSpace used to define the inner product
3575: Not collective
3577: Input Parameter:
3578: . fem - The PetscFE object
3580: Output Parameter:
3581: . sp - The PetscDualSpace object
3583: Level: intermediate
3585: .seealso: PetscFECreate()
3586: @*/
3587: PetscErrorCode PetscFEGetDualSpace(PetscFE fem, PetscDualSpace *sp)
3588: {
3592: *sp = fem->dualSpace;
3593: return(0);
3594: }
3596: /*@
3597: PetscFESetDualSpace - Sets the PetscDualSpace used to define the inner product
3599: Not collective
3601: Input Parameters:
3602: + fem - The PetscFE object
3603: - sp - The PetscDualSpace object
3605: Level: intermediate
3607: .seealso: PetscFECreate()
3608: @*/
3609: PetscErrorCode PetscFESetDualSpace(PetscFE fem, PetscDualSpace sp)
3610: {
3616: PetscDualSpaceDestroy(&fem->dualSpace);
3617: fem->dualSpace = sp;
3618: PetscObjectReference((PetscObject) fem->dualSpace);
3619: return(0);
3620: }
3622: /*@
3623: PetscFEGetQuadrature - Returns the PetscQuadrature used to calculate inner products
3625: Not collective
3627: Input Parameter:
3628: . fem - The PetscFE object
3630: Output Parameter:
3631: . q - The PetscQuadrature object
3633: Level: intermediate
3635: .seealso: PetscFECreate()
3636: @*/
3637: PetscErrorCode PetscFEGetQuadrature(PetscFE fem, PetscQuadrature *q)
3638: {
3642: *q = fem->quadrature;
3643: return(0);
3644: }
3646: /*@
3647: PetscFESetQuadrature - Sets the PetscQuadrature used to calculate inner products
3649: Not collective
3651: Input Parameters:
3652: + fem - The PetscFE object
3653: - q - The PetscQuadrature object
3655: Level: intermediate
3657: .seealso: PetscFECreate()
3658: @*/
3659: PetscErrorCode PetscFESetQuadrature(PetscFE fem, PetscQuadrature q)
3660: {
3661: PetscInt Nc, qNc;
3666: PetscFEGetNumComponents(fem, &Nc);
3667: PetscQuadratureGetNumComponents(q, &qNc);
3668: if ((qNc != 1) && (Nc != qNc)) SETERRQ2(PetscObjectComm((PetscObject) fem), PETSC_ERR_ARG_SIZ, "FE components %D != Quadrature components %D and non-scalar quadrature", Nc, qNc);
3669: PetscFERestoreTabulation(fem, 0, NULL, &fem->B, &fem->D, NULL /*&(*fem)->H*/);
3670: PetscQuadratureDestroy(&fem->quadrature);
3671: fem->quadrature = q;
3672: PetscObjectReference((PetscObject) q);
3673: return(0);
3674: }
3676: /*@
3677: PetscFEGetFaceQuadrature - Returns the PetscQuadrature used to calculate inner products on faces
3679: Not collective
3681: Input Parameter:
3682: . fem - The PetscFE object
3684: Output Parameter:
3685: . q - The PetscQuadrature object
3687: Level: intermediate
3689: .seealso: PetscFECreate()
3690: @*/
3691: PetscErrorCode PetscFEGetFaceQuadrature(PetscFE fem, PetscQuadrature *q)
3692: {
3696: *q = fem->faceQuadrature;
3697: return(0);
3698: }
3700: /*@
3701: PetscFESetFaceQuadrature - Sets the PetscQuadrature used to calculate inner products on faces
3703: Not collective
3705: Input Parameters:
3706: + fem - The PetscFE object
3707: - q - The PetscQuadrature object
3709: Level: intermediate
3711: .seealso: PetscFECreate()
3712: @*/
3713: PetscErrorCode PetscFESetFaceQuadrature(PetscFE fem, PetscQuadrature q)
3714: {
3719: PetscFERestoreTabulation(fem, 0, NULL, &fem->Bf, &fem->Df, NULL /*&(*fem)->Hf*/);
3720: PetscQuadratureDestroy(&fem->faceQuadrature);
3721: fem->faceQuadrature = q;
3722: PetscObjectReference((PetscObject) q);
3723: return(0);
3724: }
3726: /*@C
3727: PetscFEGetNumDof - Returns the number of dofs (dual basis vectors) associated to mesh points on the reference cell of a given dimension
3729: Not collective
3731: Input Parameter:
3732: . fem - The PetscFE object
3734: Output Parameter:
3735: . numDof - Array with the number of dofs per dimension
3737: Level: intermediate
3739: .seealso: PetscFECreate()
3740: @*/
3741: PetscErrorCode PetscFEGetNumDof(PetscFE fem, const PetscInt **numDof)
3742: {
3748: PetscDualSpaceGetNumDof(fem->dualSpace, numDof);
3749: return(0);
3750: }
3752: /*@C
3753: PetscFEGetDefaultTabulation - Returns the tabulation of the basis functions at the quadrature points
3755: Not collective
3757: Input Parameter:
3758: . fem - The PetscFE object
3760: Output Parameters:
3761: + B - The basis function values at quadrature points
3762: . D - The basis function derivatives at quadrature points
3763: - H - The basis function second derivatives at quadrature points
3765: Note:
3766: $ B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c
3767: $ D[((p*pdim + i)*Nc + c)*dim + d] is the derivative value at point p for basis function i, component c, in direction d
3768: $ H[(((p*pdim + i)*Nc + c)*dim + d)*dim + e] is the value at point p for basis function i, component c, in directions d and e
3770: Level: intermediate
3772: .seealso: PetscFEGetTabulation(), PetscFERestoreTabulation()
3773: @*/
3774: PetscErrorCode PetscFEGetDefaultTabulation(PetscFE fem, PetscReal **B, PetscReal **D, PetscReal **H)
3775: {
3776: PetscInt npoints;
3777: const PetscReal *points;
3778: PetscErrorCode ierr;
3785: PetscQuadratureGetData(fem->quadrature, NULL, NULL, &npoints, &points, NULL);
3786: if (!fem->B) {PetscFEGetTabulation(fem, npoints, points, &fem->B, &fem->D, NULL/*&fem->H*/);}
3787: if (B) *B = fem->B;
3788: if (D) *D = fem->D;
3789: if (H) *H = fem->H;
3790: return(0);
3791: }
3793: PetscErrorCode PetscFEGetFaceTabulation(PetscFE fem, PetscReal **Bf, PetscReal **Df, PetscReal **Hf)
3794: {
3795: PetscErrorCode ierr;
3802: if (!fem->Bf) {
3803: PetscFECellGeom cgeom;
3804: PetscQuadrature fq;
3805: PetscDualSpace sp;
3806: DM dm;
3807: const PetscInt *faces;
3808: PetscInt dim, numFaces, f, npoints, q;
3809: const PetscReal *points;
3810: PetscReal *facePoints;
3812: PetscFEGetDualSpace(fem, &sp);
3813: PetscDualSpaceGetDM(sp, &dm);
3814: DMGetDimension(dm, &dim);
3815: DMPlexGetConeSize(dm, 0, &numFaces);
3816: DMPlexGetCone(dm, 0, &faces);
3817: PetscFEGetFaceQuadrature(fem, &fq);
3818: if (fq) {
3819: PetscQuadratureGetData(fq, NULL, NULL, &npoints, &points, NULL);
3820: PetscMalloc1(numFaces*npoints*dim, &facePoints);
3821: for (f = 0; f < numFaces; ++f) {
3822: DMPlexComputeCellGeometryFEM(dm, faces[f], NULL, cgeom.v0, cgeom.J, NULL, &cgeom.detJ);
3823: for (q = 0; q < npoints; ++q) CoordinatesRefToReal(dim, dim-1, cgeom.v0, cgeom.J, &points[q*(dim-1)], &facePoints[(f*npoints+q)*dim]);
3824: }
3825: PetscFEGetTabulation(fem, numFaces*npoints, facePoints, &fem->Bf, &fem->Df, NULL/*&fem->Hf*/);
3826: PetscFree(facePoints);
3827: }
3828: }
3829: if (Bf) *Bf = fem->Bf;
3830: if (Df) *Df = fem->Df;
3831: if (Hf) *Hf = fem->Hf;
3832: return(0);
3833: }
3835: PetscErrorCode PetscFEGetFaceCentroidTabulation(PetscFE fem, PetscReal **F)
3836: {
3837: PetscErrorCode ierr;
3842: if (!fem->F) {
3843: PetscDualSpace sp;
3844: DM dm;
3845: const PetscInt *cone;
3846: PetscReal *centroids;
3847: PetscInt dim, numFaces, f;
3849: PetscFEGetDualSpace(fem, &sp);
3850: PetscDualSpaceGetDM(sp, &dm);
3851: DMGetDimension(dm, &dim);
3852: DMPlexGetConeSize(dm, 0, &numFaces);
3853: DMPlexGetCone(dm, 0, &cone);
3854: PetscMalloc1(numFaces*dim, ¢roids);
3855: for (f = 0; f < numFaces; ++f) {DMPlexComputeCellGeometryFVM(dm, cone[f], NULL, ¢roids[f*dim], NULL);}
3856: PetscFEGetTabulation(fem, numFaces, centroids, &fem->F, NULL, NULL);
3857: PetscFree(centroids);
3858: }
3859: *F = fem->F;
3860: return(0);
3861: }
3863: /*@C
3864: PetscFEGetTabulation - Tabulates the basis functions, and perhaps derivatives, at the points provided.
3866: Not collective
3868: Input Parameters:
3869: + fem - The PetscFE object
3870: . npoints - The number of tabulation points
3871: - points - The tabulation point coordinates
3873: Output Parameters:
3874: + B - The basis function values at tabulation points
3875: . D - The basis function derivatives at tabulation points
3876: - H - The basis function second derivatives at tabulation points
3878: Note:
3879: $ B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c
3880: $ D[((p*pdim + i)*Nc + c)*dim + d] is the derivative value at point p for basis function i, component c, in direction d
3881: $ H[(((p*pdim + i)*Nc + c)*dim + d)*dim + e] is the value at point p for basis function i, component c, in directions d and e
3883: Level: intermediate
3885: .seealso: PetscFERestoreTabulation(), PetscFEGetDefaultTabulation()
3886: @*/
3887: PetscErrorCode PetscFEGetTabulation(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscReal **B, PetscReal **D, PetscReal **H)
3888: {
3889: DM dm;
3890: PetscInt pdim; /* Dimension of FE space P */
3891: PetscInt dim; /* Spatial dimension */
3892: PetscInt comp; /* Field components */
3893: PetscErrorCode ierr;
3896: if (!npoints) {
3897: if (B) *B = NULL;
3898: if (D) *D = NULL;
3899: if (H) *H = NULL;
3900: return(0);
3901: }
3907: PetscDualSpaceGetDM(fem->dualSpace, &dm);
3908: DMGetDimension(dm, &dim);
3909: PetscDualSpaceGetDimension(fem->dualSpace, &pdim);
3910: PetscFEGetNumComponents(fem, &comp);
3911: if (B) {DMGetWorkArray(dm, npoints*pdim*comp, MPIU_REAL, B);}
3912: if (D) {DMGetWorkArray(dm, npoints*pdim*comp*dim, MPIU_REAL, D);}
3913: if (H) {DMGetWorkArray(dm, npoints*pdim*comp*dim*dim, MPIU_REAL, H);}
3914: (*fem->ops->gettabulation)(fem, npoints, points, B ? *B : NULL, D ? *D : NULL, H ? *H : NULL);
3915: return(0);
3916: }
3918: PetscErrorCode PetscFERestoreTabulation(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscReal **B, PetscReal **D, PetscReal **H)
3919: {
3920: DM dm;
3925: PetscDualSpaceGetDM(fem->dualSpace, &dm);
3926: if (B && *B) {DMRestoreWorkArray(dm, 0, MPIU_REAL, B);}
3927: if (D && *D) {DMRestoreWorkArray(dm, 0, MPIU_REAL, D);}
3928: if (H && *H) {DMRestoreWorkArray(dm, 0, MPIU_REAL, H);}
3929: return(0);
3930: }
3932: PetscErrorCode PetscFEDestroy_Basic(PetscFE fem)
3933: {
3934: PetscFE_Basic *b = (PetscFE_Basic *) fem->data;
3938: PetscFree(b);
3939: return(0);
3940: }
3942: PetscErrorCode PetscFEView_Basic_Ascii(PetscFE fe, PetscViewer viewer)
3943: {
3944: PetscSpace basis;
3945: PetscDualSpace dual;
3946: PetscQuadrature q = NULL;
3947: PetscInt dim, Nc, Nq;
3948: PetscViewerFormat format;
3949: PetscErrorCode ierr;
3952: PetscFEGetBasisSpace(fe, &basis);
3953: PetscFEGetDualSpace(fe, &dual);
3954: PetscFEGetQuadrature(fe, &q);
3955: PetscFEGetNumComponents(fe, &Nc);
3956: PetscQuadratureGetData(q, &dim, NULL, &Nq, NULL, NULL);
3957: PetscViewerGetFormat(viewer, &format);
3958: PetscViewerASCIIPrintf(viewer, "Basic Finite Element:\n");
3959: PetscViewerASCIIPrintf(viewer, " dimension: %d\n", dim);
3960: PetscViewerASCIIPrintf(viewer, " components: %d\n", Nc);
3961: PetscViewerASCIIPrintf(viewer, " num quad points: %d\n", Nq);
3962: if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
3963: PetscViewerASCIIPushTab(viewer);
3964: PetscQuadratureView(q, viewer);
3965: PetscViewerASCIIPopTab(viewer);
3966: }
3967: PetscViewerASCIIPushTab(viewer);
3968: PetscSpaceView(basis, viewer);
3969: PetscDualSpaceView(dual, viewer);
3970: PetscViewerASCIIPopTab(viewer);
3971: return(0);
3972: }
3974: PetscErrorCode PetscFEView_Basic(PetscFE fe, PetscViewer viewer)
3975: {
3976: PetscBool iascii;
3982: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
3983: if (iascii) {PetscFEView_Basic_Ascii(fe, viewer);}
3984: return(0);
3985: }
3987: /* Construct the change of basis from prime basis to nodal basis */
3988: PetscErrorCode PetscFESetUp_Basic(PetscFE fem)
3989: {
3990: PetscScalar *work, *invVscalar;
3991: PetscBLASInt *pivots;
3992: PetscBLASInt n, info;
3993: PetscInt pdim, j;
3997: PetscDualSpaceGetDimension(fem->dualSpace, &pdim);
3998: PetscMalloc1(pdim*pdim,&fem->invV);
3999: #if defined(PETSC_USE_COMPLEX)
4000: PetscMalloc1(pdim*pdim,&invVscalar);
4001: #else
4002: invVscalar = fem->invV;
4003: #endif
4004: for (j = 0; j < pdim; ++j) {
4005: PetscReal *Bf;
4006: PetscQuadrature f;
4007: const PetscReal *points, *weights;
4008: PetscInt Nc, Nq, q, k, c;
4010: PetscDualSpaceGetFunctional(fem->dualSpace, j, &f);
4011: PetscQuadratureGetData(f, NULL, &Nc, &Nq, &points, &weights);
4012: PetscMalloc1(Nc*Nq*pdim,&Bf);
4013: PetscSpaceEvaluate(fem->basisSpace, Nq, points, Bf, NULL, NULL);
4014: for (k = 0; k < pdim; ++k) {
4015: /* V_{jk} = n_j(\phi_k) = \int \phi_k(x) n_j(x) dx */
4016: invVscalar[j*pdim+k] = 0.0;
4018: for (q = 0; q < Nq; ++q) {
4019: for (c = 0; c < Nc; ++c) invVscalar[j*pdim+k] += Bf[(q*pdim + k)*Nc + c]*weights[q*Nc + c];
4020: }
4021: }
4022: PetscFree(Bf);
4023: }
4024: PetscMalloc2(pdim,&pivots,pdim,&work);
4025: n = pdim;
4026: PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&n, &n, invVscalar, &n, pivots, &info));
4027: PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&n, invVscalar, &n, pivots, work, &n, &info));
4028: #if defined(PETSC_USE_COMPLEX)
4029: for (j = 0; j < pdim*pdim; j++) fem->invV[j] = PetscRealPart(invVscalar[j]);
4030: PetscFree(invVscalar);
4031: #endif
4032: PetscFree2(pivots,work);
4033: return(0);
4034: }
4036: PetscErrorCode PetscFEGetDimension_Basic(PetscFE fem, PetscInt *dim)
4037: {
4041: PetscDualSpaceGetDimension(fem->dualSpace, dim);
4042: return(0);
4043: }
4045: PetscErrorCode PetscFEGetTabulation_Basic(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscReal *B, PetscReal *D, PetscReal *H)
4046: {
4047: DM dm;
4048: PetscInt pdim; /* Dimension of FE space P */
4049: PetscInt dim; /* Spatial dimension */
4050: PetscInt Nc; /* Field components */
4051: PetscReal *tmpB, *tmpD, *tmpH;
4052: PetscInt p, d, j, k, c;
4053: PetscErrorCode ierr;
4056: PetscDualSpaceGetDM(fem->dualSpace, &dm);
4057: DMGetDimension(dm, &dim);
4058: PetscDualSpaceGetDimension(fem->dualSpace, &pdim);
4059: PetscFEGetNumComponents(fem, &Nc);
4060: /* Evaluate the prime basis functions at all points */
4061: if (B) {DMGetWorkArray(dm, npoints*pdim*Nc, MPIU_REAL, &tmpB);}
4062: if (D) {DMGetWorkArray(dm, npoints*pdim*Nc*dim, MPIU_REAL, &tmpD);}
4063: if (H) {DMGetWorkArray(dm, npoints*pdim*Nc*dim*dim, MPIU_REAL, &tmpH);}
4064: PetscSpaceEvaluate(fem->basisSpace, npoints, points, B ? tmpB : NULL, D ? tmpD : NULL, H ? tmpH : NULL);
4065: /* Translate to the nodal basis */
4066: for (p = 0; p < npoints; ++p) {
4067: if (B) {
4068: /* Multiply by V^{-1} (pdim x pdim) */
4069: for (j = 0; j < pdim; ++j) {
4070: const PetscInt i = (p*pdim + j)*Nc;
4072: for (c = 0; c < Nc; ++c) {
4073: B[i+c] = 0.0;
4074: for (k = 0; k < pdim; ++k) {
4075: B[i+c] += fem->invV[k*pdim+j] * tmpB[(p*pdim + k)*Nc+c];
4076: }
4077: }
4078: }
4079: }
4080: if (D) {
4081: /* Multiply by V^{-1} (pdim x pdim) */
4082: for (j = 0; j < pdim; ++j) {
4083: for (c = 0; c < Nc; ++c) {
4084: for (d = 0; d < dim; ++d) {
4085: const PetscInt i = ((p*pdim + j)*Nc + c)*dim + d;
4087: D[i] = 0.0;
4088: for (k = 0; k < pdim; ++k) {
4089: D[i] += fem->invV[k*pdim+j] * tmpD[((p*pdim + k)*Nc + c)*dim + d];
4090: }
4091: }
4092: }
4093: }
4094: }
4095: if (H) {
4096: /* Multiply by V^{-1} (pdim x pdim) */
4097: for (j = 0; j < pdim; ++j) {
4098: for (c = 0; c < Nc; ++c) {
4099: for (d = 0; d < dim*dim; ++d) {
4100: const PetscInt i = ((p*pdim + j)*Nc + c)*dim*dim + d;
4102: H[i] = 0.0;
4103: for (k = 0; k < pdim; ++k) {
4104: H[i] += fem->invV[k*pdim+j] * tmpH[((p*pdim + k)*Nc + c)*dim*dim + d];
4105: }
4106: }
4107: }
4108: }
4109: }
4110: }
4111: if (B) {DMRestoreWorkArray(dm, npoints*pdim*Nc, MPIU_REAL, &tmpB);}
4112: if (D) {DMRestoreWorkArray(dm, npoints*pdim*Nc*dim, MPIU_REAL, &tmpD);}
4113: if (H) {DMRestoreWorkArray(dm, npoints*pdim*Nc*dim*dim, MPIU_REAL, &tmpH);}
4114: return(0);
4115: }
4117: PetscErrorCode PetscFEIntegrate_Basic(PetscFE fem, PetscDS prob, PetscInt field, PetscInt Ne, PetscFECellGeom *cgeom,
4118: const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
4119: {
4120: const PetscInt debug = 0;
4121: PetscPointFunc obj_func;
4122: PetscQuadrature quad;
4123: PetscScalar *u, *u_x, *a, *a_x, *refSpaceDer, *refSpaceDerAux;
4124: const PetscScalar *constants;
4125: PetscReal *x;
4126: PetscReal **B, **D, **BAux = NULL, **DAux = NULL;
4127: PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL;
4128: PetscInt dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, e;
4129: PetscErrorCode ierr;
4132: PetscDSGetObjective(prob, field, &obj_func);
4133: if (!obj_func) return(0);
4134: PetscFEGetSpatialDimension(fem, &dim);
4135: PetscFEGetQuadrature(fem, &quad);
4136: PetscDSGetNumFields(prob, &Nf);
4137: PetscDSGetTotalDimension(prob, &totDim);
4138: PetscDSGetDimensions(prob, &Nb);
4139: PetscDSGetComponents(prob, &Nc);
4140: PetscDSGetComponentOffsets(prob, &uOff);
4141: PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);
4142: PetscDSGetEvaluationArrays(prob, &u, NULL, &u_x);
4143: PetscDSGetRefCoordArrays(prob, &x, &refSpaceDer);
4144: PetscDSGetTabulation(prob, &B, &D);
4145: PetscDSGetConstants(prob, &numConstants, &constants);
4146: if (probAux) {
4147: PetscDSGetNumFields(probAux, &NfAux);
4148: PetscDSGetTotalDimension(probAux, &totDimAux);
4149: PetscDSGetDimensions(probAux, &NbAux);
4150: PetscDSGetComponents(probAux, &NcAux);
4151: PetscDSGetComponentOffsets(probAux, &aOff);
4152: PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);
4153: PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);
4154: PetscDSGetRefCoordArrays(probAux, NULL, &refSpaceDerAux);
4155: PetscDSGetTabulation(probAux, &BAux, &DAux);
4156: }
4157: for (e = 0; e < Ne; ++e) {
4158: const PetscReal *v0 = cgeom[e].v0;
4159: const PetscReal *J = cgeom[e].J;
4160: const PetscReal *invJ = cgeom[e].invJ;
4161: const PetscReal detJ = cgeom[e].detJ;
4162: const PetscReal *quadPoints, *quadWeights;
4163: PetscInt qNc, Nq, q;
4165: PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);
4166: if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
4167: if (debug > 1) {
4168: PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", detJ);
4169: #ifndef PETSC_USE_COMPLEX
4170: DMPrintCellMatrix(e, "invJ", dim, dim, invJ);
4171: #endif
4172: }
4173: for (q = 0; q < Nq; ++q) {
4174: PetscScalar integrand;
4176: if (debug) {PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);}
4177: CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], x);
4178: EvaluateFieldJets(dim, Nf, Nb, Nc, q, B, D, refSpaceDer, invJ, &coefficients[cOffset], NULL, u, u_x, NULL);
4179: if (probAux) EvaluateFieldJets(dim, NfAux, NbAux, NcAux, q, BAux, DAux, refSpaceDerAux, invJ, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);
4180: obj_func(dim, Nf, NfAux, uOff, uOff_x, u, NULL, u_x, aOff, aOff_x, a, NULL, a_x, 0.0, x, numConstants, constants, &integrand);
4181: integrand *= detJ*quadWeights[q];
4182: integral[e*Nf+field] += integrand;
4183: if (debug > 1) {PetscPrintf(PETSC_COMM_SELF, " int: %g %g\n", (double) PetscRealPart(integrand), (double) PetscRealPart(integral[field]));}
4184: }
4185: cOffset += totDim;
4186: cOffsetAux += totDimAux;
4187: }
4188: return(0);
4189: }
4191: PetscErrorCode PetscFEIntegrateResidual_Basic(PetscFE fem, PetscDS prob, PetscInt field, PetscInt Ne, PetscFECellGeom *cgeom,
4192: const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
4193: {
4194: const PetscInt debug = 0;
4195: PetscPointFunc f0_func;
4196: PetscPointFunc f1_func;
4197: PetscQuadrature quad;
4198: PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *refSpaceDer, *refSpaceDerAux;
4199: const PetscScalar *constants;
4200: PetscReal *x;
4201: PetscReal **B, **D, **BAux = NULL, **DAux = NULL, *BI, *DI;
4202: PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL;
4203: PetscInt dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NbI, NcI;
4204: PetscErrorCode ierr;
4207: PetscFEGetSpatialDimension(fem, &dim);
4208: PetscFEGetQuadrature(fem, &quad);
4209: PetscDSGetNumFields(prob, &Nf);
4210: PetscDSGetTotalDimension(prob, &totDim);
4211: PetscDSGetDimensions(prob, &Nb);
4212: PetscDSGetComponents(prob, &Nc);
4213: PetscDSGetComponentOffsets(prob, &uOff);
4214: PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);
4215: PetscDSGetFieldOffset(prob, field, &fOffset);
4216: PetscDSGetResidual(prob, field, &f0_func, &f1_func);
4217: PetscDSGetEvaluationArrays(prob, &u, coefficients_t ? &u_t : NULL, &u_x);
4218: PetscDSGetRefCoordArrays(prob, &x, &refSpaceDer);
4219: PetscDSGetWeakFormArrays(prob, &f0, &f1, NULL, NULL, NULL, NULL);
4220: PetscDSGetTabulation(prob, &B, &D);
4221: PetscDSGetConstants(prob, &numConstants, &constants);
4222: if (probAux) {
4223: PetscDSGetNumFields(probAux, &NfAux);
4224: PetscDSGetTotalDimension(probAux, &totDimAux);
4225: PetscDSGetDimensions(probAux, &NbAux);
4226: PetscDSGetComponents(probAux, &NcAux);
4227: PetscDSGetComponentOffsets(probAux, &aOff);
4228: PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);
4229: PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);
4230: PetscDSGetRefCoordArrays(probAux, NULL, &refSpaceDerAux);
4231: PetscDSGetTabulation(probAux, &BAux, &DAux);
4232: }
4233: NbI = Nb[field];
4234: NcI = Nc[field];
4235: BI = B[field];
4236: DI = D[field];
4237: for (e = 0; e < Ne; ++e) {
4238: const PetscReal *v0 = cgeom[e].v0;
4239: const PetscReal *J = cgeom[e].J;
4240: const PetscReal *invJ = cgeom[e].invJ;
4241: const PetscReal detJ = cgeom[e].detJ;
4242: const PetscReal *quadPoints, *quadWeights;
4243: PetscInt qNc, Nq, q;
4245: PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);
4246: if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
4247: PetscMemzero(f0, Nq*NcI* sizeof(PetscScalar));
4248: PetscMemzero(f1, Nq*NcI*dim * sizeof(PetscScalar));
4249: if (debug > 1) {
4250: PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", detJ);
4251: #ifndef PETSC_USE_COMPLEX
4252: DMPrintCellMatrix(e, "invJ", dim, dim, invJ);
4253: #endif
4254: }
4255: for (q = 0; q < Nq; ++q) {
4256: if (debug) {PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);}
4257: CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], x);
4258: EvaluateFieldJets(dim, Nf, Nb, Nc, q, B, D, refSpaceDer, invJ, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);
4259: if (probAux) EvaluateFieldJets(dim, NfAux, NbAux, NcAux, q, BAux, DAux, refSpaceDerAux, invJ, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);
4260: if (f0_func) f0_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, x, numConstants, constants, &f0[q*NcI]);
4261: if (f1_func) {
4262: PetscMemzero(refSpaceDer, NcI*dim * sizeof(PetscScalar));
4263: f1_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, x, numConstants, constants, refSpaceDer);
4264: }
4265: TransformF(dim, NcI, q, invJ, detJ, quadWeights, refSpaceDer, f0_func ? f0 : NULL, f1_func ? f1 : NULL);
4266: }
4267: UpdateElementVec(dim, Nq, NbI, NcI, BI, DI, f0, f1, &elemVec[cOffset+fOffset]);
4268: cOffset += totDim;
4269: cOffsetAux += totDimAux;
4270: }
4271: return(0);
4272: }
4274: PetscErrorCode PetscFEIntegrateBdResidual_Basic(PetscFE fem, PetscDS prob, PetscInt field, PetscInt Ne, PetscFEFaceGeom *fgeom,
4275: const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
4276: {
4277: const PetscInt debug = 0;
4278: PetscBdPointFunc f0_func;
4279: PetscBdPointFunc f1_func;
4280: PetscQuadrature quad;
4281: PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *refSpaceDer, *refSpaceDerAux;
4282: const PetscScalar *constants;
4283: PetscReal *x;
4284: PetscReal **B, **D, **BAux = NULL, **DAux = NULL, *BI, *DI;
4285: PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL;
4286: PetscInt dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NbI, NcI;
4287: PetscErrorCode ierr;
4290: PetscFEGetSpatialDimension(fem, &dim);
4291: PetscFEGetFaceQuadrature(fem, &quad);
4292: PetscDSGetNumFields(prob, &Nf);
4293: PetscDSGetTotalDimension(prob, &totDim);
4294: PetscDSGetDimensions(prob, &Nb);
4295: PetscDSGetComponents(prob, &Nc);
4296: PetscDSGetComponentOffsets(prob, &uOff);
4297: PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);
4298: PetscDSGetFieldOffset(prob, field, &fOffset);
4299: PetscDSGetBdResidual(prob, field, &f0_func, &f1_func);
4300: if (!f0_func && !f1_func) return(0);
4301: PetscDSGetEvaluationArrays(prob, &u, coefficients_t ? &u_t : NULL, &u_x);
4302: PetscDSGetRefCoordArrays(prob, &x, &refSpaceDer);
4303: PetscDSGetWeakFormArrays(prob, &f0, &f1, NULL, NULL, NULL, NULL);
4304: PetscDSGetFaceTabulation(prob, &B, &D);
4305: PetscDSGetConstants(prob, &numConstants, &constants);
4306: if (probAux) {
4307: PetscDSGetNumFields(probAux, &NfAux);
4308: PetscDSGetTotalDimension(probAux, &totDimAux);
4309: PetscDSGetDimensions(probAux, &NbAux);
4310: PetscDSGetComponents(probAux, &NcAux);
4311: PetscDSGetComponentOffsets(probAux, &aOff);
4312: PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);
4313: PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);
4314: PetscDSGetRefCoordArrays(probAux, NULL, &refSpaceDerAux);
4315: PetscDSGetFaceTabulation(probAux, &BAux, &DAux);
4316: }
4317: NbI = Nb[field];
4318: NcI = Nc[field];
4319: BI = B[field];
4320: DI = D[field];
4321: for (e = 0; e < Ne; ++e) {
4322: const PetscReal *quadPoints, *quadWeights;
4323: const PetscReal *v0 = fgeom[e].v0;
4324: const PetscReal *J = fgeom[e].J;
4325: const PetscReal *invJ = fgeom[e].invJ[0];
4326: const PetscReal detJ = fgeom[e].detJ;
4327: const PetscReal *n = fgeom[e].n;
4328: const PetscInt face = fgeom[e].face[0];
4329: PetscInt qNc, Nq, q;
4331: PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);
4332: if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
4333: PetscMemzero(f0, Nq*NcI* sizeof(PetscScalar));
4334: PetscMemzero(f1, Nq*NcI*dim * sizeof(PetscScalar));
4335: if (debug > 1) {
4336: PetscPrintf(PETSC_COMM_SELF, " detJ: %g\n", detJ);
4337: #ifndef PETSC_USE_COMPLEX
4338: DMPrintCellMatrix(e, "invJ", dim, dim, invJ);
4339: #endif
4340: }
4341: for (q = 0; q < Nq; ++q) {
4342: if (debug) {PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);}
4343: CoordinatesRefToReal(dim, dim-1, v0, J, &quadPoints[q*(dim-1)], x);
4344: EvaluateFieldJets(dim, Nf, Nb, Nc, face*Nq+q, B, D, refSpaceDer, invJ, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);
4345: if (probAux) EvaluateFieldJets(dim, NfAux, NbAux, NcAux, face*Nq+q, BAux, DAux, refSpaceDerAux, invJ, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);
4346: if (f0_func) f0_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, x, n, numConstants, constants, &f0[q*NcI]);
4347: if (f1_func) {
4348: PetscMemzero(refSpaceDer, NcI*dim * sizeof(PetscScalar));
4349: f1_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, x, n, numConstants, constants, refSpaceDer);
4350: }
4351: TransformF(dim, NcI, q, invJ, detJ, quadWeights, refSpaceDer, f0_func ? f0 : NULL, f1_func ? f1 : NULL);
4352: }
4353: UpdateElementVec(dim, Nq, NbI, NcI, &BI[face*Nq*NbI*NcI], &DI[face*Nq*NbI*NcI*dim], f0, f1, &elemVec[cOffset+fOffset]);
4354: cOffset += totDim;
4355: cOffsetAux += totDimAux;
4356: }
4357: return(0);
4358: }
4360: PetscErrorCode PetscFEIntegrateJacobian_Basic(PetscFE fem, PetscDS prob, PetscFEJacobianType jtype, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFECellGeom *geom,
4361: const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
4362: {
4363: const PetscInt debug = 0;
4364: PetscPointJac g0_func;
4365: PetscPointJac g1_func;
4366: PetscPointJac g2_func;
4367: PetscPointJac g3_func;
4368: PetscInt cOffset = 0; /* Offset into coefficients[] for element e */
4369: PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
4370: PetscInt eOffset = 0; /* Offset into elemMat[] for element e */
4371: PetscInt offsetI = 0; /* Offset into an element vector for fieldI */
4372: PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */
4373: PetscQuadrature quad;
4374: PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *refSpaceDer, *refSpaceDerAux;
4375: const PetscScalar *constants;
4376: PetscReal *x;
4377: PetscReal **B, **D, **BAux = NULL, **DAux = NULL, *BI, *DI, *BJ, *DJ;
4378: PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL;
4379: PetscInt NbI = 0, NcI = 0, NbJ = 0, NcJ = 0;
4380: PetscInt dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, e;
4381: PetscErrorCode ierr;
4384: PetscFEGetSpatialDimension(fem, &dim);
4385: PetscFEGetQuadrature(fem, &quad);
4386: PetscDSGetNumFields(prob, &Nf);
4387: PetscDSGetTotalDimension(prob, &totDim);
4388: PetscDSGetDimensions(prob, &Nb);
4389: PetscDSGetComponents(prob, &Nc);
4390: PetscDSGetComponentOffsets(prob, &uOff);
4391: PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);
4392: switch(jtype) {
4393: case PETSCFE_JACOBIAN_DYN: PetscDSGetDynamicJacobian(prob, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);break;
4394: case PETSCFE_JACOBIAN_PRE: PetscDSGetJacobianPreconditioner(prob, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);break;
4395: case PETSCFE_JACOBIAN: PetscDSGetJacobian(prob, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);break;
4396: }
4397: if (!g0_func && !g1_func && !g2_func && !g3_func) return(0);
4398: PetscDSGetEvaluationArrays(prob, &u, coefficients_t ? &u_t : NULL, &u_x);
4399: PetscDSGetRefCoordArrays(prob, &x, &refSpaceDer);
4400: PetscDSGetWeakFormArrays(prob, NULL, NULL, &g0, &g1, &g2, &g3);
4401: PetscDSGetTabulation(prob, &B, &D);
4402: PetscDSGetFieldOffset(prob, fieldI, &offsetI);
4403: PetscDSGetFieldOffset(prob, fieldJ, &offsetJ);
4404: PetscDSGetConstants(prob, &numConstants, &constants);
4405: if (probAux) {
4406: PetscDSGetNumFields(probAux, &NfAux);
4407: PetscDSGetTotalDimension(probAux, &totDimAux);
4408: PetscDSGetDimensions(probAux, &NbAux);
4409: PetscDSGetComponents(probAux, &NcAux);
4410: PetscDSGetComponentOffsets(probAux, &aOff);
4411: PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);
4412: PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);
4413: PetscDSGetRefCoordArrays(probAux, NULL, &refSpaceDerAux);
4414: PetscDSGetTabulation(probAux, &BAux, &DAux);
4415: }
4416: NbI = Nb[fieldI], NbJ = Nb[fieldJ];
4417: NcI = Nc[fieldI], NcJ = Nc[fieldJ];
4418: BI = B[fieldI], BJ = B[fieldJ];
4419: DI = D[fieldI], DJ = D[fieldJ];
4420: /* Initialize here in case the function is not defined */
4421: PetscMemzero(g0, NcI*NcJ * sizeof(PetscScalar));
4422: PetscMemzero(g1, NcI*NcJ*dim * sizeof(PetscScalar));
4423: PetscMemzero(g2, NcI*NcJ*dim * sizeof(PetscScalar));
4424: PetscMemzero(g3, NcI*NcJ*dim*dim * sizeof(PetscScalar));
4425: for (e = 0; e < Ne; ++e) {
4426: const PetscReal *v0 = geom[e].v0;
4427: const PetscReal *J = geom[e].J;
4428: const PetscReal *invJ = geom[e].invJ;
4429: const PetscReal detJ = geom[e].detJ;
4430: const PetscReal *quadPoints, *quadWeights;
4431: PetscInt qNc, Nq, q;
4433: PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);
4434: if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
4435: for (q = 0; q < Nq; ++q) {
4436: const PetscReal *BIq = &BI[q*NbI*NcI], *BJq = &BJ[q*NbJ*NcJ];
4437: const PetscReal *DIq = &DI[q*NbI*NcI*dim], *DJq = &DJ[q*NbJ*NcJ*dim];
4438: const PetscReal w = detJ*quadWeights[q];
4439: PetscInt f, g, fc, gc, c;
4441: if (debug) {PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);}
4442: CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], x);
4443: EvaluateFieldJets(dim, Nf, Nb, Nc, q, B, D, refSpaceDer, invJ, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);
4444: if (probAux) EvaluateFieldJets(dim, NfAux, NbAux, NcAux, q, BAux, DAux, refSpaceDerAux, invJ, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);
4445: if (g0_func) {
4446: PetscMemzero(g0, NcI*NcJ * sizeof(PetscScalar));
4447: g0_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, x, numConstants, constants, g0);
4448: for (c = 0; c < NcI*NcJ; ++c) g0[c] *= w;
4449: }
4450: if (g1_func) {
4451: PetscInt d, d2;
4452: PetscMemzero(refSpaceDer, NcI*NcJ*dim * sizeof(PetscScalar));
4453: g1_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, x, numConstants, constants, refSpaceDer);
4454: for (fc = 0; fc < NcI; ++fc) {
4455: for (gc = 0; gc < NcJ; ++gc) {
4456: for (d = 0; d < dim; ++d) {
4457: g1[(fc*NcJ+gc)*dim+d] = 0.0;
4458: for (d2 = 0; d2 < dim; ++d2) g1[(fc*NcJ+gc)*dim+d] += invJ[d*dim+d2]*refSpaceDer[(fc*NcJ+gc)*dim+d2];
4459: g1[(fc*NcJ+gc)*dim+d] *= w;
4460: }
4461: }
4462: }
4463: }
4464: if (g2_func) {
4465: PetscInt d, d2;
4466: PetscMemzero(refSpaceDer, NcI*NcJ*dim * sizeof(PetscScalar));
4467: g2_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, x, numConstants, constants, refSpaceDer);
4468: for (fc = 0; fc < NcI; ++fc) {
4469: for (gc = 0; gc < NcJ; ++gc) {
4470: for (d = 0; d < dim; ++d) {
4471: g2[(fc*NcJ+gc)*dim+d] = 0.0;
4472: for (d2 = 0; d2 < dim; ++d2) g2[(fc*NcJ+gc)*dim+d] += invJ[d*dim+d2]*refSpaceDer[(fc*NcJ+gc)*dim+d2];
4473: g2[(fc*NcJ+gc)*dim+d] *= w;
4474: }
4475: }
4476: }
4477: }
4478: if (g3_func) {
4479: PetscInt d, d2, dp, d3;
4480: PetscMemzero(refSpaceDer, NcI*NcJ*dim*dim * sizeof(PetscScalar));
4481: g3_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, x, numConstants, constants, refSpaceDer);
4482: for (fc = 0; fc < NcI; ++fc) {
4483: for (gc = 0; gc < NcJ; ++gc) {
4484: for (d = 0; d < dim; ++d) {
4485: for (dp = 0; dp < dim; ++dp) {
4486: g3[((fc*NcJ+gc)*dim+d)*dim+dp] = 0.0;
4487: for (d2 = 0; d2 < dim; ++d2) {
4488: for (d3 = 0; d3 < dim; ++d3) {
4489: g3[((fc*NcJ+gc)*dim+d)*dim+dp] += invJ[d*dim+d2]*refSpaceDer[((fc*NcJ+gc)*dim+d2)*dim+d3]*invJ[dp*dim+d3];
4490: }
4491: }
4492: g3[((fc*NcJ+gc)*dim+d)*dim+dp] *= w;
4493: }
4494: }
4495: }
4496: }
4497: }
4499: for (f = 0; f < NbI; ++f) {
4500: for (fc = 0; fc < NcI; ++fc) {
4501: const PetscInt fidx = f*NcI+fc; /* Test function basis index */
4502: const PetscInt i = offsetI+f; /* Element matrix row */
4503: for (g = 0; g < NbJ; ++g) {
4504: for (gc = 0; gc < NcJ; ++gc) {
4505: const PetscInt gidx = g*NcJ+gc; /* Trial function basis index */
4506: const PetscInt j = offsetJ+g; /* Element matrix column */
4507: const PetscInt fOff = eOffset+i*totDim+j;
4508: PetscInt d, d2;
4510: elemMat[fOff] += BIq[fidx]*g0[fc*NcJ+gc]*BJq[gidx];
4511: for (d = 0; d < dim; ++d) {
4512: elemMat[fOff] += BIq[fidx]*g1[(fc*NcJ+gc)*dim+d]*DJq[gidx*dim+d];
4513: elemMat[fOff] += DIq[fidx*dim+d]*g2[(fc*NcJ+gc)*dim+d]*BJq[gidx];
4514: for (d2 = 0; d2 < dim; ++d2) {
4515: elemMat[fOff] += DIq[fidx*dim+d]*g3[((fc*NcJ+gc)*dim+d)*dim+d2]*DJq[gidx*dim+d2];
4516: }
4517: }
4518: }
4519: }
4520: }
4521: }
4522: }
4523: if (debug > 1) {
4524: PetscInt fc, f, gc, g;
4526: PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);
4527: for (fc = 0; fc < NcI; ++fc) {
4528: for (f = 0; f < NbI; ++f) {
4529: const PetscInt i = offsetI + f*NcI+fc;
4530: for (gc = 0; gc < NcJ; ++gc) {
4531: for (g = 0; g < NbJ; ++g) {
4532: const PetscInt j = offsetJ + g*NcJ+gc;
4533: PetscPrintf(PETSC_COMM_SELF, " elemMat[%d,%d,%d,%d]: %g\n", f, fc, g, gc, PetscRealPart(elemMat[eOffset+i*totDim+j]));
4534: }
4535: }
4536: PetscPrintf(PETSC_COMM_SELF, "\n");
4537: }
4538: }
4539: }
4540: cOffset += totDim;
4541: cOffsetAux += totDimAux;
4542: eOffset += PetscSqr(totDim);
4543: }
4544: return(0);
4545: }
4547: PetscErrorCode PetscFEIntegrateBdJacobian_Basic(PetscFE fem, PetscDS prob, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEFaceGeom *fgeom,
4548: const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
4549: {
4550: const PetscInt debug = 0;
4551: PetscBdPointJac g0_func;
4552: PetscBdPointJac g1_func;
4553: PetscBdPointJac g2_func;
4554: PetscBdPointJac g3_func;
4555: PetscInt cOffset = 0; /* Offset into coefficients[] for element e */
4556: PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
4557: PetscInt eOffset = 0; /* Offset into elemMat[] for element e */
4558: PetscInt offsetI = 0; /* Offset into an element vector for fieldI */
4559: PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */
4560: PetscQuadrature quad;
4561: PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *refSpaceDer, *refSpaceDerAux;
4562: const PetscScalar *constants;
4563: PetscReal *x;
4564: PetscReal **B, **D, **BAux = NULL, **DAux = NULL, *BI, *DI, *BJ, *DJ;
4565: PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL;
4566: PetscInt NbI = 0, NcI = 0, NbJ = 0, NcJ = 0;
4567: PetscInt dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, e;
4568: PetscErrorCode ierr;
4571: PetscFEGetSpatialDimension(fem, &dim);
4572: PetscFEGetFaceQuadrature(fem, &quad);
4573: PetscDSGetNumFields(prob, &Nf);
4574: PetscDSGetTotalDimension(prob, &totDim);
4575: PetscDSGetDimensions(prob, &Nb);
4576: PetscDSGetComponents(prob, &Nc);
4577: PetscDSGetComponentOffsets(prob, &uOff);
4578: PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);
4579: PetscDSGetFieldOffset(prob, fieldI, &offsetI);
4580: PetscDSGetFieldOffset(prob, fieldJ, &offsetJ);
4581: PetscDSGetBdJacobian(prob, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);
4582: PetscDSGetEvaluationArrays(prob, &u, coefficients_t ? &u_t : NULL, &u_x);
4583: PetscDSGetRefCoordArrays(prob, &x, &refSpaceDer);
4584: PetscDSGetWeakFormArrays(prob, NULL, NULL, &g0, &g1, &g2, &g3);
4585: PetscDSGetFaceTabulation(prob, &B, &D);
4586: PetscDSGetConstants(prob, &numConstants, &constants);
4587: if (probAux) {
4588: PetscDSGetNumFields(probAux, &NfAux);
4589: PetscDSGetTotalDimension(probAux, &totDimAux);
4590: PetscDSGetDimensions(probAux, &NbAux);
4591: PetscDSGetComponents(probAux, &NcAux);
4592: PetscDSGetComponentOffsets(probAux, &aOff);
4593: PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);
4594: PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);
4595: PetscDSGetRefCoordArrays(probAux, NULL, &refSpaceDerAux);
4596: PetscDSGetFaceTabulation(probAux, &BAux, &DAux);
4597: }
4598: NbI = Nb[fieldI], NbJ = Nb[fieldJ];
4599: NcI = Nc[fieldI], NcJ = Nc[fieldJ];
4600: BI = B[fieldI], BJ = B[fieldJ];
4601: DI = D[fieldI], DJ = D[fieldJ];
4602: /* Initialize here in case the function is not defined */
4603: PetscMemzero(g0, NcI*NcJ * sizeof(PetscScalar));
4604: PetscMemzero(g1, NcI*NcJ*dim * sizeof(PetscScalar));
4605: PetscMemzero(g2, NcI*NcJ*dim * sizeof(PetscScalar));
4606: PetscMemzero(g3, NcI*NcJ*dim*dim * sizeof(PetscScalar));
4607: for (e = 0; e < Ne; ++e) {
4608: const PetscReal *quadPoints, *quadWeights;
4609: const PetscReal *v0 = fgeom[e].v0;
4610: const PetscReal *J = fgeom[e].J;
4611: const PetscReal *invJ = fgeom[e].invJ[0];
4612: const PetscReal detJ = fgeom[e].detJ;
4613: const PetscReal *n = fgeom[e].n;
4614: const PetscInt face = fgeom[e].face[0];
4615: PetscInt qNc, Nq, q;
4617: PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);
4618: if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
4619: for (q = 0; q < Nq; ++q) {
4620: const PetscReal *BIq = &BI[(face*Nq+q)*NbI*NcI], *BJq = &BJ[(face*Nq+q)*NbJ*NcJ];
4621: const PetscReal *DIq = &DI[(face*Nq+q)*NbI*NcI*dim], *DJq = &DJ[(face*Nq+q)*NbJ*NcJ*dim];
4622: const PetscReal w = detJ*quadWeights[q];
4623: PetscInt f, g, fc, gc, c;
4625: if (debug) {PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);}
4626: CoordinatesRefToReal(dim, dim-1, v0, J, &quadPoints[q*(dim-1)], x);
4627: EvaluateFieldJets(dim, Nf, Nb, Nc, face*Nq+q, B, D, refSpaceDer, invJ, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);
4628: if (probAux) EvaluateFieldJets(dim, NfAux, NbAux, NcAux, face*Nq+q, BAux, DAux, refSpaceDerAux, invJ, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);
4629: if (g0_func) {
4630: PetscMemzero(g0, NcI*NcJ * sizeof(PetscScalar));
4631: g0_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, x, n, numConstants, constants, g0);
4632: for (c = 0; c < NcI*NcJ; ++c) g0[c] *= w;
4633: }
4634: if (g1_func) {
4635: PetscInt d, d2;
4636: PetscMemzero(refSpaceDer, NcI*NcJ*dim * sizeof(PetscScalar));
4637: g1_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, x, n, numConstants, constants, refSpaceDer);
4638: for (fc = 0; fc < NcI; ++fc) {
4639: for (gc = 0; gc < NcJ; ++gc) {
4640: for (d = 0; d < dim; ++d) {
4641: g1[(fc*NcJ+gc)*dim+d] = 0.0;
4642: for (d2 = 0; d2 < dim; ++d2) g1[(fc*NcJ+gc)*dim+d] += invJ[d*dim+d2]*refSpaceDer[(fc*NcJ+gc)*dim+d2];
4643: g1[(fc*NcJ+gc)*dim+d] *= w;
4644: }
4645: }
4646: }
4647: }
4648: if (g2_func) {
4649: PetscInt d, d2;
4650: PetscMemzero(refSpaceDer, NcI*NcJ*dim * sizeof(PetscScalar));
4651: g2_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, x, n, numConstants, constants, refSpaceDer);
4652: for (fc = 0; fc < NcI; ++fc) {
4653: for (gc = 0; gc < NcJ; ++gc) {
4654: for (d = 0; d < dim; ++d) {
4655: g2[(fc*NcJ+gc)*dim+d] = 0.0;
4656: for (d2 = 0; d2 < dim; ++d2) g2[(fc*NcJ+gc)*dim+d] += invJ[d*dim+d2]*refSpaceDer[(fc*NcJ+gc)*dim+d2];
4657: g2[(fc*NcJ+gc)*dim+d] *= w;
4658: }
4659: }
4660: }
4661: }
4662: if (g3_func) {
4663: PetscInt d, d2, dp, d3;
4664: PetscMemzero(refSpaceDer, NcI*NcJ*dim*dim * sizeof(PetscScalar));
4665: g3_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, x, n, numConstants, constants, refSpaceDer);
4666: for (fc = 0; fc < NcI; ++fc) {
4667: for (gc = 0; gc < NcJ; ++gc) {
4668: for (d = 0; d < dim; ++d) {
4669: for (dp = 0; dp < dim; ++dp) {
4670: g3[((fc*NcJ+gc)*dim+d)*dim+dp] = 0.0;
4671: for (d2 = 0; d2 < dim; ++d2) {
4672: for (d3 = 0; d3 < dim; ++d3) {
4673: g3[((fc*NcJ+gc)*dim+d)*dim+dp] += invJ[d*dim+d2]*refSpaceDer[((fc*NcJ+gc)*dim+d2)*dim+d3]*invJ[dp*dim+d3];
4674: }
4675: }
4676: g3[((fc*NcJ+gc)*dim+d)*dim+dp] *= w;
4677: }
4678: }
4679: }
4680: }
4681: }
4683: for (f = 0; f < NbI; ++f) {
4684: for (fc = 0; fc < NcI; ++fc) {
4685: const PetscInt fidx = f*NcI+fc; /* Test function basis index */
4686: const PetscInt i = offsetI+f; /* Element matrix row */
4687: for (g = 0; g < NbJ; ++g) {
4688: for (gc = 0; gc < NcJ; ++gc) {
4689: const PetscInt gidx = g*NcJ+gc; /* Trial function basis index */
4690: const PetscInt j = offsetJ+g; /* Element matrix column */
4691: const PetscInt fOff = eOffset+i*totDim+j;
4692: PetscInt d, d2;
4694: elemMat[fOff] += BIq[fidx]*g0[fc*NcJ+gc]*BJq[gidx];
4695: for (d = 0; d < dim; ++d) {
4696: elemMat[fOff] += BIq[fidx]*g1[(fc*NcJ+gc)*dim+d]*DJq[gidx*dim+d];
4697: elemMat[fOff] += DIq[fidx*dim+d]*g2[(fc*NcJ+gc)*dim+d]*BJq[gidx];
4698: for (d2 = 0; d2 < dim; ++d2) {
4699: elemMat[fOff] += DIq[fidx*dim+d]*g3[((fc*NcJ+gc)*dim+d)*dim+d2]*DJq[gidx*dim+d2];
4700: }
4701: }
4702: }
4703: }
4704: }
4705: }
4706: }
4707: if (debug > 1) {
4708: PetscInt fc, f, gc, g;
4710: PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);
4711: for (fc = 0; fc < NcI; ++fc) {
4712: for (f = 0; f < NbI; ++f) {
4713: const PetscInt i = offsetI + f*NcI+fc;
4714: for (gc = 0; gc < NcJ; ++gc) {
4715: for (g = 0; g < NbJ; ++g) {
4716: const PetscInt j = offsetJ + g*NcJ+gc;
4717: PetscPrintf(PETSC_COMM_SELF, " elemMat[%d,%d,%d,%d]: %g\n", f, fc, g, gc, PetscRealPart(elemMat[eOffset+i*totDim+j]));
4718: }
4719: }
4720: PetscPrintf(PETSC_COMM_SELF, "\n");
4721: }
4722: }
4723: }
4724: cOffset += totDim;
4725: cOffsetAux += totDimAux;
4726: eOffset += PetscSqr(totDim);
4727: }
4728: return(0);
4729: }
4731: PetscErrorCode PetscFEInitialize_Basic(PetscFE fem)
4732: {
4734: fem->ops->setfromoptions = NULL;
4735: fem->ops->setup = PetscFESetUp_Basic;
4736: fem->ops->view = PetscFEView_Basic;
4737: fem->ops->destroy = PetscFEDestroy_Basic;
4738: fem->ops->getdimension = PetscFEGetDimension_Basic;
4739: fem->ops->gettabulation = PetscFEGetTabulation_Basic;
4740: fem->ops->integrate = PetscFEIntegrate_Basic;
4741: fem->ops->integrateresidual = PetscFEIntegrateResidual_Basic;
4742: fem->ops->integratebdresidual = PetscFEIntegrateBdResidual_Basic;
4743: fem->ops->integratejacobianaction = NULL/* PetscFEIntegrateJacobianAction_Basic */;
4744: fem->ops->integratejacobian = PetscFEIntegrateJacobian_Basic;
4745: fem->ops->integratebdjacobian = PetscFEIntegrateBdJacobian_Basic;
4746: return(0);
4747: }
4749: /*MC
4750: PETSCFEBASIC = "basic" - A PetscFE object that integrates with basic tiling and no vectorization
4752: Level: intermediate
4754: .seealso: PetscFEType, PetscFECreate(), PetscFESetType()
4755: M*/
4757: PETSC_EXTERN PetscErrorCode PetscFECreate_Basic(PetscFE fem)
4758: {
4759: PetscFE_Basic *b;
4764: PetscNewLog(fem,&b);
4765: fem->data = b;
4767: PetscFEInitialize_Basic(fem);
4768: return(0);
4769: }
4771: PetscErrorCode PetscFEDestroy_Nonaffine(PetscFE fem)
4772: {
4773: PetscFE_Nonaffine *na = (PetscFE_Nonaffine *) fem->data;
4777: PetscFree(na);
4778: return(0);
4779: }
4781: PetscErrorCode PetscFEIntegrateResidual_Nonaffine(PetscFE fem, PetscDS prob, PetscInt field, PetscInt Ne, PetscFECellGeom *cgeom,
4782: const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
4783: {
4784: const PetscInt debug = 0;
4785: PetscPointFunc f0_func;
4786: PetscPointFunc f1_func;
4787: PetscQuadrature quad;
4788: PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *refSpaceDer, *refSpaceDerAux;
4789: const PetscScalar *constants;
4790: PetscReal *x;
4791: PetscReal **B, **D, **BAux = NULL, **DAux = NULL, *BI, *DI;
4792: PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL;
4793: PetscInt dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NbI, NcI;
4794: PetscErrorCode ierr;
4797: PetscFEGetSpatialDimension(fem, &dim);
4798: PetscFEGetQuadrature(fem, &quad);
4799: PetscDSGetNumFields(prob, &Nf);
4800: PetscDSGetTotalDimension(prob, &totDim);
4801: PetscDSGetDimensions(prob, &Nb);
4802: PetscDSGetComponents(prob, &Nc);
4803: PetscDSGetComponentOffsets(prob, &uOff);
4804: PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);
4805: PetscDSGetFieldOffset(prob, field, &fOffset);
4806: PetscDSGetResidual(prob, field, &f0_func, &f1_func);
4807: PetscDSGetEvaluationArrays(prob, &u, coefficients_t ? &u_t : NULL, &u_x);
4808: PetscDSGetRefCoordArrays(prob, &x, &refSpaceDer);
4809: PetscDSGetWeakFormArrays(prob, &f0, &f1, NULL, NULL, NULL, NULL);
4810: PetscDSGetTabulation(prob, &B, &D);
4811: PetscDSGetConstants(prob, &numConstants, &constants);
4812: if (probAux) {
4813: PetscDSGetNumFields(probAux, &NfAux);
4814: PetscDSGetTotalDimension(probAux, &totDimAux);
4815: PetscDSGetDimensions(probAux, &NbAux);
4816: PetscDSGetComponents(probAux, &NcAux);
4817: PetscDSGetComponentOffsets(probAux, &aOff);
4818: PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);
4819: PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);
4820: PetscDSGetRefCoordArrays(probAux, NULL, &refSpaceDerAux);
4821: PetscDSGetTabulation(probAux, &BAux, &DAux);
4822: }
4823: NbI = Nb[field];
4824: NcI = Nc[field];
4825: BI = B[field];
4826: DI = D[field];
4827: for (e = 0; e < Ne; ++e) {
4828: const PetscReal *quadPoints, *quadWeights;
4829: PetscInt qNc, Nq, q;
4831: PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);
4832: if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
4833: PetscMemzero(f0, Nq*Nc[field]* sizeof(PetscScalar));
4834: PetscMemzero(f1, Nq*Nc[field]*dim * sizeof(PetscScalar));
4835: for (q = 0; q < Nq; ++q) {
4836: const PetscReal *v0 = cgeom[e*Nq+q].v0;
4837: const PetscReal *J = cgeom[e*Nq+q].J;
4838: const PetscReal *invJ = cgeom[e*Nq+q].invJ;
4839: const PetscReal detJ = cgeom[e*Nq+q].detJ;
4841: if (debug) {PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);}
4842: CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], x);
4843: EvaluateFieldJets(dim, Nf, Nb, Nc, q, B, D, refSpaceDer, invJ, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);
4844: if (probAux) EvaluateFieldJets(dim, NfAux, NbAux, NcAux, q, BAux, DAux, refSpaceDerAux, invJ, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);
4845: if (f0_func) f0_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, x, numConstants, constants, &f0[q*NcI]);
4846: if (f1_func) {
4847: PetscMemzero(refSpaceDer, NcI*dim * sizeof(PetscScalar));
4848: f1_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, x, numConstants, constants, refSpaceDer);
4849: }
4850: TransformF(dim, NcI, q, invJ, detJ, quadWeights, refSpaceDer, f0, f1);
4851: }
4852: UpdateElementVec(dim, Nq, NbI, NcI, BI, DI, f0, f1, &elemVec[cOffset+fOffset]);
4853: cOffset += totDim;
4854: cOffsetAux += totDimAux;
4855: }
4856: return(0);
4857: }
4859: PetscErrorCode PetscFEIntegrateBdResidual_Nonaffine(PetscFE fem, PetscDS prob, PetscInt field, PetscInt Ne, PetscFEFaceGeom *fgeom,
4860: const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
4861: {
4862: const PetscInt debug = 0;
4863: PetscBdPointFunc f0_func;
4864: PetscBdPointFunc f1_func;
4865: PetscQuadrature quad;
4866: PetscScalar *f0, *f1, *u, *u_t = NULL, *u_x, *a, *a_x, *refSpaceDer, *refSpaceDerAux;
4867: const PetscScalar *constants;
4868: PetscReal *x;
4869: PetscReal **B, **D, **BAux = NULL, **DAux = NULL, *BI, *DI;
4870: PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL;
4871: PetscInt dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, cOffset = 0, cOffsetAux = 0, fOffset, e, NbI, NcI;
4872: PetscErrorCode ierr;
4875: PetscFEGetSpatialDimension(fem, &dim);
4876: PetscFEGetFaceQuadrature(fem, &quad);
4877: PetscDSGetNumFields(prob, &Nf);
4878: PetscDSGetTotalDimension(prob, &totDim);
4879: PetscDSGetDimensions(prob, &Nb);
4880: PetscDSGetComponents(prob, &Nc);
4881: PetscDSGetComponentOffsets(prob, &uOff);
4882: PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);
4883: PetscDSGetFieldOffset(prob, field, &fOffset);
4884: PetscDSGetBdResidual(prob, field, &f0_func, &f1_func);
4885: if (!f0_func && !f1_func) return(0);
4886: PetscDSGetEvaluationArrays(prob, &u, coefficients_t ? &u_t : NULL, &u_x);
4887: PetscDSGetRefCoordArrays(prob, &x, &refSpaceDer);
4888: PetscDSGetWeakFormArrays(prob, &f0, &f1, NULL, NULL, NULL, NULL);
4889: PetscDSGetFaceTabulation(prob, &B, &D);
4890: PetscDSGetConstants(prob, &numConstants, &constants);
4891: if (probAux) {
4892: PetscDSGetNumFields(probAux, &NfAux);
4893: PetscDSGetTotalDimension(probAux, &totDimAux);
4894: PetscDSGetDimensions(probAux, &NbAux);
4895: PetscDSGetComponents(probAux, &NcAux);
4896: PetscDSGetComponentOffsets(probAux, &aOff);
4897: PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);
4898: PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);
4899: PetscDSGetRefCoordArrays(probAux, NULL, &refSpaceDerAux);
4900: PetscDSGetFaceTabulation(probAux, &BAux, &DAux);
4901: }
4902: NbI = Nb[field];
4903: NcI = Nc[field];
4904: BI = B[field];
4905: DI = D[field];
4906: for (e = 0; e < Ne; ++e) {
4907: const PetscReal *quadPoints, *quadWeights;
4908: PetscInt qNc, Nq, q, face;
4910: PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);
4911: if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
4912: face = fgeom[e*Nq].face[0];
4913: PetscMemzero(f0, Nq*NcI* sizeof(PetscScalar));
4914: PetscMemzero(f1, Nq*NcI*dim * sizeof(PetscScalar));
4915: for (q = 0; q < Nq; ++q) {
4916: const PetscReal *v0 = fgeom[e*Nq+q].v0;
4917: const PetscReal *J = fgeom[e*Nq+q].J;
4918: const PetscReal *invJ = fgeom[e*Nq+q].invJ[0];
4919: const PetscReal detJ = fgeom[e*Nq+q].detJ;
4920: const PetscReal *n = fgeom[e*Nq+q].n;
4922: if (debug) {PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);}
4923: CoordinatesRefToReal(dim, dim-1, v0, J, &quadPoints[q*(dim-1)], x);
4924: EvaluateFieldJets(dim, Nf, Nb, Nc, face*Nq+q, B, D, refSpaceDer, invJ, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);
4925: if (probAux) EvaluateFieldJets(dim, NfAux, NbAux, NcAux, face*Nq+q, BAux, DAux, refSpaceDerAux, invJ, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);
4926: if (f0_func) f0_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, x, n, numConstants, constants, &f0[q*NcI]);
4927: if (f1_func) {
4928: PetscMemzero(refSpaceDer, NcI*dim * sizeof(PetscScalar));
4929: f1_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, x, n, numConstants, constants, refSpaceDer);
4930: }
4931: TransformF(dim, NcI, q, invJ, detJ, quadWeights, refSpaceDer, f0_func ? f0 : NULL, f1_func ? f1 : NULL);
4932: }
4933: UpdateElementVec(dim, Nq, NbI, NcI, &BI[face*Nq*NbI*NcI], &DI[face*Nq*NbI*NcI*dim], f0, f1, &elemVec[cOffset+fOffset]);
4934: cOffset += totDim;
4935: cOffsetAux += totDimAux;
4936: }
4937: return(0);
4938: }
4940: PetscErrorCode PetscFEIntegrateJacobian_Nonaffine(PetscFE fem, PetscDS prob, PetscFEJacobianType jtype, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFECellGeom *cgeom,
4941: const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
4942: {
4943: const PetscInt debug = 0;
4944: PetscPointJac g0_func;
4945: PetscPointJac g1_func;
4946: PetscPointJac g2_func;
4947: PetscPointJac g3_func;
4948: PetscInt cOffset = 0; /* Offset into coefficients[] for element e */
4949: PetscInt cOffsetAux = 0; /* Offset into coefficientsAux[] for element e */
4950: PetscInt eOffset = 0; /* Offset into elemMat[] for element e */
4951: PetscInt offsetI = 0; /* Offset into an element vector for fieldI */
4952: PetscInt offsetJ = 0; /* Offset into an element vector for fieldJ */
4953: PetscQuadrature quad;
4954: PetscScalar *g0, *g1, *g2, *g3, *u, *u_t = NULL, *u_x, *a, *a_x, *refSpaceDer, *refSpaceDerAux;
4955: const PetscScalar *constants;
4956: PetscReal *x;
4957: PetscReal **B, **D, **BAux = NULL, **DAux = NULL, *BI, *DI, *BJ, *DJ;
4958: PetscInt NbI = 0, NcI = 0, NbJ = 0, NcJ = 0;
4959: PetscInt *uOff, *uOff_x, *aOff = NULL, *aOff_x = NULL, *Nb, *Nc, *NbAux = NULL, *NcAux = NULL;
4960: PetscInt dim, numConstants, Nf, NfAux = 0, totDim, totDimAux = 0, e;
4961: PetscErrorCode ierr;
4964: PetscFEGetSpatialDimension(fem, &dim);
4965: PetscFEGetQuadrature(fem, &quad);
4966: PetscDSGetNumFields(prob, &Nf);
4967: PetscDSGetTotalDimension(prob, &totDim);
4968: PetscDSGetDimensions(prob, &Nb);
4969: PetscDSGetComponents(prob, &Nc);
4970: PetscDSGetComponentOffsets(prob, &uOff);
4971: PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);
4972: switch(jtype) {
4973: case PETSCFE_JACOBIAN_DYN: PetscDSGetDynamicJacobian(prob, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);break;
4974: case PETSCFE_JACOBIAN_PRE: PetscDSGetJacobianPreconditioner(prob, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);break;
4975: case PETSCFE_JACOBIAN: PetscDSGetJacobian(prob, fieldI, fieldJ, &g0_func, &g1_func, &g2_func, &g3_func);break;
4976: }
4977: PetscDSGetEvaluationArrays(prob, &u, coefficients_t ? &u_t : NULL, &u_x);
4978: PetscDSGetRefCoordArrays(prob, &x, &refSpaceDer);
4979: PetscDSGetWeakFormArrays(prob, NULL, NULL, &g0, &g1, &g2, &g3);
4980: PetscDSGetTabulation(prob, &B, &D);
4981: PetscDSGetFieldOffset(prob, fieldI, &offsetI);
4982: PetscDSGetFieldOffset(prob, fieldJ, &offsetJ);
4983: PetscDSGetConstants(prob, &numConstants, &constants);
4984: if (probAux) {
4985: PetscDSGetNumFields(probAux, &NfAux);
4986: PetscDSGetTotalDimension(probAux, &totDimAux);
4987: PetscDSGetDimensions(probAux, &NbAux);
4988: PetscDSGetComponents(probAux, &NcAux);
4989: PetscDSGetComponentOffsets(probAux, &aOff);
4990: PetscDSGetComponentDerivativeOffsets(probAux, &aOff_x);
4991: PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);
4992: PetscDSGetRefCoordArrays(probAux, NULL, &refSpaceDerAux);
4993: PetscDSGetTabulation(probAux, &BAux, &DAux);
4994: }
4995: NbI = Nb[fieldI], NbJ = Nb[fieldJ];
4996: NcI = Nc[fieldI], NcJ = Nc[fieldJ];
4997: BI = B[fieldI], BJ = B[fieldJ];
4998: DI = D[fieldI], DJ = D[fieldJ];
4999: /* Initialize here in case the function is not defined */
5000: PetscMemzero(g0, NcI*NcJ * sizeof(PetscScalar));
5001: PetscMemzero(g1, NcI*NcJ*dim * sizeof(PetscScalar));
5002: PetscMemzero(g2, NcI*NcJ*dim * sizeof(PetscScalar));
5003: PetscMemzero(g3, NcI*NcJ*dim*dim * sizeof(PetscScalar));
5004: for (e = 0; e < Ne; ++e) {
5005: const PetscReal *quadPoints, *quadWeights;
5006: PetscInt qNc, Nq, q;
5008: PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);
5009: if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
5010: for (q = 0; q < Nq; ++q) {
5011: const PetscReal *v0 = cgeom[e*Nq+q].v0;
5012: const PetscReal *J = cgeom[e*Nq+q].J;
5013: const PetscReal *invJ = cgeom[e*Nq+q].invJ;
5014: const PetscReal detJ = cgeom[e*Nq+q].detJ;
5015: const PetscReal *BIq = &BI[q*NbI*NcI], *BJq = &BJ[q*NbJ*NcJ];
5016: const PetscReal *DIq = &DI[q*NbI*NcI*dim], *DJq = &DJ[q*NbJ*NcJ*dim];
5017: const PetscReal w = detJ*quadWeights[q];
5018: PetscInt f, g, fc, gc, c;
5020: if (debug) {PetscPrintf(PETSC_COMM_SELF, " quad point %d\n", q);}
5021: CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], x);
5022: EvaluateFieldJets(dim, Nf, Nb, Nc, q, B, D, refSpaceDer, invJ, &coefficients[cOffset], &coefficients_t[cOffset], u, u_x, u_t);
5023: if (probAux) EvaluateFieldJets(dim, NfAux, NbAux, NcAux, q, BAux, DAux, refSpaceDerAux, invJ, &coefficientsAux[cOffsetAux], NULL, a, a_x, NULL);
5024: if (g0_func) {
5025: PetscMemzero(g0, NcI*NcJ * sizeof(PetscScalar));
5026: g0_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, x, numConstants, constants, g0);
5027: for (c = 0; c < NcI*NcJ; ++c) g0[c] *= w;
5028: }
5029: if (g1_func) {
5030: PetscInt d, d2;
5031: PetscMemzero(refSpaceDer, NcI*NcJ*dim * sizeof(PetscScalar));
5032: g1_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, x, numConstants, constants, refSpaceDer);
5033: for (fc = 0; fc < NcI; ++fc) {
5034: for (gc = 0; gc < NcJ; ++gc) {
5035: for (d = 0; d < dim; ++d) {
5036: g1[(fc*NcJ+gc)*dim+d] = 0.0;
5037: for (d2 = 0; d2 < dim; ++d2) g1[(fc*NcJ+gc)*dim+d] += invJ[d*dim+d2]*refSpaceDer[(fc*NcJ+gc)*dim+d2];
5038: g1[(fc*NcJ+gc)*dim+d] *= w;
5039: }
5040: }
5041: }
5042: }
5043: if (g2_func) {
5044: PetscInt d, d2;
5045: PetscMemzero(refSpaceDer, NcI*NcJ*dim * sizeof(PetscScalar));
5046: g2_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, x, numConstants, constants, refSpaceDer);
5047: for (fc = 0; fc < NcI; ++fc) {
5048: for (gc = 0; gc < NcJ; ++gc) {
5049: for (d = 0; d < dim; ++d) {
5050: g2[(fc*NcJ+gc)*dim+d] = 0.0;
5051: for (d2 = 0; d2 < dim; ++d2) g2[(fc*NcJ+gc)*dim+d] += invJ[d*dim+d2]*refSpaceDer[(fc*NcJ+gc)*dim+d2];
5052: g2[(fc*NcJ+gc)*dim+d] *= w;
5053: }
5054: }
5055: }
5056: }
5057: if (g3_func) {
5058: PetscInt d, d2, dp, d3;
5059: PetscMemzero(refSpaceDer, NcI*NcJ*dim*dim * sizeof(PetscScalar));
5060: g3_func(dim, Nf, NfAux, uOff, uOff_x, u, u_t, u_x, aOff, aOff_x, a, NULL, a_x, t, u_tshift, x, numConstants, constants, refSpaceDer);
5061: for (fc = 0; fc < NcI; ++fc) {
5062: for (gc = 0; gc < NcJ; ++gc) {
5063: for (d = 0; d < dim; ++d) {
5064: for (dp = 0; dp < dim; ++dp) {
5065: g3[((fc*NcJ+gc)*dim+d)*dim+dp] = 0.0;
5066: for (d2 = 0; d2 < dim; ++d2) {
5067: for (d3 = 0; d3 < dim; ++d3) {
5068: g3[((fc*NcJ+gc)*dim+d)*dim+dp] += invJ[d*dim+d2]*refSpaceDer[((fc*NcJ+gc)*dim+d2)*dim+d3]*invJ[dp*dim+d3];
5069: }
5070: }
5071: g3[((fc*NcJ+gc)*dim+d)*dim+dp] *= w;
5072: }
5073: }
5074: }
5075: }
5076: }
5078: for (f = 0; f < NbI; ++f) {
5079: for (fc = 0; fc < NcI; ++fc) {
5080: const PetscInt fidx = f*NcI+fc; /* Test function basis index */
5081: const PetscInt i = offsetI+f; /* Element matrix row */
5082: for (g = 0; g < NbJ; ++g) {
5083: for (gc = 0; gc < NcJ; ++gc) {
5084: const PetscInt gidx = g*NcJ+gc; /* Trial function basis index */
5085: const PetscInt j = offsetJ+g; /* Element matrix column */
5086: const PetscInt fOff = eOffset+i*totDim+j;
5087: PetscInt d, d2;
5089: elemMat[fOff] += BIq[fidx]*g0[fc*NcJ+gc]*BJq[gidx];
5090: for (d = 0; d < dim; ++d) {
5091: elemMat[fOff] += BIq[fidx]*g1[(fc*NcJ+gc)*dim+d]*DJq[gidx*dim+d];
5092: elemMat[fOff] += DIq[fidx*dim+d]*g2[(fc*NcJ+gc)*dim+d]*BJq[gidx];
5093: for (d2 = 0; d2 < dim; ++d2) {
5094: elemMat[fOff] += DIq[fidx*dim+d]*g3[((fc*NcJ+gc)*dim+d)*dim+d2]*DJq[gidx*dim+d2];
5095: }
5096: }
5097: }
5098: }
5099: }
5100: }
5101: }
5102: if (debug > 1) {
5103: PetscInt fc, f, gc, g;
5105: PetscPrintf(PETSC_COMM_SELF, "Element matrix for fields %d and %d\n", fieldI, fieldJ);
5106: for (fc = 0; fc < NcI; ++fc) {
5107: for (f = 0; f < NbI; ++f) {
5108: const PetscInt i = offsetI + f*NcI+fc;
5109: for (gc = 0; gc < NcJ; ++gc) {
5110: for (g = 0; g < NbJ; ++g) {
5111: const PetscInt j = offsetJ + g*NcJ+gc;
5112: PetscPrintf(PETSC_COMM_SELF, " elemMat[%d,%d,%d,%d]: %g\n", f, fc, g, gc, PetscRealPart(elemMat[eOffset+i*totDim+j]));
5113: }
5114: }
5115: PetscPrintf(PETSC_COMM_SELF, "\n");
5116: }
5117: }
5118: }
5119: cOffset += totDim;
5120: cOffsetAux += totDimAux;
5121: eOffset += PetscSqr(totDim);
5122: }
5123: return(0);
5124: }
5126: PetscErrorCode PetscFEInitialize_Nonaffine(PetscFE fem)
5127: {
5129: fem->ops->setfromoptions = NULL;
5130: fem->ops->setup = PetscFESetUp_Basic;
5131: fem->ops->view = NULL;
5132: fem->ops->destroy = PetscFEDestroy_Nonaffine;
5133: fem->ops->getdimension = PetscFEGetDimension_Basic;
5134: fem->ops->gettabulation = PetscFEGetTabulation_Basic;
5135: fem->ops->integrateresidual = PetscFEIntegrateResidual_Nonaffine;
5136: fem->ops->integratebdresidual = PetscFEIntegrateBdResidual_Nonaffine;
5137: fem->ops->integratejacobianaction = NULL/* PetscFEIntegrateJacobianAction_Nonaffine */;
5138: fem->ops->integratejacobian = PetscFEIntegrateJacobian_Nonaffine;
5139: return(0);
5140: }
5142: /*MC
5143: PETSCFENONAFFINE = "nonaffine" - A PetscFE object that integrates with basic tiling and no vectorization for non-affine mappings
5145: Level: intermediate
5147: .seealso: PetscFEType, PetscFECreate(), PetscFESetType()
5148: M*/
5150: PETSC_EXTERN PetscErrorCode PetscFECreate_Nonaffine(PetscFE fem)
5151: {
5152: PetscFE_Nonaffine *na;
5153: PetscErrorCode ierr;
5157: PetscNewLog(fem, &na);
5158: fem->data = na;
5160: PetscFEInitialize_Nonaffine(fem);
5161: return(0);
5162: }
5164: #ifdef PETSC_HAVE_OPENCL
5166: PetscErrorCode PetscFEDestroy_OpenCL(PetscFE fem)
5167: {
5168: PetscFE_OpenCL *ocl = (PetscFE_OpenCL *) fem->data;
5169: PetscErrorCode ierr;
5172: clReleaseCommandQueue(ocl->queue_id);
5173: ocl->queue_id = 0;
5174: clReleaseContext(ocl->ctx_id);
5175: ocl->ctx_id = 0;
5176: PetscFree(ocl);
5177: return(0);
5178: }
5180: #define STRING_ERROR_CHECK(MSG) do { string_tail += count; if (string_tail == end_of_buffer) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, MSG);} while(0)
5181: enum {LAPLACIAN = 0, ELASTICITY = 1};
5183: /* NOTE: This is now broken for vector problems. Must redo loops to respect vector basis elements */
5184: /* dim Number of spatial dimensions: 2 */
5185: /* N_b Number of basis functions: generated */
5186: /* N_{bt} Number of total basis functions: N_b * N_{comp} */
5187: /* N_q Number of quadrature points: generated */
5188: /* N_{bs} Number of block cells LCM(N_b, N_q) */
5189: /* N_{bst} Number of block cell components LCM(N_{bt}, N_q) */
5190: /* N_{bl} Number of concurrent blocks generated */
5191: /* N_t Number of threads: N_{bl} * N_{bs} */
5192: /* N_{cbc} Number of concurrent basis cells: N_{bl} * N_q */
5193: /* N_{cqc} Number of concurrent quadrature cells: N_{bl} * N_b */
5194: /* N_{sbc} Number of serial basis cells: N_{bs} / N_q */
5195: /* N_{sqc} Number of serial quadrature cells: N_{bs} / N_b */
5196: /* N_{cb} Number of serial cell batches: input */
5197: /* N_c Number of total cells: N_{cb}*N_{t}/N_{comp} */
5198: PetscErrorCode PetscFEOpenCLGenerateIntegrationCode(PetscFE fem, char **string_buffer, PetscInt buffer_length, PetscBool useAux, PetscInt N_bl)
5199: {
5200: PetscFE_OpenCL *ocl = (PetscFE_OpenCL *) fem->data;
5201: PetscQuadrature q;
5202: char *string_tail = *string_buffer;
5203: char *end_of_buffer = *string_buffer + buffer_length;
5204: char float_str[] = "float", double_str[] = "double";
5205: char *numeric_str = &(float_str[0]);
5206: PetscInt op = ocl->op;
5207: PetscBool useField = PETSC_FALSE;
5208: PetscBool useFieldDer = PETSC_TRUE;
5209: PetscBool useFieldAux = useAux;
5210: PetscBool useFieldDerAux= PETSC_FALSE;
5211: PetscBool useF0 = PETSC_TRUE;
5212: PetscBool useF1 = PETSC_TRUE;
5213: const PetscReal *points, *weights;
5214: PetscReal *basis, *basisDer;
5215: PetscInt dim, qNc, N_b, N_c, N_q, N_t, p, d, b, c;
5216: size_t count;
5217: PetscErrorCode ierr;
5220: PetscFEGetSpatialDimension(fem, &dim);
5221: PetscFEGetDimension(fem, &N_b);
5222: PetscFEGetNumComponents(fem, &N_c);
5223: PetscFEGetQuadrature(fem, &q);
5224: PetscQuadratureGetData(q, NULL, &qNc, &N_q, &points, &weights);
5225: if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
5226: N_t = N_b * N_c * N_q * N_bl;
5227: /* Enable device extension for double precision */
5228: if (ocl->realType == PETSC_DOUBLE) {
5229: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
5230: "#if defined(cl_khr_fp64)\n"
5231: "# pragma OPENCL EXTENSION cl_khr_fp64: enable\n"
5232: "#elif defined(cl_amd_fp64)\n"
5233: "# pragma OPENCL EXTENSION cl_amd_fp64: enable\n"
5234: "#endif\n",
5235: &count);STRING_ERROR_CHECK("Message to short");
5236: numeric_str = &(double_str[0]);
5237: }
5238: /* Kernel API */
5239: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
5240: "\n"
5241: "__kernel void integrateElementQuadrature(int N_cb, __global %s *coefficients, __global %s *coefficientsAux, __global %s *jacobianInverses, __global %s *jacobianDeterminants, __global %s *elemVec)\n"
5242: "{\n",
5243: &count, numeric_str, numeric_str, numeric_str, numeric_str, numeric_str);STRING_ERROR_CHECK("Message to short");
5244: /* Quadrature */
5245: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
5246: " /* Quadrature points\n"
5247: " - (x1,y1,x2,y2,...) */\n"
5248: " const %s points[%d] = {\n",
5249: &count, numeric_str, N_q*dim);STRING_ERROR_CHECK("Message to short");
5250: for (p = 0; p < N_q; ++p) {
5251: for (d = 0; d < dim; ++d) {
5252: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "%g,\n", &count, points[p*dim+d]);STRING_ERROR_CHECK("Message to short");
5253: }
5254: }
5255: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "};\n", &count);STRING_ERROR_CHECK("Message to short");
5256: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
5257: " /* Quadrature weights\n"
5258: " - (v1,v2,...) */\n"
5259: " const %s weights[%d] = {\n",
5260: &count, numeric_str, N_q);STRING_ERROR_CHECK("Message to short");
5261: for (p = 0; p < N_q; ++p) {
5262: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "%g,\n", &count, weights[p]);STRING_ERROR_CHECK("Message to short");
5263: }
5264: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "};\n", &count);STRING_ERROR_CHECK("Message to short");
5265: /* Basis Functions */
5266: PetscFEGetDefaultTabulation(fem, &basis, &basisDer, NULL);
5267: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
5268: " /* Nodal basis function evaluations\n"
5269: " - basis component is fastest varying, the basis function, then point */\n"
5270: " const %s Basis[%d] = {\n",
5271: &count, numeric_str, N_q*N_b*N_c);STRING_ERROR_CHECK("Message to short");
5272: for (p = 0; p < N_q; ++p) {
5273: for (b = 0; b < N_b; ++b) {
5274: for (c = 0; c < N_c; ++c) {
5275: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "%g,\n", &count, basis[(p*N_b + b)*N_c + c]);STRING_ERROR_CHECK("Message to short");
5276: }
5277: }
5278: }
5279: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "};\n", &count);STRING_ERROR_CHECK("Message to short");
5280: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
5281: "\n"
5282: " /* Nodal basis function derivative evaluations,\n"
5283: " - derivative direction is fastest varying, then basis component, then basis function, then point */\n"
5284: " const %s%d BasisDerivatives[%d] = {\n",
5285: &count, numeric_str, dim, N_q*N_b*N_c);STRING_ERROR_CHECK("Message to short");
5286: for (p = 0; p < N_q; ++p) {
5287: for (b = 0; b < N_b; ++b) {
5288: for (c = 0; c < N_c; ++c) {
5289: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "(%s%d)(", &count, numeric_str, dim);STRING_ERROR_CHECK("Message to short");
5290: for (d = 0; d < dim; ++d) {
5291: if (d > 0) {
5292: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, ", %g", &count, basisDer[((p*N_b + b)*dim + d)*N_c + c]);STRING_ERROR_CHECK("Message to short");
5293: } else {
5294: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "%g", &count, basisDer[((p*N_b + b)*dim + d)*N_c + c]);STRING_ERROR_CHECK("Message to short");
5295: }
5296: }
5297: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "),\n", &count);STRING_ERROR_CHECK("Message to short");
5298: }
5299: }
5300: }
5301: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, "};\n", &count);STRING_ERROR_CHECK("Message to short");
5302: /* Sizes */
5303: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
5304: " const int dim = %d; // The spatial dimension\n"
5305: " const int N_bl = %d; // The number of concurrent blocks\n"
5306: " const int N_b = %d; // The number of basis functions\n"
5307: " const int N_comp = %d; // The number of basis function components\n"
5308: " const int N_bt = N_b*N_comp; // The total number of scalar basis functions\n"
5309: " const int N_q = %d; // The number of quadrature points\n"
5310: " const int N_bst = N_bt*N_q; // The block size, LCM(N_b*N_comp, N_q), Notice that a block is not processed simultaneously\n"
5311: " const int N_t = N_bst*N_bl; // The number of threads, N_bst * N_bl\n"
5312: " const int N_bc = N_t/N_comp; // The number of cells per batch (N_b*N_q*N_bl)\n"
5313: " const int N_sbc = N_bst / (N_q * N_comp);\n"
5314: " const int N_sqc = N_bst / N_bt;\n"
5315: " /*const int N_c = N_cb * N_bc;*/\n"
5316: "\n"
5317: " /* Calculated indices */\n"
5318: " /*const int tidx = get_local_id(0) + get_local_size(0)*get_local_id(1);*/\n"
5319: " const int tidx = get_local_id(0);\n"
5320: " const int blidx = tidx / N_bst; // Block number for this thread\n"
5321: " const int bidx = tidx %% N_bt; // Basis function mapped to this thread\n"
5322: " const int cidx = tidx %% N_comp; // Basis component mapped to this thread\n"
5323: " const int qidx = tidx %% N_q; // Quadrature point mapped to this thread\n"
5324: " const int blbidx = tidx %% N_q + blidx*N_q; // Cell mapped to this thread in the basis phase\n"
5325: " const int blqidx = tidx %% N_b + blidx*N_b; // Cell mapped to this thread in the quadrature phase\n"
5326: " const int gidx = get_group_id(1)*get_num_groups(0) + get_group_id(0);\n"
5327: " const int Goffset = gidx*N_cb*N_bc;\n",
5328: &count, dim, N_bl, N_b, N_c, N_q);STRING_ERROR_CHECK("Message to short");
5329: /* Local memory */
5330: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
5331: "\n"
5332: " /* Quadrature data */\n"
5333: " %s w; // $w_q$, Quadrature weight at $x_q$\n"
5334: " __local %s phi_i[%d]; //[N_bt*N_q]; // $\\phi_i(x_q)$, Value of the basis function $i$ at $x_q$\n"
5335: " __local %s%d phiDer_i[%d]; //[N_bt*N_q]; // $\\frac{\\partial\\phi_i(x_q)}{\\partial x_d}$, Value of the derivative of basis function $i$ in direction $x_d$ at $x_q$\n"
5336: " /* Geometric data */\n"
5337: " __local %s detJ[%d]; //[N_t]; // $|J(x_q)|$, Jacobian determinant at $x_q$\n"
5338: " __local %s invJ[%d];//[N_t*dim*dim]; // $J^{-1}(x_q)$, Jacobian inverse at $x_q$\n",
5339: &count, numeric_str, numeric_str, N_b*N_c*N_q, numeric_str, dim, N_b*N_c*N_q, numeric_str, N_t,
5340: numeric_str, N_t*dim*dim, numeric_str, N_t*N_b*N_c);STRING_ERROR_CHECK("Message to short");
5341: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
5342: " /* FEM data */\n"
5343: " __local %s u_i[%d]; //[N_t*N_bt]; // Coefficients $u_i$ of the field $u|_{\\mathcal{T}} = \\sum_i u_i \\phi_i$\n",
5344: &count, numeric_str, N_t*N_b*N_c);STRING_ERROR_CHECK("Message to short");
5345: if (useAux) {
5346: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
5347: " __local %s a_i[%d]; //[N_t]; // Coefficients $a_i$ of the auxiliary field $a|_{\\mathcal{T}} = \\sum_i a_i \\phi^R_i$\n",
5348: &count, numeric_str, N_t);STRING_ERROR_CHECK("Message to short");
5349: }
5350: if (useF0) {
5351: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
5352: " /* Intermediate calculations */\n"
5353: " __local %s f_0[%d]; //[N_t*N_sqc]; // $f_0(u(x_q), \\nabla u(x_q)) |J(x_q)| w_q$\n",
5354: &count, numeric_str, N_t*N_q);STRING_ERROR_CHECK("Message to short");
5355: }
5356: if (useF1) {
5357: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
5358: " __local %s%d f_1[%d]; //[N_t*N_sqc]; // $f_1(u(x_q), \\nabla u(x_q)) |J(x_q)| w_q$\n",
5359: &count, numeric_str, dim, N_t*N_q);STRING_ERROR_CHECK("Message to short");
5360: }
5361: /* TODO: If using elasticity, put in mu/lambda coefficients */
5362: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
5363: " /* Output data */\n"
5364: " %s e_i; // Coefficient $e_i$ of the residual\n\n",
5365: &count, numeric_str);STRING_ERROR_CHECK("Message to short");
5366: /* One-time loads */
5367: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
5368: " /* These should be generated inline */\n"
5369: " /* Load quadrature weights */\n"
5370: " w = weights[qidx];\n"
5371: " /* Load basis tabulation \\phi_i for this cell */\n"
5372: " if (tidx < N_bt*N_q) {\n"
5373: " phi_i[tidx] = Basis[tidx];\n"
5374: " phiDer_i[tidx] = BasisDerivatives[tidx];\n"
5375: " }\n\n",
5376: &count);STRING_ERROR_CHECK("Message to short");
5377: /* Batch loads */
5378: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
5379: " for (int batch = 0; batch < N_cb; ++batch) {\n"
5380: " /* Load geometry */\n"
5381: " detJ[tidx] = jacobianDeterminants[Goffset+batch*N_bc+tidx];\n"
5382: " for (int n = 0; n < dim*dim; ++n) {\n"
5383: " const int offset = n*N_t;\n"
5384: " invJ[offset+tidx] = jacobianInverses[(Goffset+batch*N_bc)*dim*dim+offset+tidx];\n"
5385: " }\n"
5386: " /* Load coefficients u_i for this cell */\n"
5387: " for (int n = 0; n < N_bt; ++n) {\n"
5388: " const int offset = n*N_t;\n"
5389: " u_i[offset+tidx] = coefficients[(Goffset*N_bt)+batch*N_t*N_b+offset+tidx];\n"
5390: " }\n",
5391: &count);STRING_ERROR_CHECK("Message to short");
5392: if (useAux) {
5393: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
5394: " /* Load coefficients a_i for this cell */\n"
5395: " /* TODO: This should not be N_t here, it should be N_bc*N_comp_aux */\n"
5396: " a_i[tidx] = coefficientsAux[Goffset+batch*N_t+tidx];\n",
5397: &count);STRING_ERROR_CHECK("Message to short");
5398: }
5399: /* Quadrature phase */
5400: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
5401: " barrier(CLK_LOCAL_MEM_FENCE);\n"
5402: "\n"
5403: " /* Map coefficients to values at quadrature points */\n"
5404: " for (int c = 0; c < N_sqc; ++c) {\n"
5405: " const int cell = c*N_bl*N_b + blqidx;\n"
5406: " const int fidx = (cell*N_q + qidx)*N_comp + cidx;\n",
5407: &count);STRING_ERROR_CHECK("Message to short");
5408: if (useField) {
5409: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
5410: " %s u[%d]; //[N_comp]; // $u(x_q)$, Value of the field at $x_q$\n",
5411: &count, numeric_str, N_c);STRING_ERROR_CHECK("Message to short");
5412: }
5413: if (useFieldDer) {
5414: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
5415: " %s%d gradU[%d]; //[N_comp]; // $\\nabla u(x_q)$, Value of the field gradient at $x_q$\n",
5416: &count, numeric_str, dim, N_c);STRING_ERROR_CHECK("Message to short");
5417: }
5418: if (useFieldAux) {
5419: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
5420: " %s a[%d]; //[1]; // $a(x_q)$, Value of the auxiliary fields at $x_q$\n",
5421: &count, numeric_str, 1);STRING_ERROR_CHECK("Message to short");
5422: }
5423: if (useFieldDerAux) {
5424: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
5425: " %s%d gradA[%d]; //[1]; // $\\nabla a(x_q)$, Value of the auxiliary field gradient at $x_q$\n",
5426: &count, numeric_str, dim, 1);STRING_ERROR_CHECK("Message to short");
5427: }
5428: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
5429: "\n"
5430: " for (int comp = 0; comp < N_comp; ++comp) {\n",
5431: &count);STRING_ERROR_CHECK("Message to short");
5432: if (useField) {PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, " u[comp] = 0.0;\n", &count);STRING_ERROR_CHECK("Message to short");}
5433: if (useFieldDer) {
5434: switch (dim) {
5435: case 1:
5436: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, " gradU[comp].x = 0.0;\n", &count);STRING_ERROR_CHECK("Message to short");break;
5437: case 2:
5438: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, " gradU[comp].x = 0.0; gradU[comp].y = 0.0;\n", &count);STRING_ERROR_CHECK("Message to short");break;
5439: case 3:
5440: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, " gradU[comp].x = 0.0; gradU[comp].y = 0.0; gradU[comp].z = 0.0;\n", &count);STRING_ERROR_CHECK("Message to short");break;
5441: }
5442: }
5443: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
5444: " }\n",
5445: &count);STRING_ERROR_CHECK("Message to short");
5446: if (useFieldAux) {
5447: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, " a[0] = 0.0;\n", &count);STRING_ERROR_CHECK("Message to short");
5448: }
5449: if (useFieldDerAux) {
5450: switch (dim) {
5451: case 1:
5452: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, " gradA[0].x = 0.0;\n", &count);STRING_ERROR_CHECK("Message to short");break;
5453: case 2:
5454: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, " gradA[0].x = 0.0; gradA[0].y = 0.0;\n", &count);STRING_ERROR_CHECK("Message to short");break;
5455: case 3:
5456: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, " gradA[0].x = 0.0; gradA[0].y = 0.0; gradA[0].z = 0.0;\n", &count);STRING_ERROR_CHECK("Message to short");break;
5457: }
5458: }
5459: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
5460: " /* Get field and derivatives at this quadrature point */\n"
5461: " for (int i = 0; i < N_b; ++i) {\n"
5462: " for (int comp = 0; comp < N_comp; ++comp) {\n"
5463: " const int b = i*N_comp+comp;\n"
5464: " const int pidx = qidx*N_bt + b;\n"
5465: " const int uidx = cell*N_bt + b;\n"
5466: " %s%d realSpaceDer;\n\n",
5467: &count, numeric_str, dim);STRING_ERROR_CHECK("Message to short");
5468: if (useField) {PetscSNPrintfCount(string_tail, end_of_buffer - string_tail," u[comp] += u_i[uidx]*phi_i[pidx];\n", &count);STRING_ERROR_CHECK("Message to short");}
5469: if (useFieldDer) {
5470: switch (dim) {
5471: case 2:
5472: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
5473: " realSpaceDer.x = invJ[cell*dim*dim+0*dim+0]*phiDer_i[pidx].x + invJ[cell*dim*dim+1*dim+0]*phiDer_i[pidx].y;\n"
5474: " gradU[comp].x += u_i[uidx]*realSpaceDer.x;\n"
5475: " realSpaceDer.y = invJ[cell*dim*dim+0*dim+1]*phiDer_i[pidx].x + invJ[cell*dim*dim+1*dim+1]*phiDer_i[pidx].y;\n"
5476: " gradU[comp].y += u_i[uidx]*realSpaceDer.y;\n",
5477: &count);STRING_ERROR_CHECK("Message to short");break;
5478: case 3:
5479: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
5480: " realSpaceDer.x = invJ[cell*dim*dim+0*dim+0]*phiDer_i[pidx].x + invJ[cell*dim*dim+1*dim+0]*phiDer_i[pidx].y + invJ[cell*dim*dim+2*dim+0]*phiDer_i[pidx].z;\n"
5481: " gradU[comp].x += u_i[uidx]*realSpaceDer.x;\n"
5482: " realSpaceDer.y = invJ[cell*dim*dim+0*dim+1]*phiDer_i[pidx].x + invJ[cell*dim*dim+1*dim+1]*phiDer_i[pidx].y + invJ[cell*dim*dim+2*dim+1]*phiDer_i[pidx].z;\n"
5483: " gradU[comp].y += u_i[uidx]*realSpaceDer.y;\n"
5484: " realSpaceDer.z = invJ[cell*dim*dim+0*dim+2]*phiDer_i[pidx].x + invJ[cell*dim*dim+1*dim+2]*phiDer_i[pidx].y + invJ[cell*dim*dim+2*dim+2]*phiDer_i[pidx].z;\n"
5485: " gradU[comp].z += u_i[uidx]*realSpaceDer.z;\n",
5486: &count);STRING_ERROR_CHECK("Message to short");break;
5487: }
5488: }
5489: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
5490: " }\n"
5491: " }\n",
5492: &count);STRING_ERROR_CHECK("Message to short");
5493: if (useFieldAux) {
5494: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail," a[0] += a_i[cell];\n", &count);STRING_ERROR_CHECK("Message to short");
5495: }
5496: /* Calculate residual at quadrature points: Should be generated by an weak form egine */
5497: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
5498: " /* Process values at quadrature points */\n",
5499: &count);STRING_ERROR_CHECK("Message to short");
5500: switch (op) {
5501: case LAPLACIAN:
5502: if (useF0) {PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, " f_0[fidx] = 4.0;\n", &count);STRING_ERROR_CHECK("Message to short");}
5503: if (useF1) {
5504: if (useAux) {PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, " f_1[fidx] = a[0]*gradU[cidx];\n", &count);STRING_ERROR_CHECK("Message to short");}
5505: else {PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, " f_1[fidx] = gradU[cidx];\n", &count);STRING_ERROR_CHECK("Message to short");}
5506: }
5507: break;
5508: case ELASTICITY:
5509: if (useF0) {PetscSNPrintfCount(string_tail, end_of_buffer - string_tail, " f_0[fidx] = 4.0;\n", &count);STRING_ERROR_CHECK("Message to short");}
5510: if (useF1) {
5511: switch (dim) {
5512: case 2:
5513: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
5514: " switch (cidx) {\n"
5515: " case 0:\n"
5516: " f_1[fidx].x = lambda*(gradU[0].x + gradU[1].y) + mu*(gradU[0].x + gradU[0].x);\n"
5517: " f_1[fidx].y = lambda*(gradU[0].x + gradU[1].y) + mu*(gradU[0].y + gradU[1].x);\n"
5518: " break;\n"
5519: " case 1:\n"
5520: " f_1[fidx].x = lambda*(gradU[0].x + gradU[1].y) + mu*(gradU[1].x + gradU[0].y);\n"
5521: " f_1[fidx].y = lambda*(gradU[0].x + gradU[1].y) + mu*(gradU[1].y + gradU[1].y);\n"
5522: " }\n",
5523: &count);STRING_ERROR_CHECK("Message to short");break;
5524: case 3:
5525: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
5526: " switch (cidx) {\n"
5527: " case 0:\n"
5528: " f_1[fidx].x = lambda*(gradU[0].x + gradU[1].y + gradU[2].z) + mu*(gradU[0].x + gradU[0].x);\n"
5529: " f_1[fidx].y = lambda*(gradU[0].x + gradU[1].y + gradU[2].z) + mu*(gradU[0].y + gradU[1].x);\n"
5530: " f_1[fidx].z = lambda*(gradU[0].x + gradU[1].y + gradU[2].z) + mu*(gradU[0].z + gradU[2].x);\n"
5531: " break;\n"
5532: " case 1:\n"
5533: " f_1[fidx].x = lambda*(gradU[0].x + gradU[1].y + gradU[2].z) + mu*(gradU[1].x + gradU[0].y);\n"
5534: " f_1[fidx].y = lambda*(gradU[0].x + gradU[1].y + gradU[2].z) + mu*(gradU[1].y + gradU[1].y);\n"
5535: " f_1[fidx].z = lambda*(gradU[0].x + gradU[1].y + gradU[2].z) + mu*(gradU[1].y + gradU[2].y);\n"
5536: " break;\n"
5537: " case 2:\n"
5538: " f_1[fidx].x = lambda*(gradU[0].x + gradU[1].y + gradU[2].z) + mu*(gradU[2].x + gradU[0].z);\n"
5539: " f_1[fidx].y = lambda*(gradU[0].x + gradU[1].y + gradU[2].z) + mu*(gradU[2].y + gradU[1].z);\n"
5540: " f_1[fidx].z = lambda*(gradU[0].x + gradU[1].y + gradU[2].z) + mu*(gradU[2].y + gradU[2].z);\n"
5541: " }\n",
5542: &count);STRING_ERROR_CHECK("Message to short");break;
5543: }}
5544: break;
5545: default:
5546: SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_SUP, "PDE operator %d is not supported", op);
5547: }
5548: if (useF0) {PetscSNPrintfCount(string_tail, end_of_buffer - string_tail," f_0[fidx] *= detJ[cell]*w;\n", &count);STRING_ERROR_CHECK("Message to short");}
5549: if (useF1) {
5550: switch (dim) {
5551: case 1:
5552: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail," f_1[fidx].x *= detJ[cell]*w;\n", &count);STRING_ERROR_CHECK("Message to short");break;
5553: case 2:
5554: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail," f_1[fidx].x *= detJ[cell]*w; f_1[fidx].y *= detJ[cell]*w;\n", &count);STRING_ERROR_CHECK("Message to short");break;
5555: case 3:
5556: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail," f_1[fidx].x *= detJ[cell]*w; f_1[fidx].y *= detJ[cell]*w; f_1[fidx].z *= detJ[cell]*w;\n", &count);STRING_ERROR_CHECK("Message to short");break;
5557: }
5558: }
5559: /* Thread transpose */
5560: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
5561: " }\n\n"
5562: " /* ==== TRANSPOSE THREADS ==== */\n"
5563: " barrier(CLK_LOCAL_MEM_FENCE);\n\n",
5564: &count);STRING_ERROR_CHECK("Message to short");
5565: /* Basis phase */
5566: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
5567: " /* Map values at quadrature points to coefficients */\n"
5568: " for (int c = 0; c < N_sbc; ++c) {\n"
5569: " const int cell = c*N_bl*N_q + blbidx; /* Cell number in batch */\n"
5570: "\n"
5571: " e_i = 0.0;\n"
5572: " for (int q = 0; q < N_q; ++q) {\n"
5573: " const int pidx = q*N_bt + bidx;\n"
5574: " const int fidx = (cell*N_q + q)*N_comp + cidx;\n"
5575: " %s%d realSpaceDer;\n\n",
5576: &count, numeric_str, dim);STRING_ERROR_CHECK("Message to short");
5578: if (useF0) {PetscSNPrintfCount(string_tail, end_of_buffer - string_tail," e_i += phi_i[pidx]*f_0[fidx];\n", &count);STRING_ERROR_CHECK("Message to short");}
5579: if (useF1) {
5580: switch (dim) {
5581: case 2:
5582: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
5583: " realSpaceDer.x = invJ[cell*dim*dim+0*dim+0]*phiDer_i[pidx].x + invJ[cell*dim*dim+1*dim+0]*phiDer_i[pidx].y;\n"
5584: " e_i += realSpaceDer.x*f_1[fidx].x;\n"
5585: " realSpaceDer.y = invJ[cell*dim*dim+0*dim+1]*phiDer_i[pidx].x + invJ[cell*dim*dim+1*dim+1]*phiDer_i[pidx].y;\n"
5586: " e_i += realSpaceDer.y*f_1[fidx].y;\n",
5587: &count);STRING_ERROR_CHECK("Message to short");break;
5588: case 3:
5589: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
5590: " realSpaceDer.x = invJ[cell*dim*dim+0*dim+0]*phiDer_i[pidx].x + invJ[cell*dim*dim+1*dim+0]*phiDer_i[pidx].y + invJ[cell*dim*dim+2*dim+0]*phiDer_i[pidx].z;\n"
5591: " e_i += realSpaceDer.x*f_1[fidx].x;\n"
5592: " realSpaceDer.y = invJ[cell*dim*dim+0*dim+1]*phiDer_i[pidx].x + invJ[cell*dim*dim+1*dim+1]*phiDer_i[pidx].y + invJ[cell*dim*dim+2*dim+1]*phiDer_i[pidx].z;\n"
5593: " e_i += realSpaceDer.y*f_1[fidx].y;\n"
5594: " realSpaceDer.z = invJ[cell*dim*dim+0*dim+2]*phiDer_i[pidx].x + invJ[cell*dim*dim+1*dim+2]*phiDer_i[pidx].y + invJ[cell*dim*dim+2*dim+2]*phiDer_i[pidx].z;\n"
5595: " e_i += realSpaceDer.z*f_1[fidx].z;\n",
5596: &count);STRING_ERROR_CHECK("Message to short");break;
5597: }
5598: }
5599: PetscSNPrintfCount(string_tail, end_of_buffer - string_tail,
5600: " }\n"
5601: " /* Write element vector for N_{cbc} cells at a time */\n"
5602: " elemVec[(Goffset + batch*N_bc + c*N_bl*N_q)*N_bt + tidx] = e_i;\n"
5603: " }\n"
5604: " /* ==== Could do one write per batch ==== */\n"
5605: " }\n"
5606: " return;\n"
5607: "}\n",
5608: &count);STRING_ERROR_CHECK("Message to short");
5609: return(0);
5610: }
5612: PetscErrorCode PetscFEOpenCLGetIntegrationKernel(PetscFE fem, PetscBool useAux, cl_program *ocl_prog, cl_kernel *ocl_kernel)
5613: {
5614: PetscFE_OpenCL *ocl = (PetscFE_OpenCL *) fem->data;
5615: PetscInt dim, N_bl;
5616: PetscBool flg;
5617: char *buffer;
5618: size_t len;
5619: char errMsg[8192];
5620: cl_int ierr2;
5621: PetscErrorCode ierr;
5624: PetscFEGetSpatialDimension(fem, &dim);
5625: PetscMalloc1(8192, &buffer);
5626: PetscFEGetTileSizes(fem, NULL, &N_bl, NULL, NULL);
5627: PetscFEOpenCLGenerateIntegrationCode(fem, &buffer, 8192, useAux, N_bl);
5628: PetscOptionsHasName(((PetscObject)fem)->options,((PetscObject)fem)->prefix, "-petscfe_opencl_kernel_print", &flg);
5629: if (flg) {PetscPrintf(PetscObjectComm((PetscObject) fem), "OpenCL FE Integration Kernel:\n%s\n", buffer);}
5630: len = strlen(buffer);
5631: *ocl_prog = clCreateProgramWithSource(ocl->ctx_id, 1, (const char **) &buffer, &len, &ierr2);CHKERRQ(ierr2);
5632: clBuildProgram(*ocl_prog, 0, NULL, NULL, NULL, NULL);
5633: if (ierr != CL_SUCCESS) {
5634: clGetProgramBuildInfo(*ocl_prog, ocl->dev_id, CL_PROGRAM_BUILD_LOG, 8192*sizeof(char), &errMsg, NULL);
5635: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Build failed! Log:\n %s", errMsg);
5636: }
5637: PetscFree(buffer);
5638: *ocl_kernel = clCreateKernel(*ocl_prog, "integrateElementQuadrature", &ierr);
5639: return(0);
5640: }
5642: PetscErrorCode PetscFEOpenCLCalculateGrid(PetscFE fem, PetscInt N, PetscInt blockSize, size_t *x, size_t *y, size_t *z)
5643: {
5644: const PetscInt Nblocks = N/blockSize;
5647: if (N % blockSize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid block size %d for %d elements", blockSize, N);
5648: *z = 1;
5649: for (*x = (size_t) (PetscSqrtReal(Nblocks) + 0.5); *x > 0; --*x) {
5650: *y = Nblocks / *x;
5651: if (*x * *y == Nblocks) break;
5652: }
5653: if (*x * *y != Nblocks) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Could not find partition for %d with block size %d", N, blockSize);
5654: return(0);
5655: }
5657: PetscErrorCode PetscFEOpenCLLogResidual(PetscFE fem, PetscLogDouble time, PetscLogDouble flops)
5658: {
5659: PetscFE_OpenCL *ocl = (PetscFE_OpenCL *) fem->data;
5660: PetscStageLog stageLog;
5661: PetscEventPerfLog eventLog = NULL;
5662: PetscInt stage;
5663: PetscErrorCode ierr;
5666: PetscLogGetStageLog(&stageLog);
5667: PetscStageLogGetCurrent(stageLog, &stage);
5668: PetscStageLogGetEventPerfLog(stageLog, stage, &eventLog);
5669: /* Log performance info */
5670: eventLog->eventInfo[ocl->residualEvent].count++;
5671: eventLog->eventInfo[ocl->residualEvent].time += time;
5672: eventLog->eventInfo[ocl->residualEvent].flops += flops;
5673: return(0);
5674: }
5676: PetscErrorCode PetscFEIntegrateResidual_OpenCL(PetscFE fem, PetscDS prob, PetscInt field, PetscInt Ne, PetscFECellGeom *cgeom,
5677: const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
5678: {
5679: /* Nbc = batchSize */
5680: PetscFE_OpenCL *ocl = (PetscFE_OpenCL *) fem->data;
5681: PetscPointFunc f0_func;
5682: PetscPointFunc f1_func;
5683: PetscQuadrature q;
5684: PetscInt dim, qNc;
5685: PetscInt N_b; /* The number of basis functions */
5686: PetscInt N_comp; /* The number of basis function components */
5687: PetscInt N_bt; /* The total number of scalar basis functions */
5688: PetscInt N_q; /* The number of quadrature points */
5689: PetscInt N_bst; /* The block size, LCM(N_bt, N_q), Notice that a block is not process simultaneously */
5690: PetscInt N_t; /* The number of threads, N_bst * N_bl */
5691: PetscInt N_bl; /* The number of blocks */
5692: PetscInt N_bc; /* The batch size, N_bl*N_q*N_b */
5693: PetscInt N_cb; /* The number of batches */
5694: PetscInt numFlops, f0Flops = 0, f1Flops = 0;
5695: PetscBool useAux = probAux ? PETSC_TRUE : PETSC_FALSE;
5696: PetscBool useField = PETSC_FALSE;
5697: PetscBool useFieldDer = PETSC_TRUE;
5698: PetscBool useF0 = PETSC_TRUE;
5699: PetscBool useF1 = PETSC_TRUE;
5700: /* OpenCL variables */
5701: cl_program ocl_prog;
5702: cl_kernel ocl_kernel;
5703: cl_event ocl_ev; /* The event for tracking kernel execution */
5704: cl_ulong ns_start; /* Nanoseconds counter on GPU at kernel start */
5705: cl_ulong ns_end; /* Nanoseconds counter on GPU at kernel stop */
5706: cl_mem o_jacobianInverses, o_jacobianDeterminants;
5707: cl_mem o_coefficients, o_coefficientsAux, o_elemVec;
5708: float *f_coeff = NULL, *f_coeffAux = NULL, *f_invJ = NULL, *f_detJ = NULL;
5709: double *d_coeff = NULL, *d_coeffAux = NULL, *d_invJ = NULL, *d_detJ = NULL;
5710: PetscReal *r_invJ = NULL, *r_detJ = NULL;
5711: void *oclCoeff, *oclCoeffAux, *oclInvJ, *oclDetJ;
5712: size_t local_work_size[3], global_work_size[3];
5713: size_t realSize, x, y, z;
5714: const PetscReal *points, *weights;
5715: PetscErrorCode ierr;
5718: if (!Ne) {PetscFEOpenCLLogResidual(fem, 0.0, 0.0); return(0);}
5719: PetscFEGetSpatialDimension(fem, &dim);
5720: PetscFEGetQuadrature(fem, &q);
5721: PetscQuadratureGetData(q, NULL, &qNc, &N_q, &points, &weights);
5722: if (qNc != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supports scalar quadrature, not %D components\n", qNc);
5723: PetscFEGetDimension(fem, &N_b);
5724: PetscFEGetNumComponents(fem, &N_comp);
5725: PetscDSGetResidual(prob, field, &f0_func, &f1_func);
5726: PetscFEGetTileSizes(fem, NULL, &N_bl, &N_bc, &N_cb);
5727: N_bt = N_b*N_comp;
5728: N_bst = N_bt*N_q;
5729: N_t = N_bst*N_bl;
5730: if (N_bc*N_comp != N_t) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of threads %d should be %d * %d", N_t, N_bc, N_comp);
5731: /* Calculate layout */
5732: if (Ne % (N_cb*N_bc)) { /* Remainder cells */
5733: PetscFEIntegrateResidual_Basic(fem, prob, field, Ne, cgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec);
5734: return(0);
5735: }
5736: PetscFEOpenCLCalculateGrid(fem, Ne, N_cb*N_bc, &x, &y, &z);
5737: local_work_size[0] = N_bc*N_comp;
5738: local_work_size[1] = 1;
5739: local_work_size[2] = 1;
5740: global_work_size[0] = x * local_work_size[0];
5741: global_work_size[1] = y * local_work_size[1];
5742: global_work_size[2] = z * local_work_size[2];
5743: PetscInfo7(fem, "GPU layout grid(%d,%d,%d) block(%d,%d,%d) with %d batches\n", x, y, z, local_work_size[0], local_work_size[1], local_work_size[2], N_cb);
5744: PetscInfo2(fem, " N_t: %d, N_cb: %d\n", N_t, N_cb);
5745: /* Generate code */
5746: if (probAux) {
5747: PetscSpace P;
5748: PetscInt NfAux, order, f;
5750: PetscDSGetNumFields(probAux, &NfAux);
5751: for (f = 0; f < NfAux; ++f) {
5752: PetscFE feAux;
5754: PetscDSGetDiscretization(probAux, f, (PetscObject *) &feAux);
5755: PetscFEGetBasisSpace(feAux, &P);
5756: PetscSpaceGetOrder(P, &order);
5757: if (order > 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Can only handle P0 coefficient fields");
5758: }
5759: }
5760: PetscFEOpenCLGetIntegrationKernel(fem, useAux, &ocl_prog, &ocl_kernel);
5761: /* Create buffers on the device and send data over */
5762: PetscDataTypeGetSize(ocl->realType, &realSize);
5763: if (sizeof(PetscReal) != realSize) {
5764: switch (ocl->realType) {
5765: case PETSC_FLOAT:
5766: {
5767: PetscInt c, b, d;
5769: PetscMalloc4(Ne*N_bt,&f_coeff,Ne,&f_coeffAux,Ne*dim*dim,&f_invJ,Ne,&f_detJ);
5770: for (c = 0; c < Ne; ++c) {
5771: f_detJ[c] = (float) cgeom[c].detJ;
5772: for (d = 0; d < dim*dim; ++d) {
5773: f_invJ[c*dim*dim+d] = (float) cgeom[c].invJ[d];
5774: }
5775: for (b = 0; b < N_bt; ++b) {
5776: f_coeff[c*N_bt+b] = (float) coefficients[c*N_bt+b];
5777: }
5778: }
5779: if (coefficientsAux) { /* Assume P0 */
5780: for (c = 0; c < Ne; ++c) {
5781: f_coeffAux[c] = (float) coefficientsAux[c];
5782: }
5783: }
5784: oclCoeff = (void *) f_coeff;
5785: if (coefficientsAux) {
5786: oclCoeffAux = (void *) f_coeffAux;
5787: } else {
5788: oclCoeffAux = NULL;
5789: }
5790: oclInvJ = (void *) f_invJ;
5791: oclDetJ = (void *) f_detJ;
5792: }
5793: break;
5794: case PETSC_DOUBLE:
5795: {
5796: PetscInt c, b, d;
5798: PetscMalloc4(Ne*N_bt,&d_coeff,Ne,&d_coeffAux,Ne*dim*dim,&d_invJ,Ne,&d_detJ);
5799: for (c = 0; c < Ne; ++c) {
5800: d_detJ[c] = (double) cgeom[c].detJ;
5801: for (d = 0; d < dim*dim; ++d) {
5802: d_invJ[c*dim*dim+d] = (double) cgeom[c].invJ[d];
5803: }
5804: for (b = 0; b < N_bt; ++b) {
5805: d_coeff[c*N_bt+b] = (double) coefficients[c*N_bt+b];
5806: }
5807: }
5808: if (coefficientsAux) { /* Assume P0 */
5809: for (c = 0; c < Ne; ++c) {
5810: d_coeffAux[c] = (double) coefficientsAux[c];
5811: }
5812: }
5813: oclCoeff = (void *) d_coeff;
5814: if (coefficientsAux) {
5815: oclCoeffAux = (void *) d_coeffAux;
5816: } else {
5817: oclCoeffAux = NULL;
5818: }
5819: oclInvJ = (void *) d_invJ;
5820: oclDetJ = (void *) d_detJ;
5821: }
5822: break;
5823: default:
5824: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported PETSc type %d", ocl->realType);
5825: }
5826: } else {
5827: PetscInt c, d;
5829: PetscMalloc2(Ne*dim*dim,&r_invJ,Ne,&r_detJ);
5830: for (c = 0; c < Ne; ++c) {
5831: r_detJ[c] = cgeom[c].detJ;
5832: for (d = 0; d < dim*dim; ++d) {
5833: r_invJ[c*dim*dim+d] = cgeom[c].invJ[d];
5834: }
5835: }
5836: oclCoeff = (void *) coefficients;
5837: oclCoeffAux = (void *) coefficientsAux;
5838: oclInvJ = (void *) r_invJ;
5839: oclDetJ = (void *) r_detJ;
5840: }
5841: o_coefficients = clCreateBuffer(ocl->ctx_id, CL_MEM_READ_ONLY | CL_MEM_COPY_HOST_PTR, Ne*N_bt * realSize, oclCoeff, &ierr);
5842: if (coefficientsAux) {
5843: o_coefficientsAux = clCreateBuffer(ocl->ctx_id, CL_MEM_READ_ONLY | CL_MEM_COPY_HOST_PTR, Ne * realSize, oclCoeffAux, &ierr);
5844: } else {
5845: o_coefficientsAux = clCreateBuffer(ocl->ctx_id, CL_MEM_READ_ONLY, Ne * realSize, oclCoeffAux, &ierr);
5846: }
5847: o_jacobianInverses = clCreateBuffer(ocl->ctx_id, CL_MEM_READ_ONLY | CL_MEM_COPY_HOST_PTR, Ne*dim*dim * realSize, oclInvJ, &ierr);
5848: o_jacobianDeterminants = clCreateBuffer(ocl->ctx_id, CL_MEM_READ_ONLY | CL_MEM_COPY_HOST_PTR, Ne * realSize, oclDetJ, &ierr);
5849: o_elemVec = clCreateBuffer(ocl->ctx_id, CL_MEM_WRITE_ONLY, Ne*N_bt * realSize, NULL, &ierr);
5850: /* Kernel launch */
5851: clSetKernelArg(ocl_kernel, 0, sizeof(cl_int), (void*) &N_cb);
5852: clSetKernelArg(ocl_kernel, 1, sizeof(cl_mem), (void*) &o_coefficients);
5853: clSetKernelArg(ocl_kernel, 2, sizeof(cl_mem), (void*) &o_coefficientsAux);
5854: clSetKernelArg(ocl_kernel, 3, sizeof(cl_mem), (void*) &o_jacobianInverses);
5855: clSetKernelArg(ocl_kernel, 4, sizeof(cl_mem), (void*) &o_jacobianDeterminants);
5856: clSetKernelArg(ocl_kernel, 5, sizeof(cl_mem), (void*) &o_elemVec);
5857: clEnqueueNDRangeKernel(ocl->queue_id, ocl_kernel, 3, NULL, global_work_size, local_work_size, 0, NULL, &ocl_ev);
5858: /* Read data back from device */
5859: if (sizeof(PetscReal) != realSize) {
5860: switch (ocl->realType) {
5861: case PETSC_FLOAT:
5862: {
5863: float *elem;
5864: PetscInt c, b;
5866: PetscFree4(f_coeff,f_coeffAux,f_invJ,f_detJ);
5867: PetscMalloc1(Ne*N_bt, &elem);
5868: clEnqueueReadBuffer(ocl->queue_id, o_elemVec, CL_TRUE, 0, Ne*N_bt * realSize, elem, 0, NULL, NULL);
5869: for (c = 0; c < Ne; ++c) {
5870: for (b = 0; b < N_bt; ++b) {
5871: elemVec[c*N_bt+b] = (PetscScalar) elem[c*N_bt+b];
5872: }
5873: }
5874: PetscFree(elem);
5875: }
5876: break;
5877: case PETSC_DOUBLE:
5878: {
5879: double *elem;
5880: PetscInt c, b;
5882: PetscFree4(d_coeff,d_coeffAux,d_invJ,d_detJ);
5883: PetscMalloc1(Ne*N_bt, &elem);
5884: clEnqueueReadBuffer(ocl->queue_id, o_elemVec, CL_TRUE, 0, Ne*N_bt * realSize, elem, 0, NULL, NULL);
5885: for (c = 0; c < Ne; ++c) {
5886: for (b = 0; b < N_bt; ++b) {
5887: elemVec[c*N_bt+b] = (PetscScalar) elem[c*N_bt+b];
5888: }
5889: }
5890: PetscFree(elem);
5891: }
5892: break;
5893: default:
5894: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported PETSc type %d", ocl->realType);
5895: }
5896: } else {
5897: PetscFree2(r_invJ,r_detJ);
5898: clEnqueueReadBuffer(ocl->queue_id, o_elemVec, CL_TRUE, 0, Ne*N_bt * realSize, elemVec, 0, NULL, NULL);
5899: }
5900: /* Log performance */
5901: clGetEventProfilingInfo(ocl_ev, CL_PROFILING_COMMAND_START, sizeof(cl_ulong), &ns_start, NULL);
5902: clGetEventProfilingInfo(ocl_ev, CL_PROFILING_COMMAND_END, sizeof(cl_ulong), &ns_end, NULL);
5903: f0Flops = 0;
5904: switch (ocl->op) {
5905: case LAPLACIAN:
5906: f1Flops = useAux ? dim : 0;break;
5907: case ELASTICITY:
5908: f1Flops = 2*dim*dim;break;
5909: }
5910: numFlops = Ne*(
5911: N_q*(
5912: N_b*N_comp*((useField ? 2 : 0) + (useFieldDer ? 2*dim*(dim + 1) : 0))
5913: /*+
5914: N_ba*N_compa*((useFieldAux ? 2 : 0) + (useFieldDerAux ? 2*dim*(dim + 1) : 0))*/
5915: +
5916: N_comp*((useF0 ? f0Flops + 2 : 0) + (useF1 ? f1Flops + 2*dim : 0)))
5917: +
5918: N_b*((useF0 ? 2 : 0) + (useF1 ? 2*dim*(dim + 1) : 0)));
5919: PetscFEOpenCLLogResidual(fem, (ns_end - ns_start)*1.0e-9, numFlops);
5920: /* Cleanup */
5921: clReleaseMemObject(o_coefficients);
5922: clReleaseMemObject(o_coefficientsAux);
5923: clReleaseMemObject(o_jacobianInverses);
5924: clReleaseMemObject(o_jacobianDeterminants);
5925: clReleaseMemObject(o_elemVec);
5926: clReleaseKernel(ocl_kernel);
5927: clReleaseProgram(ocl_prog);
5928: return(0);
5929: }
5931: PetscErrorCode PetscFEInitialize_OpenCL(PetscFE fem)
5932: {
5934: fem->ops->setfromoptions = NULL;
5935: fem->ops->setup = PetscFESetUp_Basic;
5936: fem->ops->view = NULL;
5937: fem->ops->destroy = PetscFEDestroy_OpenCL;
5938: fem->ops->getdimension = PetscFEGetDimension_Basic;
5939: fem->ops->gettabulation = PetscFEGetTabulation_Basic;
5940: fem->ops->integrateresidual = PetscFEIntegrateResidual_OpenCL;
5941: fem->ops->integratebdresidual = NULL/* PetscFEIntegrateBdResidual_OpenCL */;
5942: fem->ops->integratejacobianaction = NULL/* PetscFEIntegrateJacobianAction_OpenCL */;
5943: fem->ops->integratejacobian = PetscFEIntegrateJacobian_Basic;
5944: return(0);
5945: }
5947: /*MC
5948: PETSCFEOPENCL = "opencl" - A PetscFE object that integrates using a vectorized OpenCL implementation
5950: Level: intermediate
5952: .seealso: PetscFEType, PetscFECreate(), PetscFESetType()
5953: M*/
5955: PETSC_EXTERN PetscErrorCode PetscFECreate_OpenCL(PetscFE fem)
5956: {
5957: PetscFE_OpenCL *ocl;
5958: cl_uint num_platforms;
5959: cl_platform_id platform_ids[42];
5960: cl_uint num_devices;
5961: cl_device_id device_ids[42];
5962: cl_int ierr2;
5963: PetscErrorCode ierr;
5967: PetscNewLog(fem,&ocl);
5968: fem->data = ocl;
5970: /* Init Platform */
5971: clGetPlatformIDs(42, platform_ids, &num_platforms);
5972: if (!num_platforms) SETERRQ(PetscObjectComm((PetscObject) fem), PETSC_ERR_SUP, "No OpenCL platform found.");
5973: ocl->pf_id = platform_ids[0];
5974: /* Init Device */
5975: clGetDeviceIDs(ocl->pf_id, CL_DEVICE_TYPE_ALL, 42, device_ids, &num_devices);
5976: if (!num_devices) SETERRQ(PetscObjectComm((PetscObject) fem), PETSC_ERR_SUP, "No OpenCL device found.");
5977: ocl->dev_id = device_ids[0];
5978: /* Create context with one command queue */
5979: ocl->ctx_id = clCreateContext(0, 1, &(ocl->dev_id), NULL, NULL, &ierr2);CHKERRQ(ierr2);
5980: ocl->queue_id = clCreateCommandQueue(ocl->ctx_id, ocl->dev_id, CL_QUEUE_PROFILING_ENABLE, &ierr2);CHKERRQ(ierr2);
5981: /* Types */
5982: ocl->realType = PETSC_FLOAT;
5983: /* Register events */
5984: PetscLogEventRegister("OpenCL FEResidual", PETSCFE_CLASSID, &ocl->residualEvent);
5985: /* Equation handling */
5986: ocl->op = LAPLACIAN;
5988: PetscFEInitialize_OpenCL(fem);
5989: return(0);
5990: }
5992: PetscErrorCode PetscFEOpenCLSetRealType(PetscFE fem, PetscDataType realType)
5993: {
5994: PetscFE_OpenCL *ocl = (PetscFE_OpenCL *) fem->data;
5998: ocl->realType = realType;
5999: return(0);
6000: }
6002: PetscErrorCode PetscFEOpenCLGetRealType(PetscFE fem, PetscDataType *realType)
6003: {
6004: PetscFE_OpenCL *ocl = (PetscFE_OpenCL *) fem->data;
6009: *realType = ocl->realType;
6010: return(0);
6011: }
6013: #endif /* PETSC_HAVE_OPENCL */
6015: PetscErrorCode PetscFEDestroy_Composite(PetscFE fem)
6016: {
6017: PetscFE_Composite *cmp = (PetscFE_Composite *) fem->data;
6018: PetscErrorCode ierr;
6021: CellRefinerRestoreAffineTransforms_Internal(cmp->cellRefiner, &cmp->numSubelements, &cmp->v0, &cmp->jac, &cmp->invjac);
6022: PetscFree(cmp->embedding);
6023: PetscFree(cmp);
6024: return(0);
6025: }
6027: PetscErrorCode PetscFESetUp_Composite(PetscFE fem)
6028: {
6029: PetscFE_Composite *cmp = (PetscFE_Composite *) fem->data;
6030: DM K;
6031: PetscReal *subpoint;
6032: PetscBLASInt *pivots;
6033: PetscBLASInt n, info;
6034: PetscScalar *work, *invVscalar;
6035: PetscInt dim, pdim, spdim, j, s;
6036: PetscErrorCode ierr;
6039: /* Get affine mapping from reference cell to each subcell */
6040: PetscDualSpaceGetDM(fem->dualSpace, &K);
6041: DMGetDimension(K, &dim);
6042: DMPlexGetCellRefiner_Internal(K, &cmp->cellRefiner);
6043: CellRefinerGetAffineTransforms_Internal(cmp->cellRefiner, &cmp->numSubelements, &cmp->v0, &cmp->jac, &cmp->invjac);
6044: /* Determine dof embedding into subelements */
6045: PetscDualSpaceGetDimension(fem->dualSpace, &pdim);
6046: PetscSpaceGetDimension(fem->basisSpace, &spdim);
6047: PetscMalloc1(cmp->numSubelements*spdim,&cmp->embedding);
6048: DMGetWorkArray(K, dim, MPIU_REAL, &subpoint);
6049: for (s = 0; s < cmp->numSubelements; ++s) {
6050: PetscInt sd = 0;
6052: for (j = 0; j < pdim; ++j) {
6053: PetscBool inside;
6054: PetscQuadrature f;
6055: PetscInt d, e;
6057: PetscDualSpaceGetFunctional(fem->dualSpace, j, &f);
6058: /* Apply transform to first point, and check that point is inside subcell */
6059: for (d = 0; d < dim; ++d) {
6060: subpoint[d] = -1.0;
6061: for (e = 0; e < dim; ++e) subpoint[d] += cmp->invjac[(s*dim + d)*dim+e]*(f->points[e] - cmp->v0[s*dim+e]);
6062: }
6063: CellRefinerInCellTest_Internal(cmp->cellRefiner, subpoint, &inside);
6064: if (inside) {cmp->embedding[s*spdim+sd++] = j;}
6065: }
6066: if (sd != spdim) SETERRQ3(PetscObjectComm((PetscObject) fem), PETSC_ERR_PLIB, "Subelement %d has %d dual basis vectors != %d", s, sd, spdim);
6067: }
6068: DMRestoreWorkArray(K, dim, MPIU_REAL, &subpoint);
6069: /* Construct the change of basis from prime basis to nodal basis for each subelement */
6070: PetscMalloc1(cmp->numSubelements*spdim*spdim,&fem->invV);
6071: PetscMalloc2(spdim,&pivots,spdim,&work);
6072: #if defined(PETSC_USE_COMPLEX)
6073: PetscMalloc1(cmp->numSubelements*spdim*spdim,&invVscalar);
6074: #else
6075: invVscalar = fem->invV;
6076: #endif
6077: for (s = 0; s < cmp->numSubelements; ++s) {
6078: for (j = 0; j < spdim; ++j) {
6079: PetscReal *Bf;
6080: PetscQuadrature f;
6081: const PetscReal *points, *weights;
6082: PetscInt Nc, Nq, q, k;
6084: PetscDualSpaceGetFunctional(fem->dualSpace, cmp->embedding[s*spdim+j], &f);
6085: PetscQuadratureGetData(f, NULL, &Nc, &Nq, &points, &weights);
6086: PetscMalloc1(f->numPoints*spdim*Nc,&Bf);
6087: PetscSpaceEvaluate(fem->basisSpace, Nq, points, Bf, NULL, NULL);
6088: for (k = 0; k < spdim; ++k) {
6089: /* n_j \cdot \phi_k */
6090: invVscalar[(s*spdim + j)*spdim+k] = 0.0;
6091: for (q = 0; q < Nq; ++q) {
6092: invVscalar[(s*spdim + j)*spdim+k] += Bf[q*spdim+k]*weights[q];
6093: }
6094: }
6095: PetscFree(Bf);
6096: }
6097: n = spdim;
6098: PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&n, &n, &invVscalar[s*spdim*spdim], &n, pivots, &info));
6099: PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&n, &invVscalar[s*spdim*spdim], &n, pivots, work, &n, &info));
6100: }
6101: #if defined(PETSC_USE_COMPLEX)
6102: for (s = 0; s <cmp->numSubelements*spdim*spdim; s++) fem->invV[s] = PetscRealPart(invVscalar[s]);
6103: PetscFree(invVscalar);
6104: #endif
6105: PetscFree2(pivots,work);
6106: return(0);
6107: }
6109: PetscErrorCode PetscFEGetTabulation_Composite(PetscFE fem, PetscInt npoints, const PetscReal points[], PetscReal *B, PetscReal *D, PetscReal *H)
6110: {
6111: PetscFE_Composite *cmp = (PetscFE_Composite *) fem->data;
6112: DM dm;
6113: PetscInt pdim; /* Dimension of FE space P */
6114: PetscInt spdim; /* Dimension of subelement FE space P */
6115: PetscInt dim; /* Spatial dimension */
6116: PetscInt comp; /* Field components */
6117: PetscInt *subpoints;
6118: PetscReal *tmpB, *tmpD, *tmpH, *subpoint;
6119: PetscInt p, s, d, e, j, k;
6120: PetscErrorCode ierr;
6123: PetscDualSpaceGetDM(fem->dualSpace, &dm);
6124: DMGetDimension(dm, &dim);
6125: PetscSpaceGetDimension(fem->basisSpace, &spdim);
6126: PetscDualSpaceGetDimension(fem->dualSpace, &pdim);
6127: PetscFEGetNumComponents(fem, &comp);
6128: /* Divide points into subelements */
6129: DMGetWorkArray(dm, npoints, MPIU_INT, &subpoints);
6130: DMGetWorkArray(dm, dim, MPIU_REAL, &subpoint);
6131: for (p = 0; p < npoints; ++p) {
6132: for (s = 0; s < cmp->numSubelements; ++s) {
6133: PetscBool inside;
6135: /* Apply transform, and check that point is inside cell */
6136: for (d = 0; d < dim; ++d) {
6137: subpoint[d] = -1.0;
6138: for (e = 0; e < dim; ++e) subpoint[d] += cmp->invjac[(s*dim + d)*dim+e]*(points[p*dim+e] - cmp->v0[s*dim+e]);
6139: }
6140: CellRefinerInCellTest_Internal(cmp->cellRefiner, subpoint, &inside);
6141: if (inside) {subpoints[p] = s; break;}
6142: }
6143: if (s >= cmp->numSubelements) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Point %d was not found in any subelement", p);
6144: }
6145: DMRestoreWorkArray(dm, dim, MPIU_REAL, &subpoint);
6146: /* Evaluate the prime basis functions at all points */
6147: if (B) {DMGetWorkArray(dm, npoints*spdim, MPIU_REAL, &tmpB);}
6148: if (D) {DMGetWorkArray(dm, npoints*spdim*dim, MPIU_REAL, &tmpD);}
6149: if (H) {DMGetWorkArray(dm, npoints*spdim*dim*dim, MPIU_REAL, &tmpH);}
6150: PetscSpaceEvaluate(fem->basisSpace, npoints, points, B ? tmpB : NULL, D ? tmpD : NULL, H ? tmpH : NULL);
6151: /* Translate to the nodal basis */
6152: if (B) {PetscMemzero(B, npoints*pdim*comp * sizeof(PetscReal));}
6153: if (D) {PetscMemzero(D, npoints*pdim*comp*dim * sizeof(PetscReal));}
6154: if (H) {PetscMemzero(H, npoints*pdim*comp*dim*dim * sizeof(PetscReal));}
6155: for (p = 0; p < npoints; ++p) {
6156: const PetscInt s = subpoints[p];
6158: if (B) {
6159: /* Multiply by V^{-1} (spdim x spdim) */
6160: for (j = 0; j < spdim; ++j) {
6161: const PetscInt i = (p*pdim + cmp->embedding[s*spdim+j])*comp;
6163: B[i] = 0.0;
6164: for (k = 0; k < spdim; ++k) {
6165: B[i] += fem->invV[(s*spdim + k)*spdim+j] * tmpB[p*spdim + k];
6166: }
6167: }
6168: }
6169: if (D) {
6170: /* Multiply by V^{-1} (spdim x spdim) */
6171: for (j = 0; j < spdim; ++j) {
6172: for (d = 0; d < dim; ++d) {
6173: const PetscInt i = ((p*pdim + cmp->embedding[s*spdim+j])*comp + 0)*dim + d;
6175: D[i] = 0.0;
6176: for (k = 0; k < spdim; ++k) {
6177: D[i] += fem->invV[(s*spdim + k)*spdim+j] * tmpD[(p*spdim + k)*dim + d];
6178: }
6179: }
6180: }
6181: }
6182: if (H) {
6183: /* Multiply by V^{-1} (pdim x pdim) */
6184: for (j = 0; j < spdim; ++j) {
6185: for (d = 0; d < dim*dim; ++d) {
6186: const PetscInt i = ((p*pdim + cmp->embedding[s*spdim+j])*comp + 0)*dim*dim + d;
6188: H[i] = 0.0;
6189: for (k = 0; k < spdim; ++k) {
6190: H[i] += fem->invV[(s*spdim + k)*spdim+j] * tmpH[(p*spdim + k)*dim*dim + d];
6191: }
6192: }
6193: }
6194: }
6195: }
6196: DMRestoreWorkArray(dm, npoints, MPIU_INT, &subpoints);
6197: if (B) {DMRestoreWorkArray(dm, npoints*spdim, MPIU_REAL, &tmpB);}
6198: if (D) {DMRestoreWorkArray(dm, npoints*spdim*dim, MPIU_REAL, &tmpD);}
6199: if (H) {DMRestoreWorkArray(dm, npoints*spdim*dim*dim, MPIU_REAL, &tmpH);}
6200: return(0);
6201: }
6203: PetscErrorCode PetscFEInitialize_Composite(PetscFE fem)
6204: {
6206: fem->ops->setfromoptions = NULL;
6207: fem->ops->setup = PetscFESetUp_Composite;
6208: fem->ops->view = NULL;
6209: fem->ops->destroy = PetscFEDestroy_Composite;
6210: fem->ops->getdimension = PetscFEGetDimension_Basic;
6211: fem->ops->gettabulation = PetscFEGetTabulation_Composite;
6212: fem->ops->integrateresidual = PetscFEIntegrateResidual_Basic;
6213: fem->ops->integratebdresidual = PetscFEIntegrateBdResidual_Basic;
6214: fem->ops->integratejacobianaction = NULL/* PetscFEIntegrateJacobianAction_Basic */;
6215: fem->ops->integratejacobian = PetscFEIntegrateJacobian_Basic;
6216: return(0);
6217: }
6219: /*MC
6220: PETSCFECOMPOSITE = "composite" - A PetscFE object that represents a composite element
6222: Level: intermediate
6224: .seealso: PetscFEType, PetscFECreate(), PetscFESetType()
6225: M*/
6227: PETSC_EXTERN PetscErrorCode PetscFECreate_Composite(PetscFE fem)
6228: {
6229: PetscFE_Composite *cmp;
6230: PetscErrorCode ierr;
6234: PetscNewLog(fem, &cmp);
6235: fem->data = cmp;
6237: cmp->cellRefiner = REFINER_NOOP;
6238: cmp->numSubelements = -1;
6239: cmp->v0 = NULL;
6240: cmp->jac = NULL;
6242: PetscFEInitialize_Composite(fem);
6243: return(0);
6244: }
6246: /*@C
6247: PetscFECompositeGetMapping - Returns the mappings from the reference element to each subelement
6249: Not collective
6251: Input Parameter:
6252: . fem - The PetscFE object
6254: Output Parameters:
6255: + blockSize - The number of elements in a block
6256: . numBlocks - The number of blocks in a batch
6257: . batchSize - The number of elements in a batch
6258: - numBatches - The number of batches in a chunk
6260: Level: intermediate
6262: .seealso: PetscFECreate()
6263: @*/
6264: PetscErrorCode PetscFECompositeGetMapping(PetscFE fem, PetscInt *numSubelements, const PetscReal *v0[], const PetscReal *jac[], const PetscReal *invjac[])
6265: {
6266: PetscFE_Composite *cmp = (PetscFE_Composite *) fem->data;
6274: return(0);
6275: }
6277: /*@
6278: PetscFEGetDimension - Get the dimension of the finite element space on a cell
6280: Not collective
6282: Input Parameter:
6283: . fe - The PetscFE
6285: Output Parameter:
6286: . dim - The dimension
6288: Level: intermediate
6290: .seealso: PetscFECreate(), PetscSpaceGetDimension(), PetscDualSpaceGetDimension()
6291: @*/
6292: PetscErrorCode PetscFEGetDimension(PetscFE fem, PetscInt *dim)
6293: {
6299: if (fem->ops->getdimension) {(*fem->ops->getdimension)(fem, dim);}
6300: return(0);
6301: }
6303: /*
6304: Purpose: Compute element vector for chunk of elements
6306: Input:
6307: Sizes:
6308: Ne: number of elements
6309: Nf: number of fields
6310: PetscFE
6311: dim: spatial dimension
6312: Nb: number of basis functions
6313: Nc: number of field components
6314: PetscQuadrature
6315: Nq: number of quadrature points
6317: Geometry:
6318: PetscFECellGeom[Ne] possibly *Nq
6319: PetscReal v0s[dim]
6320: PetscReal n[dim]
6321: PetscReal jacobians[dim*dim]
6322: PetscReal jacobianInverses[dim*dim]
6323: PetscReal jacobianDeterminants
6324: FEM:
6325: PetscFE
6326: PetscQuadrature
6327: PetscReal quadPoints[Nq*dim]
6328: PetscReal quadWeights[Nq]
6329: PetscReal basis[Nq*Nb*Nc]
6330: PetscReal basisDer[Nq*Nb*Nc*dim]
6331: PetscScalar coefficients[Ne*Nb*Nc]
6332: PetscScalar elemVec[Ne*Nb*Nc]
6334: Problem:
6335: PetscInt f: the active field
6336: f0, f1
6338: Work Space:
6339: PetscFE
6340: PetscScalar f0[Nq*dim];
6341: PetscScalar f1[Nq*dim*dim];
6342: PetscScalar u[Nc];
6343: PetscScalar gradU[Nc*dim];
6344: PetscReal x[dim];
6345: PetscScalar realSpaceDer[dim];
6347: Purpose: Compute element vector for N_cb batches of elements
6349: Input:
6350: Sizes:
6351: N_cb: Number of serial cell batches
6353: Geometry:
6354: PetscReal v0s[Ne*dim]
6355: PetscReal jacobians[Ne*dim*dim] possibly *Nq
6356: PetscReal jacobianInverses[Ne*dim*dim] possibly *Nq
6357: PetscReal jacobianDeterminants[Ne] possibly *Nq
6358: FEM:
6359: static PetscReal quadPoints[Nq*dim]
6360: static PetscReal quadWeights[Nq]
6361: static PetscReal basis[Nq*Nb*Nc]
6362: static PetscReal basisDer[Nq*Nb*Nc*dim]
6363: PetscScalar coefficients[Ne*Nb*Nc]
6364: PetscScalar elemVec[Ne*Nb*Nc]
6366: ex62.c:
6367: PetscErrorCode PetscFEIntegrateResidualBatch(PetscInt Ne, PetscInt numFields, PetscInt field, PetscQuadrature quad[], const PetscScalar coefficients[],
6368: const PetscReal v0s[], const PetscReal jacobians[], const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[],
6369: void (*f0_func)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f0[]),
6370: void (*f1_func)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f1[]), PetscScalar elemVec[])
6372: ex52.c:
6373: PetscErrorCode IntegrateLaplacianBatchCPU(PetscInt Ne, PetscInt Nb, const PetscScalar coefficients[], const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscInt Nq, const PetscReal quadPoints[], const PetscReal quadWeights[], const PetscReal basisTabulation[], const PetscReal basisDerTabulation[], PetscScalar elemVec[], AppCtx *user)
6374: PetscErrorCode IntegrateElasticityBatchCPU(PetscInt Ne, PetscInt Nb, PetscInt Ncomp, const PetscScalar coefficients[], const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscInt Nq, const PetscReal quadPoints[], const PetscReal quadWeights[], const PetscReal basisTabulation[], const PetscReal basisDerTabulation[], PetscScalar elemVec[], AppCtx *user)
6376: ex52_integrateElement.cu
6377: __global__ void integrateElementQuadrature(int N_cb, realType *coefficients, realType *jacobianInverses, realType *jacobianDeterminants, realType *elemVec)
6379: PETSC_EXTERN PetscErrorCode IntegrateElementBatchGPU(PetscInt spatial_dim, PetscInt Ne, PetscInt Ncb, PetscInt Nbc, PetscInt Nbl, const PetscScalar coefficients[],
6380: const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscScalar elemVec[],
6381: PetscLogEvent event, PetscInt debug, PetscInt pde_op)
6383: ex52_integrateElementOpenCL.c:
6384: PETSC_EXTERN PetscErrorCode IntegrateElementBatchGPU(PetscInt spatial_dim, PetscInt Ne, PetscInt Ncb, PetscInt Nbc, PetscInt N_bl, const PetscScalar coefficients[],
6385: const PetscReal jacobianInverses[], const PetscReal jacobianDeterminants[], PetscScalar elemVec[],
6386: PetscLogEvent event, PetscInt debug, PetscInt pde_op)
6388: __kernel void integrateElementQuadrature(int N_cb, __global float *coefficients, __global float *jacobianInverses, __global float *jacobianDeterminants, __global float *elemVec)
6389: */
6391: /*@C
6392: PetscFEIntegrate - Produce the integral for the given field for a chunk of elements by quadrature integration
6394: Not collective
6396: Input Parameters:
6397: + fem - The PetscFE object for the field being integrated
6398: . prob - The PetscDS specifying the discretizations and continuum functions
6399: . field - The field being integrated
6400: . Ne - The number of elements in the chunk
6401: . cgeom - The cell geometry for each cell in the chunk
6402: . coefficients - The array of FEM basis coefficients for the elements
6403: . probAux - The PetscDS specifying the auxiliary discretizations
6404: - coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
6406: Output Parameter
6407: . integral - the integral for this field
6409: Level: developer
6411: .seealso: PetscFEIntegrateResidual()
6412: @*/
6413: PetscErrorCode PetscFEIntegrate(PetscFE fem, PetscDS prob, PetscInt field, PetscInt Ne, PetscFECellGeom *cgeom,
6414: const PetscScalar coefficients[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscScalar integral[])
6415: {
6421: if (fem->ops->integrate) {(*fem->ops->integrate)(fem, prob, field, Ne, cgeom, coefficients, probAux, coefficientsAux, integral);}
6422: return(0);
6423: }
6425: /*@C
6426: PetscFEIntegrateResidual - Produce the element residual vector for a chunk of elements by quadrature integration
6428: Not collective
6430: Input Parameters:
6431: + fem - The PetscFE object for the field being integrated
6432: . prob - The PetscDS specifying the discretizations and continuum functions
6433: . field - The field being integrated
6434: . Ne - The number of elements in the chunk
6435: . cgeom - The cell geometry for each cell in the chunk
6436: . coefficients - The array of FEM basis coefficients for the elements
6437: . coefficients_t - The array of FEM basis time derivative coefficients for the elements
6438: . probAux - The PetscDS specifying the auxiliary discretizations
6439: . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
6440: - t - The time
6442: Output Parameter
6443: . elemVec - the element residual vectors from each element
6445: Note:
6446: $ Loop over batch of elements (e):
6447: $ Loop over quadrature points (q):
6448: $ Make u_q and gradU_q (loops over fields,Nb,Ncomp) and x_q
6449: $ Call f_0 and f_1
6450: $ Loop over element vector entries (f,fc --> i):
6451: $ elemVec[i] += \psi^{fc}_f(q) f0_{fc}(u, \nabla u) + \nabla\psi^{fc}_f(q) \cdot f1_{fc,df}(u, \nabla u)
6453: Level: developer
6455: .seealso: PetscFEIntegrateResidual()
6456: @*/
6457: PetscErrorCode PetscFEIntegrateResidual(PetscFE fem, PetscDS prob, PetscInt field, PetscInt Ne, PetscFECellGeom *cgeom,
6458: const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
6459: {
6465: if (fem->ops->integrateresidual) {(*fem->ops->integrateresidual)(fem, prob, field, Ne, cgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec);}
6466: return(0);
6467: }
6469: /*@C
6470: PetscFEIntegrateBdResidual - Produce the element residual vector for a chunk of elements by quadrature integration over a boundary
6472: Not collective
6474: Input Parameters:
6475: + fem - The PetscFE object for the field being integrated
6476: . prob - The PetscDS specifying the discretizations and continuum functions
6477: . field - The field being integrated
6478: . Ne - The number of elements in the chunk
6479: . fgeom - The face geometry for each cell in the chunk
6480: . coefficients - The array of FEM basis coefficients for the elements
6481: . coefficients_t - The array of FEM basis time derivative coefficients for the elements
6482: . probAux - The PetscDS specifying the auxiliary discretizations
6483: . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
6484: - t - The time
6486: Output Parameter
6487: . elemVec - the element residual vectors from each element
6489: Level: developer
6491: .seealso: PetscFEIntegrateResidual()
6492: @*/
6493: PetscErrorCode PetscFEIntegrateBdResidual(PetscFE fem, PetscDS prob, PetscInt field, PetscInt Ne, PetscFEFaceGeom *fgeom,
6494: const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscScalar elemVec[])
6495: {
6500: if (fem->ops->integratebdresidual) {(*fem->ops->integratebdresidual)(fem, prob, field, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, elemVec);}
6501: return(0);
6502: }
6504: /*@C
6505: PetscFEIntegrateJacobian - Produce the element Jacobian for a chunk of elements by quadrature integration
6507: Not collective
6509: Input Parameters:
6510: + fem - The PetscFE object for the field being integrated
6511: . prob - The PetscDS specifying the discretizations and continuum functions
6512: . jtype - The type of matrix pointwise functions that should be used
6513: . fieldI - The test field being integrated
6514: . fieldJ - The basis field being integrated
6515: . Ne - The number of elements in the chunk
6516: . cgeom - The cell geometry for each cell in the chunk
6517: . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point
6518: . coefficients_t - The array of FEM basis time derivative coefficients for the elements
6519: . probAux - The PetscDS specifying the auxiliary discretizations
6520: . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
6521: . t - The time
6522: - u_tShift - A multiplier for the dF/du_t term (as opposed to the dF/du term)
6524: Output Parameter
6525: . elemMat - the element matrices for the Jacobian from each element
6527: Note:
6528: $ Loop over batch of elements (e):
6529: $ Loop over element matrix entries (f,fc,g,gc --> i,j):
6530: $ Loop over quadrature points (q):
6531: $ Make u_q and gradU_q (loops over fields,Nb,Ncomp)
6532: $ elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q)
6533: $ + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
6534: $ + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q)
6535: $ + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
6536: Level: developer
6538: .seealso: PetscFEIntegrateResidual()
6539: @*/
6540: PetscErrorCode PetscFEIntegrateJacobian(PetscFE fem, PetscDS prob, PetscFEJacobianType jtype, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFECellGeom *cgeom,
6541: const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
6542: {
6547: if (fem->ops->integratejacobian) {(*fem->ops->integratejacobian)(fem, prob, jtype, fieldI, fieldJ, Ne, cgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, u_tshift, elemMat);}
6548: return(0);
6549: }
6551: /*@C
6552: PetscFEIntegrateBdJacobian - Produce the boundary element Jacobian for a chunk of elements by quadrature integration
6554: Not collective
6556: Input Parameters:
6557: + fem = The PetscFE object for the field being integrated
6558: . prob - The PetscDS specifying the discretizations and continuum functions
6559: . fieldI - The test field being integrated
6560: . fieldJ - The basis field being integrated
6561: . Ne - The number of elements in the chunk
6562: . fgeom - The face geometry for each cell in the chunk
6563: . coefficients - The array of FEM basis coefficients for the elements for the Jacobian evaluation point
6564: . coefficients_t - The array of FEM basis time derivative coefficients for the elements
6565: . probAux - The PetscDS specifying the auxiliary discretizations
6566: . coefficientsAux - The array of FEM auxiliary basis coefficients for the elements
6567: . t - The time
6568: - u_tShift - A multiplier for the dF/du_t term (as opposed to the dF/du term)
6570: Output Parameter
6571: . elemMat - the element matrices for the Jacobian from each element
6573: Note:
6574: $ Loop over batch of elements (e):
6575: $ Loop over element matrix entries (f,fc,g,gc --> i,j):
6576: $ Loop over quadrature points (q):
6577: $ Make u_q and gradU_q (loops over fields,Nb,Ncomp)
6578: $ elemMat[i,j] += \psi^{fc}_f(q) g0_{fc,gc}(u, \nabla u) \phi^{gc}_g(q)
6579: $ + \psi^{fc}_f(q) \cdot g1_{fc,gc,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
6580: $ + \nabla\psi^{fc}_f(q) \cdot g2_{fc,gc,df}(u, \nabla u) \phi^{gc}_g(q)
6581: $ + \nabla\psi^{fc}_f(q) \cdot g3_{fc,gc,df,dg}(u, \nabla u) \nabla\phi^{gc}_g(q)
6582: Level: developer
6584: .seealso: PetscFEIntegrateJacobian(), PetscFEIntegrateResidual()
6585: @*/
6586: PetscErrorCode PetscFEIntegrateBdJacobian(PetscFE fem, PetscDS prob, PetscInt fieldI, PetscInt fieldJ, PetscInt Ne, PetscFEFaceGeom *fgeom,
6587: const PetscScalar coefficients[], const PetscScalar coefficients_t[], PetscDS probAux, const PetscScalar coefficientsAux[], PetscReal t, PetscReal u_tshift, PetscScalar elemMat[])
6588: {
6593: if (fem->ops->integratebdjacobian) {(*fem->ops->integratebdjacobian)(fem, prob, fieldI, fieldJ, Ne, fgeom, coefficients, coefficients_t, probAux, coefficientsAux, t, u_tshift, elemMat);}
6594: return(0);
6595: }
6597: PetscErrorCode PetscFEGetHeightSubspace(PetscFE fe, PetscInt height, PetscFE *subfe)
6598: {
6599: PetscSpace P, subP;
6600: PetscDualSpace Q, subQ;
6601: PetscQuadrature subq;
6602: PetscFEType fetype;
6603: PetscInt dim, Nc;
6604: PetscErrorCode ierr;
6609: if (height == 0) {
6610: *subfe = fe;
6611: return(0);
6612: }
6613: PetscFEGetBasisSpace(fe, &P);
6614: PetscFEGetDualSpace(fe, &Q);
6615: PetscFEGetNumComponents(fe, &Nc);
6616: PetscFEGetFaceQuadrature(fe, &subq);
6617: PetscDualSpaceGetDimension(Q, &dim);
6618: if (height > dim || height < 0) {SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Asked for space at height %D for dimension %D space", height, dim);}
6619: if (!fe->subspaces) {PetscCalloc1(dim, &fe->subspaces);}
6620: if (height <= dim) {
6621: if (!fe->subspaces[height-1]) {
6622: PetscFE sub;
6624: PetscSpaceGetHeightSubspace(P, height, &subP);
6625: PetscDualSpaceGetHeightSubspace(Q, height, &subQ);
6626: PetscFECreate(PetscObjectComm((PetscObject) fe), &sub);
6627: PetscFEGetType(fe, &fetype);
6628: PetscFESetType(sub, fetype);
6629: PetscFESetBasisSpace(sub, subP);
6630: PetscFESetDualSpace(sub, subQ);
6631: PetscFESetNumComponents(sub, Nc);
6632: PetscFESetUp(sub);
6633: PetscFESetQuadrature(sub, subq);
6634: fe->subspaces[height-1] = sub;
6635: }
6636: *subfe = fe->subspaces[height-1];
6637: } else {
6638: *subfe = NULL;
6639: }
6640: return(0);
6641: }
6643: /*@
6644: PetscFERefine - Create a "refined" PetscFE object that refines the reference cell into smaller copies. This is typically used
6645: to precondition a higher order method with a lower order method on a refined mesh having the same number of dofs (but more
6646: sparsity). It is also used to create an interpolation between regularly refined meshes.
6648: Collective on PetscFE
6650: Input Parameter:
6651: . fe - The initial PetscFE
6653: Output Parameter:
6654: . feRef - The refined PetscFE
6656: Level: developer
6658: .seealso: PetscFEType, PetscFECreate(), PetscFESetType()
6659: @*/
6660: PetscErrorCode PetscFERefine(PetscFE fe, PetscFE *feRef)
6661: {
6662: PetscSpace P, Pref;
6663: PetscDualSpace Q, Qref;
6664: DM K, Kref;
6665: PetscQuadrature q, qref;
6666: const PetscReal *v0, *jac;
6667: PetscInt numComp, numSubelements;
6668: PetscErrorCode ierr;
6671: PetscFEGetBasisSpace(fe, &P);
6672: PetscFEGetDualSpace(fe, &Q);
6673: PetscFEGetQuadrature(fe, &q);
6674: PetscDualSpaceGetDM(Q, &K);
6675: /* Create space */
6676: PetscObjectReference((PetscObject) P);
6677: Pref = P;
6678: /* Create dual space */
6679: PetscDualSpaceDuplicate(Q, &Qref);
6680: DMRefine(K, PetscObjectComm((PetscObject) fe), &Kref);
6681: PetscDualSpaceSetDM(Qref, Kref);
6682: DMDestroy(&Kref);
6683: PetscDualSpaceSetUp(Qref);
6684: /* Create element */
6685: PetscFECreate(PetscObjectComm((PetscObject) fe), feRef);
6686: PetscFESetType(*feRef, PETSCFECOMPOSITE);
6687: PetscFESetBasisSpace(*feRef, Pref);
6688: PetscFESetDualSpace(*feRef, Qref);
6689: PetscFEGetNumComponents(fe, &numComp);
6690: PetscFESetNumComponents(*feRef, numComp);
6691: PetscFESetUp(*feRef);
6692: PetscSpaceDestroy(&Pref);
6693: PetscDualSpaceDestroy(&Qref);
6694: /* Create quadrature */
6695: PetscFECompositeGetMapping(*feRef, &numSubelements, &v0, &jac, NULL);
6696: PetscQuadratureExpandComposite(q, numSubelements, v0, jac, &qref);
6697: PetscFESetQuadrature(*feRef, qref);
6698: PetscQuadratureDestroy(&qref);
6699: return(0);
6700: }
6702: /*@C
6703: PetscFECreateDefault - Create a PetscFE for basic FEM computation
6705: Collective on DM
6707: Input Parameters:
6708: + dm - The underlying DM for the domain
6709: . dim - The spatial dimension
6710: . Nc - The number of components
6711: . isSimplex - Flag for simplex reference cell, otherwise its a tensor product
6712: . prefix - The options prefix, or NULL
6713: - qorder - The quadrature order
6715: Output Parameter:
6716: . fem - The PetscFE object
6718: Level: beginner
6720: .keywords: PetscFE, finite element
6721: .seealso: PetscFECreate(), PetscSpaceCreate(), PetscDualSpaceCreate()
6722: @*/
6723: PetscErrorCode PetscFECreateDefault(DM dm, PetscInt dim, PetscInt Nc, PetscBool isSimplex, const char prefix[], PetscInt qorder, PetscFE *fem)
6724: {
6725: PetscQuadrature q, fq;
6726: DM K;
6727: PetscSpace P;
6728: PetscDualSpace Q;
6729: PetscInt order, quadPointsPerEdge;
6730: PetscBool tensor = isSimplex ? PETSC_FALSE : PETSC_TRUE;
6731: PetscErrorCode ierr;
6734: /* Create space */
6735: PetscSpaceCreate(PetscObjectComm((PetscObject) dm), &P);
6736: PetscObjectSetOptionsPrefix((PetscObject) P, prefix);
6737: PetscSpacePolynomialSetTensor(P, tensor);
6738: PetscSpaceSetFromOptions(P);
6739: PetscSpaceSetNumComponents(P, Nc);
6740: PetscSpacePolynomialSetNumVariables(P, dim);
6741: PetscSpaceSetUp(P);
6742: PetscSpaceGetOrder(P, &order);
6743: PetscSpacePolynomialGetTensor(P, &tensor);
6744: /* Create dual space */
6745: PetscDualSpaceCreate(PetscObjectComm((PetscObject) dm), &Q);
6746: PetscDualSpaceSetType(Q,PETSCDUALSPACELAGRANGE);
6747: PetscObjectSetOptionsPrefix((PetscObject) Q, prefix);
6748: PetscDualSpaceCreateReferenceCell(Q, dim, isSimplex, &K);
6749: PetscDualSpaceSetDM(Q, K);
6750: DMDestroy(&K);
6751: PetscDualSpaceSetNumComponents(Q, Nc);
6752: PetscDualSpaceSetOrder(Q, order);
6753: PetscDualSpaceLagrangeSetTensor(Q, tensor);
6754: PetscDualSpaceSetFromOptions(Q);
6755: PetscDualSpaceSetUp(Q);
6756: /* Create element */
6757: PetscFECreate(PetscObjectComm((PetscObject) dm), fem);
6758: PetscObjectSetOptionsPrefix((PetscObject) *fem, prefix);
6759: PetscFESetFromOptions(*fem);
6760: PetscFESetBasisSpace(*fem, P);
6761: PetscFESetDualSpace(*fem, Q);
6762: PetscFESetNumComponents(*fem, Nc);
6763: PetscFESetUp(*fem);
6764: PetscSpaceDestroy(&P);
6765: PetscDualSpaceDestroy(&Q);
6766: /* Create quadrature (with specified order if given) */
6767: qorder = qorder >= 0 ? qorder : order;
6768: PetscObjectOptionsBegin((PetscObject)*fem);
6769: PetscOptionsInt("-petscfe_default_quadrature_order","Quadrature order is one less than quadture points per edge","PetscFECreateDefault",qorder,&qorder,NULL);
6770: PetscOptionsEnd();
6771: quadPointsPerEdge = PetscMax(qorder + 1,1);
6772: if (isSimplex) {
6773: PetscDTGaussJacobiQuadrature(dim, 1, quadPointsPerEdge, -1.0, 1.0, &q);
6774: PetscDTGaussJacobiQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, &fq);
6775: }
6776: else {
6777: PetscDTGaussTensorQuadrature(dim, 1, quadPointsPerEdge, -1.0, 1.0, &q);
6778: PetscDTGaussTensorQuadrature(dim-1, 1, quadPointsPerEdge, -1.0, 1.0, &fq);
6779: }
6780: PetscFESetQuadrature(*fem, q);
6781: PetscFESetFaceQuadrature(*fem, fq);
6782: PetscQuadratureDestroy(&q);
6783: PetscQuadratureDestroy(&fq);
6784: return(0);
6785: }