Skip to content

Commit

Permalink
separate ex10.c linear and nonlinear parts of Jacobian
Browse files Browse the repository at this point in the history
inline ISLocalToGlobalMappingApply
  • Loading branch information
BarrySmith committed Sep 15, 2013
1 parent ceef4b8 commit 33b2b78
Show file tree
Hide file tree
Showing 3 changed files with 386 additions and 323 deletions.
45 changes: 44 additions & 1 deletion include/petscis.h
Expand Up @@ -138,13 +138,56 @@ typedef struct _p_ISLocalToGlobalMapping* ISLocalToGlobalMapping;
E*/ E*/
typedef enum {IS_GTOLM_MASK,IS_GTOLM_DROP} ISGlobalToLocalMappingType; typedef enum {IS_GTOLM_MASK,IS_GTOLM_DROP} ISGlobalToLocalMappingType;


#undef __FUNCT__
#define __FUNCT__ "ISLocalToGlobalMappingApply"
/*@C
ISLocalToGlobalMappingApply - Takes a list of integers in a local numbering
and converts them to the global numbering.
Not collective
Input Parameters:
+ mapping - the local to global mapping context
. N - number of integers
- in - input indices in local numbering
Output Parameter:
. out - indices in global numbering
Notes:
The in and out array parameters may be identical.
Level: advanced
.seealso: ISLocalToGlobalMappingCreate(),ISLocalToGlobalMappingDestroy(),
ISLocalToGlobalMappingApplyIS(),AOCreateBasic(),AOApplicationToPetsc(),
AOPetscToApplication(), ISGlobalToLocalMappingApply()
Concepts: mapping^local to global
@*/
PETSC_STATIC_INLINE PetscErrorCode ISLocalToGlobalMappingApply(ISLocalToGlobalMapping mapping,PetscInt N,const PetscInt in[],PetscInt out[])
{
PetscInt i,Nmax = mapping->n;
const PetscInt *idx = mapping->indices;

PetscFunctionBegin;
for (i=0; i<N; i++) {
if (in[i] < 0) {
out[i] = in[i];
continue;
}
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);
out[i] = idx[in[i]];
}
PetscFunctionReturn(0);
}

PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate(MPI_Comm,PetscInt,const PetscInt[],PetscCopyMode,ISLocalToGlobalMapping*); PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate(MPI_Comm,PetscInt,const PetscInt[],PetscCopyMode,ISLocalToGlobalMapping*);
PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreateIS(IS,ISLocalToGlobalMapping *); PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreateIS(IS,ISLocalToGlobalMapping *);
PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreateSF(PetscSF,PetscInt,ISLocalToGlobalMapping*); PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreateSF(PetscSF,PetscInt,ISLocalToGlobalMapping*);
PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingView(ISLocalToGlobalMapping,PetscViewer); PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingView(ISLocalToGlobalMapping,PetscViewer);
PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingDestroy(ISLocalToGlobalMapping*); PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingDestroy(ISLocalToGlobalMapping*);
PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingApplyIS(ISLocalToGlobalMapping,IS,IS*); PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingApplyIS(ISLocalToGlobalMapping,IS,IS*);
PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingApply(ISLocalToGlobalMapping,PetscInt,const PetscInt[],PetscInt[]);
PETSC_EXTERN PetscErrorCode ISGlobalToLocalMappingApply(ISLocalToGlobalMapping,ISGlobalToLocalMappingType,PetscInt,const PetscInt[],PetscInt*,PetscInt[]); PETSC_EXTERN PetscErrorCode ISGlobalToLocalMappingApply(ISLocalToGlobalMapping,ISGlobalToLocalMappingType,PetscInt,const PetscInt[],PetscInt*,PetscInt[]);
PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingGetSize(ISLocalToGlobalMapping,PetscInt*); PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingGetSize(ISLocalToGlobalMapping,PetscInt*);
PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt*,PetscInt*[],PetscInt*[],PetscInt**[]); PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt*,PetscInt*[],PetscInt*[],PetscInt**[]);
Expand Down
2 changes: 1 addition & 1 deletion src/mat/impls/aij/seq/aij.c
Expand Up @@ -328,7 +328,7 @@ PetscErrorCode MatSetValues_SeqAIJ(Mat A,PetscInt m,const PetscInt im[],PetscInt
} else { } else {
value = 0.; value = 0.;
} }
if (value == 0.0 && ignorezeroentries && (is == ADD_VALUES)) continue; if ((value == 0.0 && ignorezeroentries) && (is == ADD_VALUES)) continue;


if (col <= lastcol) low = 0; if (col <= lastcol) low = 0;
else high = nrow; else high = nrow;
Expand Down

0 comments on commit 33b2b78

Please sign in to comment.