Actual source code: petscis.h

petsc-3.3-p7 2013-05-11
  1: /*
  2:    An index set is a generalization of a subset of integers.  Index sets
  3:    are used for defining scatters and gathers.
  4: */
 7:  #include petscsys.h
 8:  #include petscsf.h

 10: #define IS_FILE_CLASSID 1211218
 11: PETSC_EXTERN PetscClassId IS_CLASSID;

 13: PETSC_EXTERN PetscErrorCode ISInitializePackage(const char[]);

 15: /*S
 16:      IS - Abstract PETSc object that allows indexing.

 18:    Level: beginner

 20:   Concepts: indexing, stride

 22: .seealso:  ISCreateGeneral(), ISCreateBlock(), ISCreateStride(), ISGetIndices(), ISDestroy()
 23: S*/
 24: typedef struct _p_IS* IS;

 26: /*J
 27:     ISType - String with the name of a PETSc vector or the creation function
 28:        with an optional dynamic library name, for example
 29:        http://www.mcs.anl.gov/petsc/lib.a:myveccreate()

 31:    Level: beginner

 33: .seealso: ISSetType(), IS
 34: J*/
 35: #define ISType char*
 36: #define ISGENERAL      "general"
 37: #define ISSTRIDE       "stride"
 38: #define ISBLOCK        "block"

 40: /* Dynamic creation and loading functions */
 41: PETSC_EXTERN PetscFList ISList;
 42: PETSC_EXTERN PetscBool ISRegisterAllCalled;
 43: PETSC_EXTERN PetscErrorCode ISSetType(IS, const ISType);
 44: PETSC_EXTERN PetscErrorCode ISGetType(IS, const ISType *);
 45: PETSC_EXTERN PetscErrorCode ISRegister(const char[],const char[],const char[],PetscErrorCode (*)(IS));
 46: PETSC_EXTERN PetscErrorCode ISRegisterAll(const char []);
 47: PETSC_EXTERN PetscErrorCode ISRegisterDestroy(void);
 48: PETSC_EXTERN PetscErrorCode ISCreate(MPI_Comm,IS*);

 50: /*MC
 51:   ISRegisterDynamic - Adds a new vector component implementation

 53:   Synopsis:
 54:   PetscErrorCode ISRegisterDynamic(const char *name, const char *path, const char *func_name, PetscErrorCode (*create_func)(IS))

 56:   Not Collective

 58:   Input Parameters:
 59: + name        - The name of a new user-defined creation routine
 60: . path        - The path (either absolute or relative) of the library containing this routine
 61: . func_name   - The name of routine to create method context
 62: - create_func - The creation routine itself

 64:   Notes:
 65:   ISRegisterDynamic() may be called multiple times to add several user-defined vectors

 67:   If dynamic libraries are used, then the fourth input argument (routine_create) is ignored.

 69:   Sample usage:
 70: .vb
 71:     ISRegisterDynamic("my_is_name","/home/username/my_lib/lib/libO/solaris/libmy.a", "MyISCreate", MyISCreate);
 72: .ve

 74:   Then, your vector type can be chosen with the procedural interface via
 75: .vb
 76:     ISCreate(MPI_Comm, IS *);
 77:     ISSetType(IS,"my_is_name");
 78: .ve
 79:    or at runtime via the option
 80: .vb
 81:     -is_type my_is_name
 82: .ve

 84:   Notes: $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
 85:          If your function is not being put into a shared library then use ISRegister() instead

 87:   This is no ISSetFromOptions() and the current implementations do not have a way to dynamically determine type, so
 88:   dynamic registration of custom IS types will be of limited use to users.

 90:   Level: developer

 92: .keywords: IS, register
 93: .seealso: ISRegisterAll(), ISRegisterDestroy(), ISRegister()
 94: M*/
 95: #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
 96: #define ISRegisterDynamic(a,b,c,d) ISRegister(a,b,c,0)
 97: #else
 98: #define ISRegisterDynamic(a,b,c,d) ISRegister(a,b,c,d)
 99: #endif

