1 #ifndef _PETSC_HASHMAPIJV_H 2 #define _PETSC_HASHMAPIJV_H 3 4 #include <petsc/private/hashmap.h> 5 #include <petsc/private/hashijkey.h> 6 7 /* SUBMANSEC = Sys */ 8 /* 9 Hash map from (PetscInt,PetscInt) --> PetscScalar 10 */ 11 PETSC_HASH_MAP(HMapIJV, PetscHashIJKey, PetscScalar, PetscHashIJKeyHash, PetscHashIJKeyEqual, -1) 12 13 /*MC 14 PetscHMapIJVQueryAdd - Add value to the value of a given key if the key exists, 15 otherwise, insert a new (key,value) entry in the hash table 16 17 Synopsis: 18 #include <petsc/private/hashmapijv.h> 19 PetscErrorCode PetscHMapIJVQueryAdd(PetscHMapT ht,PetscHashIJKey key,PetscScalar val,PetscBool *missing) 20 21 Input Parameters: 22 + ht - The hash table 23 . key - The key 24 - val - The value 25 26 Output Parameter: 27 . missing - `PETSC_TRUE` if the `PetscHMapIJV` did not already have the given key 28 29 Level: developer 30 31 .seealso: `PetscHMapIJVSetWithMode()`, `PetscHMapIJV`, `PetscHMapIJVGet()`, `PetscHMapIJVIterSet()`, `PetscHMapIJVSet()` 32 M*/ 33 static inline PetscErrorCode PetscHMapIJVQueryAdd(PetscHMapIJV ht, PetscHashIJKey key, PetscScalar val, PetscBool *missing) 34 { 35 int ret; 36 khiter_t iter; 37 PetscFunctionBeginHot; 38 PetscAssertPointer(ht, 1); 39 iter = kh_put(HMapIJV, ht, key, &ret); 40 PetscHashAssert(ret >= 0); 41 if (ret) kh_val(ht, iter) = val; 42 else kh_val(ht, iter) += val; 43 *missing = ret ? PETSC_TRUE : PETSC_FALSE; 44 PetscFunctionReturn(PETSC_SUCCESS); 45 } 46 #endif /* _PETSC_HASHMAPIJV_H */ 47