<html>
<head>
<meta http-equiv="Content-Type" content="text/html; charset=utf-8">
</head>
<body style="word-wrap: break-word; -webkit-nbsp-mode: space; -webkit-line-break: after-white-space;" class="">
Hi Matthew,
<div class=""><br class="">
</div>
<div class="">
<blockquote type="cite" class="">
<div dir="ltr" class="" style="font-family: Optima-Regular;">
<div class="gmail_extra">
<div class="gmail_quote">
<div class="">No. There are (k-1)^2 unknowns in the interior of the cell, so that we have</div>
<div class=""><br class="">
</div>
<div class="">  4 + 4 * (k-1) + (k-1)^2 = (k+1)^2</div>
</div>
</div>
</div>
</blockquote>
<div class=""><br class="">
</div>
Right, makes sense. With some insight from Dave, I’m getting what I expect for a simple 4-cell mesh distributed amongst 4 processors (length of 91 for the global vector, length of 25 for the local vectors, for 4th order polynomial and the GLL basis). Now I
 imagine what’s left is to figure out how these vectors are globally numbered.</div>
<div class=""><br class="">
</div>
<div class="">I usually do this by having a consistent way of numbering the dofs on the reference element (i.e. given a (u, v)-space loop through v, then u), and then getting the local->global map via finding coincident points in a kd-tree. My question is:
 in which order are the local vectors now defined? If I knew this, I could create my local->global map as before. From the way the section was created, I guess it might be something like loop vertices, edges, and then interior?</div>
<div class=""><br class="">
</div>
<div class="">Thanks for all your help,</div>
<div class="">Mike.</div>
<div class="">
<div apple-content-edited="true" class="">--<br class="">
Michael Afanasiev<br class="">
Ph.D. Candidate<br class="">
Computational Seismology<br class="">
Institut für Geophysik<br class="">
ETH Zürich<br class="">
<br class="">
Sonneggstrasse 5, NO H 39.2<br class="">
CH 8092 Zürich<br class="">
<a href="mailto:michael.afanasiev@erdw.ethz.ch" class="">michael.afanasiev@erdw.ethz.ch</a><br class="">
</div>
<br class="">
<div>
<blockquote type="cite" class="">
<div class="">On Oct 23, 2015, at 5:11 PM, Matthew Knepley <<a href="mailto:knepley@gmail.com" class="">knepley@gmail.com</a>> wrote:</div>
<br class="Apple-interchange-newline">
<div class="">
<div dir="ltr" style="font-family: Optima-Regular; font-size: 12px; font-style: normal; font-variant: normal; font-weight: normal; letter-spacing: normal; orphans: auto; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; widows: auto; word-spacing: 0px; -webkit-text-stroke-width: 0px;" class="">
<div class="gmail_extra">
<div class="gmail_quote">On Fri, Oct 23, 2015 at 9:03 AM, Afanasiev Michael<span class="Apple-converted-space"> </span><span dir="ltr" class=""><<a href="mailto:michael.afanasiev@erdw.ethz.ch" target="_blank" class="">michael.afanasiev@erdw.ethz.ch</a>></span><span class="Apple-converted-space"> </span>wrote:<br class="">
<blockquote class="gmail_quote" style="margin: 0px 0px 0px 0.8ex; border-left-width: 1px; border-left-color: rgb(204, 204, 204); border-left-style: solid; padding-left: 1ex;">
<div style="word-wrap: break-word;" class="">Hi Matthew,
<div class=""><br class="">
</div>
<div class="">
<blockquote type="cite" class="">
<div dir="ltr" style="font-family: Optima-Regular;" class="">
<div class="gmail_extra">
<div class="gmail_quote">
<div class="">Yes, however I have some questions. Starting out, I think the GLL points include the endpoints, so</div>
<div class="">that means for polynomial degree k that numDof[] should really have 4*1 dofs on each vertex,</div>
<div class="">4*(k-1) dofs on each edge, and 4*(k-1)^2 on each quad. It looks like below you have numDof[] for</div>
<div class="">a 1D mesh with DG element.</div>
</div>
</div>
</div>
</blockquote>
<div class=""><br class="">
</div>
For the GLL basis we have (k+1) points in each dimension, including the endpoints. So for example a 2D element with 4-order polynomials will have 25 locally numbered points per element. I should also mention that the velocity, displacement, and acceleration
 fields are vectors, with d=dim components at each integration point, whereas strain has (2^dim)-1 components. In what you’ve written above, does this then become (sum over the 4 fields):</div>
