trixi-framework / HOHQMesh

High Order Hex-Quad Mesh (HOHQMesh) package to automatically generate all-quadrilateral meshes with high order boundary information.

Home Page:https://trixi-framework.github.io/HOHQMesh

Geek Repo:Geek Repo

Github PK Tool:Github PK Tool

Can HOHQMesh generate split-node/discontinuous inner boundary (i.e., just a line, not a closed curve)?

thehalfspace opened this issue · comments

I want to generate a 2D rectangular mesh with an arbitrary fault discontinuity (straight line) within the domain.
It looks as follows. For this toy case, I have:

Lx, Ly = 4,4
fault_start = (2,4)
fault_end = (2,2)
creep_end = (2,0)
               (fault_start)
(0,Ly) _____________________________(Lx,Ly)
      |               |             |
      |               |             |
      |               |             |
      |    (fault_end)/(creep_start)|
      |               -             |
      |               -             |
      |               -             |
      |               -             |
 (0,0)|______________________________ | (Lx, 0)
                 (creep_end)

I want to use this type of mesh with a manually-written spectral-element code. It has Linear Elasticity/Wave propagation in the domain, a friction ODE along the fault, and fixed velocity along the creep. So I would just want to get (x, y, connectivity index) across all these interfaces.

I want to use HOHQMesh because eventually I want to have unstructured meshes with curved faults.

Creating my control file gives a segmentation fault error when I define those inner boundaries. The control file looks like this:

	open("$(@__DIR__)/../out/rectangular-structured-boundary.control", "w") do io
		println(io, raw"""

		\begin{CONTROL_INPUT}
			\begin{RUN_PARAMETERS}
				mesh file name   = rectangular-structured-boundary.mesh
				plot file name   = rectangular-structured-boundary.tec
				stats file name  = none
				mesh file format = ISM-v2
				polynomial order = 4
				plot file format = skeleton
			\end{RUN_PARAMETERS}

			\begin{BACKGROUND_GRID}
				x0 = [0.0, 0.0, 0.0]
				dx = [1.0, 1.0, 0.0]
				N  = [4,4,1]
			\end{BACKGROUND_GRID}
		\end{CONTROL_INPUT}

		\begin{MODEL}
			\begin{OUTER_BOUNDARIES}
			\begin{END_POINTS_LINE}
				name = bottom
				xStart = [0.0,0.0,0.0]
				xEnd   = [4.0,0.0,0.0]
			\end{END_POINTS_LINE}
			\begin{END_POINTS_LINE}
				name = right
				xStart = [4.0,0.0,0.0]
				xEnd   = [4.0,4.0,0.0]
			\end{END_POINTS_LINE}
			\begin{END_POINTS_LINE}
				name = top
				xStart = [4.0,4.0,0.0]
				xEnd   = [0.0,4.0,0.0]
			\end{END_POINTS_LINE}
			\begin{END_POINTS_LINE}
				name = left
				xStart = [0.0,4.0,0.0]
				xEnd   = [0.0,0.0,0.0]
			\end{END_POINTS_LINE}
			\end{OUTER_BOUNDARIES}

			\begin{INNER_BOUNDARIES}
			\begin{END_POINTS_LINE}
				name = Fault
				xStart = [ 2.0, 4.0, 0.0]
				xEnd   = [ 2.0, 2.0, 0.0]
			\end{END_POINTS_LINE}
			\begin{END_POINTS_LINE}
				name = Creep
				xStart = [ 2.0, 0.0, 0.0]
				xEnd   = [ 2.0, 2.0, 0.0]
			\end{END_POINTS_LINE}	
			
			\end{INNER_BOUNDARIES}
		\end{MODEL}
		
		\end{FILE}
		""")
	end
end
using HOHQMesh
create_mesh_control()
control_file = "$(@__DIR__)/../out/rectangular-structured-boundary.control"
generate_mesh(control_file);

The error is:

Program received signal SIGSEGV: Segmentation fault - invalid memory reference.

