1 #pragma once 2 3 #include <petsc/private/hashmap.h> 4 #include <petsc/private/hashijkey.h> 5 6 /* SUBMANSEC = Sys */ 7 /* 8 Hash map from (PetscInt,PetscInt) --> PetscScalar 9 */ 10 PETSC_HASH_MAP(HMapIJV, PetscHashIJKey, PetscScalar, PetscHashIJKeyHash, PetscHashIJKeyEqual, -1) 11 12 /*MC 13 PetscHMapIJVQueryAdd - Add value to the value of a given key if the key exists, 14 otherwise, insert a new (key,value) entry in the hash table 15 16 Synopsis: 17 #include <petsc/private/hashmapijv.h> 18 PetscErrorCode PetscHMapIJVQueryAdd(PetscHMapT ht,PetscHashIJKey key,PetscScalar val,PetscBool *missing) 19 20 Input Parameters: 21 + ht - The hash table 22 . key - The key 23 - val - The value 24 25 Output Parameter: 26 . missing - `PETSC_TRUE` if the `PetscHMapIJV` did not already have the given key 27 28 Level: developer 29 30 .seealso: `PetscHMapIJVSetWithMode()`, `PetscHMapIJV`, `PetscHMapIJVGet()`, `PetscHMapIJVIterSet()`, `PetscHMapIJVSet()` 31 M*/ 32 static inline PetscErrorCode PetscHMapIJVQueryAdd(PetscHMapIJV ht, PetscHashIJKey key, PetscScalar val, PetscBool *missing) 33 { 34 int ret; 35 khiter_t iter; 36 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