101: /*
102:     Default index set data structures that PETSc provides.
103: */
104: PETSC_EXTERN PetscErrorCode ISCreateGeneral(MPI_Comm,PetscInt,const PetscInt[],PetscCopyMode,IS *);
105: PETSC_EXTERN PetscErrorCode ISGeneralSetIndices(IS,PetscInt,const PetscInt[],PetscCopyMode);
106: PETSC_EXTERN PetscErrorCode ISCreateBlock(MPI_Comm,PetscInt,PetscInt,const PetscInt[],PetscCopyMode,IS *);
107: PETSC_EXTERN PetscErrorCode ISBlockSetIndices(IS,PetscInt,PetscInt,const PetscInt[],PetscCopyMode);
108: PETSC_EXTERN PetscErrorCode ISCreateStride(MPI_Comm,PetscInt,PetscInt,PetscInt,IS *);
109: PETSC_EXTERN PetscErrorCode ISStrideSetStride(IS,PetscInt,PetscInt,PetscInt);

111: PETSC_EXTERN PetscErrorCode ISDestroy(IS*);
112: PETSC_EXTERN PetscErrorCode ISSetPermutation(IS);
113: PETSC_EXTERN PetscErrorCode ISPermutation(IS,PetscBool *);
114: PETSC_EXTERN PetscErrorCode ISSetIdentity(IS);
115: PETSC_EXTERN PetscErrorCode ISIdentity(IS,PetscBool *);
116: PETSC_EXTERN PetscErrorCode ISContiguousLocal(IS,PetscInt,PetscInt,PetscInt*,PetscBool*);

118: PETSC_EXTERN PetscErrorCode ISGetIndices(IS,const PetscInt *[]);
119: PETSC_EXTERN PetscErrorCode ISRestoreIndices(IS,const PetscInt *[]);
120: PETSC_EXTERN PetscErrorCode ISGetTotalIndices(IS,const PetscInt *[]);
121: PETSC_EXTERN PetscErrorCode ISRestoreTotalIndices(IS,const PetscInt *[]);
122: PETSC_EXTERN PetscErrorCode ISGetNonlocalIndices(IS,const PetscInt *[]);
123: PETSC_EXTERN PetscErrorCode ISRestoreNonlocalIndices(IS,const PetscInt *[]);
124: PETSC_EXTERN PetscErrorCode ISGetNonlocalIS(IS, IS *is);
125: PETSC_EXTERN PetscErrorCode ISRestoreNonlocalIS(IS, IS *is);
126: PETSC_EXTERN PetscErrorCode ISGetSize(IS,PetscInt *);
127: PETSC_EXTERN PetscErrorCode ISGetLocalSize(IS,PetscInt *);
128: PETSC_EXTERN PetscErrorCode ISInvertPermutation(IS,PetscInt,IS*);
129: PETSC_EXTERN PetscErrorCode ISView(IS,PetscViewer);
130: PETSC_EXTERN PetscErrorCode ISEqual(IS,IS,PetscBool  *);
131: PETSC_EXTERN PetscErrorCode ISSort(IS);
132: PETSC_EXTERN PetscErrorCode ISSorted(IS,PetscBool  *);
133: PETSC_EXTERN PetscErrorCode ISDifference(IS,IS,IS*);
134: PETSC_EXTERN PetscErrorCode ISSum(IS,IS,IS*);
135: PETSC_EXTERN PetscErrorCode ISExpand(IS,IS,IS*);

137: PETSC_EXTERN PetscErrorCode ISBlockGetIndices(IS,const PetscInt *[]);
138: PETSC_EXTERN PetscErrorCode ISBlockRestoreIndices(IS,const PetscInt *[]);
139: PETSC_EXTERN PetscErrorCode ISBlockGetLocalSize(IS,PetscInt *);
140: PETSC_EXTERN PetscErrorCode ISBlockGetSize(IS,PetscInt *);
141: PETSC_EXTERN PetscErrorCode ISGetBlockSize(IS,PetscInt*);
142: PETSC_EXTERN PetscErrorCode ISSetBlockSize(IS,PetscInt);

144: PETSC_EXTERN PetscErrorCode ISStrideGetInfo(IS,PetscInt *,PetscInt*);

146: PETSC_EXTERN PetscErrorCode ISToGeneral(IS);

