Actual source code: spacetensor.c
petsc-3.10.5 2019-03-28
1: #include <petsc/private/petscfeimpl.h>
3: static PetscErrorCode PetscSpaceTensorCreateSubspace(PetscSpace space, PetscInt Nvs, PetscSpace *subspace)
4: {
5: PetscInt degree;
6: const char *prefix;
10: PetscSpaceGetDegree(space, °ree, NULL);
11: PetscObjectGetOptionsPrefix((PetscObject)space, &prefix);
12: PetscSpaceCreate(PetscObjectComm((PetscObject)space), subspace);
13: PetscSpaceSetType(*subspace, PETSCSPACEPOLYNOMIAL);
14: PetscSpaceSetNumVariables(*subspace, Nvs);
15: PetscSpaceSetNumComponents(*subspace, 1);
16: PetscSpaceSetDegree(*subspace, degree, PETSC_DETERMINE);
17: PetscObjectSetOptionsPrefix((PetscObject)*subspace, prefix);
18: PetscObjectAppendOptionsPrefix((PetscObject)*subspace, "subspace_");
19: return(0);
20: }
22: static PetscErrorCode PetscSpaceSetFromOptions_Tensor(PetscOptionItems *PetscOptionsObject,PetscSpace sp)
23: {
24: PetscSpace_Tensor *tens = (PetscSpace_Tensor *) sp->data;
25: PetscInt Ns, Nc, i, Nv, deg;
26: PetscBool uniform = PETSC_TRUE;
27: PetscErrorCode ierr;
30: PetscSpaceGetNumVariables(sp, &Nv);
31: if (!Nv) return(0);
32: PetscSpaceGetNumComponents(sp, &Nc);
33: PetscSpaceTensorGetNumSubspaces(sp, &Ns);
34: PetscSpaceGetDegree(sp, °, NULL);
35: if (Ns > 1) {
36: PetscSpace s0;
38: PetscSpaceTensorGetSubspace(sp, 0, &s0);
39: for (i = 1; i < Ns; i++) {
40: PetscSpace si;
42: PetscSpaceTensorGetSubspace(sp, i, &si);
43: if (si != s0) {uniform = PETSC_FALSE; break;}
44: }
45: }
46: Ns = (Ns == PETSC_DEFAULT) ? PetscMax(Nv,1) : Ns;
47: PetscOptionsHead(PetscOptionsObject,"PetscSpace tensor options");
48: PetscOptionsInt("-petscspace_tensor_spaces", "The number of subspaces", "PetscSpaceTensorSetNumSubspaces", Ns, &Ns, NULL);
49: PetscOptionsBool("-petscspace_tensor_uniform", "Subspaces are identical", "PetscSpaceTensorSetFromOptions", uniform, &uniform, NULL);
50: PetscOptionsTail();
51: if (Ns < 0 || (Nv > 0 && Ns == 0)) SETERRQ1(PetscObjectComm((PetscObject)sp),PETSC_ERR_ARG_OUTOFRANGE,"Cannot have a tensor space made up of %D spaces\n",Ns);
52: if (Nv > 0 && Ns > Nv) SETERRQ2(PetscObjectComm((PetscObject)sp),PETSC_ERR_ARG_OUTOFRANGE,"Cannot have a tensor space with %D subspaces over %D variables\n", Ns, Nv);
53: if (Ns != tens->numTensSpaces) {PetscSpaceTensorSetNumSubspaces(sp, Ns);}
54: if (uniform) {
55: PetscInt Nvs = Nv / Ns;
56: PetscSpace subspace;
58: if (Nv % Ns) SETERRQ2(PetscObjectComm((PetscObject)sp),PETSC_ERR_ARG_WRONG,"Cannot use %D uniform subspaces for %D variable space\n", Ns, Nv);
59: PetscSpaceTensorGetSubspace(sp, 0, &subspace);
60: if (!subspace) {PetscSpaceTensorCreateSubspace(sp, Nvs, &subspace);}
61: else {PetscObjectReference((PetscObject)subspace);}
62: PetscSpaceSetFromOptions(subspace);
63: for (i = 0; i < Ns; i++) {PetscSpaceTensorSetSubspace(sp, i, subspace);}
64: PetscSpaceDestroy(&subspace);
65: } else {
66: for (i = 0; i < Ns; i++) {
67: PetscSpace subspace;
69: PetscSpaceTensorGetSubspace(sp, i, &subspace);
70: if (!subspace) {
71: char tprefix[128];
73: PetscSpaceTensorCreateSubspace(sp, 1, &subspace);
74: PetscSNPrintf(tprefix, 128, "%d_",(int)i);
75: PetscObjectAppendOptionsPrefix((PetscObject)subspace, tprefix);
76: } else {
77: PetscObjectReference((PetscObject)subspace);
78: }
79: PetscSpaceSetFromOptions(subspace);
80: PetscSpaceTensorSetSubspace(sp, i, subspace);
81: PetscSpaceDestroy(&subspace);
82: }
83: }
84: return(0);
85: }
87: static PetscErrorCode PetscSpaceTensorView_Ascii(PetscSpace sp, PetscViewer viewer)
88: {
89: PetscSpace_Tensor *tens = (PetscSpace_Tensor *) sp->data;
90: PetscBool uniform = PETSC_TRUE;
91: PetscInt Ns = tens->numTensSpaces, i, n;
92: PetscErrorCode ierr;
95: for (i = 1; i < Ns; i++) {
96: if (tens->tensspaces[i] != tens->tensspaces[0]) {uniform = PETSC_FALSE; break;}
97: }
98: if (uniform) {PetscViewerASCIIPrintf(viewer, "Tensor space of %D subspaces (all identical)\n", Ns);
99: } else {PetscViewerASCIIPrintf(viewer, "Tensor space of %D subspaces\n", Ns);}
100: n = uniform ? 1 : Ns;
101: PetscViewerASCIIPushTab(viewer);
102: for (i = 0; i < n; i++) {
103: PetscSpaceView(tens->tensspaces[i], viewer);
104: }
105: PetscViewerASCIIPushTab(viewer);
106: return(0);
107: }
109: static PetscErrorCode PetscSpaceView_Tensor(PetscSpace sp, PetscViewer viewer)
110: {
111: PetscBool iascii;
117: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
118: if (iascii) {PetscSpaceTensorView_Ascii(sp, viewer);}
119: return(0);
120: }
122: static PetscErrorCode PetscSpaceSetUp_Tensor(PetscSpace sp)
123: {
124: PetscSpace_Tensor *tens = (PetscSpace_Tensor *) sp->data;
125: PetscInt Nv, Ns, i;
126: PetscBool uniform = PETSC_TRUE;
127: PetscInt deg, maxDeg;
128: PetscErrorCode ierr;
131: if (tens->setupCalled) return(0);
132: PetscSpaceGetNumVariables(sp, &Nv);
133: PetscSpaceTensorGetNumSubspaces(sp, &Ns);
134: if (Ns == PETSC_DEFAULT) {
135: Ns = Nv;
136: PetscSpaceTensorSetNumSubspaces(sp, Ns);
137: }
138: if (!Ns) {
139: if (Nv) SETERRQ(PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_OUTOFRANGE, "Cannot have zero subspaces");
140: } else {
141: PetscSpace s0;
143: if (Nv > 0 && Ns > Nv) SETERRQ2(PetscObjectComm((PetscObject)sp),PETSC_ERR_ARG_OUTOFRANGE,"Cannot have a tensor space with %D subspaces over %D variables\n", Ns, Nv);
144: PetscSpaceTensorGetSubspace(sp, 0, &s0);
145: for (i = 1; i < Ns; i++) {
146: PetscSpace si;
148: PetscSpaceTensorGetSubspace(sp, i, &si);
149: if (si != s0) {uniform = PETSC_FALSE; break;}
150: }
151: if (uniform) {
152: PetscInt Nvs = Nv / Ns;
154: if (Nv % Ns) SETERRQ2(PetscObjectComm((PetscObject)sp),PETSC_ERR_ARG_WRONG,"Cannot use %D uniform subspaces for %D variable space\n", Ns, Nv);
155: if (!s0) {PetscSpaceTensorCreateSubspace(sp, Nvs, &s0);}
156: else {PetscObjectReference((PetscObject) s0);}
157: PetscSpaceSetUp(s0);
158: for (i = 0; i < Ns; i++) {PetscSpaceTensorSetSubspace(sp, i, s0);}
159: PetscSpaceDestroy(&s0);
160: } else {
161: for (i = 0 ; i < Ns; i++) {
162: PetscSpace si;
164: PetscSpaceTensorGetSubspace(sp, i, &si);
165: if (!si) {PetscSpaceTensorCreateSubspace(sp, 1, &si);}
166: else {PetscObjectReference((PetscObject) si);}
167: PetscSpaceSetUp(si);
168: PetscSpaceTensorSetSubspace(sp, i, si);
169: PetscSpaceDestroy(&si);
170: }
171: }
172: }
173: deg = PETSC_MAX_INT;
174: maxDeg = 0;
175: for (i = 0; i < Ns; i++) {
176: PetscSpace si;
177: PetscInt iDeg, iMaxDeg;
179: PetscSpaceTensorGetSubspace(sp, i, &si);
180: PetscSpaceGetDegree(si, &iDeg, &iMaxDeg);
181: deg = PetscMin(deg, iDeg);
182: maxDeg += iMaxDeg;
183: }
184: sp->degree = deg;
185: sp->maxDegree = maxDeg;
186: tens->uniform = uniform;
187: tens->setupCalled = PETSC_TRUE;
188: return(0);
189: }
191: static PetscErrorCode PetscSpaceDestroy_Tensor(PetscSpace sp)
192: {
193: PetscSpace_Tensor *tens = (PetscSpace_Tensor *) sp->data;
194: PetscInt Ns, i;
195: PetscErrorCode ierr;
198: Ns = tens->numTensSpaces;
199: if (tens->heightsubspaces) {
200: PetscInt d;
202: /* sp->Nv is the spatial dimension, so it is equal to the number
203: * of subspaces on higher co-dimension points */
204: for (d = 0; d < sp->Nv; ++d) {
205: PetscSpaceDestroy(&tens->heightsubspaces[d]);
206: }
207: }
208: PetscFree(tens->heightsubspaces);
209: for (i = 0; i < Ns; i++) {PetscSpaceDestroy(&tens->tensspaces[i]);}
210: PetscObjectComposeFunction((PetscObject) sp, "PetscSpaceTensorSetSubspace_C", NULL);
211: PetscObjectComposeFunction((PetscObject) sp, "PetscSpaceTensorGetSubspace_C", NULL);
212: PetscObjectComposeFunction((PetscObject) sp, "PetscSpaceTensorSetNumSubspaces_C", NULL);
213: PetscObjectComposeFunction((PetscObject) sp, "PetscSpaceTensorGetNumSubspaces_C", NULL);
214: PetscFree(tens->tensspaces);
215: PetscFree(tens);
216: return(0);
217: }
219: static PetscErrorCode PetscSpaceGetDimension_Tensor(PetscSpace sp, PetscInt *dim)
220: {
221: PetscSpace_Tensor *tens = (PetscSpace_Tensor *) sp->data;
222: PetscInt i, Ns, Nc, d;
223: PetscErrorCode ierr;
226: PetscSpaceSetUp(sp);
227: Ns = tens->numTensSpaces;
228: Nc = sp->Nc;
229: d = 1;
230: for (i = 0; i < Ns; i++) {
231: PetscInt id;
233: PetscSpaceGetDimension(tens->tensspaces[i], &id);
234: d *= id;
235: }
236: d *= Nc;
237: *dim = d;
238: return(0);
239: }
241: static PetscErrorCode PetscSpaceEvaluate_Tensor(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[])
242: {
243: PetscSpace_Tensor *tens = (PetscSpace_Tensor *) sp->data;
244: DM dm = sp->dm;
245: PetscInt Nc = sp->Nc;
246: PetscInt Nv = sp->Nv;
247: PetscInt Ns;
248: PetscReal *lpoints, *sB = NULL, *sD = NULL, *sH = NULL;
249: PetscInt c, pdim, d, e, der, der2, i, l, si, p, s, step;
250: PetscErrorCode ierr;
253: if (!tens->setupCalled) {PetscSpaceSetUp(sp);}
254: Ns = tens->numTensSpaces;
255: PetscSpaceGetDimension(sp,&pdim);
256: pdim /= Nc;
257: DMGetWorkArray(dm, npoints*Nv, MPIU_REAL, &lpoints);
258: if (B || D || H) {DMGetWorkArray(dm, npoints*pdim, MPIU_REAL, &sB);}
259: if (D || H) {DMGetWorkArray(dm, npoints*pdim*Nv, MPIU_REAL, &sD);}
260: if (H) {DMGetWorkArray(dm, npoints*pdim*Nv*Nv, MPIU_REAL, &sH);}
261: if (B) {
262: for (i = 0; i < npoints*pdim*Nc*Nc; i++) B[i] = 0.;
263: for (i = 0; i < npoints*pdim; i++) B[i * Nc*Nc] = 1.;
264: }
265: if (D) {
266: for (i = 0; i < npoints*pdim*Nc*Nc*Nv; i++) D[i] = 0.;
267: for (i = 0; i < npoints*pdim; i++) {
268: for (l = 0; l < Nv; l++) {
269: D[i * Nc*Nc*Nv + l] = 1.;
270: }
271: }
272: }
273: if (H) {
274: for (i = 0; i < npoints*pdim*Nc*Nc*Nv*Nv; i++) D[i] = 0.;
275: for (i = 0; i < npoints*pdim; i++) {
276: for (l = 0; l < Nv*Nv; l++) {
277: H[i * Nc*Nc*Nv*Nv + l] = 1.;
278: }
279: }
280: }
281: for (s = 0, d = 0, step = 1; s < Ns; s++) {
282: PetscInt sNv, spdim;
283: PetscInt skip, j, k;
285: PetscSpaceGetNumVariables(tens->tensspaces[s], &sNv);
286: PetscSpaceGetDimension(tens->tensspaces[s], &spdim);
287: if ((pdim % step) || (pdim % spdim)) SETERRQ6(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Bad tensor loop: Nv %d, Ns %D, pdim %D, s %D, step %D, spdim %D", Nv, Ns, pdim, s, step, spdim);
288: skip = pdim / (step * spdim);
289: for (p = 0; p < npoints; ++p) {
290: for (i = 0; i < sNv; i++) {
291: lpoints[p * sNv + i] = points[p*Nv + d + i];
292: }
293: }
294: PetscSpaceEvaluate(tens->tensspaces[s], npoints, lpoints, sB, sD, sH);
295: if (B) {
296: for (p = 0; p < npoints; p++) {
297: for (k = 0; k < skip; k++) {
298: for (si = 0; si < spdim; si++) {
299: for (j = 0; j < step; j++) {
300: i = (k * spdim + si) * step + j;
301: B[(pdim * p + i) * Nc * Nc] *= sB[spdim * p + si];
302: }
303: }
304: }
305: }
306: }
307: if (D) {
308: for (p = 0; p < npoints; p++) {
309: for (k = 0; k < skip; k++) {
310: for (si = 0; si < spdim; si++) {
311: for (j = 0; j < step; j++) {
312: i = (k * spdim + si) * step + j;
313: for (der = 0; der < Nv; der++) {
314: if (der >= d && der < d + sNv) {
315: D[(pdim * p + i) * Nc*Nc*Nv + der] *= sD[(spdim * p + si) * sNv + der - d];
316: } else {
317: D[(pdim * p + i) * Nc*Nc*Nv + der] *= sB[spdim * p + si];
318: }
319: }
320: }
321: }
322: }
323: }
324: }
325: if (H) {
326: for (p = 0; p < npoints; p++) {
327: for (k = 0; k < skip; k++) {
328: for (si = 0; si < spdim; si++) {
329: for (j = 0; j < step; j++) {
330: i = (k * spdim + si) * step + j;
331: for (der = 0; der < Nv; der++) {
332: for (der2 = 0; der2 < Nv; der2++) {
333: if (der >= d && der < d + sNv && der2 >= d && der2 < d + sNv) {
334: H[((pdim * p + i) * Nc*Nc*Nv + der) * Nv + der2] *= sH[((spdim * p + si) * sNv + der - d) * sNv + der2 - d];
335: } else if (der >= d && der < d + sNv) {
336: H[((pdim * p + i) * Nc*Nc*Nv + der) * Nv + der2] *= sD[(spdim * p + si) * sNv + der - d];
337: } else if (der2 >= d && der2 < d + sNv) {
338: H[((pdim * p + i) * Nc*Nc*Nv + der) * Nv + der2] *= sD[(spdim * p + si) * sNv + der2 - d];
339: } else {
340: H[((pdim * p + i) * Nc*Nc*Nv + der) * Nv + der2] *= sB[spdim * p + si];
341: }
342: }
343: }
344: }
345: }
346: }
347: }
348: }
349: d += sNv;
350: step *= spdim;
351: }
352: if (B && Nc > 1) {
353: /* Make direct sum basis for multicomponent space */
354: for (p = 0; p < npoints; ++p) {
355: for (i = 0; i < pdim; ++i) {
356: for (c = 1; c < Nc; ++c) {
357: B[(p*pdim*Nc + i*Nc + c)*Nc + c] = B[(p*pdim + i)*Nc*Nc];
358: }
359: }
360: }
361: }
362: if (D && Nc > 1) {
363: /* Make direct sum basis for multicomponent space */
364: for (p = 0; p < npoints; ++p) {
365: for (i = 0; i < pdim; ++i) {
366: for (c = 1; c < Nc; ++c) {
367: for (d = 0; d < Nv; ++d) {
368: D[((p*pdim*Nc + i*Nc + c)*Nc + c)*Nv + d] = D[(p*pdim + i)*Nc*Nc*Nv + d];
369: }
370: }
371: }
372: }
373: }
374: if (H && Nc > 1) {
375: /* Make direct sum basis for multicomponent space */
376: for (p = 0; p < npoints; ++p) {
377: for (i = 0; i < pdim; ++i) {
378: for (c = 1; c < Nc; ++c) {
379: for (d = 0; d < Nv; ++d) {
380: for (e = 0; e < Nv; ++e) {
381: H[(((p*pdim*Nc + i*Nc + c)*Nc + c)*Nv + d)*Nv + e] = H[((p*pdim + i)*Nc*Nc*Nv + d)*Nv + e];
382: }
383: }
384: }
385: }
386: }
387: }
388: if (H) {DMRestoreWorkArray(dm, npoints*pdim*Nv*Nv, MPIU_REAL, &sH);}
389: if (D || H) {DMRestoreWorkArray(dm, npoints*pdim*Nv, MPIU_REAL, &sD);}
390: if (B || D || H) {DMRestoreWorkArray(dm, npoints*pdim, MPIU_REAL, &sB);}
391: DMRestoreWorkArray(dm, npoints*Nv, MPIU_REAL, &lpoints);
392: return(0);
393: }
395: PetscErrorCode PetscSpaceTensorSetNumSubspaces(PetscSpace sp, PetscInt numTensSpaces)
396: {
401: PetscTryMethod(sp,"PetscSpaceTensorSetNumSubspaces_C",(PetscSpace,PetscInt),(sp,numTensSpaces));
402: return(0);
403: }
405: PetscErrorCode PetscSpaceTensorGetNumSubspaces(PetscSpace sp, PetscInt *numTensSpaces)
406: {
412: PetscTryMethod(sp,"PetscSpaceTensorGetNumSubspaces_C",(PetscSpace,PetscInt*),(sp,numTensSpaces));
413: return(0);
414: }
416: PetscErrorCode PetscSpaceTensorSetSubspace(PetscSpace sp, PetscInt s, PetscSpace subsp)
417: {
423: PetscTryMethod(sp,"PetscSpaceTensorSetSubspace_C",(PetscSpace,PetscInt,PetscSpace),(sp,s,subsp));
424: return(0);
425: }
427: PetscErrorCode PetscSpaceTensorGetSubspace(PetscSpace sp, PetscInt s, PetscSpace *subsp)
428: {
434: PetscTryMethod(sp,"PetscSpaceTensorGetSubspace_C",(PetscSpace,PetscInt,PetscSpace*),(sp,s,subsp));
435: return(0);
436: }
438: static PetscErrorCode PetscSpaceTensorSetNumSubspaces_Tensor(PetscSpace space, PetscInt numTensSpaces)
439: {
440: PetscSpace_Tensor *tens = (PetscSpace_Tensor *) space->data;
441: PetscInt Ns;
442: PetscErrorCode ierr;
445: if (tens->setupCalled) SETERRQ(PetscObjectComm((PetscObject)space),PETSC_ERR_ARG_WRONGSTATE,"Cannot change number of subspaces after setup called\n");
446: Ns = tens->numTensSpaces;
447: if (numTensSpaces == Ns) return(0);
448: if (Ns >= 0) {
449: PetscInt s;
451: for (s = 0; s < Ns; s++) {PetscSpaceDestroy(&tens->tensspaces[s]);}
452: PetscFree(tens->tensspaces);
453: }
454: Ns = tens->numTensSpaces = numTensSpaces;
455: PetscCalloc1(Ns, &tens->tensspaces);
456: return(0);
457: }
459: static PetscErrorCode PetscSpaceTensorGetNumSubspaces_Tensor(PetscSpace space, PetscInt *numTensSpaces)
460: {
461: PetscSpace_Tensor *tens = (PetscSpace_Tensor *) space->data;
464: *numTensSpaces = tens->numTensSpaces;
465: return(0);
466: }
468: static PetscErrorCode PetscSpaceTensorSetSubspace_Tensor(PetscSpace space, PetscInt s, PetscSpace subspace)
469: {
470: PetscSpace_Tensor *tens = (PetscSpace_Tensor *) space->data;
471: PetscInt Ns;
472: PetscErrorCode ierr;
475: if (tens->setupCalled) SETERRQ(PetscObjectComm((PetscObject)space),PETSC_ERR_ARG_WRONGSTATE,"Cannot change subspace after setup called\n");
476: Ns = tens->numTensSpaces;
477: if (Ns < 0) SETERRQ(PetscObjectComm((PetscObject)space),PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSpaceTensorSetNumSubspaces() first\n");
478: if (s < 0 || s >= Ns) SETERRQ1(PetscObjectComm((PetscObject)space),PETSC_ERR_ARG_OUTOFRANGE,"Invalid subspace number %D\n",subspace);
479: PetscObjectReference((PetscObject)subspace);
480: PetscSpaceDestroy(&tens->tensspaces[s]);
481: tens->tensspaces[s] = subspace;
482: return(0);
483: }
485: static PetscErrorCode PetscSpaceGetHeightSubspace_Tensor(PetscSpace sp, PetscInt height, PetscSpace *subsp)
486: {
487: PetscSpace_Tensor *tens = (PetscSpace_Tensor *) sp->data;
488: PetscInt Nc, dim, order, i;
489: PetscSpace bsp;
493: if (!tens->uniform) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Can only get a generic height subspace of a uniform tensor space: this tensor space is not uniform.\n");
494: PetscSpaceGetNumComponents(sp, &Nc);
495: PetscSpaceGetNumVariables(sp, &dim);
496: PetscSpaceGetDegree(sp, &order, NULL);
497: 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);}
498: if (!tens->heightsubspaces) {PetscCalloc1(dim, &tens->heightsubspaces);}
499: if (height <= dim) {
500: if (!tens->heightsubspaces[height-1]) {
501: PetscSpace sub;
503: PetscSpaceTensorGetSubspace(sp, 0, &bsp);
504: PetscSpaceCreate(PetscObjectComm((PetscObject) sp), &sub);
505: PetscSpaceSetType(sub, PETSCSPACETENSOR);
506: PetscSpaceSetNumComponents(sub, Nc);
507: PetscSpaceSetDegree(sub, order, PETSC_DETERMINE);
508: PetscSpaceSetNumVariables(sub, dim-height);
509: PetscSpaceTensorSetNumSubspaces(sub, dim-height);
510: for (i = 0; i < dim - height; i++) {
511: PetscSpaceTensorSetSubspace(sub, i, bsp);
512: }
513: PetscSpaceSetUp(sub);
514: tens->heightsubspaces[height-1] = sub;
515: }
516: *subsp = tens->heightsubspaces[height-1];
517: } else {
518: *subsp = NULL;
519: }
520: return(0);
521: }
523: static PetscErrorCode PetscSpaceTensorGetSubspace_Tensor(PetscSpace space, PetscInt s, PetscSpace *subspace)
524: {
525: PetscSpace_Tensor *tens = (PetscSpace_Tensor *) space->data;
526: PetscInt Ns;
529: Ns = tens->numTensSpaces;
530: if (Ns < 0) SETERRQ(PetscObjectComm((PetscObject)space),PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSpaceTensorSetNumSubspaces() first\n");
531: if (s < 0 || s >= Ns) SETERRQ1(PetscObjectComm((PetscObject)space),PETSC_ERR_ARG_OUTOFRANGE,"Invalid subspace number %D\n",subspace);
532: *subspace = tens->tensspaces[s];
533: return(0);
534: }
536: static PetscErrorCode PetscSpaceInitialize_Tensor(PetscSpace sp)
537: {
541: sp->ops->setfromoptions = PetscSpaceSetFromOptions_Tensor;
542: sp->ops->setup = PetscSpaceSetUp_Tensor;
543: sp->ops->view = PetscSpaceView_Tensor;
544: sp->ops->destroy = PetscSpaceDestroy_Tensor;
545: sp->ops->getdimension = PetscSpaceGetDimension_Tensor;
546: sp->ops->evaluate = PetscSpaceEvaluate_Tensor;
547: sp->ops->getheightsubspace = PetscSpaceGetHeightSubspace_Tensor;
548: PetscObjectComposeFunction((PetscObject) sp, "PetscSpaceTensorGetNumSubspaces_C", PetscSpaceTensorGetNumSubspaces_Tensor);
549: PetscObjectComposeFunction((PetscObject) sp, "PetscSpaceTensorSetNumSubspaces_C", PetscSpaceTensorSetNumSubspaces_Tensor);
550: PetscObjectComposeFunction((PetscObject) sp, "PetscSpaceTensorGetSubspace_C", PetscSpaceTensorGetSubspace_Tensor);
551: PetscObjectComposeFunction((PetscObject) sp, "PetscSpaceTensorSetSubspace_C", PetscSpaceTensorSetSubspace_Tensor);
552: return(0);
553: }
555: /*MC
556: PETSCSPACETENSOR = "tensor" - A PetscSpace object that encapsulates a tensor product space.
557: Subspaces are scalar spaces (num of componenents = 1), so the components
558: of a vector-valued tensor space are assumed to be identical.
560: Level: intermediate
562: .seealso: PetscSpaceType, PetscSpaceCreate(), PetscSpaceSetType()
563: M*/
565: PETSC_EXTERN PetscErrorCode PetscSpaceCreate_Tensor(PetscSpace sp)
566: {
567: PetscSpace_Tensor *tens;
568: PetscErrorCode ierr;
572: PetscNewLog(sp,&tens);
573: sp->data = tens;
575: tens->numTensSpaces = PETSC_DEFAULT;
577: PetscSpaceInitialize_Tensor(sp);
578: return(0);
579: }