Dear Dr. Jonathan Guyer,


Thank you very much for your help. For the benefit of the community wishing to 
use 1D meshes in gmsh format for any of their external applications, I am 
bringing into your attention a nifty script available from the 'fluidity' CFD 
project.



For my 1D case, I have been able to use a very handy python script located 
here: https://github.com/FluidityProject/fluidity/blob/master/tools/interval.py



The usage is as follows:.  python interval.py --variable_dx=example.py 0.0 50.0 
my_variable_mesh, will create a mesh in gmsh format  called 
'my_variabel_mesh.msh' in the current directory.



The example.py can contain  any arbitrary numpy  function (with name val(X) ) 
for defining the mesh spacing.



For example something like

def val(X):

    return (5. / (X+1))



But clearly your script is much more amenable for a FiPy-centric workflow,  
since we don't have to worry about the intricacies of importing  a  1D .msh 
file.



Thank you for your help.





Krishna









-----Original Message-----
From: fipy-boun...@nist.gov [mailto:fipy-boun...@nist.gov] On Behalf Of Guyer, 
Jonathan E. Dr. (Fed)
Sent: 17 June 2016 15:28
To: FIPY <FIPY@nist.gov>
Subject: Re: casting implicit Boundary Conditions in FiPy



Gmsh does do 1D meshes, and I've got experimental code that imports them, but 
it's not ready to merge, yet.



In the meantime, this approach I've used for modeling semiconductor device 
contacts in 1D is probably better:



----



n_thickness = 1e-6 # m

p_thickness = 149e-6 # m

grid_resolution = 5e-8 # m



compression_factor = 0.8

compression_count = 10



compressed_dx = grid_resolution * 
compression_factor**numerix.arange(compression_count)

compressed_length = compressed_dx.sum()

compressed_dx = list(compressed_dx)



n_uncompressed_thickness = n_thickness - 2 * compressed_length n_dx = 
n_uncompressed_thickness / round(n_uncompressed_thickness / grid_resolution) 
n_dx = compressed_dx[::-1] + [n_dx] * int(n_uncompressed_thickness / n_dx) + 
compressed_dx



p_uncompressed_thickness = p_thickness - 2 * compressed_length p_dx = 
p_uncompressed_thickness / round(p_uncompressed_thickness / grid_resolution) 
p_dx = compressed_dx[::-1] + [p_dx] * int(p_uncompressed_thickness / p_dx + 
0.5) + compressed_dx



mesh = Grid1D(dx=n_dx + p_dx)



----







> On Jun 17, 2016, at 9:50 AM, Daniel Wheeler 
> <daniel.wheel...@gmail.com<mailto:daniel.wheel...@gmail.com>> wrote:

>

> On Thu, Jun 16, 2016 at 12:35 PM, Gopalakrishnan, Krishnakumar

> <k.gopalakrishna...@imperial.ac.uk<mailto:k.gopalakrishna...@imperial.ac.uk>> 
> wrote:

>> Thanks.

>>

>> Yes, this Is indeed only first order accurate.  I verified this by 
>> successively cutting my dx by half, running your code, and comparing against 
>> the Mathematica generated result. Each time dx is cut by half, the error is 
>> also proportionately halved.

>>

>> This code requires me to take unreasonably small dx values.  Since, neither 
>> the solution, nor the gradients are available at each implicit time-steps, I 
>> think that higher order schemes are probably ruled out.

>

> I think you're right. It's not easy.

>

>> I am thinking of using a variable mesh-sizing, let's say a log-spacing in 
>> 1D, keeping a very fine spacing (ultra-small dx) for the last cell near the 
>> boundary, and gradually taking bigger steps.  This is also physically 
>> consistent with my problem, wherein all the action takes place close to the 
>> boundary, and nothing much is happening at the middle or left edges.

>

> Maybe it's possible to make the edge cell almost infinitely small and

> have the rest of the cells evenly spaced.

>

>> This brings me to another issue.  I don't see a way to import a 1D .msh file 
>> generated by gmsh into FiPy.

>

> I'm not sure if Gmsh does 1D meshes, you can always use a 2D mesh as a

> 1D mesh of course.

>

>> Secondly,  I notice that there is an optimistic sounding  grid1DBuilder 
>> method.   I couldn't find any help on how to use this method. Is this method 
>> capable of creating the log-spaced mesh that I am considering ?

>

> Note that Grid1D takes an array for dx. So you can do

>

>>>> mesh = fp.Grid1D(dx=[0.5, 1.0, 2.0])

>>>> print(mesh.cellCenters)

>    [[ 0.25  1.   2.5 ]]

>

> As long as you can calculate the grid spacing then you don't need Gmsh.

>

> --

> Daniel Wheeler

>

> _______________________________________________

> fipy mailing list

> fipy@nist.gov<mailto:fipy@nist.gov>

> http://www.ctcms.nist.gov/fipy

>  [ NIST internal ONLY: https://email.nist.gov/mailman/listinfo/fipy ]





_______________________________________________

fipy mailing list

fipy@nist.gov<mailto:fipy@nist.gov>

http://www.ctcms.nist.gov/fipy

  [ NIST internal ONLY: https://email.nist.gov/mailman/listinfo/fipy ]
_______________________________________________
fipy mailing list
fipy@nist.gov
http://www.ctcms.nist.gov/fipy
  [ NIST internal ONLY: https://email.nist.gov/mailman/listinfo/fipy ]

Reply via email to