Actual source code: petscdmplex.h90

petsc-3.11.0 2019-03-29
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,l,nC,nD,nB,bF,bC,bP,pm,sc,e)
193:           use petscdmdef
194:           PetscSection sc
195:           PetscInt     nB
196:           DMLabel,  pointer :: l(:)
197:           PetscInt, pointer :: nC(:)
198:           PetscInt, pointer :: nD(:)
199:           PetscInt, pointer :: bF(:)
200:           IS,       pointer :: bC(:)
201:           IS,       pointer :: bP(:)
202:           IS        pm
203:           PetscErrorCode e
204:           DM m
205:         End Subroutine
206:       End Interface

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

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

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

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

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

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

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

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

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

321:       Interface
322:         Subroutine DMPlexCreateFromFile(c,str,i,m,ierr)
323:           use petscdmdef
324:           MPI_Comm :: c
325:           character(len=*) :: str
326:           PetscBool, intent(in) :: i
327:           DM, intent(out) :: m
328:           PetscErrorCode, intent(out):: ierr
329:         End Subroutine
330:       End Interface

332:       Interface
333:         Subroutine DMPlexDistribute(m,o,sf,mp,ierr)
334:           use petscdmdef
335:           DM, intent(in) :: m
336:           PetscInt, intent(in) :: o
337:           PetscSF, intent(out) :: sf
338:           DM, intent(out) :: mp
339:           PetscErrorCode, intent(out):: ierr
340:         End Subroutine
341:       End Interface

343:       Interface
344:         Subroutine DMPlexCreateDefault(m,ierr)
345:           use petscdmdef
346:           DM, intent(out) :: m
347:           PetscErrorCode, intent(out):: ierr
348:         End Subroutine
349:       End Interface