<div class=""><br class="">
</div>
<div class="">(sum(compPerField*dofPerField)) -> vertex</div>
<div class="">(sum((compPerField*dofPerField)*(k+1)) -> edge</div>
<div class="">(sum((compPerField*dofPerField)*(k+1)^2) -> quad</div>
</div>
</blockquote>
<div class=""><br class="">
</div>
<div class="">No. There are (k-1)^2 unknowns in the interior of the cell, so that we have</div>
<div class=""><br class="">
</div>
<div class="">  4 + 4 * (k-1) + (k-1)^2 = (k+1)^2</div>
<div class=""> </div>
<blockquote class="gmail_quote" style="margin: 0px 0px 0px 0.8ex; border-left-width: 1px; border-left-color: rgb(204, 204, 204); border-left-style: solid; padding-left: 1ex;">
<div style="word-wrap: break-word;" class="">
<div class="">
<blockquote type="cite" class="">
<div dir="ltr" style="font-family: Optima-Regular;" class="">
<div class="gmail_extra">
<div class="gmail_quote">
<div class="">With elements like these, it is common (I think) to eliminate the cell unknown explicitly, since the</div>
<div class="">system is dense and not connected to other cells. For this, you would retain the vertex and edge</div>
<div class="">unknowns, but not the cell unknowns. I have not tried to do this myself, so I do not know if there</div>
<div class="">are any pitfalls.</div>
</div>
</div>
</div>
</blockquote>
<br class="">
</div>
<div class="">I believe I don’t worry about this. Everything is solved in a matrix-free sense, on each element. The relevant physical quantities are calculated on each element, and then summed into the global degrees of freedom. The summed global dof are then
 brought back to the element level, where updates to the acceleration, velocity, and displacement are calculated via a Newmark time step. So no global system of equations is ever formed.</div>
</div>
</blockquote>
<div class=""><br class="">
</div>
<div class="">You accumulate into the global dofs, so you would need more storage if you do not do static condensation. It just good to know this,</div>
<div class="">but you do not have to do it.</div>
<div class=""> </div>
<blockquote class="gmail_quote" style="margin: 0px 0px 0px 0.8ex; border-left-width: 1px; border-left-color: rgb(204, 204, 204); border-left-style: solid; padding-left: 1ex;">
<div style="word-wrap: break-word;" class="">
<div class="">This being said, all the routines to a) integrate the element matrices, b) calculate the gll point locations, c) construct the local->global dof mapping, and d) do the element-wise matrix free time stepping are already written.  My problem is
 really that I do my mesh decomposition by hand (poorly), and I’d like to transfer this over to Petsc. Of course if I do this, I might as well use DM's LocalToGlobal vector routines to sum my field vectors across mesh partitions.</div>
</div>
</blockquote>
<div class=""><br class="">
</div>
<div class="">Yes.</div>
<div class=""> </div>
<blockquote class="gmail_quote" style="margin: 0px 0px 0px 0.8ex; border-left-width: 1px; border-left-color: rgb(204, 204, 204); border-left-style: solid; padding-left: 1ex;">
<div style="word-wrap: break-word;" class="">
<div class="">Perhaps a better question to ask would be in the form of a workflow:</div>
<div class=""><br class="">
</div>
<div class="">1)<span style="white-space: pre-wrap;" class=""> </span>Read in exodus mesh and partition (done)</div>
<div class="">2)<span style="white-space: pre-wrap;" class=""> </span>Set up local and global degrees of freedom on each mesh partition (done)</div>
</div>
</blockquote>
<div class=""><br class="">
</div>
<div class="">Here you just have to setup PetscSections that mirror your local and global layout. Then the LocalToGlobal and its inverse will work.</div>
<div class=""><br class="">
</div>
<div class="">   Matt</div>
<div class=""> </div>
<blockquote class="gmail_quote" style="margin: 0px 0px 0px 0.8ex; border-left-width: 1px; border-left-color: rgb(204, 204, 204); border-left-style: solid; padding-left: 1ex;">
<div style="word-wrap: break-word;" class="">
<div class="">3)<span style="white-space: pre-wrap;" class=""> </span>Integrate element matrices (done)</div>
<div class=""><br class="">
</div>
<div class="">for i 1, nTimeSteps</div>
<div class="">3)<span style="white-space: pre-wrap;" class=""> </span>Step fields forward on each partition (done)</div>
<div class="">4)<span style="white-space: pre-wrap;" class=""> </span>Sum to global degrees of freedom on each partition (done)</div>
<div class="">5)<span style="white-space: pre-wrap;" class=""> </span>Sum to global degrees of freedom across partitions (????)</div>
<div class="">6)<span style="white-space: pre-wrap;" class=""> </span>Retrieve summed global degrees of freedom across partitions (????)</div>
<div class="">7)<span style="white-space: pre-wrap;" class=""> </span>Continue</div>
<div class=""><br class="">
</div>
<div class="">So really what I hope Petsc will help with is just steps 5 and 6. I guess this involves, given a partitioned DMPlex object, which has been partitioned according to the geometry and topology defined in an exodus file, adding the degrees of freedom
 to each partition-local vector (via a DMPlexSection?). Then, ensuring that the dof added along the partition boundaries sum together properly when a LocalToGlobal vector operation is defined.</div>
