Actual source code: cartesian.c
1: #include <private/meshimpl.h> /*I "petscdmmesh.h" I*/
2: #include <CartesianSieve.hh>
6: /*@C
7: DMCartesianGetMesh - Gets the internal mesh object
9: Not collective
11: Input Parameter:
12: . dm - the mesh object
14: Output Parameter:
15: . m - the internal mesh object
17: Level: advanced
19: .seealso MeshCreate(), MeshCartesianSetMesh()
20: @*/
21: PetscErrorCode DMCartesianGetMesh(DM dm, ALE::Obj<ALE::CartesianMesh>& m)
22: {
23: DM_Cartesian *c = (DM_Cartesian *) dm->data;
27: m = c->m;
28: return(0);
29: }
33: /*@C
34: DMCartesianSetMesh - Sets the internal mesh object
36: Not collective
38: Input Parameters:
39: + mesh - the mesh object
40: - m - the internal mesh object
42: Level: advanced
44: .seealso MeshCreate(), MeshCartesianGetMesh()
45: @*/
46: PetscErrorCode DMCartesianSetMesh(DM dm, const ALE::Obj<ALE::CartesianMesh>& m)
47: {
48: DM_Cartesian *c = (DM_Cartesian *) dm->data;
52: c->m = m;
53: return(0);
54: }
58: PetscErrorCode DMDestroy_Cartesian(DM dm)
59: {
60: DM_Cartesian *c = (DM_Cartesian *) dm->data;
63: c->m = PETSC_NULL;
64: return(0);
65: }
69: PetscErrorCode DMView_Cartesian_Ascii(const ALE::Obj<ALE::CartesianMesh>& mesh, PetscViewer viewer)
70: {
71: PetscViewerFormat format;
72: PetscErrorCode ierr;
75: PetscViewerGetFormat(viewer, &format);
76: if (format == PETSC_VIEWER_ASCII_VTK) {
77: #if 0
78: VTKViewer::writeHeader(viewer);
79: VTKViewer::writeVertices(mesh, viewer);
80: VTKViewer::writeElements(mesh, viewer);
81: #endif
82: } else {
83: int dim = mesh->getDimension();
85: PetscViewerASCIIPrintf(viewer, "Mesh in %d dimensions:\n", dim);
86: /* FIX: Need to globalize */
87: PetscViewerASCIIPrintf(viewer, " %d vertices\n", mesh->getSieve()->getNumVertices());
88: PetscViewerASCIIPrintf(viewer, " %d cells\n", mesh->getSieve()->getNumCells());
89: }
90: PetscViewerFlush(viewer);
91: return(0);
92: }
94: PetscErrorCode DMView_Cartesian(DM dm, PetscViewer viewer)
95: {
96: ALE::Obj<ALE::CartesianMesh> m;
97: PetscBool iascii, isbinary, isdraw;
101: PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
102: PetscTypeCompare((PetscObject) viewer, PETSCVIEWERBINARY, &isbinary);
103: PetscTypeCompare((PetscObject) viewer, PETSCVIEWERDRAW, &isdraw);
105: DMCartesianGetMesh(dm, m);
106: if (iascii){
107: DMView_Cartesian_Ascii(m, viewer);
108: } else if (isbinary) {
109: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "Binary viewer not implemented for Cartesian Mesh");
110: } else if (isdraw){
111: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "Draw viewer not implemented for Cartesian Mesh");
112: } else {
113: SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported by this mesh object", ((PetscObject)viewer)->type_name);
114: }
115: return(0);
116: }
120: PetscErrorCode DMGetInterpolation_Cartesian(DM fineMesh, DM coarseMesh, Mat *interpolation, Vec *scaling)
121: {
122: ALE::Obj<ALE::CartesianMesh> coarse;
123: ALE::Obj<ALE::CartesianMesh> fine;
124: Mat P;
125: PetscErrorCode ierr;
128: DMCartesianGetMesh(fineMesh, fine);
129: DMCartesianGetMesh(coarseMesh, coarse);
130: #if 0
131: const ALE::Obj<ALE::Mesh::real_section_type>& coarseCoordinates = coarse->getRealSection("coordinates");
132: const ALE::Obj<ALE::Mesh::real_section_type>& fineCoordinates = fine->getRealSection("coordinates");
133: const ALE::Obj<ALE::Mesh::label_sequence>& vertices = fine->depthStratum(0);
134: const ALE::Obj<ALE::Mesh::real_section_type>& sCoarse = coarse->getRealSection("default");
135: const ALE::Obj<ALE::Mesh::real_section_type>& sFine = fine->getRealSection("default");
136: const ALE::Obj<ALE::Mesh::order_type>& coarseOrder = coarse->getFactory()->getGlobalOrder(coarse, "default", sCoarse);
137: const ALE::Obj<ALE::Mesh::order_type>& fineOrder = fine->getFactory()->getGlobalOrder(fine, "default", sFine);
138: const int dim = coarse->getDimension();
139: double *v0, *J, *invJ, detJ, *refCoords, *values;
140: #endif
142: MatCreate(fine->comm(), &P);
143: #if 0
144: MatSetSizes(P, sFine->size(), sCoarse->size(), PETSC_DETERMINE, PETSC_DETERMINE);
145: MatSetFromOptions(P);
146: PetscMalloc5(dim,double,&v0,dim*dim,double,&J,dim*dim,double,&invJ,dim,double,&refCoords,dim+1,double,&values);
147: for(ALE::Mesh::label_sequence::iterator v_iter = vertices->begin(); v_iter != vertices->end(); ++v_iter) {
148: const ALE::Mesh::real_section_type::value_type *coords = fineCoordinates->restrictPoint(*v_iter);
149: const ALE::Mesh::point_type coarseCell = coarse->locatePoint(coords);
151: coarse->computeElementGeometry(coarseCoordinates, coarseCell, v0, J, invJ, detJ);
152: for(int d = 0; d < dim; ++d) {
153: refCoords[d] = 0.0;
154: for(int e = 0; e < dim; ++e) {
155: refCoords[d] += invJ[d*dim+e]*(coords[e] - v0[e]);
156: }
157: refCoords[d] -= 1.0;
158: }
159: values[0] = 1.0/3.0 - (refCoords[0] + refCoords[1])/3.0;
160: values[1] = 0.5*(refCoords[0] + 1.0);
161: values[2] = 0.5*(refCoords[1] + 1.0);
162: updateOperatorGeneral(P, fine, sFine, fineOrder, *v_iter, sCoarse, coarseOrder, coarseCell, values, INSERT_VALUES);
163: }
164: PetscFree5(v0,J,invJ,refCoords,values);
165: #endif
166: MatAssemblyBegin(P, MAT_FINAL_ASSEMBLY);
167: MatAssemblyEnd(P, MAT_FINAL_ASSEMBLY);
168: *interpolation = P;
169: return(0);
170: }
174: PetscErrorCode DMRefine_Cartesian(DM mesh, MPI_Comm comm, DM *refinedMesh)
175: {
176: ALE::Obj<ALE::CartesianMesh> oldMesh;
177: PetscErrorCode ierr;
180: DMCartesianGetMesh(mesh, oldMesh);
181: DMCartesianCreate(comm, refinedMesh);
182: #if 0
183: ALE::Obj<ALE::Mesh> newMesh = ALE::Generator::refineMesh(oldMesh, refinementLimit, false);
184: MeshCartesianSetMesh(*refinedMesh, newMesh);
185: const ALE::Obj<ALE::CartesianMesh::real_section_type>& s = newMesh->getRealSection("default");
187: newMesh->setDiscretization(oldMesh->getDiscretization());
188: newMesh->setBoundaryCondition(oldMesh->getBoundaryCondition());
189: newMesh->setupField(s);
190: #else
191: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "Not yet implemented");
192: #endif
193: return(0);
194: }
198: PetscErrorCode DMCoarsen_Cartesian(DM mesh, MPI_Comm comm, DM *coarseMesh)
199: {
200: ALE::Obj<ALE::CartesianMesh> oldMesh;
201: PetscErrorCode ierr;
204: DMCartesianGetMesh(mesh, oldMesh);
205: DMCartesianCreate(comm, coarseMesh);
206: #if 0
207: ALE::Obj<ALE::Mesh> newMesh = ALE::Generator::refineMesh(oldMesh, refinementLimit, false);
208: MeshCartesianSetMesh(*coarseMesh, newMesh);
209: const ALE::Obj<ALE::CartesianMesh::real_section_type>& s = newMesh->getRealSection("default");
211: newMesh->setDiscretization(oldMesh->getDiscretization());
212: newMesh->setBoundaryCondition(oldMesh->getBoundaryCondition());
213: newMesh->setupField(s);
214: #else
215: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP, "Not yet implemented");
216: #endif
217: return(0);
218: }
222: PetscErrorCode DMCartesianGetSectionReal(DM dm, const char name[], SectionReal *section)
223: {
224: ALE::Obj<ALE::CartesianMesh> m;
225: PetscErrorCode ierr;
228: DMCartesianGetMesh(dm, m);
229: SectionRealCreate(m->comm(), section);
230: PetscObjectSetName((PetscObject) *section, name);
231: #if 0
232: SectionRealSetSection(*section, m->getRealSection(std::string(name)));
233: SectionRealSetBundle(*section, m);
234: #endif
235: return(0);
236: }
240: PetscErrorCode DMSetFromOptions_Cartesian(DM dm)
241: {
242: /* DM_Mesh *mesh = (DM_Mesh *) dm->data; */
243: char typeName[256];
244: PetscBool flg;
249: PetscOptionsBegin(((PetscObject) dm)->comm, ((PetscObject) dm)->prefix, "DMCartesian Options", "DMCartesian");
250: /* Handle DMCartesian refinement */
251: /* Handle associated vectors */
252: if (!VecRegisterAllCalled) {VecRegisterAll(PETSC_NULL);}
253: PetscOptionsList("-dm_vec_type", "Vector type used for created vectors", "DMSetVecType", VecList, dm->vectype, typeName, 256, &flg);
254: if (flg) {
255: DMSetVecType(dm, typeName);
256: }
257: /* process any options handlers added with PetscObjectAddOptionsHandler() */
258: PetscObjectProcessOptionsHandlers((PetscObject) dm);
259: PetscOptionsEnd();
260: return(0);
261: }
266: PetscErrorCode DMCreate_Cartesian(DM dm)
267: {
268: DM_Cartesian *mesh;
273: PetscNewLog(dm, DM_Cartesian, &mesh);
274: dm->data = mesh;
276: new(&mesh->m) ALE::Obj<ALE::CartesianMesh>(PETSC_NULL);
278: PetscStrallocpy(VECSTANDARD, &dm->vectype);
279: dm->ops->globaltolocalbegin = 0;
280: dm->ops->globaltolocalend = 0;
281: dm->ops->localtoglobalbegin = 0;
282: dm->ops->localtoglobalend = 0;
283: dm->ops->createglobalvector = 0; /* DMCreateGlobalVector_Cartesian; */
284: dm->ops->createlocalvector = 0; /* DMCreateLocalVector_Cartesian; */
285: dm->ops->getinterpolation = DMGetInterpolation_Cartesian;
286: dm->ops->getcoloring = 0;
287: dm->ops->getmatrix = 0; /* DMGetMatrix_Cartesian; */
288: dm->ops->refine = DMRefine_Cartesian;
289: dm->ops->coarsen = DMCoarsen_Cartesian;
290: dm->ops->refinehierarchy = 0;
291: dm->ops->coarsenhierarchy = 0;
292: dm->ops->getinjection = 0;
293: dm->ops->getaggregates = 0;
294: dm->ops->destroy = DMDestroy_Cartesian;
295: dm->ops->view = DMView_Cartesian;
296: dm->ops->setfromoptions = DMSetFromOptions_Cartesian;
297: dm->ops->setup = 0;
298: return(0);
299: }
304: /*@
305: DMCartesianCreate - Creates a DMCartesian object.
307: Collective on MPI_Comm
309: Input Parameter:
310: . comm - The communicator for the DMCartesian object
312: Output Parameter:
313: . mesh - The DMCartesian object
315: Level: beginner
317: .keywords: DMCartesian, create
318: @*/
319: PetscErrorCode DMCartesianCreate(MPI_Comm comm, DM *mesh)
320: {
325: DMCreate(comm, mesh);
326: DMSetType(*mesh, DMCARTESIAN);
327: return(0);
328: }