Actual source code: vsection.c
1: /*
2: This file contains routines for basic section object implementation.
3: */
5: #include <private/vecimpl.h> /*I "petscvec.h" I*/
7: #if 0
8: /* Should I protect these for C++? */
11: PetscErrorCode PetscSectionGetDof(PetscUniformSection s, PetscInt point, PetscInt *numDof)
12: {
14: if ((point < s->pStart) || (point >= s->pEnd)) {
15: SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section point %p should be in [%d, %d)", point, s->pStart, s->pEnd);
16: }
17: *numDof = s->numDof;
18: return(0);
19: }
23: PetscErrorCode PetscSectionGetOffset(PetscUniformSection s, PetscInt point, PetscInt *offset)
24: {
26: if ((point < s->pStart) || (point >= s->pEnd)) {
27: SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section point %p should be in [%d, %d)", point, s->pStart, s->pEnd);
28: }
29: *offset = s->numDof*(point - s->pStart);
30: return(0);
31: }
32: #endif
36: /*@C
37: PetscSectionCreate - Allocates PetscSection space and sets the map contents to the default.
39: Collective on MPI_Comm
41: Input Parameters:
42: + comm - the MPI communicator
43: - s - pointer to the section
45: Level: developer
47: Notes: Typical calling sequence
48: PetscLayoutCreate(MPI_Comm,PetscLayout *);
49: PetscLayoutSetBlockSize(PetscLayout,1);
50: PetscLayoutSetSize(PetscLayout,n) or PetscLayoutSetLocalSize(PetscLayout,N);
51: PetscLayoutSetUp(PetscLayout);
52: PetscLayoutGetSize(PetscLayout,PetscInt *); or PetscLayoutGetLocalSize(PetscLayout,PetscInt *;)
53: PetscLayoutDestroy(PetscLayout);
55: The PetscSection object and methods are intended to be used in the PETSc Vec and Mat implementions; it is
56: recommended they not be used in user codes unless you really gain something in their use.
58: Fortran Notes:
59: Not available from Fortran
61: .seealso: PetscSection, PetscSectionDestroy()
62: @*/
63: PetscErrorCode PetscSectionCreate(MPI_Comm comm, PetscSection *s)
64: {
68: PetscNew(struct _n_PetscSection, s);
69: (*s)->atlasLayout.comm = comm;
70: (*s)->atlasLayout.pStart = -1;
71: (*s)->atlasLayout.pEnd = -1;
72: (*s)->atlasLayout.numDof = 1;
73: (*s)->atlasDof = PETSC_NULL;
74: (*s)->atlasOff = PETSC_NULL;
75: (*s)->bc = PETSC_NULL;
76: (*s)->bcIndices = PETSC_NULL;
77: return(0);
78: }
82: PetscErrorCode PetscSectionCheckConstraints(PetscSection s)
83: {
87: if (!s->bc) {
88: PetscSectionCreate(s->atlasLayout.comm, &s->bc);
89: PetscSectionSetChart(s->bc, s->atlasLayout.pStart, s->atlasLayout.pEnd);
90: }
91: return(0);
92: }
96: PetscErrorCode PetscSectionGetChart(PetscSection s, PetscInt *pStart, PetscInt *pEnd)
97: {
99: if (pStart) {*pStart = s->atlasLayout.pStart;}
100: if (pEnd) {*pEnd = s->atlasLayout.pEnd;}
101: return(0);
102: }
106: PetscErrorCode PetscSectionSetChart(PetscSection s, PetscInt pStart, PetscInt pEnd)
107: {
111: s->atlasLayout.pStart = pStart;
112: s->atlasLayout.pEnd = pEnd;
113: PetscFree2(s->atlasDof, s->atlasOff);
114: PetscMalloc2((pEnd - pStart), PetscInt, &s->atlasDof, (pEnd - pStart), PetscInt, &s->atlasOff);
115: PetscMemzero(s->atlasDof, (pEnd - pStart)*sizeof(PetscInt));
116: return(0);
117: }
121: PetscErrorCode PetscSectionGetDof(PetscSection s, PetscInt point, PetscInt *numDof)
122: {
124: if ((point < s->atlasLayout.pStart) || (point >= s->atlasLayout.pEnd)) {
125: SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section point %p should be in [%d, %d)", point, s->atlasLayout.pStart, s->atlasLayout.pEnd);
126: }
127: *numDof = s->atlasDof[point - s->atlasLayout.pStart];
128: return(0);
129: }
133: PetscErrorCode PetscSectionSetDof(PetscSection s, PetscInt point, PetscInt numDof)
134: {
136: if ((point < s->atlasLayout.pStart) || (point >= s->atlasLayout.pEnd)) {
137: SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section point %p should be in [%d, %d)", point, s->atlasLayout.pStart, s->atlasLayout.pEnd);
138: }
139: s->atlasDof[point - s->atlasLayout.pStart] = numDof;
140: return(0);
141: }
145: PetscErrorCode PetscSectionGetConstraintDof(PetscSection s, PetscInt point, PetscInt *numDof)
146: {
150: if (s->bc) {
151: PetscSectionGetDof(s->bc, point, numDof);
152: } else {
153: *numDof = 0;
154: }
155: return(0);
156: }
160: PetscErrorCode PetscSectionSetConstraintDof(PetscSection s, PetscInt point, PetscInt numDof)
161: {
165: if (numDof) {
166: PetscSectionCheckConstraints(s);
167: PetscSectionSetDof(s->bc, point, numDof);
168: }
169: return(0);
170: }
174: PetscErrorCode PetscSectionSetUp(PetscSection s)
175: {
176: PetscInt offset = 0, p;
180: for(p = 0; p < s->atlasLayout.pEnd - s->atlasLayout.pStart; ++p) {
181: s->atlasOff[p] = offset;
182: offset += s->atlasDof[p];
183: }
184: if (s->bc) {
185: const PetscInt last = (s->bc->atlasLayout.pEnd-s->bc->atlasLayout.pStart) - 1;
187: PetscSectionSetUp(s->bc);
188: PetscMalloc((s->bc->atlasOff[last] + s->bc->atlasDof[last]) * sizeof(PetscInt), &s->bcIndices);
189: }
190: return(0);
191: }
195: PetscErrorCode PetscSectionGetStorageSize(PetscSection s, PetscInt *size)
196: {
197: PetscInt p, n = 0;
200: for(p = 0; p < s->atlasLayout.pEnd - s->atlasLayout.pStart; ++p) {
201: n += s->atlasDof[p];
202: }
203: *size = n;
204: return(0);
205: }
209: PetscErrorCode PetscSectionGetOffset(PetscSection s, PetscInt point, PetscInt *offset)
210: {
212: if ((point < s->atlasLayout.pStart) || (point >= s->atlasLayout.pEnd)) {
213: SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section point %p should be in [%d, %d)", point, s->atlasLayout.pStart, s->atlasLayout.pEnd);
214: }
215: *offset = s->atlasOff[point - s->atlasLayout.pStart];
216: return(0);
217: }
219: /*@C
220: PetscSectionDestroy - Frees a section object and frees its range if that exists.
222: Collective on MPI_Comm
224: Input Parameters:
225: . s - the PetscSection
227: Level: developer
229: The PetscSection object and methods are intended to be used in the PETSc Vec and Mat implementions; it is
230: recommended they not be used in user codes unless you really gain something in their use.
232: Fortran Notes:
233: Not available from Fortran
235: .seealso: PetscSection, PetscSectionCreate()
236: @*/
239: PetscErrorCode PetscSectionDestroy(PetscSection *s)
240: {
244: if (!*s) return(0);
245: if (!(*s)->refcnt--) {
246: PetscSectionDestroy(&(*s)->bc);
247: PetscFree((*s)->bcIndices);
248: PetscFree2((*s)->atlasDof, (*s)->atlasOff);
249: PetscFree((*s));
250: }
251: *s = PETSC_NULL;
252: return(0);
253: }
257: PetscErrorCode VecGetValuesSection(Vec v, PetscSection s, PetscInt point, PetscScalar **values)
258: {
259: PetscScalar *baseArray;
260: const PetscInt p = point - s->atlasLayout.pStart;
264: VecGetArray(v, &baseArray);
265: *values = &baseArray[s->atlasOff[p]];
266: VecRestoreArray(v, &baseArray);
267: return(0);
268: }
272: PetscErrorCode VecIntGetValuesSection(PetscInt *baseArray, PetscSection s, PetscInt point, PetscInt **values)
273: {
274: const PetscInt p = point - s->atlasLayout.pStart;
277: *values = &baseArray[s->atlasOff[p]];
278: return(0);
279: }
283: PetscErrorCode VecSetValuesSection(Vec v, PetscSection s, PetscInt point, PetscScalar values[], InsertMode mode)
284: {
285: PetscScalar *baseArray, *array;
286: const PetscBool doInsert = mode == INSERT_VALUES || mode == INSERT_ALL_VALUES ? PETSC_TRUE : PETSC_FALSE;
287: const PetscBool doBC = mode == INSERT_ALL_VALUES || mode == ADD_ALL_VALUES ? PETSC_TRUE : PETSC_FALSE;
288: const PetscInt p = point - s->atlasLayout.pStart;
289: const PetscInt orientation = 0; /* Needs to be included for use in closure operations */
290: PetscInt cDim = 0;
291: PetscErrorCode ierr;
294: PetscSectionGetConstraintDof(s, p, &cDim);
295: VecGetArray(v, &baseArray);
296: array = &baseArray[s->atlasOff[p]];
297: if (!cDim) {
298: if (orientation >= 0) {
299: const PetscInt dim = s->atlasDof[p];
300: PetscInt i;
302: if (doInsert) {
303: for(i = 0; i < dim; ++i) {
304: array[i] = values[i];
305: }
306: } else {
307: for(i = 0; i < dim; ++i) {
308: array[i] += values[i];
309: }
310: }
311: } else {
312: #if 0
313: int offset = 0;
314: int j = -1;
316: for(int space = 0; space < this->getNumSpaces(); ++space) {
317: const int& dim = this->getFiberDimension(p, space);
319: for(int i = dim-1; i >= 0; --i) {
320: array[++j] = values[i+offset];
321: }
322: offset += dim;
323: }
324: #else
325: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Fibration is not yet implemented in PetscSection");
326: #endif
327: }
328: } else {
329: if (orientation >= 0) {
330: const PetscInt dim = s->atlasDof[p];
331: PetscInt cInd = 0, i;
332: PetscInt *cDof;
334: PetscSectionGetConstraintIndices(s, point, &cDof);
335: if (doInsert) {
336: for(i = 0; i < dim; ++i) {
337: if ((cInd < cDim) && (i == cDof[cInd])) {
338: if (doBC) {array[i] = values[i];} /* Constrained update */
339: ++cInd;
340: continue;
341: }
342: array[i] = values[i]; /* Unconstrained update */
343: }
344: } else {
345: for(i = 0; i < dim; ++i) {
346: if ((cInd < cDim) && (i == cDof[cInd])) {
347: if (doBC) {array[i] += values[i];} /* Constrained update */
348: ++cInd;
349: continue;
350: }
351: array[i] += values[i]; /* Unconstrained update */
352: }
353: }
354: } else {
355: #if 0
356: const PetscInt *cDof;
357: PetscInt offset = 0;
358: PetscInt cOffset = 0;
359: PetscInt j = 0, space;
361: PetscSectionGetConstraintDof(s, point, &cDof);
362: for(space = 0; space < this->getNumSpaces(); ++space) {
363: const PetscInt dim = this->getFiberDimension(p, space);
364: const PetscInt tDim = this->getConstrainedFiberDimension(p, space);
365: const PetscInt sDim = dim - tDim;
366: PetscInt cInd = 0, i ,k;
368: for(i = 0, k = dim+offset-1; i < dim; ++i, ++j, --k) {
369: if ((cInd < sDim) && (j == cDof[cInd+cOffset])) {++cInd; continue;}
370: array[j] = values[k];
371: }
372: offset += dim;
373: cOffset += dim - tDim;
374: }
375: #else
376: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Fibration is not yet implemented in PetscSection");
377: #endif
378: }
379: }
380: VecRestoreArray(v, &baseArray);
381: return(0);
382: }
387: PetscErrorCode VecIntSetValuesSection(PetscInt *baseArray, PetscSection s, PetscInt point, PetscInt values[], InsertMode mode)
388: {
389: PetscInt *array;
390: const PetscInt p = point - s->atlasLayout.pStart;
391: const PetscInt orientation = 0; /* Needs to be included for use in closure operations */
392: PetscInt cDim = 0;
396: PetscSectionGetConstraintDof(s, p, &cDim);
397: array = &baseArray[s->atlasOff[p]];
398: if (!cDim) {
399: if (orientation >= 0) {
400: const PetscInt dim = s->atlasDof[p];
401: PetscInt i;
403: if (mode == INSERT_VALUES) {
404: for(i = 0; i < dim; ++i) {
405: array[i] = values[i];
406: }
407: } else {
408: for(i = 0; i < dim; ++i) {
409: array[i] += values[i];
410: }
411: }
412: } else {
413: #if 0
414: int offset = 0;
415: int j = -1;
417: for(int space = 0; space < this->getNumSpaces(); ++space) {
418: const int& dim = this->getFiberDimension(p, space);
420: for(int i = dim-1; i >= 0; --i) {
421: array[++j] = values[i+offset];
422: }
423: offset += dim;
424: }
425: #else
426: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Fibration is not yet implemented in PetscSection");
427: #endif
428: }
429: } else {
430: if (orientation >= 0) {
431: const PetscInt dim = s->atlasDof[p];
432: PetscInt cInd = 0, i;
433: PetscInt *cDof;
435: PetscSectionGetConstraintIndices(s, point, &cDof);
436: if (mode == INSERT_VALUES) {
437: for(i = 0; i < dim; ++i) {
438: if ((cInd < cDim) && (i == cDof[cInd])) {++cInd; continue;}
439: array[i] = values[i];
440: }
441: } else {
442: for(i = 0; i < dim; ++i) {
443: if ((cInd < cDim) && (i == cDof[cInd])) {++cInd; continue;}
444: array[i] += values[i];
445: }
446: }
447: } else {
448: #if 0
449: const PetscInt *cDof;
450: PetscInt offset = 0;
451: PetscInt cOffset = 0;
452: PetscInt j = 0, space;
454: PetscSectionGetConstraintDof(s, point, &cDof);
455: for(space = 0; space < this->getNumSpaces(); ++space) {
456: const PetscInt dim = this->getFiberDimension(p, space);
457: const PetscInt tDim = this->getConstrainedFiberDimension(p, space);
458: const PetscInt sDim = dim - tDim;
459: PetscInt cInd = 0, i ,k;
461: for(i = 0, k = dim+offset-1; i < dim; ++i, ++j, --k) {
462: if ((cInd < sDim) && (j == cDof[cInd+cOffset])) {++cInd; continue;}
463: array[j] = values[k];
464: }
465: offset += dim;
466: cOffset += dim - tDim;
467: }
468: #else
469: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Fibration is not yet implemented in PetscSection");
470: #endif
471: }
472: }
473: return(0);
474: }
478: PetscErrorCode PetscSectionGetConstraintIndices(PetscSection s, PetscInt point, PetscInt **indices)
479: {
483: if (s->bc) {
484: VecIntGetValuesSection(s->bcIndices, s->bc, point, indices);
485: } else {
486: *indices = PETSC_NULL;
487: }
488: return(0);
489: }
493: PetscErrorCode PetscSectionSetConstraintIndices(PetscSection s, PetscInt point, PetscInt indices[])
494: {
498: if (s->bc) {
499: VecIntSetValuesSection(s->bcIndices, s->bc, point, indices, INSERT_VALUES);
500: }
501: return(0);
502: }