Actual source code: petscdmplex.h90

petsc-3.9.4 2018-09-11
Report Typos and Errors

  2:       Interface
  3:         Subroutine DMPlexGetCone(m,p,cone,ierr)
  4:           use petscdmdef
  5:           PetscInt     p
  6:           PetscInt, pointer :: cone(:)
  7:           PetscErrorCode ierr
  8:           DM m
  9:         End Subroutine
 10:       End Interface

 12:       Interface
 13:         Subroutine DMPlexRestoreCone(m,p,cone,ierr)
 14:           use petscdmdef
 15:           PetscInt     p
 16:           PetscInt, pointer :: cone(:)
 17:           PetscErrorCode ierr
 18:           DM m
 19:         End Subroutine
 20:       End Interface

 22:       Interface
 23:         Subroutine DMPlexGetConeOrientation(m,p,coneOrient,ierr)
 24:           use petscdmdef
 25:           PetscInt     p
 26:           PetscInt, pointer :: coneOrient(:)
 27:           PetscErrorCode ierr
 28:           DM m
 29:         End Subroutine
 30:       End Interface

 32:       Interface
 33:         Subroutine DMPlexGetSupport(m,p,support,ierr)
 34:           use petscdmdef
 35:           PetscInt     p
 36:           PetscInt, pointer :: support(:)
 37:           PetscErrorCode ierr
 38:           DM m
 39:         End Subroutine
 40:       End Interface

 42:       Interface
 43:         Subroutine DMPlexRestoreSupport(m,p,support,ierr)
 44:           use petscdmdef
 45:           PetscInt     p
 46:           PetscInt, pointer :: support(:)
 47:           PetscErrorCode ierr
 48:           DM m
 49:         End Subroutine
 50:       End Interface

 52:       Interface
 53:         Subroutine DMPlexGetTransitiveClosure(m,p,useCone,points,ierr)
 54:           use petscdmdef
 55:           PetscInt     p
 56:           PetscBool    useCone
 57:           PetscInt, pointer :: points(:)
 58:           PetscErrorCode ierr
 59:           DM m
 60:         End Subroutine
 61:       End Interface

 63:       Interface
 64:         Subroutine DMPlexRestoreTransitiveClosure(m,p,uC,points,ierr)
 65:           use petscdmdef
 66:           PetscInt     p
 67:           PetscBool    uC
 68:           PetscInt, pointer :: points(:)
 69:           PetscErrorCode ierr
 70:           DM m
 71:         End Subroutine
 72:       End Interface

 74:       Interface
 75:         Subroutine DMPlexGetJoin(m,numPoints,points,join,ierr)
 76:           use petscdmdef
 77:           PetscInt     numPoints
 78:           PetscInt, pointer :: points(:)
 79:           PetscInt, pointer :: join(:)
 80:           PetscErrorCode ierr
 81:           DM m
 82:         End Subroutine
 83:       End Interface

 85:       Interface
 86:         Subroutine DMPlexGetFullJoin(m,numPoints,points,join,ierr)
 87:           use petscdmdef
 88:           PetscInt     numPoints
 89:           PetscInt, pointer :: points(:)
 90:           PetscInt, pointer :: join(:)
 91:           PetscErrorCode ierr
 92:           DM m
 93:         End Subroutine
 94:       End Interface

 96:       Interface
 97:         Subroutine DMPlexRestoreJoin(m,numPoints,points,join,ierr)
 98:           use petscdmdef
 99:           PetscInt     numPoints
100:           PetscInt, pointer :: points(:)
101:           PetscInt, pointer :: join(:)
102:           PetscErrorCode ierr
103:           DM m
104:         End Subroutine
105:       End Interface

107:       Interface
108:         Subroutine DMPlexGetMeet(m,numPoints,points,meet,ierr)
109:           use petscdmdef
110:           PetscInt     numPoints
111:           PetscInt, pointer :: points(:)
112:           PetscInt, pointer :: meet(:)
113:           PetscErrorCode ierr
114:           DM m
115:         End Subroutine
116:       End Interface

118:       Interface
119:         Subroutine DMPlexGetFullMeet(m,numPoints,points,meet,ierr)
120:           use petscdmdef
121:           PetscInt     numPoints
122:           PetscInt, pointer :: points(:)
123:           PetscInt, pointer :: meet(:)
124:           PetscErrorCode ierr
125:           DM m
126:         End Subroutine
127:       End Interface

129:       Interface
130:         Subroutine DMPlexRestoreMeet(m,numPoints,points,meet,ierr)
131:           use petscdmdef
132:           PetscInt     numPoints
133:           PetscInt, pointer :: points(:)
134:           PetscInt, pointer :: meet(:)
135:           PetscErrorCode ierr
136:           DM m
137:         End Subroutine
138:       End Interface

140:       Interface
141:         Subroutine DMPlexVecGetClosure(m,section,v,point,values,ierr)
142:           use petscdmdef
143:           PetscSection section
144:           Vec     v
145:           PetscInt     point
146:           PetscScalar, pointer :: values(:)
147:           PetscErrorCode ierr
148:           DM m
149:         End Subroutine
150:       End Interface

152:       Interface
153:         Subroutine DMPlexVecRestoreClosure(m,section,v,point,vs,ierr)
154:           use petscdmdef
155:           PetscSection section
156:           Vec     v
157:           PetscInt     point
158:           PetscScalar, pointer :: vs(:)
159:           PetscErrorCode ierr
160:           DM m
161:         End Subroutine
162:       End Interface

