xref: /petsc/src/sys/objects/olist.c (revision 1d0fab5e66538303e50d41f342d834e538cdb59d)
1 #define PETSC_DLL
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 "petscsys.h"
7 
8 struct _n_PetscOList {
9     char        name[256];
10     PetscObject obj;
11     PetscOList  next;
12 };
13 
14 #undef __FUNCT__
15 #define __FUNCT__ "PetscOListAdd"
16 /*@C
17      PetscOListAdd - Adds a new object to an PetscOList
18 
19     Input Parameters:
20 +     fl - the object list
21 .     name - the name to use for the object
22 -     obj - the object to attach
23 
24        Notes: Replaces item if it is already in list. Removes item if you pass in a PETSC_NULL object.
25 
26         Use PetscOListFind() or PetscOListReverseFind() to get the object back
27 
28 .seealso: PetscOListDestroy(), PetscOListFind(), PetscOListDuplicate(), PetscOListReverseFind(), PetscOListDuplicate()
29 
30 @*/
31 PetscErrorCode PETSC_DLLEXPORT PetscOListAdd(PetscOList *fl,const char name[],PetscObject obj)
32 {
33   PetscOList     olist,nlist,prev;
34   PetscErrorCode ierr;
35   PetscTruth     match;
36 
37   PetscFunctionBegin;
38 
39   if (!obj) { /* this means remove from list if it is there */
40     nlist = *fl; prev = 0;
41     while (nlist) {
42       ierr = PetscStrcmp(name,nlist->name,&match);CHKERRQ(ierr);
43       if (match) {  /* found it already in the list */
44         ierr = PetscObjectDereference(nlist->obj);CHKERRQ(ierr);
45         if (prev) prev->next = nlist->next;
46         else if (nlist->next) {
47           *fl = nlist->next;
48         } else {
49           *fl = 0;
50         }
51         ierr = PetscFree(nlist);CHKERRQ(ierr);
52         PetscFunctionReturn(0);
53       }
54       prev  = nlist;
55       nlist = nlist->next;
56     }
57     PetscFunctionReturn(0); /* did not find it to remove */
58   }
59   /* look for it already in list */
60   nlist = *fl;
61   while (nlist) {
62     ierr = PetscStrcmp(name,nlist->name,&match);CHKERRQ(ierr);
63     if (match) {  /* found it in the list */
64       ierr = PetscObjectReference(obj);CHKERRQ(ierr);
65       ierr = PetscObjectDereference(nlist->obj);CHKERRQ(ierr);
66       nlist->obj = obj;
67       PetscFunctionReturn(0);
68     }
69     nlist = nlist->next;
70   }
71 
72   /* add it to list, because it was not already there */
73 
74   ierr        = PetscNew(struct _n_PetscOList,&olist);CHKERRQ(ierr);
75   olist->next = 0;
76   olist->obj  = obj;
77   ierr = PetscObjectReference(obj);CHKERRQ(ierr);
78   ierr = PetscStrcpy(olist->name,name);CHKERRQ(ierr);
79 
80   if (!*fl) {
81     *fl = olist;
82   } else { /* go to end of list */
83     nlist = *fl;
84     while (nlist->next) {
85       nlist = nlist->next;
86     }
87     nlist->next = olist;
88   }
89   PetscFunctionReturn(0);
90 }
91 
92 #undef __FUNCT__
93 #define __FUNCT__ "PetscOListDestroy"
94 /*@C
95     PetscOListDestroy - Destroy a list of objects
96 
97     Input Parameter:
98 .   fl   - pointer to list
99 
100 .seealso: PetscOListAdd(), PetscOListFind(), PetscOListDuplicate(), PetscOListReverseFind(), PetscOListDuplicate()
101 
102 @*/
103 PetscErrorCode PETSC_DLLEXPORT PetscOListDestroy(PetscOList fl)
104 {
105   PetscOList     tmp;
106   PetscErrorCode ierr;
107 
108   PetscFunctionBegin;
109   while (fl) {
110     tmp   = fl->next;
111     ierr  = PetscObjectDereference(fl->obj);CHKERRQ(ierr);
112     ierr  = PetscFree(fl);CHKERRQ(ierr);
113     fl    = tmp;
114   }
115   PetscFunctionReturn(0);
116 }
117 
118 
119 #undef __FUNCT__
120 #define __FUNCT__ "PetscOListFind"
121 /*@C
122     PetscOListFind - givn a name, find the matching object
123 
124     Input Parameters:
125 +   fl   - pointer to list
126 -   name - name string
127 
128     Output Parameters:
129 .   ob - the PETSc object
130 
131     Notes:
132     The name must have been registered with the PetscOListAdd() before calling this
133     routine.
134 
135 .seealso: PetscOListDestroy(), PetscOListAdd(), PetscOListDuplicate(), PetscOListReverseFind(), PetscOListDuplicate()
136 
137 @*/
138 PetscErrorCode PETSC_DLLEXPORT PetscOListFind(PetscOList fl,const char name[],PetscObject *obj)
139 {
140   PetscErrorCode ierr;
141   PetscTruth     match;
142 
143   PetscFunctionBegin;
144 
145   *obj = 0;
146   while (fl) {
147     ierr = PetscStrcmp(name,fl->name,&match);CHKERRQ(ierr);
148     if (match) {
149       *obj = fl->obj;
150       break;
151     }
152     fl = fl->next;
153   }
154   PetscFunctionReturn(0);
155 }
156 
157 #undef __FUNCT__
158 #define __FUNCT__ "PetscOListReverseFind"
159 /*@C
160     PetscOListReverseFind - given a object, find the matching name if it exists
161 
162     Input Parameters:
163 +   fl   - pointer to list
164 -   ob - the PETSc object
165 
166     Output Parameters:
167 .   name - name string
168 
169     Notes:
170     The name must have been registered with the PetscOListAdd() before calling this
171     routine.
172 
173 .seealso: PetscOListDestroy(), PetscOListAdd(), PetscOListDuplicate(), PetscOListFind(), PetscOListDuplicate()
174 
175 @*/
176 PetscErrorCode PETSC_DLLEXPORT PetscOListReverseFind(PetscOList fl,PetscObject obj,char **name)
177 {
178   PetscFunctionBegin;
179 
180   *name = 0;
181   while (fl) {
182     if (fl->obj == obj) {
183       *name = fl->name;
184       break;
185     }
186     fl = fl->next;
187   }
188   PetscFunctionReturn(0);
189 }
190 
191 #undef __FUNCT__
192 #define __FUNCT__ "PetscOListDuplicate"
193 /*@C
194     PetscOListDuplicate - Creates a new list from a give object list.
195 
196     Input Parameters:
197 .   fl   - pointer to list
198 
199     Output Parameters:
200 .   nl - the new list (should point to 0 to start, otherwise appends)
201 
202 .seealso: PetscOListDestroy(), PetscOListAdd(), PetscOListReverseFind(), PetscOListFind(), PetscOListDuplicate()
203 
204 @*/
205 PetscErrorCode PETSC_DLLEXPORT PetscOListDuplicate(PetscOList fl,PetscOList *nl)
206 {
207   PetscErrorCode ierr;
208 
209   PetscFunctionBegin;
210   while (fl) {
211     ierr = PetscOListAdd(nl,fl->name,fl->obj);CHKERRQ(ierr);
212     fl = fl->next;
213   }
214   PetscFunctionReturn(0);
215 }
216 
217 
218 
219 
220 
221