<div dir="ltr"><div dir="ltr">On Wed, Oct 13, 2021 at 6:30 PM Abhishek G.S. <<a href="mailto:gsabhishek1ags@gmail.com">gsabhishek1ags@gmail.com</a>> wrote:<br></div><div class="gmail_quote"><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div>Hi,</div><div>I need some help with getting the file output working right.<br></div><div><br></div><div></div><div>
I am using a DMDACreate3D to initialize my DM.
This is my write function
</div><div><br></div><div>
<div>void write(){</div><div>PetscViewer viewer;<br>PetscViewerHDF5Open(PETSC_COMM_WORLD,filename.c_str(),FILE_MODE_WRITE,&viewer);</div><div>
DMDAVecRestoreArray(dm,global_vector,global_array)</div><div>VecView(global_vec, viewer);</div><div>
DMDAVecGetArray(dm,global_vector,global_array);</div><div>PetscViewerDestroy(&viewer);</div><div>}</div>
</div><div><br></div><div>1) I have 2 PDE's to solve. Still, I went ahead creating a single DM with dof=1 and creating two vectors using the DMCreateGlobalVector(). I want to write the file out periodically. Should I perform DMDAVecRestoreArray and DMDAVecGetArray every time is write out the global_vector? (I know that it is just indexing the pointers and there is no copying of values. But I am not sure)</div></div></blockquote><div><br></div><div>I don't think you need the Get/RestoreArray() calls here.</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div>2) I am writing out to HDF5 format. I see that the vecview is supposed to reorder the global_vector based on the DM. However, when I read the H5 files, I get an error on ViSIT and my output image becomes a 1D image rather than a 2D/3D. What might be the reason for this ?.<br></div><div>Error Msg : "In domain 0, your zonal variable "avtGhostZones" has 25600 values, but it should have 160. Some values were removed to ensure VisIt runs smoothly"</div><div>I was using a 160x160x1 DM</div></div></blockquote><div><br></div><div>I do not believe we support HDF5 <--> Visit/Paraview for DMDA. The VecView() is just writing out the vector as a linear array without mesh details. For</div><div>interfacing with the visualization, I think we use .vtu files. You should be able to get this effect using</div><div><br></div><div> VecViewFromOptions(global_vec, NULL, "-vec_view");</div><div><br></div><div>in your code, and then</div><div><br></div><div> -vec_view vtk:sol.vtu</div><div><br></div><div>on the command line.</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div>3) I tried using the "petsc_gen_xdmf.py" to generate the xdmf files for use in Paraview. Here the key ["viz/geometry"] is missing. The keys present in the output H5 file are just the two vectors I am writing and has no info about mesh. Isn't this supposed to come automatically since the vector is attached to the DM? How do I sort this out?</div></div></blockquote><div><br></div><div>This support is for unstructured grids, DMPlex and DMForest.</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div>4) Can I have multiple vectors attached to the DM by
DMCreateGlobalVector() even though I created the DMDA using dof=1.</div></div></blockquote><div><br></div><div>Yes.</div><div><br></div><div> Thanks,</div><div><br></div><div> Matt</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div dir="ltr"><div>thanks,<br></div><div>Abhishek<br></div></div>
</blockquote></div><br clear="all"><div><br></div>-- <br><div dir="ltr" class="gmail_signature"><div dir="ltr"><div><div dir="ltr"><div><div dir="ltr"><div>What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br>-- Norbert Wiener</div><div><br></div><div><a href="http://www.cse.buffalo.edu/~knepley/" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><br></div></div></div></div></div></div></div></div>