xref: /petsc/src/sys/objects/olist.c (revision 6bd3a4fd3ada34738ab3c62a7d3b0fa06ce8b47c) !
1 
2 /*
3          Provides a general mechanism to maintain a linked list of PETSc objects.
4      This is used to allow PETSc objects to carry a list of "composed" objects
5 */
6 #include <petsc/private/petscimpl.h>
7 
8 /*@C
9   PetscObjectListRemoveReference - Calls `PetscObjectDereference()` on an object in the list immediately but keeps a pointer to the object in the list.
10 
11   Input Parameters:
12 + fl   - the object list
13 - name - the name to use for the object
14 
15   Level: developer
16 
17   Notes:
18   Use `PetscObjectListAdd`(`PetscObjectList`,const char name[],NULL) to truly remove the object from the list
19 
20   Use this routine ONLY if you know that the object referenced will remain in existence until the pointing object is destroyed
21 
22   Developer Notes:
23   This is to handle some cases that otherwise would result in having circular references so reference counts never got to zero
24 
25 .seealso: `PetscObjectListDestroy()`,`PetscObjectListFind()`,`PetscObjectListDuplicate()`,`PetscObjectListReverseFind()`,
26 `PetscObject`, `PetscObjectListAdd()`
27 @*/
28 PetscErrorCode PetscObjectListRemoveReference(PetscObjectList *fl, const char name[])
29 {
30   PetscObjectList nlist;
31   PetscBool       match;
32 
33   PetscFunctionBegin;
34   PetscValidPointer(fl, 1);
35   PetscValidCharPointer(name, 2);
36   nlist = *fl;
37   while (nlist) {
38     PetscCall(PetscStrcmp(name, nlist->name, &match));
39     if (match) { /* found it in the list */
40       if (!nlist->skipdereference) PetscCall(PetscObjectDereference(nlist->obj));
41       nlist->skipdereference = PETSC_TRUE;
42       PetscFunctionReturn(PETSC_SUCCESS);
43     }
44     nlist = nlist->next;
45   }
46   PetscFunctionReturn(PETSC_SUCCESS);
47 }
48 
49 /*@C
50   PetscObjectListAdd - Adds a new object to an `PetscObjectList`
51 
52   Input Parameters:
53 + fl   - the object list
54 . name - the name to use for the object
55 - obj  - the object to attach
56 
57   Level: developer
58 
59   Notes:
60   Replaces item if it is already in list. Removes item if you pass in a `NULL` object.
61 
62   Use `PetscObjectListFind()` or `PetscObjectListReverseFind()` to get the object back
63 
64 .seealso: `PetscObjectListDestroy()`,`PetscObjectListFind()`,`PetscObjectListDuplicate()`,`PetscObjectListReverseFind()`, `PetscObject`, `PetscObjectList`
65 @*/
66 PetscErrorCode PetscObjectListAdd(PetscObjectList *fl, const char name[], PetscObject obj)
67 {
68   PetscObjectList olist, nlist, prev;
69   PetscBool       match;
70 
71   PetscFunctionBegin;
72   PetscValidPointer(fl, 1);
73   if (!obj) { /* this means remove from list if it is there */
74     nlist = *fl;
75     prev  = NULL;
76     while (nlist) {
77       PetscCall(PetscStrcmp(name, nlist->name, &match));
78       if (match) { /* found it already in the list */
79         /* Remove it first to prevent circular derefs */
80         if (prev) prev->next = nlist->next;
81         else if (nlist->next) *fl = nlist->next;
82         else *fl = NULL;
83         if (!nlist->skipdereference) PetscCall(PetscObjectDereference(nlist->obj));
84         PetscCall(PetscFree(nlist));
85         PetscFunctionReturn(PETSC_SUCCESS);
86       }
87       prev  = nlist;
88       nlist = nlist->next;
89     }
90     PetscFunctionReturn(PETSC_SUCCESS); /* did not find it to remove */
91   }
92   /* look for it already in list */
93   nlist = *fl;
94   while (nlist) {
95     PetscCall(PetscStrcmp(name, nlist->name, &match));
96     if (match) { /* found it in the list */
97       PetscCall(PetscObjectReference(obj));
98       if (!nlist->skipdereference) PetscCall(PetscObjectDereference(nlist->obj));
99       nlist->skipdereference = PETSC_FALSE;
100       nlist->obj             = obj;
101       PetscFunctionReturn(PETSC_SUCCESS);
102     }
103     nlist = nlist->next;
104   }
105 
106   /* add it to list, because it was not already there */
107   PetscCall(PetscNew(&olist));
108   olist->next = NULL;
109   olist->obj  = obj;
110 
111   PetscCall(PetscObjectReference(obj));
112   PetscCall(PetscStrncpy(olist->name, name, sizeof(olist->name)));
113 
114   if (!*fl) *fl = olist;
115   else { /* go to end of list */ nlist = *fl;
116     while (nlist->next) nlist = nlist->next;
117     nlist->next = olist;
118   }
119   PetscFunctionReturn(PETSC_SUCCESS);
120 }
121 
122 /*@C
123   PetscObjectListDestroy - Destroy a list of objects
124 
125   Input Parameter:
126 . ifl - pointer to list
127 
128   Level: developer
129 
130 .seealso: `PetscObjectList`, `PetscObject`, `PetscObjectListAdd()`, `PetscObjectListFind()`, `PetscObjectListDuplicate()`,
131           `PetscObjectListReverseFind()`
132 @*/
133 PetscErrorCode PetscObjectListDestroy(PetscObjectList *ifl)
134 {
135   PetscObjectList tmp, fl;
136 
137   PetscFunctionBegin;
138   PetscValidPointer(ifl, 1);
139   fl = *ifl;
140   while (fl) {
141     tmp = fl->next;
142     if (!fl->skipdereference) PetscCall(PetscObjectDereference(fl->obj));
143     PetscCall(PetscFree(fl));
144     fl = tmp;
145   }
146   *ifl = NULL;
147   PetscFunctionReturn(PETSC_SUCCESS);
148 }
149 
150 /*@C
151   PetscObjectListFind - given a name, find the matching object in a list
152 
153   Input Parameters:
154 + fl   - pointer to list
155 - name - name string
156 
157   Output Parameter:
158 . obj - the PETSc object
159 
160   Level: developer
161 
162   Notes:
163   The name must have been registered with the `PetscObjectListAdd()` before calling this routine.
164 
165   The reference count of the object is not increased
166 
167 .seealso: `PetscObjectListDestroy()`,`PetscObjectListAdd()`,`PetscObjectListDuplicate()`,`PetscObjectListReverseFind()`, `PetscObjectList`
168 @*/
169 PetscErrorCode PetscObjectListFind(PetscObjectList fl, const char name[], PetscObject *obj)
170 {
171   PetscFunctionBegin;
172   PetscValidPointer(obj, 3);
173   *obj = NULL;
174   while (fl) {
175     PetscBool match;
176     PetscCall(PetscStrcmp(name, fl->name, &match));
177     if (match) {
178       *obj = fl->obj;
179       break;
180     }
181     fl = fl->next;
182   }
183   PetscFunctionReturn(PETSC_SUCCESS);
184 }
185 
186 /*@C
187   PetscObjectListReverseFind - given a object, find the matching name if it exists
188 
189   Input Parameters:
190 + fl  - pointer to list
191 - obj - the PETSc object
192 
193   Output Parameters:
194 + name            - name string
195 - skipdereference - if the object is in list but does not have the increased reference count for a circular dependency
196 
197   Level: developer
198 
199   Notes:
200   The name must have been registered with the `PetscObjectListAdd()` before calling this routine.
201 
202   The reference count of the object is not increased
203 
204 .seealso: `PetscObjectListDestroy()`,`PetscObjectListAdd()`,`PetscObjectListDuplicate()`,`PetscObjectListFind()`, `PetscObjectList`
205 @*/
206 PetscErrorCode PetscObjectListReverseFind(PetscObjectList fl, PetscObject obj, char **name, PetscBool *skipdereference)
207 {
208   PetscFunctionBegin;
209   PetscValidPointer(name, 3);
210   if (skipdereference) PetscValidBoolPointer(skipdereference, 4);
211   *name = NULL;
212   while (fl) {
213     if (fl->obj == obj) {
214       *name = fl->name;
215       if (skipdereference) *skipdereference = fl->skipdereference;
216       break;
217     }
218     fl = fl->next;
219   }
220   PetscFunctionReturn(PETSC_SUCCESS);
221 }
222 
223 /*@C
224   PetscObjectListDuplicate - Creates a new list from a given object list.
225 
226   Input Parameter:
227 . fl - pointer to list
228 
229   Output Parameter:
230 . nl - the new list (should point to `NULL` to start, otherwise appends)
231 
232   Level: developer
233 
234 .seealso: `PetscObjectListDestroy()`, `PetscObjectListAdd()`, `PetscObjectListReverseFind()`,
235 `PetscObjectListFind()`, `PetscObjectListDuplicate()`, `PetscObjectList`
236 @*/
237 PetscErrorCode PetscObjectListDuplicate(PetscObjectList fl, PetscObjectList *nl)
238 {
239   PetscFunctionBegin;
240   PetscValidPointer(nl, 2);
241   while (fl) {
242     PetscCall(PetscObjectListAdd(nl, fl->name, fl->obj));
243     fl = fl->next;
244   }
245   PetscFunctionReturn(PETSC_SUCCESS);
246 }
247