<div class=""><br class="">
</div>
<div class="">Does this make sense? If so, is this something that DMPlex is designed for?<br class="">
<div class="">--<br class="">
Michael Afanasiev<br class="">
Ph.D. Candidate<br class="">
Computational Seismology<br class="">
Institut für Geophysik<br class="">
ETH Zürich<br class="">
<br class="">
Sonneggstrasse 5, NO H 39.2<br class="">
CH 8092 Zürich<br class="">
<a href="mailto:michael.afanasiev@erdw.ethz.ch" target="_blank" class="">michael.afanasiev@erdw.ethz.ch</a><br class="">
</div>
<br class="">
<div class="">
<blockquote type="cite" class="">
<div class="">On Oct 23, 2015, at 1:14 PM, Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank" class="">knepley@gmail.com</a>> wrote:</div>
<br class="">
<div class="">
<div dir="ltr" style="font-family: Optima-Regular; font-size: 12px; font-style: normal; font-variant: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px;" class="">
<div class="gmail_extra">
<div class="gmail_quote">On Fri, Oct 23, 2015 at 3:01 AM, Afanasiev Michael<span class=""> </span><span dir="ltr" class=""><<a href="mailto:michael.afanasiev@erdw.ethz.ch" target="_blank" class="">michael.afanasiev@erdw.ethz.ch</a>></span><span class=""> </span>wrote:<br class="">
<blockquote class="gmail_quote" style="margin: 0px 0px 0px 0.8ex; border-left-width: 1px; border-left-color: rgb(204, 204, 204); border-left-style: solid; padding-left: 1ex;">
<div style="word-wrap: break-word;" class="">Hi Matthew,
<div class=""><br class="">
</div>
<div class="">So I’m discretizing via a tensor product of Lagrange polynomials co-located at the Gauss-Lobatto-Legendre (GLL) points. The polynomial order might range from 4 to 10 or something like that. The current problem is solved on 2D hexes. </div>
<div class="">
<div class=""><br class="">
</div>
<div class="">I had found DMPlexCreateSection, and followed plex/ex1to get things set up. You can see my attempt below. Basically I’ve got 4 fields (displacement, velocity, acceleration, and strain) over each element. Here I’ve tried to setup a tensor product
 of 4th order (5 GLL points) Lagrange polynomials (line 11). This seemed to<span class=""> </span><i class="">somewhat</i> achieve what I wanted — I created a global vector and wrote it to a vtk file with VecView, and the numbering seemed to make sense. You
 can also see my attempt at defining a boundary condition (it looked like DMPlexCreateFromExodus labeled side sets as “Face Sets”, seems to have worked). </div>
