Actual source code: petsc-isimpl.h
petsc-3.5.4 2015-05-23
1: /*
2: Index sets for scatter-gather type operations in vectors
3: and matrices.
5: */
7: #if !defined(_IS_H)
8: #define _IS_H
10: #include <petscis.h>
11: #include <petsc-private/petscimpl.h>
13: struct _ISOps {
14: PetscErrorCode (*getsize)(IS,PetscInt*);
15: PetscErrorCode (*getlocalsize)(IS,PetscInt*);
16: PetscErrorCode (*getindices)(IS,const PetscInt*[]);
17: PetscErrorCode (*restoreindices)(IS,const PetscInt*[]);
18: PetscErrorCode (*invertpermutation)(IS,PetscInt,IS*);
19: PetscErrorCode (*sort)(IS);
20: PetscErrorCode (*sortremovedups)(IS);
21: PetscErrorCode (*sorted)(IS,PetscBool*);
22: PetscErrorCode (*duplicate)(IS,IS*);
23: PetscErrorCode (*destroy)(IS);
24: PetscErrorCode (*view)(IS,PetscViewer);
25: PetscErrorCode (*load)(IS,PetscViewer);
26: PetscErrorCode (*identity)(IS,PetscBool*);
27: PetscErrorCode (*copy)(IS,IS);
28: PetscErrorCode (*togeneral)(IS);
29: PetscErrorCode (*oncomm)(IS,MPI_Comm,PetscCopyMode,IS*);
30: PetscErrorCode (*setblocksize)(IS,PetscInt);
31: PetscErrorCode (*contiguous)(IS,PetscInt,PetscInt,PetscInt*,PetscBool*);
32: };
34: struct _p_IS {
35: PETSCHEADER(struct _ISOps);
36: PetscLayout map;
37: PetscBool isperm; /* if is a permutation */
38: PetscInt max,min; /* range of possible values */
39: void *data;
40: PetscBool isidentity;
41: PetscInt *total, *nonlocal; /* local representation of ALL indices across the comm as well as the nonlocal part. */
42: PetscInt local_offset; /* offset to the local part within the total index set */
43: IS complement; /* IS wrapping nonlocal indices. */
44: };
46: extern PetscErrorCode ISLoad_Default(IS, PetscViewer);
48: struct _p_ISLocalToGlobalMapping{
49: PETSCHEADER(int);
50: PetscInt n; /* number of local indices */
51: PetscInt bs; /* blocksize; there is one index per block */
52: PetscInt *indices; /* global index of each local index */
53: PetscInt globalstart; /* first global referenced in indices */
54: PetscInt globalend; /* last + 1 global referenced in indices */
55: PetscInt *globals; /* local index for each global index between start and end */
56: };
58: /* ----------------------------------------------------------------------------*/
59: struct _p_PetscSection {
60: PETSCHEADER(int);
61: PetscInt pStart, pEnd; /* The chart: all points are contained in [pStart, pEnd) */
62: IS perm; /* A permutation of [0, pEnd-pStart) */
63: PetscInt *atlasDof; /* Describes layout of storage, point --> # of values */
64: PetscInt *atlasOff; /* Describes layout of storage, point --> offset into storage */
65: PetscInt maxDof; /* Maximum dof on any point */
66: PetscSection bc; /* Describes constraints, point --> # local dofs which are constrained */
67: PetscInt *bcIndices; /* Local indices for constrained dofs */
68: PetscBool setup;
70: PetscInt numFields; /* The number of fields making up the degrees of freedom */
71: const char **fieldNames; /* The field names */
72: PetscInt *numFieldComponents; /* The number of components in each field */
73: PetscSection *field; /* A section describing the layout and constraints for each field */
75: PetscObject clObj; /* Key for the closure (right now we only have one) */
76: PetscSection clSection; /* Section giving the number of points in each closure */
77: IS clPoints; /* Points in each closure */
78: };
81: #endif