148: PETSC_EXTERN PetscErrorCode ISDuplicate(IS,IS*);
149: PETSC_EXTERN PetscErrorCode ISCopy(IS,IS);
150: PETSC_EXTERN PetscErrorCode ISAllGather(IS,IS*);
151: PETSC_EXTERN PetscErrorCode ISComplement(IS,PetscInt,PetscInt,IS*);
152: PETSC_EXTERN PetscErrorCode ISConcatenate(MPI_Comm,PetscInt,const IS[],IS*);
153: PETSC_EXTERN PetscErrorCode ISListToMap(MPI_Comm,PetscInt, IS[],IS*,IS*);
154: PETSC_EXTERN PetscErrorCode ISMapToList(IS,IS,PetscInt*, IS *[]);
155: PETSC_EXTERN PetscErrorCode ISMapFactorRight(IS,IS,PetscBool,IS*);
156: PETSC_EXTERN PetscErrorCode ISOnComm(IS,MPI_Comm,PetscCopyMode,IS*);

158: /* --------------------------------------------------------------------------*/
159: PETSC_EXTERN PetscClassId IS_LTOGM_CLASSID;

161: /*S
162:    ISLocalToGlobalMapping - mappings from an arbitrary
163:       local ordering from 0 to n-1 to a global PETSc ordering 
164:       used by a vector or matrix.

166:    Level: intermediate

168:    Note: mapping from Local to Global is scalable; but Global
169:   to Local may not be if the range of global values represented locally
170:   is very large.

172:    Note: the ISLocalToGlobalMapping is actually a private object; it is included
173:   here for the inline function ISLocalToGlobalMappingApply() to allow it to be inlined since
174:   it is used so often.

176: .seealso:  ISLocalToGlobalMappingCreate()
177: S*/
178: struct _p_ISLocalToGlobalMapping{
179:   PETSCHEADER(int);
180:   PetscInt n;                  /* number of local indices */
181:   PetscInt *indices;           /* global index of each local index */
182:   PetscInt globalstart;        /* first global referenced in indices */
183:   PetscInt globalend;          /* last + 1 global referenced in indices */
184:   PetscInt *globals;           /* local index for each global index between start and end */
185: };
186: typedef struct _p_ISLocalToGlobalMapping* ISLocalToGlobalMapping;

188: /*E
189:     ISGlobalToLocalMappingType - Indicates if missing global indices are 

191:    IS_GTOLM_MASK - missing global indices are replaced with -1
192:    IS_GTOLM_DROP - missing global indices are dropped

194:    Level: beginner

196: .seealso: ISGlobalToLocalMappingApply()

198: E*/
199: typedef enum {IS_GTOLM_MASK,IS_GTOLM_DROP} ISGlobalToLocalMappingType;

201: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate(MPI_Comm,PetscInt,const PetscInt[],PetscCopyMode,ISLocalToGlobalMapping*);
202: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreateIS(IS,ISLocalToGlobalMapping *);
203: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreateSF(PetscSF,PetscInt,ISLocalToGlobalMapping*);
204: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingView(ISLocalToGlobalMapping,PetscViewer);
205: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingDestroy(ISLocalToGlobalMapping*);
206: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingApplyIS(ISLocalToGlobalMapping,IS,IS*);
207: PETSC_EXTERN PetscErrorCode ISGlobalToLocalMappingApply(ISLocalToGlobalMapping,ISGlobalToLocalMappingType,PetscInt,const PetscInt[],PetscInt*,PetscInt[]);
208: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingGetSize(ISLocalToGlobalMapping,PetscInt*);
209: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt*,PetscInt*[],PetscInt*[],PetscInt**[]);
210: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingRestoreInfo(ISLocalToGlobalMapping,PetscInt*,PetscInt*[],PetscInt*[],PetscInt**[]);
211: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingGetIndices(ISLocalToGlobalMapping,const PetscInt**);
212: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingRestoreIndices(ISLocalToGlobalMapping,const PetscInt**);
213: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingBlock(ISLocalToGlobalMapping,PetscInt,ISLocalToGlobalMapping*);
214: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingUnBlock(ISLocalToGlobalMapping,PetscInt,ISLocalToGlobalMapping*);
215: PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingConcatenate(MPI_Comm,PetscInt,const ISLocalToGlobalMapping[],ISLocalToGlobalMapping*);