<div class=""><br class="">
</div>
<div class="">Does this seem to be headed in the right direction?</div>
</div>
</div>
</blockquote>
<div class=""><br class="">
</div>
<div class="">Yes, however I have some questions. Starting out, I think the GLL points include the endpoints, so</div>
<div class="">that means for polynomial degree k that numDof[] should really have 4*1 dofs on each vertex,</div>
<div class="">4*(k-1) dofs on each edge, and 4*(k-1)^2 on each quad. It looks like below you have numDof[] for</div>
<div class="">a 1D mesh with DG element.</div>
<div class=""><br class="">
</div>
<div class="">The "Face Sets" is the right label to use for boundary conditions. This will eliminate those variables</div>
<div class="">from the global system, but they will be present in the local spaces.</div>
<div class=""><br class="">
</div>
<div class="">With elements like these, it is common (I think) to eliminate the cell unknown explicitly, since the</div>
<div class="">system is dense and not connected to other cells. For this, you would retain the vertex and edge</div>
<div class="">unknowns, but not the cell unknowns. I have not tried to do this myself, so I do not know if there</div>
<div class="">are any pitfalls.</div>
<div class=""><br class="">
</div>
<div class="">You can see an example of a similar implementation specifically for the kind of spectral elements</div>
<div class="">you are considering here: <a href="https://github.com/jedbrown/dohp" target="_blank" class="">https://github.com/jedbrown/dohp</a>. It would probably be useful to understand</div>
<div class="">what is done there as you implement.</div>
<div class=""><br class="">
</div>
<div class="">  Thanks,</div>
<div class=""><br class="">
</div>
<div class="">     Matt</div>
<div class=""> </div>
<blockquote class="gmail_quote" style="margin: 0px 0px 0px 0.8ex; border-left-width: 1px; border-left-color: rgb(204, 204, 204); border-left-style: solid; padding-left: 1ex;">
<div style="word-wrap: break-word;" class="">
<div class="">
<div class="">Cheers,</div>
<div class="">Mike.</div>
<div class=""><br class="">
</div>
<div class="">
<div class="">DM</div>
<div class="">mesh::createSection(const DM &dm)</div>
<div class="">{</div>
<div class=""><br class="">
</div>
<div class="">01        // displacement, velocity, acceleration, strain</div>
<div class="">02        IS bcPointIs[1];</div>
<div class="">03        PetscInt numBc = 1;</div>
<div class="">04        PetscInt bcField[1];</div>
<div class="">05        PetscInt numFields = 4;</div>
<div class="">06        PetscInt dim; DMGetDimension(dm, &dim);</div>
<div class="">07        PetscInt numComp[numFields];</div>
<div class="">08        for (auto i=0; i<numFields; i++) {numComp[i] = dim;}</div>
<div class="">09        PetscInt numDof[numFields*(dim+1)];</div>
<div class="">10        for (auto i=0; i<numFields*(dim+1); i++) {numDof[i] = 0;}</div>
<div class="">11        for (auto i=0; i<numFields; i++) {numDof[i*(dim+1)+dim] = 5;}</div>
<div class="">12        bcField[0] = 0;</div>
<div class="">13        PetscSection section;</div>
<div class="">14        DMPlexGetStratumIS(dm, "Face Sets", 1, &bcPointIs[0]);</div>
<div class="">15        DMPlexCreateSection(dm, dim, numFields, numComp, numDof, numBc, bcField,</div>
<div class="">16                                                NULL, bcPointIs, NULL, &section);</div>
<div class="">17        ISDestroy(&bcPointIs[0]);</div>
<div class="">18        PetscSectionSetFieldName(section, 0, "u");</div>
<div class="">19        PetscSectionSetFieldName(section, 1, "v");</div>
<div class="">20        PetscSectionSetFieldName(section, 2, "a");</div>
<div class="">21        PetscSectionSetFieldName(section, 3, "e");</div>
<div class="">22        DMSetDefaultSection(dm, section);</div>
<div class="">23        return dm;</div>
<div class="">}</div>
</div>
<div class=""><span class=""><br class="">
<div class="">--<br class="">
Michael Afanasiev<br class="">
Ph.D. Candidate<br class="">
Computational Seismology<br class="">
Institut für Geophysik<br class="">
ETH Zürich<br class="">
<br class="">
Sonneggstrasse 5, NO H 39.2<br class="">
CH 8092 Zürich<br class="">
<a href="mailto:michael.afanasiev@erdw.ethz.ch" target="_blank" class="">michael.afanasiev@erdw.ethz.ch</a><br class="">
</div>
<br class="">
</span>
<div class="">
<div class="">
<div class="">
<blockquote type="cite" class="">
<div class="">On Oct 22, 2015, at 2:32 AM, Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank" class="">knepley@gmail.com</a>> wrote:</div>
<br class="">
<div class="">
<div dir="ltr" style="font-family: Optima-Regular; font-size: 12px; font-style: normal; font-variant: normal; font-weight: normal; letter-spacing: normal; text-align: start; text-indent: 0px; text-transform: none; white-space: normal; word-spacing: 0px;" class="">
<div class="gmail_extra">
<div class="gmail_quote">On Wed, Oct 21, 2015 at 3:07 PM, Dave May<span class=""> </span><span dir="ltr" class=""><<a href="mailto:dave.may@erdw.ethz.ch" target="_blank" class="">dave.may@erdw.ethz.ch</a>></span><span class=""> </span>wrote:<br class="">
<blockquote class="gmail_quote" style="margin: 0px 0px 0px 0.8ex; border-left-width: 1px; border-left-color: rgb(204, 204, 204); border-left-style: solid; padding-left: 1ex;">
<div dir="ltr" class="">Hey Mike,<br class="">
<br class="">
<br class="">
<div class="gmail_extra"><br class="">
<div class="gmail_quote">On 21 October 2015 at 18:01, Afanasiev Michael<span class=""> </span><span dir="ltr" class=""><<a href="mailto:michael.afanasiev@erdw.ethz.ch" target="_blank" class="">michael.afanasiev@erdw.ethz.ch</a>></span><span class=""> </span>wrote:<br class="">
<blockquote class="gmail_quote" style="margin: 0px 0px 0px 0.8ex; border-left-width: 1px; border-left-style: solid; border-left-color: rgb(204, 204, 204); padding-left: 1ex;">
<div style="word-wrap: break-word;" class="">Hey Dave,
<div class=""><br class="">
</div>
<div class="">So I’ve got a couple of days where there’s nothing pressing to work on… was thinking of ripping out the parallel routines (ugly) in my wave propagation code and replacing them with Petsc DM routines. I can read in my exodusii mesh with DMPLEX,
 and partition it with chaco, which gives me a nicely partitioned DM. This takes me like 5 lines of code.  That’s amazing.</div>
