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

Albany_STKNodeFieldContainer_Def.hpp

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 
00007 #include "Albany_STKNodeFieldContainer.hpp"
00008 
00009 
00010 #ifdef ALBANY_SEACAS
00011 #   include <stk_io/IossBridge.hpp>
00012 #endif
00013 
00014 
00015 #include <stk_mesh/base/FieldData.hpp>
00016 #include <stk_mesh/base/GetBuckets.hpp>
00017 #include <stk_mesh/fem/FEMMetaData.hpp>
00018 #include <stk_mesh/base/Types.hpp>
00019 #include <stk_util/parallel/Parallel.hpp>
00020 #include "Shards_Array.hpp"
00021 
00022 
00023 Teuchos::RCP<Albany::AbstractNodeFieldContainer>
00024 Albany::buildSTKNodeField(const std::string& name, const std::vector<int>& dim,
00025                           stk::mesh::fem::FEMMetaData* metaData, stk::mesh::BulkData* bulkData,
00026                           const bool output){
00027 
00028   switch(dim.size()){
00029 
00030   case 1: // scalar
00031     return Teuchos::rcp(new STKNodeField<double, 1>(name, dim, metaData, bulkData, output));
00032     break;
00033 
00034   case 2: // vector
00035     return Teuchos::rcp(new STKNodeField<double, 2>(name, dim, metaData, bulkData, output));
00036     break;
00037 
00038   case 3: // tensor
00039     return Teuchos::rcp(new STKNodeField<double, 3>(name, dim, metaData, bulkData, output));
00040     break;
00041 
00042   default:
00043     TEUCHOS_TEST_FOR_EXCEPT_MSG(true, "Error: unexpected argument for dimension");
00044   }
00045 }
00046 
00047 
00048 template<typename DataType, unsigned ArrayDim, class traits>
00049 Albany::STKNodeField<DataType, ArrayDim, traits>::STKNodeField(const std::string& name_,
00050      const std::vector<int>& dims_, stk::mesh::fem::FEMMetaData* metaData_, stk::mesh::BulkData* bulkData_,
00051      const bool output) :
00052   name(name_),
00053   dims(dims_),
00054   metaData(metaData_),
00055   bulkData(bulkData_)
00056 {
00057 
00058   node_field = traits_type::createField(name, dims, metaData_);
00059 
00060 #ifdef ALBANY_SEACAS
00061 
00062   if(output)
00063      stk::io::set_field_role(*node_field, Ioss::Field::TRANSIENT);
00064 
00065 #endif
00066 
00067 }
00068 
00069 template<typename DataType, unsigned ArrayDim, class traits>
00070 void
00071 Albany::STKNodeField<DataType, ArrayDim, traits>::saveField(const Teuchos::RCP<const Epetra_Vector>& block_mv,
00072         int offset, int blocksize){
00073 
00074  // Iterate over the processor-visible nodes
00075  stk::mesh::Selector select_owned_or_shared = metaData->locally_owned_part() | metaData->globally_shared_part();
00076 
00077  // Iterate over the overlap nodes by getting node buckets and iterating over each bucket.
00078  stk::mesh::BucketVector all_elements;
00079  stk::mesh::get_buckets(select_owned_or_shared, bulkData->buckets(metaData->node_rank()), all_elements);
00080 
00081  traits_type::saveFieldData(block_mv, all_elements, node_field, offset, blocksize);
00082 
00083 }
00084 
00085 template<typename DataType, unsigned ArrayDim, class traits>
00086 Albany::MDArray
00087 Albany::STKNodeField<DataType, ArrayDim, traits>::getMDA(const stk::mesh::Bucket& buck){
00088 
00089  stk::mesh::BucketArray<field_type> array(*node_field, buck);
00090 
00091  return array;
00092 
00093 }

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