xref: /petsc/src/snes/utils/dmplexsnes.c (revision a95d84e0230cb0d652d808aead973c41ec535430)
1 #include <petscdmplex.h> /*I "petscdmplex.h" I*/
2 #include <petscsnes.h>      /*I "petscsnes.h" I*/
3 #include <petsc-private/petscimpl.h>
4 
5 #undef __FUNCT__
6 #define __FUNCT__ "DMInterpolationCreate"
7 PetscErrorCode DMInterpolationCreate(MPI_Comm comm, DMInterpolationInfo *ctx)
8 {
9   PetscErrorCode ierr;
10 
11   PetscFunctionBegin;
12   PetscValidPointer(ctx, 2);
13   ierr = PetscMalloc(sizeof(struct _DMInterpolationInfo), ctx);CHKERRQ(ierr);
14 
15   (*ctx)->comm   = comm;
16   (*ctx)->dim    = -1;
17   (*ctx)->nInput = 0;
18   (*ctx)->points = NULL;
19   (*ctx)->cells  = NULL;
20   (*ctx)->n      = -1;
21   (*ctx)->coords = NULL;
22   PetscFunctionReturn(0);
23 }
24 
25 #undef __FUNCT__
26 #define __FUNCT__ "DMInterpolationSetDim"
27 PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo ctx, PetscInt dim)
28 {
29   PetscFunctionBegin;
30   if ((dim < 1) || (dim > 3)) SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension for points: %d", dim);
31   ctx->dim = dim;
32   PetscFunctionReturn(0);
33 }
34 
35 #undef __FUNCT__
36 #define __FUNCT__ "DMInterpolationGetDim"
37 PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo ctx, PetscInt *dim)
38 {
39   PetscFunctionBegin;
40   PetscValidIntPointer(dim, 2);
41   *dim = ctx->dim;
42   PetscFunctionReturn(0);
43 }
44 
45 #undef __FUNCT__
46 #define __FUNCT__ "DMInterpolationSetDof"
47 PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo ctx, PetscInt dof)
48 {
49   PetscFunctionBegin;
50   if (dof < 1) SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of components: %d", dof);
51   ctx->dof = dof;
52   PetscFunctionReturn(0);
53 }
54 
55 #undef __FUNCT__
56 #define __FUNCT__ "DMInterpolationGetDof"
57 PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo ctx, PetscInt *dof)
58 {
59   PetscFunctionBegin;
60   PetscValidIntPointer(dof, 2);
61   *dof = ctx->dof;
62   PetscFunctionReturn(0);
63 }
64 
65 #undef __FUNCT__
66 #define __FUNCT__ "DMInterpolationAddPoints"
67 PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo ctx, PetscInt n, PetscReal points[])
68 {
69   PetscErrorCode ierr;
70 
71   PetscFunctionBegin;
72   if (ctx->dim < 0) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
73   if (ctx->points)  SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot add points multiple times yet");
74   ctx->nInput = n;
75 
76   ierr = PetscMalloc1(n*ctx->dim, &ctx->points);CHKERRQ(ierr);
77   ierr = PetscMemcpy(ctx->points, points, n*ctx->dim * sizeof(PetscReal));CHKERRQ(ierr);
78   PetscFunctionReturn(0);
79 }
80 
81 #undef __FUNCT__
82 #define __FUNCT__ "DMInterpolationSetUp"
83 PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo ctx, DM dm, PetscBool redundantPoints)
84 {
85   MPI_Comm       comm = ctx->comm;
86   PetscScalar    *a;
87   PetscInt       p, q, i;
88   PetscMPIInt    rank, size;
89   PetscErrorCode ierr;
90   Vec            pointVec;
91   IS             cellIS;
92   PetscLayout    layout;
93   PetscReal      *globalPoints;
94   PetscScalar    *globalPointsScalar;
95   const PetscInt *ranges;
96   PetscMPIInt    *counts, *displs;
97   const PetscInt *foundCells;
98   PetscMPIInt    *foundProcs, *globalProcs;
99   PetscInt       n, N;
100 
101   PetscFunctionBegin;
102   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
103   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
104   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
105   if (ctx->dim < 0) SETERRQ(comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set");
106   /* Locate points */
107   n = ctx->nInput;
108   if (!redundantPoints) {
109     ierr = PetscLayoutCreate(comm, &layout);CHKERRQ(ierr);
110     ierr = PetscLayoutSetBlockSize(layout, 1);CHKERRQ(ierr);
111     ierr = PetscLayoutSetLocalSize(layout, n);CHKERRQ(ierr);
112     ierr = PetscLayoutSetUp(layout);CHKERRQ(ierr);
113     ierr = PetscLayoutGetSize(layout, &N);CHKERRQ(ierr);
114     /* Communicate all points to all processes */
115     ierr = PetscMalloc3(N*ctx->dim,&globalPoints,size,&counts,size,&displs);CHKERRQ(ierr);
116     ierr = PetscLayoutGetRanges(layout, &ranges);CHKERRQ(ierr);
117     for (p = 0; p < size; ++p) {
118       counts[p] = (ranges[p+1] - ranges[p])*ctx->dim;
119       displs[p] = ranges[p]*ctx->dim;
120     }
121     ierr = MPI_Allgatherv(ctx->points, n*ctx->dim, MPIU_REAL, globalPoints, counts, displs, MPIU_REAL, comm);CHKERRQ(ierr);
122   } else {
123     N = n;
124     globalPoints = ctx->points;
125     counts = displs = NULL;
126     layout = NULL;
127   }
128 #if 0
129   ierr = PetscMalloc3(N,&foundCells,N,&foundProcs,N,&globalProcs);CHKERRQ(ierr);
130   /* foundCells[p] = m->locatePoint(&globalPoints[p*ctx->dim]); */
131 #else
132 #if defined(PETSC_USE_COMPLEX)
133   ierr = PetscMalloc1(N,&globalPointsScalar);CHKERRQ(ierr);
134   for (i=0; i<N; i++) globalPointsScalar[i] = globalPoints[i];
135 #else
136   globalPointsScalar = globalPoints;
137 #endif
138   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF, ctx->dim, N*ctx->dim, globalPointsScalar, &pointVec);CHKERRQ(ierr);
139   ierr = PetscMalloc2(N,&foundProcs,N,&globalProcs);CHKERRQ(ierr);
140   ierr = DMLocatePoints(dm, pointVec, &cellIS);CHKERRQ(ierr);
141   ierr = ISGetIndices(cellIS, &foundCells);CHKERRQ(ierr);
142 #endif
143   for (p = 0; p < N; ++p) {
144     if (foundCells[p] >= 0) foundProcs[p] = rank;
145     else foundProcs[p] = size;
146   }
147   /* Let the lowest rank process own each point */
148   ierr   = MPI_Allreduce(foundProcs, globalProcs, N, MPI_INT, MPI_MIN, comm);CHKERRQ(ierr);
149   ctx->n = 0;
150   for (p = 0; p < N; ++p) {
151     if (globalProcs[p] == size) SETERRQ4(comm, PETSC_ERR_PLIB, "Point %d: %g %g %g not located in mesh", p, globalPoints[p*ctx->dim+0], ctx->dim > 1 ? globalPoints[p*ctx->dim+1] : 0.0, ctx->dim > 2 ? globalPoints[p*ctx->dim+2] : 0.0);
152     else if (globalProcs[p] == rank) ctx->n++;
153   }
154   /* Create coordinates vector and array of owned cells */
155   ierr = PetscMalloc1(ctx->n, &ctx->cells);CHKERRQ(ierr);
156   ierr = VecCreate(comm, &ctx->coords);CHKERRQ(ierr);
157   ierr = VecSetSizes(ctx->coords, ctx->n*ctx->dim, PETSC_DECIDE);CHKERRQ(ierr);
158   ierr = VecSetBlockSize(ctx->coords, ctx->dim);CHKERRQ(ierr);
159   ierr = VecSetType(ctx->coords,VECSTANDARD);CHKERRQ(ierr);
160   ierr = VecGetArray(ctx->coords, &a);CHKERRQ(ierr);
161   for (p = 0, q = 0, i = 0; p < N; ++p) {
162     if (globalProcs[p] == rank) {
163       PetscInt d;
164 
165       for (d = 0; d < ctx->dim; ++d, ++i) a[i] = globalPoints[p*ctx->dim+d];
166       ctx->cells[q++] = foundCells[p];
167     }
168   }
169   ierr = VecRestoreArray(ctx->coords, &a);CHKERRQ(ierr);
170 #if 0
171   ierr = PetscFree3(foundCells,foundProcs,globalProcs);CHKERRQ(ierr);
172 #else
173   ierr = PetscFree2(foundProcs,globalProcs);CHKERRQ(ierr);
174   ierr = ISRestoreIndices(cellIS, &foundCells);CHKERRQ(ierr);
175   ierr = ISDestroy(&cellIS);CHKERRQ(ierr);
176   ierr = VecDestroy(&pointVec);CHKERRQ(ierr);
177 #endif
178   if ((void*)globalPointsScalar != (void*)globalPoints) {ierr = PetscFree(globalPointsScalar);CHKERRQ(ierr);}
179   if (!redundantPoints) {ierr = PetscFree3(globalPoints,counts,displs);CHKERRQ(ierr);}
180   ierr = PetscLayoutDestroy(&layout);CHKERRQ(ierr);
181   PetscFunctionReturn(0);
182 }
183 
184 #undef __FUNCT__
185 #define __FUNCT__ "DMInterpolationGetCoordinates"
186 PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo ctx, Vec *coordinates)
187 {
188   PetscFunctionBegin;
189   PetscValidPointer(coordinates, 2);
190   if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
191   *coordinates = ctx->coords;
192   PetscFunctionReturn(0);
193 }
194 
195 #undef __FUNCT__
196 #define __FUNCT__ "DMInterpolationGetVector"
197 PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo ctx, Vec *v)
198 {
199   PetscErrorCode ierr;
200 
201   PetscFunctionBegin;
202   PetscValidPointer(v, 2);
203   if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
204   ierr = VecCreate(ctx->comm, v);CHKERRQ(ierr);
205   ierr = VecSetSizes(*v, ctx->n*ctx->dof, PETSC_DECIDE);CHKERRQ(ierr);
206   ierr = VecSetBlockSize(*v, ctx->dof);CHKERRQ(ierr);
207   ierr = VecSetType(*v,VECSTANDARD);CHKERRQ(ierr);
208   PetscFunctionReturn(0);
209 }
210 
211 #undef __FUNCT__
212 #define __FUNCT__ "DMInterpolationRestoreVector"
213 PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo ctx, Vec *v)
214 {
215   PetscErrorCode ierr;
216 
217   PetscFunctionBegin;
218   PetscValidPointer(v, 2);
219   if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");
220   ierr = VecDestroy(v);CHKERRQ(ierr);
221   PetscFunctionReturn(0);
222 }
223 
224 #undef __FUNCT__
225 #define __FUNCT__ "DMInterpolate_Triangle_Private"
226 PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Triangle_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
227 {
228   PetscReal      *v0, *J, *invJ, detJ;
229   PetscScalar    *a, *coords;
230   PetscInt       p;
231   PetscErrorCode ierr;
232 
233   PetscFunctionBegin;
234   ierr = PetscMalloc3(ctx->dim,&v0,ctx->dim*ctx->dim,&J,ctx->dim*ctx->dim,&invJ);CHKERRQ(ierr);
235   ierr = VecGetArray(ctx->coords, &coords);CHKERRQ(ierr);
236   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
237   for (p = 0; p < ctx->n; ++p) {
238     PetscInt     c = ctx->cells[p];
239     PetscScalar *x = NULL;
240     PetscReal    xi[4];
241     PetscInt     d, f, comp;
242 
243     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
244     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
245     ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr);
246     for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] = x[0*ctx->dof+comp];
247 
248     for (d = 0; d < ctx->dim; ++d) {
249       xi[d] = 0.0;
250       for (f = 0; f < ctx->dim; ++f) xi[d] += invJ[d*ctx->dim+f]*0.5*PetscRealPart(coords[p*ctx->dim+f] - v0[f]);
251       for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] += PetscRealPart(x[(d+1)*ctx->dof+comp] - x[0*ctx->dof+comp])*xi[d];
252     }
253     ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr);
254   }
255   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
256   ierr = VecRestoreArray(ctx->coords, &coords);CHKERRQ(ierr);
257   ierr = PetscFree3(v0, J, invJ);CHKERRQ(ierr);
258   PetscFunctionReturn(0);
259 }
260 
261 #undef __FUNCT__
262 #define __FUNCT__ "DMInterpolate_Tetrahedron_Private"
263 PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Tetrahedron_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
264 {
265   PetscReal      *v0, *J, *invJ, detJ;
266   PetscScalar    *a, *coords;
267   PetscInt       p;
268   PetscErrorCode ierr;
269 
270   PetscFunctionBegin;
271   ierr = PetscMalloc3(ctx->dim,&v0,ctx->dim*ctx->dim,&J,ctx->dim*ctx->dim,&invJ);CHKERRQ(ierr);
272   ierr = VecGetArray(ctx->coords, &coords);CHKERRQ(ierr);
273   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
274   for (p = 0; p < ctx->n; ++p) {
275     PetscInt       c = ctx->cells[p];
276     const PetscInt order[3] = {2, 1, 3};
277     PetscScalar   *x = NULL;
278     PetscReal      xi[4];
279     PetscInt       d, f, comp;
280 
281     ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr);
282     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
283     ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr);
284     for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] = x[0*ctx->dof+comp];
285 
286     for (d = 0; d < ctx->dim; ++d) {
287       xi[d] = 0.0;
288       for (f = 0; f < ctx->dim; ++f) xi[d] += invJ[d*ctx->dim+f]*0.5*PetscRealPart(coords[p*ctx->dim+f] - v0[f]);
289       for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] += PetscRealPart(x[order[d]*ctx->dof+comp] - x[0*ctx->dof+comp])*xi[d];
290     }
291     ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr);
292   }
293   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
294   ierr = VecRestoreArray(ctx->coords, &coords);CHKERRQ(ierr);
295   ierr = PetscFree3(v0, J, invJ);CHKERRQ(ierr);
296   PetscFunctionReturn(0);
297 }
298 
299 #undef __FUNCT__
300 #define __FUNCT__ "QuadMap_Private"
301 PETSC_STATIC_INLINE PetscErrorCode QuadMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx)
302 {
303   const PetscScalar *vertices = (const PetscScalar*) ctx;
304   const PetscScalar x0        = vertices[0];
305   const PetscScalar y0        = vertices[1];
306   const PetscScalar x1        = vertices[2];
307   const PetscScalar y1        = vertices[3];
308   const PetscScalar x2        = vertices[4];
309   const PetscScalar y2        = vertices[5];
310   const PetscScalar x3        = vertices[6];
311   const PetscScalar y3        = vertices[7];
312   const PetscScalar f_1       = x1 - x0;
313   const PetscScalar g_1       = y1 - y0;
314   const PetscScalar f_3       = x3 - x0;
315   const PetscScalar g_3       = y3 - y0;
316   const PetscScalar f_01      = x2 - x1 - x3 + x0;
317   const PetscScalar g_01      = y2 - y1 - y3 + y0;
318   PetscScalar       *ref, *real;
319   PetscErrorCode    ierr;
320 
321   PetscFunctionBegin;
322   ierr = VecGetArray(Xref,  &ref);CHKERRQ(ierr);
323   ierr = VecGetArray(Xreal, &real);CHKERRQ(ierr);
324   {
325     const PetscScalar p0 = ref[0];
326     const PetscScalar p1 = ref[1];
327 
328     real[0] = x0 + f_1 * p0 + f_3 * p1 + f_01 * p0 * p1;
329     real[1] = y0 + g_1 * p0 + g_3 * p1 + g_01 * p0 * p1;
330   }
331   ierr = PetscLogFlops(28);CHKERRQ(ierr);
332   ierr = VecRestoreArray(Xref,  &ref);CHKERRQ(ierr);
333   ierr = VecRestoreArray(Xreal, &real);CHKERRQ(ierr);
334   PetscFunctionReturn(0);
335 }
336 
337 #include <petsc-private/dmimpl.h>
338 #undef __FUNCT__
339 #define __FUNCT__ "QuadJacobian_Private"
340 PETSC_STATIC_INLINE PetscErrorCode QuadJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx)
341 {
342   const PetscScalar *vertices = (const PetscScalar*) ctx;
343   const PetscScalar x0        = vertices[0];
344   const PetscScalar y0        = vertices[1];
345   const PetscScalar x1        = vertices[2];
346   const PetscScalar y1        = vertices[3];
347   const PetscScalar x2        = vertices[4];
348   const PetscScalar y2        = vertices[5];
349   const PetscScalar x3        = vertices[6];
350   const PetscScalar y3        = vertices[7];
351   const PetscScalar f_01      = x2 - x1 - x3 + x0;
352   const PetscScalar g_01      = y2 - y1 - y3 + y0;
353   PetscScalar       *ref;
354   PetscErrorCode    ierr;
355 
356   PetscFunctionBegin;
357   ierr = VecGetArray(Xref,  &ref);CHKERRQ(ierr);
358   {
359     const PetscScalar x       = ref[0];
360     const PetscScalar y       = ref[1];
361     const PetscInt    rows[2] = {0, 1};
362     PetscScalar       values[4];
363 
364     values[0] = (x1 - x0 + f_01*y) * 0.5; values[1] = (x3 - x0 + f_01*x) * 0.5;
365     values[2] = (y1 - y0 + g_01*y) * 0.5; values[3] = (y3 - y0 + g_01*x) * 0.5;
366     ierr      = MatSetValues(J, 2, rows, 2, rows, values, INSERT_VALUES);CHKERRQ(ierr);
367   }
368   ierr = PetscLogFlops(30);CHKERRQ(ierr);
369   ierr = VecRestoreArray(Xref,  &ref);CHKERRQ(ierr);
370   ierr = MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
371   ierr = MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
372   PetscFunctionReturn(0);
373 }
374 
375 #undef __FUNCT__
376 #define __FUNCT__ "DMInterpolate_Quad_Private"
377 PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Quad_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
378 {
379   DM             dmCoord;
380   SNES           snes;
381   KSP            ksp;
382   PC             pc;
383   Vec            coordsLocal, r, ref, real;
384   Mat            J;
385   PetscScalar    *a, *coords;
386   PetscInt       p;
387   PetscErrorCode ierr;
388 
389   PetscFunctionBegin;
390   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
391   ierr = DMGetCoordinateDM(dm, &dmCoord);CHKERRQ(ierr);
392   ierr = SNESCreate(PETSC_COMM_SELF, &snes);CHKERRQ(ierr);
393   ierr = SNESSetOptionsPrefix(snes, "quad_interp_");CHKERRQ(ierr);
394   ierr = VecCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr);
395   ierr = VecSetSizes(r, 2, 2);CHKERRQ(ierr);
396   ierr = VecSetType(r,dm->vectype);CHKERRQ(ierr);
397   ierr = VecDuplicate(r, &ref);CHKERRQ(ierr);
398   ierr = VecDuplicate(r, &real);CHKERRQ(ierr);
399   ierr = MatCreate(PETSC_COMM_SELF, &J);CHKERRQ(ierr);
400   ierr = MatSetSizes(J, 2, 2, 2, 2);CHKERRQ(ierr);
401   ierr = MatSetType(J, MATSEQDENSE);CHKERRQ(ierr);
402   ierr = MatSetUp(J);CHKERRQ(ierr);
403   ierr = SNESSetFunction(snes, r, QuadMap_Private, NULL);CHKERRQ(ierr);
404   ierr = SNESSetJacobian(snes, J, J, QuadJacobian_Private, NULL);CHKERRQ(ierr);
405   ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr);
406   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
407   ierr = PCSetType(pc, PCLU);CHKERRQ(ierr);
408   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
409 
410   ierr = VecGetArray(ctx->coords, &coords);CHKERRQ(ierr);
411   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
412   for (p = 0; p < ctx->n; ++p) {
413     PetscScalar *x = NULL, *vertices = NULL;
414     PetscScalar *xi;
415     PetscReal    xir[2];
416     PetscInt     c = ctx->cells[p], comp, coordSize, xSize;
417 
418     /* Can make this do all points at once */
419     ierr = DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
420     if (4*2 != coordSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", coordSize, 4*2);
421     ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
422     if (4*ctx->dof != xSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", xSize, 4*ctx->dof);
423     ierr   = SNESSetFunction(snes, NULL, NULL, (void*) vertices);CHKERRQ(ierr);
424     ierr   = SNESSetJacobian(snes, NULL, NULL, NULL, (void*) vertices);CHKERRQ(ierr);
425     ierr   = VecGetArray(real, &xi);CHKERRQ(ierr);
426     xi[0]  = coords[p*ctx->dim+0];
427     xi[1]  = coords[p*ctx->dim+1];
428     ierr   = VecRestoreArray(real, &xi);CHKERRQ(ierr);
429     ierr   = SNESSolve(snes, real, ref);CHKERRQ(ierr);
430     ierr   = VecGetArray(ref, &xi);CHKERRQ(ierr);
431     xir[0] = PetscRealPart(xi[0]);
432     xir[1] = PetscRealPart(xi[1]);
433     for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] = x[0*ctx->dof+comp]*(1 - xir[0])*(1 - xir[1]) + x[1*ctx->dof+comp]*xir[0]*(1 - xir[1]) + x[2*ctx->dof+comp]*xir[0]*xir[1] + x[3*ctx->dof+comp]*(1 - xir[0])*xir[1];
434 
435     ierr = VecRestoreArray(ref, &xi);CHKERRQ(ierr);
436     ierr = DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
437     ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
438   }
439   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
440   ierr = VecRestoreArray(ctx->coords, &coords);CHKERRQ(ierr);
441 
442   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
443   ierr = VecDestroy(&r);CHKERRQ(ierr);
444   ierr = VecDestroy(&ref);CHKERRQ(ierr);
445   ierr = VecDestroy(&real);CHKERRQ(ierr);
446   ierr = MatDestroy(&J);CHKERRQ(ierr);
447   PetscFunctionReturn(0);
448 }
449 
450 #undef __FUNCT__
451 #define __FUNCT__ "HexMap_Private"
452 PETSC_STATIC_INLINE PetscErrorCode HexMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx)
453 {
454   const PetscScalar *vertices = (const PetscScalar*) ctx;
455   const PetscScalar x0        = vertices[0];
456   const PetscScalar y0        = vertices[1];
457   const PetscScalar z0        = vertices[2];
458   const PetscScalar x1        = vertices[9];
459   const PetscScalar y1        = vertices[10];
460   const PetscScalar z1        = vertices[11];
461   const PetscScalar x2        = vertices[6];
462   const PetscScalar y2        = vertices[7];
463   const PetscScalar z2        = vertices[8];
464   const PetscScalar x3        = vertices[3];
465   const PetscScalar y3        = vertices[4];
466   const PetscScalar z3        = vertices[5];
467   const PetscScalar x4        = vertices[12];
468   const PetscScalar y4        = vertices[13];
469   const PetscScalar z4        = vertices[14];
470   const PetscScalar x5        = vertices[15];
471   const PetscScalar y5        = vertices[16];
472   const PetscScalar z5        = vertices[17];
473   const PetscScalar x6        = vertices[18];
474   const PetscScalar y6        = vertices[19];
475   const PetscScalar z6        = vertices[20];
476   const PetscScalar x7        = vertices[21];
477   const PetscScalar y7        = vertices[22];
478   const PetscScalar z7        = vertices[23];
479   const PetscScalar f_1       = x1 - x0;
480   const PetscScalar g_1       = y1 - y0;
481   const PetscScalar h_1       = z1 - z0;
482   const PetscScalar f_3       = x3 - x0;
483   const PetscScalar g_3       = y3 - y0;
484   const PetscScalar h_3       = z3 - z0;
485   const PetscScalar f_4       = x4 - x0;
486   const PetscScalar g_4       = y4 - y0;
487   const PetscScalar h_4       = z4 - z0;
488   const PetscScalar f_01      = x2 - x1 - x3 + x0;
489   const PetscScalar g_01      = y2 - y1 - y3 + y0;
490   const PetscScalar h_01      = z2 - z1 - z3 + z0;
491   const PetscScalar f_12      = x7 - x3 - x4 + x0;
492   const PetscScalar g_12      = y7 - y3 - y4 + y0;
493   const PetscScalar h_12      = z7 - z3 - z4 + z0;
494   const PetscScalar f_02      = x5 - x1 - x4 + x0;
495   const PetscScalar g_02      = y5 - y1 - y4 + y0;
496   const PetscScalar h_02      = z5 - z1 - z4 + z0;
497   const PetscScalar f_012     = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
498   const PetscScalar g_012     = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
499   const PetscScalar h_012     = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
500   PetscScalar       *ref, *real;
501   PetscErrorCode    ierr;
502 
503   PetscFunctionBegin;
504   ierr = VecGetArray(Xref,  &ref);CHKERRQ(ierr);
505   ierr = VecGetArray(Xreal, &real);CHKERRQ(ierr);
506   {
507     const PetscScalar p0 = ref[0];
508     const PetscScalar p1 = ref[1];
509     const PetscScalar p2 = ref[2];
510 
511     real[0] = x0 + f_1*p0 + f_3*p1 + f_4*p2 + f_01*p0*p1 + f_12*p1*p2 + f_02*p0*p2 + f_012*p0*p1*p2;
512     real[1] = y0 + g_1*p0 + g_3*p1 + g_4*p2 + g_01*p0*p1 + g_01*p0*p1 + g_12*p1*p2 + g_02*p0*p2 + g_012*p0*p1*p2;
513     real[2] = z0 + h_1*p0 + h_3*p1 + h_4*p2 + h_01*p0*p1 + h_01*p0*p1 + h_12*p1*p2 + h_02*p0*p2 + h_012*p0*p1*p2;
514   }
515   ierr = PetscLogFlops(114);CHKERRQ(ierr);
516   ierr = VecRestoreArray(Xref,  &ref);CHKERRQ(ierr);
517   ierr = VecRestoreArray(Xreal, &real);CHKERRQ(ierr);
518   PetscFunctionReturn(0);
519 }
520 
521 #undef __FUNCT__
522 #define __FUNCT__ "HexJacobian_Private"
523 PETSC_STATIC_INLINE PetscErrorCode HexJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx)
524 {
525   const PetscScalar *vertices = (const PetscScalar*) ctx;
526   const PetscScalar x0        = vertices[0];
527   const PetscScalar y0        = vertices[1];
528   const PetscScalar z0        = vertices[2];
529   const PetscScalar x1        = vertices[9];
530   const PetscScalar y1        = vertices[10];
531   const PetscScalar z1        = vertices[11];
532   const PetscScalar x2        = vertices[6];
533   const PetscScalar y2        = vertices[7];
534   const PetscScalar z2        = vertices[8];
535   const PetscScalar x3        = vertices[3];
536   const PetscScalar y3        = vertices[4];
537   const PetscScalar z3        = vertices[5];
538   const PetscScalar x4        = vertices[12];
539   const PetscScalar y4        = vertices[13];
540   const PetscScalar z4        = vertices[14];
541   const PetscScalar x5        = vertices[15];
542   const PetscScalar y5        = vertices[16];
543   const PetscScalar z5        = vertices[17];
544   const PetscScalar x6        = vertices[18];
545   const PetscScalar y6        = vertices[19];
546   const PetscScalar z6        = vertices[20];
547   const PetscScalar x7        = vertices[21];
548   const PetscScalar y7        = vertices[22];
549   const PetscScalar z7        = vertices[23];
550   const PetscScalar f_xy      = x2 - x1 - x3 + x0;
551   const PetscScalar g_xy      = y2 - y1 - y3 + y0;
552   const PetscScalar h_xy      = z2 - z1 - z3 + z0;
553   const PetscScalar f_yz      = x7 - x3 - x4 + x0;
554   const PetscScalar g_yz      = y7 - y3 - y4 + y0;
555   const PetscScalar h_yz      = z7 - z3 - z4 + z0;
556   const PetscScalar f_xz      = x5 - x1 - x4 + x0;
557   const PetscScalar g_xz      = y5 - y1 - y4 + y0;
558   const PetscScalar h_xz      = z5 - z1 - z4 + z0;
559   const PetscScalar f_xyz     = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7;
560   const PetscScalar g_xyz     = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7;
561   const PetscScalar h_xyz     = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7;
562   PetscScalar       *ref;
563   PetscErrorCode    ierr;
564 
565   PetscFunctionBegin;
566   ierr = VecGetArray(Xref,  &ref);CHKERRQ(ierr);
567   {
568     const PetscScalar x       = ref[0];
569     const PetscScalar y       = ref[1];
570     const PetscScalar z       = ref[2];
571     const PetscInt    rows[3] = {0, 1, 2};
572     PetscScalar       values[9];
573 
574     values[0] = (x1 - x0 + f_xy*y + f_xz*z + f_xyz*y*z) / 2.0;
575     values[1] = (x3 - x0 + f_xy*x + f_yz*z + f_xyz*x*z) / 2.0;
576     values[2] = (x4 - x0 + f_yz*y + f_xz*x + f_xyz*x*y) / 2.0;
577     values[3] = (y1 - y0 + g_xy*y + g_xz*z + g_xyz*y*z) / 2.0;
578     values[4] = (y3 - y0 + g_xy*x + g_yz*z + g_xyz*x*z) / 2.0;
579     values[5] = (y4 - y0 + g_yz*y + g_xz*x + g_xyz*x*y) / 2.0;
580     values[6] = (z1 - z0 + h_xy*y + h_xz*z + h_xyz*y*z) / 2.0;
581     values[7] = (z3 - z0 + h_xy*x + h_yz*z + h_xyz*x*z) / 2.0;
582     values[8] = (z4 - z0 + h_yz*y + h_xz*x + h_xyz*x*y) / 2.0;
583 
584     ierr = MatSetValues(J, 3, rows, 3, rows, values, INSERT_VALUES);CHKERRQ(ierr);
585   }
586   ierr = PetscLogFlops(152);CHKERRQ(ierr);
587   ierr = VecRestoreArray(Xref,  &ref);CHKERRQ(ierr);
588   ierr = MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
589   ierr = MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
590   PetscFunctionReturn(0);
591 }
592 
593 #undef __FUNCT__
594 #define __FUNCT__ "DMInterpolate_Hex_Private"
595 PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Hex_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v)
596 {
597   DM             dmCoord;
598   SNES           snes;
599   KSP            ksp;
600   PC             pc;
601   Vec            coordsLocal, r, ref, real;
602   Mat            J;
603   PetscScalar    *a, *coords;
604   PetscInt       p;
605   PetscErrorCode ierr;
606 
607   PetscFunctionBegin;
608   ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr);
609   ierr = DMGetCoordinateDM(dm, &dmCoord);CHKERRQ(ierr);
610   ierr = SNESCreate(PETSC_COMM_SELF, &snes);CHKERRQ(ierr);
611   ierr = SNESSetOptionsPrefix(snes, "hex_interp_");CHKERRQ(ierr);
612   ierr = VecCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr);
613   ierr = VecSetSizes(r, 3, 3);CHKERRQ(ierr);
614   ierr = VecSetType(r,dm->vectype);CHKERRQ(ierr);
615   ierr = VecDuplicate(r, &ref);CHKERRQ(ierr);
616   ierr = VecDuplicate(r, &real);CHKERRQ(ierr);
617   ierr = MatCreate(PETSC_COMM_SELF, &J);CHKERRQ(ierr);
618   ierr = MatSetSizes(J, 3, 3, 3, 3);CHKERRQ(ierr);
619   ierr = MatSetType(J, MATSEQDENSE);CHKERRQ(ierr);
620   ierr = MatSetUp(J);CHKERRQ(ierr);
621   ierr = SNESSetFunction(snes, r, HexMap_Private, NULL);CHKERRQ(ierr);
622   ierr = SNESSetJacobian(snes, J, J, HexJacobian_Private, NULL);CHKERRQ(ierr);
623   ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr);
624   ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr);
625   ierr = PCSetType(pc, PCLU);CHKERRQ(ierr);
626   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
627 
628   ierr = VecGetArray(ctx->coords, &coords);CHKERRQ(ierr);
629   ierr = VecGetArray(v, &a);CHKERRQ(ierr);
630   for (p = 0; p < ctx->n; ++p) {
631     PetscScalar *x = NULL, *vertices = NULL;
632     PetscScalar *xi;
633     PetscReal    xir[3];
634     PetscInt     c = ctx->cells[p], comp, coordSize, xSize;
635 
636     /* Can make this do all points at once */
637     ierr = DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
638     if (8*3 != coordSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", coordSize, 8*3);
639     ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
640     if (8*ctx->dof != xSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", xSize, 8*ctx->dof);
641     ierr   = SNESSetFunction(snes, NULL, NULL, (void*) vertices);CHKERRQ(ierr);
642     ierr   = SNESSetJacobian(snes, NULL, NULL, NULL, (void*) vertices);CHKERRQ(ierr);
643     ierr   = VecGetArray(real, &xi);CHKERRQ(ierr);
644     xi[0]  = coords[p*ctx->dim+0];
645     xi[1]  = coords[p*ctx->dim+1];
646     xi[2]  = coords[p*ctx->dim+2];
647     ierr   = VecRestoreArray(real, &xi);CHKERRQ(ierr);
648     ierr   = SNESSolve(snes, real, ref);CHKERRQ(ierr);
649     ierr   = VecGetArray(ref, &xi);CHKERRQ(ierr);
650     xir[0] = PetscRealPart(xi[0]);
651     xir[1] = PetscRealPart(xi[1]);
652     xir[2] = PetscRealPart(xi[2]);
653     for (comp = 0; comp < ctx->dof; ++comp) {
654       a[p*ctx->dof+comp] =
655         x[0*ctx->dof+comp]*(1-xir[0])*(1-xir[1])*(1-xir[2]) +
656         x[3*ctx->dof+comp]*    xir[0]*(1-xir[1])*(1-xir[2]) +
657         x[2*ctx->dof+comp]*    xir[0]*    xir[1]*(1-xir[2]) +
658         x[1*ctx->dof+comp]*(1-xir[0])*    xir[1]*(1-xir[2]) +
659         x[4*ctx->dof+comp]*(1-xir[0])*(1-xir[1])*   xir[2] +
660         x[5*ctx->dof+comp]*    xir[0]*(1-xir[1])*   xir[2] +
661         x[6*ctx->dof+comp]*    xir[0]*    xir[1]*   xir[2] +
662         x[7*ctx->dof+comp]*(1-xir[0])*    xir[1]*   xir[2];
663     }
664     ierr = VecRestoreArray(ref, &xi);CHKERRQ(ierr);
665     ierr = DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr);
666     ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr);
667   }
668   ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
669   ierr = VecRestoreArray(ctx->coords, &coords);CHKERRQ(ierr);
670 
671   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
672   ierr = VecDestroy(&r);CHKERRQ(ierr);
673   ierr = VecDestroy(&ref);CHKERRQ(ierr);
674   ierr = VecDestroy(&real);CHKERRQ(ierr);
675   ierr = MatDestroy(&J);CHKERRQ(ierr);
676   PetscFunctionReturn(0);
677 }
678 
679 #undef __FUNCT__
680 #define __FUNCT__ "DMInterpolationEvaluate"
681 /*
682   Input Parameters:
683 + ctx - The DMInterpolationInfo context
684 . dm  - The DM
685 - x   - The local vector containing the field to be interpolated
686 
687   Output Parameters:
688 . v   - The vector containing the interpolated values
689 */
690 PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo ctx, DM dm, Vec x, Vec v)
691 {
692   PetscInt       dim, coneSize, n;
693   PetscErrorCode ierr;
694 
695   PetscFunctionBegin;
696   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
697   PetscValidHeaderSpecific(x, VEC_CLASSID, 3);
698   PetscValidHeaderSpecific(v, VEC_CLASSID, 4);
699   ierr = VecGetLocalSize(v, &n);CHKERRQ(ierr);
700   if (n != ctx->n*ctx->dof) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid input vector size %d should be %d", n, ctx->n*ctx->dof);
701   if (n) {
702     ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
703     ierr = DMPlexGetConeSize(dm, ctx->cells[0], &coneSize);CHKERRQ(ierr);
704     if (dim == 2) {
705       if (coneSize == 3) {
706         ierr = DMInterpolate_Triangle_Private(ctx, dm, x, v);CHKERRQ(ierr);
707       } else if (coneSize == 4) {
708         ierr = DMInterpolate_Quad_Private(ctx, dm, x, v);CHKERRQ(ierr);
709       } else SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %d for point interpolation", dim);
710     } else if (dim == 3) {
711       if (coneSize == 4) {
712         ierr = DMInterpolate_Tetrahedron_Private(ctx, dm, x, v);CHKERRQ(ierr);
713       } else {
714         ierr = DMInterpolate_Hex_Private(ctx, dm, x, v);CHKERRQ(ierr);
715       }
716     } else SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %d for point interpolation", dim);
717   }
718   PetscFunctionReturn(0);
719 }
720 
721 #undef __FUNCT__
722 #define __FUNCT__ "DMInterpolationDestroy"
723 PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *ctx)
724 {
725   PetscErrorCode ierr;
726 
727   PetscFunctionBegin;
728   PetscValidPointer(ctx, 2);
729   ierr = VecDestroy(&(*ctx)->coords);CHKERRQ(ierr);
730   ierr = PetscFree((*ctx)->points);CHKERRQ(ierr);
731   ierr = PetscFree((*ctx)->cells);CHKERRQ(ierr);
732   ierr = PetscFree(*ctx);CHKERRQ(ierr);
733   *ctx = NULL;
734   PetscFunctionReturn(0);
735 }
736 
737 #undef __FUNCT__
738 #define __FUNCT__ "SNESMonitorFields"
739 /*@C
740   SNESMonitorFields - Monitors the residual for each field separately
741 
742   Collective on SNES
743 
744   Input Parameters:
745 + snes   - the SNES context
746 . its    - iteration number
747 . fgnorm - 2-norm of residual
748 - dummy  - unused context
749 
750   Notes:
751   This routine prints the residual norm at each iteration.
752 
753   Level: intermediate
754 
755 .keywords: SNES, nonlinear, default, monitor, norm
756 .seealso: SNESMonitorSet(), SNESMonitorDefault()
757 @*/
758 PetscErrorCode SNESMonitorFields(SNES snes, PetscInt its, PetscReal fgnorm, void *dummy)
759 {
760   PetscViewer        viewer = dummy ? (PetscViewer) dummy : PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject) snes));
761   Vec                res;
762   DM                 dm;
763   PetscSection       s;
764   const PetscScalar *r;
765   PetscReal         *lnorms, *norms;
766   PetscInt           numFields, f, pStart, pEnd, p;
767   PetscErrorCode     ierr;
768 
769   PetscFunctionBegin;
770   ierr = SNESGetFunction(snes, &res, 0, 0);CHKERRQ(ierr);
771   ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr);
772   ierr = DMGetDefaultSection(dm, &s);CHKERRQ(ierr);
773   ierr = PetscSectionGetNumFields(s, &numFields);CHKERRQ(ierr);
774   ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr);
775   ierr = PetscCalloc2(numFields, &lnorms, numFields, &norms);CHKERRQ(ierr);
776   ierr = VecGetArrayRead(res, &r);CHKERRQ(ierr);
777   for (p = pStart; p < pEnd; ++p) {
778     for (f = 0; f < numFields; ++f) {
779       PetscInt fdof, foff, d;
780 
781       ierr = PetscSectionGetFieldDof(s, p, f, &fdof);CHKERRQ(ierr);
782       ierr = PetscSectionGetFieldOffset(s, p, f, &foff);CHKERRQ(ierr);
783       for (d = 0; d < fdof; ++d) lnorms[f] += PetscRealPart(PetscSqr(r[foff+d]));
784     }
785   }
786   ierr = VecRestoreArrayRead(res, &r);CHKERRQ(ierr);
787   ierr = MPI_Allreduce(lnorms, norms, numFields, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject) dm));CHKERRQ(ierr);
788   ierr = PetscViewerASCIIAddTab(viewer, ((PetscObject) snes)->tablevel);CHKERRQ(ierr);
789   ierr = PetscViewerASCIIPrintf(viewer, "%3D SNES Function norm %14.12e [", its, (double) fgnorm);CHKERRQ(ierr);
790   for (f = 0; f < numFields; ++f) {
791     if (f > 0) {ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr);}
792     ierr = PetscViewerASCIIPrintf(viewer, "%14.12e", (double) PetscSqrtReal(norms[f]));CHKERRQ(ierr);
793   }
794   ierr = PetscViewerASCIIPrintf(viewer, "]\n");CHKERRQ(ierr);
795   ierr = PetscViewerASCIISubtractTab(viewer, ((PetscObject) snes)->tablevel);CHKERRQ(ierr);
796   ierr = PetscFree2(lnorms, norms);CHKERRQ(ierr);
797   PetscFunctionReturn(0);
798 }
799