<div class=""><br class="">
</div>
<div class="">But here I’m stuck, and am having a whale of a time with the documentation. All I<span class=""> </span><i class="">think</i><span class=""> </span>I need is a way to modify the exodus-created DM, and add to it the degrees of freedom that are
 introduced by my quadrature rule. This would be really neat. I can just treat each sub-domain as its own mesh, with its own global numbering. Then whenever necessary I can scatter stuff the the<span class=""> </span><i class="">real</i> global degrees of freedom
 with something like VecLocToGlob. Most of the things I like about the code could stay the same (element-wise, matrix-free nature), just these parallel broadcasts would be infinitely nicer.</div>
<div class=""><br class="">
</div>
</div>
</blockquote>
<div class=""><br class="">
</div>
<div class="">First off - I don't use DMPLEX.</div>
</div>
</div>
</div>
</blockquote>
<div class=""><br class="">
</div>
<div class="">Dave is refreshingly candid about his shortcomings ;)</div>
<div class=""> </div>
<blockquote class="gmail_quote" style="margin: 0px 0px 0px 0.8ex; border-left-width: 1px; border-left-color: rgb(204, 204, 204); border-left-style: solid; padding-left: 1ex;">
<div dir="ltr" class="">
<div class="gmail_extra">
<div class="gmail_quote">
<div class=""> </div>
</div>
</div>
</div>
</blockquote>
<blockquote class="gmail_quote" style="margin: 0px 0px 0px 0.8ex; border-left-width: 1px; border-left-color: rgb(204, 204, 204); border-left-style: solid; padding-left: 1ex;">
<div dir="ltr" class="">
<div class="gmail_extra">
<div class="gmail_quote">
<div class=""></div>
<blockquote class="gmail_quote" style="margin: 0px 0px 0px 0.8ex; border-left-width: 1px; border-left-style: solid; border-left-color: rgb(204, 204, 204); padding-left: 1ex;">
<div style="word-wrap: break-word;" class="">
<div class=""></div>
<div class="">But I just can’t figure out how to set this up. The main problem really boils down to: what’s the best way to add my quadrature points to an already-created DM, which was constructed with an exodus file? I guess I could do this after the file
 is read, but before the partitioning. In this case though, what’s stopping the partitioner from cutting an element in half? It seems like it would be a lot cleaner to do this post-partitioning.</div>
