[cig-commits] r6719 - short/3D/PyLith/trunk/libsrc/faults

brad at geodynamics.org brad at geodynamics.org
Fri Apr 27 14:26:34 PDT 2007


Author: brad
Date: 2007-04-27 14:26:34 -0700 (Fri, 27 Apr 2007)
New Revision: 6719

Modified:
   short/3D/PyLith/trunk/libsrc/faults/FaultCohesive.cc
   short/3D/PyLith/trunk/libsrc/faults/FaultCohesive.hh
   short/3D/PyLith/trunk/libsrc/faults/FaultCohesiveKin.cc
   short/3D/PyLith/trunk/libsrc/faults/FaultCohesiveKin.hh
Log:
Worked on initialization of fault (computing orientation of vertices).

Modified: short/3D/PyLith/trunk/libsrc/faults/FaultCohesive.cc
===================================================================
--- short/3D/PyLith/trunk/libsrc/faults/FaultCohesive.cc	2007-04-27 20:49:44 UTC (rev 6718)
+++ short/3D/PyLith/trunk/libsrc/faults/FaultCohesive.cc	2007-04-27 21:26:34 UTC (rev 6719)
@@ -120,6 +120,7 @@
   const int jacobianSize = 6;
   assert(numLocs*jacobianSize == jacobian.size());
   assert(numLocs == jacobianDet.size());
+  assert(3 == upDir.size());
 
   const int cellDim = 2;
   const int spaceDim = 3;

Modified: short/3D/PyLith/trunk/libsrc/faults/FaultCohesive.hh
===================================================================
--- short/3D/PyLith/trunk/libsrc/faults/FaultCohesive.hh	2007-04-27 20:49:44 UTC (rev 6718)
+++ short/3D/PyLith/trunk/libsrc/faults/FaultCohesive.hh	2007-04-27 21:26:34 UTC (rev 6719)
@@ -197,8 +197,8 @@
   /// Not implemented
   const FaultCohesive& operator=(const FaultCohesive& m);
 
-  // PRIVATE MEMBERS ////////////////////////////////////////////////////
-private :
+// PROTECTED MEMBERS ////////////////////////////////////////////////////
+protected :
 
   ALE::Obj<ALE::Mesh>* _faultMesh;
 

Modified: short/3D/PyLith/trunk/libsrc/faults/FaultCohesiveKin.cc
===================================================================
--- short/3D/PyLith/trunk/libsrc/faults/FaultCohesiveKin.cc	2007-04-27 20:49:44 UTC (rev 6718)
+++ short/3D/PyLith/trunk/libsrc/faults/FaultCohesiveKin.cc	2007-04-27 21:26:34 UTC (rev 6719)
@@ -16,6 +16,7 @@
 
 #include "EqKinSrc.hh" // USES EqKinSrc
 
+#include "pylith/feassemble/Quadrature.hh" // USES Quadrature
 #include "pylith/utils/array.hh" // USES double_array
 
 #include <assert.h> // USES assert()
@@ -61,18 +62,51 @@
 					     const double_array& upDir)
 { // initialize
   assert(0 != _quadrature);
+  assert(0 != _faultMesh);
+  assert(!_faultMesh->isNull());
   
   if (3 != upDir.size())
     throw std::runtime_error("Up direction for fault orientation must be "
 			     "a vector with 3 components.");
+
+  // Allocate section for orientation at all fault vertices
+  ALE::Obj<real_section_type> orientation = 
+    new real_section_type((*_faultMesh)->comm(), (*_faultMesh)->debug());
+  assert(!orientation.isNull());
+  const int orientationSize = _orientationSize();
+  orientation->setFiberDimension((*_faultMesh)->depthStratum(0), 
+				 orientationSize);
+  (*_faultMesh)->allocate(orientation);
   
-  // Allocate section for orientation information (all vertices)
+  // Get section containing coordinates of vertices
+  const ALE::Obj<real_section_type>& coordinates = 
+    mesh->getRealSection("coordinates");
+  assert(!coordinates.isNull());
 
   // Loop over cells
-  //   Compute cell geometry at vertices
-  //   Compute weighted orientation of face at vertices (using geometry info)
-  //   Update weighted orientations (wt is |J|)
+  const ALE::Obj<Mesh::label_sequence>& cells = 
+    (*_faultMesh)->heightStratum(0);
+  const Mesh::label_sequence::iterator cBegin = cells->begin();
+  const Mesh::label_sequence::iterator cEnd = cells->end();
+  double_array cellOrientation(_quadrature->numBasis()*orientationSize);
+  for (Mesh::label_sequence::iterator c_iter=cBegin;
+       c_iter != cEnd;
+       ++c_iter) {
+    // Compute cell geometry at vertices
+    _quadrature->computeGeometryVert(*_faultMesh, coordinates, *c_iter);
 
+    const double_array& jacobian = _quadrature->jacobianVert();
+    const double_array& jacobianDet = _quadrature->jacobianDetVert();
+
+    // Compute weighted orientation of face at vertices (using geometry info)
+    // NEED TO CALL APPROPRIATE ORIENT FOR CELL DIM
+    //_orient1D(&cellOrientation, jacobian, jacobianDet, upDir, numVertices);
+
+    // Update weighted orientations
+    // Loop over vertices in cell
+    //   Update section with orientation for each vertex
+  } // for
+
   // Assemble orientation information
 
   // Loop over vertices
@@ -132,5 +166,29 @@
     disp->updatePoint(*v_iter, slip->restrictPoint(*v_iter));
 } // setField
 
+// ----------------------------------------------------------------------
+// Get size (fiber dimension) of orientation information.
+int 
+pylith::faults::FaultCohesiveKin::_orientationSize(void) const
+{ // _orientationSize
+  assert(0 != _quadrature);
 
+  int size = 0;
+  switch (_quadrature->cellDim()) {
+  case 1 :
+    size = 1;
+    break;
+  case 2 :
+    size = 4;
+    break;
+  case 3 :
+    size = 9;
+    break;
+  default :
+    assert(0);
+  } // switch
+  return size;
+} // _orientationSize
+
+
 // End of file 

Modified: short/3D/PyLith/trunk/libsrc/faults/FaultCohesiveKin.hh
===================================================================
--- short/3D/PyLith/trunk/libsrc/faults/FaultCohesiveKin.hh	2007-04-27 20:49:44 UTC (rev 6718)
+++ short/3D/PyLith/trunk/libsrc/faults/FaultCohesiveKin.hh	2007-04-27 21:26:34 UTC (rev 6719)
@@ -132,6 +132,15 @@
   /// Not implemented
   const FaultCohesiveKin& operator=(const FaultCohesiveKin& m);
 
+  // PRIVATE METHODS ////////////////////////////////////////////////////
+private :
+
+  /** Get size (fiber dimension) of orientation information.
+   *
+   * @returns Size of orientation information.
+   */
+  int _orientationSize(void) const;
+
   // PRIVATE MEMBERS ////////////////////////////////////////////////////
 private :
 



More information about the cig-commits mailing list