Backtrace for this error:
#0  0x10d92263c
#1  0x10d9219d3
#2  0x7ff80cd52dfc
#3  0x104f3f54e
#4  0x104f48461
#5  0x104f30033
#6  0x104f30112
#7  0x104f30cd8
#8  0x104f310af
#9  0x104f3098d
#10  0x104f310af
#11  0x104f3098d
#12  0x104f310af
#13  0x104f3098d
#14  0x104f31221
#15  0x104f4c317
#16  0x104f4b808
#17  0x104fa21d0
ERROR: LoadError: failed process: Process(setenv(`/Users/prithvi/.julia/artifacts/dcfa115e4439dcb2e24beadcb54361c4e46318e0/bin/HOHQMesh -f /var/folders/17/tlgpbfwx3p56q2gxlzt5gx2c0000gn/T/jl_o8auM1`,["_CE_M=", "LSCOLORS=Gxfxcxdxbxegedabagacad", "PATH=/Users/prithvi/.julia/artifacts/dcfa115e4439dcb2e24beadcb54361c4e46318e0/bin:/Users/prithvi/opt/chromedriver:/Users/prithvi/opt/anaconda3/bin:/opt/local/bin:/opt/local/sbin:/opt/homebrew/bin:/opt/homebrew/sbin:/usr/local/bin:/usr/bin:/bin:/usr/sbin:/sbin:/Users/prithvi/opt/chromedriver:/Users/prithvi/opt/anaconda3/bin:/Users/prithvi/opt/anaconda3/condabin:/opt/local/bin:/opt/local/sbin:/opt/homebrew/bin:/opt/homebrew/sbin", "HDF5_USE_FILE_LOCKING=FALSE", "MANPATH=/opt/local/share/man:/opt/homebrew/share/man:/usr/share/man:/usr/local/share/man:/opt/local/share/man:/opt/homebrew/share/man:", "USER=prithvi", "UCX_ERROR_SIGNALS=SIGILL,SIGBUS,SIGFPE", "LESS=-R", "CONDA_PROMPT_MODIFIER=(base) ", "VSCODE_GIT_ASKPASS_NODE=/Applications/Visual Studio Code.app/Contents/Frameworks/Code Helper.app/Contents/MacOS/Code Helper"  …  "UCX_MEMTYPE_CACHE=no", "TERM_PROGRAM_VERSION=1.68.1", "ZSH=/Users/prithvi/.oh-my-zsh", "COMMAND_MODE=unix2003", "PWD=/Users/prithvi/Library/Mobile Documents/com~apple~CloudDocs/spear-parallel", "TERM_PROGRAM=vscode", "VSCODE_GIT_ASKPASS_MAIN=/Applications/Visual Studio Code.app/Contents/Resources/app/extensions/git/dist/askpass-main.js", "ORIGINAL_XDG_CURRENT_DESKTOP=undefined", "CONDA_PREFIX=/Users/prithvi/opt/anaconda3", "OPENBLAS_MAIN_FREE=1"]), ProcessSignaled(11)) [0]

