• Main Page
  • Related Pages
  • Namespaces
  • Classes
  • Files
  • File List
  • File Members

Albany_Layouts.cpp

Go to the documentation of this file.
00001 //*****************************************************************//
00002 //    Albany 2.0:  Copyright 2012 Sandia Corporation               //
00003 //    This Software is released under the BSD license detailed     //
00004 //    in the file "license.txt" in the top-level Albany directory  //
00005 //*****************************************************************//
00006 #include "Albany_Layouts.hpp"
00007 #include "Albany_DataTypes.hpp"
00008 
00009 #include "PHAL_FactoryTraits.hpp"
00010 
00011 /*********************** Helper Functions*********************************/
00012 
00013 Albany::Layouts::Layouts (int worksetSize, int  numVertices,
00014                           int numNodes, int numQPts, int numDim, int vecDim, int numFace)
00015 // numDim is the number of spatial dimensions
00016 // vecDim is the length of a vector quantity
00017 // -- For many problems, numDim is used for both since there are 
00018 // typically numDim displacements and velocities.
00019 {
00020   using Teuchos::rcp;
00021   using PHX::MDALayout;
00022 
00023   // 
00024   if (vecDim==-1) vecDim = numDim;
00025   if (vecDim == numDim) vectorAndGradientLayoutsAreEquivalent = true;
00026   else                  vectorAndGradientLayoutsAreEquivalent = false;
00027   
00028   // Solution Fields
00029   node_scalar = rcp(new MDALayout<Cell,Node>(worksetSize,numNodes));
00030   qp_scalar   = rcp(new MDALayout<Cell,QuadPoint>(worksetSize,numQPts));
00031   cell_scalar = rcp(new MDALayout<Cell,QuadPoint>(worksetSize,1));
00032   cell_scalar2 = rcp(new MDALayout<Cell>(worksetSize));
00033   face_scalar = rcp(new MDALayout<Cell,Face>(worksetSize,numFace));
00034 
00035   node_vector = rcp(new MDALayout<Cell,Node,Dim>(worksetSize,numNodes,vecDim));
00036   qp_vector   = rcp(new MDALayout<Cell,QuadPoint,Dim>(worksetSize,numQPts,vecDim));
00037   cell_vector   = rcp(new MDALayout<Cell,Dim>(worksetSize,vecDim));
00038   face_vector   = rcp(new MDALayout<Cell,Face,Dim>(worksetSize,numFace,vecDim));
00039 
00040   node_gradient = rcp(new MDALayout<Cell,Node,Dim>(worksetSize,numNodes,numDim));
00041   qp_gradient   = rcp(new MDALayout<Cell,QuadPoint,Dim>(worksetSize,numQPts,numDim));
00042   cell_gradient   = rcp(new MDALayout<Cell,Dim>(worksetSize,numDim));
00043   face_gradient   = rcp(new MDALayout<Cell,Face,Dim>(worksetSize,numFace,numDim));
00044 
00045   node_tensor = rcp(new MDALayout<Cell,Node,Dim,Dim>(worksetSize,numNodes,numDim,numDim));
00046   qp_tensor   = rcp(new MDALayout<Cell,QuadPoint,Dim,Dim>(worksetSize,numQPts,numDim,numDim));
00047   cell_tensor   = rcp(new MDALayout<Cell,Dim,Dim>(worksetSize,numDim,numDim));
00048   face_tensor   = rcp(new MDALayout<Cell,Face,Dim,Dim>(worksetSize,numFace,numDim,numDim));
00049 
00050   node_vecgradient = rcp(new MDALayout<Cell,Node,Dim,Dim>(worksetSize,numNodes,vecDim,numDim));
00051   qp_vecgradient   = rcp(new MDALayout<Cell,QuadPoint,Dim,Dim>(worksetSize,numQPts,vecDim,numDim));
00052   cell_vecgradient = rcp(new MDALayout<Cell,Dim,Dim>(worksetSize,vecDim,numDim));
00053   face_vecgradient = rcp(new MDALayout<Cell,Face,Dim,Dim>(worksetSize,numFace,vecDim,numDim));
00054 
00055   node_tensor4   = rcp(new MDALayout<Cell,Node,Dim,Dim,Dim,Dim>(worksetSize,numNodes,numDim,numDim,numDim,numDim));
00056   qp_tensor4     = rcp(new MDALayout<Cell,QuadPoint,Dim,Dim,Dim,Dim>(worksetSize,numQPts,numDim,numDim,numDim,numDim));
00057   cell_tensor4   = rcp(new MDALayout<Cell,Dim,Dim,Dim,Dim>(worksetSize,numDim,numDim,numDim,numDim));
00058   face_tensor4   = rcp(new MDALayout<Cell,Face,Dim,Dim,Dim,Dim>(worksetSize,numFace,numDim,numDim,numDim,numDim));
00059 
00060   node_node_scalar = rcp(new MDALayout<Node>(worksetSize));
00061   node_node_vector = rcp(new MDALayout<Node,Dim>(worksetSize,vecDim));
00062   node_node_tensor = rcp(new MDALayout<Node,Dim,Dim>(worksetSize,numDim,numDim));
00063 
00064   // Coordinates:  3vector is for shells 2D topology 3 coordinates
00065   vertices_vector = rcp(new MDALayout<Cell,Vertex, Dim>(worksetSize,numVertices,numDim));
00066   node_3vector = rcp(new MDALayout<Cell,Node,Dim>(worksetSize,numNodes,3));
00067 
00068   // Basis Functions
00069   node_qp_scalar = rcp(new MDALayout<Cell,Node,QuadPoint>(worksetSize,numNodes, numQPts));
00070   node_qp_gradient = rcp(new MDALayout<Cell,Node,QuadPoint,Dim>(worksetSize,numNodes, numQPts,numDim));
00071   node_qp_vector =  node_qp_gradient;
00072 
00073   workset_scalar = rcp(new MDALayout<Dummy>(1));
00074   workset_vector = rcp(new MDALayout<Dim>(vecDim));
00075   workset_gradient = rcp(new MDALayout<Dim>(numDim));
00076   workset_tensor = rcp(new MDALayout<Dim,Dim>(numDim,numDim));
00077   workset_vecgradient = rcp(new MDALayout<Dim,Dim>(vecDim,numDim));
00078 
00079   shared_param = rcp(new MDALayout<Dim>(1));
00080   dummy = rcp(new MDALayout<Dummy>(0));
00081 
00082   // NOTE: vector data layouts here are used both for vectors fields
00083   // as well as gradients of scalar fields. This only works when
00084   // the vector fields are length numDim (which tends to be true
00085   // for displacements and velocities). Could be separated to
00086   // have separate node_vector and node_gradient layouts.
00087 }
00088 

Generated on Wed Mar 26 2014 18:36:36 for Albany: a Trilinos-based PDE code by  doxygen 1.7.1