Actual source code: grid.c
1: #ifdef PETSC_RCS_HEADER
2: static char vcid[] = "$Id: grid.c,v 1.24 2000/07/16 05:40:26 knepley Exp $";
3: #endif
5: /*
6: This file provides routines for grid vectors (vectors that are associated
7: with grids, possibly with multiple degrees of freedom per node).
9: This component is new; changes in the user interface will be occuring in the
10: near future!
11: */
13: #include "petscsles.h" /* For ALE Mesh Support *//*I "petscsles.h" I*/
14: #include "src/grid/gridimpl.h" /*I "grid.h" I*//*I "gvec.h" I*/
16: /* Logging support */
17: int GRID_COOKIE;
18: int GRID_Reform, GRID_SetUp;
20: /*---------------------------------------------- Standard Grid Functions --------------------------------------------*/
21: /*@
22: GridDestroy - Destroys a grid.
24: Collective on Grid
26: Input Parameter:
27: . grid - The grid
29: Level: beginner
31: .keywrods: grid, destroy
32: .seealso: GridView(), GridCreateGVec(), GridSetUp(), GridRefine()
33: @*/
34: int GridDestroy(Grid grid)
35: {
40: if (--grid->refct > 0) return(0);
41: (*grid->ops->destroy)(grid);
42: PetscLogObjectDestroy(grid);
43: PetscHeaderDestroy(grid);
44: return(0);
45: }
47: /*@
48: GridView - Views a grid.
50: Collective on Grid
52: Input Parameters:
53: + grid - The grid
54: - viewer - The viewer, PETSC_VIEWER_STDOUT_WORLD by default
56: Level: beginner
58: .keywords: grid, view
59: .seealso: GVecView()
60: @*/
61: int GridView(Grid grid, PetscViewer viewer)
62: {
67: /* Do not check return from setup, since we may not have a problem yet */
68: if (grid->setupcalled == PETSC_FALSE) GridSetUp(grid);
69: if (!viewer) {
70: viewer = PETSC_VIEWER_STDOUT_SELF;
71: } else {
73: }
74: (*grid->ops->view)(grid, viewer);
75: return(0);
76: }
78: EXTERN_C_BEGIN
79: int GridSerialize_Generic(MPI_Comm comm, Grid *grid, PetscViewer viewer, PetscTruth store) {
80: Grid g;
81: Mesh mesh;
82: int fd, hasParent;
83: int one = 1;
84: int zero = 0;
85: int field, len;
86: char *type_name = PETSC_NULL;
87: PetscTruth match;
88: int ierr;
94: PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_BINARY, &match);
95: if (match == PETSC_FALSE) SETERRQ(PETSC_ERR_ARG_WRONG, "Must be binary viewer");
96: PetscViewerBinaryGetDescriptor(viewer, &fd);
98: if (store) {
99: g = *grid;
101: MeshSerialize(comm, &g->mesh, viewer, store);
102: PetscStrlen(g->type_name, &len);
103: PetscBinaryWrite(fd, &len, 1, PETSC_INT, 0);
104: PetscBinaryWrite(fd, g->type_name, len, PETSC_CHAR, 0);
106: PetscBinaryWrite(fd, &g->numFields, 1, PETSC_INT, 0);
107: for(field = 0; field < g->numFields; field++) {
108: if (g->fields[field].name != PETSC_NULL) {
109: PetscStrlen(g->fields[field].name, &len);
110: } else {
111: len = 0;
112: }
113: PetscBinaryWrite(fd, &len, 1, PETSC_INT, 0);
114: PetscBinaryWrite(fd, g->fields[field].name, len, PETSC_CHAR, 0);
115: PetscBinaryWrite(fd, &g->fields[field].numComp, 1, PETSC_INT, 0);
116: PetscStrlen(g->fields[field].discType, &len);
117: PetscBinaryWrite(fd, &len, 1, PETSC_INT, 0);
118: PetscBinaryWrite(fd, g->fields[field].discType, len, PETSC_CHAR, 0);
119: DiscretizationSerialize(g->comm, &g->fields[field].disc, viewer, store);
120: PetscBinaryWrite(fd, &g->fields[field].isActive, 1, PETSC_INT, 0);
121: PetscBinaryWrite(fd, &g->fields[field].isConstrained, 1, PETSC_INT, 0);
122: PetscBinaryWrite(fd, &g->fields[field].constraintCompDiff, 1, PETSC_INT, 0);
123: }
125: if (g->gridparent != PETSC_NULL) {
126: PetscBinaryWrite(fd, &one, 1, PETSC_INT, 0);
127: GridSerialize(comm, &g->gridparent, viewer, store);
128: } else {
129: PetscBinaryWrite(fd, &zero, 1, PETSC_INT, 0);
130: }
131: PetscBinaryWrite(fd, &g->isConstrained, 1, PETSC_INT, 0);
132: #ifdef NEW_GRID_SERIALIZE
133: /* This stuff is all made in GridSetUp */
134: FieldClassSerialize(g->comm, &g->cm, viewer, store);
135: VarOrderingSerialize(g->cm, &g->order, viewer, store);
136: FieldClassSerialize(g->comm, &g->constraintCM, viewer, store);
137: VarOrderingSerialize(g->constraintCM, &g->constraintOrder, viewer, store);
138: ISSerialize(g->comm, &g->constraintOrdering, viewer, store);
139: MatSerialize(g->comm, &g->constraintMatrix, viewer, store);
140: #endif
142: PetscBinaryWrite(fd, &g->reduceAlpha, 1, PETSC_SCALAR, 0);
143: PetscBinaryWrite(fd, &g->interpolationType, 1, PETSC_INT, 0);
144: PetscBinaryWrite(fd, &g->viewField, 1, PETSC_INT, 0);
145: PetscBinaryWrite(fd, &g->viewComp, 1, PETSC_INT, 0);
146: } else {
147: MeshSerialize(comm, &mesh, viewer, store);
148: GridCreate(mesh, &g);
149: MeshDestroy(mesh);
150: PetscBinaryRead(fd, &len, 1, PETSC_INT);
151: if (len) {
152: PetscMalloc((len+1) * sizeof(char), &type_name);
153: type_name[len] = 0;
154: }
155: PetscBinaryRead(fd, type_name, len, PETSC_CHAR);
156: PetscBinaryRead(fd, &g->numFields, 1, PETSC_INT);
157: g->maxFields = g->numFields;
158: PetscFree(g->fields);
159: PetscMalloc(g->numFields * sizeof(Field), &g->fields);
160: for(field = 0; field < g->numFields; field++) {
161: PetscBinaryRead(fd, &len, 1, PETSC_INT);
162: if (len) {
163: PetscMalloc((len+1) * sizeof(char), &g->fields[field].name);
164: g->fields[field].name[len] = 0;
165: }
166: PetscBinaryRead(fd, g->fields[field].name, len, PETSC_CHAR);
167: PetscBinaryRead(fd, &g->fields[field].numComp, 1, PETSC_INT);
168: PetscBinaryRead(fd, &len, 1, PETSC_INT);
169: PetscMalloc((len+1) * sizeof(DiscretizationType), &g->fields[field].discType);
170: g->fields[field].discType[len] = 0;
171: PetscBinaryRead(fd, g->fields[field].discType, len, PETSC_CHAR);
172: DiscretizationSerialize(g->comm, &g->fields[field].disc, viewer, store);
173: PetscTypeCompare((PetscObject) g->fields[field].disc, g->fields[field].discType, &match);
174: if (match == PETSC_FALSE) SETERRQ(PETSC_ERR_ARG_CORRUPT, "Incorrect serialization of Discretization");
175: PetscBinaryRead(fd, &g->fields[field].isActive, 1, PETSC_INT);
176: PetscBinaryRead(fd, &g->fields[field].isConstrained, 1, PETSC_INT);
177: PetscBinaryRead(fd, &g->fields[field].constraintCompDiff, 1, PETSC_INT);
178: }
179: GridSetType(g, type_name);
180: if (len) {
181: PetscFree(type_name);
182: }
184: PetscBinaryRead(fd, &hasParent, 1, PETSC_INT);
185: if (hasParent) {
186: GridSerialize(comm, &g->gridparent, viewer, store);
187: }
188: PetscBinaryRead(fd, &g->isConstrained, 1, PETSC_INT);
189: #ifdef NEW_GRID_SERIALIZE
190: /* This stuff is all made in GridSetUp */
191: FieldClassSerialize(g->comm, &g->cm, viewer, store);
192: VarOrderingSerialize(g->cm, &g->order, viewer, store);
193: FieldClassSerialize(g->comm, &g->constraintCM, viewer, store);
194: VarOrderingSerialize(g->constraintCM, &g->constraintOrder, viewer, store);
195: ISSerialize(g->comm, &g->constraintOrdering, viewer, store);
196: MatSerialize(g->comm, &g->constraintMatrix, viewer, store);
197: #endif
199: PetscBinaryRead(fd, &g->reduceAlpha, 1, PETSC_SCALAR);
200: PetscBinaryRead(fd, &g->interpolationType, 1, PETSC_INT);
201: PetscBinaryRead(fd, &g->viewField, 1, PETSC_INT);
202: PetscBinaryRead(fd, &g->viewComp, 1, PETSC_INT);
204: *grid = g;
205: }
207: return(0);
208: }
209: EXTERN_C_END
211: /*@
212: GridSetUp - Set up any required internal data structures for the grid.
214: Collective on Grid
216: Input Parameter:
217: . grid - The grid
219: Level: intermediate
221: .keywords: grid, setup
222: .seealso: GridDestroy(), GridCreateGVec()
223: @*/
224: int GridSetUp(Grid grid)
225: {
226: LocalVarOrdering reduceLocOrder;
227: LocalVarOrdering locOrder;
228: int sField, tField;
229: int bd, bc, op;
230: int ierr;
234: if (grid->setupcalled == PETSC_TRUE) return(0);
235: grid->setupcalled = PETSC_TRUE;
237: PetscLogEventBegin(GRID_SetUp, grid, 0, 0, 0);
238: if (grid->numFields <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "At least one field must be defined on the Grid");
239: /* Setup boundary sizes */
240: for(bd = 0; bd < grid->numBd; bd++) {
241: if (grid->bdSize[bd] != PETSC_NULL) {
242: PetscFree(grid->bdSize[bd]);
243: }
244: PetscMalloc(grid->numFields * sizeof(int), &grid->bdSize[bd]);
245: PetscLogObjectMemory(grid, grid->numFields * sizeof(int));
246: }
247: /* Setup type-specific stuff */
248: if (grid->ops->gridsetup) {
249: (*grid->ops->gridsetup)(grid);
250: if (ierr != 0) PetscFunctionReturn(ierr);
251: }
252: /* Setup constraint structures */
253: (*grid->ops->gridsetupconstraints)(grid, grid->constraintCtx);
254: /* Setup ghost variable structures */
255: GridSetupGhostScatter(grid);
256: /* Setup boundary conditions reductions */
257: if (grid->reduceSystem == PETSC_TRUE) {
258: /* Create storage for reduction of Rhs */
259: if (grid->bdReduceVec != PETSC_NULL) {
260: GVecDestroy(grid->bdReduceVec);
261: }
262: GVecCreateRectangularGhost(grid, grid->reduceOrder, &grid->bdReduceVec);
263: grid->bdReduceVecCur = grid->bdReduceVec;
265: /* Hopefully, the user specified the right context by now */
266: GridCalcBCValues(grid, PETSC_FALSE, grid->reduceContext);
268: if (grid->reduceElement == PETSC_FALSE) {
269: /* Create storage for explicit reduction of Rhs */
270: GVecCreate(grid, &grid->reduceVec);
271: GMatCreateRectangular(grid, grid->reduceOrder, grid->order, &grid->bdReduceMat);
272: MatSetOption(grid->bdReduceMat, MAT_NEW_NONZERO_ALLOCATION_ERR);
274: /* Initialize storage */
275: MatZeroEntries(grid->bdReduceMat);
277: /* Evaluate the elimination block */
278: for(bc = 0; bc < grid->numBC; bc++) {
279: for(op = 0; op < grid->numMatOps; op++) {
280: sField = grid->matOps[op].field;
281: tField = grid->fields[sField].disc->operators[grid->matOps[op].op]->test->field;
282: if (sField == grid->bc[bc].field) {
283: LocalVarOrderingCreate(grid, 1, &tField, &locOrder);
284: LocalVarOrderingCreate(grid, 1, &sField, &reduceLocOrder);
285: GMatEvaluateOperatorGalerkin(grid->bdReduceMat, PETSC_NULL, 1, &sField, &tField, grid->matOps[op].op,
286: grid->matOps[op].alpha, MAT_FLUSH_ASSEMBLY, grid->reduceContext);
287:
288: LocalVarOrderingDestroy(locOrder);
289: LocalVarOrderingDestroy(reduceLocOrder);
290: }
291: #ifdef PETSC_USE_BOPT_g
292: #endif
293: }
294: }
295: for(bc = 0; bc < grid->numPointBC; bc++) {
296: for(op = 0; op < grid->numMatOps; op++) {
297: sField = grid->matOps[op].field;
298: tField = grid->fields[sField].disc->operators[grid->matOps[op].op]->test->field;
299: if (sField == grid->pointBC[bc].field) {
300: LocalVarOrderingCreate(grid, 1, &tField, &locOrder);
301: LocalVarOrderingCreate(grid, 1, &sField, &reduceLocOrder);
302: GMatEvaluateOperatorGalerkin(grid->bdReduceMat, PETSC_NULL, 1, &sField, &tField, grid->matOps[op].op,
303: grid->matOps[op].alpha, MAT_FLUSH_ASSEMBLY, grid->reduceContext);
304:
305: LocalVarOrderingDestroy(locOrder);
306: LocalVarOrderingDestroy(reduceLocOrder);
307: }
308: #ifdef PETSC_USE_BOPT_g
309: #endif
310: }
311: }
312: MatAssemblyBegin(grid->bdReduceMat, MAT_FINAL_ASSEMBLY);
313: MatAssemblyEnd(grid->bdReduceMat, MAT_FINAL_ASSEMBLY);
314: }
315: }
316: PetscLogEventEnd(GRID_SetUp, grid, 0, 0, 0);
317: return(0);
318: }
320: /*@
321: GridSetupGhostScatter - Set up the scatter from a global vector to local values, including ghosts.
323: Collective on Grid
325: Input Parameter:
326: . grid - The grid
328: Level: advanced
330: .keywords: grid, ghost, scatter
331: .seealso: GridDestroy(), GridCreateGVec()
332: @*/
333: int GridSetupGhostScatter(Grid grid)
334: {
335: PetscTruth isConstrained, explicitConstraints;
336: int ierr;
340: GridSetUp(grid);
341: /* Destroy old structures */
342: if (grid->ghostVec != PETSC_NULL) {
343: VecDestroy(grid->ghostVec);
344: }
345: if (grid->ghostScatter != PETSC_NULL) {
346: VecScatterDestroy(grid->ghostScatter);
347: }
348: /* Setup scatter from global vector to local ghosted vector */
349: GridIsConstrained(grid, &isConstrained);
350: GridGetExplicitConstraints(grid, &explicitConstraints);
351: if ((isConstrained == PETSC_TRUE) && (explicitConstraints == PETSC_TRUE)) {
352: (*grid->ops->gridsetupghostscatter)(grid, grid->constraintOrder, &grid->ghostVec, &grid->ghostScatter);
353:
354: } else {
355: (*grid->ops->gridsetupghostscatter)(grid, grid->order, &grid->ghostVec, &grid->ghostScatter);
356: }
357: return(0);
358: }
360: /*@
361: GridSetupBoundary - Set up any required internal data structures for the boundary.
363: Collective on Grid
365: Input Parameter:
366: . grid - The grid
368: Level: advanced
370: .keywords: grid, setup, boundary
371: .seealso: GridDestroy(), GridCreateGVec()
372: @*/
373: int GridSetupBoundary(Grid grid)
374: {
379: if (grid->bdSetupCalled == PETSC_TRUE) return(0);
380: GridSetUp(grid);
381: if (grid->ops->gridsetupboundary) {
382: (*grid->ops->gridsetupboundary)(grid);
383: }
384: return(0);
385: }
387: /*@
388: GridGetMesh - Returns the context for the underlying mesh.
390: Not collective
392: Input Parameter:
393: . grid - The initial grid
395: Output Parameter:
396: . mesh - The mesh
398: Level: intermediate
400: .keywords: grid, mesh
401: .seealso: GridDestroy()
402: @*/
403: int GridGetMesh(Grid grid, Mesh *mesh)
404: {
408: *mesh = grid->mesh;
409: return(0);
410: }
412: /*@
413: GridGetClassMap - This function returns the default field class map for the grid.
415: Not collective
417: Input Parameter:
418: . grid - The grid
420: Output Parameter:
421: . order - The default field class map
423: Level: intermediate
425: .keywords: Grid, get, class map
426: .seealso: GridGetMesh(), GridGetOrder(), FieldClassMapCreate()
427: @*/
428: int GridGetClassMap(Grid grid, FieldClassMap *map)
429: {
433: if (grid->isConstrained == PETSC_TRUE) {
434: *map = grid->constraintCM;
435: } else {
436: *map = grid->cm;
437: }
438: return(0);
439: }
441: /*@
442: GridGetOrder - This function returns the default ordering for the grid.
444: Not collective
446: Input Parameter:
447: . grid - The grid
449: Output Parameter:
450: . order - The default variable ordering
452: Level: intermediate
454: .keywords: Grid, get, order
455: .seealso: GridGetLocalOrder(), VarOrderingCreate()
456: @*/
457: int GridGetOrder(Grid grid, VarOrdering *order)
458: {
462: if (grid->isConstrained == PETSC_TRUE) {
463: *order = grid->constraintOrder;
464: } else {
465: *order = grid->order;
466: }
467: return(0);
468: }
470: /*@
471: GridGetLocalOrder - This function returns the local ordering for the grid.
473: Not collective
475: Input Parameter:
476: . grid - The grid
478: Output Parameter:
479: . order - The local variable ordering
481: Level: intermediate
483: .keywords: Grid, get, local, order
484: .seealso: GridGetOrder(), LocalVarOrderingCreate()
485: @*/
486: int GridGetLocalOrder(Grid grid, LocalVarOrdering *order)
487: {
491: *order = grid->locOrder;
492: return(0);
493: }
495: /*@
496: GridGetDiscretization - This function returns the discretization for the given field.
498: Not collective
500: Input Parameters:
501: + grid - The grid
502: - field - The field
504: Output Parameter:
505: . disc - The Discretization
507: Level: intermediate
509: .keywords: Grid, get, Discretization, field
510: .seealso: GridGetOrder(), DiscretizationCreate()
511: @*/
512: int GridGetDiscretization(Grid grid, int field, Discretization *disc) {
516: GridValidField(grid, field);
517: *disc = grid->fields[field].disc;
518: return(0);
519: }
521: /*
522: GridSetTypeFromOptions - Sets the type of grid from user options.
524: Collective on Grid
526: Input Parameter:
527: . grid - The grid
529: Level: intermediate
531: .keywords: Grid, set, options, database, type
532: .seealso: GridSetFromOptions(), GridSetType()
533: */
534: static int GridSetTypeFromOptions(Grid grid)
535: {
536: PetscTruth opt;
537: char *defaultType;
538: char typeName[256];
539: int dim = grid->dim;
540: int ierr;
543: if (grid->type_name != PETSC_NULL) {
544: defaultType = grid->type_name;
545: } else {
546: switch(dim) {
547: case 1:
548: defaultType = GRID_TRIANGULAR_1D;
549: break;
550: case 2:
551: defaultType = GRID_TRIANGULAR_2D;
552: break;
553: default:
554: SETERRQ1(PETSC_ERR_SUP, "Grid dimension %d is not supported", dim);
555: }
556: }
558: if (GridRegisterAllCalled == PETSC_FALSE) {
559: GridRegisterAll(PETSC_NULL);
560: }
561: PetscOptionsList("-grid_type", "Grid method"," GridSetType", GridList, defaultType, typeName, 256, &opt);
562: if (opt == PETSC_TRUE) {
563: GridSetType(grid, typeName);
564: } else {
565: GridSetType(grid, defaultType);
566: }
567: return(0);
568: }
570: /*
571: GridSetSerializeTypeFromOptions - Sets the type of grid serialization from user options.
573: Collective on Grid
575: Input Parameter:
576: . grid - The grid
578: Level: intermediate
580: .keywords: Grid, set, options, database, type
581: .seealso: GridSetFromOptions(), GridSetSerializeType()
582: */
583: static int GridSetSerializeTypeFromOptions(Grid grid)
584: {
585: PetscTruth opt;
586: char *defaultType;
587: char typeName[256];
588: int ierr;
591: if (grid->serialize_name != PETSC_NULL) {
592: defaultType = grid->serialize_name;
593: } else {
594: defaultType = GRID_SER_GENERIC;
595: }
597: if (GridSerializeRegisterAllCalled == PETSC_FALSE) {
598: GridSerializeRegisterAll(PETSC_NULL);
599: }
600: PetscOptionsList("-grid_serialize_type", "Grid serialization method"," GridSetSerializeType", GridSerializeList,
601: defaultType, typeName, 256, &opt);
602: if (opt == PETSC_TRUE) {
603: GridSetSerializeType(grid, typeName);
604: } else {
605: GridSetSerializeType(grid, defaultType);
606: }
607: return(0);
608: }
610: /*@
611: GridSetFromOptions - Sets various Grid parameters from user options.
613: Collective on Grid
615: Input Parameter:
616: . grid - The grid
618: Notes: To see all options, run your program with the -help option, or consult the users manual.
619: Must be called after GridCreate() but before the Grid is used.
621: Level: intermediate
623: .keywords: Grid, set, options, database
624: .seealso: GridCreate(), GridPrintHelp(), GridSetOptionsPrefix()
625: @*/
626: int GridSetFromOptions(Grid grid)
627: {
628: Mesh mesh;
629: char *intNames[NUM_INTERPOLATIONS] = {"Local", "Halo", "L2", "L2Local"};
630: char intType[256];
631: PetscTruth opt;
632: int ierr;
636: PetscOptionsBegin(grid->comm, grid->prefix, "Grid options", "Grid");
638: /* Handle generic grid options */
639: PetscOptionsHasName(PETSC_NULL, "-help", &opt);
640: if (opt == PETSC_TRUE) {
641: GridPrintHelp(grid);
642: }
644: /* Constraint options */
645: PetscOptionsName("-grid_explicit_constraints", "Explicitly form constrained system", "GridGetExplicitConstraints", &grid->explicitConstraints);
647: /* Interpolation options */
648: PetscOptionsEList("-grid_int_type", "Interpolation Method", "None", intNames, NUM_INTERPOLATIONS,
649: intNames[grid->interpolationType], intType, 256, &opt);
650: if (opt == PETSC_TRUE) {
651: PetscTruth islocal, ishalo, isl2, isl2local;
653: PetscStrcasecmp(intType, intNames[INTERPOLATION_LOCAL], &islocal);
654: PetscStrcasecmp(intType, intNames[INTERPOLATION_HALO], &ishalo);
655: PetscStrcasecmp(intType, intNames[INTERPOLATION_L2], &isl2);
656: PetscStrcasecmp(intType, intNames[INTERPOLATION_L2_LOCAL], &isl2local);
657: if (islocal == PETSC_TRUE) {
658: grid->interpolationType = INTERPOLATION_LOCAL;
659: } else if (ishalo == PETSC_TRUE) {
660: grid->interpolationType = INTERPOLATION_HALO;
661: } else if (isl2 == PETSC_TRUE) {
662: grid->interpolationType = INTERPOLATION_L2;
663: } else if (isl2local == PETSC_TRUE) {
664: grid->interpolationType = INTERPOLATION_L2_LOCAL;
665: }
666: }
668: /* Graphics options */
669: PetscOptionsInt("-grid_view_field", "View Field", "None", grid->viewField, &grid->viewField, &opt);
670: PetscOptionsInt("-grid_view_comp", "View Component", "None", grid->viewComp, &grid->viewComp, &opt);
672: /* Handle grid type options */
673: GridSetTypeFromOptions(grid);
675: /* Handle grid serialization options */
676: GridSetSerializeTypeFromOptions(grid);
678: /* Handle specific grid options */
679: if (grid->ops->setfromoptions != PETSC_NULL) {
680: (*grid->ops->setfromoptions)(grid);
681: }
683: PetscOptionsEnd();
685: /* Handle subobject options */
686: GridGetMesh(grid, &mesh);
687: MeshSetFromOptions(mesh);
689: GridViewFromOptions(grid, grid->name);
690: return(0);
691: }
693: /*@
694: GridViewFromOptions - This function visualizes the grid based upon user options.
696: Collective on Grid
698: Input Parameter:
699: . grid - The grid
701: Level: intermediate
703: .keywords: Grid, view, options, database
704: .seealso: GridSetFromOptions(), GridView()
705: @*/
706: int GridViewFromOptions(Grid grid, char *title)
707: {
708: PetscViewer viewer;
709: PetscDraw draw;
710: PetscTruth opt;
711: char *titleStr;
712: char typeName[1024];
713: char fileName[1024];
714: int len;
715: int ierr;
718: PetscOptionsHasName(grid->prefix, "-grid_view", &opt);
719: if (opt == PETSC_TRUE) {
720: PetscOptionsGetString(grid->prefix, "-grid_view", typeName, 1024, &opt);
721: PetscStrlen(typeName, &len);
722: if (len > 0) {
723: PetscViewerCreate(grid->comm, &viewer);
724: PetscViewerSetType(viewer, typeName);
725: PetscOptionsGetString(grid->prefix, "-grid_view_file", fileName, 1024, &opt);
726: if (opt == PETSC_TRUE) {
727: PetscViewerSetFilename(viewer, fileName);
728: } else {
729: PetscViewerSetFilename(viewer, grid->name);
730: }
731: GridView(grid, viewer);
732: PetscViewerFlush(viewer);
733: PetscViewerDestroy(viewer);
734: } else {
735: GridView(grid, PETSC_NULL);
736: }
737: }
738: PetscOptionsHasName(grid->prefix, "-grid_view_draw", &opt);
739: if (opt == PETSC_TRUE) {
740: PetscViewerDrawOpen(grid->comm, 0, 0, 0, 0, 300, 300, &viewer);
741: PetscViewerDrawGetDraw(viewer, 0, &draw);
742: if (title != PETSC_NULL) {
743: titleStr = title;
744: } else {
745: PetscObjectName((PetscObject) grid); CHKERRQ(ierr) ;
746: titleStr = grid->name;
747: }
748: PetscDrawSetTitle(draw, titleStr);
749: GridView(grid, viewer);
750: PetscViewerFlush(viewer);
751: PetscDrawPause(draw);
752: PetscViewerDestroy(viewer);
753: }
754: return(0);
755: }
757: /*@
758: GridPrintHelp - Prints all options for the grid.
760: Input Parameter:
761: . grid - The grid
763: Options Database Keys:
764: $ -help, -h
766: Level: intermediate
768: .keywords: Grid, help
769: .seealso: GridSetFromOptions()
770: @*/
771: int GridPrintHelp(Grid grid)
772: {
773: char p[64];
774: int ierr;
779: PetscStrcpy(p, "-");
780: if (grid->prefix != PETSC_NULL) {
781: PetscStrcat(p, grid->prefix);
782: }
784: (*PetscHelpPrintf)(grid->comm, "Grid options ------------------------------------------------n");
785: (*PetscHelpPrintf)(grid->comm," %sgrid_type <typename> : Sets the grid typen", p);
786: (*PetscHelpPrintf)(grid->comm," %sgrid_view_field <field>: Sets the field to visualizen", p);
787: (*PetscHelpPrintf)(grid->comm," %sgrid_view_comp <field>: Sets the component to visualizen", p);
788: return(0);
789: }
791: /*@C
792: GridSetOptionsPrefix - Sets the prefix used for searching for all grid options in the database.
794: Input Parameters:
795: + grid - The grid
796: - prefix - The prefix to prepend to all option names
798: Notes:
799: A hyphen (-) must NOT be given at the beginning of the prefix name.
800: The first character of all runtime options is AUTOMATICALLY the hyphen.
802: Level: intermediate
804: .keywords: grid, set, options, prefix, database
805: .seealso: GridSetFromOptions()
806: @*/
807: int GridSetOptionsPrefix(Grid grid, char *prefix)
808: {
813: PetscObjectSetOptionsPrefix((PetscObject) grid, prefix);
814: return(0);
815: }
817: /*@C
818: GridAppendOptionsPrefix - Appends to the prefix used for searching for all grid options in the database.
820: Collective on Grid
822: Input Parameters:
823: + grid - The grid
824: - prefix - The prefix to prepend to all option names
826: Notes:
827: A hyphen (-) must NOT be given at the beginning of the prefix name.
828: The first character of all runtime options is AUTOMATICALLY the hyphen.
830: Level: intermediate
832: .keywords: grid, append, options, prefix, database
833: .seealso: GridGetOptionsPrefix()
834: @*/
835: int GridAppendOptionsPrefix(Grid grid, char *prefix)
836: {
838:
841: PetscObjectAppendOptionsPrefix((PetscObject) grid, prefix);
842: return(0);
843: }
845: /*@
846: GridGetOptionsPrefix - Sets the prefix used for searching for all grid options in the database.
848: Input Parameter:
849: . grid - The grid
851: Output Parameter:
852: . prefix - A pointer to the prefix string used
854: Level: intermediate
856: .keywords: grid, get, options, prefix, database
857: .seealso: GridAppendOptionsPrefix()
858: @*/
859: int GridGetOptionsPrefix(Grid grid, char **prefix)
860: {
865: PetscObjectGetOptionsPrefix((PetscObject) grid, prefix);
866: return(0);
867: }
869: /*----------------------------------------------- Graphics Functions ------------------------------------------------*/
870: /*@C GridSetViewField
871: This function sets the field component to be used in visualization.
873: Collective on Grid
875: Input Parameters:
876: + grid - The grid
877: . field - The field
878: - comp - The component
880: Level: intermediate
882: .keywords grid, view, field
883: .seealso GridGetMeshInfo
884: @*/
885: int GridSetViewField(Grid grid, int field, int comp)
886: {
889: GridValidField(grid, field);
890: if ((comp < 0) || (comp >= grid->fields[field].numComp)) {
891: SETERRQ2(PETSC_ERR_ARG_WRONG, "Invalid field component %d should be in [0,%d)", comp, grid->fields[field].numComp);
892: }
893: grid->viewField = field;
894: grid->viewComp = comp;
895: return(0);
896: }
898: /*------------------------------------------ Interface to Mesh Functions --------------------------------------------*/
899: /*@C
900: GridGetNodeClass - This function returns the class of a given node
902: Not collective
904: Input Parameters:
905: + grid - The grid
906: - node - The node
908: Output Parameter:
909: . nclass - The node class
911: Level: intermediate
913: .keywords grid, mesh, node, class
914: .seealso GridGetMeshInfo
915: @*/
916: int GridGetNodeClass(Grid grid, int node, int *nclass)
917: {
918: FieldClassMap map;
919: int ierr;
924: GridGetClassMap(grid, &map);
925: FieldClassMapGetNodeClass(map, node, nclass);
926: return(0);
927: }
929: /*@
930: GridGetNearestNode - This function returns the node nearest to the given point on which
931: the field is defined, or -1 if the node is not contained in the local mesh.
933: Not collective
935: Input Parameters:
936: + mesh - The mesh
937: . field - The field
938: - x,y,z - The node coordinates
940: Output Parameter:
941: . node - The nearest node
943: Level: beginner
945: .keywords: grid, node, point location
946: .seealso MeshNearestNode(), MeshLocatePoint()
947: @*/
948: int GridGetNearestNode(Grid grid, int field, double x, double y, double z, int *node)
949: {
950: FieldClassMap map;
951: Mesh mesh;
952: Partition p;
953: PetscTruth defined;
954: double minDist, dist, nx, ny, nz;
955: int minNode, globalMinNode, allMinNode;
956: int numCorners, startNode, endNode;
957: int elem, corner, nearNode, nclass;
958: int ierr;
963: FieldClassMapCreateTriangular2D(grid, 1, &field, &map);
964: GridGetMesh(grid, &mesh);
965: MeshGetNumCorners(mesh, &numCorners);
966: MeshGetPartition(mesh, &p);
967: PartitionGetStartNode(p, &startNode);
968: PartitionGetEndNode(p, &endNode);
969: MeshLocatePoint(mesh, x, y, z, &elem);
970: if (elem >= 0) {
971: /* Find first node with field defined */
972: for(corner = 0, minDist = PETSC_MAX; corner < numCorners; corner++) {
973: MeshGetNodeFromElement(mesh, elem, corner, &minNode);
974: FieldClassMapGetNodeClass(map, minNode, &nclass);
975: FieldClassMapIsDefined(map, field, nclass, &defined);
976: if (defined == PETSC_TRUE) {
977: MeshGetNodeCoords(mesh, minNode, &nx, &ny, &nz);
978: minDist = PetscSqr(MeshPeriodicDiffX(mesh, nx - x)) + PetscSqr(MeshPeriodicDiffY(mesh, ny - y));
979: break;
980: }
981: }
982: if (corner == numCorners) SETERRQ2(PETSC_ERR_ARG_WRONG, "Element %d has no node with field %d defined", elem, field);
983: /* Find closest node with field defined */
984: for(corner = 1; corner < numCorners; corner++) {
985: MeshGetNodeFromElement(mesh, elem, corner, &nearNode);
986: MeshGetNodeCoords(mesh, nearNode, &nx, &ny, &nz);
987: FieldClassMapGetNodeClass(map, nearNode, &nclass);
988: FieldClassMapIsDefined(map, field, nclass, &defined);
989: dist = PetscSqr(nx - x) + PetscSqr(ny - y);
990: if ((defined == PETSC_TRUE) && (dist < minDist)) {
991: minDist = dist;
992: minNode = nearNode;
993: }
994: }
995: PartitionLocalToGlobalNodeIndex(p, minNode, &globalMinNode);
996: } else {
997: minNode = -1;
998: globalMinNode = -1;
999: }
1000: FieldClassMapDestroy(map);
1002: /* The minimum node might still be a ghost node */
1003: MPI_Allreduce(&globalMinNode, &allMinNode, 1, MPI_INT, MPI_MAX, grid->comm);
1004: if ((allMinNode >= startNode) && (allMinNode < endNode)) {
1005: *node = allMinNode - startNode;
1006: } else {
1007: *node = -1;
1008: }
1009: if (allMinNode < 0) PetscFunctionReturn(1);
1010: return(0);
1011: }
1013: /*@
1014: GridGetNearestBdNode - This function returns the boundary node nearest to the given point
1015: on which the field is defined, or -1 if the node is not contained in the local mesh.
1017: Not collective
1019: Input Parameters:
1020: + mesh - The mesh
1021: . field - The field
1022: - x,y,z - The node coordinates
1024: Output Parameter:
1025: . node - The nearest boundary node
1027: Level: beginner
1029: .keywords: grid, node, point location
1030: .seealso MeshNearestBdNode(), MeshLocatePoint()
1031: @*/
1032: int GridGetNearestBdNode(Grid grid, int field, double x, double y, double z, int *node)
1033: {
1034: FieldClassMap map;
1035: Mesh mesh;
1036: Partition p;
1037: PetscTruth defined;
1038: double minDist, dist, nx, ny, nz;
1039: int minNode, globalMinNode, allMinNode;
1040: int numCorners, startNode, endNode;
1041: int elem, corner, nearNode, nclass, bd;
1042: int ierr;
1047: FieldClassMapCreateTriangular2D(grid, 1, &field, &map);
1048: GridGetMesh(grid, &mesh);
1049: MeshGetNumCorners(mesh, &numCorners);
1050: MeshGetPartition(mesh, &p);
1051: PartitionGetStartNode(p, &startNode);
1052: PartitionGetEndNode(p, &endNode);
1053: MeshLocatePoint(mesh, x, y, z, &elem);
1054: if (elem >= 0) {
1055: /* Find first node with field defined */
1056: for(corner = 0, minDist = PETSC_MAX; corner < numCorners; corner++) {
1057: MeshGetNodeFromElement(mesh, elem, corner, &minNode);
1058: MeshGetNodeBoundary(mesh, minNode, &bd);
1059: FieldClassMapGetNodeClass(map, minNode, &nclass);
1060: FieldClassMapIsDefined(map, field, nclass, &defined);
1061: if ((bd != 0) && (defined == PETSC_TRUE)) {
1062: MeshGetNodeCoords(mesh, minNode, &nx, &ny, &nz);
1063: minDist = PetscSqr(MeshPeriodicDiffX(mesh, nx - x)) + PetscSqr(MeshPeriodicDiffY(mesh, ny - y));
1064: break;
1065: }
1066: }
1067: if (corner == numCorners) SETERRQ2(PETSC_ERR_ARG_WRONG, "Element %d has no node with field %d defined", elem, field);
1068: /* Find closest node with field defined */
1069: for(corner = 1; corner < numCorners; corner++) {
1070: MeshGetNodeFromElement(mesh, elem, corner, &nearNode);
1071: MeshGetNodeCoords(mesh, nearNode, &nx, &ny, &nz);
1072: MeshGetNodeBoundary(mesh, nearNode, &bd);
1073: FieldClassMapGetNodeClass(map, nearNode, &nclass);
1074: FieldClassMapIsDefined(map, field, nclass, &defined);
1075: dist = PetscSqr(nx - x) + PetscSqr(ny - y);
1076: if ((bd != 0) && (defined == PETSC_TRUE) && (dist < minDist)) {
1077: minDist = dist;
1078: minNode = nearNode;
1079: }
1080: }
1081: PartitionLocalToGlobalNodeIndex(p, minNode, &globalMinNode);
1082: } else {
1083: minNode = -1;
1084: globalMinNode = -1;
1085: }
1086: FieldClassMapDestroy(map);
1088: /* The minimum node might still be a ghost node */
1089: MPI_Allreduce(&globalMinNode, &allMinNode, 1, MPI_INT, MPI_MAX, grid->comm);
1090: if ((allMinNode >= startNode) && (allMinNode < endNode)) {
1091: *node = allMinNode - startNode;
1092: } else {
1093: *node = -1;
1094: }
1095: if (allMinNode < 0)
1096: PetscFunctionReturn(1);
1097: return(0);
1098: }
1100: /*@
1101: GridGetBoundarySize - Retrieves the size of the specified boundary for
1102: a given field.
1104: Not collective
1106: Input Parameters:
1107: + grid - The grid
1108: . bd - The boundary marker
1109: - field - The field
1111: Output Parameter:
1112: . size - The size of the boundary
1114: Level: intermediate
1116: .keywords: grid, boundary, size
1117: .seealso: GridGetBoundaryNext(), GridGetBoundaryStart()
1118: @*/
1119: int GridGetBoundarySize(Grid grid, int bd, int field, int *size)
1120: {
1121: int bdIndex;
1127: GridValidField(grid, field);
1128: if (grid->setupcalled == PETSC_FALSE) SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Grid must be setup first")
1129: ierr = MeshGetBoundaryIndex(grid->mesh, bd, &bdIndex);
1130: *size = grid->bdSize[bdIndex][field];
1131: return(0);
1132: }
1134: /*@
1135: GridGetBoundaryStart - Retrieves the canonical node number of the first node
1136: on the specified boundary in a class contained in a given field.
1138: Not collective
1140: Input Parameters:
1141: + grid - The grid
1142: . bd - The boundary marker
1143: . field - The field
1144: - ghost - Flag for including ghost nodes
1146: Output Parameters:
1147: + node - The canonical node number
1148: - nclass - The node class
1150: Level: intermediate
1152: .keywords: grid, boundary, start
1153: .seealso: GridGetBoundaryNext(), MeshGetBoundaryStart
1154: @*/
1155: int GridGetBoundaryStart(Grid grid, int bd, int field, PetscTruth ghost, int *node, int *nclass)
1156: {
1157: int f;
1164: if (grid->setupcalled == PETSC_FALSE) SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Grid must be setup first");
1165: /* Look for field in default class map */
1166: for(f = 0; f < grid->cm->numFields; f++) {
1167: if (grid->cm->fields[f] == field) break;
1168: }
1169: if (f == grid->cm->numFields) {
1170: SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Field %d not present in the default class map", field);
1171: }
1172: (*grid->ops->getboundarystart)(grid, bd, f, ghost, grid->cm, node, nclass);
1173: return(0);
1174: }
1176: /*@
1177: GridGetBoundaryNext - Retrieves the canonical node number of the next node
1178: on the specified boundary in a class contained in a given field, with -1
1179: indicating boundary termination.
1181: Not collective
1183: Input Parameters:
1184: + grid - The grid
1185: . bd - The boundary marker
1186: . field - The field
1187: - ghost - Flag for including ghost nodes
1189: Output Parameters:
1190: + node - The canonical node number
1191: - nclass - The node class
1193: Level: intermediate
1195: .keywords: grid, boundary, next
1196: .seealso: GridGetBoundaryStart(), MeshGetBoundaryNext
1197: @*/
1198: int GridGetBoundaryNext(Grid grid, int bd, int field, PetscTruth ghost, int *node, int *nclass)
1199: {
1200: int f;
1207: if (grid->setupcalled == PETSC_FALSE) SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Grid must be setup first");
1208: /* Look for field in default class map */
1209: for(f = 0; f < grid->cm->numFields; f++) {
1210: if (grid->cm->fields[f] == field) break;
1211: }
1212: if (f == grid->cm->numFields) {
1213: SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Field %d not present in the default class map", field);
1214: }
1215: (*grid->ops->getboundarynext)(grid, bd, field, ghost, grid->cm, node, nclass);
1216: return(0);
1217: }
1219: /*@
1220: GridRefineMesh - Refine a grid based on area constraints.
1222: Collective on Grid
1224: Input Parameters:
1225: + grid - The initial grid
1226: - area - A function which gives an area constraint when evaluated inside an element
1228: Output Parameter:
1229: . newgrid - The refined grid
1231: Level: advanced
1233: .keywords: grid, refinement
1234: .seealso: MeshRefine(), GridReformMesh()
1235: @*/
1236: int GridRefineMesh(Grid grid, PointFunction area, Grid *newgrid)
1237: {
1238: Mesh newMesh;
1239: int field;
1240: int ierr;
1245: MeshRefine(grid->mesh, area, &newMesh);
1246: GridCreate(newMesh, newgrid);
1247: MeshDestroy(newMesh);
1248: for(field = 0; field < grid->numFields; field++) {
1249: GridAddField(*newgrid, grid->fields[field].name, grid->fields[field].discType, grid->fields[field].numComp, 0);
1250: }
1251: GridSetType(*newgrid, grid->type_name);
1252: /* Copy registered operators in discretizations */
1253: /* Set parent: (*newgrid)->gridparent = grid; */
1255: return(0);
1256: }
1258: /*------------------------------------------ Distributed Data Functions ---------------------------------------------*/
1259: /*@
1260: GridGlobalToLocalGeneral - Scatters values including ghost variables
1261: from the given global vector into the given ghost vector.
1263: Collective on Grid
1265: Input Parameters:
1266: + grid - The grid
1267: . vec - The global vector
1268: . mode - The insertion mode, either INSERT_VALUES or ADD_VALUES
1269: - scatter - The scatter
1271: Output Parameters:
1272: . ghostVec - The local vector
1274: Level: intermediate
1276: .keywords: grid, scatter, global, local
1277: .seealso: GridGlobalToLocal(), GridLocalToGlobal(), GVecGlobalToLocal()
1278: @*/
1279: int GridGlobalToLocalGeneral(Grid grid, GVec vec, Vec ghostVec, InsertMode mode, VecScatter scatter)
1280: {
1284: if (mode == ADD_VALUES) SETERRQ(PETSC_ERR_SUP, "No support for adding ghost values");
1285: VecScatterBegin(vec, ghostVec, mode, SCATTER_FORWARD, scatter);
1286: VecScatterEnd(vec, ghostVec, mode, SCATTER_FORWARD, scatter);
1287: return(0);
1288: }
1290: /*@
1291: GridGlobalToLocal - Scatters values including ghost variables
1292: from the given global vector into the local grid ghost vector.
1294: Collective on Grid
1296: Input Parameters:
1297: + grid - The grid
1298: . mode - The insertion mode, either INSERT_VALUES or ADD_VALUES
1299: - vec - The grid vector
1301: Level: intermediate
1303: .keywords: grid, scatter, global, local
1304: .seealso: GridGlobalToLocalGeneral(), GridLocalToGlobal(), GVecGlobalToLocal()
1305: @*/
1306: int GridGlobalToLocal(Grid grid, InsertMode mode, GVec vec)
1307: {
1308: Grid g;
1309: int ierr;
1315: GVecGetGrid(vec, &g);
1316: if (grid != g) SETERRQ(PETSC_ERR_ARG_INCOMP, "Vector does not arise from this grid");
1317: GridGlobalToLocalGeneral(grid, vec, grid->ghostVec, mode, grid->ghostScatter);
1318: return(0);
1319: }
1321: /*@
1322: GridLocalToGlobal - Scatters local processor values including ghost variables
1323: from the local grid ghost vector into the given global vector.
1325: Collective on Grid
1327: Input Parameters:
1328: + grid - The grid
1329: . mode - The insertion mode, either SET_VALUES or ADD_VALUES
1330: - vec - The grid vector
1332: Level: intermediate
1334: .seealso: GridGlobalToLocal(), GVecLocalToGlobal()
1335: @*/
1336: int GridLocalToGlobal(Grid grid, InsertMode mode, GVec vec)
1337: {
1338: Grid g;
1339: int ierr;
1345: GVecGetGrid(vec, &g);
1346: if (grid != g) SETERRQ(PETSC_ERR_ARG_WRONG, "Vector does not arise from this grid");
1347: VecScatterBegin(grid->ghostVec, vec, mode, SCATTER_REVERSE, grid->ghostScatter);
1348: VecScatterEnd(grid->ghostVec, vec, mode, SCATTER_REVERSE, grid->ghostScatter);
1349: return(0);
1350: }
1352: /*---------------------------------------------- InterGrid Functions ------------------------------------------------*/
1353: /*@
1354: GridInterpolateElementVec - Interpolates a given element vector from one discretization to another.
1356: Not collective
1358: Input Parameters:
1359: + grid - The original grid
1360: . field - The original field
1361: . vec - The original element vector
1362: . newGrid - The grid defining the new element vector
1363: - newField - The new field
1365: Output Parameter:
1366: . newVec - The interpolated element vector
1368: Level: advanced
1370: .keywords: grid, interpolation
1371: .seealso: ElementVecCreate()
1372: @*/
1373: int GridInterpolateElementVec(Grid grid, int field, ElementVec vec, Grid newGrid, int newField, ElementVec newVec)
1374: {
1382: GridValidField(grid, field);
1383: GridValidField(newGrid, newField);
1384: if (grid->fields[field].numComp != newGrid->fields[newField].numComp) {
1385: SETERRQ2(PETSC_ERR_ARG_INCOMP, "Fields have a different number of components %d != %d",
1386: grid->fields[field].numComp, newGrid->fields[newField].numComp);
1387: }
1388: DiscretizationInterpolateElementVec(grid->fields[field].disc, vec, newGrid->fields[newField].disc, newVec);
1389: return(0);
1390: }
1392: /*@
1393: GridCreateRestriction - Generates restriction matrix between two grids.
1395: Input Parameters:
1396: + vf - The fine grid
1397: - vc - The coarse grid
1399: Output Parameter:
1400: . gmat - A sparse matrix containing restriction
1402: Level: advanced
1404: @*/
1405: int GridCreateRestriction(Grid vf, Grid vc, GMat *gmat)
1406: {
1414: if (vf->setupcalled == PETSC_FALSE) {
1415: GridSetUp(vf);
1416: }
1417: if (vc->setupcalled == PETSC_FALSE) {
1418: GridSetUp(vc);
1419: }
1420: if (!vf->ops->gridcreaterestriction) SETERRQ(PETSC_ERR_SUP, " ");
1421: (*vf->ops->gridcreaterestriction)(vf, vc, gmat);
1422: return(0);
1423: }
1425: /*---------------------------------------------- Assembly Functions -------------------------------------------------*/
1426: /*@
1427: GridLocalToElementGeneral - Scatters values including ghost variables from the given ghost vector
1428: into the given element vector.
1430: Not collective
1432: Input Parameters:
1433: + grid - The grid
1434: . ghostVec - The vector of values (including ghsot points)
1435: . reduceVec - The vector of boundary values
1436: . reduceSystem - The flag for reducing boundary conditions
1437: . reduceElement - The flag for putting boundary values in the element vector
1438: - vec - The element vector
1440: WARNING:
1441: Make sure that the indices in the element vector are local indices.
1443: Note:
1444: If reduceSystem and reduceElement are PTESC_TRUE, then boundary values are
1445: placed in vec. If reduceElement is PETSC_FALSE, then zero is used instead.
1447: Level: advanced
1449: .keywords: grid, element, scatter
1450: .seealso: GridLocalToElement(), GridGlobalToLocal()
1451: @*/
1452: int GridLocalToElementGeneral(Grid grid, Vec ghostVec, Vec reduceVec, PetscTruth reduceSystem, PetscTruth reduceElement, ElementVec vec)
1453: {
1454: PetscScalar *array, *reduceArray;
1455: PetscScalar alpha = -grid->reduceAlpha;
1456: int elemSize = vec->reduceSize;
1457: int numOverlapVars, numReduceOverlapVars;
1458: int var;
1459: int ierr;
1462: VecGetSize(ghostVec, &numOverlapVars);
1463: VecGetArray(ghostVec, &array);
1464: if (reduceSystem == PETSC_TRUE) {
1465: VecGetSize(reduceVec, &numReduceOverlapVars);
1466: VecGetArray(reduceVec, &reduceArray);
1467: for(var = 0; var < elemSize; var++) {
1468: #ifdef PETSC_USE_BOPT_g
1469: if (vec->indices[var] >= numOverlapVars) {
1470: SETERRQ2(PETSC_ERR_ARG_CORRUPT, "Invalid index %d into ghost vector should be in [0,%d)",
1471: vec->indices[var], numOverlapVars);
1472: }
1473: #endif
1474: if (vec->indices[var] < 0) {
1475: if (reduceElement == PETSC_TRUE) {
1476: #ifdef PETSC_USE_BOPT_g
1477: if (-(vec->indices[var]+1) >= numReduceOverlapVars) {
1478: SETERRQ2(PETSC_ERR_ARG_CORRUPT, "Invalid index %d into ghost vector should be in [0,%d)",
1479: -(vec->indices[var]+1), numReduceOverlapVars);
1480: }
1481: #endif
1482: vec->array[var] = alpha*reduceArray[-(vec->indices[var]+1)];
1483: } else {
1484: vec->array[var] = 0.0;
1485: }
1486: } else {
1487: vec->array[var] = array[vec->indices[var]];
1488: }
1489: }
1490: VecRestoreArray(reduceVec, &reduceArray);
1491: } else {
1492: for(var = 0; var < elemSize; var++) {
1493: #ifdef PETSC_USE_BOPT_g
1494: if ((vec->indices[var] < 0) || (vec->indices[var] >= numOverlapVars)) {
1495: SETERRQ2(PETSC_ERR_ARG_CORRUPT, "Invalid index %d into ghost vector should be in [0,%d)",
1496: vec->indices[var], numOverlapVars);
1497: }
1498: #endif
1499: vec->array[var] = array[vec->indices[var]];
1500: }
1501: }
1502: VecRestoreArray(ghostVec, &array);
1503: return(0);
1504: }
1506: /*@
1507: GridLocalToElement - Scatters values including ghost variables from the local grid ghost vector
1508: into the given element vector.
1510: Not collective
1512: Input Parameters:
1513: + grid - The grid
1514: - vec - The element vector
1516: WARNING:
1517: Make sure that the indices in the element vector are local indices.
1519: Level: advanced
1521: .keywords: grid, element, scatter
1522: .seealso: GridLocalToElementGeneral(), GridGlobalToLocal()
1523: @*/
1524: int GridLocalToElement(Grid grid, ElementVec vec)
1525: {
1531: GridLocalToElementGeneral(grid, grid->ghostVec, grid->bdReduceVecCur, grid->reduceSystem, grid->reduceElement, vec);
1532: return(0);
1533: }
1535: /*-------------------------------------------- Evaluation Functions -------------------------------------------------*/
1536: /*@C
1537: GridEvaluateRhs - This function constructs the weak form of the functions and
1538: operators specified with GridSetRhsFunction(), GridSetRhsOperator(),
1539: and GridSetRhsNonlinearOperator().
1541: Collective on Grid
1543: Input Parameters:
1544: + grid - The grid
1545: . x - The current iterate
1546: - ctx - The optional user context
1548: Output Parameters
1549: . f - The function value
1551: Level: beginner
1553: .keywords: grid, rhs
1554: .seealso: GridEvaluateSystemMatrix(), GridSetRhsFunction(), GridSetRhsOperator(), GridSetRhsNonlinearOperator()
1555: @*/
1556: int GridEvaluateRhs(Grid grid, GVec x, GVec f, PetscObject ctx)
1557: {
1558: Mesh mesh;
1559: MeshMover mover;
1560: Grid meshVelGrid;
1561: PetscScalar alpha = -grid->reduceAlpha;
1562: int ierr;
1566: /* x can be PETSC_NULL */
1568: /* Setup mesh velocity calculation for ALE variables */
1569: if (grid->ALEActive == PETSC_TRUE) {
1570: GridGetMesh(grid, &mesh);
1571: MeshGetMover(mesh, &mover);
1572: MeshMoverGetVelocityGrid(mover, &meshVelGrid);
1573: GridSetUp(meshVelGrid);
1574: }
1575: /* Calculate Rhs */
1576: (*grid->ops->gridevaluaterhs)(grid, x, f, ctx);
1577: /* Add boundary values */
1578: if ((grid->reduceSystem == PETSC_TRUE) && (grid->reduceElement == PETSC_FALSE) && (grid->bdReduceMat != PETSC_NULL)) {
1579: /* Calculate contribution of constrained variables to the Rhs */
1580: MatMult(grid->bdReduceMat, grid->bdReduceVecCur, grid->reduceVec);
1581: /* Update Rhs */
1582: VecAXPY(&alpha, grid->reduceVec, f);
1583: }
1584: return(0);
1585: }
1587: /*@C GridEvaluateRhsFunction
1588: This function constructs the weak form of the functions specified with GridSetRhsFunction().
1590: Collective on Grid
1592: Input Parameters:
1593: + grid - The grid
1594: . x - The current iterate
1595: - ctx - The optional user context
1597: Output Parameter:
1598: . f - The function value
1600: Level: beginner
1602: .keywords grid, rhs, function
1603: .seealso GridEvaluateRhsOperator(), GridEvaluateRhs(), GridEvaluateSystemMatrix()
1604: @*/
1605: int GridEvaluateRhsFunction(Grid grid, GVec x, GVec f, void *ctx)
1606: {
1607: PetscTruth oldLinear = grid->activeOpTypes[1];
1608: PetscTruth oldNonlinear = grid->activeOpTypes[2];
1609: int ierr;
1613: /* x can be PETSC_NULL */
1615: /* Turn off computation of operators */
1616: grid->activeOpTypes[1] = PETSC_FALSE;
1617: grid->activeOpTypes[2] = PETSC_FALSE;
1618: GridEvaluateRhs(grid, x, f, (PetscObject) ctx);
1619: grid->activeOpTypes[1] = oldLinear;
1620: grid->activeOpTypes[2] = oldNonlinear;
1621: return(0);
1622: }
1624: /*@C GridEvaluateRhsOperator
1625: This function constructs the weak form of the operators specified with GridSetRhsOperator().
1627: Collective on Grid
1629: Input Parameters:
1630: + grid - The grid
1631: . x - The current iterate
1632: . linear - The flag for including linear operators
1633: . nonlinear - The flag for including nonlinear operators
1634: - ctx - The optional user context
1636: Output Parameter:
1637: . f - The operator value
1639: Level: beginner
1641: .keywords grid, rhs, operator
1642: .seealso GridEvaluateRhsFunction(), GridEvaluateRhs(), GridEvaluateSystemMatrix()
1643: @*/
1644: int GridEvaluateRhsOperator(Grid grid, GVec x, GVec f, PetscTruth linear, PetscTruth nonlinear, void *ctx)
1645: {
1646: PetscTruth oldFunction = grid->activeOpTypes[0];
1647: PetscTruth oldLinear = grid->activeOpTypes[1];
1648: PetscTruth oldNonlinear = grid->activeOpTypes[2];
1649: int ierr;
1653: /* x can be PETSC_NULL */
1655: /* Turn off computation of operators */
1656: grid->activeOpTypes[0] = PETSC_FALSE;
1657: grid->activeOpTypes[1] = linear;
1658: grid->activeOpTypes[2] = nonlinear;
1659: GridEvaluateRhs(grid, x, f, (PetscObject) ctx);
1660: grid->activeOpTypes[0] = oldFunction;
1661: grid->activeOpTypes[1] = oldLinear;
1662: grid->activeOpTypes[2] = oldNonlinear;
1663: return(0);
1664: }
1666: /*@C GridEvaluateSystemMatrix
1667: This function constructs the weak form of the linear operators
1668: specified with GridSetMatOperator().
1670: Collective on Grid
1672: Input Parameters:
1673: + grid - The grid
1674: . x - The current iterate
1675: - ctx - The optional user context
1677: Output Parameter:
1678: . f - The function value
1680: Level: beginner
1682: .keywords: grid, matrix
1683: .seealso: GridEvaluateRhs()
1684: @*/
1685: int GridEvaluateSystemMatrix(Grid grid, GVec x, GMat *J, GMat *M, MatStructure *flag, PetscObject ctx)
1686: {
1687: Mesh mesh;
1688: MeshMover mover;
1689: Grid meshVelGrid;
1690: #ifdef PETSC_HAVE_MATHEMATICA
1691: PetscTruth opt;
1692: #endif
1693: int ierr;
1697: if (grid->ALEActive == PETSC_TRUE) {
1698: GridGetMesh(grid, &mesh);
1699: MeshGetMover(mesh, &mover);
1700: MeshMoverGetVelocityGrid(mover, &meshVelGrid);
1701: GridSetUp(meshVelGrid);
1702: }
1703: if (grid->isMatrixFree == PETSC_TRUE) {
1704: /* Should probably set this as the matrix context, but there is no MatShellSetContext() */
1705: grid->matrixFreeArg = x;
1706: grid->matrixFreeContext = ctx;
1707: } else {
1708: (*grid->ops->gridevaluatesystemmatrix)(grid, x, J, M, flag, ctx);
1709: }
1710: #ifdef PETSC_HAVE_MATHEMATICA
1711: /* Debugging */
1712: PetscOptionsHasName(PETSC_NULL, "-trace_grid_math", &opt);
1713: if (opt == PETSC_TRUE) {
1714: ViewerMathematicaSetName(PETSC_VIEWER_MATHEMATICA_WORLD, "jac");
1715: MatView(*J, PETSC_VIEWER_MATHEMATICA_WORLD);
1716: ViewerMathematicaClearName(PETSC_VIEWER_MATHEMATICA_WORLD);
1717: }
1718: #endif
1719: return(0);
1720: }
1722: #if 0
1723: /*----------------------------------------- Parallel Communication Functions ----------------------------------------*/
1724: /*@C
1725: GridGhostExchange - This functions transfers data between local and ghost storage without a predefined mapping.
1727: Collective on MPI_Comm
1729: Input Parameters:
1730: + comm - The communicator
1731: . numGhosts - The number of ghosts in this domain
1732: . ghostProcs - The processor from which to obtain each ghost
1733: . ghostIndices - The global index for each ghost
1734: . dataType - The type of the variables
1735: . firstVar - The first variable on each processor
1736: . addv - The insert mode, INSERT_VALUES or ADD_VALUES
1737: - mode - The direction of the transfer, SCATTER_FORWARD or SCATTER_REVERSE
1739: Output Paramters:
1740: + locVars - The local variable array
1741: - ghostVars - The ghost variables
1743: Note:
1744: The data in ghostVars is assumed contiguous and implicitly indexed by the order of
1745: ghostProcs and ghostIndices. The SCATTER_FORWARD mode will take the requested data
1746: from locVars and copy it to ghostVars in the order specified by ghostIndices. The
1747: SCATTER_REVERSE mode will take data from ghostVars and copy it to locVars.
1749: Level: advanced
1751: .keywords: ghost, exchange, grid
1752: .seealso: GridGlobalToLocal(), GridLocalToGlobal()
1753: @*/
1754: int GridGhostExchange(MPI_Comm comm, int numGhosts, int *ghostProcs, int *ghostIndices, PetscDataType dataType,
1755: int *firstVar, InsertMode addv, ScatterMode mode, void *locVars, void *ghostVars)
1756: {
1757: int *numSendGhosts; /* The number of ghosts from each domain */
1758: int *numRecvGhosts; /* The number of local variables which are ghosts in each domain */
1759: int *sumSendGhosts; /* The prefix sums of numSendGhosts */
1760: int *sumRecvGhosts; /* The prefix sums of numRecvGhosts */
1761: int *offsets; /* The offset into the send array for each domain */
1762: int totSendGhosts; /* The number of ghosts to request variables for */
1763: int totRecvGhosts; /* The number of nodes to provide class info about */
1764: int *sendIndices; /* The canonical indices of ghosts in this domain */
1765: int *recvIndices; /* The canonical indices of ghosts to return variables for */
1766: char *tempVars; /* The variables of the requested or submitted ghosts */
1767: char *locBytes = (char *) locVars;
1768: MPI_Datatype MPIType;
1769: int typeSize;
1770: #ifdef PETSC_USE_BOPT_g
1771: int numLocVars;
1772: #endif
1773: int numProcs, rank;
1774: int proc, ghost, locIndex, byte;
1775: int ierr;
1778: /* Initialize communication */
1779: MPI_Comm_size(comm, &numProcs);
1780: MPI_Comm_rank(comm, &rank);
1781: PetscMalloc(numProcs * sizeof(int), &numSendGhosts);
1782: PetscMalloc(numProcs * sizeof(int), &numRecvGhosts);
1783: PetscMalloc(numProcs * sizeof(int), &sumSendGhosts);
1784: PetscMalloc(numProcs * sizeof(int), &sumRecvGhosts);
1785: PetscMalloc(numProcs * sizeof(int), &offsets);
1786: PetscMemzero(numSendGhosts, numProcs * sizeof(int));
1787: PetscMemzero(numRecvGhosts, numProcs * sizeof(int));
1788: PetscMemzero(sumSendGhosts, numProcs * sizeof(int));
1789: PetscMemzero(sumRecvGhosts, numProcs * sizeof(int));
1790: PetscMemzero(offsets, numProcs * sizeof(int));
1791: #ifdef PETSC_USE_BOPT_g
1792: numLocVars = firstVar[rank+1] - firstVar[rank];
1793: #endif
1795: /* Get number of ghosts needed from each processor */
1796: for(ghost = 0; ghost < numGhosts; ghost++)
1797: numSendGhosts[ghostProcs[ghost]]++;
1799: /* Get number of ghosts to provide variables for */
1800: MPI_Alltoall(numSendGhosts, 1, MPI_INT, numRecvGhosts, 1, MPI_INT, comm);
1801: for(proc = 1; proc < numProcs; proc++) {
1802: sumSendGhosts[proc] = sumSendGhosts[proc-1] + numSendGhosts[proc-1];
1803: sumRecvGhosts[proc] = sumRecvGhosts[proc-1] + numRecvGhosts[proc-1];
1804: offsets[proc] = sumSendGhosts[proc];
1805: }
1806: totSendGhosts = sumSendGhosts[numProcs-1] + numSendGhosts[numProcs-1];
1807: totRecvGhosts = sumRecvGhosts[numProcs-1] + numRecvGhosts[numProcs-1];
1808: if (numGhosts != totSendGhosts) {
1809: SETERRQ2(PETSC_ERR_PLIB, "Invalid number of ghosts %d in send, should be %d", totSendGhosts, numGhosts);
1810: }
1812: PetscDataTypeGetSize(dataType, &typeSize);
1813: if (totSendGhosts) {
1814: sendIndices = (int *) PetscMalloc(totSendGhosts * sizeof(int)); CHKERRQ(sendIndices);
1815: }
1816: if (totRecvGhosts) {
1817: recvIndices = (int *) PetscMalloc(totRecvGhosts * sizeof(int)); CHKERRQ(recvIndices);
1818: tempVars = (char *) PetscMalloc(totRecvGhosts * typeSize); CHKERRQ(tempVars);
1819: }
1821: /* Must order ghosts by processor */
1822: for(ghost = 0; ghost < numGhosts; ghost++)
1823: sendIndices[offsets[ghostProcs[ghost]]++] = ghostIndices[ghost];
1825: /* Get canonical indices of ghosts to provide variables for */
1826: MPI_Alltoallv(sendIndices, numSendGhosts, sumSendGhosts, MPI_INT,
1827: recvIndices, numRecvGhosts, sumRecvGhosts, MPI_INT, comm);
1828:
1830: switch(mode)
1831: {
1832: case SCATTER_FORWARD:
1833: /* Get ghost variables */
1834: if (addv == INSERT_VALUES) {
1835: for(ghost = 0; ghost < totRecvGhosts; ghost++)
1836: {
1837: locIndex = recvIndices[ghost] - firstVar[rank];
1838: #ifdef PETSC_USE_BOPT_g
1839: if ((locIndex < 0) || (locIndex >= numLocVars)) {
1840: SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE, "Invalid ghost index %d, not in [0,%d)", locIndex, numLocVars);
1841: }
1842: #endif
1843: for(byte = 0; byte < typeSize; byte++)
1844: tempVars[ghost*typeSize+byte] = locBytes[locIndex*typeSize+byte];
1845: }
1846: } else {
1847: for(ghost = 0; ghost < totRecvGhosts; ghost++)
1848: {
1849: locIndex = recvIndices[ghost] - firstVar[rank];
1850: #ifdef PETSC_USE_BOPT_g
1851: if ((locIndex < 0) || (locIndex >= numLocVars)) {
1852: SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE, "Invalid ghost index %d, not in [0,%d)", locIndex, numLocVars);
1853: }
1854: #endif
1855: for(byte = 0; byte < typeSize; byte++)
1856: tempVars[ghost*typeSize+byte] += locBytes[locIndex*typeSize+byte];
1857: }
1858: }
1860: /* Communicate local variables to ghost storage */
1861: PetscDataTypeToMPIDataType(dataType, &MPIType);
1862: MPI_Alltoallv(tempVars, numRecvGhosts, sumRecvGhosts, MPIType,
1863: ghostVars, numSendGhosts, sumSendGhosts, MPIType, comm);
1864:
1865: break;
1866: case SCATTER_REVERSE:
1867: /* Communicate ghost variables to local storage */
1868: PetscDataTypeToMPIDataType(dataType, &MPIType);
1869: MPI_Alltoallv(ghostVars, numSendGhosts, sumSendGhosts, MPIType,
1870: tempVars, numRecvGhosts, sumRecvGhosts, MPIType, comm);
1871:
1873: /* Get ghost variables */
1874: if (addv == INSERT_VALUES) {
1875: for(ghost = 0; ghost < totRecvGhosts; ghost++)
1876: {
1877: locIndex = recvIndices[ghost] - firstVar[rank];
1878: #ifdef PETSC_USE_BOPT_g
1879: if ((locIndex < 0) || (locIndex >= numLocVars)) {
1880: SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE, "Invalid ghost index %d, not in [0,%d)", locIndex, numLocVars);
1881: }
1882: #endif
1883: for(byte = 0; byte < typeSize; byte++)
1884: locBytes[locIndex*typeSize+byte] = tempVars[ghost*typeSize+byte];
1885: }
1886: } else {
1887: /* There must be a better way to do this -- Ask Bill */
1888: if (typeSize == sizeof(int)) {
1889: int *tempInt = (int *) tempVars;
1890: int *locInt = (int *) locVars;
1892: for(ghost = 0; ghost < totRecvGhosts; ghost++) {
1893: locIndex = recvIndices[ghost] - firstVar[rank];
1894: #ifdef PETSC_USE_BOPT_g
1895: if ((locIndex < 0) || (locIndex >= numLocVars)) {
1896: SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE, "Invalid ghost index %d, not in [0,%d)", locIndex, numLocVars);
1897: }
1898: #endif
1899: locInt[locIndex] += tempInt[ghost];
1900: }
1901: } else if (typeSize == sizeof(long int)) {
1902: long int *tempLongInt = (long int *) tempVars;
1903: long int *locLongInt = (long int *) locVars;
1905: for(ghost = 0; ghost < totRecvGhosts; ghost++) {
1906: locIndex = recvIndices[ghost] - firstVar[rank];
1907: #ifdef PETSC_USE_BOPT_g
1908: if ((locIndex < 0) || (locIndex >= numLocVars)) {
1909: SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE, "Invalid ghost index %d, not in [0,%d)", locIndex, numLocVars);
1910: }
1911: #endif
1912: locLongInt[locIndex] += tempLongInt[ghost];
1913: }
1914: } else {
1915: for(ghost = 0; ghost < totRecvGhosts; ghost++) {
1916: locIndex = recvIndices[ghost] - firstVar[rank];
1917: #ifdef PETSC_USE_BOPT_g
1918: if ((locIndex < 0) || (locIndex >= numLocVars)) {
1919: SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE, "Invalid ghost index %d, not in [0,%d)", locIndex, numLocVars);
1920: }
1921: #endif
1922: for(byte = 0; byte < typeSize; byte++)
1923: locBytes[locIndex*typeSize+byte] += tempVars[ghost*typeSize+byte];
1924: }
1925: }
1926: }
1927: break;
1928: default:
1929: SETERRQ1(PETSC_ERR_ARG_WRONG, "Invalid scatter mode %d", mode);
1930: }
1932: /* Cleanup */
1933: PetscFree(numSendGhosts);
1934: PetscFree(numRecvGhosts);
1935: PetscFree(sumSendGhosts);
1936: PetscFree(sumRecvGhosts);
1937: PetscFree(offsets);
1938: if (totSendGhosts) {
1939: PetscFree(sendIndices);
1940: }
1941: if (totRecvGhosts) {
1942: PetscFree(recvIndices);
1943: PetscFree(tempVars);
1944: }
1945: return(0);
1946: }
1947: #endif