Stacktrace:
  [1] pipeline_error
    @ ./process.jl:540 [inlined]
  [2] read(cmd::Cmd)
    @ Base ./process.jl:427
  [3] read(cmd::Cmd, #unused#::Type{String})
    @ Base ./process.jl:436
  [4] readchomp
    @ ./io.jl:926 [inlined]
  [5] (::HOHQMesh.var"#17#18"{Bool, String, String, String, String})(tmppath::String, tmpio::IOStream)
    @ HOHQMesh ~/.julia/packages/HOHQMesh/wrBSY/src/HOHQMesh.jl:226
  [6] mktemp(fn::HOHQMesh.var"#17#18"{Bool, String, String, String, String}, parent::String)
    @ Base.Filesystem ./file.jl:722
  [7] mktemp
    @ ./file.jl:720 [inlined]
  [8] generate_mesh(control_file::String; output_directory::String, mesh_filename::Nothing, plot_filename::Nothing, stats_filename::Nothing, verbose::Bool)
    @ HOHQMesh ~/.julia/packages/HOHQMesh/wrBSY/src/HOHQMesh.jl:206
  [9] generate_mesh(control_file::String)
    @ HOHQMesh ~/.julia/packages/HOHQMesh/wrBSY/src/HOHQMesh.jl:174
 [10] top-level scope
    @ ~/Library/Mobile Documents/com~apple~CloudDocs/spear-parallel/src/mesh.jl:105
 [11] include(fname::String)
    @ Base.MainInclude ./client.jl:451
 [12] top-level scope
    @ REPL[29]:1
in expression starting at /Users/prithvi/Library/Mobile Documents/com~apple~CloudDocs/spear-parallel/src/mesh.jl:105

Thanks for the great package!

To answer your question, HOHQMesh does not have the capability to do what you want it to do at this time.

As noted in the documentation, all curves (interior and exterior) must be closed. As implemented the curves must not intersect, and was written with that implicit assumption. Testing for that condition was started and needs to be completed. (Answering the question "do two curves cross" requires care.) Also, the MODEL in the control file above is malformed. As noted at https://trixi-framework.github.io/HOHQMesh/the-model/, interior boundaries must be defined within a CHAIN, but that should have been caught and noted to you at execution time, so a bug fix is needed for that.

What you are trying to do is use an INTERFACE_BOUNDARY, rather than an INNER_BOUNDARY which is actually implemented in the code, but not documented/advertised because it has bugs (segmentation faults) and the algorithm doesn't generate great elements near the interface. Two things will have to be done to clear this issue: 1) Interface boundaries would have to be allowed to cross the outer boundary and 2) The segmentation fault coming from node deletion near the interface needs to be tracked down. The second was raised in an earlier Issue.

As an aside, Andrew has just completed an interactive package https://github.com/trixi-framework/HOHQMesh.jl that will construct a control file, and hence will avoid formatting errors. It does not disallow curve intersections, I don't think, so that's probably another bug to process on that side.

I think I've outlined above (for posterity) what is needed to resolve the issues.

Thanks for reporting the issue.

Thank you for the detailed response. I will eventually try using an interface boundary. For now, I am just solving on one-half of the domain since it is symmetric across the interface boundary, then in this case my right side boundary becomes an ode. It is working great so far!

As a followup, PR #50 has been made to eliminate the need to put all interior boundary curves in a CHAIN, so that leaving them out as in the control file above won't be a problem. That PR does not address the more complicated issues of defining cracks that are not closed, an possibly crossing another boundary.

As an update to this request @thehalfspace . I have made a control file that can create a mesh with this kind of "spike". As David mentioned above, HOHQMesh is not currently designed to handle this kind of geometry. So what the control file does is define a very narrow boundary spike (on the order of machine epsilon). Because of this, the mesh quality is very sensitive to the size of the background grid, overall domain size, etc. But the control file below generates a mesh in the included image.

Screenshot 2022-10-14 at 19 52 22

\begin{MODEL}
   \begin{OUTER_BOUNDARY}
      \begin{END_POINTS_LINE}
         name = right
         xEnd = [10.0,5.0,0.0]
         xStart = [10.0,-5.0,0.0]
      \end{END_POINTS_LINE}
      \begin{END_POINTS_LINE}
         name = top
         xEnd = [-10.0,5.0,0.0]
         xStart = [10.0,5.0,0.0]
      \end{END_POINTS_LINE}
      \begin{END_POINTS_LINE}
         name = left
         xEnd = [-10.0,-5.0,0.0]
         xStart = [-10.0,5.0,0.0]
      \end{END_POINTS_LINE}
      \begin{END_POINTS_LINE}
         name = bottom
         xEnd = [0.0,-5.0,0.0]
         xStart = [-10.0,-5.0,0.0]
      \end{END_POINTS_LINE}
      \begin{END_POINTS_LINE}
         name = spike
         xEnd = [0.0,-1.0,0.0]
         xStart = [0.0,-5.0,0.0]
      \end{END_POINTS_LINE}
      \begin{END_POINTS_LINE}
         name = spike
         xEnd = [1.0e-13,-5.0,0.0]
         xStart = [0.0,-1.0,0.0]
      \end{END_POINTS_LINE}
      \begin{END_POINTS_LINE}
         name = bottom
         xEnd = [10.0,-5.0,0.0]
         xStart = [1.0e-13,-5.0,0.0]
      \end{END_POINTS_LINE}
   \end{OUTER_BOUNDARY}
\end{MODEL}
\begin{CONTROL_INPUT}
   \begin{SPRING_SMOOTHER}
      smoothing type = LinearAndCrossbarSpring
      smoothing = ON
   \end{SPRING_SMOOTHER}
   \begin{BACKGROUND_GRID}
      background grid size = [2.0,2.0,0.0]
   \end{BACKGROUND_GRID}
   \begin{RUN_PARAMETERS}
      mesh file name = out/spike.mesh
      plot file format = skeleton
      plot file name = out/spike.tec
      stats file name = out/spike.txt
      mesh file format = ISM-V2
      polynomial order = 1
   \end{RUN_PARAMETERS}
\end{CONTROL_INPUT}
\end{FILE}

To follow up, in response to the issue of generating a symmetric mesh, your example can now be done by defining half of the geometry and reflecting it along the vertical line in the center. There is tutorial now at https://trixi-framework.github.io/HOHQMesh.jl/stable/tutorials/symmetric_mesh/. The trick is to define two lines in the vertical direction, one that goes from creep_end to creep_start and the other from fault_end to fault_start as the right boundary. Assign the first one the name ":symmetry" and call the other whatever you like, say "fault". The mesh will be reflected, and the extra nodes removed in the bottom half. On the top part the nodes will be doubled with two boundaries present called "fault" and "fault_R". You can then apply whatever BCs you want to those curves. Unfortunately nothing fancier, e.g. multiple faults except on a single line, can be accomplished this way. But it will solve this problem.

Thank you @DavidAKopriva . This will be really helpful for a lot of geophysics problems. I'll try the symmetric mirror mesh and let you know how it goes.

Also, feel free to close this issue since my problem was resolved.