xref: /petsc/src/dm/dt/space/impls/tensor/spacetensor.c (revision 43bfef1cf37704822d1568e92f52e31ae916812c)
1 #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
2 
3 static PetscErrorCode PetscSpaceTensorCreateSubspace(PetscSpace space, PetscInt Nvs, PetscSpace *subspace)
4 {
5   PetscInt    degree;
6   const char *prefix;
7   PetscErrorCode ierr;
8 
9   PetscFunctionBegin;
10   ierr = PetscSpaceGetDegree(space, &degree, NULL);CHKERRQ(ierr);
11   ierr = PetscObjectGetOptionsPrefix((PetscObject)space, &prefix);CHKERRQ(ierr);
12   ierr = PetscSpaceCreate(PetscObjectComm((PetscObject)space), subspace);CHKERRQ(ierr);
13   ierr = PetscSpaceSetType(*subspace, PETSCSPACEPOLYNOMIAL);CHKERRQ(ierr);
14   ierr = PetscSpaceSetNumVariables(*subspace, Nvs);CHKERRQ(ierr);
15   ierr = PetscSpaceSetNumComponents(*subspace, 1);CHKERRQ(ierr);
16   ierr = PetscSpaceSetDegree(*subspace, degree, PETSC_DETERMINE);CHKERRQ(ierr);
17   ierr = PetscObjectSetOptionsPrefix((PetscObject)*subspace, prefix);CHKERRQ(ierr);
18   ierr = PetscObjectAppendOptionsPrefix((PetscObject)*subspace, "subspace_");CHKERRQ(ierr);
19   PetscFunctionReturn(0);
20 }
21 
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;
28 
29   PetscFunctionBegin;
30   ierr = PetscSpaceGetNumVariables(sp, &Nv);CHKERRQ(ierr);
31   if (!Nv) PetscFunctionReturn(0);
32   ierr = PetscSpaceGetNumComponents(sp, &Nc);CHKERRQ(ierr);
33   ierr = PetscSpaceTensorGetNumSubspaces(sp, &Ns);CHKERRQ(ierr);
34   ierr = PetscSpaceGetDegree(sp, &deg, NULL);CHKERRQ(ierr);
35   if (Ns > 1) {
36     PetscSpace s0;
37 
38     ierr = PetscSpaceTensorGetSubspace(sp, 0, &s0);CHKERRQ(ierr);
39     for (i = 1; i < Ns; i++) {
40       PetscSpace si;
41 
42       ierr = PetscSpaceTensorGetSubspace(sp, i, &si);CHKERRQ(ierr);
43       if (si != s0) {uniform = PETSC_FALSE; break;}
44     }
45   }
46   Ns = (Ns == PETSC_DEFAULT) ? PetscMax(Nv,1) : Ns;
47   ierr = PetscOptionsHead(PetscOptionsObject,"PetscSpace tensor options");CHKERRQ(ierr);
48   ierr = PetscOptionsInt("-petscspace_tensor_spaces", "The number of subspaces", "PetscSpaceTensorSetNumSubspaces", Ns, &Ns, NULL);CHKERRQ(ierr);
49   ierr = PetscOptionsBool("-petscspace_tensor_uniform", "Subspaces are identical", "PetscSpaceTensorSetFromOptions", uniform, &uniform, NULL);CHKERRQ(ierr);
50   ierr = PetscOptionsTail();CHKERRQ(ierr);
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->numSpaces) {ierr = PetscSpaceTensorSetNumSubspaces(sp, Ns);CHKERRQ(ierr);}
54   if (uniform) {
55     PetscInt   Nvs = Nv / Ns;
56     PetscSpace subspace;
57 
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     ierr = PetscSpaceTensorGetSubspace(sp, 0, &subspace);CHKERRQ(ierr);
60     if (!subspace) {ierr = PetscSpaceTensorCreateSubspace(sp, Nvs, &subspace);CHKERRQ(ierr);}
61     else           {ierr = PetscObjectReference((PetscObject)subspace);CHKERRQ(ierr);}
62     ierr = PetscSpaceSetFromOptions(subspace);CHKERRQ(ierr);
63     for (i = 0; i < Ns; i++) {ierr = PetscSpaceTensorSetSubspace(sp, i, subspace);CHKERRQ(ierr);}
64     ierr = PetscSpaceDestroy(&subspace);CHKERRQ(ierr);
65   } else {
66     for (i = 0; i < Ns; i++) {
67       PetscSpace subspace;
68 
69       ierr = PetscSpaceTensorGetSubspace(sp, i, &subspace);CHKERRQ(ierr);
70       if (!subspace) {
71         char tprefix[128];
72 
73         ierr = PetscSpaceTensorCreateSubspace(sp, 1, &subspace);CHKERRQ(ierr);
74         ierr = PetscSNPrintf(tprefix, 128, "%d_",(int)i);CHKERRQ(ierr);
75         ierr = PetscObjectAppendOptionsPrefix((PetscObject)subspace, tprefix);CHKERRQ(ierr);
76       } else {
77         ierr = PetscObjectReference((PetscObject)subspace);CHKERRQ(ierr);
78       }
79       ierr = PetscSpaceSetFromOptions(subspace);CHKERRQ(ierr);
80       ierr = PetscSpaceTensorSetSubspace(sp, i, subspace);CHKERRQ(ierr);
81       ierr = PetscSpaceDestroy(&subspace);CHKERRQ(ierr);
82     }
83   }
84   PetscFunctionReturn(0);
85 }
86 
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->numSpaces, i, n;
92   PetscErrorCode     ierr;
93 
94   PetscFunctionBegin;
95   for (i = 1; i < Ns; i++) {
96     if (tens->spaces[i] != tens->spaces[0]) {uniform = PETSC_FALSE; break;}
97   }
98   if (uniform) {ierr = PetscViewerASCIIPrintf(viewer, "Tensor space of %D subspaces (all identical)\n", Ns);CHKERRQ(ierr);
99   } else       {ierr = PetscViewerASCIIPrintf(viewer, "Tensor space of %D subspaces\n", Ns);CHKERRQ(ierr);}
100   n = uniform ? 1 : Ns;
101   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
102   for (i = 0; i < n; i++) {
103     ierr = PetscSpaceView(tens->spaces[i], viewer);CHKERRQ(ierr);
104   }
105   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
106   PetscFunctionReturn(0);
107 }
108 
109 static PetscErrorCode PetscSpaceView_Tensor(PetscSpace sp, PetscViewer viewer)
110 {
111   PetscBool      iascii;
112   PetscErrorCode ierr;
113 
114   PetscFunctionBegin;
115   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
116   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
117   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
118   if (iascii) {ierr = PetscSpaceTensorView_Ascii(sp, viewer);CHKERRQ(ierr);}
119   PetscFunctionReturn(0);
120 }
121 
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;
129 
130   PetscFunctionBegin;
131   if (tens->setupCalled) PetscFunctionReturn(0);
132   ierr = PetscSpaceGetNumVariables(sp, &Nv);CHKERRQ(ierr);
133   ierr = PetscSpaceTensorGetNumSubspaces(sp, &Ns);CHKERRQ(ierr);
134   if (Ns == PETSC_DEFAULT) {
135     Ns = Nv;
136     ierr = PetscSpaceTensorSetNumSubspaces(sp, Ns);CHKERRQ(ierr);
137   }
138   if (!Ns) {
139     if (Nv) SETERRQ(PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_OUTOFRANGE, "Cannot have zero subspaces");
140   } else {
141     PetscSpace s0;
142 
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     ierr = PetscSpaceTensorGetSubspace(sp, 0, &s0);CHKERRQ(ierr);
145     for (i = 1; i < Ns; i++) {
146       PetscSpace si;
147 
148       ierr = PetscSpaceTensorGetSubspace(sp, i, &si);CHKERRQ(ierr);
149       if (si != s0) {uniform = PETSC_FALSE; break;}
150     }
151     if (uniform) {
152       PetscInt   Nvs = Nv / Ns;
153 
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) {ierr = PetscSpaceTensorCreateSubspace(sp, Nvs, &s0);CHKERRQ(ierr);}
156       else     {ierr = PetscObjectReference((PetscObject) s0);CHKERRQ(ierr);}
157       ierr = PetscSpaceSetUp(s0);CHKERRQ(ierr);
158       for (i = 0; i < Ns; i++) {ierr = PetscSpaceTensorSetSubspace(sp, i, s0);CHKERRQ(ierr);}
159       ierr = PetscSpaceDestroy(&s0);CHKERRQ(ierr);
160     } else {
161       for (i = 0 ; i < Ns; i++) {
162         PetscSpace si;
163 
164         ierr = PetscSpaceTensorGetSubspace(sp, i, &si);CHKERRQ(ierr);
165         if (!si) {ierr = PetscSpaceTensorCreateSubspace(sp, 1, &si);CHKERRQ(ierr);}
166         else     {ierr = PetscObjectReference((PetscObject) si);CHKERRQ(ierr);}
167         ierr = PetscSpaceSetUp(si);CHKERRQ(ierr);
168         ierr = PetscSpaceTensorSetSubspace(sp, i, si);CHKERRQ(ierr);
169         ierr = PetscSpaceDestroy(&si);CHKERRQ(ierr);
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;
178 
179     ierr = PetscSpaceTensorGetSubspace(sp, i, &si);CHKERRQ(ierr);
180     ierr = PetscSpaceGetDegree(si, &iDeg, &iMaxDeg);CHKERRQ(ierr);
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   PetscFunctionReturn(0);
189 }
190 
191 static PetscErrorCode PetscSpaceDestroy_Tensor(PetscSpace sp)
192 {
193   PetscSpace_Tensor *tens    = (PetscSpace_Tensor *) sp->data;
194   PetscInt           Ns, i;
195   PetscErrorCode     ierr;
196 
197   PetscFunctionBegin;
198   ierr = PetscSpaceTensorGetNumSubspaces(sp, &Ns);CHKERRQ(ierr);
199   for (i = 0; i < Ns; i++) {ierr = PetscSpaceDestroy(&tens->spaces[i]);CHKERRQ(ierr);}
200   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscSpaceTensorSetSubspace_C", NULL);CHKERRQ(ierr);
201   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscSpaceTensorGetSubspace_C", NULL);CHKERRQ(ierr);
202   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscSpaceTensorSetNumSubspaces_C", NULL);CHKERRQ(ierr);
203   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscSpaceTensorGetNumSubspaces_C", NULL);CHKERRQ(ierr);
204   ierr = PetscFree(tens->spaces);CHKERRQ(ierr);
205   ierr = PetscFree(tens);CHKERRQ(ierr);
206   PetscFunctionReturn(0);
207 }
208 
209 static PetscErrorCode PetscSpaceGetDimension_Tensor(PetscSpace sp, PetscInt *dim)
210 {
211   PetscSpace_Tensor *tens = (PetscSpace_Tensor *) sp->data;
212   PetscInt           i, Ns, Nc, d;
213   PetscErrorCode     ierr;
214 
215   PetscFunctionBegin;
216   ierr = PetscSpaceSetUp(sp);CHKERRQ(ierr);
217   Ns = tens->numSpaces;
218   Nc = sp->Nc;
219   d  = 1;
220   for (i = 0; i < Ns; i++) {
221     PetscInt id;
222 
223     ierr = PetscSpaceGetDimension(tens->spaces[i], &id);CHKERRQ(ierr);
224     d *= id;
225   }
226   d *= Nc;
227   *dim = d;
228   PetscFunctionReturn(0);
229 }
230 
231 static PetscErrorCode PetscSpaceEvaluate_Tensor(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[])
232 {
233   PetscSpace_Tensor *tens  = (PetscSpace_Tensor *) sp->data;
234   DM               dm      = sp->dm;
235   PetscInt         Nc      = sp->Nc;
236   PetscInt         Nv      = sp->Nv;
237   PetscInt         Ns;
238   PetscReal       *lpoints, *sB = NULL, *sD = NULL, *sH = NULL;
239   PetscInt         c, pdim, d, e, der, der2, i, l, si, p, s, step;
240   PetscErrorCode   ierr;
241 
242   PetscFunctionBegin;
243   if (!tens->setupCalled) {ierr = PetscSpaceSetUp(sp);CHKERRQ(ierr);}
244   Ns = tens->numSpaces;
245   ierr = PetscSpaceGetDimension(sp,&pdim);CHKERRQ(ierr);
246   pdim /= Nc;
247   ierr = DMGetWorkArray(dm, npoints*Nv, MPIU_REAL, &lpoints);CHKERRQ(ierr);
248   if (B || D || H) {ierr = DMGetWorkArray(dm, npoints*pdim,       MPIU_REAL, &sB);CHKERRQ(ierr);}
249   if (D || H)      {ierr = DMGetWorkArray(dm, npoints*pdim*Nv,    MPIU_REAL, &sD);CHKERRQ(ierr);}
250   if (H)           {ierr = DMGetWorkArray(dm, npoints*pdim*Nv*Nv, MPIU_REAL, &sH);CHKERRQ(ierr);}
251   if (B) {
252     for (i = 0; i < npoints*pdim; i++) B[i * Nc*Nc] = 1.;
253   }
254   if (D) {
255     for (i = 0; i < npoints*pdim; i++) {
256       for (l = 0; l < Nv; l++) {
257         D[i * Nc*Nc*Nv + l] = 1.;
258       }
259     }
260   }
261   if (H) {
262     for (i = 0; i < npoints*pdim; i++) {
263       for (l = 0; l < Nv*Nv; l++) {
264         H[i * Nc*Nc*Nv*Nv + l] = 1.;
265       }
266     }
267   }
268   for (s = 0, d = 0, step = 1; s < Ns; s++) {
269     PetscInt sNv, spdim;
270     PetscInt skip, j, k;
271 
272     ierr = PetscSpaceGetNumVariables(tens->spaces[s], &sNv);CHKERRQ(ierr);
273     ierr = PetscSpaceGetDimension(tens->spaces[s], &spdim);CHKERRQ(ierr);
274     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);
275     skip = pdim / (step * spdim);
276     for (p = 0; p < npoints; ++p) {
277       for (i = 0; i < sNv; i++) {
278         lpoints[p * sNv + i] = points[p*Nv + d + i];
279       }
280     }
281     ierr = PetscSpaceEvaluate(tens->spaces[s], npoints, lpoints, sB, sD, sH);CHKERRQ(ierr);
282     if (B) {
283       for (p = 0; p < npoints; p++) {
284         for (k = 0; k < skip; k++) {
285           for (si = 0; si < spdim; si++) {
286             for (j = 0; j < step; j++) {
287               i = (k * spdim + si) * step + j;
288               B[(pdim * p + i) * Nc * Nc] *= sB[spdim * p + si];
289             }
290           }
291         }
292       }
293     }
294     if (D) {
295       for (p = 0; p < npoints; p++) {
296         for (k = 0; k < skip; k++) {
297           for (si = 0; si < spdim; si++) {
298             for (j = 0; j < step; j++) {
299               i = (k * spdim + si) * step + j;
300               for (der = 0; der < Nv; der++) {
301                 if (der >= d && der < d + sNv) {
302                   D[(pdim * p + i) * Nc*Nc*Nv + der] *= sD[(spdim * p + si) * sNv + der - d];
303                 } else {
304                   D[(pdim * p + i) * Nc*Nc*Nv + der] *= sB[spdim * p + si];
305                 }
306               }
307             }
308           }
309         }
310       }
311     }
312     if (H) {
313       for (p = 0; p < npoints; p++) {
314         for (k = 0; k < skip; k++) {
315           for (si = 0; si < spdim; si++) {
316             for (j = 0; j < step; j++) {
317               i = (k * spdim + si) * step + j;
318               for (der = 0; der < Nv; der++) {
319                 for (der2 = 0; der2 < Nv; der2++) {
320                   if (der >= d && der < d + sNv && der2 >= d && der2 < d + sNv) {
321                     H[((pdim * p + i) * Nc*Nc*Nv + der) * Nv + der2] *= sH[((spdim * p + si) * sNv + der - d) * sNv + der2 - d];
322                   } else if (der >= d && der < d + sNv) {
323                     H[((pdim * p + i) * Nc*Nc*Nv + der) * Nv + der2] *= sD[(spdim * p + si) * sNv + der - d];
324                   } else if (der2 >= d && der2 < d + sNv) {
325                     H[((pdim * p + i) * Nc*Nc*Nv + der) * Nv + der2] *= sD[(spdim * p + si) * sNv + der2 - d];
326                   } else {
327                     H[((pdim * p + i) * Nc*Nc*Nv + der) * Nv + der2] *= sB[spdim * p + si];
328                   }
329                 }
330               }
331             }
332           }
333         }
334       }
335     }
336     d += sNv;
337     step *= spdim;
338   }
339   if (B && Nc > 1) {
340     /* Make direct sum basis for multicomponent space */
341     for (p = 0; p < npoints; ++p) {
342       for (i = 0; i < pdim; ++i) {
343         for (c = 1; c < Nc; ++c) {
344           B[(p*pdim*Nc + i*Nc + c)*Nc + c] = B[(p*pdim + i)*Nc*Nc];
345         }
346       }
347     }
348   }
349   if (D && Nc > 1) {
350     /* Make direct sum basis for multicomponent space */
351     for (p = 0; p < npoints; ++p) {
352       for (i = 0; i < pdim; ++i) {
353         for (c = 1; c < Nc; ++c) {
354           for (d = 0; d < Nv; ++d) {
355             D[((p*pdim*Nc + i*Nc + c)*Nc + c)*Nv + d] = D[(p*pdim + i)*Nc*Nc*Nv + d];
356           }
357         }
358       }
359     }
360   }
361   if (H && Nc > 1) {
362     /* Make direct sum basis for multicomponent space */
363     for (p = 0; p < npoints; ++p) {
364       for (i = 0; i < pdim; ++i) {
365         for (c = 1; c < Nc; ++c) {
366           for (d = 0; d < Nv; ++d) {
367             for (e = 0; e < Nv; ++e) {
368               H[(((p*pdim*Nc + i*Nc + c)*Nc + c)*Nv + d)*Nv + e] = H[((p*pdim + i)*Nc*Nc*Nv + d)*Nv + e];
369             }
370           }
371         }
372       }
373     }
374   }
375   if (H)           {ierr = DMRestoreWorkArray(dm, npoints*pdim*Nv*Nv, MPIU_REAL, &sH);CHKERRQ(ierr);}
376   if (D || H)      {ierr = DMRestoreWorkArray(dm, npoints*pdim*Nv,    MPIU_REAL, &sD);CHKERRQ(ierr);}
377   if (B || D || H) {ierr = DMRestoreWorkArray(dm, npoints*pdim,       MPIU_REAL, &sB);CHKERRQ(ierr);}
378   ierr = DMRestoreWorkArray(dm, npoints*Nv, MPIU_REAL, &lpoints);CHKERRQ(ierr);
379   PetscFunctionReturn(0);
380 }
381 
382 PetscErrorCode PetscSpaceTensorSetNumSubspaces(PetscSpace sp, PetscInt numSpaces)
383 {
384   PetscErrorCode ierr;
385 
386   PetscFunctionBegin;
387   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
388   ierr = PetscTryMethod(sp,"PetscSpaceTensorSetNumSubspaces_C",(PetscSpace,PetscInt),(sp,numSpaces));CHKERRQ(ierr);
389   PetscFunctionReturn(0);
390 }
391 
392 PetscErrorCode PetscSpaceTensorGetNumSubspaces(PetscSpace sp, PetscInt *numSpaces)
393 {
394   PetscErrorCode ierr;
395 
396   PetscFunctionBegin;
397   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
398   PetscValidIntPointer(numSpaces, 2);
399   ierr = PetscTryMethod(sp,"PetscSpaceTensorGetNumSubspaces_C",(PetscSpace,PetscInt*),(sp,numSpaces));CHKERRQ(ierr);
400   PetscFunctionReturn(0);
401 }
402 
403 PetscErrorCode PetscSpaceTensorSetSubspace(PetscSpace sp, PetscInt s, PetscSpace subsp)
404 {
405   PetscErrorCode ierr;
406 
407   PetscFunctionBegin;
408   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
409   if (subsp) PetscValidHeaderSpecific(subsp, PETSCSPACE_CLASSID, 3);
410   ierr = PetscTryMethod(sp,"PetscSpaceTensorSetSubspace_C",(PetscSpace,PetscInt,PetscSpace),(sp,s,subsp));CHKERRQ(ierr);
411   PetscFunctionReturn(0);
412 }
413 
414 PetscErrorCode PetscSpaceTensorGetSubspace(PetscSpace sp, PetscInt s, PetscSpace *subsp)
415 {
416   PetscErrorCode ierr;
417 
418   PetscFunctionBegin;
419   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
420   PetscValidPointer(subsp, 3);
421   ierr = PetscTryMethod(sp,"PetscSpaceTensorGetSubspace_C",(PetscSpace,PetscInt,PetscSpace*),(sp,s,subsp));CHKERRQ(ierr);
422   PetscFunctionReturn(0);
423 }
424 
425 static PetscErrorCode PetscSpaceTensorSetNumSubspaces_Tensor(PetscSpace space, PetscInt numSpaces)
426 {
427   PetscSpace_Tensor *tens = (PetscSpace_Tensor *) space->data;
428   PetscInt           Ns;
429   PetscErrorCode     ierr;
430 
431   PetscFunctionBegin;
432   if (tens->setupCalled) SETERRQ(PetscObjectComm((PetscObject)space),PETSC_ERR_ARG_WRONGSTATE,"Cannot change number of subspaces after setup called\n");
433   Ns = tens->numSpaces;
434   if (numSpaces == Ns) PetscFunctionReturn(0);
435   if (Ns >= 0) {
436     PetscInt s;
437 
438     for (s = 0; s < Ns; s++) {ierr = PetscSpaceDestroy(&tens->spaces[s]);CHKERRQ(ierr);}
439     ierr = PetscFree(tens->spaces);CHKERRQ(ierr);
440   }
441   Ns = tens->numSpaces = numSpaces;
442   ierr = PetscCalloc1(Ns, &tens->spaces);CHKERRQ(ierr);
443   PetscFunctionReturn(0);
444 }
445 
446 static PetscErrorCode PetscSpaceTensorGetNumSubspaces_Tensor(PetscSpace space, PetscInt *numSpaces)
447 {
448   PetscSpace_Tensor *tens = (PetscSpace_Tensor *) space->data;
449 
450   PetscFunctionBegin;
451   *numSpaces = tens->numSpaces;
452   PetscFunctionReturn(0);
453 }
454 
455 static PetscErrorCode PetscSpaceTensorSetSubspace_Tensor(PetscSpace space, PetscInt s, PetscSpace subspace)
456 {
457   PetscSpace_Tensor *tens = (PetscSpace_Tensor *) space->data;
458   PetscInt           Ns;
459   PetscErrorCode     ierr;
460 
461   PetscFunctionBegin;
462   if (tens->setupCalled) SETERRQ(PetscObjectComm((PetscObject)space),PETSC_ERR_ARG_WRONGSTATE,"Cannot change subspace after setup called\n");
463   Ns = tens->numSpaces;
464   if (Ns < 0) SETERRQ(PetscObjectComm((PetscObject)space),PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSpaceTensorSetNumSubspaces() first\n");
465   if (s < 0 || s >= Ns) SETERRQ1(PetscObjectComm((PetscObject)space),PETSC_ERR_ARG_OUTOFRANGE,"Invalid subspace number %D\n",subspace);
466   ierr = PetscObjectReference((PetscObject)subspace);CHKERRQ(ierr);
467   ierr = PetscSpaceDestroy(&tens->spaces[s]);CHKERRQ(ierr);
468   tens->spaces[s] = subspace;
469   PetscFunctionReturn(0);
470 }
471 
472 static PetscErrorCode PetscSpaceTensorGetSubspace_Tensor(PetscSpace space, PetscInt s, PetscSpace *subspace)
473 {
474   PetscSpace_Tensor *tens = (PetscSpace_Tensor *) space->data;
475   PetscInt           Ns;
476 
477   PetscFunctionBegin;
478   Ns = tens->numSpaces;
479   if (Ns < 0) SETERRQ(PetscObjectComm((PetscObject)space),PETSC_ERR_ARG_WRONGSTATE,"Must call PetscSpaceTensorSetNumSubspaces() first\n");
480   if (s < 0 || s >= Ns) SETERRQ1(PetscObjectComm((PetscObject)space),PETSC_ERR_ARG_OUTOFRANGE,"Invalid subspace number %D\n",subspace);
481   *subspace = tens->spaces[s];
482   PetscFunctionReturn(0);
483 }
484 
485 static PetscErrorCode PetscSpaceInitialize_Tensor(PetscSpace sp)
486 {
487   PetscErrorCode ierr;
488 
489   PetscFunctionBegin;
490   sp->ops->setfromoptions    = PetscSpaceSetFromOptions_Tensor;
491   sp->ops->setup             = PetscSpaceSetUp_Tensor;
492   sp->ops->view              = PetscSpaceView_Tensor;
493   sp->ops->destroy           = PetscSpaceDestroy_Tensor;
494   sp->ops->getdimension      = PetscSpaceGetDimension_Tensor;
495   sp->ops->evaluate          = PetscSpaceEvaluate_Tensor;
496   sp->ops->getheightsubspace = NULL;
497   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscSpaceTensorGetNumSubspaces_C", PetscSpaceTensorGetNumSubspaces_Tensor);CHKERRQ(ierr);
498   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscSpaceTensorSetNumSubspaces_C", PetscSpaceTensorSetNumSubspaces_Tensor);CHKERRQ(ierr);
499   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscSpaceTensorGetSubspace_C", PetscSpaceTensorGetSubspace_Tensor);CHKERRQ(ierr);
500   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscSpaceTensorSetSubspace_C", PetscSpaceTensorSetSubspace_Tensor);CHKERRQ(ierr);
501   PetscFunctionReturn(0);
502 }
503 
504 /*MC
505   PETSCSPACETENSOR = "tensor" - A PetscSpace object that encapsulates a tensor product space.
506                      Subspaces are scalar spaces (num of componenents = 1), vector components
507                      are assumed to be identical.
508 
509   Level: intermediate
510 
511 .seealso: PetscSpaceType, PetscSpaceCreate(), PetscSpaceSetType()
512 M*/
513 
514 PETSC_EXTERN PetscErrorCode PetscSpaceCreate_Tensor(PetscSpace sp)
515 {
516   PetscSpace_Tensor *tens;
517   PetscErrorCode   ierr;
518 
519   PetscFunctionBegin;
520   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
521   ierr     = PetscNewLog(sp,&tens);CHKERRQ(ierr);
522   sp->data = tens;
523 
524   tens->numSpaces = PETSC_DEFAULT;
525 
526   ierr = PetscSpaceInitialize_Tensor(sp);CHKERRQ(ierr);
527   PetscFunctionReturn(0);
528 }
529 
530