xref: /petsc/src/sys/objects/state.c (revision b12e77b2dafe08a3135b110b3db60a43ac3be1c7)
1 #define PETSC_DLL
2 /*
3      Provides utility routines for manulating any type of PETSc object.
4 */
5 #include "petsc.h"  /*I   "petsc.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   if (!obj) SETERRQ(PETSC_ERR_ARG_CORRUPT,"Null object");
39   *state = obj->state;
40   PetscFunctionReturn(0);
41 }
42 
43 #undef __FUNCT__
44 #define __FUNCT__ "PetscObjectSetState"
45 /*@C
46    PetscObjectSetState - Sets the state of any PetscObject,
47    regardless of the type.
48 
49    Not Collective
50 
51    Input Parameter:
52 +  obj - any PETSc object, for example a Vec, Mat or KSP. This must be
53          cast with a (PetscObject), for example,
54          PetscObjectSetState((PetscObject)mat,state);
55 -  state - the object state
56 
57    Notes: This function should be used with extreme caution. There is
58    essentially only one use for it: if the user calls Mat(Vec)GetRow(Array),
59    which increases the state, but does not alter the data, then this
60    routine can be used to reset the state.
61 
62    Level: advanced
63 
64    seealso: PetscObjectStateQuery,PetscObjectStateIncrease
65 
66    Concepts: state
67 
68 @*/
69 PetscErrorCode PETSC_DLLEXPORT PetscObjectSetState(PetscObject obj,PetscInt state)
70 {
71   PetscFunctionBegin;
72   if (!obj) SETERRQ(PETSC_ERR_ARG_CORRUPT,"Null object");
73   obj->state = state;
74   PetscFunctionReturn(0);
75 }
76 
77 PetscInt PETSC_DLLEXPORT globalcurrentstate = 0;
78 PetscInt PETSC_DLLEXPORT globalmaxstate = 10;
79 PetscErrorCode PETSC_DLLEXPORT PetscObjectComposedDataRegister(PetscInt *id)
80 {
81   PetscFunctionBegin;
82   *id = globalcurrentstate++;
83   if (globalcurrentstate > globalmaxstate) globalmaxstate += 10;
84   PetscFunctionReturn(0);
85 }
86 
87 PetscErrorCode PETSC_DLLEXPORT PetscObjectComposedDataIncreaseInt(PetscObject obj)
88 {
89   PetscInt       *ar = obj->intcomposeddata,*new_ar;
90   PetscInt       *ir = obj->intcomposedstate,*new_ir,n = obj->int_idmax,new_n,i;
91   PetscErrorCode ierr;
92 
93   PetscFunctionBegin;
94   new_n = globalmaxstate;
95   ierr = PetscMalloc(new_n*sizeof(PetscInt),&new_ar);CHKERRQ(ierr);
96   ierr = PetscMemzero(new_ar,new_n*sizeof(PetscInt));CHKERRQ(ierr);
97   ierr = PetscMalloc(new_n*sizeof(PetscInt),&new_ir);CHKERRQ(ierr);
98   ierr = PetscMemzero(new_ir,new_n*sizeof(PetscInt));CHKERRQ(ierr);
99   if (n) {
100     for (i=0; i<n; i++) {
101       new_ar[i] = ar[i]; new_ir[i] = ir[i];
102     }
103     ierr = PetscFree(ar);CHKERRQ(ierr);
104     ierr = PetscFree(ir);CHKERRQ(ierr);
105   }
106   obj->int_idmax = new_n;
107   obj->intcomposeddata = new_ar; obj->intcomposedstate = new_ir;
108   PetscFunctionReturn(0);
109 }
110 PetscErrorCode PETSC_DLLEXPORT PetscObjectComposedDataIncreaseIntstar(PetscObject obj)
111 {
112   PetscInt       **ar = obj->intstarcomposeddata,**new_ar;
113   PetscInt       *ir = obj->intstarcomposedstate,*new_ir,n = obj->intstar_idmax,new_n,i;
114   PetscErrorCode ierr;
115 
116   PetscFunctionBegin;
117   new_n = globalmaxstate;
118   ierr = PetscMalloc(new_n*sizeof(PetscInt*),&new_ar);CHKERRQ(ierr);
119   ierr = PetscMemzero(new_ar,new_n*sizeof(PetscInt*));CHKERRQ(ierr);
120   ierr = PetscMalloc(new_n*sizeof(PetscInt),&new_ir);CHKERRQ(ierr);
121   ierr = PetscMemzero(new_ir,new_n*sizeof(PetscInt));CHKERRQ(ierr);
122   if (n) {
123     for (i=0; i<n; i++) {
124       new_ar[i] = ar[i]; new_ir[i] = ir[i];
125     }
126     ierr = PetscFree(ar);CHKERRQ(ierr);
127     ierr = PetscFree(ir);CHKERRQ(ierr);
128   }
129   obj->intstar_idmax = new_n;
130   obj->intstarcomposeddata = new_ar; obj->intstarcomposedstate = new_ir;
131   PetscFunctionReturn(0);
132 }
133 
134 PetscErrorCode PETSC_DLLEXPORT PetscObjectComposedDataIncreaseReal(PetscObject obj)
135 {
136   PetscReal      *ar = obj->realcomposeddata,*new_ar;
137   PetscInt       *ir = obj->realcomposedstate,*new_ir,n = obj->real_idmax,new_n,i;
138   PetscErrorCode ierr;
139 
140   PetscFunctionBegin;
141   new_n = globalmaxstate;
142   ierr = PetscMalloc(new_n*sizeof(PetscReal),&new_ar);CHKERRQ(ierr);
143   ierr = PetscMemzero(new_ar,new_n*sizeof(PetscReal));CHKERRQ(ierr);
144   ierr = PetscMalloc(new_n*sizeof(PetscInt),&new_ir);CHKERRQ(ierr);
145   ierr = PetscMemzero(new_ir,new_n*sizeof(PetscInt));CHKERRQ(ierr);
146   if (n) {
147     for (i=0; i<n; i++) {
148       new_ar[i] = ar[i]; new_ir[i] = ir[i];
149     }
150     ierr = PetscFree(ar);CHKERRQ(ierr);
151     ierr = PetscFree(ir);CHKERRQ(ierr);
152   }
153   obj->real_idmax = new_n;
154   obj->realcomposeddata = new_ar; obj->realcomposedstate = new_ir;
155   PetscFunctionReturn(0);
156 }
157 
158 PetscErrorCode PETSC_DLLEXPORT PetscObjectComposedDataIncreaseRealstar(PetscObject obj)
159 {
160   PetscReal      **ar = obj->realstarcomposeddata,**new_ar;
161   PetscInt       *ir = obj->realstarcomposedstate,*new_ir,n = obj->realstar_idmax,new_n,i;
162   PetscErrorCode ierr;
163 
164   PetscFunctionBegin;
165   new_n = globalmaxstate;
166   ierr = PetscMalloc(new_n*sizeof(PetscReal*),&new_ar);CHKERRQ(ierr);
167   ierr = PetscMemzero(new_ar,new_n*sizeof(PetscReal*));CHKERRQ(ierr);
168   ierr = PetscMalloc(new_n*sizeof(PetscInt),&new_ir);CHKERRQ(ierr);
169   ierr = PetscMemzero(new_ir,new_n*sizeof(PetscInt));CHKERRQ(ierr);
170   if (n) {
171     for (i=0; i<n; i++) {
172       new_ar[i] = ar[i]; new_ir[i] = ir[i];
173     }
174     ierr = PetscFree(ar);CHKERRQ(ierr);
175     ierr = PetscFree(ir);CHKERRQ(ierr);
176   }
177   obj->realstar_idmax = new_n;
178   obj->realstarcomposeddata = new_ar; obj->realstarcomposedstate = new_ir;
179   PetscFunctionReturn(0);
180 }
181 
182 PetscErrorCode PETSC_DLLEXPORT PetscObjectComposedDataIncreaseScalar(PetscObject obj)
183 {
184   PetscScalar    *ar = obj->scalarcomposeddata,*new_ar;
185   PetscInt       *ir = obj->scalarcomposedstate,*new_ir,n = obj->scalar_idmax,new_n,i;
186   PetscErrorCode ierr;
187 
188   PetscFunctionBegin;
189   new_n = globalmaxstate;
190   ierr = PetscMalloc(new_n*sizeof(PetscScalar),&new_ar);CHKERRQ(ierr);
191   ierr = PetscMemzero(new_ar,new_n*sizeof(PetscScalar));CHKERRQ(ierr);
192   ierr = PetscMalloc(new_n*sizeof(PetscInt),&new_ir);CHKERRQ(ierr);
193   ierr = PetscMemzero(new_ir,new_n*sizeof(PetscInt));CHKERRQ(ierr);
194   if (n) {
195     for (i=0; i<n; i++) {
196       new_ar[i] = ar[i]; new_ir[i] = ir[i];
197     }
198     ierr = PetscFree(ar);CHKERRQ(ierr);
199     ierr = PetscFree(ir);CHKERRQ(ierr);
200   }
201   obj->scalar_idmax = new_n;
202   obj->scalarcomposeddata = new_ar; obj->scalarcomposedstate = new_ir;
203   PetscFunctionReturn(0);
204 }
205 
206 PetscErrorCode PETSC_DLLEXPORT PetscObjectComposedDataIncreaseScalarstar(PetscObject obj)
207 {
208   PetscScalar    **ar = obj->scalarstarcomposeddata,**new_ar;
209   PetscInt       *ir = obj->scalarstarcomposedstate,*new_ir,n = obj->scalarstar_idmax,new_n,i;
210   PetscErrorCode ierr;
211 
212   PetscFunctionBegin;
213   new_n = globalmaxstate;
214   ierr = PetscMalloc(new_n*sizeof(PetscScalar*),&new_ar);CHKERRQ(ierr);
215   ierr = PetscMemzero(new_ar,new_n*sizeof(PetscScalar*));CHKERRQ(ierr);
216   ierr = PetscMalloc(new_n*sizeof(PetscInt),&new_ir);CHKERRQ(ierr);
217   ierr = PetscMemzero(new_ir,new_n*sizeof(PetscInt));CHKERRQ(ierr);
218   if (n) {
219     for (i=0; i<n; i++) {
220       new_ar[i] = ar[i]; new_ir[i] = ir[i];
221     }
222     ierr = PetscFree(ar);CHKERRQ(ierr);
223     ierr = PetscFree(ir);CHKERRQ(ierr);
224   }
225   obj->scalarstar_idmax = new_n;
226   obj->scalarstarcomposeddata = new_ar; obj->scalarstarcomposedstate = new_ir;
227   PetscFunctionReturn(0);
228 }
229 
230