1 #define PETSC_DLL 2 /* 3 Provides utility routines for manulating any type of PETSc object. 4 */ 5 #include "petscsys.h" /*I "petscsys.h" I*/ 6 7 #undef __FUNCT__ 8 #define __FUNCT__ "PetscObjectStateQuery" 9 /*@C 10 PetscObjectStateQuery - Gets the state of any PetscObject, 11 regardless of the type. 12 13 Not Collective 14 15 Input Parameter: 16 . obj - any PETSc object, for example a Vec, Mat or KSP. This must be 17 cast with a (PetscObject), for example, 18 PetscObjectStateQuery((PetscObject)mat,&state); 19 20 Output Parameter: 21 . state - the object state 22 23 Notes: object state is an integer which gets increased every time 24 the object is changed. By saving and later querying the object state 25 one can determine whether information about the object is still current. 26 Currently, state is maintained for Vec and Mat objects. 27 28 Level: advanced 29 30 seealso: PetscObjectStateIncrease, PetscObjectSetState 31 32 Concepts: state 33 34 @*/ 35 PetscErrorCode PETSC_DLLEXPORT PetscObjectStateQuery(PetscObject obj,PetscInt *state) 36 { 37 PetscFunctionBegin; 38 PetscValidHeader(obj,1); 39 PetscValidIntPointer(state,2); 40 *state = obj->state; 41 PetscFunctionReturn(0); 42 } 43 44 #undef __FUNCT__ 45 #define __FUNCT__ "PetscObjectSetState" 46 /*@C 47 PetscObjectSetState - Sets the state of any PetscObject, 48 regardless of the type. 49 50 Not Collective 51 52 Input Parameter: 53 + obj - any PETSc object, for example a Vec, Mat or KSP. This must be 54 cast with a (PetscObject), for example, 55 PetscObjectSetState((PetscObject)mat,state); 56 - state - the object state 57 58 Notes: This function should be used with extreme caution. There is 59 essentially only one use for it: if the user calls Mat(Vec)GetRow(Array), 60 which increases the state, but does not alter the data, then this 61 routine can be used to reset the state. 62 63 Level: advanced 64 65 seealso: PetscObjectStateQuery,PetscObjectStateIncrease 66 67 Concepts: state 68 69 @*/ 70 PetscErrorCode PETSC_DLLEXPORT PetscObjectSetState(PetscObject obj,PetscInt state) 71 { 72 PetscFunctionBegin; 73 PetscValidHeader(obj,1); 74 obj->state = state; 75 PetscFunctionReturn(0); 76 } 77 78 PetscInt PETSC_DLLEXPORT globalcurrentstate = 0; 79 PetscInt PETSC_DLLEXPORT globalmaxstate = 10; 80 /*@C 81 PetscObjectComposedDataRegister - Get an available id for 82 composed data 83 84 Not Collective 85 86 Output parameter: 87 . id - an identifier under which data can be stored 88 89 Level: developer 90 91 seealso: PetscObjectComposedDataSetInt 92 93 @*/ 94 PetscErrorCode PETSC_DLLEXPORT PetscObjectComposedDataRegister(PetscInt *id) 95 { 96 PetscFunctionBegin; 97 *id = globalcurrentstate++; 98 if (globalcurrentstate > globalmaxstate) globalmaxstate += 10; 99 PetscFunctionReturn(0); 100 } 101 102 PetscErrorCode PETSC_DLLEXPORT PetscObjectComposedDataIncreaseInt(PetscObject obj) 103 { 104 PetscInt *ar = obj->intcomposeddata,*new_ar; 105 PetscInt *ir = obj->intcomposedstate,*new_ir,n = obj->int_idmax,new_n,i; 106 PetscErrorCode ierr; 107 108 PetscFunctionBegin; 109 new_n = globalmaxstate; 110 ierr = PetscMalloc(new_n*sizeof(PetscInt),&new_ar);CHKERRQ(ierr); 111 ierr = PetscMemzero(new_ar,new_n*sizeof(PetscInt));CHKERRQ(ierr); 112 ierr = PetscMalloc(new_n*sizeof(PetscInt),&new_ir);CHKERRQ(ierr); 113 ierr = PetscMemzero(new_ir,new_n*sizeof(PetscInt));CHKERRQ(ierr); 114 if (n) { 115 for (i=0; i<n; i++) { 116 new_ar[i] = ar[i]; new_ir[i] = ir[i]; 117 } 118 ierr = PetscFree(ar);CHKERRQ(ierr); 119 ierr = PetscFree(ir);CHKERRQ(ierr); 120 } 121 obj->int_idmax = new_n; 122 obj->intcomposeddata = new_ar; obj->intcomposedstate = new_ir; 123 PetscFunctionReturn(0); 124 } 125 PetscErrorCode PETSC_DLLEXPORT PetscObjectComposedDataIncreaseIntstar(PetscObject obj) 126 { 127 PetscInt **ar = obj->intstarcomposeddata,**new_ar; 128 PetscInt *ir = obj->intstarcomposedstate,*new_ir,n = obj->intstar_idmax,new_n,i; 129 PetscErrorCode ierr; 130 131 PetscFunctionBegin; 132 new_n = globalmaxstate; 133 ierr = PetscMalloc(new_n*sizeof(PetscInt*),&new_ar);CHKERRQ(ierr); 134 ierr = PetscMemzero(new_ar,new_n*sizeof(PetscInt*));CHKERRQ(ierr); 135 ierr = PetscMalloc(new_n*sizeof(PetscInt),&new_ir);CHKERRQ(ierr); 136 ierr = PetscMemzero(new_ir,new_n*sizeof(PetscInt));CHKERRQ(ierr); 137 if (n) { 138 for (i=0; i<n; i++) { 139 new_ar[i] = ar[i]; new_ir[i] = ir[i]; 140 } 141 ierr = PetscFree(ar);CHKERRQ(ierr); 142 ierr = PetscFree(ir);CHKERRQ(ierr); 143 } 144 obj->intstar_idmax = new_n; 145 obj->intstarcomposeddata = new_ar; obj->intstarcomposedstate = new_ir; 146 PetscFunctionReturn(0); 147 } 148 149 PetscErrorCode PETSC_DLLEXPORT PetscObjectComposedDataIncreaseReal(PetscObject obj) 150 { 151 PetscReal *ar = obj->realcomposeddata,*new_ar; 152 PetscInt *ir = obj->realcomposedstate,*new_ir,n = obj->real_idmax,new_n,i; 153 PetscErrorCode ierr; 154 155 PetscFunctionBegin; 156 new_n = globalmaxstate; 157 ierr = PetscMalloc(new_n*sizeof(PetscReal),&new_ar);CHKERRQ(ierr); 158 ierr = PetscMemzero(new_ar,new_n*sizeof(PetscReal));CHKERRQ(ierr); 159 ierr = PetscMalloc(new_n*sizeof(PetscInt),&new_ir);CHKERRQ(ierr); 160 ierr = PetscMemzero(new_ir,new_n*sizeof(PetscInt));CHKERRQ(ierr); 161 if (n) { 162 for (i=0; i<n; i++) { 163 new_ar[i] = ar[i]; new_ir[i] = ir[i]; 164 } 165 ierr = PetscFree(ar);CHKERRQ(ierr); 166 ierr = PetscFree(ir);CHKERRQ(ierr); 167 } 168 obj->real_idmax = new_n; 169 obj->realcomposeddata = new_ar; obj->realcomposedstate = new_ir; 170 PetscFunctionReturn(0); 171 } 172 173 PetscErrorCode PETSC_DLLEXPORT PetscObjectComposedDataIncreaseRealstar(PetscObject obj) 174 { 175 PetscReal **ar = obj->realstarcomposeddata,**new_ar; 176 PetscInt *ir = obj->realstarcomposedstate,*new_ir,n = obj->realstar_idmax,new_n,i; 177 PetscErrorCode ierr; 178 179 PetscFunctionBegin; 180 new_n = globalmaxstate; 181 ierr = PetscMalloc(new_n*sizeof(PetscReal*),&new_ar);CHKERRQ(ierr); 182 ierr = PetscMemzero(new_ar,new_n*sizeof(PetscReal*));CHKERRQ(ierr); 183 ierr = PetscMalloc(new_n*sizeof(PetscInt),&new_ir);CHKERRQ(ierr); 184 ierr = PetscMemzero(new_ir,new_n*sizeof(PetscInt));CHKERRQ(ierr); 185 if (n) { 186 for (i=0; i<n; i++) { 187 new_ar[i] = ar[i]; new_ir[i] = ir[i]; 188 } 189 ierr = PetscFree(ar);CHKERRQ(ierr); 190 ierr = PetscFree(ir);CHKERRQ(ierr); 191 } 192 obj->realstar_idmax = new_n; 193 obj->realstarcomposeddata = new_ar; obj->realstarcomposedstate = new_ir; 194 PetscFunctionReturn(0); 195 } 196 197 PetscErrorCode PETSC_DLLEXPORT PetscObjectComposedDataIncreaseScalar(PetscObject obj) 198 { 199 PetscScalar *ar = obj->scalarcomposeddata,*new_ar; 200 PetscInt *ir = obj->scalarcomposedstate,*new_ir,n = obj->scalar_idmax,new_n,i; 201 PetscErrorCode ierr; 202 203 PetscFunctionBegin; 204 new_n = globalmaxstate; 205 ierr = PetscMalloc(new_n*sizeof(PetscScalar),&new_ar);CHKERRQ(ierr); 206 ierr = PetscMemzero(new_ar,new_n*sizeof(PetscScalar));CHKERRQ(ierr); 207 ierr = PetscMalloc(new_n*sizeof(PetscInt),&new_ir);CHKERRQ(ierr); 208 ierr = PetscMemzero(new_ir,new_n*sizeof(PetscInt));CHKERRQ(ierr); 209 if (n) { 210 for (i=0; i<n; i++) { 211 new_ar[i] = ar[i]; new_ir[i] = ir[i]; 212 } 213 ierr = PetscFree(ar);CHKERRQ(ierr); 214 ierr = PetscFree(ir);CHKERRQ(ierr); 215 } 216 obj->scalar_idmax = new_n; 217 obj->scalarcomposeddata = new_ar; obj->scalarcomposedstate = new_ir; 218 PetscFunctionReturn(0); 219 } 220 221 PetscErrorCode PETSC_DLLEXPORT PetscObjectComposedDataIncreaseScalarstar(PetscObject obj) 222 { 223 PetscScalar **ar = obj->scalarstarcomposeddata,**new_ar; 224 PetscInt *ir = obj->scalarstarcomposedstate,*new_ir,n = obj->scalarstar_idmax,new_n,i; 225 PetscErrorCode ierr; 226 227 PetscFunctionBegin; 228 new_n = globalmaxstate; 229 ierr = PetscMalloc(new_n*sizeof(PetscScalar*),&new_ar);CHKERRQ(ierr); 230 ierr = PetscMemzero(new_ar,new_n*sizeof(PetscScalar*));CHKERRQ(ierr); 231 ierr = PetscMalloc(new_n*sizeof(PetscInt),&new_ir);CHKERRQ(ierr); 232 ierr = PetscMemzero(new_ir,new_n*sizeof(PetscInt));CHKERRQ(ierr); 233 if (n) { 234 for (i=0; i<n; i++) { 235 new_ar[i] = ar[i]; new_ir[i] = ir[i]; 236 } 237 ierr = PetscFree(ar);CHKERRQ(ierr); 238 ierr = PetscFree(ir);CHKERRQ(ierr); 239 } 240 obj->scalarstar_idmax = new_n; 241 obj->scalarstarcomposeddata = new_ar; obj->scalarstarcomposedstate = new_ir; 242 PetscFunctionReturn(0); 243 } 244 245