164:       Interface
165:         Subroutine DMPlexVecSetClosure(m,sec,v,point,vs,mode,ierr)
166:           use petscdmdef
167:           PetscSection sec
168:           Vec     v
169:           PetscInt     point
170:           InsertMode   mode
171:           PetscScalar, pointer :: vs(:)
172:           PetscErrorCode ierr
173:           DM m
174:         End Subroutine
175:       End Interface

177:       Interface
178:         Subroutine DMPlexMatSetClosure(m,s,gS,A,p,v,mode,ierr)
179:           use petscdmdef
180:           PetscSection s
181:           PetscSection gS
182:           Mat     A
183:           PetscInt     p
184:           InsertMode   mode
185:           PetscScalar, pointer :: v(:)
186:           PetscErrorCode ierr
187:           DM m
188:         End Subroutine
189:       End Interface

191:       Interface
192:         Subroutine DMPlexCreateSection(m,d,nF,nC,nD,nB,bF,bC,bP,pm,sc,e)
193:           use petscdmdef
194:           PetscSection sc
195:           PetscInt     d, nF, nB
196:           PetscInt, pointer :: nC(:)
197:           PetscInt, pointer :: nD(:)
198:           PetscInt, pointer :: bF(:)
199:           IS,  pointer :: bC(:)
200:           IS,  pointer :: bP(:)
201:           IS pm
202:           PetscErrorCode e
203:           DM m
204:         End Subroutine
205:       End Interface

207:       Interface
208:         Subroutine DMPlexComputeCellGeometryAffineFEM(m,c,v0,J,iJ,dJ,er)
209:           use petscdmdef
210:           PetscInt   c
211:           PetscReal, pointer :: v0(:)
212:           PetscReal, pointer :: J(:)
213:           PetscReal, pointer :: iJ(:)
214:           PetscReal  dJ
215:           PetscErrorCode er
216:           DM m
217:         End Subroutine
218:       End Interface

220:       Interface
221:         Subroutine DMPlexComputeCellGeometryFEM(m,c,fe,v0,J,iJ,dJ,er)
222:           use petscdmdef
223:           PetscInt   c
224:           PetscReal, pointer :: v0(:)
225:           PetscReal, pointer :: J(:)
226:           PetscReal, pointer :: iJ(:)
227:           PetscReal  dJ
228:           PetscErrorCode er
229:           PetscFE fe
230:           DM m
231:         End Subroutine
232:       End Interface

234:       Interface
235:         Subroutine DMPlexComputeCellGeometryFVM(m,cell,vol,ct,nm,ierr)
236:           use petscdmdef
237:           PetscInt   cell
238:           PetscReal  vol
239:           PetscReal, pointer :: ct(:)
240:           PetscReal, pointer :: nm(:)
241:           PetscErrorCode ierr
242:           DM m
243:         End Subroutine
244:       End Interface

246:       Interface
247:         Subroutine DMPlexGetCellFields(m,s,e,x,xt,a,u,ut,v,ierr)
248:           use petscdmdef
249:           PetscInt  s, e
250:           Vec  x, xt, a
251:           PetscScalar, pointer :: u(:)
252:           PetscScalar, pointer :: ut(:)
253:           PetscScalar, pointer :: v(:)
254:           PetscErrorCode ierr
255:           DM m
256:         End Subroutine
257:       End Interface

259:       Interface
260:         Subroutine DMPlexRestoreCellFields(m,s,e,x,xt,a,u,ut,v,ierr)
261:           use petscdmdef
262:           PetscInt  s, e
263:           Vec  x, xt, a
264:           PetscScalar, pointer :: u(:)
265:           PetscScalar, pointer :: ut(:)
266:           PetscScalar, pointer :: v(:)
267:           PetscErrorCode ierr
268:           DM m
269:         End Subroutine
270:       End Interface

272:       Interface
273:         Subroutine DMPlexGetFaceFields(m,s,e,x,xt,f,c,g,uL,uR,ierr)
274:           use petscdmdef
275:           PetscInt  s, e
276:           Vec  x, xt, f, c, g
277:           PetscScalar, pointer :: uL(:)
278:           PetscScalar, pointer :: uR(:)
279:           PetscErrorCode ierr
280:           DM m
281:         End Subroutine
282:       End Interface

284:       Interface
285:         Subroutine DMPlexRestoreFaceFields(m,s,e,x,xt,f,c,g,uL,uR,ierr)
286:           use petscdmdef
287:           PetscInt  s, e
288:           Vec  x, xt, f, c, g
289:           PetscScalar, pointer :: uL(:)
290:           PetscScalar, pointer :: uR(:)
291:           PetscErrorCode ierr
292:           DM m
293:         End Subroutine
294:       End Interface

296:       Interface
297:         Subroutine DMPlexGetFaceGeometry(m,s,e,f,c,g,v,ierr)
298:           use petscdmdef
299:           PetscInt  s, e
300:           Vec  f, c
301:           PetscScalar, pointer :: g(:)
302:           PetscScalar, pointer :: v(:)
303:           PetscErrorCode ierr
304:           DM m
305:         End Subroutine
306:       End Interface

308:       Interface
309:         Subroutine DMPlexRestoreFaceGeometry(m,s,e,f,c,g,v,ierr)
310:           use petscdmdef
311:           PetscInt  s, e
312:           Vec  f, c
313:           PetscScalar, pointer :: g(:)
314:           PetscScalar, pointer :: v(:)
315:           PetscErrorCode ierr
316:           DM m
317:         End Subroutine
318:       End Interface