219: PETSC_STATIC_INLINE PetscErrorCode ISLocalToGlobalMappingApply(ISLocalToGlobalMapping mapping,PetscInt N,const PetscInt in[],PetscInt out[])
220: {
221:   PetscInt       i,Nmax = mapping->n;
222:   const PetscInt *idx = mapping->indices;
224:   for (i=0; i<N; i++) {
225:     if (in[i] < 0) {out[i] = in[i]; continue;}
226:     if (in[i] >= Nmax) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local index %D too large %D (max) at %D",in[i],Nmax,i);
227:     out[i] = idx[in[i]];
228:   }
229:   return(0);
230: }

232: /* --------------------------------------------------------------------------*/
233: /*E
234:     ISColoringType - determines if the coloring is for the entire parallel grid/graph/matrix
235:                      or for just the local ghosted portion

237:     Level: beginner

239: $   IS_COLORING_GLOBAL - does not include the colors for ghost points, this is used when the function 
240: $                        is called synchronously in parallel. This requires generating a "parallel coloring".
241: $   IS_COLORING_GHOSTED - includes colors for ghost points, this is used when the function can be called
242: $                         seperately on individual processes with the ghost points already filled in. Does not
243: $                         require a "parallel coloring", rather each process colors its local + ghost part.
244: $                         Using this can result in much less parallel communication. In the paradigm of 
245: $                         DMGetLocalVector() and DMGetGlobalVector() this could be called IS_COLORING_LOCAL

247: .seealso: DMCreateColoring()
248: E*/
249: typedef enum {IS_COLORING_GLOBAL,IS_COLORING_GHOSTED} ISColoringType;
250: PETSC_EXTERN const char *ISColoringTypes[];
251: typedef unsigned PETSC_IS_COLOR_VALUE_TYPE ISColoringValue;
252: PETSC_EXTERN PetscErrorCode ISAllGatherColors(MPI_Comm,PetscInt,ISColoringValue*,PetscInt*,ISColoringValue*[]);

254: /*S
255:      ISColoring - sets of IS's that define a coloring
256:               of the underlying indices

258:    Level: intermediate

260:     Notes:
261:         One should not access the *is records below directly because they may not yet 
262:     have been created. One should use ISColoringGetIS() to make sure they are 
263:     created when needed.

265:     Developer Note: this is not a PetscObject

267: .seealso:  ISColoringCreate(), ISColoringGetIS(), ISColoringView(), ISColoringGetIS()
268: S*/
269: struct _n_ISColoring {
270:   PetscInt        refct;
271:   PetscInt        n;                /* number of colors */
272:   IS              *is;              /* for each color indicates columns */
273:   MPI_Comm        comm;
274:   ISColoringValue *colors;          /* for each column indicates color */
275:   PetscInt        N;                /* number of columns */
276:   ISColoringType  ctype;
277: };
278: typedef struct _n_ISColoring* ISColoring;

280: PETSC_EXTERN PetscErrorCode ISColoringCreate(MPI_Comm,PetscInt,PetscInt,const ISColoringValue[],ISColoring*);
281: PETSC_EXTERN PetscErrorCode ISColoringDestroy(ISColoring*);
282: PETSC_EXTERN PetscErrorCode ISColoringView(ISColoring,PetscViewer);
283: PETSC_EXTERN PetscErrorCode ISColoringGetIS(ISColoring,PetscInt*,IS*[]);
284: PETSC_EXTERN PetscErrorCode ISColoringRestoreIS(ISColoring,IS*[]);
285: #define ISColoringReference(coloring) ((coloring)->refct++,0)
286: #define ISColoringSetType(coloring,type) ((coloring)->ctype = type,0)

288: /* --------------------------------------------------------------------------*/

290: PETSC_EXTERN PetscErrorCode ISPartitioningToNumbering(IS,IS*);
291: PETSC_EXTERN PetscErrorCode ISPartitioningCount(IS,PetscInt,PetscInt[]);

293: PETSC_EXTERN PetscErrorCode ISCompressIndicesGeneral(PetscInt,PetscInt,PetscInt,PetscInt,const IS[],IS[]);
294: PETSC_EXTERN PetscErrorCode ISCompressIndicesSorted(PetscInt,PetscInt,PetscInt,const IS[],IS[]);
295: PETSC_EXTERN PetscErrorCode ISExpandIndicesGeneral(PetscInt,PetscInt,PetscInt,PetscInt,const IS[],IS[]);


298: /* Reset __FUNCT__ in case the user does not define it themselves */

302: #endif