<div class=""><br class="">
</div>
</div>
</blockquote>
<div class=""><br class="">
</div>
<div class="">Presumably what is read from exodus is just the vertices of the hexes, and what you want to do is define the function space (given by your GLL locations) on top of element geometry read in. Is that what you are asking about?</div>
</div>
</div>
</div>
</blockquote>
<div class=""><br class="">
</div>
<div class="">So Dave is right. We read in topology and geometry from ExodusII. Then you define a function space on top. How</div>
<div class="">exactly are you discretizing? In order to create vectors, do local to global, etc. Petsc really only need to know the</div>
<div class="">amount of data associated with each mesh piece. You can define this with PetscSection. If you give me an idea</div>
<div class="">what you want I can help you write the code easily I think.</div>
<div class=""><br class="">
</div>
<div class="">  Thanks,</div>
<div class=""><br class="">
</div>
<div class="">     Matt</div>
<div class=""> </div>
<blockquote class="gmail_quote" style="margin: 0px 0px 0px 0.8ex; border-left-width: 1px; border-left-color: rgb(204, 204, 204); border-left-style: solid; padding-left: 1ex;">
<div dir="ltr" class="">
<div class="gmail_extra">
<div class="gmail_quote">
<blockquote class="gmail_quote" style="margin: 0px 0px 0px 0.8ex; border-left-width: 1px; border-left-style: solid; border-left-color: rgb(204, 204, 204); padding-left: 1ex;">
<div style="word-wrap: break-word;" class="">
<div class=""></div>
<div class="">Any hints here?</div>
</div>
</blockquote>
<div class=""><br class="">
</div>
<div class="">Actually I have no experience with this object.<br class="">
I would just send an email to<span class=""> </span><br class="">
</div>
<div class=""> <span class=""> </span><a href="mailto:petsc-users@mcs.anl.gov" target="_blank" class="">petsc-users@mcs.anl.gov</a><br class="">
</div>
<div class="">asking for help.<br class="">
</div>
<div class=""> <br class="">
</div>
<div class="">The developer of DMPLEX (Matt Knepley) will definitely answer within in 1 day.<br class="">
</div>
<div class=""><br class="">
</div>
<div class="">Cheers,<br class="">
</div>
<div class="">  Dave<br class="">
</div>
<div class=""><br class="">
</div>
<blockquote class="gmail_quote" style="margin: 0px 0px 0px 0.8ex; border-left-width: 1px; border-left-style: solid; border-left-color: rgb(204, 204, 204); padding-left: 1ex;">
<div style="word-wrap: break-word;" class="">
<div class=""><br class="">
</div>
<div class="">Best,</div>
<div class="">Mike.<span class=""><font color="#888888" class=""><br class="">
<div class="">--<br class="">
Michael Afanasiev<br class="">
Ph.D. Candidate<br class="">
Computational Seismology<br class="">
Institut für Geophysik<br class="">
ETH Zürich<br class="">
<br class="">
Sonneggstrasse 5, NO H 39.2<br class="">
CH 8092 Zürich<br class="">
<a href="mailto:michael.afanasiev@erdw.ethz.ch" target="_blank" class="">michael.afanasiev@erdw.ethz.ch</a><br class="">
</div>
<br class="">
</font></span></div>
</div>
</blockquote>
</div>
<br class="">
</div>
</div>
</blockquote>
</div>
<br class="">
<br clear="all" class="">
<span class="HOEnZb"><font color="#888888" class="">
<div class=""><br class="">
</div>
--<span class=""> </span><br class="">
<div class="">What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br class="">
-- Norbert Wiener</div>
</font></span></div>
<span class="HOEnZb"><font color="#888888" class=""></font></span></div>
<span class="HOEnZb"><font color="#888888" class=""></font></span></div>
<span class="HOEnZb"><font color="#888888" class=""></font></span></blockquote>
<span class="HOEnZb"><font color="#888888" class=""></font></span></div>
<span class="HOEnZb"><font color="#888888" class=""><br class="">
</font></span></div>
<span class="HOEnZb"><font color="#888888" class=""></font></span></div>
<span class="HOEnZb"><font color="#888888" class=""></font></span></div>
<span class="HOEnZb"><font color="#888888" class=""></font></span></div>
<span class="HOEnZb"><font color="#888888" class=""></font></span></div>
<span class="HOEnZb"><font color="#888888" class=""></font></span></blockquote>
<span class="HOEnZb"><font color="#888888" class=""></font></span></div>
<span class="HOEnZb"><font color="#888888" class=""><br class="">
<br clear="all" class="">
<div class=""><br class="">
</div>
--<span class=""> </span><br class="">
<div class="">What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br class="">
-- Norbert Wiener</div>
</font></span></div>
</div>
</div>
</blockquote>
</div>
<br class="">
</div>
</div>
</blockquote>
</div>
<br class="">
<br clear="all" class="">
<div class=""><br class="">
</div>
--<span class="Apple-converted-space"> </span><br class="">
<div class="gmail_signature">What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br class="">
-- Norbert Wiener</div>
</div>
</div>
</div>
</blockquote>
</div>
<br class="">
</div>
</body>
</html>